d1ec5a14f540d5fa0afcceab7797c01f7d7c3962
[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, 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 "gromacs/math/vectypes.h"
43 #include "gromacs/mdlib/tgroup.h"
44
45 #define GMX_QMMM (GMX_QMMM_MOPAC || GMX_QMMM_GAMESS || GMX_QMMM_GAUSSIAN || GMX_QMMM_ORCA)
46
47 struct gmx_localtop_t;
48 struct gmx_mtop_t;
49 struct t_commrec;
50 struct t_forcerec;
51 struct t_inputrec;
52 struct t_mdatoms;
53 struct t_QMMMrec;
54
55 typedef struct {
56     int                nrQMatoms;      /* total nr of QM atoms              */
57     rvec              *xQM;            /* shifted to center of box          */
58     int               *indexQM;        /* atom i = atom indexQM[i] in mdrun */
59     int               *atomicnumberQM; /* atomic numbers of QM atoms        */
60     real              *QMcharges;      /* atomic charges of QM atoms(ONIOM) */
61     int               *shiftQM;
62     int                QMcharge;       /* charge of the QM system           */
63     int                multiplicity;   /* multipicity (no of unpaired eln)  */
64     int                QMmethod;       /* see enums.h for all methods       */
65     int                QMbasis;        /* see enums.h for all bases         */
66     int                nelectrons;     /* total number of elecs in QM region*/
67     /* Gaussian specific stuff */
68     int                nQMcpus;        /* no. of CPUs used for the QM calc. */
69     int                QMmem;          /* memory for the gaussian calc.     */
70     int                accuracy;       /* convergence criterium (E(-x))     */
71     gmx_bool           cpmcscf;        /* using cpmcscf(l1003)*/
72     char              *gauss_dir;
73     char              *gauss_exe;
74     char              *devel_dir;
75     char              *orca_basename; /* basename for I/O with orca        */
76     char              *orca_dir;      /* directory for ORCA                */
77     /* Surface hopping stuff */
78     gmx_bool           bSH;           /* surface hopping (diabatic only)   */
79     real               SAon;          /* at which energy gap the SA starts */
80     real               SAoff;         /* at which energy gap the SA stops  */
81     int                SAsteps;       /* stepwise switchinng on the SA     */
82     int                SAstep;        /* current state of SA               */
83     int                CIdim;
84     real              *CIvec1;
85     real              *CIvec2;
86     real              *CIvec1old;
87     real              *CIvec2old;
88     ivec               SHbasis;
89     int                CASelectrons;
90     int                CASorbitals;
91 } t_QMrec;
92
93 typedef struct {
94     int            nrMMatoms;   /* nr of MM atoms, updated every step*/
95     rvec          *xMM;         /* shifted to center of box          */
96     int           *indexMM;     /* atom i = atom indexMM[I] in mdrun */
97     real          *MMcharges;   /* MM point charges in std QMMM calc.*/
98     int           *shiftMM;
99     int           *MMatomtype;  /* only important for semi-emp.      */
100     real           scalefactor;
101 } t_MMrec;
102
103
104 typedef struct t_QMMMrec {
105     int             QMMMscheme; /* ONIOM (multi-layer) or normal          */
106     int             nrQMlayers; /* number of QM layers (total layers +1 (MM)) */
107     t_QMrec       **qm;         /* atoms and run params for each QM group */
108     t_MMrec        *mm;         /* there can only be one MM subsystem !   */
109 } t_QMMMrec;
110
111 void atomic_number(int nr, char ***atomtype, int *nucnum);
112
113 t_QMMMrec *mk_QMMMrec(void);
114 /* allocates memory for QMMMrec */
115
116 #if !GMX_QMMM
117 [[noreturn]]
118 #endif
119 void init_QMMMrec(const t_commrec  *cr,
120                   gmx_mtop_t       *mtop,
121                   t_inputrec       *ir,
122                   const t_forcerec *fr);
123
124 /* init_QMMMrec initializes the QMMM record. From
125  * topology->atoms.atomname and topology->atoms.atomtype the atom
126  * names and types are read; from inputrec->QMcharge
127  * resp. inputrec->QMmult the nelecs and multiplicity are determined
128  * and md->cQMMM gives numbers of the MM and QM atoms
129  */
130 #if !GMX_QMMM
131 [[noreturn]]
132 #endif
133 void update_QMMMrec(const t_commrec  *cr,
134                     const t_forcerec *fr,
135                     const rvec       *x,
136                     const t_mdatoms  *md,
137                     const matrix      box);
138
139 /* update_QMMMrec fills the MM stuff in QMMMrec. The MM atoms are
140  * taken froom the neighbourlists of the QM atoms. In a QMMM run this
141  * routine should be called at every step, since it updates the MM
142  * elements of the t_QMMMrec struct.
143  */
144 real calculate_QMMM(const t_commrec  *cr,
145                     rvec              f[],
146                     const t_forcerec *fr);
147
148 /* QMMM computes the QM forces. This routine makes either function
149  * calls to gmx QM routines (derived from MOPAC7 (semi-emp.) and MPQC
150  * (ab initio)) or generates input files for an external QM package
151  * (listed in QMMMrec.QMpackage). The binary of the QM package is
152  * called by system().
153  */
154
155 #endif