Fix segmentation fault in minimize
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_double / nb_kernel_ElecCoul_VdwNone_GeomW4P1_avx_128_fma_double.cpp
blob64311fb467523bd7e250f3eeda0de2121b797eea
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 avx_128_fma_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_avx_128_fma_double.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwNone_GeomW4P1_VF_avx_128_fma_double
51 * Electrostatics interaction: Coulomb
52 * VdW interaction: None
53 * Geometry: Water4-Particle
54 * Calculate force/pot: PotentialAndForce
56 void
57 nb_kernel_ElecCoul_VdwNone_GeomW4P1_VF_avx_128_fma_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 vdwioffset1;
80 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
81 int vdwioffset2;
82 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
83 int vdwioffset3;
84 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
85 int vdwjidx0A,vdwjidx0B;
86 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
88 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
89 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
90 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
91 real *charge;
92 __m128d dummy_mask,cutoff_mask;
93 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
94 __m128d one = _mm_set1_pd(1.0);
95 __m128d two = _mm_set1_pd(2.0);
96 x = xx[0];
97 f = ff[0];
99 nri = nlist->nri;
100 iinr = nlist->iinr;
101 jindex = nlist->jindex;
102 jjnr = nlist->jjnr;
103 shiftidx = nlist->shift;
104 gid = nlist->gid;
105 shiftvec = fr->shift_vec[0];
106 fshift = fr->fshift[0];
107 facel = _mm_set1_pd(fr->ic->epsfac);
108 charge = mdatoms->chargeA;
110 /* Setup water-specific parameters */
111 inr = nlist->iinr[0];
112 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
113 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
114 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
116 /* Avoid stupid compiler warnings */
117 jnrA = jnrB = 0;
118 j_coord_offsetA = 0;
119 j_coord_offsetB = 0;
121 outeriter = 0;
122 inneriter = 0;
124 /* Start outer loop over neighborlists */
125 for(iidx=0; iidx<nri; iidx++)
127 /* Load shift vector for this list */
128 i_shift_offset = DIM*shiftidx[iidx];
130 /* Load limits for loop over neighbors */
131 j_index_start = jindex[iidx];
132 j_index_end = jindex[iidx+1];
134 /* Get outer coordinate index */
135 inr = iinr[iidx];
136 i_coord_offset = DIM*inr;
138 /* Load i particle coords and add shift vector */
139 gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
140 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
142 fix1 = _mm_setzero_pd();
143 fiy1 = _mm_setzero_pd();
144 fiz1 = _mm_setzero_pd();
145 fix2 = _mm_setzero_pd();
146 fiy2 = _mm_setzero_pd();
147 fiz2 = _mm_setzero_pd();
148 fix3 = _mm_setzero_pd();
149 fiy3 = _mm_setzero_pd();
150 fiz3 = _mm_setzero_pd();
152 /* Reset potential sums */
153 velecsum = _mm_setzero_pd();
155 /* Start inner kernel loop */
156 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
159 /* Get j neighbor index, and coordinate index */
160 jnrA = jjnr[jidx];
161 jnrB = jjnr[jidx+1];
162 j_coord_offsetA = DIM*jnrA;
163 j_coord_offsetB = DIM*jnrB;
165 /* load j atom coordinates */
166 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
167 &jx0,&jy0,&jz0);
169 /* Calculate displacement vector */
170 dx10 = _mm_sub_pd(ix1,jx0);
171 dy10 = _mm_sub_pd(iy1,jy0);
172 dz10 = _mm_sub_pd(iz1,jz0);
173 dx20 = _mm_sub_pd(ix2,jx0);
174 dy20 = _mm_sub_pd(iy2,jy0);
175 dz20 = _mm_sub_pd(iz2,jz0);
176 dx30 = _mm_sub_pd(ix3,jx0);
177 dy30 = _mm_sub_pd(iy3,jy0);
178 dz30 = _mm_sub_pd(iz3,jz0);
180 /* Calculate squared distance and things based on it */
181 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
182 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
183 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
185 rinv10 = avx128fma_invsqrt_d(rsq10);
186 rinv20 = avx128fma_invsqrt_d(rsq20);
187 rinv30 = avx128fma_invsqrt_d(rsq30);
189 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
190 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
191 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
193 /* Load parameters for j particles */
194 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
196 fjx0 = _mm_setzero_pd();
197 fjy0 = _mm_setzero_pd();
198 fjz0 = _mm_setzero_pd();
200 /**************************
201 * CALCULATE INTERACTIONS *
202 **************************/
204 /* Compute parameters for interactions between i and j atoms */
205 qq10 = _mm_mul_pd(iq1,jq0);
207 /* COULOMB ELECTROSTATICS */
208 velec = _mm_mul_pd(qq10,rinv10);
209 felec = _mm_mul_pd(velec,rinvsq10);
211 /* Update potential sum for this i atom from the interaction with this j atom. */
212 velecsum = _mm_add_pd(velecsum,velec);
214 fscal = felec;
216 /* Update vectorial force */
217 fix1 = _mm_macc_pd(dx10,fscal,fix1);
218 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
219 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
221 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
222 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
223 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
225 /**************************
226 * CALCULATE INTERACTIONS *
227 **************************/
229 /* Compute parameters for interactions between i and j atoms */
230 qq20 = _mm_mul_pd(iq2,jq0);
232 /* COULOMB ELECTROSTATICS */
233 velec = _mm_mul_pd(qq20,rinv20);
234 felec = _mm_mul_pd(velec,rinvsq20);
236 /* Update potential sum for this i atom from the interaction with this j atom. */
237 velecsum = _mm_add_pd(velecsum,velec);
239 fscal = felec;
241 /* Update vectorial force */
242 fix2 = _mm_macc_pd(dx20,fscal,fix2);
243 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
244 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
246 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
247 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
248 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
250 /**************************
251 * CALCULATE INTERACTIONS *
252 **************************/
254 /* Compute parameters for interactions between i and j atoms */
255 qq30 = _mm_mul_pd(iq3,jq0);
257 /* COULOMB ELECTROSTATICS */
258 velec = _mm_mul_pd(qq30,rinv30);
259 felec = _mm_mul_pd(velec,rinvsq30);
261 /* Update potential sum for this i atom from the interaction with this j atom. */
262 velecsum = _mm_add_pd(velecsum,velec);
264 fscal = felec;
266 /* Update vectorial force */
267 fix3 = _mm_macc_pd(dx30,fscal,fix3);
268 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
269 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
271 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
272 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
273 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
275 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
277 /* Inner loop uses 96 flops */
280 if(jidx<j_index_end)
283 jnrA = jjnr[jidx];
284 j_coord_offsetA = DIM*jnrA;
286 /* load j atom coordinates */
287 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
288 &jx0,&jy0,&jz0);
290 /* Calculate displacement vector */
291 dx10 = _mm_sub_pd(ix1,jx0);
292 dy10 = _mm_sub_pd(iy1,jy0);
293 dz10 = _mm_sub_pd(iz1,jz0);
294 dx20 = _mm_sub_pd(ix2,jx0);
295 dy20 = _mm_sub_pd(iy2,jy0);
296 dz20 = _mm_sub_pd(iz2,jz0);
297 dx30 = _mm_sub_pd(ix3,jx0);
298 dy30 = _mm_sub_pd(iy3,jy0);
299 dz30 = _mm_sub_pd(iz3,jz0);
301 /* Calculate squared distance and things based on it */
302 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
303 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
304 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
306 rinv10 = avx128fma_invsqrt_d(rsq10);
307 rinv20 = avx128fma_invsqrt_d(rsq20);
308 rinv30 = avx128fma_invsqrt_d(rsq30);
310 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
311 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
312 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
314 /* Load parameters for j particles */
315 jq0 = _mm_load_sd(charge+jnrA+0);
317 fjx0 = _mm_setzero_pd();
318 fjy0 = _mm_setzero_pd();
319 fjz0 = _mm_setzero_pd();
321 /**************************
322 * CALCULATE INTERACTIONS *
323 **************************/
325 /* Compute parameters for interactions between i and j atoms */
326 qq10 = _mm_mul_pd(iq1,jq0);
328 /* COULOMB ELECTROSTATICS */
329 velec = _mm_mul_pd(qq10,rinv10);
330 felec = _mm_mul_pd(velec,rinvsq10);
332 /* Update potential sum for this i atom from the interaction with this j atom. */
333 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
334 velecsum = _mm_add_pd(velecsum,velec);
336 fscal = felec;
338 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
340 /* Update vectorial force */
341 fix1 = _mm_macc_pd(dx10,fscal,fix1);
342 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
343 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
345 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
346 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
347 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
349 /**************************
350 * CALCULATE INTERACTIONS *
351 **************************/
353 /* Compute parameters for interactions between i and j atoms */
354 qq20 = _mm_mul_pd(iq2,jq0);
356 /* COULOMB ELECTROSTATICS */
357 velec = _mm_mul_pd(qq20,rinv20);
358 felec = _mm_mul_pd(velec,rinvsq20);
360 /* Update potential sum for this i atom from the interaction with this j atom. */
361 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
362 velecsum = _mm_add_pd(velecsum,velec);
364 fscal = felec;
366 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
368 /* Update vectorial force */
369 fix2 = _mm_macc_pd(dx20,fscal,fix2);
370 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
371 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
373 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
374 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
375 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
377 /**************************
378 * CALCULATE INTERACTIONS *
379 **************************/
381 /* Compute parameters for interactions between i and j atoms */
382 qq30 = _mm_mul_pd(iq3,jq0);
384 /* COULOMB ELECTROSTATICS */
385 velec = _mm_mul_pd(qq30,rinv30);
386 felec = _mm_mul_pd(velec,rinvsq30);
388 /* Update potential sum for this i atom from the interaction with this j atom. */
389 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
390 velecsum = _mm_add_pd(velecsum,velec);
392 fscal = felec;
394 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
396 /* Update vectorial force */
397 fix3 = _mm_macc_pd(dx30,fscal,fix3);
398 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
399 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
401 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
402 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
403 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
405 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
407 /* Inner loop uses 96 flops */
410 /* End of innermost loop */
412 gmx_mm_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
413 f+i_coord_offset+DIM,fshift+i_shift_offset);
415 ggid = gid[iidx];
416 /* Update potential energies */
417 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
419 /* Increment number of inner iterations */
420 inneriter += j_index_end - j_index_start;
422 /* Outer loop uses 19 flops */
425 /* Increment number of outer iterations */
426 outeriter += nri;
428 /* Update outer/inner flops */
430 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*19 + inneriter*96);
433 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_avx_128_fma_double
434 * Electrostatics interaction: Coulomb
435 * VdW interaction: None
436 * Geometry: Water4-Particle
437 * Calculate force/pot: Force
439 void
440 nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_avx_128_fma_double
441 (t_nblist * gmx_restrict nlist,
442 rvec * gmx_restrict xx,
443 rvec * gmx_restrict ff,
444 struct t_forcerec * gmx_restrict fr,
445 t_mdatoms * gmx_restrict mdatoms,
446 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
447 t_nrnb * gmx_restrict nrnb)
449 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
450 * just 0 for non-waters.
451 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
452 * jnr indices corresponding to data put in the four positions in the SIMD register.
454 int i_shift_offset,i_coord_offset,outeriter,inneriter;
455 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
456 int jnrA,jnrB;
457 int j_coord_offsetA,j_coord_offsetB;
458 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
459 real rcutoff_scalar;
460 real *shiftvec,*fshift,*x,*f;
461 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
462 int vdwioffset1;
463 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
464 int vdwioffset2;
465 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
466 int vdwioffset3;
467 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
468 int vdwjidx0A,vdwjidx0B;
469 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
470 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
471 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
472 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
473 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
474 real *charge;
475 __m128d dummy_mask,cutoff_mask;
476 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
477 __m128d one = _mm_set1_pd(1.0);
478 __m128d two = _mm_set1_pd(2.0);
479 x = xx[0];
480 f = ff[0];
482 nri = nlist->nri;
483 iinr = nlist->iinr;
484 jindex = nlist->jindex;
485 jjnr = nlist->jjnr;
486 shiftidx = nlist->shift;
487 gid = nlist->gid;
488 shiftvec = fr->shift_vec[0];
489 fshift = fr->fshift[0];
490 facel = _mm_set1_pd(fr->ic->epsfac);
491 charge = mdatoms->chargeA;
493 /* Setup water-specific parameters */
494 inr = nlist->iinr[0];
495 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
496 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
497 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
499 /* Avoid stupid compiler warnings */
500 jnrA = jnrB = 0;
501 j_coord_offsetA = 0;
502 j_coord_offsetB = 0;
504 outeriter = 0;
505 inneriter = 0;
507 /* Start outer loop over neighborlists */
508 for(iidx=0; iidx<nri; iidx++)
510 /* Load shift vector for this list */
511 i_shift_offset = DIM*shiftidx[iidx];
513 /* Load limits for loop over neighbors */
514 j_index_start = jindex[iidx];
515 j_index_end = jindex[iidx+1];
517 /* Get outer coordinate index */
518 inr = iinr[iidx];
519 i_coord_offset = DIM*inr;
521 /* Load i particle coords and add shift vector */
522 gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
523 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
525 fix1 = _mm_setzero_pd();
526 fiy1 = _mm_setzero_pd();
527 fiz1 = _mm_setzero_pd();
528 fix2 = _mm_setzero_pd();
529 fiy2 = _mm_setzero_pd();
530 fiz2 = _mm_setzero_pd();
531 fix3 = _mm_setzero_pd();
532 fiy3 = _mm_setzero_pd();
533 fiz3 = _mm_setzero_pd();
535 /* Start inner kernel loop */
536 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
539 /* Get j neighbor index, and coordinate index */
540 jnrA = jjnr[jidx];
541 jnrB = jjnr[jidx+1];
542 j_coord_offsetA = DIM*jnrA;
543 j_coord_offsetB = DIM*jnrB;
545 /* load j atom coordinates */
546 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
547 &jx0,&jy0,&jz0);
549 /* Calculate displacement vector */
550 dx10 = _mm_sub_pd(ix1,jx0);
551 dy10 = _mm_sub_pd(iy1,jy0);
552 dz10 = _mm_sub_pd(iz1,jz0);
553 dx20 = _mm_sub_pd(ix2,jx0);
554 dy20 = _mm_sub_pd(iy2,jy0);
555 dz20 = _mm_sub_pd(iz2,jz0);
556 dx30 = _mm_sub_pd(ix3,jx0);
557 dy30 = _mm_sub_pd(iy3,jy0);
558 dz30 = _mm_sub_pd(iz3,jz0);
560 /* Calculate squared distance and things based on it */
561 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
562 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
563 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
565 rinv10 = avx128fma_invsqrt_d(rsq10);
566 rinv20 = avx128fma_invsqrt_d(rsq20);
567 rinv30 = avx128fma_invsqrt_d(rsq30);
569 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
570 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
571 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
573 /* Load parameters for j particles */
574 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
576 fjx0 = _mm_setzero_pd();
577 fjy0 = _mm_setzero_pd();
578 fjz0 = _mm_setzero_pd();
580 /**************************
581 * CALCULATE INTERACTIONS *
582 **************************/
584 /* Compute parameters for interactions between i and j atoms */
585 qq10 = _mm_mul_pd(iq1,jq0);
587 /* COULOMB ELECTROSTATICS */
588 velec = _mm_mul_pd(qq10,rinv10);
589 felec = _mm_mul_pd(velec,rinvsq10);
591 fscal = felec;
593 /* Update vectorial force */
594 fix1 = _mm_macc_pd(dx10,fscal,fix1);
595 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
596 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
598 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
599 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
600 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
602 /**************************
603 * CALCULATE INTERACTIONS *
604 **************************/
606 /* Compute parameters for interactions between i and j atoms */
607 qq20 = _mm_mul_pd(iq2,jq0);
609 /* COULOMB ELECTROSTATICS */
610 velec = _mm_mul_pd(qq20,rinv20);
611 felec = _mm_mul_pd(velec,rinvsq20);
613 fscal = felec;
615 /* Update vectorial force */
616 fix2 = _mm_macc_pd(dx20,fscal,fix2);
617 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
618 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
620 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
621 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
622 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
624 /**************************
625 * CALCULATE INTERACTIONS *
626 **************************/
628 /* Compute parameters for interactions between i and j atoms */
629 qq30 = _mm_mul_pd(iq3,jq0);
631 /* COULOMB ELECTROSTATICS */
632 velec = _mm_mul_pd(qq30,rinv30);
633 felec = _mm_mul_pd(velec,rinvsq30);
635 fscal = felec;
637 /* Update vectorial force */
638 fix3 = _mm_macc_pd(dx30,fscal,fix3);
639 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
640 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
642 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
643 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
644 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
646 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
648 /* Inner loop uses 93 flops */
651 if(jidx<j_index_end)
654 jnrA = jjnr[jidx];
655 j_coord_offsetA = DIM*jnrA;
657 /* load j atom coordinates */
658 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
659 &jx0,&jy0,&jz0);
661 /* Calculate displacement vector */
662 dx10 = _mm_sub_pd(ix1,jx0);
663 dy10 = _mm_sub_pd(iy1,jy0);
664 dz10 = _mm_sub_pd(iz1,jz0);
665 dx20 = _mm_sub_pd(ix2,jx0);
666 dy20 = _mm_sub_pd(iy2,jy0);
667 dz20 = _mm_sub_pd(iz2,jz0);
668 dx30 = _mm_sub_pd(ix3,jx0);
669 dy30 = _mm_sub_pd(iy3,jy0);
670 dz30 = _mm_sub_pd(iz3,jz0);
672 /* Calculate squared distance and things based on it */
673 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
674 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
675 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
677 rinv10 = avx128fma_invsqrt_d(rsq10);
678 rinv20 = avx128fma_invsqrt_d(rsq20);
679 rinv30 = avx128fma_invsqrt_d(rsq30);
681 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
682 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
683 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
685 /* Load parameters for j particles */
686 jq0 = _mm_load_sd(charge+jnrA+0);
688 fjx0 = _mm_setzero_pd();
689 fjy0 = _mm_setzero_pd();
690 fjz0 = _mm_setzero_pd();
692 /**************************
693 * CALCULATE INTERACTIONS *
694 **************************/
696 /* Compute parameters for interactions between i and j atoms */
697 qq10 = _mm_mul_pd(iq1,jq0);
699 /* COULOMB ELECTROSTATICS */
700 velec = _mm_mul_pd(qq10,rinv10);
701 felec = _mm_mul_pd(velec,rinvsq10);
703 fscal = felec;
705 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
707 /* Update vectorial force */
708 fix1 = _mm_macc_pd(dx10,fscal,fix1);
709 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
710 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
712 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
713 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
714 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
716 /**************************
717 * CALCULATE INTERACTIONS *
718 **************************/
720 /* Compute parameters for interactions between i and j atoms */
721 qq20 = _mm_mul_pd(iq2,jq0);
723 /* COULOMB ELECTROSTATICS */
724 velec = _mm_mul_pd(qq20,rinv20);
725 felec = _mm_mul_pd(velec,rinvsq20);
727 fscal = felec;
729 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
731 /* Update vectorial force */
732 fix2 = _mm_macc_pd(dx20,fscal,fix2);
733 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
734 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
736 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
737 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
738 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
740 /**************************
741 * CALCULATE INTERACTIONS *
742 **************************/
744 /* Compute parameters for interactions between i and j atoms */
745 qq30 = _mm_mul_pd(iq3,jq0);
747 /* COULOMB ELECTROSTATICS */
748 velec = _mm_mul_pd(qq30,rinv30);
749 felec = _mm_mul_pd(velec,rinvsq30);
751 fscal = felec;
753 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
755 /* Update vectorial force */
756 fix3 = _mm_macc_pd(dx30,fscal,fix3);
757 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
758 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
760 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
761 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
762 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
764 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
766 /* Inner loop uses 93 flops */
769 /* End of innermost loop */
771 gmx_mm_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
772 f+i_coord_offset+DIM,fshift+i_shift_offset);
774 /* Increment number of inner iterations */
775 inneriter += j_index_end - j_index_start;
777 /* Outer loop uses 18 flops */
780 /* Increment number of outer iterations */
781 outeriter += nri;
783 /* Update outer/inner flops */
785 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*18 + inneriter*93);