Removed simple.h from nb_kernel_sse4_1_XX
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse4_1_double / nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_sse4_1_double.c
blob7119da1ee67d0faf51a4bd869207842807014e8b
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 sse4_1_double 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_sse4_1_double.h"
49 #include "kernelutil_x86_sse4_1_double.h"
52 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_sse4_1_double
53 * Electrostatics interaction: ReactionField
54 * VdW interaction: CubicSplineTable
55 * Geometry: Particle-Particle
56 * Calculate force/pot: PotentialAndForce
58 void
59 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_sse4_1_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
76 int j_coord_offsetA,j_coord_offsetB;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
78 real rcutoff_scalar;
79 real *shiftvec,*fshift,*x,*f;
80 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
81 int vdwioffset0;
82 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
83 int vdwjidx0A,vdwjidx0B;
84 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
85 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
86 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
87 real *charge;
88 int nvdwtype;
89 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
90 int *vdwtype;
91 real *vdwparam;
92 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
93 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
94 __m128i vfitab;
95 __m128i ifour = _mm_set1_epi32(4);
96 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
97 real *vftab;
98 __m128d dummy_mask,cutoff_mask;
99 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
100 __m128d one = _mm_set1_pd(1.0);
101 __m128d two = _mm_set1_pd(2.0);
102 x = xx[0];
103 f = ff[0];
105 nri = nlist->nri;
106 iinr = nlist->iinr;
107 jindex = nlist->jindex;
108 jjnr = nlist->jjnr;
109 shiftidx = nlist->shift;
110 gid = nlist->gid;
111 shiftvec = fr->shift_vec[0];
112 fshift = fr->fshift[0];
113 facel = _mm_set1_pd(fr->epsfac);
114 charge = mdatoms->chargeA;
115 krf = _mm_set1_pd(fr->ic->k_rf);
116 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
117 crf = _mm_set1_pd(fr->ic->c_rf);
118 nvdwtype = fr->ntype;
119 vdwparam = fr->nbfp;
120 vdwtype = mdatoms->typeA;
122 vftab = kernel_data->table_vdw->data;
123 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
125 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
126 rcutoff_scalar = fr->rcoulomb;
127 rcutoff = _mm_set1_pd(rcutoff_scalar);
128 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
130 /* Avoid stupid compiler warnings */
131 jnrA = jnrB = 0;
132 j_coord_offsetA = 0;
133 j_coord_offsetB = 0;
135 outeriter = 0;
136 inneriter = 0;
138 /* Start outer loop over neighborlists */
139 for(iidx=0; iidx<nri; iidx++)
141 /* Load shift vector for this list */
142 i_shift_offset = DIM*shiftidx[iidx];
144 /* Load limits for loop over neighbors */
145 j_index_start = jindex[iidx];
146 j_index_end = jindex[iidx+1];
148 /* Get outer coordinate index */
149 inr = iinr[iidx];
150 i_coord_offset = DIM*inr;
152 /* Load i particle coords and add shift vector */
153 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
155 fix0 = _mm_setzero_pd();
156 fiy0 = _mm_setzero_pd();
157 fiz0 = _mm_setzero_pd();
159 /* Load parameters for i particles */
160 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
161 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
163 /* Reset potential sums */
164 velecsum = _mm_setzero_pd();
165 vvdwsum = _mm_setzero_pd();
167 /* Start inner kernel loop */
168 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
171 /* Get j neighbor index, and coordinate index */
172 jnrA = jjnr[jidx];
173 jnrB = jjnr[jidx+1];
174 j_coord_offsetA = DIM*jnrA;
175 j_coord_offsetB = DIM*jnrB;
177 /* load j atom coordinates */
178 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
179 &jx0,&jy0,&jz0);
181 /* Calculate displacement vector */
182 dx00 = _mm_sub_pd(ix0,jx0);
183 dy00 = _mm_sub_pd(iy0,jy0);
184 dz00 = _mm_sub_pd(iz0,jz0);
186 /* Calculate squared distance and things based on it */
187 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
189 rinv00 = gmx_mm_invsqrt_pd(rsq00);
191 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
193 /* Load parameters for j particles */
194 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
195 vdwjidx0A = 2*vdwtype[jnrA+0];
196 vdwjidx0B = 2*vdwtype[jnrB+0];
198 /**************************
199 * CALCULATE INTERACTIONS *
200 **************************/
202 if (gmx_mm_any_lt(rsq00,rcutoff2))
205 r00 = _mm_mul_pd(rsq00,rinv00);
207 /* Compute parameters for interactions between i and j atoms */
208 qq00 = _mm_mul_pd(iq0,jq0);
209 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
210 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
212 /* Calculate table index by multiplying r with table scale and truncate to integer */
213 rt = _mm_mul_pd(r00,vftabscale);
214 vfitab = _mm_cvttpd_epi32(rt);
215 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
216 vfitab = _mm_slli_epi32(vfitab,3);
218 /* REACTION-FIELD ELECTROSTATICS */
219 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
220 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
222 /* CUBIC SPLINE TABLE DISPERSION */
223 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
224 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
225 GMX_MM_TRANSPOSE2_PD(Y,F);
226 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
227 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
228 GMX_MM_TRANSPOSE2_PD(G,H);
229 Heps = _mm_mul_pd(vfeps,H);
230 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
231 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
232 vvdw6 = _mm_mul_pd(c6_00,VV);
233 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
234 fvdw6 = _mm_mul_pd(c6_00,FF);
236 /* CUBIC SPLINE TABLE REPULSION */
237 vfitab = _mm_add_epi32(vfitab,ifour);
238 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
239 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
240 GMX_MM_TRANSPOSE2_PD(Y,F);
241 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
242 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
243 GMX_MM_TRANSPOSE2_PD(G,H);
244 Heps = _mm_mul_pd(vfeps,H);
245 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
246 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
247 vvdw12 = _mm_mul_pd(c12_00,VV);
248 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
249 fvdw12 = _mm_mul_pd(c12_00,FF);
250 vvdw = _mm_add_pd(vvdw12,vvdw6);
251 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
253 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
255 /* Update potential sum for this i atom from the interaction with this j atom. */
256 velec = _mm_and_pd(velec,cutoff_mask);
257 velecsum = _mm_add_pd(velecsum,velec);
258 vvdw = _mm_and_pd(vvdw,cutoff_mask);
259 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
261 fscal = _mm_add_pd(felec,fvdw);
263 fscal = _mm_and_pd(fscal,cutoff_mask);
265 /* Calculate temporary vectorial force */
266 tx = _mm_mul_pd(fscal,dx00);
267 ty = _mm_mul_pd(fscal,dy00);
268 tz = _mm_mul_pd(fscal,dz00);
270 /* Update vectorial force */
271 fix0 = _mm_add_pd(fix0,tx);
272 fiy0 = _mm_add_pd(fiy0,ty);
273 fiz0 = _mm_add_pd(fiz0,tz);
275 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
279 /* Inner loop uses 72 flops */
282 if(jidx<j_index_end)
285 jnrA = jjnr[jidx];
286 j_coord_offsetA = DIM*jnrA;
288 /* load j atom coordinates */
289 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
290 &jx0,&jy0,&jz0);
292 /* Calculate displacement vector */
293 dx00 = _mm_sub_pd(ix0,jx0);
294 dy00 = _mm_sub_pd(iy0,jy0);
295 dz00 = _mm_sub_pd(iz0,jz0);
297 /* Calculate squared distance and things based on it */
298 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
300 rinv00 = gmx_mm_invsqrt_pd(rsq00);
302 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
304 /* Load parameters for j particles */
305 jq0 = _mm_load_sd(charge+jnrA+0);
306 vdwjidx0A = 2*vdwtype[jnrA+0];
308 /**************************
309 * CALCULATE INTERACTIONS *
310 **************************/
312 if (gmx_mm_any_lt(rsq00,rcutoff2))
315 r00 = _mm_mul_pd(rsq00,rinv00);
317 /* Compute parameters for interactions between i and j atoms */
318 qq00 = _mm_mul_pd(iq0,jq0);
319 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
321 /* Calculate table index by multiplying r with table scale and truncate to integer */
322 rt = _mm_mul_pd(r00,vftabscale);
323 vfitab = _mm_cvttpd_epi32(rt);
324 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
325 vfitab = _mm_slli_epi32(vfitab,3);
327 /* REACTION-FIELD ELECTROSTATICS */
328 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
329 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
331 /* CUBIC SPLINE TABLE DISPERSION */
332 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
333 F = _mm_setzero_pd();
334 GMX_MM_TRANSPOSE2_PD(Y,F);
335 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
336 H = _mm_setzero_pd();
337 GMX_MM_TRANSPOSE2_PD(G,H);
338 Heps = _mm_mul_pd(vfeps,H);
339 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
340 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
341 vvdw6 = _mm_mul_pd(c6_00,VV);
342 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
343 fvdw6 = _mm_mul_pd(c6_00,FF);
345 /* CUBIC SPLINE TABLE REPULSION */
346 vfitab = _mm_add_epi32(vfitab,ifour);
347 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
348 F = _mm_setzero_pd();
349 GMX_MM_TRANSPOSE2_PD(Y,F);
350 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
351 H = _mm_setzero_pd();
352 GMX_MM_TRANSPOSE2_PD(G,H);
353 Heps = _mm_mul_pd(vfeps,H);
354 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
355 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
356 vvdw12 = _mm_mul_pd(c12_00,VV);
357 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
358 fvdw12 = _mm_mul_pd(c12_00,FF);
359 vvdw = _mm_add_pd(vvdw12,vvdw6);
360 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
362 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
364 /* Update potential sum for this i atom from the interaction with this j atom. */
365 velec = _mm_and_pd(velec,cutoff_mask);
366 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
367 velecsum = _mm_add_pd(velecsum,velec);
368 vvdw = _mm_and_pd(vvdw,cutoff_mask);
369 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
370 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
372 fscal = _mm_add_pd(felec,fvdw);
374 fscal = _mm_and_pd(fscal,cutoff_mask);
376 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
378 /* Calculate temporary vectorial force */
379 tx = _mm_mul_pd(fscal,dx00);
380 ty = _mm_mul_pd(fscal,dy00);
381 tz = _mm_mul_pd(fscal,dz00);
383 /* Update vectorial force */
384 fix0 = _mm_add_pd(fix0,tx);
385 fiy0 = _mm_add_pd(fiy0,ty);
386 fiz0 = _mm_add_pd(fiz0,tz);
388 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
392 /* Inner loop uses 72 flops */
395 /* End of innermost loop */
397 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
398 f+i_coord_offset,fshift+i_shift_offset);
400 ggid = gid[iidx];
401 /* Update potential energies */
402 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
403 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
405 /* Increment number of inner iterations */
406 inneriter += j_index_end - j_index_start;
408 /* Outer loop uses 9 flops */
411 /* Increment number of outer iterations */
412 outeriter += nri;
414 /* Update outer/inner flops */
416 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*72);
419 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_sse4_1_double
420 * Electrostatics interaction: ReactionField
421 * VdW interaction: CubicSplineTable
422 * Geometry: Particle-Particle
423 * Calculate force/pot: Force
425 void
426 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_sse4_1_double
427 (t_nblist * gmx_restrict nlist,
428 rvec * gmx_restrict xx,
429 rvec * gmx_restrict ff,
430 t_forcerec * gmx_restrict fr,
431 t_mdatoms * gmx_restrict mdatoms,
432 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
433 t_nrnb * gmx_restrict nrnb)
435 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
436 * just 0 for non-waters.
437 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
438 * jnr indices corresponding to data put in the four positions in the SIMD register.
440 int i_shift_offset,i_coord_offset,outeriter,inneriter;
441 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
442 int jnrA,jnrB;
443 int j_coord_offsetA,j_coord_offsetB;
444 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
445 real rcutoff_scalar;
446 real *shiftvec,*fshift,*x,*f;
447 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
448 int vdwioffset0;
449 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
450 int vdwjidx0A,vdwjidx0B;
451 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
452 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
453 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
454 real *charge;
455 int nvdwtype;
456 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
457 int *vdwtype;
458 real *vdwparam;
459 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
460 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
461 __m128i vfitab;
462 __m128i ifour = _mm_set1_epi32(4);
463 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
464 real *vftab;
465 __m128d dummy_mask,cutoff_mask;
466 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
467 __m128d one = _mm_set1_pd(1.0);
468 __m128d two = _mm_set1_pd(2.0);
469 x = xx[0];
470 f = ff[0];
472 nri = nlist->nri;
473 iinr = nlist->iinr;
474 jindex = nlist->jindex;
475 jjnr = nlist->jjnr;
476 shiftidx = nlist->shift;
477 gid = nlist->gid;
478 shiftvec = fr->shift_vec[0];
479 fshift = fr->fshift[0];
480 facel = _mm_set1_pd(fr->epsfac);
481 charge = mdatoms->chargeA;
482 krf = _mm_set1_pd(fr->ic->k_rf);
483 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
484 crf = _mm_set1_pd(fr->ic->c_rf);
485 nvdwtype = fr->ntype;
486 vdwparam = fr->nbfp;
487 vdwtype = mdatoms->typeA;
489 vftab = kernel_data->table_vdw->data;
490 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
492 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
493 rcutoff_scalar = fr->rcoulomb;
494 rcutoff = _mm_set1_pd(rcutoff_scalar);
495 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
497 /* Avoid stupid compiler warnings */
498 jnrA = jnrB = 0;
499 j_coord_offsetA = 0;
500 j_coord_offsetB = 0;
502 outeriter = 0;
503 inneriter = 0;
505 /* Start outer loop over neighborlists */
506 for(iidx=0; iidx<nri; iidx++)
508 /* Load shift vector for this list */
509 i_shift_offset = DIM*shiftidx[iidx];
511 /* Load limits for loop over neighbors */
512 j_index_start = jindex[iidx];
513 j_index_end = jindex[iidx+1];
515 /* Get outer coordinate index */
516 inr = iinr[iidx];
517 i_coord_offset = DIM*inr;
519 /* Load i particle coords and add shift vector */
520 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
522 fix0 = _mm_setzero_pd();
523 fiy0 = _mm_setzero_pd();
524 fiz0 = _mm_setzero_pd();
526 /* Load parameters for i particles */
527 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
528 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
530 /* Start inner kernel loop */
531 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
534 /* Get j neighbor index, and coordinate index */
535 jnrA = jjnr[jidx];
536 jnrB = jjnr[jidx+1];
537 j_coord_offsetA = DIM*jnrA;
538 j_coord_offsetB = DIM*jnrB;
540 /* load j atom coordinates */
541 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
542 &jx0,&jy0,&jz0);
544 /* Calculate displacement vector */
545 dx00 = _mm_sub_pd(ix0,jx0);
546 dy00 = _mm_sub_pd(iy0,jy0);
547 dz00 = _mm_sub_pd(iz0,jz0);
549 /* Calculate squared distance and things based on it */
550 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
552 rinv00 = gmx_mm_invsqrt_pd(rsq00);
554 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
556 /* Load parameters for j particles */
557 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
558 vdwjidx0A = 2*vdwtype[jnrA+0];
559 vdwjidx0B = 2*vdwtype[jnrB+0];
561 /**************************
562 * CALCULATE INTERACTIONS *
563 **************************/
565 if (gmx_mm_any_lt(rsq00,rcutoff2))
568 r00 = _mm_mul_pd(rsq00,rinv00);
570 /* Compute parameters for interactions between i and j atoms */
571 qq00 = _mm_mul_pd(iq0,jq0);
572 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
573 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
575 /* Calculate table index by multiplying r with table scale and truncate to integer */
576 rt = _mm_mul_pd(r00,vftabscale);
577 vfitab = _mm_cvttpd_epi32(rt);
578 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
579 vfitab = _mm_slli_epi32(vfitab,3);
581 /* REACTION-FIELD ELECTROSTATICS */
582 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
584 /* CUBIC SPLINE TABLE DISPERSION */
585 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
586 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
587 GMX_MM_TRANSPOSE2_PD(Y,F);
588 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
589 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
590 GMX_MM_TRANSPOSE2_PD(G,H);
591 Heps = _mm_mul_pd(vfeps,H);
592 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
593 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
594 fvdw6 = _mm_mul_pd(c6_00,FF);
596 /* CUBIC SPLINE TABLE REPULSION */
597 vfitab = _mm_add_epi32(vfitab,ifour);
598 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
599 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
600 GMX_MM_TRANSPOSE2_PD(Y,F);
601 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
602 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
603 GMX_MM_TRANSPOSE2_PD(G,H);
604 Heps = _mm_mul_pd(vfeps,H);
605 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
606 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
607 fvdw12 = _mm_mul_pd(c12_00,FF);
608 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
610 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
612 fscal = _mm_add_pd(felec,fvdw);
614 fscal = _mm_and_pd(fscal,cutoff_mask);
616 /* Calculate temporary vectorial force */
617 tx = _mm_mul_pd(fscal,dx00);
618 ty = _mm_mul_pd(fscal,dy00);
619 tz = _mm_mul_pd(fscal,dz00);
621 /* Update vectorial force */
622 fix0 = _mm_add_pd(fix0,tx);
623 fiy0 = _mm_add_pd(fiy0,ty);
624 fiz0 = _mm_add_pd(fiz0,tz);
626 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
630 /* Inner loop uses 57 flops */
633 if(jidx<j_index_end)
636 jnrA = jjnr[jidx];
637 j_coord_offsetA = DIM*jnrA;
639 /* load j atom coordinates */
640 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
641 &jx0,&jy0,&jz0);
643 /* Calculate displacement vector */
644 dx00 = _mm_sub_pd(ix0,jx0);
645 dy00 = _mm_sub_pd(iy0,jy0);
646 dz00 = _mm_sub_pd(iz0,jz0);
648 /* Calculate squared distance and things based on it */
649 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
651 rinv00 = gmx_mm_invsqrt_pd(rsq00);
653 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
655 /* Load parameters for j particles */
656 jq0 = _mm_load_sd(charge+jnrA+0);
657 vdwjidx0A = 2*vdwtype[jnrA+0];
659 /**************************
660 * CALCULATE INTERACTIONS *
661 **************************/
663 if (gmx_mm_any_lt(rsq00,rcutoff2))
666 r00 = _mm_mul_pd(rsq00,rinv00);
668 /* Compute parameters for interactions between i and j atoms */
669 qq00 = _mm_mul_pd(iq0,jq0);
670 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
672 /* Calculate table index by multiplying r with table scale and truncate to integer */
673 rt = _mm_mul_pd(r00,vftabscale);
674 vfitab = _mm_cvttpd_epi32(rt);
675 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
676 vfitab = _mm_slli_epi32(vfitab,3);
678 /* REACTION-FIELD ELECTROSTATICS */
679 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
681 /* CUBIC SPLINE TABLE DISPERSION */
682 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
683 F = _mm_setzero_pd();
684 GMX_MM_TRANSPOSE2_PD(Y,F);
685 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
686 H = _mm_setzero_pd();
687 GMX_MM_TRANSPOSE2_PD(G,H);
688 Heps = _mm_mul_pd(vfeps,H);
689 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
690 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
691 fvdw6 = _mm_mul_pd(c6_00,FF);
693 /* CUBIC SPLINE TABLE REPULSION */
694 vfitab = _mm_add_epi32(vfitab,ifour);
695 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
696 F = _mm_setzero_pd();
697 GMX_MM_TRANSPOSE2_PD(Y,F);
698 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
699 H = _mm_setzero_pd();
700 GMX_MM_TRANSPOSE2_PD(G,H);
701 Heps = _mm_mul_pd(vfeps,H);
702 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
703 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
704 fvdw12 = _mm_mul_pd(c12_00,FF);
705 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
707 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
709 fscal = _mm_add_pd(felec,fvdw);
711 fscal = _mm_and_pd(fscal,cutoff_mask);
713 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
715 /* Calculate temporary vectorial force */
716 tx = _mm_mul_pd(fscal,dx00);
717 ty = _mm_mul_pd(fscal,dy00);
718 tz = _mm_mul_pd(fscal,dz00);
720 /* Update vectorial force */
721 fix0 = _mm_add_pd(fix0,tx);
722 fiy0 = _mm_add_pd(fiy0,ty);
723 fiz0 = _mm_add_pd(fiz0,tz);
725 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
729 /* Inner loop uses 57 flops */
732 /* End of innermost loop */
734 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
735 f+i_coord_offset,fshift+i_shift_offset);
737 /* Increment number of inner iterations */
738 inneriter += j_index_end - j_index_start;
740 /* Outer loop uses 7 flops */
743 /* Increment number of outer iterations */
744 outeriter += nri;
746 /* Update outer/inner flops */
748 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*57);