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