Sort all includes in src/gromacs
[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/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_commrec *cr, t_QMrec *qm, t_MMrec *mm)
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_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
108                 rvec f[], rvec fshift[])
109 {
110     /* do the actual QMMM calculation using directly linked mopac subroutines
111      */
112     double /* always double as the MOPAC routines are always compiled in
113               double precission! */
114     *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
115     *qmgrad, *mmgrad = NULL, energy;
116     int
117         i, j;
118     real
119         QMener = 0.0;
120     snew(qmcrd, 3*(qm->nrQMatoms));
121     snew(qmgrad, 3*(qm->nrQMatoms));
122     /* copy the data from qr into the arrays that are going to be used
123      * in the fortran routines of MOPAC
124      */
125     for (i = 0; i < qm->nrQMatoms; i++)
126     {
127         for (j = 0; j < DIM; j++)
128         {
129             qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
130         }
131     }
132     if (mm->nrMMatoms)
133     {
134         /* later we will add the point charges here. There are some
135          * conceptual problems with semi-empirical QM in combination with
136          * point charges that we need to solve first....
137          */
138         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
139                   " with MOPAC QM subroutines\n");
140     }
141     else
142     {
143         /* now compute the energy and the gradients.
144          */
145
146         snew(qmchrg, qm->nrQMatoms);
147         F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
148                                 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
149         /* add the gradients to the f[] array, and also to the fshift[].
150          * the mopac gradients are in kCal/angstrom.
151          */
152         for (i = 0; i < qm->nrQMatoms; i++)
153         {
154             for (j = 0; j < DIM; j++)
155             {
156                 f[i][j]       = (real)10*CAL2JOULE*qmgrad[3*i+j];
157                 fshift[i][j]  = (real)10*CAL2JOULE*qmgrad[3*i+j];
158             }
159         }
160         QMener = (real)CAL2JOULE*energy;
161         /* do we do something with the mulliken charges?? */
162
163         free(qmchrg);
164     }
165     free(qmgrad);
166     free(qmcrd);
167     return (QMener);
168 }
169
170 real call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
171                    rvec f[], rvec fshift[])
172 {
173     /* do the actual SH QMMM calculation using directly linked mopac
174        subroutines */
175
176     double /* always double as the MOPAC routines are always compiled in
177               double precission! */
178     *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
179     *qmgrad, *mmgrad = NULL, energy;
180     int
181         i, j;
182     real
183         QMener = 0.0;
184
185     snew(qmcrd, 3*(qm->nrQMatoms));
186     snew(qmgrad, 3*(qm->nrQMatoms));
187     /* copy the data from qr into the arrays that are going to be used
188      * in the fortran routines of MOPAC
189      */
190     for (i = 0; i < qm->nrQMatoms; i++)
191     {
192         for (j = 0; j < DIM; j++)
193         {
194             qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
195         }
196     }
197     if (mm->nrMMatoms)
198     {
199         /* later we will add the point charges here. There are some
200          * conceptual problems with semi-empirical QM in combination with
201          * point charges that we need to solve first....
202          */
203         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
204     }
205     else
206     {
207         /* now compute the energy and the gradients.
208          */
209         snew(qmchrg, qm->nrQMatoms);
210
211         F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
212                                 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
213         /* add the gradients to the f[] array, and also to the fshift[].
214          * the mopac gradients are in kCal/angstrom.
215          */
216         for (i = 0; i < qm->nrQMatoms; i++)
217         {
218             for (j = 0; j < DIM; j++)
219             {
220                 f[i][j]      = (real)10*CAL2JOULE*qmgrad[3*i+j];
221                 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
222             }
223         }
224         QMener = (real)CAL2JOULE*energy;
225     }
226     free(qmgrad);
227     free(qmcrd);
228     return (QMener);
229 } /* call_mopac_SH */
230
231 #else
232 int
233     gmx_qmmm_mopac_empty;
234 #endif