Removed include simple.h from nb_kernel_avx_128_fma_XX
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecRF_VdwCSTab_GeomW4P1_avx_128_fma_single.c
blob93a2b77218756fa836bc2a0e0dbf2e7e3954066a
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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 avx_128_fma_single kernel generator.
38 #include "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/legacyheaders/nrnb.h"
48 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
49 #include "kernelutil_x86_avx_128_fma_single.h"
52 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomW4P1_VF_avx_128_fma_single
53 * Electrostatics interaction: ReactionField
54 * VdW interaction: CubicSplineTable
55 * Geometry: Water4-Particle
56 * Calculate force/pot: PotentialAndForce
58 void
59 nb_kernel_ElecRF_VdwCSTab_GeomW4P1_VF_avx_128_fma_single
60 (t_nblist * gmx_restrict nlist,
61 rvec * gmx_restrict xx,
62 rvec * gmx_restrict ff,
63 t_forcerec * gmx_restrict fr,
64 t_mdatoms * gmx_restrict mdatoms,
65 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
66 t_nrnb * gmx_restrict nrnb)
68 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
69 * just 0 for non-waters.
70 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
71 * jnr indices corresponding to data put in the four positions in the SIMD register.
73 int i_shift_offset,i_coord_offset,outeriter,inneriter;
74 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
75 int jnrA,jnrB,jnrC,jnrD;
76 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
77 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
79 real rcutoff_scalar;
80 real *shiftvec,*fshift,*x,*f;
81 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82 real scratch[4*DIM];
83 __m128 fscal,rcutoff,rcutoff2,jidxall;
84 int vdwioffset0;
85 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
86 int vdwioffset1;
87 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
88 int vdwioffset2;
89 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
90 int vdwioffset3;
91 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
92 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
93 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
94 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
95 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
96 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
97 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
98 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
99 real *charge;
100 int nvdwtype;
101 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
102 int *vdwtype;
103 real *vdwparam;
104 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
105 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
106 __m128i vfitab;
107 __m128i ifour = _mm_set1_epi32(4);
108 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
109 real *vftab;
110 __m128 dummy_mask,cutoff_mask;
111 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
112 __m128 one = _mm_set1_ps(1.0);
113 __m128 two = _mm_set1_ps(2.0);
114 x = xx[0];
115 f = ff[0];
117 nri = nlist->nri;
118 iinr = nlist->iinr;
119 jindex = nlist->jindex;
120 jjnr = nlist->jjnr;
121 shiftidx = nlist->shift;
122 gid = nlist->gid;
123 shiftvec = fr->shift_vec[0];
124 fshift = fr->fshift[0];
125 facel = _mm_set1_ps(fr->epsfac);
126 charge = mdatoms->chargeA;
127 krf = _mm_set1_ps(fr->ic->k_rf);
128 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
129 crf = _mm_set1_ps(fr->ic->c_rf);
130 nvdwtype = fr->ntype;
131 vdwparam = fr->nbfp;
132 vdwtype = mdatoms->typeA;
134 vftab = kernel_data->table_vdw->data;
135 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
137 /* Setup water-specific parameters */
138 inr = nlist->iinr[0];
139 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
140 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
141 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
142 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
144 /* Avoid stupid compiler warnings */
145 jnrA = jnrB = jnrC = jnrD = 0;
146 j_coord_offsetA = 0;
147 j_coord_offsetB = 0;
148 j_coord_offsetC = 0;
149 j_coord_offsetD = 0;
151 outeriter = 0;
152 inneriter = 0;
154 for(iidx=0;iidx<4*DIM;iidx++)
156 scratch[iidx] = 0.0;
159 /* Start outer loop over neighborlists */
160 for(iidx=0; iidx<nri; iidx++)
162 /* Load shift vector for this list */
163 i_shift_offset = DIM*shiftidx[iidx];
165 /* Load limits for loop over neighbors */
166 j_index_start = jindex[iidx];
167 j_index_end = jindex[iidx+1];
169 /* Get outer coordinate index */
170 inr = iinr[iidx];
171 i_coord_offset = DIM*inr;
173 /* Load i particle coords and add shift vector */
174 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
175 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
177 fix0 = _mm_setzero_ps();
178 fiy0 = _mm_setzero_ps();
179 fiz0 = _mm_setzero_ps();
180 fix1 = _mm_setzero_ps();
181 fiy1 = _mm_setzero_ps();
182 fiz1 = _mm_setzero_ps();
183 fix2 = _mm_setzero_ps();
184 fiy2 = _mm_setzero_ps();
185 fiz2 = _mm_setzero_ps();
186 fix3 = _mm_setzero_ps();
187 fiy3 = _mm_setzero_ps();
188 fiz3 = _mm_setzero_ps();
190 /* Reset potential sums */
191 velecsum = _mm_setzero_ps();
192 vvdwsum = _mm_setzero_ps();
194 /* Start inner kernel loop */
195 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
198 /* Get j neighbor index, and coordinate index */
199 jnrA = jjnr[jidx];
200 jnrB = jjnr[jidx+1];
201 jnrC = jjnr[jidx+2];
202 jnrD = jjnr[jidx+3];
203 j_coord_offsetA = DIM*jnrA;
204 j_coord_offsetB = DIM*jnrB;
205 j_coord_offsetC = DIM*jnrC;
206 j_coord_offsetD = DIM*jnrD;
208 /* load j atom coordinates */
209 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
210 x+j_coord_offsetC,x+j_coord_offsetD,
211 &jx0,&jy0,&jz0);
213 /* Calculate displacement vector */
214 dx00 = _mm_sub_ps(ix0,jx0);
215 dy00 = _mm_sub_ps(iy0,jy0);
216 dz00 = _mm_sub_ps(iz0,jz0);
217 dx10 = _mm_sub_ps(ix1,jx0);
218 dy10 = _mm_sub_ps(iy1,jy0);
219 dz10 = _mm_sub_ps(iz1,jz0);
220 dx20 = _mm_sub_ps(ix2,jx0);
221 dy20 = _mm_sub_ps(iy2,jy0);
222 dz20 = _mm_sub_ps(iz2,jz0);
223 dx30 = _mm_sub_ps(ix3,jx0);
224 dy30 = _mm_sub_ps(iy3,jy0);
225 dz30 = _mm_sub_ps(iz3,jz0);
227 /* Calculate squared distance and things based on it */
228 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
229 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
230 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
231 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
233 rinv00 = gmx_mm_invsqrt_ps(rsq00);
234 rinv10 = gmx_mm_invsqrt_ps(rsq10);
235 rinv20 = gmx_mm_invsqrt_ps(rsq20);
236 rinv30 = gmx_mm_invsqrt_ps(rsq30);
238 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
239 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
240 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
242 /* Load parameters for j particles */
243 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
244 charge+jnrC+0,charge+jnrD+0);
245 vdwjidx0A = 2*vdwtype[jnrA+0];
246 vdwjidx0B = 2*vdwtype[jnrB+0];
247 vdwjidx0C = 2*vdwtype[jnrC+0];
248 vdwjidx0D = 2*vdwtype[jnrD+0];
250 fjx0 = _mm_setzero_ps();
251 fjy0 = _mm_setzero_ps();
252 fjz0 = _mm_setzero_ps();
254 /**************************
255 * CALCULATE INTERACTIONS *
256 **************************/
258 r00 = _mm_mul_ps(rsq00,rinv00);
260 /* Compute parameters for interactions between i and j atoms */
261 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
262 vdwparam+vdwioffset0+vdwjidx0B,
263 vdwparam+vdwioffset0+vdwjidx0C,
264 vdwparam+vdwioffset0+vdwjidx0D,
265 &c6_00,&c12_00);
267 /* Calculate table index by multiplying r with table scale and truncate to integer */
268 rt = _mm_mul_ps(r00,vftabscale);
269 vfitab = _mm_cvttps_epi32(rt);
270 #ifdef __XOP__
271 vfeps = _mm_frcz_ps(rt);
272 #else
273 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
274 #endif
275 twovfeps = _mm_add_ps(vfeps,vfeps);
276 vfitab = _mm_slli_epi32(vfitab,3);
278 /* CUBIC SPLINE TABLE DISPERSION */
279 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
280 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
281 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
282 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
283 _MM_TRANSPOSE4_PS(Y,F,G,H);
284 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
285 VV = _mm_macc_ps(vfeps,Fp,Y);
286 vvdw6 = _mm_mul_ps(c6_00,VV);
287 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
288 fvdw6 = _mm_mul_ps(c6_00,FF);
290 /* CUBIC SPLINE TABLE REPULSION */
291 vfitab = _mm_add_epi32(vfitab,ifour);
292 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
293 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
294 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
295 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
296 _MM_TRANSPOSE4_PS(Y,F,G,H);
297 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
298 VV = _mm_macc_ps(vfeps,Fp,Y);
299 vvdw12 = _mm_mul_ps(c12_00,VV);
300 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
301 fvdw12 = _mm_mul_ps(c12_00,FF);
302 vvdw = _mm_add_ps(vvdw12,vvdw6);
303 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
305 /* Update potential sum for this i atom from the interaction with this j atom. */
306 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
308 fscal = fvdw;
310 /* Update vectorial force */
311 fix0 = _mm_macc_ps(dx00,fscal,fix0);
312 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
313 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
315 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
316 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
317 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
319 /**************************
320 * CALCULATE INTERACTIONS *
321 **************************/
323 /* Compute parameters for interactions between i and j atoms */
324 qq10 = _mm_mul_ps(iq1,jq0);
326 /* REACTION-FIELD ELECTROSTATICS */
327 velec = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
328 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
330 /* Update potential sum for this i atom from the interaction with this j atom. */
331 velecsum = _mm_add_ps(velecsum,velec);
333 fscal = felec;
335 /* Update vectorial force */
336 fix1 = _mm_macc_ps(dx10,fscal,fix1);
337 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
338 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
340 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
341 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
342 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
344 /**************************
345 * CALCULATE INTERACTIONS *
346 **************************/
348 /* Compute parameters for interactions between i and j atoms */
349 qq20 = _mm_mul_ps(iq2,jq0);
351 /* REACTION-FIELD ELECTROSTATICS */
352 velec = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
353 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
355 /* Update potential sum for this i atom from the interaction with this j atom. */
356 velecsum = _mm_add_ps(velecsum,velec);
358 fscal = felec;
360 /* Update vectorial force */
361 fix2 = _mm_macc_ps(dx20,fscal,fix2);
362 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
363 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
365 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
366 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
367 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
369 /**************************
370 * CALCULATE INTERACTIONS *
371 **************************/
373 /* Compute parameters for interactions between i and j atoms */
374 qq30 = _mm_mul_ps(iq3,jq0);
376 /* REACTION-FIELD ELECTROSTATICS */
377 velec = _mm_mul_ps(qq30,_mm_sub_ps(_mm_macc_ps(krf,rsq30,rinv30),crf));
378 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
380 /* Update potential sum for this i atom from the interaction with this j atom. */
381 velecsum = _mm_add_ps(velecsum,velec);
383 fscal = felec;
385 /* Update vectorial force */
386 fix3 = _mm_macc_ps(dx30,fscal,fix3);
387 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
388 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
390 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
391 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
392 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
394 fjptrA = f+j_coord_offsetA;
395 fjptrB = f+j_coord_offsetB;
396 fjptrC = f+j_coord_offsetC;
397 fjptrD = f+j_coord_offsetD;
399 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
401 /* Inner loop uses 164 flops */
404 if(jidx<j_index_end)
407 /* Get j neighbor index, and coordinate index */
408 jnrlistA = jjnr[jidx];
409 jnrlistB = jjnr[jidx+1];
410 jnrlistC = jjnr[jidx+2];
411 jnrlistD = jjnr[jidx+3];
412 /* Sign of each element will be negative for non-real atoms.
413 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
414 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
416 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
417 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
418 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
419 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
420 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
421 j_coord_offsetA = DIM*jnrA;
422 j_coord_offsetB = DIM*jnrB;
423 j_coord_offsetC = DIM*jnrC;
424 j_coord_offsetD = DIM*jnrD;
426 /* load j atom coordinates */
427 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
428 x+j_coord_offsetC,x+j_coord_offsetD,
429 &jx0,&jy0,&jz0);
431 /* Calculate displacement vector */
432 dx00 = _mm_sub_ps(ix0,jx0);
433 dy00 = _mm_sub_ps(iy0,jy0);
434 dz00 = _mm_sub_ps(iz0,jz0);
435 dx10 = _mm_sub_ps(ix1,jx0);
436 dy10 = _mm_sub_ps(iy1,jy0);
437 dz10 = _mm_sub_ps(iz1,jz0);
438 dx20 = _mm_sub_ps(ix2,jx0);
439 dy20 = _mm_sub_ps(iy2,jy0);
440 dz20 = _mm_sub_ps(iz2,jz0);
441 dx30 = _mm_sub_ps(ix3,jx0);
442 dy30 = _mm_sub_ps(iy3,jy0);
443 dz30 = _mm_sub_ps(iz3,jz0);
445 /* Calculate squared distance and things based on it */
446 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
447 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
448 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
449 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
451 rinv00 = gmx_mm_invsqrt_ps(rsq00);
452 rinv10 = gmx_mm_invsqrt_ps(rsq10);
453 rinv20 = gmx_mm_invsqrt_ps(rsq20);
454 rinv30 = gmx_mm_invsqrt_ps(rsq30);
456 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
457 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
458 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
460 /* Load parameters for j particles */
461 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
462 charge+jnrC+0,charge+jnrD+0);
463 vdwjidx0A = 2*vdwtype[jnrA+0];
464 vdwjidx0B = 2*vdwtype[jnrB+0];
465 vdwjidx0C = 2*vdwtype[jnrC+0];
466 vdwjidx0D = 2*vdwtype[jnrD+0];
468 fjx0 = _mm_setzero_ps();
469 fjy0 = _mm_setzero_ps();
470 fjz0 = _mm_setzero_ps();
472 /**************************
473 * CALCULATE INTERACTIONS *
474 **************************/
476 r00 = _mm_mul_ps(rsq00,rinv00);
477 r00 = _mm_andnot_ps(dummy_mask,r00);
479 /* Compute parameters for interactions between i and j atoms */
480 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
481 vdwparam+vdwioffset0+vdwjidx0B,
482 vdwparam+vdwioffset0+vdwjidx0C,
483 vdwparam+vdwioffset0+vdwjidx0D,
484 &c6_00,&c12_00);
486 /* Calculate table index by multiplying r with table scale and truncate to integer */
487 rt = _mm_mul_ps(r00,vftabscale);
488 vfitab = _mm_cvttps_epi32(rt);
489 #ifdef __XOP__
490 vfeps = _mm_frcz_ps(rt);
491 #else
492 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
493 #endif
494 twovfeps = _mm_add_ps(vfeps,vfeps);
495 vfitab = _mm_slli_epi32(vfitab,3);
497 /* CUBIC SPLINE TABLE DISPERSION */
498 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
499 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
500 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
501 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
502 _MM_TRANSPOSE4_PS(Y,F,G,H);
503 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
504 VV = _mm_macc_ps(vfeps,Fp,Y);
505 vvdw6 = _mm_mul_ps(c6_00,VV);
506 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
507 fvdw6 = _mm_mul_ps(c6_00,FF);
509 /* CUBIC SPLINE TABLE REPULSION */
510 vfitab = _mm_add_epi32(vfitab,ifour);
511 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
512 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
513 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
514 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
515 _MM_TRANSPOSE4_PS(Y,F,G,H);
516 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
517 VV = _mm_macc_ps(vfeps,Fp,Y);
518 vvdw12 = _mm_mul_ps(c12_00,VV);
519 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
520 fvdw12 = _mm_mul_ps(c12_00,FF);
521 vvdw = _mm_add_ps(vvdw12,vvdw6);
522 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
524 /* Update potential sum for this i atom from the interaction with this j atom. */
525 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
526 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
528 fscal = fvdw;
530 fscal = _mm_andnot_ps(dummy_mask,fscal);
532 /* Update vectorial force */
533 fix0 = _mm_macc_ps(dx00,fscal,fix0);
534 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
535 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
537 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
538 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
539 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
541 /**************************
542 * CALCULATE INTERACTIONS *
543 **************************/
545 /* Compute parameters for interactions between i and j atoms */
546 qq10 = _mm_mul_ps(iq1,jq0);
548 /* REACTION-FIELD ELECTROSTATICS */
549 velec = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
550 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
552 /* Update potential sum for this i atom from the interaction with this j atom. */
553 velec = _mm_andnot_ps(dummy_mask,velec);
554 velecsum = _mm_add_ps(velecsum,velec);
556 fscal = felec;
558 fscal = _mm_andnot_ps(dummy_mask,fscal);
560 /* Update vectorial force */
561 fix1 = _mm_macc_ps(dx10,fscal,fix1);
562 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
563 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
565 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
566 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
567 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
569 /**************************
570 * CALCULATE INTERACTIONS *
571 **************************/
573 /* Compute parameters for interactions between i and j atoms */
574 qq20 = _mm_mul_ps(iq2,jq0);
576 /* REACTION-FIELD ELECTROSTATICS */
577 velec = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
578 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
580 /* Update potential sum for this i atom from the interaction with this j atom. */
581 velec = _mm_andnot_ps(dummy_mask,velec);
582 velecsum = _mm_add_ps(velecsum,velec);
584 fscal = felec;
586 fscal = _mm_andnot_ps(dummy_mask,fscal);
588 /* Update vectorial force */
589 fix2 = _mm_macc_ps(dx20,fscal,fix2);
590 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
591 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
593 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
594 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
595 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
597 /**************************
598 * CALCULATE INTERACTIONS *
599 **************************/
601 /* Compute parameters for interactions between i and j atoms */
602 qq30 = _mm_mul_ps(iq3,jq0);
604 /* REACTION-FIELD ELECTROSTATICS */
605 velec = _mm_mul_ps(qq30,_mm_sub_ps(_mm_macc_ps(krf,rsq30,rinv30),crf));
606 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
608 /* Update potential sum for this i atom from the interaction with this j atom. */
609 velec = _mm_andnot_ps(dummy_mask,velec);
610 velecsum = _mm_add_ps(velecsum,velec);
612 fscal = felec;
614 fscal = _mm_andnot_ps(dummy_mask,fscal);
616 /* Update vectorial force */
617 fix3 = _mm_macc_ps(dx30,fscal,fix3);
618 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
619 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
621 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
622 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
623 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
625 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
626 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
627 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
628 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
630 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
632 /* Inner loop uses 165 flops */
635 /* End of innermost loop */
637 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
638 f+i_coord_offset,fshift+i_shift_offset);
640 ggid = gid[iidx];
641 /* Update potential energies */
642 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
643 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
645 /* Increment number of inner iterations */
646 inneriter += j_index_end - j_index_start;
648 /* Outer loop uses 26 flops */
651 /* Increment number of outer iterations */
652 outeriter += nri;
654 /* Update outer/inner flops */
656 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*165);
659 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomW4P1_F_avx_128_fma_single
660 * Electrostatics interaction: ReactionField
661 * VdW interaction: CubicSplineTable
662 * Geometry: Water4-Particle
663 * Calculate force/pot: Force
665 void
666 nb_kernel_ElecRF_VdwCSTab_GeomW4P1_F_avx_128_fma_single
667 (t_nblist * gmx_restrict nlist,
668 rvec * gmx_restrict xx,
669 rvec * gmx_restrict ff,
670 t_forcerec * gmx_restrict fr,
671 t_mdatoms * gmx_restrict mdatoms,
672 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
673 t_nrnb * gmx_restrict nrnb)
675 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
676 * just 0 for non-waters.
677 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
678 * jnr indices corresponding to data put in the four positions in the SIMD register.
680 int i_shift_offset,i_coord_offset,outeriter,inneriter;
681 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
682 int jnrA,jnrB,jnrC,jnrD;
683 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
684 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
685 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
686 real rcutoff_scalar;
687 real *shiftvec,*fshift,*x,*f;
688 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
689 real scratch[4*DIM];
690 __m128 fscal,rcutoff,rcutoff2,jidxall;
691 int vdwioffset0;
692 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
693 int vdwioffset1;
694 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
695 int vdwioffset2;
696 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
697 int vdwioffset3;
698 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
699 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
700 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
701 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
702 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
703 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
704 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
705 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
706 real *charge;
707 int nvdwtype;
708 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
709 int *vdwtype;
710 real *vdwparam;
711 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
712 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
713 __m128i vfitab;
714 __m128i ifour = _mm_set1_epi32(4);
715 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
716 real *vftab;
717 __m128 dummy_mask,cutoff_mask;
718 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
719 __m128 one = _mm_set1_ps(1.0);
720 __m128 two = _mm_set1_ps(2.0);
721 x = xx[0];
722 f = ff[0];
724 nri = nlist->nri;
725 iinr = nlist->iinr;
726 jindex = nlist->jindex;
727 jjnr = nlist->jjnr;
728 shiftidx = nlist->shift;
729 gid = nlist->gid;
730 shiftvec = fr->shift_vec[0];
731 fshift = fr->fshift[0];
732 facel = _mm_set1_ps(fr->epsfac);
733 charge = mdatoms->chargeA;
734 krf = _mm_set1_ps(fr->ic->k_rf);
735 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
736 crf = _mm_set1_ps(fr->ic->c_rf);
737 nvdwtype = fr->ntype;
738 vdwparam = fr->nbfp;
739 vdwtype = mdatoms->typeA;
741 vftab = kernel_data->table_vdw->data;
742 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
744 /* Setup water-specific parameters */
745 inr = nlist->iinr[0];
746 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
747 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
748 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
749 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
751 /* Avoid stupid compiler warnings */
752 jnrA = jnrB = jnrC = jnrD = 0;
753 j_coord_offsetA = 0;
754 j_coord_offsetB = 0;
755 j_coord_offsetC = 0;
756 j_coord_offsetD = 0;
758 outeriter = 0;
759 inneriter = 0;
761 for(iidx=0;iidx<4*DIM;iidx++)
763 scratch[iidx] = 0.0;
766 /* Start outer loop over neighborlists */
767 for(iidx=0; iidx<nri; iidx++)
769 /* Load shift vector for this list */
770 i_shift_offset = DIM*shiftidx[iidx];
772 /* Load limits for loop over neighbors */
773 j_index_start = jindex[iidx];
774 j_index_end = jindex[iidx+1];
776 /* Get outer coordinate index */
777 inr = iinr[iidx];
778 i_coord_offset = DIM*inr;
780 /* Load i particle coords and add shift vector */
781 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
782 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
784 fix0 = _mm_setzero_ps();
785 fiy0 = _mm_setzero_ps();
786 fiz0 = _mm_setzero_ps();
787 fix1 = _mm_setzero_ps();
788 fiy1 = _mm_setzero_ps();
789 fiz1 = _mm_setzero_ps();
790 fix2 = _mm_setzero_ps();
791 fiy2 = _mm_setzero_ps();
792 fiz2 = _mm_setzero_ps();
793 fix3 = _mm_setzero_ps();
794 fiy3 = _mm_setzero_ps();
795 fiz3 = _mm_setzero_ps();
797 /* Start inner kernel loop */
798 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
801 /* Get j neighbor index, and coordinate index */
802 jnrA = jjnr[jidx];
803 jnrB = jjnr[jidx+1];
804 jnrC = jjnr[jidx+2];
805 jnrD = jjnr[jidx+3];
806 j_coord_offsetA = DIM*jnrA;
807 j_coord_offsetB = DIM*jnrB;
808 j_coord_offsetC = DIM*jnrC;
809 j_coord_offsetD = DIM*jnrD;
811 /* load j atom coordinates */
812 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
813 x+j_coord_offsetC,x+j_coord_offsetD,
814 &jx0,&jy0,&jz0);
816 /* Calculate displacement vector */
817 dx00 = _mm_sub_ps(ix0,jx0);
818 dy00 = _mm_sub_ps(iy0,jy0);
819 dz00 = _mm_sub_ps(iz0,jz0);
820 dx10 = _mm_sub_ps(ix1,jx0);
821 dy10 = _mm_sub_ps(iy1,jy0);
822 dz10 = _mm_sub_ps(iz1,jz0);
823 dx20 = _mm_sub_ps(ix2,jx0);
824 dy20 = _mm_sub_ps(iy2,jy0);
825 dz20 = _mm_sub_ps(iz2,jz0);
826 dx30 = _mm_sub_ps(ix3,jx0);
827 dy30 = _mm_sub_ps(iy3,jy0);
828 dz30 = _mm_sub_ps(iz3,jz0);
830 /* Calculate squared distance and things based on it */
831 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
832 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
833 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
834 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
836 rinv00 = gmx_mm_invsqrt_ps(rsq00);
837 rinv10 = gmx_mm_invsqrt_ps(rsq10);
838 rinv20 = gmx_mm_invsqrt_ps(rsq20);
839 rinv30 = gmx_mm_invsqrt_ps(rsq30);
841 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
842 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
843 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
845 /* Load parameters for j particles */
846 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
847 charge+jnrC+0,charge+jnrD+0);
848 vdwjidx0A = 2*vdwtype[jnrA+0];
849 vdwjidx0B = 2*vdwtype[jnrB+0];
850 vdwjidx0C = 2*vdwtype[jnrC+0];
851 vdwjidx0D = 2*vdwtype[jnrD+0];
853 fjx0 = _mm_setzero_ps();
854 fjy0 = _mm_setzero_ps();
855 fjz0 = _mm_setzero_ps();
857 /**************************
858 * CALCULATE INTERACTIONS *
859 **************************/
861 r00 = _mm_mul_ps(rsq00,rinv00);
863 /* Compute parameters for interactions between i and j atoms */
864 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
865 vdwparam+vdwioffset0+vdwjidx0B,
866 vdwparam+vdwioffset0+vdwjidx0C,
867 vdwparam+vdwioffset0+vdwjidx0D,
868 &c6_00,&c12_00);
870 /* Calculate table index by multiplying r with table scale and truncate to integer */
871 rt = _mm_mul_ps(r00,vftabscale);
872 vfitab = _mm_cvttps_epi32(rt);
873 #ifdef __XOP__
874 vfeps = _mm_frcz_ps(rt);
875 #else
876 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
877 #endif
878 twovfeps = _mm_add_ps(vfeps,vfeps);
879 vfitab = _mm_slli_epi32(vfitab,3);
881 /* CUBIC SPLINE TABLE DISPERSION */
882 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
883 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
884 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
885 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
886 _MM_TRANSPOSE4_PS(Y,F,G,H);
887 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
888 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
889 fvdw6 = _mm_mul_ps(c6_00,FF);
891 /* CUBIC SPLINE TABLE REPULSION */
892 vfitab = _mm_add_epi32(vfitab,ifour);
893 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
894 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
895 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
896 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
897 _MM_TRANSPOSE4_PS(Y,F,G,H);
898 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
899 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
900 fvdw12 = _mm_mul_ps(c12_00,FF);
901 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
903 fscal = fvdw;
905 /* Update vectorial force */
906 fix0 = _mm_macc_ps(dx00,fscal,fix0);
907 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
908 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
910 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
911 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
912 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
914 /**************************
915 * CALCULATE INTERACTIONS *
916 **************************/
918 /* Compute parameters for interactions between i and j atoms */
919 qq10 = _mm_mul_ps(iq1,jq0);
921 /* REACTION-FIELD ELECTROSTATICS */
922 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
924 fscal = felec;
926 /* Update vectorial force */
927 fix1 = _mm_macc_ps(dx10,fscal,fix1);
928 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
929 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
931 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
932 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
933 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
935 /**************************
936 * CALCULATE INTERACTIONS *
937 **************************/
939 /* Compute parameters for interactions between i and j atoms */
940 qq20 = _mm_mul_ps(iq2,jq0);
942 /* REACTION-FIELD ELECTROSTATICS */
943 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
945 fscal = felec;
947 /* Update vectorial force */
948 fix2 = _mm_macc_ps(dx20,fscal,fix2);
949 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
950 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
952 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
953 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
954 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
956 /**************************
957 * CALCULATE INTERACTIONS *
958 **************************/
960 /* Compute parameters for interactions between i and j atoms */
961 qq30 = _mm_mul_ps(iq3,jq0);
963 /* REACTION-FIELD ELECTROSTATICS */
964 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
966 fscal = felec;
968 /* Update vectorial force */
969 fix3 = _mm_macc_ps(dx30,fscal,fix3);
970 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
971 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
973 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
974 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
975 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
977 fjptrA = f+j_coord_offsetA;
978 fjptrB = f+j_coord_offsetB;
979 fjptrC = f+j_coord_offsetC;
980 fjptrD = f+j_coord_offsetD;
982 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
984 /* Inner loop uses 141 flops */
987 if(jidx<j_index_end)
990 /* Get j neighbor index, and coordinate index */
991 jnrlistA = jjnr[jidx];
992 jnrlistB = jjnr[jidx+1];
993 jnrlistC = jjnr[jidx+2];
994 jnrlistD = jjnr[jidx+3];
995 /* Sign of each element will be negative for non-real atoms.
996 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
997 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
999 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1000 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1001 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1002 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1003 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1004 j_coord_offsetA = DIM*jnrA;
1005 j_coord_offsetB = DIM*jnrB;
1006 j_coord_offsetC = DIM*jnrC;
1007 j_coord_offsetD = DIM*jnrD;
1009 /* load j atom coordinates */
1010 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1011 x+j_coord_offsetC,x+j_coord_offsetD,
1012 &jx0,&jy0,&jz0);
1014 /* Calculate displacement vector */
1015 dx00 = _mm_sub_ps(ix0,jx0);
1016 dy00 = _mm_sub_ps(iy0,jy0);
1017 dz00 = _mm_sub_ps(iz0,jz0);
1018 dx10 = _mm_sub_ps(ix1,jx0);
1019 dy10 = _mm_sub_ps(iy1,jy0);
1020 dz10 = _mm_sub_ps(iz1,jz0);
1021 dx20 = _mm_sub_ps(ix2,jx0);
1022 dy20 = _mm_sub_ps(iy2,jy0);
1023 dz20 = _mm_sub_ps(iz2,jz0);
1024 dx30 = _mm_sub_ps(ix3,jx0);
1025 dy30 = _mm_sub_ps(iy3,jy0);
1026 dz30 = _mm_sub_ps(iz3,jz0);
1028 /* Calculate squared distance and things based on it */
1029 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1030 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1031 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1032 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
1034 rinv00 = gmx_mm_invsqrt_ps(rsq00);
1035 rinv10 = gmx_mm_invsqrt_ps(rsq10);
1036 rinv20 = gmx_mm_invsqrt_ps(rsq20);
1037 rinv30 = gmx_mm_invsqrt_ps(rsq30);
1039 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
1040 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
1041 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
1043 /* Load parameters for j particles */
1044 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1045 charge+jnrC+0,charge+jnrD+0);
1046 vdwjidx0A = 2*vdwtype[jnrA+0];
1047 vdwjidx0B = 2*vdwtype[jnrB+0];
1048 vdwjidx0C = 2*vdwtype[jnrC+0];
1049 vdwjidx0D = 2*vdwtype[jnrD+0];
1051 fjx0 = _mm_setzero_ps();
1052 fjy0 = _mm_setzero_ps();
1053 fjz0 = _mm_setzero_ps();
1055 /**************************
1056 * CALCULATE INTERACTIONS *
1057 **************************/
1059 r00 = _mm_mul_ps(rsq00,rinv00);
1060 r00 = _mm_andnot_ps(dummy_mask,r00);
1062 /* Compute parameters for interactions between i and j atoms */
1063 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1064 vdwparam+vdwioffset0+vdwjidx0B,
1065 vdwparam+vdwioffset0+vdwjidx0C,
1066 vdwparam+vdwioffset0+vdwjidx0D,
1067 &c6_00,&c12_00);
1069 /* Calculate table index by multiplying r with table scale and truncate to integer */
1070 rt = _mm_mul_ps(r00,vftabscale);
1071 vfitab = _mm_cvttps_epi32(rt);
1072 #ifdef __XOP__
1073 vfeps = _mm_frcz_ps(rt);
1074 #else
1075 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
1076 #endif
1077 twovfeps = _mm_add_ps(vfeps,vfeps);
1078 vfitab = _mm_slli_epi32(vfitab,3);
1080 /* CUBIC SPLINE TABLE DISPERSION */
1081 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
1082 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
1083 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
1084 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
1085 _MM_TRANSPOSE4_PS(Y,F,G,H);
1086 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
1087 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
1088 fvdw6 = _mm_mul_ps(c6_00,FF);
1090 /* CUBIC SPLINE TABLE REPULSION */
1091 vfitab = _mm_add_epi32(vfitab,ifour);
1092 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
1093 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
1094 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
1095 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
1096 _MM_TRANSPOSE4_PS(Y,F,G,H);
1097 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
1098 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
1099 fvdw12 = _mm_mul_ps(c12_00,FF);
1100 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
1102 fscal = fvdw;
1104 fscal = _mm_andnot_ps(dummy_mask,fscal);
1106 /* Update vectorial force */
1107 fix0 = _mm_macc_ps(dx00,fscal,fix0);
1108 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
1109 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
1111 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
1112 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
1113 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
1115 /**************************
1116 * CALCULATE INTERACTIONS *
1117 **************************/
1119 /* Compute parameters for interactions between i and j atoms */
1120 qq10 = _mm_mul_ps(iq1,jq0);
1122 /* REACTION-FIELD ELECTROSTATICS */
1123 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
1125 fscal = felec;
1127 fscal = _mm_andnot_ps(dummy_mask,fscal);
1129 /* Update vectorial force */
1130 fix1 = _mm_macc_ps(dx10,fscal,fix1);
1131 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
1132 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
1134 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
1135 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
1136 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
1138 /**************************
1139 * CALCULATE INTERACTIONS *
1140 **************************/
1142 /* Compute parameters for interactions between i and j atoms */
1143 qq20 = _mm_mul_ps(iq2,jq0);
1145 /* REACTION-FIELD ELECTROSTATICS */
1146 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
1148 fscal = felec;
1150 fscal = _mm_andnot_ps(dummy_mask,fscal);
1152 /* Update vectorial force */
1153 fix2 = _mm_macc_ps(dx20,fscal,fix2);
1154 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
1155 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
1157 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
1158 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
1159 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
1161 /**************************
1162 * CALCULATE INTERACTIONS *
1163 **************************/
1165 /* Compute parameters for interactions between i and j atoms */
1166 qq30 = _mm_mul_ps(iq3,jq0);
1168 /* REACTION-FIELD ELECTROSTATICS */
1169 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
1171 fscal = felec;
1173 fscal = _mm_andnot_ps(dummy_mask,fscal);
1175 /* Update vectorial force */
1176 fix3 = _mm_macc_ps(dx30,fscal,fix3);
1177 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
1178 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
1180 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
1181 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
1182 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
1184 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1185 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1186 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1187 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1189 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1191 /* Inner loop uses 142 flops */
1194 /* End of innermost loop */
1196 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1197 f+i_coord_offset,fshift+i_shift_offset);
1199 /* Increment number of inner iterations */
1200 inneriter += j_index_end - j_index_start;
1202 /* Outer loop uses 24 flops */
1205 /* Increment number of outer iterations */
1206 outeriter += nri;
1208 /* Update outer/inner flops */
1210 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*142);