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,2016,2017, 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.
45 #include "gromacs/compat/make_unique.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/mdlib/gmx_omp_nthreads.h"
48 #include "gromacs/mdlib/qmmm.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/topology/mtop_lookup.h"
52 #include "gromacs/topology/mtop_util.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/alignedallocator.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/smalloc.h"
58 #define ALMOST_ZERO 1e-30
68 std::unique_ptr
<MDAtoms
>
69 makeMDAtoms(FILE *fp
, const gmx_mtop_t
&mtop
, const t_inputrec
&ir
)
71 auto mdAtoms
= compat::make_unique
<MDAtoms
>();
74 mdAtoms
->mdatoms_
.reset(md
);
76 md
->nenergrp
= mtop
.groups
.grps
[egcENER
].nr
;
77 md
->bVCMgrps
= (mtop
.groups
.grps
[egcVCM
].nr
> 1);
79 /* Determine the total system mass and perturbed atom counts */
80 double totalMassA
= 0.0;
81 double totalMassB
= 0.0;
83 gmx_mtop_atomloop_block_t aloop
= gmx_mtop_atomloop_block_init(&mtop
);
86 while (gmx_mtop_atomloop_block_next(aloop
, &atom
, &nmol
))
88 totalMassA
+= nmol
*atom
->m
;
89 totalMassB
+= nmol
*atom
->mB
;
91 if (ir
.efep
!= efepNO
&& PERTURBED(*atom
))
94 if (atom
->mB
!= atom
->m
)
96 md
->nMassPerturbed
+= nmol
;
98 if (atom
->qB
!= atom
->q
)
100 md
->nChargePerturbed
+= nmol
;
102 if (atom
->typeB
!= atom
->type
)
104 md
->nTypePerturbed
+= nmol
;
109 md
->tmassA
= totalMassA
;
110 md
->tmassB
= totalMassB
;
112 if (ir
.efep
!= efepNO
&& fp
)
115 "There are %d atoms and %d charges for free energy perturbation\n",
116 md
->nPerturbed
, md
->nChargePerturbed
);
119 md
->havePartiallyFrozenAtoms
= FALSE
;
120 for (int g
= 0; g
< ir
.opts
.ngfrz
; g
++)
122 for (int d
= YY
; d
< DIM
; d
++)
124 if (ir
.opts
.nFreeze
[d
] != ir
.opts
.nFreeze
[XX
])
126 md
->havePartiallyFrozenAtoms
= TRUE
;
131 md
->bOrires
= gmx_mtop_ftype_count(&mtop
, F_ORIRES
);
138 void atoms2md(const gmx_mtop_t
*mtop
, const t_inputrec
*ir
,
139 int nindex
, const int *index
,
141 gmx::MDAtoms
*mdAtoms
)
144 const t_grpopts
*opts
;
145 const gmx_groups_t
*groups
;
146 int nthreads gmx_unused
;
148 bLJPME
= EVDW_PME(ir
->vdwtype
);
152 groups
= &mtop
->groups
;
154 auto md
= mdAtoms
->mdatoms();
155 /* nindex>=0 indicates DD where we use an index */
162 md
->nr
= mtop
->natoms
;
165 if (md
->nr
> md
->nalloc
)
167 md
->nalloc
= over_alloc_dd(md
->nr
);
169 if (md
->nMassPerturbed
)
171 srenew(md
->massA
, md
->nalloc
);
172 srenew(md
->massB
, md
->nalloc
);
174 srenew(md
->massT
, md
->nalloc
);
175 /* The SIMD version of the integrator needs this aligned and padded.
176 * The padding needs to be with zeros, which we set later below.
178 gmx::AlignedAllocationPolicy::free(md
->invmass
);
179 md
->invmass
= new(gmx::AlignedAllocationPolicy::malloc((md
->nalloc
+ GMX_REAL_MAX_SIMD_WIDTH
)*sizeof(*md
->invmass
)))real
;
180 srenew(md
->invMassPerDim
, md
->nalloc
);
181 srenew(md
->chargeA
, md
->nalloc
);
182 srenew(md
->typeA
, md
->nalloc
);
185 srenew(md
->chargeB
, md
->nalloc
);
186 srenew(md
->typeB
, md
->nalloc
);
190 srenew(md
->sqrt_c6A
, md
->nalloc
);
191 srenew(md
->sigmaA
, md
->nalloc
);
192 srenew(md
->sigma3A
, md
->nalloc
);
195 srenew(md
->sqrt_c6B
, md
->nalloc
);
196 srenew(md
->sigmaB
, md
->nalloc
);
197 srenew(md
->sigma3B
, md
->nalloc
);
200 srenew(md
->ptype
, md
->nalloc
);
203 srenew(md
->cTC
, md
->nalloc
);
204 /* We always copy cTC with domain decomposition */
206 srenew(md
->cENER
, md
->nalloc
);
209 srenew(md
->cACC
, md
->nalloc
);
213 opts
->nFreeze
[0][XX
] || opts
->nFreeze
[0][YY
] || opts
->nFreeze
[0][ZZ
]))
215 srenew(md
->cFREEZE
, md
->nalloc
);
219 srenew(md
->cVCM
, md
->nalloc
);
223 srenew(md
->cORF
, md
->nalloc
);
227 srenew(md
->bPerturbed
, md
->nalloc
);
230 /* Note that these user t_mdatoms array pointers are NULL
231 * when there is only one group present.
232 * Therefore, when adding code, the user should use something like:
233 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
235 if (mtop
->groups
.grpnr
[egcUser1
] != nullptr)
237 srenew(md
->cU1
, md
->nalloc
);
239 if (mtop
->groups
.grpnr
[egcUser2
] != nullptr)
241 srenew(md
->cU2
, md
->nalloc
);
246 srenew(md
->bQM
, md
->nalloc
);
252 // cppcheck-suppress unreadVariable
253 nthreads
= gmx_omp_nthreads_get(emntDefault
);
254 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
255 for (int i
= 0; i
< md
->nr
; i
++)
263 if (index
== nullptr)
271 const t_atom
&atom
= mtopGetAtomParameters(mtop
, ag
, &molb
);
275 md
->cFREEZE
[i
] = ggrpnr(groups
, egcFREEZE
, ag
);
277 if (EI_ENERGY_MINIMIZATION(ir
->eI
))
279 /* Displacement is proportional to F, masses used for constraints */
283 else if (ir
->eI
== eiBD
)
285 /* With BD the physical masses are irrelevant.
286 * To keep the code simple we use most of the normal MD code path
287 * for BD. Thus for constraining the masses should be proportional
288 * to the friction coefficient. We set the absolute value such that
289 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
290 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
291 * correct kinetic energy and temperature using the usual code path.
292 * Thus with BD v*dt will give the displacement and the reported
293 * temperature can signal bad integration (too large time step).
297 mA
= 0.5*ir
->bd_fric
*ir
->delta_t
;
298 mB
= 0.5*ir
->bd_fric
*ir
->delta_t
;
302 /* The friction coefficient is mass/tau_t */
303 fac
= ir
->delta_t
/opts
->tau_t
[md
->cTC
? groups
->grpnr
[egcTC
][ag
] : 0];
305 mB
= 0.5*atom
.mB
*fac
;
313 if (md
->nMassPerturbed
)
323 md
->invMassPerDim
[i
][XX
] = 0;
324 md
->invMassPerDim
[i
][YY
] = 0;
325 md
->invMassPerDim
[i
][ZZ
] = 0;
327 else if (md
->cFREEZE
)
330 if (opts
->nFreeze
[g
][XX
] && opts
->nFreeze
[g
][YY
] && opts
->nFreeze
[g
][ZZ
])
332 /* Set the mass of completely frozen particles to ALMOST_ZERO
333 * iso 0 to avoid div by zero in lincs or shake.
335 md
->invmass
[i
] = ALMOST_ZERO
;
339 /* Note: Partially frozen particles use the normal invmass.
340 * If such particles are constrained, the frozen dimensions
341 * should not be updated with the constrained coordinates.
343 md
->invmass
[i
] = 1.0/mA
;
345 for (int d
= 0; d
< DIM
; d
++)
347 md
->invMassPerDim
[i
][d
] = (opts
->nFreeze
[g
][d
] ? 0 : 1.0/mA
);
352 md
->invmass
[i
] = 1.0/mA
;
353 for (int d
= 0; d
< DIM
; d
++)
355 md
->invMassPerDim
[i
][d
] = 1.0/mA
;
359 md
->chargeA
[i
] = atom
.q
;
360 md
->typeA
[i
] = atom
.type
;
363 c6
= mtop
->ffparams
.iparams
[atom
.type
*(mtop
->ffparams
.atnr
+1)].lj
.c6
;
364 c12
= mtop
->ffparams
.iparams
[atom
.type
*(mtop
->ffparams
.atnr
+1)].lj
.c12
;
365 md
->sqrt_c6A
[i
] = sqrt(c6
);
366 if (c6
== 0.0 || c12
== 0)
372 md
->sigmaA
[i
] = gmx::sixthroot(c12
/c6
);
374 md
->sigma3A
[i
] = 1/(md
->sigmaA
[i
]*md
->sigmaA
[i
]*md
->sigmaA
[i
]);
378 md
->bPerturbed
[i
] = PERTURBED(atom
);
379 md
->chargeB
[i
] = atom
.qB
;
380 md
->typeB
[i
] = atom
.typeB
;
383 c6
= mtop
->ffparams
.iparams
[atom
.typeB
*(mtop
->ffparams
.atnr
+1)].lj
.c6
;
384 c12
= mtop
->ffparams
.iparams
[atom
.typeB
*(mtop
->ffparams
.atnr
+1)].lj
.c12
;
385 md
->sqrt_c6B
[i
] = sqrt(c6
);
386 if (c6
== 0.0 || c12
== 0)
392 md
->sigmaB
[i
] = gmx::sixthroot(c12
/c6
);
394 md
->sigma3B
[i
] = 1/(md
->sigmaB
[i
]*md
->sigmaB
[i
]*md
->sigmaB
[i
]);
397 md
->ptype
[i
] = atom
.ptype
;
400 md
->cTC
[i
] = groups
->grpnr
[egcTC
][ag
];
402 md
->cENER
[i
] = ggrpnr(groups
, egcENER
, ag
);
405 md
->cACC
[i
] = groups
->grpnr
[egcACC
][ag
];
409 md
->cVCM
[i
] = groups
->grpnr
[egcVCM
][ag
];
413 md
->cORF
[i
] = ggrpnr(groups
, egcORFIT
, ag
);
418 md
->cU1
[i
] = groups
->grpnr
[egcUser1
][ag
];
422 md
->cU2
[i
] = groups
->grpnr
[egcUser2
][ag
];
427 if (groups
->grpnr
[egcQMMM
] == nullptr ||
428 groups
->grpnr
[egcQMMM
][ag
] < groups
->grps
[egcQMMM
].nr
-1)
438 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
443 /* Pad invmass with 0 so a SIMD MD update does not change v and x */
444 for (int i
= md
->nr
; i
< md
->nr
+ GMX_REAL_MAX_SIMD_WIDTH
; i
++)
451 /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
452 * For free-energy runs, these should be updated using update_mdatoms().
454 md
->tmass
= md
->tmassA
;
458 void update_mdatoms(t_mdatoms
*md
, real lambda
)
460 if (md
->nMassPerturbed
&& lambda
!= md
->lambda
)
462 real L1
= 1 - lambda
;
464 /* Update masses of perturbed atoms for the change in lambda */
465 // cppcheck-suppress unreadVariable
466 int gmx_unused nthreads
= gmx_omp_nthreads_get(emntDefault
);
467 #pragma omp parallel for num_threads(nthreads) schedule(static)
468 for (int i
= 0; i
< md
->nr
; i
++)
470 if (md
->bPerturbed
[i
])
472 md
->massT
[i
] = L1
*md
->massA
[i
] + lambda
*md
->massB
[i
];
473 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
474 * and their invmass does not depend on lambda.
476 if (md
->invmass
[i
] > 1.1*ALMOST_ZERO
)
478 md
->invmass
[i
] = 1.0/md
->massT
[i
];
479 for (int d
= 0; d
< DIM
; d
++)
481 if (md
->invMassPerDim
[i
][d
] > 1.1*ALMOST_ZERO
)
483 md
->invMassPerDim
[i
][d
] = md
->invmass
[i
];
490 /* Update the system mass for the change in lambda */
491 md
->tmass
= L1
*md
->tmassA
+ lambda
*md
->tmassB
;