2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009 Christoph Junghans, Brad Lambeth.
5 * Copyright (c) 2011,2012, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
40 #include "types/simple.h"
54 real l2 = adressr+adressw;
63 pbc_dx(pbc, (*ref), x, dx);
67 rvec_sub((*ref), x, dx);
73 /* default to explicit simulation */
76 /* constant value for weighting function = adressw */
77 return fr->adress_const_wf;
79 /* plane through center of ref, varies in x direction */
83 /* point at center of ref, assuming cubic geometry */
84 for (i = 0; i < 3; i++)
86 sqr_dl += dx[i]*dx[i];
90 /* default to explicit simulation */
96 /* molecule is coarse grained */
101 /* molecule is explicit */
102 else if (dl < adressr)
109 tmp = cos((dl-adressr)*M_PI/2/adressw);
115 update_adress_weights_com(FILE gmx_unused * fplog,
124 int icg, k, k0, k1, d;
125 real nrcg, inv_ncg, mtot, inv_mtot;
129 real adressr, adressw;
135 int n_hyb, n_ex, n_cg;
141 adresstype = fr->adress_type;
142 adressr = fr->adress_ex_width;
143 adressw = fr->adress_hy_width;
144 massT = mdatoms->massT;
146 ref = &(fr->adress_refs);
149 /* Since this is center of mass AdResS, the vsite is not guaranteed
150 * to be on the same node as the constructing atoms. Therefore we
151 * loop over the charge groups, calculate their center of mass,
152 * then use this to calculate wf for each atom. This wastes vsite
153 * construction, but it's the only way to assure that the explicit
154 * atoms have the same wf as their vsite. */
157 fprintf(fplog, "Calculating center of mass for charge groups %d to %d\n",
160 cgindex = cgs->index;
162 /* Compute the center of mass for all charge groups */
163 for (icg = cg0; (icg < cg1); icg++)
170 wf[k0] = adress_weight(x[k0], adresstype, adressr, adressw, ref, pbc, fr);
175 else if (wf[k0] == 1)
187 for (k = k0; (k < k1); k++)
196 for (k = k0; (k < k1); k++)
198 for (d = 0; (d < DIM); d++)
200 ix[d] += x[k][d]*massT[k];
203 for (d = 0; (d < DIM); d++)
208 /* Calculate the center of gravity if the charge group mtot=0 (only vsites) */
214 for (k = k0; (k < k1); k++)
216 for (d = 0; (d < DIM); d++)
221 for (d = 0; (d < DIM); d++)
227 /* Set wf of all atoms in charge group equal to wf of com */
228 wf[k0] = adress_weight(ix, adresstype, adressr, adressw, ref, pbc, fr);
234 else if (wf[k0] == 1)
243 for (k = (k0+1); (k < k1); k++)
251 void update_adress_weights_atom_per_atom(
260 int icg, k, k0, k1, d;
261 real nrcg, inv_ncg, mtot, inv_mtot;
265 real adressr, adressw;
271 int n_hyb, n_ex, n_cg;
277 adresstype = fr->adress_type;
278 adressr = fr->adress_ex_width;
279 adressw = fr->adress_hy_width;
280 massT = mdatoms->massT;
282 ref = &(fr->adress_refs);
284 cgindex = cgs->index;
286 /* Weighting function is determined for each atom individually.
287 * This is an approximation
288 * as in the theory requires an interpolation based on the center of masses.
289 * Should be used with caution */
291 for (icg = cg0; (icg < cg1); icg++)
294 k1 = cgindex[icg + 1];
297 for (k = (k0); (k < k1); k++)
299 wf[k] = adress_weight(x[k], adresstype, adressr, adressw, ref, pbc, fr);
318 update_adress_weights_cog(t_iparams ip[],
325 int i, j, k, nr, nra, inc;
326 int ftype, adresstype;
327 t_iatom avsite, ai, aj, ak, al;
329 real adressr, adressw;
332 int n_hyb, n_ex, n_cg;
334 adresstype = fr->adress_type;
335 adressr = fr->adress_ex_width;
336 adressw = fr->adress_hy_width;
338 ref = &(fr->adress_refs);
346 /* Since this is center of geometry AdResS, we know the vsite
347 * is in the same charge group node as the constructing atoms.
348 * Loop over vsite types, calculate the weight of the vsite,
349 * then assign that weight to the constructing atoms. */
351 for (ftype = 0; (ftype < F_NRE); ftype++)
353 if (interaction_function[ftype].flags & IF_VSITE)
355 nra = interaction_function[ftype].nratoms;
356 nr = ilist[ftype].nr;
357 ia = ilist[ftype].iatoms;
359 for (i = 0; (i < nr); )
361 /* The vsite and first constructing atom */
364 wf[avsite] = adress_weight(x[avsite], adresstype, adressr, adressw, ref, pbc, fr);
371 else if (wf[ai] == 1)
380 /* Assign the vsite wf to rest of constructing atoms depending on type */
429 inc = 3*ip[ia[0]].vsiten.n;
430 for (j = 3; j < inc; j += 3)
437 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
438 ftype, __FILE__, __LINE__);
441 /* Increment loop variables */
450 update_adress_weights_atom(int cg0,
461 real adressr, adressw;
466 adresstype = fr->adress_type;
467 adressr = fr->adress_ex_width;
468 adressw = fr->adress_hy_width;
469 massT = mdatoms->massT;
471 ref = &(fr->adress_refs);
472 cgindex = cgs->index;
474 /* Only use first atom in charge group.
475 * We still can't be sure that the vsite and constructing
476 * atoms are on the same processor, so we must calculate
477 * in the same way as com adress. */
479 for (icg = cg0; (icg < cg1); icg++)
483 wf[k0] = adress_weight(x[k0], adresstype, adressr, adressw, ref, pbc, fr);
485 /* Set wf of all atoms in charge group equal to wf of first atom in charge group*/
486 for (k = (k0+1); (k < k1); k++)
494 adress_thermo_force(int start,
503 int iatom, n0, nnn, nrcg, i;
505 real adressw, adressr;
507 unsigned short * ptype;
513 real w, wsq, wmin1, wmin1sq, wp, wt, rinv, sqr_dl, dl;
514 real eps, eps2, F, Geps, Heps2, Fp, dmu_dwp, dwp_dr, fscal;
516 adresstype = fr->adress_type;
517 adressw = fr->adress_hy_width;
518 adressr = fr->adress_ex_width;
519 cgindex = cgs->index;
520 ptype = mdatoms->ptype;
521 ref = &(fr->adress_refs);
524 for (iatom = start; (iatom < start+homenr); iatom++)
526 if (egp_coarsegrained(fr, mdatoms->cENER[iatom]))
528 if (ptype[iatom] == eptVSite)
531 /* is it hybrid or apply the thermodynamics force everywhere?*/
532 if (mdatoms->tf_table_index[iatom] != NO_TF_TABLE)
534 if (fr->n_adress_tf_grps > 0)
536 /* multi component tf is on, select the right table */
537 ATFtab = fr->atf_tabs[mdatoms->tf_table_index[iatom]].data;
538 tabscale = fr->atf_tabs[mdatoms->tf_table_index[iatom]].scale;
542 /* just on component*/
543 ATFtab = fr->atf_tabs[DEFAULT_TF_TABLE].data;
544 tabscale = fr->atf_tabs[DEFAULT_TF_TABLE].scale;
550 pbc_dx(pbc, (*ref), x[iatom], dr);
554 rvec_sub((*ref), x[iatom], dr);
560 /* calculate distace to adress center again */
565 /* plane through center of ref, varies in x direction */
566 sqr_dl = dr[0]*dr[0];
567 rinv = gmx_invsqrt(dr[0]*dr[0]);
570 /* point at center of ref, assuming cubic geometry */
571 for (i = 0; i < 3; i++)
573 sqr_dl += dr[i]*dr[i];
575 rinv = gmx_invsqrt(sqr_dl);
578 /* This case should not happen */
583 /* table origin is adress center */
590 Geps = eps*ATFtab[nnn+2];
591 Heps2 = eps2*ATFtab[nnn+3];
593 F = (Fp+Geps+2.0*Heps2)*tabscale;
597 f[iatom][0] += fscal*dr[0];
598 if (adresstype != eAdressXSplit)
600 f[iatom][1] += fscal*dr[1];
601 f[iatom][2] += fscal*dr[2];
609 gmx_bool egp_explicit(t_forcerec * fr, int egp_nr)
611 return fr->adress_group_explicit[egp_nr];
614 gmx_bool egp_coarsegrained(t_forcerec * fr, int egp_nr)
616 return !fr->adress_group_explicit[egp_nr];