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 #include "TVectorT.h"
00038 #include "TClass.h"
00039 #include "TMath.h"
00040 #include "TROOT.h"
00041 #include "Varargs.h"
00042
00043 #ifndef R__ALPHA
00044 templateClassImp(TVectorT)
00045 #endif
00046
00047
00048 template<class Element>
00049 void TVectorT<Element>::Delete_m(Int_t size,Element *&m)
00050 {
00051
00052
00053 if (m) {
00054 if (size > kSizeMax)
00055 delete [] m;
00056 m = 0;
00057 }
00058 }
00059
00060
00061 template<class Element>
00062 Element* TVectorT<Element>::New_m(Int_t size)
00063 {
00064
00065
00066
00067 if (size == 0) return 0;
00068 else {
00069 if ( size <= kSizeMax )
00070 return fDataStack;
00071 else {
00072 Element *heap = new Element[size];
00073 return heap;
00074 }
00075 }
00076 }
00077
00078
00079 template<class Element>
00080 void TVectorT<Element>::Add(const TVectorT<Element> &v)
00081 {
00082
00083
00084 if (gMatrixCheck && !AreCompatible(*this,v)) {
00085 Error("Add(TVectorT<Element> &)","vector's not compatible");
00086 return;
00087 }
00088
00089 const Element *sp = v.GetMatrixArray();
00090 Element *tp = this->GetMatrixArray();
00091 const Element * const tp_last = tp+fNrows;
00092 while (tp < tp_last)
00093 *tp++ += *sp++;
00094 }
00095
00096
00097 template<class Element>
00098 void TVectorT<Element>::Add(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
00099 {
00100
00101
00102 if (gMatrixCheck) {
00103 if ( !AreCompatible(*this,v1) && !AreCompatible(*this,v2)) {
00104 Error("Add(TVectorT<Element> &)","vectors not compatible");
00105 return;
00106 }
00107 }
00108
00109 const Element *sv1 = v1.GetMatrixArray();
00110 const Element *sv2 = v2.GetMatrixArray();
00111 Element *tp = this->GetMatrixArray();
00112 const Element * const tp_last = tp+fNrows;
00113 while (tp < tp_last)
00114 *tp++ = *sv1++ + *sv2++;
00115 }
00116
00117
00118 template<class Element>
00119 Int_t TVectorT<Element>::Memcpy_m(Element *newp,const Element *oldp,Int_t copySize,
00120 Int_t newSize,Int_t oldSize)
00121 {
00122
00123
00124
00125 if (copySize == 0 || oldp == newp) return 0;
00126 else {
00127 if ( newSize <= kSizeMax && oldSize <= kSizeMax ) {
00128
00129 if (newp > oldp) {
00130 for (Int_t i = copySize-1; i >= 0; i--)
00131 newp[i] = oldp[i];
00132 } else {
00133 for (Int_t i = 0; i < copySize; i++)
00134 newp[i] = oldp[i];
00135 }
00136 }
00137 else {
00138 memcpy(newp,oldp,copySize*sizeof(Element));
00139 }
00140 }
00141 return 0;
00142 }
00143
00144
00145 template<class Element>
00146 void TVectorT<Element>::Allocate(Int_t nrows,Int_t row_lwb,Int_t init)
00147 {
00148
00149
00150
00151 fIsOwner = kTRUE;
00152 fNrows = 0;
00153 fRowLwb = 0;
00154 fElements = 0;
00155
00156 if (nrows < 0)
00157 {
00158 Error("Allocate","nrows=%d",nrows);
00159 return;
00160 }
00161
00162 MakeValid();
00163 fNrows = nrows;
00164 fRowLwb = row_lwb;
00165
00166 fElements = New_m(fNrows);
00167 if (init)
00168 memset(fElements,0,fNrows*sizeof(Element));
00169 }
00170
00171
00172 template<class Element>
00173 TVectorT<Element>::TVectorT(Int_t n)
00174 {
00175
00176
00177 Allocate(n,0,1);
00178 }
00179
00180
00181 template<class Element>
00182 TVectorT<Element>::TVectorT(Int_t lwb,Int_t upb)
00183 {
00184
00185
00186 Allocate(upb-lwb+1,lwb,1);
00187 }
00188
00189
00190 template<class Element>
00191 TVectorT<Element>::TVectorT(Int_t n,const Element *elements)
00192 {
00193
00194
00195 Allocate(n,0);
00196 SetElements(elements);
00197 }
00198
00199
00200 template<class Element>
00201 TVectorT<Element>::TVectorT(Int_t lwb,Int_t upb,const Element *elements)
00202 {
00203
00204
00205 Allocate(upb-lwb+1,lwb);
00206 SetElements(elements);
00207 }
00208
00209
00210 template<class Element>
00211 TVectorT<Element>::TVectorT(const TVectorT &another) : TObject(another)
00212 {
00213
00214
00215 R__ASSERT(another.IsValid());
00216 Allocate(another.GetUpb()-another.GetLwb()+1,another.GetLwb());
00217 *this = another;
00218 }
00219
00220
00221 template<class Element>
00222 TVectorT<Element>::TVectorT(const TMatrixTRow_const<Element> &mr) : TObject()
00223 {
00224
00225
00226 const TMatrixTBase<Element> *mt = mr.GetMatrix();
00227 R__ASSERT(mt->IsValid());
00228 Allocate(mt->GetColUpb()-mt->GetColLwb()+1,mt->GetColLwb());
00229 *this = mr;
00230 }
00231
00232
00233 template<class Element>
00234 TVectorT<Element>::TVectorT(const TMatrixTColumn_const<Element> &mc) : TObject()
00235 {
00236
00237
00238 const TMatrixTBase<Element> *mt = mc.GetMatrix();
00239 R__ASSERT(mt->IsValid());
00240 Allocate(mt->GetRowUpb()-mt->GetRowLwb()+1,mt->GetRowLwb());
00241 *this = mc;
00242 }
00243
00244
00245 template<class Element>
00246 TVectorT<Element>::TVectorT(const TMatrixTDiag_const<Element> &md) : TObject()
00247 {
00248
00249
00250 const TMatrixTBase<Element> *mt = md.GetMatrix();
00251 R__ASSERT(mt->IsValid());
00252 Allocate(TMath::Min(mt->GetNrows(),mt->GetNcols()));
00253 *this = md;
00254 }
00255
00256
00257 template<class Element>
00258 TVectorT<Element>::TVectorT(Int_t lwb,Int_t upb,Element va_(iv1), ...)
00259 {
00260
00261
00262
00263
00264
00265 if (upb < lwb) {
00266 Error("TVectorT(Int_t, Int_t, ...)","upb(%d) < lwb(%d)",upb,lwb);
00267 return;
00268 }
00269
00270 Allocate(upb-lwb+1,lwb);
00271
00272 va_list args;
00273 va_start(args,va_(iv1));
00274
00275
00276 (*this)(lwb) = iv1;
00277 for (Int_t i = lwb+1; i <= upb; i++)
00278 (*this)(i) = (Element)va_arg(args,Double_t);
00279
00280 if (strcmp((char *)va_arg(args,char *),"END"))
00281 Error("TVectorT(Int_t, Int_t, ...)","argument list must be terminated by \"END\"");
00282
00283 va_end(args);
00284 }
00285
00286
00287 template<class Element>
00288 TVectorT<Element> &TVectorT<Element>::ResizeTo(Int_t lwb,Int_t upb)
00289 {
00290
00291
00292
00293
00294 R__ASSERT(IsValid());
00295 if (!fIsOwner) {
00296 Error("ResizeTo(lwb,upb)","Not owner of data array,cannot resize");
00297 return *this;
00298 }
00299
00300 const Int_t new_nrows = upb-lwb+1;
00301
00302 if (fNrows > 0) {
00303
00304 if (fNrows == new_nrows && fRowLwb == lwb)
00305 return *this;
00306 else if (new_nrows == 0) {
00307 Clear();
00308 return *this;
00309 }
00310
00311 Element *elements_old = GetMatrixArray();
00312 const Int_t nrows_old = fNrows;
00313 const Int_t rowLwb_old = fRowLwb;
00314
00315 Allocate(new_nrows,lwb);
00316 R__ASSERT(IsValid());
00317 if (fNrows > kSizeMax || nrows_old > kSizeMax)
00318 memset(GetMatrixArray(),0,fNrows*sizeof(Element));
00319 else if (fNrows > nrows_old)
00320 memset(GetMatrixArray()+nrows_old,0,(fNrows-nrows_old)*sizeof(Element));
00321
00322
00323 const Int_t rowLwb_copy = TMath::Max(fRowLwb,rowLwb_old);
00324 const Int_t rowUpb_copy = TMath::Min(fRowLwb+fNrows-1,rowLwb_old+nrows_old-1);
00325 const Int_t nrows_copy = rowUpb_copy-rowLwb_copy+1;
00326
00327 const Int_t nelems_new = fNrows;
00328 Element *elements_new = GetMatrixArray();
00329 if (nrows_copy > 0) {
00330 const Int_t rowOldOff = rowLwb_copy-rowLwb_old;
00331 const Int_t rowNewOff = rowLwb_copy-fRowLwb;
00332 Memcpy_m(elements_new+rowNewOff,elements_old+rowOldOff,nrows_copy,nelems_new,nrows_old);
00333 }
00334
00335 Delete_m(nrows_old,elements_old);
00336 } else {
00337 Allocate(upb-lwb+1,lwb,1);
00338 }
00339
00340 return *this;
00341 }
00342
00343
00344 template<class Element>
00345 TVectorT<Element> &TVectorT<Element>::Use(Int_t lwb,Int_t upb,Element *data)
00346 {
00347
00348
00349 if (upb < lwb) {
00350 Error("Use","upb(%d) < lwb(%d)",upb,lwb);
00351 return *this;
00352 }
00353
00354 Clear();
00355 fNrows = upb-lwb+1;
00356 fRowLwb = lwb;
00357 fElements = data;
00358 fIsOwner = kFALSE;
00359
00360 return *this;
00361 }
00362
00363
00364 template<class Element>
00365 TVectorT<Element> &TVectorT<Element>::GetSub(Int_t row_lwb,Int_t row_upb,TVectorT<Element> &target,Option_t *option) const
00366 {
00367
00368
00369
00370
00371
00372
00373 if (gMatrixCheck) {
00374 R__ASSERT(IsValid());
00375 if (row_lwb < fRowLwb || row_lwb > fRowLwb+fNrows-1) {
00376 Error("GetSub","row_lwb out of bounds");
00377 return target;
00378 }
00379 if (row_upb < fRowLwb || row_upb > fRowLwb+fNrows-1) {
00380 Error("GetSub","row_upb out of bounds");
00381 return target;
00382 }
00383 if (row_upb < row_lwb) {
00384 Error("GetSub","row_upb < row_lwb");
00385 return target;
00386 }
00387 }
00388
00389 TString opt(option);
00390 opt.ToUpper();
00391 const Int_t shift = (opt.Contains("S")) ? 1 : 0;
00392
00393 Int_t row_lwb_sub;
00394 Int_t row_upb_sub;
00395 if (shift) {
00396 row_lwb_sub = 0;
00397 row_upb_sub = row_upb-row_lwb;
00398 } else {
00399 row_lwb_sub = row_lwb;
00400 row_upb_sub = row_upb;
00401 }
00402
00403 target.ResizeTo(row_lwb_sub,row_upb_sub);
00404 const Int_t nrows_sub = row_upb_sub-row_lwb_sub+1;
00405
00406 const Element *ap = this->GetMatrixArray()+(row_lwb-fRowLwb);
00407 Element *bp = target.GetMatrixArray();
00408
00409 for (Int_t irow = 0; irow < nrows_sub; irow++)
00410 *bp++ = *ap++;
00411
00412 return target;
00413 }
00414
00415
00416 template<class Element>
00417 TVectorT<Element> &TVectorT<Element>::SetSub(Int_t row_lwb,const TVectorT<Element> &source)
00418 {
00419
00420
00421
00422 if (gMatrixCheck) {
00423 R__ASSERT(IsValid());
00424 R__ASSERT(source.IsValid());
00425
00426 if (row_lwb < fRowLwb && row_lwb > fRowLwb+fNrows-1) {
00427 Error("SetSub","row_lwb outof bounds");
00428 return *this;
00429 }
00430 if (row_lwb+source.GetNrows() > fRowLwb+fNrows) {
00431 Error("SetSub","source vector too large");
00432 return *this;
00433 }
00434 }
00435
00436 const Int_t nRows_source = source.GetNrows();
00437
00438 const Element *bp = source.GetMatrixArray();
00439 Element *ap = this->GetMatrixArray()+(row_lwb-fRowLwb);
00440
00441 for (Int_t irow = 0; irow < nRows_source; irow++)
00442 *ap++ = *bp++;
00443
00444 return *this;
00445 }
00446
00447
00448 template<class Element>
00449 TVectorT<Element> &TVectorT<Element>::Zero()
00450 {
00451
00452
00453 R__ASSERT(IsValid());
00454 memset(this->GetMatrixArray(),0,fNrows*sizeof(Element));
00455 return *this;
00456 }
00457
00458
00459 template<class Element>
00460 TVectorT<Element> &TVectorT<Element>::Abs()
00461 {
00462
00463
00464 R__ASSERT(IsValid());
00465
00466 Element *ep = this->GetMatrixArray();
00467 const Element * const fp = ep+fNrows;
00468 while (ep < fp) {
00469 *ep = TMath::Abs(*ep);
00470 ep++;
00471 }
00472
00473 return *this;
00474 }
00475
00476
00477 template<class Element>
00478 TVectorT<Element> &TVectorT<Element>::Sqr()
00479 {
00480
00481
00482 R__ASSERT(IsValid());
00483
00484 Element *ep = this->GetMatrixArray();
00485 const Element * const fp = ep+fNrows;
00486 while (ep < fp) {
00487 *ep = (*ep) * (*ep);
00488 ep++;
00489 }
00490
00491 return *this;
00492 }
00493
00494
00495 template<class Element>
00496 TVectorT<Element> &TVectorT<Element>::Sqrt()
00497 {
00498
00499
00500 R__ASSERT(IsValid());
00501
00502 Element *ep = this->GetMatrixArray();
00503 const Element * const fp = ep+fNrows;
00504 while (ep < fp) {
00505 R__ASSERT(*ep >= 0);
00506 if (*ep >= 0)
00507 *ep = TMath::Sqrt(*ep);
00508 else
00509 Error("Sqrt()","v(%ld) = %g < 0",Long_t(ep-this->GetMatrixArray()),(float)*ep);
00510 ep++;
00511 }
00512
00513 return *this;
00514 }
00515
00516
00517 template<class Element>
00518 TVectorT<Element> &TVectorT<Element>::Invert()
00519 {
00520
00521
00522 R__ASSERT(IsValid());
00523
00524 Element *ep = this->GetMatrixArray();
00525 const Element * const fp = ep+fNrows;
00526 while (ep < fp) {
00527 R__ASSERT(*ep != 0.0);
00528 if (*ep != 0.0)
00529 *ep = 1./ *ep;
00530 else
00531 Error("Invert()","v(%ld) = %g",Long_t(ep-this->GetMatrixArray()),(float)*ep);
00532 ep++;
00533 }
00534
00535 return *this;
00536 }
00537
00538
00539 template<class Element>
00540 TVectorT<Element> &TVectorT<Element>::SelectNonZeros(const TVectorT<Element> &select)
00541 {
00542
00543
00544 if (gMatrixCheck && !AreCompatible(*this,select)) {
00545 Error("SelectNonZeros(const TVectorT<Element> &","vector's not compatible");
00546 return *this;
00547 }
00548
00549 const Element *sp = select.GetMatrixArray();
00550 Element *ep = this->GetMatrixArray();
00551 const Element * const fp = ep+fNrows;
00552 while (ep < fp) {
00553 if (*sp == 0.0)
00554 *ep = 0.0;
00555 sp++; ep++;
00556 }
00557
00558 return *this;
00559 }
00560
00561
00562 template<class Element>
00563 Element TVectorT<Element>::Norm1() const
00564 {
00565
00566
00567 R__ASSERT(IsValid());
00568
00569 Element norm = 0;
00570 const Element *ep = this->GetMatrixArray();
00571 const Element * const fp = ep+fNrows;
00572 while (ep < fp)
00573 norm += TMath::Abs(*ep++);
00574
00575 return norm;
00576 }
00577
00578
00579 template<class Element>
00580 Element TVectorT<Element>::Norm2Sqr() const
00581 {
00582
00583
00584 R__ASSERT(IsValid());
00585
00586 Element norm = 0;
00587 const Element *ep = this->GetMatrixArray();
00588 const Element * const fp = ep+fNrows;
00589 while (ep < fp) {
00590 norm += (*ep) * (*ep);
00591 ep++;
00592 }
00593
00594 return norm;
00595 }
00596
00597
00598 template<class Element>
00599 Element TVectorT<Element>::NormInf() const
00600 {
00601
00602
00603 R__ASSERT(IsValid());
00604
00605 Element norm = 0;
00606 const Element *ep = this->GetMatrixArray();
00607 const Element * const fp = ep+fNrows;
00608 while (ep < fp)
00609 norm = TMath::Max(norm,TMath::Abs(*ep++));
00610
00611 return norm;
00612 }
00613
00614
00615 template<class Element>
00616 Int_t TVectorT<Element>::NonZeros() const
00617 {
00618
00619
00620 R__ASSERT(IsValid());
00621
00622 Int_t nr_nonzeros = 0;
00623 const Element *ep = this->GetMatrixArray();
00624 const Element * const fp = ep+fNrows;
00625 while (ep < fp)
00626 if (*ep++) nr_nonzeros++;
00627
00628 return nr_nonzeros;
00629 }
00630
00631
00632 template<class Element>
00633 Element TVectorT<Element>::Sum() const
00634 {
00635
00636
00637 R__ASSERT(IsValid());
00638
00639 Element sum = 0.0;
00640 const Element *ep = this->GetMatrixArray();
00641 const Element * const fp = ep+fNrows;
00642 while (ep < fp)
00643 sum += *ep++;
00644
00645 return sum;
00646 }
00647
00648
00649 template<class Element>
00650 Element TVectorT<Element>::Min() const
00651 {
00652
00653
00654 R__ASSERT(IsValid());
00655
00656 const Int_t index = TMath::LocMin(fNrows,fElements);
00657 return fElements[index];
00658 }
00659
00660
00661 template<class Element>
00662 Element TVectorT<Element>::Max() const
00663 {
00664
00665
00666 R__ASSERT(IsValid());
00667
00668 const Int_t index = TMath::LocMax(fNrows,fElements);
00669 return fElements[index];
00670 }
00671
00672
00673 template<class Element>
00674 TVectorT<Element> &TVectorT<Element>::operator=(const TVectorT<Element> &source)
00675 {
00676
00677
00678
00679
00680 if (gMatrixCheck && !AreCompatible(*this,source)) {
00681 Error("operator=(const TVectorT<Element> &)","vectors not compatible");
00682 return *this;
00683 }
00684
00685 if (this->GetMatrixArray() != source.GetMatrixArray()) {
00686 TObject::operator=(source);
00687 memcpy(fElements,source.GetMatrixArray(),fNrows*sizeof(Element));
00688 }
00689 return *this;
00690 }
00691
00692
00693 template<class Element>
00694 TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTRow_const<Element> &mr)
00695 {
00696
00697
00698 const TMatrixTBase<Element> *mt = mr.GetMatrix();
00699
00700 if (gMatrixCheck) {
00701 R__ASSERT(IsValid());
00702 R__ASSERT(mt->IsValid());
00703 if (mt->GetColLwb() != fRowLwb || mt->GetNcols() != fNrows) {
00704 Error("operator=(const TMatrixTRow_const &)","vector and row not compatible");
00705 return *this;
00706 }
00707 }
00708
00709 const Int_t inc = mr.GetInc();
00710 const Element *rp = mr.GetPtr();
00711 Element *ep = this->GetMatrixArray();
00712 const Element * const fp = ep+fNrows;
00713 while (ep < fp) {
00714 *ep++ = *rp;
00715 rp += inc;
00716 }
00717
00718 R__ASSERT(rp == mr.GetPtr()+mt->GetNcols());
00719
00720 return *this;
00721 }
00722
00723
00724 template<class Element>
00725 TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTColumn_const<Element> &mc)
00726 {
00727
00728
00729 const TMatrixTBase<Element> *mt = mc.GetMatrix();
00730
00731 if (gMatrixCheck) {
00732 R__ASSERT(IsValid());
00733 R__ASSERT(mt->IsValid());
00734 if (mt->GetRowLwb() != fRowLwb || mt->GetNrows() != fNrows) {
00735 Error("operator=(const TMatrixTColumn_const &)","vector and column not compatible");
00736 return *this;
00737 }
00738 }
00739
00740 const Int_t inc = mc.GetInc();
00741 const Element *cp = mc.GetPtr();
00742 Element *ep = this->GetMatrixArray();
00743 const Element * const fp = ep+fNrows;
00744 while (ep < fp) {
00745 *ep++ = *cp;
00746 cp += inc;
00747 }
00748
00749 R__ASSERT(cp == mc.GetPtr()+mt->GetNoElements());
00750
00751 return *this;
00752 }
00753
00754
00755 template<class Element>
00756 TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTDiag_const<Element> &md)
00757 {
00758
00759
00760 const TMatrixTBase<Element> *mt = md.GetMatrix();
00761
00762 if (gMatrixCheck) {
00763 R__ASSERT(IsValid());
00764 R__ASSERT(mt->IsValid());
00765 if (md.GetNdiags() != fNrows) {
00766 Error("operator=(const TMatrixTDiag_const &)","vector and matrix-diagonal not compatible");
00767 return *this;
00768 }
00769 }
00770
00771 const Int_t inc = md.GetInc();
00772 const Element *dp = md.GetPtr();
00773 Element *ep = this->GetMatrixArray();
00774 const Element * const fp = ep+fNrows;
00775 while (ep < fp) {
00776 *ep++ = *dp;
00777 dp += inc;
00778 }
00779
00780 R__ASSERT(dp < md.GetPtr()+mt->GetNoElements()+inc);
00781
00782 return *this;
00783 }
00784
00785
00786 template<class Element>
00787 TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTSparseRow_const<Element> &mr)
00788 {
00789
00790
00791
00792 const TMatrixTBase<Element> *mt = mr.GetMatrix();
00793
00794 if (gMatrixCheck) {
00795 R__ASSERT(IsValid());
00796 R__ASSERT(mt->IsValid());
00797 if (mt->GetColLwb() != fRowLwb || mt->GetNcols() != fNrows) {
00798 Error("operator=(const TMatrixTSparseRow_const &)","vector and row not compatible");
00799 return *this;
00800 }
00801 }
00802
00803 const Int_t nIndex = mr.GetNindex();
00804 const Element * const prData = mr.GetDataPtr();
00805 const Int_t * const prCol = mr.GetColPtr();
00806 Element * const pvData = this->GetMatrixArray();
00807
00808 memset(pvData,0,fNrows*sizeof(Element));
00809 for (Int_t index = 0; index < nIndex; index++) {
00810 const Int_t icol = prCol[index];
00811 pvData[icol] = prData[index];
00812 }
00813
00814 return *this;
00815 }
00816
00817
00818 template<class Element>
00819 TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTSparseDiag_const<Element> &md)
00820 {
00821
00822
00823 const TMatrixTBase<Element> *mt = md.GetMatrix();
00824
00825 if (gMatrixCheck) {
00826 R__ASSERT(IsValid());
00827 R__ASSERT(mt->IsValid());
00828 if (md.GetNdiags() != fNrows) {
00829 Error("operator=(const TMatrixTSparseDiag_const &)","vector and matrix-diagonal not compatible");
00830 return *this;
00831 }
00832 }
00833
00834 Element * const pvData = this->GetMatrixArray();
00835 for (Int_t idiag = 0; idiag < fNrows; idiag++)
00836 pvData[idiag] = md(idiag);
00837
00838 return *this;
00839 }
00840
00841
00842 template<class Element>
00843 TVectorT<Element> &TVectorT<Element>::operator=(Element val)
00844 {
00845
00846
00847 R__ASSERT(IsValid());
00848
00849 Element *ep = this->GetMatrixArray();
00850 const Element * const fp = ep+fNrows;
00851 while (ep < fp)
00852 *ep++ = val;
00853
00854 return *this;
00855 }
00856
00857
00858 template<class Element>
00859 TVectorT<Element> &TVectorT<Element>::operator+=(Element val)
00860 {
00861
00862
00863 R__ASSERT(IsValid());
00864
00865 Element *ep = this->GetMatrixArray();
00866 const Element * const fp = ep+fNrows;
00867 while (ep < fp)
00868 *ep++ += val;
00869
00870 return *this;
00871 }
00872
00873
00874 template<class Element>
00875 TVectorT<Element> &TVectorT<Element>::operator-=(Element val)
00876 {
00877
00878
00879 R__ASSERT(IsValid());
00880
00881 Element *ep = this->GetMatrixArray();
00882 const Element * const fp = ep+fNrows;
00883 while (ep < fp)
00884 *ep++ -= val;
00885
00886 return *this;
00887 }
00888
00889
00890 template<class Element>
00891 TVectorT<Element> &TVectorT<Element>::operator*=(Element val)
00892 {
00893
00894
00895 R__ASSERT(IsValid());
00896
00897 Element *ep = this->GetMatrixArray();
00898 const Element * const fp = ep+fNrows;
00899 while (ep < fp)
00900 *ep++ *= val;
00901
00902 return *this;
00903 }
00904
00905
00906 template<class Element>
00907 TVectorT<Element> &TVectorT<Element>::operator+=(const TVectorT<Element> &source)
00908 {
00909
00910
00911 if (gMatrixCheck && !AreCompatible(*this,source)) {
00912 Error("operator+=(const TVectorT<Element> &)","vector's not compatible");
00913 return *this;
00914 }
00915
00916 const Element *sp = source.GetMatrixArray();
00917 Element *tp = this->GetMatrixArray();
00918 const Element * const tp_last = tp+fNrows;
00919 while (tp < tp_last)
00920 *tp++ += *sp++;
00921
00922 return *this;
00923 }
00924
00925
00926 template<class Element>
00927 TVectorT<Element> &TVectorT<Element>::operator-=(const TVectorT<Element> &source)
00928 {
00929
00930
00931 if (gMatrixCheck && !AreCompatible(*this,source)) {
00932 Error("operator-=(const TVectorT<Element> &)","vector's not compatible");
00933 return *this;
00934 }
00935
00936 const Element *sp = source.GetMatrixArray();
00937 Element *tp = this->GetMatrixArray();
00938 const Element * const tp_last = tp+fNrows;
00939 while (tp < tp_last)
00940 *tp++ -= *sp++;
00941
00942 return *this;
00943 }
00944
00945
00946 template<class Element>
00947 TVectorT<Element> &TVectorT<Element>::operator*=(const TMatrixT<Element> &a)
00948 {
00949
00950
00951
00952 if (gMatrixCheck) {
00953 R__ASSERT(IsValid());
00954 R__ASSERT(a.IsValid());
00955 if (a.GetNcols() != fNrows || a.GetColLwb() != fRowLwb) {
00956 Error("operator*=(const TMatrixT &)","vector and matrix incompatible");
00957 return *this;
00958 }
00959 }
00960
00961 const Bool_t doResize = (fNrows != a.GetNrows() || fRowLwb != a.GetRowLwb());
00962 if (doResize && !fIsOwner) {
00963 Error("operator*=(const TMatrixT &)","vector has to be resized but not owner");
00964 return *this;
00965 }
00966
00967 Element work[kWorkMax];
00968 Bool_t isAllocated = kFALSE;
00969 Element *elements_old = work;
00970 const Int_t nrows_old = fNrows;
00971 if (nrows_old > kWorkMax) {
00972 isAllocated = kTRUE;
00973 elements_old = new Element[nrows_old];
00974 }
00975 memcpy(elements_old,fElements,nrows_old*sizeof(Element));
00976
00977 if (doResize) {
00978 const Int_t rowlwb_new = a.GetRowLwb();
00979 const Int_t nrows_new = a.GetNrows();
00980 ResizeTo(rowlwb_new,rowlwb_new+nrows_new-1);
00981 }
00982 memset(fElements,0,fNrows*sizeof(Element));
00983
00984 const Element *mp = a.GetMatrixArray();
00985 Element *tp = this->GetMatrixArray();
00986 #ifdef CBLAS
00987 if (typeid(Element) == typeid(Double_t))
00988 cblas_dgemv(CblasRowMajor,CblasNoTrans,a.GetNrows(),a.GetNcols(),1.0,mp,
00989 a.GetNcols(),elements_old,1,0.0,tp,1);
00990 else if (typeid(Element) != typeid(Float_t))
00991 cblas_sgemv(CblasRowMajor,CblasNoTrans,a.GetNrows(),a.GetNcols(),1.0,mp,
00992 a.GetNcols(),elements_old,1,0.0,tp,1);
00993 else
00994 Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
00995 #else
00996 const Element * const tp_last = tp+fNrows;
00997 while (tp < tp_last) {
00998 Element sum = 0;
00999 for (const Element *sp = elements_old; sp < elements_old+nrows_old; )
01000 sum += *sp++ * *mp++;
01001 *tp++ = sum;
01002 }
01003 R__ASSERT(mp == a.GetMatrixArray()+a.GetNoElements());
01004 #endif
01005
01006 if (isAllocated)
01007 delete [] elements_old;
01008
01009 return *this;
01010 }
01011
01012
01013 template<class Element>
01014 TVectorT<Element> &TVectorT<Element>::operator*=(const TMatrixTSparse<Element> &a)
01015 {
01016
01017
01018
01019 if (gMatrixCheck) {
01020 R__ASSERT(IsValid());
01021 R__ASSERT(a.IsValid());
01022 if (a.GetNcols() != fNrows || a.GetColLwb() != fRowLwb) {
01023 Error("operator*=(const TMatrixTSparse &)","vector and matrix incompatible");
01024 return *this;
01025 }
01026 }
01027
01028 const Bool_t doResize = (fNrows != a.GetNrows() || fRowLwb != a.GetRowLwb());
01029 if (doResize && !fIsOwner) {
01030 Error("operator*=(const TMatrixTSparse &)","vector has to be resized but not owner");
01031 return *this;
01032 }
01033
01034 Element work[kWorkMax];
01035 Bool_t isAllocated = kFALSE;
01036 Element *elements_old = work;
01037 const Int_t nrows_old = fNrows;
01038 if (nrows_old > kWorkMax) {
01039 isAllocated = kTRUE;
01040 elements_old = new Element[nrows_old];
01041 }
01042 memcpy(elements_old,fElements,nrows_old*sizeof(Element));
01043
01044 if (doResize) {
01045 const Int_t rowlwb_new = a.GetRowLwb();
01046 const Int_t nrows_new = a.GetNrows();
01047 ResizeTo(rowlwb_new,rowlwb_new+nrows_new-1);
01048 }
01049 memset(fElements,0,fNrows*sizeof(Element));
01050
01051 const Int_t * const pRowIndex = a.GetRowIndexArray();
01052 const Int_t * const pColIndex = a.GetColIndexArray();
01053 const Element * const mp = a.GetMatrixArray();
01054
01055 const Element * const sp = elements_old;
01056 Element * tp = this->GetMatrixArray();
01057
01058 for (Int_t irow = 0; irow < fNrows; irow++) {
01059 const Int_t sIndex = pRowIndex[irow];
01060 const Int_t eIndex = pRowIndex[irow+1];
01061 Element sum = 0.0;
01062 for (Int_t index = sIndex; index < eIndex; index++) {
01063 const Int_t icol = pColIndex[index];
01064 sum += mp[index]*sp[icol];
01065 }
01066 tp[irow] = sum;
01067 }
01068
01069 if (isAllocated)
01070 delete [] elements_old;
01071
01072 return *this;
01073 }
01074
01075
01076 template<class Element>
01077 TVectorT<Element> &TVectorT<Element>::operator*=(const TMatrixTSym<Element> &a)
01078 {
01079
01080
01081
01082 if (gMatrixCheck) {
01083 R__ASSERT(IsValid());
01084 R__ASSERT(a.IsValid());
01085 if (a.GetNcols() != fNrows || a.GetColLwb() != fRowLwb) {
01086 Error("operator*=(const TMatrixTSym &)","vector and matrix incompatible");
01087 return *this;
01088 }
01089 }
01090
01091 Element work[kWorkMax];
01092 Bool_t isAllocated = kFALSE;
01093 Element *elements_old = work;
01094 const Int_t nrows_old = fNrows;
01095 if (nrows_old > kWorkMax) {
01096 isAllocated = kTRUE;
01097 elements_old = new Element[nrows_old];
01098 }
01099 memcpy(elements_old,fElements,nrows_old*sizeof(Element));
01100 memset(fElements,0,fNrows*sizeof(Element));
01101
01102 const Element *mp = a.GetMatrixArray();
01103 Element *tp = this->GetMatrixArray();
01104 #ifdef CBLAS
01105 if (typeid(Element) == typeid(Double_t))
01106 cblas_dsymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
01107 fNrows,elements_old,1,0.0,tp,1);
01108 else if (typeid(Element) != typeid(Float_t))
01109 cblas_ssymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
01110 fNrows,elements_old,1,0.0,tp,1);
01111 else
01112 Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
01113 #else
01114 const Element * const tp_last = tp+fNrows;
01115 while (tp < tp_last) {
01116 Element sum = 0;
01117 for (const Element *sp = elements_old; sp < elements_old+nrows_old; )
01118 sum += *sp++ * *mp++;
01119 *tp++ = sum;
01120 }
01121 R__ASSERT(mp == a.GetMatrixArray()+a.GetNoElements());
01122 #endif
01123
01124 if (isAllocated)
01125 delete [] elements_old;
01126
01127 return *this;
01128 }
01129
01130
01131 template<class Element>
01132 Bool_t TVectorT<Element>::operator==(Element val) const
01133 {
01134
01135
01136 R__ASSERT(IsValid());
01137
01138 const Element *ep = this->GetMatrixArray();
01139 const Element * const fp = ep+fNrows;
01140 while (ep < fp)
01141 if (!(*ep++ == val))
01142 return kFALSE;
01143
01144 return kTRUE;
01145 }
01146
01147
01148 template<class Element>
01149 Bool_t TVectorT<Element>::operator!=(Element val) const
01150 {
01151
01152
01153 R__ASSERT(IsValid());
01154
01155 const Element *ep = this->GetMatrixArray();
01156 const Element * const fp = ep+fNrows;
01157 while (ep < fp)
01158 if (!(*ep++ != val))
01159 return kFALSE;
01160
01161 return kTRUE;
01162 }
01163
01164
01165 template<class Element>
01166 Bool_t TVectorT<Element>::operator<(Element val) const
01167 {
01168
01169
01170 R__ASSERT(IsValid());
01171
01172 const Element *ep = this->GetMatrixArray();
01173 const Element * const fp = ep+fNrows;
01174 while (ep < fp)
01175 if (!(*ep++ < val))
01176 return kFALSE;
01177
01178 return kTRUE;
01179 }
01180
01181
01182 template<class Element>
01183 Bool_t TVectorT<Element>::operator<=(Element val) const
01184 {
01185
01186
01187 R__ASSERT(IsValid());
01188
01189 const Element *ep = this->GetMatrixArray();
01190 const Element * const fp = ep+fNrows;
01191 while (ep < fp)
01192 if (!(*ep++ <= val))
01193 return kFALSE;
01194
01195 return kTRUE;
01196 }
01197
01198
01199 template<class Element>
01200 Bool_t TVectorT<Element>::operator>(Element val) const
01201 {
01202
01203
01204 R__ASSERT(IsValid());
01205
01206 const Element *ep = this->GetMatrixArray();
01207 const Element * const fp = ep+fNrows;
01208 while (ep < fp)
01209 if (!(*ep++ > val))
01210 return kFALSE;
01211
01212 return kTRUE;
01213 }
01214
01215
01216 template<class Element>
01217 Bool_t TVectorT<Element>::operator>=(Element val) const
01218 {
01219
01220
01221 R__ASSERT(IsValid());
01222
01223 const Element *ep = this->GetMatrixArray();
01224 const Element * const fp = ep+fNrows;
01225 while (ep < fp)
01226 if (!(*ep++ >= val))
01227 return kFALSE;
01228
01229 return kTRUE;
01230 }
01231
01232
01233 template<class Element>
01234 Bool_t TVectorT<Element>::MatchesNonZeroPattern(const TVectorT<Element> &select)
01235 {
01236
01237
01238 if (gMatrixCheck && !AreCompatible(*this,select)) {
01239 Error("MatchesNonZeroPattern(const TVectorT&)","vector's not compatible");
01240 return kFALSE;
01241 }
01242
01243 const Element *sp = select.GetMatrixArray();
01244 const Element *ep = this->GetMatrixArray();
01245 const Element * const fp = ep+fNrows;
01246 while (ep < fp) {
01247 if (*sp == 0.0 && *ep != 0.0)
01248 return kFALSE;
01249 sp++; ep++;
01250 }
01251
01252 return kTRUE;
01253 }
01254
01255
01256 template<class Element>
01257 Bool_t TVectorT<Element>::SomePositive(const TVectorT<Element> &select)
01258 {
01259
01260
01261 if (gMatrixCheck && !AreCompatible(*this,select)) {
01262 Error("SomePositive(const TVectorT&)","vector's not compatible");
01263 return kFALSE;
01264 }
01265
01266 const Element *sp = select.GetMatrixArray();
01267 const Element *ep = this->GetMatrixArray();
01268 const Element * const fp = ep+fNrows;
01269 while (ep < fp) {
01270 if (*sp != 0.0 && *ep <= 0.0)
01271 return kFALSE;
01272 sp++; ep++;
01273 }
01274
01275 return kTRUE;
01276 }
01277
01278
01279 template<class Element>
01280 void TVectorT<Element>::AddSomeConstant(Element val,const TVectorT<Element> &select)
01281 {
01282
01283
01284 if (gMatrixCheck && !AreCompatible(*this,select))
01285 Error("AddSomeConstant(Element,const TVectorT&)(const TVectorT&)","vector's not compatible");
01286
01287 const Element *sp = select.GetMatrixArray();
01288 Element *ep = this->GetMatrixArray();
01289 const Element * const fp = ep+fNrows;
01290 while (ep < fp) {
01291 if (*sp)
01292 *ep += val;
01293 sp++; ep++;
01294 }
01295 }
01296
01297 extern Double_t Drand(Double_t &ix);
01298
01299
01300 template<class Element>
01301 void TVectorT<Element>::Randomize(Element alpha,Element beta,Double_t &seed)
01302 {
01303
01304
01305 R__ASSERT(IsValid());
01306
01307 const Element scale = beta-alpha;
01308 const Element shift = alpha/scale;
01309
01310 Element * ep = GetMatrixArray();
01311 const Element * const fp = ep+fNrows;
01312 while (ep < fp)
01313 *ep++ = scale*(Drand(seed)+shift);
01314 }
01315
01316
01317 template<class Element>
01318 TVectorT<Element> &TVectorT<Element>::Apply(const TElementActionT<Element> &action)
01319 {
01320
01321
01322 R__ASSERT(IsValid());
01323 for (Element *ep = fElements; ep < fElements+fNrows; ep++)
01324 action.Operation(*ep);
01325 return *this;
01326 }
01327
01328
01329 template<class Element>
01330 TVectorT<Element> &TVectorT<Element>::Apply(const TElementPosActionT<Element> &action)
01331 {
01332
01333
01334
01335 R__ASSERT(IsValid());
01336
01337 Element *ep = fElements;
01338 for (action.fI = fRowLwb; action.fI < fRowLwb+fNrows; action.fI++)
01339 action.Operation(*ep++);
01340
01341 R__ASSERT(ep == fElements+fNrows);
01342
01343 return *this;
01344 }
01345
01346
01347 template<class Element>
01348 void TVectorT<Element>::Draw(Option_t *option)
01349 {
01350
01351
01352
01353 gROOT->ProcessLine(Form("THistPainter::PaintSpecialObjects((TObject*)0x%lx,\"%s\");",
01354 (ULong_t)this, option));
01355 }
01356
01357
01358 template<class Element>
01359 void TVectorT<Element>::Print(Option_t *flag) const
01360 {
01361
01362
01363 if (!IsValid()) {
01364 Error("Print","Vector is invalid");
01365 return;
01366 }
01367
01368 printf("\nVector (%d) %s is as follows",fNrows,flag);
01369
01370 printf("\n\n | %6d |", 1);
01371 printf("\n%s\n", "------------------");
01372 for (Int_t i = 0; i < fNrows; i++) {
01373 printf("%4d |",i+fRowLwb);
01374
01375 printf("%g \n",(*this)(i+fRowLwb));
01376 }
01377 printf("\n");
01378 }
01379
01380
01381 template<class Element>
01382 Bool_t operator==(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
01383 {
01384
01385
01386 if (!AreCompatible(v1,v2)) return kFALSE;
01387 return (memcmp(v1.GetMatrixArray(),v2.GetMatrixArray(),v1.GetNrows()*sizeof(Element)) == 0);
01388 }
01389
01390
01391 template<class Element>
01392 Element operator*(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
01393 {
01394
01395
01396 if (gMatrixCheck) {
01397 if (!AreCompatible(v1,v2)) {
01398 Error("operator*(const TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
01399 return 0.0;
01400 }
01401 }
01402
01403 return Dot(v1,v2);
01404 }
01405
01406
01407 template<class Element>
01408 TVectorT<Element> operator+(const TVectorT<Element> &source1,const TVectorT<Element> &source2)
01409 {
01410
01411
01412 TVectorT<Element> target = source1;
01413 target += source2;
01414 return target;
01415 }
01416
01417
01418 template<class Element>
01419 TVectorT<Element> operator-(const TVectorT<Element> &source1,const TVectorT<Element> &source2)
01420 {
01421
01422
01423 TVectorT<Element> target = source1;
01424 target -= source2;
01425 return target;
01426 }
01427
01428
01429 template<class Element>
01430 TVectorT<Element> operator*(const TMatrixT<Element> &a,const TVectorT<Element> &source)
01431 {
01432
01433
01434 R__ASSERT(a.IsValid());
01435 TVectorT<Element> target(a.GetRowLwb(),a.GetRowUpb());
01436 return Add(target,Element(1.0),a,source);
01437 }
01438
01439
01440 template<class Element>
01441 TVectorT<Element> operator*(const TMatrixTSym<Element> &a,const TVectorT<Element> &source)
01442 {
01443
01444
01445 R__ASSERT(a.IsValid());
01446 TVectorT<Element> target(a.GetRowLwb(),a.GetRowUpb());
01447 return Add(target,Element(1.0),a,source);
01448 }
01449
01450
01451 template<class Element>
01452 TVectorT<Element> operator*(const TMatrixTSparse<Element> &a,const TVectorT<Element> &source)
01453 {
01454
01455
01456 R__ASSERT(a.IsValid());
01457 TVectorT<Element> target(a.GetRowLwb(),a.GetRowUpb());
01458 return Add(target,Element(1.0),a,source);
01459 }
01460
01461
01462 template<class Element>
01463 TVectorT<Element> operator*(Element val,const TVectorT<Element> &source)
01464 {
01465
01466
01467 TVectorT<Element> target = source;
01468 target *= val;
01469 return target;
01470 }
01471
01472
01473 template<class Element>
01474 Element Dot(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
01475 {
01476
01477
01478 const Element *v1p = v1.GetMatrixArray();
01479 const Element *v2p = v2.GetMatrixArray();
01480 Element sum = 0.0;
01481 const Element * const fv1p = v1p+v1.GetNrows();
01482 while (v1p < fv1p)
01483 sum += *v1p++ * *v2p++;
01484
01485 return sum;
01486 }
01487
01488
01489 template <class Element1, class Element2>
01490 TMatrixT<Element1>
01491 OuterProduct(const TVectorT<Element1> &v1,const TVectorT<Element2> &v2)
01492 {
01493
01494
01495
01496
01497
01498
01499 TMatrixT<Element1> target;
01500
01501 return OuterProduct(target,v1,v2);
01502 }
01503
01504
01505 template <class Element1,class Element2,class Element3>
01506 TMatrixT<Element1>
01507 &OuterProduct(TMatrixT<Element1> &target,const TVectorT<Element2> &v1,const TVectorT<Element3> &v2)
01508 {
01509
01510
01511 target.ResizeTo(v1.GetLwb(), v1.GetUpb(), v2.GetLwb(), v2.GetUpb());
01512
01513 Element1 * mp = target.GetMatrixArray();
01514 const Element1 * const m_last = mp + target.GetNoElements();
01515
01516 const Element2 * v1p = v1.GetMatrixArray();
01517 const Element2 * const v1_last = v1p + v1.GetNrows();
01518
01519 const Element3 * const v20 = v2.GetMatrixArray();
01520 const Element3 * v2p = v20;
01521 const Element3 * const v2_last = v2p + v2.GetNrows();
01522
01523 while (v1p < v1_last) {
01524 v2p = v20;
01525 while (v2p < v2_last) {
01526 *mp++ = *v1p * *v2p++ ;
01527 }
01528 v1p++;
01529 }
01530
01531 R__ASSERT(v1p == v1_last && mp == m_last && v2p == v2_last);
01532
01533 return target;
01534 }
01535
01536
01537 template <class Element1, class Element2, class Element3>
01538 Element1 Mult(const TVectorT<Element1> &v1,const TMatrixT<Element2> &m,
01539 const TVectorT<Element3> &v2)
01540 {
01541
01542
01543 if (gMatrixCheck) {
01544 if (!AreCompatible(v1, m)) {
01545 ::Error("Mult", "Vector v1 and matrix m incompatible");
01546 return 0;
01547 }
01548 if (!AreCompatible(m, v2)) {
01549 ::Error("Mult", "Matrix m and vector v2 incompatible");
01550 return 0;
01551 }
01552 }
01553
01554 const Element1 * v1p = v1.GetMatrixArray();
01555 const Element1 * const v1_last = v1p + v1.GetNrows();
01556
01557 const Element2 * mp = m.GetMatrixArray();
01558 const Element2 * const m_last = mp + m.GetNoElements();
01559
01560 const Element3 * const v20 = v2.GetMatrixArray();
01561 const Element3 * v2p = v20;
01562 const Element3 * const v2_last = v2p + v2.GetNrows();
01563
01564 Element1 sum = 0;
01565 Element3 dot = 0;
01566
01567 while (v1p < v1_last) {
01568 v2p = v20;
01569 while (v2p < v2_last) {
01570 dot += *mp++ * *v2p++;
01571 }
01572 sum += *v1p++ * dot;
01573 dot = 0;
01574 }
01575
01576 R__ASSERT(v1p == v1_last && mp == m_last && v2p == v2_last);
01577
01578 return sum;
01579 }
01580
01581
01582 template<class Element>
01583 TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,const TVectorT<Element> &source)
01584 {
01585
01586
01587 if (gMatrixCheck && !AreCompatible(target,source)) {
01588 Error("Add(TVectorT<Element> &,Element,const TVectorT<Element> &)","vector's are incompatible");
01589 return target;
01590 }
01591
01592 const Element * sp = source.GetMatrixArray();
01593 Element * tp = target.GetMatrixArray();
01594 const Element * const ftp = tp+target.GetNrows();
01595 if (scalar == 1.0 ) {
01596 while ( tp < ftp )
01597 *tp++ += *sp++;
01598 } else if (scalar == -1.0) {
01599 while ( tp < ftp )
01600 *tp++ -= *sp++;
01601 } else {
01602 while ( tp < ftp )
01603 *tp++ += scalar * *sp++;
01604 }
01605
01606 return target;
01607 }
01608
01609
01610 template<class Element>
01611 TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,
01612 const TMatrixT<Element> &a,const TVectorT<Element> &source)
01613 {
01614
01615
01616
01617 if (gMatrixCheck) {
01618 R__ASSERT(target.IsValid());
01619 R__ASSERT(a.IsValid());
01620 if (a.GetNrows() != target.GetNrows() || a.GetRowLwb() != target.GetLwb()) {
01621 Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","target vector and matrix are incompatible");
01622 return target;
01623 }
01624
01625 R__ASSERT(source.IsValid());
01626 if (a.GetNcols() != source.GetNrows() || a.GetColLwb() != source.GetLwb()) {
01627 Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","source vector and matrix are incompatible");
01628 return target;
01629 }
01630 }
01631
01632 const Element * const sp = source.GetMatrixArray();
01633 const Element * mp = a.GetMatrixArray();
01634 Element * tp = target.GetMatrixArray();
01635 #ifdef CBLAS
01636 if (typeid(Element) == typeid(Double_t))
01637 cblas_dsymv(CblasRowMajor,CblasUpper,fNrows,scalar,mp,
01638 fNrows,sp,1,0.0,tp,1);
01639 else if (typeid(Element) != typeid(Float_t))
01640 cblas_ssymv(CblasRowMajor,CblasUpper,fNrows,scalar,mp,
01641 fNrows,sp,1,0.0,tp,1);
01642 else
01643 Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
01644 #else
01645 const Element * const sp_last = sp+source.GetNrows();
01646 const Element * const tp_last = tp+target.GetNrows();
01647 if (scalar == 1.0) {
01648 while (tp < tp_last) {
01649 const Element *sp1 = sp;
01650 Element sum = 0;
01651 while (sp1 < sp_last)
01652 sum += *sp1++ * *mp++;
01653 *tp++ += sum;
01654 }
01655 } else if (scalar == 0.0) {
01656 while (tp < tp_last) {
01657 const Element *sp1 = sp;
01658 Element sum = 0;
01659 while (sp1 < sp_last)
01660 sum += *sp1++ * *mp++;
01661 *tp++ = sum;
01662 }
01663 } else if (scalar == -1.0) {
01664 while (tp < tp_last) {
01665 const Element *sp1 = sp;
01666 Element sum = 0;
01667 while (sp1 < sp_last)
01668 sum += *sp1++ * *mp++;
01669 *tp++ -= sum;
01670 }
01671 } else {
01672 while (tp < tp_last) {
01673 const Element *sp1 = sp;
01674 Element sum = 0;
01675 while (sp1 < sp_last)
01676 sum += *sp1++ * *mp++;
01677 *tp++ += scalar * sum;
01678 }
01679 }
01680
01681 if (gMatrixCheck) R__ASSERT(mp == a.GetMatrixArray()+a.GetNoElements());
01682 #endif
01683
01684 return target;
01685 }
01686
01687
01688 template<class Element>
01689 TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,
01690 const TMatrixTSym<Element> &a,const TVectorT<Element> &source)
01691 {
01692
01693
01694
01695 if (gMatrixCheck) {
01696 R__ASSERT(target.IsValid());
01697 R__ASSERT(source.IsValid());
01698 R__ASSERT(a.IsValid());
01699 if (a.GetNrows() != target.GetNrows() || a.GetRowLwb() != target.GetLwb()) {
01700 Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","target vector and matrix are incompatible");
01701 return target;
01702 }
01703 }
01704
01705 const Element * const sp = source.GetMatrixArray();
01706 const Element * mp = a.GetMatrixArray();
01707 Element * tp = target.GetMatrixArray();
01708 #ifdef CBLAS
01709 if (typeid(Element) == typeid(Double_t))
01710 cblas_dsymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
01711 fNrows,sp,1,0.0,tp,1);
01712 else if (typeid(Element) != typeid(Float_t))
01713 cblas_ssymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
01714 fNrows,sp,1,0.0,tp,1);
01715 else
01716 Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
01717 #else
01718 const Element * const sp_last = sp+source.GetNrows();
01719 const Element * const tp_last = tp+target.GetNrows();
01720 if (scalar == 1.0) {
01721 while (tp < tp_last) {
01722 const Element *sp1 = sp;
01723 Element sum = 0;
01724 while (sp1 < sp_last)
01725 sum += *sp1++ * *mp++;
01726 *tp++ += sum;
01727 }
01728 } else if (scalar == 0.0) {
01729 while (tp < tp_last) {
01730 const Element *sp1 = sp;
01731 Element sum = 0;
01732 while (sp1 < sp_last)
01733 sum += *sp1++ * *mp++;
01734 *tp++ = sum;
01735 }
01736 } else if (scalar == -1.0) {
01737 while (tp < tp_last) {
01738 const Element *sp1 = sp;
01739 Element sum = 0;
01740 while (sp1 < sp_last)
01741 sum += *sp1++ * *mp++;
01742 *tp++ -= sum;
01743 }
01744 } else {
01745 while (tp < tp_last) {
01746 const Element *sp1 = sp;
01747 Element sum = 0;
01748 while (sp1 < sp_last)
01749 sum += *sp1++ * *mp++;
01750 *tp++ += scalar * sum;
01751 }
01752 }
01753 R__ASSERT(mp == a.GetMatrixArray()+a.GetNoElements());
01754 #endif
01755
01756 return target;
01757 }
01758
01759
01760 template<class Element>
01761 TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,
01762 const TMatrixTSparse<Element> &a,const TVectorT<Element> &source)
01763 {
01764
01765
01766
01767 if (gMatrixCheck) {
01768 R__ASSERT(target.IsValid());
01769 R__ASSERT(a.IsValid());
01770 if (a.GetNrows() != target.GetNrows() || a.GetRowLwb() != target.GetLwb()) {
01771 Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","target vector and matrix are incompatible");
01772 return target;
01773 }
01774
01775 R__ASSERT(source.IsValid());
01776 if (a.GetNcols() != source.GetNrows() || a.GetColLwb() != source.GetLwb()) {
01777 Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","source vector and matrix are incompatible");
01778 return target;
01779 }
01780 }
01781
01782 const Int_t * const pRowIndex = a.GetRowIndexArray();
01783 const Int_t * const pColIndex = a.GetColIndexArray();
01784 const Element * const mp = a.GetMatrixArray();
01785
01786 const Element * const sp = source.GetMatrixArray();
01787 Element * tp = target.GetMatrixArray();
01788
01789 if (scalar == 1.0) {
01790 for (Int_t irow = 0; irow < a.GetNrows(); irow++) {
01791 const Int_t sIndex = pRowIndex[irow];
01792 const Int_t eIndex = pRowIndex[irow+1];
01793 Element sum = 0.0;
01794 for (Int_t index = sIndex; index < eIndex; index++) {
01795 const Int_t icol = pColIndex[index];
01796 sum += mp[index]*sp[icol];
01797 }
01798 tp[irow] += sum;
01799 }
01800 } else if (scalar == 0.0) {
01801 for (Int_t irow = 0; irow < a.GetNrows(); irow++) {
01802 const Int_t sIndex = pRowIndex[irow];
01803 const Int_t eIndex = pRowIndex[irow+1];
01804 Element sum = 0.0;
01805 for (Int_t index = sIndex; index < eIndex; index++) {
01806 const Int_t icol = pColIndex[index];
01807 sum += mp[index]*sp[icol];
01808 }
01809 tp[irow] = sum;
01810 }
01811 } else if (scalar == -1.0) {
01812 for (Int_t irow = 0; irow < a.GetNrows(); irow++) {
01813 const Int_t sIndex = pRowIndex[irow];
01814 const Int_t eIndex = pRowIndex[irow+1];
01815 Element sum = 0.0;
01816 for (Int_t index = sIndex; index < eIndex; index++) {
01817 const Int_t icol = pColIndex[index];
01818 sum += mp[index]*sp[icol];
01819 }
01820 tp[irow] -= sum;
01821 }
01822 } else {
01823 for (Int_t irow = 0; irow < a.GetNrows(); irow++) {
01824 const Int_t sIndex = pRowIndex[irow];
01825 const Int_t eIndex = pRowIndex[irow+1];
01826 Element sum = 0.0;
01827 for (Int_t index = sIndex; index < eIndex; index++) {
01828 const Int_t icol = pColIndex[index];
01829 sum += mp[index]*sp[icol];
01830 }
01831 tp[irow] += scalar * sum;
01832 }
01833 }
01834
01835 return target;
01836 }
01837
01838
01839 template<class Element>
01840 TVectorT<Element> &AddElemMult(TVectorT<Element> &target,Element scalar,
01841 const TVectorT<Element> &source1,const TVectorT<Element> &source2)
01842 {
01843
01844
01845 if (gMatrixCheck && !(AreCompatible(target,source1) && AreCompatible(target,source1))) {
01846 Error("AddElemMult(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &)",
01847 "vector's are incompatible");
01848 return target;
01849 }
01850
01851 const Element * sp1 = source1.GetMatrixArray();
01852 const Element * sp2 = source2.GetMatrixArray();
01853 Element * tp = target.GetMatrixArray();
01854 const Element * const ftp = tp+target.GetNrows();
01855
01856 if (scalar == 1.0 ) {
01857 while ( tp < ftp )
01858 *tp++ += *sp1++ * *sp2++;
01859 } else if (scalar == -1.0) {
01860 while ( tp < ftp )
01861 *tp++ -= *sp1++ * *sp2++;
01862 } else {
01863 while ( tp < ftp )
01864 *tp++ += scalar * *sp1++ * *sp2++;
01865 }
01866
01867 return target;
01868 }
01869
01870
01871 template<class Element>
01872 TVectorT<Element> &AddElemMult(TVectorT<Element> &target,Element scalar,
01873 const TVectorT<Element> &source1,const TVectorT<Element> &source2,const TVectorT<Element> &select)
01874 {
01875
01876
01877
01878 if (gMatrixCheck && !( AreCompatible(target,source1) && AreCompatible(target,source1) &&
01879 AreCompatible(target,select) )) {
01880 Error("AddElemMult(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &,onst TVectorT<Element> &)",
01881 "vector's are incompatible");
01882 return target;
01883 }
01884
01885 const Element * sp1 = source1.GetMatrixArray();
01886 const Element * sp2 = source2.GetMatrixArray();
01887 const Element * mp = select.GetMatrixArray();
01888 Element * tp = target.GetMatrixArray();
01889 const Element * const ftp = tp+target.GetNrows();
01890
01891 if (scalar == 1.0 ) {
01892 while ( tp < ftp ) {
01893 if (*mp) *tp += *sp1 * *sp2;
01894 mp++; tp++; sp1++; sp2++;
01895 }
01896 } else if (scalar == -1.0) {
01897 while ( tp < ftp ) {
01898 if (*mp) *tp -= *sp1 * *sp2;
01899 mp++; tp++; sp1++; sp2++;
01900 }
01901 } else {
01902 while ( tp < ftp ) {
01903 if (*mp) *tp += scalar * *sp1 * *sp2;
01904 mp++; tp++; sp1++; sp2++;
01905 }
01906 }
01907
01908 return target;
01909 }
01910
01911
01912 template<class Element>
01913 TVectorT<Element> &AddElemDiv(TVectorT<Element> &target,Element scalar,
01914 const TVectorT<Element> &source1,const TVectorT<Element> &source2)
01915 {
01916
01917
01918 if (gMatrixCheck && !(AreCompatible(target,source1) && AreCompatible(target,source1))) {
01919 Error("AddElemDiv(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &)",
01920 "vector's are incompatible");
01921 return target;
01922 }
01923
01924 const Element * sp1 = source1.GetMatrixArray();
01925 const Element * sp2 = source2.GetMatrixArray();
01926 Element * tp = target.GetMatrixArray();
01927 const Element * const ftp = tp+target.GetNrows();
01928
01929 if (scalar == 1.0 ) {
01930 while ( tp < ftp ) {
01931 if (*sp2 != 0.0)
01932 *tp += *sp1 / *sp2;
01933 else {
01934 const Int_t irow = (sp2-source2.GetMatrixArray())/source2.GetNrows();
01935 Error("AddElemDiv","source2 (%d) is zero",irow);
01936 }
01937 tp++; sp1++; sp2++;
01938 }
01939 } else if (scalar == -1.0) {
01940 while ( tp < ftp ) {
01941 if (*sp2 != 0.0)
01942 *tp -= *sp1 / *sp2;
01943 else {
01944 const Int_t irow = (sp2-source2.GetMatrixArray())/source2.GetNrows();
01945 Error("AddElemDiv","source2 (%d) is zero",irow);
01946 }
01947 tp++; sp1++; sp2++;
01948 }
01949 } else {
01950 while ( tp < ftp ) {
01951 if (*sp2 != 0.0)
01952 *tp += scalar * *sp1 / *sp2;
01953 else {
01954 const Int_t irow = (sp2-source2.GetMatrixArray())/source2.GetNrows();
01955 Error("AddElemDiv","source2 (%d) is zero",irow);
01956 }
01957 tp++; sp1++; sp2++;
01958 }
01959 }
01960
01961 return target;
01962 }
01963
01964
01965 template<class Element>
01966 TVectorT<Element> &AddElemDiv(TVectorT<Element> &target,Element scalar,
01967 const TVectorT<Element> &source1,const TVectorT<Element> &source2,const TVectorT<Element> &select)
01968 {
01969
01970
01971
01972 if (gMatrixCheck && !( AreCompatible(target,source1) && AreCompatible(target,source1) &&
01973 AreCompatible(target,select) )) {
01974 Error("AddElemDiv(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &,onst TVectorT<Element> &)",
01975 "vector's are incompatible");
01976 return target;
01977 }
01978
01979 const Element * sp1 = source1.GetMatrixArray();
01980 const Element * sp2 = source2.GetMatrixArray();
01981 const Element * mp = select.GetMatrixArray();
01982 Element * tp = target.GetMatrixArray();
01983 const Element * const ftp = tp+target.GetNrows();
01984
01985 if (scalar == 1.0 ) {
01986 while ( tp < ftp ) {
01987 if (*mp) {
01988 if (*sp2 != 0.0)
01989 *tp += *sp1 / *sp2;
01990 else {
01991 const Int_t irow = (sp2-source2.GetMatrixArray())/source2.GetNrows();
01992 Error("AddElemDiv","source2 (%d) is zero",irow);
01993 }
01994 }
01995 mp++; tp++; sp1++; sp2++;
01996 }
01997 } else if (scalar == -1.0) {
01998 while ( tp < ftp ) {
01999 if (*mp) {
02000 if (*sp2 != 0.0)
02001 *tp -= *sp1 / *sp2;
02002 else {
02003 const Int_t irow = (sp2-source2.GetMatrixArray())/source2.GetNrows();
02004 Error("AddElemDiv","source2 (%d) is zero",irow);
02005 }
02006 }
02007 mp++; tp++; sp1++; sp2++;
02008 }
02009 } else {
02010 while ( tp < ftp ) {
02011 if (*mp) {
02012 if (*sp2 != 0.0)
02013 *tp += scalar * *sp1 / *sp2;
02014 else {
02015 const Int_t irow = (sp2-source2.GetMatrixArray())/source2.GetNrows();
02016 Error("AddElemDiv","source2 (%d) is zero",irow);
02017 }
02018 }
02019 mp++; tp++; sp1++; sp2++;
02020 }
02021 }
02022
02023 return target;
02024 }
02025
02026
02027 template<class Element>
02028 TVectorT<Element> &ElementMult(TVectorT<Element> &target,const TVectorT<Element> &source)
02029 {
02030
02031
02032 if (gMatrixCheck && !AreCompatible(target,source)) {
02033 Error("ElementMult(TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
02034 return target;
02035 }
02036
02037 const Element * sp = source.GetMatrixArray();
02038 Element * tp = target.GetMatrixArray();
02039 const Element * const ftp = tp+target.GetNrows();
02040 while ( tp < ftp )
02041 *tp++ *= *sp++;
02042
02043 return target;
02044 }
02045
02046
02047 template<class Element>
02048 TVectorT<Element> &ElementMult(TVectorT<Element> &target,const TVectorT<Element> &source,const TVectorT<Element> &select)
02049 {
02050
02051
02052 if (gMatrixCheck && !(AreCompatible(target,source) && AreCompatible(target,select))) {
02053 Error("ElementMult(TVectorT<Element> &,const TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
02054 return target;
02055 }
02056
02057 const Element * sp = source.GetMatrixArray();
02058 const Element * mp = select.GetMatrixArray();
02059 Element * tp = target.GetMatrixArray();
02060 const Element * const ftp = tp+target.GetNrows();
02061 while ( tp < ftp ) {
02062 if (*mp) *tp *= *sp;
02063 mp++; tp++; sp++;
02064 }
02065
02066 return target;
02067 }
02068
02069
02070 template<class Element>
02071 TVectorT<Element> &ElementDiv(TVectorT<Element> &target,const TVectorT<Element> &source)
02072 {
02073
02074
02075 if (gMatrixCheck && !AreCompatible(target,source)) {
02076 Error("ElementDiv(TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
02077 return target;
02078 }
02079
02080 const Element * sp = source.GetMatrixArray();
02081 Element * tp = target.GetMatrixArray();
02082 const Element * const ftp = tp+target.GetNrows();
02083 while ( tp < ftp ) {
02084 if (*sp != 0.0)
02085 *tp++ /= *sp++;
02086 else {
02087 const Int_t irow = (sp-source.GetMatrixArray())/source.GetNrows();
02088 Error("ElementDiv","source (%d) is zero",irow);
02089 }
02090 }
02091
02092 return target;
02093 }
02094
02095
02096 template<class Element>
02097 TVectorT<Element> &ElementDiv(TVectorT<Element> &target,const TVectorT<Element> &source,const TVectorT<Element> &select)
02098 {
02099
02100
02101 if (gMatrixCheck && !AreCompatible(target,source)) {
02102 Error("ElementDiv(TVectorT<Element> &,const TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
02103 return target;
02104 }
02105
02106 const Element * sp = source.GetMatrixArray();
02107 const Element * mp = select.GetMatrixArray();
02108 Element * tp = target.GetMatrixArray();
02109 const Element * const ftp = tp+target.GetNrows();
02110 while ( tp < ftp ) {
02111 if (*mp) {
02112 if (*sp != 0.0)
02113 *tp /= *sp;
02114 else {
02115 const Int_t irow = (sp-source.GetMatrixArray())/source.GetNrows();
02116 Error("ElementDiv","source (%d) is zero",irow);
02117 }
02118 }
02119 mp++; tp++; sp++;
02120 }
02121
02122 return target;
02123 }
02124
02125
02126 template<class Element1,class Element2>
02127 Bool_t AreCompatible(const TVectorT<Element1> &v1,const TVectorT<Element2> &v2,Int_t verbose)
02128 {
02129
02130
02131 if (!v1.IsValid()) {
02132 if (verbose)
02133 ::Error("AreCompatible", "vector 1 not valid");
02134 return kFALSE;
02135 }
02136 if (!v2.IsValid()) {
02137 if (verbose)
02138 ::Error("AreCompatible", "vector 2 not valid");
02139 return kFALSE;
02140 }
02141
02142 if (v1.GetNrows() != v2.GetNrows() || v1.GetLwb() != v2.GetLwb()) {
02143 if (verbose)
02144 ::Error("AreCompatible", "matrices 1 and 2 not compatible");
02145 return kFALSE;
02146 }
02147
02148 return kTRUE;
02149 }
02150
02151
02152 template<class Element1, class Element2>
02153 Bool_t AreCompatible(const TMatrixT<Element1> &m,const TVectorT<Element2> &v,Int_t verbose)
02154 {
02155
02156
02157 if (!m.IsValid()) {
02158 if (verbose)
02159 ::Error("AreCompatible", "Matrix not valid");
02160 return kFALSE;
02161 }
02162 if (!v.IsValid()) {
02163 if (verbose)
02164 ::Error("AreCompatible", "vector not valid");
02165 return kFALSE;
02166 }
02167
02168 if (m.GetNcols() != v.GetNrows() ) {
02169 if (verbose)
02170 ::Error("AreCompatible", "matrix and vector not compatible");
02171 return kFALSE;
02172 }
02173
02174 return kTRUE;
02175 }
02176
02177
02178 template<class Element1, class Element2>
02179 Bool_t AreCompatible(const TVectorT<Element1> &v,const TMatrixT<Element2> &m,Int_t verbose)
02180 {
02181
02182
02183 if (!m.IsValid()) {
02184 if (verbose)
02185 ::Error("AreCompatible", "Matrix not valid");
02186 return kFALSE;
02187 }
02188 if (!v.IsValid()) {
02189 if (verbose)
02190 ::Error("AreCompatible", "vector not valid");
02191 return kFALSE;
02192 }
02193
02194 if (v.GetNrows() != m.GetNrows() ) {
02195 if (verbose)
02196 ::Error("AreCompatible", "vector and matrix not compatible");
02197 return kFALSE;
02198 }
02199
02200 return kTRUE;
02201 }
02202
02203
02204 template<class Element>
02205 void Compare(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
02206 {
02207
02208
02209 if (!AreCompatible(v1,v2)) {
02210 Error("Compare(const TVectorT<Element> &,const TVectorT<Element> &)","vectors are incompatible");
02211 return;
02212 }
02213
02214 printf("\n\nComparison of two TVectorTs:\n");
02215
02216 Element norm1 = 0;
02217 Element norm2 = 0;
02218 Element ndiff = 0;
02219 Int_t imax = 0;
02220 Element difmax = -1;
02221 const Element *mp1 = v1.GetMatrixArray();
02222 const Element *mp2 = v2.GetMatrixArray();
02223
02224 for (Int_t i = 0; i < v1.GetNrows(); i++) {
02225 const Element mv1 = *mp1++;
02226 const Element mv2 = *mp2++;
02227 const Element diff = TMath::Abs(mv1-mv2);
02228
02229 if (diff > difmax) {
02230 difmax = diff;
02231 imax = i;
02232 }
02233 norm1 += TMath::Abs(mv1);
02234 norm2 += TMath::Abs(mv2);
02235 ndiff += TMath::Abs(diff);
02236 }
02237
02238 imax += v1.GetLwb();
02239 printf("\nMaximal discrepancy \t\t%g",difmax);
02240 printf("\n occured at the point\t\t(%d)",imax);
02241 const Element mv1 = v1(imax);
02242 const Element mv2 = v2(imax);
02243 printf("\n Vector 1 element is \t\t%g",mv1);
02244 printf("\n Vector 2 element is \t\t%g",mv2);
02245 printf("\n Absolute error v2[i]-v1[i]\t\t%g",mv2-mv1);
02246 printf("\n Relative error\t\t\t\t%g\n",
02247 (mv2-mv1)/TMath::Max(TMath::Abs(mv2+mv1)/2,(Element)1e-7));
02248
02249 printf("\n||Vector 1|| \t\t\t%g",norm1);
02250 printf("\n||Vector 2|| \t\t\t%g",norm2);
02251 printf("\n||Vector1-Vector2||\t\t\t\t%g",ndiff);
02252 printf("\n||Vector1-Vector2||/sqrt(||Vector1|| ||Vector2||)\t%g\n\n",
02253 ndiff/TMath::Max(TMath::Sqrt(norm1*norm2),1e-7));
02254 }
02255
02256
02257 template<class Element>
02258 Bool_t VerifyVectorValue(const TVectorT<Element> &v,Element val,
02259 Int_t verbose,Element maxDevAllow)
02260 {
02261
02262
02263 Int_t imax = 0;
02264 Element maxDevObs = 0;
02265
02266 if (TMath::Abs(maxDevAllow) <= 0.0)
02267 maxDevAllow = std::numeric_limits<Element>::epsilon();
02268
02269 for (Int_t i = v.GetLwb(); i <= v.GetUpb(); i++) {
02270 const Element dev = TMath::Abs(v(i)-val);
02271 if (dev > maxDevObs) {
02272 imax = i;
02273 maxDevObs = dev;
02274 }
02275 }
02276
02277 if (maxDevObs == 0)
02278 return kTRUE;
02279
02280 if (verbose) {
02281 printf("Largest dev for (%d); dev = |%g - %g| = %g\n",imax,v(imax),val,maxDevObs);
02282 if (maxDevObs > maxDevAllow)
02283 Error("VerifyVectorValue","Deviation > %g\n",maxDevAllow);
02284 }
02285
02286 if (maxDevObs > maxDevAllow)
02287 return kFALSE;
02288 return kTRUE;
02289 }
02290
02291
02292 template<class Element>
02293 Bool_t VerifyVectorIdentity(const TVectorT<Element> &v1,const TVectorT<Element> &v2,
02294 Int_t verbose, Element maxDevAllow)
02295 {
02296
02297
02298 Int_t imax = 0;
02299 Element maxDevObs = 0;
02300
02301 if (!AreCompatible(v1,v2))
02302 return kFALSE;
02303
02304 if (TMath::Abs(maxDevAllow) <= 0.0)
02305 maxDevAllow = std::numeric_limits<Element>::epsilon();
02306
02307 for (Int_t i = v1.GetLwb(); i <= v1.GetUpb(); i++) {
02308 const Element dev = TMath::Abs(v1(i)-v2(i));
02309 if (dev > maxDevObs) {
02310 imax = i;
02311 maxDevObs = dev;
02312 }
02313 }
02314
02315 if (maxDevObs == 0)
02316 return kTRUE;
02317
02318 if (verbose) {
02319 printf("Largest dev for (%d); dev = |%g - %g| = %g\n",imax,v1(imax),v2(imax),maxDevObs);
02320 if (maxDevObs > maxDevAllow)
02321 Error("VerifyVectorIdentity","Deviation > %g\n",maxDevAllow);
02322 }
02323
02324 if (maxDevObs > maxDevAllow) {
02325 return kFALSE;
02326 }
02327 return kTRUE;
02328 }
02329
02330
02331 template<class Element>
02332 void TVectorT<Element>::Streamer(TBuffer &R__b)
02333 {
02334
02335
02336 if (R__b.IsReading()) {
02337 UInt_t R__s, R__c;
02338 Version_t R__v = R__b.ReadVersion(&R__s,&R__c);
02339 if (R__v > 1) {
02340 Clear();
02341 R__b.ReadClassBuffer(TVectorT<Element>::Class(),this,R__v,R__s,R__c);
02342 } else {
02343 TObject::Streamer(R__b);
02344 R__b >> fRowLwb;
02345 fNrows = R__b.ReadArray(fElements);
02346 R__b.CheckByteCount(R__s, R__c, TVectorT<Element>::IsA());
02347 }
02348 if (fNrows > 0 && fNrows <= kSizeMax) {
02349 memcpy(fDataStack,fElements,fNrows*sizeof(Element));
02350 delete [] fElements;
02351 fElements = fDataStack;
02352 } else if (fNrows < 0)
02353 Invalidate();
02354
02355 if (R__v < 3)
02356 MakeValid();
02357 } else {
02358 R__b.WriteClassBuffer(TVectorT<Element>::Class(),this);
02359 }
02360 }
02361
02362 #ifndef ROOT_TMatrixFfwd
02363 #include "TMatrixFfwd.h"
02364 #endif
02365 #ifndef ROOT_TMatrixFSymfwd
02366 #include "TMatrixFSymfwd.h"
02367 #endif
02368 #ifndef ROOT_TMatrixFSparsefwd
02369 #include "TMatrixFSparsefwd.h"
02370 #endif
02371
02372 template class TVectorT<Float_t>;
02373
02374 template Bool_t operator== <Float_t>(const TVectorF &source1,const TVectorF &source2);
02375 template TVectorF operator+ <Float_t>(const TVectorF &source1,const TVectorF &source2);
02376 template TVectorF operator- <Float_t>(const TVectorF &source1,const TVectorF &source2);
02377 template Float_t operator* <Float_t>(const TVectorF &source1,const TVectorF &source2);
02378 template TVectorF operator* <Float_t>(const TMatrixF &a, const TVectorF &source);
02379 template TVectorF operator* <Float_t>(const TMatrixFSym &a, const TVectorF &source);
02380 template TVectorF operator* <Float_t>(const TMatrixFSparse &a, const TVectorF &source);
02381 template TVectorF operator* <Float_t>( Float_t val, const TVectorF &source);
02382
02383 template Float_t Dot <Float_t>(const TVectorF &v1, const TVectorF &v2);
02384 template TMatrixF OuterProduct <Float_t,Float_t>
02385 (const TVectorF &v1, const TVectorF &v2);
02386 template TMatrixF &OuterProduct <Float_t,Float_t,Float_t>
02387 ( TMatrixF &target, const TVectorF &v1, const TVectorF &v2);
02388 template Float_t Mult <Float_t,Float_t,Float_t>
02389 (const TVectorF &v1, const TMatrixF &m, const TVectorF &v2);
02390
02391 template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source);
02392 template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TMatrixF &a,
02393 const TVectorF &source);
02394 template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TMatrixFSym &a,
02395 const TVectorF &source);
02396 template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TMatrixFSparse &a,
02397 const TVectorF &source);
02398 template TVectorF &AddElemMult <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
02399 const TVectorF &source2);
02400 template TVectorF &AddElemMult <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
02401 const TVectorF &source2,const TVectorF &select);
02402 template TVectorF &AddElemDiv <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
02403 const TVectorF &source2);
02404 template TVectorF &AddElemDiv <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
02405 const TVectorF &source2,const TVectorF &select);
02406 template TVectorF &ElementMult <Float_t>( TVectorF &target, const TVectorF &source);
02407 template TVectorF &ElementMult <Float_t>( TVectorF &target, const TVectorF &source, const TVectorF &select);
02408 template TVectorF &ElementDiv <Float_t>( TVectorF &target, const TVectorF &source);
02409 template TVectorF &ElementDiv <Float_t>( TVectorF &target, const TVectorF &source, const TVectorF &select);
02410
02411 template Bool_t AreCompatible <Float_t,Float_t> (const TVectorF &v1,const TVectorF &v2,Int_t verbose);
02412 template Bool_t AreCompatible <Float_t,Double_t>(const TVectorF &v1,const TVectorD &v2,Int_t verbose);
02413 template Bool_t AreCompatible <Float_t,Float_t> (const TMatrixF &m, const TVectorF &v, Int_t verbose);
02414 template Bool_t AreCompatible <Float_t,Float_t> (const TVectorF &v, const TMatrixF &m, Int_t verbose);
02415
02416 template void Compare <Float_t> (const TVectorF &v1,const TVectorF &v2);
02417 template Bool_t VerifyVectorValue <Float_t> (const TVectorF &m, Float_t val,Int_t verbose,Float_t maxDevAllow);
02418 template Bool_t VerifyVectorIdentity<Float_t> (const TVectorF &m1,const TVectorF &m2, Int_t verbose,Float_t maxDevAllow);
02419
02420 #ifndef ROOT_TMatrixDfwd
02421 #include "TMatrixDfwd.h"
02422 #endif
02423 #ifndef ROOT_TMatrixDSymfwd
02424 #include "TMatrixDSymfwd.h"
02425 #endif
02426 #ifndef ROOT_TMatrixDSparsefwd
02427 #include "TMatrixDSparsefwd.h"
02428 #endif
02429
02430 template class TVectorT<Double_t>;
02431
02432 template Bool_t operator== <Double_t>(const TVectorD &source1,const TVectorD &source2);
02433 template TVectorD operator+ <Double_t>(const TVectorD &source1,const TVectorD &source2);
02434 template TVectorD operator- <Double_t>(const TVectorD &source1,const TVectorD &source2);
02435 template Double_t operator* <Double_t>(const TVectorD &source1,const TVectorD &source2);
02436 template TVectorD operator* <Double_t>(const TMatrixD &a, const TVectorD &source);
02437 template TVectorD operator* <Double_t>(const TMatrixDSym &a, const TVectorD &source);
02438 template TVectorD operator* <Double_t>(const TMatrixDSparse &a, const TVectorD &source);
02439 template TVectorD operator* <Double_t>( Double_t val, const TVectorD &source);
02440
02441 template Double_t Dot <Double_t>(const TVectorD &v1, const TVectorD &v2);
02442 template TMatrixD OuterProduct <Double_t,Double_t>
02443 (const TVectorD &v1, const TVectorD &v2);
02444 template TMatrixD &OuterProduct <Double_t,Double_t,Double_t>
02445 ( TMatrixD &target, const TVectorD &v1, const TVectorD &v2);
02446 template Double_t Mult <Double_t,Double_t,Double_t>
02447 (const TVectorD &v1, const TMatrixD &m, const TVectorD &v2);
02448
02449 template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source);
02450 template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TMatrixD &a,
02451 const TVectorD &source);
02452 template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TMatrixDSym &a
02453 , const TVectorD &source);
02454 template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TMatrixDSparse &a
02455 , const TVectorD &source);
02456 template TVectorD &AddElemMult <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
02457 const TVectorD &source2);
02458 template TVectorD &AddElemMult <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
02459 const TVectorD &source2,const TVectorD &select);
02460 template TVectorD &AddElemDiv <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
02461 const TVectorD &source2);
02462 template TVectorD &AddElemDiv <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
02463 const TVectorD &source2,const TVectorD &select);
02464 template TVectorD &ElementMult <Double_t>( TVectorD &target, const TVectorD &source);
02465 template TVectorD &ElementMult <Double_t>( TVectorD &target, const TVectorD &source, const TVectorD &select);
02466 template TVectorD &ElementDiv <Double_t>( TVectorD &target, const TVectorD &source);
02467 template TVectorD &ElementDiv <Double_t>( TVectorD &target, const TVectorD &source, const TVectorD &select);
02468
02469 template Bool_t AreCompatible <Double_t,Double_t>(const TVectorD &v1,const TVectorD &v2,Int_t verbose);
02470 template Bool_t AreCompatible <Double_t,Float_t> (const TVectorD &v1,const TVectorF &v2,Int_t verbose);
02471 template Bool_t AreCompatible <Double_t,Double_t>(const TMatrixD &m, const TVectorD &v, Int_t verbose);
02472 template Bool_t AreCompatible <Double_t,Double_t>(const TVectorD &v, const TMatrixD &m, Int_t verbose);
02473
02474 template void Compare <Double_t> (const TVectorD &v1,const TVectorD &v2);
02475 template Bool_t VerifyVectorValue <Double_t> (const TVectorD &m, Double_t val,Int_t verbose,Double_t maxDevAllow);
02476 template Bool_t VerifyVectorIdentity<Double_t> (const TVectorD &m1,const TVectorD &m2, Int_t verbose,Double_t maxDevAllow);