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