Fix segmentation fault in minimize
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse4_1_double / nb_kernel_ElecNone_VdwLJ_GeomP1P1_sse4_1_double.cpp
blob15a934fbe8687b95afb9e04a48342050291f0260
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 sse4_1_double kernel generator.
38 #include "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
47 #include "kernelutil_x86_sse4_1_double.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJ_GeomP1P1_VF_sse4_1_double
51 * Electrostatics interaction: None
52 * VdW interaction: LennardJones
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
56 void
57 nb_kernel_ElecNone_VdwLJ_GeomP1P1_VF_sse4_1_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
74 int j_coord_offsetA,j_coord_offsetB;
75 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
76 real rcutoff_scalar;
77 real *shiftvec,*fshift,*x,*f;
78 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
79 int vdwioffset0;
80 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
81 int vdwjidx0A,vdwjidx0B;
82 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
83 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
84 int nvdwtype;
85 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
86 int *vdwtype;
87 real *vdwparam;
88 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
89 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
90 __m128d dummy_mask,cutoff_mask;
91 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
92 __m128d one = _mm_set1_pd(1.0);
93 __m128d two = _mm_set1_pd(2.0);
94 x = xx[0];
95 f = ff[0];
97 nri = nlist->nri;
98 iinr = nlist->iinr;
99 jindex = nlist->jindex;
100 jjnr = nlist->jjnr;
101 shiftidx = nlist->shift;
102 gid = nlist->gid;
103 shiftvec = fr->shift_vec[0];
104 fshift = fr->fshift[0];
105 nvdwtype = fr->ntype;
106 vdwparam = fr->nbfp;
107 vdwtype = mdatoms->typeA;
109 /* Avoid stupid compiler warnings */
110 jnrA = jnrB = 0;
111 j_coord_offsetA = 0;
112 j_coord_offsetB = 0;
114 outeriter = 0;
115 inneriter = 0;
117 /* Start outer loop over neighborlists */
118 for(iidx=0; iidx<nri; iidx++)
120 /* Load shift vector for this list */
121 i_shift_offset = DIM*shiftidx[iidx];
123 /* Load limits for loop over neighbors */
124 j_index_start = jindex[iidx];
125 j_index_end = jindex[iidx+1];
127 /* Get outer coordinate index */
128 inr = iinr[iidx];
129 i_coord_offset = DIM*inr;
131 /* Load i particle coords and add shift vector */
132 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
134 fix0 = _mm_setzero_pd();
135 fiy0 = _mm_setzero_pd();
136 fiz0 = _mm_setzero_pd();
138 /* Load parameters for i particles */
139 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
141 /* Reset potential sums */
142 vvdwsum = _mm_setzero_pd();
144 /* Start inner kernel loop */
145 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
148 /* Get j neighbor index, and coordinate index */
149 jnrA = jjnr[jidx];
150 jnrB = jjnr[jidx+1];
151 j_coord_offsetA = DIM*jnrA;
152 j_coord_offsetB = DIM*jnrB;
154 /* load j atom coordinates */
155 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
156 &jx0,&jy0,&jz0);
158 /* Calculate displacement vector */
159 dx00 = _mm_sub_pd(ix0,jx0);
160 dy00 = _mm_sub_pd(iy0,jy0);
161 dz00 = _mm_sub_pd(iz0,jz0);
163 /* Calculate squared distance and things based on it */
164 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
166 rinvsq00 = sse41_inv_d(rsq00);
168 /* Load parameters for j particles */
169 vdwjidx0A = 2*vdwtype[jnrA+0];
170 vdwjidx0B = 2*vdwtype[jnrB+0];
172 /**************************
173 * CALCULATE INTERACTIONS *
174 **************************/
176 /* Compute parameters for interactions between i and j atoms */
177 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
178 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
180 /* LENNARD-JONES DISPERSION/REPULSION */
182 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
183 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
184 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
185 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
186 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
188 /* Update potential sum for this i atom from the interaction with this j atom. */
189 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
191 fscal = fvdw;
193 /* Calculate temporary vectorial force */
194 tx = _mm_mul_pd(fscal,dx00);
195 ty = _mm_mul_pd(fscal,dy00);
196 tz = _mm_mul_pd(fscal,dz00);
198 /* Update vectorial force */
199 fix0 = _mm_add_pd(fix0,tx);
200 fiy0 = _mm_add_pd(fiy0,ty);
201 fiz0 = _mm_add_pd(fiz0,tz);
203 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
205 /* Inner loop uses 32 flops */
208 if(jidx<j_index_end)
211 jnrA = jjnr[jidx];
212 j_coord_offsetA = DIM*jnrA;
214 /* load j atom coordinates */
215 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
216 &jx0,&jy0,&jz0);
218 /* Calculate displacement vector */
219 dx00 = _mm_sub_pd(ix0,jx0);
220 dy00 = _mm_sub_pd(iy0,jy0);
221 dz00 = _mm_sub_pd(iz0,jz0);
223 /* Calculate squared distance and things based on it */
224 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
226 rinvsq00 = sse41_inv_d(rsq00);
228 /* Load parameters for j particles */
229 vdwjidx0A = 2*vdwtype[jnrA+0];
231 /**************************
232 * CALCULATE INTERACTIONS *
233 **************************/
235 /* Compute parameters for interactions between i and j atoms */
236 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
238 /* LENNARD-JONES DISPERSION/REPULSION */
240 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
241 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
242 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
243 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
244 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
246 /* Update potential sum for this i atom from the interaction with this j atom. */
247 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
248 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
250 fscal = fvdw;
252 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
254 /* Calculate temporary vectorial force */
255 tx = _mm_mul_pd(fscal,dx00);
256 ty = _mm_mul_pd(fscal,dy00);
257 tz = _mm_mul_pd(fscal,dz00);
259 /* Update vectorial force */
260 fix0 = _mm_add_pd(fix0,tx);
261 fiy0 = _mm_add_pd(fiy0,ty);
262 fiz0 = _mm_add_pd(fiz0,tz);
264 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
266 /* Inner loop uses 32 flops */
269 /* End of innermost loop */
271 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
272 f+i_coord_offset,fshift+i_shift_offset);
274 ggid = gid[iidx];
275 /* Update potential energies */
276 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
278 /* Increment number of inner iterations */
279 inneriter += j_index_end - j_index_start;
281 /* Outer loop uses 7 flops */
284 /* Increment number of outer iterations */
285 outeriter += nri;
287 /* Update outer/inner flops */
289 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*32);
292 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJ_GeomP1P1_F_sse4_1_double
293 * Electrostatics interaction: None
294 * VdW interaction: LennardJones
295 * Geometry: Particle-Particle
296 * Calculate force/pot: Force
298 void
299 nb_kernel_ElecNone_VdwLJ_GeomP1P1_F_sse4_1_double
300 (t_nblist * gmx_restrict nlist,
301 rvec * gmx_restrict xx,
302 rvec * gmx_restrict ff,
303 struct t_forcerec * gmx_restrict fr,
304 t_mdatoms * gmx_restrict mdatoms,
305 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
306 t_nrnb * gmx_restrict nrnb)
308 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
309 * just 0 for non-waters.
310 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
311 * jnr indices corresponding to data put in the four positions in the SIMD register.
313 int i_shift_offset,i_coord_offset,outeriter,inneriter;
314 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
315 int jnrA,jnrB;
316 int j_coord_offsetA,j_coord_offsetB;
317 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
318 real rcutoff_scalar;
319 real *shiftvec,*fshift,*x,*f;
320 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
321 int vdwioffset0;
322 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
323 int vdwjidx0A,vdwjidx0B;
324 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
325 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
326 int nvdwtype;
327 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
328 int *vdwtype;
329 real *vdwparam;
330 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
331 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
332 __m128d dummy_mask,cutoff_mask;
333 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
334 __m128d one = _mm_set1_pd(1.0);
335 __m128d two = _mm_set1_pd(2.0);
336 x = xx[0];
337 f = ff[0];
339 nri = nlist->nri;
340 iinr = nlist->iinr;
341 jindex = nlist->jindex;
342 jjnr = nlist->jjnr;
343 shiftidx = nlist->shift;
344 gid = nlist->gid;
345 shiftvec = fr->shift_vec[0];
346 fshift = fr->fshift[0];
347 nvdwtype = fr->ntype;
348 vdwparam = fr->nbfp;
349 vdwtype = mdatoms->typeA;
351 /* Avoid stupid compiler warnings */
352 jnrA = jnrB = 0;
353 j_coord_offsetA = 0;
354 j_coord_offsetB = 0;
356 outeriter = 0;
357 inneriter = 0;
359 /* Start outer loop over neighborlists */
360 for(iidx=0; iidx<nri; iidx++)
362 /* Load shift vector for this list */
363 i_shift_offset = DIM*shiftidx[iidx];
365 /* Load limits for loop over neighbors */
366 j_index_start = jindex[iidx];
367 j_index_end = jindex[iidx+1];
369 /* Get outer coordinate index */
370 inr = iinr[iidx];
371 i_coord_offset = DIM*inr;
373 /* Load i particle coords and add shift vector */
374 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
376 fix0 = _mm_setzero_pd();
377 fiy0 = _mm_setzero_pd();
378 fiz0 = _mm_setzero_pd();
380 /* Load parameters for i particles */
381 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
383 /* Start inner kernel loop */
384 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
387 /* Get j neighbor index, and coordinate index */
388 jnrA = jjnr[jidx];
389 jnrB = jjnr[jidx+1];
390 j_coord_offsetA = DIM*jnrA;
391 j_coord_offsetB = DIM*jnrB;
393 /* load j atom coordinates */
394 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
395 &jx0,&jy0,&jz0);
397 /* Calculate displacement vector */
398 dx00 = _mm_sub_pd(ix0,jx0);
399 dy00 = _mm_sub_pd(iy0,jy0);
400 dz00 = _mm_sub_pd(iz0,jz0);
402 /* Calculate squared distance and things based on it */
403 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
405 rinvsq00 = sse41_inv_d(rsq00);
407 /* Load parameters for j particles */
408 vdwjidx0A = 2*vdwtype[jnrA+0];
409 vdwjidx0B = 2*vdwtype[jnrB+0];
411 /**************************
412 * CALCULATE INTERACTIONS *
413 **************************/
415 /* Compute parameters for interactions between i and j atoms */
416 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
417 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
419 /* LENNARD-JONES DISPERSION/REPULSION */
421 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
422 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
424 fscal = fvdw;
426 /* Calculate temporary vectorial force */
427 tx = _mm_mul_pd(fscal,dx00);
428 ty = _mm_mul_pd(fscal,dy00);
429 tz = _mm_mul_pd(fscal,dz00);
431 /* Update vectorial force */
432 fix0 = _mm_add_pd(fix0,tx);
433 fiy0 = _mm_add_pd(fiy0,ty);
434 fiz0 = _mm_add_pd(fiz0,tz);
436 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
438 /* Inner loop uses 27 flops */
441 if(jidx<j_index_end)
444 jnrA = jjnr[jidx];
445 j_coord_offsetA = DIM*jnrA;
447 /* load j atom coordinates */
448 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
449 &jx0,&jy0,&jz0);
451 /* Calculate displacement vector */
452 dx00 = _mm_sub_pd(ix0,jx0);
453 dy00 = _mm_sub_pd(iy0,jy0);
454 dz00 = _mm_sub_pd(iz0,jz0);
456 /* Calculate squared distance and things based on it */
457 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
459 rinvsq00 = sse41_inv_d(rsq00);
461 /* Load parameters for j particles */
462 vdwjidx0A = 2*vdwtype[jnrA+0];
464 /**************************
465 * CALCULATE INTERACTIONS *
466 **************************/
468 /* Compute parameters for interactions between i and j atoms */
469 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
471 /* LENNARD-JONES DISPERSION/REPULSION */
473 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
474 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
476 fscal = fvdw;
478 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
480 /* Calculate temporary vectorial force */
481 tx = _mm_mul_pd(fscal,dx00);
482 ty = _mm_mul_pd(fscal,dy00);
483 tz = _mm_mul_pd(fscal,dz00);
485 /* Update vectorial force */
486 fix0 = _mm_add_pd(fix0,tx);
487 fiy0 = _mm_add_pd(fiy0,ty);
488 fiz0 = _mm_add_pd(fiz0,tz);
490 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
492 /* Inner loop uses 27 flops */
495 /* End of innermost loop */
497 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
498 f+i_coord_offset,fshift+i_shift_offset);
500 /* Increment number of inner iterations */
501 inneriter += j_index_end - j_index_start;
503 /* Outer loop uses 6 flops */
506 /* Increment number of outer iterations */
507 outeriter += nri;
509 /* Update outer/inner flops */
511 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*27);