forked from etmc/tmLQCD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
DDalphaAMG_interface.c
1399 lines (1214 loc) · 55.5 KB
/
DDalphaAMG_interface.c
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
/***********************************************************************
*
* Copyright (C) 2016 Simone Bacchio, Jacob Finkenrath
*
* This file is part of tmLQCD.
*
* tmLQCD 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 3 of the License, or
* (at your option) any later version.
*
* tmLQCD 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 tmLQCD. If not, see <http://www.gnu.org/licenses/>.
*
* Interface for DDalphaAMG
*
*******************************************************************************/
#include "DDalphaAMG_interface.h"
#ifndef DDalphaAMG
int mg_setup_iter;
int mg_coarse_setup_iter;
int mg_update_setup_iter;
int mg_update_gauge;
int mg_omp_num_threads;
int mg_Nvec;
int mg_lvl;
int mg_blk[4];
int mg_mixed_prec;
int mg_setup_mu_set;
int mg_no_shifts = 0;
double mg_mms_mass = 0;
double mg_setup_mu;
double mg_cmu_factor;
double mg_dtau_update;
double mg_rho_update;
void MG_init(void) {
printf("ERROR: MG_init called but DDalphaAMG library not included.\n");
exit(1);
}
void MG_update_gauge(double step) {
printf("ERROR: MG_update_gauge called but DDalphaAMG library not included.\n");
exit(1);
}
void MG_update_mu(double mu_tmLQCD, double odd_tmLQCD) {
printf("ERROR: MG_update_mu called but DDalphaAMG library not included.\n");
exit(1);
}
void MG_reset(void) {
printf("ERROR: MG_reset called but DDalphaAMG library not included.\n");
exit(1);
}
void MG_finalize(void) {
printf("ERROR: MG_finalize called but DDalphaAMG library not included.\n");
exit(1);
}
int MG_solver(spinor * const phi_new, spinor * const phi_old,
const double precision, const int max_iter,const int rel_prec,
const int N, su3 **gf, matrix_mult f) {
printf("ERROR: MG_solver called but DDalphaAMG library not included.\n");
exit(1);
}
int MG_solver_eo(spinor * const Even_new, spinor * const Odd_new,
spinor * const Even, spinor * const Odd,
const double precision, const int max_iter, const int rel_prec,
const int N, su3 **gf, matrix_mult_full f_full) {
printf("ERROR: MG_solver_eo called but DDalphaAMG library not included.\n");
exit(1);
}
int MG_solver_nd(spinor * const up_new, spinor * const dn_new,
spinor * const up_old, spinor * const dn_old,
const double precision, const int max_iter, const int rel_prec,
const int N, su3 **gf, matrix_mult_nd f) {
printf("ERROR: MG_solver_nd called but DDalphaAMG library not included.\n");
exit(1);
}
#else
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include "boundary.h"
#include "linalg/convert_eo_to_lexic.h"
#include "solver/solver.h"
#include "solver/solver_field.h"
#include "gettime.h"
#include "read_input.h"
#include "DDalphaAMG.h"
#include "linalg_eo.h"
#include "phmc.h"
#include "operator/D_psi.h"
#include "operator/tm_operators.h"
#include "operator/tm_operators_nd.h"
#include "operator/clovertm_operators.h"
#include "operator/Hopping_Matrix.h"
//Enable variant for shifted operator in the ND sector.
//The variant is used in case of initial guess for the squared operator.
//It is faster and tests prove also to be safe (read Appendix A of arxiv:1801.##### by S.Bacchio et al.)
#define VARIANT_FOR_SHIFTED
DDalphaAMG_init mg_init;
DDalphaAMG_parameters mg_params;
DDalphaAMG_status mg_status;
int mg_do_setup=1; //if one do or redo the setup
int mg_update_gauge=1; //set to zero if gaugefield is up to date, set to one if it has to be updated
int mg_update_setup=0; //Number of additional setup iteration
int mg_initialized=0;
int mg_setup_iter=5;
int mg_coarse_setup_iter=3;
int mg_update_setup_iter=1;
int mg_omp_num_threads=0;
int mg_Nvec=24;
int mg_lvl=3;
int mg_blk[4] = {0, 0, 0, 0};
int mg_mixed_prec=0;
int mg_setup_mu_set = 0; //flag that enable the use of mg_setup_mu in the setup phase
int mg_no_shifts = 0; // number of shifts to invert with MG. MMS-CG is used for the others at larger mass.
double mg_mms_mass = 0.1; // mass shift value for switching from MMS-CG to MG. MMS-CG is used for larger masses than the value.
double mg_setup_mu = 0.;
double mg_cmu_factor = 1.0;
double mg_dtau_update = 0.0;
double mg_rho_update = -1.0;
double mg_tau = 0.0;
double gauge_tau = 0.0;
static int Cart_rank(MPI_Comm comm, const int coords[], int *rank)
{
int coords_l[4];
coords_l[0]=coords[0];
coords_l[1]=coords[3];
coords_l[2]=coords[2];
coords_l[3]=coords[1];
return MPI_Cart_rank(comm, coords_l, rank);
}
static int Cart_coords(MPI_Comm comm, int rank, int maxdims, int coords[])
{
int stat;
stat=MPI_Cart_coords(comm, rank, maxdims, coords);
int tmp=coords[1];
coords[1]=coords[3];
coords[3]=tmp;
return stat;
}
static int conf_index_fct(int t, int z, int y, int x, int mu)
{
int id;
id=(g_ipt[t][x][y][z])*72; //9*2*4
id+=18*((mu%2==0)?mu:((mu==1)?3:1));//9*2
return id;
}
static int vector_index_fct(int t, int z, int y, int x )
{
int id;
id=24*(g_ipt[t][x][y][z]);
return id;
}
static inline int MG_check(spinor * const phi_new, spinor * const phi_old, const int N, const double precision, matrix_mult f)
{
double differ[2], residual;
spinor ** check_vect = NULL;
double acc_factor = 2;
init_solver_field(&check_vect, VOLUMEPLUSRAND,1);
f( check_vect[0], phi_new);
diff( check_vect[0], check_vect[0], phi_old, N);
differ[0] = sqrt(square_norm(check_vect[0], N, 1));
differ[1] = sqrt(square_norm(phi_old, N, 1));
finalize_solver(check_vect, 1);
residual = differ[0]/differ[1];
if( residual > precision && residual < acc_factor*precision ) {
if(g_proc_id == 0)
printf("WARNING: solution accepted even if the residual wasn't complitely acceptable (%e > %e) \n", residual, precision);
} else if( residual > acc_factor*precision ) {
if(g_proc_id == 0) {
printf("ERROR: something bad happened... MG converged giving the wrong solution!! Trying to restart... \n");
printf("ERROR contd: || s - f_{tmLQC} * f_{DDalphaAMG}^{-1} * s || / ||s|| = %e / %e = %e > %e \n", differ[0],differ[1],differ[0]/differ[1],precision);
}
return 0;
}
if (g_debug_level > 0 && g_proc_id == 0)
printf("MGTEST: || s - f_{tmLQC} * f_{DDalphaAMG}^{-1} * s || / ||s|| = %e / %e = %e \n", differ[0],differ[1],differ[0]/differ[1]);
return 1;
}
static inline int MG_check_nd( spinor * const up_new, spinor * const dn_new, spinor * const up_old, spinor * const dn_old,
const int N, const double precision, matrix_mult_nd f)
{
double differ[2], residual;
spinor ** check_vect = NULL;
double acc_factor = 4;
#ifdef VARIANT_FOR_SHIFTED
if(( f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
&& g_shift != 0 )
acc_factor = 1/sqrt(phmc_cheb_evmin/phmc_cheb_evmax + g_shift);
#endif
init_solver_field(&check_vect, VOLUMEPLUSRAND,2);
f( check_vect[0], check_vect[1], up_new, dn_new);
diff( check_vect[0], check_vect[0], up_old, N);
diff( check_vect[1], check_vect[1], dn_old, N);
differ[0] = sqrt(square_norm(check_vect[0], N, 1)+square_norm(check_vect[1], N, 1));
differ[1] = sqrt(square_norm(up_old, N, 1)+square_norm(dn_old, N, 1));
finalize_solver(check_vect, 2);
residual = differ[0]/differ[1];
if( residual > precision && residual < acc_factor*precision ) {
if(g_proc_id == 0)
printf("WARNING: solution accepted even if the residual wasn't complitely acceptable (%e > %e). Max acc. factor %f.\n", residual, precision, acc_factor);
} else if( residual > acc_factor*precision ) {
if(g_proc_id == 0) {
printf("ERROR: something bad happened... MG converged giving the wrong solution!! Trying to restart... \n");
printf("ERROR contd: || s - f_{tmLQC} * f_{DDalphaAMG}^{-1} * s || / ||s|| = %e / %e = %e > %e \n", differ[0],differ[1],differ[0]/differ[1],precision);
}
return 0;
}
if (g_debug_level > 0 && g_proc_id == 0)
printf("MGTEST: || s - f_{tmLQC} * f_{DDalphaAMG}^{-1} * s || / ||s|| = %e / %e = %e \n", differ[0],differ[1],differ[0]/differ[1]);
return 1;
}
static inline int MG_mms_check_nd( spinor **const up_new, spinor **const dn_new,
spinor * const up_old, spinor * const dn_old,
const double * shifts, const int no_shifts,
const int N, double * precision, matrix_mult_nd f)
{
double differ[2], residual;
spinor ** check_vect = NULL;
double acc_factor = 2;
init_solver_field(&check_vect, VOLUMEPLUSRAND,2);
for( int i = 0; i < no_shifts; i++ ) {
g_shift = shifts[i]*shifts[i];
f( check_vect[0], check_vect[1], up_new[i], dn_new[i]);
diff( check_vect[0], check_vect[0], up_old, N);
diff( check_vect[1], check_vect[1], dn_old, N);
differ[0] = sqrt(square_norm(check_vect[0], N, 1)+square_norm(check_vect[1], N, 1));
differ[1] = sqrt(square_norm(up_old, N, 1)+square_norm(dn_old, N, 1));
residual = differ[0]/differ[1];
if( residual > precision[i] && residual < acc_factor*precision[i] ) {
if(g_proc_id == 0)
printf("WARNING: solution accepted even if the residual wasn't complitely acceptable (%e > %e) \n", residual, precision[i]);
} else if( residual > acc_factor*precision[i] ) {
if(g_proc_id == 0) {
printf("ERROR: something bad happened... MG converged giving the wrong solution!! Trying to restart... \n");
printf("ERROR contd: || s - f_{tmLQC} * f_{DDalphaAMG}^{-1} * s || / ||s|| = %e / %e = %e > %e \n", differ[0],differ[1],differ[0]/differ[1],precision[i]);
}
finalize_solver(check_vect, 2);
return 0;
}
if (g_debug_level > 0 && g_proc_id == 0)
printf("MGTEST: || s - f_{tmLQC} * f_{DDalphaAMG}^{-1} * s || / ||s|| = %e / %e = %e \n", differ[0],differ[1],differ[0]/differ[1]);
}
finalize_solver(check_vect, 2);
return 1;
}
static int MG_pre_solve( su3 **gf )
{
double dtau = fabs(mg_tau-gauge_tau);
// Checking if:
// mg_update_setup < mg_update_setup_iter : maybe you want to do more iteration at this run
// mg_dtau_update < dtau : regular condition for update of setup
// mg_dtau_update < -dtau : during reversability check dtau is negative!
// mg_dtau_update == 0.0 : updating at every change of configuration -> valid as well if configuration changed outside the HMC
// mg_rho_update < 0.0 : parameter ignore
// mg_rho_update == rho : updating only if this condition and the others are satisfied
if ( mg_do_setup == 0 && mg_update_setup < mg_update_setup_iter &&
( (mg_dtau_update > 0.0 && dtau > 0.0 && mg_dtau_update < dtau+1e-6) ||
(mg_dtau_update == 0.0 && mg_update_gauge == 1) ||
(mg_rho_update >= 0.0 && mg_rho_update == g_mu3) ))
mg_update_setup = mg_update_setup_iter;
if(g_debug_level > 0 && g_proc_id == 0)
printf("Tau has been increased since last MG setup update of %e\n", dtau);
if (mg_initialized==0) {
MG_init();
mg_initialized = 1;
if (g_proc_id == 0)
printf("DDalphaAMG initialized\n");
MPI_Barrier(MPI_COMM_WORLD);
}
if (mg_update_gauge==1) {
DDalphaAMG_set_configuration( (double*) &(gf[0][0]), &mg_status );
mg_update_gauge = 0;
if (mg_status.success && g_proc_id == 0)
printf("DDalphaAMG cnfg set, plaquette %e\n", mg_status.info);
else if ( g_proc_id == 0)
printf("ERROR: configuration updating did not run correctly");
}
if (mg_do_setup==1) {
if( mg_setup_mu_set ) {
if (g_proc_id == 0)
printf("DDalphaAMG using mu=%f during setup\n", mg_setup_mu);
MG_update_mu(mg_setup_mu, 0);
} else
MG_update_mu(g_mu, 0);
if (g_proc_id == 0)
printf("DDalphaAMG running setup\n");
DDalphaAMG_setup(&mg_status);
mg_do_setup = 0;
mg_tau = gauge_tau;
if (mg_status.success && g_proc_id == 0)
printf("DDalphaAMG setup ran, time %.2f sec (%.2f %% on coarse grid)\n",
mg_status.time, 100.*(mg_status.coarse_time/mg_status.time));
else if ( g_proc_id == 0)
printf("ERROR: setup procedure did not run correctly");
}
if (mg_update_setup>0) {
if( mg_setup_mu_set ) {
if (g_proc_id == 0)
printf("DDalphaAMG using mu=%f during setup\n", mg_setup_mu);
MG_update_mu(mg_setup_mu, 0);
} else
MG_update_mu(g_mu, 0);
if (g_proc_id == 0)
printf("DDalphaAMG updating setup\n");
DDalphaAMG_update_setup(mg_update_setup, &mg_status);
mg_update_setup = 0;
mg_tau = gauge_tau;
if (mg_status.success && g_proc_id == 0)
printf("DDalphaAMG setup ran, time %.2f sec (%.2f %% on coarse grid)\n",
mg_status.time, 100.*(mg_status.coarse_time/mg_status.time));
else if ( g_proc_id == 0)
printf("ERROR: setup updating did not run correctly");
}
return mg_status.success;
}
static int MG_solve(spinor * const phi_new, spinor * const phi_old, const double precision,
const int N, matrix_mult f)
{
// for rescaling convention in DDalphaAMG: (4+m)*\delta_{x,y} in tmLQCD: 1*\delta_{x,y} -> rescale by 1/4+m
double mg_scale=0.5/g_kappa;
double *old = (double*) phi_old;
double *new = (double*) phi_new;
spinor ** solver_field = NULL;
if( N != VOLUME && N != VOLUME/2 ) {
if( g_proc_id == 0 )
printf("ERROR: N = %d in MG_solve. Expettected N == VOLUME (%d) or VOLUME/2 (%d)\n", N, VOLUME, VOLUME/2);
return 0;
}
if (N==VOLUME/2) {
init_solver_field(&solver_field, VOLUMEPLUSRAND,2);
old = (double*) solver_field[0];
new = (double*) solver_field[1];
convert_odd_to_lexic( (spinor*) old, phi_old);
}
// Checking if the operator is in the list and compatible with N
if ( f == Msw_psi || // Schur complement with mu=0 on odd sites
f == Qsw_psi || // Gamma5 - Schur complement with mu=0 on odd sites
f == Mtm_plus_psi || // Schur complement with plus mu
f == Msw_plus_psi || // Schur complement with plus mu
f == Qtm_plus_psi || // Gamma5 - Schur complement with plus mu
f == Qsw_plus_psi || // Gamma5 - Schur complement with plus mu
f == Mtm_minus_psi || // Schur complement with minus mu
f == Msw_minus_psi || // Schur complement with minus mu
f == Qtm_minus_psi || // Gamma5 - Schur complement with minus mu
f == Qsw_minus_psi || // Gamma5 - Schur complement with minus mu
f == Qtm_pm_psi || // Schur complement squared
f == Qsw_pm_psi ) { // Schur complement squared
if( N != VOLUME/2 && g_proc_id == 0 )
printf("WARNING: expected N == VOLUME/2 for the required operator in MG_solve. Continuing with N == VOLUME\n");
}
else if ( f == D_psi || // Full operator with plus mu
f == Q_plus_psi || // Gamma5 - Full operator with plus mu
f == Q_minus_psi || // Gamma5 - Full operator with minus mu
f == Q_pm_psi || // Full operator squared
f == Qsw_full_plus_psi || // Gamma5 - Full operator with plus mu
f == Qsw_full_minus_psi|| //Gamma5 - Full operator with plus mu
f == Qsw_full_pm_psi || // Full operator squared
f == Msw_full_minus_psi) {// Full operator with minus mu
if( N != VOLUME && g_proc_id == 0 )
printf("WARNING: expected N == VOLUME for the required operator in MG_solve. Continuing with N == VOLUME/2\n");
}
else if( g_proc_id == 0 )
printf("WARNING: required operator unknown for MG_solve. Using standard operator: %s.\n",
N==VOLUME?"D_psi":"Msw_plus_psi");
// Setting mu
if ( f == Msw_psi || // Schur complement with mu=0 on odd sites
f == Qsw_psi ) // Gamma5 - Schur complement with mu=0 on odd sites
MG_update_mu(g_mu, -g_mu);
else if ( f == Mtm_minus_psi || // Schur complement with minus mu
f == Msw_minus_psi || // Schur complement with minus mu
f == Qtm_minus_psi || // Gamma5 - Schur complement with minus mu
f == Qsw_minus_psi || // Gamma5 - Schur complement with minus mu
f == Qsw_full_minus_psi|| //Gamma5 - Full operator with plus mu
f == Msw_full_minus_psi|| // Full operator with minus mu
f == Q_minus_psi ) // Gamma5 - Full operator with minus mu
MG_update_mu(-g_mu, -g_mu3);
else if ( f == Mtm_plus_psi || // Schur complement with plus mu
f == Msw_plus_psi || // Schur complement with plus mu
f == Qtm_plus_psi || // Gamma5 - Schur complement with plus mu
f == Qsw_plus_psi || // Gamma5 - Schur complement with plus mu
f == D_psi || // Full operator with plus mu
f == Q_plus_psi || // Gamma5 - Full operator with plus mu
f == Qtm_pm_psi || // Schur complement squared
f == Qsw_pm_psi || // Schur complement squared
f == Qsw_full_plus_psi || // Gamma5 - Full operator with plus mu
f == Qsw_full_pm_psi || // Full operator squared
f == Q_pm_psi ) // Full operator squared
MG_update_mu(g_mu, g_mu3);
else
MG_update_mu(g_mu, g_mu3);
//Solving
if ( f == Qtm_plus_psi || // Gamma5 - Schur complement with plus mu
f == Qsw_plus_psi || // Gamma5 - Schur complement with plus mu
f == Qtm_minus_psi || // Gamma5 - Schur complement with minus mu
f == Qsw_minus_psi || // Gamma5 - Schur complement with minus mu
f == Qsw_psi || // Gamma5 - Schur complement with mu=0 on odd sites
f == Q_plus_psi || // Gamma5 - Full operator with plus mu
f == Q_minus_psi || // Gamma5 - Full operator with minus mu
f == Qsw_full_plus_psi || // Gamma5 - Full operator with plus mu
f == Qsw_full_minus_psi|| //Gamma5 - Full operator with plus mu
f == Qsw_full_pm_psi ) { // Full operator squared
mul_gamma5((spinor *const) old, VOLUME);
DDalphaAMG_solve( new, old, precision, &mg_status );
if( N == VOLUME ) // in case of VOLUME/2 old is a just local vector
mul_gamma5((spinor *const) old, VOLUME);
}
else if ( f == Qtm_pm_psi || // Schur complement squared
f == Qsw_pm_psi ) { // Schur complement squared
mg_scale *= mg_scale;
DDalphaAMG_solve_squared_odd( new, old, precision, &mg_status );
}
else if ( f == Q_pm_psi ) { // Full operator squared
mg_scale *= mg_scale;
DDalphaAMG_solve_squared( new, old, precision, &mg_status );
}
else if ( f == Mtm_plus_psi || // Schur complement with plus mu
f == Msw_plus_psi || // Schur complement with plus mu
f == Mtm_minus_psi || // Schur complement with minus mu
f == Msw_minus_psi || // Schur complement with minus mu
f == Msw_psi || // Schur complement with mu=0 on odd sites
f == D_psi || // Full operator with plus mu
f == Msw_full_minus_psi) {// Full operator with minus mu
DDalphaAMG_solve( new, old, precision, &mg_status );
}
else
DDalphaAMG_solve( new, old, precision, &mg_status );
if (N==VOLUME/2) {
convert_lexic_to_odd(phi_new, (spinor*) new);
finalize_solver(solver_field, 2);
}
mul_r(phi_new ,mg_scale, phi_new, N);
if (g_proc_id == 0) {
printf("Solving time %.2f sec (%.1f %% on coarse grid)\n", mg_status.time,
100.*(mg_status.coarse_time/mg_status.time));
printf("Total iterations on fine grid %d\n", mg_status.iter_count);
printf("Total iterations on coarse grids %d\n", mg_status.coarse_iter_count);
if (!mg_status.success)
printf("ERROR: the solver did not converge!\n");
}
return mg_status.success;
}
static int MG_solve_nd( spinor * up_new, spinor * dn_new, spinor * const up_old, spinor * const dn_old,
const double precision, const int N, matrix_mult_nd f)
{
// for rescaling convention in DDalphaAMG: (4+m)*\delta_{x,y} in tmLQCD: 1*\delta_{x,y} -> rescale by 1/4+m
// moreover in the nd case, the tmLQCD is multiplied by phmc_invmaxev
double mg_scale=0.5/g_kappa/phmc_invmaxev;
double sqnorm;
int init_guess = 0;
spinor *old1 = up_old;
spinor *old2 = dn_old;
spinor *new1 = up_new, *new1tmp;
spinor *new2 = dn_new, *new2tmp;
spinor ** solver_field = NULL, ** oe_solver_field = NULL;
int no_solver_field = 0;
if( N != VOLUME && N != VOLUME/2 ) {
if( g_proc_id == 0 )
printf("ERROR: N = %d in MG_solve. Expettected N == VOLUME (%d) or VOLUME/2 (%d)\n", N, VOLUME, VOLUME/2);
return 0;
}
if (N==VOLUME/2) no_solver_field += 4;
// Checking if initial guess is given
sqnorm = square_norm(up_new, N, 1);
sqnorm += square_norm(dn_new, N, 1);
if ( sqnorm > 1e-14 ) init_guess = 1;
// In case of initial guess and squared operator, we do the inversion in two step and we need two more vectors
if ( init_guess && (
f == Qtm_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared
f == Qsw_pm_ndpsi_shift )) // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
no_solver_field += 2;
// Allocating and assigning fields
if(no_solver_field>0)
init_solver_field(&solver_field, VOLUMEPLUSRAND,no_solver_field);
int assign_solver_field = 0;
if (N==VOLUME/2) {
old1 = solver_field[assign_solver_field++];
old2 = solver_field[assign_solver_field++];
new1 = solver_field[assign_solver_field++];
new2 = solver_field[assign_solver_field++];
convert_odd_to_lexic(old1, up_old);
convert_odd_to_lexic(old2, dn_old);
set_even_to_zero(old1);
set_even_to_zero(old2);
}
if ( init_guess && (
f == Qtm_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared
f == Qsw_pm_ndpsi_shift )) {// (Gamma5 Dh tau1)^2 - Schur complement squared with shift
new1tmp = solver_field[assign_solver_field++];
new2tmp = solver_field[assign_solver_field++];
}
// Reconstructing initial guess in case of oe
if ( init_guess && N==VOLUME/2 ) {
init_solver_field(&oe_solver_field, VOLUMEPLUSRAND, 4);
spinor* tmp11 = oe_solver_field[0];
spinor* tmp21 = oe_solver_field[1];
spinor* tmp12 = oe_solver_field[2];
spinor* tmp22 = oe_solver_field[3];
if (g_debug_level > 2) {
double differ[2];
f( tmp11, tmp12, up_new, dn_new);
diff( tmp11, tmp11, up_old, N);
diff( tmp12, tmp12, dn_old, N);
differ[0] = sqrt(square_norm(tmp11, N, 1)+square_norm(tmp12, N, 1));
differ[1] = sqrt(square_norm(up_old, N, 1)+square_norm(dn_old, N, 1));
if(g_proc_id == 0)
printf("MG TEST: using initial guess. Relative residual = %e \n", differ[0]/differ[1]);
}
/* Reconstruct the even sites */
if ( f == Qtm_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
f == Qsw_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
f == Qtm_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and plus shift
f == Qtm_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and minus shift
f == Qsw_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with plus shift
f == Qsw_tau1_ndpsi_sub_Ishift ) {// Gamma5 Dh tau1 - Schur complement with minus shift
#ifdef VARIANT_FOR_SHIFTED
if(( f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
&& g_shift != 0 ) {
if( f == Qtm_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
Qtm_tau1_ndpsi_add_Ishift(tmp12, tmp11, up_new, dn_new); // tau1 exchange tmp11 <-> tmp12
} else { // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
Qsw_tau1_ndpsi_add_Ishift(tmp12, tmp11, up_new, dn_new); // tau1 exchange tmp11 <-> tmp12
}
// tau1 exchange new1tmp <-> new2tmp
convert_odd_to_lexic( new2, tmp11);
convert_odd_to_lexic( new1, tmp12);
Hopping_Matrix(EO, tmp21, tmp11);
Hopping_Matrix(EO, tmp22, tmp12);
Msw_ee_inv_ndpsi(tmp11, tmp12, tmp21, tmp22);
convert_even_to_lexic(new2, tmp11);
convert_even_to_lexic(new1, tmp12);
} else
#endif
{
// tau1 exchange tmp11 <-> tmp12
Hopping_Matrix(EO, tmp12, up_new);
Hopping_Matrix(EO, tmp11, dn_new);
Msw_ee_inv_ndpsi(tmp21, tmp22, tmp11, tmp12);
/* Assigning with plus sign for the even
* since in Hopping_Matrix the minus is missing
*/
// tau1 exchange tmp22 <-> tmp21
convert_eo_to_lexic(new1, tmp22, up_new);
convert_eo_to_lexic(new2, tmp21, dn_new);
}
} else {
Hopping_Matrix(EO, tmp11, up_new);
Hopping_Matrix(EO, tmp12, dn_new);
Msw_ee_inv_ndpsi(tmp21, tmp22, tmp11, tmp12);
/* Assigning with plus sign for the even
* since in Hopping_Matrix the minus is missing
*/
convert_eo_to_lexic(new1, tmp21, up_new);
convert_eo_to_lexic(new2, tmp22, dn_new);
}
// if squared obtaining initial guess for Gamma5 Dh
if ( f == Qtm_pm_ndpsi ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
Qtm_dagger_ndpsi(tmp11, tmp12, up_new, dn_new); // tau1 Gamma5 Dh tau1
}
else if(f == Qsw_pm_ndpsi ) { // (Gamma5 Dh tau1)^2 - Schur complement squared
Qsw_dagger_ndpsi(tmp11, tmp12, up_new, dn_new); // tau1 Gamma5 Dh tau1
}
else if(f == Qtm_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
Qtm_tau1_ndpsi_sub_Ishift(tmp12, tmp11, up_new, dn_new); // tau1 exchange tmp11 <-> tmp12
}
else if(f == Qsw_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
Qsw_tau1_ndpsi_sub_Ishift(tmp12, tmp11, up_new, dn_new); // tau1 exchange tmp11 <-> tmp12
}
if ( f == Qtm_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared
f == Qsw_pm_ndpsi_shift ){ // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
// tau1 exchange new1tmp <-> new2tmp
convert_odd_to_lexic( new2tmp, tmp11);
convert_odd_to_lexic( new1tmp, tmp12);
Hopping_Matrix(EO, tmp21, tmp11);
Hopping_Matrix(EO, tmp22, tmp12);
Msw_ee_inv_ndpsi(tmp11, tmp12, tmp21, tmp22);
convert_even_to_lexic(new2tmp, tmp11);
convert_even_to_lexic(new1tmp, tmp12);
}
finalize_solver(oe_solver_field, 4);
}
// Checking if the operator is in the list and compatible with N
if ( f == Qtm_ndpsi || // Gamma5 Dh - Schur complement with csw = 0
f == Qsw_ndpsi || // Gamma5 Dh - Schur complement
f == Qtm_dagger_ndpsi || // Gamma5 Dh - Schur complement with mu = -mubar and csw = 0
f == Qsw_dagger_ndpsi || // Gamma5 Dh - Schur complement with mu = -mubar
f == Qtm_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and plus shift
f == Qtm_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and minus shift
f == Qsw_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with plus shift
f == Qsw_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with minus shift
f == Qtm_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared
f == Qsw_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
if( N != VOLUME/2 && g_proc_id == 0 )
printf("WARNING: expected N == VOLUME/2 for the required operator in MG_solve. Continuing with N == VOLUME\n");
}
else if ( f == D_ndpsi ) { // Dh
if( N != VOLUME && g_proc_id == 0 )
printf("WARNING: expected N == VOLUME for the required operator in MG_solve. Continuing with N == VOLUME/2\n");
}
else if( g_proc_id == 0 )
printf("WARNING: required operator unknown for MG_solve. Using standard operator: %s.\n",
N==VOLUME?"":"Qsw_ndpsi");
// Setting mu and eps
if ( f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
MG_update_mubar_epsbar( g_mubar, g_epsbar, sqrt(g_shift) );
else if ( f == Qtm_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and plus shift
f == Qsw_tau1_ndpsi_add_Ishift ) // Gamma5 Dh tau1 - Schur complement with plus shift
MG_update_mubar_epsbar( g_mubar, g_epsbar, sqrt(g_shift) );
else if ( f == Qtm_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and minus shift
f == Qsw_tau1_ndpsi_sub_Ishift ) // Gamma5 Dh tau1 - Schur complement with minus shift
MG_update_mubar_epsbar( g_mubar, g_epsbar, -sqrt(g_shift) );
else if ( f == Qtm_dagger_ndpsi || // Gamma5 Dh - Schur complement with mu = -mubar csw = 0
f == Qsw_dagger_ndpsi ) // Gamma5 Dh - Schur complement with mu = -mubar
MG_update_mubar_epsbar( -g_mubar, g_epsbar, 0 );
else if ( f == Qtm_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
f == Qsw_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared
f == D_ndpsi ) // Dh
MG_update_mubar_epsbar( g_mubar, g_epsbar, 0 );
else
MG_update_mubar_epsbar( g_mubar, g_epsbar, 0 );
//Solving
if ( f == Qtm_ndpsi || // Gamma5 Dh - Schur complement with csw = 0
f == Qsw_ndpsi || // Gamma5 Dh - Schur complement
f == Qtm_dagger_ndpsi || // Gamma5 Dh - Schur complement with mu = -mubar csw = 0
f == Qsw_dagger_ndpsi ) { // Gamma5 Dh - Schur complement with mu = -mubar
mul_gamma5(old1, VOLUME);
mul_gamma5(old2, VOLUME);
if (init_guess) {
// Removing normalization from initial guess
mul_r(new1, 1/mg_scale, new1, VOLUME);
mul_r(new2, 1/mg_scale, new2, VOLUME);
DDalphaAMG_solve_doublet_with_guess( (double*) new1, (double*) old1, (double*) new2, (double*) old2,
precision, &mg_status );
} else {
DDalphaAMG_solve_doublet( (double*) new1, (double*) old1, (double*) new2, (double*) old2,
precision, &mg_status );
}
if( N == VOLUME ) { // in case of VOLUME/2 old is a just local vector
mul_gamma5(old1, VOLUME);
mul_gamma5(old2, VOLUME);
}
}
else if ( f == Qtm_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and plus shift
f == Qtm_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and minus shift
f == Qsw_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with plus shift
f == Qsw_tau1_ndpsi_sub_Ishift ) {// Gamma5 Dh tau1 - Schur complement with minus shift
mul_gamma5(old1, VOLUME);
mul_gamma5(old2, VOLUME);
// tau1 exchange new1 <-> new2
if (init_guess) {
// Removing normalization from initial guess
mul_r(new1, 1/mg_scale, new1, VOLUME);
mul_r(new2, 1/mg_scale, new2, VOLUME);
DDalphaAMG_solve_doublet_with_guess( (double*) new2, (double*) old1, (double*) new1, (double*) old2,
precision, &mg_status );
} else {
DDalphaAMG_solve_doublet( (double*) new2, (double*) old1, (double*) new1, (double*) old2,
precision, &mg_status );
}
if( N == VOLUME ) { // in case of VOLUME/2 old is a just local vector
mul_gamma5(old1, VOLUME);
mul_gamma5(old2, VOLUME);
}
}
else if ( f == Qtm_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi || // (Gamma5 Dh tau1)^2 - Schur complement squared
f == Qsw_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
// DDalphaAMG: tau1 gamma5 Dh tau1 gamma5 Dh
// tmLQCD: gamma5 Dh tau1 gamma5 Dh tau1
if (init_guess) {
mul_gamma5(old1, VOLUME);
mul_gamma5(old2, VOLUME);
// Removing normalization from initial guess
mul_r(new1tmp, 1/mg_scale, new1tmp, VOLUME);
mul_r(new2tmp, 1/mg_scale, new2tmp, VOLUME);
DDalphaAMG_solve_doublet_with_guess( (double*) new2tmp, (double*) old1, (double*) new1tmp, (double*) old2,
precision/2, &mg_status );
#ifdef VARIANT_FOR_SHIFTED
if(( f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
&& g_shift != 0 ) {
// Removing normalization from initial guess
mul_r(new1, 1/mg_scale, new1, VOLUME);
mul_r(new2, 1/mg_scale, new2, VOLUME);
MG_update_mubar_epsbar( g_mubar, g_epsbar, -sqrt(g_shift) );
DDalphaAMG_solve_doublet_with_guess( (double*) new2, (double*) old1, (double*) new1, (double*) old2,
precision/2, &mg_status );
assign_mul_add_mul(new1, -_Complex_I/2./sqrt(g_shift), new1tmp, _Complex_I/2./sqrt(g_shift), VOLUME);
assign_mul_add_mul(new2, -_Complex_I/2./sqrt(g_shift), new2tmp, _Complex_I/2./sqrt(g_shift), VOLUME);
} else
#endif
{
mul_gamma5(new1tmp, VOLUME);
mul_gamma5(new2tmp, VOLUME);
set_even_to_zero(new1tmp);
set_even_to_zero(new2tmp);
// Removing normalization from initial guess
mg_scale *= mg_scale;
mul_r(new1, 1/mg_scale, new1, VOLUME);
mul_r(new2, 1/mg_scale, new2, VOLUME);
if ( f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
MG_update_mubar_epsbar( g_mubar, g_epsbar, -sqrt(g_shift) );
DDalphaAMG_solve_doublet_with_guess( (double*) new2, (double*) new1tmp, (double*) new1, (double*) new2tmp,
precision/2, &mg_status );
}
if( N == VOLUME ) { // in case of VOLUME/2 old is a just local vector
mul_gamma5(old1, VOLUME);
mul_gamma5(old2, VOLUME);
}
} else {
mg_scale *= mg_scale;
DDalphaAMG_solve_doublet_squared_odd( (double*) new2, (double*) old2, (double*) new1, (double*) old1,
precision, &mg_status );
}
}
else if ( f == D_ndpsi ) { // Dh
if (init_guess) {
// Removing normalization from initial guess
mul_r(new1, 1/mg_scale, new1, VOLUME);
mul_r(new2, 1/mg_scale, new2, VOLUME);
DDalphaAMG_solve_doublet_with_guess( (double*) new1, (double*) old1, (double*) new2, (double*) old2,
precision, &mg_status );
} else {
DDalphaAMG_solve_doublet( (double*) new1, (double*) old1, (double*) new2, (double*) old2,
precision, &mg_status );
}
} else {
if (init_guess) {
// Removing normalization from initial guess
mul_r(new1, 1/mg_scale, new1, VOLUME);
mul_r(new2, 1/mg_scale, new2, VOLUME);
DDalphaAMG_solve_doublet_with_guess( (double*) new1, (double*) old1, (double*) new2, (double*) old2,
precision, &mg_status );
} else {
DDalphaAMG_solve_doublet( (double*) new1, (double*) old1, (double*) new2, (double*) old2,
precision, &mg_status );
}
}
if (N==VOLUME/2) {
convert_lexic_to_odd(up_new, new1);
convert_lexic_to_odd(dn_new, new2);
}
if (no_solver_field>0)
finalize_solver(solver_field, no_solver_field);
mul_r(up_new ,mg_scale, up_new, N);
mul_r(dn_new ,mg_scale, dn_new, N);
if (g_proc_id == 0) {
printf("Solving time %.2f sec (%.1f %% on coarse grid)\n", mg_status.time,
100.*(mg_status.coarse_time/mg_status.time));
printf("Total iterations on fine grid %d\n", mg_status.iter_count);
printf("Total iterations on coarse grids %d\n", mg_status.coarse_iter_count);
if (!mg_status.success)
printf("ERROR: the solver did not converge!\n");
}
return mg_status.success;
}
static int MG_mms_solve_nd( spinor **const up_new, spinor **const dn_new,
spinor * const up_old, spinor * const dn_old,
const double * shifts, const int no_shifts,
double * precision, const int N, matrix_mult_nd f)
{
// for rescaling convention in DDalphaAMG: (4+m)*\delta_{x,y} in tmLQCD: 1*\delta_{x,y} -> rescale by 1/4+m
// moreover in the nd case, the tmLQCD is multiplied by phmc_invmaxev
double mg_scale=0.5/g_kappa/phmc_invmaxev;
double *old1 = (double*) up_old;
double *old2 = (double*) dn_old;
double **new1, **new2, *mg_odd_shifts, *mg_even_shifts;
spinor ** solver_field = NULL;
// if( N != VOLUME && N != VOLUME/2 ) {
if( N != VOLUME/2 ) { // no full VOLUME functions implemented at the moment
if( g_proc_id == 0 )
printf("ERROR: N = %d in MG_solve. Expettected N == VOLUME (%d) or VOLUME/2 (%d)\n", N, VOLUME, VOLUME/2);
return 0;
}
new1 = (double**) malloc(no_shifts*sizeof(double*));
new2 = (double**) malloc(no_shifts*sizeof(double*));
mg_odd_shifts = (double*) malloc(no_shifts*sizeof(double));
mg_even_shifts = (double*) malloc(no_shifts*sizeof(double));
if( N==VOLUME/2 ) {
init_solver_field(&solver_field, VOLUMEPLUSRAND,2+2*no_shifts);
old1 = (double*) solver_field[0];
old2 = (double*) solver_field[1];
convert_odd_to_lexic( (spinor*) old1, up_old);
convert_odd_to_lexic( (spinor*) old2, dn_old);
for( int i = 0; i < no_shifts; i++ ) {
new1[i] = (double*) solver_field[2+2*i];
new2[i] = (double*) solver_field[3+2*i];
}
} else {
for( int i = 0; i < no_shifts; i++ ) {
new1[i] = (double*) up_new[i];
new2[i] = (double*) dn_new[i];
}
}
// Checking if the operator is in the list and compatible with N
if ( f == Qtm_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and plus shift
f == Qtm_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and minus shift
f == Qsw_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with plus shift
f == Qsw_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with minus shift
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
if( N != VOLUME/2 ) {
if( g_proc_id == 0 )
printf("ERROR: expected N == VOLUME/2 for the required operator in MG_mms_solve_nd.\n");
return 0;
}
} else if( g_proc_id == 0 )
printf("WARNING: required operator unknown for MG_solve. Using standard operator: %s.\n",
N==VOLUME?"":"Qsw_pm_ndpsi_shift");
// Setting mubar, epsbar and shifts
if ( f == Qtm_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and plus shift
f == Qsw_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with plus shift
f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
MG_update_mubar_epsbar( g_mubar, g_epsbar, shifts[0] );
for( int i = 0; i < no_shifts; i++ ) {
mg_odd_shifts[i] = shifts[i]*mg_scale;
mg_even_shifts[i] = 0;
}
}
else if ( f == Qtm_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and minus shift
f == Qsw_tau1_ndpsi_sub_Ishift ) {// Gamma5 Dh tau1 - Schur complement with minus shift
MG_update_mubar_epsbar( g_mubar, g_epsbar, -shifts[0] );
for( int i = 0; i < no_shifts; i++ ) {
mg_odd_shifts[i] = -shifts[i]*mg_scale;
mg_even_shifts[i] = 0;
}
}
//Solving
if ( f == Qtm_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and plus shift
f == Qtm_tau1_ndpsi_sub_Ishift || // Gamma5 Dh tau1 - Schur complement with csw = 0 and minus shift
f == Qsw_tau1_ndpsi_add_Ishift || // Gamma5 Dh tau1 - Schur complement with plus shift
f == Qsw_tau1_ndpsi_sub_Ishift ) {// Gamma5 Dh tau1 - Schur complement with minus shift
mul_gamma5((spinor *const) old1, VOLUME);
mul_gamma5((spinor *const) old2, VOLUME);
// tau1 exchange new1 <-> new2
DDalphaAMG_solve_ms_doublet( new2, old1, new1, old2, mg_even_shifts, mg_odd_shifts, no_shifts,
precision, &mg_status );
if( N == VOLUME ) { // in case of VOLUME/2 old is a just local vector
mul_gamma5((spinor *const) old1, VOLUME);
mul_gamma5((spinor *const) old2, VOLUME);
}
}
else if ( f == Qtm_pm_ndpsi_shift || // (Gamma5 Dh tau1)^2 - Schur complement squared with csw = 0 and shift
f == Qsw_pm_ndpsi_shift ) { // (Gamma5 Dh tau1)^2 - Schur complement squared with shift
mg_scale *= mg_scale;
// DDalphaAMG: tau1 gamma5 Dh tau1 gamma5 Dh
// tmLQCD: gamma5 Dh tau1 gamma5 Dh tau1
DDalphaAMG_solve_ms_doublet_squared_odd( new2, old2, new1, old1, mg_even_shifts, mg_odd_shifts, no_shifts,
precision, &mg_status );
}
else
DDalphaAMG_solve_ms_doublet( new1, old1, new2, old2, mg_even_shifts, mg_odd_shifts, no_shifts,
precision, &mg_status );
if (N==VOLUME/2) {
for( int i = 0; i < no_shifts; i++ ) {
convert_lexic_to_odd(up_new[i], (spinor*) new1[i]);
convert_lexic_to_odd(dn_new[i], (spinor*) new2[i]);
}
finalize_solver(solver_field, 2+2*no_shifts);
}
for( int i = 0; i < no_shifts; i++ ) {
mul_r(up_new[i], mg_scale, up_new[i], N);
mul_r(dn_new[i], mg_scale, dn_new[i], N);
}
if (g_proc_id == 0) {
printf("Solving time %.2f sec (%.1f %% on coarse grid)\n", mg_status.time,
100.*(mg_status.coarse_time/mg_status.time));
printf("Total iterations on fine grid %d\n", mg_status.iter_count);
printf("Total iterations on coarse grids %d\n", mg_status.coarse_iter_count);
if (!mg_status.success)