1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Groningen Machine for Chemical Simulation
62 #include "gmx_fatal.h"
66 /* ORCA interface routines */
68 void init_orca(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
73 /* ORCA settings on the system */
74 buf = getenv("BASENAME");
77 snew(qm->orca_basename, 200);
78 sscanf(buf, "%s", qm->orca_basename);
82 gmx_fatal(FARGS, "no $BASENAME\n");
85 /* ORCA directory on the system */
87 buf = getenv("ORCA_PATH");
88 fprintf(stderr, "%s", buf);
92 snew(qm->orca_dir, 200);
93 sscanf(buf, "%s", qm->orca_dir);
97 gmx_fatal(FARGS, "no $ORCA_PATH, check manual\n");
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);
109 void write_orca_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
116 *out, *pcFile, *addInputFile, *LJCoeff;
118 *buf, *orcaInput, *addInputFilename, *LJCoeffFilename,
119 *pcFilename, *exclInName, *exclOutName;
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");
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");
142 /* here we include the insertion of the additional orca-input */
144 if (addInputFile != NULL)
146 while (!feof(addInputFile))
148 if (fgets(buf, 200, addInputFile) != NULL)
156 fprintf(stderr, "No information on the calculation given in <%s>\n", addInputFilename);
157 gmx_call("qm_orca.c");
159 fclose(addInputFile);
160 if (qm->bTS || qm->bOPT)
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.
167 /* ORCA reads the exclusions from LJCoeffFilename.Excl,
168 * so we have to rename the file
171 for (i = 0; i < qm->nrQMatoms; i++)
173 if (qm->frontatoms[i])
177 fprintf(out, "%s\n", "%geom");
178 fprintf(out, " Constraints \n");
181 fprintf(out, " {C %d C}\n", i); /* counting from 0 */
186 fprintf(out, " end\n end\n");
188 /* make a file with information on the C6 and C12 coefficients */
189 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
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++)
205 fprintf(LJCoeff, "%10.7lf %10.7lf\n", qm->c6[i], qm->c12[i]);
207 fprintf(LJCoeff, "%10.7f %10.7f\n", qm->c6[i], qm->c12[i]);
210 fprintf(LJCoeff, "%d\n", mm->nrMMatoms);
211 for (i = 0; i < mm->nrMMatoms; i++)
214 fprintf(LJCoeff, "%10.7lf %10.7lf\n", mm->c6[i], mm->c12[i]);
216 fprintf(LJCoeff, "%10.7f %10.7f\n", mm->c6[i], mm->c12[i]);
222 /* write charge and multiplicity
224 fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
225 /* write the QM coordinates
227 for (i = 0; i < qm->nrQMatoms; i++)
230 if (qm->atomicnumberQM[i] == 0)
236 atomNr = qm->atomicnumberQM[i];
239 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
245 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
253 /* write the MM point charge data
255 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
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++)
266 fprintf(pcFile, "%8.4lf %10.7lf %10.7lf %10.7lf\n",
272 fprintf(pcFile, "%8.4f %10.7f %10.7f %10.7f\n",
279 fprintf(pcFile, "\n");
285 } /* write_orca_input */
287 real read_orca_output(rvec QMgrad[], rvec MMgrad[], int step, t_forcerec *fr,
288 t_QMrec *qm, t_MMrec *mm)
293 buf[300], tmp[300], orca_xyzFilename[300], orca_pcgradFilename[300], orca_engradFilename[300];
297 *xyz, *pcgrad, *engrad;
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.
307 /* we need the new xyz coordinates of the QM atoms only for separate QM-optimization
310 if (qm->bTS || qm->bOPT)
312 sprintf(orca_xyzFilename, "%s.xyz", qm->orca_basename);
313 xyz = fopen(orca_xyzFilename, "r");
314 if (fgets(buf, 300, xyz) == NULL)
316 gmx_fatal(FARGS, "Unexpected end of ORCA output");
318 if (fgets(buf, 300, xyz) == NULL)
320 gmx_fatal(FARGS, "Unexpected end of ORCA output");
322 for (i = 0; i < qm->nrQMatoms; i++)
324 if (fgets(buf, 300, xyz) == NULL)
326 gmx_fatal(FARGS, "Unexpected end of ORCA output");
329 sscanf(buf, "%s%lf%lf%lf\n",
335 sscanf(buf, "%d%f%f%f\n",
341 for (j = 0; j < DIM; j++)
343 qm->xQM[i][j] *= 0.1;
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
352 /* we can skip the first seven lines
354 for (j = 0; j < 7; j++)
356 if (fgets(buf, 300, engrad) == NULL)
358 gmx_fatal(FARGS, "Unexpected end of ORCA output");
361 /* now comes the energy
363 if (fgets(buf, 300, engrad) == NULL)
365 gmx_fatal(FARGS, "Unexpected end of ORCA output");
368 sscanf(buf, "%lf\n", &QMener);
370 sscanf(buf, "%f\n", &QMener);
372 /* we can skip the next three lines
374 for (j = 0; j < 3; j++)
376 if (fgets(buf, 300, engrad) == NULL)
378 gmx_fatal(FARGS, "Unexpected end of ORCA output");
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 ...
386 for (i = 0; i < 3*qm->nrQMatoms; i++)
389 if (fgets(buf, 300, engrad) == NULL)
391 gmx_fatal(FARGS, "Unexpected end of ORCA output");
396 sscanf(buf, "%lf\n", &QMgrad[k][XX]);
400 sscanf(buf, "%lf\n", &QMgrad[k][YY]);
404 sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
409 sscanf(buf, "%f\n", &QMgrad[k][XX]);
413 sscanf(buf, "%f\n", &QMgrad[k][YY]);
417 sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
422 /* write the MM point charge data
424 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
426 sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
427 pcgrad = fopen(orca_pcgradFilename, "r");
429 /* we read the gradient for the mm-atoms from the pcgrad file
431 /* we can skip the first line
433 if (fgets(buf, 300, pcgrad) == NULL)
435 gmx_fatal(FARGS, "Unexpected end of ORCA output");
437 for (i = 0; i < mm->nrMMatoms; i++)
439 if (fgets(buf, 300, pcgrad) == NULL)
441 gmx_fatal(FARGS, "Unexpected end of ORCA output");
444 sscanf(buf, "%lf%lf%lf\n",
449 sscanf(buf, "%f%f%f\n",
460 void do_orca(int step, char *exe, char *orca_dir, char *basename)
463 /* make the call to the orca binary through system()
464 * The location of the binary is set through the
469 sprintf(buf, "%s/%s %s.inp >> %s.out",
474 fprintf(stderr, "Calling '%s'\n", buf);
475 if (system(buf) != 0)
477 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
481 real call_orca(t_commrec *cr, t_forcerec *fr,
482 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
484 /* normal orca jobs */
497 sprintf(exe, "%s", "orca");
498 snew(QMgrad, qm->nrQMatoms);
499 snew(MMgrad, mm->nrMMatoms);
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
506 for (i = 0; i < qm->nrQMatoms; i++)
508 for (j = 0; j < DIM; j++)
510 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
511 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
514 for (i = 0; i < mm->nrMMatoms; i++)
516 for (j = 0; j < DIM; j++)
518 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
519 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
522 QMener = QMener*HARTREE2KJ*AVOGADRO;
528 /* end of orca sub routines */