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.
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/force.h"
54 #include "gromacs/mdlib/ns.h"
55 #include "gromacs/mdlib/qmmm.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/smalloc.h"
61 /* TODO: this should be made thread-safe */
63 /* Gaussian interface routines */
65 void init_gaussian(t_QMrec
*qm
)
69 basissets
[eQMbasisNR
] = {{0, 3, 0},
70 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
84 /* using the ivec above to convert the basis read form the mdp file
85 * in a human readable format into some numbers for the gaussian
86 * route. This is necessary as we are using non standard routes to
90 /* per layer we make a new subdir for integral file, checkpoint
91 * files and such. These dirs are stored in the QMrec for
96 if (!qm
->nQMcpus
) /* this we do only once per layer
97 * as we call g01 externally
101 for (i
= 0; i
< DIM
; i
++)
103 qm
->SHbasis
[i
] = basissets
[qm
->QMbasis
][i
];
106 /* init gradually switching on of the SA */
108 /* we read the number of cpus and environment from the environment
111 buf
= getenv("GMX_QM_GAUSSIAN_NCPUS");
114 sscanf(buf
, "%d", &qm
->nQMcpus
);
120 fprintf(stderr
, "number of CPUs for gaussian = %d\n", qm
->nQMcpus
);
121 buf
= getenv("GMX_QM_GAUSSIAN_MEMORY");
124 sscanf(buf
, "%d", &qm
->QMmem
);
128 qm
->QMmem
= 50000000;
130 fprintf(stderr
, "memory for gaussian = %d\n", qm
->QMmem
);
131 buf
= getenv("GMX_QM_ACCURACY");
134 sscanf(buf
, "%d", &qm
->accuracy
);
140 fprintf(stderr
, "accuracy in l510 = %d\n", qm
->accuracy
);
142 buf
= getenv("GMX_QM_CPMCSCF");
145 sscanf(buf
, "%d", &i
);
146 qm
->cpmcscf
= (i
!= 0);
154 fprintf(stderr
, "using cp-mcscf in l1003\n");
158 fprintf(stderr
, "NOT using cp-mcscf in l1003\n");
160 buf
= getenv("GMX_QM_SA_STEP");
163 sscanf(buf
, "%d", &qm
->SAstep
);
167 /* init gradually switching on of the SA */
170 /* we read the number of cpus and environment from the environment
173 fprintf(stderr
, "Level of SA at start = %d\n", qm
->SAstep
);
174 /* punch the LJ C6 and C12 coefficients to be picked up by
175 * gaussian and usd to compute the LJ interaction between the
178 if (qm
->bTS
|| qm
->bOPT
)
180 out
= fopen("LJ.dat", "w");
181 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
185 fprintf(out
, "%3d %10.7lf %10.7lf\n",
186 qm
->atomicnumberQM
[i
], qm
->c6
[i
], qm
->c12
[i
]);
188 fprintf(out
, "%3d %10.7f %10.7f\n",
189 qm
->atomicnumberQM
[i
], qm
->c6
[i
], qm
->c12
[i
]);
194 /* gaussian settings on the system */
195 buf
= getenv("GMX_QM_GAUSS_DIR");
196 fprintf(stderr
, "%s", buf
);
200 qm
->gauss_dir
= gmx_strdup(buf
);
204 gmx_fatal(FARGS
, "no $GMX_QM_GAUSS_DIR, check gaussian manual\n");
207 buf
= getenv("GMX_QM_GAUSS_EXE");
210 qm
->gauss_exe
= gmx_strdup(buf
);
214 gmx_fatal(FARGS
, "no $GMX_QM_GAUSS_EXE set, check gaussian manual\n");
216 buf
= getenv("GMX_QM_MODIFIED_LINKS_DIR");
219 qm
->devel_dir
= gmx_strdup (buf
);
223 gmx_fatal(FARGS
, "no $GMX_QM_MODIFIED_LINKS_DIR, this is were the modified links reside.\n");
227 /* reactionfield, file is needed using gaussian */
228 /* rffile=fopen("rf.dat","w");*/
229 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
233 fprintf(stderr
, "gaussian initialised...\n");
238 void write_gaussian_SH_input(int step
, gmx_bool swap
,
239 t_forcerec
*fr
, t_QMrec
*qm
, t_MMrec
*mm
)
250 bSA
= (qm
->SAstep
> 0);
252 out
= fopen("input.com", "w");
253 /* write the route */
254 fprintf(out
, "%s", "%scr=input\n");
255 fprintf(out
, "%s", "%rwf=input\n");
256 fprintf(out
, "%s", "%int=input\n");
257 fprintf(out
, "%s", "%d2e=input\n");
259 * fprintf(out,"%s","%nosave\n");
261 fprintf(out
, "%s", "%chk=input\n");
262 fprintf(out
, "%s%d\n", "%mem=", qm
->QMmem
);
263 fprintf(out
, "%s%3d\n", "%nprocshare=", qm
->nQMcpus
);
265 /* use the versions of
266 * l301 that computes the interaction between MM and QM atoms.
267 * l510 that can punch the CI coefficients
268 * l701 that can do gradients on MM atoms
272 fprintf(out
, "%s%s%s",
276 fprintf(out
, "%s%s%s",
280 fprintf(out
, "%s%s%s",
285 fprintf(out
, "%s%s%s",
289 fprintf(out
, "%s%s%s",
293 /* print the nonstandard route
296 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
298 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
302 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
305 qm
->SHbasis
[2]); /*basisset stuff */
310 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
313 qm
->SHbasis
[2]); /*basisset stuff */
316 if (step
+1) /* fetch initial guess from check point file */
317 { /* hack, to alyays read from chk file!!!!! */
318 fprintf(out
, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
320 "18=", qm
->CASorbitals
, "/1,5;\n");
322 else /* generate the first checkpoint file */
324 fprintf(out
, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
326 "18=", qm
->CASorbitals
, "/1,5;\n");
328 /* the rest of the input depends on where the system is on the PES
330 if (swap
&& bSA
) /* make a slide to the other surface */
332 if (qm
->CASorbitals
> 6) /* use direct and no full diag */
334 fprintf(out
, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
340 fprintf(out
, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
342 if (mm
->nrMMatoms
> 0)
344 fprintf(out
, " 7/7=1,16=-2,30=1/1;\n");
346 fprintf(out
, " 11/31=1,42=1,45=1/1;\n");
347 fprintf(out
, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
348 fprintf(out
, " 7/30=1/16;\n 99/10=4/99;\n");
352 fprintf(out
, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
354 fprintf(out
, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
358 else if (bSA
) /* do a "state-averaged" CAS calculation */
360 if (qm
->CASorbitals
> 6) /* no full diag */
362 fprintf(out
, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
368 fprintf(out
, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
370 if (mm
->nrMMatoms
> 0)
372 fprintf(out
, " 7/7=1,16=-2,30=1/1;\n");
374 fprintf(out
, " 11/31=1,42=1,45=1/1;\n");
375 fprintf(out
, " 10/6=1,10=700006,28=2,29=1,31=1/3;\n");
376 fprintf(out
, " 7/30=1/16;\n 99/10=4/99;\n");
380 fprintf(out
, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
382 fprintf(out
, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
386 else if (swap
) /* do a "swapped" CAS calculation */
388 if (qm
->CASorbitals
> 6)
390 fprintf(out
, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
394 fprintf(out
, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
397 fprintf(out
, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
399 else /* do a "normal" CAS calculation */
401 if (qm
->CASorbitals
> 6)
403 fprintf(out
, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
407 fprintf(out
, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
410 fprintf(out
, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
412 fprintf(out
, "\ninput-file generated by gromacs\n\n");
413 fprintf(out
, "%2d%2d\n", qm
->QMcharge
, qm
->multiplicity
);
414 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
417 fprintf(out
, "%3d %10.7lf %10.7lf %10.7lf\n",
418 qm
->atomicnumberQM
[i
],
419 qm
->xQM
[i
][XX
]/BOHR2NM
,
420 qm
->xQM
[i
][YY
]/BOHR2NM
,
421 qm
->xQM
[i
][ZZ
]/BOHR2NM
);
423 fprintf(out
, "%3d %10.7f %10.7f %10.7f\n",
424 qm
->atomicnumberQM
[i
],
425 qm
->xQM
[i
][XX
]/BOHR2NM
,
426 qm
->xQM
[i
][YY
]/BOHR2NM
,
427 qm
->xQM
[i
][ZZ
]/BOHR2NM
);
430 /* MM point charge data */
431 if (QMMMrec
->QMMMscheme
!= eQMMMschemeoniom
&& mm
->nrMMatoms
)
434 for (i
= 0; i
< mm
->nrMMatoms
; i
++)
437 fprintf(out
, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
438 mm
->xMM
[i
][XX
]/BOHR2NM
,
439 mm
->xMM
[i
][YY
]/BOHR2NM
,
440 mm
->xMM
[i
][ZZ
]/BOHR2NM
,
443 fprintf(out
, "%10.7f %10.7f %10.7f %8.4f\n",
444 mm
->xMM
[i
][XX
]/BOHR2NM
,
445 mm
->xMM
[i
][YY
]/BOHR2NM
,
446 mm
->xMM
[i
][ZZ
]/BOHR2NM
,
451 if (bSA
) /* put the SA coefficients at the end of the file */
454 fprintf(out
, "\n%10.8lf %10.8lf\n",
455 qm
->SAstep
*0.5/qm
->SAsteps
,
456 1-qm
->SAstep
*0.5/qm
->SAsteps
);
458 fprintf(out
, "\n%10.8f %10.8f\n",
459 qm
->SAstep
*0.5/qm
->SAsteps
,
460 1-qm
->SAstep
*0.5/qm
->SAsteps
);
462 fprintf(stderr
, "State Averaging level = %d/%d\n", qm
->SAstep
, qm
->SAsteps
);
466 } /* write_gaussian_SH_input */
468 void write_gaussian_input(int step
, t_forcerec
*fr
, t_QMrec
*qm
, t_MMrec
*mm
)
478 out
= fopen("input.com", "w");
479 /* write the route */
481 if (qm
->QMmethod
>= eQMmethodRHF
)
493 fprintf(out
, "%s%3d\n",
494 "%nprocshare=", qm
->nQMcpus
);
496 fprintf(out
, "%s%d\n",
498 /* use the modified links that include the LJ contribution at the QM level */
499 if (qm
->bTS
|| qm
->bOPT
)
501 fprintf(out
, "%s%s%s",
502 "%subst l701 ", qm
->devel_dir
, "/l701_LJ\n");
503 fprintf(out
, "%s%s%s",
504 "%subst l301 ", qm
->devel_dir
, "/l301_LJ\n");
508 fprintf(out
, "%s%s%s",
509 "%subst l701 ", qm
->devel_dir
, "/l701\n");
510 fprintf(out
, "%s%s%s",
511 "%subst l301 ", qm
->devel_dir
, "/l301\n");
513 fprintf(out
, "%s%s%s",
514 "%subst l9999 ", qm
->devel_dir
, "/l9999\n");
525 if (qm
->QMmethod
== eQMmethodB3LYPLAN
)
528 "B3LYP/GEN Pseudo=Read");
533 eQMmethod_names
[qm
->QMmethod
]);
535 if (qm
->QMmethod
>= eQMmethodRHF
)
537 if (qm
->QMmethod
== eQMmethodCASSCF
)
539 /* in case of cas, how many electrons and orbitals do we need?
541 fprintf(out
, "(%d,%d)",
542 qm
->CASelectrons
, qm
->CASorbitals
);
545 eQMbasis_names
[qm
->QMbasis
]);
548 if (QMMMrec
->QMMMscheme
== eQMMMschemenormal
&& mm
->nrMMatoms
)
553 if (step
|| qm
->QMmethod
== eQMmethodCASSCF
)
555 /* fetch guess from checkpoint file, always for CASSCF */
556 fprintf(out
, "%s", " guess=read");
558 fprintf(out
, "\nNosymm units=bohr\n");
562 fprintf(out
, "OPT=(Redundant,TS,noeigentest,ModRedundant) Punch=(Coord,Derivatives) ");
566 fprintf(out
, "OPT=(Redundant,ModRedundant) Punch=(Coord,Derivatives) ");
570 fprintf(out
, "FORCE Punch=(Derivatives) ");
572 fprintf(out
, "iop(3/33=1)\n\n");
573 fprintf(out
, "input-file generated by gromacs\n\n");
574 fprintf(out
, "%2d%2d\n", qm
->QMcharge
, qm
->multiplicity
);
575 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
578 fprintf(out
, "%3d %10.7lf %10.7lf %10.7lf\n",
579 qm
->atomicnumberQM
[i
],
580 qm
->xQM
[i
][XX
]/BOHR2NM
,
581 qm
->xQM
[i
][YY
]/BOHR2NM
,
582 qm
->xQM
[i
][ZZ
]/BOHR2NM
);
584 fprintf(out
, "%3d %10.7f %10.7f %10.7f\n",
585 qm
->atomicnumberQM
[i
],
586 qm
->xQM
[i
][XX
]/BOHR2NM
,
587 qm
->xQM
[i
][YY
]/BOHR2NM
,
588 qm
->xQM
[i
][ZZ
]/BOHR2NM
);
592 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
593 if (qm
->QMmethod
== eQMmethodB3LYPLAN
)
596 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
598 if (qm
->atomicnumberQM
[i
] < 21)
600 fprintf(out
, "%d ", i
+1);
603 fprintf(out
, "\n%s\n****\n", eQMbasis_names
[qm
->QMbasis
]);
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\n", "lanl2dz");
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", "lanl2dz");
626 /* MM point charge data */
627 if (QMMMrec
->QMMMscheme
!= eQMMMschemeoniom
&& mm
->nrMMatoms
)
629 fprintf(stderr
, "nr mm atoms in gaussian.c = %d\n", mm
->nrMMatoms
);
631 if (qm
->bTS
|| qm
->bOPT
)
633 /* freeze the frontier QM atoms and Link atoms. This is
634 * important only if a full QM subsystem optimization is done
635 * with a frozen MM environmeent. For dynamics, or gromacs's own
636 * optimization routines this is not important.
638 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
640 if (qm
->frontatoms
[i
])
642 fprintf(out
, "%d F\n", i
+1); /* counting from 1 */
645 /* MM point charges include LJ parameters in case of QM optimization
647 for (i
= 0; i
< mm
->nrMMatoms
; i
++)
650 fprintf(out
, "%10.7lf %10.7lf %10.7lf %8.4lf 0.0 %10.7lf %10.7lf\n",
651 mm
->xMM
[i
][XX
]/BOHR2NM
,
652 mm
->xMM
[i
][YY
]/BOHR2NM
,
653 mm
->xMM
[i
][ZZ
]/BOHR2NM
,
655 mm
->c6
[i
], mm
->c12
[i
]);
657 fprintf(out
, "%10.7f %10.7f %10.7f %8.4f 0.0 %10.7f %10.7f\n",
658 mm
->xMM
[i
][XX
]/BOHR2NM
,
659 mm
->xMM
[i
][YY
]/BOHR2NM
,
660 mm
->xMM
[i
][ZZ
]/BOHR2NM
,
662 mm
->c6
[i
], mm
->c12
[i
]);
669 for (i
= 0; i
< mm
->nrMMatoms
; i
++)
672 fprintf(out
, "%10.7lf %10.7lf %10.7lf %8.4lf\n",
673 mm
->xMM
[i
][XX
]/BOHR2NM
,
674 mm
->xMM
[i
][YY
]/BOHR2NM
,
675 mm
->xMM
[i
][ZZ
]/BOHR2NM
,
678 fprintf(out
, "%10.7f %10.7f %10.7f %8.4f\n",
679 mm
->xMM
[i
][XX
]/BOHR2NM
,
680 mm
->xMM
[i
][YY
]/BOHR2NM
,
681 mm
->xMM
[i
][ZZ
]/BOHR2NM
,
692 } /* write_gaussian_input */
694 real
read_gaussian_output(rvec QMgrad
[], rvec MMgrad
[], t_QMrec
*qm
, t_MMrec
*mm
)
705 in
= fopen("fort.7", "r");
709 /* in case of an optimization, the coordinates are printed in the
710 * fort.7 file first, followed by the energy, coordinates and (if
711 * required) the CI eigenvectors.
713 if (qm
->bTS
|| qm
->bOPT
)
715 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
717 if (NULL
== fgets(buf
, 300, in
))
719 gmx_fatal(FARGS
, "Error reading Gaussian output - not enough atom lines?");
723 sscanf(buf
, "%d %lf %lf %lf\n",
729 sscanf(buf
, "%d %f %f %f\n",
735 for (j
= 0; j
< DIM
; j
++)
737 qm
->xQM
[i
][j
] *= BOHR2NM
;
741 /* the next line is the energy and in the case of CAS, the energy
742 * difference between the two states.
744 if (NULL
== fgets(buf
, 300, in
))
746 gmx_fatal(FARGS
, "Error reading Gaussian output");
750 sscanf(buf
, "%lf\n", &QMener
);
752 sscanf(buf
, "%f\n", &QMener
);
754 /* next lines contain the gradients of the QM atoms */
755 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
757 if (NULL
== fgets(buf
, 300, in
))
759 gmx_fatal(FARGS
, "Error reading Gaussian output");
762 sscanf(buf
, "%lf %lf %lf\n",
767 sscanf(buf
, "%f %f %f\n",
773 /* the next lines are the gradients of the MM atoms */
774 if (qm
->QMmethod
>= eQMmethodRHF
)
776 for (i
= 0; i
< mm
->nrMMatoms
; i
++)
778 if (NULL
== fgets(buf
, 300, in
))
780 gmx_fatal(FARGS
, "Error reading Gaussian output");
783 sscanf(buf
, "%lf %lf %lf\n",
788 sscanf(buf
, "%f %f %f\n",
799 real
read_gaussian_SH_output(rvec QMgrad
[], rvec MMgrad
[], int step
, t_QMrec
*qm
, t_MMrec
*mm
)
810 in
= fopen("fort.7", "r");
811 /* first line is the energy and in the case of CAS, the energy
812 * difference between the two states.
814 if (NULL
== fgets(buf
, 300, in
))
816 gmx_fatal(FARGS
, "Error reading Gaussian output");
820 sscanf(buf
, "%lf %lf\n", &QMener
, &DeltaE
);
822 sscanf(buf
, "%f %f\n", &QMener
, &DeltaE
);
825 /* switch on/off the State Averaging */
827 if (DeltaE
> qm
->SAoff
)
834 else if (DeltaE
< qm
->SAon
|| (qm
->SAstep
> 0))
836 if (qm
->SAstep
< qm
->SAsteps
)
843 fprintf(stderr
, "Gap = %5f,SA = %3d\n", DeltaE
, (qm
->SAstep
> 0));
844 /* next lines contain the gradients of the QM atoms */
845 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
847 if (NULL
== fgets(buf
, 300, in
))
849 gmx_fatal(FARGS
, "Error reading Gaussian output");
853 sscanf(buf
, "%lf %lf %lf\n",
858 sscanf(buf
, "%f %f %f\n",
864 /* the next lines, are the gradients of the MM atoms */
866 for (i
= 0; i
< mm
->nrMMatoms
; i
++)
868 if (NULL
== fgets(buf
, 300, in
))
870 gmx_fatal(FARGS
, "Error reading Gaussian output");
873 sscanf(buf
, "%lf %lf %lf\n",
878 sscanf(buf
, "%f %f %f\n",
885 /* the next line contains the two CI eigenvector elements */
886 if (NULL
== fgets(buf
, 300, in
))
888 gmx_fatal(FARGS
, "Error reading Gaussian output");
892 sscanf(buf
, "%d", &qm
->CIdim
);
893 snew(qm
->CIvec1
, qm
->CIdim
);
894 snew(qm
->CIvec1old
, qm
->CIdim
);
895 snew(qm
->CIvec2
, qm
->CIdim
);
896 snew(qm
->CIvec2old
, qm
->CIdim
);
900 /* before reading in the new current CI vectors, copy the current
901 * CI vector into the old one.
903 for (i
= 0; i
< qm
->CIdim
; i
++)
905 qm
->CIvec1old
[i
] = qm
->CIvec1
[i
];
906 qm
->CIvec2old
[i
] = qm
->CIvec2
[i
];
910 for (i
= 0; i
< qm
->CIdim
; i
++)
912 if (NULL
== fgets(buf
, 300, in
))
914 gmx_fatal(FARGS
, "Error reading Gaussian output");
917 sscanf(buf
, "%lf\n", &qm
->CIvec1
[i
]);
919 sscanf(buf
, "%f\n", &qm
->CIvec1
[i
]);
923 for (i
= 0; i
< qm
->CIdim
; i
++)
925 if (NULL
== fgets(buf
, 300, in
))
927 gmx_fatal(FARGS
, "Error reading Gaussian output");
930 sscanf(buf
, "%lf\n", &qm
->CIvec2
[i
]);
932 sscanf(buf
, "%f\n", &qm
->CIvec2
[i
]);
939 real
inproduct(real
*a
, real
*b
, int n
)
946 /* computes the inner product between two vectors (a.b), both of
947 * which have length n.
949 for (i
= 0; i
< n
; i
++)
956 int hop(int step
, t_QMrec
*qm
)
961 d11
= 0.0, d12
= 0.0, d21
= 0.0, d22
= 0.0;
963 /* calculates the inproduct between the current Ci vector and the
964 * previous CI vector. A diabatic hop will be made if d12 and d21
965 * are much bigger than d11 and d22. In that case hop returns true,
966 * otherwise it returns false.
968 if (step
) /* only go on if more than one step has been done */
970 d11
= inproduct(qm
->CIvec1
, qm
->CIvec1old
, qm
->CIdim
);
971 d12
= inproduct(qm
->CIvec1
, qm
->CIvec2old
, qm
->CIdim
);
972 d21
= inproduct(qm
->CIvec2
, qm
->CIvec1old
, qm
->CIdim
);
973 d22
= inproduct(qm
->CIvec2
, qm
->CIvec2old
, qm
->CIdim
);
975 fprintf(stderr
, "-------------------\n");
976 fprintf(stderr
, "d11 = %13.8f\n", d11
);
977 fprintf(stderr
, "d12 = %13.8f\n", d12
);
978 fprintf(stderr
, "d21 = %13.8f\n", d21
);
979 fprintf(stderr
, "d22 = %13.8f\n", d22
);
980 fprintf(stderr
, "-------------------\n");
982 if ((fabs(d12
) > 0.5) && (fabs(d21
) > 0.5))
990 void do_gaussian(int step
, char *exe
)
995 /* make the call to the gaussian binary through system()
996 * The location of the binary will be picked up from the
997 * environment using getenv().
999 if (step
) /* hack to prevent long inputfiles */
1001 sprintf(buf
, "%s < %s > %s",
1008 sprintf(buf
, "%s < %s > %s",
1013 fprintf(stderr
, "Calling '%s'\n", buf
);
1014 if (system(buf
) != 0)
1016 gmx_fatal(FARGS
, "Call to '%s' failed\n", buf
);
1020 real
call_gaussian(t_forcerec
*fr
, t_QMrec
*qm
, t_MMrec
*mm
, rvec f
[], rvec fshift
[])
1022 /* normal gaussian jobs */
1035 sprintf(exe
, "%s/%s", qm
->gauss_dir
, qm
->gauss_exe
);
1036 snew(QMgrad
, qm
->nrQMatoms
);
1037 snew(MMgrad
, mm
->nrMMatoms
);
1039 write_gaussian_input(step
, fr
, qm
, mm
);
1040 do_gaussian(step
, exe
);
1041 QMener
= read_gaussian_output(QMgrad
, MMgrad
, qm
, mm
);
1042 /* put the QMMM forces in the force array and to the fshift
1044 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
1046 for (j
= 0; j
< DIM
; j
++)
1048 f
[i
][j
] = HARTREE_BOHR2MD
*QMgrad
[i
][j
];
1049 fshift
[i
][j
] = HARTREE_BOHR2MD
*QMgrad
[i
][j
];
1052 for (i
= 0; i
< mm
->nrMMatoms
; i
++)
1054 for (j
= 0; j
< DIM
; j
++)
1056 f
[i
+qm
->nrQMatoms
][j
] = HARTREE_BOHR2MD
*MMgrad
[i
][j
];
1057 fshift
[i
+qm
->nrQMatoms
][j
] = HARTREE_BOHR2MD
*MMgrad
[i
][j
];
1060 QMener
= QMener
*HARTREE2KJ
*AVOGADRO
;
1065 } /* call_gaussian */
1067 real
call_gaussian_SH(t_forcerec
*fr
, t_QMrec
*qm
, t_MMrec
*mm
, rvec f
[], rvec fshift
[])
1069 /* a gaussian call routine intended for doing diabatic surface
1070 * "sliding". See the manual for the theoretical background of this
1080 swapped
= FALSE
; /* handle for identifying the current PES */
1082 swap
= FALSE
; /* the actual swap */
1091 sprintf(exe
, "%s/%s", qm
->gauss_dir
, qm
->gauss_exe
);
1092 /* hack to do ground state simulations */
1096 buf
= getenv("GMX_QM_GROUND_STATE");
1099 sscanf(buf
, "%d", &state
);
1113 /* copy the QMMMrec pointer */
1114 snew(QMgrad
, qm
->nrQMatoms
);
1115 snew(MMgrad
, mm
->nrMMatoms
);
1116 /* at step 0 there should be no SA */
1119 /* temporray set to step + 1, since there is a chk start */
1120 write_gaussian_SH_input(step
, swapped
, fr
, qm
, mm
);
1122 do_gaussian(step
, exe
);
1123 QMener
= read_gaussian_SH_output(QMgrad
, MMgrad
, step
, qm
, mm
);
1125 /* check for a surface hop. Only possible if we were already state
1132 swap
= (step
&& hop(step
, qm
));
1135 else /* already on the other surface, so check if we go back */
1137 swap
= (step
&& hop(step
, qm
));
1138 swapped
= !swap
; /* so swapped shoud be false again */
1140 if (swap
) /* change surface, so do another call */
1142 write_gaussian_SH_input(step
, swapped
, fr
, qm
, mm
);
1143 do_gaussian(step
, exe
);
1144 QMener
= read_gaussian_SH_output(QMgrad
, MMgrad
, step
, qm
, mm
);
1147 /* add the QMMM forces to the gmx force array and fshift
1149 for (i
= 0; i
< qm
->nrQMatoms
; i
++)
1151 for (j
= 0; j
< DIM
; j
++)
1153 f
[i
][j
] = HARTREE_BOHR2MD
*QMgrad
[i
][j
];
1154 fshift
[i
][j
] = HARTREE_BOHR2MD
*QMgrad
[i
][j
];
1157 for (i
= 0; i
< mm
->nrMMatoms
; i
++)
1159 for (j
= 0; j
< DIM
; j
++)
1161 f
[i
+qm
->nrQMatoms
][j
] = HARTREE_BOHR2MD
*MMgrad
[i
][j
];
1162 fshift
[i
+qm
->nrQMatoms
][j
] = HARTREE_BOHR2MD
*MMgrad
[i
][j
];
1165 QMener
= QMener
*HARTREE2KJ
*AVOGADRO
;
1166 fprintf(stderr
, "step %5d, SA = %5d, swap = %5d\n",
1167 step
, (qm
->SAstep
> 0), swapped
);
1172 } /* call_gaussian_SH */
1174 /* end of gaussian sub routines */
1178 gmx_qmmm_gaussian_empty
;