Prepared t_mdatoms for using vector
[gromacs.git] / src / gromacs / mdlib / qm_gaussian.cpp
blobe09d641b4a97d24c41245d5f8aab6c392693377a
1 /*
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,2016,2017, 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.
37 #include "gmxpre.h"
39 #include "config.h"
41 #if GMX_QMMM_GAUSSIAN
43 #include <stdio.h>
44 #include <stdlib.h>
45 #include <string.h>
47 #include <cmath>
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/gmxlib/nrnb.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdlib/force.h"
55 #include "gromacs/mdlib/ns.h"
56 #include "gromacs/mdlib/qmmm.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
62 /* TODO: this should be made thread-safe */
64 /* Gaussian interface routines */
66 void init_gaussian(t_QMrec *qm)
68 FILE *out = NULL;
69 ivec
70 basissets[eQMbasisNR] = {{0, 3, 0},
71 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
72 {5, 0, 0},
73 {5, 0, 1},
74 {5, 0, 11},
75 {5, 6, 0},
76 {1, 6, 0},
77 {1, 6, 1},
78 {1, 6, 11},
79 {4, 6, 0}};
80 char
81 *buf = NULL;
82 int
85 /* using the ivec above to convert the basis read form the mdp file
86 * in a human readable format into some numbers for the gaussian
87 * route. This is necessary as we are using non standard routes to
88 * do SH.
91 /* per layer we make a new subdir for integral file, checkpoint
92 * files and such. These dirs are stored in the QMrec for
93 * convenience
97 if (!qm->nQMcpus) /* this we do only once per layer
98 * as we call g01 externally
102 for (i = 0; i < DIM; i++)
104 qm->SHbasis[i] = basissets[qm->QMbasis][i];
107 /* init gradually switching on of the SA */
108 qm->SAstep = 0;
109 /* we read the number of cpus and environment from the environment
110 * if set.
112 buf = getenv("GMX_QM_GAUSSIAN_NCPUS");
113 if (buf)
115 sscanf(buf, "%d", &qm->nQMcpus);
117 else
119 qm->nQMcpus = 1;
121 fprintf(stderr, "number of CPUs for gaussian = %d\n", qm->nQMcpus);
122 buf = getenv("GMX_QM_GAUSSIAN_MEMORY");
123 if (buf)
125 sscanf(buf, "%d", &qm->QMmem);
127 else
129 qm->QMmem = 50000000;
131 fprintf(stderr, "memory for gaussian = %d\n", qm->QMmem);
132 buf = getenv("GMX_QM_ACCURACY");
133 if (buf)
135 sscanf(buf, "%d", &qm->accuracy);
137 else
139 qm->accuracy = 8;
141 fprintf(stderr, "accuracy in l510 = %d\n", qm->accuracy);
143 buf = getenv("GMX_QM_CPMCSCF");
144 if (buf)
146 sscanf(buf, "%d", &i);
147 qm->cpmcscf = (i != 0);
149 else
151 qm->cpmcscf = FALSE;
153 if (qm->cpmcscf)
155 fprintf(stderr, "using cp-mcscf in l1003\n");
157 else
159 fprintf(stderr, "NOT using cp-mcscf in l1003\n");
161 buf = getenv("GMX_QM_SA_STEP");
162 if (buf)
164 sscanf(buf, "%d", &qm->SAstep);
166 else
168 /* init gradually switching on of the SA */
169 qm->SAstep = 0;
171 /* we read the number of cpus and environment from the environment
172 * if set.
174 fprintf(stderr, "Level of SA at start = %d\n", qm->SAstep);
175 /* gaussian settings on the system */
176 buf = getenv("GMX_QM_GAUSS_DIR");
178 if (buf)
180 qm->gauss_dir = gmx_strdup(buf);
182 else
184 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_DIR, check gaussian manual\n");
187 buf = getenv("GMX_QM_GAUSS_EXE");
188 if (buf)
190 qm->gauss_exe = gmx_strdup(buf);
192 else
194 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_EXE set, check gaussian manual\n");
196 buf = getenv("GMX_QM_MODIFIED_LINKS_DIR");
197 if (buf)
199 qm->devel_dir = gmx_strdup (buf);
201 else
203 gmx_fatal(FARGS, "no $GMX_QM_MODIFIED_LINKS_DIR, this is were the modified links reside.\n");
206 /* if(fr->bRF){*/
207 /* reactionfield, file is needed using gaussian */
208 /* rffile=fopen("rf.dat","w");*/
209 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
210 /* fclose(rffile);*/
211 /* }*/
213 fprintf(stderr, "gaussian initialised...\n");
218 void write_gaussian_SH_input(int step, gmx_bool swap,
219 t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
223 gmx_bool
224 bSA;
225 FILE
226 *out;
227 t_QMMMrec
228 *QMMMrec;
229 QMMMrec = fr->qr;
230 bSA = (qm->SAstep > 0);
232 out = fopen("input.com", "w");
233 /* write the route */
234 fprintf(out, "%s", "%scr=input\n");
235 fprintf(out, "%s", "%rwf=input\n");
236 fprintf(out, "%s", "%int=input\n");
237 fprintf(out, "%s", "%d2e=input\n");
238 /* if(step)
239 * fprintf(out,"%s","%nosave\n");
241 fprintf(out, "%s", "%chk=input\n");
242 fprintf(out, "%s%d\n", "%mem=", qm->QMmem);
243 fprintf(out, "%s%3d\n", "%nprocshare=", qm->nQMcpus);
245 /* use the versions of
246 * l301 that computes the interaction between MM and QM atoms.
247 * l510 that can punch the CI coefficients
248 * l701 that can do gradients on MM atoms
251 /* local version */
252 fprintf(out, "%s%s%s",
253 "%subst l510 ",
254 qm->devel_dir,
255 "/l510\n");
256 fprintf(out, "%s%s%s",
257 "%subst l301 ",
258 qm->devel_dir,
259 "/l301\n");
260 fprintf(out, "%s%s%s",
261 "%subst l701 ",
262 qm->devel_dir,
263 "/l701\n");
265 fprintf(out, "%s%s%s",
266 "%subst l1003 ",
267 qm->devel_dir,
268 "/l1003\n");
269 fprintf(out, "%s%s%s",
270 "%subst l9999 ",
271 qm->devel_dir,
272 "/l9999\n");
273 /* print the nonstandard route
275 fprintf(out, "%s",
276 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
277 fprintf(out, "%s",
278 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
279 if (mm->nrMMatoms)
281 fprintf(out,
282 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
283 qm->SHbasis[0],
284 qm->SHbasis[1],
285 qm->SHbasis[2]); /*basisset stuff */
287 else
289 fprintf(out,
290 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
291 qm->SHbasis[0],
292 qm->SHbasis[1],
293 qm->SHbasis[2]); /*basisset stuff */
295 /* development */
296 if (step+1) /* fetch initial guess from check point file */
297 { /* hack, to alyays read from chk file!!!!! */
298 fprintf(out, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
299 qm->CASelectrons,
300 "18=", qm->CASorbitals, "/1,5;\n");
302 else /* generate the first checkpoint file */
304 fprintf(out, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
305 qm->CASelectrons,
306 "18=", qm->CASorbitals, "/1,5;\n");
308 /* the rest of the input depends on where the system is on the PES
310 if (swap && bSA) /* make a slide to the other surface */
312 if (qm->CASorbitals > 6) /* use direct and no full diag */
314 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
316 else
318 if (qm->cpmcscf)
320 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
321 qm->accuracy);
322 if (mm->nrMMatoms > 0)
324 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
326 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
327 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
328 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
330 else
332 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
333 qm->accuracy);
334 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
338 else if (bSA) /* do a "state-averaged" CAS calculation */
340 if (qm->CASorbitals > 6) /* no full diag */
342 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
344 else
346 if (qm->cpmcscf)
348 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
349 qm->accuracy);
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/3;\n");
356 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
358 else
360 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
361 qm->accuracy);
362 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
366 else if (swap) /* do a "swapped" CAS calculation */
368 if (qm->CASorbitals > 6)
370 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
372 else
374 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
375 qm->accuracy);
377 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
379 else /* do a "normal" CAS calculation */
381 if (qm->CASorbitals > 6)
383 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
385 else
387 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
388 qm->accuracy);
390 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
392 fprintf(out, "\ninput-file generated by gromacs\n\n");
393 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
394 for (i = 0; i < qm->nrQMatoms; i++)
396 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
397 qm->atomicnumberQM[i],
398 qm->xQM[i][XX]/BOHR2NM,
399 qm->xQM[i][YY]/BOHR2NM,
400 qm->xQM[i][ZZ]/BOHR2NM);
402 /* MM point charge data */
403 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
405 fprintf(out, "\n");
406 for (i = 0; i < mm->nrMMatoms; i++)
408 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
409 mm->xMM[i][XX]/BOHR2NM,
410 mm->xMM[i][YY]/BOHR2NM,
411 mm->xMM[i][ZZ]/BOHR2NM,
412 mm->MMcharges[i]);
415 if (bSA) /* put the SA coefficients at the end of the file */
417 fprintf(out, "\n%10.8f %10.8f\n",
418 qm->SAstep*0.5/qm->SAsteps,
419 1-qm->SAstep*0.5/qm->SAsteps);
420 fprintf(stderr, "State Averaging level = %d/%d\n", qm->SAstep, qm->SAsteps);
422 fprintf(out, "\n");
423 fclose(out);
424 } /* write_gaussian_SH_input */
426 void write_gaussian_input(int step, t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
430 t_QMMMrec
431 *QMMMrec;
432 FILE
433 *out;
435 QMMMrec = fr->qr;
436 out = fopen("input.com", "w");
437 /* write the route */
439 if (qm->QMmethod >= eQMmethodRHF)
441 fprintf(out, "%s",
442 "%chk=input\n");
444 else
446 fprintf(out, "%s",
447 "%chk=se\n");
449 if (qm->nQMcpus > 1)
451 fprintf(out, "%s%3d\n",
452 "%nprocshare=", qm->nQMcpus);
454 fprintf(out, "%s%d\n",
455 "%mem=", qm->QMmem);
456 fprintf(out, "%s%s%s",
457 "%subst l701 ", qm->devel_dir, "/l701\n");
458 fprintf(out, "%s%s%s",
459 "%subst l301 ", qm->devel_dir, "/l301\n");
460 fprintf(out, "%s%s%s",
461 "%subst l9999 ", qm->devel_dir, "/l9999\n");
462 if (step)
464 fprintf(out, "%s",
465 "#T ");
467 else
469 fprintf(out, "%s",
470 "#P ");
472 if (qm->QMmethod == eQMmethodB3LYPLAN)
474 fprintf(out, " %s",
475 "B3LYP/GEN Pseudo=Read");
477 else
479 fprintf(out, " %s",
480 eQMmethod_names[qm->QMmethod]);
482 if (qm->QMmethod >= eQMmethodRHF)
484 if (qm->QMmethod == eQMmethodCASSCF)
486 /* in case of cas, how many electrons and orbitals do we need?
488 fprintf(out, "(%d,%d)",
489 qm->CASelectrons, qm->CASorbitals);
491 fprintf(out, "/%s",
492 eQMbasis_names[qm->QMbasis]);
495 if (QMMMrec->QMMMscheme == eQMMMschemenormal && mm->nrMMatoms)
497 fprintf(out, " %s",
498 "Charge ");
500 if (step || qm->QMmethod == eQMmethodCASSCF)
502 /* fetch guess from checkpoint file, always for CASSCF */
503 fprintf(out, "%s", " guess=read");
505 fprintf(out, "\nNosymm units=bohr\n");
507 fprintf(out, "FORCE Punch=(Derivatives) ");
508 fprintf(out, "iop(3/33=1)\n\n");
509 fprintf(out, "input-file generated by gromacs\n\n");
510 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
511 for (i = 0; i < qm->nrQMatoms; i++)
513 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
514 qm->atomicnumberQM[i],
515 qm->xQM[i][XX]/BOHR2NM,
516 qm->xQM[i][YY]/BOHR2NM,
517 qm->xQM[i][ZZ]/BOHR2NM);
520 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
521 if (qm->QMmethod == eQMmethodB3LYPLAN)
523 fprintf(out, "\n");
524 for (i = 0; i < qm->nrQMatoms; i++)
526 if (qm->atomicnumberQM[i] < 21)
528 fprintf(out, "%d ", i+1);
531 fprintf(out, "\n%s\n****\n", eQMbasis_names[qm->QMbasis]);
533 for (i = 0; i < qm->nrQMatoms; i++)
535 if (qm->atomicnumberQM[i] > 21)
537 fprintf(out, "%d ", i+1);
540 fprintf(out, "\n%s\n****\n\n", "lanl2dz");
542 for (i = 0; i < qm->nrQMatoms; i++)
544 if (qm->atomicnumberQM[i] > 21)
546 fprintf(out, "%d ", i+1);
549 fprintf(out, "\n%s\n", "lanl2dz");
554 /* MM point charge data */
555 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
557 fprintf(stderr, "nr mm atoms in gaussian.c = %d\n", mm->nrMMatoms);
558 fprintf(out, "\n");
559 for (i = 0; i < mm->nrMMatoms; i++)
561 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
562 mm->xMM[i][XX]/BOHR2NM,
563 mm->xMM[i][YY]/BOHR2NM,
564 mm->xMM[i][ZZ]/BOHR2NM,
565 mm->MMcharges[i]);
568 fprintf(out, "\n");
571 fclose(out);
573 } /* write_gaussian_input */
575 real read_gaussian_output(rvec QMgrad[], rvec MMgrad[], t_QMrec *qm, t_MMrec *mm)
578 i, j, atnum;
579 char
580 buf[300];
581 real
582 QMener;
583 FILE
584 *in;
586 in = fopen("fort.7", "r");
588 /* (There was additional content in the file in case
589 * of QM optimizations / transition state search,
590 * which was removed.
592 /* the next line is the energy and in the case of CAS, the energy
593 * difference between the two states.
595 if (NULL == fgets(buf, 300, in))
597 gmx_fatal(FARGS, "Error reading Gaussian output");
600 #if GMX_DOUBLE
601 sscanf(buf, "%lf\n", &QMener);
602 #else
603 sscanf(buf, "%f\n", &QMener);
604 #endif
605 /* next lines contain the gradients of the QM atoms */
606 for (i = 0; i < qm->nrQMatoms; i++)
608 if (NULL == fgets(buf, 300, in))
610 gmx_fatal(FARGS, "Error reading Gaussian output");
612 #if GMX_DOUBLE
613 sscanf(buf, "%lf %lf %lf\n",
614 &QMgrad[i][XX],
615 &QMgrad[i][YY],
616 &QMgrad[i][ZZ]);
617 #else
618 sscanf(buf, "%f %f %f\n",
619 &QMgrad[i][XX],
620 &QMgrad[i][YY],
621 &QMgrad[i][ZZ]);
622 #endif
624 /* the next lines are the gradients of the MM atoms */
625 if (qm->QMmethod >= eQMmethodRHF)
627 for (i = 0; i < mm->nrMMatoms; i++)
629 if (NULL == fgets(buf, 300, in))
631 gmx_fatal(FARGS, "Error reading Gaussian output");
633 #if GMX_DOUBLE
634 sscanf(buf, "%lf %lf %lf\n",
635 &MMgrad[i][XX],
636 &MMgrad[i][YY],
637 &MMgrad[i][ZZ]);
638 #else
639 sscanf(buf, "%f %f %f\n",
640 &MMgrad[i][XX],
641 &MMgrad[i][YY],
642 &MMgrad[i][ZZ]);
643 #endif
646 fclose(in);
647 return(QMener);
650 real read_gaussian_SH_output(rvec QMgrad[], rvec MMgrad[], int step, t_QMrec *qm, t_MMrec *mm)
654 char
655 buf[300];
656 real
657 QMener, DeltaE;
658 FILE
659 *in;
661 in = fopen("fort.7", "r");
662 /* first line is the energy and in the case of CAS, the energy
663 * difference between the two states.
665 if (NULL == fgets(buf, 300, in))
667 gmx_fatal(FARGS, "Error reading Gaussian output");
670 #if GMX_DOUBLE
671 sscanf(buf, "%lf %lf\n", &QMener, &DeltaE);
672 #else
673 sscanf(buf, "%f %f\n", &QMener, &DeltaE);
674 #endif
676 /* switch on/off the State Averaging */
678 if (DeltaE > qm->SAoff)
680 if (qm->SAstep > 0)
682 qm->SAstep--;
685 else if (DeltaE < qm->SAon || (qm->SAstep > 0))
687 if (qm->SAstep < qm->SAsteps)
689 qm->SAstep++;
693 /* for debugging: */
694 fprintf(stderr, "Gap = %5f,SA = %3d\n", DeltaE, (qm->SAstep > 0));
695 /* next lines contain the gradients of the QM atoms */
696 for (i = 0; i < qm->nrQMatoms; i++)
698 if (NULL == fgets(buf, 300, in))
700 gmx_fatal(FARGS, "Error reading Gaussian output");
703 #if GMX_DOUBLE
704 sscanf(buf, "%lf %lf %lf\n",
705 &QMgrad[i][XX],
706 &QMgrad[i][YY],
707 &QMgrad[i][ZZ]);
708 #else
709 sscanf(buf, "%f %f %f\n",
710 &QMgrad[i][XX],
711 &QMgrad[i][YY],
712 &QMgrad[i][ZZ]);
713 #endif
715 /* the next lines, are the gradients of the MM atoms */
717 for (i = 0; i < mm->nrMMatoms; i++)
719 if (NULL == fgets(buf, 300, in))
721 gmx_fatal(FARGS, "Error reading Gaussian output");
723 #if GMX_DOUBLE
724 sscanf(buf, "%lf %lf %lf\n",
725 &MMgrad[i][XX],
726 &MMgrad[i][YY],
727 &MMgrad[i][ZZ]);
728 #else
729 sscanf(buf, "%f %f %f\n",
730 &MMgrad[i][XX],
731 &MMgrad[i][YY],
732 &MMgrad[i][ZZ]);
733 #endif
736 /* the next line contains the two CI eigenvector elements */
737 if (NULL == fgets(buf, 300, in))
739 gmx_fatal(FARGS, "Error reading Gaussian output");
741 if (!step)
743 sscanf(buf, "%d", &qm->CIdim);
744 snew(qm->CIvec1, qm->CIdim);
745 snew(qm->CIvec1old, qm->CIdim);
746 snew(qm->CIvec2, qm->CIdim);
747 snew(qm->CIvec2old, qm->CIdim);
749 else
751 /* before reading in the new current CI vectors, copy the current
752 * CI vector into the old one.
754 for (i = 0; i < qm->CIdim; i++)
756 qm->CIvec1old[i] = qm->CIvec1[i];
757 qm->CIvec2old[i] = qm->CIvec2[i];
760 /* first vector */
761 for (i = 0; i < qm->CIdim; i++)
763 if (NULL == fgets(buf, 300, in))
765 gmx_fatal(FARGS, "Error reading Gaussian output");
767 #if GMX_DOUBLE
768 sscanf(buf, "%lf\n", &qm->CIvec1[i]);
769 #else
770 sscanf(buf, "%f\n", &qm->CIvec1[i]);
771 #endif
773 /* second vector */
774 for (i = 0; i < qm->CIdim; i++)
776 if (NULL == fgets(buf, 300, in))
778 gmx_fatal(FARGS, "Error reading Gaussian output");
780 #if GMX_DOUBLE
781 sscanf(buf, "%lf\n", &qm->CIvec2[i]);
782 #else
783 sscanf(buf, "%f\n", &qm->CIvec2[i]);
784 #endif
786 fclose(in);
787 return(QMener);
790 real inproduct(real *a, real *b, int n)
794 real
795 dot = 0.0;
797 /* computes the inner product between two vectors (a.b), both of
798 * which have length n.
800 for (i = 0; i < n; i++)
802 dot += a[i]*b[i];
804 return(dot);
807 int hop(int step, t_QMrec *qm)
810 swap = 0;
811 real
812 d11 = 0.0, d12 = 0.0, d21 = 0.0, d22 = 0.0;
814 /* calculates the inproduct between the current Ci vector and the
815 * previous CI vector. A diabatic hop will be made if d12 and d21
816 * are much bigger than d11 and d22. In that case hop returns true,
817 * otherwise it returns false.
819 if (step) /* only go on if more than one step has been done */
821 d11 = inproduct(qm->CIvec1, qm->CIvec1old, qm->CIdim);
822 d12 = inproduct(qm->CIvec1, qm->CIvec2old, qm->CIdim);
823 d21 = inproduct(qm->CIvec2, qm->CIvec1old, qm->CIdim);
824 d22 = inproduct(qm->CIvec2, qm->CIvec2old, qm->CIdim);
826 fprintf(stderr, "-------------------\n");
827 fprintf(stderr, "d11 = %13.8f\n", d11);
828 fprintf(stderr, "d12 = %13.8f\n", d12);
829 fprintf(stderr, "d21 = %13.8f\n", d21);
830 fprintf(stderr, "d22 = %13.8f\n", d22);
831 fprintf(stderr, "-------------------\n");
833 if ((fabs(d12) > 0.5) && (fabs(d21) > 0.5))
835 swap = 1;
838 return(swap);
841 void do_gaussian(int step, char *exe)
843 char
844 buf[STRLEN];
846 /* make the call to the gaussian binary through system()
847 * The location of the binary will be picked up from the
848 * environment using getenv().
850 if (step) /* hack to prevent long inputfiles */
852 sprintf(buf, "%s < %s > %s",
853 exe,
854 "input.com",
855 "input.log");
857 else
859 sprintf(buf, "%s < %s > %s",
860 exe,
861 "input.com",
862 "input.log");
864 fprintf(stderr, "Calling '%s'\n", buf);
865 if (system(buf) != 0)
867 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
871 real call_gaussian(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
873 /* normal gaussian jobs */
874 static int
875 step = 0;
877 i, j;
878 real
879 QMener = 0.0;
880 rvec
881 *QMgrad, *MMgrad;
882 char
883 *exe;
885 snew(exe, 30);
886 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
887 snew(QMgrad, qm->nrQMatoms);
888 snew(MMgrad, mm->nrMMatoms);
890 write_gaussian_input(step, fr, qm, mm);
891 do_gaussian(step, exe);
892 QMener = read_gaussian_output(QMgrad, MMgrad, qm, mm);
893 /* put the QMMM forces in the force array and to the fshift
895 for (i = 0; i < qm->nrQMatoms; i++)
897 for (j = 0; j < DIM; j++)
899 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
900 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
903 for (i = 0; i < mm->nrMMatoms; i++)
905 for (j = 0; j < DIM; j++)
907 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
908 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
911 QMener = QMener*HARTREE2KJ*AVOGADRO;
912 step++;
913 free(exe);
914 return(QMener);
916 } /* call_gaussian */
918 real call_gaussian_SH(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
920 /* a gaussian call routine intended for doing diabatic surface
921 * "sliding". See the manual for the theoretical background of this
922 * TSH method.
924 static int
925 step = 0;
927 state, i, j;
928 real
929 QMener = 0.0;
930 static gmx_bool
931 swapped = FALSE; /* handle for identifying the current PES */
932 gmx_bool
933 swap = FALSE; /* the actual swap */
934 rvec
935 *QMgrad, *MMgrad;
936 char
937 *buf;
938 char
939 *exe;
941 snew(exe, 30);
942 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
943 /* hack to do ground state simulations */
944 if (!step)
946 snew(buf, 20);
947 buf = getenv("GMX_QM_GROUND_STATE");
948 if (buf)
950 sscanf(buf, "%d", &state);
952 else
954 state = 2;
956 if (state == 1)
958 swapped = TRUE;
961 /* end of hack */
964 /* copy the QMMMrec pointer */
965 snew(QMgrad, qm->nrQMatoms);
966 snew(MMgrad, mm->nrMMatoms);
967 /* at step 0 there should be no SA */
968 /* if(!step)
969 * qr->bSA=FALSE;*/
970 /* temporray set to step + 1, since there is a chk start */
971 write_gaussian_SH_input(step, swapped, fr, qm, mm);
973 do_gaussian(step, exe);
974 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, qm, mm);
976 /* check for a surface hop. Only possible if we were already state
977 * averaging.
979 if (qm->SAstep > 0)
981 if (!swapped)
983 swap = (step && hop(step, qm));
984 swapped = swap;
986 else /* already on the other surface, so check if we go back */
988 swap = (step && hop(step, qm));
989 swapped = !swap; /* so swapped shoud be false again */
991 if (swap) /* change surface, so do another call */
993 write_gaussian_SH_input(step, swapped, fr, qm, mm);
994 do_gaussian(step, exe);
995 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, qm, mm);
998 /* add the QMMM forces to the gmx force array and fshift
1000 for (i = 0; i < qm->nrQMatoms; i++)
1002 for (j = 0; j < DIM; j++)
1004 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1005 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1008 for (i = 0; i < mm->nrMMatoms; i++)
1010 for (j = 0; j < DIM; j++)
1012 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1013 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1016 QMener = QMener*HARTREE2KJ*AVOGADRO;
1017 fprintf(stderr, "step %5d, SA = %5d, swap = %5d\n",
1018 step, (qm->SAstep > 0), swapped);
1019 step++;
1020 free(exe);
1021 return(QMener);
1023 } /* call_gaussian_SH */
1025 /* end of gaussian sub routines */
1027 #else
1029 gmx_qmmm_gaussian_empty;
1030 #endif