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 sparc64_hpc_ace_double kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
46 #include "gromacs/legacyheaders/vec.h"
49 #include "kernelutil_sparc64_hpc_ace_double.h"
52 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_VF_sparc64_hpc_ace_double
53 * Electrostatics interaction: ReactionField
54 * VdW interaction: LennardJones
55 * Geometry: Particle-Particle
56 * Calculate force/pot: PotentialAndForce
59 nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_VF_sparc64_hpc_ace_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 double precision SIMD, 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;
76 int j_coord_offsetA,j_coord_offsetB;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
79 real *shiftvec,*fshift,*x,*f;
80 _fjsp_v2r8 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
82 _fjsp_v2r8 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
83 int vdwjidx0A,vdwjidx0B;
84 _fjsp_v2r8 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
85 _fjsp_v2r8 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
86 _fjsp_v2r8 velec,felec,velecsum,facel,crf,krf,krf2;
89 _fjsp_v2r8 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
92 _fjsp_v2r8 one_sixth = gmx_fjsp_set1_v2r8(1.0/6.0);
93 _fjsp_v2r8 one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
95 _fjsp_v2r8 dummy_mask,cutoff_mask;
96 _fjsp_v2r8 one = gmx_fjsp_set1_v2r8(1.0);
97 _fjsp_v2r8 two = gmx_fjsp_set1_v2r8(2.0);
98 union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
105 jindex = nlist->jindex;
107 shiftidx = nlist->shift;
109 shiftvec = fr->shift_vec[0];
110 fshift = fr->fshift[0];
111 facel = gmx_fjsp_set1_v2r8(fr->epsfac);
112 charge = mdatoms->chargeA;
113 krf = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
114 krf2 = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
115 crf = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
116 nvdwtype = fr->ntype;
118 vdwtype = mdatoms->typeA;
120 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
121 rcutoff_scalar = fr->rcoulomb;
122 rcutoff = gmx_fjsp_set1_v2r8(rcutoff_scalar);
123 rcutoff2 = _fjsp_mul_v2r8(rcutoff,rcutoff);
125 sh_vdw_invrcut6 = gmx_fjsp_set1_v2r8(fr->ic->sh_invrc6);
126 rvdw = gmx_fjsp_set1_v2r8(fr->rvdw);
128 /* Avoid stupid compiler warnings */
136 /* Start outer loop over neighborlists */
137 for(iidx=0; iidx<nri; iidx++)
139 /* Load shift vector for this list */
140 i_shift_offset = DIM*shiftidx[iidx];
142 /* Load limits for loop over neighbors */
143 j_index_start = jindex[iidx];
144 j_index_end = jindex[iidx+1];
146 /* Get outer coordinate index */
148 i_coord_offset = DIM*inr;
150 /* Load i particle coords and add shift vector */
151 gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
153 fix0 = _fjsp_setzero_v2r8();
154 fiy0 = _fjsp_setzero_v2r8();
155 fiz0 = _fjsp_setzero_v2r8();
157 /* Load parameters for i particles */
158 iq0 = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
159 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
161 /* Reset potential sums */
162 velecsum = _fjsp_setzero_v2r8();
163 vvdwsum = _fjsp_setzero_v2r8();
165 /* Start inner kernel loop */
166 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
169 /* Get j neighbor index, and coordinate index */
172 j_coord_offsetA = DIM*jnrA;
173 j_coord_offsetB = DIM*jnrB;
175 /* load j atom coordinates */
176 gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
179 /* Calculate displacement vector */
180 dx00 = _fjsp_sub_v2r8(ix0,jx0);
181 dy00 = _fjsp_sub_v2r8(iy0,jy0);
182 dz00 = _fjsp_sub_v2r8(iz0,jz0);
184 /* Calculate squared distance and things based on it */
185 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
187 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
189 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
191 /* Load parameters for j particles */
192 jq0 = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
193 vdwjidx0A = 2*vdwtype[jnrA+0];
194 vdwjidx0B = 2*vdwtype[jnrB+0];
196 /**************************
197 * CALCULATE INTERACTIONS *
198 **************************/
200 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
203 /* Compute parameters for interactions between i and j atoms */
204 qq00 = _fjsp_mul_v2r8(iq0,jq0);
205 gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
206 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
208 /* REACTION-FIELD ELECTROSTATICS */
209 velec = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq00,rinv00),crf));
210 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
212 /* LENNARD-JONES DISPERSION/REPULSION */
214 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
215 vvdw6 = _fjsp_mul_v2r8(c6_00,rinvsix);
216 vvdw12 = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
217 vvdw = _fjsp_msub_v2r8(_fjsp_nmsub_v2r8(c12_00,_fjsp_mul_v2r8(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
218 _fjsp_mul_v2r8(_fjsp_nmsub_v2r8( c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
219 fvdw = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
221 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
223 /* Update potential sum for this i atom from the interaction with this j atom. */
224 velec = _fjsp_and_v2r8(velec,cutoff_mask);
225 velecsum = _fjsp_add_v2r8(velecsum,velec);
226 vvdw = _fjsp_and_v2r8(vvdw,cutoff_mask);
227 vvdwsum = _fjsp_add_v2r8(vvdwsum,vvdw);
229 fscal = _fjsp_add_v2r8(felec,fvdw);
231 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
233 /* Update vectorial force */
234 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
235 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
236 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
238 gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
242 /* Inner loop uses 57 flops */
249 j_coord_offsetA = DIM*jnrA;
251 /* load j atom coordinates */
252 gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
255 /* Calculate displacement vector */
256 dx00 = _fjsp_sub_v2r8(ix0,jx0);
257 dy00 = _fjsp_sub_v2r8(iy0,jy0);
258 dz00 = _fjsp_sub_v2r8(iz0,jz0);
260 /* Calculate squared distance and things based on it */
261 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
263 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
265 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
267 /* Load parameters for j particles */
268 jq0 = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
269 vdwjidx0A = 2*vdwtype[jnrA+0];
271 /**************************
272 * CALCULATE INTERACTIONS *
273 **************************/
275 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
278 /* Compute parameters for interactions between i and j atoms */
279 qq00 = _fjsp_mul_v2r8(iq0,jq0);
280 gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
281 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
283 /* REACTION-FIELD ELECTROSTATICS */
284 velec = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq00,rinv00),crf));
285 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
287 /* LENNARD-JONES DISPERSION/REPULSION */
289 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
290 vvdw6 = _fjsp_mul_v2r8(c6_00,rinvsix);
291 vvdw12 = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
292 vvdw = _fjsp_msub_v2r8(_fjsp_nmsub_v2r8(c12_00,_fjsp_mul_v2r8(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
293 _fjsp_mul_v2r8(_fjsp_nmsub_v2r8( c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
294 fvdw = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
296 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
298 /* Update potential sum for this i atom from the interaction with this j atom. */
299 velec = _fjsp_and_v2r8(velec,cutoff_mask);
300 velec = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
301 velecsum = _fjsp_add_v2r8(velecsum,velec);
302 vvdw = _fjsp_and_v2r8(vvdw,cutoff_mask);
303 vvdw = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
304 vvdwsum = _fjsp_add_v2r8(vvdwsum,vvdw);
306 fscal = _fjsp_add_v2r8(felec,fvdw);
308 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
310 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
312 /* Update vectorial force */
313 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
314 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
315 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
317 gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
321 /* Inner loop uses 57 flops */
324 /* End of innermost loop */
326 gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
327 f+i_coord_offset,fshift+i_shift_offset);
330 /* Update potential energies */
331 gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
332 gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
334 /* Increment number of inner iterations */
335 inneriter += j_index_end - j_index_start;
337 /* Outer loop uses 9 flops */
340 /* Increment number of outer iterations */
343 /* Update outer/inner flops */
345 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*57);
348 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_F_sparc64_hpc_ace_double
349 * Electrostatics interaction: ReactionField
350 * VdW interaction: LennardJones
351 * Geometry: Particle-Particle
352 * Calculate force/pot: Force
355 nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_F_sparc64_hpc_ace_double
356 (t_nblist * gmx_restrict nlist,
357 rvec * gmx_restrict xx,
358 rvec * gmx_restrict ff,
359 t_forcerec * gmx_restrict fr,
360 t_mdatoms * gmx_restrict mdatoms,
361 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
362 t_nrnb * gmx_restrict nrnb)
364 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
365 * just 0 for non-waters.
366 * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
367 * jnr indices corresponding to data put in the four positions in the SIMD register.
369 int i_shift_offset,i_coord_offset,outeriter,inneriter;
370 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
372 int j_coord_offsetA,j_coord_offsetB;
373 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
375 real *shiftvec,*fshift,*x,*f;
376 _fjsp_v2r8 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
378 _fjsp_v2r8 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
379 int vdwjidx0A,vdwjidx0B;
380 _fjsp_v2r8 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
381 _fjsp_v2r8 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
382 _fjsp_v2r8 velec,felec,velecsum,facel,crf,krf,krf2;
385 _fjsp_v2r8 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
388 _fjsp_v2r8 one_sixth = gmx_fjsp_set1_v2r8(1.0/6.0);
389 _fjsp_v2r8 one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
391 _fjsp_v2r8 dummy_mask,cutoff_mask;
392 _fjsp_v2r8 one = gmx_fjsp_set1_v2r8(1.0);
393 _fjsp_v2r8 two = gmx_fjsp_set1_v2r8(2.0);
394 union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
401 jindex = nlist->jindex;
403 shiftidx = nlist->shift;
405 shiftvec = fr->shift_vec[0];
406 fshift = fr->fshift[0];
407 facel = gmx_fjsp_set1_v2r8(fr->epsfac);
408 charge = mdatoms->chargeA;
409 krf = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
410 krf2 = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
411 crf = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
412 nvdwtype = fr->ntype;
414 vdwtype = mdatoms->typeA;
416 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
417 rcutoff_scalar = fr->rcoulomb;
418 rcutoff = gmx_fjsp_set1_v2r8(rcutoff_scalar);
419 rcutoff2 = _fjsp_mul_v2r8(rcutoff,rcutoff);
421 sh_vdw_invrcut6 = gmx_fjsp_set1_v2r8(fr->ic->sh_invrc6);
422 rvdw = gmx_fjsp_set1_v2r8(fr->rvdw);
424 /* Avoid stupid compiler warnings */
432 /* Start outer loop over neighborlists */
433 for(iidx=0; iidx<nri; iidx++)
435 /* Load shift vector for this list */
436 i_shift_offset = DIM*shiftidx[iidx];
438 /* Load limits for loop over neighbors */
439 j_index_start = jindex[iidx];
440 j_index_end = jindex[iidx+1];
442 /* Get outer coordinate index */
444 i_coord_offset = DIM*inr;
446 /* Load i particle coords and add shift vector */
447 gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
449 fix0 = _fjsp_setzero_v2r8();
450 fiy0 = _fjsp_setzero_v2r8();
451 fiz0 = _fjsp_setzero_v2r8();
453 /* Load parameters for i particles */
454 iq0 = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
455 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
457 /* Start inner kernel loop */
458 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
461 /* Get j neighbor index, and coordinate index */
464 j_coord_offsetA = DIM*jnrA;
465 j_coord_offsetB = DIM*jnrB;
467 /* load j atom coordinates */
468 gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
471 /* Calculate displacement vector */
472 dx00 = _fjsp_sub_v2r8(ix0,jx0);
473 dy00 = _fjsp_sub_v2r8(iy0,jy0);
474 dz00 = _fjsp_sub_v2r8(iz0,jz0);
476 /* Calculate squared distance and things based on it */
477 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
479 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
481 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
483 /* Load parameters for j particles */
484 jq0 = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
485 vdwjidx0A = 2*vdwtype[jnrA+0];
486 vdwjidx0B = 2*vdwtype[jnrB+0];
488 /**************************
489 * CALCULATE INTERACTIONS *
490 **************************/
492 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
495 /* Compute parameters for interactions between i and j atoms */
496 qq00 = _fjsp_mul_v2r8(iq0,jq0);
497 gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
498 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
500 /* REACTION-FIELD ELECTROSTATICS */
501 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
503 /* LENNARD-JONES DISPERSION/REPULSION */
505 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
506 fvdw = _fjsp_mul_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,c6_00),_fjsp_mul_v2r8(rinvsix,rinvsq00));
508 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
510 fscal = _fjsp_add_v2r8(felec,fvdw);
512 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
514 /* Update vectorial force */
515 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
516 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
517 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
519 gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
523 /* Inner loop uses 40 flops */
530 j_coord_offsetA = DIM*jnrA;
532 /* load j atom coordinates */
533 gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
536 /* Calculate displacement vector */
537 dx00 = _fjsp_sub_v2r8(ix0,jx0);
538 dy00 = _fjsp_sub_v2r8(iy0,jy0);
539 dz00 = _fjsp_sub_v2r8(iz0,jz0);
541 /* Calculate squared distance and things based on it */
542 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
544 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
546 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
548 /* Load parameters for j particles */
549 jq0 = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
550 vdwjidx0A = 2*vdwtype[jnrA+0];
552 /**************************
553 * CALCULATE INTERACTIONS *
554 **************************/
556 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
559 /* Compute parameters for interactions between i and j atoms */
560 qq00 = _fjsp_mul_v2r8(iq0,jq0);
561 gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
562 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
564 /* REACTION-FIELD ELECTROSTATICS */
565 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
567 /* LENNARD-JONES DISPERSION/REPULSION */
569 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
570 fvdw = _fjsp_mul_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,c6_00),_fjsp_mul_v2r8(rinvsix,rinvsq00));
572 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
574 fscal = _fjsp_add_v2r8(felec,fvdw);
576 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
578 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
580 /* Update vectorial force */
581 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
582 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
583 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
585 gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
589 /* Inner loop uses 40 flops */
592 /* End of innermost loop */
594 gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
595 f+i_coord_offset,fshift+i_shift_offset);
597 /* Increment number of inner iterations */
598 inneriter += j_index_end - j_index_start;
600 /* Outer loop uses 7 flops */
603 /* Increment number of outer iterations */
606 /* Update outer/inner flops */
608 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*40);