Fix incorrect rvdw on GPU with rvdw<rcoulomb
[gromacs.git] / src / gromacs / mdlib / qm_gaussian.cpp
blobd897e3b6d409a6dc43ea33ca0bab08f977cdd3ad
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,2018, 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 "qm_gaussian.h"
41 #include "config.h"
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/forcerec.h"
56 #include "gromacs/mdlib/ns.h"
57 #include "gromacs/mdlib/qmmm.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/smalloc.h"
63 // When not built in a configuration with QMMM support, much of this
64 // code is unreachable by design. Tell clang not to warn about it.
65 #pragma GCC diagnostic push
66 #pragma GCC diagnostic ignored "-Wmissing-noreturn"
68 /* TODO: this should be made thread-safe */
70 /* Gaussian interface routines */
72 void init_gaussian(t_QMrec *qm)
74 ivec
75 basissets[eQMbasisNR] = {{0, 3, 0},
76 {0, 3, 0}, /*added for double sto-3g entry in names.c*/
77 {5, 0, 0},
78 {5, 0, 1},
79 {5, 0, 11},
80 {5, 6, 0},
81 {1, 6, 0},
82 {1, 6, 1},
83 {1, 6, 11},
84 {4, 6, 0}};
85 char
86 *buf = nullptr;
87 int
90 if (!GMX_QMMM_GAUSSIAN)
92 gmx_fatal(FARGS, "Cannot call GAUSSIAN unless linked against it. Use cmake -DGMX_QMMM_PROGRAM=GAUSSIAN, and ensure that linking will work correctly.");
95 /* using the ivec above to convert the basis read form the mdp file
96 * in a human readable format into some numbers for the gaussian
97 * route. This is necessary as we are using non standard routes to
98 * do SH.
101 /* per layer we make a new subdir for integral file, checkpoint
102 * files and such. These dirs are stored in the QMrec for
103 * convenience
107 if (!qm->nQMcpus) /* this we do only once per layer
108 * as we call g01 externally
112 for (i = 0; i < DIM; i++)
114 qm->SHbasis[i] = basissets[qm->QMbasis][i];
117 /* init gradually switching on of the SA */
118 qm->SAstep = 0;
119 /* we read the number of cpus and environment from the environment
120 * if set.
122 buf = getenv("GMX_QM_GAUSSIAN_NCPUS");
123 if (buf)
125 sscanf(buf, "%d", &qm->nQMcpus);
127 else
129 qm->nQMcpus = 1;
131 fprintf(stderr, "number of CPUs for gaussian = %d\n", qm->nQMcpus);
132 buf = getenv("GMX_QM_GAUSSIAN_MEMORY");
133 if (buf)
135 sscanf(buf, "%d", &qm->QMmem);
137 else
139 qm->QMmem = 50000000;
141 fprintf(stderr, "memory for gaussian = %d\n", qm->QMmem);
142 buf = getenv("GMX_QM_ACCURACY");
143 if (buf)
145 sscanf(buf, "%d", &qm->accuracy);
147 else
149 qm->accuracy = 8;
151 fprintf(stderr, "accuracy in l510 = %d\n", qm->accuracy);
153 buf = getenv("GMX_QM_CPMCSCF");
154 if (buf)
156 sscanf(buf, "%d", &i);
157 qm->cpmcscf = (i != 0);
159 else
161 qm->cpmcscf = FALSE;
163 if (qm->cpmcscf)
165 fprintf(stderr, "using cp-mcscf in l1003\n");
167 else
169 fprintf(stderr, "NOT using cp-mcscf in l1003\n");
171 buf = getenv("GMX_QM_SA_STEP");
172 if (buf)
174 sscanf(buf, "%d", &qm->SAstep);
176 else
178 /* init gradually switching on of the SA */
179 qm->SAstep = 0;
181 /* we read the number of cpus and environment from the environment
182 * if set.
184 fprintf(stderr, "Level of SA at start = %d\n", qm->SAstep);
185 /* gaussian settings on the system */
186 buf = getenv("GMX_QM_GAUSS_DIR");
188 if (buf)
190 qm->gauss_dir = gmx_strdup(buf);
192 else
194 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_DIR, check gaussian manual\n");
197 buf = getenv("GMX_QM_GAUSS_EXE");
198 if (buf)
200 qm->gauss_exe = gmx_strdup(buf);
202 else
204 gmx_fatal(FARGS, "no $GMX_QM_GAUSS_EXE set, check gaussian manual\n");
206 buf = getenv("GMX_QM_MODIFIED_LINKS_DIR");
207 if (buf)
209 qm->devel_dir = gmx_strdup (buf);
211 else
213 gmx_fatal(FARGS, "no $GMX_QM_MODIFIED_LINKS_DIR, this is were the modified links reside.\n");
216 /* if(fr->bRF){*/
217 /* reactionfield, file is needed using gaussian */
218 /* rffile=fopen("rf.dat","w");*/
219 /* fprintf(rffile,"%f %f\n",fr->epsilon_r,fr->rcoulomb/BOHR2NM);*/
220 /* fclose(rffile);*/
221 /* }*/
223 fprintf(stderr, "gaussian initialised...\n");
227 static void write_gaussian_SH_input(int step, gmx_bool swap,
228 const t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
232 gmx_bool
233 bSA;
234 FILE
235 *out;
236 t_QMMMrec
237 *QMMMrec;
238 QMMMrec = fr->qr;
239 bSA = (qm->SAstep > 0);
241 out = fopen("input.com", "w");
242 /* write the route */
243 fprintf(out, "%s", "%scr=input\n");
244 fprintf(out, "%s", "%rwf=input\n");
245 fprintf(out, "%s", "%int=input\n");
246 fprintf(out, "%s", "%d2e=input\n");
247 /* if(step)
248 * fprintf(out,"%s","%nosave\n");
250 fprintf(out, "%s", "%chk=input\n");
251 fprintf(out, "%s%d\n", "%mem=", qm->QMmem);
252 fprintf(out, "%s%3d\n", "%nprocshare=", qm->nQMcpus);
254 /* use the versions of
255 * l301 that computes the interaction between MM and QM atoms.
256 * l510 that can punch the CI coefficients
257 * l701 that can do gradients on MM atoms
260 /* local version */
261 fprintf(out, "%s%s%s",
262 "%subst l510 ",
263 qm->devel_dir,
264 "/l510\n");
265 fprintf(out, "%s%s%s",
266 "%subst l301 ",
267 qm->devel_dir,
268 "/l301\n");
269 fprintf(out, "%s%s%s",
270 "%subst l701 ",
271 qm->devel_dir,
272 "/l701\n");
274 fprintf(out, "%s%s%s",
275 "%subst l1003 ",
276 qm->devel_dir,
277 "/l1003\n");
278 fprintf(out, "%s%s%s",
279 "%subst l9999 ",
280 qm->devel_dir,
281 "/l9999\n");
282 /* print the nonstandard route
284 fprintf(out, "%s",
285 "#P nonstd\n 1/18=10,20=1,38=1/1;\n");
286 fprintf(out, "%s",
287 " 2/9=110,15=1,17=6,18=5,40=1/2;\n");
288 if (mm->nrMMatoms)
290 fprintf(out,
291 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=1,94=-2/1,2,3;\n",
292 qm->SHbasis[0],
293 qm->SHbasis[1],
294 qm->SHbasis[2]); /*basisset stuff */
296 else
298 fprintf(out,
299 " 3/5=%d,6=%d,7=%d,25=1,32=1,43=0,94=-2/1,2,3;\n",
300 qm->SHbasis[0],
301 qm->SHbasis[1],
302 qm->SHbasis[2]); /*basisset stuff */
304 /* development */
305 if (step+1) /* fetch initial guess from check point file */
306 { /* hack, to alyays read from chk file!!!!! */
307 fprintf(out, "%s%d,%s%d%s", " 4/5=1,7=6,17=",
308 qm->CASelectrons,
309 "18=", qm->CASorbitals, "/1,5;\n");
311 else /* generate the first checkpoint file */
313 fprintf(out, "%s%d,%s%d%s", " 4/5=0,7=6,17=",
314 qm->CASelectrons,
315 "18=", qm->CASorbitals, "/1,5;\n");
317 /* the rest of the input depends on where the system is on the PES
319 if (swap && bSA) /* make a slide to the other surface */
321 if (qm->CASorbitals > 6) /* use direct and no full diag */
323 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6,97=100/10;\n");
325 else
327 if (qm->cpmcscf)
329 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6,97=100/10;\n",
330 qm->accuracy);
331 if (mm->nrMMatoms > 0)
333 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
335 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
336 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1,97=100/3;\n");
337 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
339 else
341 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6,97=100/10;\n",
342 qm->accuracy);
343 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
347 else if (bSA) /* do a "state-averaged" CAS calculation */
349 if (qm->CASorbitals > 6) /* no full diag */
351 fprintf(out, " 5/5=2,16=-2,17=10000000,28=2,32=2,38=6/10;\n");
353 else
355 if (qm->cpmcscf)
357 fprintf(out, " 5/5=2,6=%d,17=31000200,28=2,32=2,38=6/10;\n",
358 qm->accuracy);
359 if (mm->nrMMatoms > 0)
361 fprintf(out, " 7/7=1,16=-2,30=1/1;\n");
363 fprintf(out, " 11/31=1,42=1,45=1/1;\n");
364 fprintf(out, " 10/6=1,10=700006,28=2,29=1,31=1/3;\n");
365 fprintf(out, " 7/30=1/16;\n 99/10=4/99;\n");
367 else
369 fprintf(out, " 5/5=2,6=%d,17=11000000,28=2,32=2,38=6/10;\n",
370 qm->accuracy);
371 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
375 else if (swap) /* do a "swapped" CAS calculation */
377 if (qm->CASorbitals > 6)
379 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6,97=100/10;\n");
381 else
383 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6,97=100/10;\n",
384 qm->accuracy);
386 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
388 else /* do a "normal" CAS calculation */
390 if (qm->CASorbitals > 6)
392 fprintf(out, " 5/5=2,16=-2,17=0,28=2,32=2,38=6/10;\n");
394 else
396 fprintf(out, " 5/5=2,6=%d,17=1000000,28=2,32=2,38=6/10;\n",
397 qm->accuracy);
399 fprintf(out, " 7/7=1,16=-2,30=1/1,2,3,16;\n 99/10=4/99;\n");
401 fprintf(out, "\ninput-file generated by gromacs\n\n");
402 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
403 for (i = 0; i < qm->nrQMatoms; i++)
405 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
406 qm->atomicnumberQM[i],
407 qm->xQM[i][XX]/BOHR2NM,
408 qm->xQM[i][YY]/BOHR2NM,
409 qm->xQM[i][ZZ]/BOHR2NM);
411 /* MM point charge data */
412 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
414 fprintf(out, "\n");
415 for (i = 0; i < mm->nrMMatoms; i++)
417 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
418 mm->xMM[i][XX]/BOHR2NM,
419 mm->xMM[i][YY]/BOHR2NM,
420 mm->xMM[i][ZZ]/BOHR2NM,
421 mm->MMcharges[i]);
424 if (bSA) /* put the SA coefficients at the end of the file */
426 fprintf(out, "\n%10.8f %10.8f\n",
427 qm->SAstep*0.5/qm->SAsteps,
428 1-qm->SAstep*0.5/qm->SAsteps);
429 fprintf(stderr, "State Averaging level = %d/%d\n", qm->SAstep, qm->SAsteps);
431 fprintf(out, "\n");
432 fclose(out);
433 } /* write_gaussian_SH_input */
435 static void write_gaussian_input(int step, const t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
439 t_QMMMrec
440 *QMMMrec;
441 FILE
442 *out;
444 QMMMrec = fr->qr;
445 out = fopen("input.com", "w");
446 /* write the route */
448 if (qm->QMmethod >= eQMmethodRHF)
450 fprintf(out, "%s",
451 "%chk=input\n");
453 else
455 fprintf(out, "%s",
456 "%chk=se\n");
458 if (qm->nQMcpus > 1)
460 fprintf(out, "%s%3d\n",
461 "%nprocshare=", qm->nQMcpus);
463 fprintf(out, "%s%d\n",
464 "%mem=", qm->QMmem);
465 fprintf(out, "%s%s%s",
466 "%subst l701 ", qm->devel_dir, "/l701\n");
467 fprintf(out, "%s%s%s",
468 "%subst l301 ", qm->devel_dir, "/l301\n");
469 fprintf(out, "%s%s%s",
470 "%subst l9999 ", qm->devel_dir, "/l9999\n");
471 if (step)
473 fprintf(out, "%s",
474 "#T ");
476 else
478 fprintf(out, "%s",
479 "#P ");
481 if (qm->QMmethod == eQMmethodB3LYPLAN)
483 fprintf(out, " %s",
484 "B3LYP/GEN Pseudo=Read");
486 else
488 fprintf(out, " %s",
489 eQMmethod_names[qm->QMmethod]);
491 if (qm->QMmethod >= eQMmethodRHF)
493 if (qm->QMmethod == eQMmethodCASSCF)
495 /* in case of cas, how many electrons and orbitals do we need?
497 fprintf(out, "(%d,%d)",
498 qm->CASelectrons, qm->CASorbitals);
500 fprintf(out, "/%s",
501 eQMbasis_names[qm->QMbasis]);
504 if (QMMMrec->QMMMscheme == eQMMMschemenormal && mm->nrMMatoms)
506 fprintf(out, " %s",
507 "Charge ");
509 if (step || qm->QMmethod == eQMmethodCASSCF)
511 /* fetch guess from checkpoint file, always for CASSCF */
512 fprintf(out, "%s", " guess=read");
514 fprintf(out, "\nNosymm units=bohr\n");
516 fprintf(out, "FORCE Punch=(Derivatives) ");
517 fprintf(out, "iop(3/33=1)\n\n");
518 fprintf(out, "input-file generated by gromacs\n\n");
519 fprintf(out, "%2d%2d\n", qm->QMcharge, qm->multiplicity);
520 for (i = 0; i < qm->nrQMatoms; i++)
522 fprintf(out, "%3d %10.7f %10.7f %10.7f\n",
523 qm->atomicnumberQM[i],
524 qm->xQM[i][XX]/BOHR2NM,
525 qm->xQM[i][YY]/BOHR2NM,
526 qm->xQM[i][ZZ]/BOHR2NM);
529 /* Pseudo Potential and ECP are included here if selected (MEthod suffix LAN) */
530 if (qm->QMmethod == eQMmethodB3LYPLAN)
532 fprintf(out, "\n");
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", eQMbasis_names[qm->QMbasis]);
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****\n\n", "lanl2dz");
551 for (i = 0; i < qm->nrQMatoms; i++)
553 if (qm->atomicnumberQM[i] > 21)
555 fprintf(out, "%d ", i+1);
558 fprintf(out, "\n%s\n", "lanl2dz");
563 /* MM point charge data */
564 if (QMMMrec->QMMMscheme != eQMMMschemeoniom && mm->nrMMatoms)
566 fprintf(stderr, "nr mm atoms in gaussian.c = %d\n", mm->nrMMatoms);
567 fprintf(out, "\n");
568 for (i = 0; i < mm->nrMMatoms; i++)
570 fprintf(out, "%10.7f %10.7f %10.7f %8.4f\n",
571 mm->xMM[i][XX]/BOHR2NM,
572 mm->xMM[i][YY]/BOHR2NM,
573 mm->xMM[i][ZZ]/BOHR2NM,
574 mm->MMcharges[i]);
577 fprintf(out, "\n");
580 fclose(out);
582 } /* write_gaussian_input */
584 static real read_gaussian_output(rvec QMgrad[], rvec MMgrad[], t_QMrec *qm, t_MMrec *mm)
588 char
589 buf[300];
590 real
591 QMener;
592 FILE
593 *in;
595 in = fopen("fort.7", "r");
597 /* (There was additional content in the file in case
598 * of QM optimizations / transition state search,
599 * which was removed.
601 /* the next line is the energy and in the case of CAS, the energy
602 * difference between the two states.
604 if (nullptr == fgets(buf, 300, in))
606 gmx_fatal(FARGS, "Error reading Gaussian output");
609 #if GMX_DOUBLE
610 sscanf(buf, "%lf\n", &QMener);
611 #else
612 sscanf(buf, "%f\n", &QMener);
613 #endif
614 /* next lines contain the gradients of the QM atoms */
615 for (i = 0; i < qm->nrQMatoms; i++)
617 if (nullptr == fgets(buf, 300, in))
619 gmx_fatal(FARGS, "Error reading Gaussian output");
621 #if GMX_DOUBLE
622 sscanf(buf, "%lf %lf %lf\n",
623 &QMgrad[i][XX],
624 &QMgrad[i][YY],
625 &QMgrad[i][ZZ]);
626 #else
627 sscanf(buf, "%f %f %f\n",
628 &QMgrad[i][XX],
629 &QMgrad[i][YY],
630 &QMgrad[i][ZZ]);
631 #endif
633 /* the next lines are the gradients of the MM atoms */
634 if (qm->QMmethod >= eQMmethodRHF)
636 for (i = 0; i < mm->nrMMatoms; i++)
638 if (nullptr == fgets(buf, 300, in))
640 gmx_fatal(FARGS, "Error reading Gaussian output");
642 #if GMX_DOUBLE
643 sscanf(buf, "%lf %lf %lf\n",
644 &MMgrad[i][XX],
645 &MMgrad[i][YY],
646 &MMgrad[i][ZZ]);
647 #else
648 sscanf(buf, "%f %f %f\n",
649 &MMgrad[i][XX],
650 &MMgrad[i][YY],
651 &MMgrad[i][ZZ]);
652 #endif
655 fclose(in);
656 return(QMener);
659 static real read_gaussian_SH_output(rvec QMgrad[], rvec MMgrad[], int step, t_QMrec *qm, t_MMrec *mm)
663 char
664 buf[300];
665 real
666 QMener, DeltaE;
667 FILE
668 *in;
670 in = fopen("fort.7", "r");
671 /* first line is the energy and in the case of CAS, the energy
672 * difference between the two states.
674 if (nullptr == fgets(buf, 300, in))
676 gmx_fatal(FARGS, "Error reading Gaussian output");
679 #if GMX_DOUBLE
680 sscanf(buf, "%lf %lf\n", &QMener, &DeltaE);
681 #else
682 sscanf(buf, "%f %f\n", &QMener, &DeltaE);
683 #endif
685 /* switch on/off the State Averaging */
687 if (DeltaE > qm->SAoff)
689 if (qm->SAstep > 0)
691 qm->SAstep--;
694 else if (DeltaE < qm->SAon || (qm->SAstep > 0))
696 if (qm->SAstep < qm->SAsteps)
698 qm->SAstep++;
702 /* for debugging: */
703 fprintf(stderr, "Gap = %5f,SA = %3d\n", DeltaE, static_cast<int>(qm->SAstep > 0));
704 /* next lines contain the gradients of the QM atoms */
705 for (i = 0; i < qm->nrQMatoms; i++)
707 if (nullptr == fgets(buf, 300, in))
709 gmx_fatal(FARGS, "Error reading Gaussian output");
712 #if GMX_DOUBLE
713 sscanf(buf, "%lf %lf %lf\n",
714 &QMgrad[i][XX],
715 &QMgrad[i][YY],
716 &QMgrad[i][ZZ]);
717 #else
718 sscanf(buf, "%f %f %f\n",
719 &QMgrad[i][XX],
720 &QMgrad[i][YY],
721 &QMgrad[i][ZZ]);
722 #endif
724 /* the next lines, are the gradients of the MM atoms */
726 for (i = 0; i < mm->nrMMatoms; i++)
728 if (nullptr == fgets(buf, 300, in))
730 gmx_fatal(FARGS, "Error reading Gaussian output");
732 #if GMX_DOUBLE
733 sscanf(buf, "%lf %lf %lf\n",
734 &MMgrad[i][XX],
735 &MMgrad[i][YY],
736 &MMgrad[i][ZZ]);
737 #else
738 sscanf(buf, "%f %f %f\n",
739 &MMgrad[i][XX],
740 &MMgrad[i][YY],
741 &MMgrad[i][ZZ]);
742 #endif
745 /* the next line contains the two CI eigenvector elements */
746 if (nullptr == fgets(buf, 300, in))
748 gmx_fatal(FARGS, "Error reading Gaussian output");
750 if (!step)
752 sscanf(buf, "%d", &qm->CIdim);
753 snew(qm->CIvec1, qm->CIdim);
754 snew(qm->CIvec1old, qm->CIdim);
755 snew(qm->CIvec2, qm->CIdim);
756 snew(qm->CIvec2old, qm->CIdim);
758 else
760 /* before reading in the new current CI vectors, copy the current
761 * CI vector into the old one.
763 for (i = 0; i < qm->CIdim; i++)
765 qm->CIvec1old[i] = qm->CIvec1[i];
766 qm->CIvec2old[i] = qm->CIvec2[i];
769 /* first vector */
770 for (i = 0; i < qm->CIdim; i++)
772 if (nullptr == fgets(buf, 300, in))
774 gmx_fatal(FARGS, "Error reading Gaussian output");
776 #if GMX_DOUBLE
777 sscanf(buf, "%lf\n", &qm->CIvec1[i]);
778 #else
779 sscanf(buf, "%f\n", &qm->CIvec1[i]);
780 #endif
782 /* second vector */
783 for (i = 0; i < qm->CIdim; i++)
785 if (nullptr == fgets(buf, 300, in))
787 gmx_fatal(FARGS, "Error reading Gaussian output");
789 #if GMX_DOUBLE
790 sscanf(buf, "%lf\n", &qm->CIvec2[i]);
791 #else
792 sscanf(buf, "%f\n", &qm->CIvec2[i]);
793 #endif
795 fclose(in);
796 return(QMener);
799 static real inproduct(const real *a, const real *b, int n)
803 real
804 dot = 0.0;
806 /* computes the inner product between two vectors (a.b), both of
807 * which have length n.
809 for (i = 0; i < n; i++)
811 dot += a[i]*b[i];
813 return(dot);
816 static int hop(int step, t_QMrec *qm)
819 swap = 0;
820 real
821 d11 = 0.0, d12 = 0.0, d21 = 0.0, d22 = 0.0;
823 /* calculates the inproduct between the current Ci vector and the
824 * previous CI vector. A diabatic hop will be made if d12 and d21
825 * are much bigger than d11 and d22. In that case hop returns true,
826 * otherwise it returns false.
828 if (step) /* only go on if more than one step has been done */
830 d11 = inproduct(qm->CIvec1, qm->CIvec1old, qm->CIdim);
831 d12 = inproduct(qm->CIvec1, qm->CIvec2old, qm->CIdim);
832 d21 = inproduct(qm->CIvec2, qm->CIvec1old, qm->CIdim);
833 d22 = inproduct(qm->CIvec2, qm->CIvec2old, qm->CIdim);
835 fprintf(stderr, "-------------------\n");
836 fprintf(stderr, "d11 = %13.8f\n", d11);
837 fprintf(stderr, "d12 = %13.8f\n", d12);
838 fprintf(stderr, "d21 = %13.8f\n", d21);
839 fprintf(stderr, "d22 = %13.8f\n", d22);
840 fprintf(stderr, "-------------------\n");
842 if ((fabs(d12) > 0.5) && (fabs(d21) > 0.5))
844 swap = 1;
847 return(swap);
850 static void do_gaussian(int step, char *exe)
852 char
853 buf[STRLEN];
855 /* make the call to the gaussian binary through system()
856 * The location of the binary will be picked up from the
857 * environment using getenv().
859 if (step) /* hack to prevent long inputfiles */
861 sprintf(buf, "%s < %s > %s",
862 exe,
863 "input.com",
864 "input.log");
866 else
868 sprintf(buf, "%s < %s > %s",
869 exe,
870 "input.com",
871 "input.log");
873 fprintf(stderr, "Calling '%s'\n", buf);
874 if (system(buf) != 0)
876 gmx_fatal(FARGS, "Call to '%s' failed\n", buf);
880 real call_gaussian(const t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
882 /* normal gaussian jobs */
883 static int
884 step = 0;
886 i, j;
887 real
888 QMener = 0.0;
889 rvec
890 *QMgrad, *MMgrad;
891 char
892 *exe;
894 snew(exe, 30);
895 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
896 snew(QMgrad, qm->nrQMatoms);
897 snew(MMgrad, mm->nrMMatoms);
899 write_gaussian_input(step, fr, qm, mm);
900 do_gaussian(step, exe);
901 QMener = read_gaussian_output(QMgrad, MMgrad, qm, mm);
902 /* put the QMMM forces in the force array and to the fshift
904 for (i = 0; i < qm->nrQMatoms; i++)
906 for (j = 0; j < DIM; j++)
908 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
909 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
912 for (i = 0; i < mm->nrMMatoms; i++)
914 for (j = 0; j < DIM; j++)
916 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
917 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
920 QMener = QMener*HARTREE2KJ*AVOGADRO;
921 step++;
922 free(exe);
923 return(QMener);
925 } /* call_gaussian */
927 real call_gaussian_SH(const t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[])
929 /* a gaussian call routine intended for doing diabatic surface
930 * "sliding". See the manual for the theoretical background of this
931 * TSH method.
933 static int
934 step = 0;
936 state, i, j;
937 real
938 QMener = 0.0;
939 static gmx_bool
940 swapped = FALSE; /* handle for identifying the current PES */
941 gmx_bool
942 swap = FALSE; /* the actual swap */
943 rvec
944 *QMgrad, *MMgrad;
945 char
946 *buf;
947 char
948 *exe;
950 snew(exe, 30);
951 sprintf(exe, "%s/%s", qm->gauss_dir, qm->gauss_exe);
952 /* hack to do ground state simulations */
953 if (!step)
955 snew(buf, 20);
956 buf = getenv("GMX_QM_GROUND_STATE");
957 if (buf)
959 sscanf(buf, "%d", &state);
961 else
963 state = 2;
965 if (state == 1)
967 swapped = TRUE;
970 /* end of hack */
973 /* copy the QMMMrec pointer */
974 snew(QMgrad, qm->nrQMatoms);
975 snew(MMgrad, mm->nrMMatoms);
976 /* at step 0 there should be no SA */
977 /* if(!step)
978 * qr->bSA=FALSE;*/
979 /* temporray set to step + 1, since there is a chk start */
980 write_gaussian_SH_input(step, swapped, fr, qm, mm);
982 do_gaussian(step, exe);
983 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, qm, mm);
985 /* check for a surface hop. Only possible if we were already state
986 * averaging.
988 if (qm->SAstep > 0)
990 if (!swapped)
992 swap = ((step != 0) && (hop(step, qm) != 0));
993 swapped = swap;
995 else /* already on the other surface, so check if we go back */
997 swap = ((step != 0) && (hop(step, qm) != 0));
998 swapped = !swap; /* so swapped shoud be false again */
1000 if (swap) /* change surface, so do another call */
1002 write_gaussian_SH_input(step, swapped, fr, qm, mm);
1003 do_gaussian(step, exe);
1004 QMener = read_gaussian_SH_output(QMgrad, MMgrad, step, qm, mm);
1007 /* add the QMMM forces to the gmx force array and fshift
1009 for (i = 0; i < qm->nrQMatoms; i++)
1011 for (j = 0; j < DIM; j++)
1013 f[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1014 fshift[i][j] = HARTREE_BOHR2MD*QMgrad[i][j];
1017 for (i = 0; i < mm->nrMMatoms; i++)
1019 for (j = 0; j < DIM; j++)
1021 f[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1022 fshift[i+qm->nrQMatoms][j] = HARTREE_BOHR2MD*MMgrad[i][j];
1025 QMener = QMener*HARTREE2KJ*AVOGADRO;
1026 fprintf(stderr, "step %5d, SA = %5d, swap = %5d\n",
1027 step, static_cast<int>(qm->SAstep > 0), static_cast<int>(swapped));
1028 step++;
1029 free(exe);
1030 return(QMener);
1032 } /* call_gaussian_SH */
1034 #pragma GCC diagnostic pop