More const correctness
[gromacs.git] / src / gromacs / mdlib / qmmm.cpp
blob062c1ebf0213de8409c221f3fae0175c93412f7c
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 "qmmm.h"
41 #include "config.h"
43 #include <stdio.h>
44 #include <stdlib.h>
45 #include <string.h>
47 #include <cmath>
49 #include <algorithm>
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/fileio/confio.h"
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/gmxlib/nrnb.h"
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/force.h"
59 #include "gromacs/mdlib/ns.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/mdatom.h"
64 #include "gromacs/mdtypes/nblist.h"
65 #include "gromacs/pbcutil/ishift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/topology/mtop_lookup.h"
68 #include "gromacs/topology/mtop_util.h"
69 #include "gromacs/topology/topology.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/smalloc.h"
73 /* declarations of the interfaces to the QM packages. The _SH indicate
74 * the QM interfaces can be used for Surface Hopping simulations
76 #if GMX_QMMM_GAMESS
77 /* GAMESS interface */
79 void
80 init_gamess(const t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
82 real
83 call_gamess(const t_forcerec *fr,
84 t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
86 #elif GMX_QMMM_MOPAC
87 /* MOPAC interface */
89 void
90 init_mopac(t_QMrec *qm);
92 real
93 call_mopac(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
95 real
96 call_mopac_SH(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
98 #elif GMX_QMMM_GAUSSIAN
99 /* GAUSSIAN interface */
101 void
102 init_gaussian(t_QMrec *qm);
104 real
105 call_gaussian_SH(const t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
107 real
108 call_gaussian(const t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
110 #elif GMX_QMMM_ORCA
111 #include "gromacs/mdlib/qm_orca.h"
112 #endif
117 /* this struct and these comparison functions are needed for creating
118 * a QMMM input for the QM routines from the QMMM neighbor list.
121 typedef struct {
122 int j;
123 int shift;
124 } t_j_particle;
126 static int struct_comp(const void *a, const void *b)
129 return (int)(((t_j_particle *)a)->j)-(int)(((t_j_particle *)b)->j);
131 } /* struct_comp */
133 static real call_QMroutine(const t_commrec gmx_unused *cr, const t_forcerec gmx_unused *fr, t_QMrec gmx_unused *qm,
134 t_MMrec gmx_unused *mm, rvec gmx_unused f[], rvec gmx_unused fshift[])
136 /* makes a call to the requested QM routine (qm->QMmethod)
137 * Note that f is actually the gradient, i.e. -f
139 real
140 QMener = 0.0;
142 /* do a semi-empiprical calculation */
144 if (qm->QMmethod < eQMmethodRHF && !(mm->nrMMatoms))
146 #if GMX_QMMM_MOPAC
147 if (qm->bSH)
149 QMener = call_mopac_SH(qm, mm, f, fshift);
151 else
153 QMener = call_mopac(qm, mm, f, fshift);
155 #else
156 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
157 #endif
159 else
161 /* do an ab-initio calculation */
162 if (qm->bSH && qm->QMmethod == eQMmethodCASSCF)
164 #if GMX_QMMM_GAUSSIAN
165 QMener = call_gaussian_SH(fr, qm, mm, f, fshift);
166 #else
167 gmx_fatal(FARGS, "Ab-initio Surface-hopping only supported with Gaussian.");
168 #endif
170 else
172 #if GMX_QMMM_GAMESS
173 QMener = call_gamess(fr, qm, mm, f, fshift);
174 #elif GMX_QMMM_GAUSSIAN
175 QMener = call_gaussian(fr, qm, mm, f, fshift);
176 #elif GMX_QMMM_ORCA
177 QMener = call_orca(fr, qm, mm, f, fshift);
178 #else
179 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
180 #endif
183 return (QMener);
186 static void init_QMroutine(const t_commrec gmx_unused *cr, t_QMrec gmx_unused *qm, t_MMrec gmx_unused *mm)
188 /* makes a call to the requested QM routine (qm->QMmethod)
190 if (qm->QMmethod < eQMmethodRHF)
192 #if GMX_QMMM_MOPAC
193 /* do a semi-empiprical calculation */
194 init_mopac(qm);
195 #else
196 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
197 #endif
199 else
201 /* do an ab-initio calculation */
202 #if GMX_QMMM_GAMESS
203 init_gamess(cr, qm, mm);
204 #elif GMX_QMMM_GAUSSIAN
205 init_gaussian(qm);
206 #elif GMX_QMMM_ORCA
207 init_orca(qm);
208 #else
209 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
210 #endif
212 } /* init_QMroutine */
214 static void update_QMMM_coord(const rvec *x, const t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
216 /* shifts the QM and MM particles into the central box and stores
217 * these shifted coordinates in the coordinate arrays of the
218 * QMMMrec. These coordinates are passed on the QM subroutines.
223 /* shift the QM atoms into the central box
225 for (i = 0; i < qm->nrQMatoms; i++)
227 rvec_sub(x[qm->indexQM[i]], fr->shift_vec[qm->shiftQM[i]], qm->xQM[i]);
229 /* also shift the MM atoms into the central box, if any
231 for (i = 0; i < mm->nrMMatoms; i++)
233 rvec_sub(x[mm->indexMM[i]], fr->shift_vec[mm->shiftMM[i]], mm->xMM[i]);
235 } /* update_QMMM_coord */
237 /* end of QMMM subroutines */
239 /* QMMM core routines */
241 static t_QMrec *mk_QMrec(void)
243 t_QMrec *qm;
244 snew(qm, 1);
245 return qm;
246 } /* mk_QMrec */
248 static t_MMrec *mk_MMrec(void)
250 t_MMrec *mm;
251 snew(mm, 1);
252 return mm;
253 } /* mk_MMrec */
255 static void init_QMrec(int grpnr, t_QMrec *qm, int nr, int *atomarray,
256 gmx_mtop_t *mtop, t_inputrec *ir)
258 /* fills the t_QMrec struct of QM group grpnr
261 qm->nrQMatoms = nr;
262 snew(qm->xQM, nr);
263 snew(qm->indexQM, nr);
264 snew(qm->shiftQM, nr); /* the shifts */
265 for (int i = 0; i < nr; i++)
267 qm->indexQM[i] = atomarray[i];
270 snew(qm->atomicnumberQM, nr);
271 int molb = 0;
272 for (int i = 0; i < qm->nrQMatoms; i++)
274 const t_atom &atom = mtopGetAtomParameters(mtop, qm->indexQM[i], &molb);
275 qm->nelectrons += mtop->atomtypes.atomnumber[atom.type];
276 qm->atomicnumberQM[i] = mtop->atomtypes.atomnumber[atom.type];
279 qm->QMcharge = ir->opts.QMcharge[grpnr];
280 qm->multiplicity = ir->opts.QMmult[grpnr];
281 qm->nelectrons -= ir->opts.QMcharge[grpnr];
283 qm->QMmethod = ir->opts.QMmethod[grpnr];
284 qm->QMbasis = ir->opts.QMbasis[grpnr];
285 /* trajectory surface hopping setup (Gaussian only) */
286 qm->bSH = ir->opts.bSH[grpnr];
287 qm->CASorbitals = ir->opts.CASorbitals[grpnr];
288 qm->CASelectrons = ir->opts.CASelectrons[grpnr];
289 qm->SAsteps = ir->opts.SAsteps[grpnr];
290 qm->SAon = ir->opts.SAon[grpnr];
291 qm->SAoff = ir->opts.SAoff[grpnr];
292 /* hack to prevent gaussian from reinitializing all the time */
293 qm->nQMcpus = 0; /* number of CPU's to be used by g01, is set
294 * upon initializing gaussian
295 * (init_gaussian()
297 /* print the current layer to allow users to check their input */
298 fprintf(stderr, "Layer %d\nnr of QM atoms %d\n", grpnr, nr);
299 fprintf(stderr, "QMlevel: %s/%s\n\n",
300 eQMmethod_names[qm->QMmethod], eQMbasis_names[qm->QMbasis]);
301 } /* init_QMrec */
303 static t_QMrec *copy_QMrec(t_QMrec *qm)
305 /* copies the contents of qm into a new t_QMrec struct */
306 t_QMrec
307 *qmcopy;
311 qmcopy = mk_QMrec();
312 qmcopy->nrQMatoms = qm->nrQMatoms;
313 snew(qmcopy->xQM, qmcopy->nrQMatoms);
314 snew(qmcopy->indexQM, qmcopy->nrQMatoms);
315 snew(qmcopy->atomicnumberQM, qm->nrQMatoms);
316 snew(qmcopy->shiftQM, qmcopy->nrQMatoms); /* the shifts */
317 for (i = 0; i < qmcopy->nrQMatoms; i++)
319 qmcopy->shiftQM[i] = qm->shiftQM[i];
320 qmcopy->indexQM[i] = qm->indexQM[i];
321 qmcopy->atomicnumberQM[i] = qm->atomicnumberQM[i];
323 qmcopy->nelectrons = qm->nelectrons;
324 qmcopy->multiplicity = qm->multiplicity;
325 qmcopy->QMcharge = qm->QMcharge;
326 qmcopy->nelectrons = qm->nelectrons;
327 qmcopy->QMmethod = qm->QMmethod;
328 qmcopy->QMbasis = qm->QMbasis;
329 /* trajectory surface hopping setup (Gaussian only) */
330 qmcopy->bSH = qm->bSH;
331 qmcopy->CASorbitals = qm->CASorbitals;
332 qmcopy->CASelectrons = qm->CASelectrons;
333 qmcopy->SAsteps = qm->SAsteps;
334 qmcopy->SAon = qm->SAon;
335 qmcopy->SAoff = qm->SAoff;
337 /* Gaussian init. variables */
338 qmcopy->nQMcpus = qm->nQMcpus;
339 for (i = 0; i < DIM; i++)
341 qmcopy->SHbasis[i] = qm->SHbasis[i];
343 qmcopy->QMmem = qm->QMmem;
344 qmcopy->accuracy = qm->accuracy;
345 qmcopy->cpmcscf = qm->cpmcscf;
346 qmcopy->SAstep = qm->SAstep;
348 return(qmcopy);
350 } /*copy_QMrec */
352 t_QMMMrec *mk_QMMMrec(void)
355 t_QMMMrec *qr;
357 snew(qr, 1);
359 return qr;
361 } /* mk_QMMMrec */
363 void init_QMMMrec(const t_commrec *cr,
364 gmx_mtop_t *mtop,
365 t_inputrec *ir,
366 const t_forcerec *fr)
368 /* we put the atomsnumbers of atoms that belong to the QMMM group in
369 * an array that will be copied later to QMMMrec->indexQM[..]. Also
370 * it will be used to create an QMMMrec->bQMMM index array that
371 * simply contains true/false for QM and MM (the other) atoms.
374 gmx_groups_t *groups;
375 int *qm_arr = nullptr, vsite, ai, aj;
376 int qm_max = 0, qm_nr = 0, i, j, jmax, k, l, nrvsite2 = 0;
377 t_QMMMrec *qr;
378 t_MMrec *mm;
379 t_iatom *iatoms;
380 gmx_mtop_atomloop_all_t aloop;
381 gmx_mtop_ilistloop_all_t iloop;
382 int a_offset;
383 const t_ilist *ilist_mol;
385 if (ir->cutoff_scheme != ecutsGROUP)
387 gmx_fatal(FARGS, "QMMM is currently only supported with cutoff-scheme=group");
389 if (!EI_DYNAMICS(ir->eI))
391 gmx_fatal(FARGS, "QMMM is only supported with dynamics");
394 /* issue a fatal if the user wants to run with more than one node */
395 if (PAR(cr))
397 gmx_fatal(FARGS, "QM/MM does not work in parallel, use a single rank instead\n");
400 /* Make a local copy of the QMMMrec */
401 qr = fr->qr;
403 /* bQMMM[..] is an array containing TRUE/FALSE for atoms that are
404 * QM/not QM. We first set all elemenst at false. Afterwards we use
405 * the qm_arr (=MMrec->indexQM) to changes the elements
406 * corresponding to the QM atoms at TRUE. */
408 qr->QMMMscheme = ir->QMMMscheme;
410 /* we take the possibility into account that a user has
411 * defined more than one QM group:
413 /* an ugly work-around in case there is only one group In this case
414 * the whole system is treated as QM. Otherwise the second group is
415 * always the rest of the total system and is treated as MM.
418 /* small problem if there is only QM.... so no MM */
420 jmax = ir->opts.ngQM;
422 if (qr->QMMMscheme == eQMMMschemeoniom)
424 qr->nrQMlayers = jmax;
426 else
428 qr->nrQMlayers = 1;
431 groups = &mtop->groups;
433 /* there are jmax groups of QM atoms. In case of multiple QM groups
434 * I assume that the users wants to do ONIOM. However, maybe it
435 * should also be possible to define more than one QM subsystem with
436 * independent neighbourlists. I have to think about
437 * that.. 11-11-2003
439 snew(qr->qm, jmax);
440 for (j = 0; j < jmax; j++)
442 /* new layer */
443 aloop = gmx_mtop_atomloop_all_init(mtop);
444 const t_atom *atom;
445 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
447 if (qm_nr >= qm_max)
449 qm_max += 1000;
450 srenew(qm_arr, qm_max);
452 if (ggrpnr(groups, egcQMMM, i) == j)
454 /* hack for tip4p */
455 qm_arr[qm_nr++] = i;
458 if (qr->QMMMscheme == eQMMMschemeoniom)
460 /* add the atoms to the bQMMM array
463 /* I assume that users specify the QM groups from small to
464 * big(ger) in the mdp file
466 qr->qm[j] = mk_QMrec();
467 /* we need to throw out link atoms that in the previous layer
468 * existed to separate this QMlayer from the previous
469 * QMlayer. We use the iatoms array in the idef for that
470 * purpose. If all atoms defining the current Link Atom (Dummy2)
471 * are part of the current QM layer it needs to be removed from
472 * qm_arr[]. */
474 iloop = gmx_mtop_ilistloop_all_init(mtop);
475 while (gmx_mtop_ilistloop_all_next(iloop, &ilist_mol, &a_offset))
477 nrvsite2 = ilist_mol[F_VSITE2].nr;
478 iatoms = ilist_mol[F_VSITE2].iatoms;
480 for (k = 0; k < nrvsite2; k += 4)
482 vsite = a_offset + iatoms[k+1]; /* the vsite */
483 ai = a_offset + iatoms[k+2]; /* constructing atom */
484 aj = a_offset + iatoms[k+3]; /* constructing atom */
485 if (ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, ai)
487 ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, aj))
489 /* this dummy link atom needs to be removed from the qm_arr
490 * before making the QMrec of this layer!
492 for (i = 0; i < qm_nr; i++)
494 if (qm_arr[i] == vsite)
496 /* drop the element */
497 for (l = i; l < qm_nr; l++)
499 qm_arr[l] = qm_arr[l+1];
501 qm_nr--;
508 /* store QM atoms in this layer in the QMrec and initialise layer
510 init_QMrec(j, qr->qm[j], qm_nr, qm_arr, mtop, ir);
513 if (qr->QMMMscheme != eQMMMschemeoniom)
516 /* standard QMMM, all layers are merged together so there is one QM
517 * subsystem and one MM subsystem.
518 * Also we set the charges to zero in mtop to prevent the innerloops
519 * from doubly counting the electostatic QM MM interaction
520 * TODO: Consider doing this in grompp instead.
523 int molb = 0;
524 for (k = 0; k < qm_nr; k++)
526 int indexInMolecule;
527 mtopGetMolblockIndex(mtop, qm_arr[k], &molb, nullptr, &indexInMolecule);
528 t_atom *atom = &mtop->moltype[mtop->molblock[molb].type].atoms.atom[indexInMolecule];
529 atom->q = 0.0;
530 atom->qB = 0.0;
532 qr->qm[0] = mk_QMrec();
533 /* store QM atoms in the QMrec and initialise
535 init_QMrec(0, qr->qm[0], qm_nr, qm_arr, mtop, ir);
537 /* MM rec creation */
538 mm = mk_MMrec();
539 mm->scalefactor = ir->scalefactor;
540 mm->nrMMatoms = (mtop->natoms)-(qr->qm[0]->nrQMatoms); /* rest of the atoms */
541 qr->mm = mm;
543 else /* ONIOM */
544 { /* MM rec creation */
545 mm = mk_MMrec();
546 mm->scalefactor = ir->scalefactor;
547 mm->nrMMatoms = 0;
548 qr->mm = mm;
551 /* these variables get updated in the update QMMMrec */
553 if (qr->nrQMlayers == 1)
555 /* with only one layer there is only one initialisation
556 * needed. Multilayer is a bit more complicated as it requires
557 * re-initialisation at every step of the simulation. This is due
558 * to the use of COMMON blocks in the fortran QM subroutines.
560 if (qr->qm[0]->QMmethod < eQMmethodRHF)
562 #if GMX_QMMM_MOPAC
563 /* semi-empiprical 1-layer ONIOM calculation requested (mopac93) */
564 init_mopac(qr->qm[0]);
565 #else
566 gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
567 #endif
569 else
571 /* ab initio calculation requested (gamess/gaussian/ORCA) */
572 #if GMX_QMMM_GAMESS
573 init_gamess(cr, qr->qm[0], qr->mm);
574 #elif GMX_QMMM_GAUSSIAN
575 init_gaussian(qr->qm[0]);
576 #elif GMX_QMMM_ORCA
577 init_orca(qr->qm[0]);
578 #else
579 gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
580 #endif
583 } /* init_QMMMrec */
585 void update_QMMMrec(const t_commrec *cr,
586 const t_forcerec *fr,
587 const rvec *x,
588 const t_mdatoms *md,
589 const matrix box)
591 /* updates the coordinates of both QM atoms and MM atoms and stores
592 * them in the QMMMrec.
594 * NOTE: is NOT yet working if there are no PBC. Also in ns.c, simple
595 * ns needs to be fixed!
598 mm_max = 0, mm_nr = 0, mm_nr_new, i, j, is, k, shift;
599 t_j_particle
600 *mm_j_particles = nullptr, *qm_i_particles = nullptr;
601 t_QMMMrec
602 *qr;
603 t_nblist
604 *QMMMlist;
605 rvec
606 dx, crd;
607 t_QMrec
608 *qm;
609 t_MMrec
610 *mm;
611 t_pbc
612 pbc;
614 *parallelMMarray = nullptr;
616 /* every cpu has this array. On every processor we fill this array
617 * with 1's and 0's. 1's indicate the atoms is a QM atom on the
618 * current cpu in a later stage these arrays are all summed. indexes
619 * > 0 indicate the atom is a QM atom. Every node therefore knows
620 * whcih atoms are part of the QM subsystem.
622 /* copy some pointers */
623 qr = fr->qr;
624 mm = qr->mm;
625 QMMMlist = fr->QMMMlist;
627 /* init_pbc(box); needs to be called first, see pbc.h */
628 ivec null_ivec;
629 clear_ivec(null_ivec);
630 set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : null_ivec,
631 FALSE, box);
632 /* only in standard (normal) QMMM we need the neighbouring MM
633 * particles to provide a electric field of point charges for the QM
634 * atoms.
636 if (qr->QMMMscheme == eQMMMschemenormal) /* also implies 1 QM-layer */
638 /* we NOW create/update a number of QMMMrec entries:
640 * 1) the shiftQM, containing the shifts of the QM atoms
642 * 2) the indexMM array, containing the index of the MM atoms
644 * 3) the shiftMM, containing the shifts of the MM atoms
646 * 4) the shifted coordinates of the MM atoms
648 * the shifts are used for computing virial of the QM/MM particles.
650 qm = qr->qm[0]; /* in case of normal QMMM, there is only one group */
651 snew(qm_i_particles, QMMMlist->nri);
652 if (QMMMlist->nri)
654 qm_i_particles[0].shift = XYZ2IS(0, 0, 0);
655 for (i = 0; i < QMMMlist->nri; i++)
657 qm_i_particles[i].j = QMMMlist->iinr[i];
659 if (i)
661 qm_i_particles[i].shift = pbc_dx_aiuc(&pbc, x[QMMMlist->iinr[0]],
662 x[QMMMlist->iinr[i]], dx);
665 /* However, since nri >= nrQMatoms, we do a quicksort, and throw
666 * out double, triple, etc. entries later, as we do for the MM
667 * list too.
670 /* compute the shift for the MM j-particles with respect to
671 * the QM i-particle and store them.
674 crd[0] = IS2X(QMMMlist->shift[i]) + IS2X(qm_i_particles[i].shift);
675 crd[1] = IS2Y(QMMMlist->shift[i]) + IS2Y(qm_i_particles[i].shift);
676 crd[2] = IS2Z(QMMMlist->shift[i]) + IS2Z(qm_i_particles[i].shift);
677 is = static_cast<int>(XYZ2IS(crd[0], crd[1], crd[2]));
678 for (j = QMMMlist->jindex[i];
679 j < QMMMlist->jindex[i+1];
680 j++)
682 if (mm_nr >= mm_max)
684 mm_max += 1000;
685 srenew(mm_j_particles, mm_max);
688 mm_j_particles[mm_nr].j = QMMMlist->jjnr[j];
689 mm_j_particles[mm_nr].shift = is;
690 mm_nr++;
694 /* quicksort QM and MM shift arrays and throw away multiple entries */
698 qsort(qm_i_particles, QMMMlist->nri,
699 (size_t)sizeof(qm_i_particles[0]),
700 struct_comp);
701 /* The mm_j_particles argument to qsort is not allowed to be NULL */
702 if (mm_nr > 0)
704 qsort(mm_j_particles, mm_nr,
705 (size_t)sizeof(mm_j_particles[0]),
706 struct_comp);
708 /* remove multiples in the QM shift array, since in init_QMMM() we
709 * went through the atom numbers from 0 to md.nr, the order sorted
710 * here matches the one of QMindex already.
712 j = 0;
713 for (i = 0; i < QMMMlist->nri; i++)
715 if (i == 0 || qm_i_particles[i].j != qm_i_particles[i-1].j)
717 qm_i_particles[j++] = qm_i_particles[i];
720 mm_nr_new = 0;
721 /* Remove double entries for the MM array.
722 * Also remove mm atoms that have no charges!
723 * actually this is already done in the ns.c
725 for (i = 0; i < mm_nr; i++)
727 if ((i == 0 || mm_j_particles[i].j != mm_j_particles[i-1].j)
728 && !md->bQM[mm_j_particles[i].j]
729 && (md->chargeA[mm_j_particles[i].j]
730 || (md->chargeB && md->chargeB[mm_j_particles[i].j])))
732 mm_j_particles[mm_nr_new++] = mm_j_particles[i];
735 mm_nr = mm_nr_new;
736 /* store the data retrieved above into the QMMMrec
738 k = 0;
739 /* Keep the compiler happy,
740 * shift will always be set in the loop for i=0
742 shift = 0;
743 for (i = 0; i < qm->nrQMatoms; i++)
745 /* not all qm particles might have appeared as i
746 * particles. They might have been part of the same charge
747 * group for instance.
749 if (qm->indexQM[i] == qm_i_particles[k].j)
751 shift = qm_i_particles[k++].shift;
753 /* use previous shift, assuming they belong the same charge
754 * group anyway,
757 qm->shiftQM[i] = shift;
760 /* parallel excecution */
761 if (PAR(cr))
763 snew(parallelMMarray, 2*(md->nr));
764 /* only MM particles have a 1 at their atomnumber. The second part
765 * of the array contains the shifts. Thus:
766 * p[i]=1/0 depending on wether atomnumber i is a MM particle in the QM
767 * step or not. p[i+md->nr] is the shift of atomnumber i.
769 for (i = 0; i < 2*(md->nr); i++)
771 parallelMMarray[i] = 0;
774 for (i = 0; i < mm_nr; i++)
776 parallelMMarray[mm_j_particles[i].j] = 1;
777 parallelMMarray[mm_j_particles[i].j+(md->nr)] = mm_j_particles[i].shift;
779 gmx_sumi(md->nr, parallelMMarray, cr);
780 mm_nr = 0;
782 mm_max = 0;
783 for (i = 0; i < md->nr; i++)
785 if (parallelMMarray[i])
787 if (mm_nr >= mm_max)
789 mm_max += 1000;
790 srenew(mm->indexMM, mm_max);
791 srenew(mm->shiftMM, mm_max);
793 mm->indexMM[mm_nr] = i;
794 mm->shiftMM[mm_nr++] = parallelMMarray[i+md->nr]/parallelMMarray[i];
797 mm->nrMMatoms = mm_nr;
798 free(parallelMMarray);
800 /* serial execution */
801 else
803 mm->nrMMatoms = mm_nr;
804 srenew(mm->shiftMM, mm_nr);
805 srenew(mm->indexMM, mm_nr);
806 for (i = 0; i < mm_nr; i++)
808 mm->indexMM[i] = mm_j_particles[i].j;
809 mm->shiftMM[i] = mm_j_particles[i].shift;
813 /* (re) allocate memory for the MM coordiate array. The QM
814 * coordinate array was already allocated in init_QMMM, and is
815 * only (re)filled in the update_QMMM_coordinates routine
817 srenew(mm->xMM, mm->nrMMatoms);
818 /* now we (re) fill the array that contains the MM charges with
819 * the forcefield charges. If requested, these charges will be
820 * scaled by a factor
822 srenew(mm->MMcharges, mm->nrMMatoms);
823 for (i = 0; i < mm->nrMMatoms; i++) /* no free energy yet */
825 mm->MMcharges[i] = md->chargeA[mm->indexMM[i]]*mm->scalefactor;
827 /* the next routine fills the coordinate fields in the QMMM rec of
828 * both the qunatum atoms and the MM atoms, using the shifts
829 * calculated above.
832 update_QMMM_coord(x, fr, qr->qm[0], qr->mm);
833 free(qm_i_particles);
834 free(mm_j_particles);
836 else /* ONIOM */ /* ????? */
838 mm->nrMMatoms = 0;
839 /* do for each layer */
840 for (j = 0; j < qr->nrQMlayers; j++)
842 qm = qr->qm[j];
843 qm->shiftQM[0] = XYZ2IS(0, 0, 0);
844 for (i = 1; i < qm->nrQMatoms; i++)
846 qm->shiftQM[i] = pbc_dx_aiuc(&pbc, x[qm->indexQM[0]], x[qm->indexQM[i]],
847 dx);
849 update_QMMM_coord(x, fr, qm, mm);
852 } /* update_QMMM_rec */
855 real calculate_QMMM(const t_commrec *cr,
856 rvec f[],
857 const t_forcerec *fr)
859 real
860 QMener = 0.0;
861 /* a selection for the QM package depending on which is requested
862 * (Gaussian, GAMESS-UK, MOPAC or ORCA) needs to be implemented here. Now
863 * it works through defines.... Not so nice yet
865 t_QMMMrec
866 *qr;
867 t_QMrec
868 *qm, *qm2;
869 t_MMrec
870 *mm = nullptr;
871 rvec
872 *forces = nullptr, *fshift = nullptr,
873 *forces2 = nullptr, *fshift2 = nullptr; /* needed for multilayer ONIOM */
875 i, j, k;
876 /* make a local copy the QMMMrec pointer
878 qr = fr->qr;
879 mm = qr->mm;
881 /* now different procedures are carried out for one layer ONION and
882 * normal QMMM on one hand and multilayer oniom on the other
884 if (qr->QMMMscheme == eQMMMschemenormal || qr->nrQMlayers == 1)
886 qm = qr->qm[0];
887 snew(forces, (qm->nrQMatoms+mm->nrMMatoms));
888 snew(fshift, (qm->nrQMatoms+mm->nrMMatoms));
889 QMener = call_QMroutine(cr, fr, qm, mm, forces, fshift);
890 for (i = 0; i < qm->nrQMatoms; i++)
892 for (j = 0; j < DIM; j++)
894 f[qm->indexQM[i]][j] -= forces[i][j];
895 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
898 for (i = 0; i < mm->nrMMatoms; i++)
900 for (j = 0; j < DIM; j++)
902 f[mm->indexMM[i]][j] -= forces[qm->nrQMatoms+i][j];
903 fr->fshift[mm->shiftMM[i]][j] += fshift[qm->nrQMatoms+i][j];
907 free(forces);
908 free(fshift);
910 else /* Multi-layer ONIOM */
912 for (i = 0; i < qr->nrQMlayers-1; i++) /* last layer is special */
914 qm = qr->qm[i];
915 qm2 = copy_QMrec(qr->qm[i+1]);
917 qm2->nrQMatoms = qm->nrQMatoms;
919 for (j = 0; j < qm2->nrQMatoms; j++)
921 for (k = 0; k < DIM; k++)
923 qm2->xQM[j][k] = qm->xQM[j][k];
925 qm2->indexQM[j] = qm->indexQM[j];
926 qm2->atomicnumberQM[j] = qm->atomicnumberQM[j];
927 qm2->shiftQM[j] = qm->shiftQM[j];
930 qm2->QMcharge = qm->QMcharge;
931 /* this layer at the higher level of theory */
932 srenew(forces, qm->nrQMatoms);
933 srenew(fshift, qm->nrQMatoms);
934 /* we need to re-initialize the QMroutine every step... */
935 init_QMroutine(cr, qm, mm);
936 QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
938 /* this layer at the lower level of theory */
939 srenew(forces2, qm->nrQMatoms);
940 srenew(fshift2, qm->nrQMatoms);
941 init_QMroutine(cr, qm2, mm);
942 QMener -= call_QMroutine(cr, fr, qm2, mm, forces2, fshift2);
943 /* E = E1high-E1low The next layer includes the current layer at
944 * the lower level of theory, which provides + E2low
945 * this is similar for gradients
947 for (i = 0; i < qm->nrQMatoms; i++)
949 for (j = 0; j < DIM; j++)
951 f[qm->indexQM[i]][j] -= (forces[i][j]-forces2[i][j]);
952 fr->fshift[qm->shiftQM[i]][j] += (fshift[i][j]-fshift2[i][j]);
955 free(qm2);
957 /* now the last layer still needs to be done: */
958 qm = qr->qm[qr->nrQMlayers-1]; /* C counts from 0 */
959 init_QMroutine(cr, qm, mm);
960 srenew(forces, qm->nrQMatoms);
961 srenew(fshift, qm->nrQMatoms);
962 QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
963 for (i = 0; i < qm->nrQMatoms; i++)
965 for (j = 0; j < DIM; j++)
967 f[qm->indexQM[i]][j] -= forces[i][j];
968 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
971 free(forces);
972 free(fshift);
973 free(forces2);
974 free(fshift2);
976 return(QMener);
977 } /* calculate_QMMM */
979 /* end of QMMM core routines */