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