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