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: 2361
Committed: Wed Oct 12 21:00:50 2005 UTC (18 years, 8 months ago) by gezelter
File size: 9481 byte(s)
Log Message:
simplifying long range potential array

File Contents

# Content
1 !!
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 ! 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 #define __FORTRAN90
54 #include "UseTheForce/DarkSide/fInteractionMap.h"
55
56 logical, save :: force_globals_initialized = .false.
57
58 #ifdef IS_MPI
59 real( kind = dp ), allocatable, dimension(:,:), public :: q_Row
60 real( kind = dp ), allocatable, dimension(:,:), public :: q_Col
61 real( kind = dp ), allocatable, dimension(:,:), public :: q_group_Row
62 real( kind = dp ), allocatable, dimension(:,:), public :: q_group_Col
63 real( kind = dp ), allocatable, dimension(:,:), public :: eFrame_Row
64 real( kind = dp ), allocatable, dimension(:,:), public :: eFrame_Col
65 real( kind = dp ), allocatable, dimension(:,:), public :: A_Row
66 real( kind = dp ), allocatable, dimension(:,:), public :: A_Col
67
68 real( kind = dp ), allocatable, dimension(:,:), public :: pot_Row
69 real( kind = dp ), allocatable, dimension(:,:), public :: pot_Col
70 real( kind = dp ), allocatable, dimension(:,:), public :: pot_Temp
71 real( kind = dp ), allocatable, dimension(:,:), public :: f_Row
72 real( kind = dp ), allocatable, dimension(:,:), public :: f_Col
73 real( kind = dp ), allocatable, dimension(:,:), public :: f_Temp
74 real( kind = dp ), allocatable, dimension(:,:), public :: t_Row
75 real( kind = dp ), allocatable, dimension(:,:), public :: t_Col
76 real( kind = dp ), allocatable, dimension(:,:), public :: t_Temp
77 real( kind = dp ), allocatable, dimension(:,:), public :: rf_Row
78 real( kind = dp ), allocatable, dimension(:,:), public :: rf_Col
79 real( kind = dp ), allocatable, dimension(:,:), public :: rf_Temp
80
81 integer, allocatable, dimension(:), public :: atid_Row
82 integer, allocatable, dimension(:), public :: atid_Col
83 #endif
84
85 integer, allocatable, dimension(:), public :: atid
86
87 real( kind = dp ), allocatable, dimension(:,:), public :: rf
88 real(kind = dp), dimension(9), public :: tau_Temp = 0.0_dp
89 real(kind = dp), public :: virial_Temp = 0.0_dp
90
91 public :: InitializeForceGlobals
92
93 contains
94
95 subroutine InitializeForceGlobals(nlocal, thisStat)
96 integer, intent(out) :: thisStat
97 integer :: nAtomsInRow, nAtomsInCol
98 integer :: nGroupsInRow, nGroupsInCol
99 integer :: nlocal
100 integer :: ndim = 3
101 integer :: alloc_stat
102
103 thisStat = 0
104
105 #ifdef IS_MPI
106 nAtomsInRow = getNatomsInRow(plan_atom_row)
107 nAtomsInCol = getNatomsInCol(plan_atom_col)
108 nGroupsInRow = getNgroupsInRow(plan_group_row)
109 nGroupsInCol = getNgroupsInCol(plan_group_col)
110
111 #endif
112
113 call FreeForceGlobals()
114
115 #ifdef IS_MPI
116
117 allocate(q_Row(ndim,nAtomsInRow),stat=alloc_stat)
118 if (alloc_stat /= 0 ) then
119 thisStat = -1
120 return
121 endif
122
123 allocate(q_Col(ndim,nAtomsInCol),stat=alloc_stat)
124 if (alloc_stat /= 0 ) then
125 thisStat = -1
126 return
127 endif
128
129 allocate(q_group_Row(ndim,nGroupsInRow),stat=alloc_stat)
130 if (alloc_stat /= 0 ) then
131 thisStat = -1
132 return
133 endif
134
135 allocate(q_group_Col(ndim,nGroupsInCol),stat=alloc_stat)
136 if (alloc_stat /= 0 ) then
137 thisStat = -1
138 return
139 endif
140
141 allocate(eFrame_Row(9,nAtomsInRow),stat=alloc_stat)
142 if (alloc_stat /= 0 ) then
143 thisStat = -1
144 return
145 endif
146
147 allocate(eFrame_Col(9,nAtomsInCol),stat=alloc_stat)
148 if (alloc_stat /= 0 ) then
149 thisStat = -1
150 return
151 endif
152
153 allocate(A_row(9,nAtomsInRow),stat=alloc_stat)
154 if (alloc_stat /= 0 ) then
155 thisStat = -1
156 return
157 endif
158
159 allocate(A_Col(9,nAtomsInCol),stat=alloc_stat)
160 if (alloc_stat /= 0 ) then
161 thisStat = -1
162 return
163 endif
164
165 allocate(pot_row(LR_POT_TYPES,nAtomsInRow),stat=alloc_stat)
166 if (alloc_stat /= 0 ) then
167 thisStat = -1
168 return
169 endif
170
171 allocate(pot_Col(LR_POT_TYPES,nAtomsInCol),stat=alloc_stat)
172 if (alloc_stat /= 0 ) then
173 thisStat = -1
174 return
175 endif
176
177 allocate(pot_Temp(LR_POT_TYPES,nlocal),stat=alloc_stat)
178 if (alloc_stat /= 0 ) then
179 thisStat = -1
180 return
181 endif
182
183 allocate(f_Row(ndim,nAtomsInRow),stat=alloc_stat)
184 if (alloc_stat /= 0 ) then
185 thisStat = -1
186 return
187 endif
188
189 allocate(f_Col(ndim,nAtomsInCol),stat=alloc_stat)
190 if (alloc_stat /= 0 ) then
191 thisStat = -1
192 return
193 endif
194
195 allocate(f_Temp(ndim,nlocal),stat=alloc_stat)
196 if (alloc_stat /= 0 ) then
197 thisStat = -1
198 return
199 endif
200
201 allocate(t_Row(ndim,nAtomsInRow),stat=alloc_stat)
202 if (alloc_stat /= 0 ) then
203 thisStat = -1
204 return
205 endif
206
207 allocate(t_Col(ndim,nAtomsInCol),stat=alloc_stat)
208 if (alloc_stat /= 0 ) then
209 thisStat = -1
210 return
211 endif
212
213 allocate(t_temp(ndim,nlocal),stat=alloc_stat)
214 if (alloc_stat /= 0 ) then
215 thisStat = -1
216 return
217 endif
218
219 allocate(atid(nlocal),stat=alloc_stat)
220 if (alloc_stat /= 0 ) then
221 thisStat = -1
222 return
223 endif
224
225
226 allocate(atid_Row(nAtomsInRow),stat=alloc_stat)
227 if (alloc_stat /= 0 ) then
228 thisStat = -1
229 return
230 endif
231
232 allocate(atid_Col(nAtomsInCol),stat=alloc_stat)
233 if (alloc_stat /= 0 ) then
234 thisStat = -1
235 return
236 endif
237
238 allocate(rf_Row(ndim,nAtomsInRow),stat=alloc_stat)
239 if (alloc_stat /= 0 ) then
240 thisStat = -1
241 return
242 endif
243
244 allocate(rf_Col(ndim,nAtomsInCol),stat=alloc_stat)
245 if (alloc_stat /= 0 ) then
246 thisStat = -1
247 return
248 endif
249
250 allocate(rf_Temp(ndim,nlocal),stat=alloc_stat)
251 if (alloc_stat /= 0 ) then
252 thisStat = -1
253 return
254 endif
255
256 #else
257
258 allocate(atid(nlocal),stat=alloc_stat)
259 if (alloc_stat /= 0 ) then
260 thisStat = -1
261 return
262 end if
263
264 #endif
265
266 allocate(rf(ndim,nlocal),stat=alloc_stat)
267 if (alloc_stat /= 0 ) then
268 thisStat = -1
269 return
270 endif
271
272 force_globals_initialized = .true.
273
274 end subroutine InitializeForceGlobals
275
276 subroutine FreeForceGlobals()
277
278 !We free in the opposite order in which we allocate in.
279
280 if (allocated(rf)) deallocate(rf)
281 #ifdef IS_MPI
282 if (allocated(rf_Temp)) deallocate(rf_Temp)
283 if (allocated(rf_Col)) deallocate(rf_Col)
284 if (allocated(rf_Row)) deallocate(rf_Row)
285 if (allocated(atid_Col)) deallocate(atid_Col)
286 if (allocated(atid_Row)) deallocate(atid_Row)
287 if (allocated(atid)) deallocate(atid)
288 if (allocated(t_Temp)) deallocate(t_Temp)
289 if (allocated(t_Col)) deallocate(t_Col)
290 if (allocated(t_Row)) deallocate(t_Row)
291 if (allocated(f_Temp)) deallocate(f_Temp)
292 if (allocated(f_Col)) deallocate(f_Col)
293 if (allocated(f_Row)) deallocate(f_Row)
294 if (allocated(pot_Temp)) deallocate(pot_Temp)
295 if (allocated(pot_Col)) deallocate(pot_Col)
296 if (allocated(pot_Row)) deallocate(pot_Row)
297 if (allocated(A_Col)) deallocate(A_Col)
298 if (allocated(A_Row)) deallocate(A_Row)
299 if (allocated(eFrame_Col)) deallocate(eFrame_Col)
300 if (allocated(eFrame_Row)) deallocate(eFrame_Row)
301 if (allocated(q_group_Col)) deallocate(q_group_Col)
302 if (allocated(q_group_Row)) deallocate(q_group_Row)
303 if (allocated(q_Col)) deallocate(q_Col)
304 if (allocated(q_Row)) deallocate(q_Row)
305 #else
306 if (allocated(atid)) deallocate(atid)
307 #endif
308
309 end subroutine FreeForceGlobals
310
311 end module force_globals