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 avx_128_fma_double kernel generator.
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
49 #include "gromacs/simd/math_x86_avx_128_fma_double.h"
50 #include "kernelutil_x86_avx_128_fma_double.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW4P1_VF_avx_128_fma_double
54 * Electrostatics interaction: ReactionField
55 * VdW interaction: LennardJones
56 * Geometry: Water4-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecRF_VdwLJ_GeomW4P1_VF_avx_128_fma_double
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
77 int j_coord_offsetA,j_coord_offsetB;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
80 real *shiftvec,*fshift,*x,*f;
81 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
87 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
89 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
90 int vdwjidx0A,vdwjidx0B;
91 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
92 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
93 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
94 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
95 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
96 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
99 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
102 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
103 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
104 __m128d dummy_mask,cutoff_mask;
105 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
106 __m128d one = _mm_set1_pd(1.0);
107 __m128d two = _mm_set1_pd(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_pd(fr->epsfac);
120 charge = mdatoms->chargeA;
121 krf = _mm_set1_pd(fr->ic->k_rf);
122 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
123 crf = _mm_set1_pd(fr->ic->c_rf);
124 nvdwtype = fr->ntype;
126 vdwtype = mdatoms->typeA;
128 /* Setup water-specific parameters */
129 inr = nlist->iinr[0];
130 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
131 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
132 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
133 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
135 /* Avoid stupid compiler warnings */
143 /* Start outer loop over neighborlists */
144 for(iidx=0; iidx<nri; iidx++)
146 /* Load shift vector for this list */
147 i_shift_offset = DIM*shiftidx[iidx];
149 /* Load limits for loop over neighbors */
150 j_index_start = jindex[iidx];
151 j_index_end = jindex[iidx+1];
153 /* Get outer coordinate index */
155 i_coord_offset = DIM*inr;
157 /* Load i particle coords and add shift vector */
158 gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
159 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
161 fix0 = _mm_setzero_pd();
162 fiy0 = _mm_setzero_pd();
163 fiz0 = _mm_setzero_pd();
164 fix1 = _mm_setzero_pd();
165 fiy1 = _mm_setzero_pd();
166 fiz1 = _mm_setzero_pd();
167 fix2 = _mm_setzero_pd();
168 fiy2 = _mm_setzero_pd();
169 fiz2 = _mm_setzero_pd();
170 fix3 = _mm_setzero_pd();
171 fiy3 = _mm_setzero_pd();
172 fiz3 = _mm_setzero_pd();
174 /* Reset potential sums */
175 velecsum = _mm_setzero_pd();
176 vvdwsum = _mm_setzero_pd();
178 /* Start inner kernel loop */
179 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
182 /* Get j neighbor index, and coordinate index */
185 j_coord_offsetA = DIM*jnrA;
186 j_coord_offsetB = DIM*jnrB;
188 /* load j atom coordinates */
189 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
192 /* Calculate displacement vector */
193 dx00 = _mm_sub_pd(ix0,jx0);
194 dy00 = _mm_sub_pd(iy0,jy0);
195 dz00 = _mm_sub_pd(iz0,jz0);
196 dx10 = _mm_sub_pd(ix1,jx0);
197 dy10 = _mm_sub_pd(iy1,jy0);
198 dz10 = _mm_sub_pd(iz1,jz0);
199 dx20 = _mm_sub_pd(ix2,jx0);
200 dy20 = _mm_sub_pd(iy2,jy0);
201 dz20 = _mm_sub_pd(iz2,jz0);
202 dx30 = _mm_sub_pd(ix3,jx0);
203 dy30 = _mm_sub_pd(iy3,jy0);
204 dz30 = _mm_sub_pd(iz3,jz0);
206 /* Calculate squared distance and things based on it */
207 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
208 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
209 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
210 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
212 rinv10 = gmx_mm_invsqrt_pd(rsq10);
213 rinv20 = gmx_mm_invsqrt_pd(rsq20);
214 rinv30 = gmx_mm_invsqrt_pd(rsq30);
216 rinvsq00 = gmx_mm_inv_pd(rsq00);
217 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
218 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
219 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
221 /* Load parameters for j particles */
222 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
223 vdwjidx0A = 2*vdwtype[jnrA+0];
224 vdwjidx0B = 2*vdwtype[jnrB+0];
226 fjx0 = _mm_setzero_pd();
227 fjy0 = _mm_setzero_pd();
228 fjz0 = _mm_setzero_pd();
230 /**************************
231 * CALCULATE INTERACTIONS *
232 **************************/
234 /* Compute parameters for interactions between i and j atoms */
235 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
236 vdwparam+vdwioffset0+vdwjidx0B,&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_msub_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 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
251 /* Update vectorial force */
252 fix0 = _mm_macc_pd(dx00,fscal,fix0);
253 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
254 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
256 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
257 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
258 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
260 /**************************
261 * CALCULATE INTERACTIONS *
262 **************************/
264 /* Compute parameters for interactions between i and j atoms */
265 qq10 = _mm_mul_pd(iq1,jq0);
267 /* REACTION-FIELD ELECTROSTATICS */
268 velec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_macc_pd(krf,rsq10,rinv10),crf));
269 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
271 /* Update potential sum for this i atom from the interaction with this j atom. */
272 velecsum = _mm_add_pd(velecsum,velec);
276 /* Update vectorial force */
277 fix1 = _mm_macc_pd(dx10,fscal,fix1);
278 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
279 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
281 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
282 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
283 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
285 /**************************
286 * CALCULATE INTERACTIONS *
287 **************************/
289 /* Compute parameters for interactions between i and j atoms */
290 qq20 = _mm_mul_pd(iq2,jq0);
292 /* REACTION-FIELD ELECTROSTATICS */
293 velec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_macc_pd(krf,rsq20,rinv20),crf));
294 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
296 /* Update potential sum for this i atom from the interaction with this j atom. */
297 velecsum = _mm_add_pd(velecsum,velec);
301 /* Update vectorial force */
302 fix2 = _mm_macc_pd(dx20,fscal,fix2);
303 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
304 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
306 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
307 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
308 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
310 /**************************
311 * CALCULATE INTERACTIONS *
312 **************************/
314 /* Compute parameters for interactions between i and j atoms */
315 qq30 = _mm_mul_pd(iq3,jq0);
317 /* REACTION-FIELD ELECTROSTATICS */
318 velec = _mm_mul_pd(qq30,_mm_sub_pd(_mm_macc_pd(krf,rsq30,rinv30),crf));
319 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
321 /* Update potential sum for this i atom from the interaction with this j atom. */
322 velecsum = _mm_add_pd(velecsum,velec);
326 /* Update vectorial force */
327 fix3 = _mm_macc_pd(dx30,fscal,fix3);
328 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
329 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
331 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
332 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
333 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
335 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
337 /* Inner loop uses 143 flops */
344 j_coord_offsetA = DIM*jnrA;
346 /* load j atom coordinates */
347 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
350 /* Calculate displacement vector */
351 dx00 = _mm_sub_pd(ix0,jx0);
352 dy00 = _mm_sub_pd(iy0,jy0);
353 dz00 = _mm_sub_pd(iz0,jz0);
354 dx10 = _mm_sub_pd(ix1,jx0);
355 dy10 = _mm_sub_pd(iy1,jy0);
356 dz10 = _mm_sub_pd(iz1,jz0);
357 dx20 = _mm_sub_pd(ix2,jx0);
358 dy20 = _mm_sub_pd(iy2,jy0);
359 dz20 = _mm_sub_pd(iz2,jz0);
360 dx30 = _mm_sub_pd(ix3,jx0);
361 dy30 = _mm_sub_pd(iy3,jy0);
362 dz30 = _mm_sub_pd(iz3,jz0);
364 /* Calculate squared distance and things based on it */
365 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
366 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
367 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
368 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
370 rinv10 = gmx_mm_invsqrt_pd(rsq10);
371 rinv20 = gmx_mm_invsqrt_pd(rsq20);
372 rinv30 = gmx_mm_invsqrt_pd(rsq30);
374 rinvsq00 = gmx_mm_inv_pd(rsq00);
375 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
376 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
377 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
379 /* Load parameters for j particles */
380 jq0 = _mm_load_sd(charge+jnrA+0);
381 vdwjidx0A = 2*vdwtype[jnrA+0];
383 fjx0 = _mm_setzero_pd();
384 fjy0 = _mm_setzero_pd();
385 fjz0 = _mm_setzero_pd();
387 /**************************
388 * CALCULATE INTERACTIONS *
389 **************************/
391 /* Compute parameters for interactions between i and j atoms */
392 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
394 /* LENNARD-JONES DISPERSION/REPULSION */
396 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
397 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
398 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
399 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
400 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
402 /* Update potential sum for this i atom from the interaction with this j atom. */
403 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
404 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
408 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
410 /* Update vectorial force */
411 fix0 = _mm_macc_pd(dx00,fscal,fix0);
412 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
413 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
415 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
416 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
417 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
419 /**************************
420 * CALCULATE INTERACTIONS *
421 **************************/
423 /* Compute parameters for interactions between i and j atoms */
424 qq10 = _mm_mul_pd(iq1,jq0);
426 /* REACTION-FIELD ELECTROSTATICS */
427 velec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_macc_pd(krf,rsq10,rinv10),crf));
428 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
430 /* Update potential sum for this i atom from the interaction with this j atom. */
431 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
432 velecsum = _mm_add_pd(velecsum,velec);
436 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
438 /* Update vectorial force */
439 fix1 = _mm_macc_pd(dx10,fscal,fix1);
440 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
441 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
443 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
444 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
445 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
447 /**************************
448 * CALCULATE INTERACTIONS *
449 **************************/
451 /* Compute parameters for interactions between i and j atoms */
452 qq20 = _mm_mul_pd(iq2,jq0);
454 /* REACTION-FIELD ELECTROSTATICS */
455 velec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_macc_pd(krf,rsq20,rinv20),crf));
456 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
458 /* Update potential sum for this i atom from the interaction with this j atom. */
459 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
460 velecsum = _mm_add_pd(velecsum,velec);
464 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
466 /* Update vectorial force */
467 fix2 = _mm_macc_pd(dx20,fscal,fix2);
468 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
469 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
471 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
472 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
473 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
475 /**************************
476 * CALCULATE INTERACTIONS *
477 **************************/
479 /* Compute parameters for interactions between i and j atoms */
480 qq30 = _mm_mul_pd(iq3,jq0);
482 /* REACTION-FIELD ELECTROSTATICS */
483 velec = _mm_mul_pd(qq30,_mm_sub_pd(_mm_macc_pd(krf,rsq30,rinv30),crf));
484 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
486 /* Update potential sum for this i atom from the interaction with this j atom. */
487 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
488 velecsum = _mm_add_pd(velecsum,velec);
492 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
494 /* Update vectorial force */
495 fix3 = _mm_macc_pd(dx30,fscal,fix3);
496 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
497 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
499 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
500 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
501 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
503 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
505 /* Inner loop uses 143 flops */
508 /* End of innermost loop */
510 gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
511 f+i_coord_offset,fshift+i_shift_offset);
514 /* Update potential energies */
515 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
516 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
518 /* Increment number of inner iterations */
519 inneriter += j_index_end - j_index_start;
521 /* Outer loop uses 26 flops */
524 /* Increment number of outer iterations */
527 /* Update outer/inner flops */
529 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*143);
532 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW4P1_F_avx_128_fma_double
533 * Electrostatics interaction: ReactionField
534 * VdW interaction: LennardJones
535 * Geometry: Water4-Particle
536 * Calculate force/pot: Force
539 nb_kernel_ElecRF_VdwLJ_GeomW4P1_F_avx_128_fma_double
540 (t_nblist * gmx_restrict nlist,
541 rvec * gmx_restrict xx,
542 rvec * gmx_restrict ff,
543 t_forcerec * gmx_restrict fr,
544 t_mdatoms * gmx_restrict mdatoms,
545 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
546 t_nrnb * gmx_restrict nrnb)
548 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
549 * just 0 for non-waters.
550 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
551 * jnr indices corresponding to data put in the four positions in the SIMD register.
553 int i_shift_offset,i_coord_offset,outeriter,inneriter;
554 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
556 int j_coord_offsetA,j_coord_offsetB;
557 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
559 real *shiftvec,*fshift,*x,*f;
560 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
562 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
564 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
566 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
568 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
569 int vdwjidx0A,vdwjidx0B;
570 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
571 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
572 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
573 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
574 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
575 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
578 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
581 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
582 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
583 __m128d dummy_mask,cutoff_mask;
584 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
585 __m128d one = _mm_set1_pd(1.0);
586 __m128d two = _mm_set1_pd(2.0);
592 jindex = nlist->jindex;
594 shiftidx = nlist->shift;
596 shiftvec = fr->shift_vec[0];
597 fshift = fr->fshift[0];
598 facel = _mm_set1_pd(fr->epsfac);
599 charge = mdatoms->chargeA;
600 krf = _mm_set1_pd(fr->ic->k_rf);
601 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
602 crf = _mm_set1_pd(fr->ic->c_rf);
603 nvdwtype = fr->ntype;
605 vdwtype = mdatoms->typeA;
607 /* Setup water-specific parameters */
608 inr = nlist->iinr[0];
609 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
610 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
611 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
612 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
614 /* Avoid stupid compiler warnings */
622 /* Start outer loop over neighborlists */
623 for(iidx=0; iidx<nri; iidx++)
625 /* Load shift vector for this list */
626 i_shift_offset = DIM*shiftidx[iidx];
628 /* Load limits for loop over neighbors */
629 j_index_start = jindex[iidx];
630 j_index_end = jindex[iidx+1];
632 /* Get outer coordinate index */
634 i_coord_offset = DIM*inr;
636 /* Load i particle coords and add shift vector */
637 gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
638 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
640 fix0 = _mm_setzero_pd();
641 fiy0 = _mm_setzero_pd();
642 fiz0 = _mm_setzero_pd();
643 fix1 = _mm_setzero_pd();
644 fiy1 = _mm_setzero_pd();
645 fiz1 = _mm_setzero_pd();
646 fix2 = _mm_setzero_pd();
647 fiy2 = _mm_setzero_pd();
648 fiz2 = _mm_setzero_pd();
649 fix3 = _mm_setzero_pd();
650 fiy3 = _mm_setzero_pd();
651 fiz3 = _mm_setzero_pd();
653 /* Start inner kernel loop */
654 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
657 /* Get j neighbor index, and coordinate index */
660 j_coord_offsetA = DIM*jnrA;
661 j_coord_offsetB = DIM*jnrB;
663 /* load j atom coordinates */
664 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
667 /* Calculate displacement vector */
668 dx00 = _mm_sub_pd(ix0,jx0);
669 dy00 = _mm_sub_pd(iy0,jy0);
670 dz00 = _mm_sub_pd(iz0,jz0);
671 dx10 = _mm_sub_pd(ix1,jx0);
672 dy10 = _mm_sub_pd(iy1,jy0);
673 dz10 = _mm_sub_pd(iz1,jz0);
674 dx20 = _mm_sub_pd(ix2,jx0);
675 dy20 = _mm_sub_pd(iy2,jy0);
676 dz20 = _mm_sub_pd(iz2,jz0);
677 dx30 = _mm_sub_pd(ix3,jx0);
678 dy30 = _mm_sub_pd(iy3,jy0);
679 dz30 = _mm_sub_pd(iz3,jz0);
681 /* Calculate squared distance and things based on it */
682 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
683 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
684 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
685 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
687 rinv10 = gmx_mm_invsqrt_pd(rsq10);
688 rinv20 = gmx_mm_invsqrt_pd(rsq20);
689 rinv30 = gmx_mm_invsqrt_pd(rsq30);
691 rinvsq00 = gmx_mm_inv_pd(rsq00);
692 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
693 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
694 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
696 /* Load parameters for j particles */
697 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
698 vdwjidx0A = 2*vdwtype[jnrA+0];
699 vdwjidx0B = 2*vdwtype[jnrB+0];
701 fjx0 = _mm_setzero_pd();
702 fjy0 = _mm_setzero_pd();
703 fjz0 = _mm_setzero_pd();
705 /**************************
706 * CALCULATE INTERACTIONS *
707 **************************/
709 /* Compute parameters for interactions between i and j atoms */
710 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
711 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
713 /* LENNARD-JONES DISPERSION/REPULSION */
715 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
716 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
720 /* Update vectorial force */
721 fix0 = _mm_macc_pd(dx00,fscal,fix0);
722 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
723 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
725 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
726 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
727 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
729 /**************************
730 * CALCULATE INTERACTIONS *
731 **************************/
733 /* Compute parameters for interactions between i and j atoms */
734 qq10 = _mm_mul_pd(iq1,jq0);
736 /* REACTION-FIELD ELECTROSTATICS */
737 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
741 /* Update vectorial force */
742 fix1 = _mm_macc_pd(dx10,fscal,fix1);
743 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
744 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
746 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
747 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
748 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
750 /**************************
751 * CALCULATE INTERACTIONS *
752 **************************/
754 /* Compute parameters for interactions between i and j atoms */
755 qq20 = _mm_mul_pd(iq2,jq0);
757 /* REACTION-FIELD ELECTROSTATICS */
758 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
762 /* Update vectorial force */
763 fix2 = _mm_macc_pd(dx20,fscal,fix2);
764 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
765 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
767 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
768 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
769 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
771 /**************************
772 * CALCULATE INTERACTIONS *
773 **************************/
775 /* Compute parameters for interactions between i and j atoms */
776 qq30 = _mm_mul_pd(iq3,jq0);
778 /* REACTION-FIELD ELECTROSTATICS */
779 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
783 /* Update vectorial force */
784 fix3 = _mm_macc_pd(dx30,fscal,fix3);
785 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
786 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
788 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
789 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
790 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
792 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
794 /* Inner loop uses 123 flops */
801 j_coord_offsetA = DIM*jnrA;
803 /* load j atom coordinates */
804 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
807 /* Calculate displacement vector */
808 dx00 = _mm_sub_pd(ix0,jx0);
809 dy00 = _mm_sub_pd(iy0,jy0);
810 dz00 = _mm_sub_pd(iz0,jz0);
811 dx10 = _mm_sub_pd(ix1,jx0);
812 dy10 = _mm_sub_pd(iy1,jy0);
813 dz10 = _mm_sub_pd(iz1,jz0);
814 dx20 = _mm_sub_pd(ix2,jx0);
815 dy20 = _mm_sub_pd(iy2,jy0);
816 dz20 = _mm_sub_pd(iz2,jz0);
817 dx30 = _mm_sub_pd(ix3,jx0);
818 dy30 = _mm_sub_pd(iy3,jy0);
819 dz30 = _mm_sub_pd(iz3,jz0);
821 /* Calculate squared distance and things based on it */
822 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
823 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
824 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
825 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
827 rinv10 = gmx_mm_invsqrt_pd(rsq10);
828 rinv20 = gmx_mm_invsqrt_pd(rsq20);
829 rinv30 = gmx_mm_invsqrt_pd(rsq30);
831 rinvsq00 = gmx_mm_inv_pd(rsq00);
832 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
833 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
834 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
836 /* Load parameters for j particles */
837 jq0 = _mm_load_sd(charge+jnrA+0);
838 vdwjidx0A = 2*vdwtype[jnrA+0];
840 fjx0 = _mm_setzero_pd();
841 fjy0 = _mm_setzero_pd();
842 fjz0 = _mm_setzero_pd();
844 /**************************
845 * CALCULATE INTERACTIONS *
846 **************************/
848 /* Compute parameters for interactions between i and j atoms */
849 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
851 /* LENNARD-JONES DISPERSION/REPULSION */
853 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
854 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
858 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
860 /* Update vectorial force */
861 fix0 = _mm_macc_pd(dx00,fscal,fix0);
862 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
863 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
865 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
866 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
867 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
869 /**************************
870 * CALCULATE INTERACTIONS *
871 **************************/
873 /* Compute parameters for interactions between i and j atoms */
874 qq10 = _mm_mul_pd(iq1,jq0);
876 /* REACTION-FIELD ELECTROSTATICS */
877 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
881 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
883 /* Update vectorial force */
884 fix1 = _mm_macc_pd(dx10,fscal,fix1);
885 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
886 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
888 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
889 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
890 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
892 /**************************
893 * CALCULATE INTERACTIONS *
894 **************************/
896 /* Compute parameters for interactions between i and j atoms */
897 qq20 = _mm_mul_pd(iq2,jq0);
899 /* REACTION-FIELD ELECTROSTATICS */
900 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
904 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
906 /* Update vectorial force */
907 fix2 = _mm_macc_pd(dx20,fscal,fix2);
908 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
909 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
911 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
912 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
913 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
915 /**************************
916 * CALCULATE INTERACTIONS *
917 **************************/
919 /* Compute parameters for interactions between i and j atoms */
920 qq30 = _mm_mul_pd(iq3,jq0);
922 /* REACTION-FIELD ELECTROSTATICS */
923 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
927 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
929 /* Update vectorial force */
930 fix3 = _mm_macc_pd(dx30,fscal,fix3);
931 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
932 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
934 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
935 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
936 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
938 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
940 /* Inner loop uses 123 flops */
943 /* End of innermost loop */
945 gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
946 f+i_coord_offset,fshift+i_shift_offset);
948 /* Increment number of inner iterations */
949 inneriter += j_index_end - j_index_start;
951 /* Outer loop uses 24 flops */
954 /* Increment number of outer iterations */
957 /* Update outer/inner flops */
959 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*123);