Remove dysfunctional QMMM interface pt1
[alexxy/gromacs.git] / src / gromacs / mdlib / qm_mopac.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "qm_mopac.h"
41
42 #include "config.h"
43
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <string.h>
47
48 #include <cmath>
49
50 #include "gromacs/fileio/confio.h"
51 #include "gromacs/gmxlib/network.h"
52 #include "gromacs/gmxlib/nrnb.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.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"
59
60
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"
65
66 #if GMX_QMMM_MOPAC
67 /* mopac interface routines */
68 void F77_FUNC(domldt, DOMLDT)(int* nrqmat, int labels[], char keywords[]);
69
70 void F77_FUNC(domop, DOMOP)(int*    nrqmat,
71                             double  qmcrd[],
72                             int*    nrmmat,
73                             double  mmchrg[],
74                             double  mmcrd[],
75                             double  qmgrad[],
76                             double  mmgrad[],
77                             double* energy,
78                             double  qmcharges[]);
79
80 #else /* GMX_QMMM_MOPAC */
81 // Stub definitions to make compilation succeed when not configured
82 // for MOPAC support. In that case, the module gives a fatal error
83 // when the initialization function is called, so there is no need to
84 // issue fatal errors here, because that introduces problems with
85 // tools suggesting and prohibiting noreturn attributes.
86
87 static void F77_FUNC(domldt, DOMLDT)(int* /*unused*/, int /*unused*/[], char /*unused*/[]) {}
88
89 static void F77_FUNC(domop, DOMOP)(int* /*unused*/,
90                                    double /*unused*/[],
91                                    int* /*unused*/,
92                                    double /*unused*/[],
93                                    double /*unused*/[],
94                                    double /*unused*/[],
95                                    double /*unused*/[],
96                                    double* /*unused*/,
97                                    double /*unused*/[])
98 {
99 }
100
101 #endif
102
103
104 void init_mopac(t_QMrec* qm)
105 {
106     /* initializes the mopac routines ans sets up the semiempirical
107      * computation by calling moldat(). The inline mopac routines can
108      * only perform gradient operations. If one would like to optimize a
109      * structure or find a transition state at PM3 level, gaussian is
110      * used instead.
111      */
112     char* keywords;
113
114     if (!GMX_QMMM_MOPAC)
115     {
116         gmx_fatal(FARGS,
117                   "Cannot call MOPAC unless linked against it. Use cmake -DGMX_QMMM_PROGRAM=MOPAC, "
118                   "and ensure that linking will work correctly.");
119     }
120
121     snew(keywords, 240);
122
123     if (!qm->bSH) /* if rerun then grad should not be done! */
124     {
125         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n", qm->QMcharge,
126                 eQMmethod_names[qm->QMmethod]);
127     }
128     else
129     {
130         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
131                 qm->QMcharge, eQMmethod_names[qm->QMmethod], qm->CASorbitals, qm->CASelectrons / 2);
132     }
133     F77_FUNC(domldt, DOMLDT)(&qm->nrQMatoms, qm->atomicnumberQM, keywords);
134     fprintf(stderr, "keywords are: %s\n", keywords);
135     free(keywords);
136
137 } /* init_mopac */
138
139 real call_mopac(t_QMrec* qm, t_MMrec* mm, rvec f[], rvec fshift[])
140 {
141     /* do the actual QMMM calculation using directly linked mopac subroutines
142      */
143     double /* always double as the MOPAC routines are always compiled in
144               double precission! */
145             *qmcrd  = nullptr,
146             *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr, *qmgrad, *mmgrad = nullptr,
147             energy = 0;
148     int  i, j;
149     real QMener = 0.0;
150     snew(qmcrd, 3 * (qm->nrQMatoms));
151     snew(qmgrad, 3 * (qm->nrQMatoms));
152     /* copy the data from qr into the arrays that are going to be used
153      * in the fortran routines of MOPAC
154      */
155     for (i = 0; i < qm->nrQMatoms; i++)
156     {
157         for (j = 0; j < DIM; j++)
158         {
159             qmcrd[3 * i + j] = static_cast<double>(qm->xQM[i][j]) * 10;
160         }
161     }
162     if (mm->nrMMatoms)
163     {
164         /* later we will add the point charges here. There are some
165          * conceptual problems with semi-empirical QM in combination with
166          * point charges that we need to solve first....
167          */
168         gmx_fatal(FARGS,
169                   "At present only ONIOM is allowed in combination"
170                   " with MOPAC QM subroutines\n");
171     }
172     else
173     {
174         /* now compute the energy and the gradients.
175          */
176
177         snew(qmchrg, qm->nrQMatoms);
178         F77_FUNC(domop, DOMOP)
179         (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms, mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
180         /* add the gradients to the f[] array, and also to the fshift[].
181          * the mopac gradients are in kCal/angstrom.
182          */
183         for (i = 0; i < qm->nrQMatoms; i++)
184         {
185             for (j = 0; j < DIM; j++)
186             {
187                 f[i][j]      = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
188                 fshift[i][j] = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
189             }
190         }
191         QMener = static_cast<real> CAL2JOULE * energy;
192         /* do we do something with the mulliken charges?? */
193
194         free(qmchrg);
195     }
196     free(qmgrad);
197     free(qmcrd);
198     return (QMener);
199 }
200
201 real call_mopac_SH(t_QMrec* qm, t_MMrec* mm, rvec f[], rvec fshift[])
202 {
203     /* do the actual SH QMMM calculation using directly linked mopac
204        subroutines */
205
206     double /* always double as the MOPAC routines are always compiled in
207               double precission! */
208             *qmcrd  = nullptr,
209             *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr, *qmgrad, *mmgrad = nullptr,
210             energy = 0;
211     int  i, j;
212     real QMener = 0.0;
213
214     snew(qmcrd, 3 * (qm->nrQMatoms));
215     snew(qmgrad, 3 * (qm->nrQMatoms));
216     /* copy the data from qr into the arrays that are going to be used
217      * in the fortran routines of MOPAC
218      */
219     for (i = 0; i < qm->nrQMatoms; i++)
220     {
221         for (j = 0; j < DIM; j++)
222         {
223             qmcrd[3 * i + j] = static_cast<double>(qm->xQM[i][j]) * 10;
224         }
225     }
226     if (mm->nrMMatoms)
227     {
228         /* later we will add the point charges here. There are some
229          * conceptual problems with semi-empirical QM in combination with
230          * point charges that we need to solve first....
231          */
232         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
233     }
234     else
235     {
236         /* now compute the energy and the gradients.
237          */
238         snew(qmchrg, qm->nrQMatoms);
239
240         F77_FUNC(domop, DOMOP)
241         (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms, mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
242         /* add the gradients to the f[] array, and also to the fshift[].
243          * the mopac gradients are in kCal/angstrom.
244          */
245         for (i = 0; i < qm->nrQMatoms; i++)
246         {
247             for (j = 0; j < DIM; j++)
248             {
249                 f[i][j]      = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
250                 fshift[i][j] = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
251             }
252         }
253         QMener = static_cast<real> CAL2JOULE * energy;
254     }
255     free(qmgrad);
256     free(qmcrd);
257     return (QMener);
258 } /* call_mopac_SH */
259
260 #pragma GCC diagnostic pop