forked from milc-qcd/milc_qcd
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathlattice.h
156 lines (132 loc) · 5.19 KB
/
lattice.h
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
#ifndef _LATTICE_H
#define _LATTICE_H
/****************************** lattice.h ********************************/
/* include file for MIMD version 7
This file defines global scalars and the fields in the lattice.
Directory for dynamical improved KS action. Allow:
arbitrary paths in quark action
arbitrary paths in gauge action (eg Symanzik imp.)
If "FN" is defined,
Includes storage for Naik improvement (longlink[4], templongvec[4],
gen_pt[16], etc.
*/
#include "defines.h"
#include "../include/generic_quark_types.h"
#include "../include/random.h" /* For double_prn */
#include "../include/macros.h" /* For MAXFILENAME */
#include "../include/io_lat.h" /* For gauge_file */
#ifdef HAVE_QDP
#include <qdp.h>
#endif
/* Begin definition of site structure */
#include "../include/su3.h"
#include "../include/random.h" /* For double_prn */
/* The lattice is an array of sites. */
#define MOM_SITE /* If there is a mom member of the site struct */
typedef struct {
/* The first part is standard to all programs */
/* coordinates of this site */
short x,y,z,t;
/* is it even or odd? */
char parity;
/* my index in the array */
int index;
#ifdef SITERAND
/* The state information for a random number generator */
double_prn site_prn;
/* align to double word boundary (kludge for Intel compiler) */
int space1;
#endif
/* ------------------------------------------------------------ */
/* Now come the physical fields, program dependent */
/* ------------------------------------------------------------ */
/* gauge field */
su3_matrix link[4]; /* the fundamental field */
su3_matrix rgt; /* gauge transformation */
su3_matrix sm_link[4]; /* the smeared field */
su3_matrix stapleg;
su3_vector tempvecg;
/* antihermitian momentum matrices in each direction */
anti_hermitmat mom[4] ALIGNMENT;
/* The Kogut-Susskind phases, which have been absorbed into
the matrices. Also the antiperiodic boundary conditions. */
Real phase[4];
/* 3 element complex vectors */
su3_vector phi; /* Gaussian random source vector */
su3_vector xxx; /* solution vector = Kinverse * phi */
su3_vector resid; /* conjugate gradient residual vector */
su3_vector phi2; /* Gaussian random source vector, mass2 */
su3_vector cg_p; /* conjugate gradient change vector */
su3_vector ttt; /* temporary vector, for K*ppp */
su3_vector g_rand; /* Gaussian random vector*/
/* Use trick of combining xxx=D^adj D)^(-1) on even sites with
Dslash times this on odd sites when computing fermion force */
su3_vector tempvec[4]; /* One for each direction */
su3_matrix s_link,s_link_f,t_link_f,diag,test[4];
#ifdef FN
su3_vector templongvec[4]; /* One for each direction */
su3_vector templongv1;
#endif
} site;
/* End definition of site structure */
/* Definition of globals */
#ifdef CONTROL
#define EXTERN
#else
#define EXTERN extern
#endif
/* The following are global scalars
beta is overall gauge coupling factor
mass, mass1 and mass2 are quark masses
u0 is tadpole improvement factor, perhaps (plaq/3)^(1/4)
*/
EXTERN int nx,ny,nz,nt; /* lattice dimensions */
EXTERN int volume; /* volume of lattice = nx*ny*nz*nt */
EXTERN int iseed; /* random number seed */
EXTERN Real beta;
EXTERN int niter;
EXTERN int nflavors;
EXTERN Real u0;
EXTERN Real mass;
EXTERN Real rsqprop;
EXTERN int startflag; /* beginning lattice: CONTINUE, RELOAD, RELOAD_BINARY,
RELOAD_CHECKPOINT, FRESH */
EXTERN int saveflag; /* do with lattice: FORGET, SAVE, SAVE_BINARY,
SAVE_CHECKPOINT */
EXTERN char startfile[MAXFILENAME],savefile[MAXFILENAME];
EXTERN double g_ssplaq, g_stplaq;
EXTERN double_complex linktrsum;
EXTERN u_int32type nersc_checksum;
EXTERN char stringLFN[MAXFILENAME]; /** ILDG LFN if applicable **/
EXTERN char propfile[MAXFILENAME]; /* For blocked propagator */
EXTERN int total_iters;
EXTERN int phases_in; /* 1 if KS and BC phases absorbed into matrices */
EXTERN int source_start, source_inc, n_sources;
/* source time, increment for it, and number of source slices */
EXTERN int nrg; /* Working number of RG steps */
/* Some of these global variables are node dependent */
/* They are set in "make_lattice()" */
EXTERN size_t sites_on_node; /* number of sites on this node */
EXTERN size_t even_sites_on_node; /* number of even sites on this node */
EXTERN size_t odd_sites_on_node; /* number of odd sites on this node */
EXTERN int number_of_nodes; /* number of nodes in use */
EXTERN int this_node; /* node number of this node */
EXTERN gauge_file *startlat_p;
/* Each node maintains a structure with the pseudorandom number
generator state */
EXTERN double_prn node_prn ;
/* The lattice is a single global variable - (actually this is the
part of the lattice on this node) */
EXTERN Real boundary_phase[4];
EXTERN site *lattice;
/* Vectors for addressing */
/* Generic pointers, for gather routines */
#define N_POINTERS 16
EXTERN char ** gen_pt[N_POINTERS];
/* Storage for definition of the quark action */
EXTERN ferm_links_t fn_links;
EXTERN ferm_links_t fn_links_dmdu0;
EXTERN ks_action_paths ks_act_paths;
EXTERN ks_action_paths ks_act_paths_dmdu0;
#define smear_fac 2.5
#endif /* _LATTICE_H */