Refactor virtual site interface
[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, 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/ewald/pme.h"
42 #include "gromacs/listed_forces/manage_threading.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"
56
57 namespace gmx
58 {
59
60 /* TODO: Add a routine that collects the initial setup of the algorithms.
61  *
62  * The final solution should be an MD algorithm base class with methods
63  * for initialization and atom-data setup.
64  */
65 void mdAlgorithmsSetupAtomData(const t_commrec*        cr,
66                                const t_inputrec*       ir,
67                                const gmx_mtop_t&       top_global,
68                                gmx_localtop_t*         top,
69                                t_forcerec*             fr,
70                                PaddedHostVector<RVec>* force,
71                                MDAtoms*                mdAtoms,
72                                Constraints*            constr,
73                                VirtualSitesHandler*    vsite,
74                                gmx_shellfc_t*          shellfc)
75 {
76     bool usingDomDec = DOMAINDECOMP(cr);
77
78     int  numAtomIndex;
79     int* atomIndex;
80     int  numHomeAtoms;
81     int  numTotalAtoms;
82
83     if (usingDomDec)
84     {
85         numAtomIndex  = dd_natoms_mdatoms(cr->dd);
86         atomIndex     = cr->dd->globalAtomIndices.data();
87         numHomeAtoms  = dd_numHomeAtoms(*cr->dd);
88         numTotalAtoms = dd_natoms_mdatoms(cr->dd);
89     }
90     else
91     {
92         numAtomIndex  = -1;
93         atomIndex     = nullptr;
94         numHomeAtoms  = top_global.natoms;
95         numTotalAtoms = top_global.natoms;
96     }
97
98     if (force != nullptr)
99     {
100         /* We need to allocate one element extra, since we might use
101          * (unaligned) 4-wide SIMD loads to access rvec entries.
102          */
103         force->resizeWithPadding(numTotalAtoms);
104     }
105
106     atoms2md(&top_global, ir, numAtomIndex, atomIndex, numHomeAtoms, mdAtoms);
107
108     auto mdatoms = mdAtoms->mdatoms();
109     if (usingDomDec)
110     {
111         dd_sort_local_top(cr->dd, mdatoms, top);
112     }
113     else
114     {
115         gmx_mtop_generate_local_top(top_global, top, ir->efep != efepNO);
116     }
117
118     if (vsite)
119     {
120         vsite->setVirtualSites(top->idef.il, *mdatoms);
121     }
122
123     /* Note that with DD only flexible constraints, not shells, are supported
124      * and these don't require setup in make_local_shells().
125      *
126      * TODO: This should only happen in ShellFCElement (it is called directly by the modular
127      *       simulator ShellFCElement already, but still used here by legacy simulators)
128      */
129     if (!usingDomDec && shellfc)
130     {
131         make_local_shells(cr, mdatoms, shellfc);
132     }
133
134     setup_bonded_threading(fr->bondedThreading, fr->natoms_force, fr->gpuBonded != nullptr, top->idef);
135
136     if (EEL_PME(fr->ic->eeltype) && (cr->duty & DUTY_PME))
137     {
138         /* This handles the PP+PME rank case where fr->pmedata is valid.
139          * For PME-only ranks, gmx_pmeonly() has its own call to gmx_pme_reinit_atoms().
140          */
141         const int numPmeAtoms = numHomeAtoms - fr->n_tpi;
142         gmx_pme_reinit_atoms(fr->pmedata, numPmeAtoms, mdatoms->chargeA);
143     }
144
145     if (constr)
146     {
147         constr->setConstraints(top, *mdatoms);
148     }
149 }
150
151 } // namespace gmx