Move responsibility for bonded threading decomposition
[alexxy/gromacs.git] / src / gromacs / mdlib / mdsetup.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018, 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/mdatoms.h"
44 #include "gromacs/mdlib/shellfc.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/pbcutil/mshift.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/topology/mtop_util.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/smalloc.h"
55
56 /* TODO: Add a routine that collects the initial setup of the algorithms.
57  *
58  * The final solution should be an MD algorithm base class with methods
59  * for initialization and atom-data setup.
60  */
61
62 void mdAlgorithmsSetupAtomData(const t_commrec   *cr,
63                                const t_inputrec  *ir,
64                                const gmx_mtop_t  *top_global,
65                                gmx_localtop_t    *top,
66                                t_forcerec        *fr,
67                                t_graph          **graph,
68                                gmx::MDAtoms      *mdAtoms,
69                                gmx_vsite_t       *vsite,
70                                gmx_shellfc_t     *shellfc)
71 {
72     bool  usingDomDec = DOMAINDECOMP(cr);
73
74     int   numAtomIndex, numHomeAtoms;
75     int  *atomIndex;
76
77     if (usingDomDec)
78     {
79         numAtomIndex = dd_natoms_mdatoms(cr->dd);
80         atomIndex    = cr->dd->gatindex;
81         numHomeAtoms = cr->dd->nat_home;
82     }
83     else
84     {
85         numAtomIndex = -1;
86         atomIndex    = nullptr;
87         numHomeAtoms = top_global->natoms;
88     }
89     atoms2md(top_global, ir, numAtomIndex, atomIndex, numHomeAtoms, mdAtoms);
90
91     auto mdatoms = mdAtoms->mdatoms();
92     if (usingDomDec)
93     {
94         dd_sort_local_top(cr->dd, mdatoms, top);
95     }
96     else
97     {
98         /* Currently gmx_generate_local_top allocates and returns a pointer.
99          * We should implement a more elegant solution.
100          */
101         gmx_localtop_t *tmpTop;
102
103         tmpTop = gmx_mtop_generate_local_top(top_global, ir->efep != efepNO);
104         *top   = *tmpTop;
105         sfree(tmpTop);
106     }
107
108     if (vsite)
109     {
110         if (usingDomDec)
111         {
112             /* The vsites were already assigned by the domdec topology code.
113              * We only need to do the thread division here.
114              */
115             split_vsites_over_threads(top->idef.il, top->idef.iparams,
116                                       mdatoms, vsite);
117         }
118         else
119         {
120             set_vsite_top(vsite, top, mdatoms);
121         }
122     }
123
124     if (!usingDomDec && ir->ePBC != epbcNONE && !fr->bMolPBC)
125     {
126         GMX_ASSERT(graph != NULL, "We use a graph with PBC (no periodic mols) and without DD");
127
128         *graph = mk_graph(nullptr, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
129     }
130     else if (graph != nullptr)
131     {
132         *graph = nullptr;
133     }
134
135     /* Note that with DD only flexible constraints, not shells, are supported
136      * and these don't require setup in make_local_shells().
137      */
138     if (!usingDomDec && shellfc)
139     {
140         make_local_shells(cr, mdatoms, shellfc);
141     }
142
143     setup_bonded_threading(fr->bondedThreading,
144                            fr->natoms_force,
145                            &top->idef);
146
147     gmx_pme_reinit_atoms(fr->pmedata, numHomeAtoms, mdatoms->chargeA);
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      * TODO: this only handles the GPU logic so far, should handle CPU as well.
151      * TODO: this also does not account for TPI.
152      */
153 }