CLHEP 2.4.6.4
C++ Class Library for High Energy Physics
Loading...
Searching...
No Matches
CLHEP::HepVector Class Reference

#include <Vector.h>

+ Inheritance diagram for CLHEP::HepVector:

Public Member Functions

 HepVector ()
 
 HepVector (int p)
 
 HepVector (int p, int)
 
 HepVector (int p, HepRandom &r)
 
 HepVector (const HepVector &v)
 
 HepVector (const HepMatrix &m)
 
virtual ~HepVector ()
 
const doubleoperator() (int row) const
 
doubleoperator() (int row)
 
const doubleoperator[] (int row) const
 
doubleoperator[] (int row)
 
virtual const doubleoperator() (int row, int col) const
 
virtual doubleoperator() (int row, int col)
 
HepVectoroperator*= (double t)
 
HepVectoroperator/= (double t)
 
HepVectoroperator+= (const HepMatrix &v2)
 
HepVectoroperator+= (const HepVector &v2)
 
HepVectoroperator-= (const HepMatrix &v2)
 
HepVectoroperator-= (const HepVector &v2)
 
HepVectoroperator= (const HepVector &hm2)
 
HepVectoroperator= (const HepMatrix &)
 
HepVectoroperator= (const Hep3Vector &)
 
HepVector operator- () const
 
HepVector apply (double(*f)(double, int)) const
 
HepVector sub (int min_row, int max_row) const
 
HepVector sub (int min_row, int max_row)
 
void sub (int row, const HepVector &v1)
 
double normsq () const
 
double norm () const
 
virtual int num_row () const
 
virtual int num_col () const
 
HepMatrix T () const
 
- Public Member Functions inherited from CLHEP::HepGenMatrix
virtual ~HepGenMatrix ()
 
virtual int num_row () const =0
 
virtual int num_col () const =0
 
virtual const doubleoperator() (int row, int col) const =0
 
virtual doubleoperator() (int row, int col)=0
 
virtual void invert (int &)=0
 
HepGenMatrix_row operator[] (int)
 
const HepGenMatrix_row_const operator[] (int) const
 
virtual bool operator== (const HepGenMatrix &) const
 

Protected Member Functions

virtual int num_size () const
 
- Protected Member Functions inherited from CLHEP::HepGenMatrix
virtual int num_size () const =0
 
void delete_m (int size, double *)
 
doublenew_m (int size)
 

Friends

class HepDiagMatrix
 
class HepSymMatrix
 
class HepMatrix
 
void swap (HepVector &v1, HepVector &v2)
 
double dot (const HepVector &v1, const HepVector &v2)
 
HepVector operator+ (const HepVector &v1, const HepVector &v2)
 
HepVector operator- (const HepVector &v1, const HepVector &v2)
 
HepVector operator* (const HepSymMatrix &hm1, const HepVector &hm2)
 
HepVector operator* (const HepDiagMatrix &hm1, const HepVector &hm2)
 
HepMatrix operator* (const HepVector &hm1, const HepMatrix &hm2)
 
HepVector operator* (const HepMatrix &hm1, const HepVector &hm2)
 
HepVector solve (const HepMatrix &a, const HepVector &v)
 
void tridiagonal (HepSymMatrix *a, HepMatrix *hsm)
 
void row_house (HepMatrix *, const HepMatrix &, double, int, int, int, int)
 
void row_house (HepMatrix *, const HepVector &, double, int, int)
 
void back_solve (const HepMatrix &R, HepVector *b)
 
void col_house (HepMatrix *, const HepMatrix &, double, int, int, int, int)
 
HepVector house (const HepSymMatrix &a, int row, int col)
 
HepVector house (const HepMatrix &a, int row, int col)
 
void house_with_update (HepMatrix *a, int row, int col)
 
HepSymMatrix vT_times_v (const HepVector &v)
 
HepVector qr_solve (HepMatrix *, const HepVector &)
 

Additional Inherited Members

- Public Types inherited from CLHEP::HepGenMatrix
enum  { size_max = 25 }
 
typedef std::vector< double, Alloc< double, 25 > >::iterator mIter
 
typedef std::vector< double, Alloc< double, 25 > >::const_iterator mcIter
 
- Static Public Member Functions inherited from CLHEP::HepGenMatrix
static void swap (int &, int &)
 
static void swap (std::vector< double, Alloc< double, 25 > > &, std::vector< double, Alloc< double, 25 > > &)
 
static void error (const char *s)
 

Detailed Description

Author

Definition at line 35 of file Vector.h.

Constructor & Destructor Documentation

◆ HepVector() [1/6]

CLHEP::HepVector::HepVector ( )
inline

◆ HepVector() [2/6]

