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.
44 #include "types/commrec.h"
46 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/fileio/confio.h"
64 #include "gmx_fatal.h"
67 #include "mtop_util.h"
70 /* declarations of the interfaces to the QM packages. The _SH indicate
71 * the QM interfaces can be used for Surface Hopping simulations
73 #ifdef GMX_QMMM_GAMESS
74 /* GAMESS interface */
77 init_gamess(t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
80 call_gamess(t_commrec *cr, t_forcerec *fr,
81 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
83 #elif defined GMX_QMMM_MOPAC
87 init_mopac(t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
90 call_mopac(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
91 t_MMrec *mm, rvec f[], rvec fshift[]);
94 call_mopac_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
95 t_MMrec *mm, rvec f[], rvec fshift[]);
97 #elif defined GMX_QMMM_GAUSSIAN
98 /* GAUSSIAN interface */
101 init_gaussian(t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
104 call_gaussian_SH(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
105 t_MMrec *mm, rvec f[], rvec fshift[]);
108 call_gaussian(t_commrec *cr, t_forcerec *fr, t_QMrec *qm,
109 t_MMrec *mm, rvec f[], rvec fshift[]);
111 #elif defined GMX_QMMM_ORCA
115 init_orca(t_QMrec *qm);
118 call_orca(t_forcerec *fr, t_QMrec *qm,
119 t_MMrec *mm, rvec f[], rvec fshift[]);
126 /* this struct and these comparison functions are needed for creating
127 * a QMMM input for the QM routines from the QMMM neighbor list.
135 static int struct_comp(const void *a, const void *b)
138 return (int)(((t_j_particle *)a)->j)-(int)(((t_j_particle *)b)->j);
142 static int int_comp(const void *a, const void *b)
145 return (*(int *)a) - (*(int *)b);
149 static int QMlayer_comp(const void *a, const void *b)
152 return (int)(((t_QMrec *)a)->nrQMatoms)-(int)(((t_QMrec *)b)->nrQMatoms);
156 real call_QMroutine(t_commrec gmx_unused *cr, t_forcerec gmx_unused *fr, t_QMrec gmx_unused *qm,
157 t_MMrec gmx_unused *mm, rvec gmx_unused f[], rvec gmx_unused fshift[])
159 /* makes a call to the requested QM routine (qm->QMmethod)
160 * Note that f is actually the gradient, i.e. -f
165 /* do a semi-empiprical calculation */
167 if (qm->QMmethod < eQMmethodRHF && !(mm->nrMMatoms))
169 #ifdef GMX_QMMM_MOPAC
172 QMener = call_mopac_SH(cr, fr, qm, mm, f, fshift);
176 QMener = call_mopac(cr, fr, qm, mm, f, fshift);
179 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
184 /* do an ab-initio calculation */
185 if (qm->bSH && qm->QMmethod == eQMmethodCASSCF)
187 #ifdef GMX_QMMM_GAUSSIAN
188 QMener = call_gaussian_SH(cr, fr, qm, mm, f, fshift);
190 gmx_fatal(FARGS, "Ab-initio Surface-hopping only supported with Gaussian.");
195 #ifdef GMX_QMMM_GAMESS
196 QMener = call_gamess(cr, fr, qm, mm, f, fshift);
197 #elif defined GMX_QMMM_GAUSSIAN
198 QMener = call_gaussian(cr, fr, qm, mm, f, fshift);
199 #elif defined GMX_QMMM_ORCA
200 QMener = call_orca(fr, qm, mm, f, fshift);
202 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
209 void init_QMroutine(t_commrec gmx_unused *cr, t_QMrec gmx_unused *qm, t_MMrec gmx_unused *mm)
211 /* makes a call to the requested QM routine (qm->QMmethod)
213 if (qm->QMmethod < eQMmethodRHF)
215 #ifdef GMX_QMMM_MOPAC
216 /* do a semi-empiprical calculation */
217 init_mopac(cr, qm, mm);
219 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
224 /* do an ab-initio calculation */
225 #ifdef GMX_QMMM_GAMESS
226 init_gamess(cr, qm, mm);
227 #elif defined GMX_QMMM_GAUSSIAN
228 init_gaussian(cr, qm, mm);
229 #elif defined GMX_QMMM_ORCA
232 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
235 } /* init_QMroutine */
237 void update_QMMM_coord(rvec x[], t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
239 /* shifts the QM and MM particles into the central box and stores
240 * these shifted coordinates in the coordinate arrays of the
241 * QMMMrec. These coordinates are passed on the QM subroutines.
246 /* shift the QM atoms into the central box
248 for (i = 0; i < qm->nrQMatoms; i++)
250 rvec_sub(x[qm->indexQM[i]], fr->shift_vec[qm->shiftQM[i]], qm->xQM[i]);
252 /* also shift the MM atoms into the central box, if any
254 for (i = 0; i < mm->nrMMatoms; i++)
256 rvec_sub(x[mm->indexMM[i]], fr->shift_vec[mm->shiftMM[i]], mm->xMM[i]);
258 } /* update_QMMM_coord */
260 static void punch_QMMM_excl(t_QMrec *qm, t_MMrec *mm, t_blocka *excls)
262 /* punch a file containing the bonded interactions of each QM
263 * atom with MM atoms. These need to be excluded in the QM routines
264 * Only needed in case of QM/MM optimizations
269 i, j, k, nrexcl = 0, *excluded = NULL, max = 0;
272 out = fopen("QMMMexcl.dat", "w");
274 /* this can be done more efficiently I think
276 for (i = 0; i < qm->nrQMatoms; i++)
279 for (j = excls->index[qm->indexQM[i]];
280 j < excls->index[qm->indexQM[i]+1];
283 for (k = 0; k < mm->nrMMatoms; k++)
285 if (mm->indexMM[k] == excls->a[j]) /* the excluded MM atom */
290 srenew(excluded, max);
292 excluded[nrexcl++] = k;
298 fprintf(out, "%5d %5d\n", i+1, nrexcl);
299 for (j = 0; j < nrexcl; j++)
301 fprintf(out, "%5d ", excluded[j]);
307 } /* punch_QMMM_excl */
310 /* end of QMMM subroutines */
312 /* QMMM core routines */
314 t_QMrec *mk_QMrec(void)
321 t_MMrec *mk_MMrec(void)
328 static void init_QMrec(int grpnr, t_QMrec *qm, int nr, int *atomarray,
329 gmx_mtop_t *mtop, t_inputrec *ir)
331 /* fills the t_QMrec struct of QM group grpnr
334 gmx_mtop_atomlookup_t alook;
340 snew(qm->indexQM, nr);
341 snew(qm->shiftQM, nr); /* the shifts */
342 for (i = 0; i < nr; i++)
344 qm->indexQM[i] = atomarray[i];
347 alook = gmx_mtop_atomlookup_init(mtop);
349 snew(qm->atomicnumberQM, nr);
350 for (i = 0; i < qm->nrQMatoms; i++)
352 gmx_mtop_atomnr_to_atom(alook, qm->indexQM[i], &atom);
353 qm->nelectrons += mtop->atomtypes.atomnumber[atom->type];
354 qm->atomicnumberQM[i] = mtop->atomtypes.atomnumber[atom->type];
357 gmx_mtop_atomlookup_destroy(alook);
359 qm->QMcharge = ir->opts.QMcharge[grpnr];
360 qm->multiplicity = ir->opts.QMmult[grpnr];
361 qm->nelectrons -= ir->opts.QMcharge[grpnr];
363 qm->QMmethod = ir->opts.QMmethod[grpnr];
364 qm->QMbasis = ir->opts.QMbasis[grpnr];
365 /* trajectory surface hopping setup (Gaussian only) */
366 qm->bSH = ir->opts.bSH[grpnr];
367 qm->CASorbitals = ir->opts.CASorbitals[grpnr];
368 qm->CASelectrons = ir->opts.CASelectrons[grpnr];
369 qm->SAsteps = ir->opts.SAsteps[grpnr];
370 qm->SAon = ir->opts.SAon[grpnr];
371 qm->SAoff = ir->opts.SAoff[grpnr];
372 /* hack to prevent gaussian from reinitializing all the time */
373 qm->nQMcpus = 0; /* number of CPU's to be used by g01, is set
374 * upon initializing gaussian
377 /* print the current layer to allow users to check their input */
378 fprintf(stderr, "Layer %d\nnr of QM atoms %d\n", grpnr, nr);
379 fprintf(stderr, "QMlevel: %s/%s\n\n",
380 eQMmethod_names[qm->QMmethod], eQMbasis_names[qm->QMbasis]);
383 snew(qm->frontatoms, nr);
384 /* Lennard-Jones coefficients */
387 /* do we optimize the QM separately using the algorithms of the QM program??
389 qm->bTS = ir->opts.bTS[grpnr];
390 qm->bOPT = ir->opts.bOPT[grpnr];
394 t_QMrec *copy_QMrec(t_QMrec *qm)
396 /* copies the contents of qm into a new t_QMrec struct */
403 qmcopy->nrQMatoms = qm->nrQMatoms;
404 snew(qmcopy->xQM, qmcopy->nrQMatoms);
405 snew(qmcopy->indexQM, qmcopy->nrQMatoms);
406 snew(qmcopy->atomicnumberQM, qm->nrQMatoms);
407 snew(qmcopy->shiftQM, qmcopy->nrQMatoms); /* the shifts */
408 for (i = 0; i < qmcopy->nrQMatoms; i++)
410 qmcopy->shiftQM[i] = qm->shiftQM[i];
411 qmcopy->indexQM[i] = qm->indexQM[i];
412 qmcopy->atomicnumberQM[i] = qm->atomicnumberQM[i];
414 qmcopy->nelectrons = qm->nelectrons;
415 qmcopy->multiplicity = qm->multiplicity;
416 qmcopy->QMcharge = qm->QMcharge;
417 qmcopy->nelectrons = qm->nelectrons;
418 qmcopy->QMmethod = qm->QMmethod;
419 qmcopy->QMbasis = qm->QMbasis;
420 /* trajectory surface hopping setup (Gaussian only) */
421 qmcopy->bSH = qm->bSH;
422 qmcopy->CASorbitals = qm->CASorbitals;
423 qmcopy->CASelectrons = qm->CASelectrons;
424 qmcopy->SAsteps = qm->SAsteps;
425 qmcopy->SAon = qm->SAon;
426 qmcopy->SAoff = qm->SAoff;
427 qmcopy->bOPT = qm->bOPT;
429 /* Gaussian init. variables */
430 qmcopy->nQMcpus = qm->nQMcpus;
431 for (i = 0; i < DIM; i++)
433 qmcopy->SHbasis[i] = qm->SHbasis[i];
435 qmcopy->QMmem = qm->QMmem;
436 qmcopy->accuracy = qm->accuracy;
437 qmcopy->cpmcscf = qm->cpmcscf;
438 qmcopy->SAstep = qm->SAstep;
439 snew(qmcopy->frontatoms, qm->nrQMatoms);
440 snew(qmcopy->c12, qmcopy->nrQMatoms);
441 snew(qmcopy->c6, qmcopy->nrQMatoms);
442 if (qmcopy->bTS || qmcopy->bOPT)
444 for (i = 1; i < qmcopy->nrQMatoms; i++)
446 qmcopy->frontatoms[i] = qm->frontatoms[i];
447 qmcopy->c12[i] = qm->c12[i];
448 qmcopy->c6[i] = qm->c6[i];
456 t_QMMMrec *mk_QMMMrec(void)
467 void init_QMMMrec(t_commrec *cr,
472 /* we put the atomsnumbers of atoms that belong to the QMMM group in
473 * an array that will be copied later to QMMMrec->indexQM[..]. Also
474 * it will be used to create an QMMMrec->bQMMM index array that
475 * simply contains true/false for QM and MM (the other) atoms.
478 gmx_groups_t *groups;
479 atom_id *qm_arr = NULL, vsite, ai, aj;
480 int qm_max = 0, qm_nr = 0, i, j, jmax, k, l, nrvsite2 = 0;
485 gmx_mtop_atomloop_all_t aloop;
487 gmx_mtop_ilistloop_all_t iloop;
490 gmx_mtop_atomlookup_t alook;
492 c6au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 6));
493 c12au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 12));
494 /* issue a fatal if the user wants to run with more than one node */
497 gmx_fatal(FARGS, "QM/MM does not work in parallel, use a single rank instead\n");
500 /* Make a local copy of the QMMMrec */
503 /* bQMMM[..] is an array containing TRUE/FALSE for atoms that are
504 * QM/not QM. We first set all elemenst at false. Afterwards we use
505 * the qm_arr (=MMrec->indexQM) to changes the elements
506 * corresponding to the QM atoms at TRUE. */
508 qr->QMMMscheme = ir->QMMMscheme;
510 /* we take the possibility into account that a user has
511 * defined more than one QM group:
513 /* an ugly work-around in case there is only one group In this case
514 * the whole system is treated as QM. Otherwise the second group is
515 * always the rest of the total system and is treated as MM.
518 /* small problem if there is only QM.... so no MM */
520 jmax = ir->opts.ngQM;
522 if (qr->QMMMscheme == eQMMMschemeoniom)
524 qr->nrQMlayers = jmax;
531 groups = &mtop->groups;
533 /* there are jmax groups of QM atoms. In case of multiple QM groups
534 * I assume that the users wants to do ONIOM. However, maybe it
535 * should also be possible to define more than one QM subsystem with
536 * independent neighbourlists. I have to think about
540 for (j = 0; j < jmax; j++)
543 aloop = gmx_mtop_atomloop_all_init(mtop);
544 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
549 srenew(qm_arr, qm_max);
551 if (ggrpnr(groups, egcQMMM, i) == j)
557 if (qr->QMMMscheme == eQMMMschemeoniom)
559 /* add the atoms to the bQMMM array
562 /* I assume that users specify the QM groups from small to
563 * big(ger) in the mdp file
565 qr->qm[j] = mk_QMrec();
566 /* we need to throw out link atoms that in the previous layer
567 * existed to separate this QMlayer from the previous
568 * QMlayer. We use the iatoms array in the idef for that
569 * purpose. If all atoms defining the current Link Atom (Dummy2)
570 * are part of the current QM layer it needs to be removed from
573 iloop = gmx_mtop_ilistloop_all_init(mtop);
574 while (gmx_mtop_ilistloop_all_next(iloop, &ilist_mol, &a_offset))
576 nrvsite2 = ilist_mol[F_VSITE2].nr;
577 iatoms = ilist_mol[F_VSITE2].iatoms;
579 for (k = 0; k < nrvsite2; k += 4)
581 vsite = a_offset + iatoms[k+1]; /* the vsite */
582 ai = a_offset + iatoms[k+2]; /* constructing atom */
583 aj = a_offset + iatoms[k+3]; /* constructing atom */
584 if (ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, ai)
586 ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, aj))
588 /* this dummy link atom needs to be removed from the qm_arr
589 * before making the QMrec of this layer!
591 for (i = 0; i < qm_nr; i++)
593 if (qm_arr[i] == vsite)
595 /* drop the element */
596 for (l = i; l < qm_nr; l++)
598 qm_arr[l] = qm_arr[l+1];
607 /* store QM atoms in this layer in the QMrec and initialise layer
609 init_QMrec(j, qr->qm[j], qm_nr, qm_arr, mtop, ir);
611 /* we now store the LJ C6 and C12 parameters in QM rec in case
612 * we need to do an optimization
614 if (qr->qm[j]->bOPT || qr->qm[j]->bTS)
616 for (i = 0; i < qm_nr; i++)
618 /* nbfp now includes the 6.0/12.0 derivative prefactors */
619 qr->qm[j]->c6[i] = C6(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c6au/6.0;
620 qr->qm[j]->c12[i] = C12(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c12au/12.0;
623 /* now we check for frontier QM atoms. These occur in pairs that
624 * construct the vsite
626 iloop = gmx_mtop_ilistloop_all_init(mtop);
627 while (gmx_mtop_ilistloop_all_next(iloop, &ilist_mol, &a_offset))
629 nrvsite2 = ilist_mol[F_VSITE2].nr;
630 iatoms = ilist_mol[F_VSITE2].iatoms;
632 for (k = 0; k < nrvsite2; k += 4)
634 vsite = a_offset + iatoms[k+1]; /* the vsite */
635 ai = a_offset + iatoms[k+2]; /* constructing atom */
636 aj = a_offset + iatoms[k+3]; /* constructing atom */
637 if (ggrpnr(groups, egcQMMM, ai) < (groups->grps[egcQMMM].nr-1) &&
638 (ggrpnr(groups, egcQMMM, aj) >= (groups->grps[egcQMMM].nr-1)))
640 /* mark ai as frontier atom */
641 for (i = 0; i < qm_nr; i++)
643 if ( (qm_arr[i] == ai) || (qm_arr[i] == vsite) )
645 qr->qm[j]->frontatoms[i] = TRUE;
649 else if (ggrpnr(groups, egcQMMM, aj) < (groups->grps[egcQMMM].nr-1) &&
650 (ggrpnr(groups, egcQMMM, ai) >= (groups->grps[egcQMMM].nr-1)))
652 /* mark aj as frontier atom */
653 for (i = 0; i < qm_nr; i++)
655 if ( (qm_arr[i] == aj) || (qm_arr[i] == vsite))
657 qr->qm[j]->frontatoms[i] = TRUE;
665 if (qr->QMMMscheme != eQMMMschemeoniom)
668 /* standard QMMM, all layers are merged together so there is one QM
669 * subsystem and one MM subsystem.
670 * Also we set the charges to zero in the md->charge arrays to prevent
671 * the innerloops from doubly counting the electostatic QM MM interaction
674 alook = gmx_mtop_atomlookup_init(mtop);
676 for (k = 0; k < qm_nr; k++)
678 gmx_mtop_atomnr_to_atom(alook, qm_arr[k], &atom);
682 qr->qm[0] = mk_QMrec();
683 /* store QM atoms in the QMrec and initialise
685 init_QMrec(0, qr->qm[0], qm_nr, qm_arr, mtop, ir);
686 if (qr->qm[0]->bOPT || qr->qm[0]->bTS)
688 for (i = 0; i < qm_nr; i++)
690 gmx_mtop_atomnr_to_atom(alook, qm_arr[i], &atom);
691 /* nbfp now includes the 6.0/12.0 derivative prefactors */
692 qr->qm[0]->c6[i] = C6(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c6au/6.0;
693 qr->qm[0]->c12[i] = C12(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c12au/12.0;
697 /* find frontier atoms and mark them true in the frontieratoms array.
699 for (i = 0; i < qm_nr; i++)
701 gmx_mtop_atomnr_to_ilist(alook, qm_arr[i], &ilist_mol, &a_offset);
702 nrvsite2 = ilist_mol[F_VSITE2].nr;
703 iatoms = ilist_mol[F_VSITE2].iatoms;
705 for (k = 0; k < nrvsite2; k += 4)
707 vsite = a_offset + iatoms[k+1]; /* the vsite */
708 ai = a_offset + iatoms[k+2]; /* constructing atom */
709 aj = a_offset + iatoms[k+3]; /* constructing atom */
710 if (ggrpnr(groups, egcQMMM, ai) < (groups->grps[egcQMMM].nr-1) &&
711 (ggrpnr(groups, egcQMMM, aj) >= (groups->grps[egcQMMM].nr-1)))
713 /* mark ai as frontier atom */
714 if ( (qm_arr[i] == ai) || (qm_arr[i] == vsite) )
716 qr->qm[0]->frontatoms[i] = TRUE;
719 else if (ggrpnr(groups, egcQMMM, aj) < (groups->grps[egcQMMM].nr-1) &&
720 (ggrpnr(groups, egcQMMM, ai) >= (groups->grps[egcQMMM].nr-1)))
722 /* mark aj as frontier atom */
723 if ( (qm_arr[i] == aj) || (qm_arr[i] == vsite) )
725 qr->qm[0]->frontatoms[i] = TRUE;
731 gmx_mtop_atomlookup_destroy(alook);
733 /* MM rec creation */
735 mm->scalefactor = ir->scalefactor;
736 mm->nrMMatoms = (mtop->natoms)-(qr->qm[0]->nrQMatoms); /* rest of the atoms */
740 { /* MM rec creation */
742 mm->scalefactor = ir->scalefactor;
747 /* these variables get updated in the update QMMMrec */
749 if (qr->nrQMlayers == 1)
751 /* with only one layer there is only one initialisation
752 * needed. Multilayer is a bit more complicated as it requires
753 * re-initialisation at every step of the simulation. This is due
754 * to the use of COMMON blocks in the fortran QM subroutines.
756 if (qr->qm[0]->QMmethod < eQMmethodRHF)
758 #ifdef GMX_QMMM_MOPAC
759 /* semi-empiprical 1-layer ONIOM calculation requested (mopac93) */
760 init_mopac(cr, qr->qm[0], qr->mm);
762 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
767 /* ab initio calculation requested (gamess/gaussian/ORCA) */
768 #ifdef GMX_QMMM_GAMESS
769 init_gamess(cr, qr->qm[0], qr->mm);
770 #elif defined GMX_QMMM_GAUSSIAN
771 init_gaussian(cr, qr->qm[0], qr->mm);
772 #elif defined GMX_QMMM_ORCA
773 init_orca(qr->qm[0]);
775 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
781 void update_QMMMrec(t_commrec *cr,
788 /* updates the coordinates of both QM atoms and MM atoms and stores
789 * them in the QMMMrec.
791 * NOTE: is NOT yet working if there are no PBC. Also in ns.c, simple
792 * ns needs to be fixed!
795 mm_max = 0, mm_nr = 0, mm_nr_new, i, j, is, k, shift;
797 *mm_j_particles = NULL, *qm_i_particles = NULL;
813 *parallelMMarray = NULL;
817 c6au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 6));
818 c12au = (HARTREE2KJ*AVOGADRO*pow(BOHR2NM, 12));
820 /* every cpu has this array. On every processor we fill this array
821 * with 1's and 0's. 1's indicate the atoms is a QM atom on the
822 * current cpu in a later stage these arrays are all summed. indexes
823 * > 0 indicate the atom is a QM atom. Every node therefore knows
824 * whcih atoms are part of the QM subsystem.
826 /* copy some pointers */
829 QMMMlist = fr->QMMMlist;
833 /* init_pbc(box); needs to be called first, see pbc.h */
834 set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd : NULL, FALSE, box);
835 /* only in standard (normal) QMMM we need the neighbouring MM
836 * particles to provide a electric field of point charges for the QM
839 if (qr->QMMMscheme == eQMMMschemenormal) /* also implies 1 QM-layer */
841 /* we NOW create/update a number of QMMMrec entries:
843 * 1) the shiftQM, containing the shifts of the QM atoms
845 * 2) the indexMM array, containing the index of the MM atoms
847 * 3) the shiftMM, containing the shifts of the MM atoms
849 * 4) the shifted coordinates of the MM atoms
851 * the shifts are used for computing virial of the QM/MM particles.
853 qm = qr->qm[0]; /* in case of normal QMMM, there is only one group */
854 snew(qm_i_particles, QMMMlist.nri);
857 qm_i_particles[0].shift = XYZ2IS(0, 0, 0);
858 for (i = 0; i < QMMMlist.nri; i++)
860 qm_i_particles[i].j = QMMMlist.iinr[i];
864 qm_i_particles[i].shift = pbc_dx_aiuc(&pbc, x[QMMMlist.iinr[0]],
865 x[QMMMlist.iinr[i]], dx);
868 /* However, since nri >= nrQMatoms, we do a quicksort, and throw
869 * out double, triple, etc. entries later, as we do for the MM
873 /* compute the shift for the MM j-particles with respect to
874 * the QM i-particle and store them.
877 crd[0] = IS2X(QMMMlist.shift[i]) + IS2X(qm_i_particles[i].shift);
878 crd[1] = IS2Y(QMMMlist.shift[i]) + IS2Y(qm_i_particles[i].shift);
879 crd[2] = IS2Z(QMMMlist.shift[i]) + IS2Z(qm_i_particles[i].shift);
880 is = XYZ2IS(crd[0], crd[1], crd[2]);
881 for (j = QMMMlist.jindex[i];
882 j < QMMMlist.jindex[i+1];
888 srenew(mm_j_particles, mm_max);
891 mm_j_particles[mm_nr].j = QMMMlist.jjnr[j];
892 mm_j_particles[mm_nr].shift = is;
897 /* quicksort QM and MM shift arrays and throw away multiple entries */
901 qsort(qm_i_particles, QMMMlist.nri,
902 (size_t)sizeof(qm_i_particles[0]),
904 qsort(mm_j_particles, mm_nr,
905 (size_t)sizeof(mm_j_particles[0]),
907 /* remove multiples in the QM shift array, since in init_QMMM() we
908 * went through the atom numbers from 0 to md.nr, the order sorted
909 * here matches the one of QMindex already.
912 for (i = 0; i < QMMMlist.nri; i++)
914 if (i == 0 || qm_i_particles[i].j != qm_i_particles[i-1].j)
916 qm_i_particles[j++] = qm_i_particles[i];
920 if (qm->bTS || qm->bOPT)
922 /* only remove double entries for the MM array */
923 for (i = 0; i < mm_nr; i++)
925 if ((i == 0 || mm_j_particles[i].j != mm_j_particles[i-1].j)
926 && !md->bQM[mm_j_particles[i].j])
928 mm_j_particles[mm_nr_new++] = mm_j_particles[i];
932 /* we also remove mm atoms that have no charges!
933 * actually this is already done in the ns.c
937 for (i = 0; i < mm_nr; i++)
939 if ((i == 0 || mm_j_particles[i].j != mm_j_particles[i-1].j)
940 && !md->bQM[mm_j_particles[i].j]
941 && (md->chargeA[mm_j_particles[i].j]
942 || (md->chargeB && md->chargeB[mm_j_particles[i].j])))
944 mm_j_particles[mm_nr_new++] = mm_j_particles[i];
949 /* store the data retrieved above into the QMMMrec
952 /* Keep the compiler happy,
953 * shift will always be set in the loop for i=0
956 for (i = 0; i < qm->nrQMatoms; i++)
958 /* not all qm particles might have appeared as i
959 * particles. They might have been part of the same charge
960 * group for instance.
962 if (qm->indexQM[i] == qm_i_particles[k].j)
964 shift = qm_i_particles[k++].shift;
966 /* use previous shift, assuming they belong the same charge
970 qm->shiftQM[i] = shift;
973 /* parallel excecution */
976 snew(parallelMMarray, 2*(md->nr));
977 /* only MM particles have a 1 at their atomnumber. The second part
978 * of the array contains the shifts. Thus:
979 * p[i]=1/0 depending on wether atomnumber i is a MM particle in the QM
980 * step or not. p[i+md->nr] is the shift of atomnumber i.
982 for (i = 0; i < 2*(md->nr); i++)
984 parallelMMarray[i] = 0;
987 for (i = 0; i < mm_nr; i++)
989 parallelMMarray[mm_j_particles[i].j] = 1;
990 parallelMMarray[mm_j_particles[i].j+(md->nr)] = mm_j_particles[i].shift;
992 gmx_sumi(md->nr, parallelMMarray, cr);
996 for (i = 0; i < md->nr; i++)
998 if (parallelMMarray[i])
1000 if (mm_nr >= mm_max)
1003 srenew(mm->indexMM, mm_max);
1004 srenew(mm->shiftMM, mm_max);
1006 mm->indexMM[mm_nr] = i;
1007 mm->shiftMM[mm_nr++] = parallelMMarray[i+md->nr]/parallelMMarray[i];
1010 mm->nrMMatoms = mm_nr;
1011 free(parallelMMarray);
1013 /* serial execution */
1016 mm->nrMMatoms = mm_nr;
1017 srenew(mm->shiftMM, mm_nr);
1018 srenew(mm->indexMM, mm_nr);
1019 for (i = 0; i < mm_nr; i++)
1021 mm->indexMM[i] = mm_j_particles[i].j;
1022 mm->shiftMM[i] = mm_j_particles[i].shift;
1026 /* (re) allocate memory for the MM coordiate array. The QM
1027 * coordinate array was already allocated in init_QMMM, and is
1028 * only (re)filled in the update_QMMM_coordinates routine
1030 srenew(mm->xMM, mm->nrMMatoms);
1031 /* now we (re) fill the array that contains the MM charges with
1032 * the forcefield charges. If requested, these charges will be
1033 * scaled by a factor
1035 srenew(mm->MMcharges, mm->nrMMatoms);
1036 for (i = 0; i < mm->nrMMatoms; i++) /* no free energy yet */
1038 mm->MMcharges[i] = md->chargeA[mm->indexMM[i]]*mm->scalefactor;
1040 if (qm->bTS || qm->bOPT)
1042 /* store (copy) the c6 and c12 parameters into the MMrec struct
1044 srenew(mm->c6, mm->nrMMatoms);
1045 srenew(mm->c12, mm->nrMMatoms);
1046 for (i = 0; i < mm->nrMMatoms; i++)
1048 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1049 mm->c6[i] = C6(fr->nbfp, top->idef.atnr, md->typeA[mm->indexMM[i]], md->typeA[mm->indexMM[i]])/c6au/6.0;
1050 mm->c12[i] = C12(fr->nbfp, top->idef.atnr, md->typeA[mm->indexMM[i]], md->typeA[mm->indexMM[i]])/c12au/12.0;
1052 punch_QMMM_excl(qr->qm[0], mm, &(top->excls));
1054 /* the next routine fills the coordinate fields in the QMMM rec of
1055 * both the qunatum atoms and the MM atoms, using the shifts
1059 update_QMMM_coord(x, fr, qr->qm[0], qr->mm);
1060 free(qm_i_particles);
1061 free(mm_j_particles);
1063 else /* ONIOM */ /* ????? */
1066 /* do for each layer */
1067 for (j = 0; j < qr->nrQMlayers; j++)
1070 qm->shiftQM[0] = XYZ2IS(0, 0, 0);
1071 for (i = 1; i < qm->nrQMatoms; i++)
1073 qm->shiftQM[i] = pbc_dx_aiuc(&pbc, x[qm->indexQM[0]], x[qm->indexQM[i]],
1076 update_QMMM_coord(x, fr, qm, mm);
1079 } /* update_QMMM_rec */
1082 real calculate_QMMM(t_commrec *cr,
1088 /* a selection for the QM package depending on which is requested
1089 * (Gaussian, GAMESS-UK, MOPAC or ORCA) needs to be implemented here. Now
1090 * it works through defines.... Not so nice yet
1099 *forces = NULL, *fshift = NULL,
1100 *forces2 = NULL, *fshift2 = NULL; /* needed for multilayer ONIOM */
1103 /* make a local copy the QMMMrec pointer
1108 /* now different procedures are carried out for one layer ONION and
1109 * normal QMMM on one hand and multilayer oniom on the other
1111 if (qr->QMMMscheme == eQMMMschemenormal || qr->nrQMlayers == 1)
1114 snew(forces, (qm->nrQMatoms+mm->nrMMatoms));
1115 snew(fshift, (qm->nrQMatoms+mm->nrMMatoms));
1116 QMener = call_QMroutine(cr, fr, qm, mm, forces, fshift);
1117 for (i = 0; i < qm->nrQMatoms; i++)
1119 for (j = 0; j < DIM; j++)
1121 f[qm->indexQM[i]][j] -= forces[i][j];
1122 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
1125 for (i = 0; i < mm->nrMMatoms; i++)
1127 for (j = 0; j < DIM; j++)
1129 f[mm->indexMM[i]][j] -= forces[qm->nrQMatoms+i][j];
1130 fr->fshift[mm->shiftMM[i]][j] += fshift[qm->nrQMatoms+i][j];
1137 else /* Multi-layer ONIOM */
1139 for (i = 0; i < qr->nrQMlayers-1; i++) /* last layer is special */
1142 qm2 = copy_QMrec(qr->qm[i+1]);
1144 qm2->nrQMatoms = qm->nrQMatoms;
1146 for (j = 0; j < qm2->nrQMatoms; j++)
1148 for (k = 0; k < DIM; k++)
1150 qm2->xQM[j][k] = qm->xQM[j][k];
1152 qm2->indexQM[j] = qm->indexQM[j];
1153 qm2->atomicnumberQM[j] = qm->atomicnumberQM[j];
1154 qm2->shiftQM[j] = qm->shiftQM[j];
1157 qm2->QMcharge = qm->QMcharge;
1158 /* this layer at the higher level of theory */
1159 srenew(forces, qm->nrQMatoms);
1160 srenew(fshift, qm->nrQMatoms);
1161 /* we need to re-initialize the QMroutine every step... */
1162 init_QMroutine(cr, qm, mm);
1163 QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
1165 /* this layer at the lower level of theory */
1166 srenew(forces2, qm->nrQMatoms);
1167 srenew(fshift2, qm->nrQMatoms);
1168 init_QMroutine(cr, qm2, mm);
1169 QMener -= call_QMroutine(cr, fr, qm2, mm, forces2, fshift2);
1170 /* E = E1high-E1low The next layer includes the current layer at
1171 * the lower level of theory, which provides + E2low
1172 * this is similar for gradients
1174 for (i = 0; i < qm->nrQMatoms; i++)
1176 for (j = 0; j < DIM; j++)
1178 f[qm->indexQM[i]][j] -= (forces[i][j]-forces2[i][j]);
1179 fr->fshift[qm->shiftQM[i]][j] += (fshift[i][j]-fshift2[i][j]);
1184 /* now the last layer still needs to be done: */
1185 qm = qr->qm[qr->nrQMlayers-1]; /* C counts from 0 */
1186 init_QMroutine(cr, qm, mm);
1187 srenew(forces, qm->nrQMatoms);
1188 srenew(fshift, qm->nrQMatoms);
1189 QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
1190 for (i = 0; i < qm->nrQMatoms; i++)
1192 for (j = 0; j < DIM; j++)
1194 f[qm->indexQM[i]][j] -= forces[i][j];
1195 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
1203 if (qm->bTS || qm->bOPT)
1205 /* qm[0] still contains the largest ONIOM QM subsystem
1206 * we take the optimized coordiates and put the in x[]
1208 for (i = 0; i < qm->nrQMatoms; i++)
1210 for (j = 0; j < DIM; j++)
1212 x[qm->indexQM[i]][j] = qm->xQM[i][j];
1217 } /* calculate_QMMM */
1219 /* end of QMMM core routines */