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