2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019,2020, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "gromacs/domdec/domdec.h"
40 #include "gromacs/domdec/domdec_struct.h"
41 #include "gromacs/ewald/pme.h"
42 #include "gromacs/listed_forces/listed_forces.h"
43 #include "gromacs/mdlib/constr.h"
44 #include "gromacs/mdlib/mdatoms.h"
45 #include "gromacs/mdlib/vsite.h"
46 #include "gromacs/mdtypes/commrec.h"
47 #include "gromacs/mdtypes/forcerec.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/interaction_const.h"
50 #include "gromacs/mdtypes/mdatom.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/topology/mtop_util.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/smalloc.h"
60 /* TODO: Add a routine that collects the initial setup of the algorithms.
62 * The final solution should be an MD algorithm base class with methods
63 * for initialization and atom-data setup.
65 void mdAlgorithmsSetupAtomData(const t_commrec* cr,
67 const gmx_mtop_t& top_global,
70 PaddedHostVector<RVec>* force,
73 VirtualSitesHandler* vsite,
74 gmx_shellfc_t* shellfc)
76 bool usingDomDec = DOMAINDECOMP(cr);
84 numAtomIndex = dd_natoms_mdatoms(cr->dd);
85 numHomeAtoms = dd_numHomeAtoms(*cr->dd);
86 numTotalAtoms = dd_natoms_mdatoms(cr->dd);
91 numHomeAtoms = top_global.natoms;
92 numTotalAtoms = top_global.natoms;
97 /* We need to allocate one element extra, since we might use
98 * (unaligned) 4-wide SIMD loads to access rvec entries.
100 force->resizeWithPadding(numTotalAtoms);
103 atoms2md(&top_global, ir, numAtomIndex,
104 usingDomDec ? cr->dd->globalAtomIndices : std::vector<int>(), numHomeAtoms, mdAtoms);
106 auto mdatoms = mdAtoms->mdatoms();
109 dd_sort_local_top(cr->dd, mdatoms, top);
113 gmx_mtop_generate_local_top(top_global, top, ir->efep != efepNO);
118 vsite->setVirtualSites(top->idef.il, *mdatoms);
121 /* Note that with DD only flexible constraints, not shells, are supported
122 * and these don't require setup in make_local_shells().
124 * TODO: This should only happen in ShellFCElement (it is called directly by the modular
125 * simulator ShellFCElement already, but still used here by legacy simulators)
127 if (!usingDomDec && shellfc)
129 make_local_shells(cr, mdatoms, shellfc);
132 fr->listedForces->setup(top->idef, fr->natoms_force, fr->gpuBonded != nullptr);
134 if (EEL_PME(fr->ic->eeltype) && (cr->duty & DUTY_PME))
136 /* This handles the PP+PME rank case where fr->pmedata is valid.
137 * For PME-only ranks, gmx_pmeonly() has its own call to gmx_pme_reinit_atoms().
139 const int numPmeAtoms = numHomeAtoms - fr->n_tpi;
140 gmx_pme_reinit_atoms(fr->pmedata, numPmeAtoms, mdatoms->chargeA);
145 constr->setConstraints(top, mdatoms->nr, mdatoms->homenr, mdatoms->massT, mdatoms->invmass,
146 mdatoms->nMassPerturbed != 0, mdatoms->lambda, mdatoms->cFREEZE);