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, 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/fileio/copyrite.h"
42 #include "gromacs/legacyheaders/names.h"
43 #include "gromacs/legacyheaders/types/mdatom.h"
44 #include "gromacs/math/units.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/mdlib/force.h"
47 #include "gromacs/pbcutil/ishift.h"
48 #include "gromacs/pbcutil/mshift.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/utility/fatalerror.h"
52 real
RF_excl_correction(const t_forcerec
*fr
, t_graph
*g
,
53 const t_mdatoms
*mdatoms
, const t_blocka
*excl
,
54 rvec x
[], rvec f
[], rvec
*fshift
, const t_pbc
*pbc
,
55 real lambda
, real
*dvdlambda
)
57 /* Calculate the reaction-field energy correction for this node:
58 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
59 * and force correction for all excluded pairs, including self pairs.
61 int i
, j
, j1
, j2
, k
, ki
;
62 double q2sumA
, q2sumB
, ener
;
63 const real
*chargeA
, *chargeB
;
64 real ek
, ec
, L1
, qiA
, qiB
, qqA
, qqB
, qqL
, v
;
69 int end
= mdatoms
->homenr
;
71 gmx_bool bMolPBC
= fr
->bMolPBC
;
75 /* For test particle insertion we only correct for the test molecule */
76 start
= mdatoms
->nr
- fr
->n_tpi
;
79 ek
= fr
->epsfac
*fr
->k_rf
;
80 ec
= fr
->epsfac
*fr
->c_rf
;
81 chargeA
= mdatoms
->chargeA
;
82 chargeB
= mdatoms
->chargeB
;
98 if (mdatoms
->nChargePerturbed
== 0)
100 for (i
= start
; i
< niat
; i
++)
107 /* Do the exclusions */
109 j2
= excl
->index
[i
+1];
110 for (j
= j1
; j
< j2
; j
++)
115 qqA
= qiA
*chargeA
[k
];
120 rvec_sub(x
[i
], x
[k
], dx
);
121 ivec_sub(SHIFT_IVEC(g
, i
), SHIFT_IVEC(g
, k
), dt
);
126 ki
= pbc_dx_aiuc(pbc
, x
[i
], x
[k
], dx
);
130 rvec_sub(x
[i
], x
[k
], dx
);
132 ener
+= qqA
*(ek
*norm2(dx
) - ec
);
133 svmul(-2*qqA
*ek
, dx
, df
);
136 rvec_inc(fshift
[ki
], df
);
137 rvec_dec(fshift
[CENTRAL
], df
);
142 ener
+= -0.5*ec
*q2sumA
;
147 for (i
= start
; i
< niat
; i
++)
156 /* Do the exclusions */
158 j2
= excl
->index
[i
+1];
159 for (j
= j1
; j
< j2
; j
++)
164 qqA
= qiA
*chargeA
[k
];
165 qqB
= qiB
*chargeB
[k
];
166 if (qqA
!= 0 || qqB
!= 0)
168 qqL
= L1
*qqA
+ lambda
*qqB
;
171 rvec_sub(x
[i
], x
[k
], dx
);
172 ivec_sub(SHIFT_IVEC(g
, i
), SHIFT_IVEC(g
, k
), dt
);
177 ki
= pbc_dx_aiuc(pbc
, x
[i
], x
[k
], dx
);
181 rvec_sub(x
[i
], x
[k
], dx
);
183 v
= ek
*norm2(dx
) - ec
;
185 svmul(-2*qqL
*ek
, dx
, df
);
188 rvec_inc(fshift
[ki
], df
);
189 rvec_dec(fshift
[CENTRAL
], df
);
190 *dvdlambda
+= (qqB
- qqA
)*v
;
195 ener
+= -0.5*ec
*(L1
*q2sumA
+ lambda
*q2sumB
);
196 *dvdlambda
+= -0.5*ec
*(q2sumB
- q2sumA
);
201 fprintf(debug
, "RF exclusion energy: %g\n", ener
);
207 void calc_rffac(FILE *fplog
, int eel
, real eps_r
, real eps_rf
, real Rc
, real Temp
,
208 real zsq
, matrix box
,
209 real
*kappa
, real
*krf
, real
*crf
)
211 /* Compute constants for Generalized reaction field */
212 real k1
, k2
, I
, vol
, rmin
;
219 /* Consistency check */
222 gmx_fatal(FARGS
, "Temperature is %f while using"
223 " Generalized Reaction Field\n", Temp
);
225 /* Ionic strength (only needed for eelGRF */
227 *kappa
= std::sqrt(2*I
/(EPSILON0
*eps_rf
*BOLTZ
*Temp
));
235 /* eps == 0 signals infinite dielectric */
238 *krf
= 1/(2*Rc
*Rc
*Rc
);
243 k2
= eps_rf
*sqr((real
)(*kappa
*Rc
));
245 *krf
= ((eps_rf
- eps_r
)*k1
+ 0.5*k2
)/((2*eps_rf
+ eps_r
)*k1
+ k2
)/(Rc
*Rc
*Rc
);
247 *crf
= 1/Rc
+ *krf
*Rc
*Rc
;
248 // Make sure we don't lose resolution in pow() by casting real arg to double
249 rmin
= std::pow(static_cast<double>(*krf
*2.0), -1.0/3.0);
255 please_cite(fplog
, "Tironi95a");
256 fprintf(fplog
, "%s:\n"
257 "epsRF = %10g, I = %10g, volume = %10g, kappa = %10g\n"
258 "rc = %10g, krf = %10g, crf = %10g, epsfac = %10g\n",
259 eel_names
[eel
], eps_rf
, I
, vol
, *kappa
, Rc
, *krf
, *crf
,
264 fprintf(fplog
, "%s:\n"
265 "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
266 eel_names
[eel
], eps_rf
, Rc
, *krf
, *crf
, ONE_4PI_EPS0
/eps_r
);
269 "The electrostatics potential has its minimum at r = %g\n",
275 void init_generalized_rf(FILE *fplog
,
276 const gmx_mtop_t
*mtop
, const t_inputrec
*ir
,
280 real q
, zsq
, nrdf
, T
;
281 const gmx_moltype_t
*molt
;
284 if (ir
->efep
!= efepNO
&& fplog
)
286 fprintf(fplog
, "\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
289 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
291 molt
= &mtop
->moltype
[mtop
->molblock
[mb
].type
];
293 for (i
= 0; (i
< cgs
->nr
); i
++)
296 for (j
= cgs
->index
[i
]; (j
< cgs
->index
[i
+1]); j
++)
298 q
+= molt
->atoms
.atom
[j
].q
;
300 zsq
+= mtop
->molblock
[mb
].nmol
*q
*q
;
307 for (i
= 0; (i
< ir
->opts
.ngtc
); i
++)
309 nrdf
+= ir
->opts
.nrdf
[i
];
310 T
+= (ir
->opts
.nrdf
[i
] * ir
->opts
.ref_t
[i
]);
314 gmx_fatal(FARGS
, "No degrees of freedom!");