CLHEP::HepVector::HepVector ( int  p)
explicit

Definition at line 53 of file Vector.cc.

54 : m(p), nrow(p)
55{
56}

◆ HepVector() [3/6]

CLHEP::HepVector::HepVector ( int  p,
int  init 
)

Definition at line 58 of file Vector.cc.

59 : m(p), nrow(p)
60{
61 switch (init)
62 {
63 case 0:
64 m.assign(p,0);
65 break;
66
67 case 1:
68 {
69 mIter e = m.begin() + nrow;
70 for (mIter i=m.begin(); i<e; i++) *i = 1.0;
71 break;
72 }
73
74 default:
75 error("Vector: initialization must be either 0 or 1.");
76 }
77}
std::vector< double, Alloc< double, 25 > >::iterator mIter
Definition: GenMatrix.h:73
static void error(const char *s)
Definition: GenMatrix.cc:70
void init()
Definition: testRandom.cc:29

◆ HepVector() [4/6]

CLHEP::HepVector::HepVector ( int  p,
HepRandom r 
)

Definition at line 79 of file Vector.cc.

80 : m(p), nrow(p)
81{
82 HepGenMatrix::mIter a = m.begin();
83 HepGenMatrix::mIter b = m.begin() + nrow;
84 for(;a<b;a++) *a = r();
85}

◆ HepVector() [5/6]

CLHEP::HepVector::HepVector ( const HepVector v)

Definition at line 94 of file Vector.cc.

95 : HepGenMatrix(hm1), m(hm1.nrow), nrow(hm1.nrow)
96{
97 m = hm1.m;
98}

◆ HepVector() [6/6]

CLHEP::HepVector::HepVector ( const HepMatrix m)

Definition at line 105 of file Vector.cc.

106 : m(hm1.nrow), nrow(hm1.nrow)
107{
108 if (hm1.num_col() != 1)
109 error("Vector::Vector(Matrix) : Matrix is not Nx1");
110
111 m = hm1.m;
112}

◆ ~HepVector()

CLHEP::HepVector::~HepVector ( )
virtual

Definition at line 91 of file Vector.cc.

91 {
92}

Member Function Documentation

◆ apply()

HepVector CLHEP::HepVector::apply ( double(*)(double, int)  f) const

Definition at line 554 of file Vector.cc.

