1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
45 #include "mtop_util.h"
46 #include "gmx_omp_nthreads.h"
48 #define ALMOST_ZERO 1e-30
50 t_mdatoms *init_mdatoms(FILE *fp, gmx_mtop_t *mtop, gmx_bool bFreeEnergy)
56 gmx_mtop_atomloop_all_t aloop;
61 md->nenergrp = mtop->groups.grps[egcENER].nr;
66 aloop = gmx_mtop_atomloop_all_init(mtop);
67 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
69 if (ggrpnr(&mtop->groups, egcVCM, a) > 0)
74 if (bFreeEnergy && PERTURBED(*atom))
77 if (atom->mB != atom->m)
81 if (atom->qB != atom->q)
83 md->nChargePerturbed++;
94 if (bFreeEnergy && fp)
97 "There are %d atoms and %d charges for free energy perturbation\n",
98 md->nPerturbed, md->nChargePerturbed);
101 md->bOrires = gmx_mtop_ftype_count(mtop, F_ORIRES);
106 void atoms2md(gmx_mtop_t *mtop, t_inputrec *ir,
107 int nindex, int *index,
108 int start, int homenr,
111 gmx_mtop_atomlookup_t alook;
114 gmx_groups_t *groups;
115 gmx_molblock_t *molblock;
119 groups = &mtop->groups;
121 molblock = mtop->molblock;
123 /* Index==NULL indicates particle decomposition,
124 * unless we have an empty DD node, so also check for homenr and start.
125 * This should be signaled properly with an extra parameter or nindex==-1.
127 if (index == NULL && (homenr > 0 || start > 0))
129 md->nr = mtop->natoms;
136 if (md->nr > md->nalloc)
138 md->nalloc = over_alloc_dd(md->nr);
140 if (md->nMassPerturbed)
142 srenew(md->massA, md->nalloc);
143 srenew(md->massB, md->nalloc);
145 srenew(md->massT, md->nalloc);
146 srenew(md->invmass, md->nalloc);
147 srenew(md->chargeA, md->nalloc);
150 srenew(md->chargeB, md->nalloc);
152 srenew(md->typeA, md->nalloc);
155 srenew(md->typeB, md->nalloc);
157 srenew(md->ptype, md->nalloc);
160 srenew(md->cTC, md->nalloc);
161 /* We always copy cTC with domain decomposition */
163 srenew(md->cENER, md->nalloc);
166 srenew(md->cACC, md->nalloc);
170 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
172 srenew(md->cFREEZE, md->nalloc);
176 srenew(md->cVCM, md->nalloc);
180 srenew(md->cORF, md->nalloc);
184 srenew(md->bPerturbed, md->nalloc);
187 /* Note that these user t_mdatoms array pointers are NULL
188 * when there is only one group present.
189 * Therefore, when adding code, the user should use something like:
190 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
192 if (mtop->groups.grpnr[egcUser1] != NULL)
194 srenew(md->cU1, md->nalloc);
196 if (mtop->groups.grpnr[egcUser2] != NULL)
198 srenew(md->cU2, md->nalloc);
203 srenew(md->bQM, md->nalloc);
207 srenew(md->wf, md->nalloc);
208 srenew(md->tf_table_index, md->nalloc);
212 alook = gmx_mtop_atomlookup_init(mtop);
214 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntDefault)) schedule(static)
215 for (i = 0; i < md->nr; i++)
229 gmx_mtop_atomnr_to_atom(alook, ag, &atom);
233 md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
235 if (EI_ENERGY_MINIMIZATION(ir->eI))
237 /* Displacement is proportional to F, masses used for constraints */
241 else if (ir->eI == eiBD)
243 /* With BD the physical masses are irrelevant.
244 * To keep the code simple we use most of the normal MD code path
245 * for BD. Thus for constraining the masses should be proportional
246 * to the friction coefficient. We set the absolute value such that
247 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
248 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
249 * correct kinetic energy and temperature using the usual code path.
250 * Thus with BD v*dt will give the displacement and the reported
251 * temperature can signal bad integration (too large time step).
255 mA = 0.5*ir->bd_fric*ir->delta_t;
256 mB = 0.5*ir->bd_fric*ir->delta_t;
260 /* The friction coefficient is mass/tau_t */
261 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
262 mA = 0.5*atom->m*fac;
263 mB = 0.5*atom->mB*fac;
271 if (md->nMassPerturbed)
281 else if (md->cFREEZE)
284 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
286 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
287 * to avoid div by zero in lincs or shake.
288 * Note that constraints can still move a partially frozen particle.
290 md->invmass[i] = ALMOST_ZERO;
294 md->invmass[i] = 1.0/mA;
299 md->invmass[i] = 1.0/mA;
301 md->chargeA[i] = atom->q;
302 md->typeA[i] = atom->type;
305 md->chargeB[i] = atom->qB;
306 md->typeB[i] = atom->typeB;
307 md->bPerturbed[i] = PERTURBED(*atom);
309 md->ptype[i] = atom->ptype;
312 md->cTC[i] = groups->grpnr[egcTC][ag];
315 (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
318 md->cACC[i] = groups->grpnr[egcACC][ag];
322 md->cVCM[i] = groups->grpnr[egcVCM][ag];
326 md->cORF[i] = groups->grpnr[egcORFIT][ag];
331 md->cU1[i] = groups->grpnr[egcUser1][ag];
335 md->cU2[i] = groups->grpnr[egcUser2][ag];
340 if (groups->grpnr[egcQMMM] == 0 ||
341 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
350 /* Initialize AdResS weighting functions to adressw */
354 /* if no tf table groups specified, use default table */
355 md->tf_table_index[i] = DEFAULT_TF_TABLE;
356 if (ir->adress->n_tf_grps > 0)
358 /* if tf table groups specified, tf is only applied to thoose energy groups*/
359 md->tf_table_index[i] = NO_TF_TABLE;
360 /* check wether atom is in one of the relevant energy groups and assign a table index */
361 for (g = 0; g < ir->adress->n_tf_grps; g++)
363 if (md->cENER[i] == ir->adress->tf_table_index[g])
365 md->tf_table_index[i] = g;
372 gmx_mtop_atomlookup_destroy(alook);
379 void update_mdatoms(t_mdatoms *md, real lambda)
382 real L1 = 1.0-lambda;
386 if (md->nMassPerturbed)
388 for (al = 0; (al < end); al++)
390 if (md->bPerturbed[al])
392 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
393 if (md->invmass[al] > 1.1*ALMOST_ZERO)
395 md->invmass[al] = 1.0/md->massT[al];
399 md->tmass = L1*md->tmassA + lambda*md->tmassB;
403 md->tmass = md->tmassA;