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