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,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
41 #include "gromacs/math/functions.h"
42 #include "gromacs/math/units.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/mdlib/force.h"
45 #include "gromacs/mdtypes/inputrec.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/mdtypes/mdatom.h"
48 #include "gromacs/pbcutil/ishift.h"
49 #include "gromacs/pbcutil/mshift.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/topology/topology.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/pleasecite.h"
55 real
RF_excl_correction(const t_forcerec
*fr
, t_graph
*g
,
56 const t_mdatoms
*mdatoms
, const t_blocka
*excl
,
57 rvec x
[], rvec f
[], rvec
*fshift
, const t_pbc
*pbc
,
58 real lambda
, real
*dvdlambda
)
60 /* Calculate the reaction-field energy correction for this node:
61 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
62 * and force correction for all excluded pairs, including self pairs.
64 int i
, j
, j1
, j2
, k
, ki
;
65 double q2sumA
, q2sumB
, ener
;
66 const real
*chargeA
, *chargeB
;
67 real ek
, ec
, L1
, qiA
, qiB
, qqA
, qqB
, qqL
, v
;
72 int end
= mdatoms
->homenr
;
74 gmx_bool bMolPBC
= fr
->bMolPBC
;
78 /* For test particle insertion we only correct for the test molecule */
79 start
= mdatoms
->nr
- fr
->n_tpi
;
82 ek
= fr
->ic
->epsfac
*fr
->ic
->k_rf
;
83 ec
= fr
->ic
->epsfac
*fr
->ic
->c_rf
;
84 chargeA
= mdatoms
->chargeA
;
85 chargeB
= mdatoms
->chargeB
;
101 if (mdatoms
->nChargePerturbed
== 0)
103 for (i
= start
; i
< niat
; i
++)
110 /* Do the exclusions */
112 j2
= excl
->index
[i
+1];
113 for (j
= j1
; j
< j2
; j
++)
118 qqA
= qiA
*chargeA
[k
];
123 rvec_sub(x
[i
], x
[k
], dx
);
124 ivec_sub(SHIFT_IVEC(g
, i
), SHIFT_IVEC(g
, k
), dt
);
129 ki
= pbc_dx_aiuc(pbc
, x
[i
], x
[k
], dx
);
133 rvec_sub(x
[i
], x
[k
], dx
);
135 ener
+= qqA
*(ek
*norm2(dx
) - ec
);
136 svmul(-2*qqA
*ek
, dx
, df
);
139 rvec_inc(fshift
[ki
], df
);
140 rvec_dec(fshift
[CENTRAL
], df
);
145 ener
+= -0.5*ec
*q2sumA
;
150 for (i
= start
; i
< niat
; i
++)
159 /* Do the exclusions */
161 j2
= excl
->index
[i
+1];
162 for (j
= j1
; j
< j2
; j
++)
167 qqA
= qiA
*chargeA
[k
];
168 qqB
= qiB
*chargeB
[k
];
169 if (qqA
!= 0 || qqB
!= 0)
171 qqL
= L1
*qqA
+ lambda
*qqB
;
174 rvec_sub(x
[i
], x
[k
], dx
);
175 ivec_sub(SHIFT_IVEC(g
, i
), SHIFT_IVEC(g
, k
), dt
);
180 ki
= pbc_dx_aiuc(pbc
, x
[i
], x
[k
], dx
);
184 rvec_sub(x
[i
], x
[k
], dx
);
186 v
= ek
*norm2(dx
) - ec
;
188 svmul(-2*qqL
*ek
, dx
, df
);
191 rvec_inc(fshift
[ki
], df
);
192 rvec_dec(fshift
[CENTRAL
], df
);
193 *dvdlambda
+= (qqB
- qqA
)*v
;
198 ener
+= -0.5*ec
*(L1
*q2sumA
+ lambda
*q2sumB
);
199 *dvdlambda
+= -0.5*ec
*(q2sumB
- q2sumA
);
204 fprintf(debug
, "RF exclusion energy: %g\n", ener
);
210 void calc_rffac(FILE *fplog
, int eel
, real eps_r
, real eps_rf
, real Rc
, real Temp
,
211 real zsq
, matrix box
,
212 real
*krf
, real
*crf
)
214 /* Compute constants for Generalized reaction field */
215 real kappa
, k1
, k2
, I
, rmin
;
222 /* Consistency check */
225 gmx_fatal(FARGS
, "Temperature is %f while using"
226 " Generalized Reaction Field\n", Temp
);
228 /* Ionic strength (only needed for eelGRF */
231 kappa
= std::sqrt(2*I
/(EPSILON0
*eps_rf
*BOLTZ
*Temp
));
239 /* eps == 0 signals infinite dielectric */
242 *krf
= 1/(2*Rc
*Rc
*Rc
);
247 k2
= eps_rf
*gmx::square((real
)(kappa
*Rc
));
249 *krf
= ((eps_rf
- eps_r
)*k1
+ 0.5*k2
)/((2*eps_rf
+ eps_r
)*k1
+ k2
)/(Rc
*Rc
*Rc
);
251 *crf
= 1/Rc
+ *krf
*Rc
*Rc
;
252 // Make sure we don't lose resolution in pow() by casting real arg to double
253 rmin
= gmx::invcbrt(static_cast<double>(*krf
*2.0));
259 please_cite(fplog
, "Tironi95a");
260 fprintf(fplog
, "%s:\n"
261 "epsRF = %10g, I = %10g, volume = %10g, kappa = %10g\n"
262 "rc = %10g, krf = %10g, crf = %10g, epsfac = %10g\n",
263 eel_names
[eel
], eps_rf
, I
, vol
, kappa
, Rc
, *krf
, *crf
,
268 fprintf(fplog
, "%s:\n"
269 "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
270 eel_names
[eel
], eps_rf
, Rc
, *krf
, *crf
, ONE_4PI_EPS0
/eps_r
);
273 "The electrostatics potential has its minimum at r = %g\n",
279 void init_generalized_rf(FILE *fplog
,
280 const gmx_mtop_t
*mtop
, const t_inputrec
*ir
,
284 real q
, zsq
, nrdf
, T
;
285 const gmx_moltype_t
*molt
;
288 if (ir
->efep
!= efepNO
&& fplog
)
290 fprintf(fplog
, "\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
293 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
295 molt
= &mtop
->moltype
[molb
.type
];
297 for (i
= 0; (i
< cgs
->nr
); i
++)
300 for (j
= cgs
->index
[i
]; (j
< cgs
->index
[i
+1]); j
++)
302 q
+= molt
->atoms
.atom
[j
].q
;
304 zsq
+= molb
.nmol
*q
*q
;
311 for (i
= 0; (i
< ir
->opts
.ngtc
); i
++)
313 nrdf
+= ir
->opts
.nrdf
[i
];
314 T
+= (ir
->opts
.nrdf
[i
] * ir
->opts
.ref_t
[i
]);
318 gmx_fatal(FARGS
, "No degrees of freedom!");