2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS sse2_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gromacs/simd/math_x86_sse2_single.h"
50 #include "kernelutil_x86_sse2_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_VF_sse2_single
54 * Electrostatics interaction: GeneralizedBorn
55 * VdW interaction: LennardJones
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecGB_VdwLJ_GeomP1P1_VF_sse2_single
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
81 real *shiftvec,*fshift,*x,*f;
82 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
84 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
86 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
88 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
93 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
94 __m128 minushalf = _mm_set1_ps(-0.5);
95 real *invsqrta,*dvda,*gbtab;
97 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
100 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
101 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
103 __m128i ifour = _mm_set1_epi32(4);
104 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
106 __m128 dummy_mask,cutoff_mask;
107 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
108 __m128 one = _mm_set1_ps(1.0);
109 __m128 two = _mm_set1_ps(2.0);
115 jindex = nlist->jindex;
117 shiftidx = nlist->shift;
119 shiftvec = fr->shift_vec[0];
120 fshift = fr->fshift[0];
121 facel = _mm_set1_ps(fr->epsfac);
122 charge = mdatoms->chargeA;
123 nvdwtype = fr->ntype;
125 vdwtype = mdatoms->typeA;
127 invsqrta = fr->invsqrta;
129 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
130 gbtab = fr->gbtab.data;
131 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
133 /* Avoid stupid compiler warnings */
134 jnrA = jnrB = jnrC = jnrD = 0;
143 for(iidx=0;iidx<4*DIM;iidx++)
148 /* Start outer loop over neighborlists */
149 for(iidx=0; iidx<nri; iidx++)
151 /* Load shift vector for this list */
152 i_shift_offset = DIM*shiftidx[iidx];
154 /* Load limits for loop over neighbors */
155 j_index_start = jindex[iidx];
156 j_index_end = jindex[iidx+1];
158 /* Get outer coordinate index */
160 i_coord_offset = DIM*inr;
162 /* Load i particle coords and add shift vector */
163 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
165 fix0 = _mm_setzero_ps();
166 fiy0 = _mm_setzero_ps();
167 fiz0 = _mm_setzero_ps();
169 /* Load parameters for i particles */
170 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
171 isai0 = _mm_load1_ps(invsqrta+inr+0);
172 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
174 /* Reset potential sums */
175 velecsum = _mm_setzero_ps();
176 vgbsum = _mm_setzero_ps();
177 vvdwsum = _mm_setzero_ps();
178 dvdasum = _mm_setzero_ps();
180 /* Start inner kernel loop */
181 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
184 /* Get j neighbor index, and coordinate index */
189 j_coord_offsetA = DIM*jnrA;
190 j_coord_offsetB = DIM*jnrB;
191 j_coord_offsetC = DIM*jnrC;
192 j_coord_offsetD = DIM*jnrD;
194 /* load j atom coordinates */
195 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
196 x+j_coord_offsetC,x+j_coord_offsetD,
199 /* Calculate displacement vector */
200 dx00 = _mm_sub_ps(ix0,jx0);
201 dy00 = _mm_sub_ps(iy0,jy0);
202 dz00 = _mm_sub_ps(iz0,jz0);
204 /* Calculate squared distance and things based on it */
205 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
207 rinv00 = gmx_mm_invsqrt_ps(rsq00);
209 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
211 /* Load parameters for j particles */
212 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
213 charge+jnrC+0,charge+jnrD+0);
214 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
215 invsqrta+jnrC+0,invsqrta+jnrD+0);
216 vdwjidx0A = 2*vdwtype[jnrA+0];
217 vdwjidx0B = 2*vdwtype[jnrB+0];
218 vdwjidx0C = 2*vdwtype[jnrC+0];
219 vdwjidx0D = 2*vdwtype[jnrD+0];
221 /**************************
222 * CALCULATE INTERACTIONS *
223 **************************/
225 r00 = _mm_mul_ps(rsq00,rinv00);
227 /* Compute parameters for interactions between i and j atoms */
228 qq00 = _mm_mul_ps(iq0,jq0);
229 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
230 vdwparam+vdwioffset0+vdwjidx0B,
231 vdwparam+vdwioffset0+vdwjidx0C,
232 vdwparam+vdwioffset0+vdwjidx0D,
235 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
236 isaprod = _mm_mul_ps(isai0,isaj0);
237 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
238 gbscale = _mm_mul_ps(isaprod,gbtabscale);
240 /* Calculate generalized born table index - this is a separate table from the normal one,
241 * but we use the same procedure by multiplying r with scale and truncating to integer.
243 rt = _mm_mul_ps(r00,gbscale);
244 gbitab = _mm_cvttps_epi32(rt);
245 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
246 gbitab = _mm_slli_epi32(gbitab,2);
248 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
249 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
250 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
251 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
252 _MM_TRANSPOSE4_PS(Y,F,G,H);
253 Heps = _mm_mul_ps(gbeps,H);
254 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
255 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
256 vgb = _mm_mul_ps(gbqqfactor,VV);
258 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
259 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
260 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
261 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
266 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
267 velec = _mm_mul_ps(qq00,rinv00);
268 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
270 /* LENNARD-JONES DISPERSION/REPULSION */
272 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
273 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
274 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
275 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
276 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
278 /* Update potential sum for this i atom from the interaction with this j atom. */
279 velecsum = _mm_add_ps(velecsum,velec);
280 vgbsum = _mm_add_ps(vgbsum,vgb);
281 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
283 fscal = _mm_add_ps(felec,fvdw);
285 /* Calculate temporary vectorial force */
286 tx = _mm_mul_ps(fscal,dx00);
287 ty = _mm_mul_ps(fscal,dy00);
288 tz = _mm_mul_ps(fscal,dz00);
290 /* Update vectorial force */
291 fix0 = _mm_add_ps(fix0,tx);
292 fiy0 = _mm_add_ps(fiy0,ty);
293 fiz0 = _mm_add_ps(fiz0,tz);
295 fjptrA = f+j_coord_offsetA;
296 fjptrB = f+j_coord_offsetB;
297 fjptrC = f+j_coord_offsetC;
298 fjptrD = f+j_coord_offsetD;
299 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
301 /* Inner loop uses 71 flops */
307 /* Get j neighbor index, and coordinate index */
308 jnrlistA = jjnr[jidx];
309 jnrlistB = jjnr[jidx+1];
310 jnrlistC = jjnr[jidx+2];
311 jnrlistD = jjnr[jidx+3];
312 /* Sign of each element will be negative for non-real atoms.
313 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
314 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
316 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
317 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
318 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
319 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
320 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
321 j_coord_offsetA = DIM*jnrA;
322 j_coord_offsetB = DIM*jnrB;
323 j_coord_offsetC = DIM*jnrC;
324 j_coord_offsetD = DIM*jnrD;
326 /* load j atom coordinates */
327 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
328 x+j_coord_offsetC,x+j_coord_offsetD,
331 /* Calculate displacement vector */
332 dx00 = _mm_sub_ps(ix0,jx0);
333 dy00 = _mm_sub_ps(iy0,jy0);
334 dz00 = _mm_sub_ps(iz0,jz0);
336 /* Calculate squared distance and things based on it */
337 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
339 rinv00 = gmx_mm_invsqrt_ps(rsq00);
341 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
343 /* Load parameters for j particles */
344 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
345 charge+jnrC+0,charge+jnrD+0);
346 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
347 invsqrta+jnrC+0,invsqrta+jnrD+0);
348 vdwjidx0A = 2*vdwtype[jnrA+0];
349 vdwjidx0B = 2*vdwtype[jnrB+0];
350 vdwjidx0C = 2*vdwtype[jnrC+0];
351 vdwjidx0D = 2*vdwtype[jnrD+0];
353 /**************************
354 * CALCULATE INTERACTIONS *
355 **************************/
357 r00 = _mm_mul_ps(rsq00,rinv00);
358 r00 = _mm_andnot_ps(dummy_mask,r00);
360 /* Compute parameters for interactions between i and j atoms */
361 qq00 = _mm_mul_ps(iq0,jq0);
362 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
363 vdwparam+vdwioffset0+vdwjidx0B,
364 vdwparam+vdwioffset0+vdwjidx0C,
365 vdwparam+vdwioffset0+vdwjidx0D,
368 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
369 isaprod = _mm_mul_ps(isai0,isaj0);
370 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
371 gbscale = _mm_mul_ps(isaprod,gbtabscale);
373 /* Calculate generalized born table index - this is a separate table from the normal one,
374 * but we use the same procedure by multiplying r with scale and truncating to integer.
376 rt = _mm_mul_ps(r00,gbscale);
377 gbitab = _mm_cvttps_epi32(rt);
378 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
379 gbitab = _mm_slli_epi32(gbitab,2);
381 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
382 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
383 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
384 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
385 _MM_TRANSPOSE4_PS(Y,F,G,H);
386 Heps = _mm_mul_ps(gbeps,H);
387 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
388 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
389 vgb = _mm_mul_ps(gbqqfactor,VV);
391 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
392 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
393 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
394 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
395 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
396 /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
397 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
398 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
399 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
400 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
401 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
402 velec = _mm_mul_ps(qq00,rinv00);
403 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
405 /* LENNARD-JONES DISPERSION/REPULSION */
407 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
408 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
409 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
410 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
411 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
413 /* Update potential sum for this i atom from the interaction with this j atom. */
414 velec = _mm_andnot_ps(dummy_mask,velec);
415 velecsum = _mm_add_ps(velecsum,velec);
416 vgb = _mm_andnot_ps(dummy_mask,vgb);
417 vgbsum = _mm_add_ps(vgbsum,vgb);
418 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
419 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
421 fscal = _mm_add_ps(felec,fvdw);
423 fscal = _mm_andnot_ps(dummy_mask,fscal);
425 /* Calculate temporary vectorial force */
426 tx = _mm_mul_ps(fscal,dx00);
427 ty = _mm_mul_ps(fscal,dy00);
428 tz = _mm_mul_ps(fscal,dz00);
430 /* Update vectorial force */
431 fix0 = _mm_add_ps(fix0,tx);
432 fiy0 = _mm_add_ps(fiy0,ty);
433 fiz0 = _mm_add_ps(fiz0,tz);
435 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
436 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
437 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
438 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
439 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
441 /* Inner loop uses 72 flops */
444 /* End of innermost loop */
446 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
447 f+i_coord_offset,fshift+i_shift_offset);
450 /* Update potential energies */
451 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
452 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
453 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
454 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
455 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
457 /* Increment number of inner iterations */
458 inneriter += j_index_end - j_index_start;
460 /* Outer loop uses 10 flops */
463 /* Increment number of outer iterations */
466 /* Update outer/inner flops */
468 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*72);
471 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse2_single
472 * Electrostatics interaction: GeneralizedBorn
473 * VdW interaction: LennardJones
474 * Geometry: Particle-Particle
475 * Calculate force/pot: Force
478 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse2_single
479 (t_nblist * gmx_restrict nlist,
480 rvec * gmx_restrict xx,
481 rvec * gmx_restrict ff,
482 t_forcerec * gmx_restrict fr,
483 t_mdatoms * gmx_restrict mdatoms,
484 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
485 t_nrnb * gmx_restrict nrnb)
487 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
488 * just 0 for non-waters.
489 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
490 * jnr indices corresponding to data put in the four positions in the SIMD register.
492 int i_shift_offset,i_coord_offset,outeriter,inneriter;
493 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
494 int jnrA,jnrB,jnrC,jnrD;
495 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
496 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
497 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
499 real *shiftvec,*fshift,*x,*f;
500 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
502 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
504 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
505 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
506 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
507 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
508 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
511 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
512 __m128 minushalf = _mm_set1_ps(-0.5);
513 real *invsqrta,*dvda,*gbtab;
515 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
518 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
519 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
521 __m128i ifour = _mm_set1_epi32(4);
522 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
524 __m128 dummy_mask,cutoff_mask;
525 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
526 __m128 one = _mm_set1_ps(1.0);
527 __m128 two = _mm_set1_ps(2.0);
533 jindex = nlist->jindex;
535 shiftidx = nlist->shift;
537 shiftvec = fr->shift_vec[0];
538 fshift = fr->fshift[0];
539 facel = _mm_set1_ps(fr->epsfac);
540 charge = mdatoms->chargeA;
541 nvdwtype = fr->ntype;
543 vdwtype = mdatoms->typeA;
545 invsqrta = fr->invsqrta;
547 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
548 gbtab = fr->gbtab.data;
549 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
551 /* Avoid stupid compiler warnings */
552 jnrA = jnrB = jnrC = jnrD = 0;
561 for(iidx=0;iidx<4*DIM;iidx++)
566 /* Start outer loop over neighborlists */
567 for(iidx=0; iidx<nri; iidx++)
569 /* Load shift vector for this list */
570 i_shift_offset = DIM*shiftidx[iidx];
572 /* Load limits for loop over neighbors */
573 j_index_start = jindex[iidx];
574 j_index_end = jindex[iidx+1];
576 /* Get outer coordinate index */
578 i_coord_offset = DIM*inr;
580 /* Load i particle coords and add shift vector */
581 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
583 fix0 = _mm_setzero_ps();
584 fiy0 = _mm_setzero_ps();
585 fiz0 = _mm_setzero_ps();
587 /* Load parameters for i particles */
588 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
589 isai0 = _mm_load1_ps(invsqrta+inr+0);
590 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
592 dvdasum = _mm_setzero_ps();
594 /* Start inner kernel loop */
595 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
598 /* Get j neighbor index, and coordinate index */
603 j_coord_offsetA = DIM*jnrA;
604 j_coord_offsetB = DIM*jnrB;
605 j_coord_offsetC = DIM*jnrC;
606 j_coord_offsetD = DIM*jnrD;
608 /* load j atom coordinates */
609 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
610 x+j_coord_offsetC,x+j_coord_offsetD,
613 /* Calculate displacement vector */
614 dx00 = _mm_sub_ps(ix0,jx0);
615 dy00 = _mm_sub_ps(iy0,jy0);
616 dz00 = _mm_sub_ps(iz0,jz0);
618 /* Calculate squared distance and things based on it */
619 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
621 rinv00 = gmx_mm_invsqrt_ps(rsq00);
623 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
625 /* Load parameters for j particles */
626 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
627 charge+jnrC+0,charge+jnrD+0);
628 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
629 invsqrta+jnrC+0,invsqrta+jnrD+0);
630 vdwjidx0A = 2*vdwtype[jnrA+0];
631 vdwjidx0B = 2*vdwtype[jnrB+0];
632 vdwjidx0C = 2*vdwtype[jnrC+0];
633 vdwjidx0D = 2*vdwtype[jnrD+0];
635 /**************************
636 * CALCULATE INTERACTIONS *
637 **************************/
639 r00 = _mm_mul_ps(rsq00,rinv00);
641 /* Compute parameters for interactions between i and j atoms */
642 qq00 = _mm_mul_ps(iq0,jq0);
643 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
644 vdwparam+vdwioffset0+vdwjidx0B,
645 vdwparam+vdwioffset0+vdwjidx0C,
646 vdwparam+vdwioffset0+vdwjidx0D,
649 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
650 isaprod = _mm_mul_ps(isai0,isaj0);
651 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
652 gbscale = _mm_mul_ps(isaprod,gbtabscale);
654 /* Calculate generalized born table index - this is a separate table from the normal one,
655 * but we use the same procedure by multiplying r with scale and truncating to integer.
657 rt = _mm_mul_ps(r00,gbscale);
658 gbitab = _mm_cvttps_epi32(rt);
659 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
660 gbitab = _mm_slli_epi32(gbitab,2);
662 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
663 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
664 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
665 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
666 _MM_TRANSPOSE4_PS(Y,F,G,H);
667 Heps = _mm_mul_ps(gbeps,H);
668 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
669 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
670 vgb = _mm_mul_ps(gbqqfactor,VV);
672 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
673 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
674 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
675 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
680 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
681 velec = _mm_mul_ps(qq00,rinv00);
682 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
684 /* LENNARD-JONES DISPERSION/REPULSION */
686 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
687 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
689 fscal = _mm_add_ps(felec,fvdw);
691 /* Calculate temporary vectorial force */
692 tx = _mm_mul_ps(fscal,dx00);
693 ty = _mm_mul_ps(fscal,dy00);
694 tz = _mm_mul_ps(fscal,dz00);
696 /* Update vectorial force */
697 fix0 = _mm_add_ps(fix0,tx);
698 fiy0 = _mm_add_ps(fiy0,ty);
699 fiz0 = _mm_add_ps(fiz0,tz);
701 fjptrA = f+j_coord_offsetA;
702 fjptrB = f+j_coord_offsetB;
703 fjptrC = f+j_coord_offsetC;
704 fjptrD = f+j_coord_offsetD;
705 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
707 /* Inner loop uses 64 flops */
713 /* Get j neighbor index, and coordinate index */
714 jnrlistA = jjnr[jidx];
715 jnrlistB = jjnr[jidx+1];
716 jnrlistC = jjnr[jidx+2];
717 jnrlistD = jjnr[jidx+3];
718 /* Sign of each element will be negative for non-real atoms.
719 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
720 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
722 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
723 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
724 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
725 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
726 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
727 j_coord_offsetA = DIM*jnrA;
728 j_coord_offsetB = DIM*jnrB;
729 j_coord_offsetC = DIM*jnrC;
730 j_coord_offsetD = DIM*jnrD;
732 /* load j atom coordinates */
733 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
734 x+j_coord_offsetC,x+j_coord_offsetD,
737 /* Calculate displacement vector */
738 dx00 = _mm_sub_ps(ix0,jx0);
739 dy00 = _mm_sub_ps(iy0,jy0);
740 dz00 = _mm_sub_ps(iz0,jz0);
742 /* Calculate squared distance and things based on it */
743 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
745 rinv00 = gmx_mm_invsqrt_ps(rsq00);
747 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
749 /* Load parameters for j particles */
750 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
751 charge+jnrC+0,charge+jnrD+0);
752 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
753 invsqrta+jnrC+0,invsqrta+jnrD+0);
754 vdwjidx0A = 2*vdwtype[jnrA+0];
755 vdwjidx0B = 2*vdwtype[jnrB+0];
756 vdwjidx0C = 2*vdwtype[jnrC+0];
757 vdwjidx0D = 2*vdwtype[jnrD+0];
759 /**************************
760 * CALCULATE INTERACTIONS *
761 **************************/
763 r00 = _mm_mul_ps(rsq00,rinv00);
764 r00 = _mm_andnot_ps(dummy_mask,r00);
766 /* Compute parameters for interactions between i and j atoms */
767 qq00 = _mm_mul_ps(iq0,jq0);
768 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
769 vdwparam+vdwioffset0+vdwjidx0B,
770 vdwparam+vdwioffset0+vdwjidx0C,
771 vdwparam+vdwioffset0+vdwjidx0D,
774 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
775 isaprod = _mm_mul_ps(isai0,isaj0);
776 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
777 gbscale = _mm_mul_ps(isaprod,gbtabscale);
779 /* Calculate generalized born table index - this is a separate table from the normal one,
780 * but we use the same procedure by multiplying r with scale and truncating to integer.
782 rt = _mm_mul_ps(r00,gbscale);
783 gbitab = _mm_cvttps_epi32(rt);
784 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
785 gbitab = _mm_slli_epi32(gbitab,2);
787 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
788 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
789 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
790 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
791 _MM_TRANSPOSE4_PS(Y,F,G,H);
792 Heps = _mm_mul_ps(gbeps,H);
793 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
794 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
795 vgb = _mm_mul_ps(gbqqfactor,VV);
797 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
798 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
799 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
800 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
801 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
802 /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
803 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
804 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
805 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
806 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
807 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
808 velec = _mm_mul_ps(qq00,rinv00);
809 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
811 /* LENNARD-JONES DISPERSION/REPULSION */
813 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
814 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
816 fscal = _mm_add_ps(felec,fvdw);
818 fscal = _mm_andnot_ps(dummy_mask,fscal);
820 /* Calculate temporary vectorial force */
821 tx = _mm_mul_ps(fscal,dx00);
822 ty = _mm_mul_ps(fscal,dy00);
823 tz = _mm_mul_ps(fscal,dz00);
825 /* Update vectorial force */
826 fix0 = _mm_add_ps(fix0,tx);
827 fiy0 = _mm_add_ps(fiy0,ty);
828 fiz0 = _mm_add_ps(fiz0,tz);
830 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
831 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
832 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
833 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
834 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
836 /* Inner loop uses 65 flops */
839 /* End of innermost loop */
841 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
842 f+i_coord_offset,fshift+i_shift_offset);
844 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
845 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
847 /* Increment number of inner iterations */
848 inneriter += j_index_end - j_index_start;
850 /* Outer loop uses 7 flops */
853 /* Increment number of outer iterations */
856 /* Update outer/inner flops */
858 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*65);