ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/force_globals.F90
Revision: 648
Committed: Wed Jul 23 22:13:59 2003 UTC (20 years, 11 months ago) by chuckv
File size: 6540 byte(s)
Log Message:
Finished most code for eam....

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