Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_avx_128_fma_single.c
blob70a01ff562563072631aff7506c79a224917b320
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS avx_128_fma_single kernel generator.
38 #include "config.h"
40 #include <math.h>
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
47 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
48 #include "kernelutil_x86_avx_128_fma_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_VF_avx_128_fma_single
52 * Electrostatics interaction: ReactionField
53 * VdW interaction: LennardJones
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
57 void
58 nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_VF_avx_128_fma_single
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,C,D refer to j loop unrolling done with AVX_128, e.g. for the four 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;
74 int jnrA,jnrB,jnrC,jnrD;
75 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
78 real rcutoff_scalar;
79 real *shiftvec,*fshift,*x,*f;
80 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
81 real scratch[4*DIM];
82 __m128 fscal,rcutoff,rcutoff2,jidxall;
83 int vdwioffset0;
84 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
86 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
88 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
89 real *charge;
90 int nvdwtype;
91 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
92 int *vdwtype;
93 real *vdwparam;
94 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
95 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
96 __m128 dummy_mask,cutoff_mask;
97 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
98 __m128 one = _mm_set1_ps(1.0);
99 __m128 two = _mm_set1_ps(2.0);
100 x = xx[0];
101 f = ff[0];
103 nri = nlist->nri;
104 iinr = nlist->iinr;
105 jindex = nlist->jindex;
106 jjnr = nlist->jjnr;
107 shiftidx = nlist->shift;
108 gid = nlist->gid;
109 shiftvec = fr->shift_vec[0];
110 fshift = fr->fshift[0];
111 facel = _mm_set1_ps(fr->epsfac);
112 charge = mdatoms->chargeA;
113 krf = _mm_set1_ps(fr->ic->k_rf);
114 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
115 crf = _mm_set1_ps(fr->ic->c_rf);
116 nvdwtype = fr->ntype;
117 vdwparam = fr->nbfp;
118 vdwtype = mdatoms->typeA;
120 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
121 rcutoff_scalar = fr->rcoulomb;
122 rcutoff = _mm_set1_ps(rcutoff_scalar);
123 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
125 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
126 rvdw = _mm_set1_ps(fr->rvdw);
128 /* Avoid stupid compiler warnings */
129 jnrA = jnrB = jnrC = jnrD = 0;
130 j_coord_offsetA = 0;
131 j_coord_offsetB = 0;
132 j_coord_offsetC = 0;
133 j_coord_offsetD = 0;
135 outeriter = 0;
136 inneriter = 0;
138 for(iidx=0;iidx<4*DIM;iidx++)
140 scratch[iidx] = 0.0;
143 /* Start outer loop over neighborlists */
144 for(iidx=0; iidx<nri; iidx++)
146 /* Load shift vector for this list */
147 i_shift_offset = DIM*shiftidx[iidx];
149 /* Load limits for loop over neighbors */
150 j_index_start = jindex[iidx];
151 j_index_end = jindex[iidx+1];
153 /* Get outer coordinate index */
154 inr = iinr[iidx];
155 i_coord_offset = DIM*inr;
157 /* Load i particle coords and add shift vector */
158 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
160 fix0 = _mm_setzero_ps();
161 fiy0 = _mm_setzero_ps();
162 fiz0 = _mm_setzero_ps();
164 /* Load parameters for i particles */
165 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
166 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
168 /* Reset potential sums */
169 velecsum = _mm_setzero_ps();
170 vvdwsum = _mm_setzero_ps();
172 /* Start inner kernel loop */
173 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
176 /* Get j neighbor index, and coordinate index */
177 jnrA = jjnr[jidx];
178 jnrB = jjnr[jidx+1];
179 jnrC = jjnr[jidx+2];
180 jnrD = jjnr[jidx+3];
181 j_coord_offsetA = DIM*jnrA;
182 j_coord_offsetB = DIM*jnrB;
183 j_coord_offsetC = DIM*jnrC;
184 j_coord_offsetD = DIM*jnrD;
186 /* load j atom coordinates */
187 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
188 x+j_coord_offsetC,x+j_coord_offsetD,
189 &jx0,&jy0,&jz0);
191 /* Calculate displacement vector */
192 dx00 = _mm_sub_ps(ix0,jx0);
193 dy00 = _mm_sub_ps(iy0,jy0);
194 dz00 = _mm_sub_ps(iz0,jz0);
196 /* Calculate squared distance and things based on it */
197 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
199 rinv00 = gmx_mm_invsqrt_ps(rsq00);
201 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
203 /* Load parameters for j particles */
204 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
205 charge+jnrC+0,charge+jnrD+0);
206 vdwjidx0A = 2*vdwtype[jnrA+0];
207 vdwjidx0B = 2*vdwtype[jnrB+0];
208 vdwjidx0C = 2*vdwtype[jnrC+0];
209 vdwjidx0D = 2*vdwtype[jnrD+0];
211 /**************************
212 * CALCULATE INTERACTIONS *
213 **************************/
215 if (gmx_mm_any_lt(rsq00,rcutoff2))
218 /* Compute parameters for interactions between i and j atoms */
219 qq00 = _mm_mul_ps(iq0,jq0);
220 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
221 vdwparam+vdwioffset0+vdwjidx0B,
222 vdwparam+vdwioffset0+vdwjidx0C,
223 vdwparam+vdwioffset0+vdwjidx0D,
224 &c6_00,&c12_00);
226 /* REACTION-FIELD ELECTROSTATICS */
227 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_macc_ps(krf,rsq00,rinv00),crf));
228 felec = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
230 /* LENNARD-JONES DISPERSION/REPULSION */
232 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
233 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
234 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
235 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
236 _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
237 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
239 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
241 /* Update potential sum for this i atom from the interaction with this j atom. */
242 velec = _mm_and_ps(velec,cutoff_mask);
243 velecsum = _mm_add_ps(velecsum,velec);
244 vvdw = _mm_and_ps(vvdw,cutoff_mask);
245 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
247 fscal = _mm_add_ps(felec,fvdw);
249 fscal = _mm_and_ps(fscal,cutoff_mask);
251 /* Update vectorial force */
252 fix0 = _mm_macc_ps(dx00,fscal,fix0);
253 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
254 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
256 fjptrA = f+j_coord_offsetA;
257 fjptrB = f+j_coord_offsetB;
258 fjptrC = f+j_coord_offsetC;
259 fjptrD = f+j_coord_offsetD;
260 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
261 _mm_mul_ps(dx00,fscal),
262 _mm_mul_ps(dy00,fscal),
263 _mm_mul_ps(dz00,fscal));
267 /* Inner loop uses 57 flops */
270 if(jidx<j_index_end)
273 /* Get j neighbor index, and coordinate index */
274 jnrlistA = jjnr[jidx];
275 jnrlistB = jjnr[jidx+1];
276 jnrlistC = jjnr[jidx+2];
277 jnrlistD = jjnr[jidx+3];
278 /* Sign of each element will be negative for non-real atoms.
279 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
280 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
282 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
283 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
284 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
285 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
286 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
287 j_coord_offsetA = DIM*jnrA;
288 j_coord_offsetB = DIM*jnrB;
289 j_coord_offsetC = DIM*jnrC;
290 j_coord_offsetD = DIM*jnrD;
292 /* load j atom coordinates */
293 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
294 x+j_coord_offsetC,x+j_coord_offsetD,
295 &jx0,&jy0,&jz0);
297 /* Calculate displacement vector */
298 dx00 = _mm_sub_ps(ix0,jx0);
299 dy00 = _mm_sub_ps(iy0,jy0);
300 dz00 = _mm_sub_ps(iz0,jz0);
302 /* Calculate squared distance and things based on it */
303 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
305 rinv00 = gmx_mm_invsqrt_ps(rsq00);
307 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
309 /* Load parameters for j particles */
310 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
311 charge+jnrC+0,charge+jnrD+0);
312 vdwjidx0A = 2*vdwtype[jnrA+0];
313 vdwjidx0B = 2*vdwtype[jnrB+0];
314 vdwjidx0C = 2*vdwtype[jnrC+0];
315 vdwjidx0D = 2*vdwtype[jnrD+0];
317 /**************************
318 * CALCULATE INTERACTIONS *
319 **************************/
321 if (gmx_mm_any_lt(rsq00,rcutoff2))
324 /* Compute parameters for interactions between i and j atoms */
325 qq00 = _mm_mul_ps(iq0,jq0);
326 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
327 vdwparam+vdwioffset0+vdwjidx0B,
328 vdwparam+vdwioffset0+vdwjidx0C,
329 vdwparam+vdwioffset0+vdwjidx0D,
330 &c6_00,&c12_00);
332 /* REACTION-FIELD ELECTROSTATICS */
333 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_macc_ps(krf,rsq00,rinv00),crf));
334 felec = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
336 /* LENNARD-JONES DISPERSION/REPULSION */
338 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
339 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
340 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
341 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
342 _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
343 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
345 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
347 /* Update potential sum for this i atom from the interaction with this j atom. */
348 velec = _mm_and_ps(velec,cutoff_mask);
349 velec = _mm_andnot_ps(dummy_mask,velec);
350 velecsum = _mm_add_ps(velecsum,velec);
351 vvdw = _mm_and_ps(vvdw,cutoff_mask);
352 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
353 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
355 fscal = _mm_add_ps(felec,fvdw);
357 fscal = _mm_and_ps(fscal,cutoff_mask);
359 fscal = _mm_andnot_ps(dummy_mask,fscal);
361 /* Update vectorial force */
362 fix0 = _mm_macc_ps(dx00,fscal,fix0);
363 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
364 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
366 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
367 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
368 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
369 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
370 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
371 _mm_mul_ps(dx00,fscal),
372 _mm_mul_ps(dy00,fscal),
373 _mm_mul_ps(dz00,fscal));
377 /* Inner loop uses 57 flops */
380 /* End of innermost loop */
382 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
383 f+i_coord_offset,fshift+i_shift_offset);
385 ggid = gid[iidx];
386 /* Update potential energies */
387 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
388 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
390 /* Increment number of inner iterations */
391 inneriter += j_index_end - j_index_start;
393 /* Outer loop uses 9 flops */
396 /* Increment number of outer iterations */
397 outeriter += nri;
399 /* Update outer/inner flops */
401 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*57);
404 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_F_avx_128_fma_single
405 * Electrostatics interaction: ReactionField
406 * VdW interaction: LennardJones
407 * Geometry: Particle-Particle
408 * Calculate force/pot: Force
410 void
411 nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_F_avx_128_fma_single
412 (t_nblist * gmx_restrict nlist,
413 rvec * gmx_restrict xx,
414 rvec * gmx_restrict ff,
415 t_forcerec * gmx_restrict fr,
416 t_mdatoms * gmx_restrict mdatoms,
417 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
418 t_nrnb * gmx_restrict nrnb)
420 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
421 * just 0 for non-waters.
422 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
423 * jnr indices corresponding to data put in the four positions in the SIMD register.
425 int i_shift_offset,i_coord_offset,outeriter,inneriter;
426 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
427 int jnrA,jnrB,jnrC,jnrD;
428 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
429 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
430 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
431 real rcutoff_scalar;
432 real *shiftvec,*fshift,*x,*f;
433 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
434 real scratch[4*DIM];
435 __m128 fscal,rcutoff,rcutoff2,jidxall;
436 int vdwioffset0;
437 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
438 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
439 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
440 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
441 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
442 real *charge;
443 int nvdwtype;
444 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
445 int *vdwtype;
446 real *vdwparam;
447 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
448 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
449 __m128 dummy_mask,cutoff_mask;
450 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
451 __m128 one = _mm_set1_ps(1.0);
452 __m128 two = _mm_set1_ps(2.0);
453 x = xx[0];
454 f = ff[0];
456 nri = nlist->nri;
457 iinr = nlist->iinr;
458 jindex = nlist->jindex;
459 jjnr = nlist->jjnr;
460 shiftidx = nlist->shift;
461 gid = nlist->gid;
462 shiftvec = fr->shift_vec[0];
463 fshift = fr->fshift[0];
464 facel = _mm_set1_ps(fr->epsfac);
465 charge = mdatoms->chargeA;
466 krf = _mm_set1_ps(fr->ic->k_rf);
467 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
468 crf = _mm_set1_ps(fr->ic->c_rf);
469 nvdwtype = fr->ntype;
470 vdwparam = fr->nbfp;
471 vdwtype = mdatoms->typeA;
473 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
474 rcutoff_scalar = fr->rcoulomb;
475 rcutoff = _mm_set1_ps(rcutoff_scalar);
476 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
478 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
479 rvdw = _mm_set1_ps(fr->rvdw);
481 /* Avoid stupid compiler warnings */
482 jnrA = jnrB = jnrC = jnrD = 0;
483 j_coord_offsetA = 0;
484 j_coord_offsetB = 0;
485 j_coord_offsetC = 0;
486 j_coord_offsetD = 0;
488 outeriter = 0;
489 inneriter = 0;
491 for(iidx=0;iidx<4*DIM;iidx++)
493 scratch[iidx] = 0.0;
496 /* Start outer loop over neighborlists */
497 for(iidx=0; iidx<nri; iidx++)
499 /* Load shift vector for this list */
500 i_shift_offset = DIM*shiftidx[iidx];
502 /* Load limits for loop over neighbors */
503 j_index_start = jindex[iidx];
504 j_index_end = jindex[iidx+1];
506 /* Get outer coordinate index */
507 inr = iinr[iidx];
508 i_coord_offset = DIM*inr;
510 /* Load i particle coords and add shift vector */
511 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
513 fix0 = _mm_setzero_ps();
514 fiy0 = _mm_setzero_ps();
515 fiz0 = _mm_setzero_ps();
517 /* Load parameters for i particles */
518 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
519 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
521 /* Start inner kernel loop */
522 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
525 /* Get j neighbor index, and coordinate index */
526 jnrA = jjnr[jidx];
527 jnrB = jjnr[jidx+1];
528 jnrC = jjnr[jidx+2];
529 jnrD = jjnr[jidx+3];
530 j_coord_offsetA = DIM*jnrA;
531 j_coord_offsetB = DIM*jnrB;
532 j_coord_offsetC = DIM*jnrC;
533 j_coord_offsetD = DIM*jnrD;
535 /* load j atom coordinates */
536 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
537 x+j_coord_offsetC,x+j_coord_offsetD,
538 &jx0,&jy0,&jz0);
540 /* Calculate displacement vector */
541 dx00 = _mm_sub_ps(ix0,jx0);
542 dy00 = _mm_sub_ps(iy0,jy0);
543 dz00 = _mm_sub_ps(iz0,jz0);
545 /* Calculate squared distance and things based on it */
546 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
548 rinv00 = gmx_mm_invsqrt_ps(rsq00);
550 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
552 /* Load parameters for j particles */
553 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
554 charge+jnrC+0,charge+jnrD+0);
555 vdwjidx0A = 2*vdwtype[jnrA+0];
556 vdwjidx0B = 2*vdwtype[jnrB+0];
557 vdwjidx0C = 2*vdwtype[jnrC+0];
558 vdwjidx0D = 2*vdwtype[jnrD+0];
560 /**************************
561 * CALCULATE INTERACTIONS *
562 **************************/
564 if (gmx_mm_any_lt(rsq00,rcutoff2))
567 /* Compute parameters for interactions between i and j atoms */
568 qq00 = _mm_mul_ps(iq0,jq0);
569 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
570 vdwparam+vdwioffset0+vdwjidx0B,
571 vdwparam+vdwioffset0+vdwjidx0C,
572 vdwparam+vdwioffset0+vdwjidx0D,
573 &c6_00,&c12_00);
575 /* REACTION-FIELD ELECTROSTATICS */
576 felec = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
578 /* LENNARD-JONES DISPERSION/REPULSION */
580 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
581 fvdw = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
583 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
585 fscal = _mm_add_ps(felec,fvdw);
587 fscal = _mm_and_ps(fscal,cutoff_mask);
589 /* Update vectorial force */
590 fix0 = _mm_macc_ps(dx00,fscal,fix0);
591 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
592 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
594 fjptrA = f+j_coord_offsetA;
595 fjptrB = f+j_coord_offsetB;
596 fjptrC = f+j_coord_offsetC;
597 fjptrD = f+j_coord_offsetD;
598 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
599 _mm_mul_ps(dx00,fscal),
600 _mm_mul_ps(dy00,fscal),
601 _mm_mul_ps(dz00,fscal));
605 /* Inner loop uses 40 flops */
608 if(jidx<j_index_end)
611 /* Get j neighbor index, and coordinate index */
612 jnrlistA = jjnr[jidx];
613 jnrlistB = jjnr[jidx+1];
614 jnrlistC = jjnr[jidx+2];
615 jnrlistD = jjnr[jidx+3];
616 /* Sign of each element will be negative for non-real atoms.
617 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
618 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
620 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
621 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
622 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
623 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
624 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
625 j_coord_offsetA = DIM*jnrA;
626 j_coord_offsetB = DIM*jnrB;
627 j_coord_offsetC = DIM*jnrC;
628 j_coord_offsetD = DIM*jnrD;
630 /* load j atom coordinates */
631 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
632 x+j_coord_offsetC,x+j_coord_offsetD,
633 &jx0,&jy0,&jz0);
635 /* Calculate displacement vector */
636 dx00 = _mm_sub_ps(ix0,jx0);
637 dy00 = _mm_sub_ps(iy0,jy0);
638 dz00 = _mm_sub_ps(iz0,jz0);
640 /* Calculate squared distance and things based on it */
641 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
643 rinv00 = gmx_mm_invsqrt_ps(rsq00);
645 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
647 /* Load parameters for j particles */
648 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
649 charge+jnrC+0,charge+jnrD+0);
650 vdwjidx0A = 2*vdwtype[jnrA+0];
651 vdwjidx0B = 2*vdwtype[jnrB+0];
652 vdwjidx0C = 2*vdwtype[jnrC+0];
653 vdwjidx0D = 2*vdwtype[jnrD+0];
655 /**************************
656 * CALCULATE INTERACTIONS *
657 **************************/
659 if (gmx_mm_any_lt(rsq00,rcutoff2))
662 /* Compute parameters for interactions between i and j atoms */
663 qq00 = _mm_mul_ps(iq0,jq0);
664 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
665 vdwparam+vdwioffset0+vdwjidx0B,
666 vdwparam+vdwioffset0+vdwjidx0C,
667 vdwparam+vdwioffset0+vdwjidx0D,
668 &c6_00,&c12_00);
670 /* REACTION-FIELD ELECTROSTATICS */
671 felec = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
673 /* LENNARD-JONES DISPERSION/REPULSION */
675 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
676 fvdw = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
678 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
680 fscal = _mm_add_ps(felec,fvdw);
682 fscal = _mm_and_ps(fscal,cutoff_mask);
684 fscal = _mm_andnot_ps(dummy_mask,fscal);
686 /* Update vectorial force */
687 fix0 = _mm_macc_ps(dx00,fscal,fix0);
688 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
689 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
691 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
692 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
693 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
694 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
695 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
696 _mm_mul_ps(dx00,fscal),
697 _mm_mul_ps(dy00,fscal),
698 _mm_mul_ps(dz00,fscal));
702 /* Inner loop uses 40 flops */
705 /* End of innermost loop */
707 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
708 f+i_coord_offset,fshift+i_shift_offset);
710 /* Increment number of inner iterations */
711 inneriter += j_index_end - j_index_start;
713 /* Outer loop uses 7 flops */
716 /* Increment number of outer iterations */
717 outeriter += nri;
719 /* Update outer/inner flops */
721 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*40);