-
Notifications
You must be signed in to change notification settings - Fork 108
/
Copy pathmultigrid.cpp
2019 lines (1660 loc) · 86.8 KB
/
multigrid.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#include <cstring>
#include <multigrid.h>
#include <tune_quda.h>
#include <random_quda.h>
#include <vector_io.h>
#include <solver.hpp>
// for building the KD inverse op
#include <staggered_kd_build_xinv.h>
namespace quda
{
using namespace blas;
static bool debug = false;
MG::MG(MGParam ¶m, TimeProfile &profile_global) :
Solver(*param.matResidual, *param.matSmooth, *param.matSmoothSloppy, *param.matSmoothSloppy, param, profile),
param(param),
transfer(0),
resetTransfer(false),
presmoother(nullptr),
postsmoother(nullptr),
profile_global(profile_global),
profile("MG level " + std::to_string(param.level), false),
coarse(nullptr),
coarse_solver(nullptr),
param_coarse(nullptr),
param_presmooth(nullptr),
param_postsmooth(nullptr),
param_coarse_solver(nullptr),
B_coarse(),
r(nullptr),
b_tilde(nullptr),
r_coarse(nullptr),
x_coarse(nullptr),
tmp_coarse(nullptr),
tmp2_coarse(nullptr),
tmp_coarse_sloppy(nullptr),
tmp2_coarse_sloppy(nullptr),
xInvKD(nullptr),
xInvKD_sloppy(nullptr),
diracResidual(param.matResidual->Expose()),
diracSmoother(param.matSmooth->Expose()),
diracSmootherSloppy(param.matSmoothSloppy->Expose()),
diracCoarseResidual(nullptr),
diracCoarseSmoother(nullptr),
diracCoarseSmootherSloppy(nullptr),
matCoarseResidual(nullptr),
matCoarseSmoother(nullptr),
matCoarseSmootherSloppy(nullptr),
rng(nullptr)
{
sprintf(prefix, "MG level %d (%s): ", param.level, param.location == QUDA_CUDA_FIELD_LOCATION ? "GPU" : "CPU");
pushLevel(param.level);
if (param.level >= QUDA_MAX_MG_LEVEL)
errorQuda("Level=%d is greater than limit of multigrid recursion depth", param.level);
if (param.coarse_grid_solution_type == QUDA_MATPC_SOLUTION && param.smoother_solve_type != QUDA_DIRECT_PC_SOLVE)
errorQuda("Cannot use preconditioned coarse grid solution without preconditioned smoother solve");
// allocating vectors
{
// create residual vectors
ColorSpinorParam csParam(*(param.B[0]));
csParam.create = QUDA_NULL_FIELD_CREATE;
csParam.location = param.location;
csParam.setPrecision(param.mg_global.invert_param->cuda_prec_sloppy, QUDA_INVALID_PRECISION,
csParam.location == QUDA_CUDA_FIELD_LOCATION ? true : false);
if (csParam.location==QUDA_CUDA_FIELD_LOCATION) {
csParam.gammaBasis = param.level > 0 ? QUDA_DEGRAND_ROSSI_GAMMA_BASIS: QUDA_UKQCD_GAMMA_BASIS;
}
if (param.B[0]->Nspin() == 1) csParam.gammaBasis = param.B[0]->GammaBasis(); // hack for staggered to avoid unnecessary basis checks
r = new ColorSpinorField(csParam);
// if we're using preconditioning then allocate storage for the preconditioned source vector
if (param.smoother_solve_type == QUDA_DIRECT_PC_SOLVE) {
csParam.x[0] /= 2;
csParam.siteSubset = QUDA_PARITY_SITE_SUBSET;
b_tilde = new ColorSpinorField(csParam);
}
}
rng = new RNG(*param.B[0], 1234);
if (!is_coarsest_grid()) {
// allocate coarse temporary vectors
tmp_coarse = param.B[0]->CreateCoarse(param.geoBlockSize, param.spinBlockSize, param.Nvec, r->Precision(),
param.mg_global.location[param.level + 1]);
tmp2_coarse = param.B[0]->CreateCoarse(param.geoBlockSize, param.spinBlockSize, param.Nvec, r->Precision(),
param.mg_global.location[param.level + 1]);
r_coarse = param.B[0]->CreateCoarse(param.geoBlockSize, param.spinBlockSize, param.Nvec, r->Precision(),
param.mg_global.location[param.level + 1]);
x_coarse = param.B[0]->CreateCoarse(param.geoBlockSize, param.spinBlockSize, param.Nvec, r->Precision(),
param.mg_global.location[param.level + 1]);
}
if (param.transfer_type == QUDA_TRANSFER_AGGREGATE && !is_coarsest_grid()) {
constructNearNulls();
// only save if outfile is defined
if (strcmp(param.mg_global.vec_outfile[param.level], "") != 0)
saveVectors(param.B);
}
// in case of iterative setup with MG the coarse level may be already built
if (!transfer) reset();
popLevel();
}
void MG::reset(bool refresh) {
pushLevel(param.level);
logQuda(QUDA_VERBOSE, "%s level %d\n", transfer ? "Resetting" : "Creating", param.level);
destroySmoother();
destroyCoarseSolver();
// reset the Dirac operator pointers since these may have changed
diracResidual = param.matResidual->Expose();
diracSmoother = param.matSmooth->Expose();
diracSmootherSloppy = param.matSmoothSloppy->Expose();
// Only refresh if we needed to generate near-nulls, that is,
// if we aren't doing a staggered KD solve
if (param.level != 0 || param.transfer_type == QUDA_TRANSFER_AGGREGATE) {
// Refresh the null-space vectors if we need to
// FIXME: can we meaningfully refresh Chebyshev filter null vecs?
if (refresh && !is_coarsest_grid()) {
if (param.mg_global.setup_maxiter_refresh[param.level]) constructNearNulls(refresh);
}
}
// if not on the coarsest level, update next
if (!is_coarsest_grid()) {
if (transfer) {
// restoring FULL parity in Transfer changed at the end of this procedure
transfer->setSiteSubset(QUDA_FULL_SITE_SUBSET, QUDA_INVALID_PARITY);
if (resetTransfer || refresh) {
transfer->reset();
resetTransfer = false;
}
} else {
// create transfer operator
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Creating transfer operator\n");
transfer = new Transfer(param.B, param.Nvec, param.NblockOrtho, param.blockOrthoTwoPass, param.geoBlockSize,
param.spinBlockSize, param.mg_global.precision_null[param.level],
param.mg_global.transfer_type[param.level], profile);
for (int i=0; i<QUDA_MAX_MG_LEVEL; i++) param.mg_global.geo_block_size[param.level][i] = param.geoBlockSize[i];
int n_vec_coarse = param.mg_global.n_vec[param.level + 1];
B_coarse.resize(n_vec_coarse);
// only have single precision B vectors on the coarse grid
QudaPrecision B_coarse_precision = std::max(param.mg_global.precision_null[param.level+1], QUDA_SINGLE_PRECISION);
for (int i = 0; i < n_vec_coarse; i++)
B_coarse[i] = param.B[0]->CreateCoarse(param.geoBlockSize, param.spinBlockSize, param.Nvec, B_coarse_precision, param.mg_global.setup_location[param.level+1]);
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Transfer operator done\n");
}
// we no longer need the B fields for this level, can evict them to host memory
// (only if using managed memory and prefetching is enabled, otherwise no-op)
for (int i = 0; i < param.Nvec; i++) { param.B[i]->prefetch(QUDA_CPU_FIELD_LOCATION); }
createCoarseDirac();
}
// delay allocating smoother until after coarse-links have been created
createSmoother();
if (!is_coarsest_grid()) {
// If enabled, verify the coarse links and fine solvers were correctly built
if (param.mg_global.run_verify) verify();
// creating or resetting the coarse level temporaries and solvers
if (coarse) {
coarse->param.updateInvertParam(*param.mg_global.invert_param);
coarse->param.delta = 1e-20;
coarse->param.precision = param.mg_global.invert_param->cuda_prec_precondition;
coarse->param.matResidual = matCoarseResidual;
coarse->param.matSmooth = matCoarseSmoother;
coarse->param.matSmoothSloppy = matCoarseSmootherSloppy;
coarse->reset(refresh);
} else {
// create the next multigrid level
param_coarse = new MGParam(param, B_coarse, matCoarseResidual, matCoarseSmoother, matCoarseSmootherSloppy,
param.level + 1);
param_coarse->fine = this;
param_coarse->delta = 1e-20;
param_coarse->precision = param.mg_global.invert_param->cuda_prec_precondition;
coarse = new MG(*param_coarse, profile_global);
}
setOutputPrefix(prefix); // restore since we just popped back from coarse grid
createCoarseSolver();
}
// We're going back up the coarse construct stack now, prefetch the gauge fields on
// this level back to device memory.
diracResidual->prefetch(QUDA_CUDA_FIELD_LOCATION);
diracSmoother->prefetch(QUDA_CUDA_FIELD_LOCATION);
diracSmootherSloppy->prefetch(QUDA_CUDA_FIELD_LOCATION);
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Setup of level %d done\n", param.level);
popLevel();
}
void MG::resetStaggeredKD(cudaGaugeField *gauge_in, cudaGaugeField *fat_gauge_in, cudaGaugeField *long_gauge_in,
cudaGaugeField *gauge_sloppy_in, cudaGaugeField *fat_gauge_sloppy_in,
cudaGaugeField *long_gauge_sloppy_in, double mass)
{
if (param.level != 0) errorQuda("The staggered KD operator can only be updated from level 0");
if (param.transfer_type != QUDA_TRANSFER_OPTIMIZED_KD && param.transfer_type != QUDA_TRANSFER_OPTIMIZED_KD_DROP_LONG)
errorQuda("Attempting to update fine gauge fields of a \"coarse\" but non-KD operator");
// Need to be careful here: if we're preconditioning an ASQTAD op with
// a StaggeredKD op, we need to pass the StaggeredKD op the fat links
auto dirac_type = diracSmoother->getDiracType();
if ((dirac_type == QUDA_ASQTAD_DIRAC || dirac_type == QUDA_ASQTADPC_DIRAC)
&& param.transfer_type == QUDA_TRANSFER_OPTIMIZED_KD_DROP_LONG) {
// last nullptr is for the clover field
diracCoarseResidual->updateFields(fat_gauge_in, fat_gauge_in, long_gauge_in, nullptr);
diracCoarseSmoother->updateFields(fat_gauge_in, fat_gauge_in, long_gauge_in, nullptr);
diracCoarseSmootherSloppy->updateFields(fat_gauge_sloppy_in, fat_gauge_sloppy_in, long_gauge_sloppy_in, nullptr);
} else {
// last nullptr is for the clover field
diracCoarseResidual->updateFields(gauge_in, fat_gauge_in, long_gauge_in, nullptr);
diracCoarseSmoother->updateFields(gauge_in, fat_gauge_in, long_gauge_in, nullptr);
diracCoarseSmootherSloppy->updateFields(gauge_sloppy_in, fat_gauge_sloppy_in, long_gauge_sloppy_in, nullptr);
}
diracCoarseResidual->setMass(mass);
diracCoarseSmoother->setMass(mass);
diracCoarseSmootherSloppy->setMass(mass);
// to-do: think about updating Xinv
}
void MG::pushLevel(int level) const
{
postTrace();
pushVerbosity(param.mg_global.verbosity[level]);
pushOutputPrefix(prefix);
}
void MG::popLevel() const
{
popVerbosity();
popOutputPrefix();
postTrace();
}
void MG::destroySmoother()
{
pushLevel(param.level);
if (presmoother) {
delete presmoother;
presmoother = nullptr;
}
if (param_presmooth) {
delete param_presmooth;
param_presmooth = nullptr;
}
if (postsmoother) {
delete postsmoother;
postsmoother = nullptr;
}
if (param_postsmooth) {
delete param_postsmooth;
param_postsmooth = nullptr;
}
popLevel();
}
void MG::createSmoother() {
pushLevel(param.level);
// create the smoother for this level
logQuda(QUDA_VERBOSE, "Creating smoother\n");
destroySmoother();
param_presmooth = new SolverParam(param);
param_presmooth->is_preconditioner = false;
param_presmooth->return_residual = true; // pre-smoother returns the residual vector for subsequent coarsening
param_presmooth->use_init_guess = QUDA_USE_INIT_GUESS_NO;
param_presmooth->precision = param.mg_global.invert_param->cuda_prec_sloppy;
param_presmooth->precision_sloppy = (is_fine_grid()) ? param.mg_global.invert_param->cuda_prec_precondition :
param.mg_global.invert_param->cuda_prec_sloppy;
param_presmooth->precision_precondition = (is_fine_grid()) ? param.mg_global.invert_param->cuda_prec_precondition :
param.mg_global.invert_param->cuda_prec_sloppy;
param_presmooth->inv_type = param.smoother;
param_presmooth->inv_type_precondition = QUDA_INVALID_INVERTER;
param_presmooth->residual_type = (param_presmooth->inv_type == QUDA_MR_INVERTER) ? QUDA_INVALID_RESIDUAL : QUDA_L2_RELATIVE_RESIDUAL;
param_presmooth->Nsteps = param.mg_global.smoother_schwarz_cycle[param.level];
param_presmooth->maxiter = (is_coarsest_grid()) ? (param.nu_pre + param.nu_post) : param.nu_pre;
param_presmooth->Nkrylov = param_presmooth->maxiter;
param_presmooth->pipeline = param_presmooth->maxiter;
if (is_ca_solver(param_presmooth->inv_type)) {
param_presmooth->ca_basis = param.mg_global.smoother_solver_ca_basis[param.level];
param_presmooth->ca_lambda_min = param.mg_global.smoother_solver_ca_lambda_min[param.level];
param_presmooth->ca_lambda_max = param.mg_global.smoother_solver_ca_lambda_max[param.level];
}
param_presmooth->tol = param.smoother_tol;
param_presmooth->global_reduction = param.global_reduction;
param_presmooth->sloppy_converge = true; // this means we don't check the true residual before declaring convergence
param_presmooth->schwarz_type = param.mg_global.smoother_schwarz_type[param.level];
// inner solver should recompute the true residual after each cycle if using Schwarz preconditioning
param_presmooth->compute_true_res = (param_presmooth->schwarz_type != QUDA_INVALID_SCHWARZ) ? true : false;
presmoother = ((!is_coarsest_grid() || param_presmooth->schwarz_type != QUDA_INVALID_SCHWARZ)
&& param_presmooth->inv_type != QUDA_INVALID_INVERTER && param_presmooth->maxiter > 0) ?
Solver::create(*param_presmooth, *param.matSmooth, *param.matSmoothSloppy, *param.matSmoothSloppy,
*param.matSmoothSloppy, profile) :
nullptr;
if (!is_coarsest_grid()) { // Create the post smoother
param_postsmooth = new SolverParam(*param_presmooth);
param_postsmooth->return_residual = false; // post smoother does not need to return the residual vector
param_postsmooth->use_init_guess = QUDA_USE_INIT_GUESS_YES;
param_postsmooth->maxiter = param.nu_post;
param_postsmooth->Nkrylov = param_postsmooth->maxiter;
param_postsmooth->pipeline = param_postsmooth->maxiter;
// we never need to compute the true residual for a post smoother
param_postsmooth->compute_true_res = false;
postsmoother = (param_postsmooth->inv_type != QUDA_INVALID_INVERTER && param_postsmooth->maxiter > 0) ?
Solver::create(*param_postsmooth, *param.matSmooth, *param.matSmoothSloppy, *param.matSmoothSloppy,
*param.matSmoothSloppy, profile) :
nullptr;
}
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Smoother done\n");
popLevel();
}
void MG::createCoarseDirac() {
pushLevel(param.level);
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Creating coarse Dirac operator\n");
// check if we are coarsening the preconditioned system then
bool preconditioned_coarsen = (param.coarse_grid_solution_type == QUDA_MATPC_SOLUTION && param.smoother_solve_type == QUDA_DIRECT_PC_SOLVE);
QudaMatPCType matpc_type = param.mg_global.invert_param->matpc_type;
// use even-odd preconditioning for the coarse grid solver
if (diracCoarseResidual) delete diracCoarseResidual;
if (diracCoarseSmoother) delete diracCoarseSmoother;
if (diracCoarseSmootherSloppy) delete diracCoarseSmootherSloppy;
// custom setup for the staggered KD ops
if (param.level == 0
&& (param.mg_global.transfer_type[param.level] == QUDA_TRANSFER_OPTIMIZED_KD
|| param.mg_global.transfer_type[param.level] == QUDA_TRANSFER_OPTIMIZED_KD_DROP_LONG)) {
createOptimizedKdDirac();
} else {
// create coarse grid operator
DiracParam diracParam;
diracParam.transfer = transfer;
// Parameters that matter for coarse construction and application
diracParam.dirac = preconditioned_coarsen ? const_cast<Dirac *>(diracSmoother) : const_cast<Dirac *>(diracResidual);
diracParam.kappa = (param.B[0]->Nspin() == 1) ?
-1.0 :
diracParam.dirac->Kappa(); // -1 cancels automatic kappa in application of Y fields
diracParam.mass = diracParam.dirac->Mass();
diracParam.mu = diracParam.dirac->Mu();
diracParam.mu_factor = param.mg_global.mu_factor[param.level + 1] - param.mg_global.mu_factor[param.level];
// Need to figure out if we need to force bi-directional build. If any previous level (incl this one) was
// preconditioned, or a KD op, we have to force bi-directional builds.
diracParam.need_bidirectional = QUDA_BOOLEAN_FALSE;
for (int i = 0; i <= param.level; i++) {
if ((param.mg_global.coarse_grid_solution_type[i] == QUDA_MATPC_SOLUTION
&& param.mg_global.smoother_solve_type[i] == QUDA_DIRECT_PC_SOLVE)
|| (param.mg_global.transfer_type[i] == QUDA_TRANSFER_OPTIMIZED_KD)) {
diracParam.need_bidirectional = QUDA_BOOLEAN_TRUE;
}
}
diracParam.dagger = QUDA_DAG_NO;
diracParam.matpcType = matpc_type;
diracParam.type = QUDA_COARSE_DIRAC;
diracParam.tmp1 = tmp_coarse;
diracParam.tmp2 = tmp2_coarse;
diracParam.halo_precision = param.mg_global.precision_null[param.level];
diracParam.use_mma = param.use_mma;
diracParam.allow_truncation = (param.mg_global.allow_truncation == QUDA_BOOLEAN_TRUE) ? true : false;
diracCoarseResidual = new DiracCoarse(diracParam, param.setup_location == QUDA_CUDA_FIELD_LOCATION ? true : false,
param.mg_global.setup_minimize_memory == QUDA_BOOLEAN_TRUE ? true : false);
// create smoothing operators
diracParam.dirac = const_cast<Dirac *>(param.matSmooth->Expose());
diracParam.halo_precision = param.mg_global.smoother_halo_precision[param.level + 1];
if (param.mg_global.smoother_solve_type[param.level + 1] == QUDA_DIRECT_PC_SOLVE) {
diracParam.type = QUDA_COARSEPC_DIRAC;
diracParam.tmp1 = &(tmp_coarse->Even());
diracParam.tmp2 = &(tmp2_coarse->Even());
diracCoarseSmoother = new DiracCoarsePC(static_cast<DiracCoarse &>(*diracCoarseResidual), diracParam);
{
bool schwarz = param.mg_global.smoother_schwarz_type[param.level + 1] != QUDA_INVALID_SCHWARZ;
for (int i = 0; i < 4; i++) diracParam.commDim[i] = schwarz ? 0 : 1;
}
diracCoarseSmootherSloppy = new DiracCoarsePC(static_cast<DiracCoarse &>(*diracCoarseSmoother), diracParam);
} else {
diracParam.type = QUDA_COARSE_DIRAC;
diracParam.tmp1 = tmp_coarse;
diracParam.tmp2 = tmp2_coarse;
diracCoarseSmoother = new DiracCoarse(static_cast<DiracCoarse &>(*diracCoarseResidual), diracParam);
{
bool schwarz = param.mg_global.smoother_schwarz_type[param.level + 1] != QUDA_INVALID_SCHWARZ;
for (int i = 0; i < 4; i++) diracParam.commDim[i] = schwarz ? 0 : 1;
}
diracCoarseSmootherSloppy = new DiracCoarse(static_cast<DiracCoarse &>(*diracCoarseSmoother), diracParam);
}
}
if (matCoarseResidual) delete matCoarseResidual;
if (matCoarseSmoother) delete matCoarseSmoother;
if (matCoarseSmootherSloppy) delete matCoarseSmootherSloppy;
matCoarseResidual = new DiracM(*diracCoarseResidual);
matCoarseSmoother = new DiracM(*diracCoarseSmoother);
matCoarseSmootherSloppy = new DiracM(*diracCoarseSmootherSloppy);
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Coarse Dirac operator done\n");
popLevel();
}
void MG::createOptimizedKdDirac()
{
pushLevel(param.level);
auto dirac_type = diracSmoother->getDiracType();
auto smoother_solve_type = param.mg_global.smoother_solve_type[param.level + 1];
if (smoother_solve_type != QUDA_DIRECT_SOLVE) {
errorQuda("Invalid solve type %d for optimized KD operator", smoother_solve_type);
}
// Determine if we're doing a mixed precision solve for setup or not
bool mixed_precision_setup
= (param.mg_global.invert_param->cuda_prec_precondition != param.mg_global.invert_param->cuda_prec_sloppy);
// Allocate and build the KD inverse block (inverse coarse clover)
auto fine_dirac_type = diracSmoother->getDiracType();
if (fine_dirac_type != dirac_type)
errorQuda("Input dirac type %d does not match smoother type %d\n", dirac_type, fine_dirac_type);
// Determine if the dirac_type is naive staggered
bool is_naive_staggered = (dirac_type == QUDA_STAGGERED_DIRAC || dirac_type == QUDA_STAGGEREDPC_DIRAC);
bool is_improved_staggered = (dirac_type == QUDA_ASQTAD_DIRAC || dirac_type == QUDA_ASQTADPC_DIRAC);
bool is_coarse_naive_staggered = is_naive_staggered
|| (is_improved_staggered && param.mg_global.transfer_type[param.level] == QUDA_TRANSFER_OPTIMIZED_KD_DROP_LONG);
cudaGaugeField *fine_gauge = diracSmoother->getStaggeredShortLinkField();
cudaGaugeField *sloppy_gauge = mixed_precision_setup ? diracSmootherSloppy->getStaggeredShortLinkField() : fine_gauge;
xInvKD = AllocateAndBuildStaggeredKahlerDiracInverse(
*fine_gauge, diracSmoother->Mass(), param.mg_global.staggered_kd_dagger_approximation == QUDA_BOOLEAN_TRUE);
// Unique to the KD operator as a "coarse level", we can do a mixed-precision
// near null generation.
if (mixed_precision_setup) {
GaugeFieldParam xinv_param(*xInvKD);
// true is to force FLOAT2
xinv_param.setPrecision(param.mg_global.invert_param->cuda_prec_precondition, true);
xInvKD_sloppy = std::shared_ptr<GaugeField>(reinterpret_cast<GaugeField *>(new cudaGaugeField(xinv_param)));
xInvKD_sloppy->copy(*xInvKD);
ColorSpinorParam sloppy_tmp_param(*tmp_coarse);
sloppy_tmp_param.setPrecision(param.mg_global.invert_param->cuda_prec_precondition);
tmp_coarse_sloppy = new ColorSpinorField(sloppy_tmp_param);
tmp2_coarse_sloppy = new ColorSpinorField(sloppy_tmp_param);
} else {
// We can just alias fields
xInvKD_sloppy = xInvKD;
}
DiracParam diracParamKD;
diracParamKD.kappa
= -1.0; // Cancels automatic kappa in Y field application, which may be relevant if it propagates down
diracParamKD.mass = diracSmoother->Mass();
diracParamKD.mu = diracSmoother->Mu(); // doesn't matter
diracParamKD.mu_factor = 1.0; // doesn't matter
diracParamKD.dagger = QUDA_DAG_NO;
diracParamKD.matpcType = QUDA_MATPC_EVEN_EVEN; // We can use this to track left vs right block jacobi in the future
diracParamKD.gauge = const_cast<cudaGaugeField *>(fine_gauge);
diracParamKD.xInvKD = xInvKD.get(); // FIXME: pulling a raw unmanaged pointer out of a unique_ptr...
diracParamKD.dirac
= const_cast<Dirac *>(diracSmoother); // used to determine if the outer solve is preconditioned or not
diracParamKD.tmp1 = tmp_coarse;
diracParamKD.tmp2 = tmp2_coarse;
if (is_coarse_naive_staggered) {
diracParamKD.type = QUDA_STAGGEREDKD_DIRAC;
diracCoarseResidual = new DiracStaggeredKD(diracParamKD);
diracCoarseSmoother = new DiracStaggeredKD(diracParamKD);
if (mixed_precision_setup) {
diracParamKD.gauge = sloppy_gauge;
diracParamKD.xInvKD = xInvKD_sloppy.get();
diracParamKD.dirac = nullptr;
diracParamKD.tmp1 = tmp_coarse_sloppy;
diracParamKD.tmp2 = tmp2_coarse_sloppy;
}
diracCoarseSmootherSloppy = new DiracStaggeredKD(diracParamKD);
} else if (is_improved_staggered) {
diracParamKD.type = QUDA_ASQTADKD_DIRAC;
diracParamKD.fatGauge = fine_gauge;
diracParamKD.longGauge = diracSmoother->getStaggeredLongLinkField();
diracCoarseResidual = new DiracImprovedStaggeredKD(diracParamKD);
diracCoarseSmoother = new DiracImprovedStaggeredKD(diracParamKD);
if (mixed_precision_setup) {
diracParamKD.fatGauge = sloppy_gauge;
diracParamKD.longGauge = diracSmootherSloppy->getStaggeredLongLinkField();
diracParamKD.xInvKD = xInvKD_sloppy.get();
diracParamKD.dirac = nullptr;
diracParamKD.tmp1 = tmp_coarse_sloppy;
diracParamKD.tmp2 = tmp2_coarse_sloppy;
}
diracCoarseSmootherSloppy = new DiracImprovedStaggeredKD(diracParamKD);
} else {
errorQuda("Invalid dirac_type %d", dirac_type);
}
popLevel();
}
void MG::destroyCoarseSolver() {
pushLevel(param.level);
if (param.cycle_type == QUDA_MG_CYCLE_VCYCLE && param.level < param.Nlevel-2) {
// nothing to do
} else if (param.cycle_type == QUDA_MG_CYCLE_RECURSIVE || param.level == param.Nlevel-2) {
if (coarse_solver) {
auto &coarse_solver_inner = reinterpret_cast<PreconditionedSolver *>(coarse_solver)->ExposeSolver();
// int defl_size = coarse_solver_inner.evecs.size();
int defl_size = coarse_solver_inner.deflationSpaceSize();
if (defl_size > 0 && transfer && param.mg_global.preserve_deflation) {
// Deflation space exists and we are going to create a new solver. Extract deflation space.
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Extracting deflation space size %d to MG\n", defl_size);
coarse_solver_inner.extractDeflationSpace(evecs);
}
delete coarse_solver;
coarse_solver = nullptr;
}
if (param_coarse_solver) {
delete param_coarse_solver;
param_coarse_solver = nullptr;
}
} else {
errorQuda("Multigrid cycle type %d not supported", param.cycle_type);
}
popLevel();
}
void MG::createCoarseSolver() {
pushLevel(param.level);
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Creating coarse solver wrapper\n");
destroyCoarseSolver();
if (param.cycle_type == QUDA_MG_CYCLE_VCYCLE && param.level < param.Nlevel-2) {
// if coarse solver is not a bottom solver and on the second to bottom level then we can just use the coarse solver as is
coarse_solver = coarse;
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("Assigned coarse solver to coarse MG operator\n");
} else if (param.cycle_type == QUDA_MG_CYCLE_RECURSIVE || param.level == param.Nlevel-2) {
param_coarse_solver = new SolverParam(param);
param_coarse_solver->inv_type = param.mg_global.coarse_solver[param.level + 1];
param_coarse_solver->is_preconditioner = false;
param_coarse_solver->sloppy_converge = true; // this means we don't check the true residual before declaring convergence
param_coarse_solver->return_residual = false; // coarse solver does need to return residual vector
param_coarse_solver->use_init_guess = QUDA_USE_INIT_GUESS_NO;
// Coarse level deflation is triggered if the eig param structure exists
// on the coarsest level, and we are on the next to coarsest level.
if (param.mg_global.use_eig_solver[param.Nlevel - 1] && (param.level == param.Nlevel - 2)) {
param_coarse_solver->eig_param = *param.mg_global.eig_param[param.Nlevel - 1];
param_coarse_solver->deflate = QUDA_BOOLEAN_TRUE;
// Due to coherence between these levels, an initial guess
// might be beneficial.
if (param.mg_global.coarse_guess == QUDA_BOOLEAN_TRUE) {
param_coarse_solver->use_init_guess = QUDA_USE_INIT_GUESS_YES;
}
// Deflation is not supported for all solvers
if (!is_deflatable_solver(param_coarse_solver->inv_type)) {
errorQuda("Coarse grid deflation not supported with coarse solver %d", param_coarse_solver->inv_type);
}
if (strcmp(param_coarse_solver->eig_param.vec_infile, "") == 0 && // check that input file not already set
param.mg_global.vec_load[param.level + 1] == QUDA_BOOLEAN_TRUE
&& (strcmp(param.mg_global.vec_infile[param.level + 1], "") != 0)) {
std::string vec_infile(param.mg_global.vec_infile[param.level + 1]);
vec_infile += "_level_";
vec_infile += std::to_string(param.level + 1);
vec_infile += "_defl_";
vec_infile += std::to_string(param.mg_global.n_vec[param.level + 1]);
strcpy(param_coarse_solver->eig_param.vec_infile, vec_infile.c_str());
}
if (strcmp(param_coarse_solver->eig_param.vec_outfile, "") == 0 && // check that output file not already set
param.mg_global.vec_store[param.level + 1] == QUDA_BOOLEAN_TRUE
&& (strcmp(param.mg_global.vec_outfile[param.level + 1], "") != 0)) {
std::string vec_outfile(param.mg_global.vec_outfile[param.level + 1]);
vec_outfile += "_level_";
vec_outfile += std::to_string(param.level + 1);
vec_outfile += "_defl_";
vec_outfile += std::to_string(param.mg_global.n_vec[param.level + 1]);
strcpy(param_coarse_solver->eig_param.vec_outfile, vec_outfile.c_str());
}
}
param_coarse_solver->tol = param.mg_global.coarse_solver_tol[param.level+1];
param_coarse_solver->global_reduction = true;
param_coarse_solver->compute_true_res = false;
param_coarse_solver->delta = 1e-8;
param_coarse_solver->pipeline = 8;
param_coarse_solver->maxiter = param.mg_global.coarse_solver_maxiter[param.level+1];
param_coarse_solver->Nkrylov = param_coarse_solver->maxiter < param_coarse_solver->Nkrylov ?
param_coarse_solver->maxiter :
param_coarse_solver->Nkrylov;
if (is_ca_solver(param_coarse_solver->inv_type)) {
param_coarse_solver->ca_basis = param.mg_global.coarse_solver_ca_basis[param.level+1];
param_coarse_solver->ca_lambda_min = param.mg_global.coarse_solver_ca_lambda_min[param.level+1];
param_coarse_solver->ca_lambda_max = param.mg_global.coarse_solver_ca_lambda_max[param.level+1];
param_coarse_solver->Nkrylov = param.mg_global.coarse_solver_ca_basis_size[param.level+1];
} else if (param_coarse_solver->inv_type == QUDA_BICGSTABL_INVERTER) {
param_coarse_solver->Nkrylov = param.mg_global.coarse_solver_ca_basis_size[param.level + 1];
}
param_coarse_solver->inv_type_precondition = (param.level<param.Nlevel-2 || coarse->presmoother) ? QUDA_MG_INVERTER : QUDA_INVALID_INVERTER;
param_coarse_solver->preconditioner = (param.level<param.Nlevel-2 || coarse->presmoother) ? coarse : nullptr;
param_coarse_solver->mg_instance = true;
param_coarse_solver->verbosity_precondition = param.mg_global.verbosity[param.level+1];
// preconditioned solver wrapper is uniform precision
param_coarse_solver->precision = r_coarse->Precision();
param_coarse_solver->precision_sloppy = param_coarse_solver->precision;
param_coarse_solver->precision_precondition = param_coarse_solver->precision_sloppy;
if (param.mg_global.coarse_grid_solution_type[param.level + 1] == QUDA_MATPC_SOLUTION) {
Solver *solver = Solver::create(*param_coarse_solver, *matCoarseSmoother, *matCoarseSmoother,
*matCoarseSmoother, *matCoarseSmoother, profile);
sprintf(coarse_prefix, "MG level %d (%s): ", param.level + 1,
param.mg_global.location[param.level + 1] == QUDA_CUDA_FIELD_LOCATION ? "GPU" : "CPU");
coarse_solver = new PreconditionedSolver(*solver, *matCoarseSmoother->Expose(), *param_coarse_solver, profile,
coarse_prefix);
} else {
Solver *solver = Solver::create(*param_coarse_solver, *matCoarseResidual, *matCoarseResidual,
*matCoarseResidual, *matCoarseResidual, profile);
sprintf(coarse_prefix, "MG level %d (%s): ", param.level + 1,
param.mg_global.location[param.level + 1] == QUDA_CUDA_FIELD_LOCATION ? "GPU" : "CPU");
coarse_solver = new PreconditionedSolver(*solver, *matCoarseResidual->Expose(), *param_coarse_solver, profile, coarse_prefix);
}
setOutputPrefix(prefix); // restore since we just popped back from coarse grid
if (param.level == param.Nlevel - 2 && param.mg_global.use_eig_solver[param.level + 1]) {
// Test if a coarse grid deflation space needs to be transferred to the coarse solver to prevent recomputation
int defl_size = evecs.size();
auto &coarse_solver_inner = reinterpret_cast<PreconditionedSolver *>(coarse_solver)->ExposeSolver();
if (defl_size > 0 && transfer && param.mg_global.preserve_deflation) {
// We shall not recompute the deflation space, we shall transfer
// vectors stored in the parent MG instead
coarse_solver_inner.setDeflateCompute(false);
coarse_solver_inner.setRecomputeEvals(true);
logQuda(QUDA_VERBOSE,"Transferring deflation space size %d to coarse solver\n", defl_size);
// Create space in coarse solver to hold deflation space, destroy space in MG.
coarse_solver_inner.injectDeflationSpace(evecs);
}
// Run a dummy solve so that the deflation space is constructed and computed if needed during the MG setup,
// or the eigenvalues are recomputed during transfer.
spinorNoise(*r_coarse, *coarse->rng, QUDA_NOISE_UNIFORM);
param_coarse_solver->maxiter = 1; // do a single iteration on the dummy solve
(*coarse_solver)(*x_coarse, *r_coarse);
setOutputPrefix(prefix); // restore since we just popped back from coarse grid
param_coarse_solver->maxiter = param.mg_global.coarse_solver_maxiter[param.level + 1];
}
logQuda(QUDA_VERBOSE, "Assigned coarse solver to preconditioned GCR solver\n");
} else {
errorQuda("Multigrid cycle type %d not supported", param.cycle_type);
}
logQuda(QUDA_VERBOSE, "Coarse solver wrapper done\n");
popLevel();
}
MG::~MG()
{
pushLevel(param.level);
if (param.level < param.Nlevel - 1) {
if (coarse) delete coarse;
if (param.level == param.Nlevel-1 || param.cycle_type == QUDA_MG_CYCLE_RECURSIVE) {
if (coarse_solver) delete coarse_solver;
if (param_coarse_solver) delete param_coarse_solver;
}
int n_vec_coarse = param.mg_global.n_vec[param.level + 1];
for (int i = 0; i < n_vec_coarse; i++)
if (B_coarse[i]) delete B_coarse[i];
if (transfer) delete transfer;
if (matCoarseSmootherSloppy) delete matCoarseSmootherSloppy;
if (diracCoarseSmootherSloppy) delete diracCoarseSmootherSloppy;
if (matCoarseSmoother) delete matCoarseSmoother;
if (diracCoarseSmoother) delete diracCoarseSmoother;
if (matCoarseResidual) delete matCoarseResidual;
if (diracCoarseResidual) delete diracCoarseResidual;
if (postsmoother) delete postsmoother;
if (param_postsmooth) delete param_postsmooth;
}
if (rng) {
delete rng;
}
if (presmoother) delete presmoother;
if (param_presmooth) delete param_presmooth;
if (b_tilde && param.smoother_solve_type == QUDA_DIRECT_PC_SOLVE) delete b_tilde;
if (r) delete r;
if (r_coarse) delete r_coarse;
if (x_coarse) delete x_coarse;
if (tmp_coarse) delete tmp_coarse;
if (tmp2_coarse) delete tmp2_coarse;
if (tmp_coarse_sloppy) delete tmp_coarse_sloppy;
if (tmp2_coarse_sloppy) delete tmp2_coarse_sloppy;
if (param_coarse) delete param_coarse;
if (getVerbosity() >= QUDA_VERBOSE) profile.Print();
popLevel();
}
void MG::orthonormalize_vectors(const std::vector<ColorSpinorField*>& vecs, int count) const {
int num_vec = (count == -1) ? static_cast<int>(vecs.size()) : count;
if (num_vec > (int)vecs.size())
errorQuda("Trying to orthonormalize %d vectors which is larger than std::vector size %ld", num_vec, vecs.size());
for (int i = 0; i < num_vec; i++) {
for (int j = 0; j < i ; j++) {
Complex alpha = cDotProduct(*vecs[j], *vecs[i]);// <j,i>
caxpy(-alpha, *vecs[j], *vecs[i]); // i-<j,i>j
}
double nrm2 = norm2(*vecs[i]);
if (nrm2 > 1e-16) ax(1.0 / sqrt(nrm2), *vecs[i]);// i/<i,i>
else errorQuda("Cannot normalize %u vector", i);
}
}
// FIXME need to make this more robust (implement Solver::flops() for all solvers)
double MG::flops() const {
double flops = 0;
if (param_coarse_solver) {
flops += param_coarse_solver->gflops * 1e9;
param_coarse_solver->gflops = 0;
} else if (param.level < param.Nlevel-1) {
flops += coarse->flops();
}
if (param_presmooth) {
flops += param_presmooth->gflops * 1e9;
param_presmooth->gflops = 0;
}
if (param_postsmooth) {
flops += param_postsmooth->gflops * 1e9;
param_postsmooth->gflops = 0;
}
if (transfer) {
flops += transfer->flops();
}
return flops;
}
/**
Verification that the constructed multigrid operator is valid
*/
void MG::verify(bool recursively)
{
pushLevel(param.level);
// temporary fields used for verification
ColorSpinorParam csParam(*r);
csParam.create = QUDA_NULL_FIELD_CREATE;
ColorSpinorField tmp1(csParam);
ColorSpinorField tmp2(csParam);
double deviation;
QudaPrecision prec = (param.mg_global.precision_null[param.level] < csParam.Precision()) ?
param.mg_global.precision_null[param.level] :
csParam.Precision();
// may want to revisit this---these were relaxed for cases where ghost_precision < precision
// these were set while hacking in tests of quarter precision ghosts
// moreover, we can improve the precision of block ortho with a tighter max than 1.0
double tol;
switch (prec) {
case QUDA_QUARTER_PRECISION: tol = 5e-2; break;
case QUDA_HALF_PRECISION: tol = 5e-2; break;
case QUDA_SINGLE_PRECISION: tol = 1e-3; break;
default: tol = 1e-8;
}
// No need to check (projector) v_k for staggered case
if (param.transfer_type == QUDA_TRANSFER_AGGREGATE) {
logQuda(QUDA_SUMMARIZE, "Checking 0 = (1 - P P^\\dagger) v_k for %d vectors\n", param.Nvec);
for (int i = 0; i < param.Nvec; i++) {
// as well as copying to the correct location this also changes basis if necessary
tmp1 = *param.B[i];
transfer->R(*r_coarse, tmp1);
transfer->P(tmp2, *r_coarse);
deviation = sqrt(xmyNorm(tmp1, tmp2) / norm2(tmp1));
if (getVerbosity() >= QUDA_VERBOSE)
printfQuda(
"Vector %d: norms v_k = %e P^\\dagger v_k = %e (1 - P P^\\dagger) v_k = %e, L2 relative deviation = %e\n",
i, norm2(tmp1), norm2(*r_coarse), norm2(tmp2), deviation);
if (deviation > tol) errorQuda("L2 relative deviation for k=%d failed, %e > %e", i, deviation, tol);
}
// the oblique check
if (param.mg_global.run_oblique_proj_check) {
sprintf(prefix, "MG level %d (%s): Null vector Oblique Projections : ", param.level + 1,
param.location == QUDA_CUDA_FIELD_LOCATION ? "GPU" : "CPU");
setOutputPrefix(prefix);
// Oblique projections
if (getVerbosity() >= QUDA_SUMMARIZE)
printfQuda("Checking 1 > || (1 - DP(P^dagDP)P^dag) v_k || / || v_k || for %d vectors\n", param.Nvec);
for (int i = 0; i < param.Nvec; i++) {
transfer->R(*r_coarse, *(param.B[i]));
(*coarse_solver)(*x_coarse, *r_coarse); // this needs to be an exact solve to pass
setOutputPrefix(prefix); // restore prefix after return from coarse grid
transfer->P(tmp2, *x_coarse);
(*param.matResidual)(tmp1, tmp2);
tmp2 = *(param.B[i]);
if (getVerbosity() >= QUDA_SUMMARIZE) {
printfQuda("Vector %d: norms %e %e\n", i, norm2(*param.B[i]), norm2(tmp1));
printfQuda("relative residual = %e\n", sqrt(xmyNorm(tmp2, tmp1) / norm2(*param.B[i])));
}
}
sprintf(prefix, "MG level %d (%s): ", param.level + 1,
param.location == QUDA_CUDA_FIELD_LOCATION ? "GPU" : "CPU");
setOutputPrefix(prefix);
}
}
#if 0
if (getVerbosity() >= QUDA_SUMMARIZE)
printfQuda("Checking 1 > || (1 - D P (P^\\dagger D P) P^\\dagger v_k || / || v_k || for %d vectors\n",
param.Nvec);
for (int i=0; i<param.Nvec; i++) {
transfer->R(*r_coarse, *(param.B[i]));
(*coarse)(*x_coarse, *r_coarse); // this needs to be an exact solve to pass
setOutputPrefix(prefix); // restore output prefix
transfer->P(tmp2, *x_coarse);
param.matResidual(tmp1, tmp2);
tmp2 = *(param.B[i]);
if (getVerbosity() >= QUDA_VERBOSE) {
printfQuda("Vector %d: norms %e %e ", i, norm2(*param.B[i]), norm2(tmp1));
printfQuda("relative residual = %e\n", sqrt(xmyNorm(tmp2, tmp1) / norm2(*param.B[i])) );
}
}
#endif
if (getVerbosity() >= QUDA_SUMMARIZE) printfQuda("Checking 0 = (1 - P^\\dagger P) eta_c\n");
spinorNoise(*x_coarse, *rng, QUDA_NOISE_UNIFORM);
transfer->P(tmp2, *x_coarse);
transfer->R(*r_coarse, tmp2);
if (getVerbosity() >= QUDA_VERBOSE)
printfQuda("L2 norms %e %e (fine tmp %e) ", norm2(*x_coarse), norm2(*r_coarse), norm2(tmp2));
deviation = sqrt( xmyNorm(*x_coarse, *r_coarse) / norm2(*x_coarse) );
if (getVerbosity() >= QUDA_VERBOSE) printfQuda("relative deviation = %e\n", deviation);
if (deviation > tol ) errorQuda("L2 relative deviation = %e > %e failed", deviation, tol);
if (getVerbosity() >= QUDA_SUMMARIZE) printfQuda("Checking 0 = (D_c - P^\\dagger D P) (native coarse operator to emulated operator)\n");
zero(*tmp_coarse);
zero(*r_coarse);
#if 0 // debugging trick: point source matrix elements
tmp_coarse->Source(QUDA_POINT_SOURCE, 0, 0, 0);
#else
spinorNoise(*tmp_coarse, *rng, QUDA_NOISE_UNIFORM);
#endif
// put a non-trivial vector on the fine level as well
transfer->P(tmp1, *tmp_coarse);
// the three-hop terms in ASQTAD can break the verification depending on how we're coarsening the operator
// and if the aggregate size is too small in a direction
bool can_verify = true;
if ((param.transfer_type == QUDA_TRANSFER_OPTIMIZED_KD || param.transfer_type == QUDA_TRANSFER_OPTIMIZED_KD_DROP_LONG)
&& (diracSmoother->getDiracType() == QUDA_STAGGERED_DIRAC
|| diracSmoother->getDiracType() == QUDA_STAGGEREDPC_DIRAC || diracSmoother->getDiracType() == QUDA_ASQTAD_DIRAC
|| diracSmoother->getDiracType() == QUDA_ASQTADPC_DIRAC)) {
// If we're doing an optimized build with the staggered operator, we need to skip the verify on level 0
can_verify = false;
if (getVerbosity() >= QUDA_VERBOSE)
printfQuda("Intentionally skipping staggered -> staggered KD verify because it's not a \"real\" coarsen\n");
} else if (diracSmoother->getDiracType() == QUDA_ASQTAD_DIRAC || diracSmoother->getDiracType() == QUDA_ASQTADKD_DIRAC
|| diracSmoother->getDiracType() == QUDA_ASQTADPC_DIRAC) {
// If we're doing anything with the asqtad operator, the long links can make verification difficult
if (param.transfer_type == QUDA_TRANSFER_COARSE_KD || param.transfer_type == QUDA_TRANSFER_OPTIMIZED_KD_DROP_LONG) {
can_verify = false;
if (getVerbosity() >= QUDA_VERBOSE)
printfQuda("Using the naively coarsened KD operator with asqtad long links, skipping verify...\n");
} else if (param.transfer_type == QUDA_TRANSFER_AGGREGATE || param.transfer_type == QUDA_TRANSFER_OPTIMIZED_KD) {
// need to see if the aggregate is smaller than 3 in any direction
for (int d = 0; d < 4; d++) {
if (param.mg_global.geo_block_size[param.level][d] < 3) {
can_verify = false;
if (getVerbosity() >= QUDA_VERBOSE)
printfQuda("Aggregation geo_block_size[%d] = %d is less than 3, skipping verify for asqtad coarsen...\n",
d, param.mg_global.geo_block_size[param.level][d]);
}
}
}
}
if (can_verify) {
if (param.coarse_grid_solution_type == QUDA_MATPC_SOLUTION && param.smoother_solve_type == QUDA_DIRECT_PC_SOLVE) {
double kappa = diracResidual->Kappa();
double mass = diracResidual->Mass();
if (param.level == 0) {
if (tmp1.Nspin() == 4) {
diracSmoother->DslashXpay(tmp2.Even(), tmp1.Odd(), QUDA_EVEN_PARITY, tmp1.Even(), -kappa);
diracSmoother->DslashXpay(tmp2.Odd(), tmp1.Even(), QUDA_ODD_PARITY, tmp1.Odd(), -kappa);
} else if (tmp1.Nspin() == 2) { // if the coarse op is on top
diracSmoother->DslashXpay(tmp2.Even(), tmp1.Odd(), QUDA_EVEN_PARITY, tmp1.Even(), 1.0);
diracSmoother->DslashXpay(tmp2.Odd(), tmp1.Even(), QUDA_ODD_PARITY, tmp1.Odd(), 1.0);
} else { // staggered
diracSmoother->DslashXpay(tmp2.Even(), tmp1.Odd(), QUDA_EVEN_PARITY, tmp1.Even(),
2.0 * mass); // stag convention
diracSmoother->DslashXpay(tmp2.Odd(), tmp1.Even(), QUDA_ODD_PARITY, tmp1.Odd(),
2.0 * mass); // stag convention
}