ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/force_globals.F90
Revision: 2204
Committed: Fri Apr 15 22:04:00 2005 UTC (19 years, 2 months ago) by gezelter
File size: 9366 byte(s)
Log Message:
xemacs has been drafted to perform our indentation services

File Contents

# User Rev Content
1 gezelter 1930 !!
2     !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     !!
4     !! The University of Notre Dame grants you ("Licensee") a
5     !! non-exclusive, royalty free, license to use, modify and
6     !! redistribute this software in source and binary code form, provided
7     !! that the following conditions are met:
8     !!
9     !! 1. Acknowledgement of the program authors must be made in any
10     !! publication of scientific results based in part on use of the
11     !! program. An acceptable form of acknowledgement is citation of
12     !! the article in which the program was described (Matthew
13     !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     !! Parallel Simulation Engine for Molecular Dynamics,"
16     !! J. Comput. Chem. 26, pp. 252-271 (2005))
17     !!
18     !! 2. Redistributions of source code must retain the above copyright
19     !! notice, this list of conditions and the following disclaimer.
20     !!
21     !! 3. Redistributions in binary form must reproduce the above copyright
22     !! notice, this list of conditions and the following disclaimer in the
23     !! documentation and/or other materials provided with the
24     !! distribution.
25     !!
26     !! This software is provided "AS IS," without a warranty of any
27     !! kind. All express or implied conditions, representations and
28     !! warranties, including any implied warranty of merchantability,
29     !! fitness for a particular purpose or non-infringement, are hereby
30     !! excluded. The University of Notre Dame and its licensors shall not
31     !! be liable for any damages suffered by licensee as a result of
32     !! using, modifying or distributing the software or its
33     !! derivatives. In no event will the University of Notre Dame or its
34     !! licensors be liable for any lost revenue, profit or data, or for
35     !! direct, indirect, special, consequential, incidental or punitive
36     !! damages, however caused and regardless of the theory of liability,
37     !! arising out of the use of or inability to use software, even if the
38     !! University of Notre Dame has been advised of the possibility of
39     !! such damages.
40     !!
41    
42    
43 gezelter 1490 ! Fortran interface to C entry plug.
44    
45     module force_globals
46     use definitions
47     #ifdef IS_MPI
48     use mpiSimulation
49     #endif
50    
51     implicit none
52     PRIVATE
53    
54     logical, save :: force_globals_initialized = .false.
55    
56     #ifdef IS_MPI
57     real( kind = dp ), allocatable, dimension(:,:), public :: q_Row
58     real( kind = dp ), allocatable, dimension(:,:), public :: q_Col
59     real( kind = dp ), allocatable, dimension(:,:), public :: q_group_Row
60     real( kind = dp ), allocatable, dimension(:,:), public :: q_group_Col
61 gezelter 1930 real( kind = dp ), allocatable, dimension(:,:), public :: eFrame_Row
62     real( kind = dp ), allocatable, dimension(:,:), public :: eFrame_Col
63 gezelter 1490 real( kind = dp ), allocatable, dimension(:,:), public :: A_Row
64     real( kind = dp ), allocatable, dimension(:,:), public :: A_Col
65 gezelter 2204
66 gezelter 1490 real( kind = dp ), allocatable, dimension(:), public :: pot_Row
67     real( kind = dp ), allocatable, dimension(:), public :: pot_Col
68     real( kind = dp ), allocatable, dimension(:), public :: pot_Temp
69     real( kind = dp ), allocatable, dimension(:,:), public :: f_Row
70     real( kind = dp ), allocatable, dimension(:,:), public :: f_Col
71     real( kind = dp ), allocatable, dimension(:,:), public :: f_Temp
72     real( kind = dp ), allocatable, dimension(:,:), public :: t_Row
73     real( kind = dp ), allocatable, dimension(:,:), public :: t_Col
74     real( kind = dp ), allocatable, dimension(:,:), public :: t_Temp
75     real( kind = dp ), allocatable, dimension(:,:), public :: rf_Row
76     real( kind = dp ), allocatable, dimension(:,:), public :: rf_Col
77     real( kind = dp ), allocatable, dimension(:,:), public :: rf_Temp
78    
79     integer, allocatable, dimension(:), public :: atid_Row
80     integer, allocatable, dimension(:), public :: atid_Col
81     #endif
82    
83     integer, allocatable, dimension(:), public :: atid
84    
85     real( kind = dp ), allocatable, dimension(:,:), public :: rf
86     real(kind = dp), dimension(9), public :: tau_Temp = 0.0_dp
87     real(kind = dp), public :: virial_Temp = 0.0_dp
88 gezelter 2204
89 gezelter 1490 public :: InitializeForceGlobals
90 gezelter 2204
91 gezelter 1490 contains
92 gezelter 2204
93 gezelter 1490 subroutine InitializeForceGlobals(nlocal, thisStat)
94     integer, intent(out) :: thisStat
95     integer :: nAtomsInRow, nAtomsInCol
96     integer :: nGroupsInRow, nGroupsInCol
97     integer :: nlocal
98     integer :: ndim = 3
99     integer :: alloc_stat
100 gezelter 2204
101 gezelter 1490 thisStat = 0
102 gezelter 2204
103 gezelter 1490 #ifdef IS_MPI
104     nAtomsInRow = getNatomsInRow(plan_atom_row)
105     nAtomsInCol = getNatomsInCol(plan_atom_col)
106     nGroupsInRow = getNgroupsInRow(plan_group_row)
107     nGroupsInCol = getNgroupsInCol(plan_group_col)
108 gezelter 2204
109 gezelter 1490 #endif
110 gezelter 2204
111 gezelter 1490 call FreeForceGlobals()
112 gezelter 2204
113 gezelter 1490 #ifdef IS_MPI
114    
115     allocate(q_Row(ndim,nAtomsInRow),stat=alloc_stat)
116     if (alloc_stat /= 0 ) then
117     thisStat = -1
118     return
119     endif
120 gezelter 2204
121 gezelter 1490 allocate(q_Col(ndim,nAtomsInCol),stat=alloc_stat)
122     if (alloc_stat /= 0 ) then
123     thisStat = -1
124     return
125     endif
126    
127     allocate(q_group_Row(ndim,nGroupsInRow),stat=alloc_stat)
128     if (alloc_stat /= 0 ) then
129     thisStat = -1
130     return
131     endif
132 gezelter 2204
133 gezelter 1490 allocate(q_group_Col(ndim,nGroupsInCol),stat=alloc_stat)
134     if (alloc_stat /= 0 ) then
135     thisStat = -1
136     return
137     endif
138 gezelter 2204
139 gezelter 1930 allocate(eFrame_Row(9,nAtomsInRow),stat=alloc_stat)
140 gezelter 1490 if (alloc_stat /= 0 ) then
141     thisStat = -1
142     return
143     endif
144 gezelter 2204
145 gezelter 1930 allocate(eFrame_Col(9,nAtomsInCol),stat=alloc_stat)
146 gezelter 1490 if (alloc_stat /= 0 ) then
147     thisStat = -1
148     return
149     endif
150 gezelter 2204
151 gezelter 1490 allocate(A_row(9,nAtomsInRow),stat=alloc_stat)
152     if (alloc_stat /= 0 ) then
153     thisStat = -1
154     return
155     endif
156 gezelter 2204
157 gezelter 1490 allocate(A_Col(9,nAtomsInCol),stat=alloc_stat)
158     if (alloc_stat /= 0 ) then
159     thisStat = -1
160     return
161     endif
162 gezelter 2204
163 gezelter 1490 allocate(pot_row(nAtomsInRow),stat=alloc_stat)
164     if (alloc_stat /= 0 ) then
165     thisStat = -1
166     return
167     endif
168 gezelter 2204
169 gezelter 1490 allocate(pot_Col(nAtomsInCol),stat=alloc_stat)
170     if (alloc_stat /= 0 ) then
171     thisStat = -1
172     return
173     endif
174    
175     allocate(pot_Temp(nlocal),stat=alloc_stat)
176     if (alloc_stat /= 0 ) then
177     thisStat = -1
178     return
179     endif
180 gezelter 2204
181 gezelter 1490 allocate(f_Row(ndim,nAtomsInRow),stat=alloc_stat)
182     if (alloc_stat /= 0 ) then
183     thisStat = -1
184     return
185     endif
186 gezelter 2204
187 gezelter 1490 allocate(f_Col(ndim,nAtomsInCol),stat=alloc_stat)
188     if (alloc_stat /= 0 ) then
189     thisStat = -1
190     return
191     endif
192 gezelter 2204
193 gezelter 1490 allocate(f_Temp(ndim,nlocal),stat=alloc_stat)
194     if (alloc_stat /= 0 ) then
195     thisStat = -1
196     return
197     endif
198 gezelter 2204
199 gezelter 1490 allocate(t_Row(ndim,nAtomsInRow),stat=alloc_stat)
200     if (alloc_stat /= 0 ) then
201     thisStat = -1
202     return
203     endif
204 gezelter 2204
205 gezelter 1490 allocate(t_Col(ndim,nAtomsInCol),stat=alloc_stat)
206     if (alloc_stat /= 0 ) then
207     thisStat = -1
208     return
209     endif
210    
211     allocate(t_temp(ndim,nlocal),stat=alloc_stat)
212     if (alloc_stat /= 0 ) then
213     thisStat = -1
214     return
215     endif
216    
217     allocate(atid(nlocal),stat=alloc_stat)
218     if (alloc_stat /= 0 ) then
219     thisStat = -1
220     return
221     endif
222    
223    
224     allocate(atid_Row(nAtomsInRow),stat=alloc_stat)
225     if (alloc_stat /= 0 ) then
226     thisStat = -1
227     return
228     endif
229    
230     allocate(atid_Col(nAtomsInCol),stat=alloc_stat)
231     if (alloc_stat /= 0 ) then
232     thisStat = -1
233     return
234     endif
235    
236     allocate(rf_Row(ndim,nAtomsInRow),stat=alloc_stat)
237     if (alloc_stat /= 0 ) then
238     thisStat = -1
239     return
240     endif
241    
242     allocate(rf_Col(ndim,nAtomsInCol),stat=alloc_stat)
243     if (alloc_stat /= 0 ) then
244     thisStat = -1
245     return
246     endif
247    
248     allocate(rf_Temp(ndim,nlocal),stat=alloc_stat)
249     if (alloc_stat /= 0 ) then
250     thisStat = -1
251     return
252     endif
253    
254     #else
255    
256     allocate(atid(nlocal),stat=alloc_stat)
257     if (alloc_stat /= 0 ) then
258     thisStat = -1
259     return
260     end if
261    
262     #endif
263 gezelter 2204
264 gezelter 1490 allocate(rf(ndim,nlocal),stat=alloc_stat)
265     if (alloc_stat /= 0 ) then
266     thisStat = -1
267     return
268     endif
269    
270     force_globals_initialized = .true.
271    
272     end subroutine InitializeForceGlobals
273 gezelter 2204
274 gezelter 1490 subroutine FreeForceGlobals()
275 gezelter 2204
276 gezelter 1490 !We free in the opposite order in which we allocate in.
277 gezelter 2204
278 gezelter 1490 if (allocated(rf)) deallocate(rf)
279     #ifdef IS_MPI
280     if (allocated(rf_Temp)) deallocate(rf_Temp)
281     if (allocated(rf_Col)) deallocate(rf_Col)
282     if (allocated(rf_Row)) deallocate(rf_Row)
283     if (allocated(atid_Col)) deallocate(atid_Col)
284     if (allocated(atid_Row)) deallocate(atid_Row)
285     if (allocated(atid)) deallocate(atid)
286     if (allocated(t_Temp)) deallocate(t_Temp)
287     if (allocated(t_Col)) deallocate(t_Col)
288     if (allocated(t_Row)) deallocate(t_Row)
289     if (allocated(f_Temp)) deallocate(f_Temp)
290     if (allocated(f_Col)) deallocate(f_Col)
291     if (allocated(f_Row)) deallocate(f_Row)
292     if (allocated(pot_Temp)) deallocate(pot_Temp)
293     if (allocated(pot_Col)) deallocate(pot_Col)
294     if (allocated(pot_Row)) deallocate(pot_Row)
295     if (allocated(A_Col)) deallocate(A_Col)
296     if (allocated(A_Row)) deallocate(A_Row)
297 gezelter 1930 if (allocated(eFrame_Col)) deallocate(eFrame_Col)
298     if (allocated(eFrame_Row)) deallocate(eFrame_Row)
299 gezelter 1490 if (allocated(q_group_Col)) deallocate(q_group_Col)
300     if (allocated(q_group_Row)) deallocate(q_group_Row)
301 gezelter 1930 if (allocated(q_Col)) deallocate(q_Col)
302     if (allocated(q_Row)) deallocate(q_Row)
303 gezelter 1490 #else
304     if (allocated(atid)) deallocate(atid)
305     #endif
306 gezelter 2204
307 gezelter 1490 end subroutine FreeForceGlobals
308 gezelter 2204
309 gezelter 1490 end module force_globals