Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
/* miscmaths.cc
Mark Jenkinson & Mark Woolrich & Christian Beckmann, FMRIB Image Analysis Group
Copyright (C) 1999-2000 University of Oxford */
/* CCOPYRIGHT */
// Miscellaneous maths functions
#include "miscmaths.h"
#include "stdlib.h"
using namespace std;
namespace MISCMATHS {
// The following lines are ignored by the current SGI compiler
// (version egcs-2.91.57)
// A temporary fix of including the std:: in front of all abs() etc
// has been done for now
using std::abs;
using std::sqrt;
using std::exp;
using std::log;
// using std::pow;
using std::atan2;
// General string/IO functions
bool isnum(const string& str)
{
// assumes that initial whitespace has been removed
if (isdigit(str[0])) return true;
if ( (str[0]=='-') || (str[0]=='+') || (str[0]=='.') ) return true;
return false;
}
string skip_alpha(ifstream& fs)
{
string cline;
while (!fs.eof()) {
getline(fs,cline);
cline += " "; // force extra entry in parsing
istrstream ss(cline.c_str());
string cc="";
ss >> cc;
if (isnum(cc)) {
fs.seekg(-((int)cline.size()),ios::cur);
return cline;
}
}
return "";
}
ReturnMatrix read_ascii_matrix(int nrows, int ncols, const string& filename)
{
return read_ascii_matrix(filename,nrows,ncols);
}
ReturnMatrix read_ascii_matrix(const string& filename, int nrows, int ncols)
{
Matrix mat(nrows,ncols);
mat = 0.0;
if ( filename.size()<1 ) return mat;
ifstream fs(filename.c_str());
if (!fs) {
cerr << "Could not open matrix file " << filename << endl;
return mat;
}
mat = read_ascii_matrix(fs,nrows,ncols);
fs.close();
mat.Release();
return mat;
}
ReturnMatrix read_ascii_matrix(int nrows, int ncols, ifstream& fs)
{
return read_ascii_matrix(fs, nrows, ncols);
}
ReturnMatrix read_ascii_matrix(ifstream& fs, int nrows, int ncols)
{
Matrix mat(nrows,ncols);
mat = 0.0;
string ss="";
ss = skip_alpha(fs);
for (int r=1; r<=nrows; r++) {
for (int c=1; c<=ncols; c++) {
if (!fs.eof()) {
fs >> ss;
while ( !isnum(ss) && !fs.eof() ) {
fs >> ss;
}
Mark Jenkinson
committed
mat(r,c) = atof(ss.c_str());
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
}
}
}
mat.Release();
return mat;
}
ReturnMatrix read_ascii_matrix(const string& filename)
{
Matrix mat;
if ( filename.size()<1 ) return mat;
ifstream fs(filename.c_str());
if (!fs) {
cerr << "Could not open matrix file " << filename << endl;
mat.Release();
return mat;
}
mat = read_ascii_matrix(fs);
fs.close();
mat.Release();
return mat;
}
ReturnMatrix read_ascii_matrix(ifstream& fs)
{
Matrix mat;
int rcount=0, cmax=0;
double val;
string cline;
// skip initial non-numeric lines
// and count the number of columns in the first numeric line
cline = skip_alpha(fs);
cline += " ";
{
istrstream ss(cline.c_str());
string cc="";
while (!ss.eof()) {
cmax++;
ss >> cc;
}
}
cmax--;
RowVector newrow(cmax);
while (!fs.eof()) {
getline(fs,cline);
cline += " "; // force extra entry in parsing
istrstream ss(cline.c_str());
string cc="";
ss >> cc;
if (!isnum(cc)) break; // stop processing when non-numeric line found
rcount++; // add new row to matrix
newrow = 0.0;
int ccount = 1;
do {
val = atof(cc.c_str());
if (ccount<=cmax) newrow(ccount++) = val;
ss >> cc;
} while (!ss.eof());
if (rcount>1) mat = mat & newrow;
else mat = newrow;
}
mat.Release();
return mat;
}
#define BINFLAG 42
ReturnMatrix read_binary_matrix(const string& filename)
{
Matrix mat;
if ( filename.size()<1 ) return mat;
ifstream fs(filename.c_str(), ios::in | ios::binary);
if (!fs) {
cerr << "Could not open matrix file " << filename << endl;
return mat;
}
mat = read_binary_matrix(fs);
fs.close();
mat.Release();
return mat;
}
ReturnMatrix read_binary_matrix(ifstream& fs)
{
long int testval;
// test for byte swapping
fs.read((char*)&testval,sizeof(testval));
if (testval!=BINFLAG) {
cerr << "Different endian format - byte swapping not currently implemented" << endl;
Matrix mres;
mres.Release();
return mres;
}
// read matrix dimensions (num rows x num cols)
long int ival,nx,ny;
fs.read((char*)&ival,sizeof(ival));
nx = ival;
fs.read((char*)&ival,sizeof(ival));
ny = ival;
// set up and read matrix (rows fast, cols slow)
Real val;
Matrix mres(nx,ny);
for (int y=1; y<=ny; y++) {
for (int x=1; x<=nx; x++) {
fs.read((char*)&val,sizeof(val));
mres(x,y)=val;
}
}
mres.Release();
return mres;
}
// WRITE FUNCTIONS //
int write_ascii_matrix(const string& filename, const Matrix& mat,
int precision)
{
return write_ascii_matrix(mat, filename, precision);
}
int write_ascii_matrix(const Matrix& mat, const string& filename,
int precision)
{
Tracer tr("write_ascii_matrix");
if ( (filename.size()<1) ) return -1;
ofstream fs(filename.c_str());
if (!fs) {
cerr << "Could not open file " << filename << " for writing" << endl;
return -1;
}
int retval = write_ascii_matrix(mat,fs,precision);
fs.close();
return retval;
}
int write_ascii_matrix(ofstream& fs, const Matrix& mat,
int precision)
{
return write_ascii_matrix(mat, fs, precision);
}
int write_ascii_matrix(const Matrix& mat, ofstream& fs, int precision)
{
for (int i=1; i<=mat.Nrows(); i++) {
for (int j=1; j<=mat.Ncols(); j++) {
if (precision>0) {
fs.setf(ios::scientific | ios::showpos);
fs.precision(precision+7);
fs.precision(precision); }
fs << mat(i,j) << " ";
}
fs << endl;
}
return 0;
}
int write_vest(string p_fname, const Matrix& x, int precision)
{ return write_vest(x,p_fname,precision); }
int write_vest(const Matrix& x, string p_fname, int precision)
{
ofstream out;
out.open(p_fname.c_str(), ios::out);
if(!out)
{
cerr << "Unable to open " << p_fname << endl;
return -1;
}
out << "! VEST-Waveform File" << endl;
out << "/NumWaves\t" << x.Ncols() << endl;
out << "/NumPoints\t" << x.Nrows() << endl;
out << "/Skip" << endl;
out << endl << "/Matrix" << endl;
int retval = write_ascii_matrix(x, out, precision);
out.close();
return retval;
}
int write_binary_matrix(const Matrix& mat, const string& filename)
{
Tracer tr("write_binary_matrix");
if ( (filename.size()<1) ) return -1;
ofstream fs(filename.c_str(), ios::out | ios::binary);
if (!fs) {
cerr << "Could not open file " << filename << " for writing" << endl;
return -1;
}
int retval = write_binary_matrix(mat,fs);
fs.close();
return retval;
}
int write_binary_matrix(const Matrix& mat, ofstream& fs)
{
long int ival, nx, ny;
ival = BINFLAG;
fs.write((char*)&ival,sizeof(ival));
ival = mat.Nrows();
fs.write((char*)&ival,sizeof(ival));
ival = mat.Ncols();
fs.write((char*)&ival,sizeof(ival));
nx = mat.Nrows();
ny = mat.Ncols();
Real val;
for (int y=1; y<=ny; y++) {
for (int x=1; x<=nx; x++) {
val = mat(x,y);
fs.write((char*)&val,sizeof(val));
}
}
return 0;
}
// General mathematical functions
int round(int x) { return x; }
int round(float x)
{
if (x>0.0) return ((int) (x+0.5));
else return ((int) (x-0.5));
}
int round(double x)
{
if (x>0.0) return ((int) (x+0.5));
else return ((int) (x-0.5));
}
int sign(int x) { return (x / abs(x)); }
int sign(float x) { return (int) (x / abs(x)); }
int sign(double x) { return (int) (x / abs(x)); }
int periodicclamp(int x, int x1, int x2)
{
if (x2<x1) return periodicclamp(x,x2,x1);
int d = x2-x1+1;
int xp = x-x1;
if (xp>=0) {
return (xp % d) + x1;
} else {
xp = xp + d + std::abs(xp/d)*d;
assert(xp>0);
return periodicclamp(xp + d + std::abs(xp/d)*d,x1,x2);
}
}
ColumnVector cross(const ColumnVector& a, const ColumnVector& b)
{
Tracer tr("cross");
ColumnVector ans(3);
ans(1) = a(2)*b(3) - a(3)*b(2);
ans(2) = a(3)*b(1) - a(1)*b(3);
ans(3) = a(1)*b(2) - a(2)*b(1);
return ans;
}
ColumnVector cross(const Real *a, const Real *b)
{
Tracer tr("cross");
ColumnVector a1(3), b1(3);
a1 << a;
b1 << b;
return cross(a1,b1);
}
double norm2(const ColumnVector& x)
{
return std::sqrt(x.SumSquare());
}
int Identity(Matrix& m)
{
Tracer tr("Identity");
m=0.0;
for (int j=1; j<=m.Nrows(); j++)
m(j,j)=1.0;
return 0;
}
ReturnMatrix Identity(int num)
{
Tracer tr("Identity");
Matrix eye(num,num);
Identity(eye);
eye.Release();
return eye;
}
int diag(Matrix& m, const float diagvals[])
{
Tracer tr("diag");
m=0.0;
for (int j=1; j<=m.Nrows(); j++)
m(j,j)=diagvals[j-1];
return 0;
}
int diag(Matrix& m, const ColumnVector& diagvals)
{
Tracer tr("diag");
m=0.0;
for (int j=1; j<=m.Nrows(); j++)
m(j,j)=diagvals(j);
return 0;
}
ReturnMatrix diag(const Matrix& Mat)
{
Tracer tr("diag");
if(Mat.Ncols()==1){
Matrix retmat(Mat.Nrows(),Mat.Nrows());
diag(retmat,Mat);
retmat.Release();
return retmat;}
else{
int mindim = Min(Mat.Ncols(),Mat.Nrows());
Matrix retmat(mindim,1);
for(int ctr=1; ctr<=mindim;ctr++){
retmat(ctr,1)=Mat(ctr,ctr);
}
retmat.Release();
return retmat;
}
}
ReturnMatrix pinv(const Matrix& mat)
{
// calculates the psuedo-inverse using SVD
Tracer tr("pinv");
DiagonalMatrix D;
Matrix U, V;
SVD(mat,D,U,V);
float tol;
tol = MaximumAbsoluteValue(D) * Max(mat.Nrows(),mat.Ncols()) * 1e-16;
for (int n=1; n<=D.Nrows(); n++) {
if (fabs(D(n,n))>tol) D(n,n) = 1.0/D(n,n);
}
Matrix pinv = V * D * U.t();
pinv.Release();
return pinv;
}
ReturnMatrix sqrtaff(const Matrix& mat)
{
Tracer tr("sqrtaff");
Matrix matnew(4,4), rot(4,4), id4(4,4);
Identity(rot);
Identity(id4);
ColumnVector params(12), centre(3), trans(4);
centre = 0.0;
// Quaternion decomposition -> params(1..3) = sin(theta/2)*(unit_axis_vec)
// Want a new quaternion : q = sin(theta/4)*(unit_axis_vec)
// Therefore factor of conversion is: factor = sin(theta/4)/sin(theta/2)
// = 1/(2 * cos(theta/4)) which is calculated below
// NB: t = theta/2
decompose_aff(params,mat,centre,rotmat2quat);
double sint;
sint = std::sqrt(params(1)*params(1) + params(2)*params(2) +
params(3)*params(3));
double t = asin(sint);
double factor = 1.0/(2.0*cos(0.5*t));
params(1) = factor * params(1);
params(2) = factor * params(2);
params(3) = factor * params(3);
params(7) = std::sqrt(params(7));
params(8) = std::sqrt(params(8));
params(9) = std::sqrt(params(9));
params(10) = 0.5*params(10);
params(11) = 0.5*params(11);
params(12) = 0.5*params(12);
construct_rotmat_quat(params,3,rot,centre);
rot(1,4) = 0.0;
rot(2,4) = 0.0;
rot(3,4) = 0.0;
Matrix scale(4,4);
Identity(scale);
scale(1,1)=params(7);
scale(2,2)=params(8);
scale(3,3)=params(9);
Matrix skew(4,4);
Identity(skew);
skew(1,2)=params(10);
skew(1,3)=params(11);
skew(2,3)=params(12);
trans(1) = params(4);
trans(2) = params(5);
trans(3) = params(6);
trans(4) = 1.0;
// The translation, being independent of the 3x3 submatrix, is
// calculated so that it will be equal for each of the two
// halves of the approximate square root
// (i.e. matnew and mat*matnew.i() have exactly the same translation)
ColumnVector th(4);
th = (mat*scale.i()*skew.i()*rot.i() + id4).SubMatrix(1,3,1,3).i()
* trans.SubMatrix(1,3,1,1);
matnew = rot*skew*scale;
matnew(1,4) = th(1);
matnew(2,4) = th(2);
matnew(3,4) = th(3);
matnew.Release();
return matnew;
}
//------------------------------------------------------------------------//
// Handy MATLAB-like functions
void reshape(Matrix& r, const Matrix& m, int nrows, int ncols)
{
Tracer tr("reshape");
if (nrows*ncols != m.Nrows() * m.Ncols() ) {
cerr << "WARNING: cannot reshape " << m.Nrows() << "x"
<< m.Ncols() << " matrix into " << nrows << "x"
<< ncols << endl;
cerr << " Returning original matrix instead" << endl;
r = m;
return;
}
r.ReSize(nrows,ncols);
int rr = 1, rc = 1;
for (int mc=1; mc<=m.Ncols(); mc++) {
for (int mr=1; mr<=m.Nrows(); mr++) {
r(rr,rc) = m(mr,mc);
rr++;
if (rr>nrows) {
rc++;
rr=1;
}
}
}
}
ReturnMatrix reshape(const Matrix& m, int nrows, int ncols)
{
Tracer tr("reshape");
Matrix r;
reshape(r,m,nrows,ncols);
r.Release();
return r;
}
//------------------------------------------------------------------------//
// Spatial transformation functions (rotations and affine transforms)
int construct_rotmat_euler(const ColumnVector& params, int n, Matrix& aff,
const ColumnVector& centre)
{
Tracer tr("construct_rotmat_euler");
ColumnVector angl(3);
Matrix newaff(4,4);
Identity(aff);
if (n<=0) return 0;
// order of parameters is 3 rotation + 3 translation
// angles are in radians
// order of parameters is (Rx,Ry,Rz) and R = Rx.Ry.Rz
angl=0.0;
angl(1)=params(1);
make_rot(angl,centre,newaff);
aff = aff * newaff;
if (n==1) return 0;
angl=0.0;
angl(2)=params(2);
make_rot(angl,centre,newaff);
aff = aff * newaff;
if (n==2) return 0;
angl=0.0;
angl(3)=params(3);
make_rot(angl,centre,newaff);
aff = aff * newaff;
if (n==3) return 0;
aff(1,4)+=params(4);
if (n==4) return 0;
aff(2,4)+=params(5);
if (n==5) return 0;
aff(3,4)+=params(6);
if (n==6) return 0;
return 1;
}
int construct_rotmat_euler(const ColumnVector& params, int n, Matrix& aff)
{
Tracer tr("construct_rotmat_euler");
ColumnVector centre(3);
centre = 0.0;
return construct_rotmat_euler(params,n,aff,centre);
}
int construct_rotmat_quat(const ColumnVector& params, int n, Matrix& aff,
const ColumnVector& centre)
{
Tracer tr("construct_rotmat_quat");
Identity(aff);
if (n<=0) return 0;
// order of parameters is 3 rotation (last 3 quaternion components)
// + 3 translation
if ((n>=1) && (n<3)) { cerr<<"Can only do 3 or more, not "<< n <<endl; }
float w, w2 = 1.0 - Sqr(params(1)) - Sqr(params(2)) - Sqr(params(3));
if (w2 < 0.0) {
cerr << "Parameters do not form a valid axis - greater than unity\n";
return -1;
}
w = std::sqrt(w2);
float x=params(1), y=params(2), z=params(3);
aff(1,1) = 1 - 2*y*y - 2*z*z;
aff(2,2) = 1 - 2*x*x - 2*z*z;
aff(3,3) = 1 - 2*x*x - 2*y*y;
aff(1,2) = 2*x*y - 2*w*z;
aff(2,1) = 2*x*y + 2*w*z;
aff(1,3) = 2*x*z + 2*w*y;
aff(3,1) = 2*x*z - 2*w*y;
aff(2,3) = 2*y*z - 2*w*x;
aff(3,2) = 2*y*z + 2*w*x;
// Given Rotation matrix R: x' = Rx + (I-R)*centre
ColumnVector trans(3);
trans = aff.SubMatrix(1,3,1,3)*centre;
aff.SubMatrix(1,3,4,4) = centre - trans;
aff(1,4) += params(4);
if (n==4) return 0;
aff(2,4) += params(5);
if (n==5) return 0;
aff(3,4) += params(6);
if (n==6) return 0;
return 1;
}
int construct_rotmat_quat(const ColumnVector& params, int n, Matrix& aff)
{
Tracer tr("construct_rotmat_quat");
ColumnVector centre(3);
centre = 0.0;
return construct_rotmat_quat(params,n,aff,centre);
}
int make_rot(const ColumnVector& angl, const ColumnVector& centre,
Matrix& rot)
{
// Matrix rot must be 4x4; angl and orig must be length 3
Tracer tr("make_rot");
Identity(rot); // default return value
float theta;
theta = norm2(angl);
if (theta<1e-8) { // avoid round-off errors and return Identity
return 0;
}
ColumnVector axis = angl/theta;
ColumnVector x1(3), x2(3), x3(3);
x1 = axis;
x2(1) = -axis(2); x2(2) = axis(1); x2(3) = 0.0;
if (norm2(x2)<=0.0) {
x2(1) = 1.0; x2(2) = 0.0; x2(3) = 0.0;
}
x2 = x2/norm2(x2);
x3 = cross(x1,x2);
x3 = x3/norm2(x3);
Matrix basischange(3,3);
basischange.SubMatrix(1,3,1,1) = x2;
basischange.SubMatrix(1,3,2,2) = x3;
basischange.SubMatrix(1,3,3,3) = x1;
Matrix rotcore(3,3);
Identity(rotcore);
rotcore(1,1)=cos(theta);
rotcore(2,2)=cos(theta);
rotcore(1,2)=sin(theta);
rotcore(2,1)=-sin(theta);
rot.SubMatrix(1,3,1,3) = basischange * rotcore * basischange.t();
Matrix ident3(3,3);
Identity(ident3);
ColumnVector trans(3);
trans = (ident3 - rot.SubMatrix(1,3,1,3))*centre;
rot.SubMatrix(1,3,4,4)=trans;
return 0;
}
int getrotaxis(ColumnVector& axis, const Matrix& rotmat)
{
Tracer tr("getrotaxis");
Matrix residuals(3,3);
residuals = rotmat*rotmat.t() - Identity(3);
if (residuals.SumSquare() > 1e-4)
{ cerr << "Failed orthogonality check!" << endl; return -1; }
Matrix u(3,3), v(3,3);
DiagonalMatrix d(3);
SVD(rotmat-Identity(3),d,u,v);
// return column of V corresponding to minimum value of |S|
for (int i=1; i<=3; i++) {
if (fabs(d(i))<1e-4) axis = v.SubMatrix(1,3,i,i);
}
return 0;
}
int rotmat2euler(ColumnVector& angles, const Matrix& rotmat)
{
// uses the convention that R = Rx.Ry.Rz
Tracer tr("rotmat2euler");
float cz, sz, cy, sy, cx, sx;
cy = std::sqrt(Sqr(rotmat(1,1)) + Sqr(rotmat(1,2)));
if (cy < 1e-4) {
//cerr << "Cos y is too small - Gimbal lock condition..." << endl;
cx = rotmat(2,2);
sx = -rotmat(3,2);
sy = -rotmat(1,3);
angles(1) = atan2(sx,cx);
angles(2) = atan2(sy,(float)0.0);
angles(3) = 0.0;
} else {
// choose by convention that cy > 0
// get the same rotation if: sy stays same & all other values swap sign
cz = rotmat(1,1)/cy;
sz = rotmat(1,2)/cy;
cx = rotmat(3,3)/cy;
sx = rotmat(2,3)/cy;
sy = -rotmat(1,3);
//atan2(sin,cos) (defined as atan2(y,x))
angles(1) = atan2(sx,cx);
angles(2) = atan2(sy,cy);
angles(3) = atan2(sz,cz);
}
return 0;
}
int rotmat2quat(ColumnVector& quaternion, const Matrix& rotmat)
{
Tracer tr("rotmat2quat");
float trace = rotmat.SubMatrix(1,3,1,3).Trace();
if (trace > 0) {
float w = std::sqrt((trace + 1.0)/4.0);
quaternion(1) = (rotmat(3,2) - rotmat(2,3))/(4.0*w);
quaternion(2) = (rotmat(1,3) - rotmat(3,1))/(4.0*w);
quaternion(3) = (rotmat(2,1) - rotmat(1,2))/(4.0*w);
} else if ((rotmat(1,1) > rotmat(2,2)) && (rotmat(1,1) > rotmat(3,3))) {
// first col case
float s = std::sqrt(1.0 + rotmat(1,1) - rotmat(2,2) - rotmat(3,3)) * 2.0;
quaternion(1) = 0.5 / s;
quaternion(2) = (-rotmat(1,2) - rotmat(1,2)) / s;
quaternion(3) = (-rotmat(1,3) - rotmat(3,1)) / s;
} else if ((rotmat(2,2) > rotmat(1,1)) && (rotmat(2,2) > rotmat(3,3))) {
// 2nd col case
float s = std::sqrt(1.0 + rotmat(2,2) - rotmat(1,1) - rotmat(3,3)) * 2.0;
quaternion(1) = (-rotmat(1,2) - rotmat(2,1)) / s;
quaternion(2) = 0.5 / s;
quaternion(3) = (-rotmat(2,3) - rotmat(3,2)) / s;
} else if ((rotmat(3,3) > rotmat(1,1)) && (rotmat(3,3) > rotmat(2,2))) {
// 3rd col case
float s = std::sqrt(1.0 + rotmat(3,3) - rotmat(1,1) - rotmat(2,2)) * 2.0;
quaternion(1) = (-rotmat(1,3) - rotmat(3,1)) / s;
quaternion(2) = (-rotmat(2,3) - rotmat(3,2)) / s;
quaternion(3) = 0.5 / s;
}
return 0;
}
int decompose_aff(ColumnVector& params, const Matrix& affmat,
const ColumnVector& centre,
int (*rotmat2params)(ColumnVector& , const Matrix& ))
{
// decomposes using the convention: mat = rotmat * skew * scale
// order of parameters is 3 rotation + 3 translation + 3 scales + 3 skews
// angles are in radians
Tracer tr("decompose_aff");
if (params. Nrows() < 12)
params.ReSize(12);
if (rotmat2params==0)
{
cerr << "No rotmat2params function specified" << endl;
return -1;
}
ColumnVector x(3), y(3), z(3);
Matrix aff3(3,3);
aff3 = affmat.SubMatrix(1,3,1,3);
x = affmat.SubMatrix(1,3,1,1);
y = affmat.SubMatrix(1,3,2,2);
z = affmat.SubMatrix(1,3,3,3);
float sx, sy, sz, a, b, c;
sx = norm2(x);
sy = std::sqrt( dot(y,y) - (Sqr(dot(x,y)) / Sqr(sx)) );
a = dot(x,y)/(sx*sy);
ColumnVector x0(3), y0(3);
x0 = x/sx;
y0 = y/sy - a*x0;
sz = std::sqrt(dot(z,z) - Sqr(dot(x0,z)) - Sqr(dot(y0,z)));
b = dot(x0,z)/sz;
c = dot(y0,z)/sz;
params(7) = sx; params(8) = sy; params(9) = sz;
Matrix scales(3,3);
float diagvals[] = {sx,sy,sz};
diag(scales,diagvals);
Real skewvals[] = {1,a,b,0 , 0,1,c,0 , 0,0,1,0 , 0,0,0,1};
Matrix skew(4,4);
skew << skewvals;
params(10) = a; params(11) = b; params(12) = c;
Matrix rotmat(3,3);
rotmat = aff3 * scales.i() * (skew.SubMatrix(1,3,1,3)).i();
ColumnVector transl(3);
//transl = affmat.SubMatrix(1,3,4,4);
//transl = transl - (Identity(3) - rotmat)*centre;
transl = affmat.SubMatrix(1,3,1,3)*centre + affmat.SubMatrix(1,3,4,4)
- centre;
for (int i=1; i<=3; i++) { params(i+3) = transl(i); }
ColumnVector rotparams(3);
(*rotmat2params)(rotparams,rotmat);
for (int i=1; i<=3; i++) { params(i) = rotparams(i); }
return 0;
}
int decompose_aff(ColumnVector& params, const Matrix& affmat,
int (*rotmat2params)(ColumnVector& , const Matrix& ))
{
Tracer tr("decompose_aff");
ColumnVector centre(3);
centre = 0.0;
return decompose_aff(params,affmat,centre,rotmat2params);
}
int compose_aff(const ColumnVector& params, int n, const ColumnVector& centre,
Matrix& aff,
int (*params2rotmat)(const ColumnVector& , int , Matrix& ,
const ColumnVector& ) )
{
Tracer tr("compose_aff");
if (n<=0) return 0;
// order of parameters is 3 rotation + 3 translation + 3 scales + 3 skews
// angles are in radians
(*params2rotmat)(params,n,aff,centre);
if (n<=6) return 0;
Matrix scale(4,4);
Identity(scale);
if (n>=7) {
scale(1,1)=params(7);
if (n>=8) scale(2,2)=params(8);
else scale(2,2)=params(7);
if (n>=9) scale(3,3)=params(9);
else scale(3,3)=params(7);
}
// fix the translation so that the centre is not moved
ColumnVector strans(3);
strans = centre - scale.SubMatrix(1,3,1,3)*centre;
scale.SubMatrix(1,3,4,4) = strans;
Matrix skew(4,4);
Identity(skew);
if (n>=10) {
if (n>=10) skew(1,2)=params(10);
if (n>=11) skew(1,3)=params(11);
if (n>=12) skew(2,3)=params(12);
}
// fix the translation so that the centre is not moved
ColumnVector ktrans(3);
ktrans = centre - skew.SubMatrix(1,3,1,3)*centre;
skew.SubMatrix(1,3,4,4) = ktrans;
aff = aff * skew * scale;
return 0;
}
float rms_deviation(const Matrix& affmat1, const Matrix& affmat2,
const ColumnVector& centre, const float rmax)
{
Tracer trcr("rms_deviation");
Matrix isodiff(4,4);
try {
isodiff = affmat1*affmat2.i() - Identity(4);
} catch(...) {
cerr << "RMS_DEVIATION ERROR:: Could not invert matrix" << endl;
exit(-5);
}
Matrix adiff(3,3);
adiff = isodiff.SubMatrix(1,3,1,3);
ColumnVector tr(3);
tr = isodiff.SubMatrix(1,3,4,4) + adiff*centre;
float rms = std::sqrt( (tr.t() * tr).AsScalar() +
(rmax*rmax/5.0)*Trace(adiff.t()*adiff) );
return rms;
}
float rms_deviation(const Matrix& affmat1, const Matrix& affmat2,
const float rmax)
{
ColumnVector centre(3);
centre = 0;
return rms_deviation(affmat1,affmat2,centre,rmax);
}
// Added by MWW
// int getdiag(ColumnVector& diagvals, const Matrix& m)
// {
// Tracer ts("MiscMaths::diag");
// int num = m.Nrows();
// diagvals.ReSize(num);
// for (int j=1; j<=num; j++)
// diagvals(j)=m(j,j);
// return 0;
// }
// float var(const ColumnVector& x)
// {
// float m = mean(x);
// float ssq = (x-m).SumSquare()/(x.Nrows()-1);
// return ssq;
// }
// float mean(const ColumnVector& x)
// {
// float m = x.Sum()/x.Nrows();
// return m;
// }
float median(const ColumnVector& x)
{
ColumnVector y = x;
SortAscending(y);
float m = y((int)(y.Nrows()/2));
return m;
}
// Added by CFB --- Matlab style Matrix functions
ReturnMatrix ones(const int dim1, const int dim2)
{
int tdim = dim2;