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 * 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.
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.
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.
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.
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.
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.
66 #include "gmx_fatal.h"
70 /* ORCA interface routines */
72 void init_orca(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
77 /* ORCA settings on the system */
78 buf = getenv("BASENAME");
81 snew(qm->orca_basename, 200);
82 sscanf(buf, "%s", qm->orca_basename);
86 gmx_fatal(FARGS, "no $BASENAME\n");
89 /* ORCA directory on the system */
91 buf = getenv("ORCA_PATH");
92 fprintf(stderr, "%s", buf);
96 snew(qm->orca_dir, 200);
97 sscanf(buf, "%s", qm->orca_dir);
101 gmx_fatal(FARGS, "no $ORCA_PATH, check manual\n");
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);
113 void write_orca_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
120 *out, *pcFile, *addInputFile, *LJCoeff;
122 *buf, *orcaInput, *addInputFilename, *LJCoeffFilename,
123 *pcFilename, *exclInName, *exclOutName;
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");
135 fprintf(out, "!QMMMOpt TightSCF\n");
136 fprintf(out, "%s\n", "%geom TS_Search EF end");
140 fprintf(out, "!QMMMOpt TightSCF\n");
144 fprintf(out, "!EnGrad TightSCF\n");
146 /* here we include the insertion of the additional orca-input */
148 if (addInputFile != NULL)
150 while (!feof(addInputFile))
152 if (fgets(buf, 200, addInputFile) != NULL)
160 fprintf(stderr, "No information on the calculation given in <%s>\n", addInputFilename);
161 gmx_call("qm_orca.c");
163 fclose(addInputFile);
164 if (qm->bTS || qm->bOPT)
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.
171 /* ORCA reads the exclusions from LJCoeffFilename.Excl,
172 * so we have to rename the file
175 for (i = 0; i < qm->nrQMatoms; i++)
177 if (qm->frontatoms[i])
181 fprintf(out, "%s\n", "%geom");
182 fprintf(out, " Constraints \n");
185 fprintf(out, " {C %d C}\n", i); /* counting from 0 */
190 fprintf(out, " end\n end\n");
192 /* make a file with information on the C6 and C12 coefficients */
193 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
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++)
209 fprintf(LJCoeff, "%10.7lf %10.7lf\n", qm->c6[i], qm->c12[i]);
211 fprintf(LJCoeff, "%10.7f %10.7f\n", qm->c6[i], qm->c12[i]);
214 fprintf(LJCoeff, "%d\n", mm->nrMMatoms);
215 for (i = 0; i < mm->nrMMatoms; i++)
218 fprintf(LJCoeff, "%10.7lf %10.7lf\n", mm->c6[i], mm->c12[i]);
220 fprintf(LJCoeff, "%10.7f %10.7f\n", mm->c6[i], mm->c12[i]);
226 /* write charge and multiplicity
228 fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
229 /* write the QM coordinates
231 for (i = 0; i < qm->nrQMatoms; i++)
234 if (qm->atomicnumberQM[i] == 0)
240 atomNr = qm->atomicnumberQM[i];
243 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
249 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
257 /* write the MM point charge data
259 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
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++)
270 fprintf(pcFile, "%8.4lf %10.7lf %10.7lf %10.7lf\n",
276 fprintf(pcFile, "%8.4f %10.7f %10.7f %10.7f\n",
283 fprintf(pcFile, "\n");
289 } /* write_orca_input */
291 real read_orca_output(rvec QMgrad[], rvec MMgrad[], int step, t_forcerec *fr,
292 t_QMrec *qm, t_MMrec *mm)
297 buf[300], tmp[300], orca_xyzFilename[300], orca_pcgradFilename[300], orca_engradFilename[300];
301 *xyz, *pcgrad, *engrad;
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.
311 /* we need the new xyz coordinates of the QM atoms only for separate QM-optimization
314 if (qm->bTS || qm->bOPT)
316 sprintf(orca_xyzFilename, "%s.xyz", qm->orca_basename);
317 xyz = fopen(orca_xyzFilename, "r");
318 if (fgets(buf, 300, xyz) == NULL)
320 gmx_fatal(FARGS, "Unexpected end of ORCA output");
322 if (fgets(buf, 300, xyz) == NULL)
324 gmx_fatal(FARGS, "Unexpected end of ORCA output");
326 for (i = 0; i < qm->nrQMatoms; i++)
328 if (fgets(buf, 300, xyz) == NULL)
330 gmx_fatal(FARGS, "Unexpected end of ORCA output");
333 sscanf(buf, "%s%lf%lf%lf\n",
339 sscanf(buf, "%d%f%f%f\n",
345 for (j = 0; j < DIM; j++)
347 qm->xQM[i][j] *= 0.1;
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
356 /* we can skip the first seven lines
358 for (j = 0; j < 7; j++)
360 if (fgets(buf, 300, engrad) == NULL)
362 gmx_fatal(FARGS, "Unexpected end of ORCA output");
365 /* now comes the energy
367 if (fgets(buf, 300, engrad) == NULL)
369 gmx_fatal(FARGS, "Unexpected end of ORCA output");
372 sscanf(buf, "%lf\n", &QMener);
374 sscanf(buf, "%f\n", &QMener);
376 /* we can skip the next three lines
378 for (j = 0; j < 3; j++)
380 if (fgets(buf, 300, engrad) == NULL)
382 gmx_fatal(FARGS, "Unexpected end of ORCA output");
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 ...
390 for (i = 0; i < 3*qm->nrQMatoms; i++)
393 if (fgets(buf, 300, engrad) == NULL)
395 gmx_fatal(FARGS, "Unexpected end of ORCA output");
400 sscanf(buf, "%lf\n", &QMgrad[k][XX]);
404 sscanf(buf, "%lf\n", &QMgrad[k][YY]);
408 sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
413 sscanf(buf, "%f\n", &QMgrad[k][XX]);
417 sscanf(buf, "%f\n", &QMgrad[k][YY]);
421 sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
426 /* write the MM point charge data
428 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
430 sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
431 pcgrad = fopen(orca_pcgradFilename, "r");
433 /* we read the gradient for the mm-atoms from the pcgrad file
435 /* we can skip the first line
437 if (fgets(buf, 300, pcgrad) == NULL)
439 gmx_fatal(FARGS, "Unexpected end of ORCA output");
441 for (i = 0; i < mm->nrMMatoms; i++)
443 if (fgets(buf, 300, pcgrad) == NULL)
445 gmx_fatal(FARGS, "Unexpected end of ORCA output");
448 sscanf(buf, "%lf%lf%lf\n",
453 sscanf(buf, "%f%f%f\n",
464 void do_orca(int step, char *exe, char *orca_dir, char *basename)
467 /* make the call to the orca binary through system()
468 * The location of the binary is set through the
473 sprintf(buf, "%s/%s %s.inp >> %s.out",
478 fprintf(stderr, "Calling '%s'\n", buf);
479 if (system(buf) != 0)
481 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
485 real call_orca(t_commrec *cr, t_forcerec *fr,
486 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
488 /* normal orca jobs */
501 sprintf(exe, "%s", "orca");
502 snew(QMgrad, qm->nrQMatoms);
503 snew(MMgrad, mm->nrMMatoms);
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
510 for (i = 0; i < qm->nrQMatoms; i++)
512 for (j = 0; j < DIM; j++)
514 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
515 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
518 for (i = 0; i < mm->nrMMatoms; i++)
520 for (j = 0; j < DIM; j++)
522 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
523 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
526 QMener = QMener*HARTREE2KJ*AVOGADRO;
532 /* end of orca sub routines */