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-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, 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.
41 #ifdef GMX_QMMM_GAUSSIAN
47 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/fileio/confio.h"
65 #include "gmx_fatal.h"
70 /* TODO: this should be made thread-safe */
72 /* Gaussian interface routines */
74 void init_gaussian(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
77 *rffile = NULL, *out = NULL;
79 basissets[eQMbasisNR] = {{0, 3, 0},
80 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
94 /* using the ivec above to convert the basis read form the mdp file
95 * in a human readable format into some numbers for the gaussian
96 * route. This is necessary as we are using non standard routes to
100 /* per layer we make a new subdir for integral file, checkpoint
101 * files and such. These dirs are stored in the QMrec for
106 if (!qm->nQMcpus) /* this we do only once per layer
107 * as we call g01 externally
111 for (i = 0; i < DIM; i++)
113 qm->SHbasis[i] = basissets[qm->QMbasis][i];
116 /* init gradually switching on of the SA */
118 /* we read the number of cpus and environment from the environment
121 buf = getenv("GMX_QM_GAUSSIAN_NCPUS");
124 sscanf(buf, "%d", &qm->nQMcpus);
130 fprintf(stderr, "number of CPUs for gaussian = %d\n", qm->nQMcpus);
131 buf = getenv("GMX_QM_GAUSSIAN_MEMORY");
134 sscanf(buf, "%d", &qm->QMmem);
138 qm->QMmem = 50000000;
140 fprintf(stderr, "memory for gaussian = %d\n", qm->QMmem);
141 buf = getenv("GMX_QM_ACCURACY");
144 sscanf(buf, "%d", &qm->accuracy);
150 fprintf(stderr, "accuracy in l510 = %d\n", qm->accuracy);
152 buf = getenv("GMX_QM_CPMCSCF");
155 sscanf(buf, "%d", &i);
156 qm->cpmcscf = (i != 0);
164 fprintf(stderr, "using cp-mcscf in l1003\n");
168 fprintf(stderr, "NOT using cp-mcscf in l1003\n");
170 buf = getenv("GMX_QM_SA_STEP");
173 sscanf(buf, "%d", &qm->SAstep);
177 /* init gradually switching on of the SA */
180 /* we read the number of cpus and environment from the environment
183 fprintf(stderr, "Level of SA at start = %d\n", qm->SAstep);
184 /* punch the LJ C6 and C12 coefficients to be picked up by
185 * gaussian and usd to compute the LJ interaction between the
188 if (qm->bTS || qm->bOPT)
190 out = fopen("LJ.dat", "w");
191 for (i = 0; i < qm->nrQMatoms; i++)
195 fprintf(out, "%3d %10.7lf %10.7lf\n",
196 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
198 fprintf(out, "%3d %10.7f %10.7f\n",
199 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
204 /* gaussian settings on the system */
205 buf = getenv("GMX_QM_GAUSS_DIR");
206 fprintf(stderr, "%s", buf);
210 qm->gauss_dir = strdup(buf);
214 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_DIR, check gaussian manual\n");
217 buf = getenv("GMX_QM_GAUSS_EXE");
220 qm->gauss_exe = strdup(buf);
224 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_EXE set, check gaussian manual\n");
226 buf = getenv("GMX_QM_MODIFIED_LINKS_DIR");
229 qm->devel_dir = strdup (buf);
233 gmx_fatal(FARGS, "no $GMX_QM_MODIFIED_LINKS_DIR, this is were the modified links reside.\n");
237 /* reactionfield, file is needed using gaussian */
238 /* rffile=fopen("rf.dat","w");*/
239 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
243 fprintf(stderr, "gaussian initialised...\n");
248 void write_gaussian_SH_input(int step, gmx_bool swap,
249 t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
260 bSA = (qm->SAstep > 0);
262 out = fopen("input.com", "w");
263 /* write the route */
264 fprintf(out, "%s", "%scr=input\n");
265 fprintf(out, "%s", "%rwf=input\n");
266 fprintf(out, "%s", "%int=input\n");
267 fprintf(out, "%s", "%d2e=input\n");
269 * fprintf(out,"%s","%nosave\n");
271 fprintf(out, "%s", "%chk=input\n");
272 fprintf(out, "%s%d\n", "%mem=", qm->QMmem);
273 fprintf(out, "%s%3d\n", "%nprocshare=", qm->nQMcpus);
275 /* use the versions of
276 * l301 that computes the interaction between MM and QM atoms.
277 * l510 that can punch the CI coefficients
278 * l701 that can do gradients on MM atoms
282 fprintf(out, "%s%s%s",
286 fprintf(out, "%s%s%s",
290 fprintf(out, "%s%s%s",
295 fprintf(out, "%s%s%s",
299 fprintf(out, "%s%s%s",
303 /* print the nonstandard route
306 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
308 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
312 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
315 qm->SHbasis[2]); /*basisset stuff */
320 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
323 qm->SHbasis[2]); /*basisset stuff */
326 if (step+1) /* fetch initial guess from check point file */
327 { /* hack, to alyays read from chk file!!!!! */
328 fprintf(out, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
330 "18=", qm->CASorbitals, "/1,5;\n");
332 else /* generate the first checkpoint file */
334 fprintf(out, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
336 "18=", qm->CASorbitals, "/1,5;\n");
338 /* the rest of the input depends on where the system is on the PES
340 if (swap && bSA) /* make a slide to the other surface */
342 if (qm->CASorbitals > 6) /* use direct and no full diag */
344 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
350 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
352 if (mm->nrMMatoms > 0)
354 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
356 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
357 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
358 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
362 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
364 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
368 else if (bSA) /* do a "state-averaged" CAS calculation */
370 if (qm->CASorbitals > 6) /* no full diag */
372 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
378 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
380 if (mm->nrMMatoms > 0)
382 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
384 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
385 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1/3;\n");
386 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
390 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
392 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
396 else if (swap) /* do a "swapped" CAS calculation */
398 if (qm->CASorbitals > 6)
400 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
404 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
407 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
409 else /* do a "normal" CAS calculation */
411 if (qm->CASorbitals > 6)
413 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
417 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
420 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
422 fprintf(out, "\ninput-file generated by gromacs\n\n");
423 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
424 for (i = 0; i < qm->nrQMatoms; i++)
427 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
428 qm->atomicnumberQM[i],
429 qm->xQM[i][XX]/BOHR2NM,
430 qm->xQM[i][YY]/BOHR2NM,
431 qm->xQM[i][ZZ]/BOHR2NM);
433 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
434 qm->atomicnumberQM[i],
435 qm->xQM[i][XX]/BOHR2NM,
436 qm->xQM[i][YY]/BOHR2NM,
437 qm->xQM[i][ZZ]/BOHR2NM);
440 /* MM point charge data */
441 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
444 for (i = 0; i < mm->nrMMatoms; i++)
447 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
448 mm->xMM[i][XX]/BOHR2NM,
449 mm->xMM[i][YY]/BOHR2NM,
450 mm->xMM[i][ZZ]/BOHR2NM,
453 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
454 mm->xMM[i][XX]/BOHR2NM,
455 mm->xMM[i][YY]/BOHR2NM,
456 mm->xMM[i][ZZ]/BOHR2NM,
461 if (bSA) /* put the SA coefficients at the end of the file */
464 fprintf(out, "\n%10.8lf %10.8lf\n",
465 qm->SAstep*0.5/qm->SAsteps,
466 1-qm->SAstep*0.5/qm->SAsteps);
468 fprintf(out, "\n%10.8f %10.8f\n",
469 qm->SAstep*0.5/qm->SAsteps,
470 1-qm->SAstep*0.5/qm->SAsteps);
472 fprintf(stderr, "State Averaging level = %d/%d\n", qm->SAstep, qm->SAsteps);
476 } /* write_gaussian_SH_input */
478 void write_gaussian_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
488 out = fopen("input.com", "w");
489 /* write the route */
491 if (qm->QMmethod >= eQMmethodRHF)
503 fprintf(out, "%s%3d\n",
504 "%nprocshare=", qm->nQMcpus);
506 fprintf(out, "%s%d\n",
508 /* use the modified links that include the LJ contribution at the QM level */
509 if (qm->bTS || qm->bOPT)
511 fprintf(out, "%s%s%s",
512 "%subst l701 ", qm->devel_dir, "/l701_LJ\n");
513 fprintf(out, "%s%s%s",
514 "%subst l301 ", qm->devel_dir, "/l301_LJ\n");
518 fprintf(out, "%s%s%s",
519 "%subst l701 ", qm->devel_dir, "/l701\n");
520 fprintf(out, "%s%s%s",
521 "%subst l301 ", qm->devel_dir, "/l301\n");
523 fprintf(out, "%s%s%s",
524 "%subst l9999 ", qm->devel_dir, "/l9999\n");
535 if (qm->QMmethod == eQMmethodB3LYPLAN)
538 "B3LYP/GEN Pseudo=Read");
543 eQMmethod_names[qm->QMmethod]);
545 if (qm->QMmethod >= eQMmethodRHF)
547 if (qm->QMmethod == eQMmethodCASSCF)
549 /* in case of cas, how many electrons and orbitals do we need?
551 fprintf(out, "(%d,%d)",
552 qm->CASelectrons, qm->CASorbitals);
555 eQMbasis_names[qm->QMbasis]);
558 if (QMMMrec->QMMMscheme == eQMMMschemenormal && mm->nrMMatoms)
563 if (step || qm->QMmethod == eQMmethodCASSCF)
565 /* fetch guess from checkpoint file, always for CASSCF */
566 fprintf(out, "%s", " guess=read");
568 fprintf(out, "\nNosymm units=bohr\n");
572 fprintf(out, "OPT=(Redundant,TS,noeigentest,ModRedundant) Punch=(Coord,Derivatives) ");
576 fprintf(out, "OPT=(Redundant,ModRedundant) Punch=(Coord,Derivatives) ");
580 fprintf(out, "FORCE Punch=(Derivatives) ");
582 fprintf(out, "iop(3/33=1)\n\n");
583 fprintf(out, "input-file generated by gromacs\n\n");
584 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
585 for (i = 0; i < qm->nrQMatoms; i++)
588 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
589 qm->atomicnumberQM[i],
590 qm->xQM[i][XX]/BOHR2NM,
591 qm->xQM[i][YY]/BOHR2NM,
592 qm->xQM[i][ZZ]/BOHR2NM);
594 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
595 qm->atomicnumberQM[i],
596 qm->xQM[i][XX]/BOHR2NM,
597 qm->xQM[i][YY]/BOHR2NM,
598 qm->xQM[i][ZZ]/BOHR2NM);
602 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
603 if (qm->QMmethod == eQMmethodB3LYPLAN)
606 for (i = 0; i < qm->nrQMatoms; i++)
608 if (qm->atomicnumberQM[i] < 21)
610 fprintf(out, "%d ", i+1);
613 fprintf(out, "\n%s\n****\n", eQMbasis_names[qm->QMbasis]);
615 for (i = 0; i < qm->nrQMatoms; i++)
617 if (qm->atomicnumberQM[i] > 21)
619 fprintf(out, "%d ", i+1);
622 fprintf(out, "\n%s\n****\n\n", "lanl2dz");
624 for (i = 0; i < qm->nrQMatoms; i++)
626 if (qm->atomicnumberQM[i] > 21)
628 fprintf(out, "%d ", i+1);
631 fprintf(out, "\n%s\n", "lanl2dz");
636 /* MM point charge data */
637 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
639 fprintf(stderr, "nr mm atoms in gaussian.c = %d\n", mm->nrMMatoms);
641 if (qm->bTS || qm->bOPT)
643 /* freeze the frontier QM atoms and Link atoms. This is
644 * important only if a full QM subsystem optimization is done
645 * with a frozen MM environmeent. For dynamics, or gromacs's own
646 * optimization routines this is not important.
648 for (i = 0; i < qm->nrQMatoms; i++)
650 if (qm->frontatoms[i])
652 fprintf(out, "%d F\n", i+1); /* counting from 1 */
655 /* MM point charges include LJ parameters in case of QM optimization
657 for (i = 0; i < mm->nrMMatoms; i++)
660 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf 0.0 %10.7lf %10.7lf\n",
661 mm->xMM[i][XX]/BOHR2NM,
662 mm->xMM[i][YY]/BOHR2NM,
663 mm->xMM[i][ZZ]/BOHR2NM,
665 mm->c6[i], mm->c12[i]);
667 fprintf(out, "%10.7f %10.7f %10.7f %8.4f 0.0 %10.7f %10.7f\n",
668 mm->xMM[i][XX]/BOHR2NM,
669 mm->xMM[i][YY]/BOHR2NM,
670 mm->xMM[i][ZZ]/BOHR2NM,
672 mm->c6[i], mm->c12[i]);
679 for (i = 0; i < mm->nrMMatoms; i++)
682 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
683 mm->xMM[i][XX]/BOHR2NM,
684 mm->xMM[i][YY]/BOHR2NM,
685 mm->xMM[i][ZZ]/BOHR2NM,
688 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
689 mm->xMM[i][XX]/BOHR2NM,
690 mm->xMM[i][YY]/BOHR2NM,
691 mm->xMM[i][ZZ]/BOHR2NM,
702 } /* write_gaussian_input */
704 real read_gaussian_output(rvec QMgrad[], rvec MMgrad[], int step,
705 t_QMrec *qm, t_MMrec *mm)
716 in = fopen("fort.7", "r");
720 /* in case of an optimization, the coordinates are printed in the
721 * fort.7 file first, followed by the energy, coordinates and (if
722 * required) the CI eigenvectors.
724 if (qm->bTS || qm->bOPT)
726 for (i = 0; i < qm->nrQMatoms; i++)
728 if (NULL == fgets(buf, 300, in))
730 gmx_fatal(FARGS, "Error reading Gaussian output - not enough atom lines?");
734 sscanf(buf, "%d %lf %lf %lf\n",
740 sscanf(buf, "%d %f %f %f\n",
746 for (j = 0; j < DIM; j++)
748 qm->xQM[i][j] *= BOHR2NM;
752 /* the next line is the energy and in the case of CAS, the energy
753 * difference between the two states.
755 if (NULL == fgets(buf, 300, in))
757 gmx_fatal(FARGS, "Error reading Gaussian output");
761 sscanf(buf, "%lf\n", &QMener);
763 sscanf(buf, "%f\n", &QMener);
765 /* next lines contain the gradients of the QM atoms */
766 for (i = 0; i < qm->nrQMatoms; i++)
768 if (NULL == fgets(buf, 300, in))
770 gmx_fatal(FARGS, "Error reading Gaussian output");
773 sscanf(buf, "%lf %lf %lf\n",
778 sscanf(buf, "%f %f %f\n",
784 /* the next lines are the gradients of the MM atoms */
785 if (qm->QMmethod >= eQMmethodRHF)
787 for (i = 0; i < mm->nrMMatoms; i++)
789 if (NULL == fgets(buf, 300, in))
791 gmx_fatal(FARGS, "Error reading Gaussian output");
794 sscanf(buf, "%lf %lf %lf\n",
799 sscanf(buf, "%f %f %f\n",
810 real read_gaussian_SH_output(rvec QMgrad[], rvec MMgrad[], int step,
811 gmx_bool swapped, t_QMrec *qm, t_MMrec *mm)
822 in = fopen("fort.7", "r");
823 /* first line is the energy and in the case of CAS, the energy
824 * difference between the two states.
826 if (NULL == fgets(buf, 300, in))
828 gmx_fatal(FARGS, "Error reading Gaussian output");
832 sscanf(buf, "%lf %lf\n", &QMener, &DeltaE);
834 sscanf(buf, "%f %f\n", &QMener, &DeltaE);
837 /* switch on/off the State Averaging */
839 if (DeltaE > qm->SAoff)
846 else if (DeltaE < qm->SAon || (qm->SAstep > 0))
848 if (qm->SAstep < qm->SAsteps)
855 fprintf(stderr, "Gap = %5f,SA = %3d\n", DeltaE, (qm->SAstep > 0));
856 /* next lines contain the gradients of the QM atoms */
857 for (i = 0; i < qm->nrQMatoms; i++)
859 if (NULL == fgets(buf, 300, in))
861 gmx_fatal(FARGS, "Error reading Gaussian output");
865 sscanf(buf, "%lf %lf %lf\n",
870 sscanf(buf, "%f %f %f\n",
876 /* the next lines, are the gradients of the MM atoms */
878 for (i = 0; i < mm->nrMMatoms; i++)
880 if (NULL == fgets(buf, 300, in))
882 gmx_fatal(FARGS, "Error reading Gaussian output");
885 sscanf(buf, "%lf %lf %lf\n",
890 sscanf(buf, "%f %f %f\n",
897 /* the next line contains the two CI eigenvector elements */
898 if (NULL == fgets(buf, 300, in))
900 gmx_fatal(FARGS, "Error reading Gaussian output");
904 sscanf(buf, "%d", &qm->CIdim);
905 snew(qm->CIvec1, qm->CIdim);
906 snew(qm->CIvec1old, qm->CIdim);
907 snew(qm->CIvec2, qm->CIdim);
908 snew(qm->CIvec2old, qm->CIdim);
912 /* before reading in the new current CI vectors, copy the current
913 * CI vector into the old one.
915 for (i = 0; i < qm->CIdim; i++)
917 qm->CIvec1old[i] = qm->CIvec1[i];
918 qm->CIvec2old[i] = qm->CIvec2[i];
922 for (i = 0; i < qm->CIdim; i++)
924 if (NULL == fgets(buf, 300, in))
926 gmx_fatal(FARGS, "Error reading Gaussian output");
929 sscanf(buf, "%lf\n", &qm->CIvec1[i]);
931 sscanf(buf, "%f\n", &qm->CIvec1[i]);
935 for (i = 0; i < qm->CIdim; i++)
937 if (NULL == fgets(buf, 300, in))
939 gmx_fatal(FARGS, "Error reading Gaussian output");
942 sscanf(buf, "%lf\n", &qm->CIvec2[i]);
944 sscanf(buf, "%f\n", &qm->CIvec2[i]);
951 real inproduct(real *a, real *b, int n)
958 /* computes the inner product between two vectors (a.b), both of
959 * which have length n.
961 for (i = 0; i < n; i++)
968 int hop(int step, t_QMrec *qm)
973 d11 = 0.0, d12 = 0.0, d21 = 0.0, d22 = 0.0;
975 /* calculates the inproduct between the current Ci vector and the
976 * previous CI vector. A diabatic hop will be made if d12 and d21
977 * are much bigger than d11 and d22. In that case hop returns true,
978 * otherwise it returns false.
980 if (step) /* only go on if more than one step has been done */
982 d11 = inproduct(qm->CIvec1, qm->CIvec1old, qm->CIdim);
983 d12 = inproduct(qm->CIvec1, qm->CIvec2old, qm->CIdim);
984 d21 = inproduct(qm->CIvec2, qm->CIvec1old, qm->CIdim);
985 d22 = inproduct(qm->CIvec2, qm->CIvec2old, qm->CIdim);
987 fprintf(stderr, "-------------------\n");
988 fprintf(stderr, "d11 = %13.8f\n", d11);
989 fprintf(stderr, "d12 = %13.8f\n", d12);
990 fprintf(stderr, "d21 = %13.8f\n", d21);
991 fprintf(stderr, "d22 = %13.8f\n", d22);
992 fprintf(stderr, "-------------------\n");
994 if ((fabs(d12) > 0.5) && (fabs(d21) > 0.5))
1002 void do_gaussian(int step, char *exe)
1007 /* make the call to the gaussian binary through system()
1008 * The location of the binary will be picked up from the
1009 * environment using getenv().
1011 if (step) /* hack to prevent long inputfiles */
1013 sprintf(buf, "%s < %s > %s",
1020 sprintf(buf, "%s < %s > %s",
1025 fprintf(stderr, "Calling '%s'\n", buf);
1026 #ifdef GMX_NO_SYSTEM
1027 printf("Warning-- No calls to system(3) supported on this platform.");
1028 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1030 if (system(buf) != 0)
1032 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1037 real call_gaussian(t_commrec *cr, t_forcerec *fr,
1038 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
1040 /* normal gaussian jobs */
1053 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1054 snew(QMgrad, qm->nrQMatoms);
1055 snew(MMgrad, mm->nrMMatoms);
1057 write_gaussian_input(step, fr, qm, mm);
1058 do_gaussian(step, exe);
1059 QMener = read_gaussian_output(QMgrad, MMgrad, step, qm, mm);
1060 /* put the QMMM forces in the force array and to the fshift
1062 for (i = 0; i < qm->nrQMatoms; i++)
1064 for (j = 0; j < DIM; j++)
1066 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1067 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1070 for (i = 0; i < mm->nrMMatoms; i++)
1072 for (j = 0; j < DIM; j++)
1074 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1075 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1078 QMener = QMener*HARTREE2KJ*AVOGADRO;
1083 } /* call_gaussian */
1085 real call_gaussian_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
1086 rvec f[], rvec fshift[])
1088 /* a gaussian call routine intended for doing diabatic surface
1089 * "sliding". See the manual for the theoretical background of this
1099 swapped = FALSE; /* handle for identifying the current PES */
1101 swap = FALSE; /* the actual swap */
1110 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1111 /* hack to do ground state simulations */
1115 buf = getenv("GMX_QM_GROUND_STATE");
1118 sscanf(buf, "%d", &state);
1132 /* copy the QMMMrec pointer */
1133 snew(QMgrad, qm->nrQMatoms);
1134 snew(MMgrad, mm->nrMMatoms);
1135 /* at step 0 there should be no SA */
1138 /* temporray set to step + 1, since there is a chk start */
1139 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1141 do_gaussian(step, exe);
1142 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1144 /* check for a surface hop. Only possible if we were already state
1151 swap = (step && hop(step, qm));
1154 else /* already on the other surface, so check if we go back */
1156 swap = (step && hop(step, qm));
1157 swapped = !swap; /* so swapped shoud be false again */
1159 if (swap) /* change surface, so do another call */
1161 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1162 do_gaussian(step, exe);
1163 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1166 /* add the QMMM forces to the gmx force array and fshift
1168 for (i = 0; i < qm->nrQMatoms; i++)
1170 for (j = 0; j < DIM; j++)
1172 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1173 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1176 for (i = 0; i < mm->nrMMatoms; i++)
1178 for (j = 0; j < DIM; j++)
1180 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1181 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1184 QMener = QMener*HARTREE2KJ*AVOGADRO;
1185 fprintf(stderr, "step %5d, SA = %5d, swap = %5d\n",
1186 step, (qm->SAstep > 0), swapped);
1191 } /* call_gaussian_SH */
1193 /* end of gaussian sub routines */
1197 gmx_qmmm_gaussian_empty;