00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205 #include "TMatrixTBase.h"
00206 #include "TVectorT.h"
00207 #include "TROOT.h"
00208 #include "TClass.h"
00209 #include "TMath.h"
00210 #include <limits.h>
00211
00212 Int_t gMatrixCheck = 1;
00213
00214 #ifndef R__ALPHA
00215 templateClassImp(TMatrixTBase)
00216 #endif
00217
00218
00219 template<class Element>
00220 void TMatrixTBase<Element>::DoubleLexSort(Int_t n,Int_t *first,Int_t *second,Element *data)
00221 {
00222
00223
00224 const int incs[] = {1,5,19,41,109,209,505,929,2161,3905,8929,16001,INT_MAX};
00225
00226 Int_t kinc = 0;
00227 while (incs[kinc] <= n/2)
00228 kinc++;
00229 kinc -= 1;
00230
00231
00232
00233
00234 for( ; kinc >= 0; kinc--) {
00235 const Int_t inc = incs[kinc];
00236
00237 for (Int_t k = inc; k < n; k++) {
00238 const Element tmp = data[k];
00239 const Int_t fi = first [k];
00240 const Int_t se = second[k];
00241 Int_t j;
00242 for (j = k; j >= inc; j -= inc) {
00243 if ( fi < first[j-inc] || (fi == first[j-inc] && se < second[j-inc]) ) {
00244 data [j] = data [j-inc];
00245 first [j] = first [j-inc];
00246 second[j] = second[j-inc];
00247 } else
00248 break;
00249 }
00250 data [j] = tmp;
00251 first [j] = fi;
00252 second[j] = se;
00253 }
00254 }
00255 }
00256
00257
00258 template<class Element>
00259 void TMatrixTBase<Element>::IndexedLexSort(Int_t n,Int_t *first,Int_t swapFirst,
00260 Int_t *second,Int_t swapSecond,Int_t *index)
00261 {
00262
00263
00264 const int incs[] = {1,5,19,41,109,209,505,929,2161,3905,8929,16001,INT_MAX};
00265
00266 Int_t kinc = 0;
00267 while (incs[kinc] <= n/2)
00268 kinc++;
00269 kinc -= 1;
00270
00271
00272
00273
00274 for( ; kinc >= 0; kinc--) {
00275 const Int_t inc = incs[kinc];
00276
00277 if ( !swapFirst && !swapSecond ) {
00278 for (Int_t k = inc; k < n; k++) {
00279
00280 const Int_t ktemp = index[k];
00281 const Int_t fi = first [ktemp];
00282 const Int_t se = second[ktemp];
00283
00284 Int_t j;
00285 for (j = k; j >= inc; j -= inc) {
00286
00287 if (fi < first[index[j-inc]] || (fi == first[index[j-inc]] && se < second[index[j-inc]])) {
00288
00289
00290
00291 index[j] = index[j-inc];
00292 } else {
00293
00294
00295 break;
00296 }
00297 }
00298
00299 index[j] = ktemp;
00300
00301 }
00302 } else if ( swapSecond && !swapFirst ) {
00303 for (Int_t k = inc; k < n; k++) {
00304 const Int_t ktemp = index[k];
00305 const Int_t fi = first [ktemp];
00306 const Int_t se = second[k];
00307 Int_t j;
00308 for (j = k; j >= inc; j -= inc) {
00309 if (fi < first[index[j-inc]] || (fi == first[index[j-inc]] && se < second[j-inc])) {
00310 index [j] = index[j-inc];
00311 second[j] = second[j-inc];
00312 } else {
00313 break;
00314 }
00315 }
00316 index[j] = ktemp;
00317 second[j] = se;
00318 }
00319 } else if (swapFirst && !swapSecond) {
00320 for (Int_t k = inc; k < n; k++ ) {
00321 const Int_t ktemp = index[k];
00322 const Int_t fi = first[k];
00323 const Int_t se = second[ktemp];
00324 Int_t j;
00325 for (j = k; j >= inc; j -= inc) {
00326 if ( fi < first[j-inc] || (fi == first[j-inc] && se < second[ index[j-inc]])) {
00327 index[j] = index[j-inc];
00328 first[j] = first[j-inc];
00329 } else {
00330 break;
00331 }
00332 }
00333 index[j] = ktemp;
00334 first[j] = fi;
00335 }
00336 } else {
00337 for (Int_t k = inc; k < n; k++ ) {
00338 const Int_t ktemp = index[k];
00339 const Int_t fi = first [k];
00340 const Int_t se = second[k];
00341 Int_t j;
00342 for (j = k; j >= inc; j -= inc) {
00343 if ( fi < first[j-inc] || (fi == first[j-inc] && se < second[j-inc])) {
00344 index [j] = index [j-inc];
00345 first [j] = first [j-inc];
00346 second[j] = second[j-inc];
00347 } else {
00348 break;
00349 }
00350 }
00351 index[j] = ktemp;
00352 first[j] = fi;
00353 second[j] = se;
00354 }
00355 }
00356 }
00357 }
00358
00359
00360 template<class Element>
00361 TMatrixTBase<Element> &TMatrixTBase<Element>::SetMatrixArray(const Element *data,Option_t *option)
00362 {
00363
00364
00365
00366
00367
00368
00369
00370 R__ASSERT(IsValid());
00371
00372 TString opt = option;
00373 opt.ToUpper();
00374
00375 Element *elem = GetMatrixArray();
00376 if (opt.Contains("F")) {
00377 for (Int_t irow = 0; irow < fNrows; irow++) {
00378 const Int_t off1 = irow*fNcols;
00379 Int_t off2 = 0;
00380 for (Int_t icol = 0; icol < fNcols; icol++) {
00381 elem[off1+icol] = data[off2+irow];
00382 off2 += fNrows;
00383 }
00384 }
00385 }
00386 else
00387 memcpy(elem,data,fNelems*sizeof(Element));
00388
00389 return *this;
00390 }
00391
00392
00393 template<class Element>
00394 Bool_t TMatrixTBase<Element>::IsSymmetric() const
00395 {
00396
00397
00398 R__ASSERT(IsValid());
00399
00400 if ((fNrows != fNcols) || (fRowLwb != fColLwb))
00401 return kFALSE;
00402
00403 const Element * const elem = GetMatrixArray();
00404 for (Int_t irow = 0; irow < fNrows; irow++) {
00405 const Int_t rowOff = irow*fNcols;
00406 Int_t colOff = 0;
00407 for (Int_t icol = 0; icol < irow; icol++) {
00408 if (elem[rowOff+icol] != elem[colOff+irow])
00409 return kFALSE;
00410 colOff += fNrows;
00411 }
00412 }
00413 return kTRUE;
00414 }
00415
00416
00417 template<class Element>
00418 void TMatrixTBase<Element>::GetMatrix2Array(Element *data,Option_t *option) const
00419 {
00420
00421
00422
00423
00424
00425
00426
00427 R__ASSERT(IsValid());
00428
00429 TString opt = option;
00430 opt.ToUpper();
00431
00432 const Element * const elem = GetMatrixArray();
00433 if (opt.Contains("F")) {
00434 for (Int_t irow = 0; irow < fNrows; irow++) {
00435 const Int_t off1 = irow*fNcols;
00436 Int_t off2 = 0;
00437 for (Int_t icol = 0; icol < fNcols; icol++)
00438 data[off2+irow] = elem[off1+icol];
00439 off2 += fNrows;
00440 }
00441 }
00442 else
00443 memcpy(data,elem,fNelems*sizeof(Element));
00444 }
00445
00446
00447 template<class Element>
00448 TMatrixTBase<Element> &TMatrixTBase<Element>::InsertRow(Int_t rown,Int_t coln,const Element *v,Int_t n)
00449 {
00450
00451
00452 const Int_t arown = rown-fRowLwb;
00453 const Int_t acoln = coln-fColLwb;
00454 const Int_t nr = (n > 0) ? n : fNcols;
00455
00456 if (gMatrixCheck) {
00457 if (arown >= fNrows || arown < 0) {
00458 Error("InsertRow","row %d out of matrix range",rown);
00459 return *this;
00460 }
00461
00462 if (acoln >= fNcols || acoln < 0) {
00463 Error("InsertRow","column %d out of matrix range",coln);
00464 return *this;
00465 }
00466
00467 if (acoln+nr > fNcols || nr < 0) {
00468 Error("InsertRow","row length %d out of range",nr);
00469 return *this;
00470 }
00471 }
00472
00473 const Int_t off = arown*fNcols+acoln;
00474 Element * const elem = GetMatrixArray()+off;
00475 memcpy(elem,v,nr*sizeof(Element));
00476
00477 return *this;
00478 }
00479
00480
00481 template<class Element>
00482 void TMatrixTBase<Element>::ExtractRow(Int_t rown,Int_t coln,Element *v,Int_t n) const
00483 {
00484
00485
00486 const Int_t arown = rown-fRowLwb;
00487 const Int_t acoln = coln-fColLwb;
00488 const Int_t nr = (n > 0) ? n : fNcols;
00489
00490 if (gMatrixCheck) {
00491 if (arown >= fNrows || arown < 0) {
00492 Error("ExtractRow","row %d out of matrix range",rown);
00493 return;
00494 }
00495
00496 if (acoln >= fNcols || acoln < 0) {
00497 Error("ExtractRow","column %d out of matrix range",coln);
00498 return;
00499 }
00500
00501 if (acoln+n >= fNcols || nr < 0) {
00502 Error("ExtractRow","row length %d out of range",nr);
00503 return;
00504 }
00505 }
00506
00507 const Int_t off = arown*fNcols+acoln;
00508 const Element * const elem = GetMatrixArray()+off;
00509 memcpy(v,elem,nr*sizeof(Element));
00510 }
00511
00512
00513 template<class Element>
00514 TMatrixTBase<Element> &TMatrixTBase<Element>::Shift(Int_t row_shift,Int_t col_shift)
00515 {
00516
00517
00518
00519
00520 fRowLwb += row_shift;
00521 fColLwb += col_shift;
00522
00523 return *this;
00524 }
00525
00526
00527 template<class Element>
00528 TMatrixTBase<Element> &TMatrixTBase<Element>::Zero()
00529 {
00530
00531
00532 R__ASSERT(IsValid());
00533 memset(this->GetMatrixArray(),0,fNelems*sizeof(Element));
00534
00535 return *this;
00536 }
00537
00538
00539 template<class Element>
00540 TMatrixTBase<Element> &TMatrixTBase<Element>::Abs()
00541 {
00542
00543
00544 R__ASSERT(IsValid());
00545
00546 Element *ep = this->GetMatrixArray();
00547 const Element * const fp = ep+fNelems;
00548 while (ep < fp) {
00549 *ep = TMath::Abs(*ep);
00550 ep++;
00551 }
00552
00553 return *this;
00554 }
00555
00556
00557 template<class Element>
00558 TMatrixTBase<Element> &TMatrixTBase<Element>::Sqr()
00559 {
00560
00561
00562 R__ASSERT(IsValid());
00563
00564 Element *ep = this->GetMatrixArray();
00565 const Element * const fp = ep+fNelems;
00566 while (ep < fp) {
00567 *ep = (*ep) * (*ep);
00568 ep++;
00569 }
00570
00571 return *this;
00572 }
00573
00574
00575 template<class Element>
00576 TMatrixTBase<Element> &TMatrixTBase<Element>::Sqrt()
00577 {
00578
00579
00580 R__ASSERT(IsValid());
00581
00582 Element *ep = this->GetMatrixArray();
00583 const Element * const fp = ep+fNelems;
00584 while (ep < fp) {
00585 *ep = TMath::Sqrt(*ep);
00586 ep++;
00587 }
00588
00589 return *this;
00590 }
00591
00592
00593 template<class Element>
00594 TMatrixTBase<Element> &TMatrixTBase<Element>::UnitMatrix()
00595 {
00596
00597
00598 R__ASSERT(IsValid());
00599
00600 Element *ep = this->GetMatrixArray();
00601 memset(ep,0,fNelems*sizeof(Element));
00602 for (Int_t i = fRowLwb; i <= fRowLwb+fNrows-1; i++)
00603 for (Int_t j = fColLwb; j <= fColLwb+fNcols-1; j++)
00604 *ep++ = (i==j ? 1.0 : 0.0);
00605
00606 return *this;
00607 }
00608
00609
00610 template<class Element>
00611 TMatrixTBase<Element> &TMatrixTBase<Element>::NormByDiag(const TVectorT<Element> &v,Option_t *option)
00612 {
00613
00614
00615
00616
00617 R__ASSERT(IsValid());
00618 R__ASSERT(v.IsValid());
00619
00620 if (gMatrixCheck) {
00621 const Int_t nMax = TMath::Max(fNrows,fNcols);
00622 if (v.GetNoElements() < nMax) {
00623 Error("NormByDiag","vector shorter than matrix diagonal");
00624 return *this;
00625 }
00626 }
00627
00628 TString opt(option);
00629 opt.ToUpper();
00630 const Int_t divide = (opt.Contains("D")) ? 1 : 0;
00631
00632 const Element *pV = v.GetMatrixArray();
00633 Element *mp = this->GetMatrixArray();
00634
00635 if (divide) {
00636 for (Int_t irow = 0; irow < fNrows; irow++) {
00637 if (pV[irow] != 0.0) {
00638 for (Int_t icol = 0; icol < fNcols; icol++) {
00639 if (pV[icol] != 0.0) {
00640 const Element val = TMath::Sqrt(TMath::Abs(pV[irow]*pV[icol]));
00641 *mp++ /= val;
00642 } else {
00643 Error("NormbyDiag","vector element %d is zero",icol);
00644 mp++;
00645 }
00646 }
00647 } else {
00648 Error("NormbyDiag","vector element %d is zero",irow);
00649 mp += fNcols;
00650 }
00651 }
00652 } else {
00653 for (Int_t irow = 0; irow < fNrows; irow++) {
00654 for (Int_t icol = 0; icol < fNcols; icol++) {
00655 const Element val = TMath::Sqrt(TMath::Abs(pV[irow]*pV[icol]));
00656 *mp++ *= val;
00657 }
00658 }
00659 }
00660
00661 return *this;
00662 }
00663
00664
00665 template<class Element>
00666 Element TMatrixTBase<Element>::RowNorm() const
00667 {
00668
00669
00670
00671 R__ASSERT(IsValid());
00672
00673 const Element * ep = GetMatrixArray();
00674 const Element * const fp = ep+fNelems;
00675 Element norm = 0;
00676
00677
00678 while (ep < fp) {
00679 Element sum = 0;
00680
00681 for (Int_t j = 0; j < fNcols; j++)
00682 sum += TMath::Abs(*ep++);
00683 norm = TMath::Max(norm,sum);
00684 }
00685
00686 R__ASSERT(ep == fp);
00687
00688 return norm;
00689 }
00690
00691
00692 template<class Element>
00693 Element TMatrixTBase<Element>::ColNorm() const
00694 {
00695
00696
00697
00698 R__ASSERT(IsValid());
00699
00700 const Element * ep = GetMatrixArray();
00701 const Element * const fp = ep+fNcols;
00702 Element norm = 0;
00703
00704
00705 while (ep < fp) {
00706 Element sum = 0;
00707
00708 for (Int_t i = 0; i < fNrows; i++,ep += fNcols)
00709 sum += TMath::Abs(*ep);
00710 ep -= fNelems-1;
00711 norm = TMath::Max(norm,sum);
00712 }
00713
00714 R__ASSERT(ep == fp);
00715
00716 return norm;
00717 }
00718
00719
00720 template<class Element>
00721 Element TMatrixTBase<Element>::E2Norm() const
00722 {
00723
00724
00725 R__ASSERT(IsValid());
00726
00727 const Element * ep = GetMatrixArray();
00728 const Element * const fp = ep+fNelems;
00729 Element sum = 0;
00730
00731 for ( ; ep < fp; ep++)
00732 sum += (*ep) * (*ep);
00733
00734 return sum;
00735 }
00736
00737
00738 template<class Element>
00739 Int_t TMatrixTBase<Element>::NonZeros() const
00740 {
00741
00742
00743 R__ASSERT(IsValid());
00744
00745 Int_t nr_nonzeros = 0;
00746 const Element *ep = this->GetMatrixArray();
00747 const Element * const fp = ep+fNelems;
00748 while (ep < fp)
00749 if (*ep++ != 0.0) nr_nonzeros++;
00750
00751 return nr_nonzeros;
00752 }
00753
00754
00755 template<class Element>
00756 Element TMatrixTBase<Element>::Sum() const
00757 {
00758
00759
00760 R__ASSERT(IsValid());
00761
00762 Element sum = 0.0;
00763 const Element *ep = this->GetMatrixArray();
00764 const Element * const fp = ep+fNelems;
00765 while (ep < fp)
00766 sum += *ep++;
00767
00768 return sum;
00769 }
00770
00771
00772 template<class Element>
00773 Element TMatrixTBase<Element>::Min() const
00774 {
00775
00776
00777 R__ASSERT(IsValid());
00778
00779 const Element * const ep = this->GetMatrixArray();
00780 const Int_t index = TMath::LocMin(fNelems,ep);
00781 return ep[index];
00782 }
00783
00784
00785 template<class Element>
00786 Element TMatrixTBase<Element>::Max() const
00787 {
00788
00789
00790 R__ASSERT(IsValid());
00791
00792 const Element * const ep = this->GetMatrixArray();
00793 const Int_t index = TMath::LocMax(fNelems,ep);
00794 return ep[index];
00795 }
00796
00797
00798 template<class Element>
00799 void TMatrixTBase<Element>::Draw(Option_t *option)
00800 {
00801
00802
00803
00804 gROOT->ProcessLine(Form("THistPainter::PaintSpecialObjects((TObject*)0x%lx,\"%s\");",
00805 (ULong_t)this, option));
00806 }
00807
00808
00809 template<class Element>
00810 void TMatrixTBase<Element>::Print(Option_t *option) const
00811 {
00812
00813
00814
00815
00816
00817 if (!IsValid()) {
00818 Error("Print","Matrix is invalid");
00819 return;
00820 }
00821
00822
00823 const char *format = "%11.4g ";
00824 if (option) {
00825 const char *f = strstr(option,"f=");
00826 if (f) format = f+2;
00827 }
00828 char topbar[100];
00829 snprintf(topbar,100,format,123.456789);
00830 Int_t nch = strlen(topbar)+1;
00831 if (nch > 18) nch = 18;
00832 char ftopbar[20];
00833 for (Int_t i = 0; i < nch; i++) ftopbar[i] = ' ';
00834 Int_t nk = 1 + Int_t(TMath::Log10(fNcols));
00835 snprintf(ftopbar+nch/2,20-nch/2,"%s%dd","%",nk);
00836 Int_t nch2 = strlen(ftopbar);
00837 for (Int_t i = nch2; i < nch; i++) ftopbar[i] = ' ';
00838 ftopbar[nch] = '|';
00839 ftopbar[nch+1] = 0;
00840
00841 printf("\n%dx%d matrix is as follows",fNrows,fNcols);
00842
00843 Int_t cols_per_sheet = 5;
00844 if (nch <= 8) cols_per_sheet =10;
00845 const Int_t ncols = fNcols;
00846 const Int_t nrows = fNrows;
00847 const Int_t collwb = fColLwb;
00848 const Int_t rowlwb = fRowLwb;
00849 nk = 5+nch*TMath::Min(cols_per_sheet,fNcols);
00850 for (Int_t i = 0; i < nk; i++) topbar[i] = '-';
00851 topbar[nk] = 0;
00852 for (Int_t sheet_counter = 1; sheet_counter <= ncols; sheet_counter += cols_per_sheet) {
00853 printf("\n\n |");
00854 for (Int_t j = sheet_counter; j < sheet_counter+cols_per_sheet && j <= ncols; j++)
00855 printf(ftopbar,j+collwb-1);
00856 printf("\n%s\n",topbar);
00857 if (fNelems <= 0) continue;
00858 for (Int_t i = 1; i <= nrows; i++) {
00859 printf("%4d |",i+rowlwb-1);
00860 for (Int_t j = sheet_counter; j < sheet_counter+cols_per_sheet && j <= ncols; j++)
00861 printf(format,(*this)(i+rowlwb-1,j+collwb-1));
00862 printf("\n");
00863 }
00864 }
00865 printf("\n");
00866 }
00867
00868
00869 template<class Element>
00870 Bool_t TMatrixTBase<Element>::operator==(Element val) const
00871 {
00872
00873
00874 R__ASSERT(IsValid());
00875
00876 if (val == 0. && fNelems == 0)
00877 return kTRUE;
00878
00879 const Element * ep = GetMatrixArray();
00880 const Element * const fp = ep+fNelems;
00881 for (; ep < fp; ep++)
00882 if (!(*ep == val))
00883 return kFALSE;
00884
00885 return kTRUE;
00886 }
00887
00888
00889 template<class Element>
00890 Bool_t TMatrixTBase<Element>::operator!=(Element val) const
00891 {
00892
00893
00894 R__ASSERT(IsValid());
00895
00896 if (val == 0. && fNelems == 0)
00897 return kFALSE;
00898
00899 const Element * ep = GetMatrixArray();
00900 const Element * const fp = ep+fNelems;
00901 for (; ep < fp; ep++)
00902 if (!(*ep != val))
00903 return kFALSE;
00904
00905 return kTRUE;
00906 }
00907
00908
00909 template<class Element>
00910 Bool_t TMatrixTBase<Element>::operator<(Element val) const
00911 {
00912
00913
00914 R__ASSERT(IsValid());
00915
00916 const Element * ep = GetMatrixArray();
00917 const Element * const fp = ep+fNelems;
00918 for (; ep < fp; ep++)
00919 if (!(*ep < val))
00920 return kFALSE;
00921
00922 return kTRUE;
00923 }
00924
00925
00926 template<class Element>
00927 Bool_t TMatrixTBase<Element>::operator<=(Element val) const
00928 {
00929
00930
00931 R__ASSERT(IsValid());
00932
00933 const Element * ep = GetMatrixArray();
00934 const Element * const fp = ep+fNelems;
00935 for (; ep < fp; ep++)
00936 if (!(*ep <= val))
00937 return kFALSE;
00938
00939 return kTRUE;
00940 }
00941
00942
00943 template<class Element>
00944 Bool_t TMatrixTBase<Element>::operator>(Element val) const
00945 {
00946
00947
00948 R__ASSERT(IsValid());
00949
00950 const Element * ep = GetMatrixArray();
00951 const Element * const fp = ep+fNelems;
00952 for (; ep < fp; ep++)
00953 if (!(*ep > val))
00954 return kFALSE;
00955
00956 return kTRUE;
00957 }
00958
00959
00960 template<class Element>
00961 Bool_t TMatrixTBase<Element>::operator>=(Element val) const
00962 {
00963
00964
00965 R__ASSERT(IsValid());
00966
00967 const Element * ep = GetMatrixArray();
00968 const Element * const fp = ep+fNelems;
00969 for (; ep < fp; ep++)
00970 if (!(*ep >= val))
00971 return kFALSE;
00972
00973 return kTRUE;
00974 }
00975
00976
00977 template<class Element>
00978 TMatrixTBase<Element> &TMatrixTBase<Element>::Apply(const TElementActionT<Element> &action)
00979 {
00980
00981
00982 R__ASSERT(IsValid());
00983
00984 Element *ep = this->GetMatrixArray();
00985 const Element * const ep_last = ep+fNelems;
00986 while (ep < ep_last)
00987 action.Operation(*ep++);
00988
00989 return *this;
00990 }
00991
00992
00993 template<class Element>
00994 TMatrixTBase<Element> &TMatrixTBase<Element>::Apply(const TElementPosActionT<Element> &action)
00995 {
00996
00997
00998
00999 R__ASSERT(IsValid());
01000
01001 Element *ep = this->GetMatrixArray();
01002 for (action.fI = fRowLwb; action.fI < fRowLwb+fNrows; action.fI++)
01003 for (action.fJ = fColLwb; action.fJ < fColLwb+fNcols; action.fJ++)
01004 action.Operation(*ep++);
01005
01006 R__ASSERT(ep == this->GetMatrixArray()+fNelems);
01007
01008 return *this;
01009 }
01010
01011
01012 template<class Element>
01013 TMatrixTBase<Element> &TMatrixTBase<Element>::Randomize(Element alpha,Element beta,Double_t &seed)
01014 {
01015
01016
01017 R__ASSERT(IsValid());
01018
01019 const Element scale = beta-alpha;
01020 const Element shift = alpha/scale;
01021
01022 Element * ep = GetMatrixArray();
01023 const Element * const fp = ep+fNelems;
01024 while (ep < fp)
01025 *ep++ = scale*(Drand(seed)+shift);
01026
01027 return *this;
01028 }
01029
01030
01031 template<class Element>
01032 Bool_t operator==(const TMatrixTBase<Element> &m1,const TMatrixTBase<Element> &m2)
01033 {
01034
01035
01036 if (!AreCompatible(m1,m2)) return kFALSE;
01037 return (memcmp(m1.GetMatrixArray(),m2.GetMatrixArray(),
01038 m1.GetNoElements()*sizeof(Element)) == 0);
01039 }
01040
01041
01042 template<class Element>
01043 Element E2Norm(const TMatrixTBase<Element> &m1,const TMatrixTBase<Element> &m2)
01044 {
01045
01046
01047 if (gMatrixCheck && !AreCompatible(m1,m2)) {
01048 ::Error("E2Norm","matrices not compatible");
01049 return -1.0;
01050 }
01051
01052 const Element * mp1 = m1.GetMatrixArray();
01053 const Element * mp2 = m2.GetMatrixArray();
01054 const Element * const fmp1 = mp1+m1.GetNoElements();
01055
01056 Element sum = 0.0;
01057 for (; mp1 < fmp1; mp1++, mp2++)
01058 sum += (*mp1 - *mp2)*(*mp1 - *mp2);
01059
01060 return sum;
01061 }
01062
01063
01064 template<class Element1,class Element2>
01065 Bool_t AreCompatible(const TMatrixTBase<Element1> &m1,const TMatrixTBase<Element2> &m2,Int_t verbose)
01066 {
01067
01068 if (!m1.IsValid()) {
01069 if (verbose)
01070 ::Error("AreCompatible", "matrix 1 not valid");
01071 return kFALSE;
01072 }
01073 if (!m2.IsValid()) {
01074 if (verbose)
01075 ::Error("AreCompatible", "matrix 2 not valid");
01076 return kFALSE;
01077 }
01078
01079 if (m1.GetNrows() != m2.GetNrows() || m1.GetNcols() != m2.GetNcols() ||
01080 m1.GetRowLwb() != m2.GetRowLwb() || m1.GetColLwb() != m2.GetColLwb()) {
01081 if (verbose)
01082 ::Error("AreCompatible", "matrices 1 and 2 not compatible");
01083 return kFALSE;
01084 }
01085
01086 return kTRUE;
01087 }
01088
01089
01090 template<class Element>
01091 void Compare(const TMatrixTBase<Element> &m1,const TMatrixTBase<Element> &m2)
01092 {
01093
01094
01095 if (!AreCompatible(m1,m2)) {
01096 Error("Compare(const TMatrixTBase<Element> &,const TMatrixTBase<Element> &)","matrices are incompatible");
01097 return;
01098 }
01099
01100 printf("\n\nComparison of two TMatrices:\n");
01101
01102 Element norm1 = 0;
01103 Element norm2 = 0;
01104 Element ndiff = 0;
01105 Int_t imax = 0;
01106 Int_t jmax = 0;
01107 Element difmax = -1;
01108
01109 for (Int_t i = m1.GetRowLwb(); i <= m1.GetRowUpb(); i++) {
01110 for (Int_t j = m1.GetColLwb(); j < m1.GetColUpb(); j++) {
01111 const Element mv1 = m1(i,j);
01112 const Element mv2 = m2(i,j);
01113 const Element diff = TMath::Abs(mv1-mv2);
01114
01115 if (diff > difmax) {
01116 difmax = diff;
01117 imax = i;
01118 jmax = j;
01119 }
01120 norm1 += TMath::Abs(mv1);
01121 norm2 += TMath::Abs(mv2);
01122 ndiff += TMath::Abs(diff);
01123 }
01124 }
01125
01126 printf("\nMaximal discrepancy \t\t%g", difmax);
01127 printf("\n occured at the point\t\t(%d,%d)",imax,jmax);
01128 const Element mv1 = m1(imax,jmax);
01129 const Element mv2 = m2(imax,jmax);
01130 printf("\n Matrix 1 element is \t\t%g", mv1);
01131 printf("\n Matrix 2 element is \t\t%g", mv2);
01132 printf("\n Absolute error v2[i]-v1[i]\t\t%g", mv2-mv1);
01133 printf("\n Relative error\t\t\t\t%g\n",
01134 (mv2-mv1)/TMath::Max(TMath::Abs(mv2+mv1)/2,(Element)1e-7));
01135
01136 printf("\n||Matrix 1|| \t\t\t%g", norm1);
01137 printf("\n||Matrix 2|| \t\t\t%g", norm2);
01138 printf("\n||Matrix1-Matrix2||\t\t\t\t%g", ndiff);
01139 printf("\n||Matrix1-Matrix2||/sqrt(||Matrix1|| ||Matrix2||)\t%g\n\n",
01140 ndiff/TMath::Max(TMath::Sqrt(norm1*norm2),1e-7));
01141 }
01142
01143
01144 template<class Element>
01145 Bool_t VerifyMatrixValue(const TMatrixTBase<Element> &m,Element val,Int_t verbose,Element maxDevAllow)
01146 {
01147
01148
01149 R__ASSERT(m.IsValid());
01150
01151 if (m == 0)
01152 return kTRUE;
01153
01154 Int_t imax = 0;
01155 Int_t jmax = 0;
01156 Element maxDevObs = 0;
01157
01158 if (TMath::Abs(maxDevAllow) <= 0.0)
01159 maxDevAllow = std::numeric_limits<Element>::epsilon();
01160
01161 for (Int_t i = m.GetRowLwb(); i <= m.GetRowUpb(); i++) {
01162 for (Int_t j = m.GetColLwb(); j <= m.GetColUpb(); j++) {
01163 const Element dev = TMath::Abs(m(i,j)-val);
01164 if (dev > maxDevObs) {
01165 imax = i;
01166 jmax = j;
01167 maxDevObs = dev;
01168 }
01169 }
01170 }
01171
01172 if (maxDevObs == 0)
01173 return kTRUE;
01174
01175 if (verbose) {
01176 printf("Largest dev for (%d,%d); dev = |%g - %g| = %g\n",imax,jmax,m(imax,jmax),val,maxDevObs);
01177 if(maxDevObs > maxDevAllow)
01178 Error("VerifyElementValue","Deviation > %g\n",maxDevAllow);
01179 }
01180
01181 if(maxDevObs > maxDevAllow)
01182 return kFALSE;
01183 return kTRUE;
01184 }
01185
01186
01187 template<class Element>
01188 Bool_t VerifyMatrixIdentity(const TMatrixTBase<Element> &m1,const TMatrixTBase<Element> &m2,Int_t verbose,
01189 Element maxDevAllow)
01190 {
01191
01192
01193 if (!AreCompatible(m1,m2,verbose))
01194 return kFALSE;
01195
01196 if (m1 == 0 && m2 == 0)
01197 return kTRUE;
01198
01199 Int_t imax = 0;
01200 Int_t jmax = 0;
01201 Element maxDevObs = 0;
01202
01203 if (TMath::Abs(maxDevAllow) <= 0.0)
01204 maxDevAllow = std::numeric_limits<Element>::epsilon();
01205
01206 for (Int_t i = m1.GetRowLwb(); i <= m1.GetRowUpb(); i++) {
01207 for (Int_t j = m1.GetColLwb(); j <= m1.GetColUpb(); j++) {
01208 const Element dev = TMath::Abs(m1(i,j)-m2(i,j));
01209 if (dev > maxDevObs) {
01210 imax = i;
01211 jmax = j;
01212 maxDevObs = dev;
01213 }
01214 }
01215 }
01216
01217 if (maxDevObs == 0)
01218 return kTRUE;
01219
01220 if (verbose) {
01221 printf("Largest dev for (%d,%d); dev = |%g - %g| = %g\n",
01222 imax,jmax,m1(imax,jmax),m2(imax,jmax),maxDevObs);
01223 if (maxDevObs > maxDevAllow)
01224 Error("VerifyMatrixValue","Deviation > %g\n",maxDevAllow);
01225 }
01226
01227 if (maxDevObs > maxDevAllow)
01228 return kFALSE;
01229 return kTRUE;
01230 }
01231
01232
01233 template<class Element>
01234 void TMatrixTBase<Element>::Streamer(TBuffer &R__b)
01235 {
01236
01237
01238 if (R__b.IsReading()) {
01239 UInt_t R__s, R__c;
01240 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
01241 if (R__v > 1) {
01242 R__b.ReadClassBuffer(TMatrixTBase<Element>::Class(),this,R__v,R__s,R__c);
01243 } else {
01244 Error("TMatrixTBase<Element>::Streamer","Unknown version number: %d",R__v);
01245 R__ASSERT(0);
01246 }
01247 if (R__v < 4) MakeValid();
01248 } else {
01249 R__b.WriteClassBuffer(TMatrixTBase<Element>::Class(),this);
01250 }
01251 }
01252
01253 template class TMatrixTBase<Float_t>;
01254
01255 template Bool_t operator== <Float_t>(const TMatrixFBase &m1,const TMatrixFBase &m2);
01256 template Float_t E2Norm <Float_t>(const TMatrixFBase &m1,const TMatrixFBase &m2);
01257 template Bool_t AreCompatible<Float_t,Float_t>
01258 (const TMatrixFBase &m1,const TMatrixFBase &m2,Int_t verbose);
01259 template Bool_t AreCompatible<Float_t,Double_t>
01260 (const TMatrixFBase &m1,const TMatrixDBase &m2,Int_t verbose);
01261 template void Compare <Float_t>(const TMatrixFBase &m1,const TMatrixFBase &m2);
01262 template Bool_t VerifyMatrixValue <Float_t>(const TMatrixFBase &m,Float_t val,Int_t verbose,Float_t maxDevAllow);
01263 template Bool_t VerifyMatrixValue <Float_t>(const TMatrixFBase &m,Float_t val);
01264 template Bool_t VerifyMatrixIdentity<Float_t>(const TMatrixFBase &m1,const TMatrixFBase &m2,
01265 Int_t verbose,Float_t maxDevAllowN);
01266 template Bool_t VerifyMatrixIdentity<Float_t>(const TMatrixFBase &m1,const TMatrixFBase &m2);
01267
01268 template class TMatrixTBase<Double_t>;
01269
01270 template Bool_t operator== <Double_t>(const TMatrixDBase &m1,const TMatrixDBase &m2);
01271 template Double_t E2Norm <Double_t>(const TMatrixDBase &m1,const TMatrixDBase &m2);
01272 template Bool_t AreCompatible<Double_t,Double_t>
01273 (const TMatrixDBase &m1,const TMatrixDBase &m2,Int_t verbose);
01274 template Bool_t AreCompatible<Double_t,Float_t>
01275 (const TMatrixDBase &m1,const TMatrixFBase &m2,Int_t verbose);
01276 template void Compare <Double_t>(const TMatrixDBase &m1,const TMatrixDBase &m2);
01277 template Bool_t VerifyMatrixValue <Double_t>(const TMatrixDBase &m,Double_t val,Int_t verbose,Double_t maxDevAllow);
01278 template Bool_t VerifyMatrixValue <Double_t>(const TMatrixDBase &m,Double_t val);
01279 template Bool_t VerifyMatrixIdentity<Double_t>(const TMatrixDBase &m1,const TMatrixDBase &m2,
01280 Int_t verbose,Double_t maxDevAllow);
01281 template Bool_t VerifyMatrixIdentity<Double_t>(const TMatrixDBase &m1,const TMatrixDBase &m2);