2 C This source code is part of
6 C Copyright (c) 1991-2000, University of Groningen, The Netherlands.
7 C Copyright (c) 2001-2009, The GROMACS Development Team
9 C Gromacs is a library for molecular simulation and trajectory analysis,
10 C written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
11 C a full list of developers and information, check out http://www.gromacs.org
13 C This program is free software; you can redistribute it and/or modify it under
14 C the terms of the GNU Lesser General Public License as published by the Free
15 C Software Foundation; either version 2 of the License, or (at your option) any
17 C As a special exception, you may use this file as part of a free software
18 C library without restriction. Specifically, if other files instantiate
19 C templates or use macros or inline functions from this file, or you compile
20 C this file and link it with other files to produce an executable, this
21 C file does not by itself cause the resulting executable to be covered by
22 C the GNU Lesser General Public License.
24 C In plain-speak: do not worry about classes/macros/templates either - only
25 C changes to the library have to be LGPL, not an application linking with it.
27 C To help fund GROMACS development, we humbly ask that you cite
28 C the papers people have written on it - you can find them on the website!
36 # define gmxreal real*8
38 # define gmxreal real*4
44 C Gromacs nonbonded kernel pwr6kernel020
45 C Coulomb interaction: Not calculated
46 C VdW interaction: Buckingham
47 C water optimization: No
48 C Calculate forces: yes
50 subroutine pwr6kernel020(
83 integer*4 nri,iinr(*),jindex(*),jjnr(*),shift(*)
84 gmxreal shiftvec(*),fshift(*),pos(*),faction(*)
85 integer*4 gid(*),type(*),ntype
86 gmxreal charge(*),facel,krf,crf,Vc(*),vdwparam(*)
87 gmxreal Vvdw(*),tabscale,VFtab(*)
88 gmxreal invsqrta(*),dvda(*),gbtabscale,GBtab(*)
89 integer*4 nthreads,count,mtx,outeriter,inneriter
92 integer*4 n,ii,is3,ii3,k,nj0,nj1,jnr,j3,ggid
93 integer*4 nn0,nn1,nouter,ninner
95 gmxreal fscal,tx,ty,tz
100 gmxreal Vvdw6,Vvdwtot
102 gmxreal ix1,iy1,iz1,fix1,fiy1,fiz1
104 gmxreal dx11,dy11,dz11,rsq11,rinv11
105 gmxreal c6,cexp1,cexp2
108 C Reset outer and inner iteration counters
112 C Loop over thread workunits
113 10 call pwr6kernelsync(mtx,count,nri,nthreads,nn0,nn1)
114 if(nn1.gt.nri) nn1=nri
116 C Start outer loop over neighborlists
120 C Load shift vector for this list
123 shY = shiftvec(is3+1)
124 shZ = shiftvec(is3+2)
126 C Load limits for loop over neighbors
130 C Get outer coordinate index
134 C Load i atom data, add shift vector
135 ix1 = shX + pos(ii3+0)
136 iy1 = shY + pos(ii3+1)
137 iz1 = shZ + pos(ii3+2)
139 C Load parameters for i atom
140 nti = 3*ntype*type(ii)
142 C Zero the potential energy for this list
145 C Clear i atom forces
152 C Get j neighbor index, and coordinate index
156 C load j atom coordinates
165 rsq11 = dx11*dx11+dy11*dy11+dz11*dz11
167 C Calculate 1/r and 1/r2
169 C PowerPC intrinsics 1/sqrt lookup table
171 rinv11 = frsqrtes(rsq11)
173 rinv11 = frsqrte(dble(rsq11))
175 rinv11 = (0.5*rinv11*(3.0-((rsq11*rinv11)
178 rinv11 = (0.5*rinv11*(3.0-((rsq11*rinv11)
182 C Load parameters for j atom
183 tj = nti+3*type(jnr)+1
185 cexp1 = vdwparam(tj+1)
186 cexp2 = vdwparam(tj+2)
187 rinvsq = rinv11*rinv11
189 C Buckingham interaction
190 rinvsix = rinvsq*rinvsq*rinvsq
192 br = cexp2*rsq11*rinv11
193 Vvdwexp = cexp1*exp(-br)
194 Vvdwtot = Vvdwtot+Vvdwexp-Vvdw6
195 fscal = (br*Vvdwexp-6.0*Vvdw6)*rinvsq
197 C Calculate temporary vectorial force
202 C Increment i atom force
207 C Decrement j atom force
208 faction(j3+0) = faction(j3+0) - tx
209 faction(j3+1) = faction(j3+1) - ty
210 faction(j3+2) = faction(j3+2) - tz
212 C Inner loop uses 62 flops/iteration
216 C Add i forces to mem and shifted force list
217 faction(ii3+0) = faction(ii3+0) + fix1
218 faction(ii3+1) = faction(ii3+1) + fiy1
219 faction(ii3+2) = faction(ii3+2) + fiz1
220 fshift(is3) = fshift(is3)+fix1
221 fshift(is3+1) = fshift(is3+1)+fiy1
222 fshift(is3+2) = fshift(is3+2)+fiz1
224 C Add potential energies to the group for this list
226 Vvdw(ggid) = Vvdw(ggid) + Vvdwtot
228 C Increment number of inner iterations
229 ninner = ninner + nj1 - nj0
231 C Outer loop uses 10 flops/iteration
235 C Increment number of outer iterations
236 nouter = nouter + nn1 - nn0
237 if(nn1.lt.nri) goto 10
239 C Write outer/inner iteration count to pointers
251 C Gromacs nonbonded kernel pwr6kernel020nf
252 C Coulomb interaction: Not calculated
253 C VdW interaction: Buckingham
254 C water optimization: No
255 C Calculate forces: no
257 subroutine pwr6kernel020nf(
290 integer*4 nri,iinr(*),jindex(*),jjnr(*),shift(*)
291 gmxreal shiftvec(*),fshift(*),pos(*),faction(*)
292 integer*4 gid(*),type(*),ntype
293 gmxreal charge(*),facel,krf,crf,Vc(*),vdwparam(*)
294 gmxreal Vvdw(*),tabscale,VFtab(*)
295 gmxreal invsqrta(*),dvda(*),gbtabscale,GBtab(*)
296 integer*4 nthreads,count,mtx,outeriter,inneriter
299 integer*4 n,ii,is3,ii3,k,nj0,nj1,jnr,j3,ggid
300 integer*4 nn0,nn1,nouter,ninner
306 gmxreal Vvdw6,Vvdwtot
310 gmxreal dx11,dy11,dz11,rsq11,rinv11
311 gmxreal c6,cexp1,cexp2
314 C Reset outer and inner iteration counters
318 C Loop over thread workunits
319 10 call pwr6kernelsync(mtx,count,nri,nthreads,nn0,nn1)
320 if(nn1.gt.nri) nn1=nri
322 C Start outer loop over neighborlists
326 C Load shift vector for this list
329 shY = shiftvec(is3+1)
330 shZ = shiftvec(is3+2)
332 C Load limits for loop over neighbors
336 C Get outer coordinate index
340 C Load i atom data, add shift vector
341 ix1 = shX + pos(ii3+0)
342 iy1 = shY + pos(ii3+1)
343 iz1 = shZ + pos(ii3+2)
345 C Load parameters for i atom
346 nti = 3*ntype*type(ii)
348 C Zero the potential energy for this list
351 C Clear i atom forces
355 C Get j neighbor index, and coordinate index
359 C load j atom coordinates
368 rsq11 = dx11*dx11+dy11*dy11+dz11*dz11
370 C Calculate 1/r and 1/r2
372 C PowerPC intrinsics 1/sqrt lookup table
374 rinv11 = frsqrtes(rsq11)
376 rinv11 = frsqrte(dble(rsq11))
378 rinv11 = (0.5*rinv11*(3.0-((rsq11*rinv11)
381 rinv11 = (0.5*rinv11*(3.0-((rsq11*rinv11)
385 C Load parameters for j atom
386 tj = nti+3*type(jnr)+1
388 cexp1 = vdwparam(tj+1)
389 cexp2 = vdwparam(tj+2)
390 rinvsq = rinv11*rinv11
392 C Buckingham interaction
393 rinvsix = rinvsq*rinvsq*rinvsq
395 br = cexp2*rsq11*rinv11
396 Vvdwexp = cexp1*exp(-br)
397 Vvdwtot = Vvdwtot+Vvdwexp-Vvdw6
399 C Inner loop uses 49 flops/iteration
403 C Add i forces to mem and shifted force list
405 C Add potential energies to the group for this list
407 Vvdw(ggid) = Vvdw(ggid) + Vvdwtot
409 C Increment number of inner iterations
410 ninner = ninner + nj1 - nj0
412 C Outer loop uses 4 flops/iteration
416 C Increment number of outer iterations
417 nouter = nouter + nn1 - nn0
418 if(nn1.lt.nri) goto 10
420 C Write outer/inner iteration count to pointers