532be50cd912a0a8bc9cc6d57eccb5d6d1de7681
[alexxy/gromacs.git] / src / gromacs / ewald / long_range_correction.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "long_range_correction.h"
41
42 #include <cmath>
43
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/md_enums.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/gmxassert.h"
55
56 /* There's nothing special to do here if just masses are perturbed,
57  * but if either charge or type is perturbed then the implementation
58  * requires that B states are defined for both charge and type, and
59  * does not optimize for the cases where only one changes.
60  *
61  * The parameter vectors for B states are left undefined in atoms2md()
62  * when either FEP is inactive, or when there are no mass/charge/type
63  * perturbations. The parameter vectors for LJ-PME are likewise
64  * undefined when LJ-PME is not active. This works because
65  * bHaveChargeOrTypePerturbed handles the control flow. */
66 void ewald_LRcorrection(const int         numAtomsLocal,
67                         const t_commrec*  cr,
68                         int               numThreads,
69                         int               thread,
70                         const t_forcerec& fr,
71                         const t_inputrec& ir,
72                         const real*       chargeA,
73                         const real*       chargeB,
74                         gmx_bool          bHaveChargePerturbed,
75                         const rvec        x[],
76                         const matrix      box,
77                         const rvec        mu_tot[],
78                         rvec*             f,
79                         real*             Vcorr_q,
80                         real              lambda_q,
81                         real*             dvdlambda_q)
82 {
83     /* We need to correct only self interactions */
84     const int start = (numAtomsLocal * thread) / numThreads;
85     const int end   = (numAtomsLocal * (thread + 1)) / numThreads;
86
87     int    i, j, q;
88     double Vexcl_q, dvdl_excl_q; /* Necessary for precision */
89     real   one_4pi_eps;
90     real   Vself_q[2], Vdipole[2];
91     rvec   mutot[2], dipcorrA, dipcorrB;
92     real   L1_q, dipole_coeff;
93     real   chargecorr[2] = { 0, 0 };
94
95     /* Scale the Ewald unit cell when dimension z is not periodic */
96     matrix          scaledBox;
97     EwaldBoxZScaler boxScaler(ir);
98     boxScaler.scaleBox(box, scaledBox);
99
100     one_4pi_eps = ONE_4PI_EPS0 / fr.ic->epsilon_r;
101     Vexcl_q     = 0;
102     dvdl_excl_q = 0;
103     Vdipole[0]  = 0;
104     Vdipole[1]  = 0;
105     L1_q        = 1.0 - lambda_q;
106     /* Note that we have to transform back to gromacs units, since
107      * mu_tot contains the dipole in debye units (for output).
108      */
109     for (i = 0; (i < DIM); i++)
110     {
111         mutot[0][i] = mu_tot[0][i] * DEBYE2ENM;
112         mutot[1][i] = mu_tot[1][i] * DEBYE2ENM;
113         dipcorrA[i] = 0;
114         dipcorrB[i] = 0;
115     }
116     dipole_coeff = 0;
117
118     real boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
119     switch (ir.ewald_geometry)
120     {
121         case eewg3D:
122             if (ir.epsilon_surface != 0)
123             {
124                 dipole_coeff = 2 * M_PI * ONE_4PI_EPS0
125                                / ((2 * ir.epsilon_surface + fr.ic->epsilon_r) * boxVolume);
126                 for (i = 0; (i < DIM); i++)
127                 {
128                     dipcorrA[i] = 2 * dipole_coeff * mutot[0][i];
129                     dipcorrB[i] = 2 * dipole_coeff * mutot[1][i];
130                 }
131             }
132             break;
133         case eewg3DC:
134             dipole_coeff = 2 * M_PI * one_4pi_eps / boxVolume;
135             dipcorrA[ZZ] = 2 * dipole_coeff * mutot[0][ZZ];
136             dipcorrB[ZZ] = 2 * dipole_coeff * mutot[1][ZZ];
137             for (int q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
138             {
139                 /* Avoid charge corrections with near-zero net charge */
140                 if (fabs(fr.qsum[q]) > 1e-4)
141                 {
142                     chargecorr[q] = 2 * dipole_coeff * fr.qsum[q];
143                 }
144             }
145             break;
146         default: gmx_incons("Unsupported Ewald geometry");
147     }
148     if (debug)
149     {
150         fprintf(debug, "dipcorr = %8.3f  %8.3f  %8.3f\n", dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
151         fprintf(debug, "mutot   = %8.3f  %8.3f  %8.3f\n", mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
152     }
153     const bool bNeedLongRangeCorrection = (dipole_coeff != 0);
154     if (bNeedLongRangeCorrection && !bHaveChargePerturbed)
155     {
156         for (i = start; (i < end); i++)
157         {
158             for (j = 0; (j < DIM); j++)
159             {
160                 f[i][j] -= dipcorrA[j] * chargeA[i];
161             }
162             if (chargecorr[0] != 0)
163             {
164                 f[i][ZZ] += chargecorr[0] * chargeA[i] * x[i][ZZ];
165             }
166         }
167     }
168     else if (bNeedLongRangeCorrection)
169     {
170         for (i = start; (i < end); i++)
171         {
172             for (j = 0; (j < DIM); j++)
173             {
174                 f[i][j] -= L1_q * dipcorrA[j] * chargeA[i] + lambda_q * dipcorrB[j] * chargeB[i];
175             }
176             if (chargecorr[0] != 0 || chargecorr[1] != 0)
177             {
178                 f[i][ZZ] += (L1_q * chargecorr[0] * chargeA[i] + lambda_q * chargecorr[1]) * x[i][ZZ];
179             }
180         }
181     }
182
183     Vself_q[0] = 0;
184     Vself_q[1] = 0;
185
186     /* Global corrections only on master process */
187     if (MASTER(cr) && thread == 0)
188     {
189         for (q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
190         {
191             /* Apply surface and charged surface dipole correction:
192              * correction = dipole_coeff * ( (dipole)^2
193              *              - qsum*sum_i q_i z_i^2 - qsum^2 * box_z^2 / 12 )
194              */
195             if (dipole_coeff != 0)
196             {
197                 if (ir.ewald_geometry == eewg3D)
198                 {
199                     Vdipole[q] = dipole_coeff * iprod(mutot[q], mutot[q]);
200                 }
201                 else if (ir.ewald_geometry == eewg3DC)
202                 {
203                     Vdipole[q] = dipole_coeff * mutot[q][ZZ] * mutot[q][ZZ];
204
205                     if (chargecorr[q] != 0)
206                     {
207                         /* Here we use a non thread-parallelized loop,
208                          * because this is the only loop over atoms for
209                          * energies and they need reduction (unlike forces).
210                          * We could implement a reduction over threads,
211                          * but this case is rarely used.
212                          */
213                         const real* qPtr   = (q == 0 ? chargeA : chargeB);
214                         real        sumQZ2 = 0;
215                         for (int i = 0; i < numAtomsLocal; i++)
216                         {
217                             sumQZ2 += qPtr[i] * x[i][ZZ] * x[i][ZZ];
218                         }
219                         Vdipole[q] -= dipole_coeff * fr.qsum[q]
220                                       * (sumQZ2 + fr.qsum[q] * box[ZZ][ZZ] * box[ZZ][ZZ] / 12);
221                     }
222                 }
223             }
224         }
225     }
226     if (!bHaveChargePerturbed)
227     {
228         *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
229     }
230     else
231     {
232         *Vcorr_q = L1_q * (Vdipole[0] - Vself_q[0]) + lambda_q * (Vdipole[1] - Vself_q[1]) - Vexcl_q;
233         *dvdlambda_q += Vdipole[1] - Vself_q[1] - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
234     }
235
236     if (debug)
237     {
238         fprintf(debug, "Long Range corrections for Ewald interactions:\n");
239         fprintf(debug, "q2sum = %g, Vself_q=%g\n", L1_q * fr.q2sum[0] + lambda_q * fr.q2sum[1],
240                 L1_q * Vself_q[0] + lambda_q * Vself_q[1]);
241         fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
242         if (MASTER(cr) && thread == 0)
243         {
244             if (ir.epsilon_surface > 0 || ir.ewald_geometry == eewg3DC)
245             {
246                 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
247                         L1_q * Vdipole[0] + lambda_q * Vdipole[1]);
248             }
249         }
250     }
251 }