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