631d167227c536d881cb7ddeae714296a9337830
[alexxy/gromacs.git] / src / gromacs / domdec / mdsetup.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "mdsetup.h"
38
39 #include "gromacs/domdec/domdec.h"
40 #include "gromacs/domdec/domdec_struct.h"
41 #include "gromacs/domdec/localtopology.h"
42 #include "gromacs/ewald/pme.h"
43 #include "gromacs/listed_forces/listed_forces.h"
44 #include "gromacs/mdlib/constr.h"
45 #include "gromacs/mdlib/mdatoms.h"
46 #include "gromacs/mdlib/vsite.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/forcebuffers.h"
49 #include "gromacs/mdtypes/forcerec.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/interaction_const.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/topology/mtop_util.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/arrayref.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/smalloc.h"
59
60 namespace gmx
61 {
62
63 /* TODO: Add a routine that collects the initial setup of the algorithms.
64  *
65  * The final solution should be an MD algorithm base class with methods
66  * for initialization and atom-data setup.
67  */
68 void mdAlgorithmsSetupAtomData(const t_commrec*     cr,
69                                const t_inputrec&    inputrec,
70                                const gmx_mtop_t&    top_global,
71                                gmx_localtop_t*      top,
72                                t_forcerec*          fr,
73                                ForceBuffers*        force,
74                                MDAtoms*             mdAtoms,
75                                Constraints*         constr,
76                                VirtualSitesHandler* vsite,
77                                gmx_shellfc_t*       shellfc)
78 {
79     bool usingDomDec = DOMAINDECOMP(cr);
80
81     int numAtomIndex;
82     int numHomeAtoms;
83     int numTotalAtoms;
84
85     if (usingDomDec)
86     {
87         numAtomIndex  = dd_natoms_mdatoms(*cr->dd);
88         numHomeAtoms  = dd_numHomeAtoms(*cr->dd);
89         numTotalAtoms = dd_natoms_mdatoms(*cr->dd);
90     }
91     else
92     {
93         numAtomIndex  = -1;
94         numHomeAtoms  = top_global.natoms;
95         numTotalAtoms = top_global.natoms;
96     }
97
98     if (force != nullptr)
99     {
100         force->resize(numTotalAtoms);
101     }
102
103     atoms2md(top_global,
104              inputrec,
105              numAtomIndex,
106              usingDomDec ? cr->dd->globalAtomIndices : std::vector<int>(),
107              numHomeAtoms,
108              mdAtoms);
109
110     t_mdatoms* mdatoms = mdAtoms->mdatoms();
111     if (usingDomDec)
112     {
113         dd_sort_local_top(*cr->dd, mdatoms, top);
114     }
115     else
116     {
117         gmx_mtop_generate_local_top(top_global, top, inputrec.efep != FreeEnergyPerturbationType::No);
118     }
119
120     if (vsite)
121     {
122         vsite->setVirtualSites(top->idef.il,
123                                mdatoms->nr,
124                                mdatoms->homenr,
125                                gmx::arrayRefFromArray(mdatoms->ptype, mdatoms->nr));
126     }
127
128     /* Note that with DD only flexible constraints, not shells, are supported
129      * and these don't require setup in make_local_shells().
130      *
131      * TODO: This should only happen in ShellFCElement (it is called directly by the modular
132      *       simulator ShellFCElement already, but still used here by legacy simulators)
133      */
134     if (!usingDomDec && shellfc)
135     {
136         make_local_shells(cr, *mdatoms, shellfc);
137     }
138
139     for (auto& listedForces : fr->listedForces)
140     {
141         listedForces.setup(top->idef, fr->natoms_force, fr->listedForcesGpu != nullptr);
142     }
143
144     if (EEL_PME(fr->ic->eeltype) && (cr->duty & DUTY_PME))
145     {
146         /* This handles the PP+PME rank case where fr->pmedata is valid.
147          * For PME-only ranks, gmx_pmeonly() has its own call to gmx_pme_reinit_atoms().
148          */
149         const int numPmeAtoms = numHomeAtoms - fr->n_tpi;
150         gmx_pme_reinit_atoms(fr->pmedata,
151                              numPmeAtoms,
152                              mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
153                                               : gmx::ArrayRef<real>{},
154                              mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
155                                               : gmx::ArrayRef<real>{});
156     }
157
158     if (constr)
159     {
160         constr->setConstraints(top,
161                                mdatoms->nr,
162                                mdatoms->homenr,
163                                gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
164                                gmx::arrayRefFromArray(mdatoms->invmass, mdatoms->nr),
165                                mdatoms->nMassPerturbed != 0,
166                                mdatoms->lambda,
167                                mdatoms->cFREEZE ? gmx::arrayRefFromArray(mdatoms->cFREEZE, mdatoms->nr)
168                                                 : gmx::ArrayRef<const unsigned short>());
169     }
170 }
171
172 } // namespace gmx