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