forked from etmc/tmLQCD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathupdate_momenta.c
76 lines (68 loc) · 2.15 KB
/
update_momenta.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
/***********************************************************************
*
* Copyright (C) 2001 Martin Hasebusch
* 2002,2003,2004,2005,2006,2007,2008,2012 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<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 "monomial/monomial.h"
#include "xchange/xchange.h"
#include "operator/clover_leaf.h"
#include "read_input.h"
#include "hamiltonian_field.h"
#include "update_momenta.h"
#include "gettime.h"
/* Updates the momenta: equation 16 of Gottlieb */
void update_momenta(int * mnllist, double step, const int no,
hamiltonian_field_t * const hf) {
#ifdef OMP
#pragma omp parallel for
#endif
for(int i = 0; i < (VOLUMEPLUSRAND + g_dbw2rand);i++) {
for(int mu=0;mu<4;mu++) {
_zero_su3adj(hf->derivative[i][mu]);
}
}
for(int k = 0; k < no; k++) {
if(monomial_list[ mnllist[k] ].derivativefunction != NULL) {
monomial_list[ mnllist[k] ].derivativefunction(mnllist[k], hf);
}
}
#ifdef MPI
xchange_deri(hf->derivative);
#endif
#ifdef OMP
#pragma omp parallel for
#endif
for(int i = 0; i < VOLUME; i++) {
for(int mu = 0; mu < 4; mu++) {
/* the minus comes from an extra minus in trace_lambda */
_su3adj_minus_const_times_su3adj(hf->momenta[i][mu], step, hf->derivative[i][mu]);
}
}
return;
}