Apply clang-format to source tree
[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,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "qm_mopac.h"
40
41 #include "config.h"
42
43 #include <stdio.h>
44 #include <stdlib.h>
45 #include <string.h>
46
47 #include <cmath>
48
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"
58
59
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"
64
65 #if GMX_QMMM_MOPAC
66 /* mopac interface routines */
67 void F77_FUNC(domldt, DOMLDT)(int* nrqmat, int labels[], char keywords[]);
68
69 void F77_FUNC(domop, DOMOP)(int*    nrqmat,
70                             double  qmcrd[],
71                             int*    nrmmat,
72                             double  mmchrg[],
73                             double  mmcrd[],
74                             double  qmgrad[],
75                             double  mmgrad[],
76                             double* energy,
77                             double  qmcharges[]);
78
79 #else /* GMX_QMMM_MOPAC */
80 // Stub definitions to make compilation succeed when not configured
81 // for MOPAC support. In that case, the module gives a fatal error
82 // when the initialization function is called, so there is no need to
83 // issue fatal errors here, because that introduces problems with
84 // tools suggesting and prohibiting noreturn attributes.
85
86 static void F77_FUNC(domldt, DOMLDT)(int* /*unused*/, int /*unused*/[], char /*unused*/[]) {}
87
88 static void F77_FUNC(domop, DOMOP)(int* /*unused*/,
89                                    double /*unused*/[],
90                                    int* /*unused*/,
91                                    double /*unused*/[],
92                                    double /*unused*/[],
93                                    double /*unused*/[],
94                                    double /*unused*/[],
95                                    double* /*unused*/,
96                                    double /*unused*/[])
97 {
98 }
99
100 #endif
101
102
103 void init_mopac(t_QMrec* qm)
104 {
105     /* initializes the mopac routines ans sets up the semiempirical
106      * computation by calling moldat(). The inline mopac routines can
107      * only perform gradient operations. If one would like to optimize a
108      * structure or find a transition state at PM3 level, gaussian is
109      * used instead.
110      */
111     char* keywords;
112
113     if (!GMX_QMMM_MOPAC)
114     {
115         gmx_fatal(FARGS,
116                   "Cannot call MOPAC unless linked against it. Use cmake -DGMX_QMMM_PROGRAM=MOPAC, "
117                   "and ensure that linking will work correctly.");
118     }
119
120     snew(keywords, 240);
121
122     if (!qm->bSH) /* if rerun then grad should not be done! */
123     {
124         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n", qm->QMcharge,
125                 eQMmethod_names[qm->QMmethod]);
126     }
127     else
128     {
129         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
130                 qm->QMcharge, eQMmethod_names[qm->QMmethod], qm->CASorbitals, qm->CASelectrons / 2);
131     }
132     F77_FUNC(domldt, DOMLDT)(&qm->nrQMatoms, qm->atomicnumberQM, keywords);
133     fprintf(stderr, "keywords are: %s\n", keywords);
134     free(keywords);
135
136 } /* init_mopac */
137
138 real call_mopac(t_QMrec* qm, t_MMrec* mm, rvec f[], rvec fshift[])
139 {
140     /* do the actual QMMM calculation using directly linked mopac subroutines
141      */
142     double /* always double as the MOPAC routines are always compiled in
143               double precission! */
144             *qmcrd  = nullptr,
145             *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr, *qmgrad, *mmgrad = nullptr,
146             energy = 0;
147     int  i, j;
148     real QMener = 0.0;
149     snew(qmcrd, 3 * (qm->nrQMatoms));
150     snew(qmgrad, 3 * (qm->nrQMatoms));
151     /* copy the data from qr into the arrays that are going to be used
152      * in the fortran routines of MOPAC
153      */
154     for (i = 0; i < qm->nrQMatoms; i++)
155     {
156         for (j = 0; j < DIM; j++)
157         {
158             qmcrd[3 * i + j] = static_cast<double>(qm->xQM[i][j]) * 10;
159         }
160     }
161     if (mm->nrMMatoms)
162     {
163         /* later we will add the point charges here. There are some
164          * conceptual problems with semi-empirical QM in combination with
165          * point charges that we need to solve first....
166          */
167         gmx_fatal(FARGS,
168                   "At present only ONIOM is allowed in combination"
169                   " with MOPAC QM subroutines\n");
170     }
171     else
172     {
173         /* now compute the energy and the gradients.
174          */
175
176         snew(qmchrg, qm->nrQMatoms);
177         F77_FUNC(domop, DOMOP)
178         (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms, mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
179         /* add the gradients to the f[] array, and also to the fshift[].
180          * the mopac gradients are in kCal/angstrom.
181          */
182         for (i = 0; i < qm->nrQMatoms; i++)
183         {
184             for (j = 0; j < DIM; j++)
185             {
186                 f[i][j]      = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
187                 fshift[i][j] = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
188             }
189         }
190         QMener = static_cast<real> CAL2JOULE * energy;
191         /* do we do something with the mulliken charges?? */
192
193         free(qmchrg);
194     }
195     free(qmgrad);
196     free(qmcrd);
197     return (QMener);
198 }
199
200 real call_mopac_SH(t_QMrec* qm, t_MMrec* mm, rvec f[], rvec fshift[])
201 {
202     /* do the actual SH QMMM calculation using directly linked mopac
203        subroutines */
204
205     double /* always double as the MOPAC routines are always compiled in
206               double precission! */
207             *qmcrd  = nullptr,
208             *qmchrg = nullptr, *mmcrd = nullptr, *mmchrg = nullptr, *qmgrad, *mmgrad = nullptr,
209             energy = 0;
210     int  i, j;
211     real QMener = 0.0;
212
213     snew(qmcrd, 3 * (qm->nrQMatoms));
214     snew(qmgrad, 3 * (qm->nrQMatoms));
215     /* copy the data from qr into the arrays that are going to be used
216      * in the fortran routines of MOPAC
217      */
218     for (i = 0; i < qm->nrQMatoms; i++)
219     {
220         for (j = 0; j < DIM; j++)
221         {
222             qmcrd[3 * i + j] = static_cast<double>(qm->xQM[i][j]) * 10;
223         }
224     }
225     if (mm->nrMMatoms)
226     {
227         /* later we will add the point charges here. There are some
228          * conceptual problems with semi-empirical QM in combination with
229          * point charges that we need to solve first....
230          */
231         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
232     }
233     else
234     {
235         /* now compute the energy and the gradients.
236          */
237         snew(qmchrg, qm->nrQMatoms);
238
239         F77_FUNC(domop, DOMOP)
240         (&qm->nrQMatoms, qmcrd, &mm->nrMMatoms, mmchrg, mmcrd, qmgrad, mmgrad, &energy, qmchrg);
241         /* add the gradients to the f[] array, and also to the fshift[].
242          * the mopac gradients are in kCal/angstrom.
243          */
244         for (i = 0; i < qm->nrQMatoms; i++)
245         {
246             for (j = 0; j < DIM; j++)
247             {
248                 f[i][j]      = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
249                 fshift[i][j] = static_cast<real>(10) * CAL2JOULE * qmgrad[3 * i + j];
250             }
251         }
252         QMener = static_cast<real> CAL2JOULE * energy;
253     }
254     free(qmgrad);
255     free(qmcrd);
256     return (QMener);
257 } /* call_mopac_SH */
258
259 #pragma GCC diagnostic pop