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
51 #include "gromacs/fileio/confio.h"
63 #include "gmx_fatal.h"
68 /* TODO: this should be made thread-safe */
70 /* Gaussian interface routines */
72 void init_gaussian(t_commrec *cr, t_QMrec *qm, t_MMrec *mm)
75 *rffile = NULL, *out = NULL;
77 basissets[eQMbasisNR] = {{0, 3, 0},
78 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
92 /* using the ivec above to convert the basis read form the mdp file
93 * in a human readable format into some numbers for the gaussian
94 * route. This is necessary as we are using non standard routes to
98 /* per layer we make a new subdir for integral file, checkpoint
99 * files and such. These dirs are stored in the QMrec for
104 if (!qm->nQMcpus) /* this we do only once per layer
105 * as we call g01 externally
109 for (i = 0; i < DIM; i++)
111 qm->SHbasis[i] = basissets[qm->QMbasis][i];
114 /* init gradually switching on of the SA */
116 /* we read the number of cpus and environment from the environment
119 buf = getenv("NCPUS");
122 sscanf(buf, "%d", &qm->nQMcpus);
128 fprintf(stderr, "number of CPUs for gaussian = %d\n", qm->nQMcpus);
132 sscanf(buf, "%d", &qm->QMmem);
136 qm->QMmem = 50000000;
138 fprintf(stderr, "memory for gaussian = %d\n", qm->QMmem);
142 sscanf(buf, "%d", &qm->accuracy);
148 fprintf(stderr, "accuracy in l510 = %d\n", qm->accuracy);
150 buf = getenv("CPMCSCF");
153 sscanf(buf, "%d", &i);
154 qm->cpmcscf = (i != 0);
162 fprintf(stderr, "using cp-mcscf in l1003\n");
166 fprintf(stderr, "NOT using cp-mcscf in l1003\n");
168 buf = getenv("SASTEP");
171 sscanf(buf, "%d", &qm->SAstep);
175 /* init gradually switching on of the SA */
178 /* we read the number of cpus and environment from the environment
181 fprintf(stderr, "Level of SA at start = %d\n", qm->SAstep);
182 /* punch the LJ C6 and C12 coefficients to be picked up by
183 * gaussian and usd to compute the LJ interaction between the
186 if (qm->bTS || qm->bOPT)
188 out = fopen("LJ.dat", "w");
189 for (i = 0; i < qm->nrQMatoms; i++)
193 fprintf(out, "%3d %10.7lf %10.7lf\n",
194 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
196 fprintf(out, "%3d %10.7f %10.7f\n",
197 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
202 /* gaussian settings on the system */
203 buf = getenv("GAUSS_DIR");
204 fprintf(stderr, "%s", buf);
208 qm->gauss_dir = strdup(buf);
212 gmx_fatal(FARGS, "no $GAUSS_DIR, check gaussian manual\n");
215 buf = getenv("GAUSS_EXE");
218 qm->gauss_exe = strdup(buf);
222 gmx_fatal(FARGS, "no $GAUSS_EXE, check gaussian manual\n");
224 buf = getenv("DEVEL_DIR");
227 qm->devel_dir = strdup (buf);
231 gmx_fatal(FARGS, "no $DEVEL_DIR, this is were the modified links reside.\n");
235 /* reactionfield, file is needed using gaussian */
236 /* rffile=fopen("rf.dat","w");*/
237 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
241 fprintf(stderr, "gaussian initialised...\n");
246 void write_gaussian_SH_input(int step, gmx_bool swap,
247 t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
258 bSA = (qm->SAstep > 0);
260 out = fopen("input.com", "w");
261 /* write the route */
262 fprintf(out, "%s", "%scr=input\n");
263 fprintf(out, "%s", "%rwf=input\n");
264 fprintf(out, "%s", "%int=input\n");
265 fprintf(out, "%s", "%d2e=input\n");
267 * fprintf(out,"%s","%nosave\n");
269 fprintf(out, "%s", "%chk=input\n");
270 fprintf(out, "%s%d\n", "%mem=", qm->QMmem);
271 fprintf(out, "%s%3d\n", "%nprocshare=", qm->nQMcpus);
273 /* use the versions of
274 * l301 that computes the interaction between MM and QM atoms.
275 * l510 that can punch the CI coefficients
276 * l701 that can do gradients on MM atoms
280 fprintf(out, "%s%s%s",
284 fprintf(out, "%s%s%s",
288 fprintf(out, "%s%s%s",
293 fprintf(out, "%s%s%s",
297 fprintf(out, "%s%s%s",
301 /* print the nonstandard route
304 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
306 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
310 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
313 qm->SHbasis[2]); /*basisset stuff */
318 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
321 qm->SHbasis[2]); /*basisset stuff */
324 if (step+1) /* fetch initial guess from check point file */
325 { /* hack, to alyays read from chk file!!!!! */
326 fprintf(out, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
328 "18=", qm->CASorbitals, "/1,5;\n");
330 else /* generate the first checkpoint file */
332 fprintf(out, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
334 "18=", qm->CASorbitals, "/1,5;\n");
336 /* the rest of the input depends on where the system is on the PES
338 if (swap && bSA) /* make a slide to the other surface */
340 if (qm->CASorbitals > 6) /* use direct and no full diag */
342 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
348 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
350 if (mm->nrMMatoms > 0)
352 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
354 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
355 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
356 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
360 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
362 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
366 else if (bSA) /* do a "state-averaged" CAS calculation */
368 if (qm->CASorbitals > 6) /* no full diag */
370 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
376 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
378 if (mm->nrMMatoms > 0)
380 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
382 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
383 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1/3;\n");
384 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
388 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
390 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
394 else if (swap) /* do a "swapped" CAS calculation */
396 if (qm->CASorbitals > 6)
398 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
402 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
405 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
407 else /* do a "normal" CAS calculation */
409 if (qm->CASorbitals > 6)
411 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
415 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
418 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
420 fprintf(out, "\ninput-file generated by gromacs\n\n");
421 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
422 for (i = 0; i < qm->nrQMatoms; i++)
425 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
426 qm->atomicnumberQM[i],
427 qm->xQM[i][XX]/BOHR2NM,
428 qm->xQM[i][YY]/BOHR2NM,
429 qm->xQM[i][ZZ]/BOHR2NM);
431 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
432 qm->atomicnumberQM[i],
433 qm->xQM[i][XX]/BOHR2NM,
434 qm->xQM[i][YY]/BOHR2NM,
435 qm->xQM[i][ZZ]/BOHR2NM);
438 /* MM point charge data */
439 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
442 for (i = 0; i < mm->nrMMatoms; i++)
445 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
446 mm->xMM[i][XX]/BOHR2NM,
447 mm->xMM[i][YY]/BOHR2NM,
448 mm->xMM[i][ZZ]/BOHR2NM,
451 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
452 mm->xMM[i][XX]/BOHR2NM,
453 mm->xMM[i][YY]/BOHR2NM,
454 mm->xMM[i][ZZ]/BOHR2NM,
459 if (bSA) /* put the SA coefficients at the end of the file */
462 fprintf(out, "\n%10.8lf %10.8lf\n",
463 qm->SAstep*0.5/qm->SAsteps,
464 1-qm->SAstep*0.5/qm->SAsteps);
466 fprintf(out, "\n%10.8f %10.8f\n",
467 qm->SAstep*0.5/qm->SAsteps,
468 1-qm->SAstep*0.5/qm->SAsteps);
470 fprintf(stderr, "State Averaging level = %d/%d\n", qm->SAstep, qm->SAsteps);
474 } /* write_gaussian_SH_input */
476 void write_gaussian_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
486 out = fopen("input.com", "w");
487 /* write the route */
489 if (qm->QMmethod >= eQMmethodRHF)
501 fprintf(out, "%s%3d\n",
502 "%nprocshare=", qm->nQMcpus);
504 fprintf(out, "%s%d\n",
506 /* use the modified links that include the LJ contribution at the QM level */
507 if (qm->bTS || qm->bOPT)
509 fprintf(out, "%s%s%s",
510 "%subst l701 ", qm->devel_dir, "/l701_LJ\n");
511 fprintf(out, "%s%s%s",
512 "%subst l301 ", qm->devel_dir, "/l301_LJ\n");
516 fprintf(out, "%s%s%s",
517 "%subst l701 ", qm->devel_dir, "/l701\n");
518 fprintf(out, "%s%s%s",
519 "%subst l301 ", qm->devel_dir, "/l301\n");
521 fprintf(out, "%s%s%s",
522 "%subst l9999 ", qm->devel_dir, "/l9999\n");
533 if (qm->QMmethod == eQMmethodB3LYPLAN)
536 "B3LYP/GEN Pseudo=Read");
541 eQMmethod_names[qm->QMmethod]);
543 if (qm->QMmethod >= eQMmethodRHF)
545 if (qm->QMmethod == eQMmethodCASSCF)
547 /* in case of cas, how many electrons and orbitals do we need?
549 fprintf(out, "(%d,%d)",
550 qm->CASelectrons, qm->CASorbitals);
553 eQMbasis_names[qm->QMbasis]);
556 if (QMMMrec->QMMMscheme == eQMMMschemenormal && mm->nrMMatoms)
561 if (step || qm->QMmethod == eQMmethodCASSCF)
563 /* fetch guess from checkpoint file, always for CASSCF */
564 fprintf(out, "%s", " guess=read");
566 fprintf(out, "\nNosymm units=bohr\n");
570 fprintf(out, "OPT=(Redundant,TS,noeigentest,ModRedundant) Punch=(Coord,Derivatives) ");
574 fprintf(out, "OPT=(Redundant,ModRedundant) Punch=(Coord,Derivatives) ");
578 fprintf(out, "FORCE Punch=(Derivatives) ");
580 fprintf(out, "iop(3/33=1)\n\n");
581 fprintf(out, "input-file generated by gromacs\n\n");
582 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
583 for (i = 0; i < qm->nrQMatoms; i++)
586 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
587 qm->atomicnumberQM[i],
588 qm->xQM[i][XX]/BOHR2NM,
589 qm->xQM[i][YY]/BOHR2NM,
590 qm->xQM[i][ZZ]/BOHR2NM);
592 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
593 qm->atomicnumberQM[i],
594 qm->xQM[i][XX]/BOHR2NM,
595 qm->xQM[i][YY]/BOHR2NM,
596 qm->xQM[i][ZZ]/BOHR2NM);
600 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
601 if (qm->QMmethod == eQMmethodB3LYPLAN)
604 for (i = 0; i < qm->nrQMatoms; i++)
606 if (qm->atomicnumberQM[i] < 21)
608 fprintf(out, "%d ", i+1);
611 fprintf(out, "\n%s\n****\n", eQMbasis_names[qm->QMbasis]);
613 for (i = 0; i < qm->nrQMatoms; i++)
615 if (qm->atomicnumberQM[i] > 21)
617 fprintf(out, "%d ", i+1);
620 fprintf(out, "\n%s\n****\n\n", "lanl2dz");
622 for (i = 0; i < qm->nrQMatoms; i++)
624 if (qm->atomicnumberQM[i] > 21)
626 fprintf(out, "%d ", i+1);
629 fprintf(out, "\n%s\n", "lanl2dz");
634 /* MM point charge data */
635 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
637 fprintf(stderr, "nr mm atoms in gaussian.c = %d\n", mm->nrMMatoms);
639 if (qm->bTS || qm->bOPT)
641 /* freeze the frontier QM atoms and Link atoms. This is
642 * important only if a full QM subsystem optimization is done
643 * with a frozen MM environmeent. For dynamics, or gromacs's own
644 * optimization routines this is not important.
646 for (i = 0; i < qm->nrQMatoms; i++)
648 if (qm->frontatoms[i])
650 fprintf(out, "%d F\n", i+1); /* counting from 1 */
653 /* MM point charges include LJ parameters in case of QM optimization
655 for (i = 0; i < mm->nrMMatoms; i++)
658 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf 0.0 %10.7lf %10.7lf\n",
659 mm->xMM[i][XX]/BOHR2NM,
660 mm->xMM[i][YY]/BOHR2NM,
661 mm->xMM[i][ZZ]/BOHR2NM,
663 mm->c6[i], mm->c12[i]);
665 fprintf(out, "%10.7f %10.7f %10.7f %8.4f 0.0 %10.7f %10.7f\n",
666 mm->xMM[i][XX]/BOHR2NM,
667 mm->xMM[i][YY]/BOHR2NM,
668 mm->xMM[i][ZZ]/BOHR2NM,
670 mm->c6[i], mm->c12[i]);
677 for (i = 0; i < mm->nrMMatoms; i++)
680 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
681 mm->xMM[i][XX]/BOHR2NM,
682 mm->xMM[i][YY]/BOHR2NM,
683 mm->xMM[i][ZZ]/BOHR2NM,
686 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
687 mm->xMM[i][XX]/BOHR2NM,
688 mm->xMM[i][YY]/BOHR2NM,
689 mm->xMM[i][ZZ]/BOHR2NM,
700 } /* write_gaussian_input */
702 real read_gaussian_output(rvec QMgrad[], rvec MMgrad[], int step,
703 t_QMrec *qm, t_MMrec *mm)
714 in = fopen("fort.7", "r");
718 /* in case of an optimization, the coordinates are printed in the
719 * fort.7 file first, followed by the energy, coordinates and (if
720 * required) the CI eigenvectors.
722 if (qm->bTS || qm->bOPT)
724 for (i = 0; i < qm->nrQMatoms; i++)
726 if (NULL == fgets(buf, 300, in))
728 gmx_fatal(FARGS, "Error reading Gaussian output - not enough atom lines?");
732 sscanf(buf, "%d %lf %lf %lf\n",
738 sscanf(buf, "%d %f %f %f\n",
744 for (j = 0; j < DIM; j++)
746 qm->xQM[i][j] *= BOHR2NM;
750 /* the next line is the energy and in the case of CAS, the energy
751 * difference between the two states.
753 if (NULL == fgets(buf, 300, in))
755 gmx_fatal(FARGS, "Error reading Gaussian output");
759 sscanf(buf, "%lf\n", &QMener);
761 sscanf(buf, "%f\n", &QMener);
763 /* next lines contain the gradients of the QM atoms */
764 for (i = 0; i < qm->nrQMatoms; i++)
766 if (NULL == fgets(buf, 300, in))
768 gmx_fatal(FARGS, "Error reading Gaussian output");
771 sscanf(buf, "%lf %lf %lf\n",
776 sscanf(buf, "%f %f %f\n",
782 /* the next lines are the gradients of the MM atoms */
783 if (qm->QMmethod >= eQMmethodRHF)
785 for (i = 0; i < mm->nrMMatoms; i++)
787 if (NULL == fgets(buf, 300, in))
789 gmx_fatal(FARGS, "Error reading Gaussian output");
792 sscanf(buf, "%lf %lf %lf\n",
797 sscanf(buf, "%f %f %f\n",
808 real read_gaussian_SH_output(rvec QMgrad[], rvec MMgrad[], int step,
809 gmx_bool swapped, t_QMrec *qm, t_MMrec *mm)
820 in = fopen("fort.7", "r");
821 /* first line is the energy and in the case of CAS, the energy
822 * difference between the two states.
824 if (NULL == fgets(buf, 300, in))
826 gmx_fatal(FARGS, "Error reading Gaussian output");
830 sscanf(buf, "%lf %lf\n", &QMener, &DeltaE);
832 sscanf(buf, "%f %f\n", &QMener, &DeltaE);
835 /* switch on/off the State Averaging */
837 if (DeltaE > qm->SAoff)
844 else if (DeltaE < qm->SAon || (qm->SAstep > 0))
846 if (qm->SAstep < qm->SAsteps)
853 fprintf(stderr, "Gap = %5f,SA = %3d\n", DeltaE, (qm->SAstep > 0));
854 /* next lines contain the gradients of the QM atoms */
855 for (i = 0; i < qm->nrQMatoms; i++)
857 if (NULL == fgets(buf, 300, in))
859 gmx_fatal(FARGS, "Error reading Gaussian output");
863 sscanf(buf, "%lf %lf %lf\n",
868 sscanf(buf, "%f %f %f\n",
874 /* the next lines, are the gradients of the MM atoms */
876 for (i = 0; i < mm->nrMMatoms; i++)
878 if (NULL == fgets(buf, 300, in))
880 gmx_fatal(FARGS, "Error reading Gaussian output");
883 sscanf(buf, "%lf %lf %lf\n",
888 sscanf(buf, "%f %f %f\n",
895 /* the next line contains the two CI eigenvector elements */
896 if (NULL == fgets(buf, 300, in))
898 gmx_fatal(FARGS, "Error reading Gaussian output");
902 sscanf(buf, "%d", &qm->CIdim);
903 snew(qm->CIvec1, qm->CIdim);
904 snew(qm->CIvec1old, qm->CIdim);
905 snew(qm->CIvec2, qm->CIdim);
906 snew(qm->CIvec2old, qm->CIdim);
910 /* before reading in the new current CI vectors, copy the current
911 * CI vector into the old one.
913 for (i = 0; i < qm->CIdim; i++)
915 qm->CIvec1old[i] = qm->CIvec1[i];
916 qm->CIvec2old[i] = qm->CIvec2[i];
920 for (i = 0; i < qm->CIdim; i++)
922 if (NULL == fgets(buf, 300, in))
924 gmx_fatal(FARGS, "Error reading Gaussian output");
927 sscanf(buf, "%lf\n", &qm->CIvec1[i]);
929 sscanf(buf, "%f\n", &qm->CIvec1[i]);
933 for (i = 0; i < qm->CIdim; i++)
935 if (NULL == fgets(buf, 300, in))
937 gmx_fatal(FARGS, "Error reading Gaussian output");
940 sscanf(buf, "%lf\n", &qm->CIvec2[i]);
942 sscanf(buf, "%f\n", &qm->CIvec2[i]);
949 real inproduct(real *a, real *b, int n)
956 /* computes the inner product between two vectors (a.b), both of
957 * which have length n.
959 for (i = 0; i < n; i++)
966 int hop(int step, t_QMrec *qm)
971 d11 = 0.0, d12 = 0.0, d21 = 0.0, d22 = 0.0;
973 /* calculates the inproduct between the current Ci vector and the
974 * previous CI vector. A diabatic hop will be made if d12 and d21
975 * are much bigger than d11 and d22. In that case hop returns true,
976 * otherwise it returns false.
978 if (step) /* only go on if more than one step has been done */
980 d11 = inproduct(qm->CIvec1, qm->CIvec1old, qm->CIdim);
981 d12 = inproduct(qm->CIvec1, qm->CIvec2old, qm->CIdim);
982 d21 = inproduct(qm->CIvec2, qm->CIvec1old, qm->CIdim);
983 d22 = inproduct(qm->CIvec2, qm->CIvec2old, qm->CIdim);
985 fprintf(stderr, "-------------------\n");
986 fprintf(stderr, "d11 = %13.8f\n", d11);
987 fprintf(stderr, "d12 = %13.8f\n", d12);
988 fprintf(stderr, "d21 = %13.8f\n", d21);
989 fprintf(stderr, "d22 = %13.8f\n", d22);
990 fprintf(stderr, "-------------------\n");
992 if ((fabs(d12) > 0.5) && (fabs(d21) > 0.5))
1000 void do_gaussian(int step, char *exe)
1005 /* make the call to the gaussian binary through system()
1006 * The location of the binary will be picked up from the
1007 * environment using getenv().
1009 if (step) /* hack to prevent long inputfiles */
1011 sprintf(buf, "%s < %s > %s",
1018 sprintf(buf, "%s < %s > %s",
1023 fprintf(stderr, "Calling '%s'\n", buf);
1024 #ifdef GMX_NO_SYSTEM
1025 printf("Warning-- No calls to system(3) supported on this platform.");
1026 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1028 if (system(buf) != 0)
1030 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1035 real call_gaussian(t_commrec *cr, t_forcerec *fr,
1036 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
1038 /* normal gaussian jobs */
1051 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1052 snew(QMgrad, qm->nrQMatoms);
1053 snew(MMgrad, mm->nrMMatoms);
1055 write_gaussian_input(step, fr, qm, mm);
1056 do_gaussian(step, exe);
1057 QMener = read_gaussian_output(QMgrad, MMgrad, step, qm, mm);
1058 /* put the QMMM forces in the force array and to the fshift
1060 for (i = 0; i < qm->nrQMatoms; i++)
1062 for (j = 0; j < DIM; j++)
1064 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1065 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1068 for (i = 0; i < mm->nrMMatoms; i++)
1070 for (j = 0; j < DIM; j++)
1072 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1073 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1076 QMener = QMener*HARTREE2KJ*AVOGADRO;
1081 } /* call_gaussian */
1083 real call_gaussian_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm,
1084 rvec f[], rvec fshift[])
1086 /* a gaussian call routine intended for doing diabatic surface
1087 * "sliding". See the manual for the theoretical background of this
1097 swapped = FALSE; /* handle for identifying the current PES */
1099 swap = FALSE; /* the actual swap */
1108 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1109 /* hack to do ground state simulations */
1113 buf = getenv("STATE");
1116 sscanf(buf, "%d", &state);
1130 /* copy the QMMMrec pointer */
1131 snew(QMgrad, qm->nrQMatoms);
1132 snew(MMgrad, mm->nrMMatoms);
1133 /* at step 0 there should be no SA */
1136 /* temporray set to step + 1, since there is a chk start */
1137 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1139 do_gaussian(step, exe);
1140 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1142 /* check for a surface hop. Only possible if we were already state
1149 swap = (step && hop(step, qm));
1152 else /* already on the other surface, so check if we go back */
1154 swap = (step && hop(step, qm));
1155 swapped = !swap; /* so swapped shoud be false again */
1157 if (swap) /* change surface, so do another call */
1159 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1160 do_gaussian(step, exe);
1161 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, swapped, qm, mm);
1164 /* add the QMMM forces to the gmx force array and fshift
1166 for (i = 0; i < qm->nrQMatoms; i++)
1168 for (j = 0; j < DIM; j++)
1170 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1171 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1174 for (i = 0; i < mm->nrMMatoms; i++)
1176 for (j = 0; j < DIM; j++)
1178 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1179 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1182 QMener = QMener*HARTREE2KJ*AVOGADRO;
1183 fprintf(stderr, "step %5d, SA = %5d, swap = %5d\n",
1184 step, (qm->SAstep > 0), swapped);
1189 } /* call_gaussian_SH */
1191 /* end of gaussian sub routines */
1195 gmx_qmmm_gaussian_empty;