2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 /*! \inpublicapi \file
37 * Implements nblib tpr reading
39 * \author Victor Holanda <victor.holanda@cscs.ch>
40 * \author Joe Jordan <ejjordan@kth.se>
41 * \author Prashanth Kanduri <kanduri@cscs.ch>
42 * \author Sebastian Keller <keller@cscs.ch>
45 #include "nblib/tpr.h"
47 #include "gromacs/fileio/tpxio.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/listed_forces/listed_forces.h"
50 #include "gromacs/mdtypes/forcerec.h"
51 #include "gromacs/mdlib/mdatoms.h"
52 #include "gromacs/mdlib/forcerec.h"
53 #include "gromacs/mdtypes/commrec.h"
54 #include "gromacs/mdtypes/iforceprovider.h"
55 #include "gromacs/mdtypes/inputrec.h"
56 #include "gromacs/mdtypes/mdatom.h"
57 #include "gromacs/mdtypes/simulation_workload.h"
58 #include "gromacs/mdtypes/state.h"
59 #include "gromacs/topology/mtop_util.h"
60 #include "gromacs/topology/topology.h"
61 #include "gromacs/utility/logger.h"
62 #include "nblib/box.h"
67 TprReader::TprReader(std::string filename)
69 t_inputrec inputRecord;
71 gmx_mtop_t molecularTopology;
72 gmx::SimulationWorkload simulationWorkload;
74 // If the file does not exist, this function will throw
75 PartialDeserializedTprFile partialDeserializedTpr =
76 read_tpx_state(filename.c_str(), &inputRecord, &globalState, &molecularTopology);
79 MPI_Comm simulationCommunicator = MPI_COMM_WORLD;
80 CommrecHandle crHandle = init_commrec(simulationCommunicator);
81 t_commrec* commrec = crHandle.get();
82 assert((commrec != nullptr) && "Must have valid commrec");
85 t_forcerec forceRecord;
86 gmx::ForceProviders forceProviders;
87 forceRecord.forceProviders = &forceProviders;
88 init_forcerec(nullptr,
98 gmx::ArrayRef<const std::string>{},
101 nonbondedParameters_ = makeNonBondedParameterLists(molecularTopology.ffparams.atnr,
102 molecularTopology.ffparams.iparams,
103 forceRecord.haveBuckingham);
105 gmx_localtop_t localtop(molecularTopology.ffparams);
106 gmx_mtop_generate_local_top(molecularTopology, &localtop, false);
107 exclusionListElements_ = std::vector<int>(localtop.excls.elementsView().begin(),
108 localtop.excls.elementsView().end());
109 exclusionListRanges_ = std::vector<int>(localtop.excls.listRangesView().begin(),
110 localtop.excls.listRangesView().end());
112 int ntopatoms = molecularTopology.natoms;
113 std::unique_ptr<gmx::MDAtoms> mdAtoms =
114 gmx::makeMDAtoms(nullptr, molecularTopology, inputRecord, false);
115 atoms2md(molecularTopology, inputRecord, -1, {}, ntopatoms, mdAtoms.get());
116 update_mdatoms(mdAtoms->mdatoms(), inputRecord.fepvals->init_lambda);
117 auto numParticles = mdAtoms->mdatoms()->nr;
118 charges_.resize(numParticles);
119 particleTypeIdOfAllParticles_.resize(numParticles);
120 inverseMasses_.resize(numParticles);
121 for (int i = 0; i < numParticles; i++)
123 charges_[i] = mdAtoms->mdatoms()->chargeA[i];
124 particleTypeIdOfAllParticles_[i] = mdAtoms->mdatoms()->typeA[i];
125 inverseMasses_[i] = mdAtoms->mdatoms()->invmass[i];
127 particleInteractionFlags_ = forceRecord.atomInfo;
129 if (TRICLINIC(globalState.box))
131 throw InputException("Only rectangular unit-cells are supported here");
133 boxX_ = globalState.box[XX][XX];
134 boxY_ = globalState.box[YY][YY];
135 boxZ_ = globalState.box[ZZ][ZZ];
136 coordinates_.assign(globalState.x.begin(), globalState.x.end());
137 velocities_.assign(globalState.v.begin(), globalState.v.end());
140 Box TprReader::getBox() const
142 return Box(boxX_, boxY_, boxZ_);