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 222 by chuckv, Thu Jan 2 21:45:45 2003 UTC vs.
Revision 249 by chuckv, Mon Jan 27 21:28:19 2003 UTC

# Line 1 | Line 1
1   module simulation
2    use definitions, ONLY :dp
3 <  use force_wrappers, ONLY : alloc_force_wrappers
3 > #ifdef IS_MPI
4 >  use mpiSimulation
5 > #endif
6  
7    implicit none
8    PRIVATE
9  
10  
9  integer :: nLR
11  
12 +  type, public :: simtype
13 +     PRIVATE
14 + !     SEQUENCE
15 + !! Number of particles on this processor
16 +     integer :: nLRparticles
17 + !! Periodic Box    
18 +     real ( kind = dp ), dimension(3) :: box
19 + !! List Cutoff    
20 +     real ( kind = dp ) :: rlist = 0.0_dp
21 + !! Radial cutoff
22 +     real ( kind = dp ) :: rcut  = 0.0_dp
23 + !! List cutoff squared
24 +     real ( kind = dp ) :: rlistsq = 0.0_dp
25 + !! Radial Cutoff squared
26 +     real ( kind = dp ) :: rcutsq  = 0.0_dp
27 + !! Radial Cutoff^6
28 +     real ( kind = dp ) :: rcut6  = 0.0_dp
29  
30 <  real ( kind = dp ), dimension(3) :: box
30 >  end type simtype
31  
32 +  type (simtype), public :: thisSim
33 + !! Tag for MPI calculations  
34 +  integer, allocatable, dimension(:) :: tag
35  
36 <  real ( kind = dp ), public :: rlist
37 <  real ( kind = dp ), public :: rcut
38 <  real ( kind = dp ), public :: rlistsq
39 <  real ( kind = dp ), public :: rcutsq
36 > #ifdef IS_MPI
37 >  integer, allocatable, dimension(:) :: tag_row
38 >  integer, allocatable, dimension(:) :: tag_column
39 > #endif
40  
41 <  integer,public, allocatable, dimension(:)   :: point
42 <  integer,public, allocatable, dimension(:)   :: list
41 > !! WARNING: use_pbc hardcoded, fixme
42 >  logical :: use_pbc = .true.
43 >  logical :: setSim = .false.
44  
45 + !! array for saving previous positions for neighbor lists.  
46 +  real( kind = dp ), allocatable,dimension(:,:),save :: q0
47  
48  
25
26
49    public :: check
50    public :: save_nlist
51    public :: wrap
52    public :: getBox
53 +  public :: getRcut
54 +  public :: getRlist
55 +  public :: getNlocal
56 +  public :: setSimulation
57 + !  public :: setRcut
58  
59    interface wrap
60       module procedure wrap_1d
# Line 40 | Line 67 | module simulation
67    end interface
68  
69  
43 !MPI dependent routines
70  
45 #ifdef IS_MPI
46 ! Universal routines: All types of force calculations will need these arrays
47 ! Arrays specific to a type of force calculation should be declared in that module.
48  real( kind = dp ), allocatable, dimension(:,:) :: qRow
49  real( kind = dp ), allocatable, dimension(:,:) :: qColumn
71  
51  real( kind = dp ), allocatable, dimension(:,:) :: fRow
52  real( kind = dp ), allocatable, dimension(:,:) :: fColumn
72  
54  real( kind = dp ), allocatable, dimension(:,:) :: tRow
55  real( kind = dp ), allocatable, dimension(:,:) :: tColumn
73  
57 #endif
74  
75  
76  
77  
62
63
78   contains
79  
80 < #ifdef MPI
81 < ! Allocated work arrays for MPI
82 <  subroutine allocate_mpi_arrays(nDimensions,numComponents)
83 <    integer, intent(in) :: nDimensions
84 <    integer, intent(in) :: numComponents
80 >  subroutine setSimulation(nLRParticles,box,rlist,rcut)
81 >    integer, intent(in) :: nLRParticles
82 >    real(kind = dp ), intent(in), dimension(3) :: box
83 >    real(kind = dp ), intent(in) :: rlist
84 >    real(kind = dp ), intent(in) :: rcut
85  
86 +    if( setsim ) return  ! simulation is already initialized
87 +    setSim = .true.
88  
89 +    thisSim%nLRParticles = nLRParticles
90 +    thisSim%box          = box
91 +    thisSim%rlist        = rlist
92 +    thisSim%rcut         = rcut
93 +    thisSim%rcutsq       = rcut * rcut
94 +    thisSim%rcut6        = thisSim%rcutsq * thisSim%rcutsq * thisSim%rcutsq
95  
96 +  end subroutine setSimulation
97  
98 +  function getNparticles() result(nparticles)
99 +    integer :: nparticles
100 +    nparticles = thisSim%nLRparticles
101 +  end function getNparticles
102  
76  end subroutine allocate_mpi_arrays
77 #endif
103  
79  subroutine set_simulation(box,rlist,rcut)
80    
81
82
83  end subroutine set_simulation
84
85
86
104    subroutine change_box_size(new_box_size)
105      real(kind=dp), dimension(3) :: new_box_size
106  
107 <    box = new_box_size
107 >    thisSim%box = new_box_size
108  
109    end subroutine change_box_size
110  
111  
112 <  elemental function getBox_3d() result(thisBox)
112 >  function getBox_3d() result(thisBox)
113      real( kind = dp ), dimension(3) :: thisBox
114 <    thisBox = box
114 >    thisBox = thisSim%box
115    end function getBox_3d
116  
117    function getBox_dim(dim) result(thisBox)
118      integer, intent(in) :: dim
119      real( kind = dp ) :: thisBox
120  
121 <    thisBox = box(dim)
121 >    thisBox = thisSim%box(dim)
122    end function getBox_dim
123      
124 <  subroutine check(update_nlist)
125 <  
124 >  subroutine check(q,update_nlist)
125 >    real( kind = dp ), dimension(:,:)  :: q
126      integer :: i
127      real( kind = DP ) :: dispmx
128      logical, intent(out) :: update_nlist
129      real( kind = DP ) :: dispmx_tmp
130 <    
130 >    real( kind = dp ) :: skin_thickness
131 >    integer :: natoms
132 >
133 >    natoms = thisSim%nLRparticles
134 >    skin_thickness = thisSim%rlist
135      dispmx = 0.0E0_DP
115    
136      !! calculate the largest displacement of any atom in any direction
137      
138   #ifdef MPI
# Line 139 | Line 159 | contains
159      
160    end subroutine check
161    
162 <  subroutine save_nlist()
163 <    integer :: i
164 < #ifdef MPI
165 <    do i = 1, nlocal
146 < #else
147 <       do i = 1, natoms
148 < #endif
149 <          q0(1,i) = q(1,i)
150 <          q0(2,i) = q(2,i)
151 <          q0(3,i) = q(3,i)
152 <       end do
162 >  subroutine save_nlist(q)
163 >    real(kind = dp ), dimension(:,:), intent(in)  :: q
164 >    integer :: list_size
165 >    
166  
167 +    list_size = size(q)
168 +
169 +    if (.not. allocated(q0)) then
170 +       allocate(q0(3,list_size))
171 +    else if( list_size > size(q0)) then
172 +       deallocate(q0)
173 +       allocate(q0(3,list_size))
174 +    endif
175 +
176 +    q0 = q
177 +
178    end subroutine save_nlist
179    
180  
# Line 163 | Line 187 | contains
187      
188      if (use_pbc) then
189         !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
190 <       this_wrap = r - box(dim)*nint(r/box(dim))
190 >       this_wrap = r - thisSim%box(dim)*nint(r/thisSim%box(dim))
191      else
192         this_wrap = r
193      endif
# Line 171 | Line 195 | contains
195      return
196    end function wrap_1d
197  
198 <  elemental function wrap_3d(r) result(this_wrap)
198 >  function wrap_3d(r) result(this_wrap)
199      real( kind = dp ), dimension(3), intent(in) :: r
200      real( kind = dp ), dimension(3) :: this_wrap
201  
202      
203      if (use_pbc) then
204         !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
205 <       this_wrap(1:3) = r(1:3) - box(1:3)*nint(r(1:3)/box(1:3))
205 >       this_wrap = r - thisSim%box*nint(r/thisSim%box)
206      else
207         this_wrap = r
208      endif
209    end function wrap_3d
210  
211    
212 +
213 +  subroutine getRcut(thisrcut,rcut2,rcut6,status)
214 +    real( kind = dp ), intent(out) :: thisrcut
215 +    real( kind = dp ), intent(out), optional :: rcut2
216 +    real( kind = dp ), intent(out), optional :: rcut6
217 +    integer, optional :: status
218 +
219 +    if (present(status)) status = 0
220 +
221 +    if (.not.setSim ) then
222 +       if (present(status)) status = -1
223 +       return
224 +    end if
225 +    
226 +    thisrcut = thisSim%rcut
227 +    if(present(rcut2)) rcut2 = thisSim%rcutsq
228 +    if(present(rcut2)) rcut6 = thisSim%rcut6
229 +
230 +  end subroutine getRcut
231 +  
232 +  
233 +  
234 +
235 +  subroutine getRlist(thisrlist,rlist2,status)
236 +    real( kind = dp ), intent(out) :: thisrlist
237 +    real( kind = dp ), intent(out), optional :: rlist2
238 +
239 +    integer, optional :: status
240 +
241 +    if (present(status)) status = 0
242 +
243 +    if (.not.setSim ) then
244 +       if (present(status)) status = -1
245 +       return
246 +    end if
247 +    
248 +    thisrlist = thisSim%rlist
249 +    if(present(rlist2)) rlist2 = thisSim%rlistsq
250 +
251 +  end subroutine getRlist
252 +  
253 +  
254 +
255 + pure function getNlocal() result(nlocal)
256 +    integer :: nlocal
257 +    nlocal = thisSim%nLRparticles
258 +  end function getNlocal
259 +
260 +
261 +
262   end module simulation

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines