Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEwSh_VdwLJSh_GeomW3P1_c.c
blobed948dcc430560232c9e319eab7caf17f2ab8ccc
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS c kernel generator.
38 #include "config.h"
40 #include <math.h>
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
48 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJSh_GeomW3P1_VF_c
49 * Electrostatics interaction: Ewald
50 * VdW interaction: LennardJones
51 * Geometry: Water3-Particle
52 * Calculate force/pot: PotentialAndForce
54 void
55 nb_kernel_ElecEwSh_VdwLJSh_GeomW3P1_VF_c
56 (t_nblist * gmx_restrict nlist,
57 rvec * gmx_restrict xx,
58 rvec * gmx_restrict ff,
59 t_forcerec * gmx_restrict fr,
60 t_mdatoms * gmx_restrict mdatoms,
61 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
62 t_nrnb * gmx_restrict nrnb)
64 int i_shift_offset,i_coord_offset,j_coord_offset;
65 int j_index_start,j_index_end;
66 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
67 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
68 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
69 real *shiftvec,*fshift,*x,*f;
70 int vdwioffset0;
71 real ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
72 int vdwioffset1;
73 real ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
74 int vdwioffset2;
75 real ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
76 int vdwjidx0;
77 real jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
78 real dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
79 real dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
80 real dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
81 real velec,felec,velecsum,facel,crf,krf,krf2;
82 real *charge;
83 int nvdwtype;
84 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
85 int *vdwtype;
86 real *vdwparam;
87 int ewitab;
88 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
89 real *ewtab;
91 x = xx[0];
92 f = ff[0];
94 nri = nlist->nri;
95 iinr = nlist->iinr;
96 jindex = nlist->jindex;
97 jjnr = nlist->jjnr;
98 shiftidx = nlist->shift;
99 gid = nlist->gid;
100 shiftvec = fr->shift_vec[0];
101 fshift = fr->fshift[0];
102 facel = fr->epsfac;
103 charge = mdatoms->chargeA;
104 nvdwtype = fr->ntype;
105 vdwparam = fr->nbfp;
106 vdwtype = mdatoms->typeA;
108 sh_ewald = fr->ic->sh_ewald;
109 ewtab = fr->ic->tabq_coul_FDV0;
110 ewtabscale = fr->ic->tabq_scale;
111 ewtabhalfspace = 0.5/ewtabscale;
113 /* Setup water-specific parameters */
114 inr = nlist->iinr[0];
115 iq0 = facel*charge[inr+0];
116 iq1 = facel*charge[inr+1];
117 iq2 = facel*charge[inr+2];
118 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
120 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
121 rcutoff = fr->rcoulomb;
122 rcutoff2 = rcutoff*rcutoff;
124 sh_vdw_invrcut6 = fr->ic->sh_invrc6;
125 rvdw = fr->rvdw;
127 outeriter = 0;
128 inneriter = 0;
130 /* Start outer loop over neighborlists */
131 for(iidx=0; iidx<nri; iidx++)
133 /* Load shift vector for this list */
134 i_shift_offset = DIM*shiftidx[iidx];
135 shX = shiftvec[i_shift_offset+XX];
136 shY = shiftvec[i_shift_offset+YY];
137 shZ = shiftvec[i_shift_offset+ZZ];
139 /* Load limits for loop over neighbors */
140 j_index_start = jindex[iidx];
141 j_index_end = jindex[iidx+1];
143 /* Get outer coordinate index */
144 inr = iinr[iidx];
145 i_coord_offset = DIM*inr;
147 /* Load i particle coords and add shift vector */
148 ix0 = shX + x[i_coord_offset+DIM*0+XX];
149 iy0 = shY + x[i_coord_offset+DIM*0+YY];
150 iz0 = shZ + x[i_coord_offset+DIM*0+ZZ];
151 ix1 = shX + x[i_coord_offset+DIM*1+XX];
152 iy1 = shY + x[i_coord_offset+DIM*1+YY];
153 iz1 = shZ + x[i_coord_offset+DIM*1+ZZ];
154 ix2 = shX + x[i_coord_offset+DIM*2+XX];
155 iy2 = shY + x[i_coord_offset+DIM*2+YY];
156 iz2 = shZ + x[i_coord_offset+DIM*2+ZZ];
158 fix0 = 0.0;
159 fiy0 = 0.0;
160 fiz0 = 0.0;
161 fix1 = 0.0;
162 fiy1 = 0.0;
163 fiz1 = 0.0;
164 fix2 = 0.0;
165 fiy2 = 0.0;
166 fiz2 = 0.0;
168 /* Reset potential sums */
169 velecsum = 0.0;
170 vvdwsum = 0.0;
172 /* Start inner kernel loop */
173 for(jidx=j_index_start; jidx<j_index_end; jidx++)
175 /* Get j neighbor index, and coordinate index */
176 jnr = jjnr[jidx];
177 j_coord_offset = DIM*jnr;
179 /* load j atom coordinates */
180 jx0 = x[j_coord_offset+DIM*0+XX];
181 jy0 = x[j_coord_offset+DIM*0+YY];
182 jz0 = x[j_coord_offset+DIM*0+ZZ];
184 /* Calculate displacement vector */
185 dx00 = ix0 - jx0;
186 dy00 = iy0 - jy0;
187 dz00 = iz0 - jz0;
188 dx10 = ix1 - jx0;
189 dy10 = iy1 - jy0;
190 dz10 = iz1 - jz0;
191 dx20 = ix2 - jx0;
192 dy20 = iy2 - jy0;
193 dz20 = iz2 - jz0;
195 /* Calculate squared distance and things based on it */
196 rsq00 = dx00*dx00+dy00*dy00+dz00*dz00;
197 rsq10 = dx10*dx10+dy10*dy10+dz10*dz10;
198 rsq20 = dx20*dx20+dy20*dy20+dz20*dz20;
200 rinv00 = gmx_invsqrt(rsq00);
201 rinv10 = gmx_invsqrt(rsq10);
202 rinv20 = gmx_invsqrt(rsq20);
204 rinvsq00 = rinv00*rinv00;
205 rinvsq10 = rinv10*rinv10;
206 rinvsq20 = rinv20*rinv20;
208 /* Load parameters for j particles */
209 jq0 = charge[jnr+0];
210 vdwjidx0 = 2*vdwtype[jnr+0];
212 /**************************
213 * CALCULATE INTERACTIONS *
214 **************************/
216 if (rsq00<rcutoff2)
219 r00 = rsq00*rinv00;
221 qq00 = iq0*jq0;
222 c6_00 = vdwparam[vdwioffset0+vdwjidx0];
223 c12_00 = vdwparam[vdwioffset0+vdwjidx0+1];
225 /* EWALD ELECTROSTATICS */
227 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
228 ewrt = r00*ewtabscale;
229 ewitab = ewrt;
230 eweps = ewrt-ewitab;
231 ewitab = 4*ewitab;
232 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
233 velec = qq00*((rinv00-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
234 felec = qq00*rinv00*(rinvsq00-felec);
236 /* LENNARD-JONES DISPERSION/REPULSION */
238 rinvsix = rinvsq00*rinvsq00*rinvsq00;
239 vvdw6 = c6_00*rinvsix;
240 vvdw12 = c12_00*rinvsix*rinvsix;
241 vvdw = (vvdw12 - c12_00*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_00*sh_vdw_invrcut6)*(1.0/6.0);
242 fvdw = (vvdw12-vvdw6)*rinvsq00;
244 /* Update potential sums from outer loop */
245 velecsum += velec;
246 vvdwsum += vvdw;
248 fscal = felec+fvdw;
250 /* Calculate temporary vectorial force */
251 tx = fscal*dx00;
252 ty = fscal*dy00;
253 tz = fscal*dz00;
255 /* Update vectorial force */
256 fix0 += tx;
257 fiy0 += ty;
258 fiz0 += tz;
259 f[j_coord_offset+DIM*0+XX] -= tx;
260 f[j_coord_offset+DIM*0+YY] -= ty;
261 f[j_coord_offset+DIM*0+ZZ] -= tz;
265 /**************************
266 * CALCULATE INTERACTIONS *
267 **************************/
269 if (rsq10<rcutoff2)
272 r10 = rsq10*rinv10;
274 qq10 = iq1*jq0;
276 /* EWALD ELECTROSTATICS */
278 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
279 ewrt = r10*ewtabscale;
280 ewitab = ewrt;
281 eweps = ewrt-ewitab;
282 ewitab = 4*ewitab;
283 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
284 velec = qq10*((rinv10-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
285 felec = qq10*rinv10*(rinvsq10-felec);
287 /* Update potential sums from outer loop */
288 velecsum += velec;
290 fscal = felec;
292 /* Calculate temporary vectorial force */
293 tx = fscal*dx10;
294 ty = fscal*dy10;
295 tz = fscal*dz10;
297 /* Update vectorial force */
298 fix1 += tx;
299 fiy1 += ty;
300 fiz1 += tz;
301 f[j_coord_offset+DIM*0+XX] -= tx;
302 f[j_coord_offset+DIM*0+YY] -= ty;
303 f[j_coord_offset+DIM*0+ZZ] -= tz;
307 /**************************
308 * CALCULATE INTERACTIONS *
309 **************************/
311 if (rsq20<rcutoff2)
314 r20 = rsq20*rinv20;
316 qq20 = iq2*jq0;
318 /* EWALD ELECTROSTATICS */
320 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
321 ewrt = r20*ewtabscale;
322 ewitab = ewrt;
323 eweps = ewrt-ewitab;
324 ewitab = 4*ewitab;
325 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
326 velec = qq20*((rinv20-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
327 felec = qq20*rinv20*(rinvsq20-felec);
329 /* Update potential sums from outer loop */
330 velecsum += velec;
332 fscal = felec;
334 /* Calculate temporary vectorial force */
335 tx = fscal*dx20;
336 ty = fscal*dy20;
337 tz = fscal*dz20;
339 /* Update vectorial force */
340 fix2 += tx;
341 fiy2 += ty;
342 fiz2 += tz;
343 f[j_coord_offset+DIM*0+XX] -= tx;
344 f[j_coord_offset+DIM*0+YY] -= ty;
345 f[j_coord_offset+DIM*0+ZZ] -= tz;
349 /* Inner loop uses 143 flops */
351 /* End of innermost loop */
353 tx = ty = tz = 0;
354 f[i_coord_offset+DIM*0+XX] += fix0;
355 f[i_coord_offset+DIM*0+YY] += fiy0;
356 f[i_coord_offset+DIM*0+ZZ] += fiz0;
357 tx += fix0;
358 ty += fiy0;
359 tz += fiz0;
360 f[i_coord_offset+DIM*1+XX] += fix1;
361 f[i_coord_offset+DIM*1+YY] += fiy1;
362 f[i_coord_offset+DIM*1+ZZ] += fiz1;
363 tx += fix1;
364 ty += fiy1;
365 tz += fiz1;
366 f[i_coord_offset+DIM*2+XX] += fix2;
367 f[i_coord_offset+DIM*2+YY] += fiy2;
368 f[i_coord_offset+DIM*2+ZZ] += fiz2;
369 tx += fix2;
370 ty += fiy2;
371 tz += fiz2;
372 fshift[i_shift_offset+XX] += tx;
373 fshift[i_shift_offset+YY] += ty;
374 fshift[i_shift_offset+ZZ] += tz;
376 ggid = gid[iidx];
377 /* Update potential energies */
378 kernel_data->energygrp_elec[ggid] += velecsum;
379 kernel_data->energygrp_vdw[ggid] += vvdwsum;
381 /* Increment number of inner iterations */
382 inneriter += j_index_end - j_index_start;
384 /* Outer loop uses 32 flops */
387 /* Increment number of outer iterations */
388 outeriter += nri;
390 /* Update outer/inner flops */
392 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*32 + inneriter*143);
395 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJSh_GeomW3P1_F_c
396 * Electrostatics interaction: Ewald
397 * VdW interaction: LennardJones
398 * Geometry: Water3-Particle
399 * Calculate force/pot: Force
401 void
402 nb_kernel_ElecEwSh_VdwLJSh_GeomW3P1_F_c
403 (t_nblist * gmx_restrict nlist,
404 rvec * gmx_restrict xx,
405 rvec * gmx_restrict ff,
406 t_forcerec * gmx_restrict fr,
407 t_mdatoms * gmx_restrict mdatoms,
408 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
409 t_nrnb * gmx_restrict nrnb)
411 int i_shift_offset,i_coord_offset,j_coord_offset;
412 int j_index_start,j_index_end;
413 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
414 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
415 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
416 real *shiftvec,*fshift,*x,*f;
417 int vdwioffset0;
418 real ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
419 int vdwioffset1;
420 real ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
421 int vdwioffset2;
422 real ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
423 int vdwjidx0;
424 real jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
425 real dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
426 real dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
427 real dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
428 real velec,felec,velecsum,facel,crf,krf,krf2;
429 real *charge;
430 int nvdwtype;
431 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
432 int *vdwtype;
433 real *vdwparam;
434 int ewitab;
435 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
436 real *ewtab;
438 x = xx[0];
439 f = ff[0];
441 nri = nlist->nri;
442 iinr = nlist->iinr;
443 jindex = nlist->jindex;
444 jjnr = nlist->jjnr;
445 shiftidx = nlist->shift;
446 gid = nlist->gid;
447 shiftvec = fr->shift_vec[0];
448 fshift = fr->fshift[0];
449 facel = fr->epsfac;
450 charge = mdatoms->chargeA;
451 nvdwtype = fr->ntype;
452 vdwparam = fr->nbfp;
453 vdwtype = mdatoms->typeA;
455 sh_ewald = fr->ic->sh_ewald;
456 ewtab = fr->ic->tabq_coul_F;
457 ewtabscale = fr->ic->tabq_scale;
458 ewtabhalfspace = 0.5/ewtabscale;
460 /* Setup water-specific parameters */
461 inr = nlist->iinr[0];
462 iq0 = facel*charge[inr+0];
463 iq1 = facel*charge[inr+1];
464 iq2 = facel*charge[inr+2];
465 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
467 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
468 rcutoff = fr->rcoulomb;
469 rcutoff2 = rcutoff*rcutoff;
471 sh_vdw_invrcut6 = fr->ic->sh_invrc6;
472 rvdw = fr->rvdw;
474 outeriter = 0;
475 inneriter = 0;
477 /* Start outer loop over neighborlists */
478 for(iidx=0; iidx<nri; iidx++)
480 /* Load shift vector for this list */
481 i_shift_offset = DIM*shiftidx[iidx];
482 shX = shiftvec[i_shift_offset+XX];
483 shY = shiftvec[i_shift_offset+YY];
484 shZ = shiftvec[i_shift_offset+ZZ];
486 /* Load limits for loop over neighbors */
487 j_index_start = jindex[iidx];
488 j_index_end = jindex[iidx+1];
490 /* Get outer coordinate index */
491 inr = iinr[iidx];
492 i_coord_offset = DIM*inr;
494 /* Load i particle coords and add shift vector */
495 ix0 = shX + x[i_coord_offset+DIM*0+XX];
496 iy0 = shY + x[i_coord_offset+DIM*0+YY];
497 iz0 = shZ + x[i_coord_offset+DIM*0+ZZ];
498 ix1 = shX + x[i_coord_offset+DIM*1+XX];
499 iy1 = shY + x[i_coord_offset+DIM*1+YY];
500 iz1 = shZ + x[i_coord_offset+DIM*1+ZZ];
501 ix2 = shX + x[i_coord_offset+DIM*2+XX];
502 iy2 = shY + x[i_coord_offset+DIM*2+YY];
503 iz2 = shZ + x[i_coord_offset+DIM*2+ZZ];
505 fix0 = 0.0;
506 fiy0 = 0.0;
507 fiz0 = 0.0;
508 fix1 = 0.0;
509 fiy1 = 0.0;
510 fiz1 = 0.0;
511 fix2 = 0.0;
512 fiy2 = 0.0;
513 fiz2 = 0.0;
515 /* Start inner kernel loop */
516 for(jidx=j_index_start; jidx<j_index_end; jidx++)
518 /* Get j neighbor index, and coordinate index */
519 jnr = jjnr[jidx];
520 j_coord_offset = DIM*jnr;
522 /* load j atom coordinates */
523 jx0 = x[j_coord_offset+DIM*0+XX];
524 jy0 = x[j_coord_offset+DIM*0+YY];
525 jz0 = x[j_coord_offset+DIM*0+ZZ];
527 /* Calculate displacement vector */
528 dx00 = ix0 - jx0;
529 dy00 = iy0 - jy0;
530 dz00 = iz0 - jz0;
531 dx10 = ix1 - jx0;
532 dy10 = iy1 - jy0;
533 dz10 = iz1 - jz0;
534 dx20 = ix2 - jx0;
535 dy20 = iy2 - jy0;
536 dz20 = iz2 - jz0;
538 /* Calculate squared distance and things based on it */
539 rsq00 = dx00*dx00+dy00*dy00+dz00*dz00;
540 rsq10 = dx10*dx10+dy10*dy10+dz10*dz10;
541 rsq20 = dx20*dx20+dy20*dy20+dz20*dz20;
543 rinv00 = gmx_invsqrt(rsq00);
544 rinv10 = gmx_invsqrt(rsq10);
545 rinv20 = gmx_invsqrt(rsq20);
547 rinvsq00 = rinv00*rinv00;
548 rinvsq10 = rinv10*rinv10;
549 rinvsq20 = rinv20*rinv20;
551 /* Load parameters for j particles */
552 jq0 = charge[jnr+0];
553 vdwjidx0 = 2*vdwtype[jnr+0];
555 /**************************
556 * CALCULATE INTERACTIONS *
557 **************************/
559 if (rsq00<rcutoff2)
562 r00 = rsq00*rinv00;
564 qq00 = iq0*jq0;
565 c6_00 = vdwparam[vdwioffset0+vdwjidx0];
566 c12_00 = vdwparam[vdwioffset0+vdwjidx0+1];
568 /* EWALD ELECTROSTATICS */
570 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
571 ewrt = r00*ewtabscale;
572 ewitab = ewrt;
573 eweps = ewrt-ewitab;
574 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
575 felec = qq00*rinv00*(rinvsq00-felec);
577 /* LENNARD-JONES DISPERSION/REPULSION */
579 rinvsix = rinvsq00*rinvsq00*rinvsq00;
580 fvdw = (c12_00*rinvsix-c6_00)*rinvsix*rinvsq00;
582 fscal = felec+fvdw;
584 /* Calculate temporary vectorial force */
585 tx = fscal*dx00;
586 ty = fscal*dy00;
587 tz = fscal*dz00;
589 /* Update vectorial force */
590 fix0 += tx;
591 fiy0 += ty;
592 fiz0 += tz;
593 f[j_coord_offset+DIM*0+XX] -= tx;
594 f[j_coord_offset+DIM*0+YY] -= ty;
595 f[j_coord_offset+DIM*0+ZZ] -= tz;
599 /**************************
600 * CALCULATE INTERACTIONS *
601 **************************/
603 if (rsq10<rcutoff2)
606 r10 = rsq10*rinv10;
608 qq10 = iq1*jq0;
610 /* EWALD ELECTROSTATICS */
612 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
613 ewrt = r10*ewtabscale;
614 ewitab = ewrt;
615 eweps = ewrt-ewitab;
616 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
617 felec = qq10*rinv10*(rinvsq10-felec);
619 fscal = felec;
621 /* Calculate temporary vectorial force */
622 tx = fscal*dx10;
623 ty = fscal*dy10;
624 tz = fscal*dz10;
626 /* Update vectorial force */
627 fix1 += tx;
628 fiy1 += ty;
629 fiz1 += tz;
630 f[j_coord_offset+DIM*0+XX] -= tx;
631 f[j_coord_offset+DIM*0+YY] -= ty;
632 f[j_coord_offset+DIM*0+ZZ] -= tz;
636 /**************************
637 * CALCULATE INTERACTIONS *
638 **************************/
640 if (rsq20<rcutoff2)
643 r20 = rsq20*rinv20;
645 qq20 = iq2*jq0;
647 /* EWALD ELECTROSTATICS */
649 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
650 ewrt = r20*ewtabscale;
651 ewitab = ewrt;
652 eweps = ewrt-ewitab;
653 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
654 felec = qq20*rinv20*(rinvsq20-felec);
656 fscal = felec;
658 /* Calculate temporary vectorial force */
659 tx = fscal*dx20;
660 ty = fscal*dy20;
661 tz = fscal*dz20;
663 /* Update vectorial force */
664 fix2 += tx;
665 fiy2 += ty;
666 fiz2 += tz;
667 f[j_coord_offset+DIM*0+XX] -= tx;
668 f[j_coord_offset+DIM*0+YY] -= ty;
669 f[j_coord_offset+DIM*0+ZZ] -= tz;
673 /* Inner loop uses 109 flops */
675 /* End of innermost loop */
677 tx = ty = tz = 0;
678 f[i_coord_offset+DIM*0+XX] += fix0;
679 f[i_coord_offset+DIM*0+YY] += fiy0;
680 f[i_coord_offset+DIM*0+ZZ] += fiz0;
681 tx += fix0;
682 ty += fiy0;
683 tz += fiz0;
684 f[i_coord_offset+DIM*1+XX] += fix1;
685 f[i_coord_offset+DIM*1+YY] += fiy1;
686 f[i_coord_offset+DIM*1+ZZ] += fiz1;
687 tx += fix1;
688 ty += fiy1;
689 tz += fiz1;
690 f[i_coord_offset+DIM*2+XX] += fix2;
691 f[i_coord_offset+DIM*2+YY] += fiy2;
692 f[i_coord_offset+DIM*2+ZZ] += fiz2;
693 tx += fix2;
694 ty += fiy2;
695 tz += fiz2;
696 fshift[i_shift_offset+XX] += tx;
697 fshift[i_shift_offset+YY] += ty;
698 fshift[i_shift_offset+ZZ] += tz;
700 /* Increment number of inner iterations */
701 inneriter += j_index_end - j_index_start;
703 /* Outer loop uses 30 flops */
706 /* Increment number of outer iterations */
707 outeriter += nri;
709 /* Update outer/inner flops */
711 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*30 + inneriter*109);