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,2019, 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/qmmm.h"
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
60 // When not built in a configuration with QMMM support, much of this
61 // code is unreachable by design. Tell clang not to warn about it.
62 #pragma GCC diagnostic push
63 #pragma GCC diagnostic ignored "-Wmissing-noreturn"
66 /* mopac interface routines */
68 F77_FUNC(domldt, DOMLDT) (int *nrqmat, int labels[], char keywords[]);
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[]);
75 #else /* GMX_QMMM_MOPAC */
76 // Stub definitions to make compilation succeed when not configured
77 // for MOPAC support. In that case, the module gives a fatal error
78 // when the initialization function is called, so there is no need to
79 // issue fatal errors here, because that introduces problems with
80 // tools suggesting and prohibiting noreturn attributes.
82 static void F77_FUNC(domldt, DOMLDT) (int * /*unused*/, int /*unused*/[], char /*unused*/[])
86 static void F77_FUNC(domop, DOMOP) (int * /*unused*/, double /*unused*/[], int * /*unused*/,
87 double /*unused*/[], double /*unused*/[], double /*unused*/[],
88 double /*unused*/[], double * /*unused*/, double /*unused*/[])
95 void init_mopac(t_QMrec *qm)
97 /* initializes the mopac routines ans sets up the semiempirical
98 * computation by calling moldat(). The inline mopac routines can
99 * only perform gradient operations. If one would like to optimize a
100 * structure or find a transition state at PM3 level, gaussian is
108 gmx_fatal(FARGS, "Cannot call MOPAC unless linked against it. Use cmake -DGMX_QMMM_PROGRAM=MOPAC, and ensure that linking will work correctly.");
113 if (!qm->bSH) /* if rerun then grad should not be done! */
115 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
117 eQMmethod_names[qm->QMmethod]);
121 sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
123 eQMmethod_names[qm->QMmethod],
124 qm->CASorbitals, qm->CASelectrons/2);
126 F77_FUNC(domldt, DOMLDT) (&qm->nrQMatoms, qm->atomicnumberQM, keywords);
127 fprintf(stderr, "keywords are: %s\n", keywords);
132 real call_mopac(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
134 /* do the actual QMMM calculation using directly linked mopac subroutines
136 double /* always double as the MOPAC routines are always compiled in
137 double precission! */
138 *qmcrd = nullptr, *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr,
139 *qmgrad, *mmgrad = nullptr, energy = 0;
144 snew(qmcrd, 3*(qm->nrQMatoms));
145 snew(qmgrad, 3*(qm->nrQMatoms));
146 /* copy the data from qr into the arrays that are going to be used
147 * in the fortran routines of MOPAC
149 for (i = 0; i < qm->nrQMatoms; i++)
151 for (j = 0; j < DIM; j++)
153 qmcrd[3*i+j] = static_cast<double>(qm->xQM[i][j])*10;
158 /* later we will add the point charges here. There are some
159 * conceptual problems with semi-empirical QM in combination with
160 * point charges that we need to solve first....
162 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
163 " with MOPAC QM subroutines\n");
167 /* now compute the energy and the gradients.
170 snew(qmchrg, qm->nrQMatoms);
171 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
172 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
173 /* add the gradients to the f[] array, and also to the fshift[].
174 * the mopac gradients are in kCal/angstrom.
176 for (i = 0; i < qm->nrQMatoms; i++)
178 for (j = 0; j < DIM; j++)
180 f[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
181 fshift[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
184 QMener = static_cast<real>CAL2JOULE*energy;
185 /* do we do something with the mulliken charges?? */
194 real call_mopac_SH(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
196 /* do the actual SH QMMM calculation using directly linked mopac
199 double /* always double as the MOPAC routines are always compiled in
200 double precission! */
201 *qmcrd = nullptr, *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr,
202 *qmgrad, *mmgrad = nullptr, energy = 0;
208 snew(qmcrd, 3*(qm->nrQMatoms));
209 snew(qmgrad, 3*(qm->nrQMatoms));
210 /* copy the data from qr into the arrays that are going to be used
211 * in the fortran routines of MOPAC
213 for (i = 0; i < qm->nrQMatoms; i++)
215 for (j = 0; j < DIM; j++)
217 qmcrd[3*i+j] = static_cast<double>(qm->xQM[i][j])*10;
222 /* later we will add the point charges here. There are some
223 * conceptual problems with semi-empirical QM in combination with
224 * point charges that we need to solve first....
226 gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
230 /* now compute the energy and the gradients.
232 snew(qmchrg, qm->nrQMatoms);
234 F77_FUNC(domop, DOMOP) (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms,
235 mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
236 /* add the gradients to the f[] array, and also to the fshift[].
237 * the mopac gradients are in kCal/angstrom.
239 for (i = 0; i < qm->nrQMatoms; i++)
241 for (j = 0; j < DIM; j++)
243 f[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
244 fshift[i][j] = static_cast<real>(10)*CAL2JOULE*qmgrad[3*i+j];
247 QMener = static_cast<real>CAL2JOULE*energy;
252 } /* call_mopac_SH */
254 #pragma GCC diagnostic pop