-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathupdate_gauge.c
106 lines (97 loc) · 2.51 KB
/
update_gauge.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
/***********************************************************************
*
* Copyright (C) 2001 Martin Hasebusch
*
* some changes by C. Urbach 2002-2008,2012
*
* 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<config.h>
#endif
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <time.h>
#include "global.h"
#include "su3.h"
#include "su3adj.h"
#include "su3spinor.h"
#include "expo.h"
#include "sse.h"
#include "xchange.h"
#include "hamiltonian_field.h"
#include "update_gauge.h"
/*******************************************************
*
* Updates the gauge field corresponding to the momenta
*
*******************************************************/
void update_gauge(const double step, hamiltonian_field_t * const hf) {
#ifdef OMP
#define static
#pragma omp parallel
{
#endif
int i,mu;
static su3 v,w;
su3 *z;
static su3adj deriv;
su3adj *xm;
#ifdef _KOJAK_INST
#pragma pomp inst begin(updategauge)
#endif
#ifdef OMP
#undef static
#endif
#ifdef OMP
#pragma omp for
#endif
for(i = 0; i < VOLUME; i++) {
for(mu = 0; mu < 4; mu++){
/* moment[i][mu] = h_{i,mu}^{alpha} */
xm = &hf->momenta[i][mu];
z = &hf->gaugefield[i][mu];
_su3adj_assign_const_times_su3adj(deriv, step, *xm);
exposu3(&w,&deriv);
restoresu3(&v,&w);
_su3_times_su3(w, v, *z);
_su3_assign(*z, w);
}
}
#ifdef OMP
} /* OpenMP parallel closing brace */
#endif
#ifdef MPI
/* for parallelization */
xchange_gauge(hf->gaugefield);
#endif
/*
* The backward copy of the gauge field
* is not updated here!
*/
hf->update_gauge_copy = 1;
g_update_gauge_copy = 1;
hf->update_gauge_energy = 1;
g_update_gauge_energy = 1;
hf->update_rectangle_energy = 1;
g_update_rectangle_energy = 1;
return;
#ifdef _KOJAK_INST
#pragma pomp inst end(updategauge)
#endif
}