ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 480
Committed: Tue Apr 8 17:16:22 2003 UTC (21 years, 3 months ago) by chuckv
File size: 9242 byte(s)
Log Message:
Moved expand neighborlist to init_FF.

File Contents

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