Moved mdatom.h from legacyheader/types to mdtypes.
[gromacs.git] / src / gromacs / mdlib / qm_orca.cpp
blob9d2ff860cb7e93da577fbc73e809c8c7adaa6cb2
1 /*
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,2014,2015, 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.
37 #include "gmxpre.h"
39 #include <math.h>
40 #include <stdio.h>
41 #include <stdlib.h>
42 #include <string.h>
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/txtdump.h"
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/legacyheaders/names.h"
48 #include "gromacs/legacyheaders/nrnb.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdlib/force.h"
52 #include "gromacs/mdlib/ns.h"
53 #include "gromacs/mdlib/qmmm.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
57 /* ORCA interface routines */
59 void init_orca(t_QMrec *qm)
61 char *buf;
62 snew(buf, 200);
64 /* ORCA settings on the system */
65 buf = getenv("GMX_QM_ORCA_BASENAME");
66 if (buf)
68 snew(qm->orca_basename, 200);
69 sscanf(buf, "%s", qm->orca_basename);
71 else
73 gmx_fatal(FARGS, "$GMX_QM_ORCA_BASENAME is not set\n");
76 /* ORCA directory on the system */
77 snew(buf, 200);
78 buf = getenv("GMX_ORCA_PATH");
80 if (buf)
82 snew(qm->orca_dir, 200);
83 sscanf(buf, "%s", qm->orca_dir);
85 else
87 gmx_fatal(FARGS, "$GMX_ORCA_PATH not set, check manual\n");
90 fprintf(stderr, "Setting ORCA path to: %s...\n", qm->orca_dir);
91 fprintf(stderr, "ORCA initialised...\n\n");
92 /* since we append the output to the BASENAME.out file,
93 we should delete an existent old out-file here. */
94 sprintf(buf, "%s.out", qm->orca_basename);
95 remove(buf);
99 void write_orca_input(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
101 int i;
102 t_QMMMrec *QMMMrec;
103 FILE *out, *pcFile, *addInputFile, *LJCoeff;
104 char *buf, *orcaInput, *addInputFilename, *LJCoeffFilename, *pcFilename, *exclInName, *exclOutName;
106 QMMMrec = fr->qr;
108 /* write the first part of the input-file */
109 snew(orcaInput, 200);
110 sprintf(orcaInput, "%s.inp", qm->orca_basename);
111 out = fopen(orcaInput, "w");
113 snew(addInputFilename, 200);
114 sprintf(addInputFilename, "%s.ORCAINFO", qm->orca_basename);
115 addInputFile = fopen(addInputFilename, "r");
117 fprintf(out, "#input-file generated by GROMACS\n");
119 if (qm->bTS)
121 fprintf(out, "!QMMMOpt TightSCF\n");
122 fprintf(out, "%s\n", "%geom TS_Search EF end");
124 else if (qm->bOPT)
126 fprintf(out, "!QMMMOpt TightSCF\n");
128 else
130 fprintf(out, "!EnGrad TightSCF\n");
133 /* here we include the insertion of the additional orca-input */
134 snew(buf, 200);
135 if (addInputFile != NULL)
137 while (!feof(addInputFile))
139 if (fgets(buf, 200, addInputFile) != NULL)
141 fputs(buf, out);
145 else
147 gmx_fatal(FARGS, "No information on the calculation given in %s\n", addInputFilename);
150 fclose(addInputFile);
152 if (qm->bTS || qm->bOPT)
154 /* freeze the frontier QM atoms and Link atoms. This is
155 * important only if a full QM subsystem optimization is done
156 * with a frozen MM environmeent. For dynamics, or gromacs's own
157 * optimization routines this is not important.
159 /* ORCA reads the exclusions from LJCoeffFilename.Excl,
160 * so we have to rename the file
162 int didStart = 0;
163 for (i = 0; i < qm->nrQMatoms; i++)
165 if (qm->frontatoms[i])
167 if (!didStart)
169 fprintf(out, "%s\n", "%geom");
170 fprintf(out, " Constraints \n");
171 didStart = 1;
173 fprintf(out, " {C %d C}\n", i); /* counting from 0 */
176 if (didStart)
178 fprintf(out, " end\n end\n");
180 /* make a file with information on the C6 and C12 coefficients */
181 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
183 snew(exclInName, 200);
184 snew(exclOutName, 200);
185 sprintf(exclInName, "QMMMexcl.dat");
186 sprintf(exclOutName, "%s.LJ.Excl", qm->orca_basename);
187 rename(exclInName, exclOutName);
188 snew(LJCoeffFilename, 200);
189 sprintf(LJCoeffFilename, "%s.LJ", qm->orca_basename);
190 fprintf(out, "%s%s%s\n", "%LJCOEFFICIENTS \"", LJCoeffFilename, "\"");
191 /* make a file with information on the C6 and C12 coefficients */
192 LJCoeff = fopen(LJCoeffFilename, "w");
193 fprintf(LJCoeff, "%d\n", qm->nrQMatoms);
194 for (i = 0; i < qm->nrQMatoms; i++)
196 #ifdef GMX_DOUBLE
197 fprintf(LJCoeff, "%10.7lf %10.7lf\n", qm->c6[i], qm->c12[i]);
198 #else
199 fprintf(LJCoeff, "%10.7f %10.7f\n", qm->c6[i], qm->c12[i]);
200 #endif
202 fprintf(LJCoeff, "%d\n", mm->nrMMatoms);
203 for (i = 0; i < mm->nrMMatoms; i++)
205 #ifdef GMX_DOUBLE
206 fprintf(LJCoeff, "%10.7lf %10.7lf\n", mm->c6[i], mm->c12[i]);
207 #else
208 fprintf(LJCoeff, "%10.7f %10.7f\n", mm->c6[i], mm->c12[i]);
209 #endif
211 fclose(LJCoeff);
215 /* write charge and multiplicity */
216 fprintf(out, "*xyz %2d%2d\n", qm->QMcharge, qm->multiplicity);
218 /* write the QM coordinates */
219 for (i = 0; i < qm->nrQMatoms; i++)
221 int atomNr;
222 if (qm->atomicnumberQM[i] == 0)
224 atomNr = 1;
226 else
228 atomNr = qm->atomicnumberQM[i];
230 #ifdef GMX_DOUBLE
231 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
232 atomNr,
233 qm->xQM[i][XX]/0.1,
234 qm->xQM[i][YY]/0.1,
235 qm->xQM[i][ZZ]/0.1);
236 #else
237 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
238 atomNr,
239 qm->xQM[i][XX]/0.1,
240 qm->xQM[i][YY]/0.1,
241 qm->xQM[i][ZZ]/0.1);
242 #endif
244 fprintf(out, "*\n");
246 /* write the MM point charge data */
247 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
249 /* name of the point charge file */
250 snew(pcFilename, 200);
251 sprintf(pcFilename, "%s.pc", qm->orca_basename);
252 fprintf(out, "%s%s%s\n", "%pointcharges \"", pcFilename, "\"");
253 pcFile = fopen(pcFilename, "w");
254 fprintf(pcFile, "%d\n", mm->nrMMatoms);
255 for (i = 0; i < mm->nrMMatoms; i++)
257 #ifdef GMX_DOUBLE
258 fprintf(pcFile, "%8.4lf %10.7lf %10.7lf %10.7lf\n",
259 mm->MMcharges[i],
260 mm->xMM[i][XX]/0.1,
261 mm->xMM[i][YY]/0.1,
262 mm->xMM[i][ZZ]/0.1);
263 #else
264 fprintf(pcFile, "%8.4f %10.7f %10.7f %10.7f\n",
265 mm->MMcharges[i],
266 mm->xMM[i][XX]/0.1,
267 mm->xMM[i][YY]/0.1,
268 mm->xMM[i][ZZ]/0.1);
269 #endif
271 fprintf(pcFile, "\n");
272 fclose(pcFile);
274 fprintf(out, "\n");
276 fclose(out);
277 } /* write_orca_input */
279 real read_orca_output(rvec QMgrad[], rvec MMgrad[], t_forcerec *fr,
280 t_QMrec *qm, t_MMrec *mm)
283 i, j, atnum;
284 char
285 buf[300], orca_xyzFilename[300], orca_pcgradFilename[300], orca_engradFilename[300];
286 real
287 QMener;
288 FILE
289 *xyz, *pcgrad, *engrad;
290 int k;
291 t_QMMMrec
292 *QMMMrec;
293 QMMMrec = fr->qr;
294 /* in case of an optimization, the coordinates are printed in the
295 * xyz file, the energy and gradients for the QM part are stored in the engrad file
296 * and the gradients for the point charges are stored in the pc file.
299 /* we need the new xyz coordinates of the QM atoms only for separate QM-optimization
302 if (qm->bTS || qm->bOPT)
304 sprintf(orca_xyzFilename, "%s.xyz", qm->orca_basename);
305 xyz = fopen(orca_xyzFilename, "r");
306 if (fgets(buf, 300, xyz) == NULL)
308 gmx_fatal(FARGS, "Unexpected end of ORCA output");
310 if (fgets(buf, 300, xyz) == NULL)
312 gmx_fatal(FARGS, "Unexpected end of ORCA output");
314 for (i = 0; i < qm->nrQMatoms; i++)
316 if (fgets(buf, 300, xyz) == NULL)
318 gmx_fatal(FARGS, "Unexpected end of ORCA output");
320 #ifdef GMX_DOUBLE
321 sscanf(buf, "%d%lf%lf%lf\n",
322 &atnum,
323 &qm->xQM[i][XX],
324 &qm->xQM[i][YY],
325 &qm->xQM[i][ZZ]);
326 #else
327 sscanf(buf, "%d%f%f%f\n",
328 &atnum,
329 &qm->xQM[i][XX],
330 &qm->xQM[i][YY],
331 &qm->xQM[i][ZZ]);
332 #endif
333 for (j = 0; j < DIM; j++)
335 qm->xQM[i][j] *= 0.1;
338 fclose(xyz);
340 sprintf(orca_engradFilename, "%s.engrad", qm->orca_basename);
341 engrad = fopen(orca_engradFilename, "r");
342 /* we read the energy and the gradient for the qm-atoms from the engrad file
344 /* we can skip the first seven lines
346 for (j = 0; j < 7; j++)
348 if (fgets(buf, 300, engrad) == NULL)
350 gmx_fatal(FARGS, "Unexpected end of ORCA output");
353 /* now comes the energy
355 if (fgets(buf, 300, engrad) == NULL)
357 gmx_fatal(FARGS, "Unexpected end of ORCA output");
359 #ifdef GMX_DOUBLE
360 sscanf(buf, "%lf\n", &QMener);
361 #else
362 sscanf(buf, "%f\n", &QMener);
363 #endif
364 /* we can skip the next three lines
366 for (j = 0; j < 3; j++)
368 if (fgets(buf, 300, engrad) == NULL)
370 gmx_fatal(FARGS, "Unexpected end of ORCA output");
373 /* next lines contain the gradients of the QM atoms
374 * now comes the gradient, one value per line:
375 * (atom1 x \n atom1 y \n atom1 z \n atom2 x ...
378 for (i = 0; i < 3*qm->nrQMatoms; i++)
380 k = i/3;
381 if (fgets(buf, 300, engrad) == NULL)
383 gmx_fatal(FARGS, "Unexpected end of ORCA output");
385 #ifdef GMX_DOUBLE
386 if (i%3 == 0)
388 sscanf(buf, "%lf\n", &QMgrad[k][XX]);
390 else if (i%3 == 1)
392 sscanf(buf, "%lf\n", &QMgrad[k][YY]);
394 else if (i%3 == 2)
396 sscanf(buf, "%lf\n", &QMgrad[k][ZZ]);
398 #else
399 if (i%3 == 0)
401 sscanf(buf, "%f\n", &QMgrad[k][XX]);
403 else if (i%3 == 1)
405 sscanf(buf, "%f\n", &QMgrad[k][YY]);
407 else if (i%3 == 2)
409 sscanf(buf, "%f\n", &QMgrad[k][ZZ]);
411 #endif
413 fclose(engrad);
414 /* write the MM point charge data
416 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
418 sprintf(orca_pcgradFilename, "%s.pcgrad", qm->orca_basename);
419 pcgrad = fopen(orca_pcgradFilename, "r");
421 /* we read the gradient for the mm-atoms from the pcgrad file
423 /* we can skip the first line
425 if (fgets(buf, 300, pcgrad) == NULL)
427 gmx_fatal(FARGS, "Unexpected end of ORCA output");
429 for (i = 0; i < mm->nrMMatoms; i++)
431 if (fgets(buf, 300, pcgrad) == NULL)
433 gmx_fatal(FARGS, "Unexpected end of ORCA output");
435 #ifdef GMX_DOUBLE
436 sscanf(buf, "%lf%lf%lf\n",
437 &MMgrad[i][XX],
438 &MMgrad[i][YY],
439 &MMgrad[i][ZZ]);
440 #else
441 sscanf(buf, "%f%f%f\n",
442 &MMgrad[i][XX],
443 &MMgrad[i][YY],
444 &MMgrad[i][ZZ]);
445 #endif
447 fclose(pcgrad);
449 return(QMener);
452 void do_orca(char *orca_dir, char *basename)
455 /* make the call to the orca binary through system()
456 * The location of the binary is set through the
457 * environment.
459 char
460 buf[100];
461 sprintf(buf, "%s/%s %s.inp >> %s.out",
462 orca_dir,
463 "orca",
464 basename,
465 basename);
466 fprintf(stderr, "Calling '%s'\n", buf);
467 if (system(buf) != 0)
469 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
473 real call_orca(t_forcerec *fr,
474 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
476 /* normal orca jobs */
477 static int
478 step = 0;
480 i, j;
481 real
482 QMener;
483 rvec
484 *QMgrad, *MMgrad;
485 char
486 *exe;
488 snew(exe, 30);
489 sprintf(exe, "%s", "orca");
490 snew(QMgrad, qm->nrQMatoms);
491 snew(MMgrad, mm->nrMMatoms);
493 write_orca_input(fr, qm, mm);
494 do_orca(qm->orca_dir, qm->orca_basename);
495 QMener = read_orca_output(QMgrad, MMgrad, fr, qm, mm);
496 /* put the QMMM forces in the force array and to the fshift
498 for (i = 0; i < qm->nrQMatoms; i++)
500 for (j = 0; j < DIM; j++)
502 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
503 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
506 for (i = 0; i < mm->nrMMatoms; i++)
508 for (j = 0; j < DIM; j++)
510 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
511 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
514 QMener = QMener*HARTREE2KJ*AVOGADRO;
515 step++;
516 free(exe);
517 return(QMener);
518 } /* call_orca */
520 /* end of orca sub routines */