Move ns.h and nsgrid.h to mdlib/
[gromacs.git] / src / gromacs / mdlib / qm_mopac.cpp
blob178cfc7e1090a83b2ec49579968b8811eb98fefe
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "config.h"
41 #if GMX_QMMM_MOPAC
43 #include <math.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <string.h>
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/network.h"
51 #include "gromacs/legacyheaders/nrnb.h"
52 #include "gromacs/legacyheaders/txtdump.h"
53 #include "gromacs/legacyheaders/typedefs.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdlib/force.h"
57 #include "gromacs/mdlib/ns.h"
58 #include "gromacs/mdlib/qmmm.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
63 /* mopac interface routines */
64 void
65 F77_FUNC(domldt, DOMLDT) (int *nrqmat, int labels[], char keywords[]);
67 void
68 F77_FUNC(domop, DOMOP) (int *nrqmat, double qmcrd[], int *nrmmat,
69 double mmchrg[], double mmcrd[], double qmgrad[],
70 double mmgrad[], double *energy, double qmcharges[]);
74 void init_mopac(t_QMrec *qm)
76 /* initializes the mopac routines ans sets up the semiempirical
77 * computation by calling moldat(). The inline mopac routines can
78 * only perform gradient operations. If one would like to optimize a
79 * structure or find a transition state at PM3 level, gaussian is
80 * used instead.
82 char
83 *keywords;
85 snew(keywords, 240);
87 if (!qm->bSH) /* if rerun then grad should not be done! */
89 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
90 qm->QMcharge,
91 eQMmethod_names[qm->QMmethod]);
93 else
95 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
96 qm->QMcharge,
97 eQMmethod_names[qm->QMmethod],
98 qm->CASorbitals, qm->CASelectrons/2);
100 F77_FUNC(domldt, DOMLDT) (&qm->nrQMatoms, qm->atomicnumberQM, keywords);
101 fprintf(stderr, "keywords are: %s\n", keywords);
102 free(keywords);
104 } /* init_mopac */
106 real call_mopac(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
108 /* do the actual QMMM calculation using directly linked mopac subroutines
110 double /* always double as the MOPAC routines are always compiled in
111 double precission! */
112 *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
113 *qmgrad, *mmgrad = NULL, energy;
115 i, j;
116 real
117 QMener = 0.0;
118 snew(qmcrd, 3*(qm->nrQMatoms));
119 snew(qmgrad, 3*(qm->nrQMatoms));
120 /* copy the data from qr into the arrays that are going to be used
121 * in the fortran routines of MOPAC
123 for (i = 0; i < qm->nrQMatoms; i++)
125 for (j = 0; j < DIM; j++)
127 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
130 if (mm->nrMMatoms)
132 /* later we will add the point charges here. There are some
133 * conceptual problems with semi-empirical QM in combination with
134 * point charges that we need to solve first....
136 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
137 " with MOPAC QM subroutines\n");
139 else
141 /* now compute the energy and the gradients.
144 snew(qmchrg, qm->nrQMatoms);
145 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
146 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
147 /* add the gradients to the f[] array, and also to the fshift[].
148 * the mopac gradients are in kCal/angstrom.
150 for (i = 0; i < qm->nrQMatoms; i++)
152 for (j = 0; j < DIM; j++)
154 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
155 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
158 QMener = (real)CAL2JOULE*energy;
159 /* do we do something with the mulliken charges?? */
161 free(qmchrg);
163 free(qmgrad);
164 free(qmcrd);
165 return (QMener);
168 real call_mopac_SH(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
170 /* do the actual SH QMMM calculation using directly linked mopac
171 subroutines */
173 double /* always double as the MOPAC routines are always compiled in
174 double precission! */
175 *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
176 *qmgrad, *mmgrad = NULL, energy;
178 i, j;
179 real
180 QMener = 0.0;
182 snew(qmcrd, 3*(qm->nrQMatoms));
183 snew(qmgrad, 3*(qm->nrQMatoms));
184 /* copy the data from qr into the arrays that are going to be used
185 * in the fortran routines of MOPAC
187 for (i = 0; i < qm->nrQMatoms; i++)
189 for (j = 0; j < DIM; j++)
191 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
194 if (mm->nrMMatoms)
196 /* later we will add the point charges here. There are some
197 * conceptual problems with semi-empirical QM in combination with
198 * point charges that we need to solve first....
200 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
202 else
204 /* now compute the energy and the gradients.
206 snew(qmchrg, qm->nrQMatoms);
208 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
209 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
210 /* add the gradients to the f[] array, and also to the fshift[].
211 * the mopac gradients are in kCal/angstrom.
213 for (i = 0; i < qm->nrQMatoms; i++)
215 for (j = 0; j < DIM; j++)
217 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
218 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
221 QMener = (real)CAL2JOULE*energy;
223 free(qmgrad);
224 free(qmcrd);
225 return (QMener);
226 } /* call_mopac_SH */
228 #else
230 gmx_qmmm_mopac_empty;
231 #endif