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, 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.
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/legacyheaders/force.h"
54 #include "gromacs/fileio/confio.h"
55 #include "gromacs/legacyheaders/names.h"
56 #include "gromacs/legacyheaders/network.h"
57 #include "gromacs/legacyheaders/ns.h"
58 #include "gromacs/legacyheaders/nrnb.h"
59 #include "gromacs/legacyheaders/txtdump.h"
60 #include "gromacs/legacyheaders/qmmm.h"
61 #include "gromacs/utility/fatalerror.h"
64 /* mopac interface routines */
66 F77_FUNC(domldt, DOMLDT) (int *nrqmat, int labels[], char keywords[]);
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[]);
75 void init_mopac(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
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
88 if (!qm->bSH) /* if rerun then grad should not be done! */
90 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
92 eQMmethod_names[qm->QMmethod]);
96 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
98 eQMmethod_names[qm->QMmethod],
99 qm->CASorbitals, qm->CASelectrons/2);
101 F77_FUNC(domldt, DOMLDT) (&qm->nrQMatoms, qm->atomicnumberQM, keywords);
102 fprintf(stderr, "keywords are: %s\n", keywords);
107 real call_mopac(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
108 rvec f[], rvec fshift[])
110 /* do the actual QMMM calculation using directly linked mopac subroutines
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;
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
125 for (i = 0; i < qm->nrQMatoms; i++)
127 for (j = 0; j < DIM; j++)
129 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
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....
138 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
139 " with MOPAC QM subroutines\n");
143 /* now compute the energy and the gradients.
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.
152 for (i = 0; i < qm->nrQMatoms; i++)
154 for (j = 0; j < DIM; j++)
156 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
157 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
160 QMener = (real)CAL2JOULE*energy;
161 /* do we do something with the mulliken charges?? */
170 real call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
171 rvec f[], rvec fshift[])
173 /* do the actual SH QMMM calculation using directly linked mopac
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;
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
190 for (i = 0; i < qm->nrQMatoms; i++)
192 for (j = 0; j < DIM; j++)
194 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
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....
203 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
207 /* now compute the energy and the gradients.
209 snew(qmchrg, qm->nrQMatoms);
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.
216 for (i = 0; i < qm->nrQMatoms; i++)
218 for (j = 0; j < DIM; j++)
220 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
221 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
224 QMener = (real)CAL2JOULE*energy;
229 } /* call_mopac_SH */
233 gmx_qmmm_mopac_empty;