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