forked from etmc/tmLQCD
-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathmeasure_rectangles.c
140 lines (129 loc) · 3.06 KB
/
measure_rectangles.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
/***********************************************************************
* 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/>.
***********************************************************************/
/*******************************************************************
*
* Here the 1x2 rectangles are implemented
* for renormalization group improved gauge
* actions like the DBW2 or the Iwasaki
* gauge action.
*
* 1/3 \sum_{\mu\leq\nu;\mu,nu=1}^4 Tr U^{1x2}
*
* author: Carsten Urbach
*
*******************************************************************/
#ifdef HAVE_CONFIG_H
#include "tmlqcd_config.h"
#endif
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#ifdef TM_USE_OMP
#include <omp.h>
#endif
#include "global.h"
#include "sse.h"
#include "su3.h"
#include "su3adj.h"
#include "geometry_eo.h"
#include "measure_rectangles.h"
double measure_rectangles(const su3 ** const gf) {
static double res;
#ifdef TM_USE_MPI
double ALIGN mres;
#endif
#ifdef TM_USE_OMP
#pragma omp parallel
{
int thread_num = omp_get_thread_num();
#endif
int i, j, k, mu, nu;
su3 ALIGN pr1, pr2, tmp;
const su3 *v = NULL , *w = NULL;
double ALIGN ac, ks, kc, tr, ts, tt;
kc = 0.0;
ks = 0.0;
#ifdef TM_USE_OMP
#pragma omp for
#endif
for (i = 0; i < VOLUME; i++) {
for (mu = 0; mu < 4; mu++) {
for (nu = 0; nu < 4; nu++) {
if(nu != mu) {
/*
^
|
^
|
->
*/
j = g_iup[i][mu];
k = g_iup[j][nu];
v = &gf[i][mu];
w = &gf[j][nu];
_su3_times_su3(tmp, *v, *w);
v = &gf[k][nu];
_su3_times_su3(pr1, tmp, *v);
/*
->
^
|
^
|
*/
j = g_iup[i][nu];
k = g_iup[j][nu];
v = &gf[i][nu];
w = &gf[j][nu];
_su3_times_su3(tmp, *v, *w);
v = &gf[k][mu];
_su3_times_su3(pr2, tmp, *v);
/* Trace it */
_trace_su3_times_su3d(ac,pr1,pr2);
/* printf("i mu nu: %d %d %d, ac = %e\n", i, mu, nu, ac); */
/* Kahan summation */
tr=ac+kc;
ts=tr+ks;
tt=ts-ks;
ks=ts;
kc=tr-tt;
}
}
}
}
kc=(kc+ks)/3.0;
#ifdef TM_USE_OMP
g_omp_acc_re[thread_num] = kc;
#else
res = kc;
#endif
#ifdef TM_USE_OMP
} /* OpenMP parallel closing brace */
res = 0.0;
for(int i = 0; i < omp_num_threads; ++i)
res += g_omp_acc_re[i];
#else
#endif
#ifdef TM_USE_MPI
MPI_Allreduce(&res, &mres, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
res = mres;
#endif
return res;
}