Adjust test range and fix singleaccuracy tests
[gromacs.git] / src / gromacs / domdec / mdsetup.cpp
blob22a084ed46ec6878eec4eb0ca27248691eb08bbe
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019, 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.
35 #include "gmxpre.h"
37 #include "mdsetup.h"
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/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"
56 /* TODO: Add a routine that collects the initial setup of the algorithms.
58 * The final solution should be an MD algorithm base class with methods
59 * for initialization and atom-data setup.
61 void mdAlgorithmsSetupAtomData(const t_commrec *cr,
62 const t_inputrec *ir,
63 const gmx_mtop_t &top_global,
64 gmx_localtop_t *top,
65 t_forcerec *fr,
66 t_graph **graph,
67 gmx::MDAtoms *mdAtoms,
68 gmx::Constraints *constr,
69 gmx_vsite_t *vsite,
70 gmx_shellfc_t *shellfc)
72 bool usingDomDec = DOMAINDECOMP(cr);
74 int numAtomIndex, numHomeAtoms;
75 int *atomIndex;
77 if (usingDomDec)
79 numAtomIndex = dd_natoms_mdatoms(cr->dd);
80 atomIndex = cr->dd->globalAtomIndices.data();
81 numHomeAtoms = dd_numHomeAtoms(*cr->dd);
83 else
85 numAtomIndex = -1;
86 atomIndex = nullptr;
87 numHomeAtoms = top_global.natoms;
89 atoms2md(&top_global, ir, numAtomIndex, atomIndex, numHomeAtoms, mdAtoms);
91 auto mdatoms = mdAtoms->mdatoms();
92 if (usingDomDec)
94 dd_sort_local_top(cr->dd, mdatoms, top);
96 else
98 gmx_mtop_generate_local_top(top_global, top, ir->efep != efepNO);
101 if (vsite)
103 if (usingDomDec)
105 /* The vsites were already assigned by the domdec topology code.
106 * We only need to do the thread division here.
108 split_vsites_over_threads(top->idef.il, top->idef.iparams,
109 mdatoms, vsite);
111 else
113 set_vsite_top(vsite, top, mdatoms);
117 if (!usingDomDec && ir->ePBC != epbcNONE && !fr->bMolPBC)
119 GMX_ASSERT(graph != nullptr, "We use a graph with PBC (no periodic mols) and without DD");
121 *graph = mk_graph(nullptr, &(top->idef), 0, top_global.natoms, FALSE, FALSE);
123 else if (graph != nullptr)
125 *graph = nullptr;
128 /* Note that with DD only flexible constraints, not shells, are supported
129 * and these don't require setup in make_local_shells().
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)
134 if (!usingDomDec && shellfc)
136 make_local_shells(cr, mdatoms, shellfc);
139 setup_bonded_threading(fr->bondedThreading,
140 fr->natoms_force,
141 fr->gpuBonded != nullptr,
142 top->idef);
144 if (EEL_PME(fr->ic->eeltype) && (cr->duty & DUTY_PME))
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().
149 const int numPmeAtoms = numHomeAtoms - fr->n_tpi;
150 gmx_pme_reinit_atoms(fr->pmedata, numPmeAtoms, mdatoms->chargeA);
153 if (constr)
155 constr->setConstraints(*top, *mdatoms);