Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / mdlib / qm_orca.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-2008, 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_orca.h"
40
41 #include "config.h"
42
43 #include <cmath>
44 #include <cstdio>
45 #include <cstdlib>
46 #include <cstring>
47
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/qmmm.h"
54 #include "gromacs/mdtypes/forcerec.h"
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
58
59 // When not built in a configuration with QMMM support, much of this
60 // code is unreachable by design. Tell clang not to warn about it.
61 #pragma GCC diagnostic push
62 #pragma GCC diagnostic ignored "-Wmissing-noreturn"
63
64 /* ORCA interface routines */
65
66 void init_orca(t_QMrec* qm)
67 {
68     char* buf;
69     snew(buf, 200);
70
71     if (!GMX_QMMM_ORCA)
72     {
73         gmx_fatal(FARGS,
74                   "Cannot call ORCA unless linked against it. Use cmake -DGMX_QMMM_PROGRAM=ORCA, "
75                   "and ensure that linking will work correctly.");
76     }
77
78     /* ORCA settings on the system */
79     buf = getenv("GMX_QM_ORCA_BASENAME");
80     if (buf)
81     {
82         snew(qm->orca_basename, 200);
83         sscanf(buf, "%s", qm->orca_basename);
84     }
85     else
86     {
87         gmx_fatal(FARGS, "$GMX_QM_ORCA_BASENAME is not set\n");
88     }
89
90     /* ORCA directory on the system */
91     snew(buf, 200);
92     buf = getenv("GMX_ORCA_PATH");
93
94     if (buf)
95     {
96         snew(qm->orca_dir, 200);
97         sscanf(buf, "%s", qm->orca_dir);
98     }
99     else
100     {
101         gmx_fatal(FARGS, "$GMX_ORCA_PATH not set, check manual\n");
102     }
103
104     fprintf(stderr, "Setting ORCA path to: %s...\n", qm->orca_dir);
105     fprintf(stderr, "ORCA initialised...\n\n");
106     /* since we append the output to the BASENAME.out file,
107        we should delete an existent old out-file here. */
108     sprintf(buf, "%s.out", qm->orca_basename);
109     remove(buf);
110 }
111
112
113 static void write_orca_input(const t_forcerec* fr, t_QMrec* qm, t_MMrec* mm)
114 {
115     int        i;
116     t_QMMMrec* QMMMrec;
117     FILE *     out, *pcFile, *addInputFile;
118     char *     buf, *orcaInput, *addInputFilename, *pcFilename;
119
120     QMMMrec = fr->qr;
121
122     /* write the first part of the input-file */
123     snew(orcaInput, 200);
124     sprintf(orcaInput, "%s.inp", qm->orca_basename);
125     out = fopen(orcaInput, "w");
126
127     snew(addInputFilename, 200);
128     sprintf(addInputFilename, "%s.ORCAINFO", qm->orca_basename);
129     addInputFile = fopen(addInputFilename, "r");
130
131     fprintf(out, "#input-file generated by GROMACS\n");
132
133     fprintf(out, "!EnGrad TightSCF\n");
134
135     /* here we include the insertion of the additional orca-input */
136     snew(buf, 200);
137     if (addInputFile != nullptr)
138     {
139         while (!feof(addInputFile))
140         {
141             if (fgets(buf, 200, addInputFile) != nullptr)
142             {
143                 fputs(buf, out);
144             }
145         }
146     }
147     else
148     {
149         gmx_fatal(FARGS, "No information on the calculation given in %s\n", addInputFilename);
150     }
151
152     fclose(addInputFile);
153
154     /* write charge and multiplicity */
155     fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
156
157     /* write the QM coordinates */
158     for (i = 0; i < qm->nrQMatoms; i++)
159     {
160         int atomNr;
161         if (qm->atomicnumberQM[i] == 0)
162         {
163             atomNr = 1;
164         }
165         else
166         {
167             atomNr = qm->atomicnumberQM[i];
168         }
169         fprintf(out, "%3d %10.7f  %10.7f  %10.7f\n", atomNr, qm->xQM[i][XX] / 0.1,
170                 qm->xQM[i][YY] / 0.1, qm->xQM[i][ZZ] / 0.1);
171     }
172     fprintf(out, "*\n");
173
174     /* write the MM point charge data */
175     if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
176     {
177         /* name of the point charge file */
178         snew(pcFilename, 200);
179         sprintf(pcFilename, "%s.pc", qm->orca_basename);
180         fprintf(out, "%s%s%s\n", "%pointcharges \"", pcFilename, "\"");
181         pcFile = fopen(pcFilename, "w");
182         fprintf(pcFile, "%d\n", mm->nrMMatoms);
183         for (i = 0; i < mm->nrMMatoms; i++)
184         {
185             fprintf(pcFile, "%8.4f %10.7f  %10.7f  %10.7f\n", mm->MMcharges[i],
186                     mm->xMM[i][XX] / 0.1, mm->xMM[i][YY] / 0.1, mm->xMM[i][ZZ] / 0.1);
187         }
188         fprintf(pcFile, "\n");
189         fclose(pcFile);
190     }
191     fprintf(out, "\n");
192
193     fclose(out);
194 } /* write_orca_input */
195
196 static real read_orca_output(rvec QMgrad[], rvec MMgrad[], const t_forcerec* fr, t_QMrec* qm, t_MMrec* mm)
197 {
198     int        i, j;
199     char       buf[300], orca_pcgradFilename[300], orca_engradFilename[300];
200     real       QMener;
201     FILE *     pcgrad, *engrad;
202     int        k;
203     t_QMMMrec* QMMMrec;
204     QMMMrec = fr->qr;
205
206     /* the energy and gradients for the QM part are stored in the engrad file
207      * and the gradients for the point charges are stored in the pc file.
208      */
209     sprintf(orca_engradFilename, "%s.engrad", qm->orca_basename);
210     engrad = fopen(orca_engradFilename, "r");
211     /* we read the energy and the gradient for the qm-atoms from the engrad file
212      */
213     /* we can skip the first seven lines
214      */
215     for (j = 0; j < 7; j++)
216     {
217         if (fgets(buf, 300, engrad) == nullptr)
218         {
219             gmx_fatal(FARGS, "Unexpected end of ORCA output");
220         }
221     }
222     /* now comes the energy
223      */
224     if (fgets(buf, 300, engrad) == nullptr)
225     {
226         gmx_fatal(FARGS, "Unexpected end of ORCA output");
227     }
228 #if GMX_DOUBLE
229     sscanf(buf, "%lf\n", &QMener);
230 #else
231     sscanf(buf, "%f\n", &QMener);
232 #endif
233     /* we can skip the next three lines
234      */
235     for (j = 0; j < 3; j++)
236     {
237         if (fgets(buf, 300, engrad) == nullptr)
238         {
239             gmx_fatal(FARGS, "Unexpected end of ORCA output");
240         }
241     }
242     /* next lines contain the gradients of the QM atoms
243      * now comes the gradient, one value per line:
244      * (atom1 x \n atom1 y \n atom1 z \n atom2 x ...
245      */
246
247     for (i = 0; i < 3 * qm->nrQMatoms; i++)
248     {
249         k = i / 3;
250         if (fgets(buf, 300, engrad) == nullptr)
251         {
252             gmx_fatal(FARGS, "Unexpected end of ORCA output");
253         }
254 #if GMX_DOUBLE
255         if (i % 3 == 0)
256         {
257             sscanf(buf, "%lf\n", &QMgrad[k][XX]);
258         }
259         else if (i % 3 == 1)
260         {
261             sscanf(buf, "%lf\n", &QMgrad[k][YY]);
262         }
263         else if (i % 3 == 2)
264         {
265             sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
266         }
267 #else
268         if (i % 3 == 0)
269         {
270             sscanf(buf, "%f\n", &QMgrad[k][XX]);
271         }
272         else if (i % 3 == 1)
273         {
274             sscanf(buf, "%f\n", &QMgrad[k][YY]);
275         }
276         else if (i % 3 == 2)
277         {
278             sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
279         }
280 #endif
281     }
282     fclose(engrad);
283     /* write the MM point charge data
284      */
285     if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
286     {
287         sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
288         pcgrad = fopen(orca_pcgradFilename, "r");
289
290         /* we read the gradient for the mm-atoms from the pcgrad file
291          */
292         /* we can skip the first line
293          */
294         if (fgets(buf, 300, pcgrad) == nullptr)
295         {
296             gmx_fatal(FARGS, "Unexpected end of ORCA output");
297         }
298         for (i = 0; i < mm->nrMMatoms; i++)
299         {
300             if (fgets(buf, 300, pcgrad) == nullptr)
301             {
302                 gmx_fatal(FARGS, "Unexpected end of ORCA output");
303             }
304 #if GMX_DOUBLE
305             sscanf(buf, "%lf%lf%lf\n", &MMgrad[i][XX], &MMgrad[i][YY], &MMgrad[i][ZZ]);
306 #else
307             sscanf(buf, "%f%f%f\n", &MMgrad[i][XX], &MMgrad[i][YY], &MMgrad[i][ZZ]);
308 #endif
309         }
310         fclose(pcgrad);
311     }
312     return (QMener);
313 }
314
315 static void do_orca(char* orca_dir, char* basename)
316 {
317
318     /* make the call to the orca binary through system()
319      * The location of the binary is set through the
320      * environment.
321      */
322     char buf[100];
323     sprintf(buf, "%s/%s %s.inp >> %s.out", orca_dir, "orca", basename, basename);
324     fprintf(stderr, "Calling '%s'\n", buf);
325     if (system(buf) != 0)
326     {
327         gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
328     }
329 }
330
331 real call_orca(const t_forcerec* fr, t_QMrec* qm, t_MMrec* mm, rvec f[], rvec fshift[])
332 {
333     /* normal orca jobs */
334     static int step = 0;
335     int        i, j;
336     real       QMener;
337     rvec *     QMgrad, *MMgrad;
338     char*      exe;
339
340     snew(exe, 30);
341     sprintf(exe, "%s", "orca");
342     snew(QMgrad, qm->nrQMatoms);
343     snew(MMgrad, mm->nrMMatoms);
344
345     write_orca_input(fr, qm, mm);
346     do_orca(qm->orca_dir, qm->orca_basename);
347     QMener = read_orca_output(QMgrad, MMgrad, fr, qm, mm);
348     /* put the QMMM forces in the force array and to the fshift
349      */
350     for (i = 0; i < qm->nrQMatoms; i++)
351     {
352         for (j = 0; j < DIM; j++)
353         {
354             f[i][j]      = HARTREE_BOHR2MD * QMgrad[i][j];
355             fshift[i][j] = HARTREE_BOHR2MD * QMgrad[i][j];
356         }
357     }
358     for (i = 0; i < mm->nrMMatoms; i++)
359     {
360         for (j = 0; j < DIM; j++)
361         {
362             f[i + qm->nrQMatoms][j]      = HARTREE_BOHR2MD * MMgrad[i][j];
363             fshift[i + qm->nrQMatoms][j] = HARTREE_BOHR2MD * MMgrad[i][j];
364         }
365     }
366     QMener = QMener * HARTREE2KJ * AVOGADRO;
367     step++;
368     free(exe);
369     return (QMener);
370 } /* call_orca */
371
372 /* end of orca sub routines */
373
374 #pragma GCC diagnostic pop