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.
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
47 #include "kernelutil_sparc64_hpc_ace_double.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW3P1_VF_sparc64_hpc_ace_double
51 * Electrostatics interaction: ReactionField
52 * VdW interaction: LennardJones
53 * Geometry: Water3-Particle
54 * Calculate force/pot: PotentialAndForce
57 nb_kernel_ElecRF_VdwLJ_GeomW3P1_VF_sparc64_hpc_ace_double
58 (t_nblist * gmx_restrict nlist,
59 rvec * gmx_restrict xx,
60 rvec * gmx_restrict ff,
61 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 double precision SIMD, 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 _fjsp_v2r8 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80 _fjsp_v2r8 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
82 _fjsp_v2r8 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
84 _fjsp_v2r8 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
85 int vdwjidx0A,vdwjidx0B;
86 _fjsp_v2r8 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87 _fjsp_v2r8 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
88 _fjsp_v2r8 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
89 _fjsp_v2r8 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
90 _fjsp_v2r8 velec,felec,velecsum,facel,crf,krf,krf2;
93 _fjsp_v2r8 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
96 _fjsp_v2r8 one_sixth = gmx_fjsp_set1_v2r8(1.0/6.0);
97 _fjsp_v2r8 one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
99 _fjsp_v2r8 dummy_mask,cutoff_mask;
100 _fjsp_v2r8 one = gmx_fjsp_set1_v2r8(1.0);
101 _fjsp_v2r8 two = gmx_fjsp_set1_v2r8(2.0);
102 union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
109 jindex = nlist->jindex;
111 shiftidx = nlist->shift;
113 shiftvec = fr->shift_vec[0];
114 fshift = fr->fshift[0];
115 facel = gmx_fjsp_set1_v2r8(fr->epsfac);
116 charge = mdatoms->chargeA;
117 krf = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
118 krf2 = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
119 crf = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
120 nvdwtype = fr->ntype;
122 vdwtype = mdatoms->typeA;
124 /* Setup water-specific parameters */
125 inr = nlist->iinr[0];
126 iq0 = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+0]));
127 iq1 = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
128 iq2 = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
129 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
131 /* Avoid stupid compiler warnings */
139 /* Start outer loop over neighborlists */
140 for(iidx=0; iidx<nri; iidx++)
142 /* Load shift vector for this list */
143 i_shift_offset = DIM*shiftidx[iidx];
145 /* Load limits for loop over neighbors */
146 j_index_start = jindex[iidx];
147 j_index_end = jindex[iidx+1];
149 /* Get outer coordinate index */
151 i_coord_offset = DIM*inr;
153 /* Load i particle coords and add shift vector */
154 gmx_fjsp_load_shift_and_3rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
155 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
157 fix0 = _fjsp_setzero_v2r8();
158 fiy0 = _fjsp_setzero_v2r8();
159 fiz0 = _fjsp_setzero_v2r8();
160 fix1 = _fjsp_setzero_v2r8();
161 fiy1 = _fjsp_setzero_v2r8();
162 fiz1 = _fjsp_setzero_v2r8();
163 fix2 = _fjsp_setzero_v2r8();
164 fiy2 = _fjsp_setzero_v2r8();
165 fiz2 = _fjsp_setzero_v2r8();
167 /* Reset potential sums */
168 velecsum = _fjsp_setzero_v2r8();
169 vvdwsum = _fjsp_setzero_v2r8();
171 /* Start inner kernel loop */
172 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
175 /* Get j neighbor index, and coordinate index */
178 j_coord_offsetA = DIM*jnrA;
179 j_coord_offsetB = DIM*jnrB;
181 /* load j atom coordinates */
182 gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
185 /* Calculate displacement vector */
186 dx00 = _fjsp_sub_v2r8(ix0,jx0);
187 dy00 = _fjsp_sub_v2r8(iy0,jy0);
188 dz00 = _fjsp_sub_v2r8(iz0,jz0);
189 dx10 = _fjsp_sub_v2r8(ix1,jx0);
190 dy10 = _fjsp_sub_v2r8(iy1,jy0);
191 dz10 = _fjsp_sub_v2r8(iz1,jz0);
192 dx20 = _fjsp_sub_v2r8(ix2,jx0);
193 dy20 = _fjsp_sub_v2r8(iy2,jy0);
194 dz20 = _fjsp_sub_v2r8(iz2,jz0);
196 /* Calculate squared distance and things based on it */
197 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
198 rsq10 = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
199 rsq20 = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
201 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
202 rinv10 = gmx_fjsp_invsqrt_v2r8(rsq10);
203 rinv20 = gmx_fjsp_invsqrt_v2r8(rsq20);
205 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
206 rinvsq10 = _fjsp_mul_v2r8(rinv10,rinv10);
207 rinvsq20 = _fjsp_mul_v2r8(rinv20,rinv20);
209 /* Load parameters for j particles */
210 jq0 = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
211 vdwjidx0A = 2*vdwtype[jnrA+0];
212 vdwjidx0B = 2*vdwtype[jnrB+0];
214 fjx0 = _fjsp_setzero_v2r8();
215 fjy0 = _fjsp_setzero_v2r8();
216 fjz0 = _fjsp_setzero_v2r8();
218 /**************************
219 * CALCULATE INTERACTIONS *
220 **************************/
222 /* Compute parameters for interactions between i and j atoms */
223 qq00 = _fjsp_mul_v2r8(iq0,jq0);
224 gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
225 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
227 /* REACTION-FIELD ELECTROSTATICS */
228 velec = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq00,rinv00),crf));
229 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
231 /* LENNARD-JONES DISPERSION/REPULSION */
233 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
234 vvdw6 = _fjsp_mul_v2r8(c6_00,rinvsix);
235 vvdw12 = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
236 vvdw = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
237 fvdw = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
239 /* Update potential sum for this i atom from the interaction with this j atom. */
240 velecsum = _fjsp_add_v2r8(velecsum,velec);
241 vvdwsum = _fjsp_add_v2r8(vvdwsum,vvdw);
243 fscal = _fjsp_add_v2r8(felec,fvdw);
245 /* Update vectorial force */
246 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
247 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
248 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
250 fjx0 = _fjsp_madd_v2r8(dx00,fscal,fjx0);
251 fjy0 = _fjsp_madd_v2r8(dy00,fscal,fjy0);
252 fjz0 = _fjsp_madd_v2r8(dz00,fscal,fjz0);
254 /**************************
255 * CALCULATE INTERACTIONS *
256 **************************/
258 /* Compute parameters for interactions between i and j atoms */
259 qq10 = _fjsp_mul_v2r8(iq1,jq0);
261 /* REACTION-FIELD ELECTROSTATICS */
262 velec = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq10,rinv10),crf));
263 felec = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
265 /* Update potential sum for this i atom from the interaction with this j atom. */
266 velecsum = _fjsp_add_v2r8(velecsum,velec);
270 /* Update vectorial force */
271 fix1 = _fjsp_madd_v2r8(dx10,fscal,fix1);
272 fiy1 = _fjsp_madd_v2r8(dy10,fscal,fiy1);
273 fiz1 = _fjsp_madd_v2r8(dz10,fscal,fiz1);
275 fjx0 = _fjsp_madd_v2r8(dx10,fscal,fjx0);
276 fjy0 = _fjsp_madd_v2r8(dy10,fscal,fjy0);
277 fjz0 = _fjsp_madd_v2r8(dz10,fscal,fjz0);
279 /**************************
280 * CALCULATE INTERACTIONS *
281 **************************/
283 /* Compute parameters for interactions between i and j atoms */
284 qq20 = _fjsp_mul_v2r8(iq2,jq0);
286 /* REACTION-FIELD ELECTROSTATICS */
287 velec = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq20,rinv20),crf));
288 felec = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,krf2));
290 /* Update potential sum for this i atom from the interaction with this j atom. */
291 velecsum = _fjsp_add_v2r8(velecsum,velec);
295 /* Update vectorial force */
296 fix2 = _fjsp_madd_v2r8(dx20,fscal,fix2);
297 fiy2 = _fjsp_madd_v2r8(dy20,fscal,fiy2);
298 fiz2 = _fjsp_madd_v2r8(dz20,fscal,fiz2);
300 fjx0 = _fjsp_madd_v2r8(dx20,fscal,fjx0);
301 fjy0 = _fjsp_madd_v2r8(dy20,fscal,fjy0);
302 fjz0 = _fjsp_madd_v2r8(dz20,fscal,fjz0);
304 gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
306 /* Inner loop uses 120 flops */
313 j_coord_offsetA = DIM*jnrA;
315 /* load j atom coordinates */
316 gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
319 /* Calculate displacement vector */
320 dx00 = _fjsp_sub_v2r8(ix0,jx0);
321 dy00 = _fjsp_sub_v2r8(iy0,jy0);
322 dz00 = _fjsp_sub_v2r8(iz0,jz0);
323 dx10 = _fjsp_sub_v2r8(ix1,jx0);
324 dy10 = _fjsp_sub_v2r8(iy1,jy0);
325 dz10 = _fjsp_sub_v2r8(iz1,jz0);
326 dx20 = _fjsp_sub_v2r8(ix2,jx0);
327 dy20 = _fjsp_sub_v2r8(iy2,jy0);
328 dz20 = _fjsp_sub_v2r8(iz2,jz0);
330 /* Calculate squared distance and things based on it */
331 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
332 rsq10 = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
333 rsq20 = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
335 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
336 rinv10 = gmx_fjsp_invsqrt_v2r8(rsq10);
337 rinv20 = gmx_fjsp_invsqrt_v2r8(rsq20);
339 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
340 rinvsq10 = _fjsp_mul_v2r8(rinv10,rinv10);
341 rinvsq20 = _fjsp_mul_v2r8(rinv20,rinv20);
343 /* Load parameters for j particles */
344 jq0 = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
345 vdwjidx0A = 2*vdwtype[jnrA+0];
347 fjx0 = _fjsp_setzero_v2r8();
348 fjy0 = _fjsp_setzero_v2r8();
349 fjz0 = _fjsp_setzero_v2r8();
351 /**************************
352 * CALCULATE INTERACTIONS *
353 **************************/
355 /* Compute parameters for interactions between i and j atoms */
356 qq00 = _fjsp_mul_v2r8(iq0,jq0);
357 gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
359 /* REACTION-FIELD ELECTROSTATICS */
360 velec = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq00,rinv00),crf));
361 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
363 /* LENNARD-JONES DISPERSION/REPULSION */
365 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
366 vvdw6 = _fjsp_mul_v2r8(c6_00,rinvsix);
367 vvdw12 = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
368 vvdw = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
369 fvdw = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
371 /* Update potential sum for this i atom from the interaction with this j atom. */
372 velec = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
373 velecsum = _fjsp_add_v2r8(velecsum,velec);
374 vvdw = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
375 vvdwsum = _fjsp_add_v2r8(vvdwsum,vvdw);
377 fscal = _fjsp_add_v2r8(felec,fvdw);
379 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
381 /* Update vectorial force */
382 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
383 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
384 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
386 fjx0 = _fjsp_madd_v2r8(dx00,fscal,fjx0);
387 fjy0 = _fjsp_madd_v2r8(dy00,fscal,fjy0);
388 fjz0 = _fjsp_madd_v2r8(dz00,fscal,fjz0);
390 /**************************
391 * CALCULATE INTERACTIONS *
392 **************************/
394 /* Compute parameters for interactions between i and j atoms */
395 qq10 = _fjsp_mul_v2r8(iq1,jq0);
397 /* REACTION-FIELD ELECTROSTATICS */
398 velec = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq10,rinv10),crf));
399 felec = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
401 /* Update potential sum for this i atom from the interaction with this j atom. */
402 velec = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
403 velecsum = _fjsp_add_v2r8(velecsum,velec);
407 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
409 /* Update vectorial force */
410 fix1 = _fjsp_madd_v2r8(dx10,fscal,fix1);
411 fiy1 = _fjsp_madd_v2r8(dy10,fscal,fiy1);
412 fiz1 = _fjsp_madd_v2r8(dz10,fscal,fiz1);
414 fjx0 = _fjsp_madd_v2r8(dx10,fscal,fjx0);
415 fjy0 = _fjsp_madd_v2r8(dy10,fscal,fjy0);
416 fjz0 = _fjsp_madd_v2r8(dz10,fscal,fjz0);
418 /**************************
419 * CALCULATE INTERACTIONS *
420 **************************/
422 /* Compute parameters for interactions between i and j atoms */
423 qq20 = _fjsp_mul_v2r8(iq2,jq0);
425 /* REACTION-FIELD ELECTROSTATICS */
426 velec = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq20,rinv20),crf));
427 felec = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,krf2));
429 /* Update potential sum for this i atom from the interaction with this j atom. */
430 velec = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
431 velecsum = _fjsp_add_v2r8(velecsum,velec);
435 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
437 /* Update vectorial force */
438 fix2 = _fjsp_madd_v2r8(dx20,fscal,fix2);
439 fiy2 = _fjsp_madd_v2r8(dy20,fscal,fiy2);
440 fiz2 = _fjsp_madd_v2r8(dz20,fscal,fiz2);
442 fjx0 = _fjsp_madd_v2r8(dx20,fscal,fjx0);
443 fjy0 = _fjsp_madd_v2r8(dy20,fscal,fjy0);
444 fjz0 = _fjsp_madd_v2r8(dz20,fscal,fjz0);
446 gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
448 /* Inner loop uses 120 flops */
451 /* End of innermost loop */
453 gmx_fjsp_update_iforce_3atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
454 f+i_coord_offset,fshift+i_shift_offset);
457 /* Update potential energies */
458 gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
459 gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
461 /* Increment number of inner iterations */
462 inneriter += j_index_end - j_index_start;
464 /* Outer loop uses 20 flops */
467 /* Increment number of outer iterations */
470 /* Update outer/inner flops */
472 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*20 + inneriter*120);
475 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW3P1_F_sparc64_hpc_ace_double
476 * Electrostatics interaction: ReactionField
477 * VdW interaction: LennardJones
478 * Geometry: Water3-Particle
479 * Calculate force/pot: Force
482 nb_kernel_ElecRF_VdwLJ_GeomW3P1_F_sparc64_hpc_ace_double
483 (t_nblist * gmx_restrict nlist,
484 rvec * gmx_restrict xx,
485 rvec * gmx_restrict ff,
486 t_forcerec * gmx_restrict fr,
487 t_mdatoms * gmx_restrict mdatoms,
488 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
489 t_nrnb * gmx_restrict nrnb)
491 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
492 * just 0 for non-waters.
493 * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
494 * jnr indices corresponding to data put in the four positions in the SIMD register.
496 int i_shift_offset,i_coord_offset,outeriter,inneriter;
497 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
499 int j_coord_offsetA,j_coord_offsetB;
500 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
502 real *shiftvec,*fshift,*x,*f;
503 _fjsp_v2r8 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
505 _fjsp_v2r8 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
507 _fjsp_v2r8 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
509 _fjsp_v2r8 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
510 int vdwjidx0A,vdwjidx0B;
511 _fjsp_v2r8 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
512 _fjsp_v2r8 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
513 _fjsp_v2r8 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
514 _fjsp_v2r8 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
515 _fjsp_v2r8 velec,felec,velecsum,facel,crf,krf,krf2;
518 _fjsp_v2r8 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
521 _fjsp_v2r8 one_sixth = gmx_fjsp_set1_v2r8(1.0/6.0);
522 _fjsp_v2r8 one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
524 _fjsp_v2r8 dummy_mask,cutoff_mask;
525 _fjsp_v2r8 one = gmx_fjsp_set1_v2r8(1.0);
526 _fjsp_v2r8 two = gmx_fjsp_set1_v2r8(2.0);
527 union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
534 jindex = nlist->jindex;
536 shiftidx = nlist->shift;
538 shiftvec = fr->shift_vec[0];
539 fshift = fr->fshift[0];
540 facel = gmx_fjsp_set1_v2r8(fr->epsfac);
541 charge = mdatoms->chargeA;
542 krf = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
543 krf2 = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
544 crf = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
545 nvdwtype = fr->ntype;
547 vdwtype = mdatoms->typeA;
549 /* Setup water-specific parameters */
550 inr = nlist->iinr[0];
551 iq0 = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+0]));
552 iq1 = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
553 iq2 = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
554 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
556 /* Avoid stupid compiler warnings */
564 /* Start outer loop over neighborlists */
565 for(iidx=0; iidx<nri; iidx++)
567 /* Load shift vector for this list */
568 i_shift_offset = DIM*shiftidx[iidx];
570 /* Load limits for loop over neighbors */
571 j_index_start = jindex[iidx];
572 j_index_end = jindex[iidx+1];
574 /* Get outer coordinate index */
576 i_coord_offset = DIM*inr;
578 /* Load i particle coords and add shift vector */
579 gmx_fjsp_load_shift_and_3rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
580 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
582 fix0 = _fjsp_setzero_v2r8();
583 fiy0 = _fjsp_setzero_v2r8();
584 fiz0 = _fjsp_setzero_v2r8();
585 fix1 = _fjsp_setzero_v2r8();
586 fiy1 = _fjsp_setzero_v2r8();
587 fiz1 = _fjsp_setzero_v2r8();
588 fix2 = _fjsp_setzero_v2r8();
589 fiy2 = _fjsp_setzero_v2r8();
590 fiz2 = _fjsp_setzero_v2r8();
592 /* Start inner kernel loop */
593 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
596 /* Get j neighbor index, and coordinate index */
599 j_coord_offsetA = DIM*jnrA;
600 j_coord_offsetB = DIM*jnrB;
602 /* load j atom coordinates */
603 gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
606 /* Calculate displacement vector */
607 dx00 = _fjsp_sub_v2r8(ix0,jx0);
608 dy00 = _fjsp_sub_v2r8(iy0,jy0);
609 dz00 = _fjsp_sub_v2r8(iz0,jz0);
610 dx10 = _fjsp_sub_v2r8(ix1,jx0);
611 dy10 = _fjsp_sub_v2r8(iy1,jy0);
612 dz10 = _fjsp_sub_v2r8(iz1,jz0);
613 dx20 = _fjsp_sub_v2r8(ix2,jx0);
614 dy20 = _fjsp_sub_v2r8(iy2,jy0);
615 dz20 = _fjsp_sub_v2r8(iz2,jz0);
617 /* Calculate squared distance and things based on it */
618 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
619 rsq10 = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
620 rsq20 = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
622 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
623 rinv10 = gmx_fjsp_invsqrt_v2r8(rsq10);
624 rinv20 = gmx_fjsp_invsqrt_v2r8(rsq20);
626 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
627 rinvsq10 = _fjsp_mul_v2r8(rinv10,rinv10);
628 rinvsq20 = _fjsp_mul_v2r8(rinv20,rinv20);
630 /* Load parameters for j particles */
631 jq0 = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
632 vdwjidx0A = 2*vdwtype[jnrA+0];
633 vdwjidx0B = 2*vdwtype[jnrB+0];
635 fjx0 = _fjsp_setzero_v2r8();
636 fjy0 = _fjsp_setzero_v2r8();
637 fjz0 = _fjsp_setzero_v2r8();
639 /**************************
640 * CALCULATE INTERACTIONS *
641 **************************/
643 /* Compute parameters for interactions between i and j atoms */
644 qq00 = _fjsp_mul_v2r8(iq0,jq0);
645 gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
646 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
648 /* REACTION-FIELD ELECTROSTATICS */
649 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
651 /* LENNARD-JONES DISPERSION/REPULSION */
653 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
654 fvdw = _fjsp_mul_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,c6_00),_fjsp_mul_v2r8(rinvsix,rinvsq00));
656 fscal = _fjsp_add_v2r8(felec,fvdw);
658 /* Update vectorial force */
659 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
660 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
661 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
663 fjx0 = _fjsp_madd_v2r8(dx00,fscal,fjx0);
664 fjy0 = _fjsp_madd_v2r8(dy00,fscal,fjy0);
665 fjz0 = _fjsp_madd_v2r8(dz00,fscal,fjz0);
667 /**************************
668 * CALCULATE INTERACTIONS *
669 **************************/
671 /* Compute parameters for interactions between i and j atoms */
672 qq10 = _fjsp_mul_v2r8(iq1,jq0);
674 /* REACTION-FIELD ELECTROSTATICS */
675 felec = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
679 /* Update vectorial force */
680 fix1 = _fjsp_madd_v2r8(dx10,fscal,fix1);
681 fiy1 = _fjsp_madd_v2r8(dy10,fscal,fiy1);
682 fiz1 = _fjsp_madd_v2r8(dz10,fscal,fiz1);
684 fjx0 = _fjsp_madd_v2r8(dx10,fscal,fjx0);
685 fjy0 = _fjsp_madd_v2r8(dy10,fscal,fjy0);
686 fjz0 = _fjsp_madd_v2r8(dz10,fscal,fjz0);
688 /**************************
689 * CALCULATE INTERACTIONS *
690 **************************/
692 /* Compute parameters for interactions between i and j atoms */
693 qq20 = _fjsp_mul_v2r8(iq2,jq0);
695 /* REACTION-FIELD ELECTROSTATICS */
696 felec = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,krf2));
700 /* Update vectorial force */
701 fix2 = _fjsp_madd_v2r8(dx20,fscal,fix2);
702 fiy2 = _fjsp_madd_v2r8(dy20,fscal,fiy2);
703 fiz2 = _fjsp_madd_v2r8(dz20,fscal,fiz2);
705 fjx0 = _fjsp_madd_v2r8(dx20,fscal,fjx0);
706 fjy0 = _fjsp_madd_v2r8(dy20,fscal,fjy0);
707 fjz0 = _fjsp_madd_v2r8(dz20,fscal,fjz0);
709 gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
711 /* Inner loop uses 100 flops */
718 j_coord_offsetA = DIM*jnrA;
720 /* load j atom coordinates */
721 gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
724 /* Calculate displacement vector */
725 dx00 = _fjsp_sub_v2r8(ix0,jx0);
726 dy00 = _fjsp_sub_v2r8(iy0,jy0);
727 dz00 = _fjsp_sub_v2r8(iz0,jz0);
728 dx10 = _fjsp_sub_v2r8(ix1,jx0);
729 dy10 = _fjsp_sub_v2r8(iy1,jy0);
730 dz10 = _fjsp_sub_v2r8(iz1,jz0);
731 dx20 = _fjsp_sub_v2r8(ix2,jx0);
732 dy20 = _fjsp_sub_v2r8(iy2,jy0);
733 dz20 = _fjsp_sub_v2r8(iz2,jz0);
735 /* Calculate squared distance and things based on it */
736 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
737 rsq10 = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
738 rsq20 = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
740 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
741 rinv10 = gmx_fjsp_invsqrt_v2r8(rsq10);
742 rinv20 = gmx_fjsp_invsqrt_v2r8(rsq20);
744 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
745 rinvsq10 = _fjsp_mul_v2r8(rinv10,rinv10);
746 rinvsq20 = _fjsp_mul_v2r8(rinv20,rinv20);
748 /* Load parameters for j particles */
749 jq0 = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
750 vdwjidx0A = 2*vdwtype[jnrA+0];
752 fjx0 = _fjsp_setzero_v2r8();
753 fjy0 = _fjsp_setzero_v2r8();
754 fjz0 = _fjsp_setzero_v2r8();
756 /**************************
757 * CALCULATE INTERACTIONS *
758 **************************/
760 /* Compute parameters for interactions between i and j atoms */
761 qq00 = _fjsp_mul_v2r8(iq0,jq0);
762 gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
764 /* REACTION-FIELD ELECTROSTATICS */
765 felec = _fjsp_mul_v2r8(qq00,_fjsp_msub_v2r8(rinv00,rinvsq00,krf2));
767 /* LENNARD-JONES DISPERSION/REPULSION */
769 rinvsix = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
770 fvdw = _fjsp_mul_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,c6_00),_fjsp_mul_v2r8(rinvsix,rinvsq00));
772 fscal = _fjsp_add_v2r8(felec,fvdw);
774 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
776 /* Update vectorial force */
777 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
778 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
779 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
781 fjx0 = _fjsp_madd_v2r8(dx00,fscal,fjx0);
782 fjy0 = _fjsp_madd_v2r8(dy00,fscal,fjy0);
783 fjz0 = _fjsp_madd_v2r8(dz00,fscal,fjz0);
785 /**************************
786 * CALCULATE INTERACTIONS *
787 **************************/
789 /* Compute parameters for interactions between i and j atoms */
790 qq10 = _fjsp_mul_v2r8(iq1,jq0);
792 /* REACTION-FIELD ELECTROSTATICS */
793 felec = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
797 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
799 /* Update vectorial force */
800 fix1 = _fjsp_madd_v2r8(dx10,fscal,fix1);
801 fiy1 = _fjsp_madd_v2r8(dy10,fscal,fiy1);
802 fiz1 = _fjsp_madd_v2r8(dz10,fscal,fiz1);
804 fjx0 = _fjsp_madd_v2r8(dx10,fscal,fjx0);
805 fjy0 = _fjsp_madd_v2r8(dy10,fscal,fjy0);
806 fjz0 = _fjsp_madd_v2r8(dz10,fscal,fjz0);
808 /**************************
809 * CALCULATE INTERACTIONS *
810 **************************/
812 /* Compute parameters for interactions between i and j atoms */
813 qq20 = _fjsp_mul_v2r8(iq2,jq0);
815 /* REACTION-FIELD ELECTROSTATICS */
816 felec = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,krf2));
820 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
822 /* Update vectorial force */
823 fix2 = _fjsp_madd_v2r8(dx20,fscal,fix2);
824 fiy2 = _fjsp_madd_v2r8(dy20,fscal,fiy2);
825 fiz2 = _fjsp_madd_v2r8(dz20,fscal,fiz2);
827 fjx0 = _fjsp_madd_v2r8(dx20,fscal,fjx0);
828 fjy0 = _fjsp_madd_v2r8(dy20,fscal,fjy0);
829 fjz0 = _fjsp_madd_v2r8(dz20,fscal,fjz0);
831 gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
833 /* Inner loop uses 100 flops */
836 /* End of innermost loop */
838 gmx_fjsp_update_iforce_3atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
839 f+i_coord_offset,fshift+i_shift_offset);
841 /* Increment number of inner iterations */
842 inneriter += j_index_end - j_index_start;
844 /* Outer loop uses 18 flops */
847 /* Increment number of outer iterations */
850 /* Update outer/inner flops */
852 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*18 + inneriter*100);