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 sse2_double kernel generator.
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
47 #include "kernelutil_x86_sse2_double.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_sse2_double
51 * Electrostatics interaction: ReactionField
52 * VdW interaction: LennardJones
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
57 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_sse2_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;
74 int j_coord_offsetA,j_coord_offsetB;
75 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
77 real *shiftvec,*fshift,*x,*f;
78 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
81 int vdwjidx0A,vdwjidx0B;
82 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
83 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
84 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
87 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
90 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
91 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
92 __m128d rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
93 real rswitch_scalar,d_scalar;
94 __m128d dummy_mask,cutoff_mask;
95 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
96 __m128d one = _mm_set1_pd(1.0);
97 __m128d two = _mm_set1_pd(2.0);
103 jindex = nlist->jindex;
105 shiftidx = nlist->shift;
107 shiftvec = fr->shift_vec[0];
108 fshift = fr->fshift[0];
109 facel = _mm_set1_pd(fr->ic->epsfac);
110 charge = mdatoms->chargeA;
111 krf = _mm_set1_pd(fr->ic->k_rf);
112 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
113 crf = _mm_set1_pd(fr->ic->c_rf);
114 nvdwtype = fr->ntype;
116 vdwtype = mdatoms->typeA;
118 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
119 rcutoff_scalar = fr->ic->rcoulomb;
120 rcutoff = _mm_set1_pd(rcutoff_scalar);
121 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
123 rswitch_scalar = fr->ic->rvdw_switch;
124 rswitch = _mm_set1_pd(rswitch_scalar);
125 /* Setup switch parameters */
126 d_scalar = rcutoff_scalar-rswitch_scalar;
127 d = _mm_set1_pd(d_scalar);
128 swV3 = _mm_set1_pd(-10.0/(d_scalar*d_scalar*d_scalar));
129 swV4 = _mm_set1_pd( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
130 swV5 = _mm_set1_pd( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
131 swF2 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar));
132 swF3 = _mm_set1_pd( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
133 swF4 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
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_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
160 fix0 = _mm_setzero_pd();
161 fiy0 = _mm_setzero_pd();
162 fiz0 = _mm_setzero_pd();
164 /* Load parameters for i particles */
165 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
166 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
168 /* Reset potential sums */
169 velecsum = _mm_setzero_pd();
170 vvdwsum = _mm_setzero_pd();
172 /* Start inner kernel loop */
173 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
176 /* Get j neighbor index, and coordinate index */
179 j_coord_offsetA = DIM*jnrA;
180 j_coord_offsetB = DIM*jnrB;
182 /* load j atom coordinates */
183 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
186 /* Calculate displacement vector */
187 dx00 = _mm_sub_pd(ix0,jx0);
188 dy00 = _mm_sub_pd(iy0,jy0);
189 dz00 = _mm_sub_pd(iz0,jz0);
191 /* Calculate squared distance and things based on it */
192 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
194 rinv00 = sse2_invsqrt_d(rsq00);
196 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
198 /* Load parameters for j particles */
199 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
200 vdwjidx0A = 2*vdwtype[jnrA+0];
201 vdwjidx0B = 2*vdwtype[jnrB+0];
203 /**************************
204 * CALCULATE INTERACTIONS *
205 **************************/
207 if (gmx_mm_any_lt(rsq00,rcutoff2))
210 r00 = _mm_mul_pd(rsq00,rinv00);
212 /* Compute parameters for interactions between i and j atoms */
213 qq00 = _mm_mul_pd(iq0,jq0);
214 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
215 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
217 /* REACTION-FIELD ELECTROSTATICS */
218 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
219 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
221 /* LENNARD-JONES DISPERSION/REPULSION */
223 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
224 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
225 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
226 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
227 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
229 d = _mm_sub_pd(r00,rswitch);
230 d = _mm_max_pd(d,_mm_setzero_pd());
231 d2 = _mm_mul_pd(d,d);
232 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_add_pd(swV3,_mm_mul_pd(d,_mm_add_pd(swV4,_mm_mul_pd(d,swV5)))))));
234 dsw = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
236 /* Evaluate switch function */
237 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
238 fvdw = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
239 vvdw = _mm_mul_pd(vvdw,sw);
240 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
242 /* Update potential sum for this i atom from the interaction with this j atom. */
243 velec = _mm_and_pd(velec,cutoff_mask);
244 velecsum = _mm_add_pd(velecsum,velec);
245 vvdw = _mm_and_pd(vvdw,cutoff_mask);
246 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
248 fscal = _mm_add_pd(felec,fvdw);
250 fscal = _mm_and_pd(fscal,cutoff_mask);
252 /* Calculate temporary vectorial force */
253 tx = _mm_mul_pd(fscal,dx00);
254 ty = _mm_mul_pd(fscal,dy00);
255 tz = _mm_mul_pd(fscal,dz00);
257 /* Update vectorial force */
258 fix0 = _mm_add_pd(fix0,tx);
259 fiy0 = _mm_add_pd(fiy0,ty);
260 fiz0 = _mm_add_pd(fiz0,tz);
262 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
266 /* Inner loop uses 70 flops */
273 j_coord_offsetA = DIM*jnrA;
275 /* load j atom coordinates */
276 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
279 /* Calculate displacement vector */
280 dx00 = _mm_sub_pd(ix0,jx0);
281 dy00 = _mm_sub_pd(iy0,jy0);
282 dz00 = _mm_sub_pd(iz0,jz0);
284 /* Calculate squared distance and things based on it */
285 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
287 rinv00 = sse2_invsqrt_d(rsq00);
289 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
291 /* Load parameters for j particles */
292 jq0 = _mm_load_sd(charge+jnrA+0);
293 vdwjidx0A = 2*vdwtype[jnrA+0];
295 /**************************
296 * CALCULATE INTERACTIONS *
297 **************************/
299 if (gmx_mm_any_lt(rsq00,rcutoff2))
302 r00 = _mm_mul_pd(rsq00,rinv00);
304 /* Compute parameters for interactions between i and j atoms */
305 qq00 = _mm_mul_pd(iq0,jq0);
306 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
308 /* REACTION-FIELD ELECTROSTATICS */
309 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
310 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
312 /* LENNARD-JONES DISPERSION/REPULSION */
314 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
315 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
316 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
317 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
318 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
320 d = _mm_sub_pd(r00,rswitch);
321 d = _mm_max_pd(d,_mm_setzero_pd());
322 d2 = _mm_mul_pd(d,d);
323 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_add_pd(swV3,_mm_mul_pd(d,_mm_add_pd(swV4,_mm_mul_pd(d,swV5)))))));
325 dsw = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
327 /* Evaluate switch function */
328 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
329 fvdw = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
330 vvdw = _mm_mul_pd(vvdw,sw);
331 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
333 /* Update potential sum for this i atom from the interaction with this j atom. */
334 velec = _mm_and_pd(velec,cutoff_mask);
335 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
336 velecsum = _mm_add_pd(velecsum,velec);
337 vvdw = _mm_and_pd(vvdw,cutoff_mask);
338 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
339 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
341 fscal = _mm_add_pd(felec,fvdw);
343 fscal = _mm_and_pd(fscal,cutoff_mask);
345 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
347 /* Calculate temporary vectorial force */
348 tx = _mm_mul_pd(fscal,dx00);
349 ty = _mm_mul_pd(fscal,dy00);
350 tz = _mm_mul_pd(fscal,dz00);
352 /* Update vectorial force */
353 fix0 = _mm_add_pd(fix0,tx);
354 fiy0 = _mm_add_pd(fiy0,ty);
355 fiz0 = _mm_add_pd(fiz0,tz);
357 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
361 /* Inner loop uses 70 flops */
364 /* End of innermost loop */
366 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
367 f+i_coord_offset,fshift+i_shift_offset);
370 /* Update potential energies */
371 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
372 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
374 /* Increment number of inner iterations */
375 inneriter += j_index_end - j_index_start;
377 /* Outer loop uses 9 flops */
380 /* Increment number of outer iterations */
383 /* Update outer/inner flops */
385 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*70);
388 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_sse2_double
389 * Electrostatics interaction: ReactionField
390 * VdW interaction: LennardJones
391 * Geometry: Particle-Particle
392 * Calculate force/pot: Force
395 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_sse2_double
396 (t_nblist * gmx_restrict nlist,
397 rvec * gmx_restrict xx,
398 rvec * gmx_restrict ff,
399 struct t_forcerec * gmx_restrict fr,
400 t_mdatoms * gmx_restrict mdatoms,
401 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
402 t_nrnb * gmx_restrict nrnb)
404 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
405 * just 0 for non-waters.
406 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
407 * jnr indices corresponding to data put in the four positions in the SIMD register.
409 int i_shift_offset,i_coord_offset,outeriter,inneriter;
410 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
412 int j_coord_offsetA,j_coord_offsetB;
413 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
415 real *shiftvec,*fshift,*x,*f;
416 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
418 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
419 int vdwjidx0A,vdwjidx0B;
420 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
421 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
422 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
425 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
428 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
429 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
430 __m128d rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
431 real rswitch_scalar,d_scalar;
432 __m128d dummy_mask,cutoff_mask;
433 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
434 __m128d one = _mm_set1_pd(1.0);
435 __m128d two = _mm_set1_pd(2.0);
441 jindex = nlist->jindex;
443 shiftidx = nlist->shift;
445 shiftvec = fr->shift_vec[0];
446 fshift = fr->fshift[0];
447 facel = _mm_set1_pd(fr->ic->epsfac);
448 charge = mdatoms->chargeA;
449 krf = _mm_set1_pd(fr->ic->k_rf);
450 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
451 crf = _mm_set1_pd(fr->ic->c_rf);
452 nvdwtype = fr->ntype;
454 vdwtype = mdatoms->typeA;
456 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
457 rcutoff_scalar = fr->ic->rcoulomb;
458 rcutoff = _mm_set1_pd(rcutoff_scalar);
459 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
461 rswitch_scalar = fr->ic->rvdw_switch;
462 rswitch = _mm_set1_pd(rswitch_scalar);
463 /* Setup switch parameters */
464 d_scalar = rcutoff_scalar-rswitch_scalar;
465 d = _mm_set1_pd(d_scalar);
466 swV3 = _mm_set1_pd(-10.0/(d_scalar*d_scalar*d_scalar));
467 swV4 = _mm_set1_pd( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
468 swV5 = _mm_set1_pd( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
469 swF2 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar));
470 swF3 = _mm_set1_pd( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
471 swF4 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
473 /* Avoid stupid compiler warnings */
481 /* Start outer loop over neighborlists */
482 for(iidx=0; iidx<nri; iidx++)
484 /* Load shift vector for this list */
485 i_shift_offset = DIM*shiftidx[iidx];
487 /* Load limits for loop over neighbors */
488 j_index_start = jindex[iidx];
489 j_index_end = jindex[iidx+1];
491 /* Get outer coordinate index */
493 i_coord_offset = DIM*inr;
495 /* Load i particle coords and add shift vector */
496 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
498 fix0 = _mm_setzero_pd();
499 fiy0 = _mm_setzero_pd();
500 fiz0 = _mm_setzero_pd();
502 /* Load parameters for i particles */
503 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
504 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
506 /* Start inner kernel loop */
507 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
510 /* Get j neighbor index, and coordinate index */
513 j_coord_offsetA = DIM*jnrA;
514 j_coord_offsetB = DIM*jnrB;
516 /* load j atom coordinates */
517 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
520 /* Calculate displacement vector */
521 dx00 = _mm_sub_pd(ix0,jx0);
522 dy00 = _mm_sub_pd(iy0,jy0);
523 dz00 = _mm_sub_pd(iz0,jz0);
525 /* Calculate squared distance and things based on it */
526 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
528 rinv00 = sse2_invsqrt_d(rsq00);
530 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
532 /* Load parameters for j particles */
533 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
534 vdwjidx0A = 2*vdwtype[jnrA+0];
535 vdwjidx0B = 2*vdwtype[jnrB+0];
537 /**************************
538 * CALCULATE INTERACTIONS *
539 **************************/
541 if (gmx_mm_any_lt(rsq00,rcutoff2))
544 r00 = _mm_mul_pd(rsq00,rinv00);
546 /* Compute parameters for interactions between i and j atoms */
547 qq00 = _mm_mul_pd(iq0,jq0);
548 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
549 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
551 /* REACTION-FIELD ELECTROSTATICS */
552 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
554 /* LENNARD-JONES DISPERSION/REPULSION */
556 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
557 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
558 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
559 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
560 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
562 d = _mm_sub_pd(r00,rswitch);
563 d = _mm_max_pd(d,_mm_setzero_pd());
564 d2 = _mm_mul_pd(d,d);
565 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_add_pd(swV3,_mm_mul_pd(d,_mm_add_pd(swV4,_mm_mul_pd(d,swV5)))))));
567 dsw = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
569 /* Evaluate switch function */
570 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
571 fvdw = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
572 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
574 fscal = _mm_add_pd(felec,fvdw);
576 fscal = _mm_and_pd(fscal,cutoff_mask);
578 /* Calculate temporary vectorial force */
579 tx = _mm_mul_pd(fscal,dx00);
580 ty = _mm_mul_pd(fscal,dy00);
581 tz = _mm_mul_pd(fscal,dz00);
583 /* Update vectorial force */
584 fix0 = _mm_add_pd(fix0,tx);
585 fiy0 = _mm_add_pd(fiy0,ty);
586 fiz0 = _mm_add_pd(fiz0,tz);
588 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
592 /* Inner loop uses 61 flops */
599 j_coord_offsetA = DIM*jnrA;
601 /* load j atom coordinates */
602 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
605 /* Calculate displacement vector */
606 dx00 = _mm_sub_pd(ix0,jx0);
607 dy00 = _mm_sub_pd(iy0,jy0);
608 dz00 = _mm_sub_pd(iz0,jz0);
610 /* Calculate squared distance and things based on it */
611 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
613 rinv00 = sse2_invsqrt_d(rsq00);
615 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
617 /* Load parameters for j particles */
618 jq0 = _mm_load_sd(charge+jnrA+0);
619 vdwjidx0A = 2*vdwtype[jnrA+0];
621 /**************************
622 * CALCULATE INTERACTIONS *
623 **************************/
625 if (gmx_mm_any_lt(rsq00,rcutoff2))
628 r00 = _mm_mul_pd(rsq00,rinv00);
630 /* Compute parameters for interactions between i and j atoms */
631 qq00 = _mm_mul_pd(iq0,jq0);
632 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
634 /* REACTION-FIELD ELECTROSTATICS */
635 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
637 /* LENNARD-JONES DISPERSION/REPULSION */
639 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
640 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
641 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
642 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
643 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
645 d = _mm_sub_pd(r00,rswitch);
646 d = _mm_max_pd(d,_mm_setzero_pd());
647 d2 = _mm_mul_pd(d,d);
648 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_add_pd(swV3,_mm_mul_pd(d,_mm_add_pd(swV4,_mm_mul_pd(d,swV5)))))));
650 dsw = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
652 /* Evaluate switch function */
653 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
654 fvdw = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
655 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
657 fscal = _mm_add_pd(felec,fvdw);
659 fscal = _mm_and_pd(fscal,cutoff_mask);
661 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
663 /* Calculate temporary vectorial force */
664 tx = _mm_mul_pd(fscal,dx00);
665 ty = _mm_mul_pd(fscal,dy00);
666 tz = _mm_mul_pd(fscal,dz00);
668 /* Update vectorial force */
669 fix0 = _mm_add_pd(fix0,tx);
670 fiy0 = _mm_add_pd(fiy0,ty);
671 fiz0 = _mm_add_pd(fiz0,tz);
673 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
677 /* Inner loop uses 61 flops */
680 /* End of innermost loop */
682 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
683 f+i_coord_offset,fshift+i_shift_offset);
685 /* Increment number of inner iterations */
686 inneriter += j_index_end - j_index_start;
688 /* Outer loop uses 7 flops */
691 /* Increment number of outer iterations */
694 /* Update outer/inner flops */
696 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*61);