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 185 by chuckv, Fri Nov 22 20:39:33 2002 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 + #ifdef IS_MPI
4 +  use mpiSimulation
5 + #endif
6 +
7    implicit none
8    PRIVATE
9  
10  
11  
12 <  real ( kind = dp ), dimension(3) :: box
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 +  end type simtype
31  
32 +  type (simtype), public :: thisSim
33 + !! Tag for MPI calculations  
34 +  integer, allocatable, dimension(:) :: tag
35  
36 + #ifdef IS_MPI
37 +  integer, allocatable, dimension(:) :: tag_row
38 +  integer, allocatable, dimension(:) :: tag_column
39 + #endif
40  
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  
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
60 +     module procedure wrap_3d
61 +  end interface
62 +
63 +  interface getBox
64 +     module procedure getBox_3d
65 +     module procedure getBox_dim
66 +  end interface
67 +
68 +
69 +
70 +
71 +
72 +
73 +
74 +
75 +
76 +
77   contains
78  
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 <  subroutine set_simulation(box,rlist,rcut)
86 <    
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 set_simulation
95 >  end subroutine setSimulation
96  
97 +  function getNparticles() result(nparticles)
98 +    integer :: nparticles
99 +    nparticles = thisSim%nLRparticles
100 +  end function getNparticles
101  
102  
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 +  function getBox_3d() result(thisBox)
112 +    real( kind = dp ), dimension(3) :: thisBox
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 = thisSim%box(dim)
121 +  end function getBox_dim
122 +    
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 +    real( kind = dp ) :: skin_thickness
130 +    integer :: natoms
131 +
132 +    natoms = thisSim%nLRparticles
133 +    skin_thickness = thisSim%rlist
134 +    dispmx = 0.0E0_DP
135 +    !! calculate the largest displacement of any atom in any direction
136 +    
137 + #ifdef MPI
138 +    dispmx_tmp = 0.0E0_DP
139 +    do i = 1, nlocal
140 +       dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx )
141 +       dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx )
142 +       dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx )
143 +    end do
144 +    call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_double_precision, &
145 +       mpi_max,mpi_comm_world,mpi_err)
146 + #else
147 +    do i = 1, natoms
148 +       dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx )
149 +       dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx )
150 +       dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx )
151 +    end do
152 + #endif
153 +  
154 +    !! a conservative test of list skin crossings
155 +    dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx)
156 +    
157 +    update_nlist = (dispmx.gt.(skin_thickness))
158 +    
159 +  end subroutine check
160 +  
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 +
180 +  function wrap_1d(r,dim) result(this_wrap)
181 +    
182 +    
183 +    real( kind = DP ) :: r
184 +    real( kind = DP ) :: this_wrap
185 +    integer           :: dim
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 - thisSim%box(dim)*nint(r/thisSim%box(dim))
190 +    else
191 +       this_wrap = r
192 +    endif
193 +    
194 +    return
195 +  end function wrap_1d
196 +
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 = 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