Simplify PME GPU constants
[gromacs.git] / src / gromacs / ewald / long_range_correction.cpp
blobb62a25e58ccb22529fffbab3f8846b1212730747
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include "long_range_correction.h"
42 #include <cmath>
44 #include "gromacs/ewald/ewald_utils.h"
45 #include "gromacs/math/functions.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/forcerec.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/interaction_const.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/gmxassert.h"
57 /* There's nothing special to do here if just masses are perturbed,
58 * but if either charge or type is perturbed then the implementation
59 * requires that B states are defined for both charge and type, and
60 * does not optimize for the cases where only one changes.
62 * The parameter vectors for B states are left undefined in atoms2md()
63 * when either FEP is inactive, or when there are no mass/charge/type
64 * perturbations. The parameter vectors for LJ-PME are likewise
65 * undefined when LJ-PME is not active. This works because
66 * bHaveChargeOrTypePerturbed handles the control flow. */
67 void ewald_LRcorrection(const int numAtomsLocal,
68 const t_commrec* cr,
69 int numThreads,
70 int thread,
71 const t_forcerec& fr,
72 const t_inputrec& ir,
73 const real* chargeA,
74 const real* chargeB,
75 gmx_bool bHaveChargePerturbed,
76 const rvec x[],
77 const matrix box,
78 const rvec mu_tot[],
79 rvec* f,
80 real* Vcorr_q,
81 real lambda_q,
82 real* dvdlambda_q)
84 /* We need to correct only self interactions */
85 const int start = (numAtomsLocal * thread) / numThreads;
86 const int end = (numAtomsLocal * (thread + 1)) / numThreads;
88 int i, j, q;
89 double Vexcl_q, dvdl_excl_q; /* Necessary for precision */
90 real one_4pi_eps;
91 real Vself_q[2], Vdipole[2];
92 rvec mutot[2], dipcorrA, dipcorrB;
93 real L1_q, dipole_coeff;
94 real chargecorr[2] = { 0, 0 };
96 /* Scale the Ewald unit cell when dimension z is not periodic */
97 matrix scaledBox;
98 EwaldBoxZScaler boxScaler(ir);
99 boxScaler.scaleBox(box, scaledBox);
101 one_4pi_eps = ONE_4PI_EPS0 / fr.ic->epsilon_r;
102 Vexcl_q = 0;
103 dvdl_excl_q = 0;
104 Vdipole[0] = 0;
105 Vdipole[1] = 0;
106 L1_q = 1.0 - lambda_q;
107 /* Note that we have to transform back to gromacs units, since
108 * mu_tot contains the dipole in debye units (for output).
110 for (i = 0; (i < DIM); i++)
112 mutot[0][i] = mu_tot[0][i] * DEBYE2ENM;
113 mutot[1][i] = mu_tot[1][i] * DEBYE2ENM;
114 dipcorrA[i] = 0;
115 dipcorrB[i] = 0;
117 dipole_coeff = 0;
119 real boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
120 switch (ir.ewald_geometry)
122 case eewg3D:
123 if (ir.epsilon_surface != 0)
125 dipole_coeff = 2 * M_PI * ONE_4PI_EPS0
126 / ((2 * ir.epsilon_surface + fr.ic->epsilon_r) * boxVolume);
127 for (i = 0; (i < DIM); i++)
129 dipcorrA[i] = 2 * dipole_coeff * mutot[0][i];
130 dipcorrB[i] = 2 * dipole_coeff * mutot[1][i];
133 break;
134 case eewg3DC:
135 dipole_coeff = 2 * M_PI * one_4pi_eps / boxVolume;
136 dipcorrA[ZZ] = 2 * dipole_coeff * mutot[0][ZZ];
137 dipcorrB[ZZ] = 2 * dipole_coeff * mutot[1][ZZ];
138 for (int q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
140 /* Avoid charge corrections with near-zero net charge */
141 if (fabs(fr.qsum[q]) > 1e-4)
143 chargecorr[q] = 2 * dipole_coeff * fr.qsum[q];
146 break;
147 default: gmx_incons("Unsupported Ewald geometry");
149 if (debug)
151 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n", dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
152 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n", mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
154 const bool bNeedLongRangeCorrection = (dipole_coeff != 0);
155 if (bNeedLongRangeCorrection && !bHaveChargePerturbed)
157 for (i = start; (i < end); i++)
159 for (j = 0; (j < DIM); j++)
161 f[i][j] -= dipcorrA[j] * chargeA[i];
163 if (chargecorr[0] != 0)
165 f[i][ZZ] += chargecorr[0] * chargeA[i] * x[i][ZZ];
169 else if (bNeedLongRangeCorrection)
171 for (i = start; (i < end); i++)
173 for (j = 0; (j < DIM); j++)
175 f[i][j] -= L1_q * dipcorrA[j] * chargeA[i] + lambda_q * dipcorrB[j] * chargeB[i];
177 if (chargecorr[0] != 0 || chargecorr[1] != 0)
179 f[i][ZZ] += (L1_q * chargecorr[0] * chargeA[i] + lambda_q * chargecorr[1]) * x[i][ZZ];
184 Vself_q[0] = 0;
185 Vself_q[1] = 0;
187 /* Global corrections only on master process */
188 if (MASTER(cr) && thread == 0)
190 for (q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
192 /* Apply surface and charged surface dipole correction:
193 * correction = dipole_coeff * ( (dipole)^2
194 * - qsum*sum_i q_i z_i^2 - qsum^2 * box_z^2 / 12 )
196 if (dipole_coeff != 0)
198 if (ir.ewald_geometry == eewg3D)
200 Vdipole[q] = dipole_coeff * iprod(mutot[q], mutot[q]);
202 else if (ir.ewald_geometry == eewg3DC)
204 Vdipole[q] = dipole_coeff * mutot[q][ZZ] * mutot[q][ZZ];
206 if (chargecorr[q] != 0)
208 /* Here we use a non thread-parallelized loop,
209 * because this is the only loop over atoms for
210 * energies and they need reduction (unlike forces).
211 * We could implement a reduction over threads,
212 * but this case is rarely used.
214 const real* qPtr = (q == 0 ? chargeA : chargeB);
215 real sumQZ2 = 0;
216 for (int i = 0; i < numAtomsLocal; i++)
218 sumQZ2 += qPtr[i] * x[i][ZZ] * x[i][ZZ];
220 Vdipole[q] -= dipole_coeff * fr.qsum[q]
221 * (sumQZ2 + fr.qsum[q] * box[ZZ][ZZ] * box[ZZ][ZZ] / 12);
227 if (!bHaveChargePerturbed)
229 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
231 else
233 *Vcorr_q = L1_q * (Vdipole[0] - Vself_q[0]) + lambda_q * (Vdipole[1] - Vself_q[1]) - Vexcl_q;
234 *dvdlambda_q += Vdipole[1] - Vself_q[1] - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
237 if (debug)
239 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
240 fprintf(debug, "q2sum = %g, Vself_q=%g\n", L1_q * fr.q2sum[0] + lambda_q * fr.q2sum[1],
241 L1_q * Vself_q[0] + lambda_q * Vself_q[1]);
242 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
243 if (MASTER(cr) && thread == 0)
245 if (ir.epsilon_surface > 0 || ir.ewald_geometry == eewg3DC)
247 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
248 L1_q * Vdipole[0] + lambda_q * Vdipole[1]);