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
50 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
54 #include "gromacs/fileio/confio.h"
62 #include "gromacs/utility/fatalerror.h"
65 /* TODO: this should be made thread-safe */
67 /* Gaussian interface routines */
69 void init_gaussian(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
72 *rffile = NULL, *out = NULL;
74 basissets[eQMbasisNR] = {{0, 3, 0},
75 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
89 /* using the ivec above to convert the basis read form the mdp file
90 * in a human readable format into some numbers for the gaussian
91 * route. This is necessary as we are using non standard routes to
95 /* per layer we make a new subdir for integral file, checkpoint
96 * files and such. These dirs are stored in the QMrec for
101 if (!qm->nQMcpus) /* this we do only once per layer
102 * as we call g01 externally
106 for (i = 0; i < DIM; i++)
108 qm->SHbasis[i] = basissets[qm->QMbasis][i];
111 /* init gradually switching on of the SA */
113 /* we read the number of cpus and environment from the environment
116 buf = getenv("GMX_QM_GAUSSIAN_NCPUS");
119 sscanf(buf, "%d", &qm->nQMcpus);
125 fprintf(stderr, "number of CPUs for gaussian = %d\n", qm->nQMcpus);
126 buf = getenv("GMX_QM_GAUSSIAN_MEMORY");
129 sscanf(buf, "%d", &qm->QMmem);
133 qm->QMmem = 50000000;
135 fprintf(stderr, "memory for gaussian = %d\n", qm->QMmem);
136 buf = getenv("GMX_QM_ACCURACY");
139 sscanf(buf, "%d", &qm->accuracy);
145 fprintf(stderr, "accuracy in l510 = %d\n", qm->accuracy);
147 buf = getenv("GMX_QM_CPMCSCF");
150 sscanf(buf, "%d", &i);
151 qm->cpmcscf = (i != 0);
159 fprintf(stderr, "using cp-mcscf in l1003\n");
163 fprintf(stderr, "NOT using cp-mcscf in l1003\n");
165 buf = getenv("GMX_QM_SA_STEP");
168 sscanf(buf, "%d", &qm->SAstep);
172 /* init gradually switching on of the SA */
175 /* we read the number of cpus and environment from the environment
178 fprintf(stderr, "Level of SA at start = %d\n", qm->SAstep);
179 /* punch the LJ C6 and C12 coefficients to be picked up by
180 * gaussian and usd to compute the LJ interaction between the
183 if (qm->bTS || qm->bOPT)
185 out = fopen("LJ.dat", "w");
186 for (i = 0; i < qm->nrQMatoms; i++)
190 fprintf(out, "%3d %10.7lf %10.7lf\n",
191 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
193 fprintf(out, "%3d %10.7f %10.7f\n",
194 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
199 /* gaussian settings on the system */
200 buf = getenv("GMX_QM_GAUSS_DIR");
201 fprintf(stderr, "%s", buf);
205 qm->gauss_dir = strdup(buf);
209 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_DIR, check gaussian manual\n");
212 buf = getenv("GMX_QM_GAUSS_EXE");
215 qm->gauss_exe = strdup(buf);
219 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_EXE set, check gaussian manual\n");
221 buf = getenv("GMX_QM_MODIFIED_LINKS_DIR");
224 qm->devel_dir = strdup (buf);
228 gmx_fatal(FARGS, "no $GMX_QM_MODIFIED_LINKS_DIR, this is were the modified links reside.\n");
232 /* reactionfield, file is needed using gaussian */
233 /* rffile=fopen("rf.dat","w");*/
234 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
238 fprintf(stderr, "gaussian initialised...\n");
243 void write_gaussian_SH_input(int step, gmx_bool swap,
244 t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
255 bSA = (qm->SAstep > 0);
257 out = fopen("input.com", "w");
258 /* write the route */
259 fprintf(out, "%s", "%scr=input\n");
260 fprintf(out, "%s", "%rwf=input\n");
261 fprintf(out, "%s", "%int=input\n");
262 fprintf(out, "%s", "%d2e=input\n");
264 * fprintf(out,"%s","%nosave\n");
266 fprintf(out, "%s", "%chk=input\n");
267 fprintf(out, "%s%d\n", "%mem=", qm->QMmem);
268 fprintf(out, "%s%3d\n", "%nprocshare=", qm->nQMcpus);
270 /* use the versions of
271 * l301 that computes the interaction between MM and QM atoms.
272 * l510 that can punch the CI coefficients
273 * l701 that can do gradients on MM atoms
277 fprintf(out, "%s%s%s",
281 fprintf(out, "%s%s%s",
285 fprintf(out, "%s%s%s",
290 fprintf(out, "%s%s%s",
294 fprintf(out, "%s%s%s",
298 /* print the nonstandard route
301 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
303 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
307 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
310 qm->SHbasis[2]); /*basisset stuff */
315 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
318 qm->SHbasis[2]); /*basisset stuff */
321 if (step+1) /* fetch initial guess from check point file */
322 { /* hack, to alyays read from chk file!!!!! */
323 fprintf(out, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
325 "18=", qm->CASorbitals, "/1,5;\n");
327 else /* generate the first checkpoint file */
329 fprintf(out, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
331 "18=", qm->CASorbitals, "/1,5;\n");
333 /* the rest of the input depends on where the system is on the PES
335 if (swap && bSA) /* make a slide to the other surface */
337 if (qm->CASorbitals > 6) /* use direct and no full diag */
339 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
345 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
347 if (mm->nrMMatoms > 0)
349 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
351 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
352 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
353 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
357 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
359 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
363 else if (bSA) /* do a "state-averaged" CAS calculation */
365 if (qm->CASorbitals > 6) /* no full diag */
367 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
373 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
375 if (mm->nrMMatoms > 0)
377 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
379 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
380 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1/3;\n");
381 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
385 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
387 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
391 else if (swap) /* do a "swapped" CAS calculation */
393 if (qm->CASorbitals > 6)
395 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
399 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
402 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
404 else /* do a "normal" CAS calculation */
406 if (qm->CASorbitals > 6)
408 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
412 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
415 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
417 fprintf(out, "\ninput-file generated by gromacs\n\n");
418 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
419 for (i = 0; i < qm->nrQMatoms; i++)
422 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
423 qm->atomicnumberQM[i],
424 qm->xQM[i][XX]/BOHR2NM,
425 qm->xQM[i][YY]/BOHR2NM,
426 qm->xQM[i][ZZ]/BOHR2NM);
428 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
429 qm->atomicnumberQM[i],
430 qm->xQM[i][XX]/BOHR2NM,
431 qm->xQM[i][YY]/BOHR2NM,
432 qm->xQM[i][ZZ]/BOHR2NM);
435 /* MM point charge data */
436 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
439 for (i = 0; i < mm->nrMMatoms; i++)
442 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
443 mm->xMM[i][XX]/BOHR2NM,
444 mm->xMM[i][YY]/BOHR2NM,
445 mm->xMM[i][ZZ]/BOHR2NM,
448 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
449 mm->xMM[i][XX]/BOHR2NM,
450 mm->xMM[i][YY]/BOHR2NM,
451 mm->xMM[i][ZZ]/BOHR2NM,
456 if (bSA) /* put the SA coefficients at the end of the file */
459 fprintf(out, "\n%10.8lf %10.8lf\n",
460 qm->SAstep*0.5/qm->SAsteps,
461 1-qm->SAstep*0.5/qm->SAsteps);
463 fprintf(out, "\n%10.8f %10.8f\n",
464 qm->SAstep*0.5/qm->SAsteps,
465 1-qm->SAstep*0.5/qm->SAsteps);
467 fprintf(stderr, "State Averaging level = %d/%d\n", qm->SAstep, qm->SAsteps);
471 } /* write_gaussian_SH_input */
473 void write_gaussian_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
483 out = fopen("input.com", "w");
484 /* write the route */
486 if (qm->QMmethod >= eQMmethodRHF)
498 fprintf(out, "%s%3d\n",
499 "%nprocshare=", qm->nQMcpus);
501 fprintf(out, "%s%d\n",
503 /* use the modified links that include the LJ contribution at the QM level */
504 if (qm->bTS || qm->bOPT)
506 fprintf(out, "%s%s%s",
507 "%subst l701 ", qm->devel_dir, "/l701_LJ\n");
508 fprintf(out, "%s%s%s",
509 "%subst l301 ", qm->devel_dir, "/l301_LJ\n");
513 fprintf(out, "%s%s%s",
514 "%subst l701 ", qm->devel_dir, "/l701\n");
515 fprintf(out, "%s%s%s",
516 "%subst l301 ", qm->devel_dir, "/l301\n");
518 fprintf(out, "%s%s%s",
519 "%subst l9999 ", qm->devel_dir, "/l9999\n");
530 if (qm->QMmethod == eQMmethodB3LYPLAN)
533 "B3LYP/GEN Pseudo=Read");
538 eQMmethod_names[qm->QMmethod]);
540 if (qm->QMmethod >= eQMmethodRHF)
542 if (qm->QMmethod == eQMmethodCASSCF)
544 /* in case of cas, how many electrons and orbitals do we need?
546 fprintf(out, "(%d,%d)",
547 qm->CASelectrons, qm->CASorbitals);
550 eQMbasis_names[qm->QMbasis]);
553 if (QMMMrec->QMMMscheme == eQMMMschemenormal && mm->nrMMatoms)
558 if (step || qm->QMmethod == eQMmethodCASSCF)
560 /* fetch guess from checkpoint file, always for CASSCF */
561 fprintf(out, "%s", " guess=read");
563 fprintf(out, "\nNosymm units=bohr\n");
567 fprintf(out, "OPT=(Redundant,TS,noeigentest,ModRedundant) Punch=(Coord,Derivatives) ");
571 fprintf(out, "OPT=(Redundant,ModRedundant) Punch=(Coord,Derivatives) ");
575 fprintf(out, "FORCE Punch=(Derivatives) ");
577 fprintf(out, "iop(3/33=1)\n\n");
578 fprintf(out, "input-file generated by gromacs\n\n");
579 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
580 for (i = 0; i < qm->nrQMatoms; i++)
583 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
584 qm->atomicnumberQM[i],
585 qm->xQM[i][XX]/BOHR2NM,
586 qm->xQM[i][YY]/BOHR2NM,
587 qm->xQM[i][ZZ]/BOHR2NM);
589 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
590 qm->atomicnumberQM[i],
591 qm->xQM[i][XX]/BOHR2NM,
592 qm->xQM[i][YY]/BOHR2NM,
593 qm->xQM[i][ZZ]/BOHR2NM);
597 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
598 if (qm->QMmethod == eQMmethodB3LYPLAN)
601 for (i = 0; i < qm->nrQMatoms; i++)
603 if (qm->atomicnumberQM[i] < 21)
605 fprintf(out, "%d ", i+1);
608 fprintf(out, "\n%s\n****\n", eQMbasis_names[qm->QMbasis]);
610 for (i = 0; i < qm->nrQMatoms; i++)
612 if (qm->atomicnumberQM[i] > 21)
614 fprintf(out, "%d ", i+1);
617 fprintf(out, "\n%s\n****\n\n", "lanl2dz");
619 for (i = 0; i < qm->nrQMatoms; i++)
621 if (qm->atomicnumberQM[i] > 21)
623 fprintf(out, "%d ", i+1);
626 fprintf(out, "\n%s\n", "lanl2dz");
631 /* MM point charge data */
632 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
634 fprintf(stderr, "nr mm atoms in gaussian.c = %d\n", mm->nrMMatoms);
636 if (qm->bTS || qm->bOPT)
638 /* freeze the frontier QM atoms and Link atoms. This is
639 * important only if a full QM subsystem optimization is done
640 * with a frozen MM environmeent. For dynamics, or gromacs's own
641 * optimization routines this is not important.
643 for (i = 0; i < qm->nrQMatoms; i++)
645 if (qm->frontatoms[i])
647 fprintf(out, "%d F\n", i+1); /* counting from 1 */
650 /* MM point charges include LJ parameters in case of QM optimization
652 for (i = 0; i < mm->nrMMatoms; i++)
655 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf 0.0 %10.7lf %10.7lf\n",
656 mm->xMM[i][XX]/BOHR2NM,
657 mm->xMM[i][YY]/BOHR2NM,
658 mm->xMM[i][ZZ]/BOHR2NM,
660 mm->c6[i], mm->c12[i]);
662 fprintf(out, "%10.7f %10.7f %10.7f %8.4f 0.0 %10.7f %10.7f\n",
663 mm->xMM[i][XX]/BOHR2NM,
664 mm->xMM[i][YY]/BOHR2NM,
665 mm->xMM[i][ZZ]/BOHR2NM,
667 mm->c6[i], mm->c12[i]);
674 for (i = 0; i < mm->nrMMatoms; i++)
677 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
678 mm->xMM[i][XX]/BOHR2NM,
679 mm->xMM[i][YY]/BOHR2NM,
680 mm->xMM[i][ZZ]/BOHR2NM,
683 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
684 mm->xMM[i][XX]/BOHR2NM,
685 mm->xMM[i][YY]/BOHR2NM,
686 mm->xMM[i][ZZ]/BOHR2NM,
697 } /* write_gaussian_input */
699 real read_gaussian_output(rvec QMgrad[], rvec MMgrad[], int step,
700 t_QMrec *qm, t_MMrec *mm)
711 in = fopen("fort.7", "r");
715 /* in case of an optimization, the coordinates are printed in the
716 * fort.7 file first, followed by the energy, coordinates and (if
717 * required) the CI eigenvectors.
719 if (qm->bTS || qm->bOPT)
721 for (i = 0; i < qm->nrQMatoms; i++)
723 if (NULL == fgets(buf, 300, in))
725 gmx_fatal(FARGS, "Error reading Gaussian output - not enough atom lines?");
729 sscanf(buf, "%d %lf %lf %lf\n",
735 sscanf(buf, "%d %f %f %f\n",
741 for (j = 0; j < DIM; j++)
743 qm->xQM[i][j] *= BOHR2NM;
747 /* the next line is the energy and in the case of CAS, the energy
748 * difference between the two states.
750 if (NULL == fgets(buf, 300, in))
752 gmx_fatal(FARGS, "Error reading Gaussian output");
756 sscanf(buf, "%lf\n", &QMener);
758 sscanf(buf, "%f\n", &QMener);
760 /* next lines contain the gradients of the QM atoms */
761 for (i = 0; i < qm->nrQMatoms; i++)
763 if (NULL == fgets(buf, 300, in))
765 gmx_fatal(FARGS, "Error reading Gaussian output");
768 sscanf(buf, "%lf %lf %lf\n",
773 sscanf(buf, "%f %f %f\n",
779 /* the next lines are the gradients of the MM atoms */
780 if (qm->QMmethod >= eQMmethodRHF)
782 for (i = 0; i < mm->nrMMatoms; i++)
784 if (NULL == fgets(buf, 300, in))
786 gmx_fatal(FARGS, "Error reading Gaussian output");
789 sscanf(buf, "%lf %lf %lf\n",
794 sscanf(buf, "%f %f %f\n",
805 real read_gaussian_SH_output(rvec QMgrad[], rvec MMgrad[], int step,
806 gmx_bool swapped, t_QMrec *qm, t_MMrec *mm)
817 in = fopen("fort.7", "r");
818 /* first line is the energy and in the case of CAS, the energy
819 * difference between the two states.
821 if (NULL == fgets(buf, 300, in))
823 gmx_fatal(FARGS, "Error reading Gaussian output");
827 sscanf(buf, "%lf %lf\n", &QMener, &DeltaE);
829 sscanf(buf, "%f %f\n", &QMener, &DeltaE);
832 /* switch on/off the State Averaging */
834 if (DeltaE > qm->SAoff)
841 else if (DeltaE < qm->SAon || (qm->SAstep > 0))
843 if (qm->SAstep < qm->SAsteps)
850 fprintf(stderr, "Gap = %5f,SA = %3d\n", DeltaE, (qm->SAstep > 0));
851 /* next lines contain the gradients of the QM atoms */
852 for (i = 0; i < qm->nrQMatoms; i++)
854 if (NULL == fgets(buf, 300, in))
856 gmx_fatal(FARGS, "Error reading Gaussian output");
860 sscanf(buf, "%lf %lf %lf\n",
865 sscanf(buf, "%f %f %f\n",
871 /* the next lines, are the gradients of the MM atoms */
873 for (i = 0; i < mm->nrMMatoms; i++)
875 if (NULL == fgets(buf, 300, in))
877 gmx_fatal(FARGS, "Error reading Gaussian output");
880 sscanf(buf, "%lf %lf %lf\n",
885 sscanf(buf, "%f %f %f\n",
892 /* the next line contains the two CI eigenvector elements */
893 if (NULL == fgets(buf, 300, in))
895 gmx_fatal(FARGS, "Error reading Gaussian output");
899 sscanf(buf, "%d", &qm->CIdim);
900 snew(qm->CIvec1, qm->CIdim);
901 snew(qm->CIvec1old, qm->CIdim);
902 snew(qm->CIvec2, qm->CIdim);
903 snew(qm->CIvec2old, qm->CIdim);
907 /* before reading in the new current CI vectors, copy the current
908 * CI vector into the old one.
910 for (i = 0; i < qm->CIdim; i++)
912 qm->CIvec1old[i] = qm->CIvec1[i];
913 qm->CIvec2old[i] = qm->CIvec2[i];
917 for (i = 0; i < qm->CIdim; i++)
919 if (NULL == fgets(buf, 300, in))
921 gmx_fatal(FARGS, "Error reading Gaussian output");
924 sscanf(buf, "%lf\n", &qm->CIvec1[i]);
926 sscanf(buf, "%f\n", &qm->CIvec1[i]);
930 for (i = 0; i < qm->CIdim; i++)
932 if (NULL == fgets(buf, 300, in))
934 gmx_fatal(FARGS, "Error reading Gaussian output");
937 sscanf(buf, "%lf\n", &qm->CIvec2[i]);
939 sscanf(buf, "%f\n", &qm->CIvec2[i]);
946 real inproduct(real *a, real *b, int n)
953 /* computes the inner product between two vectors (a.b), both of
954 * which have length n.
956 for (i = 0; i < n; i++)
963 int hop(int step, t_QMrec *qm)
968 d11 = 0.0, d12 = 0.0, d21 = 0.0, d22 = 0.0;
970 /* calculates the inproduct between the current Ci vector and the
971 * previous CI vector. A diabatic hop will be made if d12 and d21
972 * are much bigger than d11 and d22. In that case hop returns true,
973 * otherwise it returns false.
975 if (step) /* only go on if more than one step has been done */
977 d11 = inproduct(qm->CIvec1, qm->CIvec1old, qm->CIdim);
978 d12 = inproduct(qm->CIvec1, qm->CIvec2old, qm->CIdim);
979 d21 = inproduct(qm->CIvec2, qm->CIvec1old, qm->CIdim);
980 d22 = inproduct(qm->CIvec2, qm->CIvec2old, qm->CIdim);
982 fprintf(stderr, "-------------------\n");
983 fprintf(stderr, "d11 = %13.8f\n", d11);
984 fprintf(stderr, "d12 = %13.8f\n", d12);
985 fprintf(stderr, "d21 = %13.8f\n", d21);
986 fprintf(stderr, "d22 = %13.8f\n", d22);
987 fprintf(stderr, "-------------------\n");
989 if ((fabs(d12) > 0.5) && (fabs(d21) > 0.5))
997 void do_gaussian(int step, char *exe)
1002 /* make the call to the gaussian binary through system()
1003 * The location of the binary will be picked up from the
1004 * environment using getenv().
1006 if (step) /* hack to prevent long inputfiles */
1008 sprintf(buf, "%s < %s > %s",
1015 sprintf(buf, "%s < %s > %s",
1020 fprintf(stderr, "Calling '%s'\n", buf);
1021 #ifdef GMX_NO_SYSTEM
1022 printf("Warning-- No calls to system(3) supported on this platform.");
1023 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1025 if (system(buf) != 0)
1027 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1032 real call_gaussian(t_commrec *cr, t_forcerec *fr,
1033 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
1035 /* normal gaussian jobs */
1048 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1049 snew(QMgrad, qm->nrQMatoms);
1050 snew(MMgrad, mm->nrMMatoms);
1052 write_gaussian_input(step, fr, qm, mm);
1053 do_gaussian(step, exe);
1054 QMener = read_gaussian_output(QMgrad, MMgrad, step, qm, mm);
1055 /* put the QMMM forces in the force array and to the fshift
1057 for (i = 0; i < qm->nrQMatoms; i++)
1059 for (j = 0; j < DIM; j++)
1061 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1062 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1065 for (i = 0; i < mm->nrMMatoms; i++)
1067 for (j = 0; j < DIM; j++)
1069 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1070 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1073 QMener = QMener*HARTREE2KJ*AVOGADRO;
1078 } /* call_gaussian */
1080 real call_gaussian_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
1081 rvec f[], rvec fshift[])
1083 /* a gaussian call routine intended for doing diabatic surface
1084 * "sliding". See the manual for the theoretical background of this
1094 swapped = FALSE; /* handle for identifying the current PES */
1096 swap = FALSE; /* the actual swap */
1105 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1106 /* hack to do ground state simulations */
1110 buf = getenv("GMX_QM_GROUND_STATE");
1113 sscanf(buf, "%d", &state);
1127 /* copy the QMMMrec pointer */
1128 snew(QMgrad, qm->nrQMatoms);
1129 snew(MMgrad, mm->nrMMatoms);
1130 /* at step 0 there should be no SA */
1133 /* temporray set to step + 1, since there is a chk start */
1134 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1136 do_gaussian(step, exe);
1137 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1139 /* check for a surface hop. Only possible if we were already state
1146 swap = (step && hop(step, qm));
1149 else /* already on the other surface, so check if we go back */
1151 swap = (step && hop(step, qm));
1152 swapped = !swap; /* so swapped shoud be false again */
1154 if (swap) /* change surface, so do another call */
1156 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1157 do_gaussian(step, exe);
1158 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1161 /* add the QMMM forces to the gmx force array and fshift
1163 for (i = 0; i < qm->nrQMatoms; i++)
1165 for (j = 0; j < DIM; j++)
1167 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1168 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1171 for (i = 0; i < mm->nrMMatoms; i++)
1173 for (j = 0; j < DIM; j++)
1175 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1176 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1179 QMener = QMener*HARTREE2KJ*AVOGADRO;
1180 fprintf(stderr, "step %5d, SA = %5d, swap = %5d\n",
1181 step, (qm->SAstep > 0), swapped);
1186 } /* call_gaussian_SH */
1188 /* end of gaussian sub routines */
1192 gmx_qmmm_gaussian_empty;