Code beautification with uncrustify
[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         {
231             continue;
232         }
233
234         for (k = nj0; (k < nj1); k++)
235         {
236             jnr              = nlist->jjnr[k];
237             weight_cg2       = wf[jnr];
238             weight_product   = weight_cg1*weight_cg2;
239
240             if (weight_product < ALMOST_ZERO)
241             {
242                 /* if it's a explicit loop, skip this atom */
243                 if (!bCG)
244                 {
245                     continue;
246                 }
247                 else /* if it's a coarse grained loop, include this atom */
248                 {
249                     hybscal = 1.0;
250                 }
251             }
252             else if (weight_product >= ALMOST_ONE)
253             {
254
255                 /* if it's a explicit loop, include this atom */
256                 if (!bCG)
257                 {
258                     hybscal = 1.0;
259                 }
260                 else  /* if it's a coarse grained loop, skip this atom */
261                 {
262                     continue;
263                 }
264             }
265             /* both have double identity, get hybrid scaling factor */
266             else
267             {
268                 hybscal = weight_product;
269
270                 if (bCG)
271                 {
272                     hybscal = 1.0 - hybscal;
273                 }
274             }
275
276             j3               = 3*jnr;
277             jx               = x[j3+0];
278             jy               = x[j3+1];
279             jz               = x[j3+2];
280             dx               = ix - jx;
281             dy               = iy - jy;
282             dz               = iz - jz;
283             rsq              = dx*dx+dy*dy+dz*dz;
284             rinv             = gmx_invsqrt(rsq);
285             rinvsq           = rinv*rinv;
286             felec            = 0;
287             fvdw             = 0;
288             velec            = 0;
289             vvdw             = 0;
290
291             if (bExactCutoff && rsq > rcutoff2)
292             {
293                 continue;
294             }
295
296             if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
297             {
298                 r                = rsq*rinv;
299                 rt               = r*tabscale;
300                 n0               = rt;
301                 eps              = rt-n0;
302                 eps2             = eps*eps;
303                 nnn              = table_nelements*n0;
304             }
305
306             /* Coulomb interaction. ielec==0 means no interaction */
307             if (ielec != GMX_NBKERNEL_ELEC_NONE)
308             {
309                 qq               = iq*charge[jnr];
310
311                 switch (ielec)
312                 {
313                     case GMX_NBKERNEL_ELEC_NONE:
314                         break;
315
316                     case GMX_NBKERNEL_ELEC_COULOMB:
317                         /* Vanilla cutoff coulomb */
318                         velec            = qq*rinv;
319                         felec            = velec*rinvsq;
320                         break;
321
322                     case GMX_NBKERNEL_ELEC_REACTIONFIELD:
323                         /* Reaction-field */
324                         velec            = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
325                         felec            = qq*(rinv*rinvsq-2.0*fr->k_rf);
326                         break;
327
328                     case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
329                         /* Tabulated coulomb */
330                         Y                = VFtab[nnn];
331                         F                = VFtab[nnn+1];
332                         Geps             = eps*VFtab[nnn+2];
333                         Heps2            = eps2*VFtab[nnn+3];
334                         Fp               = F+Geps+Heps2;
335                         VV               = Y+eps*Fp;
336                         FF               = Fp+Geps+2.0*Heps2;
337                         velec            = qq*VV;
338                         felec            = -qq*FF*tabscale*rinv;
339                         break;
340
341                     case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
342                         /* GB */
343                         gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
344                         break;
345
346                     case GMX_NBKERNEL_ELEC_EWALD:
347                         ewrt             = rsq*rinv*ewtabscale;
348                         ewitab           = ewrt;
349                         eweps            = ewrt-ewitab;
350                         ewitab           = 4*ewitab;
351                         felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
352                         rinvcorr         = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
353                         velec            = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
354                         felec            = qq*rinv*(rinvsq-felec);
355                         break;
356
357                     default:
358                         gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
359                         break;
360                 }
361                 if (fr->coulomb_modifier == eintmodPOTSWITCH)
362                 {
363                     d                = rsq*rinv-fr->rcoulomb_switch;
364                     d                = (d > 0.0) ? d : 0.0;
365                     d2               = d*d;
366                     sw               = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
367                     dsw              = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
368                     /* Apply switch function. Note that felec=f/r since it will be multiplied
369                      * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
370                      * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
371                      */
372                     felec            = felec*sw - rinv*velec*dsw;
373                     /* Once we have used velec to update felec we can modify velec too */
374                     velec           *= sw;
375                 }
376                 if (bExactElecCutoff)
377                 {
378                     felec            = (rsq <= rcoulomb2) ? felec : 0.0;
379                     velec            = (rsq <= rcoulomb2) ? velec : 0.0;
380                 }
381                 vctot           += velec;
382             } /* End of coulomb interactions */
383
384
385             /* VdW interaction. ivdw==0 means no interaction */
386             if (ivdw != GMX_NBKERNEL_VDW_NONE)
387             {
388                 tj               = nti+nvdwparam*type[jnr];
389
390                 switch (ivdw)
391                 {
392                     case GMX_NBKERNEL_VDW_NONE:
393                         break;
394
395                     case GMX_NBKERNEL_VDW_LENNARDJONES:
396                         /* Vanilla Lennard-Jones cutoff */
397                         c6               = vdwparam[tj];
398                         c12              = vdwparam[tj+1];
399                         rinvsix          = rinvsq*rinvsq*rinvsq;
400                         vvdw_disp        = c6*rinvsix;
401                         vvdw_rep         = c12*rinvsix*rinvsix;
402                         fvdw             = (vvdw_rep-vvdw_disp)*rinvsq;
403                         if (fr->vdw_modifier == eintmodPOTSHIFT)
404                         {
405                             vvdw             = (vvdw_rep-c12*sh_invrc6*sh_invrc6)*(1.0/12.0)-(vvdw_disp-c6*sh_invrc6)*(1.0/6.0);
406                         }
407                         else
408                         {
409                             vvdw             = vvdw_rep/12.0-vvdw_disp/6.0;
410                         }
411                         break;
412
413                     case GMX_NBKERNEL_VDW_BUCKINGHAM:
414                         /* Buckingham */
415                         c6               = vdwparam[tj];
416                         cexp1            = vdwparam[tj+1];
417                         cexp2            = vdwparam[tj+2];
418
419                         rinvsix          = rinvsq*rinvsq*rinvsq;
420                         vvdw_disp        = c6*rinvsix;
421                         br               = cexp2*rsq*rinv;
422                         vvdw_rep         = cexp1*exp(-br);
423                         fvdw             = (br*vvdw_rep-vvdw_disp)*rinvsq;
424                         if (fr->vdw_modifier == eintmodPOTSHIFT)
425                         {
426                             vvdw             = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp-c6*sh_invrc6)/6.0;
427                         }
428                         else
429                         {
430                             vvdw             = vvdw_rep-vvdw_disp/6.0;
431                         }
432                         break;
433
434                     case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
435                         /* Tabulated VdW */
436                         c6               = vdwparam[tj];
437                         c12              = vdwparam[tj+1];
438                         Y                = VFtab[nnn+4];
439                         F                = VFtab[nnn+5];
440                         Geps             = eps*VFtab[nnn+6];
441                         Heps2            = eps2*VFtab[nnn+7];
442                         Fp               = F+Geps+Heps2;
443                         VV               = Y+eps*Fp;
444                         FF               = Fp+Geps+2.0*Heps2;
445                         vvdw_disp        = c6*VV;
446                         fijD             = c6*FF;
447                         Y                = VFtab[nnn+8];
448                         F                = VFtab[nnn+9];
449                         Geps             = eps*VFtab[nnn+10];
450                         Heps2            = eps2*VFtab[nnn+11];
451                         Fp               = F+Geps+Heps2;
452                         VV               = Y+eps*Fp;
453                         FF               = Fp+Geps+2.0*Heps2;
454                         vvdw_rep         = c12*VV;
455                         fijR             = c12*FF;
456                         fvdw             = -(fijD+fijR)*tabscale*rinv;
457                         vvdw             = vvdw_disp + vvdw_rep;
458                         break;
459
460                     default:
461                         gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
462                         break;
463                 }
464                 if (fr->vdw_modifier == eintmodPOTSWITCH)
465                 {
466                     d                = rsq*rinv-fr->rvdw_switch;
467                     d                = (d > 0.0) ? d : 0.0;
468                     d2               = d*d;
469                     sw               = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
470                     dsw              = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
471                     /* See coulomb interaction for the force-switch formula */
472                     fvdw             = fvdw*sw - rinv*vvdw*dsw;
473                     vvdw            *= sw;
474                 }
475                 if (bExactVdwCutoff)
476                 {
477                     fvdw             = (rsq <= rvdw2) ? fvdw : 0.0;
478                     vvdw             = (rsq <= rvdw2) ? vvdw : 0.0;
479                 }
480                 vvdwtot         += vvdw;
481             } /* end VdW interactions */
482
483             fscal            = felec+fvdw;
484
485             if (!bCG && force_cap > 0 && (fabs(fscal) > force_cap))
486             {
487                 fscal = force_cap*fscal/fabs(fscal);
488             }
489
490             fscal           *= hybscal;
491
492             tx               = fscal*dx;
493             ty               = fscal*dy;
494             tz               = fscal*dz;
495             fix              = fix + tx;
496             fiy              = fiy + ty;
497             fiz              = fiz + tz;
498             f[j3+0]          = f[j3+0] - tx;
499             f[j3+1]          = f[j3+1] - ty;
500             f[j3+2]          = f[j3+2] - tz;
501         }
502
503         f[ii3+0]         = f[ii3+0] + fix;
504         f[ii3+1]         = f[ii3+1] + fiy;
505         f[ii3+2]         = f[ii3+2] + fiz;
506         fshift[is3]      = fshift[is3]+fix;
507         fshift[is3+1]    = fshift[is3+1]+fiy;
508         fshift[is3+2]    = fshift[is3+2]+fiz;
509         ggid             = nlist->gid[n];
510         velecgrp[ggid]  += vctot;
511         vvdwgrp[ggid]   += vvdwtot;
512     }
513     /* Estimate flops, average for generic kernel:
514      * 12 flops per outer iteration
515      * 50 flops per inner iteration
516      */
517     inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);
518 }