Tidy: modernize-use-nullptr
[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, 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 <math.h>
40 #include <stdio.h>
41 #include <stdlib.h>
42 #include <string.h>
43
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/gmxlib/network.h"
46 #include "gromacs/gmxlib/nrnb.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/force.h"
50 #include "gromacs/mdlib/ns.h"
51 #include "gromacs/mdlib/qmmm.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
55
56 /* ORCA interface routines */
57
58 void init_orca(t_QMrec *qm)
59 {
60     char *buf;
61     snew(buf, 200);
62
63     /* ORCA settings on the system */
64     buf = getenv("GMX_QM_ORCA_BASENAME");
65     if (buf)
66     {
67         snew(qm->orca_basename, 200);
68         sscanf(buf, "%s", qm->orca_basename);
69     }
70     else
71     {
72         gmx_fatal(FARGS, "$GMX_QM_ORCA_BASENAME is not set\n");
73     }
74
75     /* ORCA directory on the system */
76     snew(buf, 200);
77     buf = getenv("GMX_ORCA_PATH");
78
79     if (buf)
80     {
81         snew(qm->orca_dir, 200);
82         sscanf(buf, "%s", qm->orca_dir);
83     }
84     else
85     {
86         gmx_fatal(FARGS, "$GMX_ORCA_PATH not set, check manual\n");
87     }
88
89     fprintf(stderr, "Setting ORCA path to: %s...\n", qm->orca_dir);
90     fprintf(stderr, "ORCA initialised...\n\n");
91     /* since we append the output to the BASENAME.out file,
92        we should delete an existent old out-file here. */
93     sprintf(buf, "%s.out", qm->orca_basename);
94     remove(buf);
95 }
96
97
98 void write_orca_input(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
99 {
100     int        i;
101     t_QMMMrec *QMMMrec;
102     FILE      *out, *pcFile, *addInputFile, *LJCoeff;
103     char      *buf, *orcaInput, *addInputFilename, *LJCoeffFilename, *pcFilename, *exclInName, *exclOutName;
104
105     QMMMrec = fr->qr;
106
107     /* write the first part of the input-file */
108     snew(orcaInput, 200);
109     sprintf(orcaInput, "%s.inp", qm->orca_basename);
110     out = fopen(orcaInput, "w");
111
112     snew(addInputFilename, 200);
113     sprintf(addInputFilename, "%s.ORCAINFO", qm->orca_basename);
114     addInputFile = fopen(addInputFilename, "r");
115
116     fprintf(out, "#input-file generated by GROMACS\n");
117
118     if (qm->bTS)
119     {
120         fprintf(out, "!QMMMOpt TightSCF\n");
121         fprintf(out, "%s\n", "%geom TS_Search EF end");
122     }
123     else if (qm->bOPT)
124     {
125         fprintf(out, "!QMMMOpt TightSCF\n");
126     }
127     else
128     {
129         fprintf(out, "!EnGrad TightSCF\n");
130     }
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     if (qm->bTS || qm->bOPT)
152     {
153         /* freeze the frontier QM atoms and Link atoms. This is
154          * important only if a full QM subsystem optimization is done
155          * with a frozen MM environmeent. For dynamics, or gromacs's own
156          * optimization routines this is not important.
157          */
158         /* ORCA reads the exclusions from LJCoeffFilename.Excl,
159          * so we have to rename the file
160          */
161         int didStart = 0;
162         for (i = 0; i < qm->nrQMatoms; i++)
163         {
164             if (qm->frontatoms[i])
165             {
166                 if (!didStart)
167                 {
168                     fprintf(out, "%s\n", "%geom");
169                     fprintf(out, "   Constraints \n");
170                     didStart = 1;
171                 }
172                 fprintf(out, "        {C %d C}\n", i); /* counting from 0 */
173             }
174         }
175         if (didStart)
176         {
177             fprintf(out, "     end\n   end\n");
178         }
179         /* make a file with information on the C6 and C12 coefficients */
180         if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
181         {
182             snew(exclInName, 200);
183             snew(exclOutName, 200);
184             sprintf(exclInName, "QMMMexcl.dat");
185             sprintf(exclOutName, "%s.LJ.Excl", qm->orca_basename);
186             rename(exclInName, exclOutName);
187             snew(LJCoeffFilename, 200);
188             sprintf(LJCoeffFilename, "%s.LJ", qm->orca_basename);
189             fprintf(out, "%s%s%s\n", "%LJCOEFFICIENTS \"", LJCoeffFilename, "\"");
190             /* make a file with information on the C6 and C12 coefficients */
191             LJCoeff = fopen(LJCoeffFilename, "w");
192             fprintf(LJCoeff, "%d\n", qm->nrQMatoms);
193             for (i = 0; i < qm->nrQMatoms; i++)
194             {
195                 fprintf(LJCoeff, "%10.7f  %10.7f\n", qm->c6[i], qm->c12[i]);
196             }
197             fprintf(LJCoeff, "%d\n", mm->nrMMatoms);
198             for (i = 0; i < mm->nrMMatoms; i++)
199             {
200                 fprintf(LJCoeff, "%10.7f  %10.7f\n", mm->c6[i], mm->c12[i]);
201             }
202             fclose(LJCoeff);
203         }
204     }
205
206     /* write charge and multiplicity */
207     fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
208
209     /* write the QM coordinates */
210     for (i = 0; i < qm->nrQMatoms; i++)
211     {
212         int atomNr;
213         if (qm->atomicnumberQM[i] == 0)
214         {
215             atomNr = 1;
216         }
217         else
218         {
219             atomNr = qm->atomicnumberQM[i];
220         }
221         fprintf(out, "%3d %10.7f  %10.7f  %10.7f\n",
222                 atomNr,
223                 qm->xQM[i][XX]/0.1,
224                 qm->xQM[i][YY]/0.1,
225                 qm->xQM[i][ZZ]/0.1);
226     }
227     fprintf(out, "*\n");
228
229     /* write the MM point charge data */
230     if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
231     {
232         /* name of the point charge file */
233         snew(pcFilename, 200);
234         sprintf(pcFilename, "%s.pc", qm->orca_basename);
235         fprintf(out, "%s%s%s\n", "%pointcharges \"", pcFilename, "\"");
236         pcFile = fopen(pcFilename, "w");
237         fprintf(pcFile, "%d\n", mm->nrMMatoms);
238         for (i = 0; i < mm->nrMMatoms; i++)
239         {
240             fprintf(pcFile, "%8.4f %10.7f  %10.7f  %10.7f\n",
241                     mm->MMcharges[i],
242                     mm->xMM[i][XX]/0.1,
243                     mm->xMM[i][YY]/0.1,
244                     mm->xMM[i][ZZ]/0.1);
245         }
246         fprintf(pcFile, "\n");
247         fclose(pcFile);
248     }
249     fprintf(out, "\n");
250
251     fclose(out);
252 }  /* write_orca_input */
253
254 real read_orca_output(rvec QMgrad[], rvec MMgrad[], t_forcerec *fr,
255                       t_QMrec *qm, t_MMrec *mm)
256 {
257     int
258         i, j, atnum;
259     char
260         buf[300], orca_xyzFilename[300], orca_pcgradFilename[300], orca_engradFilename[300];
261     real
262         QMener;
263     FILE
264        *xyz, *pcgrad, *engrad;
265     int k;
266     t_QMMMrec
267        *QMMMrec;
268     QMMMrec = fr->qr;
269     /* in case of an optimization, the coordinates are printed in the
270      * xyz file, the energy and gradients for the QM part are stored in the engrad file
271      * and the gradients for the point charges are stored in the pc file.
272      */
273
274     /* we need the new xyz coordinates of the QM atoms only for separate QM-optimization
275      */
276
277     if (qm->bTS || qm->bOPT)
278     {
279         sprintf(orca_xyzFilename, "%s.xyz", qm->orca_basename);
280         xyz = fopen(orca_xyzFilename, "r");
281         if (fgets(buf, 300, xyz) == nullptr)
282         {
283             gmx_fatal(FARGS, "Unexpected end of ORCA output");
284         }
285         if (fgets(buf, 300, xyz) == nullptr)
286         {
287             gmx_fatal(FARGS, "Unexpected end of ORCA output");
288         }
289         for (i = 0; i < qm->nrQMatoms; i++)
290         {
291             if (fgets(buf, 300, xyz) == nullptr)
292             {
293                 gmx_fatal(FARGS, "Unexpected end of ORCA output");
294             }
295 #if GMX_DOUBLE
296             sscanf(buf, "%d%lf%lf%lf\n",
297                    &atnum,
298                    &qm->xQM[i][XX],
299                    &qm->xQM[i][YY],
300                    &qm->xQM[i][ZZ]);
301 #else
302             sscanf(buf, "%d%f%f%f\n",
303                    &atnum,
304                    &qm->xQM[i][XX],
305                    &qm->xQM[i][YY],
306                    &qm->xQM[i][ZZ]);
307 #endif
308             for (j = 0; j < DIM; j++)
309             {
310                 qm->xQM[i][j] *= 0.1;
311             }
312         }
313         fclose(xyz);
314     }
315     sprintf(orca_engradFilename, "%s.engrad", qm->orca_basename);
316     engrad = fopen(orca_engradFilename, "r");
317     /* we read the energy and the gradient for the qm-atoms from the engrad file
318      */
319     /* we can skip the first seven lines
320      */
321     for (j = 0; j < 7; j++)
322     {
323         if (fgets(buf, 300, engrad) == nullptr)
324         {
325             gmx_fatal(FARGS, "Unexpected end of ORCA output");
326         }
327     }
328     /* now comes the energy
329      */
330     if (fgets(buf, 300, engrad) == nullptr)
331     {
332         gmx_fatal(FARGS, "Unexpected end of ORCA output");
333     }
334 #if GMX_DOUBLE
335     sscanf(buf, "%lf\n", &QMener);
336 #else
337     sscanf(buf, "%f\n", &QMener);
338 #endif
339     /* we can skip the next three lines
340      */
341     for (j = 0; j < 3; j++)
342     {
343         if (fgets(buf, 300, engrad) == nullptr)
344         {
345             gmx_fatal(FARGS, "Unexpected end of ORCA output");
346         }
347     }
348     /* next lines contain the gradients of the QM atoms
349      * now comes the gradient, one value per line:
350      * (atom1 x \n atom1 y \n atom1 z \n atom2 x ...
351      */
352
353     for (i = 0; i < 3*qm->nrQMatoms; i++)
354     {
355         k = i/3;
356         if (fgets(buf, 300, engrad) == nullptr)
357         {
358             gmx_fatal(FARGS, "Unexpected end of ORCA output");
359         }
360 #if GMX_DOUBLE
361         if (i%3 == 0)
362         {
363             sscanf(buf, "%lf\n", &QMgrad[k][XX]);
364         }
365         else if (i%3 == 1)
366         {
367             sscanf(buf, "%lf\n", &QMgrad[k][YY]);
368         }
369         else if (i%3 == 2)
370         {
371             sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
372         }
373 #else
374         if (i%3 == 0)
375         {
376             sscanf(buf, "%f\n", &QMgrad[k][XX]);
377         }
378         else if (i%3 == 1)
379         {
380             sscanf(buf, "%f\n", &QMgrad[k][YY]);
381         }
382         else if (i%3 == 2)
383         {
384             sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
385         }
386 #endif
387     }
388     fclose(engrad);
389     /* write the MM point charge data
390      */
391     if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
392     {
393         sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
394         pcgrad = fopen(orca_pcgradFilename, "r");
395
396         /* we read the gradient for the mm-atoms from the pcgrad file
397          */
398         /* we can skip the first line
399          */
400         if (fgets(buf, 300, pcgrad) == nullptr)
401         {
402             gmx_fatal(FARGS, "Unexpected end of ORCA output");
403         }
404         for (i = 0; i < mm->nrMMatoms; i++)
405         {
406             if (fgets(buf, 300, pcgrad) == nullptr)
407             {
408                 gmx_fatal(FARGS, "Unexpected end of ORCA output");
409             }
410     #if GMX_DOUBLE
411             sscanf(buf, "%lf%lf%lf\n",
412                    &MMgrad[i][XX],
413                    &MMgrad[i][YY],
414                    &MMgrad[i][ZZ]);
415     #else
416             sscanf(buf, "%f%f%f\n",
417                    &MMgrad[i][XX],
418                    &MMgrad[i][YY],
419                    &MMgrad[i][ZZ]);
420     #endif
421         }
422         fclose(pcgrad);
423     }
424     return(QMener);
425 }
426
427 void do_orca(char *orca_dir, char *basename)
428 {
429
430     /* make the call to the orca binary through system()
431      * The location of the binary is set through the
432      * environment.
433      */
434     char
435         buf[100];
436     sprintf(buf, "%s/%s %s.inp >> %s.out",
437             orca_dir,
438             "orca",
439             basename,
440             basename);
441     fprintf(stderr, "Calling '%s'\n", buf);
442     if (system(buf) != 0)
443     {
444         gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
445     }
446 }
447
448 real call_orca(t_forcerec *fr,
449                t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
450 {
451     /* normal orca jobs */
452     static int
453         step = 0;
454     int
455         i, j;
456     real
457         QMener;
458     rvec
459        *QMgrad, *MMgrad;
460     char
461        *exe;
462
463     snew(exe, 30);
464     sprintf(exe, "%s", "orca");
465     snew(QMgrad, qm->nrQMatoms);
466     snew(MMgrad, mm->nrMMatoms);
467
468     write_orca_input(fr, qm, mm);
469     do_orca(qm->orca_dir, qm->orca_basename);
470     QMener = read_orca_output(QMgrad, MMgrad, fr, qm, mm);
471     /* put the QMMM forces in the force array and to the fshift
472      */
473     for (i = 0; i < qm->nrQMatoms; i++)
474     {
475         for (j = 0; j < DIM; j++)
476         {
477             f[i][j]      = HARTREE_BOHR2MD*QMgrad[i][j];
478             fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
479         }
480     }
481     for (i = 0; i < mm->nrMMatoms; i++)
482     {
483         for (j = 0; j < DIM; j++)
484         {
485             f[i+qm->nrQMatoms][j]      = HARTREE_BOHR2MD*MMgrad[i][j];
486             fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
487         }
488     }
489     QMener = QMener*HARTREE2KJ*AVOGADRO;
490     step++;
491     free(exe);
492     return(QMener);
493 } /* call_orca */
494
495 /* end of orca sub routines */