ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/simulation_module.F90
Revision: 367
Committed: Wed Mar 19 17:29:49 2003 UTC (21 years, 5 months ago) by mmeineke
File size: 9756 byte(s)
Log Message:
*** empty log message ***

File Contents

# User Rev Content
1 chuckv 290 !! Fortran interface to C entry plug.
2    
3 mmeineke 270 module simulation
4 gezelter 325 use definitions
5     use neighborLists
6 gezelter 360 use force_globals
7 gezelter 358 use vector_class
8 gezelter 328 use atype_module
9 gezelter 360 use lj
10 mmeineke 270 #ifdef IS_MPI
11     use mpiSimulation
12     #endif
13    
14     implicit none
15     PRIVATE
16    
17 mmeineke 285 #define __FORTRAN90
18 gezelter 309 #include "fSimulation.h"
19 mmeineke 270
20     type (simtype), public :: thisSim
21    
22 gezelter 325 logical, save :: simulation_setup_complete = .false.
23 mmeineke 270
24 gezelter 361 integer, public, save :: natoms
25 gezelter 330 integer, public, save :: nExcludes_Global = 0
26     integer, public, save :: nExcludes_Local = 0
27 gezelter 360 integer, allocatable, dimension(:,:), public :: excludesLocal
28     integer, allocatable, dimension(:), public :: excludesGlobal
29 gezelter 332
30 gezelter 325 real(kind=dp), save :: rcut2 = 0.0_DP
31     real(kind=dp), save :: rcut6 = 0.0_DP
32     real(kind=dp), save :: rlist2 = 0.0_DP
33 gezelter 330 real(kind=dp), public, dimension(3), save :: box
34 mmeineke 270
35 gezelter 325
36     public :: SimulationSetup
37 gezelter 361 public :: getNlocal
38 gezelter 360 public :: setBox
39 mmeineke 367 public :: setBox_3d
40 mmeineke 270 public :: getBox
41 gezelter 360 public :: setRcut
42 mmeineke 270 public :: getRcut
43     public :: getRlist
44 gezelter 312 public :: getRrf
45     public :: getRt
46 gezelter 329 public :: getDielect
47 gezelter 325 public :: SimUsesPBC
48     public :: SimUsesLJ
49     public :: SimUsesDipoles
50     public :: SimUsesSticky
51     public :: SimUsesRF
52     public :: SimUsesGB
53     public :: SimUsesEAM
54     public :: SimRequiresPrepairCalc
55     public :: SimRequiresPostpairCalc
56 gezelter 332 public :: SimUsesDirectionalAtoms
57 mmeineke 285
58 mmeineke 270 interface getBox
59     module procedure getBox_3d
60 gezelter 360 module procedure getBox_1d
61 mmeineke 270 end interface
62 gezelter 325
63 gezelter 360 interface setBox
64     module procedure setBox_3d
65     module procedure setBox_1d
66     end interface
67    
68 mmeineke 270 contains
69 gezelter 325
70     subroutine SimulationSetup(setThisSim, nComponents, c_idents, &
71 gezelter 332 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
72     status)
73 chuckv 290 type (simtype) :: setThisSim
74 gezelter 325 integer, intent(inout) :: nComponents
75     integer, dimension(nComponents),intent(inout) :: c_idents
76 gezelter 332
77     integer :: CnLocalExcludes
78     integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
79     integer :: CnGlobalExcludes
80     integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
81 gezelter 325 !! Result status, success = 0, status = -1
82     integer, intent(out) :: status
83     integer :: i, me, thisStat, alloc_stat, myNode
84     #ifdef IS_MPI
85     integer, allocatable, dimension(:) :: c_idents_Row
86     integer, allocatable, dimension(:) :: c_idents_Col
87     integer :: nrow
88     integer :: ncol
89 gezelter 332 #endif
90 gezelter 325
91     simulation_setup_complete = .false.
92     status = 0
93    
94 gezelter 312 ! copy C struct into fortran type
95 chuckv 290 thisSim = setThisSim
96 gezelter 325 natoms = nComponents
97     rcut2 = thisSim%rcut * thisSim%rcut
98     rcut6 = rcut2 * rcut2 * rcut2
99     rlist2 = thisSim%rlist * thisSim%rlist
100 gezelter 330 box = thisSim%box
101 gezelter 332 nExcludes_Global = CnGlobalExcludes
102     nExcludes_Local = CnLocalExcludes
103 mmeineke 270
104 gezelter 360 call InitializeForceGlobals(natoms, thisStat)
105 gezelter 332 if (thisStat /= 0) then
106     status = -1
107     return
108     endif
109    
110 gezelter 360 call InitializeSimGlobals(thisStat)
111     if (thisStat /= 0) then
112     status = -1
113     return
114     endif
115    
116 gezelter 325 #ifdef IS_MPI
117     ! We can only set up forces if mpiSimulation has been setup.
118     if (.not. isMPISimSet()) then
119     write(default_error,*) "MPI is not set"
120     status = -1
121     return
122     endif
123     nrow = getNrow(plan_row)
124     ncol = getNcol(plan_col)
125     mynode = getMyNode()
126    
127     allocate(c_idents_Row(nrow),stat=alloc_stat)
128     if (alloc_stat /= 0 ) then
129     status = -1
130     return
131     endif
132 mmeineke 270
133 gezelter 325 allocate(c_idents_Col(ncol),stat=alloc_stat)
134     if (alloc_stat /= 0 ) then
135     status = -1
136     return
137     endif
138    
139     call gather(c_idents, c_idents_Row, plan_row)
140     call gather(c_idents, c_idents_Col, plan_col)
141    
142     do i = 1, nrow
143     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
144     atid_Row(i) = me
145     enddo
146    
147     do i = 1, ncol
148     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
149     atid_Col(i) = me
150     enddo
151    
152     !! free temporary ident arrays
153     if (allocated(c_idents_Col)) then
154     deallocate(c_idents_Col)
155     end if
156     if (allocated(c_idents_Row)) then
157     deallocate(c_idents_Row)
158     endif
159    
160     #else
161     do i = 1, nComponents
162    
163     me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
164     atid(i) = me
165    
166     enddo
167     #endif
168    
169     !! Create neighbor lists
170     call expandNeighborList(nComponents, thisStat)
171     if (thisStat /= 0) then
172     status = -1
173     return
174     endif
175 gezelter 332
176     do i = 1, nExcludes_Local
177     excludesLocal(1,i) = CexcludesLocal(1,i)
178     excludesLocal(2,i) = CexcludesLocal(2,i)
179     enddo
180 gezelter 325
181 gezelter 332 do i = 1, nExcludes_Global
182     excludesGlobal(i) = CexcludesGlobal(i)
183     enddo
184    
185 gezelter 325 if (status == 0) simulation_setup_complete = .true.
186 gezelter 332
187 gezelter 325 end subroutine SimulationSetup
188    
189 gezelter 360 subroutine setBox_3d(new_box_size)
190 mmeineke 270 real(kind=dp), dimension(3) :: new_box_size
191 gezelter 360 integer :: smallest, status, i
192    
193 mmeineke 270 thisSim%box = new_box_size
194 gezelter 330 box = thisSim%box
195 mmeineke 270
196 gezelter 360 smallest = 1
197     do i = 2, 3
198     if (new_box_size(i) .lt. new_box_size(smallest)) smallest = i
199     end do
200     if (thisSim%rcut .gt. 0.5_dp * new_box_size(smallest)) &
201     call setRcut(0.5_dp * new_box_size(smallest), status)
202     return
203     end subroutine setBox_3d
204    
205     subroutine setBox_1d(dim, new_box_size)
206     integer :: dim, status
207     real(kind=dp) :: new_box_size
208     thisSim%box(dim) = new_box_size
209     box(dim) = thisSim%box(dim)
210     if (thisSim%rcut .gt. 0.5_dp * new_box_size) &
211     call setRcut(0.5_dp * new_box_size, status)
212     end subroutine setBox_1d
213    
214     subroutine setRcut(new_rcut, status)
215     real(kind = dp) :: new_rcut
216     integer :: myStatus, status
217     thisSim%rcut = new_rcut
218     rcut2 = thisSim%rcut * thisSim%rcut
219     rcut6 = rcut2 * rcut2 * rcut2
220     myStatus = 0
221     call LJ_new_rcut(new_rcut, myStatus)
222     if (myStatus .ne. 0) then
223     write(default_error, *) 'LJ module refused our rcut!'
224     status = -1
225     return
226     endif
227     status = 0
228     return
229     end subroutine setRcut
230    
231 mmeineke 270 function getBox_3d() result(thisBox)
232     real( kind = dp ), dimension(3) :: thisBox
233     thisBox = thisSim%box
234     end function getBox_3d
235 gezelter 360
236     function getBox_1d(dim) result(thisBox)
237 mmeineke 270 integer, intent(in) :: dim
238     real( kind = dp ) :: thisBox
239 gezelter 312
240 mmeineke 270 thisBox = thisSim%box(dim)
241 gezelter 360 end function getBox_1d
242 gezelter 312
243 gezelter 325 subroutine getRcut(thisrcut,rc2,rc6,status)
244 mmeineke 270 real( kind = dp ), intent(out) :: thisrcut
245 gezelter 325 real( kind = dp ), intent(out), optional :: rc2
246     real( kind = dp ), intent(out), optional :: rc6
247 mmeineke 270 integer, optional :: status
248    
249     if (present(status)) status = 0
250 gezelter 325
251     if (.not.simulation_setup_complete ) then
252 mmeineke 270 if (present(status)) status = -1
253     return
254     end if
255    
256     thisrcut = thisSim%rcut
257 gezelter 325 if(present(rc2)) rc2 = rcut2
258     if(present(rc6)) rc6 = rcut6
259 mmeineke 270 end subroutine getRcut
260    
261 gezelter 325 subroutine getRlist(thisrlist,rl2,status)
262 mmeineke 270 real( kind = dp ), intent(out) :: thisrlist
263 gezelter 325 real( kind = dp ), intent(out), optional :: rl2
264 mmeineke 270
265     integer, optional :: status
266    
267     if (present(status)) status = 0
268    
269 gezelter 325 if (.not.simulation_setup_complete ) then
270 mmeineke 270 if (present(status)) status = -1
271     return
272     end if
273    
274     thisrlist = thisSim%rlist
275 gezelter 325 if(present(rl2)) rl2 = rlist2
276 gezelter 312 end subroutine getRlist
277 mmeineke 270
278 gezelter 312 function getRrf() result(rrf)
279     real( kind = dp ) :: rrf
280     rrf = thisSim%rrf
281     end function getRrf
282 mmeineke 270
283 gezelter 312 function getRt() result(rt)
284     real( kind = dp ) :: rt
285     rt = thisSim%rt
286     end function getRt
287 gezelter 329
288     function getDielect() result(dielect)
289     real( kind = dp ) :: dielect
290     dielect = thisSim%dielect
291     end function getDielect
292 gezelter 360
293 gezelter 325 function SimUsesPBC() result(doesit)
294     logical :: doesit
295     doesit = thisSim%SIM_uses_PBC
296     end function SimUsesPBC
297    
298     function SimUsesLJ() result(doesit)
299     logical :: doesit
300     doesit = thisSim%SIM_uses_LJ
301     end function SimUsesLJ
302    
303     function SimUsesSticky() result(doesit)
304     logical :: doesit
305     doesit = thisSim%SIM_uses_sticky
306     end function SimUsesSticky
307    
308     function SimUsesDipoles() result(doesit)
309     logical :: doesit
310     doesit = thisSim%SIM_uses_dipoles
311     end function SimUsesDipoles
312    
313     function SimUsesRF() result(doesit)
314     logical :: doesit
315     doesit = thisSim%SIM_uses_RF
316     end function SimUsesRF
317    
318     function SimUsesGB() result(doesit)
319     logical :: doesit
320     doesit = thisSim%SIM_uses_GB
321     end function SimUsesGB
322    
323     function SimUsesEAM() result(doesit)
324     logical :: doesit
325     doesit = thisSim%SIM_uses_EAM
326     end function SimUsesEAM
327    
328 gezelter 329 function SimUsesDirectionalAtoms() result(doesit)
329     logical :: doesit
330     doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
331     thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
332 gezelter 330 end function SimUsesDirectionalAtoms
333 gezelter 329
334 gezelter 325 function SimRequiresPrepairCalc() result(doesit)
335     logical :: doesit
336     doesit = thisSim%SIM_uses_EAM
337     end function SimRequiresPrepairCalc
338    
339     function SimRequiresPostpairCalc() result(doesit)
340     logical :: doesit
341     doesit = thisSim%SIM_uses_RF
342     end function SimRequiresPostpairCalc
343 gezelter 312
344 gezelter 360 subroutine InitializeSimGlobals(thisStat)
345 gezelter 325 integer, intent(out) :: thisStat
346     integer :: alloc_stat
347    
348     thisStat = 0
349    
350 gezelter 360 call FreeSimGlobals()
351 gezelter 325
352 gezelter 360 allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
353 gezelter 325 if (alloc_stat /= 0 ) then
354 gezelter 332 thisStat = -1
355 gezelter 325 return
356     endif
357    
358 gezelter 360 allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
359 gezelter 325 if (alloc_stat /= 0 ) then
360 gezelter 332 thisStat = -1
361 gezelter 325 return
362     endif
363    
364 gezelter 360 end subroutine InitializeSimGlobals
365 gezelter 312
366 gezelter 360 subroutine FreeSimGlobals()
367 gezelter 325
368     !We free in the opposite order in which we allocate in.
369 gezelter 332
370     if (allocated(excludesGlobal)) deallocate(excludesGlobal)
371     if (allocated(excludesLocal)) deallocate(excludesLocal)
372 gezelter 360
373     end subroutine FreeSimGlobals
374 gezelter 361
375     pure function getNlocal() result(nlocal)
376     integer :: nlocal
377     nlocal = natoms
378     end function getNlocal
379    
380 gezelter 309
381 mmeineke 270 end module simulation