mirror of
https://github.com/paboyle/Grid.git
synced 2025-06-18 07:47:06 +01:00
Imported changes from feature/gparity_HMC branch:
Added storage of final true residual in mixed-prec CG and enhanced log output Fixed const correctness of multi-shift constructor Added a mixed precision variant of the multi-shift algorithm that uses a single precision operator and applies periodic reliable update to the residual Added tests/solver/Test_dwf_multishift_mixedprec to test the above Fixed local coherence lanczos using the (large!) max approx to the chebyshev eval as the scale from which to judge the quality of convergence, resulting a test that always passes Added a method to local coherence lanczos class that returns the fine eval/evec pair Added iterative log output to power method Added optional disabling of the plaquette check in Nerscio to support loading old G-parity configs which have a factor of 2 error in the plaquette G-parity Dirac op no longer allows GPBC in the time direction; instead we toggle between periodic and antiperiodic Replaced thread_for G-parity 5D force insertion implementation with accelerator_for version capable of running on GPUs Generalized tests/lanczos/Test_dwf_lanczos to support regular DWF as well as Gparity, with the action chosen by a command line option Modified tests/forces/Test_dwf_gpforce,Test_gpdwf_force,Test_gpwilson_force to use GPBC a spatial direction rather than the t-direction, and antiperiodic BCs for time direction tests/core/Test_gparity now supports using APBC in time direction using command line toggle
This commit is contained in:
@ -55,13 +55,17 @@ static_assert(same_vComplex == 1, "Dirac Operators must have same underlying SIM
|
||||
int main (int argc, char ** argv)
|
||||
{
|
||||
int nu = 0;
|
||||
|
||||
int tbc_aprd = 0; //use antiperiodic BCs in the time direction?
|
||||
|
||||
Grid_init(&argc,&argv);
|
||||
|
||||
for(int i=1;i<argc;i++){
|
||||
if(std::string(argv[i]) == "--Gparity-dir"){
|
||||
std::stringstream ss; ss << argv[i+1]; ss >> nu;
|
||||
std::cout << GridLogMessage << "Set Gparity direction to " << nu << std::endl;
|
||||
}else if(std::string(argv[i]) == "--Tbc-APRD"){
|
||||
tbc_aprd = 1;
|
||||
std::cout << GridLogMessage << "Using antiperiodic BCs in the time direction" << std::endl;
|
||||
}
|
||||
}
|
||||
|
||||
@ -155,13 +159,18 @@ int main (int argc, char ** argv)
|
||||
|
||||
//Coordinate grid for reference
|
||||
LatticeInteger xcoor_1f5(FGrid_1f);
|
||||
LatticeCoordinate(xcoor_1f5,1+nu);
|
||||
LatticeCoordinate(xcoor_1f5,1+nu); //note '1+nu'! This is because for 5D fields the s-direction is direction 0
|
||||
Replicate(src,src_1f);
|
||||
src_1f = where( xcoor_1f5 >= Integer(L), 2.0*src_1f,src_1f );
|
||||
|
||||
RealD mass=0.0;
|
||||
RealD M5=1.8;
|
||||
StandardDiracOp Ddwf(Umu_1f,*FGrid_1f,*FrbGrid_1f,*UGrid_1f,*UrbGrid_1f,mass,M5 DOP_PARAMS);
|
||||
|
||||
//Standard Dirac op
|
||||
AcceleratorVector<Complex,4> bc_std(Nd, 1.0);
|
||||
if(tbc_aprd) bc_std[Nd-1] = -1.; //antiperiodic time BC
|
||||
StandardDiracOp::ImplParams std_params(bc_std);
|
||||
StandardDiracOp Ddwf(Umu_1f,*FGrid_1f,*FrbGrid_1f,*UGrid_1f,*UrbGrid_1f,mass,M5 DOP_PARAMS, std_params);
|
||||
|
||||
StandardFermionField src_o_1f(FrbGrid_1f);
|
||||
StandardFermionField result_o_1f(FrbGrid_1f);
|
||||
@ -172,9 +181,11 @@ int main (int argc, char ** argv)
|
||||
ConjugateGradient<StandardFermionField> CG(1.0e-8,10000);
|
||||
CG(HermOpEO,src_o_1f,result_o_1f);
|
||||
|
||||
// const int nu = 3;
|
||||
//Gparity Dirac op
|
||||
std::vector<int> twists(Nd,0);
|
||||
twists[nu] = 1;
|
||||
if(tbc_aprd) twists[Nd-1] = 1;
|
||||
|
||||
GparityDiracOp::ImplParams params;
|
||||
params.twists = twists;
|
||||
GparityDiracOp GPDdwf(Umu_2f,*FGrid_2f,*FrbGrid_2f,*UGrid_2f,*UrbGrid_2f,mass,M5 DOP_PARAMS,params);
|
||||
@ -271,8 +282,11 @@ int main (int argc, char ** argv)
|
||||
std::cout << "2f cb "<<result_o_2f.Checkerboard()<<std::endl;
|
||||
std::cout << "1f cb "<<result_o_1f.Checkerboard()<<std::endl;
|
||||
|
||||
std::cout << " result norms " <<norm2(result_o_2f)<<" " <<norm2(result_o_1f)<<std::endl;
|
||||
//Compare norms
|
||||
std::cout << " result norms 2f: " <<norm2(result_o_2f)<<" 1f: " <<norm2(result_o_1f)<<std::endl;
|
||||
|
||||
|
||||
//Take the 2f solution and convert into the corresponding 1f solution (odd cb only)
|
||||
StandardFermionField res0o (FrbGrid_2f);
|
||||
StandardFermionField res1o (FrbGrid_2f);
|
||||
StandardFermionField res0 (FGrid_2f);
|
||||
@ -281,14 +295,15 @@ int main (int argc, char ** argv)
|
||||
res0=Zero();
|
||||
res1=Zero();
|
||||
|
||||
res0o = PeekIndex<0>(result_o_2f,0);
|
||||
res1o = PeekIndex<0>(result_o_2f,1);
|
||||
res0o = PeekIndex<0>(result_o_2f,0); //flavor 0, odd cb
|
||||
res1o = PeekIndex<0>(result_o_2f,1); //flavor 1, odd cb
|
||||
|
||||
std::cout << "res cb "<<res0o.Checkerboard()<<std::endl;
|
||||
std::cout << "res cb "<<res1o.Checkerboard()<<std::endl;
|
||||
|
||||
setCheckerboard(res0,res0o);
|
||||
setCheckerboard(res1,res1o);
|
||||
//poke odd onto non-cb field
|
||||
setCheckerboard(res0,res0o);
|
||||
setCheckerboard(res1,res1o);
|
||||
|
||||
StandardFermionField replica (FGrid_1f);
|
||||
StandardFermionField replica0(FGrid_1f);
|
||||
@ -296,12 +311,13 @@ int main (int argc, char ** argv)
|
||||
Replicate(res0,replica0);
|
||||
Replicate(res1,replica1);
|
||||
|
||||
//2nd half of doubled lattice has f=1
|
||||
replica = where( xcoor_1f5 >= Integer(L), replica1,replica0 );
|
||||
|
||||
replica0 = Zero();
|
||||
setCheckerboard(replica0,result_o_1f);
|
||||
|
||||
std::cout << "Norm2 solutions is " <<norm2(replica)<<" "<< norm2(replica0)<<std::endl;
|
||||
std::cout << "Norm2 solutions 1f reconstructed from 2f: " <<norm2(replica)<<" Actual 1f: "<< norm2(replica0)<<std::endl;
|
||||
|
||||
replica = replica - replica0;
|
||||
|
||||
|
@ -71,26 +71,14 @@ int main (int argc, char ** argv)
|
||||
////////////////////////////////////
|
||||
RealD mass=0.2; //kills the diagonal term
|
||||
RealD M5=1.8;
|
||||
// const int nu = 3;
|
||||
// std::vector<int> twists(Nd,0); // twists[nu] = 1;
|
||||
// GparityDomainWallFermionR::ImplParams params; params.twists = twists;
|
||||
// GparityDomainWallFermionR Ddwf(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,params);
|
||||
|
||||
// DomainWallFermionR Dw (U, Grid,RBGrid,mass,M5);
|
||||
|
||||
const int nu = 3;
|
||||
const int nu = 0; //gparity direction
|
||||
std::vector<int> twists(Nd,0);
|
||||
twists[nu] = 1;
|
||||
twists[Nd-1] = 1; //antiperiodic in time
|
||||
GparityDomainWallFermionR::ImplParams params;
|
||||
params.twists = twists;
|
||||
|
||||
/*
|
||||
params.boundary_phases[0] = 1.0;
|
||||
params.boundary_phases[1] = 1.0;
|
||||
params.boundary_phases[2] = 1.0;
|
||||
params.boundary_phases[3] =- 1.0;
|
||||
*/
|
||||
|
||||
|
||||
GparityDomainWallFermionR Dw(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,params);
|
||||
|
||||
Dw.M (phi,Mphi);
|
||||
|
@ -71,8 +71,10 @@ int main (int argc, char ** argv)
|
||||
RealD mass=0.01;
|
||||
RealD M5=1.8;
|
||||
|
||||
const int nu = 3;
|
||||
std::vector<int> twists(Nd,0); twists[nu] = 1;
|
||||
const int nu = 1;
|
||||
std::vector<int> twists(Nd,0);
|
||||
twists[nu] = 1;
|
||||
twists[3] = 1;
|
||||
GparityDomainWallFermionR::ImplParams params; params.twists = twists;
|
||||
GparityDomainWallFermionR Ddwf(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,params);
|
||||
Ddwf.M (phi,Mphi);
|
||||
|
@ -64,8 +64,12 @@ int main (int argc, char ** argv)
|
||||
////////////////////////////////////
|
||||
RealD mass=0.01;
|
||||
|
||||
const int nu = 3;
|
||||
std::vector<int> twists(Nd,0); twists[nu] = 1;
|
||||
const int nu = 1;
|
||||
const int Lnu=latt_size[nu];
|
||||
|
||||
std::vector<int> twists(Nd,0);
|
||||
twists[nu] = 1;
|
||||
twists[3]=1;
|
||||
GparityWilsonFermionR::ImplParams params; params.twists = twists;
|
||||
GparityWilsonFermionR Wil(U,*UGrid,*UrbGrid,mass,params);
|
||||
Wil.M (phi,Mphi);
|
||||
|
@ -31,14 +31,38 @@ using namespace std;
|
||||
using namespace Grid;
|
||||
;
|
||||
|
||||
typedef typename GparityDomainWallFermionR::FermionField FermionField;
|
||||
template<typename Action>
|
||||
struct Setup{};
|
||||
|
||||
RealD AllZero(RealD x){ return 0.;}
|
||||
template<>
|
||||
struct Setup<GparityMobiusFermionR>{
|
||||
static GparityMobiusFermionR* getAction(LatticeGaugeField &Umu,
|
||||
GridCartesian* FGrid, GridRedBlackCartesian* FrbGrid, GridCartesian* UGrid, GridRedBlackCartesian* UrbGrid){
|
||||
RealD mass=0.01;
|
||||
RealD M5=1.8;
|
||||
RealD mob_b=1.5;
|
||||
GparityMobiusFermionD ::ImplParams params;
|
||||
std::vector<int> twists({1,1,1,0});
|
||||
params.twists = twists;
|
||||
return new GparityMobiusFermionR(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,mob_b,mob_b-1.,params);
|
||||
}
|
||||
};
|
||||
|
||||
int main (int argc, char ** argv)
|
||||
{
|
||||
Grid_init(&argc,&argv);
|
||||
template<>
|
||||
struct Setup<DomainWallFermionR>{
|
||||
static DomainWallFermionR* getAction(LatticeGaugeField &Umu,
|
||||
GridCartesian* FGrid, GridRedBlackCartesian* FrbGrid, GridCartesian* UGrid, GridRedBlackCartesian* UrbGrid){
|
||||
RealD mass=0.01;
|
||||
RealD M5=1.8;
|
||||
return new DomainWallFermionR(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5);
|
||||
}
|
||||
};
|
||||
|
||||
|
||||
|
||||
template<typename Action>
|
||||
void run(){
|
||||
typedef typename Action::FermionField FermionField;
|
||||
const int Ls=8;
|
||||
|
||||
GridCartesian * UGrid = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi());
|
||||
@ -56,24 +80,10 @@ int main (int argc, char ** argv)
|
||||
LatticeGaugeField Umu(UGrid);
|
||||
SU<Nc>::HotConfiguration(RNG4, Umu);
|
||||
|
||||
std::vector<LatticeColourMatrix> U(4,UGrid);
|
||||
for(int mu=0;mu<Nd;mu++){
|
||||
U[mu] = PeekIndex<LorentzIndex>(Umu,mu);
|
||||
}
|
||||
|
||||
RealD mass=0.01;
|
||||
RealD M5=1.8;
|
||||
RealD mob_b=1.5;
|
||||
// DomainWallFermionR Ddwf(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5);
|
||||
GparityMobiusFermionD ::ImplParams params;
|
||||
std::vector<int> twists({1,1,1,0});
|
||||
params.twists = twists;
|
||||
GparityMobiusFermionR Ddwf(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,mob_b,mob_b-1.,params);
|
||||
|
||||
// MdagMLinearOperator<DomainWallFermionR,LatticeFermion> HermOp(Ddwf);
|
||||
// SchurDiagTwoOperator<DomainWallFermionR,LatticeFermion> HermOp(Ddwf);
|
||||
SchurDiagTwoOperator<GparityMobiusFermionR,FermionField> HermOp(Ddwf);
|
||||
// SchurDiagMooeeOperator<DomainWallFermionR,LatticeFermion> HermOp(Ddwf);
|
||||
Action *action = Setup<Action>::getAction(Umu,FGrid,FrbGrid,UGrid,UrbGrid);
|
||||
|
||||
//MdagMLinearOperator<Action,FermionField> HermOp(Ddwf);
|
||||
SchurDiagTwoOperator<Action,FermionField> HermOp(*action);
|
||||
|
||||
const int Nstop = 30;
|
||||
const int Nk = 40;
|
||||
@ -90,8 +100,7 @@ int main (int argc, char ** argv)
|
||||
PlainHermOp<FermionField> Op (HermOp);
|
||||
|
||||
ImplicitlyRestartedLanczos<FermionField> IRL(OpCheby,Op,Nstop,Nk,Nm,resid,MaxIt);
|
||||
|
||||
|
||||
|
||||
std::vector<RealD> eval(Nm);
|
||||
FermionField src(FrbGrid);
|
||||
gaussian(RNG5rb,src);
|
||||
@ -103,6 +112,28 @@ int main (int argc, char ** argv)
|
||||
int Nconv;
|
||||
IRL.calc(eval,evec,src,Nconv);
|
||||
|
||||
delete action;
|
||||
}
|
||||
|
||||
int main (int argc, char ** argv)
|
||||
{
|
||||
Grid_init(&argc,&argv);
|
||||
|
||||
std::string action = "GparityMobius";
|
||||
for(int i=1;i<argc;i++){
|
||||
if(std::string(argv[i]) == "-action"){
|
||||
action = argv[i+1];
|
||||
}
|
||||
}
|
||||
|
||||
if(action == "GparityMobius"){
|
||||
run<GparityMobiusFermionR>();
|
||||
}else if(action == "DWF"){
|
||||
run<DomainWallFermionR>();
|
||||
}else{
|
||||
std::cout << "Unknown action" << std::endl;
|
||||
exit(1);
|
||||
}
|
||||
|
||||
Grid_finalize();
|
||||
}
|
||||
|
184
tests/solver/Test_dwf_multishift_mixedprec.cc
Normal file
184
tests/solver/Test_dwf_multishift_mixedprec.cc
Normal file
@ -0,0 +1,184 @@
|
||||
/*************************************************************************************
|
||||
|
||||
Grid physics library, www.github.com/paboyle/Grid
|
||||
|
||||
Source file: ./tests/Test_dwf_multishift_mixedprec.cc
|
||||
|
||||
Copyright (C) 2015
|
||||
|
||||
Author: Christopher Kelly <ckelly@bnl.gov>
|
||||
|
||||
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 */
|
||||
#include <Grid/Grid.h>
|
||||
|
||||
using namespace Grid;
|
||||
|
||||
template<typename SpeciesD, typename SpeciesF, typename GaugeStatisticsType>
|
||||
void run_test(int argc, char ** argv, const typename SpeciesD::ImplParams ¶ms){
|
||||
const int Ls = 16;
|
||||
GridCartesian* UGrid_d = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd, vComplexD::Nsimd()), GridDefaultMpi());
|
||||
GridRedBlackCartesian* UrbGrid_d = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid_d);
|
||||
GridCartesian* FGrid_d = SpaceTimeGrid::makeFiveDimGrid(Ls, UGrid_d);
|
||||
GridRedBlackCartesian* FrbGrid_d = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls, UGrid_d);
|
||||
|
||||
GridCartesian* UGrid_f = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd, vComplexF::Nsimd()), GridDefaultMpi());
|
||||
GridRedBlackCartesian* UrbGrid_f = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid_f);
|
||||
GridCartesian* FGrid_f = SpaceTimeGrid::makeFiveDimGrid(Ls, UGrid_f);
|
||||
GridRedBlackCartesian* FrbGrid_f = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls, UGrid_f);
|
||||
|
||||
typedef typename SpeciesD::FermionField FermionFieldD;
|
||||
typedef typename SpeciesF::FermionField FermionFieldF;
|
||||
|
||||
std::vector<int> seeds4({1, 2, 3, 4});
|
||||
std::vector<int> seeds5({5, 6, 7, 8});
|
||||
GridParallelRNG RNG5(FGrid_d);
|
||||
RNG5.SeedFixedIntegers(seeds5);
|
||||
GridParallelRNG RNG4(UGrid_d);
|
||||
RNG4.SeedFixedIntegers(seeds4);
|
||||
|
||||
FermionFieldD src_d(FGrid_d);
|
||||
random(RNG5, src_d);
|
||||
|
||||
LatticeGaugeFieldD Umu_d(UGrid_d);
|
||||
|
||||
//CPS-created G-parity ensembles have a factor of 2 error in the plaquette that causes the read to fail unless we workaround it
|
||||
bool gparity_plaquette_fix = false;
|
||||
for(int i=1;i<argc;i++){
|
||||
if(std::string(argv[i]) == "--gparity_plaquette_fix"){
|
||||
gparity_plaquette_fix=true;
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
bool cfg_loaded=false;
|
||||
for(int i=1;i<argc;i++){
|
||||
if(std::string(argv[i]) == "--load_config"){
|
||||
assert(i != argc-1);
|
||||
std::string file = argv[i+1];
|
||||
NerscIO io;
|
||||
FieldMetaData metadata;
|
||||
|
||||
if(gparity_plaquette_fix) NerscIO::exitOnReadPlaquetteMismatch() = false;
|
||||
|
||||
io.readConfiguration<GaugeStatisticsType>(Umu_d, metadata, file);
|
||||
|
||||
if(gparity_plaquette_fix){
|
||||
metadata.plaquette *= 2.; //correct header value
|
||||
|
||||
//Get the true plaquette
|
||||
FieldMetaData tmp;
|
||||
GaugeStatisticsType gs; gs(Umu_d, tmp);
|
||||
|
||||
std::cout << "After correction: plaqs " << tmp.plaquette << " " << metadata.plaquette << std::endl;
|
||||
assert(fabs(tmp.plaquette -metadata.plaquette ) < 1.0e-5 );
|
||||
}
|
||||
|
||||
cfg_loaded=true;
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if(!cfg_loaded)
|
||||
SU<Nc>::HotConfiguration(RNG4, Umu_d);
|
||||
|
||||
LatticeGaugeFieldF Umu_f(UGrid_f);
|
||||
precisionChange(Umu_f, Umu_d);
|
||||
|
||||
std::cout << GridLogMessage << "Lattice dimensions: " << GridDefaultLatt() << " Ls: " << Ls << std::endl;
|
||||
|
||||
RealD mass = 0.01;
|
||||
RealD M5 = 1.8;
|
||||
SpeciesD Ddwf_d(Umu_d, *FGrid_d, *FrbGrid_d, *UGrid_d, *UrbGrid_d, mass, M5, params);
|
||||
SpeciesF Ddwf_f(Umu_f, *FGrid_f, *FrbGrid_f, *UGrid_f, *UrbGrid_f, mass, M5, params);
|
||||
|
||||
FermionFieldD src_o_d(FrbGrid_d);
|
||||
pickCheckerboard(Odd, src_o_d, src_d);
|
||||
|
||||
SchurDiagMooeeOperator<SpeciesD, FermionFieldD> HermOpEO_d(Ddwf_d);
|
||||
SchurDiagMooeeOperator<SpeciesF, FermionFieldF> HermOpEO_f(Ddwf_f);
|
||||
|
||||
AlgRemez remez(1e-4, 64, 50);
|
||||
int order = 15;
|
||||
remez.generateApprox(order, 1, 2); //sqrt
|
||||
|
||||
MultiShiftFunction shifts(remez, 1e-10, false);
|
||||
|
||||
int relup_freq = 50;
|
||||
double t1=usecond();
|
||||
ConjugateGradientMultiShiftMixedPrec<FermionFieldD,FermionFieldF> mcg(10000, shifts, FrbGrid_f, HermOpEO_f, relup_freq);
|
||||
|
||||
std::vector<FermionFieldD> results_o_d(order, FrbGrid_d);
|
||||
mcg(HermOpEO_d, src_o_d, results_o_d);
|
||||
double t2=usecond();
|
||||
|
||||
//Crosscheck double and mixed prec results
|
||||
ConjugateGradientMultiShift<FermionFieldD> dmcg(10000, shifts);
|
||||
std::vector<FermionFieldD> results_o_d_2(order, FrbGrid_d);
|
||||
dmcg(HermOpEO_d, src_o_d, results_o_d_2);
|
||||
double t3=usecond();
|
||||
|
||||
std::cout << GridLogMessage << "Comparison of mixed prec results to double prec results |mixed - double|^2 :" << std::endl;
|
||||
FermionFieldD tmp(FrbGrid_d);
|
||||
for(int i=0;i<order;i++){
|
||||
RealD ndiff = axpy_norm(tmp, -1., results_o_d[i], results_o_d_2[i]);
|
||||
std::cout << i << " " << ndiff << std::endl;
|
||||
}
|
||||
|
||||
std::cout<<GridLogMessage << "Mixed precision algorithm: Total usec = "<< (t2-t1)<<std::endl;
|
||||
std::cout<<GridLogMessage << "Double precision algorithm: Total usec = "<< (t3-t2)<<std::endl;
|
||||
}
|
||||
|
||||
|
||||
|
||||
|
||||
|
||||
int main (int argc, char ** argv)
|
||||
{
|
||||
Grid_init(&argc, &argv);
|
||||
|
||||
bool gparity = false;
|
||||
int gpdir;
|
||||
|
||||
for(int i=1;i<argc;i++){
|
||||
std::string arg(argv[i]);
|
||||
if(arg == "--Gparity"){
|
||||
assert(i!=argc-1);
|
||||
gpdir = std::stoi(argv[i+1]);
|
||||
assert(gpdir >= 0 && gpdir <= 2); //spatial!
|
||||
gparity = true;
|
||||
}
|
||||
}
|
||||
if(gparity){
|
||||
std::cout << "Running test with G-parity BCs in " << gpdir << " direction" << std::endl;
|
||||
GparityWilsonImplParams params;
|
||||
params.twists[gpdir] = 1;
|
||||
|
||||
std::vector<int> conj_dirs(Nd,0);
|
||||
conj_dirs[gpdir] = 1;
|
||||
ConjugateGimplD::setDirections(conj_dirs);
|
||||
|
||||
run_test<GparityDomainWallFermionD, GparityDomainWallFermionF, ConjugateGaugeStatistics>(argc,argv,params);
|
||||
}else{
|
||||
std::cout << "Running test with periodic BCs" << std::endl;
|
||||
WilsonImplParams params;
|
||||
run_test<DomainWallFermionD, DomainWallFermionF, PeriodicGaugeStatistics>(argc,argv,params);
|
||||
}
|
||||
|
||||
Grid_finalize();
|
||||
}
|
Reference in New Issue
Block a user