5174d0af5c885fa7b9aab2e1f62443c53636462c
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_double / nb_kernel_ElecCoul_VdwNone_GeomW4P1_avx_256_double.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_double kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46
47 #include "gromacs/simd/math_x86_avx_256_double.h"
48 #include "kernelutil_x86_avx_256_double.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW4P1_VF_avx_256_double
52  * Electrostatics interaction: Coulomb
53  * VdW interaction:            None
54  * Geometry:                   Water4-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCoul_VdwNone_GeomW4P1_VF_avx_256_double
59                     (t_nblist                    * gmx_restrict       nlist,
60                      rvec                        * gmx_restrict          xx,
61                      rvec                        * gmx_restrict          ff,
62                      t_forcerec                  * gmx_restrict          fr,
63                      t_mdatoms                   * gmx_restrict     mdatoms,
64                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65                      t_nrnb                      * gmx_restrict        nrnb)
66 {
67     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
68      * just 0 for non-waters.
69      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
70      * jnr indices corresponding to data put in the four positions in the SIMD register.
71      */
72     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
73     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74     int              jnrA,jnrB,jnrC,jnrD;
75     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
77     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
79     real             rcutoff_scalar;
80     real             *shiftvec,*fshift,*x,*f;
81     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82     real             scratch[4*DIM];
83     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
84     real *           vdwioffsetptr1;
85     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
86     real *           vdwioffsetptr2;
87     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
88     real *           vdwioffsetptr3;
89     __m256d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
90     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
91     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
92     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
93     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
94     __m256d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
95     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
96     real             *charge;
97     __m256d          dummy_mask,cutoff_mask;
98     __m128           tmpmask0,tmpmask1;
99     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
100     __m256d          one     = _mm256_set1_pd(1.0);
101     __m256d          two     = _mm256_set1_pd(2.0);
102     x                = xx[0];
103     f                = ff[0];
104
105     nri              = nlist->nri;
106     iinr             = nlist->iinr;
107     jindex           = nlist->jindex;
108     jjnr             = nlist->jjnr;
109     shiftidx         = nlist->shift;
110     gid              = nlist->gid;
111     shiftvec         = fr->shift_vec[0];
112     fshift           = fr->fshift[0];
113     facel            = _mm256_set1_pd(fr->epsfac);
114     charge           = mdatoms->chargeA;
115
116     /* Setup water-specific parameters */
117     inr              = nlist->iinr[0];
118     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
119     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
120     iq3              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+3]));
121
122     /* Avoid stupid compiler warnings */
123     jnrA = jnrB = jnrC = jnrD = 0;
124     j_coord_offsetA = 0;
125     j_coord_offsetB = 0;
126     j_coord_offsetC = 0;
127     j_coord_offsetD = 0;
128
129     outeriter        = 0;
130     inneriter        = 0;
131
132     for(iidx=0;iidx<4*DIM;iidx++)
133     {
134         scratch[iidx] = 0.0;
135     }
136
137     /* Start outer loop over neighborlists */
138     for(iidx=0; iidx<nri; iidx++)
139     {
140         /* Load shift vector for this list */
141         i_shift_offset   = DIM*shiftidx[iidx];
142
143         /* Load limits for loop over neighbors */
144         j_index_start    = jindex[iidx];
145         j_index_end      = jindex[iidx+1];
146
147         /* Get outer coordinate index */
148         inr              = iinr[iidx];
149         i_coord_offset   = DIM*inr;
150
151         /* Load i particle coords and add shift vector */
152         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
153                                                     &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
154
155         fix1             = _mm256_setzero_pd();
156         fiy1             = _mm256_setzero_pd();
157         fiz1             = _mm256_setzero_pd();
158         fix2             = _mm256_setzero_pd();
159         fiy2             = _mm256_setzero_pd();
160         fiz2             = _mm256_setzero_pd();
161         fix3             = _mm256_setzero_pd();
162         fiy3             = _mm256_setzero_pd();
163         fiz3             = _mm256_setzero_pd();
164
165         /* Reset potential sums */
166         velecsum         = _mm256_setzero_pd();
167
168         /* Start inner kernel loop */
169         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
170         {
171
172             /* Get j neighbor index, and coordinate index */
173             jnrA             = jjnr[jidx];
174             jnrB             = jjnr[jidx+1];
175             jnrC             = jjnr[jidx+2];
176             jnrD             = jjnr[jidx+3];
177             j_coord_offsetA  = DIM*jnrA;
178             j_coord_offsetB  = DIM*jnrB;
179             j_coord_offsetC  = DIM*jnrC;
180             j_coord_offsetD  = DIM*jnrD;
181
182             /* load j atom coordinates */
183             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
184                                                  x+j_coord_offsetC,x+j_coord_offsetD,
185                                                  &jx0,&jy0,&jz0);
186
187             /* Calculate displacement vector */
188             dx10             = _mm256_sub_pd(ix1,jx0);
189             dy10             = _mm256_sub_pd(iy1,jy0);
190             dz10             = _mm256_sub_pd(iz1,jz0);
191             dx20             = _mm256_sub_pd(ix2,jx0);
192             dy20             = _mm256_sub_pd(iy2,jy0);
193             dz20             = _mm256_sub_pd(iz2,jz0);
194             dx30             = _mm256_sub_pd(ix3,jx0);
195             dy30             = _mm256_sub_pd(iy3,jy0);
196             dz30             = _mm256_sub_pd(iz3,jz0);
197
198             /* Calculate squared distance and things based on it */
199             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
200             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
201             rsq30            = gmx_mm256_calc_rsq_pd(dx30,dy30,dz30);
202
203             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
204             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
205             rinv30           = gmx_mm256_invsqrt_pd(rsq30);
206
207             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
208             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
209             rinvsq30         = _mm256_mul_pd(rinv30,rinv30);
210
211             /* Load parameters for j particles */
212             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
213                                                                  charge+jnrC+0,charge+jnrD+0);
214
215             fjx0             = _mm256_setzero_pd();
216             fjy0             = _mm256_setzero_pd();
217             fjz0             = _mm256_setzero_pd();
218
219             /**************************
220              * CALCULATE INTERACTIONS *
221              **************************/
222
223             /* Compute parameters for interactions between i and j atoms */
224             qq10             = _mm256_mul_pd(iq1,jq0);
225
226             /* COULOMB ELECTROSTATICS */
227             velec            = _mm256_mul_pd(qq10,rinv10);
228             felec            = _mm256_mul_pd(velec,rinvsq10);
229
230             /* Update potential sum for this i atom from the interaction with this j atom. */
231             velecsum         = _mm256_add_pd(velecsum,velec);
232
233             fscal            = felec;
234
235             /* Calculate temporary vectorial force */
236             tx               = _mm256_mul_pd(fscal,dx10);
237             ty               = _mm256_mul_pd(fscal,dy10);
238             tz               = _mm256_mul_pd(fscal,dz10);
239
240             /* Update vectorial force */
241             fix1             = _mm256_add_pd(fix1,tx);
242             fiy1             = _mm256_add_pd(fiy1,ty);
243             fiz1             = _mm256_add_pd(fiz1,tz);
244
245             fjx0             = _mm256_add_pd(fjx0,tx);
246             fjy0             = _mm256_add_pd(fjy0,ty);
247             fjz0             = _mm256_add_pd(fjz0,tz);
248
249             /**************************
250              * CALCULATE INTERACTIONS *
251              **************************/
252
253             /* Compute parameters for interactions between i and j atoms */
254             qq20             = _mm256_mul_pd(iq2,jq0);
255
256             /* COULOMB ELECTROSTATICS */
257             velec            = _mm256_mul_pd(qq20,rinv20);
258             felec            = _mm256_mul_pd(velec,rinvsq20);
259
260             /* Update potential sum for this i atom from the interaction with this j atom. */
261             velecsum         = _mm256_add_pd(velecsum,velec);
262
263             fscal            = felec;
264
265             /* Calculate temporary vectorial force */
266             tx               = _mm256_mul_pd(fscal,dx20);
267             ty               = _mm256_mul_pd(fscal,dy20);
268             tz               = _mm256_mul_pd(fscal,dz20);
269
270             /* Update vectorial force */
271             fix2             = _mm256_add_pd(fix2,tx);
272             fiy2             = _mm256_add_pd(fiy2,ty);
273             fiz2             = _mm256_add_pd(fiz2,tz);
274
275             fjx0             = _mm256_add_pd(fjx0,tx);
276             fjy0             = _mm256_add_pd(fjy0,ty);
277             fjz0             = _mm256_add_pd(fjz0,tz);
278
279             /**************************
280              * CALCULATE INTERACTIONS *
281              **************************/
282
283             /* Compute parameters for interactions between i and j atoms */
284             qq30             = _mm256_mul_pd(iq3,jq0);
285
286             /* COULOMB ELECTROSTATICS */
287             velec            = _mm256_mul_pd(qq30,rinv30);
288             felec            = _mm256_mul_pd(velec,rinvsq30);
289
290             /* Update potential sum for this i atom from the interaction with this j atom. */
291             velecsum         = _mm256_add_pd(velecsum,velec);
292
293             fscal            = felec;
294
295             /* Calculate temporary vectorial force */
296             tx               = _mm256_mul_pd(fscal,dx30);
297             ty               = _mm256_mul_pd(fscal,dy30);
298             tz               = _mm256_mul_pd(fscal,dz30);
299
300             /* Update vectorial force */
301             fix3             = _mm256_add_pd(fix3,tx);
302             fiy3             = _mm256_add_pd(fiy3,ty);
303             fiz3             = _mm256_add_pd(fiz3,tz);
304
305             fjx0             = _mm256_add_pd(fjx0,tx);
306             fjy0             = _mm256_add_pd(fjy0,ty);
307             fjz0             = _mm256_add_pd(fjz0,tz);
308
309             fjptrA             = f+j_coord_offsetA;
310             fjptrB             = f+j_coord_offsetB;
311             fjptrC             = f+j_coord_offsetC;
312             fjptrD             = f+j_coord_offsetD;
313
314             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
315
316             /* Inner loop uses 84 flops */
317         }
318
319         if(jidx<j_index_end)
320         {
321
322             /* Get j neighbor index, and coordinate index */
323             jnrlistA         = jjnr[jidx];
324             jnrlistB         = jjnr[jidx+1];
325             jnrlistC         = jjnr[jidx+2];
326             jnrlistD         = jjnr[jidx+3];
327             /* Sign of each element will be negative for non-real atoms.
328              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
329              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
330              */
331             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
332
333             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
334             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
335             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
336
337             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
338             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
339             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
340             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
341             j_coord_offsetA  = DIM*jnrA;
342             j_coord_offsetB  = DIM*jnrB;
343             j_coord_offsetC  = DIM*jnrC;
344             j_coord_offsetD  = DIM*jnrD;
345
346             /* load j atom coordinates */
347             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
348                                                  x+j_coord_offsetC,x+j_coord_offsetD,
349                                                  &jx0,&jy0,&jz0);
350
351             /* Calculate displacement vector */
352             dx10             = _mm256_sub_pd(ix1,jx0);
353             dy10             = _mm256_sub_pd(iy1,jy0);
354             dz10             = _mm256_sub_pd(iz1,jz0);
355             dx20             = _mm256_sub_pd(ix2,jx0);
356             dy20             = _mm256_sub_pd(iy2,jy0);
357             dz20             = _mm256_sub_pd(iz2,jz0);
358             dx30             = _mm256_sub_pd(ix3,jx0);
359             dy30             = _mm256_sub_pd(iy3,jy0);
360             dz30             = _mm256_sub_pd(iz3,jz0);
361
362             /* Calculate squared distance and things based on it */
363             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
364             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
365             rsq30            = gmx_mm256_calc_rsq_pd(dx30,dy30,dz30);
366
367             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
368             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
369             rinv30           = gmx_mm256_invsqrt_pd(rsq30);
370
371             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
372             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
373             rinvsq30         = _mm256_mul_pd(rinv30,rinv30);
374
375             /* Load parameters for j particles */
376             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
377                                                                  charge+jnrC+0,charge+jnrD+0);
378
379             fjx0             = _mm256_setzero_pd();
380             fjy0             = _mm256_setzero_pd();
381             fjz0             = _mm256_setzero_pd();
382
383             /**************************
384              * CALCULATE INTERACTIONS *
385              **************************/
386
387             /* Compute parameters for interactions between i and j atoms */
388             qq10             = _mm256_mul_pd(iq1,jq0);
389
390             /* COULOMB ELECTROSTATICS */
391             velec            = _mm256_mul_pd(qq10,rinv10);
392             felec            = _mm256_mul_pd(velec,rinvsq10);
393
394             /* Update potential sum for this i atom from the interaction with this j atom. */
395             velec            = _mm256_andnot_pd(dummy_mask,velec);
396             velecsum         = _mm256_add_pd(velecsum,velec);
397
398             fscal            = felec;
399
400             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
401
402             /* Calculate temporary vectorial force */
403             tx               = _mm256_mul_pd(fscal,dx10);
404             ty               = _mm256_mul_pd(fscal,dy10);
405             tz               = _mm256_mul_pd(fscal,dz10);
406
407             /* Update vectorial force */
408             fix1             = _mm256_add_pd(fix1,tx);
409             fiy1             = _mm256_add_pd(fiy1,ty);
410             fiz1             = _mm256_add_pd(fiz1,tz);
411
412             fjx0             = _mm256_add_pd(fjx0,tx);
413             fjy0             = _mm256_add_pd(fjy0,ty);
414             fjz0             = _mm256_add_pd(fjz0,tz);
415
416             /**************************
417              * CALCULATE INTERACTIONS *
418              **************************/
419
420             /* Compute parameters for interactions between i and j atoms */
421             qq20             = _mm256_mul_pd(iq2,jq0);
422
423             /* COULOMB ELECTROSTATICS */
424             velec            = _mm256_mul_pd(qq20,rinv20);
425             felec            = _mm256_mul_pd(velec,rinvsq20);
426
427             /* Update potential sum for this i atom from the interaction with this j atom. */
428             velec            = _mm256_andnot_pd(dummy_mask,velec);
429             velecsum         = _mm256_add_pd(velecsum,velec);
430
431             fscal            = felec;
432
433             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
434
435             /* Calculate temporary vectorial force */
436             tx               = _mm256_mul_pd(fscal,dx20);
437             ty               = _mm256_mul_pd(fscal,dy20);
438             tz               = _mm256_mul_pd(fscal,dz20);
439
440             /* Update vectorial force */
441             fix2             = _mm256_add_pd(fix2,tx);
442             fiy2             = _mm256_add_pd(fiy2,ty);
443             fiz2             = _mm256_add_pd(fiz2,tz);
444
445             fjx0             = _mm256_add_pd(fjx0,tx);
446             fjy0             = _mm256_add_pd(fjy0,ty);
447             fjz0             = _mm256_add_pd(fjz0,tz);
448
449             /**************************
450              * CALCULATE INTERACTIONS *
451              **************************/
452
453             /* Compute parameters for interactions between i and j atoms */
454             qq30             = _mm256_mul_pd(iq3,jq0);
455
456             /* COULOMB ELECTROSTATICS */
457             velec            = _mm256_mul_pd(qq30,rinv30);
458             felec            = _mm256_mul_pd(velec,rinvsq30);
459
460             /* Update potential sum for this i atom from the interaction with this j atom. */
461             velec            = _mm256_andnot_pd(dummy_mask,velec);
462             velecsum         = _mm256_add_pd(velecsum,velec);
463
464             fscal            = felec;
465
466             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
467
468             /* Calculate temporary vectorial force */
469             tx               = _mm256_mul_pd(fscal,dx30);
470             ty               = _mm256_mul_pd(fscal,dy30);
471             tz               = _mm256_mul_pd(fscal,dz30);
472
473             /* Update vectorial force */
474             fix3             = _mm256_add_pd(fix3,tx);
475             fiy3             = _mm256_add_pd(fiy3,ty);
476             fiz3             = _mm256_add_pd(fiz3,tz);
477
478             fjx0             = _mm256_add_pd(fjx0,tx);
479             fjy0             = _mm256_add_pd(fjy0,ty);
480             fjz0             = _mm256_add_pd(fjz0,tz);
481
482             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
483             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
484             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
485             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
486
487             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
488
489             /* Inner loop uses 84 flops */
490         }
491
492         /* End of innermost loop */
493
494         gmx_mm256_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
495                                                  f+i_coord_offset+DIM,fshift+i_shift_offset);
496
497         ggid                        = gid[iidx];
498         /* Update potential energies */
499         gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
500
501         /* Increment number of inner iterations */
502         inneriter                  += j_index_end - j_index_start;
503
504         /* Outer loop uses 19 flops */
505     }
506
507     /* Increment number of outer iterations */
508     outeriter        += nri;
509
510     /* Update outer/inner flops */
511
512     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*19 + inneriter*84);
513 }
514 /*
515  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_avx_256_double
516  * Electrostatics interaction: Coulomb
517  * VdW interaction:            None
518  * Geometry:                   Water4-Particle
519  * Calculate force/pot:        Force
520  */
521 void
522 nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_avx_256_double
523                     (t_nblist                    * gmx_restrict       nlist,
524                      rvec                        * gmx_restrict          xx,
525                      rvec                        * gmx_restrict          ff,
526                      t_forcerec                  * gmx_restrict          fr,
527                      t_mdatoms                   * gmx_restrict     mdatoms,
528                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
529                      t_nrnb                      * gmx_restrict        nrnb)
530 {
531     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
532      * just 0 for non-waters.
533      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
534      * jnr indices corresponding to data put in the four positions in the SIMD register.
535      */
536     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
537     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
538     int              jnrA,jnrB,jnrC,jnrD;
539     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
540     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
541     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
542     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
543     real             rcutoff_scalar;
544     real             *shiftvec,*fshift,*x,*f;
545     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
546     real             scratch[4*DIM];
547     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
548     real *           vdwioffsetptr1;
549     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
550     real *           vdwioffsetptr2;
551     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
552     real *           vdwioffsetptr3;
553     __m256d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
554     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
555     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
556     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
557     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
558     __m256d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
559     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
560     real             *charge;
561     __m256d          dummy_mask,cutoff_mask;
562     __m128           tmpmask0,tmpmask1;
563     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
564     __m256d          one     = _mm256_set1_pd(1.0);
565     __m256d          two     = _mm256_set1_pd(2.0);
566     x                = xx[0];
567     f                = ff[0];
568
569     nri              = nlist->nri;
570     iinr             = nlist->iinr;
571     jindex           = nlist->jindex;
572     jjnr             = nlist->jjnr;
573     shiftidx         = nlist->shift;
574     gid              = nlist->gid;
575     shiftvec         = fr->shift_vec[0];
576     fshift           = fr->fshift[0];
577     facel            = _mm256_set1_pd(fr->epsfac);
578     charge           = mdatoms->chargeA;
579
580     /* Setup water-specific parameters */
581     inr              = nlist->iinr[0];
582     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
583     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
584     iq3              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+3]));
585
586     /* Avoid stupid compiler warnings */
587     jnrA = jnrB = jnrC = jnrD = 0;
588     j_coord_offsetA = 0;
589     j_coord_offsetB = 0;
590     j_coord_offsetC = 0;
591     j_coord_offsetD = 0;
592
593     outeriter        = 0;
594     inneriter        = 0;
595
596     for(iidx=0;iidx<4*DIM;iidx++)
597     {
598         scratch[iidx] = 0.0;
599     }
600
601     /* Start outer loop over neighborlists */
602     for(iidx=0; iidx<nri; iidx++)
603     {
604         /* Load shift vector for this list */
605         i_shift_offset   = DIM*shiftidx[iidx];
606
607         /* Load limits for loop over neighbors */
608         j_index_start    = jindex[iidx];
609         j_index_end      = jindex[iidx+1];
610
611         /* Get outer coordinate index */
612         inr              = iinr[iidx];
613         i_coord_offset   = DIM*inr;
614
615         /* Load i particle coords and add shift vector */
616         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
617                                                     &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
618
619         fix1             = _mm256_setzero_pd();
620         fiy1             = _mm256_setzero_pd();
621         fiz1             = _mm256_setzero_pd();
622         fix2             = _mm256_setzero_pd();
623         fiy2             = _mm256_setzero_pd();
624         fiz2             = _mm256_setzero_pd();
625         fix3             = _mm256_setzero_pd();
626         fiy3             = _mm256_setzero_pd();
627         fiz3             = _mm256_setzero_pd();
628
629         /* Start inner kernel loop */
630         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
631         {
632
633             /* Get j neighbor index, and coordinate index */
634             jnrA             = jjnr[jidx];
635             jnrB             = jjnr[jidx+1];
636             jnrC             = jjnr[jidx+2];
637             jnrD             = jjnr[jidx+3];
638             j_coord_offsetA  = DIM*jnrA;
639             j_coord_offsetB  = DIM*jnrB;
640             j_coord_offsetC  = DIM*jnrC;
641             j_coord_offsetD  = DIM*jnrD;
642
643             /* load j atom coordinates */
644             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
645                                                  x+j_coord_offsetC,x+j_coord_offsetD,
646                                                  &jx0,&jy0,&jz0);
647
648             /* Calculate displacement vector */
649             dx10             = _mm256_sub_pd(ix1,jx0);
650             dy10             = _mm256_sub_pd(iy1,jy0);
651             dz10             = _mm256_sub_pd(iz1,jz0);
652             dx20             = _mm256_sub_pd(ix2,jx0);
653             dy20             = _mm256_sub_pd(iy2,jy0);
654             dz20             = _mm256_sub_pd(iz2,jz0);
655             dx30             = _mm256_sub_pd(ix3,jx0);
656             dy30             = _mm256_sub_pd(iy3,jy0);
657             dz30             = _mm256_sub_pd(iz3,jz0);
658
659             /* Calculate squared distance and things based on it */
660             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
661             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
662             rsq30            = gmx_mm256_calc_rsq_pd(dx30,dy30,dz30);
663
664             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
665             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
666             rinv30           = gmx_mm256_invsqrt_pd(rsq30);
667
668             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
669             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
670             rinvsq30         = _mm256_mul_pd(rinv30,rinv30);
671
672             /* Load parameters for j particles */
673             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
674                                                                  charge+jnrC+0,charge+jnrD+0);
675
676             fjx0             = _mm256_setzero_pd();
677             fjy0             = _mm256_setzero_pd();
678             fjz0             = _mm256_setzero_pd();
679
680             /**************************
681              * CALCULATE INTERACTIONS *
682              **************************/
683
684             /* Compute parameters for interactions between i and j atoms */
685             qq10             = _mm256_mul_pd(iq1,jq0);
686
687             /* COULOMB ELECTROSTATICS */
688             velec            = _mm256_mul_pd(qq10,rinv10);
689             felec            = _mm256_mul_pd(velec,rinvsq10);
690
691             fscal            = felec;
692
693             /* Calculate temporary vectorial force */
694             tx               = _mm256_mul_pd(fscal,dx10);
695             ty               = _mm256_mul_pd(fscal,dy10);
696             tz               = _mm256_mul_pd(fscal,dz10);
697
698             /* Update vectorial force */
699             fix1             = _mm256_add_pd(fix1,tx);
700             fiy1             = _mm256_add_pd(fiy1,ty);
701             fiz1             = _mm256_add_pd(fiz1,tz);
702
703             fjx0             = _mm256_add_pd(fjx0,tx);
704             fjy0             = _mm256_add_pd(fjy0,ty);
705             fjz0             = _mm256_add_pd(fjz0,tz);
706
707             /**************************
708              * CALCULATE INTERACTIONS *
709              **************************/
710
711             /* Compute parameters for interactions between i and j atoms */
712             qq20             = _mm256_mul_pd(iq2,jq0);
713
714             /* COULOMB ELECTROSTATICS */
715             velec            = _mm256_mul_pd(qq20,rinv20);
716             felec            = _mm256_mul_pd(velec,rinvsq20);
717
718             fscal            = felec;
719
720             /* Calculate temporary vectorial force */
721             tx               = _mm256_mul_pd(fscal,dx20);
722             ty               = _mm256_mul_pd(fscal,dy20);
723             tz               = _mm256_mul_pd(fscal,dz20);
724
725             /* Update vectorial force */
726             fix2             = _mm256_add_pd(fix2,tx);
727             fiy2             = _mm256_add_pd(fiy2,ty);
728             fiz2             = _mm256_add_pd(fiz2,tz);
729
730             fjx0             = _mm256_add_pd(fjx0,tx);
731             fjy0             = _mm256_add_pd(fjy0,ty);
732             fjz0             = _mm256_add_pd(fjz0,tz);
733
734             /**************************
735              * CALCULATE INTERACTIONS *
736              **************************/
737
738             /* Compute parameters for interactions between i and j atoms */
739             qq30             = _mm256_mul_pd(iq3,jq0);
740
741             /* COULOMB ELECTROSTATICS */
742             velec            = _mm256_mul_pd(qq30,rinv30);
743             felec            = _mm256_mul_pd(velec,rinvsq30);
744
745             fscal            = felec;
746
747             /* Calculate temporary vectorial force */
748             tx               = _mm256_mul_pd(fscal,dx30);
749             ty               = _mm256_mul_pd(fscal,dy30);
750             tz               = _mm256_mul_pd(fscal,dz30);
751
752             /* Update vectorial force */
753             fix3             = _mm256_add_pd(fix3,tx);
754             fiy3             = _mm256_add_pd(fiy3,ty);
755             fiz3             = _mm256_add_pd(fiz3,tz);
756
757             fjx0             = _mm256_add_pd(fjx0,tx);
758             fjy0             = _mm256_add_pd(fjy0,ty);
759             fjz0             = _mm256_add_pd(fjz0,tz);
760
761             fjptrA             = f+j_coord_offsetA;
762             fjptrB             = f+j_coord_offsetB;
763             fjptrC             = f+j_coord_offsetC;
764             fjptrD             = f+j_coord_offsetD;
765
766             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
767
768             /* Inner loop uses 81 flops */
769         }
770
771         if(jidx<j_index_end)
772         {
773
774             /* Get j neighbor index, and coordinate index */
775             jnrlistA         = jjnr[jidx];
776             jnrlistB         = jjnr[jidx+1];
777             jnrlistC         = jjnr[jidx+2];
778             jnrlistD         = jjnr[jidx+3];
779             /* Sign of each element will be negative for non-real atoms.
780              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
781              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
782              */
783             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
784
785             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
786             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
787             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
788
789             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
790             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
791             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
792             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
793             j_coord_offsetA  = DIM*jnrA;
794             j_coord_offsetB  = DIM*jnrB;
795             j_coord_offsetC  = DIM*jnrC;
796             j_coord_offsetD  = DIM*jnrD;
797
798             /* load j atom coordinates */
799             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
800                                                  x+j_coord_offsetC,x+j_coord_offsetD,
801                                                  &jx0,&jy0,&jz0);
802
803             /* Calculate displacement vector */
804             dx10             = _mm256_sub_pd(ix1,jx0);
805             dy10             = _mm256_sub_pd(iy1,jy0);
806             dz10             = _mm256_sub_pd(iz1,jz0);
807             dx20             = _mm256_sub_pd(ix2,jx0);
808             dy20             = _mm256_sub_pd(iy2,jy0);
809             dz20             = _mm256_sub_pd(iz2,jz0);
810             dx30             = _mm256_sub_pd(ix3,jx0);
811             dy30             = _mm256_sub_pd(iy3,jy0);
812             dz30             = _mm256_sub_pd(iz3,jz0);
813
814             /* Calculate squared distance and things based on it */
815             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
816             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
817             rsq30            = gmx_mm256_calc_rsq_pd(dx30,dy30,dz30);
818
819             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
820             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
821             rinv30           = gmx_mm256_invsqrt_pd(rsq30);
822
823             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
824             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
825             rinvsq30         = _mm256_mul_pd(rinv30,rinv30);
826
827             /* Load parameters for j particles */
828             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
829                                                                  charge+jnrC+0,charge+jnrD+0);
830
831             fjx0             = _mm256_setzero_pd();
832             fjy0             = _mm256_setzero_pd();
833             fjz0             = _mm256_setzero_pd();
834
835             /**************************
836              * CALCULATE INTERACTIONS *
837              **************************/
838
839             /* Compute parameters for interactions between i and j atoms */
840             qq10             = _mm256_mul_pd(iq1,jq0);
841
842             /* COULOMB ELECTROSTATICS */
843             velec            = _mm256_mul_pd(qq10,rinv10);
844             felec            = _mm256_mul_pd(velec,rinvsq10);
845
846             fscal            = felec;
847
848             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
849
850             /* Calculate temporary vectorial force */
851             tx               = _mm256_mul_pd(fscal,dx10);
852             ty               = _mm256_mul_pd(fscal,dy10);
853             tz               = _mm256_mul_pd(fscal,dz10);
854
855             /* Update vectorial force */
856             fix1             = _mm256_add_pd(fix1,tx);
857             fiy1             = _mm256_add_pd(fiy1,ty);
858             fiz1             = _mm256_add_pd(fiz1,tz);
859
860             fjx0             = _mm256_add_pd(fjx0,tx);
861             fjy0             = _mm256_add_pd(fjy0,ty);
862             fjz0             = _mm256_add_pd(fjz0,tz);
863
864             /**************************
865              * CALCULATE INTERACTIONS *
866              **************************/
867
868             /* Compute parameters for interactions between i and j atoms */
869             qq20             = _mm256_mul_pd(iq2,jq0);
870
871             /* COULOMB ELECTROSTATICS */
872             velec            = _mm256_mul_pd(qq20,rinv20);
873             felec            = _mm256_mul_pd(velec,rinvsq20);
874
875             fscal            = felec;
876
877             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
878
879             /* Calculate temporary vectorial force */
880             tx               = _mm256_mul_pd(fscal,dx20);
881             ty               = _mm256_mul_pd(fscal,dy20);
882             tz               = _mm256_mul_pd(fscal,dz20);
883
884             /* Update vectorial force */
885             fix2             = _mm256_add_pd(fix2,tx);
886             fiy2             = _mm256_add_pd(fiy2,ty);
887             fiz2             = _mm256_add_pd(fiz2,tz);
888
889             fjx0             = _mm256_add_pd(fjx0,tx);
890             fjy0             = _mm256_add_pd(fjy0,ty);
891             fjz0             = _mm256_add_pd(fjz0,tz);
892
893             /**************************
894              * CALCULATE INTERACTIONS *
895              **************************/
896
897             /* Compute parameters for interactions between i and j atoms */
898             qq30             = _mm256_mul_pd(iq3,jq0);
899
900             /* COULOMB ELECTROSTATICS */
901             velec            = _mm256_mul_pd(qq30,rinv30);
902             felec            = _mm256_mul_pd(velec,rinvsq30);
903
904             fscal            = felec;
905
906             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
907
908             /* Calculate temporary vectorial force */
909             tx               = _mm256_mul_pd(fscal,dx30);
910             ty               = _mm256_mul_pd(fscal,dy30);
911             tz               = _mm256_mul_pd(fscal,dz30);
912
913             /* Update vectorial force */
914             fix3             = _mm256_add_pd(fix3,tx);
915             fiy3             = _mm256_add_pd(fiy3,ty);
916             fiz3             = _mm256_add_pd(fiz3,tz);
917
918             fjx0             = _mm256_add_pd(fjx0,tx);
919             fjy0             = _mm256_add_pd(fjy0,ty);
920             fjz0             = _mm256_add_pd(fjz0,tz);
921
922             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
923             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
924             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
925             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
926
927             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
928
929             /* Inner loop uses 81 flops */
930         }
931
932         /* End of innermost loop */
933
934         gmx_mm256_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
935                                                  f+i_coord_offset+DIM,fshift+i_shift_offset);
936
937         /* Increment number of inner iterations */
938         inneriter                  += j_index_end - j_index_start;
939
940         /* Outer loop uses 18 flops */
941     }
942
943     /* Increment number of outer iterations */
944     outeriter        += nri;
945
946     /* Update outer/inner flops */
947
948     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*18 + inneriter*81);
949 }