ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 460
Committed: Fri Apr 4 22:22:30 2003 UTC (21 years, 3 months ago) by chuckv
File size: 9253 byte(s)
Log Message:
Breaking c and fortran, c gets smarter, fortran gets dumber...

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     status = -1
108     return
109     endif
110    
111     call InitializeSimGlobals(thisStat)
112     if (thisStat /= 0) then
113     status = -1
114     return
115     endif
116    
117     #ifdef IS_MPI
118     ! We can only set up forces if mpiSimulation has been setup.
119     if (.not. isMPISimSet()) then
120     write(default_error,*) "MPI is not set"
121     status = -1
122     return
123     endif
124     nrow = getNrow(plan_row)
125     ncol = getNcol(plan_col)
126     mynode = getMyNode()
127    
128     allocate(c_idents_Row(nrow),stat=alloc_stat)
129     if (alloc_stat /= 0 ) then
130     status = -1
131     return
132     endif
133    
134     allocate(c_idents_Col(ncol),stat=alloc_stat)
135     if (alloc_stat /= 0 ) then
136     status = -1
137     return
138     endif
139    
140     call gather(c_idents, c_idents_Row, plan_row)
141     call gather(c_idents, c_idents_Col, plan_col)
142    
143     do i = 1, nrow
144     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
145     atid_Row(i) = me
146     enddo
147    
148     do i = 1, ncol
149     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
150     atid_Col(i) = me
151     enddo
152    
153     !! free temporary ident arrays
154     if (allocated(c_idents_Col)) then
155     deallocate(c_idents_Col)
156     end if
157     if (allocated(c_idents_Row)) then
158     deallocate(c_idents_Row)
159     endif
160    
161     #else
162     do i = 1, nComponents
163    
164     me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
165     atid(i) = me
166    
167     enddo
168     #endif
169    
170     !! Create neighbor lists
171     call expandNeighborList(nComponents, thisStat)
172     if (thisStat /= 0) then
173     status = -1
174     return
175     endif
176    
177 chuckv 388
178 mmeineke 377 do i = 1, nExcludes_Local
179     excludesLocal(1,i) = CexcludesLocal(1,i)
180     excludesLocal(2,i) = CexcludesLocal(2,i)
181     enddo
182    
183     do i = 1, nExcludes_Global
184     excludesGlobal(i) = CexcludesGlobal(i)
185     enddo
186 mmeineke 435
187 mmeineke 377 if (status == 0) simulation_setup_complete = .true.
188    
189     end subroutine SimulationSetup
190    
191     subroutine setBox_3d(new_box_size)
192     real(kind=dp), dimension(3) :: new_box_size
193     integer :: smallest, status, i
194    
195     thisSim%box = new_box_size
196     box = thisSim%box
197    
198     return
199     end subroutine setBox_3d
200    
201     subroutine setBox_1d(dim, new_box_size)
202     integer :: dim, status
203     real(kind=dp) :: new_box_size
204     thisSim%box(dim) = new_box_size
205     box(dim) = thisSim%box(dim)
206     end subroutine setBox_1d
207    
208     subroutine setRcut(new_rcut, status)
209     real(kind = dp) :: new_rcut
210     integer :: myStatus, status
211     thisSim%rcut = new_rcut
212     rcut2 = thisSim%rcut * thisSim%rcut
213     rcut6 = rcut2 * rcut2 * rcut2
214     status = 0
215     return
216     end subroutine setRcut
217    
218     function getBox_3d() result(thisBox)
219     real( kind = dp ), dimension(3) :: thisBox
220     thisBox = thisSim%box
221     end function getBox_3d
222    
223     function getBox_1d(dim) result(thisBox)
224     integer, intent(in) :: dim
225     real( kind = dp ) :: thisBox
226    
227     thisBox = thisSim%box(dim)
228     end function getBox_1d
229    
230     subroutine getRcut(thisrcut,rc2,rc6,status)
231     real( kind = dp ), intent(out) :: thisrcut
232     real( kind = dp ), intent(out), optional :: rc2
233     real( kind = dp ), intent(out), optional :: rc6
234     integer, optional :: status
235    
236     if (present(status)) status = 0
237    
238     if (.not.simulation_setup_complete ) then
239     if (present(status)) status = -1
240     return
241     end if
242    
243     thisrcut = thisSim%rcut
244     if(present(rc2)) rc2 = rcut2
245     if(present(rc6)) rc6 = rcut6
246     end subroutine getRcut
247    
248     subroutine getRlist(thisrlist,rl2,status)
249     real( kind = dp ), intent(out) :: thisrlist
250     real( kind = dp ), intent(out), optional :: rl2
251    
252     integer, optional :: status
253    
254     if (present(status)) status = 0
255    
256     if (.not.simulation_setup_complete ) then
257     if (present(status)) status = -1
258     return
259     end if
260    
261     thisrlist = thisSim%rlist
262     if(present(rl2)) rl2 = rlist2
263     end subroutine getRlist
264    
265     function getRrf() result(rrf)
266     real( kind = dp ) :: rrf
267     rrf = thisSim%rrf
268 gezelter 394 write(*,*) 'getRrf = ', rrf, thisSim%rrf
269 mmeineke 377 end function getRrf
270    
271     function getRt() result(rt)
272     real( kind = dp ) :: rt
273     rt = thisSim%rt
274     end function getRt
275    
276     function getDielect() result(dielect)
277     real( kind = dp ) :: dielect
278     dielect = thisSim%dielect
279     end function getDielect
280    
281     function SimUsesPBC() result(doesit)
282     logical :: doesit
283     doesit = thisSim%SIM_uses_PBC
284     end function SimUsesPBC
285    
286     function SimUsesLJ() result(doesit)
287     logical :: doesit
288     doesit = thisSim%SIM_uses_LJ
289     end function SimUsesLJ
290    
291     function SimUsesSticky() result(doesit)
292     logical :: doesit
293     doesit = thisSim%SIM_uses_sticky
294     end function SimUsesSticky
295    
296     function SimUsesDipoles() result(doesit)
297     logical :: doesit
298     doesit = thisSim%SIM_uses_dipoles
299     end function SimUsesDipoles
300    
301     function SimUsesRF() result(doesit)
302     logical :: doesit
303     doesit = thisSim%SIM_uses_RF
304     end function SimUsesRF
305    
306     function SimUsesGB() result(doesit)
307     logical :: doesit
308     doesit = thisSim%SIM_uses_GB
309     end function SimUsesGB
310    
311     function SimUsesEAM() result(doesit)
312     logical :: doesit
313     doesit = thisSim%SIM_uses_EAM
314     end function SimUsesEAM
315    
316     function SimUsesDirectionalAtoms() result(doesit)
317     logical :: doesit
318     doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
319     thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
320     end function SimUsesDirectionalAtoms
321    
322     function SimRequiresPrepairCalc() result(doesit)
323     logical :: doesit
324     doesit = thisSim%SIM_uses_EAM
325     end function SimRequiresPrepairCalc
326    
327     function SimRequiresPostpairCalc() result(doesit)
328     logical :: doesit
329     doesit = thisSim%SIM_uses_RF
330     end function SimRequiresPostpairCalc
331    
332     subroutine InitializeSimGlobals(thisStat)
333     integer, intent(out) :: thisStat
334     integer :: alloc_stat
335    
336     thisStat = 0
337    
338     call FreeSimGlobals()
339    
340     allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
341     if (alloc_stat /= 0 ) then
342     thisStat = -1
343     return
344     endif
345    
346     allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
347     if (alloc_stat /= 0 ) then
348     thisStat = -1
349     return
350     endif
351    
352     end subroutine InitializeSimGlobals
353    
354     subroutine FreeSimGlobals()
355    
356     !We free in the opposite order in which we allocate in.
357    
358     if (allocated(excludesGlobal)) deallocate(excludesGlobal)
359     if (allocated(excludesLocal)) deallocate(excludesLocal)
360    
361     end subroutine FreeSimGlobals
362    
363     pure function getNlocal() result(nlocal)
364     integer :: nlocal
365     nlocal = natoms
366     end function getNlocal
367    
368    
369     end module simulation