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) 2012,2013, 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.
47 #include "mtop_util.h"
48 #include "gmx_omp_nthreads.h"
50 #define ALMOST_ZERO 1e-30
52 t_mdatoms *init_mdatoms(FILE *fp, gmx_mtop_t *mtop, gmx_bool bFreeEnergy)
58 gmx_mtop_atomloop_all_t aloop;
63 md->nenergrp = mtop->groups.grps[egcENER].nr;
68 aloop = gmx_mtop_atomloop_all_init(mtop);
69 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
71 if (ggrpnr(&mtop->groups, egcVCM, a) > 0)
76 if (bFreeEnergy && PERTURBED(*atom))
79 if (atom->mB != atom->m)
83 if (atom->qB != atom->q)
85 md->nChargePerturbed++;
87 if (atom->typeB != atom->type)
100 if (bFreeEnergy && fp)
103 "There are %d atoms and %d charges for free energy perturbation\n",
104 md->nPerturbed, md->nChargePerturbed);
107 md->bOrires = gmx_mtop_ftype_count(mtop, F_ORIRES);
112 void atoms2md(gmx_mtop_t *mtop, t_inputrec *ir,
113 int nindex, int *index,
114 int start, int homenr,
117 gmx_mtop_atomlookup_t alook;
121 gmx_groups_t *groups;
122 gmx_molblock_t *molblock;
126 groups = &mtop->groups;
128 molblock = mtop->molblock;
130 /* Index==NULL indicates particle decomposition,
131 * unless we have an empty DD node, so also check for homenr and start.
132 * This should be signaled properly with an extra parameter or nindex==-1.
134 if (index == NULL && (homenr > 0 || start > 0))
136 md->nr = mtop->natoms;
143 if (md->nr > md->nalloc)
145 md->nalloc = over_alloc_dd(md->nr);
147 if (md->nMassPerturbed)
149 srenew(md->massA, md->nalloc);
150 srenew(md->massB, md->nalloc);
152 srenew(md->massT, md->nalloc);
153 srenew(md->invmass, md->nalloc);
154 srenew(md->chargeA, md->nalloc);
155 srenew(md->c6A, md->nalloc);
156 srenew(md->sigmaA, md->nalloc);
157 srenew(md->sigma3A, md->nalloc);
160 srenew(md->chargeB, md->nalloc);
161 srenew(md->c6B, md->nalloc);
162 srenew(md->sigmaB, md->nalloc);
163 srenew(md->sigma3B, md->nalloc);
165 srenew(md->typeA, md->nalloc);
168 srenew(md->typeB, md->nalloc);
170 srenew(md->ptype, md->nalloc);
173 srenew(md->cTC, md->nalloc);
174 /* We always copy cTC with domain decomposition */
176 srenew(md->cENER, md->nalloc);
179 srenew(md->cACC, md->nalloc);
183 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
185 srenew(md->cFREEZE, md->nalloc);
189 srenew(md->cVCM, md->nalloc);
193 srenew(md->cORF, md->nalloc);
197 srenew(md->bPerturbed, md->nalloc);
200 /* Note that these user t_mdatoms array pointers are NULL
201 * when there is only one group present.
202 * Therefore, when adding code, the user should use something like:
203 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
205 if (mtop->groups.grpnr[egcUser1] != NULL)
207 srenew(md->cU1, md->nalloc);
209 if (mtop->groups.grpnr[egcUser2] != NULL)
211 srenew(md->cU2, md->nalloc);
216 srenew(md->bQM, md->nalloc);
220 srenew(md->wf, md->nalloc);
221 srenew(md->tf_table_index, md->nalloc);
225 alook = gmx_mtop_atomlookup_init(mtop);
227 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntDefault)) schedule(static)
228 for (i = 0; i < md->nr; i++)
242 gmx_mtop_atomnr_to_atom(alook, ag, &atom);
246 md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
248 if (EI_ENERGY_MINIMIZATION(ir->eI))
250 /* Displacement is proportional to F, masses used for constraints */
254 else if (ir->eI == eiBD)
256 /* With BD the physical masses are irrelevant.
257 * To keep the code simple we use most of the normal MD code path
258 * for BD. Thus for constraining the masses should be proportional
259 * to the friction coefficient. We set the absolute value such that
260 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
261 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
262 * correct kinetic energy and temperature using the usual code path.
263 * Thus with BD v*dt will give the displacement and the reported
264 * temperature can signal bad integration (too large time step).
268 mA = 0.5*ir->bd_fric*ir->delta_t;
269 mB = 0.5*ir->bd_fric*ir->delta_t;
273 /* The friction coefficient is mass/tau_t */
274 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
275 mA = 0.5*atom->m*fac;
276 mB = 0.5*atom->mB*fac;
284 if (md->nMassPerturbed)
294 else if (md->cFREEZE)
297 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
299 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
300 * to avoid div by zero in lincs or shake.
301 * Note that constraints can still move a partially frozen particle.
303 md->invmass[i] = ALMOST_ZERO;
307 md->invmass[i] = 1.0/mA;
312 md->invmass[i] = 1.0/mA;
314 md->chargeA[i] = atom->q;
315 md->typeA[i] = atom->type;
316 c6 = mtop->ffparams.iparams[atom->type*(mtop->ffparams.atnr+1)].lj.c6;
317 c12 = mtop->ffparams.iparams[atom->type*(mtop->ffparams.atnr+1)].lj.c12;
318 md->c6A[i] = sqrt(c6);
319 if (c6 == 0.0 || c12 == 0.0)
325 md->sigmaA[i] = pow(c12/c6, 1.0/6.0);
327 md->sigma3A[i] = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
331 md->chargeB[i] = atom->qB;
332 md->typeB[i] = atom->typeB;
333 c6 = mtop->ffparams.iparams[atom->typeB*(mtop->ffparams.atnr+1)].lj.c6;
334 c12 = mtop->ffparams.iparams[atom->typeB*(mtop->ffparams.atnr+1)].lj.c12;
335 md->c6B[i] = sqrt(c6);
336 if (c6 == 0.0 || c12 == 0.0)
342 md->sigmaB[i] = pow(c12/c6, 1.0/6.0);
344 md->bPerturbed[i] = PERTURBED(*atom);
345 md->sigma3B[i] = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
347 md->ptype[i] = atom->ptype;
350 md->cTC[i] = groups->grpnr[egcTC][ag];
353 (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
356 md->cACC[i] = groups->grpnr[egcACC][ag];
360 md->cVCM[i] = groups->grpnr[egcVCM][ag];
364 md->cORF[i] = groups->grpnr[egcORFIT][ag];
369 md->cU1[i] = groups->grpnr[egcUser1][ag];
373 md->cU2[i] = groups->grpnr[egcUser2][ag];
378 if (groups->grpnr[egcQMMM] == 0 ||
379 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
388 /* Initialize AdResS weighting functions to adressw */
392 /* if no tf table groups specified, use default table */
393 md->tf_table_index[i] = DEFAULT_TF_TABLE;
394 if (ir->adress->n_tf_grps > 0)
396 /* if tf table groups specified, tf is only applied to thoose energy groups*/
397 md->tf_table_index[i] = NO_TF_TABLE;
398 /* check wether atom is in one of the relevant energy groups and assign a table index */
399 for (g = 0; g < ir->adress->n_tf_grps; g++)
401 if (md->cENER[i] == ir->adress->tf_table_index[g])
403 md->tf_table_index[i] = g;
410 gmx_mtop_atomlookup_destroy(alook);
417 void update_mdatoms(t_mdatoms *md, real lambda)
420 real L1 = 1.0-lambda;
424 if (md->nMassPerturbed)
426 for (al = 0; (al < end); al++)
428 if (md->bPerturbed[al])
430 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
431 if (md->invmass[al] > 1.1*ALMOST_ZERO)
433 md->invmass[al] = 1.0/md->massT[al];
437 md->tmass = L1*md->tmassA + lambda*md->tmassB;
441 md->tmass = md->tmassA;