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