/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./lib/algorithms/CoarsenedMatrix.h Copyright (C) 2015 Author: Azusa Yamaguchi Author: Peter Boyle Author: Peter Boyle Author: paboyle 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 2 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, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. See the full license in the file "LICENSE" in the top level distribution directory *************************************************************************************/ /* END LEGAL */ #ifndef GRID_ALGORITHM_COARSENED_MATRIX_H #define GRID_ALGORITHM_COARSENED_MATRIX_H #include namespace Grid { class Geometry { // int dimension; public: int npoint; std::vector directions ; std::vector displacements; Geometry(int _d) { int base = (_d==5) ? 1:0; // make coarse grid stencil for 4d , not 5d if ( _d==5 ) _d=4; npoint = 2*_d+1; directions.resize(npoint); displacements.resize(npoint); for(int d=0;d<_d;d++){ directions[2*d ] = d+base; directions[2*d+1] = d+base; displacements[2*d ] = +1; displacements[2*d+1] = -1; } directions [2*_d]=0; displacements[2*_d]=0; //// report back std::cout< GetDelta(int point) { std::vector delta(dimension,0); delta[directions[point]] = displacements[point]; return delta; }; */ }; template class Aggregation { public: typedef iVector siteVector; typedef Lattice CoarseVector; typedef Lattice > CoarseMatrix; typedef Lattice< CComplex > CoarseScalar; // used for inner products on fine field typedef Lattice FineField; GridBase *CoarseGrid; GridBase *FineGrid; std::vector > subspace; Aggregation(GridBase *_CoarseGrid,GridBase *_FineGrid) : CoarseGrid(_CoarseGrid), FineGrid(_FineGrid), subspace(nbasis,_FineGrid) { }; void Orthogonalise(void){ CoarseScalar InnerProd(CoarseGrid); blockOrthogonalise(InnerProd,subspace); } void CheckOrthogonal(void){ CoarseVector iProj(CoarseGrid); CoarseVector eProj(CoarseGrid); Lattice pokey(CoarseGrid); for(int i=0;ioSites();ss++){ eProj._odata[ss](i)=CComplex(1.0); } eProj=eProj - iProj; std::cout< &hermop,int nn=nbasis) { RealD scale; ConjugateGradient CG(1.0e-2,10000); FineField noise(FineGrid); FineField Mn(FineGrid); for(int b=0;b "< "< class CoarsenedMatrix : public SparseMatrixBase > > { public: typedef iVector siteVector; typedef Lattice CoarseVector; typedef Lattice > CoarseMatrix; typedef Lattice< CComplex > CoarseScalar; // used for inner products on fine field typedef Lattice FineField; //////////////////// // Data members //////////////////// Geometry geom; GridBase * _grid; CartesianStencil > Stencil; std::vector A; std::vector > comm_buf; /////////////////////// // Interface /////////////////////// GridBase * Grid(void) { return _grid; }; // this is all the linalg routines need to know RealD M (const CoarseVector &in, CoarseVector &out){ conformable(_grid,in._grid); conformable(in._grid,out._grid); SimpleCompressor compressor; Stencil.HaloExchange(in,comm_buf,compressor); PARALLEL_FOR_LOOP for(int ss=0;ssoSites();ss++){ siteVector res = zero; siteVector nbr; int ptype; StencilEntry *SE; for(int point=0;point_is_local&&SE->_permute) { permute(nbr,in._odata[SE->_offset],ptype); } else if(SE->_is_local) { nbr = in._odata[SE->_offset]; } else { nbr = comm_buf[SE->_offset]; } res = res + A[point]._odata[ss]*nbr; } vstream(out._odata[ss],res); } return norm2(out); }; RealD Mdag (const CoarseVector &in, CoarseVector &out){ return M(in,out); }; // Defer support for further coarsening for now void Mdiag (const CoarseVector &in, CoarseVector &out){}; void Mdir (const CoarseVector &in, CoarseVector &out,int dir, int disp){}; CoarsenedMatrix(GridCartesian &CoarseGrid) : _grid(&CoarseGrid), geom(CoarseGrid._ndimension), Stencil(&CoarseGrid,geom.npoint,Even,geom.directions,geom.displacements), A(geom.npoint,&CoarseGrid) { comm_buf.resize(Stencil._unified_buffer_size); }; void CoarsenOperator(GridBase *FineGrid,LinearOperatorBase > &linop, Aggregation & Subspace){ FineField iblock(FineGrid); // contributions from within this block FineField oblock(FineGrid); // contributions from outwith this block FineField phi(FineGrid); FineField tmp(FineGrid); FineField zz(FineGrid); zz=zero; FineField Mphi(FineGrid); Lattice > coor(FineGrid); CoarseVector iProj(Grid()); CoarseVector oProj(Grid()); CoarseScalar InnerProd(Grid()); // Orthogonalise the subblocks over the basis blockOrthogonalise(InnerProd,Subspace.subspace); // Compute the matrix elements of linop between this orthonormal // set of vectors. int self_stencil=-1; for(int p=0;p_rdimensions[dir])/(Grid()->_rdimensions[dir]); LatticeCoordinate(coor,dir); if ( disp==0 ){ linop.OpDiag(phi,Mphi); } else { linop.OpDir(phi,Mphi,dir,disp); } //////////////////////////////////////////////////////////////////////// // Pick out contributions coming from this cell and neighbour cell //////////////////////////////////////////////////////////////////////// if ( disp==0 ) { iblock = Mphi; oblock = zero; } else if ( disp==1 ) { oblock = where(mod(coor,block)==(block-1),Mphi,zz); iblock = where(mod(coor,block)!=(block-1),Mphi,zz); } else if ( disp==-1 ) { oblock = where(mod(coor,block)==(Integer)0,Mphi,zz); iblock = where(mod(coor,block)!=(Integer)0,Mphi,zz); } else { assert(0); } Subspace.ProjectToSubspace(iProj,iblock); Subspace.ProjectToSubspace(oProj,oblock); // blockProject(iProj,iblock,Subspace.subspace); // blockProject(oProj,oblock,Subspace.subspace); PARALLEL_FOR_LOOP for(int ss=0;ssoSites();ss++){ for(int j=0;j bc(FineGrid->_ndimension,0); blockPick(Grid(),phi,tmp,bc); // Pick out a block linop.Op(tmp,Mphi); // Apply big dop blockProject(iProj,Mphi,Subspace.subspace); // project it and print it std::cout< > val(Grid()); random(RNG,val); Complex one(1.0); iMatrix ident; ident=one; val = val*adj(val); val = val + 1.0; A[8] = val*ident; // for(int s=0;soSites();s++) { // A[8]._odata[s]=val._odata[s]; // } } void ForceHermitian(void) { for(int d=0;d<4;d++){ int dd=d+1; A[2*d] = adj(Cshift(A[2*d+1],dd,1)); } // A[8] = 0.5*(A[8] + adj(A[8])); } void AssertHermitian(void) { CoarseMatrix AA (Grid()); CoarseMatrix AAc (Grid()); CoarseMatrix Diff (Grid()); for(int d=0;d<4;d++){ int dd=d+1; AAc = Cshift(A[2*d+1],dd,1); AA = A[2*d]; Diff = AA - adj(AAc); std::cout<