/* LA: linear algebra C++ interface library Copyright (C) 2024 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 . */ #include #include "tensor.h" #include "laerror.h" #include "qsort.h" #include "miscfunc.h" #include #include "bitvector.h" namespace LA { template int Tensor:: calcrank() { int r=0; for(int i=0; i *>(&shape))[i]; if(sh->number==0) laerror("empty index group"); r+=sh->number; } myrank=r; return r; } template LA_largeindex Tensor::calcsize() { groupsizes.resize(shape.size()); cumsizes.resize(shape.size()); LA_largeindex s=1; for(int i=0; i *>(&shape))[i]; if(sh->number==0) laerror("empty index group"); if(sh->range==0) return 0; cumsizes[i]=s; switch(sh->symmetry) { case 0: s *= groupsizes[i] = longpow(sh->range,sh->number); break; case 1: s *= groupsizes[i] = simplicial(sh->number,sh->range); break; case -1: s *= groupsizes[i] = simplicial(sh->number,sh->range-sh->number+1); break; default: laerror("illegal index group symmetry"); break; } } return s; } LA_largeindex subindex(int *sign, const INDEXGROUP &g, const NRVec &I) //index of one subgroup { #ifdef DEBUG if(I.size()<=0) laerror("empty index group in subindex"); if(g.number!=I.size()) laerror("mismatch in the number of indices in a group"); for(int i=0; i= g.offset+g.range) laerror("index out of range in tensor subindex"); #endif switch(I.size()) //a few special cases for efficiency { case 0: *sign=0; return 0; break; case 1: *sign=1; return I[0]-g.offset; break; case 2: { *sign=1; if(g.symmetry==0) return (I[1]-g.offset)*g.range+I[0]-g.offset; LA_index i0,i1; if(I[0]>I[1]) {i1=I[0]; i0=I[1]; if(g.symmetry<0) *sign = -1;} else {i1=I[1]; i0=I[0];} i0 -= g.offset; i1 -= g.offset; if(g.symmetry<0) { if(i0==i1) {*sign=0; return -1;} return i1*(i1-1)/2+i0; } else { return i1*(i1+1)/2+i0; } } break; default: //general case { *sign=1; if(g.symmetry==0) //rectangular case { LA_largeindex r=0; for(int i=I.size()-1; i>=0; --i) { r*= g.range; r+= I[i]-g.offset; } return r; } } //compressed storage case NRVec II(I); II.copyonwrite(); if(g.offset!=0) II -= g.offset; int parity=netsort(II.size(),&II[0]); if(g.symmetry<0 && (parity&1)) *sign= -1; if(g.symmetry<0) //antisymmetric { for(int i=0; i inverse_subindex(const INDEXGROUP &g, LA_largeindex s) { NRVec I(g.number); if(g.number==1) {I[0]=s+g.offset; return I;} switch(g.symmetry) { case 0: for(int i=0; i0; --i) { I[i-1] = inverse_simplicial(i,s); s -= simplicial(i,I[i-1]); } break; case -1: for(int i=g.number-1; i>=0; --i) { I[i] = i + inverse_simplicial(i+1,s); s -= simplicial(i+1,I[i]-i); } break; default: laerror("illegal index symmetry"); } if(g.offset!=0) I += g.offset; return I; } template SUPERINDEX Tensor::inverse_index(LA_largeindex s) const { SUPERINDEX I(shape.size()); for(int g=shape.size()-1; g>=0; --g) { LA_largeindex groupindex; if(g>0) { groupindex = s/cumsizes[g]; s %= cumsizes[g]; } else groupindex=s; I[g] = inverse_subindex(shape[g],groupindex); } return I; } template LA_largeindex Tensor::index(int *sign, const SUPERINDEX &I) const { //check index structure and ranges #ifdef DEBUG if(I.size()!=shape.size()) laerror("mismatch in the number of tensor index groups"); for(int i=0; i= shape[i].offset+shape[i].range) { std::cerr<<"error in index group no. "< LA_largeindex Tensor::index(int *sign, const FLATINDEX &I) const { #ifdef DEBUG if(rank()!=I.size()) laerror("tensor rank mismatch in index"); #endif LA_largeindex r=0; *sign=1; int gstart=0; for(int g=0; g subI = I.subvector(gstart,gend); gstart=gend+1; LA_largeindex groupindex = subindex(&gsign,shape[g],subI); //std::cout <<"FLATINDEX TEST group "< LA_largeindex Tensor::vindex(int *sign, LA_index i1, va_list args) const { NRVec I(rank()); I[0]=i1; for(int i=1; i void Tensor::put(int fd) const { shape.put(fd,true); groupsizes.put(fd,true); cumsizes.put(fd,true); data.put(fd,true); } template void Tensor::get(int fd) { shape.get(fd,true); myrank=calcrank(); //is not stored but recomputed groupsizes.put(fd,true); cumsizes.get(fd,true); data.get(fd,true); } template Tensor::Tensor(const NRVec &x) : data(x) { myrank=1; shape.resize(1); shape[0].number=1; shape[0].symmetry=0; #ifndef LA_TENSOR_ZERO_OFFSET shape[0].offset=0; #endif shape[0].range=x.size(); calcsize(); } template Tensor::Tensor(const NRMat &x) : data(&x(0,0),x.nrows()*x.ncols()) { myrank=2; if(x.nrows()==x.ncols()) { shape.resize(1); shape[0].number=2; shape[0].symmetry=0; #ifndef LA_TENSOR_ZERO_OFFSET shape[0].offset=0; #endif shape[0].range=x.nrows(); } else { shape.resize(2); shape[0].number=1; shape[1].number=1; shape[0].symmetry=0; shape[1].symmetry=0; #ifndef LA_TENSOR_ZERO_OFFSET shape[0].offset=0; shape[1].offset=0; #endif shape[0].range=x.ncols(); shape[1].range=x.nrows(); } calcsize(); } template Tensor::Tensor(const NRSMat &x) : data(NRVec(x)) { myrank=2; shape.resize(1); shape[0].number=2; shape[0].symmetry=1; #ifndef LA_TENSOR_ZERO_OFFSET shape[0].offset=0; #endif shape[0].range=x.nrows(); calcsize(); } template void loopingroups(Tensor &t, int ngroup, int igroup, T **p, SUPERINDEX &I, void (*callback)(const SUPERINDEX &, T *)) { LA_index istart,iend; const indexgroup *sh = &(* const_cast *>(&t.shape))[ngroup]; switch(sh->symmetry) { case 0: istart= sh->offset; iend= sh->offset+sh->range-1; break; case 1: istart= sh->offset; if(igroup==sh->number-1) iend= sh->offset+sh->range-1; else iend = I[ngroup][igroup+1]; break; case -1: istart= sh->offset + igroup; if(igroup==sh->number-1) iend= sh->offset+sh->range-1; else iend = I[ngroup][igroup+1]-1; break; } for(LA_index i = istart; i<=iend; ++i) { I[ngroup][igroup]=i; if(ngroup==0 && igroup==0) { int sign; //std::cout <<"TEST "< *>(&t.shape))[newngroup]; newigroup=sh2->number-1; } loopingroups(t,newngroup,newigroup,p,I,callback); } } } template void Tensor::loopover(void (*callback)(const SUPERINDEX &, T *)) { SUPERINDEX I(shape.size()); for(int i=0; i *>(&shape))[i]; I[i].resize(sh->number); I[i] = sh->offset; } T *pp=&data[0]; int ss=shape.size()-1; const indexgroup *sh = &(* const_cast *>(&shape))[ss]; loopingroups(*this,ss,sh->number-1,&pp,I,callback); } static std::ostream *sout; template static void outputcallback(const SUPERINDEX &I, T *v) { //print indices flat for(int i=0; i>(std::istream &s, INDEXGROUP &x) { s>>x.number>>x.symmetry; #ifndef LA_TENSOR_ZERO_OFFSET s>>x.offset; #endif s>>x.range; return s; } template std::ostream & operator<<(std::ostream &s, const Tensor &x) { s< *>(&x)->loopover(&outputcallback); return s; } template std::istream & operator>>(std::istream &s, Tensor &x) { s>>x.shape; x.data.resize(x.calcsize()); x.calcrank(); FLATINDEX I(x.rank()); for(LA_largeindex i=0; i>I[j]; T val; s>>val; x.lhs(I) = val; } return s; } template void loopovergroups(Tensor &t, int ngroup, T **p, GROUPINDEX &I, void (*callback)(const GROUPINDEX &, T *)) { for(LA_largeindex i = 0; i void Tensor::grouploopover(void (*callback)(const GROUPINDEX &, T *)) { GROUPINDEX I(shape.size()); T *pp=&data[0]; loopovergroups(*this,shape.size()-1,&pp,I,callback); } const NRPerm *help_p; template Tensor *help_t; template const Tensor *help_tt; template static void permutecallback(const GROUPINDEX &I, T *v) { LA_largeindex target=0; for(int i=0; i< help_t->shape.size(); ++i) { target += I[(*help_p)[i+1]-1] * help_t->cumsizes[i]; } help_t->data[target] = *v; } template Tensor Tensor::permute_index_groups(const NRPerm &p) const { NRVec newshape=shape.permuted(p); Tensor r(newshape); //prepare statics for the callback help_p = &p; help_t = &r; //now rearrange the data const_cast *>(this)->grouploopover(permutecallback); return r; } FLATINDEX superindex2flat(const SUPERINDEX &I) { int rank=0; for(int i=0; i &shape) { int ii=0; for(int g=0; g &shape) { int ii=0; for(int g=0; g static void unwind_callback(const SUPERINDEX &I, T *v) { FLATINDEX J = superindex2flat(I); FLATINDEX JP = J.permuted(*help_p,true); //std::cout <<"TEST unwind_callback: from "<)(JP); //rhs operator() generates the redundant elements for the unwinded lhs tensor } template Tensor Tensor::unwind_index(int group, int index) const { if(group<0||group>=shape.size()) laerror("wrong group number in unwind_index"); if(index<0||index>=shape[group].number) laerror("wrong index number in unwind_index"); if(shape[group].number==1) //single index in the group { if(group==0) return *this; //is already the least significant group NRPerm p(shape.size()); p[1]= 1+group; int ii=1; if(ii==1+group) ii++; //skip this for(int i=2; i<=shape.size(); ++i) { p[i]=ii++; if(ii==1+group) ii++; //skip this } if(!p.is_valid()) laerror("internal error in unwind_index"); return permute_index_groups(p); } //general case - recalculate the shape and allocate the new tensor NRVec newshape(shape.size()+1); newshape[0].number=1; newshape[0].symmetry=0; newshape[0].range=shape[group].range; #ifndef LA_TENSOR_ZERO_OFFSET newshape[0].offset = shape[group].offset; #endif int flatindex=0; //(group,index) in flat form for(int i=0; i r(newshape); if(r.rank()!=rank()) laerror("internal error 2 in unwind_index"); //compute the corresponding permutation of FLATINDEX for use in the callback NRPerm indexperm(rank()); indexperm[1]=flatindex+1; int ii=1; if(ii==flatindex+1) ii++; for(int i=2; i<=rank(); ++i) { indexperm[i] = ii++; if(ii==flatindex+1) ii++; //skip this } if(!indexperm.is_valid()) { std::cout << "indexperm = "< = this; help_p = &indexperm; r.loopover(unwind_callback); return r; } template Tensor Tensor::unwind_indices(const INDEXLIST &il) const { if(il.size()==0) return *this; if(il.size()==1) return unwind_index(il[0].group,il[0].index); for(int i=0; i=shape.size()) laerror("wrong group number in unwind_indices"); if(il[i].index<0||il[i].index>=shape[il[i].group].number) laerror("wrong index number in unwind_indices"); } //all indices are solo in their groups - permute groups bool sologroups=true; int nonsolo=0; for(int i=0; i p(shape.size()); bitvector waslisted(shape.size()); waslisted.clear(); for(int i=0; i oldshape(shape); oldshape.copyonwrite(); NRVec newshape(shape.size()+nonsolo); //first the unwound indices as solo groups for(int i=0; i0) { newshape[ii++] = oldshape[i]; } else ++emptied_groups; if(emptied_groups) newshape.resize(newshape.size()-emptied_groups,true); Tensor r(newshape); if(r.rank()!=rank()) laerror("internal error 2 in unwind_indces"); //compute the corresponding permutation of FLATINDEX for use in the callback NRPerm indexperm(rank()); bitvector waslisted(rank()); waslisted.clear(); //first unwound indices ii=0; for(int i=0; i = this; help_p = &indexperm; r.loopover(unwind_callback); return r; } template static void auxmatmult(int nn, int mm, int kk, T *r, T *a, T *b, T alpha=1, T beta=0, bool conjugate=false) //R(nn,mm) = A(nn,kk) * B^T(mm,kk) { for(int i=0; i::conjugate(b[j*kk+k]) : b[j*kk+k]; } } template<> void auxmatmult(int nn, int mm, int kk, double *r, double *a, double *b, double alpha, double beta, bool conjugate) { cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasTrans, nn, mm, kk, alpha, a, kk, b, kk, beta, r, mm); } template<> void auxmatmult >(int nn, int mm, int kk, std::complex *r, std::complex *a, std::complex *b, std::complex alpha, std::complex beta, bool conjugate) { cblas_zgemm(CblasRowMajor, CblasNoTrans, (conjugate?CblasConjTrans:CblasTrans), nn, mm, kk, &alpha, a, kk, b, kk, &beta, r, mm); } //Conntraction could be implemented without the temporary storage for unwinding, but then we would need //double recursion over indices of both tensors. Hopefully using the matrix multiplication here //makes it also more efficient, even for (anti)symmetric indices //The index unwinding is unfortunately a big burden, and in principle could be eliminated in case of non-symmetric indices // template void Tensor::addcontraction(const Tensor &rhs1, int group, int index, const Tensor &rhs, int rhsgroup, int rhsindex, T alpha, T beta, bool doresize, bool conjugate1, bool conjugate) { if(group<0||group>=rhs1.shape.size()) laerror("wrong group number in contraction"); if(rhsgroup<0||rhsgroup>=rhs.shape.size()) laerror("wrong rhsgroup number in contraction"); if(index<0||index>=rhs1.shape[group].number) laerror("wrong index number in conntraction"); if(rhsindex<0||rhsindex>=rhs.shape[rhsgroup].number) laerror("wrong index number in conntraction"); if(rhs1.shape[group].offset != rhs.shape[rhsgroup].offset) laerror("incompatible index offset in contraction"); if(rhs1.shape[group].range != rhs.shape[rhsgroup].range) laerror("incompatible index range in contraction"); Tensor u = rhs1.unwind_index(group,index); if(conjugate1) u.conjugateme(); Tensor rhsu = rhs.unwind_index(rhsgroup,rhsindex); NRVec newshape(u.shape.size()+rhsu.shape.size()-2); int ii=0; for(int i=1; i(nn,mm,kk,&data[0],&u.data[0], &rhsu.data[0],alpha,beta,conjugate); } template void Tensor::addcontractions(const Tensor &rhs1, const INDEXLIST &il1, const Tensor &rhs2, const INDEXLIST &il2, T alpha, T beta, bool doresize, bool conjugate1, bool conjugate2) { if(il1.size()==0) laerror("empty contraction - outer product not implemented"); if(il1.size()!=il2.size()) laerror("mismatch in index lists in addcontractions"); for(int i=0; i=rhs1.shape.size()) laerror("wrong group1 number in contractions"); if(il2[i].group<0||il2[i].group>=rhs2.shape.size()) laerror("wrong group2 number in contractions"); if(il1[i].index<0||il1[i].index>=rhs1.shape[il1[i].group].number) laerror("wrong index1 number in conntractions"); if(il2[i].index<0||il2[i].index>=rhs2.shape[il2[i].group].number) laerror("wrong index2 number in conntractions"); if(rhs1.shape[il1[i].group].offset != rhs2.shape[il2[i].group].offset) laerror("incompatible index offset in contractions"); if(rhs1.shape[il1[i].group].range != rhs2.shape[il2[i].group].range) laerror("incompatible index range in contractions"); } Tensor u = rhs1.unwind_indices(il1); if(conjugate1) u.conjugateme(); Tensor rhsu = rhs2.unwind_indices(il2); NRVec newshape(u.shape.size()+rhsu.shape.size()-2*il1.size()); int ii=0; for(int i=il1.size(); i(nn,mm,kk,&data[0],&u.data[0], &rhsu.data[0],alpha,beta,conjugate2); } template static const PermutationAlgebra *help_pa; static bool help_inverse; template static T help_alpha; template static void permutationalgebra_callback(const SUPERINDEX &I, T *v) { FLATINDEX J = superindex2flat(I); for(int p=0; p->size(); ++p) { FLATINDEX Jp = J.permuted((*help_pa)[p].perm,help_inverse); *v += help_alpha * (*help_pa)[p].weight * (*help_t)(Jp); } } static int help_tn; template static Tensor *help_tv; template static void permutationalgebra_callback2(const SUPERINDEX &I, T *v) { FLATINDEX J = superindex2flat(I); for(int p=0; p->size(); ++p) { FLATINDEX Jp = J.permuted((*help_pa)[p].perm,help_inverse); T product; //build product from the vector of tensors using the Jp index int rankshift=0; for(int i=0; i[i].rank(); FLATINDEX indi = Jp.subvector(rankshift,rankshift+rank-1); if(i==0) product= help_tv[i](indi); else product *= help_tv[i](indi); rankshift += rank; } *v += help_alpha * (*help_pa)[p].weight * product; } } template void Tensor::apply_permutation_algebra(const Tensor &rhs, const PermutationAlgebra &pa, bool inverse, T alpha, T beta) { if(beta!=(T)0) {if(beta!=(T)1) *this *= beta;} else clear(); if(alpha==(T)0) return; copyonwrite(); if(rank()!=rhs.rank()) laerror("rank mismatch in apply_permutation_algebra"); help_t = const_cast *>(&rhs); help_pa = &pa; help_inverse = inverse; help_alpha = alpha; loopover(permutationalgebra_callback); } template void Tensor::apply_permutation_algebra(const NRVec > &rhsvec, const PermutationAlgebra &pa, bool inverse, T alpha, T beta) { if(beta!=(T)0) {if(beta!=(T)1) *this *= beta;} else clear(); if(alpha==(T)0) return; copyonwrite(); int totrank=0; for(int i=0; i = const_cast *>(&rhsvec[0]); help_tn = rhsvec.size(); help_pa = &pa; help_inverse = inverse; help_alpha = alpha; loopover(permutationalgebra_callback2); } template void Tensor::split_index_group(int group) { if(group<0||group >= shape.size()) laerror("illegal index group number"); if(shape[group].number==1) return; //nothing to split if(shape[group].symmetry!=0) laerror("only non-symmetric index group can be splitted"); NRVec newshape(shape.size()+shape[group].number-1); int gg=0; for(int g=0; g void Tensor:: merge_adjacent_index_groups(int groupfrom, int groupto) { if(groupfrom<0||groupfrom>= shape.size()) laerror("illegal index group number"); if(groupto<0||groupto>= shape.size()) laerror("illegal index group number"); if(groupfrom==groupto) return; if(groupfrom>groupto) {int t=groupfrom; groupfrom=groupto; groupto=t;} int newnumber=0; for(int g=groupfrom; g<=groupto; ++g) { if(shape[g].symmetry!=0) laerror("only non-symmetric index groups can be merged"); if(shape[g].offset!=shape[groupfrom].offset) laerror("incompatible offset in merge_adjacent_index_groups"); if(shape[g].range!=shape[groupfrom].range) laerror("incompatible range in merge_adjacent_index_groups"); newnumber += shape[g].number; } NRVec newshape(shape.size()-(groupto-groupfrom+1)+1); for(int g=0; g<=groupfrom; ++g) newshape[g]=shape[g]; newshape[groupfrom].number=newnumber; for(int g=groupfrom+1; g Tensor Tensor::merge_index_groups(const NRVec &groups) const { if(groups.size()<=1) return *this; NRPerm p(shape.size()); int gg=0; bitvector selected(shape.size()); selected.clear(); for(int g=0; g=shape.size()) laerror("illegal group number in merge_index_groups"); if(selected[g]) laerror("repeated group number in merge_index_groups"); selected.set(g); p[gg++] = 1+groups[g]; } for(int g=0; g r = permute_index_groups(p); r.merge_adjacent_index_groups(0,groups.size()-1); return r; } template NRVec > Tensor::Tucker(typename LA_traits::normtype thr) { int r=rank(); NRVec > ret(r); if(r<2) return ret; int rr=0; for(int i=0; i um; NRVec ushape; { Tensor u=unwind_index(i,j); ushape=u.shape; um=u.matrix(); } int mini=um.nrows(); if(um.ncols() u(um.nrows(),mini),vt(mini,um.ncols()); NRVec::normtype> w(mini); singular_decomposition(um,&u,w,&vt,0); um.resize(0,0); //deallocate int preserve=mini; for(int k=0; k umnew; if(preserve newdata(umnew); *this = Tensor(ushape,newdata); } return ret; } template class Tensor; template class Tensor >; template std::ostream & operator<<(std::ostream &s, const Tensor &x); template std::ostream & operator<<(std::ostream &s, const Tensor > &x); template std::istream & operator>>(std::istream &s, Tensor &x); template std::istream & operator>>(std::istream &s, Tensor > &x); }//namespace