Introduce gmxpre.h for truly global definitions
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_double / nb_kernel_ElecRFCut_VdwLJSw_GeomW3W3_sse2_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 sse2_double kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48
49 #include "gromacs/simd/math_x86_sse2_double.h"
50 #include "kernelutil_x86_sse2_double.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwLJSw_GeomW3W3_VF_sse2_double
54  * Electrostatics interaction: ReactionField
55  * VdW interaction:            LennardJones
56  * Geometry:                   Water3-Water3
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecRFCut_VdwLJSw_GeomW3W3_VF_sse2_double
61                     (t_nblist                    * gmx_restrict       nlist,
62                      rvec                        * gmx_restrict          xx,
63                      rvec                        * gmx_restrict          ff,
64                      t_forcerec                  * gmx_restrict          fr,
65                      t_mdatoms                   * gmx_restrict     mdatoms,
66                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67                      t_nrnb                      * gmx_restrict        nrnb)
68 {
69     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70      * just 0 for non-waters.
71      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
72      * jnr indices corresponding to data put in the four positions in the SIMD register.
73      */
74     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
75     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76     int              jnrA,jnrB;
77     int              j_coord_offsetA,j_coord_offsetB;
78     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
79     real             rcutoff_scalar;
80     real             *shiftvec,*fshift,*x,*f;
81     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
82     int              vdwioffset0;
83     __m128d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84     int              vdwioffset1;
85     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
86     int              vdwioffset2;
87     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
88     int              vdwjidx0A,vdwjidx0B;
89     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
90     int              vdwjidx1A,vdwjidx1B;
91     __m128d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
92     int              vdwjidx2A,vdwjidx2B;
93     __m128d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
94     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
95     __m128d          dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
96     __m128d          dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
97     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
98     __m128d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
99     __m128d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
100     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
101     __m128d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
102     __m128d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
103     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
104     real             *charge;
105     int              nvdwtype;
106     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
107     int              *vdwtype;
108     real             *vdwparam;
109     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
110     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
111     __m128d          rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
112     real             rswitch_scalar,d_scalar;
113     __m128d          dummy_mask,cutoff_mask;
114     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
115     __m128d          one     = _mm_set1_pd(1.0);
116     __m128d          two     = _mm_set1_pd(2.0);
117     x                = xx[0];
118     f                = ff[0];
119
120     nri              = nlist->nri;
121     iinr             = nlist->iinr;
122     jindex           = nlist->jindex;
123     jjnr             = nlist->jjnr;
124     shiftidx         = nlist->shift;
125     gid              = nlist->gid;
126     shiftvec         = fr->shift_vec[0];
127     fshift           = fr->fshift[0];
128     facel            = _mm_set1_pd(fr->epsfac);
129     charge           = mdatoms->chargeA;
130     krf              = _mm_set1_pd(fr->ic->k_rf);
131     krf2             = _mm_set1_pd(fr->ic->k_rf*2.0);
132     crf              = _mm_set1_pd(fr->ic->c_rf);
133     nvdwtype         = fr->ntype;
134     vdwparam         = fr->nbfp;
135     vdwtype          = mdatoms->typeA;
136
137     /* Setup water-specific parameters */
138     inr              = nlist->iinr[0];
139     iq0              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+0]));
140     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
141     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
142     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
143
144     jq0              = _mm_set1_pd(charge[inr+0]);
145     jq1              = _mm_set1_pd(charge[inr+1]);
146     jq2              = _mm_set1_pd(charge[inr+2]);
147     vdwjidx0A        = 2*vdwtype[inr+0];
148     qq00             = _mm_mul_pd(iq0,jq0);
149     c6_00            = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A]);
150     c12_00           = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A+1]);
151     qq01             = _mm_mul_pd(iq0,jq1);
152     qq02             = _mm_mul_pd(iq0,jq2);
153     qq10             = _mm_mul_pd(iq1,jq0);
154     qq11             = _mm_mul_pd(iq1,jq1);
155     qq12             = _mm_mul_pd(iq1,jq2);
156     qq20             = _mm_mul_pd(iq2,jq0);
157     qq21             = _mm_mul_pd(iq2,jq1);
158     qq22             = _mm_mul_pd(iq2,jq2);
159
160     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
161     rcutoff_scalar   = fr->rcoulomb;
162     rcutoff          = _mm_set1_pd(rcutoff_scalar);
163     rcutoff2         = _mm_mul_pd(rcutoff,rcutoff);
164
165     rswitch_scalar   = fr->rvdw_switch;
166     rswitch          = _mm_set1_pd(rswitch_scalar);
167     /* Setup switch parameters */
168     d_scalar         = rcutoff_scalar-rswitch_scalar;
169     d                = _mm_set1_pd(d_scalar);
170     swV3             = _mm_set1_pd(-10.0/(d_scalar*d_scalar*d_scalar));
171     swV4             = _mm_set1_pd( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
172     swV5             = _mm_set1_pd( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
173     swF2             = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar));
174     swF3             = _mm_set1_pd( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
175     swF4             = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
176
177     /* Avoid stupid compiler warnings */
178     jnrA = jnrB = 0;
179     j_coord_offsetA = 0;
180     j_coord_offsetB = 0;
181
182     outeriter        = 0;
183     inneriter        = 0;
184
185     /* Start outer loop over neighborlists */
186     for(iidx=0; iidx<nri; iidx++)
187     {
188         /* Load shift vector for this list */
189         i_shift_offset   = DIM*shiftidx[iidx];
190
191         /* Load limits for loop over neighbors */
192         j_index_start    = jindex[iidx];
193         j_index_end      = jindex[iidx+1];
194
195         /* Get outer coordinate index */
196         inr              = iinr[iidx];
197         i_coord_offset   = DIM*inr;
198
199         /* Load i particle coords and add shift vector */
200         gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
201                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
202
203         fix0             = _mm_setzero_pd();
204         fiy0             = _mm_setzero_pd();
205         fiz0             = _mm_setzero_pd();
206         fix1             = _mm_setzero_pd();
207         fiy1             = _mm_setzero_pd();
208         fiz1             = _mm_setzero_pd();
209         fix2             = _mm_setzero_pd();
210         fiy2             = _mm_setzero_pd();
211         fiz2             = _mm_setzero_pd();
212
213         /* Reset potential sums */
214         velecsum         = _mm_setzero_pd();
215         vvdwsum          = _mm_setzero_pd();
216
217         /* Start inner kernel loop */
218         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
219         {
220
221             /* Get j neighbor index, and coordinate index */
222             jnrA             = jjnr[jidx];
223             jnrB             = jjnr[jidx+1];
224             j_coord_offsetA  = DIM*jnrA;
225             j_coord_offsetB  = DIM*jnrB;
226
227             /* load j atom coordinates */
228             gmx_mm_load_3rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
229                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
230
231             /* Calculate displacement vector */
232             dx00             = _mm_sub_pd(ix0,jx0);
233             dy00             = _mm_sub_pd(iy0,jy0);
234             dz00             = _mm_sub_pd(iz0,jz0);
235             dx01             = _mm_sub_pd(ix0,jx1);
236             dy01             = _mm_sub_pd(iy0,jy1);
237             dz01             = _mm_sub_pd(iz0,jz1);
238             dx02             = _mm_sub_pd(ix0,jx2);
239             dy02             = _mm_sub_pd(iy0,jy2);
240             dz02             = _mm_sub_pd(iz0,jz2);
241             dx10             = _mm_sub_pd(ix1,jx0);
242             dy10             = _mm_sub_pd(iy1,jy0);
243             dz10             = _mm_sub_pd(iz1,jz0);
244             dx11             = _mm_sub_pd(ix1,jx1);
245             dy11             = _mm_sub_pd(iy1,jy1);
246             dz11             = _mm_sub_pd(iz1,jz1);
247             dx12             = _mm_sub_pd(ix1,jx2);
248             dy12             = _mm_sub_pd(iy1,jy2);
249             dz12             = _mm_sub_pd(iz1,jz2);
250             dx20             = _mm_sub_pd(ix2,jx0);
251             dy20             = _mm_sub_pd(iy2,jy0);
252             dz20             = _mm_sub_pd(iz2,jz0);
253             dx21             = _mm_sub_pd(ix2,jx1);
254             dy21             = _mm_sub_pd(iy2,jy1);
255             dz21             = _mm_sub_pd(iz2,jz1);
256             dx22             = _mm_sub_pd(ix2,jx2);
257             dy22             = _mm_sub_pd(iy2,jy2);
258             dz22             = _mm_sub_pd(iz2,jz2);
259
260             /* Calculate squared distance and things based on it */
261             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
262             rsq01            = gmx_mm_calc_rsq_pd(dx01,dy01,dz01);
263             rsq02            = gmx_mm_calc_rsq_pd(dx02,dy02,dz02);
264             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
265             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
266             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
267             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
268             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
269             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
270
271             rinv00           = gmx_mm_invsqrt_pd(rsq00);
272             rinv01           = gmx_mm_invsqrt_pd(rsq01);
273             rinv02           = gmx_mm_invsqrt_pd(rsq02);
274             rinv10           = gmx_mm_invsqrt_pd(rsq10);
275             rinv11           = gmx_mm_invsqrt_pd(rsq11);
276             rinv12           = gmx_mm_invsqrt_pd(rsq12);
277             rinv20           = gmx_mm_invsqrt_pd(rsq20);
278             rinv21           = gmx_mm_invsqrt_pd(rsq21);
279             rinv22           = gmx_mm_invsqrt_pd(rsq22);
280
281             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
282             rinvsq01         = _mm_mul_pd(rinv01,rinv01);
283             rinvsq02         = _mm_mul_pd(rinv02,rinv02);
284             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
285             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
286             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
287             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
288             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
289             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
290
291             fjx0             = _mm_setzero_pd();
292             fjy0             = _mm_setzero_pd();
293             fjz0             = _mm_setzero_pd();
294             fjx1             = _mm_setzero_pd();
295             fjy1             = _mm_setzero_pd();
296             fjz1             = _mm_setzero_pd();
297             fjx2             = _mm_setzero_pd();
298             fjy2             = _mm_setzero_pd();
299             fjz2             = _mm_setzero_pd();
300
301             /**************************
302              * CALCULATE INTERACTIONS *
303              **************************/
304
305             if (gmx_mm_any_lt(rsq00,rcutoff2))
306             {
307
308             r00              = _mm_mul_pd(rsq00,rinv00);
309
310             /* REACTION-FIELD ELECTROSTATICS */
311             velec            = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
312             felec            = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
313
314             /* LENNARD-JONES DISPERSION/REPULSION */
315
316             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
317             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
318             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
319             vvdw             = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
320             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
321
322             d                = _mm_sub_pd(r00,rswitch);
323             d                = _mm_max_pd(d,_mm_setzero_pd());
324             d2               = _mm_mul_pd(d,d);
325             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)))))));
326
327             dsw              = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
328
329             /* Evaluate switch function */
330             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
331             fvdw             = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
332             vvdw             = _mm_mul_pd(vvdw,sw);
333             cutoff_mask      = _mm_cmplt_pd(rsq00,rcutoff2);
334
335             /* Update potential sum for this i atom from the interaction with this j atom. */
336             velec            = _mm_and_pd(velec,cutoff_mask);
337             velecsum         = _mm_add_pd(velecsum,velec);
338             vvdw             = _mm_and_pd(vvdw,cutoff_mask);
339             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
340
341             fscal            = _mm_add_pd(felec,fvdw);
342
343             fscal            = _mm_and_pd(fscal,cutoff_mask);
344
345             /* Calculate temporary vectorial force */
346             tx               = _mm_mul_pd(fscal,dx00);
347             ty               = _mm_mul_pd(fscal,dy00);
348             tz               = _mm_mul_pd(fscal,dz00);
349
350             /* Update vectorial force */
351             fix0             = _mm_add_pd(fix0,tx);
352             fiy0             = _mm_add_pd(fiy0,ty);
353             fiz0             = _mm_add_pd(fiz0,tz);
354
355             fjx0             = _mm_add_pd(fjx0,tx);
356             fjy0             = _mm_add_pd(fjy0,ty);
357             fjz0             = _mm_add_pd(fjz0,tz);
358
359             }
360
361             /**************************
362              * CALCULATE INTERACTIONS *
363              **************************/
364
365             if (gmx_mm_any_lt(rsq01,rcutoff2))
366             {
367
368             /* REACTION-FIELD ELECTROSTATICS */
369             velec            = _mm_mul_pd(qq01,_mm_sub_pd(_mm_add_pd(rinv01,_mm_mul_pd(krf,rsq01)),crf));
370             felec            = _mm_mul_pd(qq01,_mm_sub_pd(_mm_mul_pd(rinv01,rinvsq01),krf2));
371
372             cutoff_mask      = _mm_cmplt_pd(rsq01,rcutoff2);
373
374             /* Update potential sum for this i atom from the interaction with this j atom. */
375             velec            = _mm_and_pd(velec,cutoff_mask);
376             velecsum         = _mm_add_pd(velecsum,velec);
377
378             fscal            = felec;
379
380             fscal            = _mm_and_pd(fscal,cutoff_mask);
381
382             /* Calculate temporary vectorial force */
383             tx               = _mm_mul_pd(fscal,dx01);
384             ty               = _mm_mul_pd(fscal,dy01);
385             tz               = _mm_mul_pd(fscal,dz01);
386
387             /* Update vectorial force */
388             fix0             = _mm_add_pd(fix0,tx);
389             fiy0             = _mm_add_pd(fiy0,ty);
390             fiz0             = _mm_add_pd(fiz0,tz);
391
392             fjx1             = _mm_add_pd(fjx1,tx);
393             fjy1             = _mm_add_pd(fjy1,ty);
394             fjz1             = _mm_add_pd(fjz1,tz);
395
396             }
397
398             /**************************
399              * CALCULATE INTERACTIONS *
400              **************************/
401
402             if (gmx_mm_any_lt(rsq02,rcutoff2))
403             {
404
405             /* REACTION-FIELD ELECTROSTATICS */
406             velec            = _mm_mul_pd(qq02,_mm_sub_pd(_mm_add_pd(rinv02,_mm_mul_pd(krf,rsq02)),crf));
407             felec            = _mm_mul_pd(qq02,_mm_sub_pd(_mm_mul_pd(rinv02,rinvsq02),krf2));
408
409             cutoff_mask      = _mm_cmplt_pd(rsq02,rcutoff2);
410
411             /* Update potential sum for this i atom from the interaction with this j atom. */
412             velec            = _mm_and_pd(velec,cutoff_mask);
413             velecsum         = _mm_add_pd(velecsum,velec);
414
415             fscal            = felec;
416
417             fscal            = _mm_and_pd(fscal,cutoff_mask);
418
419             /* Calculate temporary vectorial force */
420             tx               = _mm_mul_pd(fscal,dx02);
421             ty               = _mm_mul_pd(fscal,dy02);
422             tz               = _mm_mul_pd(fscal,dz02);
423
424             /* Update vectorial force */
425             fix0             = _mm_add_pd(fix0,tx);
426             fiy0             = _mm_add_pd(fiy0,ty);
427             fiz0             = _mm_add_pd(fiz0,tz);
428
429             fjx2             = _mm_add_pd(fjx2,tx);
430             fjy2             = _mm_add_pd(fjy2,ty);
431             fjz2             = _mm_add_pd(fjz2,tz);
432
433             }
434
435             /**************************
436              * CALCULATE INTERACTIONS *
437              **************************/
438
439             if (gmx_mm_any_lt(rsq10,rcutoff2))
440             {
441
442             /* REACTION-FIELD ELECTROSTATICS */
443             velec            = _mm_mul_pd(qq10,_mm_sub_pd(_mm_add_pd(rinv10,_mm_mul_pd(krf,rsq10)),crf));
444             felec            = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
445
446             cutoff_mask      = _mm_cmplt_pd(rsq10,rcutoff2);
447
448             /* Update potential sum for this i atom from the interaction with this j atom. */
449             velec            = _mm_and_pd(velec,cutoff_mask);
450             velecsum         = _mm_add_pd(velecsum,velec);
451
452             fscal            = felec;
453
454             fscal            = _mm_and_pd(fscal,cutoff_mask);
455
456             /* Calculate temporary vectorial force */
457             tx               = _mm_mul_pd(fscal,dx10);
458             ty               = _mm_mul_pd(fscal,dy10);
459             tz               = _mm_mul_pd(fscal,dz10);
460
461             /* Update vectorial force */
462             fix1             = _mm_add_pd(fix1,tx);
463             fiy1             = _mm_add_pd(fiy1,ty);
464             fiz1             = _mm_add_pd(fiz1,tz);
465
466             fjx0             = _mm_add_pd(fjx0,tx);
467             fjy0             = _mm_add_pd(fjy0,ty);
468             fjz0             = _mm_add_pd(fjz0,tz);
469
470             }
471
472             /**************************
473              * CALCULATE INTERACTIONS *
474              **************************/
475
476             if (gmx_mm_any_lt(rsq11,rcutoff2))
477             {
478
479             /* REACTION-FIELD ELECTROSTATICS */
480             velec            = _mm_mul_pd(qq11,_mm_sub_pd(_mm_add_pd(rinv11,_mm_mul_pd(krf,rsq11)),crf));
481             felec            = _mm_mul_pd(qq11,_mm_sub_pd(_mm_mul_pd(rinv11,rinvsq11),krf2));
482
483             cutoff_mask      = _mm_cmplt_pd(rsq11,rcutoff2);
484
485             /* Update potential sum for this i atom from the interaction with this j atom. */
486             velec            = _mm_and_pd(velec,cutoff_mask);
487             velecsum         = _mm_add_pd(velecsum,velec);
488
489             fscal            = felec;
490
491             fscal            = _mm_and_pd(fscal,cutoff_mask);
492
493             /* Calculate temporary vectorial force */
494             tx               = _mm_mul_pd(fscal,dx11);
495             ty               = _mm_mul_pd(fscal,dy11);
496             tz               = _mm_mul_pd(fscal,dz11);
497
498             /* Update vectorial force */
499             fix1             = _mm_add_pd(fix1,tx);
500             fiy1             = _mm_add_pd(fiy1,ty);
501             fiz1             = _mm_add_pd(fiz1,tz);
502
503             fjx1             = _mm_add_pd(fjx1,tx);
504             fjy1             = _mm_add_pd(fjy1,ty);
505             fjz1             = _mm_add_pd(fjz1,tz);
506
507             }
508
509             /**************************
510              * CALCULATE INTERACTIONS *
511              **************************/
512
513             if (gmx_mm_any_lt(rsq12,rcutoff2))
514             {
515
516             /* REACTION-FIELD ELECTROSTATICS */
517             velec            = _mm_mul_pd(qq12,_mm_sub_pd(_mm_add_pd(rinv12,_mm_mul_pd(krf,rsq12)),crf));
518             felec            = _mm_mul_pd(qq12,_mm_sub_pd(_mm_mul_pd(rinv12,rinvsq12),krf2));
519
520             cutoff_mask      = _mm_cmplt_pd(rsq12,rcutoff2);
521
522             /* Update potential sum for this i atom from the interaction with this j atom. */
523             velec            = _mm_and_pd(velec,cutoff_mask);
524             velecsum         = _mm_add_pd(velecsum,velec);
525
526             fscal            = felec;
527
528             fscal            = _mm_and_pd(fscal,cutoff_mask);
529
530             /* Calculate temporary vectorial force */
531             tx               = _mm_mul_pd(fscal,dx12);
532             ty               = _mm_mul_pd(fscal,dy12);
533             tz               = _mm_mul_pd(fscal,dz12);
534
535             /* Update vectorial force */
536             fix1             = _mm_add_pd(fix1,tx);
537             fiy1             = _mm_add_pd(fiy1,ty);
538             fiz1             = _mm_add_pd(fiz1,tz);
539
540             fjx2             = _mm_add_pd(fjx2,tx);
541             fjy2             = _mm_add_pd(fjy2,ty);
542             fjz2             = _mm_add_pd(fjz2,tz);
543
544             }
545
546             /**************************
547              * CALCULATE INTERACTIONS *
548              **************************/
549
550             if (gmx_mm_any_lt(rsq20,rcutoff2))
551             {
552
553             /* REACTION-FIELD ELECTROSTATICS */
554             velec            = _mm_mul_pd(qq20,_mm_sub_pd(_mm_add_pd(rinv20,_mm_mul_pd(krf,rsq20)),crf));
555             felec            = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
556
557             cutoff_mask      = _mm_cmplt_pd(rsq20,rcutoff2);
558
559             /* Update potential sum for this i atom from the interaction with this j atom. */
560             velec            = _mm_and_pd(velec,cutoff_mask);
561             velecsum         = _mm_add_pd(velecsum,velec);
562
563             fscal            = felec;
564
565             fscal            = _mm_and_pd(fscal,cutoff_mask);
566
567             /* Calculate temporary vectorial force */
568             tx               = _mm_mul_pd(fscal,dx20);
569             ty               = _mm_mul_pd(fscal,dy20);
570             tz               = _mm_mul_pd(fscal,dz20);
571
572             /* Update vectorial force */
573             fix2             = _mm_add_pd(fix2,tx);
574             fiy2             = _mm_add_pd(fiy2,ty);
575             fiz2             = _mm_add_pd(fiz2,tz);
576
577             fjx0             = _mm_add_pd(fjx0,tx);
578             fjy0             = _mm_add_pd(fjy0,ty);
579             fjz0             = _mm_add_pd(fjz0,tz);
580
581             }
582
583             /**************************
584              * CALCULATE INTERACTIONS *
585              **************************/
586
587             if (gmx_mm_any_lt(rsq21,rcutoff2))
588             {
589
590             /* REACTION-FIELD ELECTROSTATICS */
591             velec            = _mm_mul_pd(qq21,_mm_sub_pd(_mm_add_pd(rinv21,_mm_mul_pd(krf,rsq21)),crf));
592             felec            = _mm_mul_pd(qq21,_mm_sub_pd(_mm_mul_pd(rinv21,rinvsq21),krf2));
593
594             cutoff_mask      = _mm_cmplt_pd(rsq21,rcutoff2);
595
596             /* Update potential sum for this i atom from the interaction with this j atom. */
597             velec            = _mm_and_pd(velec,cutoff_mask);
598             velecsum         = _mm_add_pd(velecsum,velec);
599
600             fscal            = felec;
601
602             fscal            = _mm_and_pd(fscal,cutoff_mask);
603
604             /* Calculate temporary vectorial force */
605             tx               = _mm_mul_pd(fscal,dx21);
606             ty               = _mm_mul_pd(fscal,dy21);
607             tz               = _mm_mul_pd(fscal,dz21);
608
609             /* Update vectorial force */
610             fix2             = _mm_add_pd(fix2,tx);
611             fiy2             = _mm_add_pd(fiy2,ty);
612             fiz2             = _mm_add_pd(fiz2,tz);
613
614             fjx1             = _mm_add_pd(fjx1,tx);
615             fjy1             = _mm_add_pd(fjy1,ty);
616             fjz1             = _mm_add_pd(fjz1,tz);
617
618             }
619
620             /**************************
621              * CALCULATE INTERACTIONS *
622              **************************/
623
624             if (gmx_mm_any_lt(rsq22,rcutoff2))
625             {
626
627             /* REACTION-FIELD ELECTROSTATICS */
628             velec            = _mm_mul_pd(qq22,_mm_sub_pd(_mm_add_pd(rinv22,_mm_mul_pd(krf,rsq22)),crf));
629             felec            = _mm_mul_pd(qq22,_mm_sub_pd(_mm_mul_pd(rinv22,rinvsq22),krf2));
630
631             cutoff_mask      = _mm_cmplt_pd(rsq22,rcutoff2);
632
633             /* Update potential sum for this i atom from the interaction with this j atom. */
634             velec            = _mm_and_pd(velec,cutoff_mask);
635             velecsum         = _mm_add_pd(velecsum,velec);
636
637             fscal            = felec;
638
639             fscal            = _mm_and_pd(fscal,cutoff_mask);
640
641             /* Calculate temporary vectorial force */
642             tx               = _mm_mul_pd(fscal,dx22);
643             ty               = _mm_mul_pd(fscal,dy22);
644             tz               = _mm_mul_pd(fscal,dz22);
645
646             /* Update vectorial force */
647             fix2             = _mm_add_pd(fix2,tx);
648             fiy2             = _mm_add_pd(fiy2,ty);
649             fiz2             = _mm_add_pd(fiz2,tz);
650
651             fjx2             = _mm_add_pd(fjx2,tx);
652             fjy2             = _mm_add_pd(fjy2,ty);
653             fjz2             = _mm_add_pd(fjz2,tz);
654
655             }
656
657             gmx_mm_decrement_3rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
658
659             /* Inner loop uses 358 flops */
660         }
661
662         if(jidx<j_index_end)
663         {
664
665             jnrA             = jjnr[jidx];
666             j_coord_offsetA  = DIM*jnrA;
667
668             /* load j atom coordinates */
669             gmx_mm_load_3rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
670                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
671
672             /* Calculate displacement vector */
673             dx00             = _mm_sub_pd(ix0,jx0);
674             dy00             = _mm_sub_pd(iy0,jy0);
675             dz00             = _mm_sub_pd(iz0,jz0);
676             dx01             = _mm_sub_pd(ix0,jx1);
677             dy01             = _mm_sub_pd(iy0,jy1);
678             dz01             = _mm_sub_pd(iz0,jz1);
679             dx02             = _mm_sub_pd(ix0,jx2);
680             dy02             = _mm_sub_pd(iy0,jy2);
681             dz02             = _mm_sub_pd(iz0,jz2);
682             dx10             = _mm_sub_pd(ix1,jx0);
683             dy10             = _mm_sub_pd(iy1,jy0);
684             dz10             = _mm_sub_pd(iz1,jz0);
685             dx11             = _mm_sub_pd(ix1,jx1);
686             dy11             = _mm_sub_pd(iy1,jy1);
687             dz11             = _mm_sub_pd(iz1,jz1);
688             dx12             = _mm_sub_pd(ix1,jx2);
689             dy12             = _mm_sub_pd(iy1,jy2);
690             dz12             = _mm_sub_pd(iz1,jz2);
691             dx20             = _mm_sub_pd(ix2,jx0);
692             dy20             = _mm_sub_pd(iy2,jy0);
693             dz20             = _mm_sub_pd(iz2,jz0);
694             dx21             = _mm_sub_pd(ix2,jx1);
695             dy21             = _mm_sub_pd(iy2,jy1);
696             dz21             = _mm_sub_pd(iz2,jz1);
697             dx22             = _mm_sub_pd(ix2,jx2);
698             dy22             = _mm_sub_pd(iy2,jy2);
699             dz22             = _mm_sub_pd(iz2,jz2);
700
701             /* Calculate squared distance and things based on it */
702             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
703             rsq01            = gmx_mm_calc_rsq_pd(dx01,dy01,dz01);
704             rsq02            = gmx_mm_calc_rsq_pd(dx02,dy02,dz02);
705             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
706             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
707             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
708             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
709             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
710             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
711
712             rinv00           = gmx_mm_invsqrt_pd(rsq00);
713             rinv01           = gmx_mm_invsqrt_pd(rsq01);
714             rinv02           = gmx_mm_invsqrt_pd(rsq02);
715             rinv10           = gmx_mm_invsqrt_pd(rsq10);
716             rinv11           = gmx_mm_invsqrt_pd(rsq11);
717             rinv12           = gmx_mm_invsqrt_pd(rsq12);
718             rinv20           = gmx_mm_invsqrt_pd(rsq20);
719             rinv21           = gmx_mm_invsqrt_pd(rsq21);
720             rinv22           = gmx_mm_invsqrt_pd(rsq22);
721
722             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
723             rinvsq01         = _mm_mul_pd(rinv01,rinv01);
724             rinvsq02         = _mm_mul_pd(rinv02,rinv02);
725             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
726             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
727             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
728             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
729             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
730             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
731
732             fjx0             = _mm_setzero_pd();
733             fjy0             = _mm_setzero_pd();
734             fjz0             = _mm_setzero_pd();
735             fjx1             = _mm_setzero_pd();
736             fjy1             = _mm_setzero_pd();
737             fjz1             = _mm_setzero_pd();
738             fjx2             = _mm_setzero_pd();
739             fjy2             = _mm_setzero_pd();
740             fjz2             = _mm_setzero_pd();
741
742             /**************************
743              * CALCULATE INTERACTIONS *
744              **************************/
745
746             if (gmx_mm_any_lt(rsq00,rcutoff2))
747             {
748
749             r00              = _mm_mul_pd(rsq00,rinv00);
750
751             /* REACTION-FIELD ELECTROSTATICS */
752             velec            = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
753             felec            = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
754
755             /* LENNARD-JONES DISPERSION/REPULSION */
756
757             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
758             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
759             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
760             vvdw             = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
761             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
762
763             d                = _mm_sub_pd(r00,rswitch);
764             d                = _mm_max_pd(d,_mm_setzero_pd());
765             d2               = _mm_mul_pd(d,d);
766             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)))))));
767
768             dsw              = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
769
770             /* Evaluate switch function */
771             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
772             fvdw             = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
773             vvdw             = _mm_mul_pd(vvdw,sw);
774             cutoff_mask      = _mm_cmplt_pd(rsq00,rcutoff2);
775
776             /* Update potential sum for this i atom from the interaction with this j atom. */
777             velec            = _mm_and_pd(velec,cutoff_mask);
778             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
779             velecsum         = _mm_add_pd(velecsum,velec);
780             vvdw             = _mm_and_pd(vvdw,cutoff_mask);
781             vvdw             = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
782             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
783
784             fscal            = _mm_add_pd(felec,fvdw);
785
786             fscal            = _mm_and_pd(fscal,cutoff_mask);
787
788             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
789
790             /* Calculate temporary vectorial force */
791             tx               = _mm_mul_pd(fscal,dx00);
792             ty               = _mm_mul_pd(fscal,dy00);
793             tz               = _mm_mul_pd(fscal,dz00);
794
795             /* Update vectorial force */
796             fix0             = _mm_add_pd(fix0,tx);
797             fiy0             = _mm_add_pd(fiy0,ty);
798             fiz0             = _mm_add_pd(fiz0,tz);
799
800             fjx0             = _mm_add_pd(fjx0,tx);
801             fjy0             = _mm_add_pd(fjy0,ty);
802             fjz0             = _mm_add_pd(fjz0,tz);
803
804             }
805
806             /**************************
807              * CALCULATE INTERACTIONS *
808              **************************/
809
810             if (gmx_mm_any_lt(rsq01,rcutoff2))
811             {
812
813             /* REACTION-FIELD ELECTROSTATICS */
814             velec            = _mm_mul_pd(qq01,_mm_sub_pd(_mm_add_pd(rinv01,_mm_mul_pd(krf,rsq01)),crf));
815             felec            = _mm_mul_pd(qq01,_mm_sub_pd(_mm_mul_pd(rinv01,rinvsq01),krf2));
816
817             cutoff_mask      = _mm_cmplt_pd(rsq01,rcutoff2);
818
819             /* Update potential sum for this i atom from the interaction with this j atom. */
820             velec            = _mm_and_pd(velec,cutoff_mask);
821             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
822             velecsum         = _mm_add_pd(velecsum,velec);
823
824             fscal            = felec;
825
826             fscal            = _mm_and_pd(fscal,cutoff_mask);
827
828             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
829
830             /* Calculate temporary vectorial force */
831             tx               = _mm_mul_pd(fscal,dx01);
832             ty               = _mm_mul_pd(fscal,dy01);
833             tz               = _mm_mul_pd(fscal,dz01);
834
835             /* Update vectorial force */
836             fix0             = _mm_add_pd(fix0,tx);
837             fiy0             = _mm_add_pd(fiy0,ty);
838             fiz0             = _mm_add_pd(fiz0,tz);
839
840             fjx1             = _mm_add_pd(fjx1,tx);
841             fjy1             = _mm_add_pd(fjy1,ty);
842             fjz1             = _mm_add_pd(fjz1,tz);
843
844             }
845
846             /**************************
847              * CALCULATE INTERACTIONS *
848              **************************/
849
850             if (gmx_mm_any_lt(rsq02,rcutoff2))
851             {
852
853             /* REACTION-FIELD ELECTROSTATICS */
854             velec            = _mm_mul_pd(qq02,_mm_sub_pd(_mm_add_pd(rinv02,_mm_mul_pd(krf,rsq02)),crf));
855             felec            = _mm_mul_pd(qq02,_mm_sub_pd(_mm_mul_pd(rinv02,rinvsq02),krf2));
856
857             cutoff_mask      = _mm_cmplt_pd(rsq02,rcutoff2);
858
859             /* Update potential sum for this i atom from the interaction with this j atom. */
860             velec            = _mm_and_pd(velec,cutoff_mask);
861             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
862             velecsum         = _mm_add_pd(velecsum,velec);
863
864             fscal            = felec;
865
866             fscal            = _mm_and_pd(fscal,cutoff_mask);
867
868             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
869
870             /* Calculate temporary vectorial force */
871             tx               = _mm_mul_pd(fscal,dx02);
872             ty               = _mm_mul_pd(fscal,dy02);
873             tz               = _mm_mul_pd(fscal,dz02);
874
875             /* Update vectorial force */
876             fix0             = _mm_add_pd(fix0,tx);
877             fiy0             = _mm_add_pd(fiy0,ty);
878             fiz0             = _mm_add_pd(fiz0,tz);
879
880             fjx2             = _mm_add_pd(fjx2,tx);
881             fjy2             = _mm_add_pd(fjy2,ty);
882             fjz2             = _mm_add_pd(fjz2,tz);
883
884             }
885
886             /**************************
887              * CALCULATE INTERACTIONS *
888              **************************/
889
890             if (gmx_mm_any_lt(rsq10,rcutoff2))
891             {
892
893             /* REACTION-FIELD ELECTROSTATICS */
894             velec            = _mm_mul_pd(qq10,_mm_sub_pd(_mm_add_pd(rinv10,_mm_mul_pd(krf,rsq10)),crf));
895             felec            = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
896
897             cutoff_mask      = _mm_cmplt_pd(rsq10,rcutoff2);
898
899             /* Update potential sum for this i atom from the interaction with this j atom. */
900             velec            = _mm_and_pd(velec,cutoff_mask);
901             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
902             velecsum         = _mm_add_pd(velecsum,velec);
903
904             fscal            = felec;
905
906             fscal            = _mm_and_pd(fscal,cutoff_mask);
907
908             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
909
910             /* Calculate temporary vectorial force */
911             tx               = _mm_mul_pd(fscal,dx10);
912             ty               = _mm_mul_pd(fscal,dy10);
913             tz               = _mm_mul_pd(fscal,dz10);
914
915             /* Update vectorial force */
916             fix1             = _mm_add_pd(fix1,tx);
917             fiy1             = _mm_add_pd(fiy1,ty);
918             fiz1             = _mm_add_pd(fiz1,tz);
919
920             fjx0             = _mm_add_pd(fjx0,tx);
921             fjy0             = _mm_add_pd(fjy0,ty);
922             fjz0             = _mm_add_pd(fjz0,tz);
923
924             }
925
926             /**************************
927              * CALCULATE INTERACTIONS *
928              **************************/
929
930             if (gmx_mm_any_lt(rsq11,rcutoff2))
931             {
932
933             /* REACTION-FIELD ELECTROSTATICS */
934             velec            = _mm_mul_pd(qq11,_mm_sub_pd(_mm_add_pd(rinv11,_mm_mul_pd(krf,rsq11)),crf));
935             felec            = _mm_mul_pd(qq11,_mm_sub_pd(_mm_mul_pd(rinv11,rinvsq11),krf2));
936
937             cutoff_mask      = _mm_cmplt_pd(rsq11,rcutoff2);
938
939             /* Update potential sum for this i atom from the interaction with this j atom. */
940             velec            = _mm_and_pd(velec,cutoff_mask);
941             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
942             velecsum         = _mm_add_pd(velecsum,velec);
943
944             fscal            = felec;
945
946             fscal            = _mm_and_pd(fscal,cutoff_mask);
947
948             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
949
950             /* Calculate temporary vectorial force */
951             tx               = _mm_mul_pd(fscal,dx11);
952             ty               = _mm_mul_pd(fscal,dy11);
953             tz               = _mm_mul_pd(fscal,dz11);
954
955             /* Update vectorial force */
956             fix1             = _mm_add_pd(fix1,tx);
957             fiy1             = _mm_add_pd(fiy1,ty);
958             fiz1             = _mm_add_pd(fiz1,tz);
959
960             fjx1             = _mm_add_pd(fjx1,tx);
961             fjy1             = _mm_add_pd(fjy1,ty);
962             fjz1             = _mm_add_pd(fjz1,tz);
963
964             }
965
966             /**************************
967              * CALCULATE INTERACTIONS *
968              **************************/
969
970             if (gmx_mm_any_lt(rsq12,rcutoff2))
971             {
972
973             /* REACTION-FIELD ELECTROSTATICS */
974             velec            = _mm_mul_pd(qq12,_mm_sub_pd(_mm_add_pd(rinv12,_mm_mul_pd(krf,rsq12)),crf));
975             felec            = _mm_mul_pd(qq12,_mm_sub_pd(_mm_mul_pd(rinv12,rinvsq12),krf2));
976
977             cutoff_mask      = _mm_cmplt_pd(rsq12,rcutoff2);
978
979             /* Update potential sum for this i atom from the interaction with this j atom. */
980             velec            = _mm_and_pd(velec,cutoff_mask);
981             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
982             velecsum         = _mm_add_pd(velecsum,velec);
983
984             fscal            = felec;
985
986             fscal            = _mm_and_pd(fscal,cutoff_mask);
987
988             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
989
990             /* Calculate temporary vectorial force */
991             tx               = _mm_mul_pd(fscal,dx12);
992             ty               = _mm_mul_pd(fscal,dy12);
993             tz               = _mm_mul_pd(fscal,dz12);
994
995             /* Update vectorial force */
996             fix1             = _mm_add_pd(fix1,tx);
997             fiy1             = _mm_add_pd(fiy1,ty);
998             fiz1             = _mm_add_pd(fiz1,tz);
999
1000             fjx2             = _mm_add_pd(fjx2,tx);
1001             fjy2             = _mm_add_pd(fjy2,ty);
1002             fjz2             = _mm_add_pd(fjz2,tz);
1003
1004             }
1005
1006             /**************************
1007              * CALCULATE INTERACTIONS *
1008              **************************/
1009
1010             if (gmx_mm_any_lt(rsq20,rcutoff2))
1011             {
1012
1013             /* REACTION-FIELD ELECTROSTATICS */
1014             velec            = _mm_mul_pd(qq20,_mm_sub_pd(_mm_add_pd(rinv20,_mm_mul_pd(krf,rsq20)),crf));
1015             felec            = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
1016
1017             cutoff_mask      = _mm_cmplt_pd(rsq20,rcutoff2);
1018
1019             /* Update potential sum for this i atom from the interaction with this j atom. */
1020             velec            = _mm_and_pd(velec,cutoff_mask);
1021             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1022             velecsum         = _mm_add_pd(velecsum,velec);
1023
1024             fscal            = felec;
1025
1026             fscal            = _mm_and_pd(fscal,cutoff_mask);
1027
1028             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1029
1030             /* Calculate temporary vectorial force */
1031             tx               = _mm_mul_pd(fscal,dx20);
1032             ty               = _mm_mul_pd(fscal,dy20);
1033             tz               = _mm_mul_pd(fscal,dz20);
1034
1035             /* Update vectorial force */
1036             fix2             = _mm_add_pd(fix2,tx);
1037             fiy2             = _mm_add_pd(fiy2,ty);
1038             fiz2             = _mm_add_pd(fiz2,tz);
1039
1040             fjx0             = _mm_add_pd(fjx0,tx);
1041             fjy0             = _mm_add_pd(fjy0,ty);
1042             fjz0             = _mm_add_pd(fjz0,tz);
1043
1044             }
1045
1046             /**************************
1047              * CALCULATE INTERACTIONS *
1048              **************************/
1049
1050             if (gmx_mm_any_lt(rsq21,rcutoff2))
1051             {
1052
1053             /* REACTION-FIELD ELECTROSTATICS */
1054             velec            = _mm_mul_pd(qq21,_mm_sub_pd(_mm_add_pd(rinv21,_mm_mul_pd(krf,rsq21)),crf));
1055             felec            = _mm_mul_pd(qq21,_mm_sub_pd(_mm_mul_pd(rinv21,rinvsq21),krf2));
1056
1057             cutoff_mask      = _mm_cmplt_pd(rsq21,rcutoff2);
1058
1059             /* Update potential sum for this i atom from the interaction with this j atom. */
1060             velec            = _mm_and_pd(velec,cutoff_mask);
1061             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1062             velecsum         = _mm_add_pd(velecsum,velec);
1063
1064             fscal            = felec;
1065
1066             fscal            = _mm_and_pd(fscal,cutoff_mask);
1067
1068             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1069
1070             /* Calculate temporary vectorial force */
1071             tx               = _mm_mul_pd(fscal,dx21);
1072             ty               = _mm_mul_pd(fscal,dy21);
1073             tz               = _mm_mul_pd(fscal,dz21);
1074
1075             /* Update vectorial force */
1076             fix2             = _mm_add_pd(fix2,tx);
1077             fiy2             = _mm_add_pd(fiy2,ty);
1078             fiz2             = _mm_add_pd(fiz2,tz);
1079
1080             fjx1             = _mm_add_pd(fjx1,tx);
1081             fjy1             = _mm_add_pd(fjy1,ty);
1082             fjz1             = _mm_add_pd(fjz1,tz);
1083
1084             }
1085
1086             /**************************
1087              * CALCULATE INTERACTIONS *
1088              **************************/
1089
1090             if (gmx_mm_any_lt(rsq22,rcutoff2))
1091             {
1092
1093             /* REACTION-FIELD ELECTROSTATICS */
1094             velec            = _mm_mul_pd(qq22,_mm_sub_pd(_mm_add_pd(rinv22,_mm_mul_pd(krf,rsq22)),crf));
1095             felec            = _mm_mul_pd(qq22,_mm_sub_pd(_mm_mul_pd(rinv22,rinvsq22),krf2));
1096
1097             cutoff_mask      = _mm_cmplt_pd(rsq22,rcutoff2);
1098
1099             /* Update potential sum for this i atom from the interaction with this j atom. */
1100             velec            = _mm_and_pd(velec,cutoff_mask);
1101             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1102             velecsum         = _mm_add_pd(velecsum,velec);
1103
1104             fscal            = felec;
1105
1106             fscal            = _mm_and_pd(fscal,cutoff_mask);
1107
1108             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1109
1110             /* Calculate temporary vectorial force */
1111             tx               = _mm_mul_pd(fscal,dx22);
1112             ty               = _mm_mul_pd(fscal,dy22);
1113             tz               = _mm_mul_pd(fscal,dz22);
1114
1115             /* Update vectorial force */
1116             fix2             = _mm_add_pd(fix2,tx);
1117             fiy2             = _mm_add_pd(fiy2,ty);
1118             fiz2             = _mm_add_pd(fiz2,tz);
1119
1120             fjx2             = _mm_add_pd(fjx2,tx);
1121             fjy2             = _mm_add_pd(fjy2,ty);
1122             fjz2             = _mm_add_pd(fjz2,tz);
1123
1124             }
1125
1126             gmx_mm_decrement_3rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1127
1128             /* Inner loop uses 358 flops */
1129         }
1130
1131         /* End of innermost loop */
1132
1133         gmx_mm_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1134                                               f+i_coord_offset,fshift+i_shift_offset);
1135
1136         ggid                        = gid[iidx];
1137         /* Update potential energies */
1138         gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
1139         gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
1140
1141         /* Increment number of inner iterations */
1142         inneriter                  += j_index_end - j_index_start;
1143
1144         /* Outer loop uses 20 flops */
1145     }
1146
1147     /* Increment number of outer iterations */
1148     outeriter        += nri;
1149
1150     /* Update outer/inner flops */
1151
1152     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_VF,outeriter*20 + inneriter*358);
1153 }
1154 /*
1155  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwLJSw_GeomW3W3_F_sse2_double
1156  * Electrostatics interaction: ReactionField
1157  * VdW interaction:            LennardJones
1158  * Geometry:                   Water3-Water3
1159  * Calculate force/pot:        Force
1160  */
1161 void
1162 nb_kernel_ElecRFCut_VdwLJSw_GeomW3W3_F_sse2_double
1163                     (t_nblist                    * gmx_restrict       nlist,
1164                      rvec                        * gmx_restrict          xx,
1165                      rvec                        * gmx_restrict          ff,
1166                      t_forcerec                  * gmx_restrict          fr,
1167                      t_mdatoms                   * gmx_restrict     mdatoms,
1168                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
1169                      t_nrnb                      * gmx_restrict        nrnb)
1170 {
1171     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
1172      * just 0 for non-waters.
1173      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
1174      * jnr indices corresponding to data put in the four positions in the SIMD register.
1175      */
1176     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
1177     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1178     int              jnrA,jnrB;
1179     int              j_coord_offsetA,j_coord_offsetB;
1180     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
1181     real             rcutoff_scalar;
1182     real             *shiftvec,*fshift,*x,*f;
1183     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
1184     int              vdwioffset0;
1185     __m128d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
1186     int              vdwioffset1;
1187     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1188     int              vdwioffset2;
1189     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1190     int              vdwjidx0A,vdwjidx0B;
1191     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1192     int              vdwjidx1A,vdwjidx1B;
1193     __m128d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1194     int              vdwjidx2A,vdwjidx2B;
1195     __m128d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1196     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1197     __m128d          dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
1198     __m128d          dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
1199     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
1200     __m128d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1201     __m128d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1202     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
1203     __m128d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1204     __m128d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1205     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
1206     real             *charge;
1207     int              nvdwtype;
1208     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1209     int              *vdwtype;
1210     real             *vdwparam;
1211     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
1212     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
1213     __m128d          rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
1214     real             rswitch_scalar,d_scalar;
1215     __m128d          dummy_mask,cutoff_mask;
1216     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
1217     __m128d          one     = _mm_set1_pd(1.0);
1218     __m128d          two     = _mm_set1_pd(2.0);
1219     x                = xx[0];
1220     f                = ff[0];
1221
1222     nri              = nlist->nri;
1223     iinr             = nlist->iinr;
1224     jindex           = nlist->jindex;
1225     jjnr             = nlist->jjnr;
1226     shiftidx         = nlist->shift;
1227     gid              = nlist->gid;
1228     shiftvec         = fr->shift_vec[0];
1229     fshift           = fr->fshift[0];
1230     facel            = _mm_set1_pd(fr->epsfac);
1231     charge           = mdatoms->chargeA;
1232     krf              = _mm_set1_pd(fr->ic->k_rf);
1233     krf2             = _mm_set1_pd(fr->ic->k_rf*2.0);
1234     crf              = _mm_set1_pd(fr->ic->c_rf);
1235     nvdwtype         = fr->ntype;
1236     vdwparam         = fr->nbfp;
1237     vdwtype          = mdatoms->typeA;
1238
1239     /* Setup water-specific parameters */
1240     inr              = nlist->iinr[0];
1241     iq0              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+0]));
1242     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
1243     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
1244     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
1245
1246     jq0              = _mm_set1_pd(charge[inr+0]);
1247     jq1              = _mm_set1_pd(charge[inr+1]);
1248     jq2              = _mm_set1_pd(charge[inr+2]);
1249     vdwjidx0A        = 2*vdwtype[inr+0];
1250     qq00             = _mm_mul_pd(iq0,jq0);
1251     c6_00            = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A]);
1252     c12_00           = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A+1]);
1253     qq01             = _mm_mul_pd(iq0,jq1);
1254     qq02             = _mm_mul_pd(iq0,jq2);
1255     qq10             = _mm_mul_pd(iq1,jq0);
1256     qq11             = _mm_mul_pd(iq1,jq1);
1257     qq12             = _mm_mul_pd(iq1,jq2);
1258     qq20             = _mm_mul_pd(iq2,jq0);
1259     qq21             = _mm_mul_pd(iq2,jq1);
1260     qq22             = _mm_mul_pd(iq2,jq2);
1261
1262     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
1263     rcutoff_scalar   = fr->rcoulomb;
1264     rcutoff          = _mm_set1_pd(rcutoff_scalar);
1265     rcutoff2         = _mm_mul_pd(rcutoff,rcutoff);
1266
1267     rswitch_scalar   = fr->rvdw_switch;
1268     rswitch          = _mm_set1_pd(rswitch_scalar);
1269     /* Setup switch parameters */
1270     d_scalar         = rcutoff_scalar-rswitch_scalar;
1271     d                = _mm_set1_pd(d_scalar);
1272     swV3             = _mm_set1_pd(-10.0/(d_scalar*d_scalar*d_scalar));
1273     swV4             = _mm_set1_pd( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
1274     swV5             = _mm_set1_pd( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
1275     swF2             = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar));
1276     swF3             = _mm_set1_pd( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
1277     swF4             = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
1278
1279     /* Avoid stupid compiler warnings */
1280     jnrA = jnrB = 0;
1281     j_coord_offsetA = 0;
1282     j_coord_offsetB = 0;
1283
1284     outeriter        = 0;
1285     inneriter        = 0;
1286
1287     /* Start outer loop over neighborlists */
1288     for(iidx=0; iidx<nri; iidx++)
1289     {
1290         /* Load shift vector for this list */
1291         i_shift_offset   = DIM*shiftidx[iidx];
1292
1293         /* Load limits for loop over neighbors */
1294         j_index_start    = jindex[iidx];
1295         j_index_end      = jindex[iidx+1];
1296
1297         /* Get outer coordinate index */
1298         inr              = iinr[iidx];
1299         i_coord_offset   = DIM*inr;
1300
1301         /* Load i particle coords and add shift vector */
1302         gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
1303                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
1304
1305         fix0             = _mm_setzero_pd();
1306         fiy0             = _mm_setzero_pd();
1307         fiz0             = _mm_setzero_pd();
1308         fix1             = _mm_setzero_pd();
1309         fiy1             = _mm_setzero_pd();
1310         fiz1             = _mm_setzero_pd();
1311         fix2             = _mm_setzero_pd();
1312         fiy2             = _mm_setzero_pd();
1313         fiz2             = _mm_setzero_pd();
1314
1315         /* Start inner kernel loop */
1316         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
1317         {
1318
1319             /* Get j neighbor index, and coordinate index */
1320             jnrA             = jjnr[jidx];
1321             jnrB             = jjnr[jidx+1];
1322             j_coord_offsetA  = DIM*jnrA;
1323             j_coord_offsetB  = DIM*jnrB;
1324
1325             /* load j atom coordinates */
1326             gmx_mm_load_3rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1327                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1328
1329             /* Calculate displacement vector */
1330             dx00             = _mm_sub_pd(ix0,jx0);
1331             dy00             = _mm_sub_pd(iy0,jy0);
1332             dz00             = _mm_sub_pd(iz0,jz0);
1333             dx01             = _mm_sub_pd(ix0,jx1);
1334             dy01             = _mm_sub_pd(iy0,jy1);
1335             dz01             = _mm_sub_pd(iz0,jz1);
1336             dx02             = _mm_sub_pd(ix0,jx2);
1337             dy02             = _mm_sub_pd(iy0,jy2);
1338             dz02             = _mm_sub_pd(iz0,jz2);
1339             dx10             = _mm_sub_pd(ix1,jx0);
1340             dy10             = _mm_sub_pd(iy1,jy0);
1341             dz10             = _mm_sub_pd(iz1,jz0);
1342             dx11             = _mm_sub_pd(ix1,jx1);
1343             dy11             = _mm_sub_pd(iy1,jy1);
1344             dz11             = _mm_sub_pd(iz1,jz1);
1345             dx12             = _mm_sub_pd(ix1,jx2);
1346             dy12             = _mm_sub_pd(iy1,jy2);
1347             dz12             = _mm_sub_pd(iz1,jz2);
1348             dx20             = _mm_sub_pd(ix2,jx0);
1349             dy20             = _mm_sub_pd(iy2,jy0);
1350             dz20             = _mm_sub_pd(iz2,jz0);
1351             dx21             = _mm_sub_pd(ix2,jx1);
1352             dy21             = _mm_sub_pd(iy2,jy1);
1353             dz21             = _mm_sub_pd(iz2,jz1);
1354             dx22             = _mm_sub_pd(ix2,jx2);
1355             dy22             = _mm_sub_pd(iy2,jy2);
1356             dz22             = _mm_sub_pd(iz2,jz2);
1357
1358             /* Calculate squared distance and things based on it */
1359             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
1360             rsq01            = gmx_mm_calc_rsq_pd(dx01,dy01,dz01);
1361             rsq02            = gmx_mm_calc_rsq_pd(dx02,dy02,dz02);
1362             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
1363             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
1364             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
1365             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
1366             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
1367             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
1368
1369             rinv00           = gmx_mm_invsqrt_pd(rsq00);
1370             rinv01           = gmx_mm_invsqrt_pd(rsq01);
1371             rinv02           = gmx_mm_invsqrt_pd(rsq02);
1372             rinv10           = gmx_mm_invsqrt_pd(rsq10);
1373             rinv11           = gmx_mm_invsqrt_pd(rsq11);
1374             rinv12           = gmx_mm_invsqrt_pd(rsq12);
1375             rinv20           = gmx_mm_invsqrt_pd(rsq20);
1376             rinv21           = gmx_mm_invsqrt_pd(rsq21);
1377             rinv22           = gmx_mm_invsqrt_pd(rsq22);
1378
1379             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
1380             rinvsq01         = _mm_mul_pd(rinv01,rinv01);
1381             rinvsq02         = _mm_mul_pd(rinv02,rinv02);
1382             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
1383             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
1384             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
1385             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
1386             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
1387             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
1388
1389             fjx0             = _mm_setzero_pd();
1390             fjy0             = _mm_setzero_pd();
1391             fjz0             = _mm_setzero_pd();
1392             fjx1             = _mm_setzero_pd();
1393             fjy1             = _mm_setzero_pd();
1394             fjz1             = _mm_setzero_pd();
1395             fjx2             = _mm_setzero_pd();
1396             fjy2             = _mm_setzero_pd();
1397             fjz2             = _mm_setzero_pd();
1398
1399             /**************************
1400              * CALCULATE INTERACTIONS *
1401              **************************/
1402
1403             if (gmx_mm_any_lt(rsq00,rcutoff2))
1404             {
1405
1406             r00              = _mm_mul_pd(rsq00,rinv00);
1407
1408             /* REACTION-FIELD ELECTROSTATICS */
1409             felec            = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
1410
1411             /* LENNARD-JONES DISPERSION/REPULSION */
1412
1413             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1414             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
1415             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
1416             vvdw             = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
1417             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
1418
1419             d                = _mm_sub_pd(r00,rswitch);
1420             d                = _mm_max_pd(d,_mm_setzero_pd());
1421             d2               = _mm_mul_pd(d,d);
1422             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)))))));
1423
1424             dsw              = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
1425
1426             /* Evaluate switch function */
1427             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1428             fvdw             = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
1429             cutoff_mask      = _mm_cmplt_pd(rsq00,rcutoff2);
1430
1431             fscal            = _mm_add_pd(felec,fvdw);
1432
1433             fscal            = _mm_and_pd(fscal,cutoff_mask);
1434
1435             /* Calculate temporary vectorial force */
1436             tx               = _mm_mul_pd(fscal,dx00);
1437             ty               = _mm_mul_pd(fscal,dy00);
1438             tz               = _mm_mul_pd(fscal,dz00);
1439
1440             /* Update vectorial force */
1441             fix0             = _mm_add_pd(fix0,tx);
1442             fiy0             = _mm_add_pd(fiy0,ty);
1443             fiz0             = _mm_add_pd(fiz0,tz);
1444
1445             fjx0             = _mm_add_pd(fjx0,tx);
1446             fjy0             = _mm_add_pd(fjy0,ty);
1447             fjz0             = _mm_add_pd(fjz0,tz);
1448
1449             }
1450
1451             /**************************
1452              * CALCULATE INTERACTIONS *
1453              **************************/
1454
1455             if (gmx_mm_any_lt(rsq01,rcutoff2))
1456             {
1457
1458             /* REACTION-FIELD ELECTROSTATICS */
1459             felec            = _mm_mul_pd(qq01,_mm_sub_pd(_mm_mul_pd(rinv01,rinvsq01),krf2));
1460
1461             cutoff_mask      = _mm_cmplt_pd(rsq01,rcutoff2);
1462
1463             fscal            = felec;
1464
1465             fscal            = _mm_and_pd(fscal,cutoff_mask);
1466
1467             /* Calculate temporary vectorial force */
1468             tx               = _mm_mul_pd(fscal,dx01);
1469             ty               = _mm_mul_pd(fscal,dy01);
1470             tz               = _mm_mul_pd(fscal,dz01);
1471
1472             /* Update vectorial force */
1473             fix0             = _mm_add_pd(fix0,tx);
1474             fiy0             = _mm_add_pd(fiy0,ty);
1475             fiz0             = _mm_add_pd(fiz0,tz);
1476
1477             fjx1             = _mm_add_pd(fjx1,tx);
1478             fjy1             = _mm_add_pd(fjy1,ty);
1479             fjz1             = _mm_add_pd(fjz1,tz);
1480
1481             }
1482
1483             /**************************
1484              * CALCULATE INTERACTIONS *
1485              **************************/
1486
1487             if (gmx_mm_any_lt(rsq02,rcutoff2))
1488             {
1489
1490             /* REACTION-FIELD ELECTROSTATICS */
1491             felec            = _mm_mul_pd(qq02,_mm_sub_pd(_mm_mul_pd(rinv02,rinvsq02),krf2));
1492
1493             cutoff_mask      = _mm_cmplt_pd(rsq02,rcutoff2);
1494
1495             fscal            = felec;
1496
1497             fscal            = _mm_and_pd(fscal,cutoff_mask);
1498
1499             /* Calculate temporary vectorial force */
1500             tx               = _mm_mul_pd(fscal,dx02);
1501             ty               = _mm_mul_pd(fscal,dy02);
1502             tz               = _mm_mul_pd(fscal,dz02);
1503
1504             /* Update vectorial force */
1505             fix0             = _mm_add_pd(fix0,tx);
1506             fiy0             = _mm_add_pd(fiy0,ty);
1507             fiz0             = _mm_add_pd(fiz0,tz);
1508
1509             fjx2             = _mm_add_pd(fjx2,tx);
1510             fjy2             = _mm_add_pd(fjy2,ty);
1511             fjz2             = _mm_add_pd(fjz2,tz);
1512
1513             }
1514
1515             /**************************
1516              * CALCULATE INTERACTIONS *
1517              **************************/
1518
1519             if (gmx_mm_any_lt(rsq10,rcutoff2))
1520             {
1521
1522             /* REACTION-FIELD ELECTROSTATICS */
1523             felec            = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
1524
1525             cutoff_mask      = _mm_cmplt_pd(rsq10,rcutoff2);
1526
1527             fscal            = felec;
1528
1529             fscal            = _mm_and_pd(fscal,cutoff_mask);
1530
1531             /* Calculate temporary vectorial force */
1532             tx               = _mm_mul_pd(fscal,dx10);
1533             ty               = _mm_mul_pd(fscal,dy10);
1534             tz               = _mm_mul_pd(fscal,dz10);
1535
1536             /* Update vectorial force */
1537             fix1             = _mm_add_pd(fix1,tx);
1538             fiy1             = _mm_add_pd(fiy1,ty);
1539             fiz1             = _mm_add_pd(fiz1,tz);
1540
1541             fjx0             = _mm_add_pd(fjx0,tx);
1542             fjy0             = _mm_add_pd(fjy0,ty);
1543             fjz0             = _mm_add_pd(fjz0,tz);
1544
1545             }
1546
1547             /**************************
1548              * CALCULATE INTERACTIONS *
1549              **************************/
1550
1551             if (gmx_mm_any_lt(rsq11,rcutoff2))
1552             {
1553
1554             /* REACTION-FIELD ELECTROSTATICS */
1555             felec            = _mm_mul_pd(qq11,_mm_sub_pd(_mm_mul_pd(rinv11,rinvsq11),krf2));
1556
1557             cutoff_mask      = _mm_cmplt_pd(rsq11,rcutoff2);
1558
1559             fscal            = felec;
1560
1561             fscal            = _mm_and_pd(fscal,cutoff_mask);
1562
1563             /* Calculate temporary vectorial force */
1564             tx               = _mm_mul_pd(fscal,dx11);
1565             ty               = _mm_mul_pd(fscal,dy11);
1566             tz               = _mm_mul_pd(fscal,dz11);
1567
1568             /* Update vectorial force */
1569             fix1             = _mm_add_pd(fix1,tx);
1570             fiy1             = _mm_add_pd(fiy1,ty);
1571             fiz1             = _mm_add_pd(fiz1,tz);
1572
1573             fjx1             = _mm_add_pd(fjx1,tx);
1574             fjy1             = _mm_add_pd(fjy1,ty);
1575             fjz1             = _mm_add_pd(fjz1,tz);
1576
1577             }
1578
1579             /**************************
1580              * CALCULATE INTERACTIONS *
1581              **************************/
1582
1583             if (gmx_mm_any_lt(rsq12,rcutoff2))
1584             {
1585
1586             /* REACTION-FIELD ELECTROSTATICS */
1587             felec            = _mm_mul_pd(qq12,_mm_sub_pd(_mm_mul_pd(rinv12,rinvsq12),krf2));
1588
1589             cutoff_mask      = _mm_cmplt_pd(rsq12,rcutoff2);
1590
1591             fscal            = felec;
1592
1593             fscal            = _mm_and_pd(fscal,cutoff_mask);
1594
1595             /* Calculate temporary vectorial force */
1596             tx               = _mm_mul_pd(fscal,dx12);
1597             ty               = _mm_mul_pd(fscal,dy12);
1598             tz               = _mm_mul_pd(fscal,dz12);
1599
1600             /* Update vectorial force */
1601             fix1             = _mm_add_pd(fix1,tx);
1602             fiy1             = _mm_add_pd(fiy1,ty);
1603             fiz1             = _mm_add_pd(fiz1,tz);
1604
1605             fjx2             = _mm_add_pd(fjx2,tx);
1606             fjy2             = _mm_add_pd(fjy2,ty);
1607             fjz2             = _mm_add_pd(fjz2,tz);
1608
1609             }
1610
1611             /**************************
1612              * CALCULATE INTERACTIONS *
1613              **************************/
1614
1615             if (gmx_mm_any_lt(rsq20,rcutoff2))
1616             {
1617
1618             /* REACTION-FIELD ELECTROSTATICS */
1619             felec            = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
1620
1621             cutoff_mask      = _mm_cmplt_pd(rsq20,rcutoff2);
1622
1623             fscal            = felec;
1624
1625             fscal            = _mm_and_pd(fscal,cutoff_mask);
1626
1627             /* Calculate temporary vectorial force */
1628             tx               = _mm_mul_pd(fscal,dx20);
1629             ty               = _mm_mul_pd(fscal,dy20);
1630             tz               = _mm_mul_pd(fscal,dz20);
1631
1632             /* Update vectorial force */
1633             fix2             = _mm_add_pd(fix2,tx);
1634             fiy2             = _mm_add_pd(fiy2,ty);
1635             fiz2             = _mm_add_pd(fiz2,tz);
1636
1637             fjx0             = _mm_add_pd(fjx0,tx);
1638             fjy0             = _mm_add_pd(fjy0,ty);
1639             fjz0             = _mm_add_pd(fjz0,tz);
1640
1641             }
1642
1643             /**************************
1644              * CALCULATE INTERACTIONS *
1645              **************************/
1646
1647             if (gmx_mm_any_lt(rsq21,rcutoff2))
1648             {
1649
1650             /* REACTION-FIELD ELECTROSTATICS */
1651             felec            = _mm_mul_pd(qq21,_mm_sub_pd(_mm_mul_pd(rinv21,rinvsq21),krf2));
1652
1653             cutoff_mask      = _mm_cmplt_pd(rsq21,rcutoff2);
1654
1655             fscal            = felec;
1656
1657             fscal            = _mm_and_pd(fscal,cutoff_mask);
1658
1659             /* Calculate temporary vectorial force */
1660             tx               = _mm_mul_pd(fscal,dx21);
1661             ty               = _mm_mul_pd(fscal,dy21);
1662             tz               = _mm_mul_pd(fscal,dz21);
1663
1664             /* Update vectorial force */
1665             fix2             = _mm_add_pd(fix2,tx);
1666             fiy2             = _mm_add_pd(fiy2,ty);
1667             fiz2             = _mm_add_pd(fiz2,tz);
1668
1669             fjx1             = _mm_add_pd(fjx1,tx);
1670             fjy1             = _mm_add_pd(fjy1,ty);
1671             fjz1             = _mm_add_pd(fjz1,tz);
1672
1673             }
1674
1675             /**************************
1676              * CALCULATE INTERACTIONS *
1677              **************************/
1678
1679             if (gmx_mm_any_lt(rsq22,rcutoff2))
1680             {
1681
1682             /* REACTION-FIELD ELECTROSTATICS */
1683             felec            = _mm_mul_pd(qq22,_mm_sub_pd(_mm_mul_pd(rinv22,rinvsq22),krf2));
1684
1685             cutoff_mask      = _mm_cmplt_pd(rsq22,rcutoff2);
1686
1687             fscal            = felec;
1688
1689             fscal            = _mm_and_pd(fscal,cutoff_mask);
1690
1691             /* Calculate temporary vectorial force */
1692             tx               = _mm_mul_pd(fscal,dx22);
1693             ty               = _mm_mul_pd(fscal,dy22);
1694             tz               = _mm_mul_pd(fscal,dz22);
1695
1696             /* Update vectorial force */
1697             fix2             = _mm_add_pd(fix2,tx);
1698             fiy2             = _mm_add_pd(fiy2,ty);
1699             fiz2             = _mm_add_pd(fiz2,tz);
1700
1701             fjx2             = _mm_add_pd(fjx2,tx);
1702             fjy2             = _mm_add_pd(fjy2,ty);
1703             fjz2             = _mm_add_pd(fjz2,tz);
1704
1705             }
1706
1707             gmx_mm_decrement_3rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1708
1709             /* Inner loop uses 301 flops */
1710         }
1711
1712         if(jidx<j_index_end)
1713         {
1714
1715             jnrA             = jjnr[jidx];
1716             j_coord_offsetA  = DIM*jnrA;
1717
1718             /* load j atom coordinates */
1719             gmx_mm_load_3rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
1720                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1721
1722             /* Calculate displacement vector */
1723             dx00             = _mm_sub_pd(ix0,jx0);
1724             dy00             = _mm_sub_pd(iy0,jy0);
1725             dz00             = _mm_sub_pd(iz0,jz0);
1726             dx01             = _mm_sub_pd(ix0,jx1);
1727             dy01             = _mm_sub_pd(iy0,jy1);
1728             dz01             = _mm_sub_pd(iz0,jz1);
1729             dx02             = _mm_sub_pd(ix0,jx2);
1730             dy02             = _mm_sub_pd(iy0,jy2);
1731             dz02             = _mm_sub_pd(iz0,jz2);
1732             dx10             = _mm_sub_pd(ix1,jx0);
1733             dy10             = _mm_sub_pd(iy1,jy0);
1734             dz10             = _mm_sub_pd(iz1,jz0);
1735             dx11             = _mm_sub_pd(ix1,jx1);
1736             dy11             = _mm_sub_pd(iy1,jy1);
1737             dz11             = _mm_sub_pd(iz1,jz1);
1738             dx12             = _mm_sub_pd(ix1,jx2);
1739             dy12             = _mm_sub_pd(iy1,jy2);
1740             dz12             = _mm_sub_pd(iz1,jz2);
1741             dx20             = _mm_sub_pd(ix2,jx0);
1742             dy20             = _mm_sub_pd(iy2,jy0);
1743             dz20             = _mm_sub_pd(iz2,jz0);
1744             dx21             = _mm_sub_pd(ix2,jx1);
1745             dy21             = _mm_sub_pd(iy2,jy1);
1746             dz21             = _mm_sub_pd(iz2,jz1);
1747             dx22             = _mm_sub_pd(ix2,jx2);
1748             dy22             = _mm_sub_pd(iy2,jy2);
1749             dz22             = _mm_sub_pd(iz2,jz2);
1750
1751             /* Calculate squared distance and things based on it */
1752             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
1753             rsq01            = gmx_mm_calc_rsq_pd(dx01,dy01,dz01);
1754             rsq02            = gmx_mm_calc_rsq_pd(dx02,dy02,dz02);
1755             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
1756             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
1757             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
1758             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
1759             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
1760             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
1761
1762             rinv00           = gmx_mm_invsqrt_pd(rsq00);
1763             rinv01           = gmx_mm_invsqrt_pd(rsq01);
1764             rinv02           = gmx_mm_invsqrt_pd(rsq02);
1765             rinv10           = gmx_mm_invsqrt_pd(rsq10);
1766             rinv11           = gmx_mm_invsqrt_pd(rsq11);
1767             rinv12           = gmx_mm_invsqrt_pd(rsq12);
1768             rinv20           = gmx_mm_invsqrt_pd(rsq20);
1769             rinv21           = gmx_mm_invsqrt_pd(rsq21);
1770             rinv22           = gmx_mm_invsqrt_pd(rsq22);
1771
1772             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
1773             rinvsq01         = _mm_mul_pd(rinv01,rinv01);
1774             rinvsq02         = _mm_mul_pd(rinv02,rinv02);
1775             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
1776             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
1777             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
1778             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
1779             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
1780             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
1781
1782             fjx0             = _mm_setzero_pd();
1783             fjy0             = _mm_setzero_pd();
1784             fjz0             = _mm_setzero_pd();
1785             fjx1             = _mm_setzero_pd();
1786             fjy1             = _mm_setzero_pd();
1787             fjz1             = _mm_setzero_pd();
1788             fjx2             = _mm_setzero_pd();
1789             fjy2             = _mm_setzero_pd();
1790             fjz2             = _mm_setzero_pd();
1791
1792             /**************************
1793              * CALCULATE INTERACTIONS *
1794              **************************/
1795
1796             if (gmx_mm_any_lt(rsq00,rcutoff2))
1797             {
1798
1799             r00              = _mm_mul_pd(rsq00,rinv00);
1800
1801             /* REACTION-FIELD ELECTROSTATICS */
1802             felec            = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
1803
1804             /* LENNARD-JONES DISPERSION/REPULSION */
1805
1806             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1807             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
1808             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
1809             vvdw             = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
1810             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
1811
1812             d                = _mm_sub_pd(r00,rswitch);
1813             d                = _mm_max_pd(d,_mm_setzero_pd());
1814             d2               = _mm_mul_pd(d,d);
1815             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)))))));
1816
1817             dsw              = _mm_mul_pd(d2,_mm_add_pd(swF2,_mm_mul_pd(d,_mm_add_pd(swF3,_mm_mul_pd(d,swF4)))));
1818
1819             /* Evaluate switch function */
1820             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1821             fvdw             = _mm_sub_pd( _mm_mul_pd(fvdw,sw) , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
1822             cutoff_mask      = _mm_cmplt_pd(rsq00,rcutoff2);
1823
1824             fscal            = _mm_add_pd(felec,fvdw);
1825
1826             fscal            = _mm_and_pd(fscal,cutoff_mask);
1827
1828             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1829
1830             /* Calculate temporary vectorial force */
1831             tx               = _mm_mul_pd(fscal,dx00);
1832             ty               = _mm_mul_pd(fscal,dy00);
1833             tz               = _mm_mul_pd(fscal,dz00);
1834
1835             /* Update vectorial force */
1836             fix0             = _mm_add_pd(fix0,tx);
1837             fiy0             = _mm_add_pd(fiy0,ty);
1838             fiz0             = _mm_add_pd(fiz0,tz);
1839
1840             fjx0             = _mm_add_pd(fjx0,tx);
1841             fjy0             = _mm_add_pd(fjy0,ty);
1842             fjz0             = _mm_add_pd(fjz0,tz);
1843
1844             }
1845
1846             /**************************
1847              * CALCULATE INTERACTIONS *
1848              **************************/
1849
1850             if (gmx_mm_any_lt(rsq01,rcutoff2))
1851             {
1852
1853             /* REACTION-FIELD ELECTROSTATICS */
1854             felec            = _mm_mul_pd(qq01,_mm_sub_pd(_mm_mul_pd(rinv01,rinvsq01),krf2));
1855
1856             cutoff_mask      = _mm_cmplt_pd(rsq01,rcutoff2);
1857
1858             fscal            = felec;
1859
1860             fscal            = _mm_and_pd(fscal,cutoff_mask);
1861
1862             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1863
1864             /* Calculate temporary vectorial force */
1865             tx               = _mm_mul_pd(fscal,dx01);
1866             ty               = _mm_mul_pd(fscal,dy01);
1867             tz               = _mm_mul_pd(fscal,dz01);
1868
1869             /* Update vectorial force */
1870             fix0             = _mm_add_pd(fix0,tx);
1871             fiy0             = _mm_add_pd(fiy0,ty);
1872             fiz0             = _mm_add_pd(fiz0,tz);
1873
1874             fjx1             = _mm_add_pd(fjx1,tx);
1875             fjy1             = _mm_add_pd(fjy1,ty);
1876             fjz1             = _mm_add_pd(fjz1,tz);
1877
1878             }
1879
1880             /**************************
1881              * CALCULATE INTERACTIONS *
1882              **************************/
1883
1884             if (gmx_mm_any_lt(rsq02,rcutoff2))
1885             {
1886
1887             /* REACTION-FIELD ELECTROSTATICS */
1888             felec            = _mm_mul_pd(qq02,_mm_sub_pd(_mm_mul_pd(rinv02,rinvsq02),krf2));
1889
1890             cutoff_mask      = _mm_cmplt_pd(rsq02,rcutoff2);
1891
1892             fscal            = felec;
1893
1894             fscal            = _mm_and_pd(fscal,cutoff_mask);
1895
1896             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1897
1898             /* Calculate temporary vectorial force */
1899             tx               = _mm_mul_pd(fscal,dx02);
1900             ty               = _mm_mul_pd(fscal,dy02);
1901             tz               = _mm_mul_pd(fscal,dz02);
1902
1903             /* Update vectorial force */
1904             fix0             = _mm_add_pd(fix0,tx);
1905             fiy0             = _mm_add_pd(fiy0,ty);
1906             fiz0             = _mm_add_pd(fiz0,tz);
1907
1908             fjx2             = _mm_add_pd(fjx2,tx);
1909             fjy2             = _mm_add_pd(fjy2,ty);
1910             fjz2             = _mm_add_pd(fjz2,tz);
1911
1912             }
1913
1914             /**************************
1915              * CALCULATE INTERACTIONS *
1916              **************************/
1917
1918             if (gmx_mm_any_lt(rsq10,rcutoff2))
1919             {
1920
1921             /* REACTION-FIELD ELECTROSTATICS */
1922             felec            = _mm_mul_pd(qq10,_mm_sub_pd(_mm_mul_pd(rinv10,rinvsq10),krf2));
1923
1924             cutoff_mask      = _mm_cmplt_pd(rsq10,rcutoff2);
1925
1926             fscal            = felec;
1927
1928             fscal            = _mm_and_pd(fscal,cutoff_mask);
1929
1930             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1931
1932             /* Calculate temporary vectorial force */
1933             tx               = _mm_mul_pd(fscal,dx10);
1934             ty               = _mm_mul_pd(fscal,dy10);
1935             tz               = _mm_mul_pd(fscal,dz10);
1936
1937             /* Update vectorial force */
1938             fix1             = _mm_add_pd(fix1,tx);
1939             fiy1             = _mm_add_pd(fiy1,ty);
1940             fiz1             = _mm_add_pd(fiz1,tz);
1941
1942             fjx0             = _mm_add_pd(fjx0,tx);
1943             fjy0             = _mm_add_pd(fjy0,ty);
1944             fjz0             = _mm_add_pd(fjz0,tz);
1945
1946             }
1947
1948             /**************************
1949              * CALCULATE INTERACTIONS *
1950              **************************/
1951
1952             if (gmx_mm_any_lt(rsq11,rcutoff2))
1953             {
1954
1955             /* REACTION-FIELD ELECTROSTATICS */
1956             felec            = _mm_mul_pd(qq11,_mm_sub_pd(_mm_mul_pd(rinv11,rinvsq11),krf2));
1957
1958             cutoff_mask      = _mm_cmplt_pd(rsq11,rcutoff2);
1959
1960             fscal            = felec;
1961
1962             fscal            = _mm_and_pd(fscal,cutoff_mask);
1963
1964             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1965
1966             /* Calculate temporary vectorial force */
1967             tx               = _mm_mul_pd(fscal,dx11);
1968             ty               = _mm_mul_pd(fscal,dy11);
1969             tz               = _mm_mul_pd(fscal,dz11);
1970
1971             /* Update vectorial force */
1972             fix1             = _mm_add_pd(fix1,tx);
1973             fiy1             = _mm_add_pd(fiy1,ty);
1974             fiz1             = _mm_add_pd(fiz1,tz);
1975
1976             fjx1             = _mm_add_pd(fjx1,tx);
1977             fjy1             = _mm_add_pd(fjy1,ty);
1978             fjz1             = _mm_add_pd(fjz1,tz);
1979
1980             }
1981
1982             /**************************
1983              * CALCULATE INTERACTIONS *
1984              **************************/
1985
1986             if (gmx_mm_any_lt(rsq12,rcutoff2))
1987             {
1988
1989             /* REACTION-FIELD ELECTROSTATICS */
1990             felec            = _mm_mul_pd(qq12,_mm_sub_pd(_mm_mul_pd(rinv12,rinvsq12),krf2));
1991
1992             cutoff_mask      = _mm_cmplt_pd(rsq12,rcutoff2);
1993
1994             fscal            = felec;
1995
1996             fscal            = _mm_and_pd(fscal,cutoff_mask);
1997
1998             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1999
2000             /* Calculate temporary vectorial force */
2001             tx               = _mm_mul_pd(fscal,dx12);
2002             ty               = _mm_mul_pd(fscal,dy12);
2003             tz               = _mm_mul_pd(fscal,dz12);
2004
2005             /* Update vectorial force */
2006             fix1             = _mm_add_pd(fix1,tx);
2007             fiy1             = _mm_add_pd(fiy1,ty);
2008             fiz1             = _mm_add_pd(fiz1,tz);
2009
2010             fjx2             = _mm_add_pd(fjx2,tx);
2011             fjy2             = _mm_add_pd(fjy2,ty);
2012             fjz2             = _mm_add_pd(fjz2,tz);
2013
2014             }
2015
2016             /**************************
2017              * CALCULATE INTERACTIONS *
2018              **************************/
2019
2020             if (gmx_mm_any_lt(rsq20,rcutoff2))
2021             {
2022
2023             /* REACTION-FIELD ELECTROSTATICS */
2024             felec            = _mm_mul_pd(qq20,_mm_sub_pd(_mm_mul_pd(rinv20,rinvsq20),krf2));
2025
2026             cutoff_mask      = _mm_cmplt_pd(rsq20,rcutoff2);
2027
2028             fscal            = felec;
2029
2030             fscal            = _mm_and_pd(fscal,cutoff_mask);
2031
2032             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2033
2034             /* Calculate temporary vectorial force */
2035             tx               = _mm_mul_pd(fscal,dx20);
2036             ty               = _mm_mul_pd(fscal,dy20);
2037             tz               = _mm_mul_pd(fscal,dz20);
2038
2039             /* Update vectorial force */
2040             fix2             = _mm_add_pd(fix2,tx);
2041             fiy2             = _mm_add_pd(fiy2,ty);
2042             fiz2             = _mm_add_pd(fiz2,tz);
2043
2044             fjx0             = _mm_add_pd(fjx0,tx);
2045             fjy0             = _mm_add_pd(fjy0,ty);
2046             fjz0             = _mm_add_pd(fjz0,tz);
2047
2048             }
2049
2050             /**************************
2051              * CALCULATE INTERACTIONS *
2052              **************************/
2053
2054             if (gmx_mm_any_lt(rsq21,rcutoff2))
2055             {
2056
2057             /* REACTION-FIELD ELECTROSTATICS */
2058             felec            = _mm_mul_pd(qq21,_mm_sub_pd(_mm_mul_pd(rinv21,rinvsq21),krf2));
2059
2060             cutoff_mask      = _mm_cmplt_pd(rsq21,rcutoff2);
2061
2062             fscal            = felec;
2063
2064             fscal            = _mm_and_pd(fscal,cutoff_mask);
2065
2066             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2067
2068             /* Calculate temporary vectorial force */
2069             tx               = _mm_mul_pd(fscal,dx21);
2070             ty               = _mm_mul_pd(fscal,dy21);
2071             tz               = _mm_mul_pd(fscal,dz21);
2072
2073             /* Update vectorial force */
2074             fix2             = _mm_add_pd(fix2,tx);
2075             fiy2             = _mm_add_pd(fiy2,ty);
2076             fiz2             = _mm_add_pd(fiz2,tz);
2077
2078             fjx1             = _mm_add_pd(fjx1,tx);
2079             fjy1             = _mm_add_pd(fjy1,ty);
2080             fjz1             = _mm_add_pd(fjz1,tz);
2081
2082             }
2083
2084             /**************************
2085              * CALCULATE INTERACTIONS *
2086              **************************/
2087
2088             if (gmx_mm_any_lt(rsq22,rcutoff2))
2089             {
2090
2091             /* REACTION-FIELD ELECTROSTATICS */
2092             felec            = _mm_mul_pd(qq22,_mm_sub_pd(_mm_mul_pd(rinv22,rinvsq22),krf2));
2093
2094             cutoff_mask      = _mm_cmplt_pd(rsq22,rcutoff2);
2095
2096             fscal            = felec;
2097
2098             fscal            = _mm_and_pd(fscal,cutoff_mask);
2099
2100             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2101
2102             /* Calculate temporary vectorial force */
2103             tx               = _mm_mul_pd(fscal,dx22);
2104             ty               = _mm_mul_pd(fscal,dy22);
2105             tz               = _mm_mul_pd(fscal,dz22);
2106
2107             /* Update vectorial force */
2108             fix2             = _mm_add_pd(fix2,tx);
2109             fiy2             = _mm_add_pd(fiy2,ty);
2110             fiz2             = _mm_add_pd(fiz2,tz);
2111
2112             fjx2             = _mm_add_pd(fjx2,tx);
2113             fjy2             = _mm_add_pd(fjy2,ty);
2114             fjz2             = _mm_add_pd(fjz2,tz);
2115
2116             }
2117
2118             gmx_mm_decrement_3rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
2119
2120             /* Inner loop uses 301 flops */
2121         }
2122
2123         /* End of innermost loop */
2124
2125         gmx_mm_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
2126                                               f+i_coord_offset,fshift+i_shift_offset);
2127
2128         /* Increment number of inner iterations */
2129         inneriter                  += j_index_end - j_index_start;
2130
2131         /* Outer loop uses 18 flops */
2132     }
2133
2134     /* Increment number of outer iterations */
2135     outeriter        += nri;
2136
2137     /* Update outer/inner flops */
2138
2139     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_F,outeriter*18 + inneriter*301);
2140 }