287c97b64ac0ce5f7a5625a75b55e2bf55e1c874
[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 <stdio.h>
42 #include <stdlib.h>
43 #include <string.h>
44
45 #include <cmath>
46
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/ns.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 /* ORCA interface routines */
60
61 void init_orca(t_QMrec *qm)
62 {
63     char *buf;
64     snew(buf, 200);
65
66     /* ORCA settings on the system */
67     buf = getenv("GMX_QM_ORCA_BASENAME");
68     if (buf)
69     {
70         snew(qm->orca_basename, 200);
71         sscanf(buf, "%s", qm->orca_basename);
72     }
73     else
74     {
75         gmx_fatal(FARGS, "$GMX_QM_ORCA_BASENAME is not set\n");
76     }
77
78     /* ORCA directory on the system */
79     snew(buf, 200);
80     buf = getenv("GMX_ORCA_PATH");
81
82     if (buf)
83     {
84         snew(qm->orca_dir, 200);
85         sscanf(buf, "%s", qm->orca_dir);
86     }
87     else
88     {
89         gmx_fatal(FARGS, "$GMX_ORCA_PATH not set, check manual\n");
90     }
91
92     fprintf(stderr, "Setting ORCA path to: %s...\n", qm->orca_dir);
93     fprintf(stderr, "ORCA initialised...\n\n");
94     /* since we append the output to the BASENAME.out file,
95        we should delete an existent old out-file here. */
96     sprintf(buf, "%s.out", qm->orca_basename);
97     remove(buf);
98 }
99
100
101 static void write_orca_input(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
102 {
103     int        i;
104     t_QMMMrec *QMMMrec;
105     FILE      *out, *pcFile, *addInputFile;
106     char      *buf, *orcaInput, *addInputFilename, *pcFilename;
107
108     QMMMrec = fr->qr;
109
110     /* write the first part of the input-file */
111     snew(orcaInput, 200);
112     sprintf(orcaInput, "%s.inp", qm->orca_basename);
113     out = fopen(orcaInput, "w");
114
115     snew(addInputFilename, 200);
116     sprintf(addInputFilename, "%s.ORCAINFO", qm->orca_basename);
117     addInputFile = fopen(addInputFilename, "r");
118
119     fprintf(out, "#input-file generated by GROMACS\n");
120
121     fprintf(out, "!EnGrad TightSCF\n");
122
123     /* here we include the insertion of the additional orca-input */
124     snew(buf, 200);
125     if (addInputFile != nullptr)
126     {
127         while (!feof(addInputFile))
128         {
129             if (fgets(buf, 200, addInputFile) != nullptr)
130             {
131                 fputs(buf, out);
132             }
133         }
134     }
135     else
136     {
137         gmx_fatal(FARGS, "No information on the calculation given in %s\n", addInputFilename);
138     }
139
140     fclose(addInputFile);
141
142     /* write charge and multiplicity */
143     fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
144
145     /* write the QM coordinates */
146     for (i = 0; i < qm->nrQMatoms; i++)
147     {
148         int atomNr;
149         if (qm->atomicnumberQM[i] == 0)
150         {
151             atomNr = 1;
152         }
153         else
154         {
155             atomNr = qm->atomicnumberQM[i];
156         }
157         fprintf(out, "%3d %10.7f  %10.7f  %10.7f\n",
158                 atomNr,
159                 qm->xQM[i][XX]/0.1,
160                 qm->xQM[i][YY]/0.1,
161                 qm->xQM[i][ZZ]/0.1);
162     }
163     fprintf(out, "*\n");
164
165     /* write the MM point charge data */
166     if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
167     {
168         /* name of the point charge file */
169         snew(pcFilename, 200);
170         sprintf(pcFilename, "%s.pc", qm->orca_basename);
171         fprintf(out, "%s%s%s\n", "%pointcharges \"", pcFilename, "\"");
172         pcFile = fopen(pcFilename, "w");
173         fprintf(pcFile, "%d\n", mm->nrMMatoms);
174         for (i = 0; i < mm->nrMMatoms; i++)
175         {
176             fprintf(pcFile, "%8.4f %10.7f  %10.7f  %10.7f\n",
177                     mm->MMcharges[i],
178                     mm->xMM[i][XX]/0.1,
179                     mm->xMM[i][YY]/0.1,
180                     mm->xMM[i][ZZ]/0.1);
181         }
182         fprintf(pcFile, "\n");
183         fclose(pcFile);
184     }
185     fprintf(out, "\n");
186
187     fclose(out);
188 }  /* write_orca_input */
189
190 static real read_orca_output(rvec QMgrad[], rvec MMgrad[], t_forcerec *fr,
191                              t_QMrec *qm, t_MMrec *mm)
192 {
193     int
194         i, j;
195     char
196         buf[300], orca_pcgradFilename[300], orca_engradFilename[300];
197     real
198         QMener;
199     FILE
200        *pcgrad, *engrad;
201     int k;
202     t_QMMMrec
203        *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",
306                    &MMgrad[i][XX],
307                    &MMgrad[i][YY],
308                    &MMgrad[i][ZZ]);
309     #else
310             sscanf(buf, "%f%f%f\n",
311                    &MMgrad[i][XX],
312                    &MMgrad[i][YY],
313                    &MMgrad[i][ZZ]);
314     #endif
315         }
316         fclose(pcgrad);
317     }
318     return(QMener);
319 }
320
321 static void do_orca(char *orca_dir, char *basename)
322 {
323
324     /* make the call to the orca binary through system()
325      * The location of the binary is set through the
326      * environment.
327      */
328     char
329         buf[100];
330     sprintf(buf, "%s/%s %s.inp >> %s.out",
331             orca_dir,
332             "orca",
333             basename,
334             basename);
335     fprintf(stderr, "Calling '%s'\n", buf);
336     if (system(buf) != 0)
337     {
338         gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
339     }
340 }
341
342 real call_orca(t_forcerec *fr,
343                t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
344 {
345     /* normal orca jobs */
346     static int
347         step = 0;
348     int
349         i, j;
350     real
351         QMener;
352     rvec
353        *QMgrad, *MMgrad;
354     char
355        *exe;
356
357     snew(exe, 30);
358     sprintf(exe, "%s", "orca");
359     snew(QMgrad, qm->nrQMatoms);
360     snew(MMgrad, mm->nrMMatoms);
361
362     write_orca_input(fr, qm, mm);
363     do_orca(qm->orca_dir, qm->orca_basename);
364     QMener = read_orca_output(QMgrad, MMgrad, fr, qm, mm);
365     /* put the QMMM forces in the force array and to the fshift
366      */
367     for (i = 0; i < qm->nrQMatoms; i++)
368     {
369         for (j = 0; j < DIM; j++)
370         {
371             f[i][j]      = HARTREE_BOHR2MD*QMgrad[i][j];
372             fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
373         }
374     }
375     for (i = 0; i < mm->nrMMatoms; i++)
376     {
377         for (j = 0; j < DIM; j++)
378         {
379             f[i+qm->nrQMatoms][j]      = HARTREE_BOHR2MD*MMgrad[i][j];
380             fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
381         }
382     }
383     QMener = QMener*HARTREE2KJ*AVOGADRO;
384     step++;
385     free(exe);
386     return(QMener);
387 } /* call_orca */
388
389 /* end of orca sub routines */