Merge branch release-5-1
[alexxy/gromacs.git] / src / gromacs / mdlib / qmmm.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "qmmm.h"
40
41 #include "config.h"
42
43 #include <math.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <string.h>
47
48 #include <cmath>
49
50 #include <algorithm>
51
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/fileio/confio.h"
54 #include "gromacs/gmxlib/network.h"
55 #include "gromacs/gmxlib/nrnb.h"
56 #include "gromacs/math/functions.h"
57 #include "gromacs/math/units.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/mdlib/force.h"
60 #include "gromacs/mdlib/ns.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/mdtypes/inputrec.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/mdtypes/mdatom.h"
65 #include "gromacs/mdtypes/nblist.h"
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/pbcutil/pbc.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"
72
73 /* declarations of the interfaces to the QM packages. The _SH indicate
74  * the QM interfaces can be used for Surface Hopping simulations
75  */
76 #if GMX_QMMM_GAMESS
77 /* GAMESS interface */
78
79 void
80 init_gamess(t_commrec *cr, t_QMrec *qm, t_MMrec *mm);
81
82 real
83 call_gamess(t_forcerec *fr,
84             t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
85
86 #elif GMX_QMMM_MOPAC
87 /* MOPAC interface */
88
89 void
90 init_mopac(t_QMrec *qm);
91
92 real
93 call_mopac(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
94
95 real
96 call_mopac_SH(t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
97
98 #elif GMX_QMMM_GAUSSIAN
99 /* GAUSSIAN interface */
100
101 void
102 init_gaussian(t_QMrec *qm);
103
104 real
105 call_gaussian_SH(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
106
107 real
108 call_gaussian(t_forcerec *fr, t_QMrec *qm, t_MMrec *mm, rvec f[], rvec fshift[]);
109
110 #elif GMX_QMMM_ORCA
111 /* ORCA interface */
112
113 void
114 init_orca(t_QMrec *qm);
115
116 real
117 call_orca(t_forcerec *fr, t_QMrec *qm,
118           t_MMrec *mm, rvec f[], rvec fshift[]);
119
120 #endif
121
122
123
124
125 /* this struct and these comparison functions are needed for creating
126  * a QMMM input for the QM routines from the QMMM neighbor list.
127  */
128
129 typedef struct {
130     int      j;
131     int      shift;
132 } t_j_particle;
133
134 static int struct_comp(const void *a, const void *b)
135 {
136
137     return (int)(((t_j_particle *)a)->j)-(int)(((t_j_particle *)b)->j);
138
139 } /* struct_comp */
140
141 real call_QMroutine(t_commrec gmx_unused *cr, t_forcerec gmx_unused *fr, t_QMrec gmx_unused *qm,
142                     t_MMrec gmx_unused *mm, rvec gmx_unused f[], rvec gmx_unused fshift[])
143 {
144     /* makes a call to the requested QM routine (qm->QMmethod)
145      * Note that f is actually the gradient, i.e. -f
146      */
147     real
148         QMener = 0.0;
149
150     /* do a semi-empiprical calculation */
151
152     if (qm->QMmethod < eQMmethodRHF && !(mm->nrMMatoms))
153     {
154 #if GMX_QMMM_MOPAC
155         if (qm->bSH)
156         {
157             QMener = call_mopac_SH(qm, mm, f, fshift);
158         }
159         else
160         {
161             QMener = call_mopac(qm, mm, f, fshift);
162         }
163 #else
164         gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
165 #endif
166     }
167     else
168     {
169         /* do an ab-initio calculation */
170         if (qm->bSH && qm->QMmethod == eQMmethodCASSCF)
171         {
172 #if GMX_QMMM_GAUSSIAN
173             QMener = call_gaussian_SH(fr, qm, mm, f, fshift);
174 #else
175             gmx_fatal(FARGS, "Ab-initio Surface-hopping only supported with Gaussian.");
176 #endif
177         }
178         else
179         {
180 #if GMX_QMMM_GAMESS
181             QMener = call_gamess(fr, qm, mm, f, fshift);
182 #elif GMX_QMMM_GAUSSIAN
183             QMener = call_gaussian(fr, qm, mm, f, fshift);
184 #elif GMX_QMMM_ORCA
185             QMener = call_orca(fr, qm, mm, f, fshift);
186 #else
187             gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
188 #endif
189         }
190     }
191     return (QMener);
192 }
193
194 void init_QMroutine(t_commrec gmx_unused *cr, t_QMrec gmx_unused *qm, t_MMrec gmx_unused *mm)
195 {
196     /* makes a call to the requested QM routine (qm->QMmethod)
197      */
198     if (qm->QMmethod < eQMmethodRHF)
199     {
200 #if GMX_QMMM_MOPAC
201         /* do a semi-empiprical calculation */
202         init_mopac(qm);
203 #else
204         gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
205 #endif
206     }
207     else
208     {
209         /* do an ab-initio calculation */
210 #if GMX_QMMM_GAMESS
211         init_gamess(cr, qm, mm);
212 #elif GMX_QMMM_GAUSSIAN
213         init_gaussian(qm);
214 #elif GMX_QMMM_ORCA
215         init_orca(qm);
216 #else
217         gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
218 #endif
219     }
220 } /* init_QMroutine */
221
222 void update_QMMM_coord(rvec x[], t_forcerec *fr, t_QMrec *qm, t_MMrec *mm)
223 {
224     /* shifts the QM and MM particles into the central box and stores
225      * these shifted coordinates in the coordinate arrays of the
226      * QMMMrec. These coordinates are passed on the QM subroutines.
227      */
228     int
229         i;
230
231     /* shift the QM atoms into the central box
232      */
233     for (i = 0; i < qm->nrQMatoms; i++)
234     {
235         rvec_sub(x[qm->indexQM[i]], fr->shift_vec[qm->shiftQM[i]], qm->xQM[i]);
236     }
237     /* also shift the MM atoms into the central box, if any
238      */
239     for (i = 0; i < mm->nrMMatoms; i++)
240     {
241         rvec_sub(x[mm->indexMM[i]], fr->shift_vec[mm->shiftMM[i]], mm->xMM[i]);
242     }
243 } /* update_QMMM_coord */
244
245 static void punch_QMMM_excl(t_QMrec *qm, t_MMrec *mm, t_blocka *excls)
246 {
247     /* punch a file containing the bonded interactions of each QM
248      * atom with MM atoms. These need to be excluded in the QM routines
249      * Only needed in case of QM/MM optimizations
250      */
251     FILE
252        *out = NULL;
253     int
254         i, j, k, nrexcl = 0, *excluded = NULL, max_excl = 0;
255
256
257     out = fopen("QMMMexcl.dat", "w");
258
259     /* this can be done more efficiently I think
260      */
261     for (i = 0; i < qm->nrQMatoms; i++)
262     {
263         nrexcl = 0;
264         for (j = excls->index[qm->indexQM[i]];
265              j < excls->index[qm->indexQM[i]+1];
266              j++)
267         {
268             for (k = 0; k < mm->nrMMatoms; k++)
269             {
270                 if (mm->indexMM[k] == excls->a[j]) /* the excluded MM atom */
271                 {
272                     if (nrexcl >= max_excl)
273                     {
274                         max_excl += 1000;
275                         srenew(excluded, max_excl);
276                     }
277                     excluded[nrexcl++] = k;
278                     continue;
279                 }
280             }
281         }
282         /* write to file: */
283         fprintf(out, "%5d %5d\n", i+1, nrexcl);
284         for (j = 0; j < nrexcl; j++)
285         {
286             fprintf(out, "%5d ", excluded[j]);
287         }
288         fprintf(out, "\n");
289     }
290     free(excluded);
291     fclose(out);
292 } /* punch_QMMM_excl */
293
294
295 /* end of QMMM subroutines */
296
297 /* QMMM core routines */
298
299 t_QMrec *mk_QMrec(void)
300 {
301     t_QMrec *qm;
302     snew(qm, 1);
303     return qm;
304 } /* mk_QMrec */
305
306 t_MMrec *mk_MMrec(void)
307 {
308     t_MMrec *mm;
309     snew(mm, 1);
310     return mm;
311 } /* mk_MMrec */
312
313 static void init_QMrec(int grpnr, t_QMrec *qm, int nr, int *atomarray,
314                        gmx_mtop_t *mtop, t_inputrec *ir)
315 {
316     /* fills the t_QMrec struct of QM group grpnr
317      */
318     int                   i;
319     gmx_mtop_atomlookup_t alook;
320     t_atom               *atom;
321
322
323     qm->nrQMatoms = nr;
324     snew(qm->xQM, nr);
325     snew(qm->indexQM, nr);
326     snew(qm->shiftQM, nr); /* the shifts */
327     for (i = 0; i < nr; i++)
328     {
329         qm->indexQM[i] = atomarray[i];
330     }
331
332     alook = gmx_mtop_atomlookup_init(mtop);
333
334     snew(qm->atomicnumberQM, nr);
335     for (i = 0; i < qm->nrQMatoms; i++)
336     {
337         gmx_mtop_atomnr_to_atom(alook, qm->indexQM[i], &atom);
338         qm->nelectrons       += mtop->atomtypes.atomnumber[atom->type];
339         qm->atomicnumberQM[i] = mtop->atomtypes.atomnumber[atom->type];
340     }
341
342     gmx_mtop_atomlookup_destroy(alook);
343
344     qm->QMcharge       = ir->opts.QMcharge[grpnr];
345     qm->multiplicity   = ir->opts.QMmult[grpnr];
346     qm->nelectrons    -= ir->opts.QMcharge[grpnr];
347
348     qm->QMmethod       = ir->opts.QMmethod[grpnr];
349     qm->QMbasis        = ir->opts.QMbasis[grpnr];
350     /* trajectory surface hopping setup (Gaussian only) */
351     qm->bSH            = ir->opts.bSH[grpnr];
352     qm->CASorbitals    = ir->opts.CASorbitals[grpnr];
353     qm->CASelectrons   = ir->opts.CASelectrons[grpnr];
354     qm->SAsteps        = ir->opts.SAsteps[grpnr];
355     qm->SAon           = ir->opts.SAon[grpnr];
356     qm->SAoff          = ir->opts.SAoff[grpnr];
357     /* hack to prevent gaussian from reinitializing all the time */
358     qm->nQMcpus        = 0; /* number of CPU's to be used by g01, is set
359                              * upon initializing gaussian
360                              * (init_gaussian()
361                              */
362     /* print the current layer to allow users to check their input */
363     fprintf(stderr, "Layer %d\nnr of QM atoms %d\n", grpnr, nr);
364     fprintf(stderr, "QMlevel: %s/%s\n\n",
365             eQMmethod_names[qm->QMmethod], eQMbasis_names[qm->QMbasis]);
366
367     /* frontier atoms */
368     snew(qm->frontatoms, nr);
369     /* Lennard-Jones coefficients */
370     snew(qm->c6, nr);
371     snew(qm->c12, nr);
372     /* do we optimize the QM separately using the algorithms of the QM program??
373      */
374     qm->bTS      = ir->opts.bTS[grpnr];
375     qm->bOPT     = ir->opts.bOPT[grpnr];
376
377 } /* init_QMrec */
378
379 t_QMrec *copy_QMrec(t_QMrec *qm)
380 {
381     /* copies the contents of qm into a new t_QMrec struct */
382     t_QMrec
383        *qmcopy;
384     int
385         i;
386
387     qmcopy            = mk_QMrec();
388     qmcopy->nrQMatoms = qm->nrQMatoms;
389     snew(qmcopy->xQM, qmcopy->nrQMatoms);
390     snew(qmcopy->indexQM, qmcopy->nrQMatoms);
391     snew(qmcopy->atomicnumberQM, qm->nrQMatoms);
392     snew(qmcopy->shiftQM, qmcopy->nrQMatoms); /* the shifts */
393     for (i = 0; i < qmcopy->nrQMatoms; i++)
394     {
395         qmcopy->shiftQM[i]        = qm->shiftQM[i];
396         qmcopy->indexQM[i]        = qm->indexQM[i];
397         qmcopy->atomicnumberQM[i] = qm->atomicnumberQM[i];
398     }
399     qmcopy->nelectrons   = qm->nelectrons;
400     qmcopy->multiplicity = qm->multiplicity;
401     qmcopy->QMcharge     = qm->QMcharge;
402     qmcopy->nelectrons   = qm->nelectrons;
403     qmcopy->QMmethod     = qm->QMmethod;
404     qmcopy->QMbasis      = qm->QMbasis;
405     /* trajectory surface hopping setup (Gaussian only) */
406     qmcopy->bSH          = qm->bSH;
407     qmcopy->CASorbitals  = qm->CASorbitals;
408     qmcopy->CASelectrons = qm->CASelectrons;
409     qmcopy->SAsteps      = qm->SAsteps;
410     qmcopy->SAon         = qm->SAon;
411     qmcopy->SAoff        = qm->SAoff;
412     qmcopy->bOPT         = qm->bOPT;
413
414     /* Gaussian init. variables */
415     qmcopy->nQMcpus      = qm->nQMcpus;
416     for (i = 0; i < DIM; i++)
417     {
418         qmcopy->SHbasis[i] = qm->SHbasis[i];
419     }
420     qmcopy->QMmem        = qm->QMmem;
421     qmcopy->accuracy     = qm->accuracy;
422     qmcopy->cpmcscf      = qm->cpmcscf;
423     qmcopy->SAstep       = qm->SAstep;
424     snew(qmcopy->frontatoms, qm->nrQMatoms);
425     snew(qmcopy->c12, qmcopy->nrQMatoms);
426     snew(qmcopy->c6, qmcopy->nrQMatoms);
427     if (qmcopy->bTS || qmcopy->bOPT)
428     {
429         for (i = 1; i < qmcopy->nrQMatoms; i++)
430         {
431             qmcopy->frontatoms[i] = qm->frontatoms[i];
432             qmcopy->c12[i]        = qm->c12[i];
433             qmcopy->c6[i]         = qm->c6[i];
434         }
435     }
436
437     return(qmcopy);
438
439 } /*copy_QMrec */
440
441 t_QMMMrec *mk_QMMMrec(void)
442 {
443
444     t_QMMMrec *qr;
445
446     snew(qr, 1);
447
448     return qr;
449
450 } /* mk_QMMMrec */
451
452 void init_QMMMrec(t_commrec  *cr,
453                   gmx_mtop_t *mtop,
454                   t_inputrec *ir,
455                   t_forcerec *fr)
456 {
457     /* we put the atomsnumbers of atoms that belong to the QMMM group in
458      * an array that will be copied later to QMMMrec->indexQM[..]. Also
459      * it will be used to create an QMMMrec->bQMMM index array that
460      * simply contains true/false for QM and MM (the other) atoms.
461      */
462
463     gmx_groups_t            *groups;
464     int                     *qm_arr = NULL, vsite, ai, aj;
465     int                      qm_max = 0, qm_nr = 0, i, j, jmax, k, l, nrvsite2 = 0;
466     t_QMMMrec               *qr;
467     t_MMrec                 *mm;
468     t_iatom                 *iatoms;
469     real                     c12au, c6au;
470     gmx_mtop_atomloop_all_t  aloop;
471     t_atom                  *atom;
472     gmx_mtop_ilistloop_all_t iloop;
473     int                      a_offset;
474     t_ilist                 *ilist_mol;
475     gmx_mtop_atomlookup_t    alook;
476
477     if (ir->cutoff_scheme != ecutsGROUP)
478     {
479         gmx_fatal(FARGS, "QMMM is currently only supported with cutoff-scheme=group");
480     }
481     if (!EI_DYNAMICS(ir->eI))
482     {
483         gmx_fatal(FARGS, "QMMM is only supported with dynamics");
484     }
485
486     c6au  = (HARTREE2KJ*AVOGADRO*gmx::power6(BOHR2NM));
487     c12au = (HARTREE2KJ*AVOGADRO*gmx::power12(BOHR2NM));
488     /* issue a fatal if the user wants to run with more than one node */
489     if (PAR(cr))
490     {
491         gmx_fatal(FARGS, "QM/MM does not work in parallel, use a single rank instead\n");
492     }
493
494     /* Make a local copy of the QMMMrec */
495     qr = fr->qr;
496
497     /* bQMMM[..] is an array containing TRUE/FALSE for atoms that are
498      * QM/not QM. We first set all elemenst at false. Afterwards we use
499      * the qm_arr (=MMrec->indexQM) to changes the elements
500      * corresponding to the QM atoms at TRUE.  */
501
502     qr->QMMMscheme     = ir->QMMMscheme;
503
504     /* we take the possibility into account that a user has
505      * defined more than one QM group:
506      */
507     /* an ugly work-around in case there is only one group In this case
508      * the whole system is treated as QM. Otherwise the second group is
509      * always the rest of the total system and is treated as MM.
510      */
511
512     /* small problem if there is only QM.... so no MM */
513
514     jmax = ir->opts.ngQM;
515
516     if (qr->QMMMscheme == eQMMMschemeoniom)
517     {
518         qr->nrQMlayers = jmax;
519     }
520     else
521     {
522         qr->nrQMlayers = 1;
523     }
524
525     groups = &mtop->groups;
526
527     /* there are jmax groups of QM atoms. In case of multiple QM groups
528      * I assume that the users wants to do ONIOM. However, maybe it
529      * should also be possible to define more than one QM subsystem with
530      * independent neighbourlists. I have to think about
531      * that.. 11-11-2003
532      */
533     snew(qr->qm, jmax);
534     for (j = 0; j < jmax; j++)
535     {
536         /* new layer */
537         aloop = gmx_mtop_atomloop_all_init(mtop);
538         while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
539         {
540             if (qm_nr >= qm_max)
541             {
542                 qm_max += 1000;
543                 srenew(qm_arr, qm_max);
544             }
545             if (ggrpnr(groups, egcQMMM, i) == j)
546             {
547                 /* hack for tip4p */
548                 qm_arr[qm_nr++] = i;
549             }
550         }
551         if (qr->QMMMscheme == eQMMMschemeoniom)
552         {
553             /* add the atoms to the bQMMM array
554              */
555
556             /* I assume that users specify the QM groups from small to
557              * big(ger) in the mdp file
558              */
559             qr->qm[j] = mk_QMrec();
560             /* we need to throw out link atoms that in the previous layer
561              * existed to separate this QMlayer from the previous
562              * QMlayer. We use the iatoms array in the idef for that
563              * purpose. If all atoms defining the current Link Atom (Dummy2)
564              * are part of the current QM layer it needs to be removed from
565              * qm_arr[].  */
566
567             iloop = gmx_mtop_ilistloop_all_init(mtop);
568             while (gmx_mtop_ilistloop_all_next(iloop, &ilist_mol, &a_offset))
569             {
570                 nrvsite2 = ilist_mol[F_VSITE2].nr;
571                 iatoms   = ilist_mol[F_VSITE2].iatoms;
572
573                 for (k = 0; k < nrvsite2; k += 4)
574                 {
575                     vsite = a_offset + iatoms[k+1]; /* the vsite         */
576                     ai    = a_offset + iatoms[k+2]; /* constructing atom */
577                     aj    = a_offset + iatoms[k+3]; /* constructing atom */
578                     if (ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, ai)
579                         &&
580                         ggrpnr(groups, egcQMMM, vsite) == ggrpnr(groups, egcQMMM, aj))
581                     {
582                         /* this dummy link atom needs to be removed from the qm_arr
583                          * before making the QMrec of this layer!
584                          */
585                         for (i = 0; i < qm_nr; i++)
586                         {
587                             if (qm_arr[i] == vsite)
588                             {
589                                 /* drop the element */
590                                 for (l = i; l < qm_nr; l++)
591                                 {
592                                     qm_arr[l] = qm_arr[l+1];
593                                 }
594                                 qm_nr--;
595                             }
596                         }
597                     }
598                 }
599             }
600
601             /* store QM atoms in this layer in the QMrec and initialise layer
602              */
603             init_QMrec(j, qr->qm[j], qm_nr, qm_arr, mtop, ir);
604
605             /* we now store the LJ C6 and C12 parameters in QM rec in case
606              * we need to do an optimization
607              */
608             if (qr->qm[j]->bOPT || qr->qm[j]->bTS)
609             {
610                 for (i = 0; i < qm_nr; i++)
611                 {
612                     /* nbfp now includes the 6.0/12.0 derivative prefactors */
613                     qr->qm[j]->c6[i]  =  C6(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c6au/6.0;
614                     qr->qm[j]->c12[i] = C12(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c12au/12.0;
615                 }
616             }
617             /* now we check for frontier QM atoms. These occur in pairs that
618              * construct the vsite
619              */
620             iloop = gmx_mtop_ilistloop_all_init(mtop);
621             while (gmx_mtop_ilistloop_all_next(iloop, &ilist_mol, &a_offset))
622             {
623                 nrvsite2 = ilist_mol[F_VSITE2].nr;
624                 iatoms   = ilist_mol[F_VSITE2].iatoms;
625
626                 for (k = 0; k < nrvsite2; k += 4)
627                 {
628                     vsite = a_offset + iatoms[k+1]; /* the vsite         */
629                     ai    = a_offset + iatoms[k+2]; /* constructing atom */
630                     aj    = a_offset + iatoms[k+3]; /* constructing atom */
631                     if (ggrpnr(groups, egcQMMM, ai) < (groups->grps[egcQMMM].nr-1) &&
632                         (ggrpnr(groups, egcQMMM, aj) >= (groups->grps[egcQMMM].nr-1)))
633                     {
634                         /* mark ai as frontier atom */
635                         for (i = 0; i < qm_nr; i++)
636                         {
637                             if ( (qm_arr[i] == ai) || (qm_arr[i] == vsite) )
638                             {
639                                 qr->qm[j]->frontatoms[i] = TRUE;
640                             }
641                         }
642                     }
643                     else if (ggrpnr(groups, egcQMMM, aj) < (groups->grps[egcQMMM].nr-1) &&
644                              (ggrpnr(groups, egcQMMM, ai) >= (groups->grps[egcQMMM].nr-1)))
645                     {
646                         /* mark aj as frontier atom */
647                         for (i = 0; i < qm_nr; i++)
648                         {
649                             if ( (qm_arr[i] == aj) || (qm_arr[i] == vsite))
650                             {
651                                 qr->qm[j]->frontatoms[i] = TRUE;
652                             }
653                         }
654                     }
655                 }
656             }
657         }
658     }
659     if (qr->QMMMscheme != eQMMMschemeoniom)
660     {
661
662         /* standard QMMM, all layers are merged together so there is one QM
663          * subsystem and one MM subsystem.
664          * Also we set the charges to zero in the md->charge arrays to prevent
665          * the innerloops from doubly counting the electostatic QM MM interaction
666          */
667
668         alook = gmx_mtop_atomlookup_init(mtop);
669
670         for (k = 0; k < qm_nr; k++)
671         {
672             gmx_mtop_atomnr_to_atom(alook, qm_arr[k], &atom);
673             atom->q  = 0.0;
674             atom->qB = 0.0;
675         }
676         qr->qm[0] = mk_QMrec();
677         /* store QM atoms in the QMrec and initialise
678          */
679         init_QMrec(0, qr->qm[0], qm_nr, qm_arr, mtop, ir);
680         if (qr->qm[0]->bOPT || qr->qm[0]->bTS)
681         {
682             for (i = 0; i < qm_nr; i++)
683             {
684                 gmx_mtop_atomnr_to_atom(alook, qm_arr[i], &atom);
685                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
686                 qr->qm[0]->c6[i]  =  C6(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c6au/6.0;
687                 qr->qm[0]->c12[i] = C12(fr->nbfp, mtop->ffparams.atnr, atom->type, atom->type)/c12au/12.0;
688             }
689         }
690
691         /* find frontier atoms and mark them true in the frontieratoms array.
692          */
693         for (i = 0; i < qm_nr; i++)
694         {
695             gmx_mtop_atomnr_to_ilist(alook, qm_arr[i], &ilist_mol, &a_offset);
696             nrvsite2 = ilist_mol[F_VSITE2].nr;
697             iatoms   = ilist_mol[F_VSITE2].iatoms;
698
699             for (k = 0; k < nrvsite2; k += 4)
700             {
701                 vsite = a_offset + iatoms[k+1]; /* the vsite         */
702                 ai    = a_offset + iatoms[k+2]; /* constructing atom */
703                 aj    = a_offset + iatoms[k+3]; /* constructing atom */
704                 if (ggrpnr(groups, egcQMMM, ai) < (groups->grps[egcQMMM].nr-1) &&
705                     (ggrpnr(groups, egcQMMM, aj) >= (groups->grps[egcQMMM].nr-1)))
706                 {
707                     /* mark ai as frontier atom */
708                     if ( (qm_arr[i] == ai) || (qm_arr[i] == vsite) )
709                     {
710                         qr->qm[0]->frontatoms[i] = TRUE;
711                     }
712                 }
713                 else if (ggrpnr(groups, egcQMMM, aj) < (groups->grps[egcQMMM].nr-1) &&
714                          (ggrpnr(groups, egcQMMM, ai) >= (groups->grps[egcQMMM].nr-1)))
715                 {
716                     /* mark aj as frontier atom */
717                     if ( (qm_arr[i] == aj) || (qm_arr[i] == vsite) )
718                     {
719                         qr->qm[0]->frontatoms[i] = TRUE;
720                     }
721                 }
722             }
723         }
724
725         gmx_mtop_atomlookup_destroy(alook);
726
727         /* MM rec creation */
728         mm               = mk_MMrec();
729         mm->scalefactor  = ir->scalefactor;
730         mm->nrMMatoms    = (mtop->natoms)-(qr->qm[0]->nrQMatoms); /* rest of the atoms */
731         qr->mm           = mm;
732     }
733     else /* ONIOM */
734     {    /* MM rec creation */
735         mm               = mk_MMrec();
736         mm->scalefactor  = ir->scalefactor;
737         mm->nrMMatoms    = 0;
738         qr->mm           = mm;
739     }
740
741     /* these variables get updated in the update QMMMrec */
742
743     if (qr->nrQMlayers == 1)
744     {
745         /* with only one layer there is only one initialisation
746          * needed. Multilayer is a bit more complicated as it requires
747          * re-initialisation at every step of the simulation. This is due
748          * to the use of COMMON blocks in the fortran QM subroutines.
749          */
750         if (qr->qm[0]->QMmethod < eQMmethodRHF)
751         {
752 #if GMX_QMMM_MOPAC
753             /* semi-empiprical 1-layer ONIOM calculation requested (mopac93) */
754             init_mopac(qr->qm[0]);
755 #else
756             gmx_fatal(FARGS, "Semi-empirical QM only supported with Mopac.");
757 #endif
758         }
759         else
760         {
761             /* ab initio calculation requested (gamess/gaussian/ORCA) */
762 #if GMX_QMMM_GAMESS
763             init_gamess(cr, qr->qm[0], qr->mm);
764 #elif GMX_QMMM_GAUSSIAN
765             init_gaussian(qr->qm[0]);
766 #elif GMX_QMMM_ORCA
767             init_orca(qr->qm[0]);
768 #else
769             gmx_fatal(FARGS, "Ab-initio calculation only supported with Gamess, Gaussian or ORCA.");
770 #endif
771         }
772     }
773 } /* init_QMMMrec */
774
775 void update_QMMMrec(t_commrec      *cr,
776                     t_forcerec     *fr,
777                     rvec            x[],
778                     t_mdatoms      *md,
779                     matrix          box,
780                     gmx_localtop_t *top)
781 {
782     /* updates the coordinates of both QM atoms and MM atoms and stores
783      * them in the QMMMrec.
784      *
785      * NOTE: is NOT yet working if there are no PBC. Also in ns.c, simple
786      * ns needs to be fixed!
787      */
788     int
789         mm_max = 0, mm_nr = 0, mm_nr_new, i, j, is, k, shift;
790     t_j_particle
791        *mm_j_particles = NULL, *qm_i_particles = NULL;
792     t_QMMMrec
793        *qr;
794     t_nblist
795        *QMMMlist;
796     rvec
797         dx, crd;
798     t_QMrec
799        *qm;
800     t_MMrec
801        *mm;
802     t_pbc
803         pbc;
804     int
805        *parallelMMarray = NULL;
806     real
807         c12au, c6au;
808
809     c6au  = (HARTREE2KJ*AVOGADRO*gmx::power6(BOHR2NM));
810     c12au = (HARTREE2KJ*AVOGADRO*gmx::power12(BOHR2NM));
811
812     /* every cpu has this array. On every processor we fill this array
813      * with 1's and 0's. 1's indicate the atoms is a QM atom on the
814      * current cpu in a later stage these arrays are all summed. indexes
815      * > 0 indicate the atom is a QM atom. Every node therefore knows
816      * whcih atoms are part of the QM subsystem.
817      */
818     /* copy some pointers */
819     qr          = fr->qr;
820     mm          = qr->mm;
821     QMMMlist    = fr->QMMMlist;
822
823
824
825     /*  init_pbc(box);  needs to be called first, see pbc.h */
826     ivec null_ivec;
827     clear_ivec(null_ivec);
828     set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : null_ivec,
829                FALSE, box);
830     /* only in standard (normal) QMMM we need the neighbouring MM
831      * particles to provide a electric field of point charges for the QM
832      * atoms.
833      */
834     if (qr->QMMMscheme == eQMMMschemenormal) /* also implies 1 QM-layer */
835     {
836         /* we NOW create/update a number of QMMMrec entries:
837          *
838          * 1) the shiftQM, containing the shifts of the QM atoms
839          *
840          * 2) the indexMM array, containing the index of the MM atoms
841          *
842          * 3) the shiftMM, containing the shifts of the MM atoms
843          *
844          * 4) the shifted coordinates of the MM atoms
845          *
846          * the shifts are used for computing virial of the QM/MM particles.
847          */
848         qm = qr->qm[0]; /* in case of normal QMMM, there is only one group */
849         snew(qm_i_particles, QMMMlist->nri);
850         if (QMMMlist->nri)
851         {
852             qm_i_particles[0].shift = XYZ2IS(0, 0, 0);
853             for (i = 0; i < QMMMlist->nri; i++)
854             {
855                 qm_i_particles[i].j     = QMMMlist->iinr[i];
856
857                 if (i)
858                 {
859                     qm_i_particles[i].shift = pbc_dx_aiuc(&pbc, x[QMMMlist->iinr[0]],
860                                                           x[QMMMlist->iinr[i]], dx);
861
862                 }
863                 /* However, since nri >= nrQMatoms, we do a quicksort, and throw
864                  * out double, triple, etc. entries later, as we do for the MM
865                  * list too.
866                  */
867
868                 /* compute the shift for the MM j-particles with respect to
869                  * the QM i-particle and store them.
870                  */
871
872                 crd[0] = IS2X(QMMMlist->shift[i]) + IS2X(qm_i_particles[i].shift);
873                 crd[1] = IS2Y(QMMMlist->shift[i]) + IS2Y(qm_i_particles[i].shift);
874                 crd[2] = IS2Z(QMMMlist->shift[i]) + IS2Z(qm_i_particles[i].shift);
875                 is     = static_cast<int>(XYZ2IS(crd[0], crd[1], crd[2]));
876                 for (j = QMMMlist->jindex[i];
877                      j < QMMMlist->jindex[i+1];
878                      j++)
879                 {
880                     if (mm_nr >= mm_max)
881                     {
882                         mm_max += 1000;
883                         srenew(mm_j_particles, mm_max);
884                     }
885
886                     mm_j_particles[mm_nr].j     = QMMMlist->jjnr[j];
887                     mm_j_particles[mm_nr].shift = is;
888                     mm_nr++;
889                 }
890             }
891
892             /* quicksort QM and MM shift arrays and throw away multiple entries */
893
894
895
896             qsort(qm_i_particles, QMMMlist->nri,
897                   (size_t)sizeof(qm_i_particles[0]),
898                   struct_comp);
899             /* The mm_j_particles argument to qsort is not allowed to be NULL */
900             if (mm_nr > 0)
901             {
902                 qsort(mm_j_particles, mm_nr,
903                       (size_t)sizeof(mm_j_particles[0]),
904                       struct_comp);
905             }
906             /* remove multiples in the QM shift array, since in init_QMMM() we
907              * went through the atom numbers from 0 to md.nr, the order sorted
908              * here matches the one of QMindex already.
909              */
910             j = 0;
911             for (i = 0; i < QMMMlist->nri; i++)
912             {
913                 if (i == 0 || qm_i_particles[i].j != qm_i_particles[i-1].j)
914                 {
915                     qm_i_particles[j++] = qm_i_particles[i];
916                 }
917             }
918             mm_nr_new = 0;
919             if (qm->bTS || qm->bOPT)
920             {
921                 /* only remove double entries for the MM array */
922                 for (i = 0; i < mm_nr; i++)
923                 {
924                     if ((i == 0 || mm_j_particles[i].j != mm_j_particles[i-1].j)
925                         && !md->bQM[mm_j_particles[i].j])
926                     {
927                         mm_j_particles[mm_nr_new++] = mm_j_particles[i];
928                     }
929                 }
930             }
931             /* we also remove mm atoms that have no charges!
932              * actually this is already done in the ns.c
933              */
934             else
935             {
936                 for (i = 0; i < mm_nr; i++)
937                 {
938                     if ((i == 0 || mm_j_particles[i].j != mm_j_particles[i-1].j)
939                         && !md->bQM[mm_j_particles[i].j]
940                         && (md->chargeA[mm_j_particles[i].j]
941                             || (md->chargeB && md->chargeB[mm_j_particles[i].j])))
942                     {
943                         mm_j_particles[mm_nr_new++] = mm_j_particles[i];
944                     }
945                 }
946             }
947             mm_nr = mm_nr_new;
948             /* store the data retrieved above into the QMMMrec
949              */
950             k = 0;
951             /* Keep the compiler happy,
952              * shift will always be set in the loop for i=0
953              */
954             shift = 0;
955             for (i = 0; i < qm->nrQMatoms; i++)
956             {
957                 /* not all qm particles might have appeared as i
958                  * particles. They might have been part of the same charge
959                  * group for instance.
960                  */
961                 if (qm->indexQM[i] == qm_i_particles[k].j)
962                 {
963                     shift = qm_i_particles[k++].shift;
964                 }
965                 /* use previous shift, assuming they belong the same charge
966                  * group anyway,
967                  */
968
969                 qm->shiftQM[i] = shift;
970             }
971         }
972         /* parallel excecution */
973         if (PAR(cr))
974         {
975             snew(parallelMMarray, 2*(md->nr));
976             /* only MM particles have a 1 at their atomnumber. The second part
977              * of the array contains the shifts. Thus:
978              * p[i]=1/0 depending on wether atomnumber i is a MM particle in the QM
979              * step or not. p[i+md->nr] is the shift of atomnumber i.
980              */
981             for (i = 0; i < 2*(md->nr); i++)
982             {
983                 parallelMMarray[i] = 0;
984             }
985
986             for (i = 0; i < mm_nr; i++)
987             {
988                 parallelMMarray[mm_j_particles[i].j]          = 1;
989                 parallelMMarray[mm_j_particles[i].j+(md->nr)] = mm_j_particles[i].shift;
990             }
991             gmx_sumi(md->nr, parallelMMarray, cr);
992             mm_nr = 0;
993
994             mm_max = 0;
995             for (i = 0; i < md->nr; i++)
996             {
997                 if (parallelMMarray[i])
998                 {
999                     if (mm_nr >= mm_max)
1000                     {
1001                         mm_max += 1000;
1002                         srenew(mm->indexMM, mm_max);
1003                         srenew(mm->shiftMM, mm_max);
1004                     }
1005                     mm->indexMM[mm_nr]   = i;
1006                     mm->shiftMM[mm_nr++] = parallelMMarray[i+md->nr]/parallelMMarray[i];
1007                 }
1008             }
1009             mm->nrMMatoms = mm_nr;
1010             free(parallelMMarray);
1011         }
1012         /* serial execution */
1013         else
1014         {
1015             mm->nrMMatoms = mm_nr;
1016             srenew(mm->shiftMM, mm_nr);
1017             srenew(mm->indexMM, mm_nr);
1018             for (i = 0; i < mm_nr; i++)
1019             {
1020                 mm->indexMM[i] = mm_j_particles[i].j;
1021                 mm->shiftMM[i] = mm_j_particles[i].shift;
1022             }
1023
1024         }
1025         /* (re) allocate memory for the MM coordiate array. The QM
1026          * coordinate array was already allocated in init_QMMM, and is
1027          * only (re)filled in the update_QMMM_coordinates routine
1028          */
1029         srenew(mm->xMM, mm->nrMMatoms);
1030         /* now we (re) fill the array that contains the MM charges with
1031          * the forcefield charges. If requested, these charges will be
1032          * scaled by a factor
1033          */
1034         srenew(mm->MMcharges, mm->nrMMatoms);
1035         for (i = 0; i < mm->nrMMatoms; i++) /* no free energy yet */
1036         {
1037             mm->MMcharges[i] = md->chargeA[mm->indexMM[i]]*mm->scalefactor;
1038         }
1039         if (qm->bTS || qm->bOPT)
1040         {
1041             /* store (copy) the c6 and c12 parameters into the MMrec struct
1042              */
1043             srenew(mm->c6, mm->nrMMatoms);
1044             srenew(mm->c12, mm->nrMMatoms);
1045             for (i = 0; i < mm->nrMMatoms; i++)
1046             {
1047                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1048                 mm->c6[i]  = C6(fr->nbfp, top->idef.atnr, md->typeA[mm->indexMM[i]], md->typeA[mm->indexMM[i]])/c6au/6.0;
1049                 mm->c12[i] = C12(fr->nbfp, top->idef.atnr, md->typeA[mm->indexMM[i]], md->typeA[mm->indexMM[i]])/c12au/12.0;
1050             }
1051             punch_QMMM_excl(qr->qm[0], mm, &(top->excls));
1052         }
1053         /* the next routine fills the coordinate fields in the QMMM rec of
1054          * both the qunatum atoms and the MM atoms, using the shifts
1055          * calculated above.
1056          */
1057
1058         update_QMMM_coord(x, fr, qr->qm[0], qr->mm);
1059         free(qm_i_particles);
1060         free(mm_j_particles);
1061     }
1062     else /* ONIOM */ /* ????? */
1063     {
1064         mm->nrMMatoms = 0;
1065         /* do for each layer */
1066         for (j = 0; j < qr->nrQMlayers; j++)
1067         {
1068             qm             = qr->qm[j];
1069             qm->shiftQM[0] = XYZ2IS(0, 0, 0);
1070             for (i = 1; i < qm->nrQMatoms; i++)
1071             {
1072                 qm->shiftQM[i] = pbc_dx_aiuc(&pbc, x[qm->indexQM[0]], x[qm->indexQM[i]],
1073                                              dx);
1074             }
1075             update_QMMM_coord(x, fr, qm, mm);
1076         }
1077     }
1078 } /* update_QMMM_rec */
1079
1080
1081 real calculate_QMMM(t_commrec *cr,
1082                     rvec x[], rvec f[],
1083                     t_forcerec *fr)
1084 {
1085     real
1086         QMener = 0.0;
1087     /* a selection for the QM package depending on which is requested
1088      * (Gaussian, GAMESS-UK, MOPAC or ORCA) needs to be implemented here. Now
1089      * it works through defines.... Not so nice yet
1090      */
1091     t_QMMMrec
1092     *qr;
1093     t_QMrec
1094     *qm, *qm2;
1095     t_MMrec
1096     *mm = NULL;
1097     rvec
1098     *forces  = NULL, *fshift = NULL,
1099     *forces2 = NULL, *fshift2 = NULL; /* needed for multilayer ONIOM */
1100     int
1101         i, j, k;
1102     /* make a local copy the QMMMrec pointer
1103      */
1104     qr = fr->qr;
1105     mm = qr->mm;
1106
1107     /* now different procedures are carried out for one layer ONION and
1108      * normal QMMM on one hand and multilayer oniom on the other
1109      */
1110     if (qr->QMMMscheme == eQMMMschemenormal || qr->nrQMlayers == 1)
1111     {
1112         qm = qr->qm[0];
1113         snew(forces, (qm->nrQMatoms+mm->nrMMatoms));
1114         snew(fshift, (qm->nrQMatoms+mm->nrMMatoms));
1115         QMener = call_QMroutine(cr, fr, qm, mm, forces, fshift);
1116         for (i = 0; i < qm->nrQMatoms; i++)
1117         {
1118             for (j = 0; j < DIM; j++)
1119             {
1120                 f[qm->indexQM[i]][j]          -= forces[i][j];
1121                 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
1122             }
1123         }
1124         for (i = 0; i < mm->nrMMatoms; i++)
1125         {
1126             for (j = 0; j < DIM; j++)
1127             {
1128                 f[mm->indexMM[i]][j]          -= forces[qm->nrQMatoms+i][j];
1129                 fr->fshift[mm->shiftMM[i]][j] += fshift[qm->nrQMatoms+i][j];
1130             }
1131
1132         }
1133         free(forces);
1134         free(fshift);
1135     }
1136     else                                       /* Multi-layer ONIOM */
1137     {
1138         for (i = 0; i < qr->nrQMlayers-1; i++) /* last layer is special */
1139         {
1140             qm  = qr->qm[i];
1141             qm2 = copy_QMrec(qr->qm[i+1]);
1142
1143             qm2->nrQMatoms = qm->nrQMatoms;
1144
1145             for (j = 0; j < qm2->nrQMatoms; j++)
1146             {
1147                 for (k = 0; k < DIM; k++)
1148                 {
1149                     qm2->xQM[j][k]       = qm->xQM[j][k];
1150                 }
1151                 qm2->indexQM[j]        = qm->indexQM[j];
1152                 qm2->atomicnumberQM[j] = qm->atomicnumberQM[j];
1153                 qm2->shiftQM[j]        = qm->shiftQM[j];
1154             }
1155
1156             qm2->QMcharge = qm->QMcharge;
1157             /* this layer at the higher level of theory */
1158             srenew(forces, qm->nrQMatoms);
1159             srenew(fshift, qm->nrQMatoms);
1160             /* we need to re-initialize the QMroutine every step... */
1161             init_QMroutine(cr, qm, mm);
1162             QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
1163
1164             /* this layer at the lower level of theory */
1165             srenew(forces2, qm->nrQMatoms);
1166             srenew(fshift2, qm->nrQMatoms);
1167             init_QMroutine(cr, qm2, mm);
1168             QMener -= call_QMroutine(cr, fr, qm2, mm, forces2, fshift2);
1169             /* E = E1high-E1low The next layer includes the current layer at
1170              * the lower level of theory, which provides + E2low
1171              * this is similar for gradients
1172              */
1173             for (i = 0; i < qm->nrQMatoms; i++)
1174             {
1175                 for (j = 0; j < DIM; j++)
1176                 {
1177                     f[qm->indexQM[i]][j]          -= (forces[i][j]-forces2[i][j]);
1178                     fr->fshift[qm->shiftQM[i]][j] += (fshift[i][j]-fshift2[i][j]);
1179                 }
1180             }
1181             free(qm2);
1182         }
1183         /* now the last layer still needs to be done: */
1184         qm      = qr->qm[qr->nrQMlayers-1]; /* C counts from 0 */
1185         init_QMroutine(cr, qm, mm);
1186         srenew(forces, qm->nrQMatoms);
1187         srenew(fshift, qm->nrQMatoms);
1188         QMener += call_QMroutine(cr, fr, qm, mm, forces, fshift);
1189         for (i = 0; i < qm->nrQMatoms; i++)
1190         {
1191             for (j = 0; j < DIM; j++)
1192             {
1193                 f[qm->indexQM[i]][j]          -= forces[i][j];
1194                 fr->fshift[qm->shiftQM[i]][j] += fshift[i][j];
1195             }
1196         }
1197         free(forces);
1198         free(fshift);
1199         free(forces2);
1200         free(fshift2);
1201     }
1202     if (qm->bTS || qm->bOPT)
1203     {
1204         /* qm[0] still contains the largest ONIOM QM subsystem
1205          * we take the optimized coordiates and put the in x[]
1206          */
1207         for (i = 0; i < qm->nrQMatoms; i++)
1208         {
1209             for (j = 0; j < DIM; j++)
1210             {
1211                 x[qm->indexQM[i]][j] = qm->xQM[i][j];
1212             }
1213         }
1214     }
1215     return(QMener);
1216 } /* calculate_QMMM */
1217
1218 /* end of QMMM core routines */