Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / mimic / communicator.cpp
blobd83ccc3cee5c9b47c0a79fd6346bdf79aad16a5d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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 "communicator.h"
39 #include "config.h"
41 #include <unordered_set>
43 #include "gromacs/math/units.h"
44 #include "gromacs/utility/fatalerror.h"
46 #if GMX_MIMIC
47 #include <DataTypes.h>
48 #include <MessageApi.h>
49 #endif
51 // When not built in a configuration with QMMM support, much of this
52 // code is unreachable by design. Tell clang not to warn about it.
53 #pragma GCC diagnostic push
54 #pragma GCC diagnostic ignored "-Wmissing-noreturn"
56 #if !GMX_MIMIC
57 //! \brief Definitions to stub the ones defined in DataTypes.h
58 constexpr int TYPE_INT = 0, TYPE_DOUBLE = 0;
60 /*! \brief Stub communication library function to call in case if
61 * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
63 static void MCL_init_client(const char *) // NOLINT(readability-named-parameter)
65 GMX_RELEASE_ASSERT(GMX_MIMIC, "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
68 /*! \brief Stub communication library function to call in case if
69 * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
71 static void MCL_send(void *, int, int, int) // NOLINT(readability-named-parameter)
73 GMX_RELEASE_ASSERT(GMX_MIMIC, "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
76 /*! \brief Stub communication library function to call in case if
77 * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
79 static void MCL_receive(void *, int, int, int) // NOLINT(readability-named-parameter)
81 GMX_RELEASE_ASSERT(GMX_MIMIC, "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
84 /*! \brief Stub communication library function to call in case if
85 * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
87 static void MCL_destroy()
89 GMX_RELEASE_ASSERT(GMX_MIMIC, "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
91 #endif
93 void gmx::MimicCommunicator::init()
95 char path[GMX_PATH_MAX];
96 gmx_getcwd(path, GMX_PATH_MAX);
97 return MCL_init_client(path);
100 void gmx::MimicCommunicator::sendInitData(gmx_mtop_t *mtop,
101 PaddedHostVector<gmx::RVec> coords)
103 MCL_send(&mtop->natoms, 1, TYPE_INT, 0);
104 MCL_send(&mtop->atomtypes.nr, 1, TYPE_INT, 0);
106 std::vector<int> atomTypes;
107 std::vector<int> nAtomsMol;
108 std::vector<int> idOrder;
109 std::vector<double> charges;
110 std::vector<double> masses(mtop->atomtypes.nr, -1);
111 std::vector<int> elements(mtop->atomtypes.nr, -1);
112 std::vector<int> bonds;
113 std::vector<double> bondLengths;
114 std::unordered_set<int> existingTypes;
116 atomTypes.reserve(static_cast<size_t>(mtop->natoms));
117 charges.reserve(static_cast<size_t>(mtop->natoms));
119 int offset = 0;
120 for (const gmx_molblock_t &molblock : mtop->molblock)
122 gmx_moltype_t *type = &mtop->moltype[molblock.type];
123 for (int mol = 0; mol < molblock.nmol; ++mol)
125 int nconstr = type->ilist[F_CONSTR].size() / 3;
126 int nconstrc = type->ilist[F_CONSTRNC].size() / 3;
127 int nsettle = type->ilist[F_SETTLE].size() / 4;
129 for (int ncon = 0; ncon < nconstr + nconstrc; ++ncon)
131 int contype = type->ilist[F_CONSTR].iatoms[0];
132 int at1 = type->ilist[F_CONSTR].iatoms[1];
133 int at2 = type->ilist[F_CONSTR].iatoms[2];
134 bonds.push_back(offset + at1 + 1);
135 bonds.push_back(offset + at2 + 1);
136 bondLengths.push_back(static_cast<double>(mtop->ffparams.iparams[contype].constr.dA) / BOHR2NM);
139 for (int ncon = 0; ncon < nsettle; ++ncon)
141 t_iatom ox;
142 t_iatom h1;
143 t_iatom h2;
145 int contype = type->ilist[F_SETTLE].iatoms[0];
147 ox = type->ilist[F_SETTLE].iatoms[1];
148 h1 = type->ilist[F_SETTLE].iatoms[2];
149 h2 = type->ilist[F_SETTLE].iatoms[3];
151 bonds.push_back(offset + ox + 1);
152 bonds.push_back(offset + h1 + 1);
154 bonds.push_back(offset + ox + 1);
155 bonds.push_back(offset + h2 + 1);
157 bonds.push_back(offset + h1 + 1);
158 bonds.push_back(offset + h2 + 1);
159 bondLengths.push_back(static_cast<double>(mtop->ffparams.iparams[contype].constr.dA) / BOHR2NM);
160 bondLengths.push_back(static_cast<double>(mtop->ffparams.iparams[contype].constr.dA) / BOHR2NM);
161 bondLengths.push_back(static_cast<double>(mtop->ffparams.iparams[contype].constr.dB) / BOHR2NM);
164 nAtomsMol.push_back(type->atoms.nr);
165 for (int at = 0; at < type->atoms.nr; ++at)
167 int atomtype = type->atoms.atom[at].type;
168 auto charge = static_cast<double>(type->atoms.atom[at].q);
169 idOrder.push_back(offset + 1);
170 offset++;
171 atomTypes.push_back(atomtype + 1);
172 charges.push_back(charge);
173 if (existingTypes.insert(atomtype).second)
175 masses[atomtype] = type->atoms.atom[at].m;
176 elements[atomtype] = type->atoms.atom[at].atomnumber;
181 // sending atom types
182 MCL_send(&*atomTypes.begin(), mtop->natoms, TYPE_INT, 0);
184 int max_multipole_order = 0;
185 // sending multipole orders
186 MCL_send(&max_multipole_order, 1, TYPE_INT, 0);
188 int nMolecules = nAtomsMol.size();
189 // sending molecule number
190 MCL_send(&nMolecules, 1, TYPE_INT, 0);
192 // sending number of atoms per molecules
193 MCL_send(&*nAtomsMol.begin(), nAtomsMol.size(), TYPE_INT, 0);
195 int nBonds = bonds.size() / 2;
196 // sending number of bond constraints
197 MCL_send(&nBonds, 1, TYPE_INT, 0);
199 // sending number of angle constraints
200 MCL_send(&max_multipole_order, 1, TYPE_INT, 0);
202 if (nBonds > 0)
204 // sending bonded atoms indices
205 MCL_send(&*bonds.begin(), bonds.size(), TYPE_INT, 0);
207 // sending bond lengths
208 MCL_send(&*bondLengths.begin(), bondLengths.size(), TYPE_DOUBLE, 0);
211 // sending array of atomic charges
212 MCL_send(&*charges.begin(), mtop->natoms, TYPE_DOUBLE, 0);
214 // sending array of atomic masses
215 MCL_send(&*masses.begin(), mtop->atomtypes.nr, TYPE_DOUBLE, 0);
217 // sending ids of atoms per molecule
218 MCL_send(&*idOrder.begin(), idOrder.size(), TYPE_INT, 0);
220 // sending list of elements
221 MCL_send(&*elements.begin(), mtop->atomtypes.nr, TYPE_INT, 0);
223 std::vector<double> convertedCoords;
224 for (auto &coord : coords)
226 convertedCoords.push_back(static_cast<double>(coord[0]) / BOHR2NM);
227 convertedCoords.push_back(static_cast<double>(coord[1]) / BOHR2NM);
228 convertedCoords.push_back(static_cast<double>(coord[2]) / BOHR2NM);
231 // sending array of coordinates
232 MCL_send(&*convertedCoords.begin(), 3 * mtop->natoms, TYPE_DOUBLE, 0);
235 int64_t gmx::MimicCommunicator::getStepNumber()
237 int steps;
238 MCL_receive(&steps, 1, TYPE_INT, 0);
239 return steps;
242 void gmx::MimicCommunicator::getCoords(PaddedHostVector<RVec> *x, const int natoms)
244 std::vector<double> coords(natoms * 3);
245 MCL_receive(&*coords.begin(), 3 * natoms, TYPE_DOUBLE, 0);
246 for (int j = 0; j < natoms; ++j)
248 (*x)[j][0] = static_cast<real>(coords[j * 3] * BOHR2NM);
249 (*x)[j][1] = static_cast<real>(coords[j * 3 + 1] * BOHR2NM);
250 (*x)[j][2] = static_cast<real>(coords[j * 3 + 2] * BOHR2NM);
254 void gmx::MimicCommunicator::sendEnergies(real energy)
256 double convertedEnergy = energy / (HARTREE2KJ * AVOGADRO);
257 MCL_send(&convertedEnergy, 1, TYPE_DOUBLE, 0);
260 void gmx::MimicCommunicator::sendForces(gmx::ArrayRef<gmx::RVec> forces, int natoms)
262 std::vector<double> convertedForce;
263 for (int j = 0; j < natoms; ++j)
265 convertedForce.push_back(static_cast<real>(forces[j][0]) / HARTREE_BOHR2MD);
266 convertedForce.push_back(static_cast<real>(forces[j][1]) / HARTREE_BOHR2MD);
267 convertedForce.push_back(static_cast<real>(forces[j][2]) / HARTREE_BOHR2MD);
269 MCL_send(&*convertedForce.begin(), convertedForce.size(), TYPE_DOUBLE, 0);
272 void gmx::MimicCommunicator::finalize()
274 MCL_destroy();
277 #pragma GCC diagnostic pop