Remove nb-parameters from t_forcerec
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEw_VdwNone_GeomW3P1_c.c
blob15b5761e142a44039044f512d88c01f55b35efdf
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014.2015,2017, 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 "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
48 * Gromacs nonbonded kernel: nb_kernel_ElecEw_VdwNone_GeomW3P1_VF_c
49 * Electrostatics interaction: Ewald
50 * VdW interaction: None
51 * Geometry: Water3-Particle
52 * Calculate force/pot: PotentialAndForce
54 void
55 nb_kernel_ElecEw_VdwNone_GeomW3P1_VF_c
56 (t_nblist * gmx_restrict nlist,
57 rvec * gmx_restrict xx,
58 rvec * gmx_restrict ff,
59 struct 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 ewitab;
84 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
85 real *ewtab;
87 x = xx[0];
88 f = ff[0];
90 nri = nlist->nri;
91 iinr = nlist->iinr;
92 jindex = nlist->jindex;
93 jjnr = nlist->jjnr;
94 shiftidx = nlist->shift;
95 gid = nlist->gid;
96 shiftvec = fr->shift_vec[0];
97 fshift = fr->fshift[0];
98 facel = fr->ic->epsfac;
99 charge = mdatoms->chargeA;
101 sh_ewald = fr->ic->sh_ewald;
102 ewtab = fr->ic->tabq_coul_FDV0;
103 ewtabscale = fr->ic->tabq_scale;
104 ewtabhalfspace = 0.5/ewtabscale;
106 /* Setup water-specific parameters */
107 inr = nlist->iinr[0];
108 iq0 = facel*charge[inr+0];
109 iq1 = facel*charge[inr+1];
110 iq2 = facel*charge[inr+2];
112 outeriter = 0;
113 inneriter = 0;
115 /* Start outer loop over neighborlists */
116 for(iidx=0; iidx<nri; iidx++)
118 /* Load shift vector for this list */
119 i_shift_offset = DIM*shiftidx[iidx];
120 shX = shiftvec[i_shift_offset+XX];
121 shY = shiftvec[i_shift_offset+YY];
122 shZ = shiftvec[i_shift_offset+ZZ];
124 /* Load limits for loop over neighbors */
125 j_index_start = jindex[iidx];
126 j_index_end = jindex[iidx+1];
128 /* Get outer coordinate index */
129 inr = iinr[iidx];
130 i_coord_offset = DIM*inr;
132 /* Load i particle coords and add shift vector */
133 ix0 = shX + x[i_coord_offset+DIM*0+XX];
134 iy0 = shY + x[i_coord_offset+DIM*0+YY];
135 iz0 = shZ + x[i_coord_offset+DIM*0+ZZ];
136 ix1 = shX + x[i_coord_offset+DIM*1+XX];
137 iy1 = shY + x[i_coord_offset+DIM*1+YY];
138 iz1 = shZ + x[i_coord_offset+DIM*1+ZZ];
139 ix2 = shX + x[i_coord_offset+DIM*2+XX];
140 iy2 = shY + x[i_coord_offset+DIM*2+YY];
141 iz2 = shZ + x[i_coord_offset+DIM*2+ZZ];
143 fix0 = 0.0;
144 fiy0 = 0.0;
145 fiz0 = 0.0;
146 fix1 = 0.0;
147 fiy1 = 0.0;
148 fiz1 = 0.0;
149 fix2 = 0.0;
150 fiy2 = 0.0;
151 fiz2 = 0.0;
153 /* Reset potential sums */
154 velecsum = 0.0;
156 /* Start inner kernel loop */
157 for(jidx=j_index_start; jidx<j_index_end; jidx++)
159 /* Get j neighbor index, and coordinate index */
160 jnr = jjnr[jidx];
161 j_coord_offset = DIM*jnr;
163 /* load j atom coordinates */
164 jx0 = x[j_coord_offset+DIM*0+XX];
165 jy0 = x[j_coord_offset+DIM*0+YY];
166 jz0 = x[j_coord_offset+DIM*0+ZZ];
168 /* Calculate displacement vector */
169 dx00 = ix0 - jx0;
170 dy00 = iy0 - jy0;
171 dz00 = iz0 - jz0;
172 dx10 = ix1 - jx0;
173 dy10 = iy1 - jy0;
174 dz10 = iz1 - jz0;
175 dx20 = ix2 - jx0;
176 dy20 = iy2 - jy0;
177 dz20 = iz2 - jz0;
179 /* Calculate squared distance and things based on it */
180 rsq00 = dx00*dx00+dy00*dy00+dz00*dz00;
181 rsq10 = dx10*dx10+dy10*dy10+dz10*dz10;
182 rsq20 = dx20*dx20+dy20*dy20+dz20*dz20;
184 rinv00 = 1.0/sqrt(rsq00);
185 rinv10 = 1.0/sqrt(rsq10);
186 rinv20 = 1.0/sqrt(rsq20);
188 rinvsq00 = rinv00*rinv00;
189 rinvsq10 = rinv10*rinv10;
190 rinvsq20 = rinv20*rinv20;
192 /* Load parameters for j particles */
193 jq0 = charge[jnr+0];
195 /**************************
196 * CALCULATE INTERACTIONS *
197 **************************/
199 r00 = rsq00*rinv00;
201 qq00 = iq0*jq0;
203 /* EWALD ELECTROSTATICS */
205 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
206 ewrt = r00*ewtabscale;
207 ewitab = ewrt;
208 eweps = ewrt-ewitab;
209 ewitab = 4*ewitab;
210 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
211 velec = qq00*(rinv00-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
212 felec = qq00*rinv00*(rinvsq00-felec);
214 /* Update potential sums from outer loop */
215 velecsum += velec;
217 fscal = felec;
219 /* Calculate temporary vectorial force */
220 tx = fscal*dx00;
221 ty = fscal*dy00;
222 tz = fscal*dz00;
224 /* Update vectorial force */
225 fix0 += tx;
226 fiy0 += ty;
227 fiz0 += tz;
228 f[j_coord_offset+DIM*0+XX] -= tx;
229 f[j_coord_offset+DIM*0+YY] -= ty;
230 f[j_coord_offset+DIM*0+ZZ] -= tz;
232 /**************************
233 * CALCULATE INTERACTIONS *
234 **************************/
236 r10 = rsq10*rinv10;
238 qq10 = iq1*jq0;
240 /* EWALD ELECTROSTATICS */
242 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
243 ewrt = r10*ewtabscale;
244 ewitab = ewrt;
245 eweps = ewrt-ewitab;
246 ewitab = 4*ewitab;
247 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
248 velec = qq10*(rinv10-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
249 felec = qq10*rinv10*(rinvsq10-felec);
251 /* Update potential sums from outer loop */
252 velecsum += velec;
254 fscal = felec;
256 /* Calculate temporary vectorial force */
257 tx = fscal*dx10;
258 ty = fscal*dy10;
259 tz = fscal*dz10;
261 /* Update vectorial force */
262 fix1 += tx;
263 fiy1 += ty;
264 fiz1 += tz;
265 f[j_coord_offset+DIM*0+XX] -= tx;
266 f[j_coord_offset+DIM*0+YY] -= ty;
267 f[j_coord_offset+DIM*0+ZZ] -= tz;
269 /**************************
270 * CALCULATE INTERACTIONS *
271 **************************/
273 r20 = rsq20*rinv20;
275 qq20 = iq2*jq0;
277 /* EWALD ELECTROSTATICS */
279 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
280 ewrt = r20*ewtabscale;
281 ewitab = ewrt;
282 eweps = ewrt-ewitab;
283 ewitab = 4*ewitab;
284 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
285 velec = qq20*(rinv20-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
286 felec = qq20*rinv20*(rinvsq20-felec);
288 /* Update potential sums from outer loop */
289 velecsum += velec;
291 fscal = felec;
293 /* Calculate temporary vectorial force */
294 tx = fscal*dx20;
295 ty = fscal*dy20;
296 tz = fscal*dz20;
298 /* Update vectorial force */
299 fix2 += tx;
300 fiy2 += ty;
301 fiz2 += tz;
302 f[j_coord_offset+DIM*0+XX] -= tx;
303 f[j_coord_offset+DIM*0+YY] -= ty;
304 f[j_coord_offset+DIM*0+ZZ] -= tz;
306 /* Inner loop uses 123 flops */
308 /* End of innermost loop */
310 tx = ty = tz = 0;
311 f[i_coord_offset+DIM*0+XX] += fix0;
312 f[i_coord_offset+DIM*0+YY] += fiy0;
313 f[i_coord_offset+DIM*0+ZZ] += fiz0;
314 tx += fix0;
315 ty += fiy0;
316 tz += fiz0;
317 f[i_coord_offset+DIM*1+XX] += fix1;
318 f[i_coord_offset+DIM*1+YY] += fiy1;
319 f[i_coord_offset+DIM*1+ZZ] += fiz1;
320 tx += fix1;
321 ty += fiy1;
322 tz += fiz1;
323 f[i_coord_offset+DIM*2+XX] += fix2;
324 f[i_coord_offset+DIM*2+YY] += fiy2;
325 f[i_coord_offset+DIM*2+ZZ] += fiz2;
326 tx += fix2;
327 ty += fiy2;
328 tz += fiz2;
329 fshift[i_shift_offset+XX] += tx;
330 fshift[i_shift_offset+YY] += ty;
331 fshift[i_shift_offset+ZZ] += tz;
333 ggid = gid[iidx];
334 /* Update potential energies */
335 kernel_data->energygrp_elec[ggid] += velecsum;
337 /* Increment number of inner iterations */
338 inneriter += j_index_end - j_index_start;
340 /* Outer loop uses 31 flops */
343 /* Increment number of outer iterations */
344 outeriter += nri;
346 /* Update outer/inner flops */
348 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_VF,outeriter*31 + inneriter*123);
351 * Gromacs nonbonded kernel: nb_kernel_ElecEw_VdwNone_GeomW3P1_F_c
352 * Electrostatics interaction: Ewald
353 * VdW interaction: None
354 * Geometry: Water3-Particle
355 * Calculate force/pot: Force
357 void
358 nb_kernel_ElecEw_VdwNone_GeomW3P1_F_c
359 (t_nblist * gmx_restrict nlist,
360 rvec * gmx_restrict xx,
361 rvec * gmx_restrict ff,
362 struct t_forcerec * gmx_restrict fr,
363 t_mdatoms * gmx_restrict mdatoms,
364 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
365 t_nrnb * gmx_restrict nrnb)
367 int i_shift_offset,i_coord_offset,j_coord_offset;
368 int j_index_start,j_index_end;
369 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
370 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
371 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
372 real *shiftvec,*fshift,*x,*f;
373 int vdwioffset0;
374 real ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
375 int vdwioffset1;
376 real ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
377 int vdwioffset2;
378 real ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
379 int vdwjidx0;
380 real jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
381 real dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
382 real dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
383 real dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
384 real velec,felec,velecsum,facel,crf,krf,krf2;
385 real *charge;
386 int ewitab;
387 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
388 real *ewtab;
390 x = xx[0];
391 f = ff[0];
393 nri = nlist->nri;
394 iinr = nlist->iinr;
395 jindex = nlist->jindex;
396 jjnr = nlist->jjnr;
397 shiftidx = nlist->shift;
398 gid = nlist->gid;
399 shiftvec = fr->shift_vec[0];
400 fshift = fr->fshift[0];
401 facel = fr->ic->epsfac;
402 charge = mdatoms->chargeA;
404 sh_ewald = fr->ic->sh_ewald;
405 ewtab = fr->ic->tabq_coul_F;
406 ewtabscale = fr->ic->tabq_scale;
407 ewtabhalfspace = 0.5/ewtabscale;
409 /* Setup water-specific parameters */
410 inr = nlist->iinr[0];
411 iq0 = facel*charge[inr+0];
412 iq1 = facel*charge[inr+1];
413 iq2 = facel*charge[inr+2];
415 outeriter = 0;
416 inneriter = 0;
418 /* Start outer loop over neighborlists */
419 for(iidx=0; iidx<nri; iidx++)
421 /* Load shift vector for this list */
422 i_shift_offset = DIM*shiftidx[iidx];
423 shX = shiftvec[i_shift_offset+XX];
424 shY = shiftvec[i_shift_offset+YY];
425 shZ = shiftvec[i_shift_offset+ZZ];
427 /* Load limits for loop over neighbors */
428 j_index_start = jindex[iidx];
429 j_index_end = jindex[iidx+1];
431 /* Get outer coordinate index */
432 inr = iinr[iidx];
433 i_coord_offset = DIM*inr;
435 /* Load i particle coords and add shift vector */
436 ix0 = shX + x[i_coord_offset+DIM*0+XX];
437 iy0 = shY + x[i_coord_offset+DIM*0+YY];
438 iz0 = shZ + x[i_coord_offset+DIM*0+ZZ];
439 ix1 = shX + x[i_coord_offset+DIM*1+XX];
440 iy1 = shY + x[i_coord_offset+DIM*1+YY];
441 iz1 = shZ + x[i_coord_offset+DIM*1+ZZ];
442 ix2 = shX + x[i_coord_offset+DIM*2+XX];
443 iy2 = shY + x[i_coord_offset+DIM*2+YY];
444 iz2 = shZ + x[i_coord_offset+DIM*2+ZZ];
446 fix0 = 0.0;
447 fiy0 = 0.0;
448 fiz0 = 0.0;
449 fix1 = 0.0;
450 fiy1 = 0.0;
451 fiz1 = 0.0;
452 fix2 = 0.0;
453 fiy2 = 0.0;
454 fiz2 = 0.0;
456 /* Start inner kernel loop */
457 for(jidx=j_index_start; jidx<j_index_end; jidx++)
459 /* Get j neighbor index, and coordinate index */
460 jnr = jjnr[jidx];
461 j_coord_offset = DIM*jnr;
463 /* load j atom coordinates */
464 jx0 = x[j_coord_offset+DIM*0+XX];
465 jy0 = x[j_coord_offset+DIM*0+YY];
466 jz0 = x[j_coord_offset+DIM*0+ZZ];
468 /* Calculate displacement vector */
469 dx00 = ix0 - jx0;
470 dy00 = iy0 - jy0;
471 dz00 = iz0 - jz0;
472 dx10 = ix1 - jx0;
473 dy10 = iy1 - jy0;
474 dz10 = iz1 - jz0;
475 dx20 = ix2 - jx0;
476 dy20 = iy2 - jy0;
477 dz20 = iz2 - jz0;
479 /* Calculate squared distance and things based on it */
480 rsq00 = dx00*dx00+dy00*dy00+dz00*dz00;
481 rsq10 = dx10*dx10+dy10*dy10+dz10*dz10;
482 rsq20 = dx20*dx20+dy20*dy20+dz20*dz20;
484 rinv00 = 1.0/sqrt(rsq00);
485 rinv10 = 1.0/sqrt(rsq10);
486 rinv20 = 1.0/sqrt(rsq20);
488 rinvsq00 = rinv00*rinv00;
489 rinvsq10 = rinv10*rinv10;
490 rinvsq20 = rinv20*rinv20;
492 /* Load parameters for j particles */
493 jq0 = charge[jnr+0];
495 /**************************
496 * CALCULATE INTERACTIONS *
497 **************************/
499 r00 = rsq00*rinv00;
501 qq00 = iq0*jq0;
503 /* EWALD ELECTROSTATICS */
505 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
506 ewrt = r00*ewtabscale;
507 ewitab = ewrt;
508 eweps = ewrt-ewitab;
509 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
510 felec = qq00*rinv00*(rinvsq00-felec);
512 fscal = felec;
514 /* Calculate temporary vectorial force */
515 tx = fscal*dx00;
516 ty = fscal*dy00;
517 tz = fscal*dz00;
519 /* Update vectorial force */
520 fix0 += tx;
521 fiy0 += ty;
522 fiz0 += tz;
523 f[j_coord_offset+DIM*0+XX] -= tx;
524 f[j_coord_offset+DIM*0+YY] -= ty;
525 f[j_coord_offset+DIM*0+ZZ] -= tz;
527 /**************************
528 * CALCULATE INTERACTIONS *
529 **************************/
531 r10 = rsq10*rinv10;
533 qq10 = iq1*jq0;
535 /* EWALD ELECTROSTATICS */
537 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
538 ewrt = r10*ewtabscale;
539 ewitab = ewrt;
540 eweps = ewrt-ewitab;
541 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
542 felec = qq10*rinv10*(rinvsq10-felec);
544 fscal = felec;
546 /* Calculate temporary vectorial force */
547 tx = fscal*dx10;
548 ty = fscal*dy10;
549 tz = fscal*dz10;
551 /* Update vectorial force */
552 fix1 += tx;
553 fiy1 += ty;
554 fiz1 += tz;
555 f[j_coord_offset+DIM*0+XX] -= tx;
556 f[j_coord_offset+DIM*0+YY] -= ty;
557 f[j_coord_offset+DIM*0+ZZ] -= tz;
559 /**************************
560 * CALCULATE INTERACTIONS *
561 **************************/
563 r20 = rsq20*rinv20;
565 qq20 = iq2*jq0;
567 /* EWALD ELECTROSTATICS */
569 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
570 ewrt = r20*ewtabscale;
571 ewitab = ewrt;
572 eweps = ewrt-ewitab;
573 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
574 felec = qq20*rinv20*(rinvsq20-felec);
576 fscal = felec;
578 /* Calculate temporary vectorial force */
579 tx = fscal*dx20;
580 ty = fscal*dy20;
581 tz = fscal*dz20;
583 /* Update vectorial force */
584 fix2 += tx;
585 fiy2 += ty;
586 fiz2 += tz;
587 f[j_coord_offset+DIM*0+XX] -= tx;
588 f[j_coord_offset+DIM*0+YY] -= ty;
589 f[j_coord_offset+DIM*0+ZZ] -= tz;
591 /* Inner loop uses 102 flops */
593 /* End of innermost loop */
595 tx = ty = tz = 0;
596 f[i_coord_offset+DIM*0+XX] += fix0;
597 f[i_coord_offset+DIM*0+YY] += fiy0;
598 f[i_coord_offset+DIM*0+ZZ] += fiz0;
599 tx += fix0;
600 ty += fiy0;
601 tz += fiz0;
602 f[i_coord_offset+DIM*1+XX] += fix1;
603 f[i_coord_offset+DIM*1+YY] += fiy1;
604 f[i_coord_offset+DIM*1+ZZ] += fiz1;
605 tx += fix1;
606 ty += fiy1;
607 tz += fiz1;
608 f[i_coord_offset+DIM*2+XX] += fix2;
609 f[i_coord_offset+DIM*2+YY] += fiy2;
610 f[i_coord_offset+DIM*2+ZZ] += fiz2;
611 tx += fix2;
612 ty += fiy2;
613 tz += fiz2;
614 fshift[i_shift_offset+XX] += tx;
615 fshift[i_shift_offset+YY] += ty;
616 fshift[i_shift_offset+ZZ] += tz;
618 /* Increment number of inner iterations */
619 inneriter += j_index_end - j_index_start;
621 /* Outer loop uses 30 flops */
624 /* Increment number of outer iterations */
625 outeriter += nri;
627 /* Update outer/inner flops */
629 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_F,outeriter*30 + inneriter*102);