558{
559#else
560{
561 HepVector mret(num_row());
562#endif
563 HepGenMatrix::mcIter a = m.begin();
564 HepGenMatrix::mIter b = mret.m.begin();
565 for(int ir=1;ir<=num_row();ir++) {
566 *(b++) = (*f)(*(a++), ir);
567 }
568 return mret;
569}
std::vector< double, Alloc< double, 25 > >::const_iterator mcIter
Definition: GenMatrix.h:74
virtual int num_row() const
Definition: Vector.cc:116
void f(void g())
Definition: excDblThrow.cc:38

Referenced by main().

◆ norm()

double CLHEP::HepVector::norm ( ) const
inline

Referenced by main().

◆ normsq()

double CLHEP::HepVector::normsq ( ) const
inline

Referenced by main().

◆ num_col()

int CLHEP::HepVector::num_col ( ) const
virtual

Implements CLHEP::HepGenMatrix.

Definition at line 118 of file Vector.cc.

118{ return 1; }

Referenced by CLHEP::HepMatrix::operator+=(), and CLHEP::HepMatrix::operator-=().

◆ num_row()

◆ num_size()

int CLHEP::HepVector::num_size ( ) const
protectedvirtual

Implements CLHEP::HepGenMatrix.

Definition at line 117 of file Vector.cc.

117{return nrow;}

◆ operator()() [1/4]

double & CLHEP::HepVector::operator() ( int  row)
inline

◆ operator()() [2/4]

const double & CLHEP::HepVector::operator() ( int  row) const
inline

◆ operator()() [3/4]

double & CLHEP::HepVector::operator() ( int  row,
int  col 
)
virtual

Implements CLHEP::HepGenMatrix.

Definition at line 128 of file Vector.cc.

129{
130#endif
131
132 return *(m.begin()+(row-1));
133}

◆ operator()() [4/4]

const double & CLHEP::HepVector::operator() ( int  row,
int  col 
) const
virtual

Implements CLHEP::HepGenMatrix.

Definition at line 141 of file Vector.cc.

142{
143#endif
144
145 return *(m.begin()+(row-1));
146}

◆ operator*=()

HepVector & CLHEP::HepVector::operator*= ( double  t)

Definition at line 449 of file Vector.cc.

450{
451 SIMPLE_UOP(*=)
452 return (*this);
453}
#define SIMPLE_UOP(OPER)
Definition: DiagMatrix.cc:26

◆ operator+=() [1/2]

HepVector & CLHEP::HepVector::operator+= ( const HepMatrix v2)

Definition at line 408 of file Vector.cc.

409{
410 CHK_DIM_2(num_row(),hm2.num_row(),1,hm2.num_col(),+=);
411 SIMPLE_BOP(+=)
412 return (*this);
413}
#define CHK_DIM_2(r1, r2, c1, c2, fun)
Definition: DiagMatrix.cc:44
#define SIMPLE_BOP(OPER)
Definition: DiagMatrix.cc:31

◆ operator+=() [2/2]

HepVector & CLHEP::HepVector::operator+= ( const HepVector v2)

Definition at line 415 of file Vector.cc.

416{
417 CHK_DIM_1(num_row(),hm2.num_row(),+=);
418 SIMPLE_BOP(+=)
419 return (*this);
420}
#define CHK_DIM_1(c1, r2, fun)
Definition: DiagMatrix.cc:49

◆ operator-()

HepVector CLHEP::HepVector::operator- ( ) const

Definition at line 212 of file Vector.cc.

215{
216#else
217{
218 HepVector hm2(nrow);
219#endif
220 HepGenMatrix::mcIter a=m.begin();
221 HepGenMatrix::mIter b=hm2.m.begin();
222 HepGenMatrix::mcIter e=m.begin()+num_size();
223 for(;a<e; a++, b++) (*b) = -(*a);
224 return hm2;
225}
virtual int num_size() const
Definition: Vector.cc:117

◆ operator-=() [1/2]

HepVector & CLHEP::HepVector::operator-= ( const HepMatrix v2)

Definition at line 429 of file Vector.cc.

430{
431 CHK_DIM_2(num_row(),hm2.num_row(),1,hm2.num_col(),-=);
432 SIMPLE_BOP(-=)
433 return (*this);
434}

◆ operator-=() [2/2]

HepVector & CLHEP::HepVector::operator-= ( const HepVector v2)

Definition at line 436 of file Vector.cc.

437{
438 CHK_DIM_1(num_row(),hm2.num_row(),-=);
439 SIMPLE_BOP(-=)
440 return (*this);
441}

◆ operator/=()

HepVector & CLHEP::HepVector::operator/= ( double  t)

Definition at line 443 of file Vector.cc.

444{
445 SIMPLE_UOP(/=)
446 return (*this);
447}

◆ operator=() [1/3]

HepVector & CLHEP::HepVector::operator= ( const Hep3Vector v)

Definition at line 493 of file Vector.cc.

494{
495 if(nrow != 3)
496 {
497 nrow = 3;
498 m.resize(nrow);
499 }
500 m[0] = v.x();
501 m[1] = v.y();
502 m[2] = v.z();
503 return (*this);
504}

◆ operator=() [2/3]

HepVector & CLHEP::HepVector::operator= ( const HepMatrix hm1)

Definition at line 479 of file Vector.cc.

480{
481 if (hm1.num_col() != 1)
482 error("Vector::operator=(Matrix) : Matrix is not Nx1");
483
484 if(hm1.nrow != nrow)
485 {
486 nrow = hm1.nrow;
487 m.resize(nrow);
488 }
489 m = hm1.m;
490 return (*this);
491}

◆ operator=() [3/3]

HepVector & CLHEP::HepVector::operator= ( const HepVector hm2)

Definition at line 468 of file Vector.cc.

469{
470 if(hm1.nrow != nrow)
471 {
472 nrow = hm1.nrow;
473 m.resize(nrow);
474 }
475 m = hm1.m;
476 return (*this);
477}

◆ operator[]() [1/2]

double & CLHEP::HepVector::operator[] ( int  row)
inline

◆ operator[]() [2/2]

const double & CLHEP::HepVector::operator[] ( int  row) const
inline

◆ sub() [1/3]

HepVector CLHEP::HepVector::sub ( int  min_row,
int  max_row 
)

Definition at line 167 of file Vector.cc.

168{
169 HepVector vret(max_row-min_row+1);
170 if(max_row > num_row())
171 error("HepVector::sub: Index out of range");
172 HepGenMatrix::mIter a = vret.m.begin();
173 HepGenMatrix::mIter b = m.begin() + min_row - 1;
174 HepGenMatrix::mIter e = vret.m.begin() + vret.num_row();
175 for(;a<e;) *(a++) = *(b++);
176 return vret;
177}

◆ sub() [2/3]

HepVector CLHEP::HepVector::sub ( int  min_row,
int  max_row 
) const

Definition at line 150 of file Vector.cc.

153{
154#else
155{
156 HepVector vret(max_row-min_row+1);
157#endif
158 if(max_row > num_row())
159 error("HepVector::sub: Index out of range");
160 HepGenMatrix::mIter a = vret.m.begin();
161 HepGenMatrix::mcIter b = m.begin() + min_row - 1;
162 HepGenMatrix::mIter e = vret.m.begin() + vret.num_row();
163 for(;a<e;) *(a++) = *(b++);
164 return vret;
165}

Referenced by main(), and vector_test().

◆ sub() [3/3]

void CLHEP::HepVector::sub ( int  row,
const HepVector v1 
)

Definition at line 179 of file Vector.cc.

180{
181 if(row <1 || row+v1.num_row()-1 > num_row())
182 error("HepVector::sub: Index out of range");
183 HepGenMatrix::mcIter a = v1.m.begin();
184 HepGenMatrix::mIter b = m.begin() + row - 1;
185 HepGenMatrix::mcIter e = v1.m.begin() + v1.num_row();
186 for(;a<e;) *(b++) = *(a++);
187}

◆ T()

HepMatrix CLHEP::HepVector::T ( ) const

Definition at line 530 of file Vector.cc.

533{
534#else
535{
536 HepMatrix mret(1,num_row());
537#endif
538 mret.m = m;
539 return mret;
540}
friend class HepMatrix
Definition: Vector.h:134

Referenced by main(), and symmatrix_test().

Friends And Related Function Documentation

◆ back_solve

void back_solve ( const HepMatrix R,
HepVector b 
)
friend

Definition at line 63 of file MatrixLinear.cc.

64{
65 (*b)(b->num_row()) /= R(b->num_row(),b->num_row());
66 int n = R.num_col();
67 int nb = b->num_row();
68 HepMatrix::mIter br = b->m.begin() + b->num_row() - 2;
69 HepMatrix::mcIter Rrr = R.m.begin() + (nb-2) * (n+1);
70 for (int r=b->num_row()-1;r>=1;--r) {
71 HepMatrix::mIter bc = br+1;
72 HepMatrix::mcIter Rrc = Rrr+1;
73 for (int c=r+1;c<=b->num_row();c++) {
74 (*br)-=(*(Rrc++))*(*(bc++));
75 }
76 (*br)/=(*Rrr);
77 if(r>1) {
78 br--;
79 Rrr -= n+1;
80 }
81 }
82}

◆ col_house

void col_house ( HepMatrix a,
const HepMatrix v,
double  vnormsq,
int  row,
int  col,
int  row_start,
int  col_start 
)
friend

Definition at line 154 of file MatrixLinear.cc.

155 {
156 double beta=-2/vnormsq;
157
158 // This is a fast way of calculating w=beta*A.sub(row,n,col,n)*v
159
160 HepVector w(a->num_col()-col+1,0);
161/* not tested */
162 HepMatrix::mIter wptr = w.m.begin();
163 int n = a->num_col();
164 int nv = v.num_col();
165 HepMatrix::mIter acrb = a->m.begin() + (col-1) * n + (row-1);
166 int c;
167 for (c=col;c<=a->num_col();c++) {
168 HepMatrix::mcIter vp = v.m.begin() + (row_start-1) * nv + (col_start-1);
169 HepMatrix::mcIter acr = acrb;
170 for (int r=row;r<=a->num_row();r++) {
171 (*wptr)+=(*(acr++))*(*vp);
172 vp += nv;
173 }
174 wptr++;
175 if(c<a->num_col()) acrb += n;
176 }
177 w*=beta;
178
179 // Fast way of calculating A.sub=A.sub+w*v.T()
180
181 HepMatrix::mIter arcb = a->m.begin() + (row-1) * n + col-1;
182 wptr = w.m.begin();
183 for (int r=row; r<=a->num_row();r++) {
184 HepMatrix::mIter arc = arcb;
185 HepMatrix::mcIter vp = v.m.begin() + (row_start-1) * nv + col_start;
186 for (c=col;c<=a->num_col();c++) {
187 (*(arc++))+=(*vp)*(*wptr);
188 vp += nv;
189 }
190 wptr++;
191 if(r<a->num_row()) arcb += n;
192 }
193}
virtual int num_col() const
Definition: Vector.cc:118

◆ dot

double dot ( const HepVector v1,
const HepVector v2 
)
friend

Definition at line 542 of file Vector.cc.

543{
544 if(v1.num_row()!=v2.num_row())
545 HepGenMatrix::error("v1 and v2 need to be the same size in dot(HepVector, HepVector)");
546 double d= 0;
547 HepGenMatrix::mcIter a = v1.m.begin();
548 HepGenMatrix::mcIter b = v2.m.begin();
549 HepGenMatrix::mcIter e = a + v1.num_size();
550 for(;a<e;) d += (*(a++)) * (*(b++));
551 return d;
552}

◆ HepDiagMatrix

friend class HepDiagMatrix
friend

Definition at line 132 of file Vector.h.

◆ HepMatrix

friend class HepMatrix
friend

Definition at line 134 of file Vector.h.

◆ HepSymMatrix

friend class HepSymMatrix
friend

Definition at line 133 of file Vector.h.

◆ house [1/2]

HepVector house ( const HepMatrix a,
int  row = 1,
int  col = 1 
)
friend

Definition at line 371 of file MatrixLinear.cc.

372{
373 HepVector v(a.num_row()-row+1);
374/* not tested */
375 int n = a.num_col();
376 HepMatrix::mcIter aic = a.m.begin() + (row - 1) * n + (col - 1) ;
377 HepMatrix::mIter vp = v.m.begin();
378 for (int i=row;i<=a.num_row();i++) {
379 (*(vp++))=(*aic);
380 aic += n;
381 }
382 v(1)+=sign(a(row,col))*v.norm();
383 return v;
384}

◆ house [2/2]

HepVector house ( const HepSymMatrix a,
int  row = 1,
int  col = 1 
)
friend

Definition at line 353 of file MatrixLinear.cc.

354{
355 HepVector v(a.num_row()-row+1);
356/* not tested */
357 HepMatrix::mIter vp = v.m.begin();
358 HepMatrix::mcIter aci = a.m.begin() + col * (col - 1) / 2 + row - 1;
359 int i;
360 for (i=row;i<=col;i++) {
361 (*(vp++))=(*(aci++));
362 }
363 for (;i<=a.num_row();i++) {
364 (*(vp++))=(*aci);
365 aci += i;
366 }
367 v(1)+=sign(a(row,col))*v.norm();
368 return v;
369}

◆ house_with_update

void house_with_update ( HepMatrix a,
int  row = 1,
int  col = 1 
)
friend

Definition at line 396 of file MatrixLinear.cc.

397{
398 HepVector v(a->num_row()-row+1);
399/* not tested */
400 HepMatrix::mIter vp = v.m.begin();
401 int n = a->num_col();
402 HepMatrix::mIter arc = a->m.begin() + (row-1) * n + (col-1);
403 int r;
404 for (r=row;r<=a->num_row();r++) {
405 (*(vp++))=(*arc);
406 if(r<a->num_row()) arc += n;
407 }
408 double normsq=v.normsq();
409 double norm=sqrt(normsq);
410 normsq-=v(1)*v(1);
411 v(1)+=sign((*a)(row,col))*norm;
412 normsq+=v(1)*v(1);
413 (*a)(row,col)=-sign((*a)(row,col))*norm;
414 if (row<a->num_row()) {
415 arc = a->m.begin() + row * n + (col-1);
416 for (r=row+1;r<=a->num_row();r++) {
417 (*arc)=0;
418 if(r<a->num_row()) arc += n;
419 }
420 row_house(a,v,normsq,row,col+1);
421 }
422}
double norm() const
friend void row_house(HepMatrix *, const HepMatrix &, double, int, int, int, int)
double normsq() const

◆ operator* [1/4]

HepVector operator* ( const HepDiagMatrix hm1,
const HepVector hm2 
)
friend

Definition at line 430 of file DiagMatrix.cc.

433{
434#else
435{
436 HepVector mret(hm1.num_row());
437#endif
438 CHK_DIM_1(hm1.num_col(),hm2.num_row(),*);
439 HepGenMatrix::mIter mir=mret.m.begin();
440 HepGenMatrix::mcIter mi1 = hm1.m.begin(), mi2 = hm2.m.begin();
441 for(int icol=1;icol<=hm1.num_col();icol++) {
442 *(mir++) = *(mi1++) * *(mi2++);
443 }
444 return mret;
445}

◆ operator* [2/4]

HepVector operator* ( const HepMatrix hm1,
const HepVector hm2 
)
friend

Definition at line 354 of file Vector.cc.

357{
358#else
359{
360 HepVector mret(hm1.num_row());
361#endif
362 CHK_DIM_1(hm1.num_col(),hm2.num_row(),*);
363 HepGenMatrix::mcIter hm1p,hm2p,vp;
365 double temp;
366 m3p=mret.m.begin();
367 for(hm1p=hm1.m.begin();hm1p<hm1.m.begin()+hm1.num_row()*hm1.num_col();hm1p=hm2p)
368 {
369 temp=0;
370 vp=hm2.m.begin();
371 hm2p=hm1p;
372 while(hm2p<hm1p+hm1.num_col())
373 temp+=(*(hm2p++))*(*(vp++));
374 *(m3p++)=temp;
375 }
376 return mret;
377}

◆ operator* [3/4]

HepVector operator* ( const HepSymMatrix hm1,
const HepVector hm2 
)
friend

Definition at line 507 of file SymMatrix.cc.

510{
511#else
512{
513 HepVector mret(hm1.num_row());
514#endif
515 CHK_DIM_1(hm1.num_col(),hm2.num_row(),*);
516 HepMatrix::mcIter sp,snp,vpt;
517 double temp;
518 int step,stept;
519 HepMatrix::mIter vrp=mret.m.begin();
520 for(step=1,snp=hm1.m.begin();step<=hm1.num_row();++step)
521 {
522 sp=snp;
523 vpt=hm2.m.begin();
524 snp+=step;
525 temp=0;
526 while(sp<snp)
527 temp+=*(sp++)*(*(vpt++));
528 if(step<hm1.num_row()) sp+=step-1;
529 for(stept=step+1;stept<=hm1.num_row();stept++)
530 {
531 temp+=*sp*(*(vpt++));
532 if(stept<hm1.num_row()) sp+=stept;
533 }
534 *(vrp++)=temp;
535 } // for(step
536 return mret;
537}

◆ operator* [4/4]

HepMatrix operator* ( const HepVector hm1,
const HepMatrix hm2 
)
friend

Definition at line 379 of file Vector.cc.

382{
383#else
384{
385 HepMatrix mret(hm1.num_row(),hm2.num_col());
386#endif
387 CHK_DIM_1(1,hm2.num_row(),*);
390 HepMatrix::mIter mrp=mret.m.begin();
391 for(hm1p=hm1.m.begin();hm1p<hm1.m.begin()+hm1.num_row();hm1p++)
392 for(hm2p=hm2.m.begin();hm2p<hm2.m.begin()+hm2.num_col();hm2p++)
393 *(mrp++)=*hm1p*(*hm2p);
394 return mret;
395}

◆ operator+

HepVector operator+ ( const HepVector v1,
const HepVector v2 
)
friend

Definition at line 255 of file Vector.cc.

258{
259#else
260{
261 HepVector mret(hm1.num_row());
262#endif
263 CHK_DIM_1(hm1.num_row(),hm2.num_row(),+);
264 SIMPLE_TOP(+)
265 return mret;
266}
#define SIMPLE_TOP(OPER)
Definition: DiagMatrix.cc:37

◆ operator-

HepVector operator- ( const HepVector v1,
const HepVector v2 
)
friend

Definition at line 299 of file Vector.cc.

302{
303#else
304{
305 HepVector mret(hm1.num_row());
306#endif
307 CHK_DIM_1(hm1.num_row(),hm2.num_row(),-);
308 SIMPLE_TOP(-)
309 return mret;
310}

◆ qr_solve

HepVector qr_solve ( HepMatrix A,
const HepVector b 
)
friend

Definition at line 710 of file MatrixLinear.cc.

711{
712 HepMatrix Q=qr_decomp(A);
713 // Quick way to to Q.T*b.
714 HepVector b2(Q.num_col(),0);
715 HepMatrix::mIter b2r = b2.m.begin();
716 HepMatrix::mIter Qr = Q.m.begin();
717 int n = Q.num_col();
718 for (int r=1;r<=b2.num_row();r++) {
719 HepMatrix::mcIter bc = b.m.begin();
720 HepMatrix::mIter Qcr = Qr;
721 for (int c=1;c<=b.num_row();c++) {
722 *b2r += (*Qcr) * (*(bc++));
723 if(c<b.num_row()) Qcr += n;
724 }
725 b2r++;
726 Qr++;
727 }
728 back_solve(*A,&b2);
729 return b2;
730}
friend void back_solve(const HepMatrix &R, HepVector *b)
Definition: MatrixLinear.cc:63
void qr_decomp(HepMatrix *A, HepMatrix *hsm)

◆ row_house [1/2]

void row_house ( HepMatrix a,
const HepMatrix v,
double  vnormsq,
int  row,
int  col,
int  row_start,
int  col_start 
)
friend

Definition at line 652 of file MatrixLinear.cc.

653 {
654 double beta=-2/vnormsq;
655
656 // This is a fast way of calculating w=beta*A.sub(row,n,col,n).T()*v
657 HepVector w(a->num_col()-col+1,0);
658 int na = a->num_col();
659 int nv = v.num_col();
660 HepMatrix::mIter wptr = w.m.begin();
661 HepMatrix::mIter arcb = a->m.begin() + (row-1) * na + (col-1);
662 HepMatrix::mcIter vpcb = v.m.begin() + (row_start-1) * nv + (col_start-1);
663 int c;
664 for (c=col;c<=a->num_col();c++) {
665 HepMatrix::mIter arc = arcb;
666 HepMatrix::mcIter vpc = vpcb;
667 for (int r=row;r<=a->num_row();r++) {
668 (*wptr)+=(*arc)*(*vpc);
669 if(r<a->num_row()) {
670 arc += na;
671 vpc += nv;
672 }
673 }
674 wptr++;
675 arcb++;
676 }
677 w*=beta;
678
679 arcb = a->m.begin() + (row-1) * na + (col-1);
680 HepMatrix::mcIter vpc = v.m.begin() + (row_start-1) * nv + (col_start-1);
681 for (int r=row; r<=a->num_row();r++) {
682 HepMatrix::mIter arc = arcb;
683 HepMatrix::mIter wptr2 = w.m.begin();
684 for (c=col;c<=a->num_col();c++) {
685 (*(arc++))+=(*vpc)*(*(wptr2++));
686 }
687 if(r<a->num_row()) {
688 arcb += na;
689 vpc += nv;
690 }
691 }
692}

◆ row_house [2/2]

void row_house ( HepMatrix a,
const HepVector v,
double  vnormsq,
int  row = 1,
int  col = 1 
)
friend

Definition at line 613 of file MatrixLinear.cc.

614 {
615 double beta=-2/vnormsq;
616
617 // This is a fast way of calculating w=beta*A.sub(row,n,col,n).T()*v
618
619 HepVector w(a->num_col()-col+1,0);
620/* not tested */
621 int na = a->num_col();
622 HepMatrix::mIter wptr = w.m.begin();
623 HepMatrix::mIter arcb = a->m.begin() + (row-1) * na + (col-1);
624 int c;
625 for (c=col;c<=a->num_col();c++) {
626 HepMatrix::mcIter vp = v.m.begin();
627 HepMatrix::mIter arc = arcb;
628 for (int r=row;r<=a->num_row();r++) {
629 (*wptr)+=(*arc)*(*(vp++));
630 if(r<a->num_row()) arc += na;
631 }
632 wptr++;
633 arcb++;
634 }
635 w*=beta;
636
637 // Fast way of calculating A.sub=A.sub+v*w.T()
638
639 arcb = a->m.begin() + (row-1) * na + (col-1);
640 HepMatrix::mcIter vp = v.m.begin();
641 for (int r=row; r<=a->num_row();r++) {
642 HepMatrix::mIter wptr2 = w.m.begin();
643 HepMatrix::mIter arc = arcb;
644 for (c=col;c<=a->num_col();c++) {
645 (*(arc++))+=(*vp)*(*(wptr2++));
646 }
647 vp++;
648 if(r<a->num_row()) arcb += na;
649 }
650}

◆ solve

HepVector solve ( const HepMatrix a,
const HepVector v 
)
friend

Definition at line 575 of file Vector.cc.

578{
579#else
580{
581 HepVector vret(v);
582#endif
583 static CLHEP_THREAD_LOCAL int max_array = 20;
584 static CLHEP_THREAD_LOCAL int *ir = new int [max_array+1];
585
586 if(a.ncol != a.nrow)
587 HepGenMatrix::error("Matrix::solve Matrix is not NxN");
588 if(a.ncol != v.nrow)
589 HepGenMatrix::error("Matrix::solve Vector has wrong number of rows");
590
591 int n = a.ncol;
592 if (n > max_array) {
593 delete [] ir;
594 max_array = n;
595 ir = new int [max_array+1];
596 }
597 double det;
598 HepMatrix mt(a);
599 int i = mt.dfact_matrix(det, ir);
600 if (i!=0) {
601 for (i=1;i<=n;i++) vret(i) = 0;
602 return vret;
603 }
604 double s21, s22;
605 int nxch = ir[n];
606 if (nxch!=0) {
607 for (int hmm=1;hmm<=nxch;hmm++) {
608 int ij = ir[hmm];
609 i = ij >> 12;
610 int j = ij%4096;
611 double te = vret(i);
612 vret(i) = vret(j);
613 vret(j) = te;
614 }
615 }
616 vret(1) = mt(1,1) * vret(1);
617 if (n!=1) {
618 for (i=2;i<=n;i++) {
619 s21 = -vret(i);
620 for (int j=1;j<i;j++) {
621 s21 += mt(i,j) * vret(j);
622 }
623 vret(i) = -mt(i,i)*s21;
624 }
625 for (i=1;i<n;i++) {
626 int nmi = n-i;
627 s22 = -vret(nmi);
628 for (int j=1;j<=i;j++) {
629 s22 += mt(nmi,n-j+1) * vret(n-j+1);
630 }
631 vret(nmi) = -s22;
632 }
633 }
634 return vret;
635}

◆ swap

void swap ( HepVector v1,
HepVector v2 
)
friend

◆ tridiagonal

void tridiagonal ( HepSymMatrix a,
HepMatrix hsm 
)
friend

Definition at line 777 of file MatrixLinear.cc.

778{
779 int nh = hsm->num_col();
780 for (int k=1;k<=a->num_col()-2;k++) {
781
782 // If this row is already in tridiagonal form, we can skip the
783 // transformation.
784
785 double scale=0;
786 HepMatrix::mIter ajk = a->m.begin() + k * (k+5) / 2;
787 int j;
788 for (j=k+2;j<=a->num_row();j++) {
789 scale+=fabs(*ajk);
790 if(j<a->num_row()) ajk += j;
791 }
792 if (scale ==0) {
793 HepMatrix::mIter hsmjkp = hsm->m.begin() + k * (nh+1) - 1;
794 for (j=k+1;j<=hsm->num_row();j++) {
795 *hsmjkp = 0;
796 if(j<hsm->num_row()) hsmjkp += nh;
797 }
798 } else {
799 house_with_update2(a,hsm,k+1,k);
800 double normsq=0;
801 HepMatrix::mIter hsmrptrkp = hsm->m.begin() + k * (nh+1) - 1;
802 int rptr;
803 for (rptr=k+1;rptr<=hsm->num_row();rptr++) {
804 normsq+=(*hsmrptrkp)*(*hsmrptrkp);
805 if(rptr<hsm->num_row()) hsmrptrkp += nh;
806 }
807 HepVector p(a->num_row()-k,0);
808 rptr=k+1;
809 HepMatrix::mIter pr = p.m.begin();
810 int r;
811 for (r=1;r<=p.num_row();r++) {
812 HepMatrix::mIter hsmcptrkp = hsm->m.begin() + k * (nh+1) - 1;
813 int cptr;
814 for (cptr=k+1;cptr<=rptr;cptr++) {
815 (*pr)+=a->fast(rptr,cptr)*(*hsmcptrkp);
816 if(cptr<a->num_col()) hsmcptrkp += nh;
817 }
818 for (;cptr<=a->num_col();cptr++) {
819 (*pr)+=a->fast(cptr,rptr)*(*hsmcptrkp);
820 if(cptr<a->num_col()) hsmcptrkp += nh;
821 }
822 (*pr)*=2/normsq;
823 rptr++;
824 pr++;
825 }
826 double pdotv=0;
827 pr = p.m.begin();
828 hsmrptrkp = hsm->m.begin() + k * (nh+1) - 1;
829 for (r=1;r<=p.num_row();r++) {
830 pdotv+=(*(pr++))*(*hsmrptrkp);
831 if(r<p.num_row()) hsmrptrkp += nh;
832 }
833 pr = p.m.begin();
834 hsmrptrkp = hsm->m.begin() + k * (nh+1) - 1;
835 for (r=1;r<=p.num_row();r++) {
836 (*(pr++))-=pdotv*(*hsmrptrkp)/normsq;
837 if(r<p.num_row()) hsmrptrkp += nh;
838 }
839 rptr=k+1;
840 pr = p.m.begin();
841 hsmrptrkp = hsm->m.begin() + k * (nh+1) - 1;
842 for (r=1;r<=p.num_row();r++) {
843 int cptr=k+1;
844 HepMatrix::mIter mpc = p.m.begin();
845 HepMatrix::mIter hsmcptrkp = hsm->m.begin() + k * (nh+1) - 1;
846 for (int c=1;c<=r;c++) {
847 a->fast(rptr,cptr)-= (*hsmrptrkp)*(*(mpc++))+(*pr)*(*hsmcptrkp);
848 cptr++;
849 if(c<r) hsmcptrkp += nh;
850 }
851 pr++;
852 rptr++;
853 if(r<p.num_row()) hsmrptrkp += nh;
854 }
855 }
856 }
857}
void house_with_update2(HepSymMatrix *a, HepMatrix *v, int row=1, int col=1)

◆ vT_times_v

HepSymMatrix vT_times_v ( const HepVector v)
friend

Definition at line 539 of file SymMatrix.cc.

542{
543#else
544{
545 HepSymMatrix mret(v.num_row());
546#endif
547 HepMatrix::mIter mr=mret.m.begin();
548 HepMatrix::mcIter vt1,vt2;
549 for(vt1=v.m.begin();vt1<v.m.begin()+v.num_row();vt1++)
550 for(vt2=v.m.begin();vt2<=vt1;vt2++)
551 *(mr++)=(*vt1)*(*vt2);
552 return mret;
553}
friend class HepSymMatrix
Definition: Vector.h:133

The documentation for this class was generated from the following files: