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 247 by chuckv, Mon Jan 27 18:28:11 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 :: setRcut
57  
58    interface wrap
59       module procedure wrap_1d
# Line 40 | Line 66 | module simulation
66    end interface
67  
68  
43 !MPI dependent routines
69  
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
70  
51  real( kind = dp ), allocatable, dimension(:,:) :: fRow
52  real( kind = dp ), allocatable, dimension(:,:) :: fColumn
71  
54  real( kind = dp ), allocatable, dimension(:,:) :: tRow
55  real( kind = dp ), allocatable, dimension(:,:) :: tColumn
72  
57 #endif
73  
74  
75  
76  
62
63
77   contains
78  
79 < #ifdef MPI
80 < ! Allocated work arrays for MPI
81 <  subroutine allocate_mpi_arrays(nDimensions,numComponents)
82 <    integer, intent(in) :: nDimensions
83 <    integer, intent(in) :: numComponents
79 >  subroutine setSimulation(nLRParticles,box,rlist,rcut)
80 >    integer, intent(in) :: nLRParticles
81 >    real(kind = dp ), intent(in), dimension(3) :: box
82 >    real(kind = dp ), intent(in) :: rlist
83 >    real(kind = dp ), intent(in) :: rcut
84  
85 +    if( setsim ) return  ! simulation is already initialized
86 +    setSim = .true.
87  
88 +    thisSim%nLRParticles = nLRParticles
89 +    thisSim%box          = box
90 +    thisSim%rlist        = rlist
91 +    thisSim%rcut         = rcut
92 +    thisSim%rcutsq       = rcut * rcut
93 +    thisSim%rcut6        = thisSim%rcutsq * thisSim%rcutsq * thisSim%rcutsq
94  
95 +  end subroutine setSimulation
96  
97 +  function getNparticles() result(nparticles)
98 +    integer :: nparticles
99 +    nparticles = thisSim%nLRparticles
100 +  end function getNparticles
101  
76  end subroutine allocate_mpi_arrays
77 #endif
102  
79  subroutine set_simulation(box,rlist,rcut)
80    
81
82
83  end subroutine set_simulation
84
85
86
103    subroutine change_box_size(new_box_size)
104      real(kind=dp), dimension(3) :: new_box_size
105  
106 <    box = new_box_size
106 >    thisSim%box = new_box_size
107  
108    end subroutine change_box_size
109  
110  
111 <  elemental function getBox_3d() result(thisBox)
111 >  function getBox_3d() result(thisBox)
112      real( kind = dp ), dimension(3) :: thisBox
113 <    thisBox = box
113 >    thisBox = thisSim%box
114    end function getBox_3d
115  
116    function getBox_dim(dim) result(thisBox)
117      integer, intent(in) :: dim
118      real( kind = dp ) :: thisBox
119  
120 <    thisBox = box(dim)
120 >    thisBox = thisSim%box(dim)
121    end function getBox_dim
122      
123 <  subroutine check(update_nlist)
124 <  
123 >  subroutine check(q,update_nlist)
124 >    real( kind = dp ), dimension(:,:)  :: q
125      integer :: i
126      real( kind = DP ) :: dispmx
127      logical, intent(out) :: update_nlist
128      real( kind = DP ) :: dispmx_tmp
129 <    
129 >    real( kind = dp ) :: skin_thickness
130 >    integer :: natoms
131 >
132 >    natoms = thisSim%nLRparticles
133 >    skin_thickness = thisSim%rlist
134      dispmx = 0.0E0_DP
115    
135      !! calculate the largest displacement of any atom in any direction
136      
137   #ifdef MPI
# Line 139 | Line 158 | contains
158      
159    end subroutine check
160    
161 <  subroutine save_nlist()
162 <    integer :: i
163 < #ifdef MPI
164 <    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
161 >  subroutine save_nlist(q)
162 >    real(kind = dp ), dimension(:,:), intent(in)  :: q
163 >    integer :: list_size
164 >    
165  
166 +    list_size = size(q)
167 +
168 +    if (.not. allocated(q0)) then
169 +       allocate(q0(3,list_size))
170 +    else if( list_size > size(q0)) then
171 +       deallocate(q0)
172 +       allocate(q0(3,list_size))
173 +    endif
174 +
175 +    q0 = q
176 +
177    end subroutine save_nlist
178    
179  
# Line 163 | Line 186 | contains
186      
187      if (use_pbc) then
188         !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
189 <       this_wrap = r - box(dim)*nint(r/box(dim))
189 >       this_wrap = r - thisSim%box(dim)*nint(r/thisSim%box(dim))
190      else
191         this_wrap = r
192      endif
# Line 171 | Line 194 | contains
194      return
195    end function wrap_1d
196  
197 <  elemental function wrap_3d(r) result(this_wrap)
197 >  function wrap_3d(r) result(this_wrap)
198      real( kind = dp ), dimension(3), intent(in) :: r
199      real( kind = dp ), dimension(3) :: this_wrap
200  
201      
202      if (use_pbc) then
203         !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
204 <       this_wrap(1:3) = r(1:3) - box(1:3)*nint(r(1:3)/box(1:3))
204 >       this_wrap = r - thisSim%box*nint(r/thisSim%box)
205      else
206         this_wrap = r
207      endif
208    end function wrap_3d
209  
210    
211 +
212 +  subroutine getRcut(thisrcut,rcut2,rcut6,status)
213 +    real( kind = dp ), intent(out) :: thisrcut
214 +    real( kind = dp ), intent(out), optional :: rcut2
215 +    real( kind = dp ), intent(out), optional :: rcut6
216 +    integer, optional :: status
217 +
218 +    if (present(status)) status = 0
219 +
220 +    if (.not.setSim ) then
221 +       if (present(status)) status = -1
222 +       return
223 +    end if
224 +    
225 +    thisrcut = thisSim%rcut
226 +    if(present(rcut2)) rcut2 = thisSim%rcutsq
227 +    if(present(rcut2)) rcut6 = thisSim%rcut6
228 +
229 +  end subroutine getRcut
230 +  
231 +  
232 +  
233 +
234 +  subroutine getRlist(thisrlist,rlist2,status)
235 +    real( kind = dp ), intent(out) :: thisrlist
236 +    real( kind = dp ), intent(out), optional :: rlist2
237 +
238 +    integer, optional :: status
239 +
240 +    if (present(status)) status = 0
241 +
242 +    if (.not.setSim ) then
243 +       if (present(status)) status = -1
244 +       return
245 +    end if
246 +    
247 +    thisrlist = thisSim%rlist
248 +    if(present(rlist2)) rlist2 = thisSim%rlistsq
249 +
250 +  end subroutine getRlist
251 +  
252 +  
253 +
254 + pure function getNlocal() result(nlocal)
255 +    integer :: nlocal
256 +    nlocal = thisSim%nLRparticles
257 +  end function getNlocal
258 +
259 +
260 +
261   end module simulation

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines