470b9a5b1bb9a99cc5ece89031dc4a01cbd784d8
[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 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_orca.h"
41
42 #include "config.h"
43
44 #include <cmath>
45 #include <cstdio>
46 #include <cstdlib>
47 #include <cstring>
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 // 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_QMMMrec* QMMMrec, t_QMrec* qm, t_MMrec* mm)
114 {
115     int   i;
116     FILE *pcFile, *addInputFile;
117     char *buf, *orcaInput, *addInputFilename, *pcFilename;
118
119     /* write the first part of the input-file */
120     snew(orcaInput, 200);
121     sprintf(orcaInput, "%s.inp", qm->orca_basename);
122     FILE* out = fopen(orcaInput, "w");
123
124     snew(addInputFilename, 200);
125     sprintf(addInputFilename, "%s.ORCAINFO", qm->orca_basename);
126     addInputFile = fopen(addInputFilename, "r");
127
128     fprintf(out, "#input-file generated by GROMACS\n");
129
130     fprintf(out, "!EnGrad TightSCF\n");
131
132     /* here we include the insertion of the additional orca-input */
133     snew(buf, 200);
134     if (addInputFile != nullptr)
135     {
136         while (!feof(addInputFile))
137         {
138             if (fgets(buf, 200, addInputFile) != nullptr)
139             {
140                 fputs(buf, out);
141             }
142         }
143     }
144     else
145     {
146         gmx_fatal(FARGS, "No information on the calculation given in %s\n", addInputFilename);
147     }
148
149     fclose(addInputFile);
150
151     /* write charge and multiplicity */
152     fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
153
154     /* write the QM coordinates */
155     for (i = 0; i < qm->nrQMatoms; i++)
156     {
157         int atomNr;
158         if (qm->atomicnumberQM[i] == 0)
159         {
160             atomNr = 1;
161         }
162         else
163         {
164             atomNr = qm->atomicnumberQM[i];
165         }
166         fprintf(out, "%3d %10.7f  %10.7f  %10.7f\n", atomNr, qm->xQM[i][XX] / 0.1,
167                 qm->xQM[i][YY] / 0.1, qm->xQM[i][ZZ] / 0.1);
168     }
169     fprintf(out, "*\n");
170
171     /* write the MM point charge data */
172     if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
173     {
174         /* name of the point charge file */
175         snew(pcFilename, 200);
176         sprintf(pcFilename, "%s.pc", qm->orca_basename);
177         fprintf(out, "%s%s%s\n", "%pointcharges \"", pcFilename, "\"");
178         pcFile = fopen(pcFilename, "w");
179         fprintf(pcFile, "%d\n", mm->nrMMatoms);
180         for (i = 0; i < mm->nrMMatoms; i++)
181         {
182             fprintf(pcFile, "%8.4f %10.7f  %10.7f  %10.7f\n", mm->MMcharges[i],
183                     mm->xMM[i][XX] / 0.1, mm->xMM[i][YY] / 0.1, mm->xMM[i][ZZ] / 0.1);
184         }
185         fprintf(pcFile, "\n");
186         fclose(pcFile);
187     }
188     fprintf(out, "\n");
189
190     fclose(out);
191 } /* write_orca_input */
192
193 static real read_orca_output(rvec QMgrad[], rvec MMgrad[], const t_QMMMrec* QMMMrec, t_QMrec* qm, t_MMrec* mm)
194 {
195     int   i, j;
196     char  buf[300], orca_pcgradFilename[300], orca_engradFilename[300];
197     real  QMener;
198     FILE *pcgrad, *engrad;
199     int   k;
200
201     /* the energy and gradients for the QM part are stored in the engrad file
202      * and the gradients for the point charges are stored in the pc file.
203      */
204     sprintf(orca_engradFilename, "%s.engrad", qm->orca_basename);
205     engrad = fopen(orca_engradFilename, "r");
206     /* we read the energy and the gradient for the qm-atoms from the engrad file
207      */
208     /* we can skip the first seven lines
209      */
210     for (j = 0; j < 7; j++)
211     {
212         if (fgets(buf, 300, engrad) == nullptr)
213         {
214             gmx_fatal(FARGS, "Unexpected end of ORCA output");
215         }
216     }
217     /* now comes the energy
218      */
219     if (fgets(buf, 300, engrad) == nullptr)
220     {
221         gmx_fatal(FARGS, "Unexpected end of ORCA output");
222     }
223 #if GMX_DOUBLE
224     sscanf(buf, "%lf\n", &QMener);
225 #else
226     sscanf(buf, "%f\n", &QMener);
227 #endif
228     /* we can skip the next three lines
229      */
230     for (j = 0; j < 3; j++)
231     {
232         if (fgets(buf, 300, engrad) == nullptr)
233         {
234             gmx_fatal(FARGS, "Unexpected end of ORCA output");
235         }
236     }
237     /* next lines contain the gradients of the QM atoms
238      * now comes the gradient, one value per line:
239      * (atom1 x \n atom1 y \n atom1 z \n atom2 x ...
240      */
241
242     for (i = 0; i < 3 * qm->nrQMatoms; i++)
243     {
244         k = i / 3;
245         if (fgets(buf, 300, engrad) == nullptr)
246         {
247             gmx_fatal(FARGS, "Unexpected end of ORCA output");
248         }
249 #if GMX_DOUBLE
250         if (i % 3 == 0)
251         {
252             sscanf(buf, "%lf\n", &QMgrad[k][XX]);
253         }
254         else if (i % 3 == 1)
255         {
256             sscanf(buf, "%lf\n", &QMgrad[k][YY]);
257         }
258         else if (i % 3 == 2)
259         {
260             sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
261         }
262 #else
263         if (i % 3 == 0)
264         {
265             sscanf(buf, "%f\n", &QMgrad[k][XX]);
266         }
267         else if (i % 3 == 1)
268         {
269             sscanf(buf, "%f\n", &QMgrad[k][YY]);
270         }
271         else if (i % 3 == 2)
272         {
273             sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
274         }
275 #endif
276     }
277     fclose(engrad);
278     /* write the MM point charge data
279      */
280     if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
281     {
282         sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
283         pcgrad = fopen(orca_pcgradFilename, "r");
284
285         /* we read the gradient for the mm-atoms from the pcgrad file
286          */
287         /* we can skip the first line
288          */
289         if (fgets(buf, 300, pcgrad) == nullptr)
290         {
291             gmx_fatal(FARGS, "Unexpected end of ORCA output");
292         }
293         for (i = 0; i < mm->nrMMatoms; i++)
294         {
295             if (fgets(buf, 300, pcgrad) == nullptr)
296             {
297                 gmx_fatal(FARGS, "Unexpected end of ORCA output");
298             }
299 #if GMX_DOUBLE
300             sscanf(buf, "%lf%lf%lf\n", &MMgrad[i][XX], &MMgrad[i][YY], &MMgrad[i][ZZ]);
301 #else
302             sscanf(buf, "%f%f%f\n", &MMgrad[i][XX], &MMgrad[i][YY], &MMgrad[i][ZZ]);
303 #endif
304         }
305         fclose(pcgrad);
306     }
307     return (QMener);
308 }
309
310 static void do_orca(char* orca_dir, char* basename)
311 {
312
313     /* make the call to the orca binary through system()
314      * The location of the binary is set through the
315      * environment.
316      */
317     char buf[100];
318     sprintf(buf, "%s/%s %s.inp >> %s.out", orca_dir, "orca", basename, basename);
319     fprintf(stderr, "Calling '%s'\n", buf);
320     if (system(buf) != 0)
321     {
322         gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
323     }
324 }
325
326 real call_orca(const t_QMMMrec* qmmm, t_QMrec* qm, t_MMrec* mm, rvec f[], rvec fshift[])
327 {
328     /* normal orca jobs */
329     static int step = 0;
330     int        i, j;
331     real       QMener;
332     rvec *     QMgrad, *MMgrad;
333     char*      exe;
334
335     snew(exe, 30);
336     sprintf(exe, "%s", "orca");
337     snew(QMgrad, qm->nrQMatoms);
338     snew(MMgrad, mm->nrMMatoms);
339
340     write_orca_input(qmmm, qm, mm);
341     do_orca(qm->orca_dir, qm->orca_basename);
342     QMener = read_orca_output(QMgrad, MMgrad, qmmm, qm, mm);
343     /* put the QMMM forces in the force array and to the fshift
344      */
345     for (i = 0; i < qm->nrQMatoms; i++)
346     {
347         for (j = 0; j < DIM; j++)
348         {
349             f[i][j]      = HARTREE_BOHR2MD * QMgrad[i][j];
350             fshift[i][j] = HARTREE_BOHR2MD * QMgrad[i][j];
351         }
352     }
353     for (i = 0; i < mm->nrMMatoms; i++)
354     {
355         for (j = 0; j < DIM; j++)
356         {
357             f[i + qm->nrQMatoms][j]      = HARTREE_BOHR2MD * MMgrad[i][j];
358             fshift[i + qm->nrQMatoms][j] = HARTREE_BOHR2MD * MMgrad[i][j];
359         }
360     }
361     QMener = QMener * HARTREE2KJ * AVOGADRO;
362     step++;
363     free(exe);
364     return (QMener);
365 } /* call_orca */
366
367 /* end of orca sub routines */
368
369 #pragma GCC diagnostic pop