Removed simple.h from nb_kernel_sse2_XX
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_double / nb_kernel_ElecRF_VdwCSTab_GeomW3P1_sse2_double.c
blob8082cf5ffef3b4c343a91900761f23e86fc4159c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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 sse2_double kernel generator.
38 #include "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/legacyheaders/nrnb.h"
48 #include "gromacs/simd/math_x86_sse2_double.h"
49 #include "kernelutil_x86_sse2_double.h"
52 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomW3P1_VF_sse2_double
53 * Electrostatics interaction: ReactionField
54 * VdW interaction: CubicSplineTable
55 * Geometry: Water3-Particle
56 * Calculate force/pot: PotentialAndForce
58 void
59 nb_kernel_ElecRF_VdwCSTab_GeomW3P1_VF_sse2_double
60 (t_nblist * gmx_restrict nlist,
61 rvec * gmx_restrict xx,
62 rvec * gmx_restrict ff,
63 t_forcerec * gmx_restrict fr,
64 t_mdatoms * gmx_restrict mdatoms,
65 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
66 t_nrnb * gmx_restrict nrnb)
68 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
69 * just 0 for non-waters.
70 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
71 * jnr indices corresponding to data put in the four positions in the SIMD register.
73 int i_shift_offset,i_coord_offset,outeriter,inneriter;
74 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
75 int jnrA,jnrB;
76 int j_coord_offsetA,j_coord_offsetB;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
78 real rcutoff_scalar;
79 real *shiftvec,*fshift,*x,*f;
80 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
81 int vdwioffset0;
82 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
83 int vdwioffset1;
84 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
85 int vdwioffset2;
86 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
87 int vdwjidx0A,vdwjidx0B;
88 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
91 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
92 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
93 real *charge;
94 int nvdwtype;
95 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
96 int *vdwtype;
97 real *vdwparam;
98 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
99 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
100 __m128i vfitab;
101 __m128i ifour = _mm_set1_epi32(4);
102 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
103 real *vftab;
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);
108 x = xx[0];
109 f = ff[0];
111 nri = nlist->nri;
112 iinr = nlist->iinr;
113 jindex = nlist->jindex;
114 jjnr = nlist->jjnr;
115 shiftidx = nlist->shift;
116 gid = nlist->gid;
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;
125 vdwparam = fr->nbfp;
126 vdwtype = mdatoms->typeA;
128 vftab = kernel_data->table_vdw->data;
129 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
131 /* Setup water-specific parameters */
132 inr = nlist->iinr[0];
133 iq0 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+0]));
134 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
135 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
136 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
138 /* Avoid stupid compiler warnings */
139 jnrA = jnrB = 0;
140 j_coord_offsetA = 0;
141 j_coord_offsetB = 0;
143 outeriter = 0;
144 inneriter = 0;
146 /* Start outer loop over neighborlists */
147 for(iidx=0; iidx<nri; iidx++)
149 /* Load shift vector for this list */
150 i_shift_offset = DIM*shiftidx[iidx];
152 /* Load limits for loop over neighbors */
153 j_index_start = jindex[iidx];
154 j_index_end = jindex[iidx+1];
156 /* Get outer coordinate index */
157 inr = iinr[iidx];
158 i_coord_offset = DIM*inr;
160 /* Load i particle coords and add shift vector */
161 gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
162 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
164 fix0 = _mm_setzero_pd();
165 fiy0 = _mm_setzero_pd();
166 fiz0 = _mm_setzero_pd();
167 fix1 = _mm_setzero_pd();
168 fiy1 = _mm_setzero_pd();
169 fiz1 = _mm_setzero_pd();
170 fix2 = _mm_setzero_pd();
171 fiy2 = _mm_setzero_pd();
172 fiz2 = _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 */
183 jnrA = jjnr[jidx];
184 jnrB = jjnr[jidx+1];
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,
190 &jx0,&jy0,&jz0);
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);
203 /* Calculate squared distance and things based on it */
204 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
205 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
206 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
208 rinv00 = gmx_mm_invsqrt_pd(rsq00);
209 rinv10 = gmx_mm_invsqrt_pd(rsq10);
210 rinv20 = gmx_mm_invsqrt_pd(rsq20);
212 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
213 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
214 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
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 r00 = _mm_mul_pd(rsq00,rinv00);
231 /* Compute parameters for interactions between i and j atoms */
232 qq00 = _mm_mul_pd(iq0,jq0);
233 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
234 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
236 /* Calculate table index by multiplying r with table scale and truncate to integer */
237 rt = _mm_mul_pd(r00,vftabscale);
238 vfitab = _mm_cvttpd_epi32(rt);
239 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
240 vfitab = _mm_slli_epi32(vfitab,3);
242 /* REACTION-FIELD ELECTROSTATICS */
243 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
244 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
246 /* CUBIC SPLINE TABLE DISPERSION */
247 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
248 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
249 GMX_MM_TRANSPOSE2_PD(Y,F);
250 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
251 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
252 GMX_MM_TRANSPOSE2_PD(G,H);
253 Heps = _mm_mul_pd(vfeps,H);
254 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
255 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
256 vvdw6 = _mm_mul_pd(c6_00,VV);
257 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
258 fvdw6 = _mm_mul_pd(c6_00,FF);
260 /* CUBIC SPLINE TABLE REPULSION */
261 vfitab = _mm_add_epi32(vfitab,ifour);
262 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
263 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
264 GMX_MM_TRANSPOSE2_PD(Y,F);
265 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
266 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
267 GMX_MM_TRANSPOSE2_PD(G,H);
268 Heps = _mm_mul_pd(vfeps,H);
269 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
270 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
271 vvdw12 = _mm_mul_pd(c12_00,VV);
272 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
273 fvdw12 = _mm_mul_pd(c12_00,FF);
274 vvdw = _mm_add_pd(vvdw12,vvdw6);
275 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
277 /* Update potential sum for this i atom from the interaction with this j atom. */
278 velecsum = _mm_add_pd(velecsum,velec);
279 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
281 fscal = _mm_add_pd(felec,fvdw);
283 /* Calculate temporary vectorial force */
284 tx = _mm_mul_pd(fscal,dx00);
285 ty = _mm_mul_pd(fscal,dy00);
286 tz = _mm_mul_pd(fscal,dz00);
288 /* Update vectorial force */
289 fix0 = _mm_add_pd(fix0,tx);
290 fiy0 = _mm_add_pd(fiy0,ty);
291 fiz0 = _mm_add_pd(fiz0,tz);
293 fjx0 = _mm_add_pd(fjx0,tx);
294 fjy0 = _mm_add_pd(fjy0,ty);
295 fjz0 = _mm_add_pd(fjz0,tz);
297 /**************************
298 * CALCULATE INTERACTIONS *
299 **************************/
301 /* Compute parameters for interactions between i and j atoms */
302 qq10 = _mm_mul_pd(iq1,jq0);
304 /* REACTION-FIELD ELECTROSTATICS */
305 velec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_add_pd(rinv10,_mm_mul_pd(krf,rsq10)),crf));
306 felec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
308 /* Update potential sum for this i atom from the interaction with this j atom. */
309 velecsum = _mm_add_pd(velecsum,velec);
311 fscal = felec;
313 /* Calculate temporary vectorial force */
314 tx = _mm_mul_pd(fscal,dx10);
315 ty = _mm_mul_pd(fscal,dy10);
316 tz = _mm_mul_pd(fscal,dz10);
318 /* Update vectorial force */
319 fix1 = _mm_add_pd(fix1,tx);
320 fiy1 = _mm_add_pd(fiy1,ty);
321 fiz1 = _mm_add_pd(fiz1,tz);
323 fjx0 = _mm_add_pd(fjx0,tx);
324 fjy0 = _mm_add_pd(fjy0,ty);
325 fjz0 = _mm_add_pd(fjz0,tz);
327 /**************************
328 * CALCULATE INTERACTIONS *
329 **************************/
331 /* Compute parameters for interactions between i and j atoms */
332 qq20 = _mm_mul_pd(iq2,jq0);
334 /* REACTION-FIELD ELECTROSTATICS */
335 velec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_add_pd(rinv20,_mm_mul_pd(krf,rsq20)),crf));
336 felec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
338 /* Update potential sum for this i atom from the interaction with this j atom. */
339 velecsum = _mm_add_pd(velecsum,velec);
341 fscal = felec;
343 /* Calculate temporary vectorial force */
344 tx = _mm_mul_pd(fscal,dx20);
345 ty = _mm_mul_pd(fscal,dy20);
346 tz = _mm_mul_pd(fscal,dz20);
348 /* Update vectorial force */
349 fix2 = _mm_add_pd(fix2,tx);
350 fiy2 = _mm_add_pd(fiy2,ty);
351 fiz2 = _mm_add_pd(fiz2,tz);
353 fjx0 = _mm_add_pd(fjx0,tx);
354 fjy0 = _mm_add_pd(fjy0,ty);
355 fjz0 = _mm_add_pd(fjz0,tz);
357 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
359 /* Inner loop uses 134 flops */
362 if(jidx<j_index_end)
365 jnrA = jjnr[jidx];
366 j_coord_offsetA = DIM*jnrA;
368 /* load j atom coordinates */
369 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
370 &jx0,&jy0,&jz0);
372 /* Calculate displacement vector */
373 dx00 = _mm_sub_pd(ix0,jx0);
374 dy00 = _mm_sub_pd(iy0,jy0);
375 dz00 = _mm_sub_pd(iz0,jz0);
376 dx10 = _mm_sub_pd(ix1,jx0);
377 dy10 = _mm_sub_pd(iy1,jy0);
378 dz10 = _mm_sub_pd(iz1,jz0);
379 dx20 = _mm_sub_pd(ix2,jx0);
380 dy20 = _mm_sub_pd(iy2,jy0);
381 dz20 = _mm_sub_pd(iz2,jz0);
383 /* Calculate squared distance and things based on it */
384 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
385 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
386 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
388 rinv00 = gmx_mm_invsqrt_pd(rsq00);
389 rinv10 = gmx_mm_invsqrt_pd(rsq10);
390 rinv20 = gmx_mm_invsqrt_pd(rsq20);
392 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
393 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
394 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
396 /* Load parameters for j particles */
397 jq0 = _mm_load_sd(charge+jnrA+0);
398 vdwjidx0A = 2*vdwtype[jnrA+0];
400 fjx0 = _mm_setzero_pd();
401 fjy0 = _mm_setzero_pd();
402 fjz0 = _mm_setzero_pd();
404 /**************************
405 * CALCULATE INTERACTIONS *
406 **************************/
408 r00 = _mm_mul_pd(rsq00,rinv00);
410 /* Compute parameters for interactions between i and j atoms */
411 qq00 = _mm_mul_pd(iq0,jq0);
412 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
414 /* Calculate table index by multiplying r with table scale and truncate to integer */
415 rt = _mm_mul_pd(r00,vftabscale);
416 vfitab = _mm_cvttpd_epi32(rt);
417 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
418 vfitab = _mm_slli_epi32(vfitab,3);
420 /* REACTION-FIELD ELECTROSTATICS */
421 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
422 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
424 /* CUBIC SPLINE TABLE DISPERSION */
425 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
426 F = _mm_setzero_pd();
427 GMX_MM_TRANSPOSE2_PD(Y,F);
428 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
429 H = _mm_setzero_pd();
430 GMX_MM_TRANSPOSE2_PD(G,H);
431 Heps = _mm_mul_pd(vfeps,H);
432 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
433 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
434 vvdw6 = _mm_mul_pd(c6_00,VV);
435 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
436 fvdw6 = _mm_mul_pd(c6_00,FF);
438 /* CUBIC SPLINE TABLE REPULSION */
439 vfitab = _mm_add_epi32(vfitab,ifour);
440 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
441 F = _mm_setzero_pd();
442 GMX_MM_TRANSPOSE2_PD(Y,F);
443 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
444 H = _mm_setzero_pd();
445 GMX_MM_TRANSPOSE2_PD(G,H);
446 Heps = _mm_mul_pd(vfeps,H);
447 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
448 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
449 vvdw12 = _mm_mul_pd(c12_00,VV);
450 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
451 fvdw12 = _mm_mul_pd(c12_00,FF);
452 vvdw = _mm_add_pd(vvdw12,vvdw6);
453 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
455 /* Update potential sum for this i atom from the interaction with this j atom. */
456 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
457 velecsum = _mm_add_pd(velecsum,velec);
458 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
459 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
461 fscal = _mm_add_pd(felec,fvdw);
463 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
465 /* Calculate temporary vectorial force */
466 tx = _mm_mul_pd(fscal,dx00);
467 ty = _mm_mul_pd(fscal,dy00);
468 tz = _mm_mul_pd(fscal,dz00);
470 /* Update vectorial force */
471 fix0 = _mm_add_pd(fix0,tx);
472 fiy0 = _mm_add_pd(fiy0,ty);
473 fiz0 = _mm_add_pd(fiz0,tz);
475 fjx0 = _mm_add_pd(fjx0,tx);
476 fjy0 = _mm_add_pd(fjy0,ty);
477 fjz0 = _mm_add_pd(fjz0,tz);
479 /**************************
480 * CALCULATE INTERACTIONS *
481 **************************/
483 /* Compute parameters for interactions between i and j atoms */
484 qq10 = _mm_mul_pd(iq1,jq0);
486 /* REACTION-FIELD ELECTROSTATICS */
487 velec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_add_pd(rinv10,_mm_mul_pd(krf,rsq10)),crf));
488 felec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
490 /* Update potential sum for this i atom from the interaction with this j atom. */
491 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
492 velecsum = _mm_add_pd(velecsum,velec);
494 fscal = felec;
496 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
498 /* Calculate temporary vectorial force */
499 tx = _mm_mul_pd(fscal,dx10);
500 ty = _mm_mul_pd(fscal,dy10);
501 tz = _mm_mul_pd(fscal,dz10);
503 /* Update vectorial force */
504 fix1 = _mm_add_pd(fix1,tx);
505 fiy1 = _mm_add_pd(fiy1,ty);
506 fiz1 = _mm_add_pd(fiz1,tz);
508 fjx0 = _mm_add_pd(fjx0,tx);
509 fjy0 = _mm_add_pd(fjy0,ty);
510 fjz0 = _mm_add_pd(fjz0,tz);
512 /**************************
513 * CALCULATE INTERACTIONS *
514 **************************/
516 /* Compute parameters for interactions between i and j atoms */
517 qq20 = _mm_mul_pd(iq2,jq0);
519 /* REACTION-FIELD ELECTROSTATICS */
520 velec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_add_pd(rinv20,_mm_mul_pd(krf,rsq20)),crf));
521 felec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
523 /* Update potential sum for this i atom from the interaction with this j atom. */
524 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
525 velecsum = _mm_add_pd(velecsum,velec);
527 fscal = felec;
529 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
531 /* Calculate temporary vectorial force */
532 tx = _mm_mul_pd(fscal,dx20);
533 ty = _mm_mul_pd(fscal,dy20);
534 tz = _mm_mul_pd(fscal,dz20);
536 /* Update vectorial force */
537 fix2 = _mm_add_pd(fix2,tx);
538 fiy2 = _mm_add_pd(fiy2,ty);
539 fiz2 = _mm_add_pd(fiz2,tz);
541 fjx0 = _mm_add_pd(fjx0,tx);
542 fjy0 = _mm_add_pd(fjy0,ty);
543 fjz0 = _mm_add_pd(fjz0,tz);
545 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
547 /* Inner loop uses 134 flops */
550 /* End of innermost loop */
552 gmx_mm_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
553 f+i_coord_offset,fshift+i_shift_offset);
555 ggid = gid[iidx];
556 /* Update potential energies */
557 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
558 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
560 /* Increment number of inner iterations */
561 inneriter += j_index_end - j_index_start;
563 /* Outer loop uses 20 flops */
566 /* Increment number of outer iterations */
567 outeriter += nri;
569 /* Update outer/inner flops */
571 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*20 + inneriter*134);
574 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomW3P1_F_sse2_double
575 * Electrostatics interaction: ReactionField
576 * VdW interaction: CubicSplineTable
577 * Geometry: Water3-Particle
578 * Calculate force/pot: Force
580 void
581 nb_kernel_ElecRF_VdwCSTab_GeomW3P1_F_sse2_double
582 (t_nblist * gmx_restrict nlist,
583 rvec * gmx_restrict xx,
584 rvec * gmx_restrict ff,
585 t_forcerec * gmx_restrict fr,
586 t_mdatoms * gmx_restrict mdatoms,
587 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
588 t_nrnb * gmx_restrict nrnb)
590 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
591 * just 0 for non-waters.
592 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
593 * jnr indices corresponding to data put in the four positions in the SIMD register.
595 int i_shift_offset,i_coord_offset,outeriter,inneriter;
596 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
597 int jnrA,jnrB;
598 int j_coord_offsetA,j_coord_offsetB;
599 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
600 real rcutoff_scalar;
601 real *shiftvec,*fshift,*x,*f;
602 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
603 int vdwioffset0;
604 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
605 int vdwioffset1;
606 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
607 int vdwioffset2;
608 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
609 int vdwjidx0A,vdwjidx0B;
610 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
611 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
612 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
613 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
614 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
615 real *charge;
616 int nvdwtype;
617 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
618 int *vdwtype;
619 real *vdwparam;
620 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
621 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
622 __m128i vfitab;
623 __m128i ifour = _mm_set1_epi32(4);
624 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
625 real *vftab;
626 __m128d dummy_mask,cutoff_mask;
627 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
628 __m128d one = _mm_set1_pd(1.0);
629 __m128d two = _mm_set1_pd(2.0);
630 x = xx[0];
631 f = ff[0];
633 nri = nlist->nri;
634 iinr = nlist->iinr;
635 jindex = nlist->jindex;
636 jjnr = nlist->jjnr;
637 shiftidx = nlist->shift;
638 gid = nlist->gid;
639 shiftvec = fr->shift_vec[0];
640 fshift = fr->fshift[0];
641 facel = _mm_set1_pd(fr->epsfac);
642 charge = mdatoms->chargeA;
643 krf = _mm_set1_pd(fr->ic->k_rf);
644 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
645 crf = _mm_set1_pd(fr->ic->c_rf);
646 nvdwtype = fr->ntype;
647 vdwparam = fr->nbfp;
648 vdwtype = mdatoms->typeA;
650 vftab = kernel_data->table_vdw->data;
651 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
653 /* Setup water-specific parameters */
654 inr = nlist->iinr[0];
655 iq0 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+0]));
656 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
657 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
658 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
660 /* Avoid stupid compiler warnings */
661 jnrA = jnrB = 0;
662 j_coord_offsetA = 0;
663 j_coord_offsetB = 0;
665 outeriter = 0;
666 inneriter = 0;
668 /* Start outer loop over neighborlists */
669 for(iidx=0; iidx<nri; iidx++)
671 /* Load shift vector for this list */
672 i_shift_offset = DIM*shiftidx[iidx];
674 /* Load limits for loop over neighbors */
675 j_index_start = jindex[iidx];
676 j_index_end = jindex[iidx+1];
678 /* Get outer coordinate index */
679 inr = iinr[iidx];
680 i_coord_offset = DIM*inr;
682 /* Load i particle coords and add shift vector */
683 gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
684 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
686 fix0 = _mm_setzero_pd();
687 fiy0 = _mm_setzero_pd();
688 fiz0 = _mm_setzero_pd();
689 fix1 = _mm_setzero_pd();
690 fiy1 = _mm_setzero_pd();
691 fiz1 = _mm_setzero_pd();
692 fix2 = _mm_setzero_pd();
693 fiy2 = _mm_setzero_pd();
694 fiz2 = _mm_setzero_pd();
696 /* Start inner kernel loop */
697 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
700 /* Get j neighbor index, and coordinate index */
701 jnrA = jjnr[jidx];
702 jnrB = jjnr[jidx+1];
703 j_coord_offsetA = DIM*jnrA;
704 j_coord_offsetB = DIM*jnrB;
706 /* load j atom coordinates */
707 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
708 &jx0,&jy0,&jz0);
710 /* Calculate displacement vector */
711 dx00 = _mm_sub_pd(ix0,jx0);
712 dy00 = _mm_sub_pd(iy0,jy0);
713 dz00 = _mm_sub_pd(iz0,jz0);
714 dx10 = _mm_sub_pd(ix1,jx0);
715 dy10 = _mm_sub_pd(iy1,jy0);
716 dz10 = _mm_sub_pd(iz1,jz0);
717 dx20 = _mm_sub_pd(ix2,jx0);
718 dy20 = _mm_sub_pd(iy2,jy0);
719 dz20 = _mm_sub_pd(iz2,jz0);
721 /* Calculate squared distance and things based on it */
722 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
723 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
724 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
726 rinv00 = gmx_mm_invsqrt_pd(rsq00);
727 rinv10 = gmx_mm_invsqrt_pd(rsq10);
728 rinv20 = gmx_mm_invsqrt_pd(rsq20);
730 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
731 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
732 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
734 /* Load parameters for j particles */
735 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
736 vdwjidx0A = 2*vdwtype[jnrA+0];
737 vdwjidx0B = 2*vdwtype[jnrB+0];
739 fjx0 = _mm_setzero_pd();
740 fjy0 = _mm_setzero_pd();
741 fjz0 = _mm_setzero_pd();
743 /**************************
744 * CALCULATE INTERACTIONS *
745 **************************/
747 r00 = _mm_mul_pd(rsq00,rinv00);
749 /* Compute parameters for interactions between i and j atoms */
750 qq00 = _mm_mul_pd(iq0,jq0);
751 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
752 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
754 /* Calculate table index by multiplying r with table scale and truncate to integer */
755 rt = _mm_mul_pd(r00,vftabscale);
756 vfitab = _mm_cvttpd_epi32(rt);
757 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
758 vfitab = _mm_slli_epi32(vfitab,3);
760 /* REACTION-FIELD ELECTROSTATICS */
761 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
763 /* CUBIC SPLINE TABLE DISPERSION */
764 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
765 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
766 GMX_MM_TRANSPOSE2_PD(Y,F);
767 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
768 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
769 GMX_MM_TRANSPOSE2_PD(G,H);
770 Heps = _mm_mul_pd(vfeps,H);
771 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
772 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
773 fvdw6 = _mm_mul_pd(c6_00,FF);
775 /* CUBIC SPLINE TABLE REPULSION */
776 vfitab = _mm_add_epi32(vfitab,ifour);
777 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
778 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
779 GMX_MM_TRANSPOSE2_PD(Y,F);
780 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
781 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
782 GMX_MM_TRANSPOSE2_PD(G,H);
783 Heps = _mm_mul_pd(vfeps,H);
784 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
785 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
786 fvdw12 = _mm_mul_pd(c12_00,FF);
787 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
789 fscal = _mm_add_pd(felec,fvdw);
791 /* Calculate temporary vectorial force */
792 tx = _mm_mul_pd(fscal,dx00);
793 ty = _mm_mul_pd(fscal,dy00);
794 tz = _mm_mul_pd(fscal,dz00);
796 /* Update vectorial force */
797 fix0 = _mm_add_pd(fix0,tx);
798 fiy0 = _mm_add_pd(fiy0,ty);
799 fiz0 = _mm_add_pd(fiz0,tz);
801 fjx0 = _mm_add_pd(fjx0,tx);
802 fjy0 = _mm_add_pd(fjy0,ty);
803 fjz0 = _mm_add_pd(fjz0,tz);
805 /**************************
806 * CALCULATE INTERACTIONS *
807 **************************/
809 /* Compute parameters for interactions between i and j atoms */
810 qq10 = _mm_mul_pd(iq1,jq0);
812 /* REACTION-FIELD ELECTROSTATICS */
813 felec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
815 fscal = felec;
817 /* Calculate temporary vectorial force */
818 tx = _mm_mul_pd(fscal,dx10);
819 ty = _mm_mul_pd(fscal,dy10);
820 tz = _mm_mul_pd(fscal,dz10);
822 /* Update vectorial force */
823 fix1 = _mm_add_pd(fix1,tx);
824 fiy1 = _mm_add_pd(fiy1,ty);
825 fiz1 = _mm_add_pd(fiz1,tz);
827 fjx0 = _mm_add_pd(fjx0,tx);
828 fjy0 = _mm_add_pd(fjy0,ty);
829 fjz0 = _mm_add_pd(fjz0,tz);
831 /**************************
832 * CALCULATE INTERACTIONS *
833 **************************/
835 /* Compute parameters for interactions between i and j atoms */
836 qq20 = _mm_mul_pd(iq2,jq0);
838 /* REACTION-FIELD ELECTROSTATICS */
839 felec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
841 fscal = felec;
843 /* Calculate temporary vectorial force */
844 tx = _mm_mul_pd(fscal,dx20);
845 ty = _mm_mul_pd(fscal,dy20);
846 tz = _mm_mul_pd(fscal,dz20);
848 /* Update vectorial force */
849 fix2 = _mm_add_pd(fix2,tx);
850 fiy2 = _mm_add_pd(fiy2,ty);
851 fiz2 = _mm_add_pd(fiz2,tz);
853 fjx0 = _mm_add_pd(fjx0,tx);
854 fjy0 = _mm_add_pd(fjy0,ty);
855 fjz0 = _mm_add_pd(fjz0,tz);
857 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
859 /* Inner loop uses 111 flops */
862 if(jidx<j_index_end)
865 jnrA = jjnr[jidx];
866 j_coord_offsetA = DIM*jnrA;
868 /* load j atom coordinates */
869 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
870 &jx0,&jy0,&jz0);
872 /* Calculate displacement vector */
873 dx00 = _mm_sub_pd(ix0,jx0);
874 dy00 = _mm_sub_pd(iy0,jy0);
875 dz00 = _mm_sub_pd(iz0,jz0);
876 dx10 = _mm_sub_pd(ix1,jx0);
877 dy10 = _mm_sub_pd(iy1,jy0);
878 dz10 = _mm_sub_pd(iz1,jz0);
879 dx20 = _mm_sub_pd(ix2,jx0);
880 dy20 = _mm_sub_pd(iy2,jy0);
881 dz20 = _mm_sub_pd(iz2,jz0);
883 /* Calculate squared distance and things based on it */
884 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
885 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
886 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
888 rinv00 = gmx_mm_invsqrt_pd(rsq00);
889 rinv10 = gmx_mm_invsqrt_pd(rsq10);
890 rinv20 = gmx_mm_invsqrt_pd(rsq20);
892 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
893 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
894 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
896 /* Load parameters for j particles */
897 jq0 = _mm_load_sd(charge+jnrA+0);
898 vdwjidx0A = 2*vdwtype[jnrA+0];
900 fjx0 = _mm_setzero_pd();
901 fjy0 = _mm_setzero_pd();
902 fjz0 = _mm_setzero_pd();
904 /**************************
905 * CALCULATE INTERACTIONS *
906 **************************/
908 r00 = _mm_mul_pd(rsq00,rinv00);
910 /* Compute parameters for interactions between i and j atoms */
911 qq00 = _mm_mul_pd(iq0,jq0);
912 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
914 /* Calculate table index by multiplying r with table scale and truncate to integer */
915 rt = _mm_mul_pd(r00,vftabscale);
916 vfitab = _mm_cvttpd_epi32(rt);
917 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
918 vfitab = _mm_slli_epi32(vfitab,3);
920 /* REACTION-FIELD ELECTROSTATICS */
921 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
923 /* CUBIC SPLINE TABLE DISPERSION */
924 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
925 F = _mm_setzero_pd();
926 GMX_MM_TRANSPOSE2_PD(Y,F);
927 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
928 H = _mm_setzero_pd();
929 GMX_MM_TRANSPOSE2_PD(G,H);
930 Heps = _mm_mul_pd(vfeps,H);
931 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
932 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
933 fvdw6 = _mm_mul_pd(c6_00,FF);
935 /* CUBIC SPLINE TABLE REPULSION */
936 vfitab = _mm_add_epi32(vfitab,ifour);
937 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
938 F = _mm_setzero_pd();
939 GMX_MM_TRANSPOSE2_PD(Y,F);
940 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
941 H = _mm_setzero_pd();
942 GMX_MM_TRANSPOSE2_PD(G,H);
943 Heps = _mm_mul_pd(vfeps,H);
944 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
945 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
946 fvdw12 = _mm_mul_pd(c12_00,FF);
947 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
949 fscal = _mm_add_pd(felec,fvdw);
951 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
953 /* Calculate temporary vectorial force */
954 tx = _mm_mul_pd(fscal,dx00);
955 ty = _mm_mul_pd(fscal,dy00);
956 tz = _mm_mul_pd(fscal,dz00);
958 /* Update vectorial force */
959 fix0 = _mm_add_pd(fix0,tx);
960 fiy0 = _mm_add_pd(fiy0,ty);
961 fiz0 = _mm_add_pd(fiz0,tz);
963 fjx0 = _mm_add_pd(fjx0,tx);
964 fjy0 = _mm_add_pd(fjy0,ty);
965 fjz0 = _mm_add_pd(fjz0,tz);
967 /**************************
968 * CALCULATE INTERACTIONS *
969 **************************/
971 /* Compute parameters for interactions between i and j atoms */
972 qq10 = _mm_mul_pd(iq1,jq0);
974 /* REACTION-FIELD ELECTROSTATICS */
975 felec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
977 fscal = felec;
979 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
981 /* Calculate temporary vectorial force */
982 tx = _mm_mul_pd(fscal,dx10);
983 ty = _mm_mul_pd(fscal,dy10);
984 tz = _mm_mul_pd(fscal,dz10);
986 /* Update vectorial force */
987 fix1 = _mm_add_pd(fix1,tx);
988 fiy1 = _mm_add_pd(fiy1,ty);
989 fiz1 = _mm_add_pd(fiz1,tz);
991 fjx0 = _mm_add_pd(fjx0,tx);
992 fjy0 = _mm_add_pd(fjy0,ty);
993 fjz0 = _mm_add_pd(fjz0,tz);
995 /**************************
996 * CALCULATE INTERACTIONS *
997 **************************/
999 /* Compute parameters for interactions between i and j atoms */
1000 qq20 = _mm_mul_pd(iq2,jq0);
1002 /* REACTION-FIELD ELECTROSTATICS */
1003 felec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
1005 fscal = felec;
1007 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1009 /* Calculate temporary vectorial force */
1010 tx = _mm_mul_pd(fscal,dx20);
1011 ty = _mm_mul_pd(fscal,dy20);
1012 tz = _mm_mul_pd(fscal,dz20);
1014 /* Update vectorial force */
1015 fix2 = _mm_add_pd(fix2,tx);
1016 fiy2 = _mm_add_pd(fiy2,ty);
1017 fiz2 = _mm_add_pd(fiz2,tz);
1019 fjx0 = _mm_add_pd(fjx0,tx);
1020 fjy0 = _mm_add_pd(fjy0,ty);
1021 fjz0 = _mm_add_pd(fjz0,tz);
1023 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
1025 /* Inner loop uses 111 flops */
1028 /* End of innermost loop */
1030 gmx_mm_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1031 f+i_coord_offset,fshift+i_shift_offset);
1033 /* Increment number of inner iterations */
1034 inneriter += j_index_end - j_index_start;
1036 /* Outer loop uses 18 flops */
1039 /* Increment number of outer iterations */
1040 outeriter += nri;
1042 /* Update outer/inner flops */
1044 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*18 + inneriter*111);