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, 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.
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/legacyheaders/types/commrec.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/legacyheaders/force.h"
53 #include "gromacs/fileio/confio.h"
54 #include "gromacs/legacyheaders/names.h"
55 #include "gromacs/legacyheaders/network.h"
56 #include "gromacs/legacyheaders/ns.h"
57 #include "gromacs/legacyheaders/nrnb.h"
58 #include "gromacs/legacyheaders/txtdump.h"
59 #include "gromacs/legacyheaders/qmmm.h"
60 #include "gromacs/legacyheaders/typedefs.h"
61 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/pbcutil/ishift.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/smalloc.h"
68 /* declarations of the interfaces to the QM packages. The _SH indicate
69 * the QM interfaces can be used for Surface Hopping simulations
71 #ifdef GMX_QMMM_GAMESS
72 /* GAMESS interface */
75 init_gamess(t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
78 call_gamess(t_commrec *cr, t_forcerec *fr,
79 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
81 #elif defined GMX_QMMM_MOPAC
85 init_mopac(t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
88 call_mopac(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
89 t_MMrec *mm, rvec f[], rvec fshift[]);
92 call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
93 t_MMrec *mm, rvec f[], rvec fshift[]);
95 #elif defined GMX_QMMM_GAUSSIAN
96 /* GAUSSIAN interface */
99 init_gaussian(t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
102 call_gaussian_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
103 t_MMrec *mm, rvec f[], rvec fshift[]);
106 call_gaussian(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
107 t_MMrec *mm, rvec f[], rvec fshift[]);
109 #elif defined GMX_QMMM_ORCA
113 init_orca(t_QMrec *qm);
116 call_orca(t_forcerec *fr, t_QMrec *qm,
117 t_MMrec *mm, rvec f[], rvec fshift[]);
124 /* this struct and these comparison functions are needed for creating
125 * a QMMM input for the QM routines from the QMMM neighbor list.
133 static int struct_comp(const void *a, const void *b)
136 return (int)(((t_j_particle *)a)->j)-(int)(((t_j_particle *)b)->j);
140 static int int_comp(const void *a, const void *b)
143 return (*(int *)a) - (*(int *)b);
147 static int QMlayer_comp(const void *a, const void *b)
150 return (int)(((t_QMrec *)a)->nrQMatoms)-(int)(((t_QMrec *)b)->nrQMatoms);
154 real call_QMroutine(t_commrec gmx_unused *cr, t_forcerec gmx_unused *fr, t_QMrec gmx_unused *qm,
155 t_MMrec gmx_unused *mm, rvec gmx_unused f[], rvec gmx_unused fshift[])
157 /* makes a call to the requested QM routine (qm->QMmethod)
158 * Note that f is actually the gradient, i.e. -f
163 /* do a semi-empiprical calculation */
165 if (qm->QMmethod < eQMmethodRHF && !(mm->nrMMatoms))
167 #ifdef GMX_QMMM_MOPAC
170 QMener = call_mopac_SH(cr, fr, qm, mm, f, fshift);
174 QMener = call_mopac(cr, fr, qm, mm, f, fshift);
177 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
182 /* do an ab-initio calculation */
183 if (qm->bSH && qm->QMmethod == eQMmethodCASSCF)
185 #ifdef GMX_QMMM_GAUSSIAN
186 QMener = call_gaussian_SH(cr, fr, qm, mm, f, fshift);
188 gmx_fatal(FARGS, "Ab-initio Surface-hopping only supported with Gaussian.");
193 #ifdef GMX_QMMM_GAMESS
194 QMener = call_gamess(cr, fr, qm, mm, f, fshift);
195 #elif defined GMX_QMMM_GAUSSIAN
196 QMener = call_gaussian(cr, fr, qm, mm, f, fshift);
197 #elif defined GMX_QMMM_ORCA
198 QMener = call_orca(fr, qm, mm, f, fshift);
200 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
207 void init_QMroutine(t_commrec gmx_unused *cr, t_QMrec gmx_unused *qm, t_MMrec gmx_unused *mm)
209 /* makes a call to the requested QM routine (qm->QMmethod)
211 if (qm->QMmethod < eQMmethodRHF)
213 #ifdef GMX_QMMM_MOPAC
214 /* do a semi-empiprical calculation */
215 init_mopac(cr, qm, mm);
217 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
222 /* do an ab-initio calculation */
223 #ifdef GMX_QMMM_GAMESS
224 init_gamess(cr, qm, mm);
225 #elif defined GMX_QMMM_GAUSSIAN
226 init_gaussian(cr, qm, mm);
227 #elif defined GMX_QMMM_ORCA
230 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
233 } /* init_QMroutine */
235 void update_QMMM_coord(rvec x[], t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
237 /* shifts the QM and MM particles into the central box and stores
238 * these shifted coordinates in the coordinate arrays of the
239 * QMMMrec. These coordinates are passed on the QM subroutines.
244 /* shift the QM atoms into the central box
246 for (i = 0; i < qm->nrQMatoms; i++)
248 rvec_sub(x[qm->indexQM[i]], fr->shift_vec[qm->shiftQM[i]], qm->xQM[i]);
250 /* also shift the MM atoms into the central box, if any
252 for (i = 0; i < mm->nrMMatoms; i++)
254 rvec_sub(x[mm->indexMM[i]], fr->shift_vec[mm->shiftMM[i]], mm->xMM[i]);
256 } /* update_QMMM_coord */
258 static void punch_QMMM_excl(t_QMrec *qm, t_MMrec *mm, t_blocka *excls)
260 /* punch a file containing the bonded interactions of each QM
261 * atom with MM atoms. These need to be excluded in the QM routines
262 * Only needed in case of QM/MM optimizations
267 i, j, k, nrexcl = 0, *excluded = NULL, max = 0;
270 out = fopen("QMMMexcl.dat", "w");
272 /* this can be done more efficiently I think
274 for (i = 0; i < qm->nrQMatoms; i++)
277 for (j = excls->index[qm->indexQM[i]];
278 j < excls->index[qm->indexQM[i]+1];
281 for (k = 0; k < mm->nrMMatoms; k++)
283 if (mm->indexMM[k] == excls->a[j]) /* the excluded MM atom */
288 srenew(excluded, max);
290 excluded[nrexcl++] = k;
296 fprintf(out, "%5d %5d\n", i+1, nrexcl);
297 for (j = 0; j < nrexcl; j++)
299 fprintf(out, "%5d ", excluded[j]);
305 } /* punch_QMMM_excl */
308 /* end of QMMM subroutines */
310 /* QMMM core routines */
312 t_QMrec *mk_QMrec(void)
319 t_MMrec *mk_MMrec(void)
326 static void init_QMrec(int grpnr, t_QMrec *qm, int nr, int *atomarray,
327 gmx_mtop_t *mtop, t_inputrec *ir)
329 /* fills the t_QMrec struct of QM group grpnr
332 gmx_mtop_atomlookup_t alook;
338 snew(qm->indexQM, nr);
339 snew(qm->shiftQM, nr); /* the shifts */
340 for (i = 0; i < nr; i++)
342 qm->indexQM[i] = atomarray[i];
345 alook = gmx_mtop_atomlookup_init(mtop);
347 snew(qm->atomicnumberQM, nr);
348 for (i = 0; i < qm->nrQMatoms; i++)
350 gmx_mtop_atomnr_to_atom(alook, qm->indexQM[i], &atom);
351 qm->nelectrons += mtop->atomtypes.atomnumber[atom->type];
352 qm->atomicnumberQM[i] = mtop->atomtypes.atomnumber[atom->type];
355 gmx_mtop_atomlookup_destroy(alook);
357 qm->QMcharge = ir->opts.QMcharge[grpnr];
358 qm->multiplicity = ir->opts.QMmult[grpnr];
359 qm->nelectrons -= ir->opts.QMcharge[grpnr];
361 qm->QMmethod = ir->opts.QMmethod[grpnr];
362 qm->QMbasis = ir->opts.QMbasis[grpnr];
363 /* trajectory surface hopping setup (Gaussian only) */
364 qm->bSH = ir->opts.bSH[grpnr];
365 qm->CASorbitals = ir->opts.CASorbitals[grpnr];
366 qm->CASelectrons = ir->opts.CASelectrons[grpnr];
367 qm->SAsteps = ir->opts.SAsteps[grpnr];
368 qm->SAon = ir->opts.SAon[grpnr];
369 qm->SAoff = ir->opts.SAoff[grpnr];
370 /* hack to prevent gaussian from reinitializing all the time */
371 qm->nQMcpus = 0; /* number of CPU's to be used by g01, is set
372 * upon initializing gaussian
375 /* print the current layer to allow users to check their input */
376 fprintf(stderr, "Layer %d\nnr of QM atoms %d\n", grpnr, nr);
377 fprintf(stderr, "QMlevel: %s/%s\n\n",
378 eQMmethod_names[qm->QMmethod], eQMbasis_names[qm->QMbasis]);
381 snew(qm->frontatoms, nr);
382 /* Lennard-Jones coefficients */
385 /* do we optimize the QM separately using the algorithms of the QM program??
387 qm->bTS = ir->opts.bTS[grpnr];
388 qm->bOPT = ir->opts.bOPT[grpnr];
392 t_QMrec *copy_QMrec(t_QMrec *qm)
394 /* copies the contents of qm into a new t_QMrec struct */
401 qmcopy->nrQMatoms = qm->nrQMatoms;
402 snew(qmcopy->xQM, qmcopy->nrQMatoms);
403 snew(qmcopy->indexQM, qmcopy->nrQMatoms);
404 snew(qmcopy->atomicnumberQM, qm->nrQMatoms);
405 snew(qmcopy->shiftQM, qmcopy->nrQMatoms); /* the shifts */
406 for (i = 0; i < qmcopy->nrQMatoms; i++)
408 qmcopy->shiftQM[i] = qm->shiftQM[i];
409 qmcopy->indexQM[i] = qm->indexQM[i];
410 qmcopy->atomicnumberQM[i] = qm->atomicnumberQM[i];
412 qmcopy->nelectrons = qm->nelectrons;
413 qmcopy->multiplicity = qm->multiplicity;
414 qmcopy->QMcharge = qm->QMcharge;
415 qmcopy->nelectrons = qm->nelectrons;
416 qmcopy->QMmethod = qm->QMmethod;
417 qmcopy->QMbasis = qm->QMbasis;
418 /* trajectory surface hopping setup (Gaussian only) */
419 qmcopy->bSH = qm->bSH;
420 qmcopy->CASorbitals = qm->CASorbitals;
421 qmcopy->CASelectrons = qm->CASelectrons;
422 qmcopy->SAsteps = qm->SAsteps;
423 qmcopy->SAon = qm->SAon;
424 qmcopy->SAoff = qm->SAoff;
425 qmcopy->bOPT = qm->bOPT;
427 /* Gaussian init. variables */
428 qmcopy->nQMcpus = qm->nQMcpus;
429 for (i = 0; i < DIM; i++)
431 qmcopy->SHbasis[i] = qm->SHbasis[i];
433 qmcopy->QMmem = qm->QMmem;
434 qmcopy->accuracy = qm->accuracy;
435 qmcopy->cpmcscf = qm->cpmcscf;
436 qmcopy->SAstep = qm->SAstep;
437 snew(qmcopy->frontatoms, qm->nrQMatoms);
438 snew(qmcopy->c12, qmcopy->nrQMatoms);
439 snew(qmcopy->c6, qmcopy->nrQMatoms);
440 if (qmcopy->bTS || qmcopy->bOPT)
442 for (i = 1; i < qmcopy->nrQMatoms; i++)
444 qmcopy->frontatoms[i] = qm->frontatoms[i];
445 qmcopy->c12[i] = qm->c12[i];
446 qmcopy->c6[i] = qm->c6[i];
454 t_QMMMrec *mk_QMMMrec(void)
465 void init_QMMMrec(t_commrec *cr,
470 /* we put the atomsnumbers of atoms that belong to the QMMM group in
471 * an array that will be copied later to QMMMrec->indexQM[..]. Also
472 * it will be used to create an QMMMrec->bQMMM index array that
473 * simply contains true/false for QM and MM (the other) atoms.
476 gmx_groups_t *groups;
477 atom_id *qm_arr = NULL, vsite, ai, aj;
478 int qm_max = 0, qm_nr = 0, i, j, jmax, k, l, nrvsite2 = 0;
483 gmx_mtop_atomloop_all_t aloop;
485 gmx_mtop_ilistloop_all_t iloop;
488 gmx_mtop_atomlookup_t alook;
490 c6au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 6));
491 c12au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 12));
492 /* issue a fatal if the user wants to run with more than one node */
495 gmx_fatal(FARGS, "QM/MM does not work in parallel, use a single rank instead\n");
498 /* Make a local copy of the QMMMrec */
501 /* bQMMM[..] is an array containing TRUE/FALSE for atoms that are
502 * QM/not QM. We first set all elemenst at false. Afterwards we use
503 * the qm_arr (=MMrec->indexQM) to changes the elements
504 * corresponding to the QM atoms at TRUE. */
506 qr->QMMMscheme = ir->QMMMscheme;
508 /* we take the possibility into account that a user has
509 * defined more than one QM group:
511 /* an ugly work-around in case there is only one group In this case
512 * the whole system is treated as QM. Otherwise the second group is
513 * always the rest of the total system and is treated as MM.
516 /* small problem if there is only QM.... so no MM */
518 jmax = ir->opts.ngQM;
520 if (qr->QMMMscheme == eQMMMschemeoniom)
522 qr->nrQMlayers = jmax;
529 groups = &mtop->groups;
531 /* there are jmax groups of QM atoms. In case of multiple QM groups
532 * I assume that the users wants to do ONIOM. However, maybe it
533 * should also be possible to define more than one QM subsystem with
534 * independent neighbourlists. I have to think about
538 for (j = 0; j < jmax; j++)
541 aloop = gmx_mtop_atomloop_all_init(mtop);
542 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
547 srenew(qm_arr, qm_max);
549 if (ggrpnr(groups, egcQMMM, i) == j)
555 if (qr->QMMMscheme == eQMMMschemeoniom)
557 /* add the atoms to the bQMMM array
560 /* I assume that users specify the QM groups from small to
561 * big(ger) in the mdp file
563 qr->qm[j] = mk_QMrec();
564 /* we need to throw out link atoms that in the previous layer
565 * existed to separate this QMlayer from the previous
566 * QMlayer. We use the iatoms array in the idef for that
567 * purpose. If all atoms defining the current Link Atom (Dummy2)
568 * are part of the current QM layer it needs to be removed from
571 iloop = gmx_mtop_ilistloop_all_init(mtop);
572 while (gmx_mtop_ilistloop_all_next(iloop, &ilist_mol, &a_offset))
574 nrvsite2 = ilist_mol[F_VSITE2].nr;
575 iatoms = ilist_mol[F_VSITE2].iatoms;
577 for (k = 0; k < nrvsite2; k += 4)
579 vsite = a_offset + iatoms[k+1]; /* the vsite */
580 ai = a_offset + iatoms[k+2]; /* constructing atom */
581 aj = a_offset + iatoms[k+3]; /* constructing atom */
582 if (ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, ai)
584 ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, aj))
586 /* this dummy link atom needs to be removed from the qm_arr
587 * before making the QMrec of this layer!
589 for (i = 0; i < qm_nr; i++)
591 if (qm_arr[i] == vsite)
593 /* drop the element */
594 for (l = i; l < qm_nr; l++)
596 qm_arr[l] = qm_arr[l+1];
605 /* store QM atoms in this layer in the QMrec and initialise layer
607 init_QMrec(j, qr->qm[j], qm_nr, qm_arr, mtop, ir);
609 /* we now store the LJ C6 and C12 parameters in QM rec in case
610 * we need to do an optimization
612 if (qr->qm[j]->bOPT || qr->qm[j]->bTS)
614 for (i = 0; i < qm_nr; i++)
616 /* nbfp now includes the 6.0/12.0 derivative prefactors */
617 qr->qm[j]->c6[i] = C6(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c6au/6.0;
618 qr->qm[j]->c12[i] = C12(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c12au/12.0;
621 /* now we check for frontier QM atoms. These occur in pairs that
622 * construct the vsite
624 iloop = gmx_mtop_ilistloop_all_init(mtop);
625 while (gmx_mtop_ilistloop_all_next(iloop, &ilist_mol, &a_offset))
627 nrvsite2 = ilist_mol[F_VSITE2].nr;
628 iatoms = ilist_mol[F_VSITE2].iatoms;
630 for (k = 0; k < nrvsite2; k += 4)
632 vsite = a_offset + iatoms[k+1]; /* the vsite */
633 ai = a_offset + iatoms[k+2]; /* constructing atom */
634 aj = a_offset + iatoms[k+3]; /* constructing atom */
635 if (ggrpnr(groups, egcQMMM, ai) < (groups->grps[egcQMMM].nr-1) &&
636 (ggrpnr(groups, egcQMMM, aj) >= (groups->grps[egcQMMM].nr-1)))
638 /* mark ai as frontier atom */
639 for (i = 0; i < qm_nr; i++)
641 if ( (qm_arr[i] == ai) || (qm_arr[i] == vsite) )
643 qr->qm[j]->frontatoms[i] = TRUE;
647 else if (ggrpnr(groups, egcQMMM, aj) < (groups->grps[egcQMMM].nr-1) &&
648 (ggrpnr(groups, egcQMMM, ai) >= (groups->grps[egcQMMM].nr-1)))
650 /* mark aj as frontier atom */
651 for (i = 0; i < qm_nr; i++)
653 if ( (qm_arr[i] == aj) || (qm_arr[i] == vsite))
655 qr->qm[j]->frontatoms[i] = TRUE;
663 if (qr->QMMMscheme != eQMMMschemeoniom)
666 /* standard QMMM, all layers are merged together so there is one QM
667 * subsystem and one MM subsystem.
668 * Also we set the charges to zero in the md->charge arrays to prevent
669 * the innerloops from doubly counting the electostatic QM MM interaction
672 alook = gmx_mtop_atomlookup_init(mtop);
674 for (k = 0; k < qm_nr; k++)
676 gmx_mtop_atomnr_to_atom(alook, qm_arr[k], &atom);
680 qr->qm[0] = mk_QMrec();
681 /* store QM atoms in the QMrec and initialise
683 init_QMrec(0, qr->qm[0], qm_nr, qm_arr, mtop, ir);
684 if (qr->qm[0]->bOPT || qr->qm[0]->bTS)
686 for (i = 0; i < qm_nr; i++)
688 gmx_mtop_atomnr_to_atom(alook, qm_arr[i], &atom);
689 /* nbfp now includes the 6.0/12.0 derivative prefactors */
690 qr->qm[0]->c6[i] = C6(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c6au/6.0;
691 qr->qm[0]->c12[i] = C12(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c12au/12.0;
695 /* find frontier atoms and mark them true in the frontieratoms array.
697 for (i = 0; i < qm_nr; i++)
699 gmx_mtop_atomnr_to_ilist(alook, qm_arr[i], &ilist_mol, &a_offset);
700 nrvsite2 = ilist_mol[F_VSITE2].nr;
701 iatoms = ilist_mol[F_VSITE2].iatoms;
703 for (k = 0; k < nrvsite2; k += 4)
705 vsite = a_offset + iatoms[k+1]; /* the vsite */
706 ai = a_offset + iatoms[k+2]; /* constructing atom */
707 aj = a_offset + iatoms[k+3]; /* constructing atom */
708 if (ggrpnr(groups, egcQMMM, ai) < (groups->grps[egcQMMM].nr-1) &&
709 (ggrpnr(groups, egcQMMM, aj) >= (groups->grps[egcQMMM].nr-1)))
711 /* mark ai as frontier atom */
712 if ( (qm_arr[i] == ai) || (qm_arr[i] == vsite) )
714 qr->qm[0]->frontatoms[i] = TRUE;
717 else if (ggrpnr(groups, egcQMMM, aj) < (groups->grps[egcQMMM].nr-1) &&
718 (ggrpnr(groups, egcQMMM, ai) >= (groups->grps[egcQMMM].nr-1)))
720 /* mark aj as frontier atom */
721 if ( (qm_arr[i] == aj) || (qm_arr[i] == vsite) )
723 qr->qm[0]->frontatoms[i] = TRUE;
729 gmx_mtop_atomlookup_destroy(alook);
731 /* MM rec creation */
733 mm->scalefactor = ir->scalefactor;
734 mm->nrMMatoms = (mtop->natoms)-(qr->qm[0]->nrQMatoms); /* rest of the atoms */
738 { /* MM rec creation */
740 mm->scalefactor = ir->scalefactor;
745 /* these variables get updated in the update QMMMrec */
747 if (qr->nrQMlayers == 1)
749 /* with only one layer there is only one initialisation
750 * needed. Multilayer is a bit more complicated as it requires
751 * re-initialisation at every step of the simulation. This is due
752 * to the use of COMMON blocks in the fortran QM subroutines.
754 if (qr->qm[0]->QMmethod < eQMmethodRHF)
756 #ifdef GMX_QMMM_MOPAC
757 /* semi-empiprical 1-layer ONIOM calculation requested (mopac93) */
758 init_mopac(cr, qr->qm[0], qr->mm);
760 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
765 /* ab initio calculation requested (gamess/gaussian/ORCA) */
766 #ifdef GMX_QMMM_GAMESS
767 init_gamess(cr, qr->qm[0], qr->mm);
768 #elif defined GMX_QMMM_GAUSSIAN
769 init_gaussian(cr, qr->qm[0], qr->mm);
770 #elif defined GMX_QMMM_ORCA
771 init_orca(qr->qm[0]);
773 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
779 void update_QMMMrec(t_commrec *cr,
786 /* updates the coordinates of both QM atoms and MM atoms and stores
787 * them in the QMMMrec.
789 * NOTE: is NOT yet working if there are no PBC. Also in ns.c, simple
790 * ns needs to be fixed!
793 mm_max = 0, mm_nr = 0, mm_nr_new, i, j, is, k, shift;
795 *mm_j_particles = NULL, *qm_i_particles = NULL;
811 *parallelMMarray = NULL;
815 c6au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 6));
816 c12au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 12));
818 /* every cpu has this array. On every processor we fill this array
819 * with 1's and 0's. 1's indicate the atoms is a QM atom on the
820 * current cpu in a later stage these arrays are all summed. indexes
821 * > 0 indicate the atom is a QM atom. Every node therefore knows
822 * whcih atoms are part of the QM subsystem.
824 /* copy some pointers */
827 QMMMlist = fr->QMMMlist;
831 /* init_pbc(box); needs to be called first, see pbc.h */
832 set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd : NULL, FALSE, box);
833 /* only in standard (normal) QMMM we need the neighbouring MM
834 * particles to provide a electric field of point charges for the QM
837 if (qr->QMMMscheme == eQMMMschemenormal) /* also implies 1 QM-layer */
839 /* we NOW create/update a number of QMMMrec entries:
841 * 1) the shiftQM, containing the shifts of the QM atoms
843 * 2) the indexMM array, containing the index of the MM atoms
845 * 3) the shiftMM, containing the shifts of the MM atoms
847 * 4) the shifted coordinates of the MM atoms
849 * the shifts are used for computing virial of the QM/MM particles.
851 qm = qr->qm[0]; /* in case of normal QMMM, there is only one group */
852 snew(qm_i_particles, QMMMlist.nri);
855 qm_i_particles[0].shift = XYZ2IS(0, 0, 0);
856 for (i = 0; i < QMMMlist.nri; i++)
858 qm_i_particles[i].j = QMMMlist.iinr[i];
862 qm_i_particles[i].shift = pbc_dx_aiuc(&pbc, x[QMMMlist.iinr[0]],
863 x[QMMMlist.iinr[i]], dx);
866 /* However, since nri >= nrQMatoms, we do a quicksort, and throw
867 * out double, triple, etc. entries later, as we do for the MM
871 /* compute the shift for the MM j-particles with respect to
872 * the QM i-particle and store them.
875 crd[0] = IS2X(QMMMlist.shift[i]) + IS2X(qm_i_particles[i].shift);
876 crd[1] = IS2Y(QMMMlist.shift[i]) + IS2Y(qm_i_particles[i].shift);
877 crd[2] = IS2Z(QMMMlist.shift[i]) + IS2Z(qm_i_particles[i].shift);
878 is = XYZ2IS(crd[0], crd[1], crd[2]);
879 for (j = QMMMlist.jindex[i];
880 j < QMMMlist.jindex[i+1];
886 srenew(mm_j_particles, mm_max);
889 mm_j_particles[mm_nr].j = QMMMlist.jjnr[j];
890 mm_j_particles[mm_nr].shift = is;
895 /* quicksort QM and MM shift arrays and throw away multiple entries */
899 qsort(qm_i_particles, QMMMlist.nri,
900 (size_t)sizeof(qm_i_particles[0]),
902 qsort(mm_j_particles, mm_nr,
903 (size_t)sizeof(mm_j_particles[0]),
905 /* remove multiples in the QM shift array, since in init_QMMM() we
906 * went through the atom numbers from 0 to md.nr, the order sorted
907 * here matches the one of QMindex already.
910 for (i = 0; i < QMMMlist.nri; i++)
912 if (i == 0 || qm_i_particles[i].j != qm_i_particles[i-1].j)
914 qm_i_particles[j++] = qm_i_particles[i];
918 if (qm->bTS || qm->bOPT)
920 /* only remove double entries for the MM array */
921 for (i = 0; i < mm_nr; i++)
923 if ((i == 0 || mm_j_particles[i].j != mm_j_particles[i-1].j)
924 && !md->bQM[mm_j_particles[i].j])
926 mm_j_particles[mm_nr_new++] = mm_j_particles[i];
930 /* we also remove mm atoms that have no charges!
931 * actually this is already done in the ns.c
935 for (i = 0; i < mm_nr; i++)
937 if ((i == 0 || mm_j_particles[i].j != mm_j_particles[i-1].j)
938 && !md->bQM[mm_j_particles[i].j]
939 && (md->chargeA[mm_j_particles[i].j]
940 || (md->chargeB && md->chargeB[mm_j_particles[i].j])))
942 mm_j_particles[mm_nr_new++] = mm_j_particles[i];
947 /* store the data retrieved above into the QMMMrec
950 /* Keep the compiler happy,
951 * shift will always be set in the loop for i=0
954 for (i = 0; i < qm->nrQMatoms; i++)
956 /* not all qm particles might have appeared as i
957 * particles. They might have been part of the same charge
958 * group for instance.
960 if (qm->indexQM[i] == qm_i_particles[k].j)
962 shift = qm_i_particles[k++].shift;
964 /* use previous shift, assuming they belong the same charge
968 qm->shiftQM[i] = shift;
971 /* parallel excecution */
974 snew(parallelMMarray, 2*(md->nr));
975 /* only MM particles have a 1 at their atomnumber. The second part
976 * of the array contains the shifts. Thus:
977 * p[i]=1/0 depending on wether atomnumber i is a MM particle in the QM
978 * step or not. p[i+md->nr] is the shift of atomnumber i.
980 for (i = 0; i < 2*(md->nr); i++)
982 parallelMMarray[i] = 0;
985 for (i = 0; i < mm_nr; i++)
987 parallelMMarray[mm_j_particles[i].j] = 1;
988 parallelMMarray[mm_j_particles[i].j+(md->nr)] = mm_j_particles[i].shift;
990 gmx_sumi(md->nr, parallelMMarray, cr);
994 for (i = 0; i < md->nr; i++)
996 if (parallelMMarray[i])
1001 srenew(mm->indexMM, mm_max);
1002 srenew(mm->shiftMM, mm_max);
1004 mm->indexMM[mm_nr] = i;
1005 mm->shiftMM[mm_nr++] = parallelMMarray[i+md->nr]/parallelMMarray[i];
1008 mm->nrMMatoms = mm_nr;
1009 free(parallelMMarray);
1011 /* serial execution */
1014 mm->nrMMatoms = mm_nr;
1015 srenew(mm->shiftMM, mm_nr);
1016 srenew(mm->indexMM, mm_nr);
1017 for (i = 0; i < mm_nr; i++)
1019 mm->indexMM[i] = mm_j_particles[i].j;
1020 mm->shiftMM[i] = mm_j_particles[i].shift;
1024 /* (re) allocate memory for the MM coordiate array. The QM
1025 * coordinate array was already allocated in init_QMMM, and is
1026 * only (re)filled in the update_QMMM_coordinates routine
1028 srenew(mm->xMM, mm->nrMMatoms);
1029 /* now we (re) fill the array that contains the MM charges with
1030 * the forcefield charges. If requested, these charges will be
1031 * scaled by a factor
1033 srenew(mm->MMcharges, mm->nrMMatoms);
1034 for (i = 0; i < mm->nrMMatoms; i++) /* no free energy yet */
1036 mm->MMcharges[i] = md->chargeA[mm->indexMM[i]]*mm->scalefactor;
1038 if (qm->bTS || qm->bOPT)
1040 /* store (copy) the c6 and c12 parameters into the MMrec struct
1042 srenew(mm->c6, mm->nrMMatoms);
1043 srenew(mm->c12, mm->nrMMatoms);
1044 for (i = 0; i < mm->nrMMatoms; i++)
1046 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1047 mm->c6[i] = C6(fr->nbfp, top->idef.atnr, md->typeA[mm->indexMM[i]], md->typeA[mm->indexMM[i]])/c6au/6.0;
1048 mm->c12[i] = C12(fr->nbfp, top->idef.atnr, md->typeA[mm->indexMM[i]], md->typeA[mm->indexMM[i]])/c12au/12.0;
1050 punch_QMMM_excl(qr->qm[0], mm, &(top->excls));
1052 /* the next routine fills the coordinate fields in the QMMM rec of
1053 * both the qunatum atoms and the MM atoms, using the shifts
1057 update_QMMM_coord(x, fr, qr->qm[0], qr->mm);
1058 free(qm_i_particles);
1059 free(mm_j_particles);
1061 else /* ONIOM */ /* ????? */
1064 /* do for each layer */
1065 for (j = 0; j < qr->nrQMlayers; j++)
1068 qm->shiftQM[0] = XYZ2IS(0, 0, 0);
1069 for (i = 1; i < qm->nrQMatoms; i++)
1071 qm->shiftQM[i] = pbc_dx_aiuc(&pbc, x[qm->indexQM[0]], x[qm->indexQM[i]],
1074 update_QMMM_coord(x, fr, qm, mm);
1077 } /* update_QMMM_rec */
1080 real calculate_QMMM(t_commrec *cr,
1086 /* a selection for the QM package depending on which is requested
1087 * (Gaussian, GAMESS-UK, MOPAC or ORCA) needs to be implemented here. Now
1088 * it works through defines.... Not so nice yet
1097 *forces = NULL, *fshift = NULL,
1098 *forces2 = NULL, *fshift2 = NULL; /* needed for multilayer ONIOM */
1101 /* make a local copy the QMMMrec pointer
1106 /* now different procedures are carried out for one layer ONION and
1107 * normal QMMM on one hand and multilayer oniom on the other
1109 if (qr->QMMMscheme == eQMMMschemenormal || qr->nrQMlayers == 1)
1112 snew(forces, (qm->nrQMatoms+mm->nrMMatoms));
1113 snew(fshift, (qm->nrQMatoms+mm->nrMMatoms));
1114 QMener = call_QMroutine(cr, fr, qm, mm, forces, fshift);
1115 for (i = 0; i < qm->nrQMatoms; i++)
1117 for (j = 0; j < DIM; j++)
1119 f[qm->indexQM[i]][j] -= forces[i][j];
1120 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
1123 for (i = 0; i < mm->nrMMatoms; i++)
1125 for (j = 0; j < DIM; j++)
1127 f[mm->indexMM[i]][j] -= forces[qm->nrQMatoms+i][j];
1128 fr->fshift[mm->shiftMM[i]][j] += fshift[qm->nrQMatoms+i][j];
1135 else /* Multi-layer ONIOM */
1137 for (i = 0; i < qr->nrQMlayers-1; i++) /* last layer is special */
1140 qm2 = copy_QMrec(qr->qm[i+1]);
1142 qm2->nrQMatoms = qm->nrQMatoms;
1144 for (j = 0; j < qm2->nrQMatoms; j++)
1146 for (k = 0; k < DIM; k++)
1148 qm2->xQM[j][k] = qm->xQM[j][k];
1150 qm2->indexQM[j] = qm->indexQM[j];
1151 qm2->atomicnumberQM[j] = qm->atomicnumberQM[j];
1152 qm2->shiftQM[j] = qm->shiftQM[j];
1155 qm2->QMcharge = qm->QMcharge;
1156 /* this layer at the higher level of theory */
1157 srenew(forces, qm->nrQMatoms);
1158 srenew(fshift, qm->nrQMatoms);
1159 /* we need to re-initialize the QMroutine every step... */
1160 init_QMroutine(cr, qm, mm);
1161 QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
1163 /* this layer at the lower level of theory */
1164 srenew(forces2, qm->nrQMatoms);
1165 srenew(fshift2, qm->nrQMatoms);
1166 init_QMroutine(cr, qm2, mm);
1167 QMener -= call_QMroutine(cr, fr, qm2, mm, forces2, fshift2);
1168 /* E = E1high-E1low The next layer includes the current layer at
1169 * the lower level of theory, which provides + E2low
1170 * this is similar for gradients
1172 for (i = 0; i < qm->nrQMatoms; i++)
1174 for (j = 0; j < DIM; j++)
1176 f[qm->indexQM[i]][j] -= (forces[i][j]-forces2[i][j]);
1177 fr->fshift[qm->shiftQM[i]][j] += (fshift[i][j]-fshift2[i][j]);
1182 /* now the last layer still needs to be done: */
1183 qm = qr->qm[qr->nrQMlayers-1]; /* C counts from 0 */
1184 init_QMroutine(cr, qm, mm);
1185 srenew(forces, qm->nrQMatoms);
1186 srenew(fshift, qm->nrQMatoms);
1187 QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
1188 for (i = 0; i < qm->nrQMatoms; i++)
1190 for (j = 0; j < DIM; j++)
1192 f[qm->indexQM[i]][j] -= forces[i][j];
1193 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
1201 if (qm->bTS || qm->bOPT)
1203 /* qm[0] still contains the largest ONIOM QM subsystem
1204 * we take the optimized coordiates and put the in x[]
1206 for (i = 0; i < qm->nrQMatoms; i++)
1208 for (j = 0; j < DIM; j++)
1210 x[qm->indexQM[i]][j] = qm->xQM[i][j];
1215 } /* calculate_QMMM */
1217 /* end of QMMM core routines */