1
0
mirror of https://github.com/paboyle/Grid.git synced 2025-04-09 05:30:46 +01:00

WilsonMG: Revert CoarsenedMatrix.h and Lattice_transfer.h back to state of develop branch

This commit is contained in:
Daniel Richtmann 2018-05-16 16:02:17 +02:00
parent 61812ab7f1
commit df8c208f5c
2 changed files with 40 additions and 91 deletions

View File

@ -93,16 +93,12 @@ namespace Grid {
template<class Fobj,class CComplex,int nbasis> template<class Fobj,class CComplex,int nbasis>
class Aggregation { class Aggregation {
public: public:
typedef iVector<CComplex,nbasis > siteVector;
typedef Lattice<siteVector> CoarseVector;
typedef Lattice<iMatrix<CComplex,nbasis > > CoarseMatrix;
typedef typename CComplex::vector_type innerType; typedef Lattice< CComplex > CoarseScalar; // used for inner products on fine field
typedef iScalar<iScalar<iScalar<innerType > > > siteScalar; // used for inner products on fine field typedef Lattice<Fobj > FineField;
typedef iScalar<iVector<iVector<innerType, nbasis >, 1 > > siteVector;
typedef iScalar<iMatrix<iMatrix<innerType, nbasis >, 1 > > siteMatrix;
typedef Lattice<siteScalar> CoarseScalar; // used for inner products on fine field
typedef Lattice<siteVector> CoarseVector;
typedef Lattice<siteMatrix> CoarseMatrix;
typedef Lattice<Fobj> FineField;
GridBase *CoarseGrid; GridBase *CoarseGrid;
GridBase *FineGrid; GridBase *FineGrid;
@ -119,12 +115,12 @@ namespace Grid {
void Orthogonalise(void){ void Orthogonalise(void){
CoarseScalar InnerProd(CoarseGrid); CoarseScalar InnerProd(CoarseGrid);
std::cout << GridLogMessage <<"Gram-Schmidt pass 1"<<std::endl; std::cout << GridLogMessage <<" Gramm-Schmidt pass 1"<<std::endl;
blockOrthogonalise(InnerProd,subspace); blockOrthogonalise(InnerProd,subspace);
std::cout << GridLogMessage <<"Gram-Schmidt pass 2"<<std::endl; std::cout << GridLogMessage <<" Gramm-Schmidt pass 2"<<std::endl;
blockOrthogonalise(InnerProd,subspace); blockOrthogonalise(InnerProd,subspace);
std::cout << GridLogMessage <<"Gram-Schmidt checking orthogonality"<<std::endl; // std::cout << GridLogMessage <<" Gramm-Schmidt checking orthogonality"<<std::endl;
CheckOrthogonal(); // CheckOrthogonal();
} }
void CheckOrthogonal(void){ void CheckOrthogonal(void){
CoarseVector iProj(CoarseGrid); CoarseVector iProj(CoarseGrid);
@ -133,7 +129,7 @@ namespace Grid {
blockProject(iProj,subspace[i],subspace); blockProject(iProj,subspace[i],subspace);
eProj=zero; eProj=zero;
parallel_for(int ss=0;ss<CoarseGrid->oSites();ss++){ parallel_for(int ss=0;ss<CoarseGrid->oSites();ss++){
eProj._odata[ss]()(0)(i)=innerType(1.0); eProj._odata[ss](i)=CComplex(1.0);
} }
eProj=eProj - iProj; eProj=eProj - iProj;
std::cout<<GridLogMessage<<"Orthog check error "<<i<<" " << norm2(eProj)<<std::endl; std::cout<<GridLogMessage<<"Orthog check error "<<i<<" " << norm2(eProj)<<std::endl;
@ -209,8 +205,7 @@ namespace Grid {
RealD scale; RealD scale;
TrivialPrecon<FineField> TrivialPrec; ConjugateGradient<FineField> CG(1.0e-2,10000);
FlexibleGeneralisedMinimalResidual<FineField> FGMRES(1.0e-14,1,TrivialPrec,1,false); // TODO: need to use GMRES as long as Mdag doesn't work on coarser levels (i.e., MdagM isn't hermitian)
FineField noise(FineGrid); FineField noise(FineGrid);
FineField Mn(FineGrid); FineField Mn(FineGrid);
@ -222,9 +217,9 @@ namespace Grid {
hermop.Op(noise,Mn); std::cout<<GridLogMessage << "noise ["<<b<<"] <n|MdagM|n> "<<norm2(Mn)<<std::endl; hermop.Op(noise,Mn); std::cout<<GridLogMessage << "noise ["<<b<<"] <n|MdagM|n> "<<norm2(Mn)<<std::endl;
for(int i=0;i<3;i++){ for(int i=0;i<1;i++){
FGMRES(hermop,noise,subspace[b]); CG(hermop,noise,subspace[b]);
noise = subspace[b]; noise = subspace[b];
scale = std::pow(norm2(noise),-0.5); scale = std::pow(norm2(noise),-0.5);
@ -244,18 +239,15 @@ namespace Grid {
// Fine Object == (per site) type of fine field // Fine Object == (per site) type of fine field
// nbasis == number of deflation vectors // nbasis == number of deflation vectors
template<class Fobj,class CComplex,int nbasis> template<class Fobj,class CComplex,int nbasis>
class CoarsenedMatrix : public SparseMatrixBase<Lattice<iScalar<iVector<iVector<typename CComplex::vector_type, nbasis >, 1 > > > > { class CoarsenedMatrix : public SparseMatrixBase<Lattice<iVector<CComplex,nbasis > > > {
public: public:
typedef iVector<CComplex,nbasis > siteVector;
typedef Lattice<siteVector> CoarseVector;
typedef Lattice<iMatrix<CComplex,nbasis > > CoarseMatrix;
typedef typename CComplex::vector_type innerType; typedef Lattice< CComplex > CoarseScalar; // used for inner products on fine field
typedef iScalar<iScalar<iScalar<innerType > > > siteScalar; typedef Lattice<Fobj > FineField;
typedef iScalar<iVector<iVector<innerType, nbasis >, 1 > > siteVector;
typedef iScalar<iMatrix<iMatrix<innerType, nbasis >, 1 > > siteMatrix;
typedef Lattice<siteScalar> CoarseScalar; // used for inner products on fine field
typedef Lattice<siteVector> CoarseVector;
typedef Lattice<siteMatrix> CoarseMatrix;
typedef Lattice<Fobj> FineField;
//////////////////// ////////////////////
// Data members // Data members
@ -303,50 +295,13 @@ namespace Grid {
return norm2(out); return norm2(out);
}; };
RealD Mdag (const CoarseVector &in, CoarseVector &out){ // TODO: get this correct RealD Mdag (const CoarseVector &in, CoarseVector &out){
return M(in,out); return M(in,out);
}; };
void Mdir(const CoarseVector &in, CoarseVector &out, int dir, int disp) { // Defer support for further coarsening for now
void Mdiag (const CoarseVector &in, CoarseVector &out){};
conformable(_grid,in._grid); void Mdir (const CoarseVector &in, CoarseVector &out,int dir, int disp){};
conformable(in._grid,out._grid);
SimpleCompressor<siteVector> compressor;
Stencil.HaloExchange(in,compressor);
auto point = [dir, disp](){
if(dir == 0 and disp == 0)
return 8;
else
return (4 * dir + 1 - disp) / 2;
}();
parallel_for(int ss=0;ss<Grid()->oSites();ss++){
siteVector res = zero;
siteVector nbr;
int ptype;
StencilEntry *SE;
SE=Stencil.GetEntry(ptype,point,ss);
if(SE->_is_local&&SE->_permute) {
permute(nbr,in._odata[SE->_offset],ptype);
} else if(SE->_is_local) {
nbr = in._odata[SE->_offset];
} else {
nbr = Stencil.CommBuf()[SE->_offset];
}
res = res + A[point]._odata[ss]*nbr;
vstream(out._odata[ss],res);
}
};
void Mdiag(const CoarseVector &in, CoarseVector &out) {
Mdir(in, out, 0, 0); // use the self coupling (= last) point of the stencil
};
CoarsenedMatrix(GridCartesian &CoarseGrid) : CoarsenedMatrix(GridCartesian &CoarseGrid) :
@ -372,9 +327,10 @@ namespace Grid {
CoarseVector iProj(Grid()); CoarseVector iProj(Grid());
CoarseVector oProj(Grid()); CoarseVector oProj(Grid());
CoarseScalar InnerProd(Grid());
// Orthogonalise the subblocks over the basis // Orthogonalise the subblocks over the basis
Subspace.Orthogonalise(); blockOrthogonalise(InnerProd,Subspace.subspace);
// Compute the matrix elements of linop between this orthonormal // Compute the matrix elements of linop between this orthonormal
// set of vectors. // set of vectors.
@ -431,9 +387,9 @@ namespace Grid {
parallel_for(int ss=0;ss<Grid()->oSites();ss++){ parallel_for(int ss=0;ss<Grid()->oSites();ss++){
for(int j=0;j<nbasis;j++){ for(int j=0;j<nbasis;j++){
if( disp!= 0 ) { if( disp!= 0 ) {
A[p]._odata[ss]()(0,0)(j,i) = oProj._odata[ss]()(0)(j); A[p]._odata[ss](j,i) = oProj._odata[ss](j);
} }
A[self_stencil]._odata[ss]()(0,0)(j,i) = A[self_stencil]._odata[ss]()(0,0)(j,i) + iProj._odata[ss]()(0)(j); A[self_stencil]._odata[ss](j,i) = A[self_stencil]._odata[ss](j,i) + iProj._odata[ss](j);
} }
} }
} }
@ -461,7 +417,7 @@ namespace Grid {
std::cout<<GridLogMessage<<"Computed Coarse Operator"<<std::endl; std::cout<<GridLogMessage<<"Computed Coarse Operator"<<std::endl;
#endif #endif
// ForceHermitian(); // ForceHermitian();
// AssertHermitian(); AssertHermitian();
// ForceDiagonal(); // ForceDiagonal();
} }
void ForceDiagonal(void) { void ForceDiagonal(void) {

View File

@ -80,8 +80,8 @@ inline void subdivides(GridBase *coarse,GridBase *fine)
} }
template<class vobj,class vobjC> template<class vobj,class CComplex,int nbasis>
inline void blockProject(Lattice<vobjC> &coarseData, inline void blockProject(Lattice<iVector<CComplex,nbasis > > &coarseData,
const Lattice<vobj> &fineData, const Lattice<vobj> &fineData,
const std::vector<Lattice<vobj> > &Basis) const std::vector<Lattice<vobj> > &Basis)
{ {
@ -90,8 +90,7 @@ inline void blockProject(Lattice<vobjC> &coarseData,
int _ndimension = coarse->_ndimension; int _ndimension = coarse->_ndimension;
// checks // checks
assert((Basis.size() != 0) && ((Basis.size() & 0x1) == 0)); assert( nbasis == Basis.size() );
auto nbasis = Basis.size();
subdivides(coarse,fine); subdivides(coarse,fine);
for(int i=0;i<nbasis;i++){ for(int i=0;i<nbasis;i++){
conformable(Basis[i],fineData); conformable(Basis[i],fineData);
@ -119,8 +118,8 @@ inline void blockProject(Lattice<vobjC> &coarseData,
PARALLEL_CRITICAL PARALLEL_CRITICAL
for(int i=0;i<nbasis;i++) { for(int i=0;i<nbasis;i++) {
coarseData._odata[sc]()(0)(i)=coarseData._odata[sc]()(0)(i) coarseData._odata[sc](i)=coarseData._odata[sc](i)
+ TensorRemove(innerProduct(Basis[i]._odata[sf],fineData._odata[sf])); + innerProduct(Basis[i]._odata[sf],fineData._odata[sf]);
} }
} }
@ -286,9 +285,9 @@ inline void blockOrthogonalise(Lattice<CComplex> &ip,std::vector<Lattice<vobj> >
} }
} }
template<class vobj,class vobjC> template<class vobj,class CComplex,int nbasis>
inline void blockPromote(const Lattice<vobjC> &coarseData, inline void blockPromote(const Lattice<iVector<CComplex,nbasis > > &coarseData,
Lattice<vobj> &fineData, Lattice<vobj> &fineData,
const std::vector<Lattice<vobj> > &Basis) const std::vector<Lattice<vobj> > &Basis)
{ {
GridBase * fine = fineData._grid; GridBase * fine = fineData._grid;
@ -296,9 +295,7 @@ inline void blockPromote(const Lattice<vobjC> &coarseData,
int _ndimension = coarse->_ndimension; int _ndimension = coarse->_ndimension;
// checks // checks
assert((Basis.size() != 0) && ((Basis.size() & 0x1) == 0)); assert( nbasis == Basis.size() );
auto nbasis = Basis.size();
subdivides(coarse,fine); subdivides(coarse,fine);
for(int i=0;i<nbasis;i++){ for(int i=0;i<nbasis;i++){
conformable(Basis[i]._grid,fine); conformable(Basis[i]._grid,fine);
@ -322,13 +319,9 @@ inline void blockPromote(const Lattice<vobjC> &coarseData,
for(int d=0;d<_ndimension;d++) coor_c[d]=coor_f[d]/block_r[d]; for(int d=0;d<_ndimension;d++) coor_c[d]=coor_f[d]/block_r[d];
Lexicographic::IndexFromCoor(coor_c,sc,coarse->_rdimensions); Lexicographic::IndexFromCoor(coor_c,sc,coarse->_rdimensions);
// The temporary is necessary, since a pure instance of Grid::simd<...> is
// not a valid argument to operator+ with an iVector, we need an an iScalar
typename vobjC::tensor_reduced tmp; // iScalar<iVector<iVector<...>>> -> iScalar<iScalar<iScalar<...>>>
for(int i=0;i<nbasis;i++) { for(int i=0;i<nbasis;i++) {
tmp = coarseData._odata[sc]()(0)(i); if(i==0) fineData._odata[sf]=coarseData._odata[sc](i) * Basis[i]._odata[sf];
if(i==0) fineData._odata[sf] = tmp * Basis[i]._odata[sf]; else fineData._odata[sf]=fineData._odata[sf]+coarseData._odata[sc](i)*Basis[i]._odata[sf];
else fineData._odata[sf]=fineData._odata[sf]+tmp*Basis[i]._odata[sf];
} }
} }
} }