Remove group scheme search code
[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
68     F77_FUNC(domldt, DOMLDT) (int *nrqmat, int labels[], char keywords[]);
69
70 void
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[]);
74
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.
81
82 static void F77_FUNC(domldt, DOMLDT) (int * /*unused*/, int  /*unused*/[], char  /*unused*/[])
83 {
84 }
85
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*/[])
89 {
90 }
91
92 #endif
93
94
95 void init_mopac(t_QMrec *qm)
96 {
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
101      * used instead.
102      */
103     char
104     *keywords;
105
106     if (!GMX_QMMM_MOPAC)
107     {
108         gmx_fatal(FARGS, "Cannot call MOPAC unless linked against it. Use cmake -DGMX_QMMM_PROGRAM=MOPAC, and ensure that linking will work correctly.");
109     }
110
111     snew(keywords, 240);
112
113     if (!qm->bSH)  /* if rerun then grad should not be done! */
114     {
115         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d GRAD MMOK ANALYT %s\n",
116                 qm->QMcharge,
117                 eQMmethod_names[qm->QMmethod]);
118     }
119     else
120     {
121         sprintf(keywords, "PRECISE GEO-OK CHARGE=%d SINGLET GRAD %s C.I.=(%d,%d) root=2 MECI \n",
122                 qm->QMcharge,
123                 eQMmethod_names[qm->QMmethod],
124                 qm->CASorbitals, qm->CASelectrons/2);
125     }
126     F77_FUNC(domldt, DOMLDT) (&qm->nrQMatoms, qm->atomicnumberQM, keywords);
127     fprintf(stderr, "keywords are: %s\n", keywords);
128     free(keywords);
129
130 } /* init_mopac */
131
132 real call_mopac(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
133 {
134     /* do the actual QMMM calculation using directly linked mopac subroutines
135      */
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;
140     int
141         i, j;
142     real
143         QMener = 0.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
148      */
149     for (i = 0; i < qm->nrQMatoms; i++)
150     {
151         for (j = 0; j < DIM; j++)
152         {
153             qmcrd[3*i+j] = static_cast<double>(qm->xQM[i][j])*10;
154         }
155     }
156     if (mm->nrMMatoms)
157     {
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....
161          */
162         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination"
163                   " with MOPAC QM subroutines\n");
164     }
165     else
166     {
167         /* now compute the energy and the gradients.
168          */
169
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.
175          */
176         for (i = 0; i < qm->nrQMatoms; i++)
177         {
178             for (j = 0; j < DIM; j++)
179             {
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];
182             }
183         }
184         QMener = static_cast<real>CAL2JOULE*energy;
185         /* do we do something with the mulliken charges?? */
186
187         free(qmchrg);
188     }
189     free(qmgrad);
190     free(qmcrd);
191     return (QMener);
192 }
193
194 real call_mopac_SH(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
195 {
196     /* do the actual SH QMMM calculation using directly linked mopac
197        subroutines */
198
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;
203     int
204         i, j;
205     real
206         QMener = 0.0;
207
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
212      */
213     for (i = 0; i < qm->nrQMatoms; i++)
214     {
215         for (j = 0; j < DIM; j++)
216         {
217             qmcrd[3*i+j] = static_cast<double>(qm->xQM[i][j])*10;
218         }
219     }
220     if (mm->nrMMatoms)
221     {
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....
225          */
226         gmx_fatal(FARGS, "At present only ONIOM is allowed in combination with MOPAC\n");
227     }
228     else
229     {
230         /* now compute the energy and the gradients.
231          */
232         snew(qmchrg, qm->nrQMatoms);
233
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.
238          */
239         for (i = 0; i < qm->nrQMatoms; i++)
240         {
241             for (j = 0; j < DIM; j++)
242             {
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];
245             }
246         }
247         QMener = static_cast<real>CAL2JOULE*energy;
248     }
249     free(qmgrad);
250     free(qmcrd);
251     return (QMener);
252 } /* call_mopac_SH */
253
254 #pragma GCC diagnostic pop