2 * Note: this file was generated by the Gromacs sse4_1_double kernel generator.
4 * This source code is part of
8 * Copyright (c) 2001-2012, The GROMACS Development Team
10 * Gromacs is a library for molecular simulation and trajectory analysis,
11 * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12 * a full list of developers and information, check out http://www.gromacs.org
14 * This program is free software; you can redistribute it and/or modify it under
15 * the terms of the GNU Lesser General Public License as published by the Free
16 * Software Foundation; either version 2 of the License, or (at your option) any
19 * To help fund GROMACS development, we humbly ask that you cite
20 * the papers people have written on it - you can find them on the website.
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
33 #include "gmx_math_x86_sse4_1_double.h"
34 #include "kernelutil_x86_sse4_1_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwLJ_GeomP1P1_VF_sse4_1_double
38 * Electrostatics interaction: Coulomb
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecCoul_VdwLJ_GeomP1P1_VF_sse4_1_double
45 (t_nblist * gmx_restrict nlist,
46 rvec * gmx_restrict xx,
47 rvec * gmx_restrict ff,
48 t_forcerec * gmx_restrict fr,
49 t_mdatoms * gmx_restrict mdatoms,
50 nb_kernel_data_t * gmx_restrict kernel_data,
51 t_nrnb * gmx_restrict nrnb)
53 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
54 * just 0 for non-waters.
55 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
56 * jnr indices corresponding to data put in the four positions in the SIMD register.
58 int i_shift_offset,i_coord_offset,outeriter,inneriter;
59 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
61 int j_coord_offsetA,j_coord_offsetB;
62 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
64 real *shiftvec,*fshift,*x,*f;
65 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
67 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
68 int vdwjidx0A,vdwjidx0B;
69 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
70 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
71 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
74 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
77 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
78 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
79 __m128d dummy_mask,cutoff_mask;
80 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
81 __m128d one = _mm_set1_pd(1.0);
82 __m128d two = _mm_set1_pd(2.0);
88 jindex = nlist->jindex;
90 shiftidx = nlist->shift;
92 shiftvec = fr->shift_vec[0];
93 fshift = fr->fshift[0];
94 facel = _mm_set1_pd(fr->epsfac);
95 charge = mdatoms->chargeA;
98 vdwtype = mdatoms->typeA;
100 /* Avoid stupid compiler warnings */
108 /* Start outer loop over neighborlists */
109 for(iidx=0; iidx<nri; iidx++)
111 /* Load shift vector for this list */
112 i_shift_offset = DIM*shiftidx[iidx];
114 /* Load limits for loop over neighbors */
115 j_index_start = jindex[iidx];
116 j_index_end = jindex[iidx+1];
118 /* Get outer coordinate index */
120 i_coord_offset = DIM*inr;
122 /* Load i particle coords and add shift vector */
123 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
125 fix0 = _mm_setzero_pd();
126 fiy0 = _mm_setzero_pd();
127 fiz0 = _mm_setzero_pd();
129 /* Load parameters for i particles */
130 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
131 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
133 /* Reset potential sums */
134 velecsum = _mm_setzero_pd();
135 vvdwsum = _mm_setzero_pd();
137 /* Start inner kernel loop */
138 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
141 /* Get j neighbor index, and coordinate index */
144 j_coord_offsetA = DIM*jnrA;
145 j_coord_offsetB = DIM*jnrB;
147 /* load j atom coordinates */
148 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
151 /* Calculate displacement vector */
152 dx00 = _mm_sub_pd(ix0,jx0);
153 dy00 = _mm_sub_pd(iy0,jy0);
154 dz00 = _mm_sub_pd(iz0,jz0);
156 /* Calculate squared distance and things based on it */
157 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
159 rinv00 = gmx_mm_invsqrt_pd(rsq00);
161 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
163 /* Load parameters for j particles */
164 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
165 vdwjidx0A = 2*vdwtype[jnrA+0];
166 vdwjidx0B = 2*vdwtype[jnrB+0];
168 /**************************
169 * CALCULATE INTERACTIONS *
170 **************************/
172 /* Compute parameters for interactions between i and j atoms */
173 qq00 = _mm_mul_pd(iq0,jq0);
174 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
175 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
177 /* COULOMB ELECTROSTATICS */
178 velec = _mm_mul_pd(qq00,rinv00);
179 felec = _mm_mul_pd(velec,rinvsq00);
181 /* LENNARD-JONES DISPERSION/REPULSION */
183 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
184 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
185 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
186 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
187 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
189 /* Update potential sum for this i atom from the interaction with this j atom. */
190 velecsum = _mm_add_pd(velecsum,velec);
191 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
193 fscal = _mm_add_pd(felec,fvdw);
195 /* Calculate temporary vectorial force */
196 tx = _mm_mul_pd(fscal,dx00);
197 ty = _mm_mul_pd(fscal,dy00);
198 tz = _mm_mul_pd(fscal,dz00);
200 /* Update vectorial force */
201 fix0 = _mm_add_pd(fix0,tx);
202 fiy0 = _mm_add_pd(fiy0,ty);
203 fiz0 = _mm_add_pd(fiz0,tz);
205 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
207 /* Inner loop uses 40 flops */
214 j_coord_offsetA = DIM*jnrA;
216 /* load j atom coordinates */
217 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
220 /* Calculate displacement vector */
221 dx00 = _mm_sub_pd(ix0,jx0);
222 dy00 = _mm_sub_pd(iy0,jy0);
223 dz00 = _mm_sub_pd(iz0,jz0);
225 /* Calculate squared distance and things based on it */
226 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
228 rinv00 = gmx_mm_invsqrt_pd(rsq00);
230 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
232 /* Load parameters for j particles */
233 jq0 = _mm_load_sd(charge+jnrA+0);
234 vdwjidx0A = 2*vdwtype[jnrA+0];
236 /**************************
237 * CALCULATE INTERACTIONS *
238 **************************/
240 /* Compute parameters for interactions between i and j atoms */
241 qq00 = _mm_mul_pd(iq0,jq0);
242 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
244 /* COULOMB ELECTROSTATICS */
245 velec = _mm_mul_pd(qq00,rinv00);
246 felec = _mm_mul_pd(velec,rinvsq00);
248 /* LENNARD-JONES DISPERSION/REPULSION */
250 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
251 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
252 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
253 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
254 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
256 /* Update potential sum for this i atom from the interaction with this j atom. */
257 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
258 velecsum = _mm_add_pd(velecsum,velec);
259 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
260 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
262 fscal = _mm_add_pd(felec,fvdw);
264 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
266 /* Calculate temporary vectorial force */
267 tx = _mm_mul_pd(fscal,dx00);
268 ty = _mm_mul_pd(fscal,dy00);
269 tz = _mm_mul_pd(fscal,dz00);
271 /* Update vectorial force */
272 fix0 = _mm_add_pd(fix0,tx);
273 fiy0 = _mm_add_pd(fiy0,ty);
274 fiz0 = _mm_add_pd(fiz0,tz);
276 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
278 /* Inner loop uses 40 flops */
281 /* End of innermost loop */
283 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
284 f+i_coord_offset,fshift+i_shift_offset);
287 /* Update potential energies */
288 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
289 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
291 /* Increment number of inner iterations */
292 inneriter += j_index_end - j_index_start;
294 /* Outer loop uses 9 flops */
297 /* Increment number of outer iterations */
300 /* Update outer/inner flops */
302 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*40);
305 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwLJ_GeomP1P1_F_sse4_1_double
306 * Electrostatics interaction: Coulomb
307 * VdW interaction: LennardJones
308 * Geometry: Particle-Particle
309 * Calculate force/pot: Force
312 nb_kernel_ElecCoul_VdwLJ_GeomP1P1_F_sse4_1_double
313 (t_nblist * gmx_restrict nlist,
314 rvec * gmx_restrict xx,
315 rvec * gmx_restrict ff,
316 t_forcerec * gmx_restrict fr,
317 t_mdatoms * gmx_restrict mdatoms,
318 nb_kernel_data_t * gmx_restrict kernel_data,
319 t_nrnb * gmx_restrict nrnb)
321 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
322 * just 0 for non-waters.
323 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
324 * jnr indices corresponding to data put in the four positions in the SIMD register.
326 int i_shift_offset,i_coord_offset,outeriter,inneriter;
327 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
329 int j_coord_offsetA,j_coord_offsetB;
330 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
332 real *shiftvec,*fshift,*x,*f;
333 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
335 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
336 int vdwjidx0A,vdwjidx0B;
337 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
338 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
339 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
342 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
345 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
346 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
347 __m128d dummy_mask,cutoff_mask;
348 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
349 __m128d one = _mm_set1_pd(1.0);
350 __m128d two = _mm_set1_pd(2.0);
356 jindex = nlist->jindex;
358 shiftidx = nlist->shift;
360 shiftvec = fr->shift_vec[0];
361 fshift = fr->fshift[0];
362 facel = _mm_set1_pd(fr->epsfac);
363 charge = mdatoms->chargeA;
364 nvdwtype = fr->ntype;
366 vdwtype = mdatoms->typeA;
368 /* Avoid stupid compiler warnings */
376 /* Start outer loop over neighborlists */
377 for(iidx=0; iidx<nri; iidx++)
379 /* Load shift vector for this list */
380 i_shift_offset = DIM*shiftidx[iidx];
382 /* Load limits for loop over neighbors */
383 j_index_start = jindex[iidx];
384 j_index_end = jindex[iidx+1];
386 /* Get outer coordinate index */
388 i_coord_offset = DIM*inr;
390 /* Load i particle coords and add shift vector */
391 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
393 fix0 = _mm_setzero_pd();
394 fiy0 = _mm_setzero_pd();
395 fiz0 = _mm_setzero_pd();
397 /* Load parameters for i particles */
398 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
399 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
401 /* Start inner kernel loop */
402 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
405 /* Get j neighbor index, and coordinate index */
408 j_coord_offsetA = DIM*jnrA;
409 j_coord_offsetB = DIM*jnrB;
411 /* load j atom coordinates */
412 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
415 /* Calculate displacement vector */
416 dx00 = _mm_sub_pd(ix0,jx0);
417 dy00 = _mm_sub_pd(iy0,jy0);
418 dz00 = _mm_sub_pd(iz0,jz0);
420 /* Calculate squared distance and things based on it */
421 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
423 rinv00 = gmx_mm_invsqrt_pd(rsq00);
425 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
427 /* Load parameters for j particles */
428 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
429 vdwjidx0A = 2*vdwtype[jnrA+0];
430 vdwjidx0B = 2*vdwtype[jnrB+0];
432 /**************************
433 * CALCULATE INTERACTIONS *
434 **************************/
436 /* Compute parameters for interactions between i and j atoms */
437 qq00 = _mm_mul_pd(iq0,jq0);
438 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
439 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
441 /* COULOMB ELECTROSTATICS */
442 velec = _mm_mul_pd(qq00,rinv00);
443 felec = _mm_mul_pd(velec,rinvsq00);
445 /* LENNARD-JONES DISPERSION/REPULSION */
447 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
448 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
450 fscal = _mm_add_pd(felec,fvdw);
452 /* Calculate temporary vectorial force */
453 tx = _mm_mul_pd(fscal,dx00);
454 ty = _mm_mul_pd(fscal,dy00);
455 tz = _mm_mul_pd(fscal,dz00);
457 /* Update vectorial force */
458 fix0 = _mm_add_pd(fix0,tx);
459 fiy0 = _mm_add_pd(fiy0,ty);
460 fiz0 = _mm_add_pd(fiz0,tz);
462 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
464 /* Inner loop uses 34 flops */
471 j_coord_offsetA = DIM*jnrA;
473 /* load j atom coordinates */
474 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
477 /* Calculate displacement vector */
478 dx00 = _mm_sub_pd(ix0,jx0);
479 dy00 = _mm_sub_pd(iy0,jy0);
480 dz00 = _mm_sub_pd(iz0,jz0);
482 /* Calculate squared distance and things based on it */
483 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
485 rinv00 = gmx_mm_invsqrt_pd(rsq00);
487 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
489 /* Load parameters for j particles */
490 jq0 = _mm_load_sd(charge+jnrA+0);
491 vdwjidx0A = 2*vdwtype[jnrA+0];
493 /**************************
494 * CALCULATE INTERACTIONS *
495 **************************/
497 /* Compute parameters for interactions between i and j atoms */
498 qq00 = _mm_mul_pd(iq0,jq0);
499 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
501 /* COULOMB ELECTROSTATICS */
502 velec = _mm_mul_pd(qq00,rinv00);
503 felec = _mm_mul_pd(velec,rinvsq00);
505 /* LENNARD-JONES DISPERSION/REPULSION */
507 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
508 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
510 fscal = _mm_add_pd(felec,fvdw);
512 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
514 /* Calculate temporary vectorial force */
515 tx = _mm_mul_pd(fscal,dx00);
516 ty = _mm_mul_pd(fscal,dy00);
517 tz = _mm_mul_pd(fscal,dz00);
519 /* Update vectorial force */
520 fix0 = _mm_add_pd(fix0,tx);
521 fiy0 = _mm_add_pd(fiy0,ty);
522 fiz0 = _mm_add_pd(fiz0,tz);
524 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
526 /* Inner loop uses 34 flops */
529 /* End of innermost loop */
531 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
532 f+i_coord_offset,fshift+i_shift_offset);
534 /* Increment number of inner iterations */
535 inneriter += j_index_end - j_index_start;
537 /* Outer loop uses 7 flops */
540 /* Increment number of outer iterations */
543 /* Update outer/inner flops */
545 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*34);