Merge "Extended genbox to insert molecules at given positions"
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_generic_adress.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009 Christoph Junghans, Brad Lambeth.
5  * Copyright (c) 2011 Christoph Junghans, Sebastian Fritsch
6  * Copyright (c) 2012, The GROMACS development team,
7  * check out http://www.gromacs.org for more information.
8  * Copyright (c) 2012, by the GROMACS development team, led by
9  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
10  * others, as listed in the AUTHORS file in the top-level source
11  * directory and at http://www.gromacs.org.
12  *
13  * GROMACS is free software; you can redistribute it and/or
14  * modify it under the terms of the GNU Lesser General Public License
15  * as published by the Free Software Foundation; either version 2.1
16  * of the License, or (at your option) any later version.
17  *
18  * GROMACS is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
21  * Lesser General Public License for more details.
22  *
23  * You should have received a copy of the GNU Lesser General Public
24  * License along with GROMACS; if not, see
25  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
26  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
27  *
28  * If you want to redistribute modifications to GROMACS, please
29  * consider that scientific software is very special. Version
30  * control is crucial - bugs must be traceable. We will be happy to
31  * consider code for inclusion in the official distribution, but
32  * derived work must not be called official GROMACS. Details are found
33  * in the README & COPYING files - if they are missing, get the
34  * official version at http://www.gromacs.org.
35  *
36  * To help us fund GROMACS development, we humbly ask that you cite
37  * the research papers on the package. Check out http://www.gromacs.org.
38  */
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <math.h>
44
45 #include "types/simple.h"
46 #include "vec.h"
47 #include "typedefs.h"
48 #include "nb_generic_adress.h"
49 #include "nrnb.h"
50
51 #include "nonbonded.h"
52 #include "nb_kernel.h"
53
54 #define ALMOST_ZERO 1e-30
55 #define ALMOST_ONE 1-(1e-30)
56 void
57 gmx_nb_generic_adress_kernel(t_nblist *                nlist,
58                       rvec *                    xx,
59                       rvec *                    ff,
60                       t_forcerec *              fr,
61                       t_mdatoms *               mdatoms,
62                       nb_kernel_data_t *        kernel_data,
63                       t_nrnb *                  nrnb)
64 {
65     int           nri,ntype,table_nelements,ielec,ivdw;
66     real          facel,gbtabscale;
67     int           n,ii,is3,ii3,k,nj0,nj1,jnr,j3,ggid,nnn,n0;
68     real          shX,shY,shZ;
69     real          fscal,felec,fvdw,velec,vvdw,tx,ty,tz;
70     real          rinvsq;
71     real          iq;
72     real          qq,vctot;
73     int           nti,nvdwparam;
74     int           tj;
75     real          rt,r,eps,eps2,Y,F,Geps,Heps2,VV,FF,Fp,fijD,fijR;
76     real          rinvsix;
77     real          vvdwtot;
78     real          vvdw_rep,vvdw_disp;
79     real          ix,iy,iz,fix,fiy,fiz;
80     real          jx,jy,jz;
81     real          dx,dy,dz,rsq,rinv;
82     real          c6,c12,cexp1,cexp2,br;
83     real *        charge;
84     real *        shiftvec;
85     real *        vdwparam;
86     int *         shift;
87     int *         type;
88     real *        fshift;
89     real *        velecgrp;
90     real *        vvdwgrp;
91     real          tabscale;
92     real *        VFtab;
93     real *        x;
94     real *        f;
95     int           ewitab;
96     real          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
97     real *        ewtab;
98     real          rcoulomb2,rvdw,rvdw2,sh_invrc6;
99     real          rcutoff,rcutoff2;
100     real          rswitch_elec,rswitch_vdw,d,d2,sw,dsw,rinvcorr;
101     real          elec_swV3,elec_swV4,elec_swV5,elec_swF2,elec_swF3,elec_swF4;
102     real          vdw_swV3,vdw_swV4,vdw_swV5,vdw_swF2,vdw_swF3,vdw_swF4;
103     gmx_bool      bExactElecCutoff,bExactVdwCutoff,bExactCutoff;
104
105     real *     wf;
106     real       weight_cg1;
107     real       weight_cg2;
108     real       weight_product;
109     real       hybscal; /* the multiplicator to the force for hybrid interactions*/
110     real       force_cap;
111     gmx_bool   bCG;
112     int egp_nr;
113
114     wf                  = mdatoms->wf;
115
116     force_cap           = fr->adress_ex_forcecap;
117
118     x                   = xx[0];
119     f                   = ff[0];
120     ielec               = nlist->ielec;
121     ivdw                = nlist->ivdw;
122
123     fshift              = fr->fshift[0];
124     velecgrp            = kernel_data->energygrp_elec;
125     vvdwgrp             = kernel_data->energygrp_vdw;
126     tabscale            = kernel_data->table_elec_vdw->scale;
127     VFtab               = kernel_data->table_elec_vdw->data;
128
129     sh_ewald            = fr->ic->sh_ewald;
130     ewtab               = fr->ic->tabq_coul_FDV0;
131     ewtabscale          = fr->ic->tabq_scale;
132     ewtabhalfspace      = 0.5/ewtabscale;
133
134     rcoulomb2           = fr->rcoulomb*fr->rcoulomb;
135     rvdw                = fr->rvdw;
136     rvdw2               = rvdw*rvdw;
137     sh_invrc6           = fr->ic->sh_invrc6;
138
139     if(fr->coulomb_modifier==eintmodPOTSWITCH)
140     {
141         d               = fr->rcoulomb-fr->rcoulomb_switch;
142         elec_swV3       = -10.0/(d*d*d);
143         elec_swV4       =  15.0/(d*d*d*d);
144         elec_swV5       =  -6.0/(d*d*d*d*d);
145         elec_swF2       = -30.0/(d*d*d);
146         elec_swF3       =  60.0/(d*d*d*d);
147         elec_swF4       = -30.0/(d*d*d*d*d);
148     }
149     else
150     {
151         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
152         elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
153     }
154     if(fr->vdw_modifier==eintmodPOTSWITCH)
155     {
156         d               = fr->rvdw-fr->rvdw_switch;
157         vdw_swV3        = -10.0/(d*d*d);
158         vdw_swV4        =  15.0/(d*d*d*d);
159         vdw_swV5        =  -6.0/(d*d*d*d*d);
160         vdw_swF2        = -30.0/(d*d*d);
161         vdw_swF3        =  60.0/(d*d*d*d);
162         vdw_swF4        = -30.0/(d*d*d*d*d);
163     }
164     else
165     {
166         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
167         vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
168     }
169
170     bExactElecCutoff    = (fr->coulomb_modifier!=eintmodNONE) || fr->eeltype == eelRF_ZERO;
171     bExactVdwCutoff     = (fr->vdw_modifier!=eintmodNONE);
172     bExactCutoff        = bExactElecCutoff || bExactVdwCutoff;
173
174     if(bExactCutoff)
175     {
176         rcutoff  = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
177         rcutoff2 = rcutoff*rcutoff;
178     }
179     else
180     {
181         /* Fix warnings for stupid compilers */
182         rcutoff = rcutoff2 = 1e30;
183     }
184
185     /* avoid compiler warnings for cases that cannot happen */
186     nnn                 = 0;
187     eps                 = 0.0;
188     eps2                = 0.0;
189
190     /* 3 VdW parameters for buckingham, otherwise 2 */
191     nvdwparam           = (ivdw==GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
192     table_nelements     = 12;
193
194     charge              = mdatoms->chargeA;
195     type                = mdatoms->typeA;
196     facel               = fr->epsfac;
197     shiftvec            = fr->shift_vec[0];
198     vdwparam            = fr->nbfp;
199     ntype               = fr->ntype;
200
201     for(n=0; (n<nlist->nri); n++)
202     {
203         is3              = 3*nlist->shift[n];
204         shX              = shiftvec[is3];  
205         shY              = shiftvec[is3+1];
206         shZ              = shiftvec[is3+2];
207         nj0              = nlist->jindex[n];      
208         nj1              = nlist->jindex[n+1];    
209         ii               = nlist->iinr[n];        
210         ii3              = 3*ii;           
211         ix               = shX + x[ii3+0];
212         iy               = shY + x[ii3+1];
213         iz               = shZ + x[ii3+2];
214         iq               = facel*charge[ii];
215         nti              = nvdwparam*ntype*type[ii];
216         vctot            = 0;              
217         vvdwtot          = 0;
218         fix              = 0;
219         fiy              = 0;
220         fiz              = 0;              
221
222         /* We need to find out if this i atom is part of an
223            all-atom or CG energy group  */
224         egp_nr=mdatoms->cENER[ii];
225         bCG= !fr->adress_group_explicit[egp_nr];
226         
227         weight_cg1       = wf[ii];
228
229         if ((!bCG) && weight_cg1 < ALMOST_ZERO){
230             continue;
231         }
232
233         for(k=nj0; (k<nj1); k++)
234         {
235             jnr              = nlist->jjnr[k];        
236             weight_cg2       = wf[jnr];
237             weight_product   = weight_cg1*weight_cg2;
238
239             if (weight_product < ALMOST_ZERO)
240             {                
241                 /* if it's a explicit loop, skip this atom */
242                 if (!bCG)
243                 {
244                     continue;
245                 }
246                 else /* if it's a coarse grained loop, include this atom */
247                 {
248                     hybscal = 1.0;
249                 }
250             }
251             else if (weight_product >= ALMOST_ONE)
252             {
253                 
254                 /* if it's a explicit loop, include this atom */
255                 if(!bCG)
256                 {
257                     hybscal = 1.0;
258                 }             
259                 else  /* if it's a coarse grained loop, skip this atom */
260                 {
261                     continue;
262                 }
263             }
264             /* both have double identity, get hybrid scaling factor */
265             else
266             {
267                 hybscal = weight_product;
268
269                 if(bCG)
270                 {
271                     hybscal = 1.0 - hybscal;
272                 }
273             }
274
275             j3               = 3*jnr;          
276             jx               = x[j3+0];      
277             jy               = x[j3+1];      
278             jz               = x[j3+2];      
279             dx               = ix - jx;      
280             dy               = iy - jy;      
281             dz               = iz - jz;      
282             rsq              = dx*dx+dy*dy+dz*dz;
283             rinv             = gmx_invsqrt(rsq);
284             rinvsq           = rinv*rinv;  
285             felec            = 0;
286             fvdw             = 0;
287             velec            = 0;
288             vvdw             = 0;
289
290             if(bExactCutoff && rsq>rcutoff2)
291             {
292                 continue;
293             }
294
295             if(ielec==GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw==GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
296             {
297                 r                = rsq*rinv;
298                 rt               = r*tabscale;
299                 n0               = rt;
300                 eps              = rt-n0;
301                 eps2             = eps*eps;
302                 nnn              = table_nelements*n0;
303             }
304
305             /* Coulomb interaction. ielec==0 means no interaction */
306             if(ielec!=GMX_NBKERNEL_ELEC_NONE)
307             {
308                 qq               = iq*charge[jnr];
309
310                 switch(ielec)
311                 {
312                     case GMX_NBKERNEL_ELEC_NONE:
313                         break;
314
315                     case GMX_NBKERNEL_ELEC_COULOMB:
316                         /* Vanilla cutoff coulomb */
317                         velec            = qq*rinv;
318                         felec            = velec*rinvsq;
319                         break;
320
321                     case GMX_NBKERNEL_ELEC_REACTIONFIELD:
322                         /* Reaction-field */
323                         velec            = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
324                         felec            = qq*(rinv*rinvsq-2.0*fr->k_rf);
325                         break;
326
327                     case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
328                         /* Tabulated coulomb */
329                         Y                = VFtab[nnn];
330                         F                = VFtab[nnn+1];
331                         Geps             = eps*VFtab[nnn+2];
332                         Heps2            = eps2*VFtab[nnn+3];
333                         Fp               = F+Geps+Heps2;
334                         VV               = Y+eps*Fp;
335                         FF               = Fp+Geps+2.0*Heps2;
336                         velec            = qq*VV;
337                         felec            = -qq*FF*tabscale*rinv;
338                         break;
339
340                     case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
341                         /* GB */
342                         gmx_fatal(FARGS,"Death & horror! GB generic interaction not implemented.\n");
343                         break;
344
345                     case GMX_NBKERNEL_ELEC_EWALD:
346                         ewrt             = rsq*rinv*ewtabscale;
347                         ewitab           = ewrt;
348                         eweps            = ewrt-ewitab;
349                         ewitab           = 4*ewitab;
350                         felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
351                         rinvcorr         = (fr->coulomb_modifier==eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
352                         velec            = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
353                         felec            = qq*rinv*(rinvsq-felec);
354                         break;
355
356                     default:
357                         gmx_fatal(FARGS,"Death & horror! No generic coulomb interaction for ielec=%d.\n",ielec);
358                         break;
359                 }
360                 if(fr->coulomb_modifier==eintmodPOTSWITCH)
361                 {
362                     d                = rsq*rinv-fr->rcoulomb_switch;
363                     d                = (d>0.0) ? d : 0.0;
364                     d2               = d*d;
365                     sw               = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
366                     dsw              = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
367                     /* Apply switch function. Note that felec=f/r since it will be multiplied
368                      * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
369                      * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
370                      */
371                     felec            = felec*sw - rinv*velec*dsw;
372                     /* Once we have used velec to update felec we can modify velec too */
373                     velec           *= sw;
374                 }
375                 if(bExactElecCutoff)
376                 {
377                     felec            = (rsq<=rcoulomb2) ? felec : 0.0;
378                     velec            = (rsq<=rcoulomb2) ? velec : 0.0;
379                 }
380                 vctot           += velec;
381             } /* End of coulomb interactions */
382
383
384             /* VdW interaction. ivdw==0 means no interaction */
385             if(ivdw!=GMX_NBKERNEL_VDW_NONE)
386             {
387                 tj               = nti+nvdwparam*type[jnr];
388
389                 switch(ivdw)
390                 {
391                     case GMX_NBKERNEL_VDW_NONE:
392                         break;
393
394                     case GMX_NBKERNEL_VDW_LENNARDJONES:
395                         /* Vanilla Lennard-Jones cutoff */
396                         c6               = vdwparam[tj];
397                         c12              = vdwparam[tj+1];
398                         rinvsix          = rinvsq*rinvsq*rinvsq;
399                         vvdw_disp        = c6*rinvsix;
400                         vvdw_rep         = c12*rinvsix*rinvsix;
401                         fvdw             = (vvdw_rep-vvdw_disp)*rinvsq;
402                         if(fr->vdw_modifier==eintmodPOTSHIFT)
403                         {
404                             vvdw             = (vvdw_rep-c12*sh_invrc6*sh_invrc6)*(1.0/12.0)-(vvdw_disp-c6*sh_invrc6)*(1.0/6.0);
405                         }
406                         else
407                         {
408                             vvdw             = vvdw_rep/12.0-vvdw_disp/6.0;
409                         }
410                         break;
411
412                     case GMX_NBKERNEL_VDW_BUCKINGHAM:
413                         /* Buckingham */
414                         c6               = vdwparam[tj];
415                         cexp1            = vdwparam[tj+1];
416                         cexp2            = vdwparam[tj+2];
417
418                         rinvsix          = rinvsq*rinvsq*rinvsq;
419                         vvdw_disp        = c6*rinvsix;
420                         br               = cexp2*rsq*rinv;
421                         vvdw_rep         = cexp1*exp(-br);
422                         fvdw             = (br*vvdw_rep-vvdw_disp)*rinvsq;
423                         if(fr->vdw_modifier==eintmodPOTSHIFT)
424                         {
425                             vvdw             = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp-c6*sh_invrc6)/6.0;
426                         }
427                         else
428                         {
429                             vvdw             = vvdw_rep-vvdw_disp/6.0;
430                         }
431                         break;
432
433                     case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
434                         /* Tabulated VdW */
435                         c6               = vdwparam[tj];
436                         c12              = vdwparam[tj+1];
437                         Y                = VFtab[nnn+4];
438                         F                = VFtab[nnn+5];
439                         Geps             = eps*VFtab[nnn+6];
440                         Heps2            = eps2*VFtab[nnn+7];
441                         Fp               = F+Geps+Heps2;
442                         VV               = Y+eps*Fp;
443                         FF               = Fp+Geps+2.0*Heps2;
444                         vvdw_disp        = c6*VV;
445                         fijD             = c6*FF;
446                         Y                = VFtab[nnn+8];
447                         F                = VFtab[nnn+9];
448                         Geps             = eps*VFtab[nnn+10];
449                         Heps2            = eps2*VFtab[nnn+11];
450                         Fp               = F+Geps+Heps2;
451                         VV               = Y+eps*Fp;
452                         FF               = Fp+Geps+2.0*Heps2;
453                         vvdw_rep         = c12*VV;
454                         fijR             = c12*FF;
455                         fvdw             = -(fijD+fijR)*tabscale*rinv;
456                         vvdw             = vvdw_disp + vvdw_rep;
457                         break;
458
459                     default:
460                         gmx_fatal(FARGS,"Death & horror! No generic VdW interaction for ivdw=%d.\n",ivdw);
461                         break;
462                 }
463                 if(fr->vdw_modifier==eintmodPOTSWITCH)
464                 {
465                     d                = rsq*rinv-fr->rvdw_switch;
466                     d                = (d>0.0) ? d : 0.0;
467                     d2               = d*d;
468                     sw               = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
469                     dsw              = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
470                     /* See coulomb interaction for the force-switch formula */
471                     fvdw             = fvdw*sw - rinv*vvdw*dsw;
472                     vvdw            *= sw;
473                 }
474                 if(bExactVdwCutoff)
475                 {
476                     fvdw             = (rsq<=rvdw2) ? fvdw : 0.0;
477                     vvdw             = (rsq<=rvdw2) ? vvdw : 0.0;
478                 }
479                 vvdwtot         += vvdw;
480             } /* end VdW interactions */
481
482             fscal            = felec+fvdw;
483
484             if(!bCG && force_cap>0 && (fabs(fscal)> force_cap))
485             {
486                  fscal=force_cap*fscal/fabs(fscal);
487             }
488
489             fscal           *= hybscal;
490
491             tx               = fscal*dx;
492             ty               = fscal*dy;
493             tz               = fscal*dz;
494             fix              = fix + tx;
495             fiy              = fiy + ty;
496             fiz              = fiz + tz;
497             f[j3+0]          = f[j3+0] - tx;
498             f[j3+1]          = f[j3+1] - ty;
499             f[j3+2]          = f[j3+2] - tz;
500         }
501
502         f[ii3+0]         = f[ii3+0] + fix;
503         f[ii3+1]         = f[ii3+1] + fiy;
504         f[ii3+2]         = f[ii3+2] + fiz;
505         fshift[is3]      = fshift[is3]+fix;
506         fshift[is3+1]    = fshift[is3+1]+fiy;
507         fshift[is3+2]    = fshift[is3+2]+fiz;
508         ggid             = nlist->gid[n];
509         velecgrp[ggid]  += vctot;
510         vvdwgrp[ggid]   += vvdwtot;
511     }
512     /* Estimate flops, average for generic kernel:
513      * 12 flops per outer iteration
514      * 50 flops per inner iteration
515      */
516     inc_nrnb(nrnb,eNR_NBKERNEL_GENERIC,nlist->nri*12 + nlist->jindex[n]*50);
517 }