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