ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 1138
Committed: Wed Apr 28 21:39:22 2004 UTC (20 years, 2 months ago) by gezelter
File size: 7795 byte(s)
Log Message:
work on molecular cutoffs

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 gezelter 747 type (simtype), public, save :: thisSim
20 mmeineke 377
21     logical, save :: simulation_setup_complete = .false.
22    
23 mmeineke 491 integer, public, save :: nLocal, nGlobal
24 mmeineke 377 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 chuckv 482 integer, allocatable, dimension(:), public :: molMembershipList
29 mmeineke 377
30 mmeineke 569 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
31 gezelter 570 logical, public, save :: boxIsOrthorhombic
32 mmeineke 377
33     public :: SimulationSetup
34     public :: getNlocal
35     public :: setBox
36     public :: getDielect
37     public :: SimUsesPBC
38     public :: SimUsesLJ
39 gezelter 941 public :: SimUsesCharges
40 mmeineke 377 public :: SimUsesDipoles
41     public :: SimUsesSticky
42     public :: SimUsesRF
43     public :: SimUsesGB
44     public :: SimUsesEAM
45     public :: SimRequiresPrepairCalc
46     public :: SimRequiresPostpairCalc
47     public :: SimUsesDirectionalAtoms
48 gezelter 1138 public :: SimUsesMolecularCutoffs
49 mmeineke 377
50     contains
51    
52 mmeineke 491 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
53 gezelter 483 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
54     CmolMembership, &
55 mmeineke 377 status)
56    
57     type (simtype) :: setThisSim
58 mmeineke 491 integer, intent(inout) :: CnGlobal, CnLocal
59     integer, dimension(CnLocal),intent(inout) :: c_idents
60 mmeineke 377
61     integer :: CnLocalExcludes
62     integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
63     integer :: CnGlobalExcludes
64     integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
65 mmeineke 491 integer, dimension(CnGlobal),intent(in) :: CmolMembership
66 mmeineke 377 !! Result status, success = 0, status = -1
67     integer, intent(out) :: status
68     integer :: i, me, thisStat, alloc_stat, myNode
69     #ifdef IS_MPI
70     integer, allocatable, dimension(:) :: c_idents_Row
71     integer, allocatable, dimension(:) :: c_idents_Col
72     integer :: nrow
73     integer :: ncol
74     #endif
75    
76     simulation_setup_complete = .false.
77     status = 0
78    
79     ! copy C struct into fortran type
80 mmeineke 491
81     nLocal = CnLocal
82     nGlobal = CnGlobal
83    
84 mmeineke 377 thisSim = setThisSim
85 mmeineke 620
86 mmeineke 377 nExcludes_Global = CnGlobalExcludes
87     nExcludes_Local = CnLocalExcludes
88    
89 mmeineke 491 call InitializeForceGlobals(nLocal, thisStat)
90 mmeineke 377 if (thisStat /= 0) then
91 chuckv 480 write(default_error,*) "SimSetup: InitializeForceGlobals error"
92 mmeineke 377 status = -1
93     return
94     endif
95    
96     call InitializeSimGlobals(thisStat)
97     if (thisStat /= 0) then
98 chuckv 480 write(default_error,*) "SimSetup: InitializeSimGlobals error"
99 mmeineke 377 status = -1
100     return
101     endif
102    
103     #ifdef IS_MPI
104     ! We can only set up forces if mpiSimulation has been setup.
105     if (.not. isMPISimSet()) then
106     write(default_error,*) "MPI is not set"
107     status = -1
108     return
109     endif
110     nrow = getNrow(plan_row)
111     ncol = getNcol(plan_col)
112     mynode = getMyNode()
113    
114     allocate(c_idents_Row(nrow),stat=alloc_stat)
115     if (alloc_stat /= 0 ) then
116     status = -1
117     return
118     endif
119    
120     allocate(c_idents_Col(ncol),stat=alloc_stat)
121     if (alloc_stat /= 0 ) then
122     status = -1
123     return
124     endif
125    
126     call gather(c_idents, c_idents_Row, plan_row)
127     call gather(c_idents, c_idents_Col, plan_col)
128    
129     do i = 1, nrow
130     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
131     atid_Row(i) = me
132     enddo
133    
134     do i = 1, ncol
135     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
136     atid_Col(i) = me
137     enddo
138    
139     !! free temporary ident arrays
140     if (allocated(c_idents_Col)) then
141     deallocate(c_idents_Col)
142     end if
143     if (allocated(c_idents_Row)) then
144     deallocate(c_idents_Row)
145     endif
146    
147 chuckv 648 #endif
148    
149     ! We build the local atid's for both mpi and nonmpi
150 gezelter 490 do i = 1, nLocal
151 mmeineke 377
152     me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
153     atid(i) = me
154    
155     enddo
156    
157    
158 chuckv 388
159 chuckv 648
160 mmeineke 377 do i = 1, nExcludes_Local
161     excludesLocal(1,i) = CexcludesLocal(1,i)
162     excludesLocal(2,i) = CexcludesLocal(2,i)
163     enddo
164    
165     do i = 1, nExcludes_Global
166     excludesGlobal(i) = CexcludesGlobal(i)
167     enddo
168 mmeineke 435
169 gezelter 490 do i = 1, nGlobal
170 gezelter 483 molMemberShipList(i) = CmolMembership(i)
171 mmeineke 491 enddo
172 chuckv 482
173 mmeineke 377 if (status == 0) simulation_setup_complete = .true.
174    
175     end subroutine SimulationSetup
176    
177 mmeineke 569 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
178 gezelter 570 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
179 mmeineke 569 integer :: cBoxIsOrthorhombic
180 mmeineke 377 integer :: smallest, status, i
181 gezelter 570
182 mmeineke 569 Hmat = cHmat
183     HmatInv = cHmatInv
184 gezelter 570 if (cBoxIsOrthorhombic .eq. 0 ) then
185 mmeineke 569 boxIsOrthorhombic = .false.
186 gezelter 570 else
187     boxIsOrthorhombic = .true.
188 mmeineke 569 endif
189    
190 mmeineke 377 return
191 mmeineke 569 end subroutine setBox
192 mmeineke 377
193     function getDielect() result(dielect)
194     real( kind = dp ) :: dielect
195     dielect = thisSim%dielect
196     end function getDielect
197    
198     function SimUsesPBC() result(doesit)
199     logical :: doesit
200     doesit = thisSim%SIM_uses_PBC
201     end function SimUsesPBC
202    
203     function SimUsesLJ() result(doesit)
204     logical :: doesit
205     doesit = thisSim%SIM_uses_LJ
206     end function SimUsesLJ
207    
208     function SimUsesSticky() result(doesit)
209     logical :: doesit
210     doesit = thisSim%SIM_uses_sticky
211     end function SimUsesSticky
212    
213 gezelter 941 function SimUsesCharges() result(doesit)
214     logical :: doesit
215     doesit = thisSim%SIM_uses_charges
216     end function SimUsesCharges
217    
218 mmeineke 377 function SimUsesDipoles() result(doesit)
219     logical :: doesit
220     doesit = thisSim%SIM_uses_dipoles
221     end function SimUsesDipoles
222    
223     function SimUsesRF() result(doesit)
224     logical :: doesit
225     doesit = thisSim%SIM_uses_RF
226     end function SimUsesRF
227    
228     function SimUsesGB() result(doesit)
229     logical :: doesit
230     doesit = thisSim%SIM_uses_GB
231     end function SimUsesGB
232    
233     function SimUsesEAM() result(doesit)
234     logical :: doesit
235     doesit = thisSim%SIM_uses_EAM
236     end function SimUsesEAM
237    
238     function SimUsesDirectionalAtoms() result(doesit)
239     logical :: doesit
240     doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
241     thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
242     end function SimUsesDirectionalAtoms
243    
244 gezelter 1138 function SimUsesMolecularCutoffs() result(doesit)
245     logical :: doesit
246     doesit = thisSim%SIM_uses_molecular_cutoffs
247     end function SimUsesMolecularCutoffs
248    
249 mmeineke 377 function SimRequiresPrepairCalc() result(doesit)
250     logical :: doesit
251     doesit = thisSim%SIM_uses_EAM
252     end function SimRequiresPrepairCalc
253    
254     function SimRequiresPostpairCalc() result(doesit)
255     logical :: doesit
256     doesit = thisSim%SIM_uses_RF
257     end function SimRequiresPostpairCalc
258    
259     subroutine InitializeSimGlobals(thisStat)
260     integer, intent(out) :: thisStat
261     integer :: alloc_stat
262    
263     thisStat = 0
264    
265     call FreeSimGlobals()
266    
267     allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
268     if (alloc_stat /= 0 ) then
269     thisStat = -1
270     return
271     endif
272    
273     allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
274     if (alloc_stat /= 0 ) then
275     thisStat = -1
276     return
277     endif
278 chuckv 482
279 mmeineke 491 allocate(molMembershipList(nGlobal), stat=alloc_stat)
280 chuckv 482 if (alloc_stat /= 0 ) then
281     thisStat = -1
282     return
283     endif
284 mmeineke 377
285     end subroutine InitializeSimGlobals
286    
287     subroutine FreeSimGlobals()
288    
289     !We free in the opposite order in which we allocate in.
290 gezelter 483
291     if (allocated(molMembershipList)) deallocate(molMembershipList)
292 mmeineke 377 if (allocated(excludesGlobal)) deallocate(excludesGlobal)
293     if (allocated(excludesLocal)) deallocate(excludesLocal)
294 gezelter 483
295 mmeineke 377 end subroutine FreeSimGlobals
296    
297 mmeineke 491 pure function getNlocal() result(n)
298     integer :: n
299     n = nLocal
300 mmeineke 377 end function getNlocal
301    
302    
303     end module simulation