forked from etmc/tmLQCD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
quda_interface.c
831 lines (696 loc) · 28.6 KB
/
quda_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
/***********************************************************************
*
* Copyright (C) 2015 Mario Schroeck
*
* 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/>.
*
*
***********************************************************************/
/***********************************************************************
*
* File quda_interface.h
*
* Author: Mario Schroeck <mario.schroeck@roma3.infn.it>
*
* Last changes: 06/2015
*
*
* Interface to QUDA for multi-GPU inverters
*
* The externally accessible functions are
*
* void _initQuda()
* Initializes the QUDA library. Carries over the lattice size and the
* MPI process grid and thus must be called after initializing MPI.
* Currently it is called in init_operators() if optr->use_qudainverter
* flag is set.
* Memory for the QUDA gaugefield on the host is allocated but not filled
* yet (the latter is done in _loadGaugeQuda(), see below).
* Performance critical settings are done here and can be changed.
*
* void _endQuda()
* Finalizes the QUDA library. Call before MPI_Finalize().
*
* void _loadGaugeQuda()
* Copies and reorders the gaugefield on the host and copies it to the GPU.
* Must be called between last changes on the gaugefield (smearing etc.)
* and first call of the inverter. In particular, 'boundary(const double kappa)'
* must be called before if nontrivial boundary conditions are to be used since
* those will be applied directly to the gaugefield. Currently it is called just
* before the inversion is done (might result in wasted loads...).
*
* The functions
*
* int invert_eo_quda(...);
* int invert_doublet_eo_quda(...);
* void M_full_quda(...);
* void D_psi_quda(...);
*
* mimic their tmLQCD counterparts in functionality as well as input and
* output parameters. The invert functions will check the parameters
* g_mu, g_c_sw do decide which QUDA operator to create.
*
* To activate those, set "UseQudaInverter = yes" in the operator
* declaration of the input file. For details see the documentation.
*
*
* Notes:
*
* Minimum QUDA version is 0.7.0 (see https://github.com/lattice/quda/issues/151
* and https://github.com/lattice/quda/issues/157).
*
*
**************************************************************************/
#include "quda_interface.h"
#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 "boundary.h"
#include "quda.h"
double X0, X1, X2, X3;
// define order of the spatial indices
// default is LX-LY-LZ-T, see below def. of local lattice size, this is related to
// the gamma basis transformation from tmLQCD -> UKQCD
// for details see https://github.com/lattice/quda/issues/157
#define USE_LZ_LY_LX_T 0
// TRIVIAL_BC are trivial (anti-)periodic boundary conditions,
// i.e. 1 or -1 on last timeslice
// tmLQCD uses twisted BC, i.e. phases on all timeslices.
// if using TRIVIAL_BC: can't compare inversion result to tmLQCD
// if not using TRIVIAL_BC: BC will be applied to gauge field,
// can't use 12 parameter reconstruction
#define TRIVIAL_BC 0
#define MAX(a,b) ((a)>(b)?(a):(b))
// gauge and invert paramameter structs; init. in _initQuda()
QudaGaugeParam gauge_param;
QudaInvertParam inv_param;
// pointer to the QUDA gaugefield
double *gauge_quda[4];
// pointer to a temp. spinor, used for reordering etc.
double *tempSpinor;
// function that maps coordinates in the communication grid to MPI ranks
int commsMap(const int *coords, void *fdata) {
#if USE_LZ_LY_LX_T
int n[4] = {coords[3], coords[2], coords[1], coords[0]};
#else
int n[4] = {coords[3], coords[0], coords[1], coords[2]};
#endif
int rank = 0;
#ifdef MPI
MPI_Cart_rank( g_cart_grid, n, &rank );
#endif
return rank;
}
// variable to check if quda has been initialized
static int quda_initialized = 0;
void _initQuda() {
if( quda_initialized )
return;
if( g_debug_level > 0 )
if(g_proc_id == 0)
printf("\n# QUDA: Detected QUDA version %d.%d.%d\n\n", QUDA_VERSION_MAJOR, QUDA_VERSION_MINOR, QUDA_VERSION_SUBMINOR);
if( QUDA_VERSION_MAJOR == 0 && QUDA_VERSION_MINOR < 7) {
fprintf(stderr, "Error: minimum QUDA version required is 0.7.0 (for support of chiral basis and removal of bug in mass normalization with preconditioning).\n");
exit(-2);
}
gauge_param = newQudaGaugeParam();
inv_param = newQudaInvertParam();
// *** QUDA parameters begin here (sloppy prec. will be adjusted in invert)
QudaPrecision cpu_prec = QUDA_DOUBLE_PRECISION;
QudaPrecision cuda_prec = QUDA_DOUBLE_PRECISION;
QudaPrecision cuda_prec_sloppy = QUDA_SINGLE_PRECISION;
QudaPrecision cuda_prec_precondition = QUDA_HALF_PRECISION;
QudaTune tune = QUDA_TUNE_YES;
// *** the remainder should not be changed for this application
// local lattice size
#if USE_LZ_LY_LX_T
gauge_param.X[0] = LZ;
gauge_param.X[1] = LY;
gauge_param.X[2] = LX;
gauge_param.X[3] = T;
#else
gauge_param.X[0] = LX;
gauge_param.X[1] = LY;
gauge_param.X[2] = LZ;
gauge_param.X[3] = T;
#endif
inv_param.Ls = 1;
gauge_param.anisotropy = 1.0;
gauge_param.type = QUDA_WILSON_LINKS;
gauge_param.gauge_order = QUDA_QDP_GAUGE_ORDER;
gauge_param.cpu_prec = cpu_prec;
gauge_param.cuda_prec = cuda_prec;
gauge_param.reconstruct = 18;
gauge_param.cuda_prec_sloppy = cuda_prec_sloppy;
gauge_param.reconstruct_sloppy = 18;
gauge_param.cuda_prec_precondition = cuda_prec_precondition;
gauge_param.reconstruct_precondition = 18;
gauge_param.gauge_fix = QUDA_GAUGE_FIXED_NO;
inv_param.dagger = QUDA_DAG_NO;
inv_param.mass_normalization = QUDA_KAPPA_NORMALIZATION;
inv_param.solver_normalization = QUDA_DEFAULT_NORMALIZATION;
inv_param.pipeline = 0;
inv_param.gcrNkrylov = 10;
// require both L2 relative and heavy quark residual to determine convergence
// inv_param.residual_type = (QudaResidualType)(QUDA_L2_RELATIVE_RESIDUAL | QUDA_HEAVY_QUARK_RESIDUAL);
inv_param.tol_hq = 1.0;//1e-3; // specify a tolerance for the residual for heavy quark residual
inv_param.reliable_delta = 1e-2; // ignored by multi-shift solver
// domain decomposition preconditioner parameters
inv_param.inv_type_precondition = QUDA_CG_INVERTER;
inv_param.schwarz_type = QUDA_ADDITIVE_SCHWARZ;
inv_param.precondition_cycle = 1;
inv_param.tol_precondition = 1e-1;
inv_param.maxiter_precondition = 10;
inv_param.verbosity_precondition = QUDA_SILENT;
inv_param.cuda_prec_precondition = cuda_prec_precondition;
inv_param.omega = 1.0;
inv_param.cpu_prec = cpu_prec;
inv_param.cuda_prec = cuda_prec;
inv_param.cuda_prec_sloppy = cuda_prec_sloppy;
inv_param.clover_cpu_prec = cpu_prec;
inv_param.clover_cuda_prec = cuda_prec;
inv_param.clover_cuda_prec_sloppy = cuda_prec_sloppy;
inv_param.clover_cuda_prec_precondition = cuda_prec_precondition;
inv_param.preserve_source = QUDA_PRESERVE_SOURCE_YES;
inv_param.gamma_basis = QUDA_CHIRAL_GAMMA_BASIS;
inv_param.dirac_order = QUDA_DIRAC_ORDER;
inv_param.input_location = QUDA_CPU_FIELD_LOCATION;
inv_param.output_location = QUDA_CPU_FIELD_LOCATION;
inv_param.tune = tune ? QUDA_TUNE_YES : QUDA_TUNE_NO;
gauge_param.ga_pad = 0; // 24*24*24/2;
inv_param.sp_pad = 0; // 24*24*24/2;
inv_param.cl_pad = 0; // 24*24*24/2;
// For multi-GPU, ga_pad must be large enough to store a time-slice
int x_face_size = gauge_param.X[1]*gauge_param.X[2]*gauge_param.X[3]/2;
int y_face_size = gauge_param.X[0]*gauge_param.X[2]*gauge_param.X[3]/2;
int z_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[3]/2;
int t_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[2]/2;
int pad_size =MAX(x_face_size, y_face_size);
pad_size = MAX(pad_size, z_face_size);
pad_size = MAX(pad_size, t_face_size);
gauge_param.ga_pad = pad_size;
// solver verbosity
if( g_debug_level == 0 )
inv_param.verbosity = QUDA_SILENT;
else if( g_debug_level == 1 )
inv_param.verbosity = QUDA_SUMMARIZE;
else
inv_param.verbosity = QUDA_VERBOSE;
// general verbosity
setVerbosityQuda( QUDA_SUMMARIZE, "# QUDA: ", stdout);
// declare the grid mapping used for communications in a multi-GPU grid
#if USE_LZ_LY_LX_T
int grid[4] = {g_nproc_z, g_nproc_y, g_nproc_x, g_nproc_t};
#else
int grid[4] = {g_nproc_x, g_nproc_y, g_nproc_z, g_nproc_t};
#endif
initCommsGridQuda(4, grid, commsMap, NULL);
// alloc gauge_quda
size_t gSize = (gauge_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
for (int dir = 0; dir < 4; dir++) {
gauge_quda[dir] = (double*) malloc(VOLUME*18*gSize);
if(gauge_quda[dir] == NULL) {
fprintf(stderr, "_initQuda: malloc for gauge_quda[dir] failed");
exit(-2);
}
}
// alloc space for a temp. spinor, used throughout this module
tempSpinor = (double*)malloc( 2*VOLUME*24*sizeof(double) ); /* factor 2 for doublet */
if(tempSpinor == NULL) {
fprintf(stderr, "_initQuda: malloc for tempSpinor failed");
exit(-2);
}
// initialize the QUDA library
#ifdef MPI
initQuda(-1); //sets device numbers automatically
#else
initQuda(0); //scalar build: use device 0
#endif
quda_initialized = 1;
}
// finalize the QUDA library
void _endQuda() {
if( quda_initialized ) {
freeGaugeQuda();
free((void*)tempSpinor);
endQuda();
}
}
void _loadGaugeQuda( const int compression ) {
if( inv_param.verbosity > QUDA_SILENT )
if(g_proc_id == 0)
printf("# QUDA: Called _loadGaugeQuda\n");
_Complex double tmpcplx;
size_t gSize = (gauge_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
// now copy and reorder
for( int x0=0; x0<T; x0++ )
for( int x1=0; x1<LX; x1++ )
for( int x2=0; x2<LY; x2++ )
for( int x3=0; x3<LZ; x3++ ) {
#if USE_LZ_LY_LX_T
int j = x3 + LZ*x2 + LY*LZ*x1 + LX*LY*LZ*x0;
int tm_idx = x1 + LX*x2 + LY*LX*x3 + LZ*LY*LX*x0;
#else
int j = x1 + LX*x2 + LY*LX*x3 + LZ*LY*LX*x0;
int tm_idx = x3 + LZ*x2 + LY*LZ*x1 + LX*LY*LZ*x0;
#endif
int oddBit = (x0+x1+x2+x3) & 1;
int quda_idx = 18*(oddBit*VOLUME/2+j/2);
#if USE_LZ_LY_LX_T
memcpy( &(gauge_quda[0][quda_idx]), &(g_gauge_field[tm_idx][3]), 18*gSize);
memcpy( &(gauge_quda[1][quda_idx]), &(g_gauge_field[tm_idx][2]), 18*gSize);
memcpy( &(gauge_quda[2][quda_idx]), &(g_gauge_field[tm_idx][1]), 18*gSize);
memcpy( &(gauge_quda[3][quda_idx]), &(g_gauge_field[tm_idx][0]), 18*gSize);
#else
memcpy( &(gauge_quda[0][quda_idx]), &(g_gauge_field[tm_idx][1]), 18*gSize);
memcpy( &(gauge_quda[1][quda_idx]), &(g_gauge_field[tm_idx][2]), 18*gSize);
memcpy( &(gauge_quda[2][quda_idx]), &(g_gauge_field[tm_idx][3]), 18*gSize);
memcpy( &(gauge_quda[3][quda_idx]), &(g_gauge_field[tm_idx][0]), 18*gSize);
if( compression == 18 ) {
// apply boundary conditions
for( int i=0; i<9; i++ ) {
tmpcplx = gauge_quda[0][quda_idx+2*i] + I*gauge_quda[0][quda_idx+2*i+1];
tmpcplx *= -phase_1/g_kappa;
gauge_quda[0][quda_idx+2*i] = creal(tmpcplx);
gauge_quda[0][quda_idx+2*i+1] = cimag(tmpcplx);
tmpcplx = gauge_quda[1][quda_idx+2*i] + I*gauge_quda[1][quda_idx+2*i+1];
tmpcplx *= -phase_2/g_kappa;
gauge_quda[1][quda_idx+2*i] = creal(tmpcplx);
gauge_quda[1][quda_idx+2*i+1] = cimag(tmpcplx);
tmpcplx = gauge_quda[2][quda_idx+2*i] + I*gauge_quda[2][quda_idx+2*i+1];
tmpcplx *= -phase_3/g_kappa;
gauge_quda[2][quda_idx+2*i] = creal(tmpcplx);
gauge_quda[2][quda_idx+2*i+1] = cimag(tmpcplx);
tmpcplx = gauge_quda[3][quda_idx+2*i] + I*gauge_quda[3][quda_idx+2*i+1];
tmpcplx *= -phase_0/g_kappa;
gauge_quda[3][quda_idx+2*i] = creal(tmpcplx);
gauge_quda[3][quda_idx+2*i+1] = cimag(tmpcplx);
}
}
#endif
}
loadGaugeQuda((void*)gauge_quda, &gauge_param);
}
// reorder spinor to QUDA format
void reorder_spinor_toQuda( double* spinor, QudaPrecision precision, int doublet, double* spinor2 ) {
double startTime = gettime();
if( doublet ) {
memcpy( tempSpinor, spinor, VOLUME*24*sizeof(double) );
memcpy( tempSpinor+VOLUME*24, spinor2, VOLUME*24*sizeof(double) );
}
else {
memcpy( tempSpinor, spinor, VOLUME*24*sizeof(double) );
}
// now copy and reorder from tempSpinor to spinor
for( int x0=0; x0<T; x0++ )
for( int x1=0; x1<LX; x1++ )
for( int x2=0; x2<LY; x2++ )
for( int x3=0; x3<LZ; x3++ ) {
#if USE_LZ_LY_LX_T
int j = x3 + LZ*x2 + LY*LZ*x1 + LX*LY*LZ*x0;
int tm_idx = x1 + LX*x2 + LY*LX*x3 + LZ*LY*LX*x0;
#else
int j = x1 + LX*x2 + LY*LX*x3 + LZ*LY*LX*x0;
int tm_idx = x3 + LZ*x2 + LY*LZ*x1 + LX*LY*LZ*x0;
#endif
int oddBit = (x0+x1+x2+x3) & 1;
if( doublet ) {
memcpy( &(spinor[24*(oddBit*VOLUME+j/2)]), &(tempSpinor[24* tm_idx ]), 24*sizeof(double));
memcpy( &(spinor[24*(oddBit*VOLUME+j/2+VOLUME/2)]), &(tempSpinor[24*(tm_idx+VOLUME)]), 24*sizeof(double));
}
else {
memcpy( &(spinor[24*(oddBit*VOLUME/2+j/2)]), &(tempSpinor[24*tm_idx]), 24*sizeof(double));
}
}
double endTime = gettime();
double diffTime = endTime - startTime;
if(g_proc_id == 0)
printf("# QUDA: time spent in reorder_spinor_toQuda: %f secs\n", diffTime);
}
// reorder spinor from QUDA format
void reorder_spinor_fromQuda( double* spinor, QudaPrecision precision, int doublet, double* spinor2 ) {
double startTime = gettime();
if( doublet ) {
memcpy( tempSpinor, spinor, 2*VOLUME*24*sizeof(double) );
}
else {
memcpy( tempSpinor, spinor, VOLUME*24*sizeof(double) );
}
// now copy and reorder from tempSpinor to spinor
for( int x0=0; x0<T; x0++ )
for( int x1=0; x1<LX; x1++ )
for( int x2=0; x2<LY; x2++ )
for( int x3=0; x3<LZ; x3++ ) {
#if USE_LZ_LY_LX_T
int j = x3 + LZ*x2 + LY*LZ*x1 + LX*LY*LZ*x0;
int tm_idx = x1 + LX*x2 + LY*LX*x3 + LZ*LY*LX*x0;
#else
int j = x1 + LX*x2 + LY*LX*x3 + LZ*LY*LX*x0;
int tm_idx = x3 + LZ*x2 + LY*LZ*x1 + LX*LY*LZ*x0;
#endif
int oddBit = (x0+x1+x2+x3) & 1;
if( doublet ) {
memcpy( &(spinor [24* tm_idx]), &(tempSpinor[24*(oddBit*VOLUME+j/2) ]), 24*sizeof(double));
memcpy( &(spinor2[24*(tm_idx)]), &(tempSpinor[24*(oddBit*VOLUME+j/2+VOLUME/2)]), 24*sizeof(double));
}
else {
memcpy( &(spinor[24*tm_idx]), &(tempSpinor[24*(oddBit*VOLUME/2+j/2)]), 24*sizeof(double));
}
}
double endTime = gettime();
double diffTime = endTime - startTime;
if(g_proc_id == 0)
printf("# QUDA: time spent in reorder_spinor_fromQuda: %f secs\n", diffTime);
}
void set_boundary_conditions( CompressionType* compression ) {
// we can't have compression and theta-BC
if( fabs(X1)>0.0 || fabs(X2)>0.0 || fabs(X3)>0.0 || (fabs(X0)!=0.0 && fabs(X0)!=1.0) ) {
if( *compression!=NO_COMPRESSION ) {
if(g_proc_id == 0) {
printf("\n# QUDA: WARNING you can't use compression %d with boundary conditions for fermion fields (t,x,y,z)*pi: (%f,%f,%f,%f) \n", *compression,X0,X1,X2,X3);
printf("# QUDA: disabling compression.\n\n");
*compression=NO_COMPRESSION;
}
}
}
QudaReconstructType link_recon;
QudaReconstructType link_recon_sloppy;
if( *compression==NO_COMPRESSION ) { // theta BC
gauge_param.t_boundary = QUDA_PERIODIC_T; // BC will be applied to gaugefield
link_recon = 18;
link_recon_sloppy = 18;
}
else { // trivial BC
gauge_param.t_boundary = ( fabs(X0)>0.0 ? QUDA_ANTI_PERIODIC_T : QUDA_PERIODIC_T );
link_recon = 12;
link_recon_sloppy = *compression;
if( g_debug_level > 0 )
if(g_proc_id == 0)
printf("\n# QUDA: WARNING using %d compression with trivial (A)PBC instead of theta-BC ((t,x,y,z)*pi: (%f,%f,%f,%f))! This works fine but the residual check on the host (CPU) will fail.\n",*compression,X0,X1,X2,X3);
}
gauge_param.reconstruct = link_recon;
gauge_param.reconstruct_sloppy = link_recon_sloppy;
gauge_param.reconstruct_precondition = link_recon_sloppy;
}
void set_sloppy_prec( const SloppyPrecision sloppy_precision ) {
// choose sloppy prec.
QudaPrecision cuda_prec_sloppy;
if( sloppy_precision==SLOPPY_DOUBLE ) {
cuda_prec_sloppy = QUDA_DOUBLE_PRECISION;
if(g_proc_id == 0) printf("# QUDA: Using double prec. as sloppy!\n");
}
else if( sloppy_precision==SLOPPY_HALF ) {
cuda_prec_sloppy = QUDA_HALF_PRECISION;
if(g_proc_id == 0) printf("# QUDA: Using half prec. as sloppy!\n");
}
else {
cuda_prec_sloppy = QUDA_SINGLE_PRECISION;
if(g_proc_id == 0) printf("# QUDA: Using single prec. as sloppy!\n");
}
gauge_param.cuda_prec_sloppy = cuda_prec_sloppy;
inv_param.cuda_prec_sloppy = cuda_prec_sloppy;
inv_param.clover_cuda_prec_sloppy = cuda_prec_sloppy;
}
int invert_eo_quda(spinor * const Even_new, spinor * const Odd_new,
spinor * const Even, spinor * const Odd,
const double precision, const int max_iter,
const int solver_flag, const int rel_prec,
const int even_odd_flag, solver_params_t solver_params,
SloppyPrecision sloppy_precision,
CompressionType compression) {
spinor ** solver_field = NULL;
const int nr_sf = 2;
init_solver_field(&solver_field, VOLUMEPLUSRAND, nr_sf);
convert_eo_to_lexic(solver_field[0], Even, Odd);
// convert_eo_to_lexic(solver_field[1], Even_new, Odd_new);
void *spinorIn = (void*)solver_field[0]; // source
void *spinorOut = (void*)solver_field[1]; // solution
if ( rel_prec )
inv_param.residual_type = QUDA_L2_RELATIVE_RESIDUAL;
else
inv_param.residual_type = QUDA_L2_ABSOLUTE_RESIDUAL;
inv_param.kappa = g_kappa;
// figure out which BC to use (theta, trivial...)
set_boundary_conditions(&compression);
// set the sloppy precision of the mixed prec solver
set_sloppy_prec(sloppy_precision);
// load gauge after setting precision
_loadGaugeQuda(compression);
// choose dslash type
if( g_mu != 0.0 && g_c_sw > 0.0 ) {
inv_param.dslash_type = QUDA_TWISTED_CLOVER_DSLASH;
inv_param.matpc_type = QUDA_MATPC_EVEN_EVEN;
inv_param.solution_type = QUDA_MAT_SOLUTION;
inv_param.clover_order = QUDA_PACKED_CLOVER_ORDER;
inv_param.mu = fabs(g_mu/2./g_kappa);
inv_param.clover_coeff = g_c_sw*g_kappa;
}
else if( g_mu != 0.0 ) {
inv_param.dslash_type = QUDA_TWISTED_MASS_DSLASH;
inv_param.matpc_type = QUDA_MATPC_EVEN_EVEN_ASYMMETRIC;
inv_param.solution_type = QUDA_MAT_SOLUTION;
inv_param.mu = fabs(g_mu/2./g_kappa);
}
else if( g_c_sw > 0.0 ) {
inv_param.dslash_type = QUDA_CLOVER_WILSON_DSLASH;
inv_param.matpc_type = QUDA_MATPC_EVEN_EVEN;
inv_param.solution_type = QUDA_MAT_SOLUTION;
inv_param.clover_order = QUDA_PACKED_CLOVER_ORDER;
inv_param.clover_coeff = g_c_sw*g_kappa;
}
else {
inv_param.dslash_type = QUDA_WILSON_DSLASH;
inv_param.matpc_type = QUDA_MATPC_EVEN_EVEN;
inv_param.solution_type = QUDA_MAT_SOLUTION;
}
// choose solver
if(solver_flag == BICGSTAB) {
if(g_proc_id == 0) {printf("# QUDA: Using BiCGstab!\n"); fflush(stdout);}
inv_param.inv_type = QUDA_BICGSTAB_INVERTER;
}
else {
/* Here we invert the hermitean operator squared */
inv_param.inv_type = QUDA_CG_INVERTER;
if(g_proc_id == 0) {
printf("# QUDA: Using mixed precision CG!\n");
printf("# QUDA: mu = %f, kappa = %f\n", g_mu/2./g_kappa, g_kappa);
fflush(stdout);
}
}
// direct or norm-op. solve
if( inv_param.inv_type == QUDA_CG_INVERTER ) {
if( even_odd_flag ) {
inv_param.solve_type = QUDA_NORMOP_PC_SOLVE;
if(g_proc_id == 0) printf("# QUDA: Using preconditioning!\n");
}
else {
inv_param.solve_type = QUDA_NORMOP_SOLVE;
if(g_proc_id == 0) printf("# QUDA: Not using preconditioning!\n");
}
}
else {
if( even_odd_flag ) {
inv_param.solve_type = QUDA_DIRECT_PC_SOLVE;
if(g_proc_id == 0) printf("# QUDA: Using preconditioning!\n");
}
else {
inv_param.solve_type = QUDA_DIRECT_SOLVE;
if(g_proc_id == 0) printf("# QUDA: Not using preconditioning!\n");
}
}
inv_param.tol = sqrt(precision)*0.25;
inv_param.maxiter = max_iter;
// IMPORTANT: use opposite TM flavor since gamma5 -> -gamma5 (until LXLYLZT prob. resolved)
inv_param.twist_flavor = (g_mu < 0.0 ? QUDA_TWIST_PLUS : QUDA_TWIST_MINUS);
inv_param.Ls = 1;
// NULL pointers to host fields to force
// construction instead of download of the clover field:
if( g_c_sw > 0.0 )
loadCloverQuda(NULL, NULL, &inv_param);
// reorder spinor
reorder_spinor_toQuda( (double*)spinorIn, inv_param.cpu_prec, 0, NULL );
// perform the inversion
invertQuda(spinorOut, spinorIn, &inv_param);
if( inv_param.verbosity == QUDA_VERBOSE )
if(g_proc_id == 0)
printf("# QUDA: Device memory used: Spinor: %f GiB, Gauge: %f GiB, Clover: %f GiB\n",
inv_param.spinorGiB, gauge_param.gaugeGiB, inv_param.cloverGiB);
if( inv_param.verbosity > QUDA_SILENT )
if(g_proc_id == 0)
printf("# QUDA: Done: %i iter / %g secs = %g Gflops\n",
inv_param.iter, inv_param.secs, inv_param.gflops/inv_param.secs);
// number of CG iterations
int iteration = inv_param.iter;
// reorder spinor
reorder_spinor_fromQuda( (double*)spinorIn, inv_param.cpu_prec, 0, NULL );
reorder_spinor_fromQuda( (double*)spinorOut, inv_param.cpu_prec, 0, NULL );
convert_lexic_to_eo(Even, Odd, solver_field[0]);
convert_lexic_to_eo(Even_new, Odd_new, solver_field[1]);
finalize_solver(solver_field, nr_sf);
freeGaugeQuda();
if(iteration >= max_iter)
return(-1);
return(iteration);
}
int invert_doublet_eo_quda(spinor * const Even_new_s, spinor * const Odd_new_s,
spinor * const Even_new_c, spinor * const Odd_new_c,
spinor * const Even_s, spinor * const Odd_s,
spinor * const Even_c, spinor * const Odd_c,
const double precision, const int max_iter,
const int solver_flag, const int rel_prec, const int even_odd_flag,
const SloppyPrecision sloppy_precision,
CompressionType compression) {
spinor ** solver_field = NULL;
const int nr_sf = 4;
init_solver_field(&solver_field, VOLUMEPLUSRAND, nr_sf);
convert_eo_to_lexic(solver_field[0], Even_s, Odd_s);
convert_eo_to_lexic(solver_field[1], Even_c, Odd_c);
// convert_eo_to_lexic(g_spinor_field[DUM_DERI+1], Even_new, Odd_new);
void *spinorIn = (void*)solver_field[0]; // source
void *spinorIn_c = (void*)solver_field[1]; // charme source
void *spinorOut = (void*)solver_field[2]; // solution
void *spinorOut_c = (void*)solver_field[3]; // charme solution
if ( rel_prec )
inv_param.residual_type = QUDA_L2_RELATIVE_RESIDUAL;
else
inv_param.residual_type = QUDA_L2_ABSOLUTE_RESIDUAL;
inv_param.kappa = g_kappa;
// IMPORTANT: use opposite TM mu-flavor since gamma5 -> -gamma5
inv_param.mu = -g_mubar /2./g_kappa;
inv_param.epsilon = g_epsbar/2./g_kappa;
// figure out which BC to use (theta, trivial...)
set_boundary_conditions(&compression);
// set the sloppy precision of the mixed prec solver
set_sloppy_prec(sloppy_precision);
// load gauge after setting precision
_loadGaugeQuda(compression);
// choose dslash type
if( g_c_sw > 0.0 ) {
inv_param.dslash_type = QUDA_TWISTED_CLOVER_DSLASH;
inv_param.matpc_type = QUDA_MATPC_EVEN_EVEN;
inv_param.solution_type = QUDA_MAT_SOLUTION;
inv_param.clover_order = QUDA_PACKED_CLOVER_ORDER;
inv_param.clover_coeff = g_c_sw*g_kappa;
}
else {
inv_param.dslash_type = QUDA_TWISTED_MASS_DSLASH;
inv_param.matpc_type = QUDA_MATPC_EVEN_EVEN_ASYMMETRIC;
inv_param.solution_type = QUDA_MAT_SOLUTION;
}
// choose solver
if(solver_flag == BICGSTAB) {
if(g_proc_id == 0) {printf("# QUDA: Using BiCGstab!\n"); fflush(stdout);}
inv_param.inv_type = QUDA_BICGSTAB_INVERTER;
}
else {
/* Here we invert the hermitean operator squared */
inv_param.inv_type = QUDA_CG_INVERTER;
if(g_proc_id == 0) {
printf("# QUDA: Using mixed precision CG!\n");
printf("# QUDA: mu = %f, kappa = %f\n", g_mu/2./g_kappa, g_kappa);
fflush(stdout);
}
}
if( even_odd_flag ) {
inv_param.solve_type = QUDA_NORMOP_PC_SOLVE;
if(g_proc_id == 0) printf("# QUDA: Using preconditioning!\n");
}
else {
inv_param.solve_type = QUDA_NORMOP_SOLVE;
if(g_proc_id == 0) printf("# QUDA: Not using preconditioning!\n");
}
inv_param.tol = sqrt(precision)*0.25;
inv_param.maxiter = max_iter;
inv_param.twist_flavor = QUDA_TWIST_NONDEG_DOUBLET;
inv_param.Ls = 2;
// NULL pointers to host fields to force
// construction instead of download of the clover field:
if( g_c_sw > 0.0 )
loadCloverQuda(NULL, NULL, &inv_param);
// reorder spinor
reorder_spinor_toQuda( (double*)spinorIn, inv_param.cpu_prec, 1, (double*)spinorIn_c );
// perform the inversion
invertQuda(spinorOut, spinorIn, &inv_param);
if( inv_param.verbosity == QUDA_VERBOSE )
if(g_proc_id == 0)
printf("# QUDA: Device memory used: Spinor: %f GiB, Gauge: %f GiB, Clover: %f GiB\n",
inv_param.spinorGiB, gauge_param.gaugeGiB, inv_param.cloverGiB);
if( inv_param.verbosity > QUDA_SILENT )
if(g_proc_id == 0)
printf("# QUDA: Done: %i iter / %g secs = %g Gflops\n",
inv_param.iter, inv_param.secs, inv_param.gflops/inv_param.secs);
// number of CG iterations
int iteration = inv_param.iter;
// reorder spinor
reorder_spinor_fromQuda( (double*)spinorIn, inv_param.cpu_prec, 1, (double*)spinorIn_c );
reorder_spinor_fromQuda( (double*)spinorOut, inv_param.cpu_prec, 1, (double*)spinorOut_c );
convert_lexic_to_eo(Even_s, Odd_s, solver_field[0]);
convert_lexic_to_eo(Even_c, Odd_c, solver_field[1]);
convert_lexic_to_eo(Even_new_s, Odd_new_s, solver_field[2]);
convert_lexic_to_eo(Even_new_c, Odd_new_c, solver_field[3]);
finalize_solver(solver_field, nr_sf);
freeGaugeQuda();
if(iteration >= max_iter)
return(-1);
return(iteration);
}
// if even_odd_flag set
void M_full_quda(spinor * const Even_new, spinor * const Odd_new, spinor * const Even, spinor * const Odd) {
inv_param.kappa = g_kappa;
inv_param.mu = fabs(g_mu);
inv_param.epsilon = 0.0;
// IMPORTANT: use opposite TM flavor since gamma5 -> -gamma5 (until LXLYLZT prob. resolved)
inv_param.twist_flavor = (g_mu < 0.0 ? QUDA_TWIST_PLUS : QUDA_TWIST_MINUS);
inv_param.Ls = (inv_param.twist_flavor == QUDA_TWIST_NONDEG_DOUBLET ||
inv_param.twist_flavor == QUDA_TWIST_DEG_DOUBLET ) ? 2 : 1;
void *spinorIn = (void*)g_spinor_field[DUM_DERI]; // source
void *spinorOut = (void*)g_spinor_field[DUM_DERI+1]; // solution
// reorder spinor
convert_eo_to_lexic( spinorIn, Even, Odd );
reorder_spinor_toQuda( (double*)spinorIn, inv_param.cpu_prec, 0, NULL );
// multiply
inv_param.solution_type = QUDA_MAT_SOLUTION;
MatQuda( spinorOut, spinorIn, &inv_param);
// reorder spinor
reorder_spinor_fromQuda( (double*)spinorOut, inv_param.cpu_prec, 0, NULL );
convert_lexic_to_eo( Even_new, Odd_new, spinorOut );
}
// no even-odd
void D_psi_quda(spinor * const P, spinor * const Q) {
inv_param.kappa = g_kappa;
inv_param.mu = fabs(g_mu);
inv_param.epsilon = 0.0;
// IMPORTANT: use opposite TM flavor since gamma5 -> -gamma5 (until LXLYLZT prob. resolved)
inv_param.twist_flavor = (g_mu < 0.0 ? QUDA_TWIST_PLUS : QUDA_TWIST_MINUS);
inv_param.Ls = (inv_param.twist_flavor == QUDA_TWIST_NONDEG_DOUBLET ||
inv_param.twist_flavor == QUDA_TWIST_DEG_DOUBLET ) ? 2 : 1;
void *spinorIn = (void*)Q;
void *spinorOut = (void*)P;
// reorder spinor
reorder_spinor_toQuda( (double*)spinorIn, inv_param.cpu_prec, 0, NULL );
// multiply
inv_param.solution_type = QUDA_MAT_SOLUTION;
MatQuda( spinorOut, spinorIn, &inv_param);
// reorder spinor
reorder_spinor_fromQuda( (double*)spinorIn, inv_param.cpu_prec, 0, NULL );
reorder_spinor_fromQuda( (double*)spinorOut, inv_param.cpu_prec, 0, NULL );
}