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-2004, 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 * GROwing Monsters And Cloning Shrimps
39 #ifdef GMX_QMMM_GAUSSIAN
64 #include "gmx_fatal.h"
69 /* TODO: this should be made thread-safe */
71 /* Gaussian interface routines */
73 void init_gaussian(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
76 *rffile = NULL, *out = NULL;
78 basissets[eQMbasisNR] = {{0, 3, 0},
79 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
93 /* using the ivec above to convert the basis read form the mdp file
94 * in a human readable format into some numbers for the gaussian
95 * route. This is necessary as we are using non standard routes to
99 /* per layer we make a new subdir for integral file, checkpoint
100 * files and such. These dirs are stored in the QMrec for
105 if (!qm->nQMcpus) /* this we do only once per layer
106 * as we call g01 externally
110 for (i = 0; i < DIM; i++)
112 qm->SHbasis[i] = basissets[qm->QMbasis][i];
115 /* init gradually switching on of the SA */
117 /* we read the number of cpus and environment from the environment
120 buf = getenv("NCPUS");
123 sscanf(buf, "%d", &qm->nQMcpus);
129 fprintf(stderr, "number of CPUs for gaussian = %d\n", qm->nQMcpus);
133 sscanf(buf, "%d", &qm->QMmem);
137 qm->QMmem = 50000000;
139 fprintf(stderr, "memory for gaussian = %d\n", qm->QMmem);
143 sscanf(buf, "%d", &qm->accuracy);
149 fprintf(stderr, "accuracy in l510 = %d\n", qm->accuracy);
151 buf = getenv("CPMCSCF");
154 sscanf(buf, "%d", &i);
155 qm->cpmcscf = (i != 0);
163 fprintf(stderr, "using cp-mcscf in l1003\n");
167 fprintf(stderr, "NOT using cp-mcscf in l1003\n");
169 buf = getenv("SASTEP");
172 sscanf(buf, "%d", &qm->SAstep);
176 /* init gradually switching on of the SA */
179 /* we read the number of cpus and environment from the environment
182 fprintf(stderr, "Level of SA at start = %d\n", qm->SAstep);
183 /* punch the LJ C6 and C12 coefficients to be picked up by
184 * gaussian and usd to compute the LJ interaction between the
187 if (qm->bTS || qm->bOPT)
189 out = fopen("LJ.dat", "w");
190 for (i = 0; i < qm->nrQMatoms; i++)
194 fprintf(out, "%3d %10.7lf %10.7lf\n",
195 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
197 fprintf(out, "%3d %10.7f %10.7f\n",
198 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
203 /* gaussian settings on the system */
204 buf = getenv("GAUSS_DIR");
205 fprintf(stderr, "%s", buf);
209 qm->gauss_dir = strdup(buf);
213 gmx_fatal(FARGS, "no $GAUSS_DIR, check gaussian manual\n");
216 buf = getenv("GAUSS_EXE");
219 qm->gauss_exe = strdup(buf);
223 gmx_fatal(FARGS, "no $GAUSS_EXE, check gaussian manual\n");
225 buf = getenv("DEVEL_DIR");
228 qm->devel_dir = strdup (buf);
232 gmx_fatal(FARGS, "no $DEVEL_DIR, this is were the modified links reside.\n");
236 /* reactionfield, file is needed using gaussian */
237 /* rffile=fopen("rf.dat","w");*/
238 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
242 fprintf(stderr, "gaussian initialised...\n");
247 void write_gaussian_SH_input(int step, gmx_bool swap,
248 t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
259 bSA = (qm->SAstep > 0);
261 out = fopen("input.com", "w");
262 /* write the route */
263 fprintf(out, "%s", "%scr=input\n");
264 fprintf(out, "%s", "%rwf=input\n");
265 fprintf(out, "%s", "%int=input\n");
266 fprintf(out, "%s", "%d2e=input\n");
268 * fprintf(out,"%s","%nosave\n");
270 fprintf(out, "%s", "%chk=input\n");
271 fprintf(out, "%s%d\n", "%mem=", qm->QMmem);
272 fprintf(out, "%s%3d\n", "%nprocshare=", qm->nQMcpus);
274 /* use the versions of
275 * l301 that computes the interaction between MM and QM atoms.
276 * l510 that can punch the CI coefficients
277 * l701 that can do gradients on MM atoms
281 fprintf(out, "%s%s%s",
285 fprintf(out, "%s%s%s",
289 fprintf(out, "%s%s%s",
294 fprintf(out, "%s%s%s",
298 fprintf(out, "%s%s%s",
302 /* print the nonstandard route
305 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
307 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
311 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
314 qm->SHbasis[2]); /*basisset stuff */
319 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
322 qm->SHbasis[2]); /*basisset stuff */
325 if (step+1) /* fetch initial guess from check point file */
326 { /* hack, to alyays read from chk file!!!!! */
327 fprintf(out, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
329 "18=", qm->CASorbitals, "/1,5;\n");
331 else /* generate the first checkpoint file */
333 fprintf(out, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
335 "18=", qm->CASorbitals, "/1,5;\n");
337 /* the rest of the input depends on where the system is on the PES
339 if (swap && bSA) /* make a slide to the other surface */
341 if (qm->CASorbitals > 6) /* use direct and no full diag */
343 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
349 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
351 if (mm->nrMMatoms > 0)
353 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
355 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
356 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
357 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
361 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
363 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
367 else if (bSA) /* do a "state-averaged" CAS calculation */
369 if (qm->CASorbitals > 6) /* no full diag */
371 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
377 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
379 if (mm->nrMMatoms > 0)
381 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
383 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
384 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1/3;\n");
385 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
389 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
391 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
395 else if (swap) /* do a "swapped" CAS calculation */
397 if (qm->CASorbitals > 6)
399 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
403 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
406 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
408 else /* do a "normal" CAS calculation */
410 if (qm->CASorbitals > 6)
412 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
416 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
419 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
421 fprintf(out, "\ninput-file generated by gromacs\n\n");
422 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
423 for (i = 0; i < qm->nrQMatoms; i++)
426 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
427 qm->atomicnumberQM[i],
428 qm->xQM[i][XX]/BOHR2NM,
429 qm->xQM[i][YY]/BOHR2NM,
430 qm->xQM[i][ZZ]/BOHR2NM);
432 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
433 qm->atomicnumberQM[i],
434 qm->xQM[i][XX]/BOHR2NM,
435 qm->xQM[i][YY]/BOHR2NM,
436 qm->xQM[i][ZZ]/BOHR2NM);
439 /* MM point charge data */
440 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
443 for (i = 0; i < mm->nrMMatoms; i++)
446 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
447 mm->xMM[i][XX]/BOHR2NM,
448 mm->xMM[i][YY]/BOHR2NM,
449 mm->xMM[i][ZZ]/BOHR2NM,
452 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
453 mm->xMM[i][XX]/BOHR2NM,
454 mm->xMM[i][YY]/BOHR2NM,
455 mm->xMM[i][ZZ]/BOHR2NM,
460 if (bSA) /* put the SA coefficients at the end of the file */
463 fprintf(out, "\n%10.8lf %10.8lf\n",
464 qm->SAstep*0.5/qm->SAsteps,
465 1-qm->SAstep*0.5/qm->SAsteps);
467 fprintf(out, "\n%10.8f %10.8f\n",
468 qm->SAstep*0.5/qm->SAsteps,
469 1-qm->SAstep*0.5/qm->SAsteps);
471 fprintf(stderr, "State Averaging level = %d/%d\n", qm->SAstep, qm->SAsteps);
475 } /* write_gaussian_SH_input */
477 void write_gaussian_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
487 out = fopen("input.com", "w");
488 /* write the route */
490 if (qm->QMmethod >= eQMmethodRHF)
502 fprintf(out, "%s%3d\n",
503 "%nprocshare=", qm->nQMcpus);
505 fprintf(out, "%s%d\n",
507 /* use the modified links that include the LJ contribution at the QM level */
508 if (qm->bTS || qm->bOPT)
510 fprintf(out, "%s%s%s",
511 "%subst l701 ", qm->devel_dir, "/l701_LJ\n");
512 fprintf(out, "%s%s%s",
513 "%subst l301 ", qm->devel_dir, "/l301_LJ\n");
517 fprintf(out, "%s%s%s",
518 "%subst l701 ", qm->devel_dir, "/l701\n");
519 fprintf(out, "%s%s%s",
520 "%subst l301 ", qm->devel_dir, "/l301\n");
522 fprintf(out, "%s%s%s",
523 "%subst l9999 ", qm->devel_dir, "/l9999\n");
534 if (qm->QMmethod == eQMmethodB3LYPLAN)
537 "B3LYP/GEN Pseudo=Read");
542 eQMmethod_names[qm->QMmethod]);
544 if (qm->QMmethod >= eQMmethodRHF)
546 if (qm->QMmethod == eQMmethodCASSCF)
548 /* in case of cas, how many electrons and orbitals do we need?
550 fprintf(out, "(%d,%d)",
551 qm->CASelectrons, qm->CASorbitals);
554 eQMbasis_names[qm->QMbasis]);
557 if (QMMMrec->QMMMscheme == eQMMMschemenormal && mm->nrMMatoms)
562 if (step || qm->QMmethod == eQMmethodCASSCF)
564 /* fetch guess from checkpoint file, always for CASSCF */
565 fprintf(out, "%s", " guess=read");
567 fprintf(out, "\nNosymm units=bohr\n");
571 fprintf(out, "OPT=(Redundant,TS,noeigentest,ModRedundant) Punch=(Coord,Derivatives) ");
575 fprintf(out, "OPT=(Redundant,ModRedundant) Punch=(Coord,Derivatives) ");
579 fprintf(out, "FORCE Punch=(Derivatives) ");
581 fprintf(out, "iop(3/33=1)\n\n");
582 fprintf(out, "input-file generated by gromacs\n\n");
583 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
584 for (i = 0; i < qm->nrQMatoms; i++)
587 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
588 qm->atomicnumberQM[i],
589 qm->xQM[i][XX]/BOHR2NM,
590 qm->xQM[i][YY]/BOHR2NM,
591 qm->xQM[i][ZZ]/BOHR2NM);
593 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
594 qm->atomicnumberQM[i],
595 qm->xQM[i][XX]/BOHR2NM,
596 qm->xQM[i][YY]/BOHR2NM,
597 qm->xQM[i][ZZ]/BOHR2NM);
601 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
602 if (qm->QMmethod == eQMmethodB3LYPLAN)
605 for (i = 0; i < qm->nrQMatoms; i++)
607 if (qm->atomicnumberQM[i] < 21)
609 fprintf(out, "%d ", i+1);
612 fprintf(out, "\n%s\n****\n", eQMbasis_names[qm->QMbasis]);
614 for (i = 0; i < qm->nrQMatoms; i++)
616 if (qm->atomicnumberQM[i] > 21)
618 fprintf(out, "%d ", i+1);
621 fprintf(out, "\n%s\n****\n\n", "lanl2dz");
623 for (i = 0; i < qm->nrQMatoms; i++)
625 if (qm->atomicnumberQM[i] > 21)
627 fprintf(out, "%d ", i+1);
630 fprintf(out, "\n%s\n", "lanl2dz");
635 /* MM point charge data */
636 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
638 fprintf(stderr, "nr mm atoms in gaussian.c = %d\n", mm->nrMMatoms);
640 if (qm->bTS || qm->bOPT)
642 /* freeze the frontier QM atoms and Link atoms. This is
643 * important only if a full QM subsystem optimization is done
644 * with a frozen MM environmeent. For dynamics, or gromacs's own
645 * optimization routines this is not important.
647 for (i = 0; i < qm->nrQMatoms; i++)
649 if (qm->frontatoms[i])
651 fprintf(out, "%d F\n", i+1); /* counting from 1 */
654 /* MM point charges include LJ parameters in case of QM optimization
656 for (i = 0; i < mm->nrMMatoms; i++)
659 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf 0.0 %10.7lf %10.7lf\n",
660 mm->xMM[i][XX]/BOHR2NM,
661 mm->xMM[i][YY]/BOHR2NM,
662 mm->xMM[i][ZZ]/BOHR2NM,
664 mm->c6[i], mm->c12[i]);
666 fprintf(out, "%10.7f %10.7f %10.7f %8.4f 0.0 %10.7f %10.7f\n",
667 mm->xMM[i][XX]/BOHR2NM,
668 mm->xMM[i][YY]/BOHR2NM,
669 mm->xMM[i][ZZ]/BOHR2NM,
671 mm->c6[i], mm->c12[i]);
678 for (i = 0; i < mm->nrMMatoms; i++)
681 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
682 mm->xMM[i][XX]/BOHR2NM,
683 mm->xMM[i][YY]/BOHR2NM,
684 mm->xMM[i][ZZ]/BOHR2NM,
687 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
688 mm->xMM[i][XX]/BOHR2NM,
689 mm->xMM[i][YY]/BOHR2NM,
690 mm->xMM[i][ZZ]/BOHR2NM,
701 } /* write_gaussian_input */
703 real read_gaussian_output(rvec QMgrad[], rvec MMgrad[], int step,
704 t_QMrec *qm, t_MMrec *mm)
715 in = fopen("fort.7", "r");
719 /* in case of an optimization, the coordinates are printed in the
720 * fort.7 file first, followed by the energy, coordinates and (if
721 * required) the CI eigenvectors.
723 if (qm->bTS || qm->bOPT)
725 for (i = 0; i < qm->nrQMatoms; i++)
727 if (NULL == fgets(buf, 300, in))
729 gmx_fatal(FARGS, "Error reading Gaussian output - not enough atom lines?");
733 sscanf(buf, "%d %lf %lf %lf\n",
739 sscanf(buf, "%d %f %f %f\n",
745 for (j = 0; j < DIM; j++)
747 qm->xQM[i][j] *= BOHR2NM;
751 /* the next line is the energy and in the case of CAS, the energy
752 * difference between the two states.
754 if (NULL == fgets(buf, 300, in))
756 gmx_fatal(FARGS, "Error reading Gaussian output");
760 sscanf(buf, "%lf\n", &QMener);
762 sscanf(buf, "%f\n", &QMener);
764 /* next lines contain the gradients of the QM atoms */
765 for (i = 0; i < qm->nrQMatoms; i++)
767 if (NULL == fgets(buf, 300, in))
769 gmx_fatal(FARGS, "Error reading Gaussian output");
772 sscanf(buf, "%lf %lf %lf\n",
777 sscanf(buf, "%f %f %f\n",
783 /* the next lines are the gradients of the MM atoms */
784 if (qm->QMmethod >= eQMmethodRHF)
786 for (i = 0; i < mm->nrMMatoms; i++)
788 if (NULL == fgets(buf, 300, in))
790 gmx_fatal(FARGS, "Error reading Gaussian output");
793 sscanf(buf, "%lf %lf %lf\n",
798 sscanf(buf, "%f %f %f\n",
809 real read_gaussian_SH_output(rvec QMgrad[], rvec MMgrad[], int step,
810 gmx_bool swapped, t_QMrec *qm, t_MMrec *mm)
821 in = fopen("fort.7", "r");
822 /* first line is the energy and in the case of CAS, the energy
823 * difference between the two states.
825 if (NULL == fgets(buf, 300, in))
827 gmx_fatal(FARGS, "Error reading Gaussian output");
831 sscanf(buf, "%lf %lf\n", &QMener, &DeltaE);
833 sscanf(buf, "%f %f\n", &QMener, &DeltaE);
836 /* switch on/off the State Averaging */
838 if (DeltaE > qm->SAoff)
845 else if (DeltaE < qm->SAon || (qm->SAstep > 0))
847 if (qm->SAstep < qm->SAsteps)
854 fprintf(stderr, "Gap = %5f,SA = %3d\n", DeltaE, (qm->SAstep > 0));
855 /* next lines contain the gradients of the QM atoms */
856 for (i = 0; i < qm->nrQMatoms; i++)
858 if (NULL == fgets(buf, 300, in))
860 gmx_fatal(FARGS, "Error reading Gaussian output");
864 sscanf(buf, "%lf %lf %lf\n",
869 sscanf(buf, "%f %f %f\n",
875 /* the next lines, are the gradients of the MM atoms */
877 for (i = 0; i < mm->nrMMatoms; i++)
879 if (NULL == fgets(buf, 300, in))
881 gmx_fatal(FARGS, "Error reading Gaussian output");
884 sscanf(buf, "%lf %lf %lf\n",
889 sscanf(buf, "%f %f %f\n",
896 /* the next line contains the two CI eigenvector elements */
897 if (NULL == fgets(buf, 300, in))
899 gmx_fatal(FARGS, "Error reading Gaussian output");
903 sscanf(buf, "%d", &qm->CIdim);
904 snew(qm->CIvec1, qm->CIdim);
905 snew(qm->CIvec1old, qm->CIdim);
906 snew(qm->CIvec2, qm->CIdim);
907 snew(qm->CIvec2old, qm->CIdim);
911 /* before reading in the new current CI vectors, copy the current
912 * CI vector into the old one.
914 for (i = 0; i < qm->CIdim; i++)
916 qm->CIvec1old[i] = qm->CIvec1[i];
917 qm->CIvec2old[i] = qm->CIvec2[i];
921 for (i = 0; i < qm->CIdim; i++)
923 if (NULL == fgets(buf, 300, in))
925 gmx_fatal(FARGS, "Error reading Gaussian output");
928 sscanf(buf, "%lf\n", &qm->CIvec1[i]);
930 sscanf(buf, "%f\n", &qm->CIvec1[i]);
934 for (i = 0; i < qm->CIdim; i++)
936 if (NULL == fgets(buf, 300, in))
938 gmx_fatal(FARGS, "Error reading Gaussian output");
941 sscanf(buf, "%lf\n", &qm->CIvec2[i]);
943 sscanf(buf, "%f\n", &qm->CIvec2[i]);
950 real inproduct(real *a, real *b, int n)
957 /* computes the inner product between two vectors (a.b), both of
958 * which have length n.
960 for (i = 0; i < n; i++)
967 int hop(int step, t_QMrec *qm)
972 d11 = 0.0, d12 = 0.0, d21 = 0.0, d22 = 0.0;
974 /* calculates the inproduct between the current Ci vector and the
975 * previous CI vector. A diabatic hop will be made if d12 and d21
976 * are much bigger than d11 and d22. In that case hop returns true,
977 * otherwise it returns false.
979 if (step) /* only go on if more than one step has been done */
981 d11 = inproduct(qm->CIvec1, qm->CIvec1old, qm->CIdim);
982 d12 = inproduct(qm->CIvec1, qm->CIvec2old, qm->CIdim);
983 d21 = inproduct(qm->CIvec2, qm->CIvec1old, qm->CIdim);
984 d22 = inproduct(qm->CIvec2, qm->CIvec2old, qm->CIdim);
986 fprintf(stderr, "-------------------\n");
987 fprintf(stderr, "d11 = %13.8f\n", d11);
988 fprintf(stderr, "d12 = %13.8f\n", d12);
989 fprintf(stderr, "d21 = %13.8f\n", d21);
990 fprintf(stderr, "d22 = %13.8f\n", d22);
991 fprintf(stderr, "-------------------\n");
993 if ((fabs(d12) > 0.5) && (fabs(d21) > 0.5))
1001 void do_gaussian(int step, char *exe)
1006 /* make the call to the gaussian binary through system()
1007 * The location of the binary will be picked up from the
1008 * environment using getenv().
1010 if (step) /* hack to prevent long inputfiles */
1012 sprintf(buf, "%s < %s > %s",
1019 sprintf(buf, "%s < %s > %s",
1024 fprintf(stderr, "Calling '%s'\n", buf);
1025 #ifdef GMX_NO_SYSTEM
1026 printf("Warning-- No calls to system(3) supported on this platform.");
1027 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1029 if (system(buf) != 0)
1031 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1036 real call_gaussian(t_commrec *cr, t_forcerec *fr,
1037 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
1039 /* normal gaussian jobs */
1052 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1053 snew(QMgrad, qm->nrQMatoms);
1054 snew(MMgrad, mm->nrMMatoms);
1056 write_gaussian_input(step, fr, qm, mm);
1057 do_gaussian(step, exe);
1058 QMener = read_gaussian_output(QMgrad, MMgrad, step, qm, mm);
1059 /* put the QMMM forces in the force array and to the fshift
1061 for (i = 0; i < qm->nrQMatoms; i++)
1063 for (j = 0; j < DIM; j++)
1065 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1066 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1069 for (i = 0; i < mm->nrMMatoms; i++)
1071 for (j = 0; j < DIM; j++)
1073 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1074 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1077 QMener = QMener*HARTREE2KJ*AVOGADRO;
1082 } /* call_gaussian */
1084 real call_gaussian_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
1085 rvec f[], rvec fshift[])
1087 /* a gaussian call routine intended for doing diabatic surface
1088 * "sliding". See the manual for the theoretical background of this
1098 swapped = FALSE; /* handle for identifying the current PES */
1100 swap = FALSE; /* the actual swap */
1109 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1110 /* hack to do ground state simulations */
1114 buf = getenv("STATE");
1117 sscanf(buf, "%d", &state);
1131 /* copy the QMMMrec pointer */
1132 snew(QMgrad, qm->nrQMatoms);
1133 snew(MMgrad, mm->nrMMatoms);
1134 /* at step 0 there should be no SA */
1137 /* temporray set to step + 1, since there is a chk start */
1138 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1140 do_gaussian(step, exe);
1141 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1143 /* check for a surface hop. Only possible if we were already state
1150 swap = (step && hop(step, qm));
1153 else /* already on the other surface, so check if we go back */
1155 swap = (step && hop(step, qm));
1156 swapped = !swap; /* so swapped shoud be false again */
1158 if (swap) /* change surface, so do another call */
1160 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1161 do_gaussian(step, exe);
1162 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1165 /* add the QMMM forces to the gmx force array and fshift
1167 for (i = 0; i < qm->nrQMatoms; i++)
1169 for (j = 0; j < DIM; j++)
1171 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1172 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1175 for (i = 0; i < mm->nrMMatoms; i++)
1177 for (j = 0; j < DIM; j++)
1179 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1180 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1183 QMener = QMener*HARTREE2KJ*AVOGADRO;
1184 fprintf(stderr, "step %5d, SA = %5d, swap = %5d\n",
1185 step, (qm->SAstep > 0), swapped);
1190 } /* call_gaussian_SH */
1192 /* end of gaussian sub routines */
1196 gmx_qmmm_gaussian_empty;