Grid 0.7.0
DomainDecomposedBoundaryTwoFlavourRatioPseudoFermion.h
Go to the documentation of this file.
1/*************************************************************************************
2
3 Grid physics library, www.github.com/paboyle/Grid
4
5 Source file: ./lib/qcd/action/pseudofermion/DomainDecomposedTwoFlavourBoundary.h
6
7 Copyright (C) 2021
8
9Author: Peter Boyle <paboyle@ph.ed.ac.uk>
10
11 This program is free software; you can redistribute it and/or modify
12 it under the terms of the GNU General Public License as published by
13 the Free Software Foundation; either version 2 of the License, or
14 (at your option) any later version.
15
16 This program is distributed in the hope that it will be useful,
17 but WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 GNU General Public License for more details.
20
21 You should have received a copy of the GNU General Public License along
22 with this program; if not, write to the Free Software Foundation, Inc.,
23 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24
25 See the full license in the file "LICENSE" in the top level distribution directory
26*************************************************************************************/
27/* END LEGAL */
28#pragma once
29
31
33// Two flavour ratio
35template<class ImplD,class ImplF>
36class DomainDecomposedBoundaryTwoFlavourRatioPseudoFermion : public Action<typename ImplD::GaugeField> {
37public:
39
40private:
41 SchurFactoredFermionOperator<ImplD,ImplF> & NumOp;// the basic operator
42 SchurFactoredFermionOperator<ImplD,ImplF> & DenOp;// the basic operator
43
47
48 FermionField Phi; // the pseudo fermion field for this trajectory
49
50public:
51 DomainDecomposedBoundaryTwoFlavourRatioPseudoFermion(SchurFactoredFermionOperator<ImplD,ImplF> &_NumOp,
52 SchurFactoredFermionOperator<ImplD,ImplF> &_DenOp,
53 RealD _DerivativeTol, RealD _ActionTol, RealD _InnerTol=1.0e-6)
54 : NumOp(_NumOp), DenOp(_DenOp),
55 Phi(_NumOp.PeriodicFermOpD.FermionGrid()),
56 InnerStoppingCondition(_InnerTol),
57 DerivativeStoppingCondition(_DerivativeTol),
58 ActionStoppingCondition(_ActionTol)
59 {};
60
61 virtual std::string action_name(){return "DomainDecomposedBoundaryTwoFlavourRatioPseudoFermion";}
62
63 virtual std::string LogParameters(){
64 std::stringstream sstream;
65 return sstream.str();
66 }
67
68 virtual void refresh(const GaugeField &U, GridSerialRNG& sRNG, GridParallelRNG& pRNG)
69 {
70 NumOp.ImportGauge(U);
71 DenOp.ImportGauge(U);
72
73 FermionField eta(NumOp.PeriodicFermOpD.FermionGrid());
74 FermionField tmp(NumOp.PeriodicFermOpD.FermionGrid());
75
76 // P(phi) = e^{- phi^dag P^dag Rdag^-1 R^-1 P phi}
77 //
78 // NumOp == P
79 // DenOp == R
80 //
81 // Take phi = P^{-1} R eta ; eta = R^-1 P Phi
82 //
83 // P(eta) = e^{- eta^dag eta}
84 //
85 // e^{x^2/2 sig^2} => sig^2 = 0.5.
86 //
87 // So eta should be of width sig = 1/sqrt(2) and must multiply by 0.707....
88 //
89 RealD scale = std::sqrt(0.5);
90
91 gaussian(pRNG,eta); eta=eta*scale;
92
93 NumOp.ProjectBoundaryBar(eta);
98 DenOp.R(eta,tmp);
99 NumOp.RInv(tmp,Phi);
100 DumpSliceNorm("Phi",Phi);
101
102 };
103
105 // S = phi^dag Pdag Rdag^-1 R^-1 P phi
107 virtual RealD S(const GaugeField &U) {
108
109 NumOp.ImportGauge(U);
110 DenOp.ImportGauge(U);
111
112 FermionField X(NumOp.PeriodicFermOpD.FermionGrid());
113 FermionField Y(NumOp.PeriodicFermOpD.FermionGrid());
114
119 NumOp.R(Phi,Y);
120 DenOp.RInv(Y,X);
121
122 RealD action = norm2(X);
123 // std::cout << " DD boundary action is " <<action<<std::endl;
124
125 return action;
126 };
127
128 virtual void deriv(const GaugeField &U,GaugeField & dSdU)
129 {
130 NumOp.ImportGauge(U);
131 DenOp.ImportGauge(U);
132
133 GridBase *fgrid = NumOp.PeriodicFermOpD.FermionGrid();
134 GridBase *ugrid = NumOp.PeriodicFermOpD.GaugeGrid();
135
136 FermionField X(fgrid);
137 FermionField Y(fgrid);
138 FermionField tmp(fgrid);
139
140 GaugeField force(ugrid);
141
142 FermionField DobiDdbPhi(fgrid); // Vector A in my notes
143 FermionField DoiDdDobiDdbPhi(fgrid); // Vector B in my notes
144 FermionField DiDdbP_Phi(fgrid); // Vector C in my notes
145 FermionField DidRinvP_Phi(fgrid); // Vector D in my notes
146 FermionField DdbdDidRinvP_Phi(fgrid);
147 FermionField DoidRinvDagRinvP_Phi(fgrid); // Vector E in my notes
148 FermionField DobidDddDoidRinvDagRinvP_Phi(fgrid); // Vector F in my notes
149
150 FermionField P_Phi(fgrid);
151 FermionField RinvP_Phi(fgrid);
152 FermionField RinvDagRinvP_Phi(fgrid);
153 FermionField PdagRinvDagRinvP_Phi(fgrid);
154
155 // RealD action = S(U);
160
161 // P term
162 NumOp.dBoundaryBar(Phi,tmp);
163 NumOp.dOmegaBarInv(tmp,DobiDdbPhi); // Vector A
164 NumOp.dBoundary(DobiDdbPhi,tmp);
165 NumOp.dOmegaInv(tmp,DoiDdDobiDdbPhi); // Vector B
166 P_Phi = Phi - DoiDdDobiDdbPhi;
167 NumOp.ProjectBoundaryBar(P_Phi);
168
169 // R^-1 P term
170 DenOp.dBoundaryBar(P_Phi,tmp);
171 DenOp.Dinverse(tmp,DiDdbP_Phi); // Vector C
172 RinvP_Phi = P_Phi - DiDdbP_Phi;
173 DenOp.ProjectBoundaryBar(RinvP_Phi); // Correct to here
174
175
176 // R^-dagger R^-1 P term
177 DenOp.DinverseDag(RinvP_Phi,DidRinvP_Phi); // Vector D
178 DenOp.dBoundaryBarDag(DidRinvP_Phi,DdbdDidRinvP_Phi);
179 RinvDagRinvP_Phi = RinvP_Phi - DdbdDidRinvP_Phi;
180 DenOp.ProjectBoundaryBar(RinvDagRinvP_Phi);
181
182
183 // P^dag R^-dagger R^-1 P term
184 NumOp.dOmegaDagInv(RinvDagRinvP_Phi,DoidRinvDagRinvP_Phi); // Vector E
185 NumOp.dBoundaryDag(DoidRinvDagRinvP_Phi,tmp);
186 NumOp.dOmegaBarDagInv(tmp,DobidDddDoidRinvDagRinvP_Phi); // Vector F
187 NumOp.dBoundaryBarDag(DobidDddDoidRinvDagRinvP_Phi,tmp);
188 PdagRinvDagRinvP_Phi = RinvDagRinvP_Phi- tmp;
189 NumOp.ProjectBoundaryBar(PdagRinvDagRinvP_Phi);
190
191 /*
192 std::cout << "S eval "<< action << std::endl;
193 std::cout << "S - IP1 "<< innerProduct(Phi,PdagRinvDagRinvP_Phi) << std::endl;
194 std::cout << "S - IP2 "<< norm2(RinvP_Phi) << std::endl;
195
196 NumOp.R(Phi,tmp);
197 tmp = tmp - P_Phi;
198 std::cout << "diff1 "<<norm2(tmp) <<std::endl;
199
200
201 DenOp.RInv(P_Phi,tmp);
202 tmp = tmp - RinvP_Phi;
203 std::cout << "diff2 "<<norm2(tmp) <<std::endl;
204
205 DenOp.RDagInv(RinvP_Phi,tmp);
206 tmp = tmp - RinvDagRinvP_Phi;
207 std::cout << "diff3 "<<norm2(tmp) <<std::endl;
208
209 DenOp.RDag(RinvDagRinvP_Phi,tmp);
210 tmp = tmp - PdagRinvDagRinvP_Phi;
211 std::cout << "diff4 "<<norm2(tmp) <<std::endl;
212 */
213
214 dSdU=Zero();
215
216 X = DobiDdbPhi;
217 Y = DobidDddDoidRinvDagRinvP_Phi;
218 NumOp.DirichletFermOpD.MDeriv(force,Y,X,DaggerNo); dSdU=dSdU+force;
219 NumOp.DirichletFermOpD.MDeriv(force,X,Y,DaggerYes); dSdU=dSdU+force;
220
221 X = DoiDdDobiDdbPhi;
222 Y = DoidRinvDagRinvP_Phi;
223 NumOp.DirichletFermOpD.MDeriv(force,Y,X,DaggerNo); dSdU=dSdU+force;
224 NumOp.DirichletFermOpD.MDeriv(force,X,Y,DaggerYes); dSdU=dSdU+force;
225
226 X = DiDdbP_Phi;
227 Y = DidRinvP_Phi;
228 DenOp.PeriodicFermOpD.MDeriv(force,Y,X,DaggerNo); dSdU=dSdU+force;
229 DenOp.PeriodicFermOpD.MDeriv(force,X,Y,DaggerYes); dSdU=dSdU+force;
230
231 dSdU *= -1.0;
232
233 };
234};
235
237
void DumpSliceNorm(std::string s, const Lattice< vobj > &f, int mu=-1)
Definition Lattice_crc.h:32
RealD norm2(const Lattice< vobj > &arg)
void gaussian(GridParallelRNG &rng, Lattice< vobj > &l)
#define NAMESPACE_BEGIN(A)
Definition Namespace.h:35
#define NAMESPACE_END(A)
Definition Namespace.h:36
static constexpr int DaggerYes
Definition QCD.h:70
static constexpr int DaggerNo
Definition QCD.h:69
double RealD
Definition Simd.h:61
static INTERNAL_PRECISION U
Definition Zolotarev.cc:230
Base class for all actions.
Definition ActionBase.h:64
virtual void refresh(const GaugeField &U, GridSerialRNG &sRNG, GridParallelRNG &pRNG)
Refresh pseudofermion fields.
DomainDecomposedBoundaryTwoFlavourRatioPseudoFermion(SchurFactoredFermionOperator< ImplD, ImplF > &_NumOp, SchurFactoredFermionOperator< ImplD, ImplF > &_DenOp, RealD _DerivativeTol, RealD _ActionTol, RealD _InnerTol=1.0e-6)
virtual RealD S(const GaugeField &U)
Evaluate this action with the given gauge field.
Definition Simd.h:194