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 sse2_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_sse2_single.h"
48 #include "kernelutil_x86_sse2_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_VF_sse2_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_sse2_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_cvtepi32_ps(gbitab));
244 gbitab = _mm_slli_epi32(gbitab,2);
246 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
247 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
248 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
249 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
250 _MM_TRANSPOSE4_PS(Y,F,G,H);
251 Heps = _mm_mul_ps(gbeps,H);
252 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
253 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
254 vgb = _mm_mul_ps(gbqqfactor,VV);
256 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
257 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
258 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
259 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
264 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
265 velec = _mm_mul_ps(qq00,rinv00);
266 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
268 /* LENNARD-JONES DISPERSION/REPULSION */
270 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
271 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
272 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
273 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
274 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
276 /* Update potential sum for this i atom from the interaction with this j atom. */
277 velecsum = _mm_add_ps(velecsum,velec);
278 vgbsum = _mm_add_ps(vgbsum,vgb);
279 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
281 fscal = _mm_add_ps(felec,fvdw);
283 /* Calculate temporary vectorial force */
284 tx = _mm_mul_ps(fscal,dx00);
285 ty = _mm_mul_ps(fscal,dy00);
286 tz = _mm_mul_ps(fscal,dz00);
288 /* Update vectorial force */
289 fix0 = _mm_add_ps(fix0,tx);
290 fiy0 = _mm_add_ps(fiy0,ty);
291 fiz0 = _mm_add_ps(fiz0,tz);
293 fjptrA = f+j_coord_offsetA;
294 fjptrB = f+j_coord_offsetB;
295 fjptrC = f+j_coord_offsetC;
296 fjptrD = f+j_coord_offsetD;
297 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
299 /* Inner loop uses 71 flops */
305 /* Get j neighbor index, and coordinate index */
306 jnrlistA = jjnr[jidx];
307 jnrlistB = jjnr[jidx+1];
308 jnrlistC = jjnr[jidx+2];
309 jnrlistD = jjnr[jidx+3];
310 /* Sign of each element will be negative for non-real atoms.
311 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
312 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
314 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
315 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
316 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
317 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
318 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
319 j_coord_offsetA = DIM*jnrA;
320 j_coord_offsetB = DIM*jnrB;
321 j_coord_offsetC = DIM*jnrC;
322 j_coord_offsetD = DIM*jnrD;
324 /* load j atom coordinates */
325 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
326 x+j_coord_offsetC,x+j_coord_offsetD,
329 /* Calculate displacement vector */
330 dx00 = _mm_sub_ps(ix0,jx0);
331 dy00 = _mm_sub_ps(iy0,jy0);
332 dz00 = _mm_sub_ps(iz0,jz0);
334 /* Calculate squared distance and things based on it */
335 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
337 rinv00 = gmx_mm_invsqrt_ps(rsq00);
339 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
341 /* Load parameters for j particles */
342 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
343 charge+jnrC+0,charge+jnrD+0);
344 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
345 invsqrta+jnrC+0,invsqrta+jnrD+0);
346 vdwjidx0A = 2*vdwtype[jnrA+0];
347 vdwjidx0B = 2*vdwtype[jnrB+0];
348 vdwjidx0C = 2*vdwtype[jnrC+0];
349 vdwjidx0D = 2*vdwtype[jnrD+0];
351 /**************************
352 * CALCULATE INTERACTIONS *
353 **************************/
355 r00 = _mm_mul_ps(rsq00,rinv00);
356 r00 = _mm_andnot_ps(dummy_mask,r00);
358 /* Compute parameters for interactions between i and j atoms */
359 qq00 = _mm_mul_ps(iq0,jq0);
360 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
361 vdwparam+vdwioffset0+vdwjidx0B,
362 vdwparam+vdwioffset0+vdwjidx0C,
363 vdwparam+vdwioffset0+vdwjidx0D,
366 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
367 isaprod = _mm_mul_ps(isai0,isaj0);
368 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
369 gbscale = _mm_mul_ps(isaprod,gbtabscale);
371 /* Calculate generalized born table index - this is a separate table from the normal one,
372 * but we use the same procedure by multiplying r with scale and truncating to integer.
374 rt = _mm_mul_ps(r00,gbscale);
375 gbitab = _mm_cvttps_epi32(rt);
376 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
377 gbitab = _mm_slli_epi32(gbitab,2);
379 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
380 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
381 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
382 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
383 _MM_TRANSPOSE4_PS(Y,F,G,H);
384 Heps = _mm_mul_ps(gbeps,H);
385 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
386 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
387 vgb = _mm_mul_ps(gbqqfactor,VV);
389 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
390 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
391 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
392 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
393 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
394 /* 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. */
395 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
396 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
397 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
398 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
399 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
400 velec = _mm_mul_ps(qq00,rinv00);
401 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
403 /* LENNARD-JONES DISPERSION/REPULSION */
405 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
406 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
407 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
408 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
409 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
411 /* Update potential sum for this i atom from the interaction with this j atom. */
412 velec = _mm_andnot_ps(dummy_mask,velec);
413 velecsum = _mm_add_ps(velecsum,velec);
414 vgb = _mm_andnot_ps(dummy_mask,vgb);
415 vgbsum = _mm_add_ps(vgbsum,vgb);
416 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
417 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
419 fscal = _mm_add_ps(felec,fvdw);
421 fscal = _mm_andnot_ps(dummy_mask,fscal);
423 /* Calculate temporary vectorial force */
424 tx = _mm_mul_ps(fscal,dx00);
425 ty = _mm_mul_ps(fscal,dy00);
426 tz = _mm_mul_ps(fscal,dz00);
428 /* Update vectorial force */
429 fix0 = _mm_add_ps(fix0,tx);
430 fiy0 = _mm_add_ps(fiy0,ty);
431 fiz0 = _mm_add_ps(fiz0,tz);
433 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
434 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
435 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
436 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
437 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
439 /* Inner loop uses 72 flops */
442 /* End of innermost loop */
444 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
445 f+i_coord_offset,fshift+i_shift_offset);
448 /* Update potential energies */
449 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
450 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
451 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
452 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
453 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
455 /* Increment number of inner iterations */
456 inneriter += j_index_end - j_index_start;
458 /* Outer loop uses 10 flops */
461 /* Increment number of outer iterations */
464 /* Update outer/inner flops */
466 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*72);
469 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse2_single
470 * Electrostatics interaction: GeneralizedBorn
471 * VdW interaction: LennardJones
472 * Geometry: Particle-Particle
473 * Calculate force/pot: Force
476 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse2_single
477 (t_nblist * gmx_restrict nlist,
478 rvec * gmx_restrict xx,
479 rvec * gmx_restrict ff,
480 t_forcerec * gmx_restrict fr,
481 t_mdatoms * gmx_restrict mdatoms,
482 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
483 t_nrnb * gmx_restrict nrnb)
485 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
486 * just 0 for non-waters.
487 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
488 * jnr indices corresponding to data put in the four positions in the SIMD register.
490 int i_shift_offset,i_coord_offset,outeriter,inneriter;
491 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
492 int jnrA,jnrB,jnrC,jnrD;
493 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
494 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
495 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
497 real *shiftvec,*fshift,*x,*f;
498 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
500 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
502 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
503 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
504 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
505 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
506 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
509 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
510 __m128 minushalf = _mm_set1_ps(-0.5);
511 real *invsqrta,*dvda,*gbtab;
513 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
516 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
517 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
519 __m128i ifour = _mm_set1_epi32(4);
520 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
522 __m128 dummy_mask,cutoff_mask;
523 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
524 __m128 one = _mm_set1_ps(1.0);
525 __m128 two = _mm_set1_ps(2.0);
531 jindex = nlist->jindex;
533 shiftidx = nlist->shift;
535 shiftvec = fr->shift_vec[0];
536 fshift = fr->fshift[0];
537 facel = _mm_set1_ps(fr->epsfac);
538 charge = mdatoms->chargeA;
539 nvdwtype = fr->ntype;
541 vdwtype = mdatoms->typeA;
543 invsqrta = fr->invsqrta;
545 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
546 gbtab = fr->gbtab.data;
547 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
549 /* Avoid stupid compiler warnings */
550 jnrA = jnrB = jnrC = jnrD = 0;
559 for(iidx=0;iidx<4*DIM;iidx++)
564 /* Start outer loop over neighborlists */
565 for(iidx=0; iidx<nri; iidx++)
567 /* Load shift vector for this list */
568 i_shift_offset = DIM*shiftidx[iidx];
570 /* Load limits for loop over neighbors */
571 j_index_start = jindex[iidx];
572 j_index_end = jindex[iidx+1];
574 /* Get outer coordinate index */
576 i_coord_offset = DIM*inr;
578 /* Load i particle coords and add shift vector */
579 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
581 fix0 = _mm_setzero_ps();
582 fiy0 = _mm_setzero_ps();
583 fiz0 = _mm_setzero_ps();
585 /* Load parameters for i particles */
586 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
587 isai0 = _mm_load1_ps(invsqrta+inr+0);
588 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
590 dvdasum = _mm_setzero_ps();
592 /* Start inner kernel loop */
593 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
596 /* Get j neighbor index, and coordinate index */
601 j_coord_offsetA = DIM*jnrA;
602 j_coord_offsetB = DIM*jnrB;
603 j_coord_offsetC = DIM*jnrC;
604 j_coord_offsetD = DIM*jnrD;
606 /* load j atom coordinates */
607 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
608 x+j_coord_offsetC,x+j_coord_offsetD,
611 /* Calculate displacement vector */
612 dx00 = _mm_sub_ps(ix0,jx0);
613 dy00 = _mm_sub_ps(iy0,jy0);
614 dz00 = _mm_sub_ps(iz0,jz0);
616 /* Calculate squared distance and things based on it */
617 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
619 rinv00 = gmx_mm_invsqrt_ps(rsq00);
621 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
623 /* Load parameters for j particles */
624 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
625 charge+jnrC+0,charge+jnrD+0);
626 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
627 invsqrta+jnrC+0,invsqrta+jnrD+0);
628 vdwjidx0A = 2*vdwtype[jnrA+0];
629 vdwjidx0B = 2*vdwtype[jnrB+0];
630 vdwjidx0C = 2*vdwtype[jnrC+0];
631 vdwjidx0D = 2*vdwtype[jnrD+0];
633 /**************************
634 * CALCULATE INTERACTIONS *
635 **************************/
637 r00 = _mm_mul_ps(rsq00,rinv00);
639 /* Compute parameters for interactions between i and j atoms */
640 qq00 = _mm_mul_ps(iq0,jq0);
641 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
642 vdwparam+vdwioffset0+vdwjidx0B,
643 vdwparam+vdwioffset0+vdwjidx0C,
644 vdwparam+vdwioffset0+vdwjidx0D,
647 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
648 isaprod = _mm_mul_ps(isai0,isaj0);
649 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
650 gbscale = _mm_mul_ps(isaprod,gbtabscale);
652 /* Calculate generalized born table index - this is a separate table from the normal one,
653 * but we use the same procedure by multiplying r with scale and truncating to integer.
655 rt = _mm_mul_ps(r00,gbscale);
656 gbitab = _mm_cvttps_epi32(rt);
657 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
658 gbitab = _mm_slli_epi32(gbitab,2);
660 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
661 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
662 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
663 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
664 _MM_TRANSPOSE4_PS(Y,F,G,H);
665 Heps = _mm_mul_ps(gbeps,H);
666 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
667 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
668 vgb = _mm_mul_ps(gbqqfactor,VV);
670 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
671 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
672 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
673 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
678 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
679 velec = _mm_mul_ps(qq00,rinv00);
680 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
682 /* LENNARD-JONES DISPERSION/REPULSION */
684 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
685 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
687 fscal = _mm_add_ps(felec,fvdw);
689 /* Calculate temporary vectorial force */
690 tx = _mm_mul_ps(fscal,dx00);
691 ty = _mm_mul_ps(fscal,dy00);
692 tz = _mm_mul_ps(fscal,dz00);
694 /* Update vectorial force */
695 fix0 = _mm_add_ps(fix0,tx);
696 fiy0 = _mm_add_ps(fiy0,ty);
697 fiz0 = _mm_add_ps(fiz0,tz);
699 fjptrA = f+j_coord_offsetA;
700 fjptrB = f+j_coord_offsetB;
701 fjptrC = f+j_coord_offsetC;
702 fjptrD = f+j_coord_offsetD;
703 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
705 /* Inner loop uses 64 flops */
711 /* Get j neighbor index, and coordinate index */
712 jnrlistA = jjnr[jidx];
713 jnrlistB = jjnr[jidx+1];
714 jnrlistC = jjnr[jidx+2];
715 jnrlistD = jjnr[jidx+3];
716 /* Sign of each element will be negative for non-real atoms.
717 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
718 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
720 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
721 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
722 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
723 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
724 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
725 j_coord_offsetA = DIM*jnrA;
726 j_coord_offsetB = DIM*jnrB;
727 j_coord_offsetC = DIM*jnrC;
728 j_coord_offsetD = DIM*jnrD;
730 /* load j atom coordinates */
731 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
732 x+j_coord_offsetC,x+j_coord_offsetD,
735 /* Calculate displacement vector */
736 dx00 = _mm_sub_ps(ix0,jx0);
737 dy00 = _mm_sub_ps(iy0,jy0);
738 dz00 = _mm_sub_ps(iz0,jz0);
740 /* Calculate squared distance and things based on it */
741 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
743 rinv00 = gmx_mm_invsqrt_ps(rsq00);
745 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
747 /* Load parameters for j particles */
748 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
749 charge+jnrC+0,charge+jnrD+0);
750 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
751 invsqrta+jnrC+0,invsqrta+jnrD+0);
752 vdwjidx0A = 2*vdwtype[jnrA+0];
753 vdwjidx0B = 2*vdwtype[jnrB+0];
754 vdwjidx0C = 2*vdwtype[jnrC+0];
755 vdwjidx0D = 2*vdwtype[jnrD+0];
757 /**************************
758 * CALCULATE INTERACTIONS *
759 **************************/
761 r00 = _mm_mul_ps(rsq00,rinv00);
762 r00 = _mm_andnot_ps(dummy_mask,r00);
764 /* Compute parameters for interactions between i and j atoms */
765 qq00 = _mm_mul_ps(iq0,jq0);
766 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
767 vdwparam+vdwioffset0+vdwjidx0B,
768 vdwparam+vdwioffset0+vdwjidx0C,
769 vdwparam+vdwioffset0+vdwjidx0D,
772 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
773 isaprod = _mm_mul_ps(isai0,isaj0);
774 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
775 gbscale = _mm_mul_ps(isaprod,gbtabscale);
777 /* Calculate generalized born table index - this is a separate table from the normal one,
778 * but we use the same procedure by multiplying r with scale and truncating to integer.
780 rt = _mm_mul_ps(r00,gbscale);
781 gbitab = _mm_cvttps_epi32(rt);
782 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
783 gbitab = _mm_slli_epi32(gbitab,2);
785 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
786 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
787 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
788 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
789 _MM_TRANSPOSE4_PS(Y,F,G,H);
790 Heps = _mm_mul_ps(gbeps,H);
791 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
792 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
793 vgb = _mm_mul_ps(gbqqfactor,VV);
795 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
796 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
797 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
798 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
799 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
800 /* 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. */
801 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
802 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
803 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
804 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
805 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
806 velec = _mm_mul_ps(qq00,rinv00);
807 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
809 /* LENNARD-JONES DISPERSION/REPULSION */
811 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
812 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
814 fscal = _mm_add_ps(felec,fvdw);
816 fscal = _mm_andnot_ps(dummy_mask,fscal);
818 /* Calculate temporary vectorial force */
819 tx = _mm_mul_ps(fscal,dx00);
820 ty = _mm_mul_ps(fscal,dy00);
821 tz = _mm_mul_ps(fscal,dz00);
823 /* Update vectorial force */
824 fix0 = _mm_add_ps(fix0,tx);
825 fiy0 = _mm_add_ps(fiy0,ty);
826 fiz0 = _mm_add_ps(fiz0,tz);
828 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
829 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
830 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
831 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
832 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
834 /* Inner loop uses 65 flops */
837 /* End of innermost loop */
839 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
840 f+i_coord_offset,fshift+i_shift_offset);
842 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
843 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
845 /* Increment number of inner iterations */
846 inneriter += j_index_end - j_index_start;
848 /* Outer loop uses 7 flops */
851 /* Increment number of outer iterations */
854 /* Update outer/inner flops */
856 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*65);