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 sse2_double kernel generator.
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
47 #include "gromacs/simd/math_x86_sse2_double.h"
48 #include "kernelutil_x86_sse2_double.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomP1P1_VF_sse2_double
52 * Electrostatics interaction: ReactionField
53 * VdW interaction: CubicSplineTable
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecRF_VdwCSTab_GeomP1P1_VF_sse2_double
59 (t_nblist * gmx_restrict nlist,
60 rvec * gmx_restrict xx,
61 rvec * gmx_restrict ff,
62 t_forcerec * gmx_restrict fr,
63 t_mdatoms * gmx_restrict mdatoms,
64 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65 t_nrnb * gmx_restrict nrnb)
67 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68 * just 0 for non-waters.
69 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
70 * jnr indices corresponding to data put in the four positions in the SIMD register.
72 int i_shift_offset,i_coord_offset,outeriter,inneriter;
73 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
75 int j_coord_offsetA,j_coord_offsetB;
76 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
78 real *shiftvec,*fshift,*x,*f;
79 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
81 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
82 int vdwjidx0A,vdwjidx0B;
83 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
84 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
85 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
88 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
91 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
92 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
94 __m128i ifour = _mm_set1_epi32(4);
95 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
97 __m128d dummy_mask,cutoff_mask;
98 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
99 __m128d one = _mm_set1_pd(1.0);
100 __m128d two = _mm_set1_pd(2.0);
106 jindex = nlist->jindex;
108 shiftidx = nlist->shift;
110 shiftvec = fr->shift_vec[0];
111 fshift = fr->fshift[0];
112 facel = _mm_set1_pd(fr->epsfac);
113 charge = mdatoms->chargeA;
114 krf = _mm_set1_pd(fr->ic->k_rf);
115 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
116 crf = _mm_set1_pd(fr->ic->c_rf);
117 nvdwtype = fr->ntype;
119 vdwtype = mdatoms->typeA;
121 vftab = kernel_data->table_vdw->data;
122 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
124 /* Avoid stupid compiler warnings */
132 /* Start outer loop over neighborlists */
133 for(iidx=0; iidx<nri; iidx++)
135 /* Load shift vector for this list */
136 i_shift_offset = DIM*shiftidx[iidx];
138 /* Load limits for loop over neighbors */
139 j_index_start = jindex[iidx];
140 j_index_end = jindex[iidx+1];
142 /* Get outer coordinate index */
144 i_coord_offset = DIM*inr;
146 /* Load i particle coords and add shift vector */
147 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
149 fix0 = _mm_setzero_pd();
150 fiy0 = _mm_setzero_pd();
151 fiz0 = _mm_setzero_pd();
153 /* Load parameters for i particles */
154 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
155 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
157 /* Reset potential sums */
158 velecsum = _mm_setzero_pd();
159 vvdwsum = _mm_setzero_pd();
161 /* Start inner kernel loop */
162 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
165 /* Get j neighbor index, and coordinate index */
168 j_coord_offsetA = DIM*jnrA;
169 j_coord_offsetB = DIM*jnrB;
171 /* load j atom coordinates */
172 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
175 /* Calculate displacement vector */
176 dx00 = _mm_sub_pd(ix0,jx0);
177 dy00 = _mm_sub_pd(iy0,jy0);
178 dz00 = _mm_sub_pd(iz0,jz0);
180 /* Calculate squared distance and things based on it */
181 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
183 rinv00 = gmx_mm_invsqrt_pd(rsq00);
185 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
187 /* Load parameters for j particles */
188 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
189 vdwjidx0A = 2*vdwtype[jnrA+0];
190 vdwjidx0B = 2*vdwtype[jnrB+0];
192 /**************************
193 * CALCULATE INTERACTIONS *
194 **************************/
196 r00 = _mm_mul_pd(rsq00,rinv00);
198 /* Compute parameters for interactions between i and j atoms */
199 qq00 = _mm_mul_pd(iq0,jq0);
200 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
201 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
203 /* Calculate table index by multiplying r with table scale and truncate to integer */
204 rt = _mm_mul_pd(r00,vftabscale);
205 vfitab = _mm_cvttpd_epi32(rt);
206 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
207 vfitab = _mm_slli_epi32(vfitab,3);
209 /* REACTION-FIELD ELECTROSTATICS */
210 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
211 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
213 /* CUBIC SPLINE TABLE DISPERSION */
214 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
215 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
216 GMX_MM_TRANSPOSE2_PD(Y,F);
217 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
218 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
219 GMX_MM_TRANSPOSE2_PD(G,H);
220 Heps = _mm_mul_pd(vfeps,H);
221 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
222 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
223 vvdw6 = _mm_mul_pd(c6_00,VV);
224 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
225 fvdw6 = _mm_mul_pd(c6_00,FF);
227 /* CUBIC SPLINE TABLE REPULSION */
228 vfitab = _mm_add_epi32(vfitab,ifour);
229 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
230 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
231 GMX_MM_TRANSPOSE2_PD(Y,F);
232 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
233 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
234 GMX_MM_TRANSPOSE2_PD(G,H);
235 Heps = _mm_mul_pd(vfeps,H);
236 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
237 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
238 vvdw12 = _mm_mul_pd(c12_00,VV);
239 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
240 fvdw12 = _mm_mul_pd(c12_00,FF);
241 vvdw = _mm_add_pd(vvdw12,vvdw6);
242 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
244 /* Update potential sum for this i atom from the interaction with this j atom. */
245 velecsum = _mm_add_pd(velecsum,velec);
246 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
248 fscal = _mm_add_pd(felec,fvdw);
250 /* Calculate temporary vectorial force */
251 tx = _mm_mul_pd(fscal,dx00);
252 ty = _mm_mul_pd(fscal,dy00);
253 tz = _mm_mul_pd(fscal,dz00);
255 /* Update vectorial force */
256 fix0 = _mm_add_pd(fix0,tx);
257 fiy0 = _mm_add_pd(fiy0,ty);
258 fiz0 = _mm_add_pd(fiz0,tz);
260 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
262 /* Inner loop uses 67 flops */
269 j_coord_offsetA = DIM*jnrA;
271 /* load j atom coordinates */
272 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
275 /* Calculate displacement vector */
276 dx00 = _mm_sub_pd(ix0,jx0);
277 dy00 = _mm_sub_pd(iy0,jy0);
278 dz00 = _mm_sub_pd(iz0,jz0);
280 /* Calculate squared distance and things based on it */
281 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
283 rinv00 = gmx_mm_invsqrt_pd(rsq00);
285 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
287 /* Load parameters for j particles */
288 jq0 = _mm_load_sd(charge+jnrA+0);
289 vdwjidx0A = 2*vdwtype[jnrA+0];
291 /**************************
292 * CALCULATE INTERACTIONS *
293 **************************/
295 r00 = _mm_mul_pd(rsq00,rinv00);
297 /* Compute parameters for interactions between i and j atoms */
298 qq00 = _mm_mul_pd(iq0,jq0);
299 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
301 /* Calculate table index by multiplying r with table scale and truncate to integer */
302 rt = _mm_mul_pd(r00,vftabscale);
303 vfitab = _mm_cvttpd_epi32(rt);
304 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
305 vfitab = _mm_slli_epi32(vfitab,3);
307 /* REACTION-FIELD ELECTROSTATICS */
308 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
309 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
311 /* CUBIC SPLINE TABLE DISPERSION */
312 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
313 F = _mm_setzero_pd();
314 GMX_MM_TRANSPOSE2_PD(Y,F);
315 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
316 H = _mm_setzero_pd();
317 GMX_MM_TRANSPOSE2_PD(G,H);
318 Heps = _mm_mul_pd(vfeps,H);
319 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
320 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
321 vvdw6 = _mm_mul_pd(c6_00,VV);
322 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
323 fvdw6 = _mm_mul_pd(c6_00,FF);
325 /* CUBIC SPLINE TABLE REPULSION */
326 vfitab = _mm_add_epi32(vfitab,ifour);
327 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
328 F = _mm_setzero_pd();
329 GMX_MM_TRANSPOSE2_PD(Y,F);
330 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
331 H = _mm_setzero_pd();
332 GMX_MM_TRANSPOSE2_PD(G,H);
333 Heps = _mm_mul_pd(vfeps,H);
334 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
335 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
336 vvdw12 = _mm_mul_pd(c12_00,VV);
337 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
338 fvdw12 = _mm_mul_pd(c12_00,FF);
339 vvdw = _mm_add_pd(vvdw12,vvdw6);
340 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
342 /* Update potential sum for this i atom from the interaction with this j atom. */
343 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
344 velecsum = _mm_add_pd(velecsum,velec);
345 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
346 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
348 fscal = _mm_add_pd(felec,fvdw);
350 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
352 /* Calculate temporary vectorial force */
353 tx = _mm_mul_pd(fscal,dx00);
354 ty = _mm_mul_pd(fscal,dy00);
355 tz = _mm_mul_pd(fscal,dz00);
357 /* Update vectorial force */
358 fix0 = _mm_add_pd(fix0,tx);
359 fiy0 = _mm_add_pd(fiy0,ty);
360 fiz0 = _mm_add_pd(fiz0,tz);
362 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
364 /* Inner loop uses 67 flops */
367 /* End of innermost loop */
369 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
370 f+i_coord_offset,fshift+i_shift_offset);
373 /* Update potential energies */
374 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
375 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
377 /* Increment number of inner iterations */
378 inneriter += j_index_end - j_index_start;
380 /* Outer loop uses 9 flops */
383 /* Increment number of outer iterations */
386 /* Update outer/inner flops */
388 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*67);
391 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomP1P1_F_sse2_double
392 * Electrostatics interaction: ReactionField
393 * VdW interaction: CubicSplineTable
394 * Geometry: Particle-Particle
395 * Calculate force/pot: Force
398 nb_kernel_ElecRF_VdwCSTab_GeomP1P1_F_sse2_double
399 (t_nblist * gmx_restrict nlist,
400 rvec * gmx_restrict xx,
401 rvec * gmx_restrict ff,
402 t_forcerec * gmx_restrict fr,
403 t_mdatoms * gmx_restrict mdatoms,
404 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
405 t_nrnb * gmx_restrict nrnb)
407 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
408 * just 0 for non-waters.
409 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
410 * jnr indices corresponding to data put in the four positions in the SIMD register.
412 int i_shift_offset,i_coord_offset,outeriter,inneriter;
413 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
415 int j_coord_offsetA,j_coord_offsetB;
416 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
418 real *shiftvec,*fshift,*x,*f;
419 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
421 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
422 int vdwjidx0A,vdwjidx0B;
423 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
424 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
425 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
428 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
431 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
432 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
434 __m128i ifour = _mm_set1_epi32(4);
435 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
437 __m128d dummy_mask,cutoff_mask;
438 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
439 __m128d one = _mm_set1_pd(1.0);
440 __m128d two = _mm_set1_pd(2.0);
446 jindex = nlist->jindex;
448 shiftidx = nlist->shift;
450 shiftvec = fr->shift_vec[0];
451 fshift = fr->fshift[0];
452 facel = _mm_set1_pd(fr->epsfac);
453 charge = mdatoms->chargeA;
454 krf = _mm_set1_pd(fr->ic->k_rf);
455 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
456 crf = _mm_set1_pd(fr->ic->c_rf);
457 nvdwtype = fr->ntype;
459 vdwtype = mdatoms->typeA;
461 vftab = kernel_data->table_vdw->data;
462 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
464 /* Avoid stupid compiler warnings */
472 /* Start outer loop over neighborlists */
473 for(iidx=0; iidx<nri; iidx++)
475 /* Load shift vector for this list */
476 i_shift_offset = DIM*shiftidx[iidx];
478 /* Load limits for loop over neighbors */
479 j_index_start = jindex[iidx];
480 j_index_end = jindex[iidx+1];
482 /* Get outer coordinate index */
484 i_coord_offset = DIM*inr;
486 /* Load i particle coords and add shift vector */
487 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
489 fix0 = _mm_setzero_pd();
490 fiy0 = _mm_setzero_pd();
491 fiz0 = _mm_setzero_pd();
493 /* Load parameters for i particles */
494 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
495 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
497 /* Start inner kernel loop */
498 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
501 /* Get j neighbor index, and coordinate index */
504 j_coord_offsetA = DIM*jnrA;
505 j_coord_offsetB = DIM*jnrB;
507 /* load j atom coordinates */
508 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
511 /* Calculate displacement vector */
512 dx00 = _mm_sub_pd(ix0,jx0);
513 dy00 = _mm_sub_pd(iy0,jy0);
514 dz00 = _mm_sub_pd(iz0,jz0);
516 /* Calculate squared distance and things based on it */
517 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
519 rinv00 = gmx_mm_invsqrt_pd(rsq00);
521 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
523 /* Load parameters for j particles */
524 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
525 vdwjidx0A = 2*vdwtype[jnrA+0];
526 vdwjidx0B = 2*vdwtype[jnrB+0];
528 /**************************
529 * CALCULATE INTERACTIONS *
530 **************************/
532 r00 = _mm_mul_pd(rsq00,rinv00);
534 /* Compute parameters for interactions between i and j atoms */
535 qq00 = _mm_mul_pd(iq0,jq0);
536 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
537 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
539 /* Calculate table index by multiplying r with table scale and truncate to integer */
540 rt = _mm_mul_pd(r00,vftabscale);
541 vfitab = _mm_cvttpd_epi32(rt);
542 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
543 vfitab = _mm_slli_epi32(vfitab,3);
545 /* REACTION-FIELD ELECTROSTATICS */
546 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
548 /* CUBIC SPLINE TABLE DISPERSION */
549 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
550 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
551 GMX_MM_TRANSPOSE2_PD(Y,F);
552 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
553 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
554 GMX_MM_TRANSPOSE2_PD(G,H);
555 Heps = _mm_mul_pd(vfeps,H);
556 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
557 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
558 fvdw6 = _mm_mul_pd(c6_00,FF);
560 /* CUBIC SPLINE TABLE REPULSION */
561 vfitab = _mm_add_epi32(vfitab,ifour);
562 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
563 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
564 GMX_MM_TRANSPOSE2_PD(Y,F);
565 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
566 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
567 GMX_MM_TRANSPOSE2_PD(G,H);
568 Heps = _mm_mul_pd(vfeps,H);
569 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
570 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
571 fvdw12 = _mm_mul_pd(c12_00,FF);
572 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
574 fscal = _mm_add_pd(felec,fvdw);
576 /* Calculate temporary vectorial force */
577 tx = _mm_mul_pd(fscal,dx00);
578 ty = _mm_mul_pd(fscal,dy00);
579 tz = _mm_mul_pd(fscal,dz00);
581 /* Update vectorial force */
582 fix0 = _mm_add_pd(fix0,tx);
583 fiy0 = _mm_add_pd(fiy0,ty);
584 fiz0 = _mm_add_pd(fiz0,tz);
586 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
588 /* Inner loop uses 54 flops */
595 j_coord_offsetA = DIM*jnrA;
597 /* load j atom coordinates */
598 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
601 /* Calculate displacement vector */
602 dx00 = _mm_sub_pd(ix0,jx0);
603 dy00 = _mm_sub_pd(iy0,jy0);
604 dz00 = _mm_sub_pd(iz0,jz0);
606 /* Calculate squared distance and things based on it */
607 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
609 rinv00 = gmx_mm_invsqrt_pd(rsq00);
611 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
613 /* Load parameters for j particles */
614 jq0 = _mm_load_sd(charge+jnrA+0);
615 vdwjidx0A = 2*vdwtype[jnrA+0];
617 /**************************
618 * CALCULATE INTERACTIONS *
619 **************************/
621 r00 = _mm_mul_pd(rsq00,rinv00);
623 /* Compute parameters for interactions between i and j atoms */
624 qq00 = _mm_mul_pd(iq0,jq0);
625 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
627 /* Calculate table index by multiplying r with table scale and truncate to integer */
628 rt = _mm_mul_pd(r00,vftabscale);
629 vfitab = _mm_cvttpd_epi32(rt);
630 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
631 vfitab = _mm_slli_epi32(vfitab,3);
633 /* REACTION-FIELD ELECTROSTATICS */
634 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
636 /* CUBIC SPLINE TABLE DISPERSION */
637 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
638 F = _mm_setzero_pd();
639 GMX_MM_TRANSPOSE2_PD(Y,F);
640 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
641 H = _mm_setzero_pd();
642 GMX_MM_TRANSPOSE2_PD(G,H);
643 Heps = _mm_mul_pd(vfeps,H);
644 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
645 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
646 fvdw6 = _mm_mul_pd(c6_00,FF);
648 /* CUBIC SPLINE TABLE REPULSION */
649 vfitab = _mm_add_epi32(vfitab,ifour);
650 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
651 F = _mm_setzero_pd();
652 GMX_MM_TRANSPOSE2_PD(Y,F);
653 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
654 H = _mm_setzero_pd();
655 GMX_MM_TRANSPOSE2_PD(G,H);
656 Heps = _mm_mul_pd(vfeps,H);
657 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
658 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
659 fvdw12 = _mm_mul_pd(c12_00,FF);
660 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
662 fscal = _mm_add_pd(felec,fvdw);
664 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
666 /* Calculate temporary vectorial force */
667 tx = _mm_mul_pd(fscal,dx00);
668 ty = _mm_mul_pd(fscal,dy00);
669 tz = _mm_mul_pd(fscal,dz00);
671 /* Update vectorial force */
672 fix0 = _mm_add_pd(fix0,tx);
673 fiy0 = _mm_add_pd(fiy0,ty);
674 fiz0 = _mm_add_pd(fiz0,tz);
676 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
678 /* Inner loop uses 54 flops */
681 /* End of innermost loop */
683 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
684 f+i_coord_offset,fshift+i_shift_offset);
686 /* Increment number of inner iterations */
687 inneriter += j_index_end - j_index_start;
689 /* Outer loop uses 7 flops */
692 /* Increment number of outer iterations */
695 /* Update outer/inner flops */
697 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*54);