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,2017,2018, 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.
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/gmxlib/nrnb.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdlib/ns.h"
55 #include "gromacs/mdlib/qmmm.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
61 // When not built in a configuration with QMMM support, much of this
62 // code is unreachable by design. Tell clang not to warn about it.
63 #pragma GCC diagnostic push
64 #pragma GCC diagnostic ignored "-Wmissing-noreturn"
67 /* mopac interface routines */
69 F77_FUNC(domldt, DOMLDT) (int *nrqmat, int labels[], char keywords[]);
72 F77_FUNC(domop, DOMOP) (int *nrqmat, double qmcrd[], int *nrmmat,
73 double mmchrg[], double mmcrd[], double qmgrad[],
74 double mmgrad[], double *energy, double qmcharges[]);
76 #else /* GMX_QMMM_MOPAC */
77 // Stub definitions to make compilation succeed when not configured
78 // for MOPAC support. In that case, the module gives a fatal error
79 // when the initialization function is called, so there is no need to
80 // issue fatal errors here, because that introduces problems with
81 // tools suggesting and prohibiting noreturn attributes.
83 static void F77_FUNC(domldt, DOMLDT) (int * /*unused*/, int /*unused*/[], char /*unused*/[])
87 static void F77_FUNC(domop, DOMOP) (int * /*unused*/, double /*unused*/[], int * /*unused*/,
88 double /*unused*/[], double /*unused*/[], double /*unused*/[],
89 double /*unused*/[], double * /*unused*/, double /*unused*/[])
96 void init_mopac(t_QMrec *qm)
98 /* initializes the mopac routines ans sets up the semiempirical
99 * computation by calling moldat(). The inline mopac routines can
100 * only perform gradient operations. If one would like to optimize a
101 * structure or find a transition state at PM3 level, gaussian is
109 gmx_fatal(FARGS, "Cannot call MOPAC unless linked against it. Use cmake -DGMX_QMMM_PROGRAM=MOPAC, and ensure that linking will work correctly.");
114 if (!qm->bSH) /* if rerun then grad should not be done! */
116 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
118 eQMmethod_names[qm->QMmethod]);
122 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
124 eQMmethod_names[qm->QMmethod],
125 qm->CASorbitals, qm->CASelectrons/2);
127 F77_FUNC(domldt, DOMLDT) (&qm->nrQMatoms, qm->atomicnumberQM, keywords);
128 fprintf(stderr, "keywords are: %s\n", keywords);
133 real call_mopac(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
135 /* do the actual QMMM calculation using directly linked mopac subroutines
137 double /* always double as the MOPAC routines are always compiled in
138 double precission! */
139 *qmcrd = nullptr, *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr,
140 *qmgrad, *mmgrad = nullptr, energy = 0;
145 snew(qmcrd, 3*(qm->nrQMatoms));
146 snew(qmgrad, 3*(qm->nrQMatoms));
147 /* copy the data from qr into the arrays that are going to be used
148 * in the fortran routines of MOPAC
150 for (i = 0; i < qm->nrQMatoms; i++)
152 for (j = 0; j < DIM; j++)
154 qmcrd[3*i+j] = static_cast<double>(qm->xQM[i][j])*10;
159 /* later we will add the point charges here. There are some
160 * conceptual problems with semi-empirical QM in combination with
161 * point charges that we need to solve first....
163 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
164 " with MOPAC QM subroutines\n");
168 /* now compute the energy and the gradients.
171 snew(qmchrg, qm->nrQMatoms);
172 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
173 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
174 /* add the gradients to the f[] array, and also to the fshift[].
175 * the mopac gradients are in kCal/angstrom.
177 for (i = 0; i < qm->nrQMatoms; i++)
179 for (j = 0; j < DIM; j++)
181 f[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
182 fshift[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
185 QMener = static_cast<real>CAL2JOULE*energy;
186 /* do we do something with the mulliken charges?? */
195 real call_mopac_SH(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
197 /* do the actual SH QMMM calculation using directly linked mopac
200 double /* always double as the MOPAC routines are always compiled in
201 double precission! */
202 *qmcrd = nullptr, *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr,
203 *qmgrad, *mmgrad = nullptr, energy = 0;
209 snew(qmcrd, 3*(qm->nrQMatoms));
210 snew(qmgrad, 3*(qm->nrQMatoms));
211 /* copy the data from qr into the arrays that are going to be used
212 * in the fortran routines of MOPAC
214 for (i = 0; i < qm->nrQMatoms; i++)
216 for (j = 0; j < DIM; j++)
218 qmcrd[3*i+j] = static_cast<double>(qm->xQM[i][j])*10;
223 /* later we will add the point charges here. There are some
224 * conceptual problems with semi-empirical QM in combination with
225 * point charges that we need to solve first....
227 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
231 /* now compute the energy and the gradients.
233 snew(qmchrg, qm->nrQMatoms);
235 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
236 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
237 /* add the gradients to the f[] array, and also to the fshift[].
238 * the mopac gradients are in kCal/angstrom.
240 for (i = 0; i < qm->nrQMatoms; i++)
242 for (j = 0; j < DIM; j++)
244 f[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
245 fshift[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
248 QMener = static_cast<real>CAL2JOULE*energy;
253 } /* call_mopac_SH */
255 #pragma GCC diagnostic pop