Newer
Older
Mark Jenkinson & Mark Woolrich & Christian Beckmann & Tim Behrens, FMRIB Image Analysis Group
Copyright (C) 1999-2000 University of Oxford */
/* CCOPYRIGHT */
// Miscellaneous maths functions
#include "miscmaths.h"
#include "newmatio.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;
string size(const Matrix& mat)
{
string str = num2str(mat.Nrows())+"*"+num2str(mat.Ncols());
return str;
}
float Sinc(const float x) {
if (fabs(x)<1e-9) {
return 1-x*x*M_PI*M_PI/6.0;
} else {
return sin(M_PI*x)/(M_PI*x);
}
}
double Sinc(const double x) {
if (fabs(x)<1e-9) {
return 1-x*x*M_PI*M_PI/6.0;
} else {
return sin(M_PI*x)/(M_PI*x);
}
}
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
98
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
// 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());
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
}
}
}
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)
{
int rcount=0, cmax=0;
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--;
Mark Jenkinson
committed
do {
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
Mark Jenkinson
committed
} while (!fs.eof());
Mark Jenkinson
committed
// now know the size of matrix
fs.clear();
fs.seekg(0,ios::beg);
return read_ascii_matrix(fs,rcount,cmax);
}
#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)
{
Mark Jenkinson
committed
bool swapbytes = false;
unsigned int testval;
// test for byte swapping
fs.read((char*)&testval,sizeof(testval));
if (testval!=BINFLAG) {
Mark Jenkinson
committed
swapbytes = true;
Swap_Nbytes(1,sizeof(testval),&testval);
if (testval!=BINFLAG) {
cerr << "Unrecognised binary matrix file format" << endl;
Matrix mres;
mres.Release();
return mres;
}
}
// read matrix dimensions (num rows x num cols)
Mark Jenkinson
committed
unsigned int ival,nx,ny;
// ignore the padding (reserved for future use)
fs.read((char*)&ival,sizeof(ival));
Mark Jenkinson
committed
if (swapbytes) Swap_Nbytes(1,sizeof(ival),&ival);
nx = ival;
fs.read((char*)&ival,sizeof(ival));
Mark Jenkinson
committed
if (swapbytes) Swap_Nbytes(1,sizeof(ival),&ival);
ny = ival;
// set up and read matrix (rows fast, cols slow)
Mark Jenkinson
committed
double val;
Mark Jenkinson
committed
for (unsigned int y=1; y<=ny; y++) {
for (unsigned int x=1; x<=nx; x++) {
Mark Jenkinson
committed
if (swapbytes) Swap_Nbytes(1,sizeof(val),&val);
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
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)
{
Mark Jenkinson
committed
unsigned int ival, nx, ny;
ival = BINFLAG;
fs.write((char*)&ival,sizeof(ival));
ival = 0; // padding (reserved for future use)
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();
Mark Jenkinson
committed
double val;
for (unsigned int y=1; y<=ny; y++) {
for (unsigned 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));
double rounddouble(double x){
return ( floor(x+0.5));
}
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
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());
}
double norm2(double a, double b, double c)
{
return a*a + b*b + c*c;
}
float norm2(float a, float b, float c)
{
return a*a + b*b + c*c;
}
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
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(DiagonalMatrix& m, const ColumnVector& diagvals)
{
Tracer tr("diag");
m.ReSize(diagvals.Nrows());
m=0.0;
for (int j=1; j<=diagvals.Nrows(); j++)
m(j)=diagvals(j);
return 0;
}
int diag(Matrix& m, const ColumnVector& diagvals)
{
Tracer tr("diag");
m.ReSize(diagvals.Nrows(),diagvals.Nrows());
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
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 rpinv(const Matrix& mat) // note that mat must be a fat matrix (e.g. x.t() )
{
// calculates the right psuedo-inverse using SVD
Tracer tr("rpinv");
DiagonalMatrix D;
Matrix U, V;
SVD(mat.t(),D,U,V); // feed transpose of input into SVD
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 rpinv = U * D * V.t();
rpinv.Release();
return rpinv;
}
int rank(const Matrix& X)
{
// calculates the rank of matrix X
Tracer tr("rank");
DiagonalMatrix eigenvals;
SVD(X,eigenvals);
double tolerance = Max(X.Nrows(),X.Ncols()) * eigenvals.Maximum() * 1e-16;
int therank=0;
for(int i=0; i<eigenvals.Nrows(); i++)
if (eigenvals(i+1)>tolerance)
therank++;
// cout << "tolerance = " << tolerance << "\n" << "eigenvalues = " << eigenvals << "\n" << "rank = " << therank << endl;
return therank;
}
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
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;
}
void cart2sph(const ColumnVector& dir, float& th, float& ph)
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
{
float mag=sqrt(dir(1)*dir(1)+dir(2)*dir(2)+dir(3)*dir(3));
if(mag==0){
ph=M_PI/2;
th=M_PI/2;
}
else{
if(dir(1)==0 && dir(2)>=0) ph=M_PI/2;
else if(dir(1)==0 && dir(2)<0) ph=-M_PI/2;
else if(dir(1)>0) ph=atan(dir(2)/dir(1));
else if(dir(2)>0) ph=atan(dir(2)/dir(1))+M_PI;
else ph=atan(dir(2)/dir(1))-M_PI;
if(dir(3)==0) th=M_PI/2;
else if(dir(3)>0) th=atan(sqrt(dir(1)*dir(1)+dir(2)*dir(2))/dir(3));
else th=atan(sqrt(dir(1)*dir(1)+dir(2)*dir(2))/dir(3))+M_PI;
}
}
void cart2sph(const Matrix& dir,ColumnVector& th,ColumnVector& ph)
{
if(th.Nrows()!=dir.Ncols()){
th.ReSize(dir.Ncols());
}
if(ph.Nrows()!=dir.Ncols()){
ph.ReSize(dir.Ncols());
}
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
for (int i=1;i<=dir.Ncols();i++) {
float mag=sqrt(dir(1,i)*dir(1,i)+dir(2,i)*dir(2,i)+dir(3,i)*dir(3,i));
if(mag==0){
ph(i)=M_PI/2;
th(i)=M_PI/2;
}
else{
if(dir(1,i)==0 && dir(2,i)>=0) ph(i)=M_PI/2;
else if(dir(1,i)==0 && dir(2,i)<0) ph(i)=-M_PI/2;
else if(dir(1,i)>0) ph(i)=atan(dir(2,i)/dir(1,i));
else if(dir(2,i)>0) ph(i)=atan(dir(2,i)/dir(1,i))+M_PI;
else ph(i)=atan(dir(2,i)/dir(1,i))-M_PI;
if(dir(3,i)==0) th(i)=M_PI/2;
else if(dir(3,i)>0) th(i)=atan(sqrt(dir(1,i)*dir(1,i)+dir(2,i)*dir(2,i))/dir(3,i));
else th(i)=atan(sqrt(dir(1,i)*dir(1,i)+dir(2,i)*dir(2,i))/dir(3,i))+M_PI;
}
}
}
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
// Added by CFB --- Matlab style Matrix functions
ReturnMatrix ones(const int dim1, const int dim2)
{
int tdim = dim2;
if(tdim<0){tdim=dim1;}
Matrix res(dim1,tdim); res = 1.0;
res.Release();
return res;
}
ReturnMatrix zeros(const int dim1, const int dim2)
{
int tdim = dim2;
if(tdim<0){tdim=dim1;}
Matrix res(dim1,tdim); res = 0.0;
res.Release();
return res;
}
ReturnMatrix repmat(const Matrix &mat, const int rows, const int cols)
{
Matrix res = mat;
for(int ctr = 1; ctr < cols; ctr++){res |= mat;}
Matrix tmpres = res;
for(int ctr = 1; ctr < rows; ctr++){res &= tmpres;}
res.Release();
return res;
}
ReturnMatrix dist2(const Matrix &mat1, const Matrix &mat2)
{
Matrix res(mat1.Ncols(),mat2.Ncols());
for(int ctr1 = 1; ctr1 <= mat1.Ncols(); ctr1++)
for(int ctr2 =1; ctr2 <= mat2.Ncols(); ctr2++)
{
ColumnVector tmp;
tmp=mat1.Column(ctr1)-mat2.Column(ctr2);
res(ctr1,ctr2) = std::sqrt(tmp.SumSquare());
}
res.Release();
return res;
}
ReturnMatrix abs(const Matrix& mat)
{
Matrix res = mat;
for (int mc=1; mc<=mat.Ncols(); mc++) {
for (int mr=1; mr<=mat.Nrows(); mr++) {
res(mr,mc)=std::abs(res(mr,mc));
}
}
res.Release();
return res;
}
ReturnMatrix sqrt(const Matrix& mat)
{
Matrix res = mat;
bool neg_flag = false;
for (int mc=1; mc<=mat.Ncols(); mc++) {
for (int mr=1; mr<=mat.Nrows(); mr++) {
if(res(mr,mc)<0){ neg_flag = true; }
res(mr,mc)=std::sqrt(std::abs(res(mr,mc)));
}
}
if(neg_flag){
//cerr << " Matrix contained negative elements" << endl;
//cerr << " return sqrt(abs(X)) instead" << endl;
}
res.Release();
return res;
}
ReturnMatrix log(const Matrix& mat)
{
Matrix res = mat;
bool neg_flag = false;
for (int mc=1; mc<=mat.Ncols(); mc++) {
for (int mr=1; mr<=mat.Nrows(); mr++) {
if(res(mr,mc)<0){ neg_flag = true; }
res(mr,mc)=std::log(std::abs(res(mr,mc)));
}
}
if(neg_flag){
// cerr << " Matrix contained negative elements" << endl;
// cerr << " return log(abs(X)) instead" << endl;
}
res.Release();
return res;
}
ReturnMatrix exp(const Matrix& mat)
{
Matrix res = mat;
for (int mc=1; mc<=mat.Ncols(); mc++) {
for (int mr=1; mr<=mat.Nrows(); mr++) {
res(mr,mc)=std::exp(res(mr,mc));
}
}
res.Release();
return res;
}
ReturnMatrix tanh(const Matrix& mat)
{
Matrix res = mat;
for (int mc=1; mc<=mat.Ncols(); mc++) {
for (int mr=1; mr<=mat.Nrows(); mr++) {
res(mr,mc)=std::tanh(res(mr,mc));
}
}
res.Release();
return res;
}
ReturnMatrix pow(const Matrix& mat, const double exp)
{
Matrix res = mat;
for (int mc=1; mc<=mat.Ncols(); mc++) {
for (int mr=1; mr<=mat.Nrows(); mr++) {
res(mr,mc)=std::pow(res(mr,mc),exp);
}
}
res.Release();
return res;
}
ReturnMatrix max(const Matrix& mat)
{
Matrix res;
if(mat.Nrows()>1){
res=zeros(1,mat.Ncols());
res=mat.Row(1);
for(int mc=1; mc<=mat.Ncols();mc++){
for(int mr=2; mr<=mat.Nrows();mr++){
if(mat(mr,mc)>res(1,mc)){res(1,mc)=mat(mr,mc);}
}
}
}
else{
res=zeros(1);
res=mat(1,1);
for(int mc=2; mc<=mat.Ncols(); mc++){
if(mat(1,mc)>res(1,1)){res(1,1)=mat(1,mc);}
}
}
res.Release();
return res;
}
ReturnMatrix max(const Matrix& mat,ColumnVector& index)
{
index.ReSize(mat.Nrows());
index=1;
Matrix res;
if(mat.Nrows()>1){
res=zeros(1,mat.Ncols());
res=mat.Row(1);
for(int mc=1; mc<=mat.Ncols();mc++){
for(int mr=2; mr<=mat.Nrows();mr++){
if(mat(mr,mc)>res(1,mc))
{
res(1,mc)=mat(mr,mc);
index(mr)=mc;
}
}
}
}
else{
res=zeros(1);
res=mat(1,1);
for(int mc=2; mc<=mat.Ncols(); mc++){
if(mat(1,mc)>res(1,1))
{
res(1,1)=mat(1,mc);
index(1)=mc;
}
}
}
res.Release();
return res;
}
ReturnMatrix min(const Matrix& mat)
{
Matrix res;
if(mat.Nrows()>1){
res=zeros(1,mat.Ncols());
res=mat.Row(1);
for(int mc=1; mc<=mat.Ncols();mc++){
for(int mr=2; mr<=mat.Nrows();mr++){
if(mat(mr,mc)<res(1,mc)){res(1,mc)=mat(mr,mc);}
}
}
}
else{
res=zeros(1);
res=mat(1,1);
for(int mc=2; mc<=mat.Ncols(); mc++){
if(mat(1,mc)<res(1,1)){res(1,1)=mat(1,mc);}
}
}
res.Release();
return res;
}
ReturnMatrix sum(const Matrix& mat, const int dim)
{
Matrix tmp;
if (dim == 1) {tmp=mat;}
else {tmp=mat.t();}
Matrix res(1,tmp.Ncols());
res = 0.0;
for (int mc=1; mc<=tmp.Ncols(); mc++) {
for (int mr=1; mr<=tmp.Nrows(); mr++) {
res(1,mc) += tmp(mr,mc);
}
}
if (!(dim == 1)) {res=res.t();}
res.Release();
return res;
}
ReturnMatrix mean(const Matrix& mat, const int dim)
{
Matrix tmp;
if (dim == 1) {tmp=mat;}
else {tmp=mat.t();}
int N = tmp.Nrows();
Matrix res(1,tmp.Ncols());
res = 0.0;
for (int mc=1; mc<=tmp.Ncols(); mc++) {
for (int mr=1; mr<=tmp.Nrows(); mr++) {
res(1,mc) += tmp(mr,mc)/N;
}
}
if (!(dim == 1)) {res=res.t();}
res.Release();
return res;
}
ReturnMatrix var(const Matrix& mat, const int dim)
{
Matrix tmp;
if (dim == 1) {tmp=mat;}
else {tmp=mat.t();}
int N = tmp.Nrows();
Matrix res(1,tmp.Ncols());
res = 0.0;
if(N>1){
tmp -= ones(tmp.Nrows(),1)*mean(tmp,1);
for (int mc=1; mc<=tmp.Ncols(); mc++)
for (int mr=1; mr<=tmp.Nrows(); mr++)
res(1,mc) += tmp(mr,mc) / (N-1) * tmp(mr,mc);
}
if (!(dim == 1)) {res=res.t();}
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
res.Release();
return res;
}
ReturnMatrix stdev(const Matrix& mat, const int dim)
{
return sqrt(var(mat,dim));
}
ReturnMatrix gt(const Matrix& mat1,const Matrix& mat2)
{
int ctrcol = std::min(mat1.Ncols(),mat2.Ncols());
int ctrrow = std::min(mat1.Nrows(),mat2.Nrows());
Matrix res(ctrrow,ctrcol);
res=0.0;
for (int ctr1 = 1; ctr1 <= ctrrow; ctr1++) {
for (int ctr2 =1; ctr2 <= ctrcol; ctr2++) {
if( mat1(ctr1,ctr2) > mat2(ctr1,ctr2)){
res(ctr1,ctr2) = 1.0;
}
}
}
res.Release();
return res;
}
ReturnMatrix lt(const Matrix& mat1,const Matrix& mat2)
{
int ctrcol = std::min(mat1.Ncols(),mat2.Ncols());
int ctrrow = std::min(mat1.Nrows(),mat2.Nrows());
Matrix res(ctrrow,ctrcol);
res=0.0;
for (int ctr1 = 1; ctr1 <= ctrrow; ctr1++) {
for (int ctr2 =1; ctr2 <= ctrcol; ctr2++) {
if( mat1(ctr1,ctr2) < mat2(ctr1,ctr2)){
res(ctr1,ctr2) = 1.0;
}
}
}
res.Release();
return res;
}
ReturnMatrix geqt(const Matrix& mat1,const Matrix& mat2)
{
int ctrcol = std::min(mat1.Ncols(),mat2.Ncols());
int ctrrow = std::min(mat1.Nrows(),mat2.Nrows());
Matrix res(ctrrow,ctrcol);
res=0.0;
for (int ctr1 = 1; ctr1 <= ctrrow; ctr1++) {
for (int ctr2 =1; ctr2 <= ctrcol; ctr2++) {
if( mat1(ctr1,ctr2) >= mat2(ctr1,ctr2)){
res(ctr1,ctr2) = 1.0;
}
}
}
res.Release();
return res;
}
ReturnMatrix leqt(const Matrix& mat1,const Matrix& mat2)
{
int ctrcol = std::min(mat1.Ncols(),mat2.Ncols());
int ctrrow = std::min(mat1.Nrows(),mat2.Nrows());
Matrix res(ctrrow,ctrcol);
res=0.0;
for (int ctr1 = 1; ctr1 <= ctrrow; ctr1++) {
for (int ctr2 =1; ctr2 <= ctrcol; ctr2++) {
if( mat1(ctr1,ctr2) <= mat2(ctr1,ctr2)){
res(ctr1,ctr2) = 1.0;
}
}
}
res.Release();
return res;
}
ReturnMatrix eq(const Matrix& mat1,const Matrix& mat2)
{
int ctrcol = std::min(mat1.Ncols(),mat2.Ncols());
int ctrrow = std::min(mat1.Nrows(),mat2.Nrows());
Matrix res(ctrrow,ctrcol);
res=0.0;
for (int ctr1 = 1; ctr1 <= ctrrow; ctr1++) {
for (int ctr2 =1; ctr2 <= ctrcol; ctr2++) {
if( mat1(ctr1,ctr2) == mat2(ctr1,ctr2)){
res(ctr1,ctr2) = 1.0;
}
}
}
res.Release();
return res;
}
ReturnMatrix neq(const Matrix& mat1,const Matrix& mat2)
{
int ctrcol = std::min(mat1.Ncols(),mat2.Ncols());
int ctrrow = std::min(mat1.Nrows(),mat2.Nrows());
Matrix res(ctrrow,ctrcol);
res=0.0;
for (int ctr1 = 1; ctr1 <= ctrrow; ctr1++) {
for (int ctr2 =1; ctr2 <= ctrcol; ctr2++) {
if( mat1(ctr1,ctr2) != mat2(ctr1,ctr2)){
res(ctr1,ctr2) = 1.0;
}
}
}
res.Release();
return res;
}
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
ReturnMatrix SD(const Matrix& mat1,const Matrix& mat2)
{
if((mat1.Nrows() != mat2.Nrows()) ||
(mat1.Ncols() != mat2.Ncols()) ){
cerr <<"MISCMATHS::SD - matrices are of different dimensions"<<endl;
exit(-1);
}
Matrix ret(mat1.Nrows(),mat1.Ncols());
for (int r = 1; r <= mat1.Nrows(); r++) {
for (int c =1; c <= mat1.Ncols(); c++) {
if( mat2(r,c)==0)
ret(r,c)=0;
else
ret(r,c) = mat1(r,c)/mat2(r,c);
}
}
ret.Release();
return ret;
}
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
ReturnMatrix vox_to_vox(const ColumnVector& xyz1,const ColumnVector& dims1,const ColumnVector& dims2,const Matrix& xfm){
ColumnVector xyz1_mm(4),xyz2_mm,xyz2(3);
xyz1_mm<<xyz1(1)*dims1(1)<<xyz1(2)*dims1(2)<<xyz1(3)*dims1(3)<<1;
xyz2_mm=xfm*xyz1_mm;
xyz2_mm=xyz2_mm/xyz2_mm(4);
xyz2<<xyz2_mm(1)/dims2(1)<<xyz2_mm(2)/dims2(2)<<xyz2_mm(3)/dims2(3);
xyz2.Release();
return xyz2;
}
ReturnMatrix mni_to_imgvox(const ColumnVector& mni,const ColumnVector& mni_origin,const Matrix& mni2img, const ColumnVector& img_dims){
ColumnVector mni_new_origin(4),img_mm;//homogeneous
ColumnVector img_vox(3);
mni_new_origin<<mni(1)+mni_origin(1)<<mni(2)+mni_origin(2)<<mni(3)+mni_origin(3)<<1;
img_mm=mni2img*mni_new_origin;
img_vox<<img_mm(1)/img_dims(1)<<img_mm(2)/img_dims(2)<<img_mm(3)/img_dims(3);
img_vox.Release();
return img_vox;
}
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
ReturnMatrix remmean(const Matrix& mat, const int dim)
{
Matrix res;
if (dim == 1) {res=mat;}
else {res=mat.t();}
Matrix Mean;
Mean = mean(res);
for (int ctr = 1; ctr <= res.Nrows(); ctr++) {
for (int ctr2 =1; ctr2 <= res.Ncols(); ctr2++) {
res(ctr,ctr2)-=Mean(1,ctr2);
}
}
if (dim != 1) {res=res.t();}
res.Release();
return res;
}
void remmean(const Matrix& mat, Matrix& demeanedmat, Matrix& Mean, const int dim)
{
if (dim == 1) {demeanedmat=mat;}
else {demeanedmat=mat.t();}
Mean = mean(demeanedmat);
for (int ctr = 1; ctr <= demeanedmat.Nrows(); ctr++) {
for (int ctr2 =1; ctr2 <= demeanedmat.Ncols(); ctr2++) {
demeanedmat(ctr,ctr2)-=Mean(1,ctr2);
}
}
if (dim != 1){demeanedmat = demeanedmat.t();Mean = Mean.t();}
}
ReturnMatrix cov(const Matrix& mat, const int norm)
{
SymmetricMatrix res;
Matrix tmp;
int N;
tmp=remmean(mat);
if (norm == 1) {N = mat.Nrows();}
else {N = mat.Nrows()-1;}
res << tmp.t()*tmp;
res = res/N;
res.Release();
return res;
}
ReturnMatrix corrcoef(const Matrix& mat, const int norm)
{
SymmetricMatrix res;
SymmetricMatrix C;
C = cov(mat,norm);
Matrix D;
D=diag(C);
D=pow(sqrt(D*D.t()),-1);
res << SP(C,D);
res.Release();
return res;
}
void symm_orth(Matrix &Mat)
{
SymmetricMatrix Metric;
Metric << Mat.t()*Mat;
Metric = Metric.i();
Matrix tmpE;
DiagonalMatrix tmpD;
EigenValues(Metric,tmpD,tmpE);
Mat = Mat * tmpE * sqrt(abs(tmpD)) * tmpE.t();
}
void powerspectrum(const Matrix &Mat1, Matrix &Result, bool useLog)
//calculates the powerspectrum for every column of Mat1
{
Matrix res;
for(int ctr=1; ctr <= Mat1.Ncols(); ctr++)
{
ColumnVector tmpCol;
tmpCol=Mat1.Column(ctr);
ColumnVector FtmpCol_real;
ColumnVector FtmpCol_imag;
ColumnVector tmpPow;
RealFFT(tmpCol,FtmpCol_real,FtmpCol_imag);
tmpPow = pow(FtmpCol_real,2)+pow(FtmpCol_imag,2);
tmpPow = tmpPow.Rows(2,tmpPow.Nrows());
if(useLog){tmpPow = log(tmpPow);}
if(res.Storage()==0){res= tmpPow;}else{res|=tmpPow;}
}
Result=res;
}
void element_mod_n(Matrix& Mat,double n)
{
//represent each element in modulo n (useful for wrapping phases (n=2*M_PI))
double tmp;
for( int j=1;j<=Mat.Ncols();j++){
tmp = ( Mat(i,j) - rounddouble(Mat(i,j)/n)*n );
Mat(i,j)= tmp > 0 ? tmp : tmp + n;
}
}
}
}
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
int nextpow2(int n)
{
return (int)pow(2,ceil(log(float(n))/log(float(2))));
}
void xcorr(const Matrix& p_ts, Matrix& ret, int lag, int p_zeropad)
{
Tracer tr("MISCMATHS::xcorr");
int sizeTS = p_ts.Nrows();
int numTS = p_ts.Ncols();
if(p_zeropad == 0)
p_zeropad = sizeTS;
if(lag == 0)
lag = sizeTS;
ColumnVector x(p_zeropad);
x = 0;
ColumnVector fft_real;
ColumnVector fft_im;
ColumnVector dummy(p_zeropad);
ColumnVector dummy2;
dummy = 0;
ColumnVector realifft(p_zeropad);
ret.ReSize(lag,numTS);
ret = 0;
for(int i = 1; i <= numTS; i++)
{
x.Rows(1,sizeTS) = p_ts.Column(i);
FFT(x, dummy, fft_real, fft_im);
for(int j = 1; j <= p_zeropad; j++)
{
// (x+iy)(x-iy) = x^2 + y^2
fft_real(j) = fft_real(j)*fft_real(j) + fft_im(j)*fft_im(j);
fft_im(j) = 0;
}
FFTI(fft_real, fft_im, realifft, dummy2);
float varx = var(x.Rows(1,sizeTS)).AsScalar();
ret.Column(i) = realifft.Rows(1,lag);
for(int j = 1; j <= lag-1; j++)
{
// Correction to make autocorr unbiased and normalised
ret(j,i) = ret(j,i)/((sizeTS-j)*varx);
}
}
}
ReturnMatrix xcorr(const Matrix& p_ts, int lag, int p_zeropad )
{
Matrix r;
xcorr(p_ts,r,lag,p_zeropad);
r.Release();
return r;
}
void detrend(Matrix& p_ts, int p_level)
{
Tracer trace("MISCMATHS::detrend");
int sizeTS = p_ts.Nrows();
// p_ts = b*a + e (OLS regression)
// e is detrended data
Matrix a(sizeTS, p_level+1);
// Create a
for(int t = 1; t <= sizeTS; t++)
{
for(int l = 0; l <= p_level; l++)
a(t,l+1) = pow((float)t/sizeTS,l);
}
// Form residual forming matrix R:
Matrix R = Identity(sizeTS)-a*pinv(a);
for(int t = 1; t <= sizeTS; t++)
{
p_ts.Column(t) = R*p_ts.Column(t);
}
}
ReturnMatrix read_vest(string p_fname)
{
ifstream in;
in.open(p_fname.c_str(), ios::in);
if(!in)
{
//cerr << "Unable to open " << p_fname << endl;
throw Exception("Unable to open vest file");
}
int numWaves = 0;
int numPoints = 0;
string str;
while(true)
{
if(!in.good())
{
cerr << p_fname << "is not a valid vest file" << endl;
throw Exception("Not a valid vest file");
}
in >> str;
if(str == "/Matrix")
break;
else if(str == "/NumWaves")
{
in >> numWaves;
}
else if(str == "/NumPoints" || str == "/NumContrasts")
{
in >> numPoints;
}
}
Matrix p_mat(numPoints, numWaves);
for(int i = 1; i <= numPoints; i++)
{
for(int j = 1; j <= numWaves; j++)
{
in >> p_mat(i,j);
}
}
in.close();
p_mat.Release();
return p_mat;
}
void ols(const Matrix& data,const Matrix& des,const Matrix& tc, Matrix& cope,Matrix& varcope){
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
// ols
// data is t x v
// des is t x ev (design matrix)
// tc is cons x ev (contrast matrix)
// cope and varcope will be cons x v
// but will be resized if they are wrong
// hence may be passed in uninitialised
// TB 2004
if(data.Nrows() != des.Nrows()){
cerr <<"MISCMATHS::ols - data and design have different number of time points"<<endl;
exit(-1);
}
if(des.Ncols() != tc.Ncols()){
cerr <<"MISCMATHS::ols - design and contrast matrix have different number of EVs"<<endl;
exit(-1);
}
Matrix pdes = pinv(des);
Matrix prevar=diag(tc*pdes*pdes.t()*tc.t());
Matrix R=Identity(des.Nrows())-des*pdes;
float tR=R.Trace();
Matrix pe=pdes*data;
cope=tc*pe;
Matrix res=data-des*pe;
Matrix sigsq=sum(SP(res,res))/tR;
varcope=prevar*sigsq;
}
int ols_dof(const Matrix& des){
Matrix pdes = pinv(des);
Matrix R=Identity(des.Nrows())-des*pdes;
return int(R.Trace());
}
int conjgrad(ColumnVector& x, const Matrix& A, const ColumnVector& b, int maxit,
float reltol)
{
// solves: A * x = b (for x)
// implementation of algorithm in Golub and Van Loan (3rd ed, page 527)
ColumnVector rk1, rk2, pk, apk;
double betak, alphak, rk1rk1=0, rk2rk2, r00=0;
int k=0;
rk1 = b - A*x; // a *big* calculation
for (int n=1; n<=maxit; n++) {
k++;
if (k==1) {
pk = rk1;
rk1rk1 = (rk1.t() * rk1).AsScalar();
} else {
rk2rk2 = rk1rk1; // from before
rk1rk1 = (rk1.t() * rk1).AsScalar();
if (rk2rk2<1e-10*rk1rk1) {
cerr << "WARNING:: Conj Grad - low demoninator (rk2rk2)" << endl;
if (rk2rk2<=0) {
cerr << "Aborting conj grad ..." << endl;
return 1;
}
}
betak = rk1rk1 / rk2rk2;
pk = rk1 + betak * pk; // note RHS pk is p(k-1) in algorithm
}
// stop if sufficient accuracy is achieved
if (rk1rk1<reltol*reltol*r00) return 0;
apk = A * pk; // the *big* calculation in this algorithm
ColumnVector pap = pk.t() * apk;
if (pap.AsScalar()<0) {
cerr << "ERROR:: Conj Grad - negative eigenvector found (matrix must be symmetric and positive-definite)\nAborting ... " << endl;
return 2;
} else if (pap.AsScalar()<1e-10) {
cerr << "WARNING:: Conj Grad - nearly null eigenvector found (terminating early)" << endl;
return 1;
} else {
alphak = rk1rk1 / pap.AsScalar();
}
x = x + alphak * pk; // note LHS is x(k) and RHS is x(k-1) in algorithm
rk2 = rk1; // update prior to the next step
rk1 = rk1 - alphak * apk; // note LHS is r(k) in algorithm
}
return 0;
}
int conjgrad(ColumnVector& x, const Matrix& A, const ColumnVector& b, int maxit)
{
return conjgrad(x,A,b,maxit,1e-10);
}
float csevl(const float x, const ColumnVector& cs, const int n)
{
float b0 = 0;
float b1 = 0;
float b2 = 0;
const float twox=2*x;
for(int i=1; i<=n; i++)
{
b2=b1;
b1=b0;
b0=twox*b1-b2+cs(n+1-i);
}
return 0.5*(b0-b2);
}
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
float digamma(const float x)
{
int ntapsi(16);
int ntpsi(23);
ColumnVector psics(ntpsi);
ColumnVector apsics(ntapsi);
psics << -.038057080835217922E0<<
.49141539302938713E0<<
-.056815747821244730E0<<
.008357821225914313E0<<
-.001333232857994342E0<<
.000220313287069308E0<<
-.000037040238178456E0<<
.000006283793654854E0<<
-.000001071263908506E0<<
.000000183128394654E0<<
-.000000031353509361E0<<
.000000005372808776E0<<
-.000000000921168141E0<<
.000000000157981265E0<<
-.000000000027098646E0<<
.000000000004648722E0<<
-.000000000000797527E0<<
.000000000000136827E0<<
-.000000000000023475E0<<
.000000000000004027E0<<
-.000000000000000691E0<<
.000000000000000118E0<<
-.000000000000000020E0;
apsics <<-.0204749044678185E0<<
-.0101801271534859E0<<
.0000559718725387E0<<
-.0000012917176570E0<<
.0000000572858606E0<<
-.0000000038213539E0<<
.0000000003397434E0<<
-.0000000000374838E0<<
.0000000000048990E0<<
-.0000000000007344E0<<
.0000000000001233E0<<
-.0000000000000228E0<<
.0000000000000045E0<<
-.0000000000000009E0<<
.0000000000000002E0<<
-.0000000000000000E0;
float y = fabs(x);
float psi;
if(y<2.0)
{
// do we need to deal with the following case?
// c psi(x) for -2. .lt. x .lt. 2.
int n = int(floor(x));
y = x - n;
n = n - 1;
psi = csevl(2*y-1, psics, ntpsi);
if(n==-1)
{
psi = psi - 1.0/x;
}
}
else
{
const float aux = csevl(8/(Sqr(y))-1, apsics, ntapsi);
psi = log(x) - 0.5/x + aux;
}
return psi;
}
void glm_vb(const Matrix& X, const ColumnVector& Y, ColumnVector& B, SymmetricMatrix& ilambda_B, int niters)
// Does Variational Bayes inference on GLM Y=XB+e with ARD priors on B
// design matrix X should be num_tpts*num_evs
/////////////////////
// setup
OUT("Setup");
int ntpts=Y.Nrows();
throw Exception("COCK");
OUT(nevs);
OUT(ntpts);
ColumnVector gam_m(nevs);
gam_m=1e10;
float gam_y;
ColumnVector lambdaB(nevs);
if(nevs<ntpts-10)
// initialise with OLS
B=pinv(X)*Y;
ColumnVector res=Y-X*B;
gam_y=(ntpts-nevs)/(res.t()*res).AsScalar();
ilambda_B << (X.t()*X*gam_y).i();
lambdaB=0;
for(int l=1; l <= nevs; l++)
{
lambdaB(l)=ilambda_B(l,l);
}
else
{
OUT("no ols");
B.ReSize(nevs);
B=0;
lambdaB=1;
// ColumnVector res=Y-X*B;
// gam_y=ntpts/(res.t()*res).AsScalar();
gam_y=10;
}
// OUT(B(1));
// OUT(lambdaB(1));
float trace_ilambdaZZ=1;
SymmetricMatrix ZZ;
float YY=0;
for(int t=1; t <= ntpts; t++)
YY += Sqr(Y(t));
/////////////////////
// iterate
OUT("Iterate");
int i = 1;;
for(; i<=niters; i++)
{
////////////////////
// update phim
for(int l=1; l <= nevs; l++)
{
float b_m0 = 1e10;
float b_m = 1.0/(0.5*(Sqr(B(l))+lambdaB(l))+1.0/b_m0);
////////////////////
// update B
ColumnVector beta(nevs);
beta = 0;
SymmetricMatrix lambda_B(nevs);
lambda_B = 0;
for(int l=1; l <= nevs; l++)
lambda_B(l,l)=gam_m(l);
SymmetricMatrix tmp = lambda_B + gam_y*ZZ;
lambda_B << tmp;
beta += gam_y*ZY;
ilambda_B << lambda_B.i();
lambdaB.ReSize(nevs);
lambdaB=0;
for(int l=1; l <= nevs; l++)
{
lambdaB(l)=ilambda_B(l,l);
}
////////////////////
// compute trace for noise precision phiy update
SymmetricMatrix tmp3;
tmp3 << ilambda_B;
SymmetricMatrix tmp2;
tmp2 << tmp3*ZZ;
trace_ilambdaZZ=tmp2.Trace();
/////////////////////
// update phiy
float b_y0 = 1e10;
float c_y0 = 1;
float sum = YY + (B.t()*ZZ*B).AsScalar() - 2*(B.t()*ZY).AsScalar();
float b_y = 1.0/(0.5*(sum + trace_ilambdaZZ)+1/b_y0);
gam_y = b_y*c_y;
// OUT(gam_y);
}
vector<float> ColumnVector2vector(const ColumnVector& col)
{
vector<float> vec(col.Nrows());
for(int c = 0; c < col.Nrows(); c++)
vec[c] = col(c+1);
return vec;
}
Mark Jenkinson
committed
/////////////////////////////////////////////////////////////////////////////////////////////////////
// Uninteresting byte swapping functions
Mark Jenkinson
committed
typedef struct { unsigned char a,b ; } TWObytes ;
Mark Jenkinson
committed
void Swap_2bytes( int n , void *ar ) /* 2 bytes at a time */
{
register TWObytes *tb = (TWObytes *)ar ;
Mark Jenkinson
committed
register unsigned char tt ;
for( ii=0 ; ii < n ; ii++ ){
tt = tb[ii].a ; tb[ii].a = tb[ii].b ; tb[ii].b = tt ;
}
return ;
}
/*---------------------------------------------------------------------------*/
typedef struct { unsigned char a,b,c,d ; } FOURbytes ;
Mark Jenkinson
committed
void Swap_4bytes( int n , void *ar ) /* 4 bytes at a time */
{
register int ii ;
register FOURbytes *tb = (FOURbytes *)ar ;
Mark Jenkinson
committed
register unsigned char tt ;
for( ii=0 ; ii < n ; ii++ ){
Mark Jenkinson
committed
tt = tb[ii].a ; tb[ii].a = tb[ii].d ; tb[ii].d = tt ;
tt = tb[ii].b ; tb[ii].b = tb[ii].c ; tb[ii].c = tt ;
}
return ;
}
/*---------------------------------------------------------------------------*/
typedef struct { unsigned char a,b,c,d , D,C,B,A ; } EIGHTbytes ;
Mark Jenkinson
committed
void Swap_8bytes( int n , void *ar ) /* 8 bytes at a time */
{
register int ii ;
register EIGHTbytes *tb = (EIGHTbytes *)ar ;
Mark Jenkinson
committed
register unsigned char tt ;
for( ii=0 ; ii < n ; ii++ ){
tt = tb[ii].a ; tb[ii].a = tb[ii].A ; tb[ii].A = tt ;
tt = tb[ii].b ; tb[ii].b = tb[ii].B ; tb[ii].B = tt ;
tt = tb[ii].c ; tb[ii].c = tb[ii].C ; tb[ii].C = tt ;
tt = tb[ii].d ; tb[ii].d = tb[ii].D ; tb[ii].D = tt ;
}
return ;
}
/*---------------------------------------------------------------------------*/
typedef struct { unsigned char a,b,c,d,e,f,g,h ,
H,G,F,E,D,C,B,A ; } SIXTEENbytes ;
Mark Jenkinson
committed
void Swap_16bytes( int n , void *ar ) /* 16 bytes at a time */
{
register int ii ;
register SIXTEENbytes *tb = (SIXTEENbytes *)ar ;
Mark Jenkinson
committed
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
register unsigned char tt ;
for( ii=0 ; ii < n ; ii++ ){
tt = tb[ii].a ; tb[ii].a = tb[ii].A ; tb[ii].A = tt ;
tt = tb[ii].b ; tb[ii].b = tb[ii].B ; tb[ii].B = tt ;
tt = tb[ii].c ; tb[ii].c = tb[ii].C ; tb[ii].C = tt ;
tt = tb[ii].d ; tb[ii].d = tb[ii].D ; tb[ii].D = tt ;
tt = tb[ii].e ; tb[ii].e = tb[ii].E ; tb[ii].E = tt ;
tt = tb[ii].f ; tb[ii].f = tb[ii].F ; tb[ii].F = tt ;
tt = tb[ii].g ; tb[ii].g = tb[ii].G ; tb[ii].G = tt ;
tt = tb[ii].h ; tb[ii].h = tb[ii].H ; tb[ii].H = tt ;
}
return ;
}
/*---------------------------------------------------------------------------*/
void Swap_Nbytes( int n , int siz , void *ar ) /* subsuming case */
{
switch( siz ){
case 2: Swap_2bytes ( n , ar ) ; break ;
case 4: Swap_4bytes ( n , ar ) ; break ;
case 8: Swap_8bytes ( n , ar ) ; break ;
case 16: Swap_16bytes( n , ar ) ; break ;
}
return ;
}
// end namespace MISCMATHS
}