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