2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2013, 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.
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.
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.
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.
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.
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.
51 #include "gromacs/fileio/confio.h"
63 #include "gmx_fatal.h"
67 /* ORCA interface routines */
69 void init_orca(t_QMrec *qm)
74 /* ORCA settings on the system */
75 buf = getenv("BASENAME");
78 snew(qm->orca_basename, 200);
79 sscanf(buf, "%s", qm->orca_basename);
83 gmx_fatal(FARGS, "$BASENAME not set\n");
86 /* ORCA directory on the system */
88 buf = getenv("ORCA_PATH");
92 snew(qm->orca_dir, 200);
93 sscanf(buf, "%s", qm->orca_dir);
97 gmx_fatal(FARGS, "$ORCA_PATH not set, check manual\n");
100 fprintf(stderr, "Setting ORCA path to: %s...\n", qm->orca_dir);
101 fprintf(stderr, "ORCA initialised...\n\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);
109 void write_orca_input(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
113 FILE *out, *pcFile, *addInputFile, *LJCoeff;
114 char *buf, *orcaInput, *addInputFilename, *LJCoeffFilename, *pcFilename, *exclInName, *exclOutName;
118 /* write the first part of the input-file */
119 snew(orcaInput, 200);
120 sprintf(orcaInput, "%s.inp", qm->orca_basename);
121 out = fopen(orcaInput, "w");
123 snew(addInputFilename, 200);
124 sprintf(addInputFilename, "%s.ORCAINFO", qm->orca_basename);
125 addInputFile = fopen(addInputFilename, "r");
127 fprintf(out, "#input-file generated by GROMACS\n");
131 fprintf(out, "!QMMMOpt TightSCF\n");
132 fprintf(out, "%s\n", "%geom TS_Search EF end");
136 fprintf(out, "!QMMMOpt TightSCF\n");
140 fprintf(out, "!EnGrad TightSCF\n");
143 /* here we include the insertion of the additional orca-input */
145 if (addInputFile != NULL)
147 while (!feof(addInputFile))
149 if (fgets(buf, 200, addInputFile) != NULL)
157 gmx_fatal(FARGS, "No information on the calculation given in %s\n", addInputFilename);
160 fclose(addInputFile);
162 if (qm->bTS || qm->bOPT)
164 /* freeze the frontier QM atoms and Link atoms. This is
165 * important only if a full QM subsystem optimization is done
166 * with a frozen MM environmeent. For dynamics, or gromacs's own
167 * optimization routines this is not important.
169 /* ORCA reads the exclusions from LJCoeffFilename.Excl,
170 * so we have to rename the file
173 for (i = 0; i < qm->nrQMatoms; i++)
175 if (qm->frontatoms[i])
179 fprintf(out, "%s\n", "%geom");
180 fprintf(out, " Constraints \n");
183 fprintf(out, " {C %d C}\n", i); /* counting from 0 */
188 fprintf(out, " end\n end\n");
190 /* make a file with information on the C6 and C12 coefficients */
191 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
193 snew(exclInName, 200);
194 snew(exclOutName, 200);
195 sprintf(exclInName, "QMMMexcl.dat");
196 sprintf(exclOutName, "%s.LJ.Excl", qm->orca_basename);
197 rename(exclInName, exclOutName);
198 snew(LJCoeffFilename, 200);
199 sprintf(LJCoeffFilename, "%s.LJ", qm->orca_basename);
200 fprintf(out, "%s%s%s\n", "%LJCOEFFICIENTS \"", LJCoeffFilename, "\"");
201 /* make a file with information on the C6 and C12 coefficients */
202 LJCoeff = fopen(LJCoeffFilename, "w");
203 fprintf(LJCoeff, "%d\n", qm->nrQMatoms);
204 for (i = 0; i < qm->nrQMatoms; i++)
207 fprintf(LJCoeff, "%10.7lf %10.7lf\n", qm->c6[i], qm->c12[i]);
209 fprintf(LJCoeff, "%10.7f %10.7f\n", qm->c6[i], qm->c12[i]);
212 fprintf(LJCoeff, "%d\n", mm->nrMMatoms);
213 for (i = 0; i < mm->nrMMatoms; i++)
216 fprintf(LJCoeff, "%10.7lf %10.7lf\n", mm->c6[i], mm->c12[i]);
218 fprintf(LJCoeff, "%10.7f %10.7f\n", mm->c6[i], mm->c12[i]);
225 /* write charge and multiplicity */
226 fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
228 /* write the QM coordinates */
229 for (i = 0; i < qm->nrQMatoms; i++)
232 if (qm->atomicnumberQM[i] == 0)
238 atomNr = qm->atomicnumberQM[i];
241 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
247 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
256 /* write the MM point charge data */
257 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
259 /* name of the point charge file */
260 snew(pcFilename, 200);
261 sprintf(pcFilename, "%s.pc", qm->orca_basename);
262 fprintf(out, "%s%s%s\n", "%pointcharges \"", pcFilename, "\"");
263 pcFile = fopen(pcFilename, "w");
264 fprintf(pcFile, "%d\n", mm->nrMMatoms);
265 for (i = 0; i < mm->nrMMatoms; i++)
268 fprintf(pcFile, "%8.4lf %10.7lf %10.7lf %10.7lf\n",
274 fprintf(pcFile, "%8.4f %10.7f %10.7f %10.7f\n",
281 fprintf(pcFile, "\n");
287 } /* write_orca_input */
289 real read_orca_output(rvec QMgrad[], rvec MMgrad[], t_forcerec *fr,
290 t_QMrec *qm, t_MMrec *mm)
295 buf[300], tmp[300], orca_xyzFilename[300], orca_pcgradFilename[300], orca_engradFilename[300];
299 *xyz, *pcgrad, *engrad;
304 /* in case of an optimization, the coordinates are printed in the
305 * xyz file, the energy and gradients for the QM part are stored in the engrad file
306 * and the gradients for the point charges are stored in the pc file.
309 /* we need the new xyz coordinates of the QM atoms only for separate QM-optimization
312 if (qm->bTS || qm->bOPT)
314 sprintf(orca_xyzFilename, "%s.xyz", qm->orca_basename);
315 xyz = fopen(orca_xyzFilename, "r");
316 if (fgets(buf, 300, xyz) == NULL)
318 gmx_fatal(FARGS, "Unexpected end of ORCA output");
320 if (fgets(buf, 300, xyz) == NULL)
322 gmx_fatal(FARGS, "Unexpected end of ORCA output");
324 for (i = 0; i < qm->nrQMatoms; i++)
326 if (fgets(buf, 300, xyz) == NULL)
328 gmx_fatal(FARGS, "Unexpected end of ORCA output");
331 sscanf(buf, "%s%lf%lf%lf\n",
337 sscanf(buf, "%d%f%f%f\n",
343 for (j = 0; j < DIM; j++)
345 qm->xQM[i][j] *= 0.1;
350 sprintf(orca_engradFilename, "%s.engrad", qm->orca_basename);
351 engrad = fopen(orca_engradFilename, "r");
352 /* we read the energy and the gradient for the qm-atoms from the engrad file
354 /* we can skip the first seven lines
356 for (j = 0; j < 7; j++)
358 if (fgets(buf, 300, engrad) == NULL)
360 gmx_fatal(FARGS, "Unexpected end of ORCA output");
363 /* now comes the energy
365 if (fgets(buf, 300, engrad) == NULL)
367 gmx_fatal(FARGS, "Unexpected end of ORCA output");
370 sscanf(buf, "%lf\n", &QMener);
372 sscanf(buf, "%f\n", &QMener);
374 /* we can skip the next three lines
376 for (j = 0; j < 3; j++)
378 if (fgets(buf, 300, engrad) == NULL)
380 gmx_fatal(FARGS, "Unexpected end of ORCA output");
383 /* next lines contain the gradients of the QM atoms
384 * now comes the gradient, one value per line:
385 * (atom1 x \n atom1 y \n atom1 z \n atom2 x ...
388 for (i = 0; i < 3*qm->nrQMatoms; i++)
391 if (fgets(buf, 300, engrad) == NULL)
393 gmx_fatal(FARGS, "Unexpected end of ORCA output");
398 sscanf(buf, "%lf\n", &QMgrad[k][XX]);
402 sscanf(buf, "%lf\n", &QMgrad[k][YY]);
406 sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
411 sscanf(buf, "%f\n", &QMgrad[k][XX]);
415 sscanf(buf, "%f\n", &QMgrad[k][YY]);
419 sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
424 /* write the MM point charge data
426 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
428 sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
429 pcgrad = fopen(orca_pcgradFilename, "r");
431 /* we read the gradient for the mm-atoms from the pcgrad file
433 /* we can skip the first line
435 if (fgets(buf, 300, pcgrad) == NULL)
437 gmx_fatal(FARGS, "Unexpected end of ORCA output");
439 for (i = 0; i < mm->nrMMatoms; i++)
441 if (fgets(buf, 300, pcgrad) == NULL)
443 gmx_fatal(FARGS, "Unexpected end of ORCA output");
446 sscanf(buf, "%lf%lf%lf\n",
451 sscanf(buf, "%f%f%f\n",
462 void do_orca(char *orca_dir, char *basename)
465 /* make the call to the orca binary through system()
466 * The location of the binary is set through the
471 sprintf(buf, "%s/%s %s.inp >> %s.out",
476 fprintf(stderr, "Calling '%s'\n", buf);
477 if (system(buf) != 0)
479 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
483 real call_orca(t_forcerec *fr,
484 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
486 /* normal orca jobs */
499 sprintf(exe, "%s", "orca");
500 snew(QMgrad, qm->nrQMatoms);
501 snew(MMgrad, mm->nrMMatoms);
503 write_orca_input(fr, qm, mm);
504 do_orca(qm->orca_dir, qm->orca_basename);
505 QMener = read_orca_output(QMgrad, MMgrad, fr, qm, mm);
506 /* put the QMMM forces in the force array and to the fshift
508 for (i = 0; i < qm->nrQMatoms; i++)
510 for (j = 0; j < DIM; j++)
512 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
513 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
516 for (i = 0; i < mm->nrMMatoms; i++)
518 for (j = 0; j < DIM; j++)
520 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
521 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
524 QMener = QMener*HARTREE2KJ*AVOGADRO;
530 /* end of orca sub routines */