Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_double / nb_kernel_ElecCoul_VdwLJ_GeomW4P1_avx_128_fma_double.c
blob1b9d3490d772814b9cbf13cfafaceeef2a7e3fd4
1 /*
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.
38 #include "config.h"
40 #include <math.h>
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
47 #include "gromacs/simd/math_x86_avx_128_fma_double.h"
48 #include "kernelutil_x86_avx_128_fma_double.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwLJ_GeomW4P1_VF_avx_128_fma_double
52 * Electrostatics interaction: Coulomb
53 * VdW interaction: LennardJones
54 * Geometry: Water4-Particle
55 * Calculate force/pot: PotentialAndForce
57 void
58 nb_kernel_ElecCoul_VdwLJ_GeomW4P1_VF_avx_128_fma_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
75 int j_coord_offsetA,j_coord_offsetB;
76 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
77 real rcutoff_scalar;
78 real *shiftvec,*fshift,*x,*f;
79 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80 int vdwioffset0;
81 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
82 int vdwioffset1;
83 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
84 int vdwioffset2;
85 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
86 int vdwioffset3;
87 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
88 int vdwjidx0A,vdwjidx0B;
89 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
90 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
91 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
92 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
93 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
94 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
95 real *charge;
96 int nvdwtype;
97 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
98 int *vdwtype;
99 real *vdwparam;
100 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
101 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
102 __m128d dummy_mask,cutoff_mask;
103 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
104 __m128d one = _mm_set1_pd(1.0);
105 __m128d two = _mm_set1_pd(2.0);
106 x = xx[0];
107 f = ff[0];
109 nri = nlist->nri;
110 iinr = nlist->iinr;
111 jindex = nlist->jindex;
112 jjnr = nlist->jjnr;
113 shiftidx = nlist->shift;
114 gid = nlist->gid;
115 shiftvec = fr->shift_vec[0];
116 fshift = fr->fshift[0];
117 facel = _mm_set1_pd(fr->epsfac);
118 charge = mdatoms->chargeA;
119 nvdwtype = fr->ntype;
120 vdwparam = fr->nbfp;
121 vdwtype = mdatoms->typeA;
123 /* Setup water-specific parameters */
124 inr = nlist->iinr[0];
125 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
126 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
127 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
128 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
130 /* Avoid stupid compiler warnings */
131 jnrA = jnrB = 0;
132 j_coord_offsetA = 0;
133 j_coord_offsetB = 0;
135 outeriter = 0;
136 inneriter = 0;
138 /* Start outer loop over neighborlists */
139 for(iidx=0; iidx<nri; iidx++)
141 /* Load shift vector for this list */
142 i_shift_offset = DIM*shiftidx[iidx];
144 /* Load limits for loop over neighbors */
145 j_index_start = jindex[iidx];
146 j_index_end = jindex[iidx+1];
148 /* Get outer coordinate index */
149 inr = iinr[iidx];
150 i_coord_offset = DIM*inr;
152 /* Load i particle coords and add shift vector */
153 gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
154 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
156 fix0 = _mm_setzero_pd();
157 fiy0 = _mm_setzero_pd();
158 fiz0 = _mm_setzero_pd();
159 fix1 = _mm_setzero_pd();
160 fiy1 = _mm_setzero_pd();
161 fiz1 = _mm_setzero_pd();
162 fix2 = _mm_setzero_pd();
163 fiy2 = _mm_setzero_pd();
164 fiz2 = _mm_setzero_pd();
165 fix3 = _mm_setzero_pd();
166 fiy3 = _mm_setzero_pd();
167 fiz3 = _mm_setzero_pd();
169 /* Reset potential sums */
170 velecsum = _mm_setzero_pd();
171 vvdwsum = _mm_setzero_pd();
173 /* Start inner kernel loop */
174 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
177 /* Get j neighbor index, and coordinate index */
178 jnrA = jjnr[jidx];
179 jnrB = jjnr[jidx+1];
180 j_coord_offsetA = DIM*jnrA;
181 j_coord_offsetB = DIM*jnrB;
183 /* load j atom coordinates */
184 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
185 &jx0,&jy0,&jz0);
187 /* Calculate displacement vector */
188 dx00 = _mm_sub_pd(ix0,jx0);
189 dy00 = _mm_sub_pd(iy0,jy0);
190 dz00 = _mm_sub_pd(iz0,jz0);
191 dx10 = _mm_sub_pd(ix1,jx0);
192 dy10 = _mm_sub_pd(iy1,jy0);
193 dz10 = _mm_sub_pd(iz1,jz0);
194 dx20 = _mm_sub_pd(ix2,jx0);
195 dy20 = _mm_sub_pd(iy2,jy0);
196 dz20 = _mm_sub_pd(iz2,jz0);
197 dx30 = _mm_sub_pd(ix3,jx0);
198 dy30 = _mm_sub_pd(iy3,jy0);
199 dz30 = _mm_sub_pd(iz3,jz0);
201 /* Calculate squared distance and things based on it */
202 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
203 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
204 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
205 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
207 rinv10 = gmx_mm_invsqrt_pd(rsq10);
208 rinv20 = gmx_mm_invsqrt_pd(rsq20);
209 rinv30 = gmx_mm_invsqrt_pd(rsq30);
211 rinvsq00 = gmx_mm_inv_pd(rsq00);
212 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
213 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
214 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
216 /* Load parameters for j particles */
217 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
218 vdwjidx0A = 2*vdwtype[jnrA+0];
219 vdwjidx0B = 2*vdwtype[jnrB+0];
221 fjx0 = _mm_setzero_pd();
222 fjy0 = _mm_setzero_pd();
223 fjz0 = _mm_setzero_pd();
225 /**************************
226 * CALCULATE INTERACTIONS *
227 **************************/
229 /* Compute parameters for interactions between i and j atoms */
230 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
231 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
233 /* LENNARD-JONES DISPERSION/REPULSION */
235 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
236 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
237 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
238 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
239 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
241 /* Update potential sum for this i atom from the interaction with this j atom. */
242 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
244 fscal = fvdw;
246 /* Update vectorial force */
247 fix0 = _mm_macc_pd(dx00,fscal,fix0);
248 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
249 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
251 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
252 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
253 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
255 /**************************
256 * CALCULATE INTERACTIONS *
257 **************************/
259 /* Compute parameters for interactions between i and j atoms */
260 qq10 = _mm_mul_pd(iq1,jq0);
262 /* COULOMB ELECTROSTATICS */
263 velec = _mm_mul_pd(qq10,rinv10);
264 felec = _mm_mul_pd(velec,rinvsq10);
266 /* Update potential sum for this i atom from the interaction with this j atom. */
267 velecsum = _mm_add_pd(velecsum,velec);
269 fscal = felec;
271 /* Update vectorial force */
272 fix1 = _mm_macc_pd(dx10,fscal,fix1);
273 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
274 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
276 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
277 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
278 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
280 /**************************
281 * CALCULATE INTERACTIONS *
282 **************************/
284 /* Compute parameters for interactions between i and j atoms */
285 qq20 = _mm_mul_pd(iq2,jq0);
287 /* COULOMB ELECTROSTATICS */
288 velec = _mm_mul_pd(qq20,rinv20);
289 felec = _mm_mul_pd(velec,rinvsq20);
291 /* Update potential sum for this i atom from the interaction with this j atom. */
292 velecsum = _mm_add_pd(velecsum,velec);
294 fscal = felec;
296 /* Update vectorial force */
297 fix2 = _mm_macc_pd(dx20,fscal,fix2);
298 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
299 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
301 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
302 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
303 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
305 /**************************
306 * CALCULATE INTERACTIONS *
307 **************************/
309 /* Compute parameters for interactions between i and j atoms */
310 qq30 = _mm_mul_pd(iq3,jq0);
312 /* COULOMB ELECTROSTATICS */
313 velec = _mm_mul_pd(qq30,rinv30);
314 felec = _mm_mul_pd(velec,rinvsq30);
316 /* Update potential sum for this i atom from the interaction with this j atom. */
317 velecsum = _mm_add_pd(velecsum,velec);
319 fscal = felec;
321 /* Update vectorial force */
322 fix3 = _mm_macc_pd(dx30,fscal,fix3);
323 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
324 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
326 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
327 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
328 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
330 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
332 /* Inner loop uses 131 flops */
335 if(jidx<j_index_end)
338 jnrA = jjnr[jidx];
339 j_coord_offsetA = DIM*jnrA;
341 /* load j atom coordinates */
342 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
343 &jx0,&jy0,&jz0);
345 /* Calculate displacement vector */
346 dx00 = _mm_sub_pd(ix0,jx0);
347 dy00 = _mm_sub_pd(iy0,jy0);
348 dz00 = _mm_sub_pd(iz0,jz0);
349 dx10 = _mm_sub_pd(ix1,jx0);
350 dy10 = _mm_sub_pd(iy1,jy0);
351 dz10 = _mm_sub_pd(iz1,jz0);
352 dx20 = _mm_sub_pd(ix2,jx0);
353 dy20 = _mm_sub_pd(iy2,jy0);
354 dz20 = _mm_sub_pd(iz2,jz0);
355 dx30 = _mm_sub_pd(ix3,jx0);
356 dy30 = _mm_sub_pd(iy3,jy0);
357 dz30 = _mm_sub_pd(iz3,jz0);
359 /* Calculate squared distance and things based on it */
360 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
361 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
362 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
363 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
365 rinv10 = gmx_mm_invsqrt_pd(rsq10);
366 rinv20 = gmx_mm_invsqrt_pd(rsq20);
367 rinv30 = gmx_mm_invsqrt_pd(rsq30);
369 rinvsq00 = gmx_mm_inv_pd(rsq00);
370 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
371 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
372 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
374 /* Load parameters for j particles */
375 jq0 = _mm_load_sd(charge+jnrA+0);
376 vdwjidx0A = 2*vdwtype[jnrA+0];
378 fjx0 = _mm_setzero_pd();
379 fjy0 = _mm_setzero_pd();
380 fjz0 = _mm_setzero_pd();
382 /**************************
383 * CALCULATE INTERACTIONS *
384 **************************/
386 /* Compute parameters for interactions between i and j atoms */
387 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
389 /* LENNARD-JONES DISPERSION/REPULSION */
391 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
392 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
393 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
394 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
395 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
397 /* Update potential sum for this i atom from the interaction with this j atom. */
398 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
399 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
401 fscal = fvdw;
403 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
405 /* Update vectorial force */
406 fix0 = _mm_macc_pd(dx00,fscal,fix0);
407 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
408 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
410 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
411 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
412 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
414 /**************************
415 * CALCULATE INTERACTIONS *
416 **************************/
418 /* Compute parameters for interactions between i and j atoms */
419 qq10 = _mm_mul_pd(iq1,jq0);
421 /* COULOMB ELECTROSTATICS */
422 velec = _mm_mul_pd(qq10,rinv10);
423 felec = _mm_mul_pd(velec,rinvsq10);
425 /* Update potential sum for this i atom from the interaction with this j atom. */
426 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
427 velecsum = _mm_add_pd(velecsum,velec);
429 fscal = felec;
431 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
433 /* Update vectorial force */
434 fix1 = _mm_macc_pd(dx10,fscal,fix1);
435 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
436 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
438 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
439 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
440 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
442 /**************************
443 * CALCULATE INTERACTIONS *
444 **************************/
446 /* Compute parameters for interactions between i and j atoms */
447 qq20 = _mm_mul_pd(iq2,jq0);
449 /* COULOMB ELECTROSTATICS */
450 velec = _mm_mul_pd(qq20,rinv20);
451 felec = _mm_mul_pd(velec,rinvsq20);
453 /* Update potential sum for this i atom from the interaction with this j atom. */
454 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
455 velecsum = _mm_add_pd(velecsum,velec);
457 fscal = felec;
459 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
461 /* Update vectorial force */
462 fix2 = _mm_macc_pd(dx20,fscal,fix2);
463 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
464 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
466 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
467 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
468 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
470 /**************************
471 * CALCULATE INTERACTIONS *
472 **************************/
474 /* Compute parameters for interactions between i and j atoms */
475 qq30 = _mm_mul_pd(iq3,jq0);
477 /* COULOMB ELECTROSTATICS */
478 velec = _mm_mul_pd(qq30,rinv30);
479 felec = _mm_mul_pd(velec,rinvsq30);
481 /* Update potential sum for this i atom from the interaction with this j atom. */
482 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
483 velecsum = _mm_add_pd(velecsum,velec);
485 fscal = felec;
487 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
489 /* Update vectorial force */
490 fix3 = _mm_macc_pd(dx30,fscal,fix3);
491 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
492 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
494 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
495 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
496 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
498 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
500 /* Inner loop uses 131 flops */
503 /* End of innermost loop */
505 gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
506 f+i_coord_offset,fshift+i_shift_offset);
508 ggid = gid[iidx];
509 /* Update potential energies */
510 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
511 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
513 /* Increment number of inner iterations */
514 inneriter += j_index_end - j_index_start;
516 /* Outer loop uses 26 flops */
519 /* Increment number of outer iterations */
520 outeriter += nri;
522 /* Update outer/inner flops */
524 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*131);
527 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwLJ_GeomW4P1_F_avx_128_fma_double
528 * Electrostatics interaction: Coulomb
529 * VdW interaction: LennardJones
530 * Geometry: Water4-Particle
531 * Calculate force/pot: Force
533 void
534 nb_kernel_ElecCoul_VdwLJ_GeomW4P1_F_avx_128_fma_double
535 (t_nblist * gmx_restrict nlist,
536 rvec * gmx_restrict xx,
537 rvec * gmx_restrict ff,
538 t_forcerec * gmx_restrict fr,
539 t_mdatoms * gmx_restrict mdatoms,
540 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
541 t_nrnb * gmx_restrict nrnb)
543 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
544 * just 0 for non-waters.
545 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
546 * jnr indices corresponding to data put in the four positions in the SIMD register.
548 int i_shift_offset,i_coord_offset,outeriter,inneriter;
549 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
550 int jnrA,jnrB;
551 int j_coord_offsetA,j_coord_offsetB;
552 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
553 real rcutoff_scalar;
554 real *shiftvec,*fshift,*x,*f;
555 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
556 int vdwioffset0;
557 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
558 int vdwioffset1;
559 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
560 int vdwioffset2;
561 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
562 int vdwioffset3;
563 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
564 int vdwjidx0A,vdwjidx0B;
565 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
566 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
567 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
568 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
569 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
570 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
571 real *charge;
572 int nvdwtype;
573 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
574 int *vdwtype;
575 real *vdwparam;
576 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
577 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
578 __m128d dummy_mask,cutoff_mask;
579 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
580 __m128d one = _mm_set1_pd(1.0);
581 __m128d two = _mm_set1_pd(2.0);
582 x = xx[0];
583 f = ff[0];
585 nri = nlist->nri;
586 iinr = nlist->iinr;
587 jindex = nlist->jindex;
588 jjnr = nlist->jjnr;
589 shiftidx = nlist->shift;
590 gid = nlist->gid;
591 shiftvec = fr->shift_vec[0];
592 fshift = fr->fshift[0];
593 facel = _mm_set1_pd(fr->epsfac);
594 charge = mdatoms->chargeA;
595 nvdwtype = fr->ntype;
596 vdwparam = fr->nbfp;
597 vdwtype = mdatoms->typeA;
599 /* Setup water-specific parameters */
600 inr = nlist->iinr[0];
601 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
602 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
603 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
604 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
606 /* Avoid stupid compiler warnings */
607 jnrA = jnrB = 0;
608 j_coord_offsetA = 0;
609 j_coord_offsetB = 0;
611 outeriter = 0;
612 inneriter = 0;
614 /* Start outer loop over neighborlists */
615 for(iidx=0; iidx<nri; iidx++)
617 /* Load shift vector for this list */
618 i_shift_offset = DIM*shiftidx[iidx];
620 /* Load limits for loop over neighbors */
621 j_index_start = jindex[iidx];
622 j_index_end = jindex[iidx+1];
624 /* Get outer coordinate index */
625 inr = iinr[iidx];
626 i_coord_offset = DIM*inr;
628 /* Load i particle coords and add shift vector */
629 gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
630 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
632 fix0 = _mm_setzero_pd();
633 fiy0 = _mm_setzero_pd();
634 fiz0 = _mm_setzero_pd();
635 fix1 = _mm_setzero_pd();
636 fiy1 = _mm_setzero_pd();
637 fiz1 = _mm_setzero_pd();
638 fix2 = _mm_setzero_pd();
639 fiy2 = _mm_setzero_pd();
640 fiz2 = _mm_setzero_pd();
641 fix3 = _mm_setzero_pd();
642 fiy3 = _mm_setzero_pd();
643 fiz3 = _mm_setzero_pd();
645 /* Start inner kernel loop */
646 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
649 /* Get j neighbor index, and coordinate index */
650 jnrA = jjnr[jidx];
651 jnrB = jjnr[jidx+1];
652 j_coord_offsetA = DIM*jnrA;
653 j_coord_offsetB = DIM*jnrB;
655 /* load j atom coordinates */
656 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
657 &jx0,&jy0,&jz0);
659 /* Calculate displacement vector */
660 dx00 = _mm_sub_pd(ix0,jx0);
661 dy00 = _mm_sub_pd(iy0,jy0);
662 dz00 = _mm_sub_pd(iz0,jz0);
663 dx10 = _mm_sub_pd(ix1,jx0);
664 dy10 = _mm_sub_pd(iy1,jy0);
665 dz10 = _mm_sub_pd(iz1,jz0);
666 dx20 = _mm_sub_pd(ix2,jx0);
667 dy20 = _mm_sub_pd(iy2,jy0);
668 dz20 = _mm_sub_pd(iz2,jz0);
669 dx30 = _mm_sub_pd(ix3,jx0);
670 dy30 = _mm_sub_pd(iy3,jy0);
671 dz30 = _mm_sub_pd(iz3,jz0);
673 /* Calculate squared distance and things based on it */
674 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
675 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
676 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
677 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
679 rinv10 = gmx_mm_invsqrt_pd(rsq10);
680 rinv20 = gmx_mm_invsqrt_pd(rsq20);
681 rinv30 = gmx_mm_invsqrt_pd(rsq30);
683 rinvsq00 = gmx_mm_inv_pd(rsq00);
684 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
685 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
686 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
688 /* Load parameters for j particles */
689 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
690 vdwjidx0A = 2*vdwtype[jnrA+0];
691 vdwjidx0B = 2*vdwtype[jnrB+0];
693 fjx0 = _mm_setzero_pd();
694 fjy0 = _mm_setzero_pd();
695 fjz0 = _mm_setzero_pd();
697 /**************************
698 * CALCULATE INTERACTIONS *
699 **************************/
701 /* Compute parameters for interactions between i and j atoms */
702 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
703 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
705 /* LENNARD-JONES DISPERSION/REPULSION */
707 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
708 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
710 fscal = fvdw;
712 /* Update vectorial force */
713 fix0 = _mm_macc_pd(dx00,fscal,fix0);
714 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
715 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
717 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
718 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
719 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
721 /**************************
722 * CALCULATE INTERACTIONS *
723 **************************/
725 /* Compute parameters for interactions between i and j atoms */
726 qq10 = _mm_mul_pd(iq1,jq0);
728 /* COULOMB ELECTROSTATICS */
729 velec = _mm_mul_pd(qq10,rinv10);
730 felec = _mm_mul_pd(velec,rinvsq10);
732 fscal = felec;
734 /* Update vectorial force */
735 fix1 = _mm_macc_pd(dx10,fscal,fix1);
736 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
737 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
739 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
740 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
741 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
743 /**************************
744 * CALCULATE INTERACTIONS *
745 **************************/
747 /* Compute parameters for interactions between i and j atoms */
748 qq20 = _mm_mul_pd(iq2,jq0);
750 /* COULOMB ELECTROSTATICS */
751 velec = _mm_mul_pd(qq20,rinv20);
752 felec = _mm_mul_pd(velec,rinvsq20);
754 fscal = felec;
756 /* Update vectorial force */
757 fix2 = _mm_macc_pd(dx20,fscal,fix2);
758 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
759 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
761 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
762 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
763 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
765 /**************************
766 * CALCULATE INTERACTIONS *
767 **************************/
769 /* Compute parameters for interactions between i and j atoms */
770 qq30 = _mm_mul_pd(iq3,jq0);
772 /* COULOMB ELECTROSTATICS */
773 velec = _mm_mul_pd(qq30,rinv30);
774 felec = _mm_mul_pd(velec,rinvsq30);
776 fscal = felec;
778 /* Update vectorial force */
779 fix3 = _mm_macc_pd(dx30,fscal,fix3);
780 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
781 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
783 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
784 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
785 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
787 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
789 /* Inner loop uses 123 flops */
792 if(jidx<j_index_end)
795 jnrA = jjnr[jidx];
796 j_coord_offsetA = DIM*jnrA;
798 /* load j atom coordinates */
799 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
800 &jx0,&jy0,&jz0);
802 /* Calculate displacement vector */
803 dx00 = _mm_sub_pd(ix0,jx0);
804 dy00 = _mm_sub_pd(iy0,jy0);
805 dz00 = _mm_sub_pd(iz0,jz0);
806 dx10 = _mm_sub_pd(ix1,jx0);
807 dy10 = _mm_sub_pd(iy1,jy0);
808 dz10 = _mm_sub_pd(iz1,jz0);
809 dx20 = _mm_sub_pd(ix2,jx0);
810 dy20 = _mm_sub_pd(iy2,jy0);
811 dz20 = _mm_sub_pd(iz2,jz0);
812 dx30 = _mm_sub_pd(ix3,jx0);
813 dy30 = _mm_sub_pd(iy3,jy0);
814 dz30 = _mm_sub_pd(iz3,jz0);
816 /* Calculate squared distance and things based on it */
817 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
818 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
819 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
820 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
822 rinv10 = gmx_mm_invsqrt_pd(rsq10);
823 rinv20 = gmx_mm_invsqrt_pd(rsq20);
824 rinv30 = gmx_mm_invsqrt_pd(rsq30);
826 rinvsq00 = gmx_mm_inv_pd(rsq00);
827 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
828 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
829 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
831 /* Load parameters for j particles */
832 jq0 = _mm_load_sd(charge+jnrA+0);
833 vdwjidx0A = 2*vdwtype[jnrA+0];
835 fjx0 = _mm_setzero_pd();
836 fjy0 = _mm_setzero_pd();
837 fjz0 = _mm_setzero_pd();
839 /**************************
840 * CALCULATE INTERACTIONS *
841 **************************/
843 /* Compute parameters for interactions between i and j atoms */
844 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
846 /* LENNARD-JONES DISPERSION/REPULSION */
848 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
849 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
851 fscal = fvdw;
853 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
855 /* Update vectorial force */
856 fix0 = _mm_macc_pd(dx00,fscal,fix0);
857 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
858 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
860 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
861 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
862 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
864 /**************************
865 * CALCULATE INTERACTIONS *
866 **************************/
868 /* Compute parameters for interactions between i and j atoms */
869 qq10 = _mm_mul_pd(iq1,jq0);
871 /* COULOMB ELECTROSTATICS */
872 velec = _mm_mul_pd(qq10,rinv10);
873 felec = _mm_mul_pd(velec,rinvsq10);
875 fscal = felec;
877 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
879 /* Update vectorial force */
880 fix1 = _mm_macc_pd(dx10,fscal,fix1);
881 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
882 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
884 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
885 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
886 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
888 /**************************
889 * CALCULATE INTERACTIONS *
890 **************************/
892 /* Compute parameters for interactions between i and j atoms */
893 qq20 = _mm_mul_pd(iq2,jq0);
895 /* COULOMB ELECTROSTATICS */
896 velec = _mm_mul_pd(qq20,rinv20);
897 felec = _mm_mul_pd(velec,rinvsq20);
899 fscal = felec;
901 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
903 /* Update vectorial force */
904 fix2 = _mm_macc_pd(dx20,fscal,fix2);
905 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
906 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
908 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
909 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
910 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
912 /**************************
913 * CALCULATE INTERACTIONS *
914 **************************/
916 /* Compute parameters for interactions between i and j atoms */
917 qq30 = _mm_mul_pd(iq3,jq0);
919 /* COULOMB ELECTROSTATICS */
920 velec = _mm_mul_pd(qq30,rinv30);
921 felec = _mm_mul_pd(velec,rinvsq30);
923 fscal = felec;
925 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
927 /* Update vectorial force */
928 fix3 = _mm_macc_pd(dx30,fscal,fix3);
929 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
930 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
932 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
933 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
934 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
936 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
938 /* Inner loop uses 123 flops */
941 /* End of innermost loop */
943 gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
944 f+i_coord_offset,fshift+i_shift_offset);
946 /* Increment number of inner iterations */
947 inneriter += j_index_end - j_index_start;
949 /* Outer loop uses 24 flops */
952 /* Increment number of outer iterations */
953 outeriter += nri;
955 /* Update outer/inner flops */
957 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*123);