Moved mdatom.h from legacyheader/types to mdtypes.
[gromacs.git] / src / gromacs / mdlib / mdatoms.cpp
blobf2992cf46c67669a2768da003fce7fbe24e44c98
1 /*
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,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.
37 #include "gmxpre.h"
39 #include "mdatoms.h"
41 #include <math.h>
43 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
44 #include "gromacs/mdlib/qmmm.h"
45 #include "gromacs/topology/mtop_util.h"
46 #include "gromacs/topology/topology.h"
47 #include "gromacs/utility/exceptions.h"
48 #include "gromacs/utility/smalloc.h"
50 #define ALMOST_ZERO 1e-30
52 t_mdatoms *init_mdatoms(FILE *fp, const gmx_mtop_t *mtop, gmx_bool bFreeEnergy)
54 int a;
55 double tmA, tmB;
56 t_atom *atom;
57 t_mdatoms *md;
58 gmx_mtop_atomloop_all_t aloop;
60 snew(md, 1);
62 md->nenergrp = mtop->groups.grps[egcENER].nr;
63 md->bVCMgrps = FALSE;
64 tmA = 0.0;
65 tmB = 0.0;
67 aloop = gmx_mtop_atomloop_all_init(mtop);
68 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
70 if (ggrpnr(&mtop->groups, egcVCM, a) > 0)
72 md->bVCMgrps = TRUE;
75 if (bFreeEnergy && PERTURBED(*atom))
77 md->nPerturbed++;
78 if (atom->mB != atom->m)
80 md->nMassPerturbed++;
82 if (atom->qB != atom->q)
84 md->nChargePerturbed++;
86 if (atom->typeB != atom->type)
88 md->nTypePerturbed++;
92 tmA += atom->m;
93 tmB += atom->mB;
96 md->tmassA = tmA;
97 md->tmassB = tmB;
99 if (bFreeEnergy && fp)
101 fprintf(fp,
102 "There are %d atoms and %d charges for free energy perturbation\n",
103 md->nPerturbed, md->nChargePerturbed);
106 md->bOrires = gmx_mtop_ftype_count(mtop, F_ORIRES);
108 return md;
111 void atoms2md(const gmx_mtop_t *mtop, const t_inputrec *ir,
112 int nindex, const int *index,
113 int homenr,
114 t_mdatoms *md)
116 gmx_bool bLJPME;
117 gmx_mtop_atomlookup_t alook;
118 int i;
119 const t_grpopts *opts;
120 const gmx_groups_t *groups;
121 int nthreads gmx_unused;
122 const real oneOverSix = 1.0 / 6.0;
124 bLJPME = EVDW_PME(ir->vdwtype);
126 opts = &ir->opts;
128 groups = &mtop->groups;
130 /* Index==NULL indicates no DD (unless we have a DD node with no
131 * atoms), so also check for homenr. This should be
132 * signaled properly with an extra parameter or nindex==-1.
134 if (index == NULL && (homenr > 0))
136 md->nr = mtop->natoms;
138 else
140 md->nr = nindex;
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->typeA, md->nalloc);
156 if (md->nPerturbed)
158 srenew(md->chargeB, md->nalloc);
159 srenew(md->typeB, md->nalloc);
161 if (bLJPME)
163 srenew(md->sqrt_c6A, md->nalloc);
164 srenew(md->sigmaA, md->nalloc);
165 srenew(md->sigma3A, md->nalloc);
166 if (md->nPerturbed)
168 srenew(md->sqrt_c6B, md->nalloc);
169 srenew(md->sigmaB, md->nalloc);
170 srenew(md->sigma3B, md->nalloc);
173 srenew(md->ptype, md->nalloc);
174 if (opts->ngtc > 1)
176 srenew(md->cTC, md->nalloc);
177 /* We always copy cTC with domain decomposition */
179 srenew(md->cENER, md->nalloc);
180 if (opts->ngacc > 1)
182 srenew(md->cACC, md->nalloc);
184 if (opts->nFreeze &&
185 (opts->ngfrz > 1 ||
186 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
188 srenew(md->cFREEZE, md->nalloc);
190 if (md->bVCMgrps)
192 srenew(md->cVCM, md->nalloc);
194 if (md->bOrires)
196 srenew(md->cORF, md->nalloc);
198 if (md->nPerturbed)
200 srenew(md->bPerturbed, md->nalloc);
203 /* Note that these user t_mdatoms array pointers are NULL
204 * when there is only one group present.
205 * Therefore, when adding code, the user should use something like:
206 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
208 if (mtop->groups.grpnr[egcUser1] != NULL)
210 srenew(md->cU1, md->nalloc);
212 if (mtop->groups.grpnr[egcUser2] != NULL)
214 srenew(md->cU2, md->nalloc);
217 if (ir->bQMMM)
219 srenew(md->bQM, md->nalloc);
223 alook = gmx_mtop_atomlookup_init(mtop);
225 // cppcheck-suppress unreadVariable
226 nthreads = gmx_omp_nthreads_get(emntDefault);
227 #pragma omp parallel for num_threads(nthreads) schedule(static)
228 for (i = 0; i < md->nr; i++)
232 int g, ag;
233 real mA, mB, fac;
234 real c6, c12;
235 t_atom *atom;
237 if (index == NULL)
239 ag = i;
241 else
243 ag = index[i];
245 gmx_mtop_atomnr_to_atom(alook, ag, &atom);
247 if (md->cFREEZE)
249 md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
251 if (EI_ENERGY_MINIMIZATION(ir->eI))
253 /* Displacement is proportional to F, masses used for constraints */
254 mA = 1.0;
255 mB = 1.0;
257 else if (ir->eI == eiBD)
259 /* With BD the physical masses are irrelevant.
260 * To keep the code simple we use most of the normal MD code path
261 * for BD. Thus for constraining the masses should be proportional
262 * to the friction coefficient. We set the absolute value such that
263 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
264 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
265 * correct kinetic energy and temperature using the usual code path.
266 * Thus with BD v*dt will give the displacement and the reported
267 * temperature can signal bad integration (too large time step).
269 if (ir->bd_fric > 0)
271 mA = 0.5*ir->bd_fric*ir->delta_t;
272 mB = 0.5*ir->bd_fric*ir->delta_t;
274 else
276 /* The friction coefficient is mass/tau_t */
277 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
278 mA = 0.5*atom->m*fac;
279 mB = 0.5*atom->mB*fac;
282 else
284 mA = atom->m;
285 mB = atom->mB;
287 if (md->nMassPerturbed)
289 md->massA[i] = mA;
290 md->massB[i] = mB;
292 md->massT[i] = mA;
293 if (mA == 0.0)
295 md->invmass[i] = 0;
297 else if (md->cFREEZE)
299 g = md->cFREEZE[i];
300 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
302 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
303 * to avoid div by zero in lincs or shake.
304 * Note that constraints can still move a partially frozen particle.
306 md->invmass[i] = ALMOST_ZERO;
308 else
310 md->invmass[i] = 1.0/mA;
313 else
315 md->invmass[i] = 1.0/mA;
317 md->chargeA[i] = atom->q;
318 md->typeA[i] = atom->type;
319 if (bLJPME)
321 c6 = mtop->ffparams.iparams[atom->type*(mtop->ffparams.atnr+1)].lj.c6;
322 c12 = mtop->ffparams.iparams[atom->type*(mtop->ffparams.atnr+1)].lj.c12;
323 md->sqrt_c6A[i] = sqrt(c6);
324 if (c6 == 0.0 || c12 == 0)
326 md->sigmaA[i] = 1.0;
328 else
330 md->sigmaA[i] = pow(c12/c6, oneOverSix);
332 md->sigma3A[i] = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
334 if (md->nPerturbed)
336 md->bPerturbed[i] = PERTURBED(*atom);
337 md->chargeB[i] = atom->qB;
338 md->typeB[i] = atom->typeB;
339 if (bLJPME)
341 c6 = mtop->ffparams.iparams[atom->typeB*(mtop->ffparams.atnr+1)].lj.c6;
342 c12 = mtop->ffparams.iparams[atom->typeB*(mtop->ffparams.atnr+1)].lj.c12;
343 md->sqrt_c6B[i] = sqrt(c6);
344 if (c6 == 0.0 || c12 == 0)
346 md->sigmaB[i] = 1.0;
348 else
350 md->sigmaB[i] = pow(c12/c6, oneOverSix);
352 md->sigma3B[i] = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
355 md->ptype[i] = atom->ptype;
356 if (md->cTC)
358 md->cTC[i] = groups->grpnr[egcTC][ag];
360 md->cENER[i] =
361 (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
362 if (md->cACC)
364 md->cACC[i] = groups->grpnr[egcACC][ag];
366 if (md->cVCM)
368 md->cVCM[i] = groups->grpnr[egcVCM][ag];
370 if (md->cORF)
372 md->cORF[i] = groups->grpnr[egcORFIT][ag];
375 if (md->cU1)
377 md->cU1[i] = groups->grpnr[egcUser1][ag];
379 if (md->cU2)
381 md->cU2[i] = groups->grpnr[egcUser2][ag];
384 if (ir->bQMMM)
386 if (groups->grpnr[egcQMMM] == 0 ||
387 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
389 md->bQM[i] = TRUE;
391 else
393 md->bQM[i] = FALSE;
397 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
400 gmx_mtop_atomlookup_destroy(alook);
402 md->homenr = homenr;
403 md->lambda = 0;
406 void update_mdatoms(t_mdatoms *md, real lambda)
408 int al, end;
409 real L1 = 1.0-lambda;
411 end = md->nr;
413 if (md->nMassPerturbed)
415 for (al = 0; (al < end); al++)
417 if (md->bPerturbed[al])
419 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
420 if (md->invmass[al] > 1.1*ALMOST_ZERO)
422 md->invmass[al] = 1.0/md->massT[al];
426 md->tmass = L1*md->tmassA + lambda*md->tmassB;
428 else
430 md->tmass = md->tmassA;
432 md->lambda = lambda;