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.
47 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/fileio/confio.h"
65 #include "gmx_fatal.h"
70 /* mopac interface routines */
72 F77_FUNC(domldt, DOMLDT) (int *nrqmat, int labels[], char keywords[]);
75 F77_FUNC(domop, DOMOP) (int *nrqmat, double qmcrd[], int *nrmmat,
76 double mmchrg[], double mmcrd[], double qmgrad[],
77 double mmgrad[], double *energy, double qmcharges[]);
81 void init_mopac(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
83 /* initializes the mopac routines ans sets up the semiempirical
84 * computation by calling moldat(). The inline mopac routines can
85 * only perform gradient operations. If one would like to optimize a
86 * structure or find a transition state at PM3 level, gaussian is
94 if (!qm->bSH) /* if rerun then grad should not be done! */
96 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
98 eQMmethod_names[qm->QMmethod]);
102 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
104 eQMmethod_names[qm->QMmethod],
105 qm->CASorbitals, qm->CASelectrons/2);
107 F77_FUNC(domldt, DOMLDT) (&qm->nrQMatoms, qm->atomicnumberQM, keywords);
108 fprintf(stderr, "keywords are: %s\n", keywords);
113 real call_mopac(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
114 rvec f[], rvec fshift[])
116 /* do the actual QMMM calculation using directly linked mopac subroutines
118 double /* always double as the MOPAC routines are always compiled in
119 double precission! */
120 *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
121 *qmgrad, *mmgrad = NULL, energy;
126 snew(qmcrd, 3*(qm->nrQMatoms));
127 snew(qmgrad, 3*(qm->nrQMatoms));
128 /* copy the data from qr into the arrays that are going to be used
129 * in the fortran routines of MOPAC
131 for (i = 0; i < qm->nrQMatoms; i++)
133 for (j = 0; j < DIM; j++)
135 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
140 /* later we will add the point charges here. There are some
141 * conceptual problems with semi-empirical QM in combination with
142 * point charges that we need to solve first....
144 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
145 " with MOPAC QM subroutines\n");
149 /* now compute the energy and the gradients.
152 snew(qmchrg, qm->nrQMatoms);
153 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
154 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
155 /* add the gradients to the f[] array, and also to the fshift[].
156 * the mopac gradients are in kCal/angstrom.
158 for (i = 0; i < qm->nrQMatoms; i++)
160 for (j = 0; j < DIM; j++)
162 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
163 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
166 QMener = (real)CAL2JOULE*energy;
167 /* do we do something with the mulliken charges?? */
176 real call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
177 rvec f[], rvec fshift[])
179 /* do the actual SH QMMM calculation using directly linked mopac
182 double /* always double as the MOPAC routines are always compiled in
183 double precission! */
184 *qmcrd = NULL, *qmchrg = NULL, *mmcrd = NULL, *mmchrg = NULL,
185 *qmgrad, *mmgrad = NULL, energy;
191 snew(qmcrd, 3*(qm->nrQMatoms));
192 snew(qmgrad, 3*(qm->nrQMatoms));
193 /* copy the data from qr into the arrays that are going to be used
194 * in the fortran routines of MOPAC
196 for (i = 0; i < qm->nrQMatoms; i++)
198 for (j = 0; j < DIM; j++)
200 qmcrd[3*i+j] = (double)qm->xQM[i][j]*10;
205 /* later we will add the point charges here. There are some
206 * conceptual problems with semi-empirical QM in combination with
207 * point charges that we need to solve first....
209 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
213 /* now compute the energy and the gradients.
215 snew(qmchrg, qm->nrQMatoms);
217 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
218 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
219 /* add the gradients to the f[] array, and also to the fshift[].
220 * the mopac gradients are in kCal/angstrom.
222 for (i = 0; i < qm->nrQMatoms; i++)
224 for (j = 0; j < DIM; j++)
226 f[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
227 fshift[i][j] = (real)10*CAL2JOULE*qmgrad[3*i+j];
230 QMener = (real)CAL2JOULE*energy;
235 } /* call_mopac_SH */
239 gmx_qmmm_mopac_empty;