Grid 0.7.0
AdefGeneric.h
Go to the documentation of this file.
1 /*************************************************************************************
2
3 Grid physics library, www.github.com/paboyle/Grid
4
5 Source file: ./lib/algorithms/iterative/AdefGeneric.h
6
7 Copyright (C) 2015
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#ifndef GRID_ALGORITHMS_ITERATIVE_GENERIC_PCG
29#define GRID_ALGORITHMS_ITERATIVE_GENERIC_PCG
30
31 /*
32 * Compared to Tang-2009: P=Pleft. P^T = PRight Q=MssInv.
33 * Script A = SolverMatrix
34 * Script P = Preconditioner
35 *
36 * Implement ADEF-2
37 *
38 * Vstart = P^Tx + Qb
39 * M1 = P^TM + Q
40 * M2=M3=1
41 */
43
44
45template<class Field>
46class TwoLevelCG : public LinearFunction<Field>
47{
48 public:
52
53 // Fine operator, Smoother, CoarseSolver
56
57 // more most opertor functions
59 Integer maxit,
61 LinearFunction<Field> &Smoother,
62 GridBase *fine) :
63 Tolerance(tol),
64 MaxIterations(maxit),
65 _FineLinop(FineLinop),
66 _Smoother(Smoother)
67 {
68 grid = fine;
69 };
70
71 virtual void operator() (const Field &src, Field &x)
72 {
73 std::cout << GridLogMessage<<"HDCG: fPcg starting single RHS"<<std::endl;
74 RealD f;
75 RealD rtzp,rtz,a,d,b;
76 RealD rptzp;
77
79 // Set up history vectors
81 int mmax = 5;
82 std::cout << GridLogMessage<<"HDCG: fPcg allocating"<<std::endl;
83 std::vector<Field> p(mmax,grid);
84 std::vector<Field> mmp(mmax,grid);
85 std::vector<RealD> pAp(mmax);
86 Field z(grid);
87 Field tmp(grid);
88 Field mp (grid);
89 Field r (grid);
90 Field mu (grid);
91
92 std::cout << GridLogMessage<<"HDCG: fPcg allocated"<<std::endl;
93 //Initial residual computation & set up
94 RealD guess = norm2(x);
95 std::cout << GridLogMessage<<"HDCG: fPcg guess nrm "<<guess<<std::endl;
96 RealD src_nrm = norm2(src);
97 std::cout << GridLogMessage<<"HDCG: fPcg src nrm "<<src_nrm<<std::endl;
98
99 if ( src_nrm == 0.0 ) {
100 std::cout << GridLogMessage<<"HDCG: fPcg given trivial source norm "<<src_nrm<<std::endl;
101 x=Zero();
102 }
103 RealD tn;
104
105 GridStopWatch HDCGTimer;
106 HDCGTimer.Start();
108 // x0 = Vstart -- possibly modify guess
110 Vstart(x,src);
111
112 // r0 = b -A x0
113 _FineLinop.HermOp(x,mmp[0]);
114 axpy (r, -1.0,mmp[0], src); // Recomputes r=src-Ax0
115 {
116 double n1 = norm2(x);
117 double n2 = norm2(mmp[0]);
118 double n3 = norm2(r);
119 std::cout<<GridLogMessage<<"x,vstart,r = "<<n1<<" "<<n2<<" "<<n3<<std::endl;
120 }
121
123 // Compute z = M1 x
125 PcgM1(r,z);
126 rtzp =real(innerProduct(r,z));
127
129 // Solve for Mss mu = P A z and set p = z-mu
130 // Def2 p = 1 - Q Az = Pright z
131 // Other algos M2 is trivial
133 PcgM2(z,p[0]);
134
135 RealD ssq = norm2(src);
136 RealD rsq = ssq*Tolerance*Tolerance;
137
138 std::cout << GridLogMessage<<"HDCG: k=0 residual "<<rtzp<<" rsq "<<rsq<<"\n";
139
140 Field pp(grid);
141
142 for (int k=0;k<=MaxIterations;k++){
143
144 int peri_k = k % mmax;
145 int peri_kp = (k+1) % mmax;
146
147 rtz=rtzp;
148 d= PcgM3(p[peri_k],mmp[peri_k]);
149 a = rtz/d;
150
151 // Memorise this
152 pAp[peri_k] = d;
153
154 axpy(x,a,p[peri_k],x);
155 RealD rn = axpy_norm(r,-a,mmp[peri_k],r);
156
157 // Compute z = M x
158 PcgM1(r,z);
159
160 {
161 RealD n1,n2;
162 n1=norm2(r);
163 n2=norm2(z);
164 std::cout << GridLogMessage<<"HDCG::fPcg iteration "<<k<<" : vector r,z "<<n1<<" "<<n2<<"\n";
165 }
166 rtzp =real(innerProduct(r,z));
167 std::cout << GridLogMessage<<"HDCG::fPcg iteration "<<k<<" : inner rtzp "<<rtzp<<"\n";
168
169 // PcgM2(z,p[0]);
170 PcgM2(z,mu); // ADEF-2 this is identity. Axpy possible to eliminate
171
172 p[peri_kp]=mu;
173
174 // Standard search direction p -> z + b p
175 b = (rtzp)/rtz;
176
177 int northog;
178 // k=zero <=> peri_kp=1; northog = 1
179 // k=1 <=> peri_kp=2; northog = 2
180 // ... ... ...
181 // k=mmax-2<=> peri_kp=mmax-1; northog = mmax-1
182 // k=mmax-1<=> peri_kp=0; northog = 1
183
184 // northog = (peri_kp==0)?1:peri_kp; // This is the fCG(mmax) algorithm
185 northog = (k>mmax-1)?(mmax-1):k; // This is the fCG-Tr(mmax-1) algorithm
186
187 std::cout<<GridLogMessage<<"HDCG::fPcg iteration "<<k<<" : orthogonalising to last "<<northog<<" vectors\n";
188 for(int back=0; back < northog; back++){
189 int peri_back = (k-back)%mmax;
190 RealD pbApk= real(innerProduct(mmp[peri_back],p[peri_kp]));
191 RealD beta = -pbApk/pAp[peri_back];
192 axpy(p[peri_kp],beta,p[peri_back],p[peri_kp]);
193 }
194
195 RealD rrn=sqrt(rn/ssq);
196 RealD rtn=sqrt(rtz/ssq);
197 RealD rtnp=sqrt(rtzp/ssq);
198
199 std::cout<<GridLogMessage<<"HDCG: fPcg k= "<<k<<" residual = "<<rrn<<"\n";
200
201 // Stopping condition
202 if ( rn <= rsq ) {
203
204 HDCGTimer.Stop();
205 std::cout<<GridLogMessage<<"HDCG: fPcg converged in "<<k<<" iterations and "<<HDCGTimer.Elapsed()<<std::endl;;
206
207 _FineLinop.HermOp(x,mmp[0]);
208 axpy(tmp,-1.0,src,mmp[0]);
209
210 RealD mmpnorm = sqrt(norm2(mmp[0]));
211 RealD xnorm = sqrt(norm2(x));
212 RealD srcnorm = sqrt(norm2(src));
213 RealD tmpnorm = sqrt(norm2(tmp));
214 RealD true_residual = tmpnorm/srcnorm;
215 std::cout<<GridLogMessage
216 <<"HDCG: true residual is "<<true_residual
217 <<" solution "<<xnorm
218 <<" source "<<srcnorm
219 <<" mmp "<<mmpnorm
220 <<std::endl;
221
222 return;
223 }
224
225 }
226 HDCGTimer.Stop();
227 std::cout<<GridLogMessage<<"HDCG: not converged "<<HDCGTimer.Elapsed()<<std::endl;
228 RealD xnorm = sqrt(norm2(x));
229 RealD srcnorm = sqrt(norm2(src));
230 std::cout<<GridLogMessage<<"HDCG: non-converged solution "<<xnorm<<" source "<<srcnorm<<std::endl;
231 }
232
233
234
235 virtual void operator() (std::vector<Field> &src, std::vector<Field> &x)
236 {
237 std::cout << GridLogMessage<<"HDCG: mrhs fPcg starting"<<std::endl;
238 src[0].Grid()->Barrier();
239 int nrhs = src.size();
240 std::vector<RealD> f(nrhs);
241 std::vector<RealD> rtzp(nrhs);
242 std::vector<RealD> rtz(nrhs);
243 std::vector<RealD> a(nrhs);
244 std::vector<RealD> d(nrhs);
245 std::vector<RealD> b(nrhs);
246 std::vector<RealD> rptzp(nrhs);
248 // Set up history vectors
250 int mmax = 3;
251 std::cout << GridLogMessage<<"HDCG: fPcg allocating"<<std::endl;
252 src[0].Grid()->Barrier();
253 std::vector<std::vector<Field> > p(nrhs); for(int r=0;r<nrhs;r++) p[r].resize(mmax,grid);
254 std::cout << GridLogMessage<<"HDCG: fPcg allocated p"<<std::endl;
255 src[0].Grid()->Barrier();
256 std::vector<std::vector<Field> > mmp(nrhs); for(int r=0;r<nrhs;r++) mmp[r].resize(mmax,grid);
257 std::cout << GridLogMessage<<"HDCG: fPcg allocated mmp"<<std::endl;
258 src[0].Grid()->Barrier();
259 std::vector<std::vector<RealD> > pAp(nrhs); for(int r=0;r<nrhs;r++) pAp[r].resize(mmax);
260 std::cout << GridLogMessage<<"HDCG: fPcg allocated pAp"<<std::endl;
261 src[0].Grid()->Barrier();
262 std::vector<Field> z(nrhs,grid);
263 std::vector<Field> mp (nrhs,grid);
264 std::vector<Field> r (nrhs,grid);
265 std::vector<Field> mu (nrhs,grid);
266 std::cout << GridLogMessage<<"HDCG: fPcg allocated z,mp,r,mu"<<std::endl;
267 src[0].Grid()->Barrier();
268
269 //Initial residual computation & set up
270 std::vector<RealD> src_nrm(nrhs);
271 for(int rhs=0;rhs<nrhs;rhs++) {
272 src_nrm[rhs]=norm2(src[rhs]);
273 assert(src_nrm[rhs]!=0.0);
274 }
275 std::vector<RealD> tn(nrhs);
276
277 GridStopWatch HDCGTimer;
278 HDCGTimer.Start();
280 // x0 = Vstart -- possibly modify guess
282 Vstart(x,src);
283
284 for(int rhs=0;rhs<nrhs;rhs++){
285 // r0 = b -A x0
286 _FineLinop.HermOp(x[rhs],mmp[rhs][0]);
287 axpy (r[rhs], -1.0,mmp[rhs][0], src[rhs]); // Recomputes r=src-Ax0
288 }
289
291 // Compute z = M1 x
293 // This needs a multiRHS version for acceleration
294 PcgM1(r,z);
295
296 std::vector<RealD> ssq(nrhs);
297 std::vector<RealD> rsq(nrhs);
298 std::vector<Field> pp(nrhs,grid);
299
300 for(int rhs=0;rhs<nrhs;rhs++){
301 rtzp[rhs] =real(innerProduct(r[rhs],z[rhs]));
302 p[rhs][0]=z[rhs];
303 ssq[rhs]=norm2(src[rhs]);
304 rsq[rhs]= ssq[rhs]*Tolerance*Tolerance;
305 std::cout << GridLogMessage<<"mrhs HDCG: "<<rhs<<" k=0 residual "<<rtzp[rhs]<<" rsq "<<rsq[rhs]<<"\n";
306 }
307
308 std::vector<RealD> rn(nrhs);
309 for (int k=0;k<=MaxIterations;k++){
310
311 int peri_k = k % mmax;
312 int peri_kp = (k+1) % mmax;
313
314 for(int rhs=0;rhs<nrhs;rhs++){
315 rtz[rhs]=rtzp[rhs];
316 d[rhs]= PcgM3(p[rhs][peri_k],mmp[rhs][peri_k]);
317 a[rhs] = rtz[rhs]/d[rhs];
318
319 // Memorise this
320 pAp[rhs][peri_k] = d[rhs];
321
322 axpy(x[rhs],a[rhs],p[rhs][peri_k],x[rhs]);
323 rn[rhs] = axpy_norm(r[rhs],-a[rhs],mmp[rhs][peri_k],r[rhs]);
324 }
325
326 // Compute z = M x (for *all* RHS)
327 PcgM1(r,z);
328 std::cout << GridLogMessage<<"HDCG::fPcg M1 complete"<<std::endl;
329 grid->Barrier();
330
331 RealD max_rn=0.0;
332 for(int rhs=0;rhs<nrhs;rhs++){
333
334 rtzp[rhs] =real(innerProduct(r[rhs],z[rhs]));
335
336 std::cout << GridLogMessage<<"HDCG::fPcg rhs"<<rhs<<" iteration "<<k<<" : inner rtzp "<<rtzp[rhs]<<"\n";
337
338 mu[rhs]=z[rhs];
339
340 p[rhs][peri_kp]=mu[rhs];
341
342 // Standard search direction p == z + b p
343 b[rhs] = (rtzp[rhs])/rtz[rhs];
344
345 int northog = (k>mmax-1)?(mmax-1):k; // This is the fCG-Tr(mmax-1) algorithm
346 std::cout<<GridLogMessage<<"HDCG::fPcg iteration "<<k<<" : orthogonalising to last "<<northog<<" vectors\n";
347 for(int back=0; back < northog; back++){
348 int peri_back = (k-back)%mmax;
349 RealD pbApk= real(innerProduct(mmp[rhs][peri_back],p[rhs][peri_kp]));
350 RealD beta = -pbApk/pAp[rhs][peri_back];
351 axpy(p[rhs][peri_kp],beta,p[rhs][peri_back],p[rhs][peri_kp]);
352 }
353
354 RealD rrn=sqrt(rn[rhs]/ssq[rhs]);
355 RealD rtn=sqrt(rtz[rhs]/ssq[rhs]);
356 RealD rtnp=sqrt(rtzp[rhs]/ssq[rhs]);
357
358 std::cout<<GridLogMessage<<"HDCG: rhs "<<rhs<<"fPcg k= "<<k<<" residual = "<<rrn<<"\n";
359 if ( rrn > max_rn ) max_rn = rrn;
360 }
361
362 // Stopping condition based on worst case
363 if ( max_rn <= Tolerance ) {
364
365 HDCGTimer.Stop();
366 std::cout<<GridLogMessage<<"HDCG: mrhs fPcg converged in "<<k<<" iterations and "<<HDCGTimer.Elapsed()<<std::endl;;
367
368 for(int rhs=0;rhs<nrhs;rhs++){
369 _FineLinop.HermOp(x[rhs],mmp[rhs][0]);
370 Field tmp(grid);
371 axpy(tmp,-1.0,src[rhs],mmp[rhs][0]);
372
373 RealD mmpnorm = sqrt(norm2(mmp[rhs][0]));
374 RealD xnorm = sqrt(norm2(x[rhs]));
375 RealD srcnorm = sqrt(norm2(src[rhs]));
376 RealD tmpnorm = sqrt(norm2(tmp));
377 RealD true_residual = tmpnorm/srcnorm;
378 std::cout<<GridLogMessage
379 <<"HDCG: true residual ["<<rhs<<"] is "<<true_residual
380 <<" solution "<<xnorm
381 <<" source "<<srcnorm
382 <<" mmp "<<mmpnorm
383 <<std::endl;
384 }
385 return;
386 }
387
388 }
389 HDCGTimer.Stop();
390 std::cout<<GridLogMessage<<"HDCG: not converged "<<HDCGTimer.Elapsed()<<std::endl;
391 for(int rhs=0;rhs<nrhs;rhs++){
392 RealD xnorm = sqrt(norm2(x[rhs]));
393 RealD srcnorm = sqrt(norm2(src[rhs]));
394 std::cout<<GridLogMessage<<"HDCG: non-converged solution "<<xnorm<<" source "<<srcnorm<<std::endl;
395 }
396 }
397
398
399 public:
400
401 virtual void PcgM1(std::vector<Field> & in,std::vector<Field> & out)
402 {
403 std::cout << "PcgM1 default (cheat) mrhs version"<<std::endl;
404 for(int rhs=0;rhs<in.size();rhs++){
405 this->PcgM1(in[rhs],out[rhs]);
406 }
407 }
408 virtual void PcgM1(Field & in, Field & out) =0;
409 virtual void Vstart(std::vector<Field> & x,std::vector<Field> & src)
410 {
411 std::cout << "Vstart default (cheat) mrhs version"<<std::endl;
412 for(int rhs=0;rhs<x.size();rhs++){
413 this->Vstart(x[rhs],src[rhs]);
414 }
415 }
416 virtual void Vstart(Field & x,const Field & src)=0;
417
418 virtual void PcgM2(const Field & in, Field & out) {
419 out=in;
420 }
421
422 virtual RealD PcgM3(const Field & p, Field & mmp){
423 RealD dd;
424 _FineLinop.HermOp(p,mmp);
425 ComplexD dot = innerProduct(p,mmp);
426 dd=real(dot);
427 return dd;
428 }
429
431 // Only Def1 has non-trivial Vout.
433
434};
435
436template<class Field, class CoarseField, class Aggregation>
437class TwoLevelADEF2 : public TwoLevelCG<Field>
438{
439 public:
441 // Need something that knows how to get from Coarse to fine and back again
442 // void ProjectToSubspace(CoarseVector &CoarseVec,const FineField &FineVec){
443 // void PromoteFromSubspace(const CoarseVector &CoarseVec,FineField &FineVec){
450
451 // more most opertor functions
453 Integer maxit,
454 LinearOperatorBase<Field> &FineLinop,
455 LinearFunction<Field> &Smoother,
456 LinearFunction<CoarseField> &CoarseSolver,
457 LinearFunction<CoarseField> &CoarseSolverPrecise,
458 Aggregation &Aggregates
459 ) :
460 TwoLevelCG<Field>(tol,maxit,FineLinop,Smoother,Aggregates.FineGrid),
461 _CoarseSolver(CoarseSolver),
462 _CoarseSolverPrecise(CoarseSolverPrecise),
463 _Aggregates(Aggregates)
464 {
465 coarsegrid = Aggregates.CoarseGrid;
466 };
467
468 virtual void PcgM1(Field & in, Field & out)
469 {
470 GRID_TRACE("MultiGridPreconditioner ");
471 // [PTM+Q] in = [1 - Q A] M in + Q in = Min + Q [ in -A Min]
472
473 Field tmp(this->grid);
474 Field Min(this->grid);
475 CoarseField PleftProj(this->coarsegrid);
476 CoarseField PleftMss_proj(this->coarsegrid);
477
478 GridStopWatch SmootherTimer;
479 GridStopWatch MatrixTimer;
480 SmootherTimer.Start();
481 this->_Smoother(in,Min);
482 SmootherTimer.Stop();
483
484 MatrixTimer.Start();
485 this->_FineLinop.HermOp(Min,out);
486 MatrixTimer.Stop();
487 axpy(tmp,-1.0,out,in); // tmp = in - A Min
488
489 GridStopWatch ProjTimer;
490 GridStopWatch CoarseTimer;
491 GridStopWatch PromTimer;
492 ProjTimer.Start();
493 this->_Aggregates.ProjectToSubspace(PleftProj,tmp);
494 ProjTimer.Stop();
495 CoarseTimer.Start();
496 this->_CoarseSolver(PleftProj,PleftMss_proj); // Ass^{-1} [in - A Min]_s
497 CoarseTimer.Stop();
498 PromTimer.Start();
499 this->_Aggregates.PromoteFromSubspace(PleftMss_proj,tmp);// tmp = Q[in - A Min]
500 PromTimer.Stop();
501 std::cout << GridLogPerformance << "PcgM1 breakdown "<<std::endl;
502 std::cout << GridLogPerformance << "\tSmoother " << SmootherTimer.Elapsed() <<std::endl;
503 std::cout << GridLogPerformance << "\tMatrix " << MatrixTimer.Elapsed() <<std::endl;
504 std::cout << GridLogPerformance << "\tProj " << ProjTimer.Elapsed() <<std::endl;
505 std::cout << GridLogPerformance << "\tCoarse " << CoarseTimer.Elapsed() <<std::endl;
506 std::cout << GridLogPerformance << "\tProm " << PromTimer.Elapsed() <<std::endl;
507
508 axpy(out,1.0,Min,tmp); // Min+tmp
509 }
510
511 virtual void Vstart(Field & x,const Field & src)
512 {
513 std::cout << GridLogMessage<<"HDCG: fPcg Vstart "<<std::endl;
515 // Choose x_0 such that
516 // x_0 = guess + (A_ss^inv) r_s = guess + Ass_inv [src -Aguess]
517 // = [1 - Ass_inv A] Guess + Assinv src
518 // = P^T guess + Assinv src
519 // = Vstart [Tang notation]
520 // This gives:
521 // W^T (src - A x_0) = src_s - A guess_s - r_s
522 // = src_s - (A guess)_s - src_s + (A guess)_s
523 // = 0
525 Field r(this->grid);
526 Field mmp(this->grid);
527 CoarseField PleftProj(this->coarsegrid);
528 CoarseField PleftMss_proj(this->coarsegrid);
529
530 std::cout << GridLogMessage<<"HDCG: fPcg Vstart projecting "<<std::endl;
531 this->_Aggregates.ProjectToSubspace(PleftProj,src);
532 std::cout << GridLogMessage<<"HDCG: fPcg Vstart coarse solve "<<std::endl;
533 this->_CoarseSolverPrecise(PleftProj,PleftMss_proj); // Ass^{-1} r_s
534 std::cout << GridLogMessage<<"HDCG: fPcg Vstart promote "<<std::endl;
535 this->_Aggregates.PromoteFromSubspace(PleftMss_proj,x);
536
537 }
538
539};
540
541
542template<class Field>
543class TwoLevelADEF1defl : public TwoLevelCG<Field>
544{
545public:
546 const std::vector<Field> &evec;
547 const std::vector<RealD> &eval;
548
550 Integer maxit,
551 LinearOperatorBase<Field> &FineLinop,
552 LinearFunction<Field> &Smoother,
553 std::vector<Field> &_evec,
554 std::vector<RealD> &_eval) :
555 TwoLevelCG<Field>(tol,maxit,FineLinop,Smoother,_evec[0].Grid()),
556 evec(_evec),
557 eval(_eval)
558 {};
559
560 // Can just inherit existing M2
561 // Can just inherit existing M3
562
563 // Simple vstart - do nothing
564 virtual void Vstart(Field & x,const Field & src){
565 x=src; // Could apply Q
566 };
567
568 // Override PcgM1
569 virtual void PcgM1(Field & in, Field & out)
570 {
571 GRID_TRACE("EvecPreconditioner ");
572 int N=evec.size();
573 Field Pin(this->grid);
574 Field Qin(this->grid);
575
576 //MP + Q = M(1-AQ) + Q = M
577 // // If we are eigenvector deflating in coarse space
578 // // Q = Sum_i |phi_i> 1/lambda_i <phi_i|
579 // // A Q = Sum_i |phi_i> <phi_i|
580 // // M(1-AQ) = M(1-proj) + Q
581 Qin.Checkerboard()=in.Checkerboard();
582 Qin = Zero();
583 Pin = in;
584 for (int i=0;i<N;i++) {
585 const Field& tmp = evec[i];
586 auto ip = TensorRemove(innerProduct(tmp,in));
587 axpy(Qin, ip / eval[i],tmp,Qin);
588 axpy(Pin, -ip ,tmp,Pin);
589 }
590
591 this->_Smoother(Pin,out);
592
593 out = out + Qin;
594 }
595};
596
598
599#endif
accelerator_inline Grid_simd< S, V > sqrt(const Grid_simd< S, V > &r)
RealD axpy_norm(Lattice< vobj > &ret, sobj a, const Lattice< vobj > &x, const Lattice< vobj > &y)
void axpy(Lattice< vobj > &ret, sobj a, const Lattice< vobj > &x, const Lattice< vobj > &y)
Lattice< vobj > real(const Lattice< vobj > &lhs)
ComplexD innerProduct(const Lattice< vobj > &left, const Lattice< vobj > &right)
RealD norm2(const Lattice< vobj > &arg)
GridLogger GridLogPerformance(1, "Performance", GridLogColours, "GREEN")
GridLogger GridLogMessage(1, "Message", GridLogColours, "NORMAL")
#define NAMESPACE_BEGIN(A)
Definition Namespace.h:35
#define NAMESPACE_END(A)
Definition Namespace.h:36
uint32_t Integer
Definition Simd.h:58
std::complex< RealD > ComplexD
Definition Simd.h:79
double RealD
Definition Simd.h:61
accelerator_inline std::enable_if<!isGridTensor< T >::value, T >::type TensorRemove(T arg)
#define GRID_TRACE(name)
Definition Tracing.h:68
void PromoteFromSubspace(const CoarseVector &CoarseVec, FineField &FineVec)
Definition Aggregates.h:78
GridBase * CoarseGrid
Definition Aggregates.h:56
void ProjectToSubspace(CoarseVector &CoarseVec, const FineField &FineVec)
Definition Aggregates.h:75
void Start(void)
Definition Timer.h:92
GridTime Elapsed(void) const
Definition Timer.h:113
void Stop(void)
Definition Timer.h:99
const std::vector< Field > & evec
virtual void Vstart(Field &x, const Field &src)
virtual void PcgM1(Field &in, Field &out)
const std::vector< RealD > & eval
TwoLevelADEF1defl(RealD tol, Integer maxit, LinearOperatorBase< Field > &FineLinop, LinearFunction< Field > &Smoother, std::vector< Field > &_evec, std::vector< RealD > &_eval)
LinearFunction< CoarseField > & _CoarseSolverPrecise
virtual void Vstart(Field &x, const Field &src)
Aggregation & _Aggregates
TwoLevelADEF2(RealD tol, Integer maxit, LinearOperatorBase< Field > &FineLinop, LinearFunction< Field > &Smoother, LinearFunction< CoarseField > &CoarseSolver, LinearFunction< CoarseField > &CoarseSolverPrecise, Aggregation &Aggregates)
GridBase * coarsegrid
LinearFunction< CoarseField > & _CoarseSolver
virtual void PcgM1(Field &in, Field &out)
virtual void Vstart(Field &x, const Field &src)=0
virtual void Vstart(std::vector< Field > &x, std::vector< Field > &src)
virtual void PcgM1(std::vector< Field > &in, std::vector< Field > &out)
LinearOperatorBase< Field > & _FineLinop
Definition AdefGeneric.h:54
TwoLevelCG(RealD tol, Integer maxit, LinearOperatorBase< Field > &FineLinop, LinearFunction< Field > &Smoother, GridBase *fine)
Definition AdefGeneric.h:58
LinearFunction< Field > & _Smoother
Definition AdefGeneric.h:55
RealD Tolerance
Definition AdefGeneric.h:49
Integer MaxIterations
Definition AdefGeneric.h:50
virtual RealD PcgM3(const Field &p, Field &mmp)
virtual void PcgM2(const Field &in, Field &out)
virtual void operator()(const Field &src, Field &x)
Definition AdefGeneric.h:71
GridBase * grid
Definition AdefGeneric.h:51
virtual void PcgM1(Field &in, Field &out)=0
Definition Simd.h:194