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

File Contents

# Content
1 ! Fortran interface to C entry plug.
2
3 module force_globals
4 use definitions
5 #ifdef IS_MPI
6 use mpiSimulation
7 #endif
8
9 implicit none
10 PRIVATE
11
12 logical, save :: force_globals_initialized = .false.
13
14 #ifdef IS_MPI
15 real( kind = dp ), allocatable, dimension(:,:), public :: q_Row
16 real( kind = dp ), allocatable, dimension(:,:), public :: q_Col
17 real( kind = dp ), allocatable, dimension(:,:), public :: qcom_Row
18 real( kind = dp ), allocatable, dimension(:,:), public :: qcom_Col
19 real( kind = dp ), allocatable, dimension(:,:), public :: u_l_Row
20 real( kind = dp ), allocatable, dimension(:,:), public :: u_l_Col
21 real( kind = dp ), allocatable, dimension(:,:), public :: A_Row
22 real( kind = dp ), allocatable, dimension(:,:), public :: A_Col
23
24 real( kind = dp ), allocatable, dimension(:), public :: pot_Row
25 real( kind = dp ), allocatable, dimension(:), public :: pot_Col
26 real( kind = dp ), allocatable, dimension(:), public :: pot_Temp
27 real( kind = dp ), allocatable, dimension(:,:), public :: f_Row
28 real( kind = dp ), allocatable, dimension(:,:), public :: f_Col
29 real( kind = dp ), allocatable, dimension(:,:), public :: f_Temp
30 real( kind = dp ), allocatable, dimension(:,:), public :: t_Row
31 real( kind = dp ), allocatable, dimension(:,:), public :: t_Col
32 real( kind = dp ), allocatable, dimension(:,:), public :: t_Temp
33 real( kind = dp ), allocatable, dimension(:,:), public :: rf_Row
34 real( kind = dp ), allocatable, dimension(:,:), public :: rf_Col
35 real( kind = dp ), allocatable, dimension(:,:), public :: rf_Temp
36
37 integer, allocatable, dimension(:), public :: atid_Row
38 integer, allocatable, dimension(:), public :: atid_Col
39 #endif
40
41 integer, allocatable, dimension(:), public :: atid
42
43 real( kind = dp ), allocatable, dimension(:,:), public :: rf
44 real(kind = dp), dimension(9), public :: tau_Temp = 0.0_dp
45 real(kind = dp), public :: virial_Temp = 0.0_dp
46
47 public :: InitializeForceGlobals
48
49 contains
50
51 subroutine InitializeForceGlobals(nlocal, thisStat)
52 integer, intent(out) :: thisStat
53 integer :: nrow, ncol
54 integer :: nlocal
55 integer :: ndim = 3
56 integer :: alloc_stat
57
58 thisStat = 0
59
60 #ifdef IS_MPI
61 nrow = getNrow(plan_row)
62 ncol = getNcol(plan_col)
63 #endif
64
65 call FreeForceGlobals()
66
67 #ifdef IS_MPI
68
69 allocate(q_Row(ndim,nrow),stat=alloc_stat)
70 if (alloc_stat /= 0 ) then
71 thisStat = -1
72 return
73 endif
74
75 allocate(q_Col(ndim,ncol),stat=alloc_stat)
76 if (alloc_stat /= 0 ) then
77 thisStat = -1
78 return
79 endif
80
81 allocate(qcom_Row(ndim,nrow),stat=alloc_stat)
82 if (alloc_stat /= 0 ) then
83 thisStat = -1
84 return
85 endif
86
87 allocate(qcom_Col(ndim,ncol),stat=alloc_stat)
88 if (alloc_stat /= 0 ) then
89 thisStat = -1
90 return
91 endif
92
93 allocate(u_l_Row(ndim,nrow),stat=alloc_stat)
94 if (alloc_stat /= 0 ) then
95 thisStat = -1
96 return
97 endif
98
99 allocate(u_l_Col(ndim,ncol),stat=alloc_stat)
100 if (alloc_stat /= 0 ) then
101 thisStat = -1
102 return
103 endif
104
105 allocate(A_row(9,nrow),stat=alloc_stat)
106 if (alloc_stat /= 0 ) then
107 thisStat = -1
108 return
109 endif
110
111 allocate(A_Col(9,ncol),stat=alloc_stat)
112 if (alloc_stat /= 0 ) then
113 thisStat = -1
114 return
115 endif
116
117 allocate(pot_row(nrow),stat=alloc_stat)
118 if (alloc_stat /= 0 ) then
119 thisStat = -1
120 return
121 endif
122
123 allocate(pot_Col(ncol),stat=alloc_stat)
124 if (alloc_stat /= 0 ) then
125 thisStat = -1
126 return
127 endif
128
129 allocate(pot_Temp(nlocal),stat=alloc_stat)
130 if (alloc_stat /= 0 ) then
131 thisStat = -1
132 return
133 endif
134
135 allocate(f_Row(ndim,nrow),stat=alloc_stat)
136 if (alloc_stat /= 0 ) then
137 thisStat = -1
138 return
139 endif
140
141 allocate(f_Col(ndim,ncol),stat=alloc_stat)
142 if (alloc_stat /= 0 ) then
143 thisStat = -1
144 return
145 endif
146
147 allocate(f_Temp(ndim,nlocal),stat=alloc_stat)
148 if (alloc_stat /= 0 ) then
149 thisStat = -1
150 return
151 endif
152
153 allocate(t_Row(ndim,nrow),stat=alloc_stat)
154 if (alloc_stat /= 0 ) then
155 thisStat = -1
156 return
157 endif
158
159 allocate(t_Col(ndim,ncol),stat=alloc_stat)
160 if (alloc_stat /= 0 ) then
161 thisStat = -1
162 return
163 endif
164
165 allocate(t_temp(ndim,nlocal),stat=alloc_stat)
166 if (alloc_stat /= 0 ) then
167 thisStat = -1
168 return
169 endif
170
171 allocate(atid(nlocal),stat=alloc_stat)
172 if (alloc_stat /= 0 ) then
173 thisStat = -1
174 return
175 endif
176
177
178 allocate(atid_Row(nrow),stat=alloc_stat)
179 if (alloc_stat /= 0 ) then
180 thisStat = -1
181 return
182 endif
183
184 allocate(atid_Col(ncol),stat=alloc_stat)
185 if (alloc_stat /= 0 ) then
186 thisStat = -1
187 return
188 endif
189
190 allocate(rf_Row(ndim,nrow),stat=alloc_stat)
191 if (alloc_stat /= 0 ) then
192 thisStat = -1
193 return
194 endif
195
196 allocate(rf_Col(ndim,ncol),stat=alloc_stat)
197 if (alloc_stat /= 0 ) then
198 thisStat = -1
199 return
200 endif
201
202 allocate(rf_Temp(ndim,nlocal),stat=alloc_stat)
203 if (alloc_stat /= 0 ) then
204 thisStat = -1
205 return
206 endif
207
208 #else
209
210 allocate(atid(nlocal),stat=alloc_stat)
211 if (alloc_stat /= 0 ) then
212 thisStat = -1
213 return
214 end if
215
216 #endif
217
218 allocate(rf(ndim,nlocal),stat=alloc_stat)
219 if (alloc_stat /= 0 ) then
220 thisStat = -1
221 return
222 endif
223
224 force_globals_initialized = .true.
225
226 end subroutine InitializeForceGlobals
227
228 subroutine FreeForceGlobals()
229
230 !We free in the opposite order in which we allocate in.
231
232 if (allocated(rf)) deallocate(rf)
233 #ifdef IS_MPI
234 if (allocated(rf_Temp)) deallocate(rf_Temp)
235 if (allocated(rf_Col)) deallocate(rf_Col)
236 if (allocated(rf_Row)) deallocate(rf_Row)
237 if (allocated(atid_Col)) deallocate(atid_Col)
238 if (allocated(atid_Row)) deallocate(atid_Row)
239 if (allocated(atid)) deallocate(atid)
240 if (allocated(t_Temp)) deallocate(t_Temp)
241 if (allocated(t_Col)) deallocate(t_Col)
242 if (allocated(t_Row)) deallocate(t_Row)
243 if (allocated(f_Temp)) deallocate(f_Temp)
244 if (allocated(f_Col)) deallocate(f_Col)
245 if (allocated(f_Row)) deallocate(f_Row)
246 if (allocated(pot_Temp)) deallocate(pot_Temp)
247 if (allocated(pot_Col)) deallocate(pot_Col)
248 if (allocated(pot_Row)) deallocate(pot_Row)
249 if (allocated(A_Col)) deallocate(A_Col)
250 if (allocated(A_Row)) deallocate(A_Row)
251 if (allocated(u_l_Col)) deallocate(u_l_Col)
252 if (allocated(u_l_Row)) deallocate(u_l_Row)
253 if (allocated(qcom_Col)) deallocate(qcom_Col)
254 if (allocated(qcom_Row)) deallocate(qcom_Row)
255 if (allocated(q_Col)) deallocate(q_Col)
256 if (allocated(q_Row)) deallocate(q_Row)
257 #else
258 if (allocated(atid)) deallocate(atid)
259 #endif
260
261 end subroutine FreeForceGlobals
262
263 end module force_globals