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