Fix segmentation fault in minimize
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_single / nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_sse2_single.cpp
blobb520d55179048db94df1ca05f7b0793f962da411
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017,2018, 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.
38 #include "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
47 #include "kernelutil_x86_sse2_single.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_VF_sse2_single
51 * Electrostatics interaction: Ewald
52 * VdW interaction: LennardJones
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
56 void
57 nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_VF_sse2_single
58 (t_nblist * gmx_restrict nlist,
59 rvec * gmx_restrict xx,
60 rvec * gmx_restrict ff,
61 struct t_forcerec * gmx_restrict fr,
62 t_mdatoms * gmx_restrict mdatoms,
63 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64 t_nrnb * gmx_restrict nrnb)
66 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
67 * just 0 for non-waters.
68 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
69 * jnr indices corresponding to data put in the four positions in the SIMD register.
71 int i_shift_offset,i_coord_offset,outeriter,inneriter;
72 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73 int jnrA,jnrB,jnrC,jnrD;
74 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
75 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
76 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
77 real rcutoff_scalar;
78 real *shiftvec,*fshift,*x,*f;
79 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
80 real scratch[4*DIM];
81 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
82 int vdwioffset0;
83 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
85 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
86 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
87 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
88 real *charge;
89 int nvdwtype;
90 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
91 int *vdwtype;
92 real *vdwparam;
93 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
94 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
95 __m128i ewitab;
96 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
97 real *ewtab;
98 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
99 real rswitch_scalar,d_scalar;
100 __m128 dummy_mask,cutoff_mask;
101 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
102 __m128 one = _mm_set1_ps(1.0);
103 __m128 two = _mm_set1_ps(2.0);
104 x = xx[0];
105 f = ff[0];
107 nri = nlist->nri;
108 iinr = nlist->iinr;
109 jindex = nlist->jindex;
110 jjnr = nlist->jjnr;
111 shiftidx = nlist->shift;
112 gid = nlist->gid;
113 shiftvec = fr->shift_vec[0];
114 fshift = fr->fshift[0];
115 facel = _mm_set1_ps(fr->ic->epsfac);
116 charge = mdatoms->chargeA;
117 nvdwtype = fr->ntype;
118 vdwparam = fr->nbfp;
119 vdwtype = mdatoms->typeA;
121 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
122 ewtab = fr->ic->tabq_coul_FDV0;
123 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
124 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
126 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
127 rcutoff_scalar = fr->ic->rcoulomb;
128 rcutoff = _mm_set1_ps(rcutoff_scalar);
129 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
131 rswitch_scalar = fr->ic->rcoulomb_switch;
132 rswitch = _mm_set1_ps(rswitch_scalar);
133 /* Setup switch parameters */
134 d_scalar = rcutoff_scalar-rswitch_scalar;
135 d = _mm_set1_ps(d_scalar);
136 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
137 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
138 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
139 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
140 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
141 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
143 /* Avoid stupid compiler warnings */
144 jnrA = jnrB = jnrC = jnrD = 0;
145 j_coord_offsetA = 0;
146 j_coord_offsetB = 0;
147 j_coord_offsetC = 0;
148 j_coord_offsetD = 0;
150 outeriter = 0;
151 inneriter = 0;
153 for(iidx=0;iidx<4*DIM;iidx++)
155 scratch[iidx] = 0.0;
158 /* Start outer loop over neighborlists */
159 for(iidx=0; iidx<nri; iidx++)
161 /* Load shift vector for this list */
162 i_shift_offset = DIM*shiftidx[iidx];
164 /* Load limits for loop over neighbors */
165 j_index_start = jindex[iidx];
166 j_index_end = jindex[iidx+1];
168 /* Get outer coordinate index */
169 inr = iinr[iidx];
170 i_coord_offset = DIM*inr;
172 /* Load i particle coords and add shift vector */
173 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
175 fix0 = _mm_setzero_ps();
176 fiy0 = _mm_setzero_ps();
177 fiz0 = _mm_setzero_ps();
179 /* Load parameters for i particles */
180 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
181 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
183 /* Reset potential sums */
184 velecsum = _mm_setzero_ps();
185 vvdwsum = _mm_setzero_ps();
187 /* Start inner kernel loop */
188 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
191 /* Get j neighbor index, and coordinate index */
192 jnrA = jjnr[jidx];
193 jnrB = jjnr[jidx+1];
194 jnrC = jjnr[jidx+2];
195 jnrD = jjnr[jidx+3];
196 j_coord_offsetA = DIM*jnrA;
197 j_coord_offsetB = DIM*jnrB;
198 j_coord_offsetC = DIM*jnrC;
199 j_coord_offsetD = DIM*jnrD;
201 /* load j atom coordinates */
202 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
203 x+j_coord_offsetC,x+j_coord_offsetD,
204 &jx0,&jy0,&jz0);
206 /* Calculate displacement vector */
207 dx00 = _mm_sub_ps(ix0,jx0);
208 dy00 = _mm_sub_ps(iy0,jy0);
209 dz00 = _mm_sub_ps(iz0,jz0);
211 /* Calculate squared distance and things based on it */
212 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
214 rinv00 = sse2_invsqrt_f(rsq00);
216 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
218 /* Load parameters for j particles */
219 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
220 charge+jnrC+0,charge+jnrD+0);
221 vdwjidx0A = 2*vdwtype[jnrA+0];
222 vdwjidx0B = 2*vdwtype[jnrB+0];
223 vdwjidx0C = 2*vdwtype[jnrC+0];
224 vdwjidx0D = 2*vdwtype[jnrD+0];
226 /**************************
227 * CALCULATE INTERACTIONS *
228 **************************/
230 if (gmx_mm_any_lt(rsq00,rcutoff2))
233 r00 = _mm_mul_ps(rsq00,rinv00);
235 /* Compute parameters for interactions between i and j atoms */
236 qq00 = _mm_mul_ps(iq0,jq0);
237 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
238 vdwparam+vdwioffset0+vdwjidx0B,
239 vdwparam+vdwioffset0+vdwjidx0C,
240 vdwparam+vdwioffset0+vdwjidx0D,
241 &c6_00,&c12_00);
243 /* EWALD ELECTROSTATICS */
245 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
246 ewrt = _mm_mul_ps(r00,ewtabscale);
247 ewitab = _mm_cvttps_epi32(ewrt);
248 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
249 ewitab = _mm_slli_epi32(ewitab,2);
250 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
251 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
252 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
253 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
254 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
255 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
256 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
257 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
258 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
260 /* LENNARD-JONES DISPERSION/REPULSION */
262 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
263 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
264 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
265 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
266 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
268 d = _mm_sub_ps(r00,rswitch);
269 d = _mm_max_ps(d,_mm_setzero_ps());
270 d2 = _mm_mul_ps(d,d);
271 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
273 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
275 /* Evaluate switch function */
276 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
277 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
278 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
279 velec = _mm_mul_ps(velec,sw);
280 vvdw = _mm_mul_ps(vvdw,sw);
281 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
283 /* Update potential sum for this i atom from the interaction with this j atom. */
284 velec = _mm_and_ps(velec,cutoff_mask);
285 velecsum = _mm_add_ps(velecsum,velec);
286 vvdw = _mm_and_ps(vvdw,cutoff_mask);
287 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
289 fscal = _mm_add_ps(felec,fvdw);
291 fscal = _mm_and_ps(fscal,cutoff_mask);
293 /* Calculate temporary vectorial force */
294 tx = _mm_mul_ps(fscal,dx00);
295 ty = _mm_mul_ps(fscal,dy00);
296 tz = _mm_mul_ps(fscal,dz00);
298 /* Update vectorial force */
299 fix0 = _mm_add_ps(fix0,tx);
300 fiy0 = _mm_add_ps(fiy0,ty);
301 fiz0 = _mm_add_ps(fiz0,tz);
303 fjptrA = f+j_coord_offsetA;
304 fjptrB = f+j_coord_offsetB;
305 fjptrC = f+j_coord_offsetC;
306 fjptrD = f+j_coord_offsetD;
307 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
311 /* Inner loop uses 83 flops */
314 if(jidx<j_index_end)
317 /* Get j neighbor index, and coordinate index */
318 jnrlistA = jjnr[jidx];
319 jnrlistB = jjnr[jidx+1];
320 jnrlistC = jjnr[jidx+2];
321 jnrlistD = jjnr[jidx+3];
322 /* Sign of each element will be negative for non-real atoms.
323 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
324 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
326 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
327 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
328 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
329 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
330 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
331 j_coord_offsetA = DIM*jnrA;
332 j_coord_offsetB = DIM*jnrB;
333 j_coord_offsetC = DIM*jnrC;
334 j_coord_offsetD = DIM*jnrD;
336 /* load j atom coordinates */
337 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
338 x+j_coord_offsetC,x+j_coord_offsetD,
339 &jx0,&jy0,&jz0);
341 /* Calculate displacement vector */
342 dx00 = _mm_sub_ps(ix0,jx0);
343 dy00 = _mm_sub_ps(iy0,jy0);
344 dz00 = _mm_sub_ps(iz0,jz0);
346 /* Calculate squared distance and things based on it */
347 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
349 rinv00 = sse2_invsqrt_f(rsq00);
351 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
353 /* Load parameters for j particles */
354 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
355 charge+jnrC+0,charge+jnrD+0);
356 vdwjidx0A = 2*vdwtype[jnrA+0];
357 vdwjidx0B = 2*vdwtype[jnrB+0];
358 vdwjidx0C = 2*vdwtype[jnrC+0];
359 vdwjidx0D = 2*vdwtype[jnrD+0];
361 /**************************
362 * CALCULATE INTERACTIONS *
363 **************************/
365 if (gmx_mm_any_lt(rsq00,rcutoff2))
368 r00 = _mm_mul_ps(rsq00,rinv00);
369 r00 = _mm_andnot_ps(dummy_mask,r00);
371 /* Compute parameters for interactions between i and j atoms */
372 qq00 = _mm_mul_ps(iq0,jq0);
373 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
374 vdwparam+vdwioffset0+vdwjidx0B,
375 vdwparam+vdwioffset0+vdwjidx0C,
376 vdwparam+vdwioffset0+vdwjidx0D,
377 &c6_00,&c12_00);
379 /* EWALD ELECTROSTATICS */
381 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
382 ewrt = _mm_mul_ps(r00,ewtabscale);
383 ewitab = _mm_cvttps_epi32(ewrt);
384 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
385 ewitab = _mm_slli_epi32(ewitab,2);
386 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
387 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
388 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
389 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
390 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
391 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
392 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
393 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
394 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
396 /* LENNARD-JONES DISPERSION/REPULSION */
398 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
399 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
400 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
401 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
402 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
404 d = _mm_sub_ps(r00,rswitch);
405 d = _mm_max_ps(d,_mm_setzero_ps());
406 d2 = _mm_mul_ps(d,d);
407 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
409 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
411 /* Evaluate switch function */
412 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
413 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
414 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
415 velec = _mm_mul_ps(velec,sw);
416 vvdw = _mm_mul_ps(vvdw,sw);
417 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
419 /* Update potential sum for this i atom from the interaction with this j atom. */
420 velec = _mm_and_ps(velec,cutoff_mask);
421 velec = _mm_andnot_ps(dummy_mask,velec);
422 velecsum = _mm_add_ps(velecsum,velec);
423 vvdw = _mm_and_ps(vvdw,cutoff_mask);
424 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
425 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
427 fscal = _mm_add_ps(felec,fvdw);
429 fscal = _mm_and_ps(fscal,cutoff_mask);
431 fscal = _mm_andnot_ps(dummy_mask,fscal);
433 /* Calculate temporary vectorial force */
434 tx = _mm_mul_ps(fscal,dx00);
435 ty = _mm_mul_ps(fscal,dy00);
436 tz = _mm_mul_ps(fscal,dz00);
438 /* Update vectorial force */
439 fix0 = _mm_add_ps(fix0,tx);
440 fiy0 = _mm_add_ps(fiy0,ty);
441 fiz0 = _mm_add_ps(fiz0,tz);
443 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
444 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
445 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
446 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
447 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
451 /* Inner loop uses 84 flops */
454 /* End of innermost loop */
456 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
457 f+i_coord_offset,fshift+i_shift_offset);
459 ggid = gid[iidx];
460 /* Update potential energies */
461 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
462 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
464 /* Increment number of inner iterations */
465 inneriter += j_index_end - j_index_start;
467 /* Outer loop uses 9 flops */
470 /* Increment number of outer iterations */
471 outeriter += nri;
473 /* Update outer/inner flops */
475 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*84);
478 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_F_sse2_single
479 * Electrostatics interaction: Ewald
480 * VdW interaction: LennardJones
481 * Geometry: Particle-Particle
482 * Calculate force/pot: Force
484 void
485 nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_F_sse2_single
486 (t_nblist * gmx_restrict nlist,
487 rvec * gmx_restrict xx,
488 rvec * gmx_restrict ff,
489 struct t_forcerec * gmx_restrict fr,
490 t_mdatoms * gmx_restrict mdatoms,
491 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
492 t_nrnb * gmx_restrict nrnb)
494 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
495 * just 0 for non-waters.
496 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
497 * jnr indices corresponding to data put in the four positions in the SIMD register.
499 int i_shift_offset,i_coord_offset,outeriter,inneriter;
500 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
501 int jnrA,jnrB,jnrC,jnrD;
502 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
503 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
504 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
505 real rcutoff_scalar;
506 real *shiftvec,*fshift,*x,*f;
507 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
508 real scratch[4*DIM];
509 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
510 int vdwioffset0;
511 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
512 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
513 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
514 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
515 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
516 real *charge;
517 int nvdwtype;
518 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
519 int *vdwtype;
520 real *vdwparam;
521 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
522 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
523 __m128i ewitab;
524 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
525 real *ewtab;
526 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
527 real rswitch_scalar,d_scalar;
528 __m128 dummy_mask,cutoff_mask;
529 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
530 __m128 one = _mm_set1_ps(1.0);
531 __m128 two = _mm_set1_ps(2.0);
532 x = xx[0];
533 f = ff[0];
535 nri = nlist->nri;
536 iinr = nlist->iinr;
537 jindex = nlist->jindex;
538 jjnr = nlist->jjnr;
539 shiftidx = nlist->shift;
540 gid = nlist->gid;
541 shiftvec = fr->shift_vec[0];
542 fshift = fr->fshift[0];
543 facel = _mm_set1_ps(fr->ic->epsfac);
544 charge = mdatoms->chargeA;
545 nvdwtype = fr->ntype;
546 vdwparam = fr->nbfp;
547 vdwtype = mdatoms->typeA;
549 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
550 ewtab = fr->ic->tabq_coul_FDV0;
551 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
552 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
554 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
555 rcutoff_scalar = fr->ic->rcoulomb;
556 rcutoff = _mm_set1_ps(rcutoff_scalar);
557 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
559 rswitch_scalar = fr->ic->rcoulomb_switch;
560 rswitch = _mm_set1_ps(rswitch_scalar);
561 /* Setup switch parameters */
562 d_scalar = rcutoff_scalar-rswitch_scalar;
563 d = _mm_set1_ps(d_scalar);
564 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
565 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
566 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
567 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
568 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
569 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
571 /* Avoid stupid compiler warnings */
572 jnrA = jnrB = jnrC = jnrD = 0;
573 j_coord_offsetA = 0;
574 j_coord_offsetB = 0;
575 j_coord_offsetC = 0;
576 j_coord_offsetD = 0;
578 outeriter = 0;
579 inneriter = 0;
581 for(iidx=0;iidx<4*DIM;iidx++)
583 scratch[iidx] = 0.0;
586 /* Start outer loop over neighborlists */
587 for(iidx=0; iidx<nri; iidx++)
589 /* Load shift vector for this list */
590 i_shift_offset = DIM*shiftidx[iidx];
592 /* Load limits for loop over neighbors */
593 j_index_start = jindex[iidx];
594 j_index_end = jindex[iidx+1];
596 /* Get outer coordinate index */
597 inr = iinr[iidx];
598 i_coord_offset = DIM*inr;
600 /* Load i particle coords and add shift vector */
601 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
603 fix0 = _mm_setzero_ps();
604 fiy0 = _mm_setzero_ps();
605 fiz0 = _mm_setzero_ps();
607 /* Load parameters for i particles */
608 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
609 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
611 /* Start inner kernel loop */
612 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
615 /* Get j neighbor index, and coordinate index */
616 jnrA = jjnr[jidx];
617 jnrB = jjnr[jidx+1];
618 jnrC = jjnr[jidx+2];
619 jnrD = jjnr[jidx+3];
620 j_coord_offsetA = DIM*jnrA;
621 j_coord_offsetB = DIM*jnrB;
622 j_coord_offsetC = DIM*jnrC;
623 j_coord_offsetD = DIM*jnrD;
625 /* load j atom coordinates */
626 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
627 x+j_coord_offsetC,x+j_coord_offsetD,
628 &jx0,&jy0,&jz0);
630 /* Calculate displacement vector */
631 dx00 = _mm_sub_ps(ix0,jx0);
632 dy00 = _mm_sub_ps(iy0,jy0);
633 dz00 = _mm_sub_ps(iz0,jz0);
635 /* Calculate squared distance and things based on it */
636 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
638 rinv00 = sse2_invsqrt_f(rsq00);
640 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
642 /* Load parameters for j particles */
643 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
644 charge+jnrC+0,charge+jnrD+0);
645 vdwjidx0A = 2*vdwtype[jnrA+0];
646 vdwjidx0B = 2*vdwtype[jnrB+0];
647 vdwjidx0C = 2*vdwtype[jnrC+0];
648 vdwjidx0D = 2*vdwtype[jnrD+0];
650 /**************************
651 * CALCULATE INTERACTIONS *
652 **************************/
654 if (gmx_mm_any_lt(rsq00,rcutoff2))
657 r00 = _mm_mul_ps(rsq00,rinv00);
659 /* Compute parameters for interactions between i and j atoms */
660 qq00 = _mm_mul_ps(iq0,jq0);
661 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
662 vdwparam+vdwioffset0+vdwjidx0B,
663 vdwparam+vdwioffset0+vdwjidx0C,
664 vdwparam+vdwioffset0+vdwjidx0D,
665 &c6_00,&c12_00);
667 /* EWALD ELECTROSTATICS */
669 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
670 ewrt = _mm_mul_ps(r00,ewtabscale);
671 ewitab = _mm_cvttps_epi32(ewrt);
672 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
673 ewitab = _mm_slli_epi32(ewitab,2);
674 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
675 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
676 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
677 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
678 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
679 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
680 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
681 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
682 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
684 /* LENNARD-JONES DISPERSION/REPULSION */
686 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
687 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
688 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
689 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
690 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
692 d = _mm_sub_ps(r00,rswitch);
693 d = _mm_max_ps(d,_mm_setzero_ps());
694 d2 = _mm_mul_ps(d,d);
695 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
697 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
699 /* Evaluate switch function */
700 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
701 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
702 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
703 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
705 fscal = _mm_add_ps(felec,fvdw);
707 fscal = _mm_and_ps(fscal,cutoff_mask);
709 /* Calculate temporary vectorial force */
710 tx = _mm_mul_ps(fscal,dx00);
711 ty = _mm_mul_ps(fscal,dy00);
712 tz = _mm_mul_ps(fscal,dz00);
714 /* Update vectorial force */
715 fix0 = _mm_add_ps(fix0,tx);
716 fiy0 = _mm_add_ps(fiy0,ty);
717 fiz0 = _mm_add_ps(fiz0,tz);
719 fjptrA = f+j_coord_offsetA;
720 fjptrB = f+j_coord_offsetB;
721 fjptrC = f+j_coord_offsetC;
722 fjptrD = f+j_coord_offsetD;
723 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
727 /* Inner loop uses 77 flops */
730 if(jidx<j_index_end)
733 /* Get j neighbor index, and coordinate index */
734 jnrlistA = jjnr[jidx];
735 jnrlistB = jjnr[jidx+1];
736 jnrlistC = jjnr[jidx+2];
737 jnrlistD = jjnr[jidx+3];
738 /* Sign of each element will be negative for non-real atoms.
739 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
740 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
742 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
743 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
744 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
745 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
746 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
747 j_coord_offsetA = DIM*jnrA;
748 j_coord_offsetB = DIM*jnrB;
749 j_coord_offsetC = DIM*jnrC;
750 j_coord_offsetD = DIM*jnrD;
752 /* load j atom coordinates */
753 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
754 x+j_coord_offsetC,x+j_coord_offsetD,
755 &jx0,&jy0,&jz0);
757 /* Calculate displacement vector */
758 dx00 = _mm_sub_ps(ix0,jx0);
759 dy00 = _mm_sub_ps(iy0,jy0);
760 dz00 = _mm_sub_ps(iz0,jz0);
762 /* Calculate squared distance and things based on it */
763 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
765 rinv00 = sse2_invsqrt_f(rsq00);
767 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
769 /* Load parameters for j particles */
770 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
771 charge+jnrC+0,charge+jnrD+0);
772 vdwjidx0A = 2*vdwtype[jnrA+0];
773 vdwjidx0B = 2*vdwtype[jnrB+0];
774 vdwjidx0C = 2*vdwtype[jnrC+0];
775 vdwjidx0D = 2*vdwtype[jnrD+0];
777 /**************************
778 * CALCULATE INTERACTIONS *
779 **************************/
781 if (gmx_mm_any_lt(rsq00,rcutoff2))
784 r00 = _mm_mul_ps(rsq00,rinv00);
785 r00 = _mm_andnot_ps(dummy_mask,r00);
787 /* Compute parameters for interactions between i and j atoms */
788 qq00 = _mm_mul_ps(iq0,jq0);
789 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
790 vdwparam+vdwioffset0+vdwjidx0B,
791 vdwparam+vdwioffset0+vdwjidx0C,
792 vdwparam+vdwioffset0+vdwjidx0D,
793 &c6_00,&c12_00);
795 /* EWALD ELECTROSTATICS */
797 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
798 ewrt = _mm_mul_ps(r00,ewtabscale);
799 ewitab = _mm_cvttps_epi32(ewrt);
800 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
801 ewitab = _mm_slli_epi32(ewitab,2);
802 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
803 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
804 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
805 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
806 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
807 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
808 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
809 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
810 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
812 /* LENNARD-JONES DISPERSION/REPULSION */
814 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
815 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
816 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
817 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
818 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
820 d = _mm_sub_ps(r00,rswitch);
821 d = _mm_max_ps(d,_mm_setzero_ps());
822 d2 = _mm_mul_ps(d,d);
823 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
825 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
827 /* Evaluate switch function */
828 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
829 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
830 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
831 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
833 fscal = _mm_add_ps(felec,fvdw);
835 fscal = _mm_and_ps(fscal,cutoff_mask);
837 fscal = _mm_andnot_ps(dummy_mask,fscal);
839 /* Calculate temporary vectorial force */
840 tx = _mm_mul_ps(fscal,dx00);
841 ty = _mm_mul_ps(fscal,dy00);
842 tz = _mm_mul_ps(fscal,dz00);
844 /* Update vectorial force */
845 fix0 = _mm_add_ps(fix0,tx);
846 fiy0 = _mm_add_ps(fiy0,ty);
847 fiz0 = _mm_add_ps(fiz0,tz);
849 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
850 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
851 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
852 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
853 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
857 /* Inner loop uses 78 flops */
860 /* End of innermost loop */
862 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
863 f+i_coord_offset,fshift+i_shift_offset);
865 /* Increment number of inner iterations */
866 inneriter += j_index_end - j_index_start;
868 /* Outer loop uses 7 flops */
871 /* Increment number of outer iterations */
872 outeriter += nri;
874 /* Update outer/inner flops */
876 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*78);