2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS c kernel generator.
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
48 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_VF_c
49 * Electrostatics interaction: ReactionField
50 * VdW interaction: LennardJones
51 * Geometry: Particle-Particle
52 * Calculate force/pot: PotentialAndForce
55 nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_VF_c
56 (t_nblist * gmx_restrict nlist,
57 rvec * gmx_restrict xx,
58 rvec * gmx_restrict ff,
59 t_forcerec * gmx_restrict fr,
60 t_mdatoms * gmx_restrict mdatoms,
61 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
62 t_nrnb * gmx_restrict nrnb)
64 int i_shift_offset,i_coord_offset,j_coord_offset;
65 int j_index_start,j_index_end;
66 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
67 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
68 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
69 real *shiftvec,*fshift,*x,*f;
71 real ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
73 real jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
74 real dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
75 real velec,felec,velecsum,facel,crf,krf,krf2;
78 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
87 jindex = nlist->jindex;
89 shiftidx = nlist->shift;
91 shiftvec = fr->shift_vec[0];
92 fshift = fr->fshift[0];
94 charge = mdatoms->chargeA;
100 vdwtype = mdatoms->typeA;
102 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
103 rcutoff = fr->rcoulomb;
104 rcutoff2 = rcutoff*rcutoff;
106 sh_vdw_invrcut6 = fr->ic->sh_invrc6;
112 /* Start outer loop over neighborlists */
113 for(iidx=0; iidx<nri; iidx++)
115 /* Load shift vector for this list */
116 i_shift_offset = DIM*shiftidx[iidx];
117 shX = shiftvec[i_shift_offset+XX];
118 shY = shiftvec[i_shift_offset+YY];
119 shZ = shiftvec[i_shift_offset+ZZ];
121 /* Load limits for loop over neighbors */
122 j_index_start = jindex[iidx];
123 j_index_end = jindex[iidx+1];
125 /* Get outer coordinate index */
127 i_coord_offset = DIM*inr;
129 /* Load i particle coords and add shift vector */
130 ix0 = shX + x[i_coord_offset+DIM*0+XX];
131 iy0 = shY + x[i_coord_offset+DIM*0+YY];
132 iz0 = shZ + x[i_coord_offset+DIM*0+ZZ];
138 /* Load parameters for i particles */
139 iq0 = facel*charge[inr+0];
140 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
142 /* Reset potential sums */
146 /* Start inner kernel loop */
147 for(jidx=j_index_start; jidx<j_index_end; jidx++)
149 /* Get j neighbor index, and coordinate index */
151 j_coord_offset = DIM*jnr;
153 /* load j atom coordinates */
154 jx0 = x[j_coord_offset+DIM*0+XX];
155 jy0 = x[j_coord_offset+DIM*0+YY];
156 jz0 = x[j_coord_offset+DIM*0+ZZ];
158 /* Calculate displacement vector */
163 /* Calculate squared distance and things based on it */
164 rsq00 = dx00*dx00+dy00*dy00+dz00*dz00;
166 rinv00 = gmx_invsqrt(rsq00);
168 rinvsq00 = rinv00*rinv00;
170 /* Load parameters for j particles */
172 vdwjidx0 = 2*vdwtype[jnr+0];
174 /**************************
175 * CALCULATE INTERACTIONS *
176 **************************/
182 c6_00 = vdwparam[vdwioffset0+vdwjidx0];
183 c12_00 = vdwparam[vdwioffset0+vdwjidx0+1];
185 /* REACTION-FIELD ELECTROSTATICS */
186 velec = qq00*(rinv00+krf*rsq00-crf);
187 felec = qq00*(rinv00*rinvsq00-krf2);
189 /* LENNARD-JONES DISPERSION/REPULSION */
191 rinvsix = rinvsq00*rinvsq00*rinvsq00;
192 vvdw6 = c6_00*rinvsix;
193 vvdw12 = c12_00*rinvsix*rinvsix;
194 vvdw = (vvdw12 - c12_00*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_00*sh_vdw_invrcut6)*(1.0/6.0);
195 fvdw = (vvdw12-vvdw6)*rinvsq00;
197 /* Update potential sums from outer loop */
203 /* Calculate temporary vectorial force */
208 /* Update vectorial force */
212 f[j_coord_offset+DIM*0+XX] -= tx;
213 f[j_coord_offset+DIM*0+YY] -= ty;
214 f[j_coord_offset+DIM*0+ZZ] -= tz;
218 /* Inner loop uses 49 flops */
220 /* End of innermost loop */
223 f[i_coord_offset+DIM*0+XX] += fix0;
224 f[i_coord_offset+DIM*0+YY] += fiy0;
225 f[i_coord_offset+DIM*0+ZZ] += fiz0;
229 fshift[i_shift_offset+XX] += tx;
230 fshift[i_shift_offset+YY] += ty;
231 fshift[i_shift_offset+ZZ] += tz;
234 /* Update potential energies */
235 kernel_data->energygrp_elec[ggid] += velecsum;
236 kernel_data->energygrp_vdw[ggid] += vvdwsum;
238 /* Increment number of inner iterations */
239 inneriter += j_index_end - j_index_start;
241 /* Outer loop uses 15 flops */
244 /* Increment number of outer iterations */
247 /* Update outer/inner flops */
249 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*15 + inneriter*49);
252 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_F_c
253 * Electrostatics interaction: ReactionField
254 * VdW interaction: LennardJones
255 * Geometry: Particle-Particle
256 * Calculate force/pot: Force
259 nb_kernel_ElecRFCut_VdwLJSh_GeomP1P1_F_c
260 (t_nblist * gmx_restrict nlist,
261 rvec * gmx_restrict xx,
262 rvec * gmx_restrict ff,
263 t_forcerec * gmx_restrict fr,
264 t_mdatoms * gmx_restrict mdatoms,
265 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
266 t_nrnb * gmx_restrict nrnb)
268 int i_shift_offset,i_coord_offset,j_coord_offset;
269 int j_index_start,j_index_end;
270 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
271 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
272 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
273 real *shiftvec,*fshift,*x,*f;
275 real ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
277 real jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
278 real dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
279 real velec,felec,velecsum,facel,crf,krf,krf2;
282 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
291 jindex = nlist->jindex;
293 shiftidx = nlist->shift;
295 shiftvec = fr->shift_vec[0];
296 fshift = fr->fshift[0];
298 charge = mdatoms->chargeA;
302 nvdwtype = fr->ntype;
304 vdwtype = mdatoms->typeA;
306 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
307 rcutoff = fr->rcoulomb;
308 rcutoff2 = rcutoff*rcutoff;
310 sh_vdw_invrcut6 = fr->ic->sh_invrc6;
316 /* Start outer loop over neighborlists */
317 for(iidx=0; iidx<nri; iidx++)
319 /* Load shift vector for this list */
320 i_shift_offset = DIM*shiftidx[iidx];
321 shX = shiftvec[i_shift_offset+XX];
322 shY = shiftvec[i_shift_offset+YY];
323 shZ = shiftvec[i_shift_offset+ZZ];
325 /* Load limits for loop over neighbors */
326 j_index_start = jindex[iidx];
327 j_index_end = jindex[iidx+1];
329 /* Get outer coordinate index */
331 i_coord_offset = DIM*inr;
333 /* Load i particle coords and add shift vector */
334 ix0 = shX + x[i_coord_offset+DIM*0+XX];
335 iy0 = shY + x[i_coord_offset+DIM*0+YY];
336 iz0 = shZ + x[i_coord_offset+DIM*0+ZZ];
342 /* Load parameters for i particles */
343 iq0 = facel*charge[inr+0];
344 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
346 /* Start inner kernel loop */
347 for(jidx=j_index_start; jidx<j_index_end; jidx++)
349 /* Get j neighbor index, and coordinate index */
351 j_coord_offset = DIM*jnr;
353 /* load j atom coordinates */
354 jx0 = x[j_coord_offset+DIM*0+XX];
355 jy0 = x[j_coord_offset+DIM*0+YY];
356 jz0 = x[j_coord_offset+DIM*0+ZZ];
358 /* Calculate displacement vector */
363 /* Calculate squared distance and things based on it */
364 rsq00 = dx00*dx00+dy00*dy00+dz00*dz00;
366 rinv00 = gmx_invsqrt(rsq00);
368 rinvsq00 = rinv00*rinv00;
370 /* Load parameters for j particles */
372 vdwjidx0 = 2*vdwtype[jnr+0];
374 /**************************
375 * CALCULATE INTERACTIONS *
376 **************************/
382 c6_00 = vdwparam[vdwioffset0+vdwjidx0];
383 c12_00 = vdwparam[vdwioffset0+vdwjidx0+1];
385 /* REACTION-FIELD ELECTROSTATICS */
386 felec = qq00*(rinv00*rinvsq00-krf2);
388 /* LENNARD-JONES DISPERSION/REPULSION */
390 rinvsix = rinvsq00*rinvsq00*rinvsq00;
391 fvdw = (c12_00*rinvsix-c6_00)*rinvsix*rinvsq00;
395 /* Calculate temporary vectorial force */
400 /* Update vectorial force */
404 f[j_coord_offset+DIM*0+XX] -= tx;
405 f[j_coord_offset+DIM*0+YY] -= ty;
406 f[j_coord_offset+DIM*0+ZZ] -= tz;
410 /* Inner loop uses 34 flops */
412 /* End of innermost loop */
415 f[i_coord_offset+DIM*0+XX] += fix0;
416 f[i_coord_offset+DIM*0+YY] += fiy0;
417 f[i_coord_offset+DIM*0+ZZ] += fiz0;
421 fshift[i_shift_offset+XX] += tx;
422 fshift[i_shift_offset+YY] += ty;
423 fshift[i_shift_offset+ZZ] += tz;
425 /* Increment number of inner iterations */
426 inneriter += j_index_end - j_index_start;
428 /* Outer loop uses 13 flops */
431 /* Increment number of outer iterations */
434 /* Update outer/inner flops */
436 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*13 + inneriter*34);