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