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.
37 * QMMMTopologyPrepocessor class responsible for
38 * all modificatios of the topology during input pre-processing
40 * \author Dmitry Morozov <dmitry.morozov@jyu.fi>
41 * \ingroup module_applied_forces
43 #ifndef GMX_APPLIED_FORCES_QMMMTOPOLOGYPREPROCESSOR_H
44 #define GMX_APPLIED_FORCES_QMMMTOPOLOGYPREPROCESSOR_H
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/arrayref.h"
52 #include "qmmmtypes.h"
60 * \brief Contains various information about topology modifications
61 * Used for statistics during topology pre-processing within QMMMTopologyPreprocessor class
63 struct QMMMTopologyInfo
65 //! Total number of MM atoms
67 //! Total number of QM atoms
69 //! Total remaining charge of MM part
70 real remainingMMCharge = 0.0;
71 //! Total classical charge removed from QM atoms
72 real totalClassicalChargeOfQMAtoms = 0.0;
73 //! Total number of Non-bonded (LJ) exclusions made for QM-QM interactions
74 int numExclusionsMade = 0;
75 //! Total number of removed classical Bonds between QM-QM atoms
76 int numBondsRemoved = 0;
77 //! Total number of removed classical Angles between QM-QM atoms
78 int numAnglesRemoved = 0;
79 //! Total number of removed classical Dihedrals between QM-QM atoms
80 int numDihedralsRemoved = 0;
81 //! Total number of removed F_SETTLE between QM-QM atoms
82 int numSettleRemoved = 0;
83 //! Total number of empty chemical bonds (F_CONNBONDS) added between QM-QM atoms
84 int numConnBondsAdded = 0;
85 //! Total number of virtual sites, that consisting of QM atoms only, which charge has been removed
86 int numVirtualSitesModified = 0;
87 //! Total number of constrained bonds within QM subsystem
88 int numConstrainedBondsInQMSubsystem = 0;
89 //! Total number of broken bonds between QM and MM atoms (Link Frontier)
94 * \brief Class implementing gmx_mtop_t QMMM modifications during preprocessing
95 * 1) Split QM-containing molecules from other molecules in blocks
96 * 2) Nullify charges on all virtual sites consisting of QM only atoms
97 * 3) Nullifies charges on all QM atoms
98 * 4) Excludes LJ interactions between QM atoms
99 * 5) Builds vector with atomic numbers of all atoms
100 * 6) Makes F_CONNBOND between atoms within QM region
101 * 7) Removes angles and settles containing 2 or more QM atoms
102 * 8) Removes dihedrals containing 3 or more QM atoms
103 * 9) Builds vector containing pairs of bonded QM - MM atoms (Link frontier)
105 class QMMMTopologyPreprocessor
108 /*! \brief Constructor for QMMMTopologyPreprocessor from its parameters
110 * \param[in] qmIndices Array with global indicies of QM atoms
112 QMMMTopologyPreprocessor(ArrayRef<const index> qmIndices);
114 /*! \brief Pocesses mtop topology and prepares atomNumbers_ and linkFrontier_ vectors
115 * Builds topInfo_ containing information about topology modifications
117 * \param[in,out] mtop Topology that needs to be modified
119 void preprocess(gmx_mtop_t* mtop);
121 //! \brief Returns data about modifications made via QMMMTopologyInfo
122 const QMMMTopologyInfo& topInfo() const;
124 //! \brief Returns view of atomic numbers for all atoms in the processed topology
125 ArrayRef<const int> atomNumbers() const;
127 //! \brief Returns view of point charges for all atoms in the processed topology
128 ArrayRef<const real> atomCharges() const;
130 //! \brief Returns view of the whole Link Frontier for the processed topology
131 ArrayRef<const LinkFrontier> linkFrontier() const;
134 //! Retruns true if globalAtomIndex belongs to QM region
135 bool isQMAtom(index globalAtomIndex);
137 /*! \brief Splits QM containing molecules out of MM blocks in topology
138 * Modifies blocks in topology
139 * Updates bQMBlock vector containing QM flags of all blocks in modified mtop
141 void splitQMblocks(gmx_mtop_t* mtop);
143 /*! \brief Removes classical charges from QM atoms
144 * Provides data about removed charge via topInfo_
146 void removeQMClassicalCharges(gmx_mtop_t* mtop);
148 //! \brief Build exlusion list for LJ interactions between QM atoms
149 void addQMLJExclusions(gmx_mtop_t* mtop);
151 /*! \brief Builds atomNumbers_ vector
152 * Provides data about total number of QM and MM atoms via topInfo_
154 void buildQMMMAtomNumbers(gmx_mtop_t* mtop);
156 /*! \brief Modifies pairwise bonded interactions
157 * Removes any other pairwise bonded interactions between QM-QM atoms
158 * Creates F_CONNBOND between QM atoms
159 * Any restraints and constraints will be kept
160 * Provides data about modifications via topInfo_
162 void modifyQMMMTwoCenterInteractions(gmx_mtop_t* mtop);
164 /*! \brief Builds link_ vector with pairs of atoms indicting broken QM - MM chemical bonds.
165 * Also performs search of constrained bonds within QM subsystem.
167 void buildQMMMLink(gmx_mtop_t* mtop);
169 /*! \brief Modifies three-centers interactions (i.e. Angles, Settles)
170 * Removes any other three-centers bonded interactions including 2 or more QM atoms
171 * Any restraints and constraints will be kept
172 * Any F_SETTLE containing QM atoms will be converted to the pair of F_CONNBONDS
173 * Provides data about modifications via topInfo_
175 void modifyQMMMThreeCenterInteractions(gmx_mtop_t* mtop);
177 /*! \brief Modifies four-centers interactions
178 * Removes any other four-centers bonded interactions including 3 or more QM atoms
179 * Any restraints and constraints will be kept
180 * Provides data about modifications via topInfo_
182 void modifyQMMMFourCenterInteractions(gmx_mtop_t* mtop);
184 //! \brief Removes charge from all virtual sites which are consists of only QM atoms
185 void modifyQMMMVirtualSites(gmx_mtop_t* mtop);
187 //! Vector indicating which molblocks have QM atoms
188 std::vector<bool> bQMBlock_;
189 /*! \brief Global indices of QM atoms;
190 * The dominant operation is search and we also expect the set of qm atoms to be very small
191 * relative to the rest, so set should outperform unordered set, i.e. unsorted std::vector.
193 std::set<int> qmIndices_;
194 //! Vector with atom numbers for the whole system
195 std::vector<int> atomNumbers_;
196 //! Vector with atom point charges for the whole system
197 std::vector<real> atomCharges_;
198 //! Vector with pairs of indices defining broken bonds in QMMM
199 std::vector<LinkFrontier> linkFrontier_;
200 //! Structure with information about modifications made
201 QMMMTopologyInfo topInfo_;