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