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