2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2017,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_MDLIB_QMMM_H
39 #define GMX_MDLIB_QMMM_H
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdlib/tgroup.h"
47 #include "gromacs/utility/arrayref.h"
49 #define GMX_QMMM (GMX_QMMM_MOPAC || GMX_QMMM_GAMESS || GMX_QMMM_GAUSSIAN || GMX_QMMM_ORCA)
51 struct gmx_localtop_t;
61 class ForceWithShiftForces;
66 int nrQMatoms; /* total nr of QM atoms */
67 rvec* xQM; /* shifted to center of box */
68 int* indexQM; /* atom i = atom indexQM[i] in mdrun */
69 int* atomicnumberQM; /* atomic numbers of QM atoms */
70 real* QMcharges; /* atomic charges of QM atoms(ONIOM) */
72 int QMcharge; /* charge of the QM system */
73 int multiplicity; /* multipicity (no of unpaired eln) */
74 int QMmethod; /* see enums.h for all methods */
75 int QMbasis; /* see enums.h for all bases */
76 int nelectrons; /* total number of elecs in QM region*/
77 /* Gaussian specific stuff */
78 int nQMcpus; /* no. of CPUs used for the QM calc. */
79 int QMmem; /* memory for the gaussian calc. */
80 int accuracy; /* convergence criterium (E(-x)) */
81 gmx_bool cpmcscf; /* using cpmcscf(l1003)*/
85 char* orca_basename; /* basename for I/O with orca */
86 char* orca_dir; /* directory for ORCA */
87 /* Surface hopping stuff */
88 gmx_bool bSH; /* surface hopping (diabatic only) */
89 real SAon; /* at which energy gap the SA starts */
90 real SAoff; /* at which energy gap the SA stops */
91 int SAsteps; /* stepwise switchinng on the SA */
92 int SAstep; /* current state of SA */
105 int nrMMatoms; /* nr of MM atoms, updated every step*/
106 rvec* xMM; /* shifted to center of box */
107 int* indexMM; /* atom i = atom indexMM[I] in mdrun */
108 real* MMcharges; /* MM point charges in std QMMM calc.*/
110 int* MMatomtype; /* only important for semi-emp. */
115 typedef struct t_QMMMrec
117 int QMMMscheme; /* ONIOM (multi-layer) or normal */
118 int nrQMlayers; /* number of QM layers (total layers +1 (MM)) */
119 t_QMrec** qm; /* atoms and run params for each QM group */
120 t_MMrec* mm; /* there can only be one MM subsystem ! */
123 void atomic_number(int nr, char*** atomtype, int* nucnum);
125 t_QMMMrec* mk_QMMMrec();
126 /* allocates memory for QMMMrec */
128 void init_QMMMrec(const t_commrec* cr, const gmx_mtop_t* mtop, const t_inputrec* ir, const t_forcerec* fr);
130 /* init_QMMMrec initializes the QMMM record. From
131 * topology->atoms.atomname and topology->atoms.atomtype the atom
132 * names and types are read; from inputrec->QMcharge
133 * resp. inputrec->QMmult the nelecs and multiplicity are determined
134 * and md->cQMMM gives numbers of the MM and QM atoms
136 void update_QMMMrec(const t_commrec* cr, const t_forcerec* fr, const rvec* x, const t_mdatoms* md, const matrix box);
138 /* update_QMMMrec fills the MM stuff in QMMMrec. The MM atoms are
139 * taken froom the neighbourlists of the QM atoms. In a QMMM run this
140 * routine should be called at every step, since it updates the MM
141 * elements of the t_QMMMrec struct.
143 real calculate_QMMM(const t_commrec* cr, gmx::ForceWithShiftForces* forceWithShiftForces, const t_QMMMrec* qmmm);
145 /* QMMM computes the QM forces. This routine makes either function
146 * calls to gmx QM routines (derived from MOPAC7 (semi-emp.) and MPQC
147 * (ab initio)) or generates input files for an external QM package
148 * (listed in QMMMrec.QMpackage). The binary of the QM package is
149 * called by system().
153 * Return vector of atom indices for atoms in the QMMM region.
155 * \param[in] mtop Topology to use for populating array.
156 * \param[in] ir Inputrec used in simulation.
157 * \returns Vector of atoms.
159 std::vector<int> qmmmAtomIndices(const t_inputrec& ir, const gmx_mtop_t& mtop);
162 * Remove charges from QMMM atoms.
164 * \param[in] mtop Topology used for removing atoms.
165 * \param[in] qmmmAtoms ArrayRef to vector conatining qmmm atom indices.
167 void removeQmmmAtomCharges(gmx_mtop_t* mtop, gmx::ArrayRef<const int> qmmmAtoms);