-
Notifications
You must be signed in to change notification settings - Fork 102
/
Copy pathgauge_stout.cuh
199 lines (173 loc) · 6.8 KB
/
gauge_stout.cuh
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
#pragma once
#include <gauge_field_order.h>
#include <index_helper.cuh>
#include <quda_matrix.h>
#include <su3_project.cuh>
#include <kernels/gauge_utils.cuh>
#include <kernel.h>
#include <thread_local_cache.h>
namespace quda
{
template <typename Float_, int nColor_, QudaReconstructType recon_, int stoutDim_> struct STOUTArg : kernel_param<> {
using Float = Float_;
static constexpr int nColor = nColor_;
static_assert(nColor == 3, "Only nColor=3 enabled at this time");
static constexpr QudaReconstructType recon = recon_;
static constexpr int stoutDim = stoutDim_;
typedef typename gauge_mapper<Float, recon>::type Gauge;
Gauge out;
const Gauge in;
int X[4]; // grid dimensions
int border[4];
const Float rho;
const Float staple_coeff;
const Float rectangle_coeff;
const int dir_ignore;
STOUTArg(GaugeField &out, const GaugeField &in, Float rho, Float epsilon, int dir_ignore) :
kernel_param(dim3(1, 2, stoutDim)),
out(out),
in(in),
rho(rho),
staple_coeff(rho * (5.0 - 2.0 * epsilon) / 3.0),
rectangle_coeff(rho * (1.0 - epsilon) / 12.0),
dir_ignore(dir_ignore)
{
for (int dir = 0; dir < 4; ++dir) {
border[dir] = in.R()[dir];
X[dir] = in.X()[dir] - border[dir] * 2;
this->threads.x *= X[dir];
}
this->threads.x /= 2;
}
};
template <typename Arg> struct STOUT : computeStapleOps {
using real = typename Arg::Float;
using Complex = complex<real>;
using Link = Matrix<complex<real>, Arg::nColor>;
const Arg &arg;
template <typename... OpsArgs> constexpr STOUT(const Arg &arg, const OpsArgs &...ops) : KernelOpsT(ops...), arg(arg)
{
}
static constexpr const char *filename() { return KERNEL_FILE; }
__device__ __host__ inline void operator()(int x_cb, int parity, int dir)
{
// Compute spacetime and local coords
int X[4];
for (int dr = 0; dr < 4; ++dr) X[dr] = arg.X[dr];
int x[4];
getCoords(x, x_cb, X, parity);
for (int dr = 0; dr < 4; ++dr) {
x[dr] += arg.border[dr];
X[dr] += 2 * arg.border[dr];
}
dir = dir + (dir >= arg.dir_ignore);
Link U, Stap, Q;
// This function gets stap = S_{mu,nu} i.e., the staple of length 3,
computeStaple(*this, x, X, parity, dir, Stap, arg.dir_ignore);
// Get link U
U = arg.in(dir, linkIndex(x, X), parity);
// Compute Omega_{mu}=[Sum_{mu neq nu}rho_{mu,nu}C_{mu,nu}]*U_{mu}^dag
//--------------------------------------------------------------------
// Compute \Omega = \rho * S * U^{\dagger}
Q = (arg.rho * Stap) * conj(U);
// Compute \Q_{mu} = i/2[Omega_{mu}^dag - Omega_{mu}
// - 1/3 Tr(Omega_{mu}^dag - Omega_{mu})]
makeHerm(Q);
// Q is now defined.
Link exp_iQ = exponentiate_iQ(Q);
U = exp_iQ * U;
arg.out(dir, linkIndex(x, X), parity) = U;
// Debug tools
#if 0
//Test for Traceless:
double error = getTrace(Q).real();
printf("Trace test %d %d %.15e\n", x_cb, dir, error);
//Test for hermiticity:
Link Q_diff = conj(Q) - Q; //This should be the zero matrix. Test by ReTr(Q_diff^2);
Q_diff *= Q_diff;
error = getTrace(Q_diff).real();
printf("Herm test %d %d %.15e\n", x_cb, dir, error);
//Test for expiQ unitarity:
error = ErrorSU3(exp_iQ);
printf("expiQ test %d %d %.15e\n", x_cb, dir, error);
//Test for expiQ*U unitarity:
error = ErrorSU3(U);
printf("expiQ*u test %d %d %.15e\n", x_cb, dir, error);
#endif
}
};
//------------------------//
// Over-Improved routines //
//------------------------//
template <typename Arg> struct OvrImpSTOUTOps {
using real = typename Arg::Float;
using Complex = complex<real>;
using Link = Matrix<complex<real>, Arg::nColor>;
using StapCacheT = ThreadLocalCache<Link, 0, computeStapleRectangleOps>; // offset by computeStapleRectangleOps
using RectCacheT = ThreadLocalCache<Link, 0, StapCacheT>; // offset by StapCacheT
using Ops = combineOps<computeStapleRectangleOps, KernelOps<StapCacheT, RectCacheT>>;
};
template <typename Arg> struct OvrImpSTOUT : OvrImpSTOUTOps<Arg>::Ops {
using real = typename Arg::Float;
using Complex = complex<real>;
using Link = Matrix<complex<real>, Arg::nColor>;
using typename OvrImpSTOUTOps<Arg>::Ops::KernelOpsT;
const Arg &arg;
template <typename... OpsArgs>
constexpr OvrImpSTOUT(const Arg &arg, const OpsArgs &...ops) : KernelOpsT(ops...), arg(arg)
{
}
static constexpr const char *filename() { return KERNEL_FILE; }
__device__ __host__ inline void operator()(int x_cb, int parity, int dir)
{
// Compute spacetime and local coords
int X[4];
for (int dr = 0; dr < 4; ++dr) X[dr] = arg.X[dr];
int x[4];
getCoords(x, x_cb, X, parity);
for (int dr = 0; dr < 4; ++dr) {
x[dr] += arg.border[dr];
X[dr] += 2 * arg.border[dr];
}
dir = dir + (dir >= arg.dir_ignore);
Link U, Q;
typename OvrImpSTOUTOps<Arg>::StapCacheT Stap {*this};
typename OvrImpSTOUTOps<Arg>::RectCacheT Rect {*this};
// This function gets stap = S_{mu,nu} i.e., the staple of length 3,
// and the 1x2 and 2x1 rectangles of length 5. From the following paper:
// https://arxiv.org/abs/0801.1165
computeStapleRectangle(*this, x, X, parity, dir, Stap, Rect, arg.dir_ignore);
// Get link U
U = arg.in(dir, linkIndex(x, X), parity);
// Compute Omega_{mu}=[Sum_{mu neq nu}rho_{mu,nu}C_{mu,nu}]*U_{mu}^dag
//-------------------------------------------------------------------
// Compute \rho * staple_coeff * S - \rho * rectangle_coeff * R
Q = ((arg.staple_coeff * static_cast<const Link &>(Stap)) - (arg.rectangle_coeff * static_cast<const Link &>(Rect)))
* conj(U);
// Compute \Q_{mu} = i/2[Omega_{mu}^dag - Omega_{mu}
// - 1/3 Tr(Omega_{mu}^dag - Omega_{mu})]
makeHerm(Q);
// Q is now defined.
Link exp_iQ = exponentiate_iQ(Q);
U = exp_iQ * U;
arg.out(dir, linkIndex(x, X), parity) = U;
// Debug tools
#if 0
//Test for Traceless:
double error = getTrace(Q).real();
printf("Trace test %d %d %.15e\n", x_cb, dir, error);
//Test for hermiticity:
Link Q_diff = conj(Q) - Q; //This should be the zero matrix. Test by ReTr(Q_diff^2);
Q_diff *= Q_diff;
error = getTrace(Q_diff).real();
printf("Herm test %d %d %.15e\n", x_cb, dir, error);
//Test for expiQ unitarity:
error = ErrorSU3(exp_iQ);
printf("expiQ test %d %d %.15e\n", x_cb, dir, error);
//Test for expiQ*U unitarity:
error = ErrorSU3(U);
printf("expiQ*u test %d %d %.15e\n", x_cb, dir, error);
#endif
}
};
} // namespace quda