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,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.
41 #ifdef GMX_QMMM_GAUSSIAN
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/legacyheaders/force.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/network.h"
53 #include "gromacs/legacyheaders/nrnb.h"
54 #include "gromacs/legacyheaders/ns.h"
55 #include "gromacs/legacyheaders/qmmm.h"
56 #include "gromacs/legacyheaders/txtdump.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/smalloc.h"
64 /* TODO: this should be made thread-safe */
66 /* Gaussian interface routines */
68 void init_gaussian(t_QMrec *qm)
72 basissets[eQMbasisNR] = {{0, 3, 0},
73 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
87 /* using the ivec above to convert the basis read form the mdp file
88 * in a human readable format into some numbers for the gaussian
89 * route. This is necessary as we are using non standard routes to
93 /* per layer we make a new subdir for integral file, checkpoint
94 * files and such. These dirs are stored in the QMrec for
99 if (!qm->nQMcpus) /* this we do only once per layer
100 * as we call g01 externally
104 for (i = 0; i < DIM; i++)
106 qm->SHbasis[i] = basissets[qm->QMbasis][i];
109 /* init gradually switching on of the SA */
111 /* we read the number of cpus and environment from the environment
114 buf = getenv("GMX_QM_GAUSSIAN_NCPUS");
117 sscanf(buf, "%d", &qm->nQMcpus);
123 fprintf(stderr, "number of CPUs for gaussian = %d\n", qm->nQMcpus);
124 buf = getenv("GMX_QM_GAUSSIAN_MEMORY");
127 sscanf(buf, "%d", &qm->QMmem);
131 qm->QMmem = 50000000;
133 fprintf(stderr, "memory for gaussian = %d\n", qm->QMmem);
134 buf = getenv("GMX_QM_ACCURACY");
137 sscanf(buf, "%d", &qm->accuracy);
143 fprintf(stderr, "accuracy in l510 = %d\n", qm->accuracy);
145 buf = getenv("GMX_QM_CPMCSCF");
148 sscanf(buf, "%d", &i);
149 qm->cpmcscf = (i != 0);
157 fprintf(stderr, "using cp-mcscf in l1003\n");
161 fprintf(stderr, "NOT using cp-mcscf in l1003\n");
163 buf = getenv("GMX_QM_SA_STEP");
166 sscanf(buf, "%d", &qm->SAstep);
170 /* init gradually switching on of the SA */
173 /* we read the number of cpus and environment from the environment
176 fprintf(stderr, "Level of SA at start = %d\n", qm->SAstep);
177 /* punch the LJ C6 and C12 coefficients to be picked up by
178 * gaussian and usd to compute the LJ interaction between the
181 if (qm->bTS || qm->bOPT)
183 out = fopen("LJ.dat", "w");
184 for (i = 0; i < qm->nrQMatoms; i++)
188 fprintf(out, "%3d %10.7lf %10.7lf\n",
189 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
191 fprintf(out, "%3d %10.7f %10.7f\n",
192 qm->atomicnumberQM[i], qm->c6[i], qm->c12[i]);
197 /* gaussian settings on the system */
198 buf = getenv("GMX_QM_GAUSS_DIR");
199 fprintf(stderr, "%s", buf);
203 qm->gauss_dir = gmx_strdup(buf);
207 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_DIR, check gaussian manual\n");
210 buf = getenv("GMX_QM_GAUSS_EXE");
213 qm->gauss_exe = gmx_strdup(buf);
217 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_EXE set, check gaussian manual\n");
219 buf = getenv("GMX_QM_MODIFIED_LINKS_DIR");
222 qm->devel_dir = gmx_strdup (buf);
226 gmx_fatal(FARGS, "no $GMX_QM_MODIFIED_LINKS_DIR, this is were the modified links reside.\n");
230 /* reactionfield, file is needed using gaussian */
231 /* rffile=fopen("rf.dat","w");*/
232 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
236 fprintf(stderr, "gaussian initialised...\n");
241 void write_gaussian_SH_input(int step, gmx_bool swap,
242 t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
253 bSA = (qm->SAstep > 0);
255 out = fopen("input.com", "w");
256 /* write the route */
257 fprintf(out, "%s", "%scr=input\n");
258 fprintf(out, "%s", "%rwf=input\n");
259 fprintf(out, "%s", "%int=input\n");
260 fprintf(out, "%s", "%d2e=input\n");
262 * fprintf(out,"%s","%nosave\n");
264 fprintf(out, "%s", "%chk=input\n");
265 fprintf(out, "%s%d\n", "%mem=", qm->QMmem);
266 fprintf(out, "%s%3d\n", "%nprocshare=", qm->nQMcpus);
268 /* use the versions of
269 * l301 that computes the interaction between MM and QM atoms.
270 * l510 that can punch the CI coefficients
271 * l701 that can do gradients on MM atoms
275 fprintf(out, "%s%s%s",
279 fprintf(out, "%s%s%s",
283 fprintf(out, "%s%s%s",
288 fprintf(out, "%s%s%s",
292 fprintf(out, "%s%s%s",
296 /* print the nonstandard route
299 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
301 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
305 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
308 qm->SHbasis[2]); /*basisset stuff */
313 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
316 qm->SHbasis[2]); /*basisset stuff */
319 if (step+1) /* fetch initial guess from check point file */
320 { /* hack, to alyays read from chk file!!!!! */
321 fprintf(out, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
323 "18=", qm->CASorbitals, "/1,5;\n");
325 else /* generate the first checkpoint file */
327 fprintf(out, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
329 "18=", qm->CASorbitals, "/1,5;\n");
331 /* the rest of the input depends on where the system is on the PES
333 if (swap && bSA) /* make a slide to the other surface */
335 if (qm->CASorbitals > 6) /* use direct and no full diag */
337 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
343 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
345 if (mm->nrMMatoms > 0)
347 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
349 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
350 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
351 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
355 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
357 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
361 else if (bSA) /* do a "state-averaged" CAS calculation */
363 if (qm->CASorbitals > 6) /* no full diag */
365 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
371 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
373 if (mm->nrMMatoms > 0)
375 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
377 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
378 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1/3;\n");
379 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
383 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
385 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
389 else if (swap) /* do a "swapped" CAS calculation */
391 if (qm->CASorbitals > 6)
393 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
397 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
400 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
402 else /* do a "normal" CAS calculation */
404 if (qm->CASorbitals > 6)
406 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
410 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
413 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
415 fprintf(out, "\ninput-file generated by gromacs\n\n");
416 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
417 for (i = 0; i < qm->nrQMatoms; i++)
420 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
421 qm->atomicnumberQM[i],
422 qm->xQM[i][XX]/BOHR2NM,
423 qm->xQM[i][YY]/BOHR2NM,
424 qm->xQM[i][ZZ]/BOHR2NM);
426 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
427 qm->atomicnumberQM[i],
428 qm->xQM[i][XX]/BOHR2NM,
429 qm->xQM[i][YY]/BOHR2NM,
430 qm->xQM[i][ZZ]/BOHR2NM);
433 /* MM point charge data */
434 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
437 for (i = 0; i < mm->nrMMatoms; i++)
440 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
441 mm->xMM[i][XX]/BOHR2NM,
442 mm->xMM[i][YY]/BOHR2NM,
443 mm->xMM[i][ZZ]/BOHR2NM,
446 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
447 mm->xMM[i][XX]/BOHR2NM,
448 mm->xMM[i][YY]/BOHR2NM,
449 mm->xMM[i][ZZ]/BOHR2NM,
454 if (bSA) /* put the SA coefficients at the end of the file */
457 fprintf(out, "\n%10.8lf %10.8lf\n",
458 qm->SAstep*0.5/qm->SAsteps,
459 1-qm->SAstep*0.5/qm->SAsteps);
461 fprintf(out, "\n%10.8f %10.8f\n",
462 qm->SAstep*0.5/qm->SAsteps,
463 1-qm->SAstep*0.5/qm->SAsteps);
465 fprintf(stderr, "State Averaging level = %d/%d\n", qm->SAstep, qm->SAsteps);
469 } /* write_gaussian_SH_input */
471 void write_gaussian_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
481 out = fopen("input.com", "w");
482 /* write the route */
484 if (qm->QMmethod >= eQMmethodRHF)
496 fprintf(out, "%s%3d\n",
497 "%nprocshare=", qm->nQMcpus);
499 fprintf(out, "%s%d\n",
501 /* use the modified links that include the LJ contribution at the QM level */
502 if (qm->bTS || qm->bOPT)
504 fprintf(out, "%s%s%s",
505 "%subst l701 ", qm->devel_dir, "/l701_LJ\n");
506 fprintf(out, "%s%s%s",
507 "%subst l301 ", qm->devel_dir, "/l301_LJ\n");
511 fprintf(out, "%s%s%s",
512 "%subst l701 ", qm->devel_dir, "/l701\n");
513 fprintf(out, "%s%s%s",
514 "%subst l301 ", qm->devel_dir, "/l301\n");
516 fprintf(out, "%s%s%s",
517 "%subst l9999 ", qm->devel_dir, "/l9999\n");
528 if (qm->QMmethod == eQMmethodB3LYPLAN)
531 "B3LYP/GEN Pseudo=Read");
536 eQMmethod_names[qm->QMmethod]);
538 if (qm->QMmethod >= eQMmethodRHF)
540 if (qm->QMmethod == eQMmethodCASSCF)
542 /* in case of cas, how many electrons and orbitals do we need?
544 fprintf(out, "(%d,%d)",
545 qm->CASelectrons, qm->CASorbitals);
548 eQMbasis_names[qm->QMbasis]);
551 if (QMMMrec->QMMMscheme == eQMMMschemenormal && mm->nrMMatoms)
556 if (step || qm->QMmethod == eQMmethodCASSCF)
558 /* fetch guess from checkpoint file, always for CASSCF */
559 fprintf(out, "%s", " guess=read");
561 fprintf(out, "\nNosymm units=bohr\n");
565 fprintf(out, "OPT=(Redundant,TS,noeigentest,ModRedundant) Punch=(Coord,Derivatives) ");
569 fprintf(out, "OPT=(Redundant,ModRedundant) Punch=(Coord,Derivatives) ");
573 fprintf(out, "FORCE Punch=(Derivatives) ");
575 fprintf(out, "iop(3/33=1)\n\n");
576 fprintf(out, "input-file generated by gromacs\n\n");
577 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
578 for (i = 0; i < qm->nrQMatoms; i++)
581 fprintf(out, "%3d %10.7lf %10.7lf %10.7lf\n",
582 qm->atomicnumberQM[i],
583 qm->xQM[i][XX]/BOHR2NM,
584 qm->xQM[i][YY]/BOHR2NM,
585 qm->xQM[i][ZZ]/BOHR2NM);
587 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
588 qm->atomicnumberQM[i],
589 qm->xQM[i][XX]/BOHR2NM,
590 qm->xQM[i][YY]/BOHR2NM,
591 qm->xQM[i][ZZ]/BOHR2NM);
595 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
596 if (qm->QMmethod == eQMmethodB3LYPLAN)
599 for (i = 0; i < qm->nrQMatoms; i++)
601 if (qm->atomicnumberQM[i] < 21)
603 fprintf(out, "%d ", i+1);
606 fprintf(out, "\n%s\n****\n", eQMbasis_names[qm->QMbasis]);
608 for (i = 0; i < qm->nrQMatoms; i++)
610 if (qm->atomicnumberQM[i] > 21)
612 fprintf(out, "%d ", i+1);
615 fprintf(out, "\n%s\n****\n\n", "lanl2dz");
617 for (i = 0; i < qm->nrQMatoms; i++)
619 if (qm->atomicnumberQM[i] > 21)
621 fprintf(out, "%d ", i+1);
624 fprintf(out, "\n%s\n", "lanl2dz");
629 /* MM point charge data */
630 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
632 fprintf(stderr, "nr mm atoms in gaussian.c = %d\n", mm->nrMMatoms);
634 if (qm->bTS || qm->bOPT)
636 /* freeze the frontier QM atoms and Link atoms. This is
637 * important only if a full QM subsystem optimization is done
638 * with a frozen MM environmeent. For dynamics, or gromacs's own
639 * optimization routines this is not important.
641 for (i = 0; i < qm->nrQMatoms; i++)
643 if (qm->frontatoms[i])
645 fprintf(out, "%d F\n", i+1); /* counting from 1 */
648 /* MM point charges include LJ parameters in case of QM optimization
650 for (i = 0; i < mm->nrMMatoms; i++)
653 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf 0.0 %10.7lf %10.7lf\n",
654 mm->xMM[i][XX]/BOHR2NM,
655 mm->xMM[i][YY]/BOHR2NM,
656 mm->xMM[i][ZZ]/BOHR2NM,
658 mm->c6[i], mm->c12[i]);
660 fprintf(out, "%10.7f %10.7f %10.7f %8.4f 0.0 %10.7f %10.7f\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]);
672 for (i = 0; i < mm->nrMMatoms; i++)
675 fprintf(out, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
676 mm->xMM[i][XX]/BOHR2NM,
677 mm->xMM[i][YY]/BOHR2NM,
678 mm->xMM[i][ZZ]/BOHR2NM,
681 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
682 mm->xMM[i][XX]/BOHR2NM,
683 mm->xMM[i][YY]/BOHR2NM,
684 mm->xMM[i][ZZ]/BOHR2NM,
695 } /* write_gaussian_input */
697 real read_gaussian_output(rvec QMgrad[], rvec MMgrad[], t_QMrec *qm, t_MMrec *mm)
708 in = fopen("fort.7", "r");
712 /* in case of an optimization, the coordinates are printed in the
713 * fort.7 file first, followed by the energy, coordinates and (if
714 * required) the CI eigenvectors.
716 if (qm->bTS || qm->bOPT)
718 for (i = 0; i < qm->nrQMatoms; i++)
720 if (NULL == fgets(buf, 300, in))
722 gmx_fatal(FARGS, "Error reading Gaussian output - not enough atom lines?");
726 sscanf(buf, "%d %lf %lf %lf\n",
732 sscanf(buf, "%d %f %f %f\n",
738 for (j = 0; j < DIM; j++)
740 qm->xQM[i][j] *= BOHR2NM;
744 /* the next line is the energy and in the case of CAS, the energy
745 * difference between the two states.
747 if (NULL == fgets(buf, 300, in))
749 gmx_fatal(FARGS, "Error reading Gaussian output");
753 sscanf(buf, "%lf\n", &QMener);
755 sscanf(buf, "%f\n", &QMener);
757 /* next lines contain the gradients of the QM atoms */
758 for (i = 0; i < qm->nrQMatoms; i++)
760 if (NULL == fgets(buf, 300, in))
762 gmx_fatal(FARGS, "Error reading Gaussian output");
765 sscanf(buf, "%lf %lf %lf\n",
770 sscanf(buf, "%f %f %f\n",
776 /* the next lines are the gradients of the MM atoms */
777 if (qm->QMmethod >= eQMmethodRHF)
779 for (i = 0; i < mm->nrMMatoms; i++)
781 if (NULL == fgets(buf, 300, in))
783 gmx_fatal(FARGS, "Error reading Gaussian output");
786 sscanf(buf, "%lf %lf %lf\n",
791 sscanf(buf, "%f %f %f\n",
802 real read_gaussian_SH_output(rvec QMgrad[], rvec MMgrad[], int step, t_QMrec *qm, t_MMrec *mm)
813 in = fopen("fort.7", "r");
814 /* first line is the energy and in the case of CAS, the energy
815 * difference between the two states.
817 if (NULL == fgets(buf, 300, in))
819 gmx_fatal(FARGS, "Error reading Gaussian output");
823 sscanf(buf, "%lf %lf\n", &QMener, &DeltaE);
825 sscanf(buf, "%f %f\n", &QMener, &DeltaE);
828 /* switch on/off the State Averaging */
830 if (DeltaE > qm->SAoff)
837 else if (DeltaE < qm->SAon || (qm->SAstep > 0))
839 if (qm->SAstep < qm->SAsteps)
846 fprintf(stderr, "Gap = %5f,SA = %3d\n", DeltaE, (qm->SAstep > 0));
847 /* next lines contain the gradients of the QM atoms */
848 for (i = 0; i < qm->nrQMatoms; i++)
850 if (NULL == fgets(buf, 300, in))
852 gmx_fatal(FARGS, "Error reading Gaussian output");
856 sscanf(buf, "%lf %lf %lf\n",
861 sscanf(buf, "%f %f %f\n",
867 /* the next lines, are the gradients of the MM atoms */
869 for (i = 0; i < mm->nrMMatoms; i++)
871 if (NULL == fgets(buf, 300, in))
873 gmx_fatal(FARGS, "Error reading Gaussian output");
876 sscanf(buf, "%lf %lf %lf\n",
881 sscanf(buf, "%f %f %f\n",
888 /* the next line contains the two CI eigenvector elements */
889 if (NULL == fgets(buf, 300, in))
891 gmx_fatal(FARGS, "Error reading Gaussian output");
895 sscanf(buf, "%d", &qm->CIdim);
896 snew(qm->CIvec1, qm->CIdim);
897 snew(qm->CIvec1old, qm->CIdim);
898 snew(qm->CIvec2, qm->CIdim);
899 snew(qm->CIvec2old, qm->CIdim);
903 /* before reading in the new current CI vectors, copy the current
904 * CI vector into the old one.
906 for (i = 0; i < qm->CIdim; i++)
908 qm->CIvec1old[i] = qm->CIvec1[i];
909 qm->CIvec2old[i] = qm->CIvec2[i];
913 for (i = 0; i < qm->CIdim; i++)
915 if (NULL == fgets(buf, 300, in))
917 gmx_fatal(FARGS, "Error reading Gaussian output");
920 sscanf(buf, "%lf\n", &qm->CIvec1[i]);
922 sscanf(buf, "%f\n", &qm->CIvec1[i]);
926 for (i = 0; i < qm->CIdim; i++)
928 if (NULL == fgets(buf, 300, in))
930 gmx_fatal(FARGS, "Error reading Gaussian output");
933 sscanf(buf, "%lf\n", &qm->CIvec2[i]);
935 sscanf(buf, "%f\n", &qm->CIvec2[i]);
942 real inproduct(real *a, real *b, int n)
949 /* computes the inner product between two vectors (a.b), both of
950 * which have length n.
952 for (i = 0; i < n; i++)
959 int hop(int step, t_QMrec *qm)
964 d11 = 0.0, d12 = 0.0, d21 = 0.0, d22 = 0.0;
966 /* calculates the inproduct between the current Ci vector and the
967 * previous CI vector. A diabatic hop will be made if d12 and d21
968 * are much bigger than d11 and d22. In that case hop returns true,
969 * otherwise it returns false.
971 if (step) /* only go on if more than one step has been done */
973 d11 = inproduct(qm->CIvec1, qm->CIvec1old, qm->CIdim);
974 d12 = inproduct(qm->CIvec1, qm->CIvec2old, qm->CIdim);
975 d21 = inproduct(qm->CIvec2, qm->CIvec1old, qm->CIdim);
976 d22 = inproduct(qm->CIvec2, qm->CIvec2old, qm->CIdim);
978 fprintf(stderr, "-------------------\n");
979 fprintf(stderr, "d11 = %13.8f\n", d11);
980 fprintf(stderr, "d12 = %13.8f\n", d12);
981 fprintf(stderr, "d21 = %13.8f\n", d21);
982 fprintf(stderr, "d22 = %13.8f\n", d22);
983 fprintf(stderr, "-------------------\n");
985 if ((fabs(d12) > 0.5) && (fabs(d21) > 0.5))
993 void do_gaussian(int step, char *exe)
998 /* make the call to the gaussian binary through system()
999 * The location of the binary will be picked up from the
1000 * environment using getenv().
1002 if (step) /* hack to prevent long inputfiles */
1004 sprintf(buf, "%s < %s > %s",
1011 sprintf(buf, "%s < %s > %s",
1016 fprintf(stderr, "Calling '%s'\n", buf);
1017 if (system(buf) != 0)
1019 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
1023 real call_gaussian(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
1025 /* normal gaussian jobs */
1038 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1039 snew(QMgrad, qm->nrQMatoms);
1040 snew(MMgrad, mm->nrMMatoms);
1042 write_gaussian_input(step, fr, qm, mm);
1043 do_gaussian(step, exe);
1044 QMener = read_gaussian_output(QMgrad, MMgrad, qm, mm);
1045 /* put the QMMM forces in the force array and to the fshift
1047 for (i = 0; i < qm->nrQMatoms; i++)
1049 for (j = 0; j < DIM; j++)
1051 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1052 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1055 for (i = 0; i < mm->nrMMatoms; i++)
1057 for (j = 0; j < DIM; j++)
1059 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1060 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1063 QMener = QMener*HARTREE2KJ*AVOGADRO;
1068 } /* call_gaussian */
1070 real call_gaussian_SH(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
1072 /* a gaussian call routine intended for doing diabatic surface
1073 * "sliding". See the manual for the theoretical background of this
1083 swapped = FALSE; /* handle for identifying the current PES */
1085 swap = FALSE; /* the actual swap */
1094 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
1095 /* hack to do ground state simulations */
1099 buf = getenv("GMX_QM_GROUND_STATE");
1102 sscanf(buf, "%d", &state);
1116 /* copy the QMMMrec pointer */
1117 snew(QMgrad, qm->nrQMatoms);
1118 snew(MMgrad, mm->nrMMatoms);
1119 /* at step 0 there should be no SA */
1122 /* temporray set to step + 1, since there is a chk start */
1123 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1125 do_gaussian(step, exe);
1126 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, qm, mm);
1128 /* check for a surface hop. Only possible if we were already state
1135 swap = (step && hop(step, qm));
1138 else /* already on the other surface, so check if we go back */
1140 swap = (step && hop(step, qm));
1141 swapped = !swap; /* so swapped shoud be false again */
1143 if (swap) /* change surface, so do another call */
1145 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1146 do_gaussian(step, exe);
1147 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, qm, mm);
1150 /* add the QMMM forces to the gmx force array and fshift
1152 for (i = 0; i < qm->nrQMatoms; i++)
1154 for (j = 0; j < DIM; j++)
1156 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1157 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1160 for (i = 0; i < mm->nrMMatoms; i++)
1162 for (j = 0; j < DIM; j++)
1164 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1165 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1168 QMener = QMener*HARTREE2KJ*AVOGADRO;
1169 fprintf(stderr, "step %5d, SA = %5d, swap = %5d\n",
1170 step, (qm->SAstep > 0), swapped);
1175 } /* call_gaussian_SH */
1177 /* end of gaussian sub routines */
1181 gmx_qmmm_gaussian_empty;