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