forked from etmc/tmLQCD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
get_rectangle_staples.c
187 lines (172 loc) · 4.51 KB
/
get_rectangle_staples.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
/***********************************************************************
* Copyright (C) 2002,2003,2004,2005,2006,2007,2008 Carsten Urbach
*
* 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/>.
***********************************************************************/
#ifdef HAVE_CONFIG_H
#include "tmlqcd_config.h"
#endif
#include <stdlib.h>
#include <stdio.h>
#include "global.h"
#include "su3.h"
#include "get_rectangle_staples.h"
void get_rectangle_staples(su3 * const v, const int x, const int mu) {
get_rectangle_staples_general(v,x,mu,g_gauge_field);
}
void get_rectangle_staples_general(su3 * const v, const int x, const int mu, const su3** const gf) {
su3 ALIGN tmp1, tmp2;
int y, z, nu;
su3 * a, * b, * c, * d, * e;
#ifdef _KOJAK_INST
#pragma pomp inst begin(rectstaples)
#endif
#ifdef XLC
#pragma disjoint(*v, tmp1, tmp2, *a, *b, *c, *d, *e)
#endif
_su3_zero((*v));
for(nu = 0; nu < 4; nu++) {
if(mu != nu) {
/* first contr. starting from x
* a b c e^+ d^+
* c
* _
* b| |e
* a| |d
*/
a = &gf[x][nu];
y = g_iup[x][nu];
b = &gf[y][nu];
_su3_times_su3(tmp1, *a, *b);
z = g_iup[y][nu];
c = &gf[z][mu];
_su3_times_su3(tmp2, tmp1, *c);
y = g_iup[x][mu];
d = &gf[y][nu];
z = g_iup[y][nu];
e = &gf[z][nu];
_su3_times_su3(tmp1, *d, *e);
_su3_times_su3d_acc((*v), tmp2, tmp1);
/* 1 contr. starting idn[idn[x][nu]][nu]
* e^+ d^+ a b c
*
*e| |c
*d|_|b
* a
*/
y = g_idn[x][nu];
z = g_idn[y][nu];
d = &gf[z][nu];
a = &gf[z][mu];
_su3d_times_su3(tmp1, *d, *a);
e = &gf[y][nu];
_su3d_times_su3(tmp2, *e, tmp1);
y = g_iup[z][mu];
b = &gf[y][nu];
z = g_iup[y][nu];
c = &gf[z][nu];
_su3_times_su3(tmp1, *b, *c);
_su3_times_su3_acc((*v), tmp2, tmp1);
/* second contr. starting from x
* a b c e^+ d^+
*
* bc
* __
* a| _|e
* d
*/
a = &gf[x][nu];
y = g_iup[x][nu];
b = &gf[y][mu];
_su3_times_su3(tmp1, *a, *b);
z = g_iup[y][mu];
c = &gf[z][mu];
_su3_times_su3(tmp2, tmp1, *c);
y = g_iup[x][mu];
d = &gf[y][mu];
z = g_iup[y][mu];
e = &gf[z][nu];
_su3_times_su3(tmp1, *d, *e);
_su3_times_su3d_acc((*v), tmp2, tmp1);
/* 1 contr. starting idn[x][nu]
* d^+ a b c e^+
*
* e
* _
* d|__|c
* ab
*/
y = g_idn[x][nu];
d = &gf[y][nu];
a = &gf[y][mu];
_su3d_times_su3(tmp1, *d, *a);
z = g_iup[y][mu];
b = &gf[z][mu];
_su3_times_su3(tmp2, tmp1, *b);
y = g_iup[z][mu];
c = &gf[y][nu];
z = g_iup[x][mu];
e = &gf[z][mu];
_su3_times_su3d(tmp1, *c, *e);
_su3_times_su3_acc((*v), tmp2, tmp1);
/* 1 contr. starting idn[idn[x][mu]][nu]
* e^+ d^+ a b c
*
* e
* _
*d|__|c
* ab
*/
y = g_idn[x][mu];
z = g_idn[y][nu];
d = &gf[z][nu];
a = &gf[z][mu];
_su3d_times_su3(tmp1, *d, *a);
e = &gf[y][mu];
_su3d_times_su3(tmp2, *e, tmp1);
y = g_idn[x][nu];
b = &gf[y][mu];
z = g_iup[y][mu];
c = &gf[z][nu];
_su3_times_su3(tmp1, *b, *c);
_su3_times_su3_acc((*v), tmp2, tmp1);
/* 1 contr. starting idn[x][mu]
* d^+ a b c e^+
*
* bc
* __
*a|_ |e
* d
*/
y = g_idn[x][mu];
d = &gf[y][mu];
z = g_iup[y][nu];
a = &gf[y][nu];
_su3d_times_su3(tmp1, *d, *a);
b = &gf[z][mu];
_su3_times_su3(tmp2, tmp1, *b);
y = g_iup[x][mu];
e = &gf[y][nu];
z = g_iup[x][nu];
c = &gf[z][mu];
_su3_times_su3d(tmp1, *c, *e);
_su3_times_su3_acc((*v), tmp2, tmp1);
}
}
#ifdef _KOJAK_INST
#pragma pomp inst end(rectstaples)
#endif
}