/* LA: linear algebra C++ interface library Copyright (C) 2008 Jiri Pittner or This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program. If not, see . */ // g++ -D _GLIBCPP_NO_TEMPLATE_EXPORT -g testblas.cc testblas2.cc nrutil_modif.cc -L/usr/local/lib/atlas -lstrassen -lf77blas -lcblas -latlas -ltraceback -lbfd -liberty #include #include "la.h" #include "vecmat3.h" #include "quaternion.h" #include "permutation.h" #include "polynomial.h" #include "contfrac.h" #include "simple.h" #include "graph.h" #include "numbers.h" #include "miscfunc.h" using namespace std; using namespace LA_Vecmat3; using namespace LA_Quaternion; using namespace LA_Simple; using namespace LA; extern void test(const NRVec &x); double ad; void f1(const double *c) { ad=*c; } void f2(double *c) { *c=ad; } inline int randind(const int n) { return int(RANDDOUBLE()*n); } complex mycident (const complex&x) {return x;} void printme(const NRPerm &p) { PERM_RANK_TYPE rank=p.rank(); int n=p.size(); cout< qq(n,rank); if(qq!=p) laerror("error in rank algorithm"); for(int i=0; i<4; ++i) { NRVec_from1 inv=p.inversions(i); NRPerm q(i,inv); if(q!=p) laerror("error in inversions algorithm"); } } void printme0(const NRPerm &p) { cout< &p) { cout < Snmtable; static int unitary_n; static PERM_RANK_TYPE space_dim; static NRVec > allyoung; static NRVec >allyoungmat; static NRVec >allyoungregular; static NRVec allyoung_irrep; int current_irrep; int allyoung_index; void yprintme(const YoungTableaux&y) { cout < y(p); PERM_RANK_TYPE dim=y.generate_all_standard(yprintme); Partition q=p.adjoint(); if(dim!=snd) laerror("inconsistency in standard tableaux generation"); PERM_RANK_TYPE und=p.Un_irrep_dim(unitary_n); cout<<"U("< qc(q); cout <<"("< u; u=x+y; cout <<"u= "< aa(0.,3,3); aa[0][0]=aa[1][1]=aa(2,2)=2.; NRMat bb(aa); double *p; aa.copyonwrite(); p= &aa[2][2]; *p=3.; bb.copyonwrite(); bb(0,2)=1.; cout << "aa= " < cc=aa & bb; cout << "aa o+ bb= " << cc <<"\n"; cout << cc.rsum() <<"\n"; cout << cc.csum() <<"\n"; NRVecw(3); w[0]=1; w[1]=2;w[2]=3; NRVec v(0.,3); v.gemv(0.,bb,'n',1.,w); cout << " v= " < bb(1.,n,n); for(int i=0;i amat,bmat,cmat; cin >>amat; cin >>bmat; cmat=amat*bmat; cout< amat(1.,2,2); NRMat bmat(amat); NRMat dmat(amat); //NRMat cmat; cmat=bmat*2.; NRMat cmat(bmat*2); //more efficient dmat.copyonwrite(); dmat[0][0]=0; cout< amat; NRVec avec; cin >>amat; cin >>avec; cout << amat*avec; cout << avec*amat; NRVec avec(0.,10); f1(avec); f2(avec); NRVec uu(3); uu[0]=1; uu[1]=2; uu[2]=3; cout << uu << (uu|uu) <<"\n"; NRSMat sa(0.,3); sa(0,0)=1; sa(0,2)=5; sa(2,2)=10;sa(1,0)=2;sa(1,1)=3; sa(2,1)=-1; NRSMat sb(0.,3); sb(0,0)=-2; sb(0,2)=1; sb(2,2)=2;sb(1,0)=-1;sb(1,1)=7; sb(2,1)=3; cout << "symetr\n" < m10(10.,3,3); cout << "10 + sa" << m10 + sa <<"\n"; */ /* const int dim=256; NRMat big1(dim,dim),big2(dim,dim),big3; for(int i=0;i atest, btest,ctest; { int cc,c1,c2,c3; cin >>cc>>c1>>c2>>c3; atest.s_cutoff(cc,c1,c2,c3); } cin>>atest; cin>>btest; NRMat dtest(atest.nrows(),btest.ncols()); dtest.gemm(0., atest, 't', btest, 'n', 1.); cout << dtest; NRMat etest(atest.nrows(),btest.ncols()); etest.strassen(0., atest, 't', btest, 'n', 1.); cout << etest; } if(0) { int dim; cin >>dim; NRMat big1(dim,dim),big2(dim,dim),big3,big4(dim,dim); for(int i=0;i a(3,3),b; NRVec v(3); for(int i=0;i<3;i++) for(int j=0;j<3;j++) { a[i][j]= i*i+j; v[i]=10-i;} b=a; b*= sin(1.)+1; cout << a < a(3,3),b; NRVec v(10); v[0]=2;v[1]=3;v[2]=1;v[3]=-3;v[4]=2;v[5]=-1;v[6]=3;v[7]=-2;v[8]=1;v[9]=1; for(int i=0;i<3;i++) for(int j=0;j<3;j++) { a[i][j]= (i+j)/10.; } cout < a(3,3); for(int i=0;i<3;i++) for(int j=0;j<3;j++) { a[i][j]= (i+j)/10.; } NRSMat b(a); NRMat c(b); cout < a(3,3); a[0][0]=1; a[0][1]=2;a[0][2]=3; a[1][0]=4; a[1][1]=-5;a[1][2]=7; a[2][0]=-3;a[2][1]=10;a[2][2]=2; NRMat b(2,3); b[0][0]=1;b[0][1]=2;b[0][2]=3; b[1][0]=2;b[1][1]=4;b[1][2]=6; cout < a(3,3); for(int i=0;i<3;i++) for(int j=0;j<3;j++) { a[i][j]= (i+j)/10.; } NRVec b(3); cout < a(3); NRMatv(3,3); for(int i=0;i<3;i++) for(int j=0;j<3;j++) { a(i,j)= (i+j)/10.; } NRVec b(3); cout <c=(NRMat)a; //nebo NRMatc(a); NRMatd=exp(c); diagonalize(a,b,&v); cout < a; cin >>a ; NRMat abak=a; NRMat abak2=a; NRMat u(a.nrows(),a.nrows()),v(a.ncols(),a.ncols()); NRVecs(a.ncols() sdiag(0., u.ncols(),v.nrows()); sdiag.diagonalset(s); cout < ai=calcinverse(abak2); cout <<"regular inverse "<ss(s);ss.copyonwrite(); for(int i=0; i1e-12) ss[i]=1./ss[i]; else ss[i]=0; sdiag.diagonalset(ss); NRMat aii = v.transpose() *sdiag*u.transpose(); cout <<"svd inverse "< > a; cin >>a ; NRMat > abak=a; NRMat > u(a.nrows(),a.nrows()),v(a.ncols(),a.ncols()); NRVecs(a.ncols() > sdiag(0., u.ncols(),v.nrows()); NRVec > ss = s; sdiag.diagonalset(ss); cout < a(aa,3,3); NRMat a; cin >>a; cout < u(n,n),v(n,n); NRVecwr(n),wi(n); gdiagonalize(a,wr,wi,&u,&v,0); cout <z=diagofproduct(u,v,1); for(int i=0;i a; cin >>a; cout < > u(n,n),v(n,n); NRVec >w(n); gdiagonalize(a,w,&u,&v); cout < >z=diagofproduct(u,v,1,1); //NRMat > zz=u*v.transpose(1); cout < a; cin >>a; int n=a.nrows(); NRMat > u(n,n),v(n,n); NRVec >w(n); gdiagonalize(a,w,&u,&v,0,n,0,1); cout < >z=diagofproduct(u,v,1,1); cout < > a; cin >>a; int n=a.nrows(); NRMat > u(n,n),v(n,n); NRVec >w(n); gdiagonalize(a,w,&u,&v); //gdiagonalize(a,w,&u,&v,0,n,0,1); cout < >z=diagofproduct(u,v,1,1); cout < a(4,4); NRVec v(4); v[0]=1;v[1]=2;v[2]=3;v[3]=4; a=1.; a.copyonwrite(); a.add(3,0,.5); a.add(0,2,.2); a.add(2,1,.1); a.add(3,3,1.); a.add(1,1,-1.); SparseMat c(a); c*=10.; cout <b(c); cout < a(4,4),b(4,4); a=1.; a.copyonwrite(); a.add(3,0,.5); b.add(0,2,.2); b.add(2,1,.1); b.add(3,3,1.); b.add(1,1,-1.); SparseMatc=a+b; cout < a(4,4),b(4,4); a=0.; b=2; a.add(3,0,.5); a.add(0,2,.2); a.add(1,1,1); a.add(1,0,.2); b.add(2,1,.1); b.add(3,3,1.); b.add(1,1,-1.); NRMat aa(a),bb(b); SparseMatc; NRMatcc; //cout << NRMat(c); //cout <(c); cout <<"norms2 "< aa(n,n); cout << "\n\n\ntiming for size "< a(0.,n,n); for(int i=0; i b(exp(a)); //cout <(a); } else { for(int i=0; i bb(exp(aa)); //cout <>n; SparseMat aa(n,n); for(int i=0; i bb=exp(aa); NRVec v(n); for(int i=0; i res1=bb*v; NRVec res2=exptimes(aa,v); cout <<"difference = "<<(res1-res2).norm()<>n>>k>>m; { NRMat a(n,k),b(k,m),c(n,m),d(n,m); for(int i=0;i aa(a); d.gemm(0., aa, 'n', b, 'n', .6); cout< v(10.,10); v+= 5.; cout <>n; NRMat a(n,n); for(int i=0;i b; b|=a; NRVec er(n),ei(n); NRMat vr(n,n),vl(n,n); gdiagonalize(b,er,ei,&vl,&vr,1,0,1,1); cout < u=exp(a*.1); gdiagonalize(u,er,ei,&vl,&vr,1,0,1,1); cout <>k; int n=2*k; NRMat a(n,n); //matrix with known spectrum for(int i=0;i er(n),ei(n); NRMat vr(n,n),vl(n,n); cout <<"input matrix\n"<>n; NRMat a(n,n); for(int i=0;i b=exp0(a); cout <<"b=exp(a) matrix\n"< c=log(b); //matrixfunction(a,&mycident,1); cout <<"c=log(exp(a))\n"< e=exp(c); cout <<"e=exp(c)\n"< b=exp0(a); cout <<"exp0 took "< c=exp(a,true); cout <<" horner exp took "<1e-10); } if(0) { int n=3; NRMat a(n,n); a(0,0)=1.; a(0,1)=2.; a(1,0)=2.; a(1,1)=6.; a(2,2)=-4; a(0,2)=1; cout < c=calcinverse(a,&d); cout < a(3,3); NRMat b=a; for(int i=1; i<4;i++) b=b*b; } if(0) { NRMat aa,bb,cc; cin >>aa; cc=copytest(aa); cout < > a,b,c; a=complexify(aa); c=copytest(a); cout < > a,b,c; cin>>a; c=copytest(a); cout < a,b,c; cin >>a; c=copytest(a); cout < a; NRMat b=exp(a); NRMat c=log(b); cout < a; cin >>a; NRMat c=log(a); //matrixfunction(a,&mycident,1); cout < b=exp(c); cout <<"exp(log(x))\n"< a; cin >>a; NRMat aa(a); NRMat b=exp(aa); NRMat c=matrixfunction(a,&exp); cout < h; NRMat t; cin >>h; cin >>t; NRMat r1= exp(-t) * h * exp(t); NRMat r2=BCHexpansion(h,t,30); cout <>n; SparseMat a(n,n); int spars=n*n/3; for(int i=0; i aa(a); NRVec v(aa[0],n*n); cout < amat,bmat; cin >>amat; cin >>bmat; NRVec v(amat.nrows()); diagonalize(amat,v,1,1,0,&bmat,1); cout <>n >>m; NRSMat a(n,n); NRVec rr(n); for(int i=0;i aa; NRMat vv(n,n); aa=a; diagonalize(aa,rr,&vv); NRVec r(m); NRVec *eivecs = new NRVec[m]; davidson(a,r,eivecs,NULL,m,1,1e-6,1,200); cout <<"Davidson energies " <>n >>m; NRMat a(n,n); NRVec rr(n),ii(n); double tmp=0.; for(int i=0;i aa=a; NRMat vv=aa; gdiagonalize(aa, rr, ii, NULL, &vv, 1, 0, 2, 0, NULL,NULL); NRVec r(m); NRVec *eivecs = new NRVec[m]; davidson(a,r,eivecs,NULL,m,1,1e-6,1,200); cout <<"Davidson energies " <>n>>m; SparseMat aa(n,n); aa.setsymmetric(); for(int i=0; i r(m); davidson(aa,r,(NRVec *)NULL,"eivecs",m,1,1e-5,0,300,300); cout <>n>>m; SparseMat aa(n,n); aa.setsymmetric(); for(int i=0; i r(m); NRVec r2(m); davidson(aa,r,(NRVec *)NULL,"eivecs",m,1,1e-5,1,300,300); SparseMat bb(n,n); for(int i=0; i e1,e2,cc; e1=exp(bb); e2=exp(bb*-1.); aa.setunsymmetric(); cc=e1*aa*e2; davidson(cc,r2,(NRVec *)NULL,"eivecs2",m,1,1e-5,1,300,300); cout <<"original matrix" <>n>>m; SparseMat aa(n,n); for(int i=0; i r(m); davidson(aa,r,(NRVec *)NULL,"eivecs",m,1,1e-5,0,300,300); cout <>n >>m; NRMat a(n,m); NRVec b(n); NRVec x(m); for(int i=0;i bb=b; cout <aa=a; double d; linear_solve(aa,bb,&d); cout <<"det = "<>n >>m; SparseMat aa(n,m); NRVec b(n); NRVec x(m); //tridiagonal for(int i=0; i A(3,3); A=1.; double *p = (double *)A; *p=2.; cout < > diis(5); int dim; cin>>dim; NRVec solution(dim), deviation(dim); for(i=0; i1e-8 ; ++iter) { NRVec trial=solution; trial.copyonwrite(); for(i=0; i a,b; cin >>a; b=realsqrt(a); cout < a; NRMat b; cin >>a>>b; cout < a,b; cin >>a >>b; cout < a; cin >>a; if(a.nrows()!=a.ncols()) laerror("must be square matrix"); int n=a.nrows(); NRMatvl(n,n),vr(n,n); NRVec wr(n),wi(n); NRMat awork(a); gdiagonalize(awork,wr,wi,&vl,&vr,1,0,1,biorthonormalize); for(int i=0; i eival(0.,n,n); eival.diagonalset(wr); cout <<"test biorthonormality "<< (vl.transpose() * vr - 1.).norm()< hlp; hlp = a *vr - vr * eival; cout <<"test right eivectors "< ader(n,n); for(int i=0; i vlx(n,n),vrx(n,n); NRVec wrx(n),wix(n); gdiagonalize(awork,wrx,wix,&vlx,&vrx,1,0,1,biorthonormalize); for(int i=0; i vly(n,n),vry(n,n); NRVec wry(n),wiy(n); gdiagonalize(awork,wry,wiy,&vly,&vry,1,0,1,biorthonormalize); for(int i=0; i vld,vrd; NRVec wrd; vld = (vlx-vly) * (.5/h); vrd = (vrx-vry) * (.5/h); wrd = (wrx-wry) * (.5/h); NRMat vlg,vrg; NRVec wrg(n); //compute analytic derivative NRMat tmp(n,n); tmp.gemm(0.,vl,'t', ader * vr,'n',1.); hlp |= tmp; cout <<" C~+ VH C = \n"< Y = tmp; NRMat S = vr.transpose() * vr; cout <<"test S\n"< tmp2 = S * tmp; cout <<"test tmp2\n"< vri = calcinverse(vr); NRMat numX = vri * vrd; cout <<" numerical X matrix \n"<< numX; cout <<" numerical X matrix test = "<< (vr * numX - vrd).norm()< v(3); v[0]=1; v[1]=2; v[2]=3; NRVec > vv = v; NRVecu(v); vv += u; cout < scale; cin >> scale; NRMat > h; cin >>h; NRVec > x(h.nrows()); NRVec > y,z; x.randomize(1.); y=exptimes(h*scale,x,false,1.); z=exptimes(h,x,false,scale); cout <>nocc>>nvirt; int n=nocc+nvirt; NRMat t(nocc,nvirt); t.randomize(0.5); NRSMat A(n); A=1.; for(int i=0; i B(n,n); NRVec E(n); cout < Awork=A; diagonalize(Awork,E,&B); cout < U = A*B; for(int p=0; p > h(ham); h *= complex(0,1); h.simplify(); cout <<"norms of input (should be same) "< > hexp = exp(h); cout <<"sparse exp time "< >hamexp2(hexp); cout <<"norm of results "< >hamexp3 = complexmatrix(c,s); cout <<"sincos error = "<<(hamexp1-hamexp3).norm()< hreal(hamreal); SparseMat si,co; t=clock()/((double) (CLOCKS_PER_SEC)); sincos(si,co,hreal); cout <<"sparse sincos time "< >hamexp4 = complexmatrix(NRMat(co),NRMat(si)); cout <<"sincos error 2 = "<<(hamexp1-hamexp4).norm()< > rhs(ham.ncols()); rhs.randomize(1.); t=clock()/((double) (CLOCKS_PER_SEC)); NRVec > r1 = exptimes(ham*complex(0,1),rhs); cout <<"dense exptimes "< > r2 = exptimes(h,rhs); cout <<"sparse exptimes "< > r3 = cov + siv*complex(0,1); cout <<"sincostimes error = "<<(r1-r3).norm()< > siv2,cov2; t=clock()/((double) (CLOCKS_PER_SEC)); sincostimes(hreal,siv2,cov2,rhs); cout <<"sparse sincostimes "< > r4 = cov2 + siv2*complex(0,1); cout <<"sincostimes error 2 = "<<(r1-r4).norm()< > h; cin >> h; h *= complex(0,1); h.simplify(); cout <<"length of hamiltonian "< > hexp = exp(h); cout <<"sparse exp time "< > he(hexp); NRMat > hec(he); NRMat > het(he); hec.conjugateme(); het.transposeme(); //cout <>n; NRSMat hd(n); hd.randomize(1); SparseSMat h(hd); NRMat rd = hd*hd; SparseSMat r = h*h; NRSMatrx(r); NRMat r2(rx); cout <<"Error = "<<(r2-rd).norm()< h0; cin>>h0; cout <<"matrix read\n"; cout.flush(); SparseSMat h1 = h0; //.submatrix(0,2047,0,2047); SparseSMat > h = imagmatrix(h1); double t=clock()/((double) (CLOCKS_PER_SEC)); SparseSMat > r = h*h; cout <<"SparseSMat mult time "< > h3(h); NRMat > h2(h3); NRMat >r2 = exp(h2); cout <<"error = "<<(r2-NRMat >(NRSMat >(r))).norm()< > b(bb); //cholesky(b,0); //if(n<100)cout <<"Dense Cholesky result\n"< > c; c= cc.cholesky(); NRMat > cx(c); if(n<100)cout <<"Sparse pivoted Cholesky result \n"< bb(n,n); bb.gemm(0.,bh,'c',bh,'n',1.); if(n<1000) cout <<"Input matrix\n"< br(n,n); br.gemm(0.,b,'c',b,'n',1.); if(n<1000) cout <<"Result of reconstruction\n"<cgpu; a.moveto(gpu1); b.moveto(gpu1); t=clock()/((double) (CLOCKS_PER_SEC)); for(int i=0; i h1 = h0; CSRMat > h = imagmatrix(h1); double t=clock()/((double) (CLOCKS_PER_SEC)); CSRMat > r = h*h; cout <<"CSRMat mult time "< > h2(h); NRMat >r2 = exp(h2); cout <<"error = "<<(r2-NRMat >(r)).norm()< > m; ifstream f("libormat"); f >> m; f.close(); NRVec eivals(m.nrows()); NRMat > m_aux = m; NRMat > m_test = m; cout << "hermiticity error = " <<(m_aux.transpose(true) - m_aux).norm()< > eivalsc(eivals); NRMat > m4 = m.transpose(true); cout <<"unitarity error "<< (m4*m).norm(true)< > m5(m.nrows(),m.nrows()); for(int i=0; i(m(j,i).real(), - m(j,i).imag()); cout << "conjugatetest "<<(m4-m5).norm()< > m1= m_aux * m; NRMat > m1x = m; m1x.diagmultr(eivalsc); cout << "test m1 m1x "<<(m1-m1x).norm()< > m2= m.transpose(true) * m1; //NRMat > m2b= m * m_aux * m.transpose(true); cout <<"check "< > m3= vl.transpose(true)* m_aux *vr; cout <<"check2 "<>v; r= q*v*q1; cout <<"rotated "< q(1.); Quaternion p,r(q); p=r; q+=p; q-=r; r=p+q; r+=Quaternion(0,8,8,8); r+= 2.; r/=4.; cout< s(-1,2,1,2); Quaternion t=r*s.conjugate(); cout <tt = t.inverse(); cout < rotmat; quat2rotmat(r,rotmat); cout << rotmat[0][1]< rotmat2(3,3),rotmat3(3,3); Quaternion > rotmatder; rotmatder[0].resize(3,3); rotmatder[1].resize(3,3); rotmatder[2].resize(3,3); rotmatder[3].resize(3,3); normquat2rotmatder(r,rotmatder); cout << rotmatder; normquat2rotmat(r,rotmat2); normquat2rotmat(-r,rotmat3); cout << rotmat2< rr; rotmat2normquat(rotmat2,rr); cout <<"orig "< eul; r.normalize(NULL,true); r.normquat2eulerzyx(eul); cout<< "euler "< rback; rback.eulerzyx2normquat(eul); cout <<"q from euler back "< rm; euler2rotmat((double *)eul,rm,"xyz",false,false,false); cout < eulback; rotmat2euler((double *)eulback,rm,"xyz",false,false,false); cout <<"euler back from rm error "< axis={0,1/sqrt(2),1/sqrt(2)}; Quaternion rrr; rrr.axis2normquat(axis,0.5); NRVec axis2(3); double angle; rrr.normquat2axis(&axis2[0],angle); cout <<"back angle "< vec={1,2,3}; Quaternion qvec(vec); Quaternion rvec = qvec.rotateby(rrr); Quaternion rvec2 = qvec.rotateby(rrr.conjugate()); cout < rrvec = rvec.rotateby(rrr.conjugate()); cout < rotvec; rrr.rotate(rotvec,vec); Quaternion qq={1.5,2,-2,-1.123}; cout << " test "< qq2={-.5,1,2,3.23}; Mat3 m(1,2,3,4,5,6,7,8,-9); NRMat mm(m,3,3); Vec3 v(10,20,30); NRVec vv(v,3); cout << m*v< m2 = m.inverse(); cout< mm2(3,3); mm2=1.; multiply_by_inverse(mm2,mm); cout < xeul,yeul; qq.normquat2euler((double *)xeul,"zyx"); qq.normquat2eulerzyx((double *)yeul); cout < xqq; xqq.euler2normquat((double *)xeul,"zyx"); cout <<"normquat2euler test "< a,b,c; cin >>a>>b; c=a+b; cout< p; cin >>p; int n=p.size(); NRVec_from1 v(n); int i; for(i=1; i<=n; ++i) v[i]=10.*i; cout < c; cin>>c >>n; c.simplify(); cout< p(c); cout < cc(p); cout < pp(cc); cout <<(p == pp)< cc1,cc2; cc1=cc.pow_simple(n); cc2=cc.pow(n); cout<>n; NRPerm p(n); p.randomize(); cout < cc(p); cout < pp(cc,n); cout < v(n); for(int i=0; i vv(v); v.permuteme(cc); cout < vvv= vv.permuted(pp); cout<>n; NRVec v(n); v.randomize(1.); NRVec vv(v); NRPerm p(n); vv.sort(0,p); NRVec vvv=v.permuted(p,true); NRVec v4=vv.permuted(p,false); cout<>n; NRPerm p(n); p.identity(); do{ cout <>n; NRPerm p(n); int tot=p.generate_all_lex(printme); cout <<"generated "<>n; NRPerm p(n); int seed; int f=open("/dev/random",O_RDONLY); if(sizeof(int)!=read(f,&seed,sizeof(int))) laerror("cannot read /dev/random"); close(f); srand(seed); for(int i=1; i<=n; ++i) p[i]=1+RANDINT32()%(n-1); cout <<"Initial = "<>n; NRPerm p(n); NRVec_from1 classes; cin>>classes; if(classes.size()!=p.size()) laerror("sizes do not match"); PermutationAlgebra tot = p.list_restricted(classes,-2); tot.copyonwrite(); cout <<"generated " < all= p.list_all(); all.sort(); cout <<"ALL\n"<>n >>unitary_n; Sn_characters Sn(n); cout < p(n); space_dim=0; int tot=p.generate_all(pprintme,0); cout <<"partitions generated "< r=allyoung[i]*allyoung[j]; NRMat rm(r,false,n); NRMat rm2 = allyoungmat[i]*allyoungmat[j]; NRMat rreg=RegularRepresentation(r,Snmtable); NRMat rreg2=allyoungregular[i]*allyoungregular[j]; cout <<"Product of Young "<>n ; Sn_characters Sn(n); cout < p(1); p[0]=p[1]=1; CycleIndex c(Sn); PERM_RANK_TYPE d; Polynomial pp=c.substitute(p,&d); for(int i=0; i>n >>m >>x; if(n p(n),q(m); p.randomize(1.); q.randomize(1.); NRVec qq(q); Polynomial qqq(qq); Polynomial pp(n); pp.randomize(1.); p*=10.; Polynomial a,b; p.polydiv(q,a,b); Polynomial r=p*q; Polynomial z=value(p,q); //p(q(x)) Polynomial y=value(q,p); cout < u(5,5); u.randomize(1.); cout << (value(p,u)*value(q,u) -value(r,u)).norm()<>n ; NRVec r(n); //r.randomize(1.); //wilkinson's ill-conditionel polynomial for(int i=0; i p=polyfromroots(r); cout < rr= p.realroots(1e-10); rr.resize(n,true); cout <>n ; NRVec x(n+1),y(n+1); x.randomize(1.); y.randomize(1.); Polynomial p=lagrange_interpolation(x,y); cout < yy=values(p,x); cout<<"interpolation error= "<<(y-yy).norm()<q=p.integral(2); Polynomialpp=q.derivative(2); cout<<"test deriv. "<<(pp-p).norm()<>n; NRVec rr(n); rr.randomize(1.); rr.sort(0); if(rr.size()>2) rr[1]=rr[0];//make a degenerate root NRVec pr(2*n); NRVec qr(2*n); pr.randomize(1.); qr.randomize(1.); for(int i=0; i p=polyfromroots(pr); Polynomial q=polyfromroots(qr); Polynomial g=poly_gcd(p,q,1e-8); cout <<"GCD ="< gg=svd_gcd(p,q,1e-13); cout <<"SVDGCD ="< rrr=g.realroots(1e-5); NRVec rrrr=gg.realroots(1e-5); cout <>n; Polynomial p; p.binomial(n); Polynomial q(1); q[0]=q[1]=1.; Polynomial qq=q.pow(n); cout <<"test binom "<<(p-qq).norm()<>n ; Polynomial p=hermite_polynomial(n); cout <<"Hermite = "< v({1,2,3,4}); cout < vv({1,2}); cout < m({{1,2,3},{4,5,6}}); cout< mm({{1,2,3},{4,5,6}}); cout< x({1,2,3}); cout< y({{1,2,3},{4,5,6},{7,8,9}}); cout < q({1,2,3,4}); cout< p({1,4,2,3}); cout < c({NRVec_from1({1,2}),NRVec_from1({3,4,5})}); cout < pp({1,2,3,4,5}); cout<> scale; NRMat tmp(3,3); tmp.randomize(scale); Mat3 mm(tmp); mm.symmetrize(); mm+= 1.; NRMat m(&mm[0][0],3,3); cout < w(3); diagonalize(m,w); cout << w< ww; Mat3 vv; mm.eivec_sym(ww,vv); cout <<"3\n"< www(&ww[0],3); NRMat vvv(&vv[0][0],3,3); cout<<"eival error = "<<(w-www).norm()< tmp(3,3); tmp.randomize(2.); Mat3 mm(tmp); NRMat m(&mm[0][0],3,3); cout < mmi=mm.inverse(); double dd=simple_gaussj(mm.elements(),rr); cout <<"det="<(rr); cout<<"linear solve det="< fit; double funcs[3]; for(int i=0; i<100; ++i) { double x = 10*RANDDOUBLESIGNED(); funcs[0]=1; funcs[1]=x; funcs[2]=x*x; double y = 2*funcs[2] -3*funcs[1] + funcs[0]; //cout <<"test "< m; cin >>m; int n=m.nrows(); int seed; int f=open("/dev/random",O_RDONLY); if(sizeof(int)!=read(f,&seed,sizeof(int))) laerror("cannot read /dev/random"); close(f); srand(seed); NRPerm p(n); p.randomize(); cout < mm=m.permuted_rows(p); cout < p,q; cin>>p>>q; cout <>x; ContFrac xx(x,20,1000000); cout< r({11,101}); ContFrac x(r); ContFrac y= x+Rational({2,3}); cout<(y)< z= x*Rational({2,3}); cout<(z)< x(11,101); ContFrac v(3,7); ContFrac y= x+v; cout<(y)< z= x*v; cout<(z)< h({{{12,4},{3,1}},{{0,1},{-1,0}}}); ContFrac zz=h.value(x,v); cout<(zz)<(x)+3)*(Rational(v)+4)/(Rational(x)-Rational(v))<>x; ContFrac xx(x,15,100000); cout < xm= -xx; cout<< xm< x; cin >>x; ContFrac xx(x); cout< xm= -xx; cout<< xm; ContFrac yy(-x); cout< ym= -yy; cout<< ym; cout <<"ZEROs\n"<>x; ContFrac xx(x,25,100000); ContFrac xx1(x+1e-4,25,100000); ContFrac xx2(x-1e-4,25,100000); xx.canonicalize(); xx1.canonicalize(); xx2.canonicalize(); cout<<"small "<xx2) < > uubar=expi(hbar); NRMat > ubar = NRMat >(ainv) * u * NRMat >(a); cout <<"error of similarity transformed exponential = "<<(ubar-uubar).norm()< > ubarh = ubar.transpose(true); NRMat > ubarh2 = NRMat >(atr) * u.transpose(true) * NRMat >(ainvtr); cout <<"error of ubar hermit conjugation = "<<(ubarh-ubarh2).norm()< > uubarh = uubar.transpose(true); //cout <<"H = "< adj; cin >>adj; NRVec clique = findclique(adj); cout <<"First clique is "< cover = cliquecover(adj); cout <<"Clique cover is "< p(cover.size()); cover.sort(0,p,true); cout<<"permutation to disentabgle the cliques = "< adjperm = adj.permuted(p); cout <<"resorted graph = "< a; a.randomize(10.); cout< b; b.randomize(10.); cout < x = a.linear_solve(b); cout < a; a.randomize(10.); cout< ata = a.transpose()*a; Mat3 work(ata); Vec3 w; Mat3 v; work.eivec_sym(w,v,true); //cout < aat = a*a.transpose(); Mat3 work2(aat); Vec3 w2; Mat3 v2; work2.eivec_sym(w2,v2,true); //cout < u3,v3; Vec3 w3; a.svd(u3,w3,v3); cout <<"Mat3 SVD\n"; cout < aa(a); NRMat uu(3,3),vv(3,3); NRVec ss(3); singular_decomposition(aa,&uu,ss,&vv,0); cout <<"REGULAR SVD\n"; cout < ainv= a.inverse(); Mat3 ainv2= a.svdinverse(); cout <>n; bitvector v(n); v.randomize(); do{ cout <>=1; }while(!v.iszero()); cout <>n; bitvector v(n),u(n); u.randomize(); v.randomize(); bitvector w=u*v; bitvector z=v*u; cout <>n; cout <>ir; if(!ir.is_irreducible()) laerror("input must be an irreducible polynomial"); bitvector a; cin >>a; bitvector ai = a.field_inv(ir); cout<< "inverse = "<>d; bitvector ir=find_irreducible(d,3); if(ir.is_zero()) ir=find_irreducible(d,5); if(ir.is_zero()) ir=find_irreducible(d,7); if(ir.is_zero()) {cout<<"cannot find IR polynomial\n"; exit(1);} cout <<"IR = "<>n; NRPerm p(n); PermutationAlgebra a=p.list_all(); for(int i=0; i > groups; cin >> groups; int ntot=0; for(int i=0; i indices(ntot); cin >>indices; PermutationAlgebra a=general_antisymmetrizer(groups,-2,true); for(int i=0; i pindices=applypermutation(a[i].perm,indices,false); for(int j=1; j<=pindices.size(); ++j) cout < a; cin>>a; cout <>n; int count=0; for(int i=4; i<=n; ++i) for(int j=3; j>d>>n; unsigned long long s1; s1=simplicial(d,n); cout < d({6,2,1,4,3,5}); d.copyonwrite(); netsort(d.size(),&d[0]); cout < epsilon(g); epsilon.clear(); cout < I({1,2,3}); NRVec > II(I,1); epsilon.lhs(II)=1; II.copyonwrite(); II[0][0]=3; II[0][1]=2; II[0][2]=1; epsilon *= 2.; epsilon.lhs(1,2,3) -= 1.; cout < s(4); for(int i=0; i<4; ++i) { s[i].number=1; s[i].symmetry=0; s[i].offset=0; s[i].range=n; } Tensor t(s); t.randomize(1.); cout < p({3,1,4,2}); Tensor tt=t.permute_index_groups(p); cout < e(g); e.randomize(1.); Tensor eu = e.unwind_index(0,1); for(int i=0; i a(ag); a.randomize(1.); INDEXGROUP bg; bg.number=3; bg.symmetry= 0; bg.offset=0; bg.range=n; Tensor b(bg); b.randomize(1.); Tensor cc = a.contraction(0,0,b,0,1); cout < shape({cgb,cga}); Tensor c(shape); c.clear(); for(int i=0; i1e-13) laerror("internal error in conntraction"); } //cout <