12#include "CLHEP/Matrix/Matrix.h"
13#include "CLHEP/Matrix/Vector.h"
14#include "CLHEP/Matrix/SymMatrix.h"
21static int sign(
double x) {
return (x>0 ? 1: -1);}
70 for (
int r=b->
num_row()-1;r>=1;--r) {
73 for (
int c=r+1;c<=b->
num_row();c++) {
74 (*br)-=(*(Rrc++))*(*(bc++));
96 for (
int i=1;i<=b->
num_col();i++) {
100 for (
int r=b->
num_row()-1;r>=1;--r) {
103 for (
int c=r+1;c<=b->
num_row();c++) {
104 (*bri)-= (*(Rrc++)) * (*bci);
105 if(c<b->num_row()) bci += nc;
125 int k1,
int k2,
int row_min,
int row_max) {
126 if (row_max<=0) row_max = A->
num_row();
130 for (
int j=row_min;j<=row_max;j++) {
131 double tau1=(*Ajk1);
double tau2=(*Ajk2);
132 (*Ajk1)=c*tau1-ds*tau2;(*Ajk2)=ds*tau1+c*tau2;
155 int row,
int col,
int row_start,
int col_start) {
156 double beta=-2/vnormsq;
167 for (c=col;c<=a->
num_col();c++) {
170 for (
int r=row;r<=a->
num_row();r++) {
171 (*wptr)+=(*(acr++))*(*vp);
175 if(c<a->num_col()) acrb += n;
183 for (
int r=row; r<=a->
num_row();r++) {
186 for (c=col;c<=a->
num_col();c++) {
187 (*(arc++))+=(*vp)*(*wptr);
191 if(r<a->num_row()) arcb += n;
206 max=min=fabs(mcopy(1,1));
210 for (
int i=2; i<=n; i++) {
211 if (max<fabs(*mii)) max=fabs(*mii);
212 if (min>fabs(*mii)) min=fabs(*mii);
229 double d=(t->
fast(end-1,end-1)-t->
fast(end,end))/2;
230 double mu=t->
fast(end,end)-t->
fast(end,end-1)*t->
fast(end,end-1)/
231 (d+sign(d)*sqrt(d*d+ t->
fast(end,end-1)*t->
fast(end,end-1)));
232 double x=t->
fast(begin,begin)-mu;
233 double z=t->
fast(begin+1,begin);
237 for (
int k=begin;k<=end-1;k++) {
247 *(tkk-1)= *(tkk-1)*c-(*(tkp1k-1))*ds;
252 double aq=(*tkp1k+1);
253 (*tkk)=ap*c*c-2*c*bp*ds+aq*ds*ds;
254 (*tkp1k)=c*ap*ds+bp*c*c-bp*ds*ds-ds*aq*c;
255 (*(tkp1k+1))=ap*ds*ds+2*c*bp*ds+aq*c*c;
258 double bq=(*(tkp2k+1));
266 if(k<end-2) tkp2k += k+3;
272 double d=(t->
fast(end-1,end-1)-t->
fast(end,end))/2;
273 double mu=t->
fast(end,end)-t->
fast(end,end-1)*t->
fast(end,end-1)/
274 (d+sign(d)*sqrt(d*d+ t->
fast(end,end-1)*t->
fast(end,end-1)));
275 double x=t->
fast(begin,begin)-mu;
276 double z=t->
fast(begin+1,begin);
280 for (
int k=begin;k<=end-1;k++) {
290 *(tkk-1)= (*(tkk-1))*c-(*(tkp1k-1))*ds;
295 double aq=(*(tkp1k+1));
296 (*tkk)=ap*c*c-2*c*bp*ds+aq*ds*ds;
297 (*tkp1k)=c*ap*ds+bp*c*c-bp*ds*ds-ds*aq*c;
298 (*(tkp1k+1))=ap*ds*ds+2*c*bp*ds+aq*c*c;
300 double bq=(*(tkp2k+1));
308 if(k<end-2) tkp2k += k+3;
320 const double tolerance = 1e-12;
328 for (
int i=begin;i<=end-1;i++) {
330 tolerance*(fabs(*sii)+fabs(*(sip1i+1)))) {
338 while(begin<end && hms->fast(begin+1,begin) ==0) begin++;
339 while(end>begin && hms->
fast(end,end-1) ==0) end--;
360 for (i=row;i<=col;i++) {
361 (*(vp++))=(*(aci++));
367 v(1)+=sign(a(row,col))*v.
norm();
378 for (
int i=row;i<=a.
num_row();i++) {
382 v(1)+=sign(a(row,col))*v.
norm();
404 for (r=row;r<=a->
num_row();r++) {
406 if(r<a->num_row()) arc += n;
409 double norm=sqrt(normsq);
411 v(1)+=sign((*a)(row,col))*
norm;
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++) {
418 if(r<a->num_row()) arc += n;
432 for (r=row;r<=a->
num_row();r++) {
434 normsq+=(*vrc)*(*vrc);
440 double norm=sqrt(normsq);
441 vrc = v->m.begin() + (row-1) * nv + (col-1);
442 normsq-=(*vrc)*(*vrc);
443 (*vrc)+=sign((*a)(row,col))*
norm;
444 normsq+=(*vrc)*(*vrc);
445 (*a)(row,col)=-sign((*a)(row,col))*
norm;
446 if (row<a->num_row()) {
447 arc = a->m.begin() + row * na + (col-1);
448 for (r=row+1;r<=a->
num_row();r++) {
450 if(r<a->num_row()) arc += na;
452 row_house(a,*v,normsq,row,col+1,row,col);
469 for (r=row;r<=a->
num_row();r++)
472 normsq+=(*vrc)*(*vrc);
478 double norm=sqrt(normsq);
479 vrc = v->m.begin() + (row-1) * nv + (col-1);
480 arc = a->m.begin() + (row-1) * row / 2 + (col-1);
481 (*vrc)+=sign(*arc)*
norm;
482 (*arc)=-sign(*arc)*
norm;
484 for (r=row+1;r<=a->
num_row();r++) {
486 if(r<a->num_row()) arc += r;
535 const double small = 1e-10;
539 for (
int i=0;i<n;i++)
541 if (fabs(t=A[i].normsq())<small) {
546 D +=
dot(A[i],
B[i])*(1-2/t)*A[i]+
B[i];
571 for (
int j=hsm.
num_col();j>=1;--j)
588 int k1,
int k2,
int col_min,
int col_max) {
590 if (col_max==0) col_max = A->
num_col();
594 for (
int j=col_min;j<=col_max;j++) {
595 double tau1=(*Ak1j);
double tau2=(*Ak2j);
596 (*(Ak1j++))=c*tau1-ds*tau2;(*(Ak2j++))=ds*tau1+c*tau2;
615 double beta=-2/vnormsq;
625 for (c=col;c<=a->
num_col();c++) {
628 for (
int r=row;r<=a->
num_row();r++) {
629 (*wptr)+=(*arc)*(*(vp++));
630 if(r<a->num_row()) arc += na;
639 arcb = a->m.begin() + (row-1) * na + (col-1);
641 for (
int r=row; r<=a->
num_row();r++) {
644 for (c=col;c<=a->
num_col();c++) {
645 (*(arc++))+=(*vp)*(*(wptr2++));
648 if(r<a->num_row()) arcb += na;
653 int row,
int col,
int row_start,
int col_start) {
654 double beta=-2/vnormsq;
664 for (c=col;c<=a->
num_col();c++) {
667 for (
int r=row;r<=a->
num_row();r++) {
668 (*wptr)+=(*arc)*(*vpc);
679 arcb = a->m.begin() + (row-1) * na + (col-1);
681 for (
int r=row; r<=a->
num_row();r++) {
684 for (c=col;c<=a->
num_col();c++) {
685 (*(arc++))+=(*vpc)*(*(wptr2++));
718 for (
int r=1;r<=b2.
num_row();r++) {
721 for (
int c=1;c<=b.
num_row();c++) {
722 *b2r += (*Qcr) * (*(bc++));
747 for (
int i=1;i<=b.
num_col();i++) {
750 for (
int r=1;r<=b2.
num_row();r++) {
753 for (
int c=1;c<=b.
num_row();c++) {
754 *b2ri += (*Qcr) * (*bci);
780 for (
int k=1;k<=a->
num_col()-2;k++) {
788 for (j=k+2;j<=a->
num_row();j++) {
790 if(j<a->num_row()) ajk += j;
794 for (j=k+1;j<=hsm->
num_row();j++) {
796 if(j<hsm->num_row()) hsmjkp += nh;
803 for (rptr=k+1;rptr<=hsm->
num_row();rptr++) {
804 normsq+=(*hsmrptrkp)*(*hsmrptrkp);
805 if(rptr<hsm->num_row()) hsmrptrkp += nh;
814 for (cptr=k+1;cptr<=rptr;cptr++) {
815 (*pr)+=a->
fast(rptr,cptr)*(*hsmcptrkp);
816 if(cptr<a->num_col()) hsmcptrkp += nh;
818 for (;cptr<=a->
num_col();cptr++) {
819 (*pr)+=a->
fast(cptr,rptr)*(*hsmcptrkp);
820 if(cptr<a->num_col()) hsmcptrkp += nh;
828 hsmrptrkp = hsm->m.begin() + k * (nh+1) - 1;
830 pdotv+=(*(pr++))*(*hsmrptrkp);
831 if(r<p.
num_row()) hsmrptrkp += nh;
834 hsmrptrkp = hsm->m.begin() + k * (nh+1) - 1;
836 (*(pr++))-=pdotv*(*hsmrptrkp)/normsq;
837 if(r<p.
num_row()) hsmrptrkp += nh;
841 hsmrptrkp = hsm->m.begin() + k * (nh+1) - 1;
846 for (
int c=1;c<=r;c++) {
847 a->
fast(rptr,cptr)-= (*hsmrptrkp)*(*(mpc++))+(*pr)*(*hsmcptrkp);
849 if(c<r) hsmcptrkp += nh;
853 if(r<p.
num_row()) hsmrptrkp += nh;
866 for (
int j=hsm.
num_col();j>=1;--j) {
874 int row_start,
int col_start)
877 for (
int i=row_start;i<=row_start+a->
num_row()-row;i++)
878 normsq+=v(i,col)*v(i,col);
879 col_house(a,v,normsq,row,col,row_start,col_start);
882void givens(
double a,
double b,
double *c,
double *ds)
884 if (b ==0) { *c=1; *ds=0; }
886 if (fabs(b)>fabs(a)) {
888 *ds=1/sqrt(1+tau*tau);
892 *c=1/sqrt(1+tau*tau);
900 for (
int i=1;i<=A->
num_col();i++)
905 int row_start,
int col_start)
908 int end = row_start+a->
num_row()-row;
909 for (
int i=row_start; i<=end; i++)
910 normsq += v(i,col)*v(i,col);
913 row_house(a,v,normsq,row,col,row_start,col_start);
std::vector< double, Alloc< double, 25 > >::const_iterator mcIter
std::vector< double, Alloc< double, 25 > >::iterator mIter
static void error(const char *s)
virtual int num_col() const
virtual int num_row() const
const double & fast(int row, int col) const
virtual int num_row() const
void col_house(HepMatrix *a, const HepMatrix &v, double vnormsq, int row, int col, int row_start, int col_start)
void tridiagonal(HepSymMatrix *a, HepMatrix *hsm)
HepVector qr_solve(const HepMatrix &A, const HepVector &b)
HepSymMatrix vT_times_v(const HepVector &v)
void house_with_update2(HepSymMatrix *a, HepMatrix *v, int row=1, int col=1)
double norm(const HepGenMatrix &m)
HepMatrix qr_inverse(const HepMatrix &A)
void house_with_update(HepMatrix *a, int row=1, int col=1)
void back_solve(const HepMatrix &R, HepVector *b)
void qr_decomp(HepMatrix *A, HepMatrix *hsm)
double dot(const HepVector &v1, const HepVector &v2)
HepVector house(const HepMatrix &a, int row=1, int col=1)
HepVector min_line_dist(const HepVector *const A, const HepVector *const B, int n)
void row_givens(HepMatrix *A, double c, double s, int k1, int k2, int col_min=1, int col_max=0)
void givens(double a, double b, double *c, double *s)
void col_givens(HepMatrix *A, double c, double s, int k1, int k2, int row_min=1, int row_max=0)
void row_house(HepMatrix *a, const HepVector &v, double vnormsq, int row=1, int col=1)
double condition(const HepSymMatrix &m)
HepMatrix diagonalize(HepSymMatrix *s)
void diag_step(HepSymMatrix *t, int begin, int end)