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