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