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_ElecRFCut_VdwLJSw_GeomW4P1_avx_128_fma_single.c
blobcc9e15b5723de47e55753a1d21439643771f71f4
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_ElecRFCut_VdwLJSw_GeomW4P1_VF_avx_128_fma_single
53 * Electrostatics interaction: ReactionField
54 * VdW interaction: LennardJones
55 * Geometry: Water4-Particle
56 * Calculate force/pot: PotentialAndForce
58 void
59 nb_kernel_ElecRFCut_VdwLJSw_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 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
107 real rswitch_scalar,d_scalar;
108 __m128 dummy_mask,cutoff_mask;
109 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
110 __m128 one = _mm_set1_ps(1.0);
111 __m128 two = _mm_set1_ps(2.0);
112 x = xx[0];
113 f = ff[0];
115 nri = nlist->nri;
116 iinr = nlist->iinr;
117 jindex = nlist->jindex;
118 jjnr = nlist->jjnr;
119 shiftidx = nlist->shift;
120 gid = nlist->gid;
121 shiftvec = fr->shift_vec[0];
122 fshift = fr->fshift[0];
123 facel = _mm_set1_ps(fr->epsfac);
124 charge = mdatoms->chargeA;
125 krf = _mm_set1_ps(fr->ic->k_rf);
126 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
127 crf = _mm_set1_ps(fr->ic->c_rf);
128 nvdwtype = fr->ntype;
129 vdwparam = fr->nbfp;
130 vdwtype = mdatoms->typeA;
132 /* Setup water-specific parameters */
133 inr = nlist->iinr[0];
134 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
135 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
136 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
137 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
139 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
140 rcutoff_scalar = fr->rcoulomb;
141 rcutoff = _mm_set1_ps(rcutoff_scalar);
142 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
144 rswitch_scalar = fr->rvdw_switch;
145 rswitch = _mm_set1_ps(rswitch_scalar);
146 /* Setup switch parameters */
147 d_scalar = rcutoff_scalar-rswitch_scalar;
148 d = _mm_set1_ps(d_scalar);
149 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
150 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
151 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
152 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
153 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
154 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
156 /* Avoid stupid compiler warnings */
157 jnrA = jnrB = jnrC = jnrD = 0;
158 j_coord_offsetA = 0;
159 j_coord_offsetB = 0;
160 j_coord_offsetC = 0;
161 j_coord_offsetD = 0;
163 outeriter = 0;
164 inneriter = 0;
166 for(iidx=0;iidx<4*DIM;iidx++)
168 scratch[iidx] = 0.0;
171 /* Start outer loop over neighborlists */
172 for(iidx=0; iidx<nri; iidx++)
174 /* Load shift vector for this list */
175 i_shift_offset = DIM*shiftidx[iidx];
177 /* Load limits for loop over neighbors */
178 j_index_start = jindex[iidx];
179 j_index_end = jindex[iidx+1];
181 /* Get outer coordinate index */
182 inr = iinr[iidx];
183 i_coord_offset = DIM*inr;
185 /* Load i particle coords and add shift vector */
186 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
187 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
189 fix0 = _mm_setzero_ps();
190 fiy0 = _mm_setzero_ps();
191 fiz0 = _mm_setzero_ps();
192 fix1 = _mm_setzero_ps();
193 fiy1 = _mm_setzero_ps();
194 fiz1 = _mm_setzero_ps();
195 fix2 = _mm_setzero_ps();
196 fiy2 = _mm_setzero_ps();
197 fiz2 = _mm_setzero_ps();
198 fix3 = _mm_setzero_ps();
199 fiy3 = _mm_setzero_ps();
200 fiz3 = _mm_setzero_ps();
202 /* Reset potential sums */
203 velecsum = _mm_setzero_ps();
204 vvdwsum = _mm_setzero_ps();
206 /* Start inner kernel loop */
207 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
210 /* Get j neighbor index, and coordinate index */
211 jnrA = jjnr[jidx];
212 jnrB = jjnr[jidx+1];
213 jnrC = jjnr[jidx+2];
214 jnrD = jjnr[jidx+3];
215 j_coord_offsetA = DIM*jnrA;
216 j_coord_offsetB = DIM*jnrB;
217 j_coord_offsetC = DIM*jnrC;
218 j_coord_offsetD = DIM*jnrD;
220 /* load j atom coordinates */
221 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
222 x+j_coord_offsetC,x+j_coord_offsetD,
223 &jx0,&jy0,&jz0);
225 /* Calculate displacement vector */
226 dx00 = _mm_sub_ps(ix0,jx0);
227 dy00 = _mm_sub_ps(iy0,jy0);
228 dz00 = _mm_sub_ps(iz0,jz0);
229 dx10 = _mm_sub_ps(ix1,jx0);
230 dy10 = _mm_sub_ps(iy1,jy0);
231 dz10 = _mm_sub_ps(iz1,jz0);
232 dx20 = _mm_sub_ps(ix2,jx0);
233 dy20 = _mm_sub_ps(iy2,jy0);
234 dz20 = _mm_sub_ps(iz2,jz0);
235 dx30 = _mm_sub_ps(ix3,jx0);
236 dy30 = _mm_sub_ps(iy3,jy0);
237 dz30 = _mm_sub_ps(iz3,jz0);
239 /* Calculate squared distance and things based on it */
240 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
241 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
242 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
243 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
245 rinv00 = gmx_mm_invsqrt_ps(rsq00);
246 rinv10 = gmx_mm_invsqrt_ps(rsq10);
247 rinv20 = gmx_mm_invsqrt_ps(rsq20);
248 rinv30 = gmx_mm_invsqrt_ps(rsq30);
250 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
251 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
252 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
253 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
255 /* Load parameters for j particles */
256 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
257 charge+jnrC+0,charge+jnrD+0);
258 vdwjidx0A = 2*vdwtype[jnrA+0];
259 vdwjidx0B = 2*vdwtype[jnrB+0];
260 vdwjidx0C = 2*vdwtype[jnrC+0];
261 vdwjidx0D = 2*vdwtype[jnrD+0];
263 fjx0 = _mm_setzero_ps();
264 fjy0 = _mm_setzero_ps();
265 fjz0 = _mm_setzero_ps();
267 /**************************
268 * CALCULATE INTERACTIONS *
269 **************************/
271 if (gmx_mm_any_lt(rsq00,rcutoff2))
274 r00 = _mm_mul_ps(rsq00,rinv00);
276 /* Compute parameters for interactions between i and j atoms */
277 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
278 vdwparam+vdwioffset0+vdwjidx0B,
279 vdwparam+vdwioffset0+vdwjidx0C,
280 vdwparam+vdwioffset0+vdwjidx0D,
281 &c6_00,&c12_00);
283 /* LENNARD-JONES DISPERSION/REPULSION */
285 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
286 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
287 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
288 vvdw = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
289 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
291 d = _mm_sub_ps(r00,rswitch);
292 d = _mm_max_ps(d,_mm_setzero_ps());
293 d2 = _mm_mul_ps(d,d);
294 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
296 dsw = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
298 /* Evaluate switch function */
299 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
300 fvdw = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
301 vvdw = _mm_mul_ps(vvdw,sw);
302 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
304 /* Update potential sum for this i atom from the interaction with this j atom. */
305 vvdw = _mm_and_ps(vvdw,cutoff_mask);
306 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
308 fscal = fvdw;
310 fscal = _mm_and_ps(fscal,cutoff_mask);
312 /* Update vectorial force */
313 fix0 = _mm_macc_ps(dx00,fscal,fix0);
314 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
315 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
317 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
318 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
319 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
323 /**************************
324 * CALCULATE INTERACTIONS *
325 **************************/
327 if (gmx_mm_any_lt(rsq10,rcutoff2))
330 /* Compute parameters for interactions between i and j atoms */
331 qq10 = _mm_mul_ps(iq1,jq0);
333 /* REACTION-FIELD ELECTROSTATICS */
334 velec = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
335 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
337 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
339 /* Update potential sum for this i atom from the interaction with this j atom. */
340 velec = _mm_and_ps(velec,cutoff_mask);
341 velecsum = _mm_add_ps(velecsum,velec);
343 fscal = felec;
345 fscal = _mm_and_ps(fscal,cutoff_mask);
347 /* Update vectorial force */
348 fix1 = _mm_macc_ps(dx10,fscal,fix1);
349 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
350 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
352 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
353 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
354 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
358 /**************************
359 * CALCULATE INTERACTIONS *
360 **************************/
362 if (gmx_mm_any_lt(rsq20,rcutoff2))
365 /* Compute parameters for interactions between i and j atoms */
366 qq20 = _mm_mul_ps(iq2,jq0);
368 /* REACTION-FIELD ELECTROSTATICS */
369 velec = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
370 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
372 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
374 /* Update potential sum for this i atom from the interaction with this j atom. */
375 velec = _mm_and_ps(velec,cutoff_mask);
376 velecsum = _mm_add_ps(velecsum,velec);
378 fscal = felec;
380 fscal = _mm_and_ps(fscal,cutoff_mask);
382 /* Update vectorial force */
383 fix2 = _mm_macc_ps(dx20,fscal,fix2);
384 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
385 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
387 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
388 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
389 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
393 /**************************
394 * CALCULATE INTERACTIONS *
395 **************************/
397 if (gmx_mm_any_lt(rsq30,rcutoff2))
400 /* Compute parameters for interactions between i and j atoms */
401 qq30 = _mm_mul_ps(iq3,jq0);
403 /* REACTION-FIELD ELECTROSTATICS */
404 velec = _mm_mul_ps(qq30,_mm_sub_ps(_mm_macc_ps(krf,rsq30,rinv30),crf));
405 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
407 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
409 /* Update potential sum for this i atom from the interaction with this j atom. */
410 velec = _mm_and_ps(velec,cutoff_mask);
411 velecsum = _mm_add_ps(velecsum,velec);
413 fscal = felec;
415 fscal = _mm_and_ps(fscal,cutoff_mask);
417 /* Update vectorial force */
418 fix3 = _mm_macc_ps(dx30,fscal,fix3);
419 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
420 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
422 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
423 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
424 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
428 fjptrA = f+j_coord_offsetA;
429 fjptrB = f+j_coord_offsetB;
430 fjptrC = f+j_coord_offsetC;
431 fjptrD = f+j_coord_offsetD;
433 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
435 /* Inner loop uses 179 flops */
438 if(jidx<j_index_end)
441 /* Get j neighbor index, and coordinate index */
442 jnrlistA = jjnr[jidx];
443 jnrlistB = jjnr[jidx+1];
444 jnrlistC = jjnr[jidx+2];
445 jnrlistD = jjnr[jidx+3];
446 /* Sign of each element will be negative for non-real atoms.
447 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
448 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
450 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
451 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
452 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
453 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
454 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
455 j_coord_offsetA = DIM*jnrA;
456 j_coord_offsetB = DIM*jnrB;
457 j_coord_offsetC = DIM*jnrC;
458 j_coord_offsetD = DIM*jnrD;
460 /* load j atom coordinates */
461 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
462 x+j_coord_offsetC,x+j_coord_offsetD,
463 &jx0,&jy0,&jz0);
465 /* Calculate displacement vector */
466 dx00 = _mm_sub_ps(ix0,jx0);
467 dy00 = _mm_sub_ps(iy0,jy0);
468 dz00 = _mm_sub_ps(iz0,jz0);
469 dx10 = _mm_sub_ps(ix1,jx0);
470 dy10 = _mm_sub_ps(iy1,jy0);
471 dz10 = _mm_sub_ps(iz1,jz0);
472 dx20 = _mm_sub_ps(ix2,jx0);
473 dy20 = _mm_sub_ps(iy2,jy0);
474 dz20 = _mm_sub_ps(iz2,jz0);
475 dx30 = _mm_sub_ps(ix3,jx0);
476 dy30 = _mm_sub_ps(iy3,jy0);
477 dz30 = _mm_sub_ps(iz3,jz0);
479 /* Calculate squared distance and things based on it */
480 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
481 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
482 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
483 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
485 rinv00 = gmx_mm_invsqrt_ps(rsq00);
486 rinv10 = gmx_mm_invsqrt_ps(rsq10);
487 rinv20 = gmx_mm_invsqrt_ps(rsq20);
488 rinv30 = gmx_mm_invsqrt_ps(rsq30);
490 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
491 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
492 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
493 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
495 /* Load parameters for j particles */
496 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
497 charge+jnrC+0,charge+jnrD+0);
498 vdwjidx0A = 2*vdwtype[jnrA+0];
499 vdwjidx0B = 2*vdwtype[jnrB+0];
500 vdwjidx0C = 2*vdwtype[jnrC+0];
501 vdwjidx0D = 2*vdwtype[jnrD+0];
503 fjx0 = _mm_setzero_ps();
504 fjy0 = _mm_setzero_ps();
505 fjz0 = _mm_setzero_ps();
507 /**************************
508 * CALCULATE INTERACTIONS *
509 **************************/
511 if (gmx_mm_any_lt(rsq00,rcutoff2))
514 r00 = _mm_mul_ps(rsq00,rinv00);
515 r00 = _mm_andnot_ps(dummy_mask,r00);
517 /* Compute parameters for interactions between i and j atoms */
518 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
519 vdwparam+vdwioffset0+vdwjidx0B,
520 vdwparam+vdwioffset0+vdwjidx0C,
521 vdwparam+vdwioffset0+vdwjidx0D,
522 &c6_00,&c12_00);
524 /* LENNARD-JONES DISPERSION/REPULSION */
526 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
527 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
528 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
529 vvdw = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
530 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
532 d = _mm_sub_ps(r00,rswitch);
533 d = _mm_max_ps(d,_mm_setzero_ps());
534 d2 = _mm_mul_ps(d,d);
535 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
537 dsw = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
539 /* Evaluate switch function */
540 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
541 fvdw = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
542 vvdw = _mm_mul_ps(vvdw,sw);
543 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
545 /* Update potential sum for this i atom from the interaction with this j atom. */
546 vvdw = _mm_and_ps(vvdw,cutoff_mask);
547 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
548 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
550 fscal = fvdw;
552 fscal = _mm_and_ps(fscal,cutoff_mask);
554 fscal = _mm_andnot_ps(dummy_mask,fscal);
556 /* Update vectorial force */
557 fix0 = _mm_macc_ps(dx00,fscal,fix0);
558 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
559 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
561 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
562 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
563 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
567 /**************************
568 * CALCULATE INTERACTIONS *
569 **************************/
571 if (gmx_mm_any_lt(rsq10,rcutoff2))
574 /* Compute parameters for interactions between i and j atoms */
575 qq10 = _mm_mul_ps(iq1,jq0);
577 /* REACTION-FIELD ELECTROSTATICS */
578 velec = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
579 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
581 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
583 /* Update potential sum for this i atom from the interaction with this j atom. */
584 velec = _mm_and_ps(velec,cutoff_mask);
585 velec = _mm_andnot_ps(dummy_mask,velec);
586 velecsum = _mm_add_ps(velecsum,velec);
588 fscal = felec;
590 fscal = _mm_and_ps(fscal,cutoff_mask);
592 fscal = _mm_andnot_ps(dummy_mask,fscal);
594 /* Update vectorial force */
595 fix1 = _mm_macc_ps(dx10,fscal,fix1);
596 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
597 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
599 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
600 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
601 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
605 /**************************
606 * CALCULATE INTERACTIONS *
607 **************************/
609 if (gmx_mm_any_lt(rsq20,rcutoff2))
612 /* Compute parameters for interactions between i and j atoms */
613 qq20 = _mm_mul_ps(iq2,jq0);
615 /* REACTION-FIELD ELECTROSTATICS */
616 velec = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
617 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
619 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
621 /* Update potential sum for this i atom from the interaction with this j atom. */
622 velec = _mm_and_ps(velec,cutoff_mask);
623 velec = _mm_andnot_ps(dummy_mask,velec);
624 velecsum = _mm_add_ps(velecsum,velec);
626 fscal = felec;
628 fscal = _mm_and_ps(fscal,cutoff_mask);
630 fscal = _mm_andnot_ps(dummy_mask,fscal);
632 /* Update vectorial force */
633 fix2 = _mm_macc_ps(dx20,fscal,fix2);
634 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
635 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
637 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
638 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
639 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
643 /**************************
644 * CALCULATE INTERACTIONS *
645 **************************/
647 if (gmx_mm_any_lt(rsq30,rcutoff2))
650 /* Compute parameters for interactions between i and j atoms */
651 qq30 = _mm_mul_ps(iq3,jq0);
653 /* REACTION-FIELD ELECTROSTATICS */
654 velec = _mm_mul_ps(qq30,_mm_sub_ps(_mm_macc_ps(krf,rsq30,rinv30),crf));
655 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
657 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
659 /* Update potential sum for this i atom from the interaction with this j atom. */
660 velec = _mm_and_ps(velec,cutoff_mask);
661 velec = _mm_andnot_ps(dummy_mask,velec);
662 velecsum = _mm_add_ps(velecsum,velec);
664 fscal = felec;
666 fscal = _mm_and_ps(fscal,cutoff_mask);
668 fscal = _mm_andnot_ps(dummy_mask,fscal);
670 /* Update vectorial force */
671 fix3 = _mm_macc_ps(dx30,fscal,fix3);
672 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
673 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
675 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
676 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
677 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
681 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
682 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
683 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
684 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
686 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
688 /* Inner loop uses 180 flops */
691 /* End of innermost loop */
693 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
694 f+i_coord_offset,fshift+i_shift_offset);
696 ggid = gid[iidx];
697 /* Update potential energies */
698 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
699 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
701 /* Increment number of inner iterations */
702 inneriter += j_index_end - j_index_start;
704 /* Outer loop uses 26 flops */
707 /* Increment number of outer iterations */
708 outeriter += nri;
710 /* Update outer/inner flops */
712 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*180);
715 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomW4P1_F_avx_128_fma_single
716 * Electrostatics interaction: ReactionField
717 * VdW interaction: LennardJones
718 * Geometry: Water4-Particle
719 * Calculate force/pot: Force
721 void
722 nb_kernel_ElecRFCut_VdwLJSw_GeomW4P1_F_avx_128_fma_single
723 (t_nblist * gmx_restrict nlist,
724 rvec * gmx_restrict xx,
725 rvec * gmx_restrict ff,
726 t_forcerec * gmx_restrict fr,
727 t_mdatoms * gmx_restrict mdatoms,
728 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
729 t_nrnb * gmx_restrict nrnb)
731 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
732 * just 0 for non-waters.
733 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
734 * jnr indices corresponding to data put in the four positions in the SIMD register.
736 int i_shift_offset,i_coord_offset,outeriter,inneriter;
737 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
738 int jnrA,jnrB,jnrC,jnrD;
739 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
740 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
741 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
742 real rcutoff_scalar;
743 real *shiftvec,*fshift,*x,*f;
744 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
745 real scratch[4*DIM];
746 __m128 fscal,rcutoff,rcutoff2,jidxall;
747 int vdwioffset0;
748 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
749 int vdwioffset1;
750 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
751 int vdwioffset2;
752 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
753 int vdwioffset3;
754 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
755 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
756 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
757 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
758 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
759 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
760 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
761 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
762 real *charge;
763 int nvdwtype;
764 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
765 int *vdwtype;
766 real *vdwparam;
767 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
768 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
769 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
770 real rswitch_scalar,d_scalar;
771 __m128 dummy_mask,cutoff_mask;
772 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
773 __m128 one = _mm_set1_ps(1.0);
774 __m128 two = _mm_set1_ps(2.0);
775 x = xx[0];
776 f = ff[0];
778 nri = nlist->nri;
779 iinr = nlist->iinr;
780 jindex = nlist->jindex;
781 jjnr = nlist->jjnr;
782 shiftidx = nlist->shift;
783 gid = nlist->gid;
784 shiftvec = fr->shift_vec[0];
785 fshift = fr->fshift[0];
786 facel = _mm_set1_ps(fr->epsfac);
787 charge = mdatoms->chargeA;
788 krf = _mm_set1_ps(fr->ic->k_rf);
789 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
790 crf = _mm_set1_ps(fr->ic->c_rf);
791 nvdwtype = fr->ntype;
792 vdwparam = fr->nbfp;
793 vdwtype = mdatoms->typeA;
795 /* Setup water-specific parameters */
796 inr = nlist->iinr[0];
797 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
798 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
799 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
800 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
802 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
803 rcutoff_scalar = fr->rcoulomb;
804 rcutoff = _mm_set1_ps(rcutoff_scalar);
805 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
807 rswitch_scalar = fr->rvdw_switch;
808 rswitch = _mm_set1_ps(rswitch_scalar);
809 /* Setup switch parameters */
810 d_scalar = rcutoff_scalar-rswitch_scalar;
811 d = _mm_set1_ps(d_scalar);
812 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
813 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
814 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
815 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
816 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
817 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
819 /* Avoid stupid compiler warnings */
820 jnrA = jnrB = jnrC = jnrD = 0;
821 j_coord_offsetA = 0;
822 j_coord_offsetB = 0;
823 j_coord_offsetC = 0;
824 j_coord_offsetD = 0;
826 outeriter = 0;
827 inneriter = 0;
829 for(iidx=0;iidx<4*DIM;iidx++)
831 scratch[iidx] = 0.0;
834 /* Start outer loop over neighborlists */
835 for(iidx=0; iidx<nri; iidx++)
837 /* Load shift vector for this list */
838 i_shift_offset = DIM*shiftidx[iidx];
840 /* Load limits for loop over neighbors */
841 j_index_start = jindex[iidx];
842 j_index_end = jindex[iidx+1];
844 /* Get outer coordinate index */
845 inr = iinr[iidx];
846 i_coord_offset = DIM*inr;
848 /* Load i particle coords and add shift vector */
849 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
850 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
852 fix0 = _mm_setzero_ps();
853 fiy0 = _mm_setzero_ps();
854 fiz0 = _mm_setzero_ps();
855 fix1 = _mm_setzero_ps();
856 fiy1 = _mm_setzero_ps();
857 fiz1 = _mm_setzero_ps();
858 fix2 = _mm_setzero_ps();
859 fiy2 = _mm_setzero_ps();
860 fiz2 = _mm_setzero_ps();
861 fix3 = _mm_setzero_ps();
862 fiy3 = _mm_setzero_ps();
863 fiz3 = _mm_setzero_ps();
865 /* Start inner kernel loop */
866 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
869 /* Get j neighbor index, and coordinate index */
870 jnrA = jjnr[jidx];
871 jnrB = jjnr[jidx+1];
872 jnrC = jjnr[jidx+2];
873 jnrD = jjnr[jidx+3];
874 j_coord_offsetA = DIM*jnrA;
875 j_coord_offsetB = DIM*jnrB;
876 j_coord_offsetC = DIM*jnrC;
877 j_coord_offsetD = DIM*jnrD;
879 /* load j atom coordinates */
880 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
881 x+j_coord_offsetC,x+j_coord_offsetD,
882 &jx0,&jy0,&jz0);
884 /* Calculate displacement vector */
885 dx00 = _mm_sub_ps(ix0,jx0);
886 dy00 = _mm_sub_ps(iy0,jy0);
887 dz00 = _mm_sub_ps(iz0,jz0);
888 dx10 = _mm_sub_ps(ix1,jx0);
889 dy10 = _mm_sub_ps(iy1,jy0);
890 dz10 = _mm_sub_ps(iz1,jz0);
891 dx20 = _mm_sub_ps(ix2,jx0);
892 dy20 = _mm_sub_ps(iy2,jy0);
893 dz20 = _mm_sub_ps(iz2,jz0);
894 dx30 = _mm_sub_ps(ix3,jx0);
895 dy30 = _mm_sub_ps(iy3,jy0);
896 dz30 = _mm_sub_ps(iz3,jz0);
898 /* Calculate squared distance and things based on it */
899 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
900 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
901 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
902 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
904 rinv00 = gmx_mm_invsqrt_ps(rsq00);
905 rinv10 = gmx_mm_invsqrt_ps(rsq10);
906 rinv20 = gmx_mm_invsqrt_ps(rsq20);
907 rinv30 = gmx_mm_invsqrt_ps(rsq30);
909 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
910 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
911 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
912 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
914 /* Load parameters for j particles */
915 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
916 charge+jnrC+0,charge+jnrD+0);
917 vdwjidx0A = 2*vdwtype[jnrA+0];
918 vdwjidx0B = 2*vdwtype[jnrB+0];
919 vdwjidx0C = 2*vdwtype[jnrC+0];
920 vdwjidx0D = 2*vdwtype[jnrD+0];
922 fjx0 = _mm_setzero_ps();
923 fjy0 = _mm_setzero_ps();
924 fjz0 = _mm_setzero_ps();
926 /**************************
927 * CALCULATE INTERACTIONS *
928 **************************/
930 if (gmx_mm_any_lt(rsq00,rcutoff2))
933 r00 = _mm_mul_ps(rsq00,rinv00);
935 /* Compute parameters for interactions between i and j atoms */
936 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
937 vdwparam+vdwioffset0+vdwjidx0B,
938 vdwparam+vdwioffset0+vdwjidx0C,
939 vdwparam+vdwioffset0+vdwjidx0D,
940 &c6_00,&c12_00);
942 /* LENNARD-JONES DISPERSION/REPULSION */
944 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
945 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
946 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
947 vvdw = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
948 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
950 d = _mm_sub_ps(r00,rswitch);
951 d = _mm_max_ps(d,_mm_setzero_ps());
952 d2 = _mm_mul_ps(d,d);
953 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
955 dsw = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
957 /* Evaluate switch function */
958 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
959 fvdw = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
960 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
962 fscal = fvdw;
964 fscal = _mm_and_ps(fscal,cutoff_mask);
966 /* Update vectorial force */
967 fix0 = _mm_macc_ps(dx00,fscal,fix0);
968 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
969 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
971 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
972 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
973 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
977 /**************************
978 * CALCULATE INTERACTIONS *
979 **************************/
981 if (gmx_mm_any_lt(rsq10,rcutoff2))
984 /* Compute parameters for interactions between i and j atoms */
985 qq10 = _mm_mul_ps(iq1,jq0);
987 /* REACTION-FIELD ELECTROSTATICS */
988 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
990 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
992 fscal = felec;
994 fscal = _mm_and_ps(fscal,cutoff_mask);
996 /* Update vectorial force */
997 fix1 = _mm_macc_ps(dx10,fscal,fix1);
998 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
999 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
1001 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
1002 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
1003 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
1007 /**************************
1008 * CALCULATE INTERACTIONS *
1009 **************************/
1011 if (gmx_mm_any_lt(rsq20,rcutoff2))
1014 /* Compute parameters for interactions between i and j atoms */
1015 qq20 = _mm_mul_ps(iq2,jq0);
1017 /* REACTION-FIELD ELECTROSTATICS */
1018 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
1020 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
1022 fscal = felec;
1024 fscal = _mm_and_ps(fscal,cutoff_mask);
1026 /* Update vectorial force */
1027 fix2 = _mm_macc_ps(dx20,fscal,fix2);
1028 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
1029 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
1031 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
1032 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
1033 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
1037 /**************************
1038 * CALCULATE INTERACTIONS *
1039 **************************/
1041 if (gmx_mm_any_lt(rsq30,rcutoff2))
1044 /* Compute parameters for interactions between i and j atoms */
1045 qq30 = _mm_mul_ps(iq3,jq0);
1047 /* REACTION-FIELD ELECTROSTATICS */
1048 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
1050 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
1052 fscal = felec;
1054 fscal = _mm_and_ps(fscal,cutoff_mask);
1056 /* Update vectorial force */
1057 fix3 = _mm_macc_ps(dx30,fscal,fix3);
1058 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
1059 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
1061 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
1062 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
1063 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
1067 fjptrA = f+j_coord_offsetA;
1068 fjptrB = f+j_coord_offsetB;
1069 fjptrC = f+j_coord_offsetC;
1070 fjptrD = f+j_coord_offsetD;
1072 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1074 /* Inner loop uses 158 flops */
1077 if(jidx<j_index_end)
1080 /* Get j neighbor index, and coordinate index */
1081 jnrlistA = jjnr[jidx];
1082 jnrlistB = jjnr[jidx+1];
1083 jnrlistC = jjnr[jidx+2];
1084 jnrlistD = jjnr[jidx+3];
1085 /* Sign of each element will be negative for non-real atoms.
1086 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1087 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1089 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1090 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1091 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1092 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1093 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1094 j_coord_offsetA = DIM*jnrA;
1095 j_coord_offsetB = DIM*jnrB;
1096 j_coord_offsetC = DIM*jnrC;
1097 j_coord_offsetD = DIM*jnrD;
1099 /* load j atom coordinates */
1100 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1101 x+j_coord_offsetC,x+j_coord_offsetD,
1102 &jx0,&jy0,&jz0);
1104 /* Calculate displacement vector */
1105 dx00 = _mm_sub_ps(ix0,jx0);
1106 dy00 = _mm_sub_ps(iy0,jy0);
1107 dz00 = _mm_sub_ps(iz0,jz0);
1108 dx10 = _mm_sub_ps(ix1,jx0);
1109 dy10 = _mm_sub_ps(iy1,jy0);
1110 dz10 = _mm_sub_ps(iz1,jz0);
1111 dx20 = _mm_sub_ps(ix2,jx0);
1112 dy20 = _mm_sub_ps(iy2,jy0);
1113 dz20 = _mm_sub_ps(iz2,jz0);
1114 dx30 = _mm_sub_ps(ix3,jx0);
1115 dy30 = _mm_sub_ps(iy3,jy0);
1116 dz30 = _mm_sub_ps(iz3,jz0);
1118 /* Calculate squared distance and things based on it */
1119 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1120 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1121 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1122 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
1124 rinv00 = gmx_mm_invsqrt_ps(rsq00);
1125 rinv10 = gmx_mm_invsqrt_ps(rsq10);
1126 rinv20 = gmx_mm_invsqrt_ps(rsq20);
1127 rinv30 = gmx_mm_invsqrt_ps(rsq30);
1129 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
1130 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
1131 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
1132 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
1134 /* Load parameters for j particles */
1135 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1136 charge+jnrC+0,charge+jnrD+0);
1137 vdwjidx0A = 2*vdwtype[jnrA+0];
1138 vdwjidx0B = 2*vdwtype[jnrB+0];
1139 vdwjidx0C = 2*vdwtype[jnrC+0];
1140 vdwjidx0D = 2*vdwtype[jnrD+0];
1142 fjx0 = _mm_setzero_ps();
1143 fjy0 = _mm_setzero_ps();
1144 fjz0 = _mm_setzero_ps();
1146 /**************************
1147 * CALCULATE INTERACTIONS *
1148 **************************/
1150 if (gmx_mm_any_lt(rsq00,rcutoff2))
1153 r00 = _mm_mul_ps(rsq00,rinv00);
1154 r00 = _mm_andnot_ps(dummy_mask,r00);
1156 /* Compute parameters for interactions between i and j atoms */
1157 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1158 vdwparam+vdwioffset0+vdwjidx0B,
1159 vdwparam+vdwioffset0+vdwjidx0C,
1160 vdwparam+vdwioffset0+vdwjidx0D,
1161 &c6_00,&c12_00);
1163 /* LENNARD-JONES DISPERSION/REPULSION */
1165 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1166 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
1167 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
1168 vvdw = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
1169 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
1171 d = _mm_sub_ps(r00,rswitch);
1172 d = _mm_max_ps(d,_mm_setzero_ps());
1173 d2 = _mm_mul_ps(d,d);
1174 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1176 dsw = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1178 /* Evaluate switch function */
1179 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1180 fvdw = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
1181 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
1183 fscal = fvdw;
1185 fscal = _mm_and_ps(fscal,cutoff_mask);
1187 fscal = _mm_andnot_ps(dummy_mask,fscal);
1189 /* Update vectorial force */
1190 fix0 = _mm_macc_ps(dx00,fscal,fix0);
1191 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
1192 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
1194 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
1195 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
1196 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
1200 /**************************
1201 * CALCULATE INTERACTIONS *
1202 **************************/
1204 if (gmx_mm_any_lt(rsq10,rcutoff2))
1207 /* Compute parameters for interactions between i and j atoms */
1208 qq10 = _mm_mul_ps(iq1,jq0);
1210 /* REACTION-FIELD ELECTROSTATICS */
1211 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
1213 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
1215 fscal = felec;
1217 fscal = _mm_and_ps(fscal,cutoff_mask);
1219 fscal = _mm_andnot_ps(dummy_mask,fscal);
1221 /* Update vectorial force */
1222 fix1 = _mm_macc_ps(dx10,fscal,fix1);
1223 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
1224 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
1226 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
1227 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
1228 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
1232 /**************************
1233 * CALCULATE INTERACTIONS *
1234 **************************/
1236 if (gmx_mm_any_lt(rsq20,rcutoff2))
1239 /* Compute parameters for interactions between i and j atoms */
1240 qq20 = _mm_mul_ps(iq2,jq0);
1242 /* REACTION-FIELD ELECTROSTATICS */
1243 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
1245 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
1247 fscal = felec;
1249 fscal = _mm_and_ps(fscal,cutoff_mask);
1251 fscal = _mm_andnot_ps(dummy_mask,fscal);
1253 /* Update vectorial force */
1254 fix2 = _mm_macc_ps(dx20,fscal,fix2);
1255 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
1256 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
1258 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
1259 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
1260 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
1264 /**************************
1265 * CALCULATE INTERACTIONS *
1266 **************************/
1268 if (gmx_mm_any_lt(rsq30,rcutoff2))
1271 /* Compute parameters for interactions between i and j atoms */
1272 qq30 = _mm_mul_ps(iq3,jq0);
1274 /* REACTION-FIELD ELECTROSTATICS */
1275 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
1277 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
1279 fscal = felec;
1281 fscal = _mm_and_ps(fscal,cutoff_mask);
1283 fscal = _mm_andnot_ps(dummy_mask,fscal);
1285 /* Update vectorial force */
1286 fix3 = _mm_macc_ps(dx30,fscal,fix3);
1287 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
1288 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
1290 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
1291 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
1292 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
1296 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1297 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1298 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1299 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1301 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1303 /* Inner loop uses 159 flops */
1306 /* End of innermost loop */
1308 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1309 f+i_coord_offset,fshift+i_shift_offset);
1311 /* Increment number of inner iterations */
1312 inneriter += j_index_end - j_index_start;
1314 /* Outer loop uses 24 flops */
1317 /* Increment number of outer iterations */
1318 outeriter += nri;
1320 /* Update outer/inner flops */
1322 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*159);