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
61 #include "gmx_fatal.h"
65 /* ORCA interface routines */
67 void init_orca(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
72 /* ORCA settings on the system */
73 buf = getenv("BASENAME");
76 snew(qm->orca_basename, 200);
77 sscanf(buf, "%s", qm->orca_basename);
81 gmx_fatal(FARGS, "no $BASENAME\n");
84 /* ORCA directory on the system */
86 buf = getenv("ORCA_PATH");
87 fprintf(stderr, "%s", buf);
91 snew(qm->orca_dir, 200);
92 sscanf(buf, "%s", qm->orca_dir);
96 gmx_fatal(FARGS, "no $ORCA_PATH, check manual\n");
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);
108 void write_orca_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
115 *out, *pcFile, *addInputFile, *LJCoeff;
117 *buf, *orcaInput, *addInputFilename, *LJCoeffFilename,
118 *pcFilename, *exclInName, *exclOutName;
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");
130 fprintf(out, "!QMMMOpt TightSCF\n");
131 fprintf(out, "%s\n", "%geom TS_Search EF end");
135 fprintf(out, "!QMMMOpt TightSCF\n");
139 fprintf(out, "!EnGrad TightSCF\n");
141 /* here we include the insertion of the additional orca-input */
143 if (addInputFile != NULL)
145 while (!feof(addInputFile))
147 if (fgets(buf, 200, addInputFile) != NULL)
155 fprintf(stderr, "No information on the calculation given in <%s>\n", addInputFilename);
156 gmx_call("qm_orca.c");
158 fclose(addInputFile);
159 if (qm->bTS || qm->bOPT)
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.
166 /* ORCA reads the exclusions from LJCoeffFilename.Excl,
167 * so we have to rename the file
170 for (i = 0; i < qm->nrQMatoms; i++)
172 if (qm->frontatoms[i])
176 fprintf(out, "%s\n", "%geom");
177 fprintf(out, " Constraints \n");
180 fprintf(out, " {C %d C}\n", i); /* counting from 0 */
185 fprintf(out, " end\n end\n");
187 /* make a file with information on the C6 and C12 coefficients */
188 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
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++)
204 fprintf(LJCoeff, "%10.7lf %10.7lf\n", qm->c6[i], qm->c12[i]);
206 fprintf(LJCoeff, "%10.7f %10.7f\n", qm->c6[i], qm->c12[i]);
209 fprintf(LJCoeff, "%d\n", mm->nrMMatoms);
210 for (i = 0; i < mm->nrMMatoms; i++)
213 fprintf(LJCoeff, "%10.7lf %10.7lf\n", mm->c6[i], mm->c12[i]);
215 fprintf(LJCoeff, "%10.7f %10.7f\n", mm->c6[i], mm->c12[i]);
221 /* write charge and multiplicity
223 fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
224 /* write the QM coordinates
226 for (i = 0; i < qm->nrQMatoms; i++)
229 if (qm->atomicnumberQM[i] == 0)
235 atomNr = qm->atomicnumberQM[i];
238 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
244 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
252 /* write the MM point charge data
254 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
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++)
265 fprintf(pcFile, "%8.4lf %10.7lf %10.7lf %10.7lf\n",
271 fprintf(pcFile, "%8.4f %10.7f %10.7f %10.7f\n",
278 fprintf(pcFile, "\n");
284 } /* write_orca_input */
286 real read_orca_output(rvec QMgrad[], rvec MMgrad[], int step, t_forcerec *fr,
287 t_QMrec *qm, t_MMrec *mm)
292 buf[300], tmp[300], orca_xyzFilename[300], orca_pcgradFilename[300], orca_engradFilename[300];
296 *xyz, *pcgrad, *engrad;
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.
306 /* we need the new xyz coordinates of the QM atoms only for separate QM-optimization
309 if (qm->bTS || qm->bOPT)
311 sprintf(orca_xyzFilename, "%s.xyz", qm->orca_basename);
312 xyz = fopen(orca_xyzFilename, "r");
313 if (fgets(buf, 300, xyz) == NULL)
315 gmx_fatal(FARGS, "Unexpected end of ORCA output");
317 if (fgets(buf, 300, xyz) == NULL)
319 gmx_fatal(FARGS, "Unexpected end of ORCA output");
321 for (i = 0; i < qm->nrQMatoms; i++)
323 if (fgets(buf, 300, xyz) == NULL)
325 gmx_fatal(FARGS, "Unexpected end of ORCA output");
328 sscanf(buf, "%s%lf%lf%lf\n",
334 sscanf(buf, "%d%f%f%f\n",
340 for (j = 0; j < DIM; j++)
342 qm->xQM[i][j] *= 0.1;
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
351 /* we can skip the first seven lines
353 for (j = 0; j < 7; j++)
355 if (fgets(buf, 300, engrad) == NULL)
357 gmx_fatal(FARGS, "Unexpected end of ORCA output");
360 /* now comes the energy
362 if (fgets(buf, 300, engrad) == NULL)
364 gmx_fatal(FARGS, "Unexpected end of ORCA output");
367 sscanf(buf, "%lf\n", &QMener);
369 sscanf(buf, "%f\n", &QMener);
371 /* we can skip the next three lines
373 for (j = 0; j < 3; j++)
375 if (fgets(buf, 300, engrad) == NULL)
377 gmx_fatal(FARGS, "Unexpected end of ORCA output");
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 ...
385 for (i = 0; i < 3*qm->nrQMatoms; i++)
388 if (fgets(buf, 300, engrad) == NULL)
390 gmx_fatal(FARGS, "Unexpected end of ORCA output");
395 sscanf(buf, "%lf\n", &QMgrad[k][XX]);
399 sscanf(buf, "%lf\n", &QMgrad[k][YY]);
403 sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
408 sscanf(buf, "%f\n", &QMgrad[k][XX]);
412 sscanf(buf, "%f\n", &QMgrad[k][YY]);
416 sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
421 /* write the MM point charge data
423 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
425 sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
426 pcgrad = fopen(orca_pcgradFilename, "r");
428 /* we read the gradient for the mm-atoms from the pcgrad file
430 /* we can skip the first line
432 if (fgets(buf, 300, pcgrad) == NULL)
434 gmx_fatal(FARGS, "Unexpected end of ORCA output");
436 for (i = 0; i < mm->nrMMatoms; i++)
438 if (fgets(buf, 300, pcgrad) == NULL)
440 gmx_fatal(FARGS, "Unexpected end of ORCA output");
443 sscanf(buf, "%lf%lf%lf\n",
448 sscanf(buf, "%f%f%f\n",
459 void do_orca(int step, char *exe, char *orca_dir, char *basename)
462 /* make the call to the orca binary through system()
463 * The location of the binary is set through the
468 sprintf(buf, "%s/%s %s.inp >> %s.out",
473 fprintf(stderr, "Calling '%s'\n", buf);
474 if (system(buf) != 0)
476 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
480 real call_orca(t_commrec *cr, t_forcerec *fr,
481 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
483 /* normal orca jobs */
496 sprintf(exe, "%s", "orca");
497 snew(QMgrad, qm->nrQMatoms);
498 snew(MMgrad, mm->nrMMatoms);
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
505 for (i = 0; i < qm->nrQMatoms; i++)
507 for (j = 0; j < DIM; j++)
509 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
510 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
513 for (i = 0; i < mm->nrMMatoms; i++)
515 for (j = 0; j < DIM; j++)
517 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
518 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
521 QMener = QMener*HARTREE2KJ*AVOGADRO;
527 /* end of orca sub routines */