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 sse4_1_single kernel generator.
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
47 #include "gromacs/simd/math_x86_sse4_1_single.h"
48 #include "kernelutil_x86_sse4_1_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_VF_sse4_1_single
52 * Electrostatics interaction: GeneralizedBorn
53 * VdW interaction: LennardJones
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecGB_VdwLJ_GeomP1P1_VF_sse4_1_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 SSE, 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;
79 real *shiftvec,*fshift,*x,*f;
80 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
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;
91 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
92 __m128 minushalf = _mm_set1_ps(-0.5);
93 real *invsqrta,*dvda,*gbtab;
95 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
98 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
99 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
101 __m128i ifour = _mm_set1_epi32(4);
102 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
104 __m128 dummy_mask,cutoff_mask;
105 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
106 __m128 one = _mm_set1_ps(1.0);
107 __m128 two = _mm_set1_ps(2.0);
113 jindex = nlist->jindex;
115 shiftidx = nlist->shift;
117 shiftvec = fr->shift_vec[0];
118 fshift = fr->fshift[0];
119 facel = _mm_set1_ps(fr->epsfac);
120 charge = mdatoms->chargeA;
121 nvdwtype = fr->ntype;
123 vdwtype = mdatoms->typeA;
125 invsqrta = fr->invsqrta;
127 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
128 gbtab = fr->gbtab.data;
129 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
131 /* Avoid stupid compiler warnings */
132 jnrA = jnrB = jnrC = jnrD = 0;
141 for(iidx=0;iidx<4*DIM;iidx++)
146 /* Start outer loop over neighborlists */
147 for(iidx=0; iidx<nri; iidx++)
149 /* Load shift vector for this list */
150 i_shift_offset = DIM*shiftidx[iidx];
152 /* Load limits for loop over neighbors */
153 j_index_start = jindex[iidx];
154 j_index_end = jindex[iidx+1];
156 /* Get outer coordinate index */
158 i_coord_offset = DIM*inr;
160 /* Load i particle coords and add shift vector */
161 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
163 fix0 = _mm_setzero_ps();
164 fiy0 = _mm_setzero_ps();
165 fiz0 = _mm_setzero_ps();
167 /* Load parameters for i particles */
168 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
169 isai0 = _mm_load1_ps(invsqrta+inr+0);
170 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
172 /* Reset potential sums */
173 velecsum = _mm_setzero_ps();
174 vgbsum = _mm_setzero_ps();
175 vvdwsum = _mm_setzero_ps();
176 dvdasum = _mm_setzero_ps();
178 /* Start inner kernel loop */
179 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
182 /* Get j neighbor index, and coordinate index */
187 j_coord_offsetA = DIM*jnrA;
188 j_coord_offsetB = DIM*jnrB;
189 j_coord_offsetC = DIM*jnrC;
190 j_coord_offsetD = DIM*jnrD;
192 /* load j atom coordinates */
193 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
194 x+j_coord_offsetC,x+j_coord_offsetD,
197 /* Calculate displacement vector */
198 dx00 = _mm_sub_ps(ix0,jx0);
199 dy00 = _mm_sub_ps(iy0,jy0);
200 dz00 = _mm_sub_ps(iz0,jz0);
202 /* Calculate squared distance and things based on it */
203 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
205 rinv00 = gmx_mm_invsqrt_ps(rsq00);
207 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
209 /* Load parameters for j particles */
210 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
211 charge+jnrC+0,charge+jnrD+0);
212 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
213 invsqrta+jnrC+0,invsqrta+jnrD+0);
214 vdwjidx0A = 2*vdwtype[jnrA+0];
215 vdwjidx0B = 2*vdwtype[jnrB+0];
216 vdwjidx0C = 2*vdwtype[jnrC+0];
217 vdwjidx0D = 2*vdwtype[jnrD+0];
219 /**************************
220 * CALCULATE INTERACTIONS *
221 **************************/
223 r00 = _mm_mul_ps(rsq00,rinv00);
225 /* Compute parameters for interactions between i and j atoms */
226 qq00 = _mm_mul_ps(iq0,jq0);
227 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
228 vdwparam+vdwioffset0+vdwjidx0B,
229 vdwparam+vdwioffset0+vdwjidx0C,
230 vdwparam+vdwioffset0+vdwjidx0D,
233 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
234 isaprod = _mm_mul_ps(isai0,isaj0);
235 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
236 gbscale = _mm_mul_ps(isaprod,gbtabscale);
238 /* Calculate generalized born table index - this is a separate table from the normal one,
239 * but we use the same procedure by multiplying r with scale and truncating to integer.
241 rt = _mm_mul_ps(r00,gbscale);
242 gbitab = _mm_cvttps_epi32(rt);
243 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
244 gbitab = _mm_slli_epi32(gbitab,2);
245 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
246 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
247 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
248 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
249 _MM_TRANSPOSE4_PS(Y,F,G,H);
250 Heps = _mm_mul_ps(gbeps,H);
251 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
252 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
253 vgb = _mm_mul_ps(gbqqfactor,VV);
255 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
256 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
257 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
258 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
263 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
264 velec = _mm_mul_ps(qq00,rinv00);
265 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
267 /* LENNARD-JONES DISPERSION/REPULSION */
269 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
270 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
271 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
272 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
273 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
275 /* Update potential sum for this i atom from the interaction with this j atom. */
276 velecsum = _mm_add_ps(velecsum,velec);
277 vgbsum = _mm_add_ps(vgbsum,vgb);
278 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
280 fscal = _mm_add_ps(felec,fvdw);
282 /* Calculate temporary vectorial force */
283 tx = _mm_mul_ps(fscal,dx00);
284 ty = _mm_mul_ps(fscal,dy00);
285 tz = _mm_mul_ps(fscal,dz00);
287 /* Update vectorial force */
288 fix0 = _mm_add_ps(fix0,tx);
289 fiy0 = _mm_add_ps(fiy0,ty);
290 fiz0 = _mm_add_ps(fiz0,tz);
292 fjptrA = f+j_coord_offsetA;
293 fjptrB = f+j_coord_offsetB;
294 fjptrC = f+j_coord_offsetC;
295 fjptrD = f+j_coord_offsetD;
296 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
298 /* Inner loop uses 71 flops */
304 /* Get j neighbor index, and coordinate index */
305 jnrlistA = jjnr[jidx];
306 jnrlistB = jjnr[jidx+1];
307 jnrlistC = jjnr[jidx+2];
308 jnrlistD = jjnr[jidx+3];
309 /* Sign of each element will be negative for non-real atoms.
310 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
311 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
313 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
314 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
315 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
316 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
317 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
318 j_coord_offsetA = DIM*jnrA;
319 j_coord_offsetB = DIM*jnrB;
320 j_coord_offsetC = DIM*jnrC;
321 j_coord_offsetD = DIM*jnrD;
323 /* load j atom coordinates */
324 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
325 x+j_coord_offsetC,x+j_coord_offsetD,
328 /* Calculate displacement vector */
329 dx00 = _mm_sub_ps(ix0,jx0);
330 dy00 = _mm_sub_ps(iy0,jy0);
331 dz00 = _mm_sub_ps(iz0,jz0);
333 /* Calculate squared distance and things based on it */
334 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
336 rinv00 = gmx_mm_invsqrt_ps(rsq00);
338 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
340 /* Load parameters for j particles */
341 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
342 charge+jnrC+0,charge+jnrD+0);
343 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
344 invsqrta+jnrC+0,invsqrta+jnrD+0);
345 vdwjidx0A = 2*vdwtype[jnrA+0];
346 vdwjidx0B = 2*vdwtype[jnrB+0];
347 vdwjidx0C = 2*vdwtype[jnrC+0];
348 vdwjidx0D = 2*vdwtype[jnrD+0];
350 /**************************
351 * CALCULATE INTERACTIONS *
352 **************************/
354 r00 = _mm_mul_ps(rsq00,rinv00);
355 r00 = _mm_andnot_ps(dummy_mask,r00);
357 /* Compute parameters for interactions between i and j atoms */
358 qq00 = _mm_mul_ps(iq0,jq0);
359 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
360 vdwparam+vdwioffset0+vdwjidx0B,
361 vdwparam+vdwioffset0+vdwjidx0C,
362 vdwparam+vdwioffset0+vdwjidx0D,
365 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
366 isaprod = _mm_mul_ps(isai0,isaj0);
367 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
368 gbscale = _mm_mul_ps(isaprod,gbtabscale);
370 /* Calculate generalized born table index - this is a separate table from the normal one,
371 * but we use the same procedure by multiplying r with scale and truncating to integer.
373 rt = _mm_mul_ps(r00,gbscale);
374 gbitab = _mm_cvttps_epi32(rt);
375 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
376 gbitab = _mm_slli_epi32(gbitab,2);
377 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
378 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
379 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
380 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
381 _MM_TRANSPOSE4_PS(Y,F,G,H);
382 Heps = _mm_mul_ps(gbeps,H);
383 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
384 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
385 vgb = _mm_mul_ps(gbqqfactor,VV);
387 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
388 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
389 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
390 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
391 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
392 /* 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. */
393 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
394 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
395 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
396 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
397 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
398 velec = _mm_mul_ps(qq00,rinv00);
399 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
401 /* LENNARD-JONES DISPERSION/REPULSION */
403 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
404 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
405 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
406 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
407 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
409 /* Update potential sum for this i atom from the interaction with this j atom. */
410 velec = _mm_andnot_ps(dummy_mask,velec);
411 velecsum = _mm_add_ps(velecsum,velec);
412 vgb = _mm_andnot_ps(dummy_mask,vgb);
413 vgbsum = _mm_add_ps(vgbsum,vgb);
414 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
415 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
417 fscal = _mm_add_ps(felec,fvdw);
419 fscal = _mm_andnot_ps(dummy_mask,fscal);
421 /* Calculate temporary vectorial force */
422 tx = _mm_mul_ps(fscal,dx00);
423 ty = _mm_mul_ps(fscal,dy00);
424 tz = _mm_mul_ps(fscal,dz00);
426 /* Update vectorial force */
427 fix0 = _mm_add_ps(fix0,tx);
428 fiy0 = _mm_add_ps(fiy0,ty);
429 fiz0 = _mm_add_ps(fiz0,tz);
431 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
432 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
433 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
434 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
435 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
437 /* Inner loop uses 72 flops */
440 /* End of innermost loop */
442 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
443 f+i_coord_offset,fshift+i_shift_offset);
446 /* Update potential energies */
447 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
448 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
449 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
450 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
451 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
453 /* Increment number of inner iterations */
454 inneriter += j_index_end - j_index_start;
456 /* Outer loop uses 10 flops */
459 /* Increment number of outer iterations */
462 /* Update outer/inner flops */
464 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*72);
467 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse4_1_single
468 * Electrostatics interaction: GeneralizedBorn
469 * VdW interaction: LennardJones
470 * Geometry: Particle-Particle
471 * Calculate force/pot: Force
474 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse4_1_single
475 (t_nblist * gmx_restrict nlist,
476 rvec * gmx_restrict xx,
477 rvec * gmx_restrict ff,
478 t_forcerec * gmx_restrict fr,
479 t_mdatoms * gmx_restrict mdatoms,
480 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
481 t_nrnb * gmx_restrict nrnb)
483 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
484 * just 0 for non-waters.
485 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
486 * jnr indices corresponding to data put in the four positions in the SIMD register.
488 int i_shift_offset,i_coord_offset,outeriter,inneriter;
489 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
490 int jnrA,jnrB,jnrC,jnrD;
491 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
492 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
493 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
495 real *shiftvec,*fshift,*x,*f;
496 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
498 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
500 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
501 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
502 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
503 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
504 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
507 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
508 __m128 minushalf = _mm_set1_ps(-0.5);
509 real *invsqrta,*dvda,*gbtab;
511 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
514 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
515 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
517 __m128i ifour = _mm_set1_epi32(4);
518 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
520 __m128 dummy_mask,cutoff_mask;
521 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
522 __m128 one = _mm_set1_ps(1.0);
523 __m128 two = _mm_set1_ps(2.0);
529 jindex = nlist->jindex;
531 shiftidx = nlist->shift;
533 shiftvec = fr->shift_vec[0];
534 fshift = fr->fshift[0];
535 facel = _mm_set1_ps(fr->epsfac);
536 charge = mdatoms->chargeA;
537 nvdwtype = fr->ntype;
539 vdwtype = mdatoms->typeA;
541 invsqrta = fr->invsqrta;
543 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
544 gbtab = fr->gbtab.data;
545 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
547 /* Avoid stupid compiler warnings */
548 jnrA = jnrB = jnrC = jnrD = 0;
557 for(iidx=0;iidx<4*DIM;iidx++)
562 /* Start outer loop over neighborlists */
563 for(iidx=0; iidx<nri; iidx++)
565 /* Load shift vector for this list */
566 i_shift_offset = DIM*shiftidx[iidx];
568 /* Load limits for loop over neighbors */
569 j_index_start = jindex[iidx];
570 j_index_end = jindex[iidx+1];
572 /* Get outer coordinate index */
574 i_coord_offset = DIM*inr;
576 /* Load i particle coords and add shift vector */
577 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
579 fix0 = _mm_setzero_ps();
580 fiy0 = _mm_setzero_ps();
581 fiz0 = _mm_setzero_ps();
583 /* Load parameters for i particles */
584 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
585 isai0 = _mm_load1_ps(invsqrta+inr+0);
586 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
588 dvdasum = _mm_setzero_ps();
590 /* Start inner kernel loop */
591 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
594 /* Get j neighbor index, and coordinate index */
599 j_coord_offsetA = DIM*jnrA;
600 j_coord_offsetB = DIM*jnrB;
601 j_coord_offsetC = DIM*jnrC;
602 j_coord_offsetD = DIM*jnrD;
604 /* load j atom coordinates */
605 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
606 x+j_coord_offsetC,x+j_coord_offsetD,
609 /* Calculate displacement vector */
610 dx00 = _mm_sub_ps(ix0,jx0);
611 dy00 = _mm_sub_ps(iy0,jy0);
612 dz00 = _mm_sub_ps(iz0,jz0);
614 /* Calculate squared distance and things based on it */
615 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
617 rinv00 = gmx_mm_invsqrt_ps(rsq00);
619 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
621 /* Load parameters for j particles */
622 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
623 charge+jnrC+0,charge+jnrD+0);
624 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
625 invsqrta+jnrC+0,invsqrta+jnrD+0);
626 vdwjidx0A = 2*vdwtype[jnrA+0];
627 vdwjidx0B = 2*vdwtype[jnrB+0];
628 vdwjidx0C = 2*vdwtype[jnrC+0];
629 vdwjidx0D = 2*vdwtype[jnrD+0];
631 /**************************
632 * CALCULATE INTERACTIONS *
633 **************************/
635 r00 = _mm_mul_ps(rsq00,rinv00);
637 /* Compute parameters for interactions between i and j atoms */
638 qq00 = _mm_mul_ps(iq0,jq0);
639 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
640 vdwparam+vdwioffset0+vdwjidx0B,
641 vdwparam+vdwioffset0+vdwjidx0C,
642 vdwparam+vdwioffset0+vdwjidx0D,
645 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
646 isaprod = _mm_mul_ps(isai0,isaj0);
647 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
648 gbscale = _mm_mul_ps(isaprod,gbtabscale);
650 /* Calculate generalized born table index - this is a separate table from the normal one,
651 * but we use the same procedure by multiplying r with scale and truncating to integer.
653 rt = _mm_mul_ps(r00,gbscale);
654 gbitab = _mm_cvttps_epi32(rt);
655 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
656 gbitab = _mm_slli_epi32(gbitab,2);
657 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
658 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
659 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
660 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
661 _MM_TRANSPOSE4_PS(Y,F,G,H);
662 Heps = _mm_mul_ps(gbeps,H);
663 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
664 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
665 vgb = _mm_mul_ps(gbqqfactor,VV);
667 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
668 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
669 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
670 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
675 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
676 velec = _mm_mul_ps(qq00,rinv00);
677 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
679 /* LENNARD-JONES DISPERSION/REPULSION */
681 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
682 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
684 fscal = _mm_add_ps(felec,fvdw);
686 /* Calculate temporary vectorial force */
687 tx = _mm_mul_ps(fscal,dx00);
688 ty = _mm_mul_ps(fscal,dy00);
689 tz = _mm_mul_ps(fscal,dz00);
691 /* Update vectorial force */
692 fix0 = _mm_add_ps(fix0,tx);
693 fiy0 = _mm_add_ps(fiy0,ty);
694 fiz0 = _mm_add_ps(fiz0,tz);
696 fjptrA = f+j_coord_offsetA;
697 fjptrB = f+j_coord_offsetB;
698 fjptrC = f+j_coord_offsetC;
699 fjptrD = f+j_coord_offsetD;
700 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
702 /* Inner loop uses 64 flops */
708 /* Get j neighbor index, and coordinate index */
709 jnrlistA = jjnr[jidx];
710 jnrlistB = jjnr[jidx+1];
711 jnrlistC = jjnr[jidx+2];
712 jnrlistD = jjnr[jidx+3];
713 /* Sign of each element will be negative for non-real atoms.
714 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
715 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
717 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
718 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
719 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
720 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
721 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
722 j_coord_offsetA = DIM*jnrA;
723 j_coord_offsetB = DIM*jnrB;
724 j_coord_offsetC = DIM*jnrC;
725 j_coord_offsetD = DIM*jnrD;
727 /* load j atom coordinates */
728 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
729 x+j_coord_offsetC,x+j_coord_offsetD,
732 /* Calculate displacement vector */
733 dx00 = _mm_sub_ps(ix0,jx0);
734 dy00 = _mm_sub_ps(iy0,jy0);
735 dz00 = _mm_sub_ps(iz0,jz0);
737 /* Calculate squared distance and things based on it */
738 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
740 rinv00 = gmx_mm_invsqrt_ps(rsq00);
742 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
744 /* Load parameters for j particles */
745 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
746 charge+jnrC+0,charge+jnrD+0);
747 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
748 invsqrta+jnrC+0,invsqrta+jnrD+0);
749 vdwjidx0A = 2*vdwtype[jnrA+0];
750 vdwjidx0B = 2*vdwtype[jnrB+0];
751 vdwjidx0C = 2*vdwtype[jnrC+0];
752 vdwjidx0D = 2*vdwtype[jnrD+0];
754 /**************************
755 * CALCULATE INTERACTIONS *
756 **************************/
758 r00 = _mm_mul_ps(rsq00,rinv00);
759 r00 = _mm_andnot_ps(dummy_mask,r00);
761 /* Compute parameters for interactions between i and j atoms */
762 qq00 = _mm_mul_ps(iq0,jq0);
763 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
764 vdwparam+vdwioffset0+vdwjidx0B,
765 vdwparam+vdwioffset0+vdwjidx0C,
766 vdwparam+vdwioffset0+vdwjidx0D,
769 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
770 isaprod = _mm_mul_ps(isai0,isaj0);
771 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
772 gbscale = _mm_mul_ps(isaprod,gbtabscale);
774 /* Calculate generalized born table index - this is a separate table from the normal one,
775 * but we use the same procedure by multiplying r with scale and truncating to integer.
777 rt = _mm_mul_ps(r00,gbscale);
778 gbitab = _mm_cvttps_epi32(rt);
779 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
780 gbitab = _mm_slli_epi32(gbitab,2);
781 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
782 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
783 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
784 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
785 _MM_TRANSPOSE4_PS(Y,F,G,H);
786 Heps = _mm_mul_ps(gbeps,H);
787 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
788 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
789 vgb = _mm_mul_ps(gbqqfactor,VV);
791 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
792 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
793 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
794 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
795 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
796 /* 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. */
797 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
798 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
799 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
800 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
801 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
802 velec = _mm_mul_ps(qq00,rinv00);
803 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
805 /* LENNARD-JONES DISPERSION/REPULSION */
807 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
808 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
810 fscal = _mm_add_ps(felec,fvdw);
812 fscal = _mm_andnot_ps(dummy_mask,fscal);
814 /* Calculate temporary vectorial force */
815 tx = _mm_mul_ps(fscal,dx00);
816 ty = _mm_mul_ps(fscal,dy00);
817 tz = _mm_mul_ps(fscal,dz00);
819 /* Update vectorial force */
820 fix0 = _mm_add_ps(fix0,tx);
821 fiy0 = _mm_add_ps(fiy0,ty);
822 fiz0 = _mm_add_ps(fiz0,tz);
824 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
825 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
826 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
827 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
828 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
830 /* Inner loop uses 65 flops */
833 /* End of innermost loop */
835 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
836 f+i_coord_offset,fshift+i_shift_offset);
838 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
839 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
841 /* Increment number of inner iterations */
842 inneriter += j_index_end - j_index_start;
844 /* Outer loop uses 7 flops */
847 /* Increment number of outer iterations */
850 /* Update outer/inner flops */
852 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*65);