Move QMMM source to C++
[alexxy/gromacs.git] / src / gromacs / mdlib / qm_mopac.cpp
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-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015, 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 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #ifdef GMX_QMMM_MOPAC
42
43 #include <math.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <string.h>
47
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/legacyheaders/force.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/network.h"
53 #include "gromacs/legacyheaders/nrnb.h"
54 #include "gromacs/legacyheaders/ns.h"
55 #include "gromacs/legacyheaders/qmmm.h"
56 #include "gromacs/legacyheaders/txtdump.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/smalloc.h"
62
63
64 /* mopac interface routines */
65 void
66     F77_FUNC(domldt, DOMLDT) (int *nrqmat, int labels[], char keywords[]);
67
68 void
69     F77_FUNC(domop, DOMOP) (int *nrqmat, double qmcrd[], int *nrmmat,
70                             double mmchrg[], double mmcrd[], double qmgrad[],
71                             double mmgrad[], double *energy, double qmcharges[]);
72
73
74
75 void init_mopac(t_QMrec *qm)
76 {
77     /* initializes the mopac routines ans sets up the semiempirical
78      * computation by calling moldat(). The inline mopac routines can
79      * only perform gradient operations. If one would like to optimize a
80      * structure or find a transition state at PM3 level, gaussian is
81      * used instead.
82      */
83     char
84     *keywords;
85
86     snew(keywords, 240);
87
88     if (!qm->bSH)  /* if rerun then grad should not be done! */
89     {
90         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
91                 qm->QMcharge,
92                 eQMmethod_names[qm->QMmethod]);
93     }
94     else
95     {
96         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
97                 qm->QMcharge,
98                 eQMmethod_names[qm->QMmethod],
99                 qm->CASorbitals, qm->CASelectrons/2);
100     }
101     F77_FUNC(domldt, DOMLDT) (&qm->nrQMatoms, qm->atomicnumberQM, keywords);
102     fprintf(stderr, "keywords are: %s\n", keywords);
103     free(keywords);
104
105 } /* init_mopac */
106
107 real call_mopac(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
108 {
109     /* do the actual QMMM calculation using directly linked mopac subroutines
110      */
111     double /* always double as the MOPAC routines are always compiled in
112               double precission! */
113     *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
114     *qmgrad, *mmgrad = NULL, energy;
115     int
116         i, j;
117     real
118         QMener = 0.0;
119     snew(qmcrd, 3*(qm->nrQMatoms));
120     snew(qmgrad, 3*(qm->nrQMatoms));
121     /* copy the data from qr into the arrays that are going to be used
122      * in the fortran routines of MOPAC
123      */
124     for (i = 0; i < qm->nrQMatoms; i++)
125     {
126         for (j = 0; j < DIM; j++)
127         {
128             qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
129         }
130     }
131     if (mm->nrMMatoms)
132     {
133         /* later we will add the point charges here. There are some
134          * conceptual problems with semi-empirical QM in combination with
135          * point charges that we need to solve first....
136          */
137         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
138                   " with MOPAC QM subroutines\n");
139     }
140     else
141     {
142         /* now compute the energy and the gradients.
143          */
144
145         snew(qmchrg, qm->nrQMatoms);
146         F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
147                                 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
148         /* add the gradients to the f[] array, and also to the fshift[].
149          * the mopac gradients are in kCal/angstrom.
150          */
151         for (i = 0; i < qm->nrQMatoms; i++)
152         {
153             for (j = 0; j < DIM; j++)
154             {
155                 f[i][j]       = (real)10*CAL2JOULE*qmgrad[3*i+j];
156                 fshift[i][j]  = (real)10*CAL2JOULE*qmgrad[3*i+j];
157             }
158         }
159         QMener = (real)CAL2JOULE*energy;
160         /* do we do something with the mulliken charges?? */
161
162         free(qmchrg);
163     }
164     free(qmgrad);
165     free(qmcrd);
166     return (QMener);
167 }
168
169 real call_mopac_SH(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
170 {
171     /* do the actual SH QMMM calculation using directly linked mopac
172        subroutines */
173
174     double /* always double as the MOPAC routines are always compiled in
175               double precission! */
176     *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
177     *qmgrad, *mmgrad = NULL, energy;
178     int
179         i, j;
180     real
181         QMener = 0.0;
182
183     snew(qmcrd, 3*(qm->nrQMatoms));
184     snew(qmgrad, 3*(qm->nrQMatoms));
185     /* copy the data from qr into the arrays that are going to be used
186      * in the fortran routines of MOPAC
187      */
188     for (i = 0; i < qm->nrQMatoms; i++)
189     {
190         for (j = 0; j < DIM; j++)
191         {
192             qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
193         }
194     }
195     if (mm->nrMMatoms)
196     {
197         /* later we will add the point charges here. There are some
198          * conceptual problems with semi-empirical QM in combination with
199          * point charges that we need to solve first....
200          */
201         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
202     }
203     else
204     {
205         /* now compute the energy and the gradients.
206          */
207         snew(qmchrg, qm->nrQMatoms);
208
209         F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
210                                 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
211         /* add the gradients to the f[] array, and also to the fshift[].
212          * the mopac gradients are in kCal/angstrom.
213          */
214         for (i = 0; i < qm->nrQMatoms; i++)
215         {
216             for (j = 0; j < DIM; j++)
217             {
218                 f[i][j]      = (real)10*CAL2JOULE*qmgrad[3*i+j];
219                 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
220             }
221         }
222         QMener = (real)CAL2JOULE*energy;
223     }
224     free(qmgrad);
225     free(qmcrd);
226     return (QMener);
227 } /* call_mopac_SH */
228
229 #else
230 int
231     gmx_qmmm_mopac_empty;
232 #endif