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