ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/simulation_module.F90
(Generate patch)

Comparing trunk/mdtools/md_code/simulation_module.F90 (file contents):
Revision 252 by chuckv, Tue Jan 28 22:16:55 2003 UTC vs.
Revision 264 by chuckv, Tue Feb 4 20:16:08 2003 UTC

# Line 89 | Line 89 | contains
89      thisSim%nLRParticles = nLRParticles
90      thisSim%box          = box
91      thisSim%rlist        = rlist
92 +    thisSIm%rlistsq      = rlist * rlist
93      thisSim%rcut         = rcut
94      thisSim%rcutsq       = rcut * rcut
95      thisSim%rcut6        = thisSim%rcutsq * thisSim%rcutsq * thisSim%rcutsq
# Line 134 | Line 135 | contains
135      integer :: natoms
136  
137      natoms = thisSim%nLRparticles
138 <    skin_thickness = thisSim%rlist
138 >    skin_thickness = thisSim%rlist - thisSim%rcut
139      dispmx = 0.0E0_DP
140      !! calculate the largest displacement of any atom in any direction
141      
# Line 252 | Line 253 | contains
253      thisrlist = thisSim%rlist
254      if(present(rlist2)) rlist2 = thisSim%rlistsq
255  
256 +
257    end subroutine getRlist
258    
259    

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines