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