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 241 by mmeineke, Thu Jan 23 16:19:59 2003 UTC vs.
Revision 246 by chuckv, Fri Jan 24 21:36:52 2003 UTC

# Line 11 | Line 11 | module simulation
11  
12    type, public :: simtype
13       PRIVATE
14 <     SEQUENCE
14 > !     SEQUENCE
15   !! Number of particles on this processor
16       integer :: nLRparticles
17   !! Periodic Box    
# Line 38 | Line 38 | module simulation
38    integer, allocatable, dimension(:) :: tag_column
39   #endif
40  
41 <
41 > !! WARNING: use_pbc hardcoded, fixme
42 >  logical :: use_pbc = .true.
43    logical :: setSim = .false.
44    
45 <  real( kind = dp ), allocatable(:,:),save :: q0
45 >  real( kind = dp ), allocatable,dimension(:,:),save :: q0
46  
47    public :: check
48    public :: save_nlist
49    public :: wrap
50    public :: getBox
51    public :: getRcut
52 <  public :: setRcut
52 >  public :: getRlist
53 >  public :: getNlocal
54 > !  public :: setRcut
55  
56    interface wrap
57       module procedure wrap_1d
# Line 85 | Line 88 | contains
88      thisSim%rlist        = rlist
89      thisSim%rcut         = rcut
90      thisSim%rcutsq       = rcut * rcut
91 <    thisSim%rcut6        = rcutsq * rcutsq * rcutsq
91 >    thisSim%rcut6        = thisSim%rcutsq * thisSim%rcutsq * thisSim%rcutsq
92  
93    end subroutine setSimulation
94  
# Line 103 | Line 106 | contains
106    end subroutine change_box_size
107  
108  
109 <  elemental function getBox_3d() result(thisBox)
109 >  function getBox_3d() result(thisBox)
110      real( kind = dp ), dimension(3) :: thisBox
111      thisBox = thisSim%box
112    end function getBox_3d
# Line 115 | Line 118 | contains
118      thisBox = thisSim%box(dim)
119    end function getBox_dim
120      
121 <  subroutine check(update_nlist)
122 <  
121 >  subroutine check(q,update_nlist)
122 >    real( kind = dp ), dimension(:,:)  :: q
123      integer :: i
124      real( kind = DP ) :: dispmx
125      logical, intent(out) :: update_nlist
126      real( kind = DP ) :: dispmx_tmp
127 <    
127 >    real( kind = dp ) :: skin_thickness
128 >    integer :: natoms
129 >
130 >    natoms = thisSim%nLRparticles
131 >    skin_thickness = thisSim%rlist
132      dispmx = 0.0E0_DP
126    
133      !! calculate the largest displacement of any atom in any direction
134      
135   #ifdef MPI
# Line 159 | Line 165 | contains
165  
166      if (.not. allocated(q0)) then
167         allocate(q0(3,list_size))
168 <    else if( list_size > size(q0))
168 >    else if( list_size > size(q0)) then
169         deallocate(q0)
170         allocate(q0(3,list_size))
171      endif
# Line 186 | Line 192 | contains
192      return
193    end function wrap_1d
194  
195 <  elemental function wrap_3d(r) result(this_wrap)
195 >  function wrap_3d(r) result(this_wrap)
196      real( kind = dp ), dimension(3), intent(in) :: r
197      real( kind = dp ), dimension(3) :: this_wrap
198  
199      
200      if (use_pbc) then
201         !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
202 <       this_wrap(1:3) = r(1:3) - box(1:3)*nint(r(1:3)/box(1:3))
202 >       this_wrap = r - thisSim%box*nint(r/thisSim%box)
203      else
204         this_wrap = r
205      endif
206    end function wrap_3d
207  
208 <
209 <
208 >  
209 >
210    subroutine getRcut(thisrcut,rcut2,rcut6,status)
211      real( kind = dp ), intent(out) :: thisrcut
212      real( kind = dp ), intent(out), optional :: rcut2
213 <    real( kind = dp ), intent(out), optional :: thisrcut6
213 >    real( kind = dp ), intent(out), optional :: rcut6
214      integer, optional :: status
215  
216      if (present(status)) status = 0
# Line 214 | Line 220 | contains
220         return
221      end if
222      
223 <    thisrcut = rcut
224 <    if(present(rcut2)) rcut2 = rcutsq
225 <    if(present(rcut2)) rcut6 = rcut_6
223 >    thisrcut = thisSim%rcut
224 >    if(present(rcut2)) rcut2 = thisSim%rcutsq
225 >    if(present(rcut2)) rcut6 = thisSim%rcut6
226  
227    end subroutine getRcut
228    
229 +  
230 +  
231 +
232 +  subroutine getRlist(thisrlist,rlist2,status)
233 +    real( kind = dp ), intent(out) :: thisrlist
234 +    real( kind = dp ), intent(out), optional :: rlist2
235  
236 +    integer, optional :: status
237  
238 +    if (present(status)) status = 0
239  
240 +    if (.not.setSim ) then
241 +       if (present(status)) status = -1
242 +       return
243 +    end if
244 +    
245 +    thisrlist = thisSim%rlist
246 +    if(present(rlist2)) rlist2 = thisSim%rlistsq
247 +
248 +  end subroutine getRlist
249 +  
250 +  
251 +
252 + pure function getNlocal() result(nlocal)
253 +    integer :: nlocal
254 +    nlocal = thisSim%nLRparticles
255 +  end function getNlocal
256 +
257 +
258 +
259   end module simulation

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines