MFIX  2016-1
check_data_20.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
2 ! C
3 ! Subroutine: CHECK_DATA_20 C
4 ! Purpose: C
5 ! - check whether field variables are initialized in all cells C
6 ! - check whether the sum of void and volume fractions is 1.0 C
7 ! in all fluid and mass inflow cells C
8 ! - check whether mu_gmax is specified if k_epsilon or l_scale C
9 ! C
10 ! Author: M. Syamlal Date: 30-JAN-92 C
11 ! Reviewer: P. Nicoletti, W. Rogers, S. Venkatesan Date: 31-JAN-92 C
12 ! C
13 ! Revision Number: C
14 ! Purpose: C
15 ! Author: Date: dd-mmm-yy C
16 ! Reviewer: Date: dd-mmm-yy C
17 ! C
18 ! Literature/Document References: C
19 ! C
20 ! Variables referenced: C
21 ! Variables modified: C
22 ! Local variables: C
23 ! C
24 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
25 
26  SUBROUTINE check_data_20
27 
28 !-----------------------------------------------
29 ! Modules
30 !-----------------------------------------------
31  USE param
32  USE param1
33  USE toleranc
34  USE fldvar
35  USE run
36  USE geometry
37  USE constant
38  USE physprop
39  USE indices
40  USE funits
41  USE visc_g
42  USE rxns
43  USE scalars
44  USE compar
45  USE sendrecv
46  USE discretelement
47  USE mfix_pic
48  USE functions
49 
50  use mpi_utility
51  use error_manager
52 
53  IMPLICIT NONE
54 
55 !-----------------------------------------------
56 ! Local variables
57 !-----------------------------------------------
58 ! Indices
59  INTEGER :: I, J, K, IJK, IMJK, IJMK, IJKM
60 ! Solids phase
61  INTEGER :: M
62 ! Species index
63  INTEGER :: NN
64 ! Logical variable to set, if there is an error
65  LOGICAL :: ABORT
66 ! Whether L_scale is nonzero
67  LOGICAL :: NONZERO
68 ! 1.0 - sum of all volume fractions
69  DOUBLE PRECISION DIF
70 !-----------------------------------------------
71 
72  CALL init_err_msg("CHECK_DATA_20")
73 
74  call send_recv(p_g,2)
75  call send_recv(ep_g,2)
76  call send_recv(w_s,2)
77  call send_recv(w_g,2)
78  call send_recv(u_s,2)
79  call send_recv(u_g,2)
80  call send_recv(v_s,2)
81  call send_recv(v_g,2)
82  call send_recv(ro_s,2)
83  call send_recv(rop_s,2)
84  call send_recv( p_star, 2 )
85  call send_recv( rop_g, 2 )
86  call send_recv( ro_g, 2 )
87  call send_recv( t_g, 2 )
88  call send_recv( t_s, 2 )
89  call send_recv( x_g, 2 )
90  call send_recv( x_s, 2 )
91  IF(granular_energy) call send_recv( theta_m, 2 )
92 
93  abort = .false.
94 
95 ! Check whether all field variables are initialized in all fluid cells
96 ! and flow boundary cells
97 ! ---------------------------------------------------------------->>>
98  DO k = kstart2, kend2
99  DO j = jstart2, jend2
100  DO i = istart2, iend2
101  ijk = funijk(i,j,k)
102  IF (.NOT.wall_at(ijk)) THEN
103 
104  imjk = im_of(ijk)
105  ijmk = jm_of(ijk)
106  ijkm = km_of(ijk)
107 
108 ! check gas phase fields
109  IF(ep_g(ijk) == undefined) &
110  CALL report_error(abort, i, j, k, 'EP_G')
111  IF(p_g(ijk) == undefined) &
112  CALL report_error(abort, i, j, k, 'P_G')
113  IF(p_star(ijk) == undefined) &
114  CALL report_error(abort, i, j, k, 'P_STAR')
115  IF(ro_g(ijk) == undefined) &
116  CALL report_error(abort, i, j, k, 'RO_G')
117  IF(rop_g(ijk) == undefined) &
118  CALL report_error(abort, i, j, k, 'ROP_G')
119 
120  IF(u_g(ijk) == undefined) &
121  CALL report_error(abort, i, j, k, 'U_G')
122  IF(v_g(ijk) == undefined) &
123  CALL report_error(abort, i, j, k, 'V_G')
124  IF(w_g(ijk) == undefined) &
125  CALL report_error(abort, i, j, k, 'W_G')
126 
127  IF(u_g(imjk) == undefined) &
128  CALL report_error(abort, i-1, j, k, 'U_G')
129  IF(v_g(ijmk) == undefined) &
130  CALL report_error(abort, i, j-1, k, 'V_G')
131  IF(w_g(ijkm) == undefined) &
132  CALL report_error(abort, i, j, k-1, 'W_G')
133 
134  IF(t_g(ijk) == undefined) THEN
135  IF(energy_eq .OR. ro_g0==undefined .OR. &
136  mu_g0==undefined) &
137  CALL report_error(abort, i, j, k, 'T_G')
138  ENDIF
139 
140  IF(species_eq(0) .OR. ro_g0==undefined .AND. &
141  mw_avg==undefined) THEN
142  DO nn = 1, nmax(0)
143  IF(x_g(ijk,nn) == undefined) &
144  CALL report_error(abort, i, j, k, 'X_G',nn)
145  ENDDO
146  ENDIF
147 
148  DO nn = 1, nscalar
149  IF(scalar(ijk,nn) == undefined) &
150  CALL report_error(abort, i, j, k, 'SCALAR',nn)
151  ENDDO
152 
153 ! check solids phase fields. these quantities are specified via the
154 ! subroutines set_ic and set_bc0/set_bc1 that employ the initial and
155 ! boundary conditions set in the mfix.dat.
156  DO m = 1, smax
157 
158  IF (ro_s(ijk,m) == undefined) &
159  CALL report_error(abort, i, j, k, 'RO_S',m)
160 
161  IF (rop_s(ijk,m) == undefined) &
162  CALL report_error(abort, i, j, k, 'ROP_S',m)
163 
164  IF (u_s(ijk,m)==undefined .AND. i/=imax2) &
165  CALL report_error(abort, i, j, k, 'U_S',m)
166  IF (v_s(ijk,m)==undefined .AND. j/=jmax2) &
167  CALL report_error(abort, i, j, k, 'V_S',m)
168  IF (w_s(ijk,m)==undefined .AND. k/=kmax2) &
169  CALL report_error(abort, i, j, k, 'W_S',m)
170 
171  IF (u_s(imjk,m) == undefined) &
172  CALL report_error(abort, i, j, k, 'U_S',m)
173  IF (v_s(ijmk,m) == undefined) &
174  CALL report_error(abort, i, j, k, 'V_S',m)
175  IF (w_s(ijkm,m) == undefined) &
176  CALL report_error(abort, i, j, k, 'W_S',m)
177 
178  IF (t_s(ijk,m) == undefined) THEN
179  IF (energy_eq) &
180  CALL report_error(abort, i, j, k, 'T_S',m)
181  ENDIF
182 
183  IF (species_eq(m)) THEN
184  DO nn = 1, nmax(m)
185  IF(x_s(ijk,m,nn) == undefined) &
186  CALL report_error(abort, i, j, k, 'X_S', m, nn)
187  ENDDO
188  ENDIF
189  ENDDO ! end do m=1,smax
190 
191  ENDIF ! IF (.NOT.WALL_AT(IJK)) THEN
192  ENDDO ! end do I = istart2, iend2
193  ENDDO ! end do J = jstart2, jend2
194  ENDDO ! end do K = kstart2, kend2
195 
196 
197  CALL global_all_or(abort)
198  IF(abort) THEN
199  WRITE(err_msg, 2000)
200  CALL flush_err_msg(header=.false., abort=.true.)
201  CALL finl_err_msg
202  RETURN
203  ENDIF
204 
205 
206 ! Additional check for fluid or mass inflow cells
207 ! --------------------------------------------------------------------//
208  nonzero = .false.
209 
210  DO k = kstart2, kend2
211  DO j = jstart2, jend2
212  DO i = istart2, iend2
213  ijk = funijk(i,j,k)
214 
215  IF (flag(ijk)==1 .OR. flag(ijk)==20) THEN
216 
217 ! Check whether L_scale is non-zero anywhere
218  IF (l_scale(ijk) /= zero) nonzero = .true.
219 
220 ! Ep_g must have a value > 0 and < 1
221  IF(ep_g(ijk) < small_number .OR. ep_g(ijk) > one) &
222  CALL report_unphysical(abort, i, j, k, 'EP_G', ep_g(ijk))
223 
224 ! Check the sum of volume fractions. This is skipped for DES as the
225 ! solids volume fractions may not yet be calculated.
226  IF(.NOT.discrete_element .AND. smax>0) THEN
227  dif = one - ep_g(ijk)
228  dif = dif - sum(rop_s(ijk,:mmax)/ro_s(ijk,:mmax))
229  IF (abs(dif) > small_number) CALL report_unphysical( &
230  abort, i, j, k, 'Volume Fraction SUM', 1.0-dif)
231  ENDIF
232  ENDIF ! IF (FLAG(IJK)==1 .OR. FLAG(IJK)==20) THEN
233  ENDDO ! I = istart2, iend2
234  ENDDO ! J = jstart2, jend2
235  ENDDO ! K = kstart2, kend2
236 
237 
238  CALL global_all_or(abort)
239  IF(abort) THEN
240  WRITE(err_msg, 2000)
241  CALL flush_err_msg(header=.false., abort=.true.)
242  CALL finl_err_msg
243  RETURN
244  ENDIF
245 
246 ! Check whether MU_gmax is specified
247  CALL global_all_or(nonzero)
248  IF (nonzero .AND. mu_gmax==undefined) THEN
249  WRITE(err_msg, 1300)
250  CALL flush_err_msg(abort=.true.)
251  CALL finl_err_msg
252  ENDIF
253 
254  1300 FORMAT('Error 1300: Message: Turbulent length scale is nonzero.',&
255  /'Specify MU_gmax in the mfix.dat file.')
256 
257  CALL finl_err_msg
258  RETURN
259 
260  2000 FORMAT('Please correct the mfix.dat file.')
261 
262 
263 
264  CONTAINS
265 
266 
267  SUBROUTINE report_error(ABORT, pI, pJ, pK, VAR, LC1, LC2)
269  LOGICAL, INTENT(INOUT) :: ABORT
270  INTEGER, INTENT(IN) :: pI, pJ, pK
271  CHARACTER(LEN=*), INTENT(IN) :: VAR
272  INTEGER, INTENT(IN), OPTIONAL :: LC1, LC2
273  CHARACTER(LEN=32) :: VAR_FULL
274 
275  var_full=''
276  IF(PRESENT(lc2)) THEN
277  var_full = ivar(var,lc1,lc2)
278  ELSEIF(PRESENT(lc1)) THEN
279  var_full = ivar(var,lc1)
280  ELSE
281  var_full = var
282  ENDIF
283 
284  IF(.NOT.abort) THEN
285  WRITE(err_msg,1000)
286  CALL flush_err_msg(footer=.false.)
287  abort = .true.
288  ENDIF
289 
290  1000 FORMAT('Error 1000: The following field variables are undefined')
291 
292 
293  WRITE(err_msg, 1010) i, j, k, trim(var_full)
294  CALL flush_err_msg(header=.false., footer=.false.)
295 
296  1010 FORMAT(1x,'I = ',i6,' J = ',i6,' K = ',i6,5x,a)
297 
298 
299 
300  END SUBROUTINE report_error
301 
302  SUBROUTINE report_unphysical(ABORT, pI, pJ, pK, VAR, VALUE)
304  LOGICAL, INTENT(INOUT) :: ABORT
305  INTEGER, INTENT(IN) :: pI, pJ, pK
306  CHARACTER(LEN=*), INTENT(IN) :: VAR
307  DOUBLE PRECISION, INTENT(IN) :: VALUE
308 
309  IF(.NOT.abort) THEN
310  WRITE(err_msg,1100)
311  CALL flush_err_msg(footer=.false.)
312  abort = .true.
313  ENDIF
314 
315  1100 FORMAT('Error 1100: The following field variables are ',&
316  'out of range')
317 
318  WRITE(err_msg, 1110) i, j, k, trim(var), VALUE
319  CALL flush_err_msg(header=.false., footer=.false.)
320 
321  1110 FORMAT(1x,'I = ',i6,' J = ',i6,' K = ',i6,2x,a,'Value:',g11.4)
322 
323 
324  END SUBROUTINE report_unphysical
325 
326 
327  END SUBROUTINE check_data_20
integer jend2
Definition: compar_mod.f:80
double precision, dimension(:,:), allocatable v_s
Definition: fldvar_mod.f:105
integer imax2
Definition: geometry_mod.f:61
character(len=32) function ivar(VAR, i1, i2, i3)
double precision, dimension(:), allocatable ep_g
Definition: fldvar_mod.f:15
subroutine finl_err_msg
double precision, parameter one
Definition: param1_mod.f:29
subroutine check_data_20
Definition: check_data_20.f:27
integer istart2
Definition: compar_mod.f:80
double precision, dimension(:,:), allocatable w_s
Definition: fldvar_mod.f:117
Definition: rxns_mod.f:1
logical, dimension(0:dim_m) species_eq
Definition: run_mod.f:115
integer iend2
Definition: compar_mod.f:80
double precision, dimension(:), allocatable t_g
Definition: fldvar_mod.f:63
double precision mu_g0
Definition: physprop_mod.f:62
double precision, dimension(:,:), allocatable scalar
Definition: fldvar_mod.f:155
double precision, dimension(:), allocatable p_g
Definition: fldvar_mod.f:26
double precision mu_gmax
Definition: constant_mod.f:180
double precision, parameter undefined
Definition: param1_mod.f:18
integer kend2
Definition: compar_mod.f:80
double precision, dimension(:), allocatable a
Definition: scalars_mod.f:29
integer kstart2
Definition: compar_mod.f:80
double precision, dimension(:), allocatable l_scale
Definition: visc_g_mod.f:24
double precision, dimension(:,:), allocatable u_s
Definition: fldvar_mod.f:93
subroutine init_err_msg(CALLER)
integer mmax
Definition: physprop_mod.f:19
double precision, dimension(:,:,:), allocatable x_s
Definition: fldvar_mod.f:78
double precision ro_g0
Definition: physprop_mod.f:59
double precision, parameter small_number
Definition: param1_mod.f:24
integer jmax2
Definition: geometry_mod.f:63
double precision, dimension(:,:), allocatable t_s
Definition: fldvar_mod.f:66
double precision, dimension(:,:), allocatable x_g
Definition: fldvar_mod.f:75
integer jstart2
Definition: compar_mod.f:80
double precision, dimension(:,:), allocatable theta_m
Definition: fldvar_mod.f:149
double precision, dimension(:), allocatable v_g
Definition: fldvar_mod.f:99
integer kmax2
Definition: geometry_mod.f:65
double precision, dimension(:), allocatable w_g
Definition: fldvar_mod.f:111
Definition: run_mod.f:13
subroutine report_unphysical(ABORT, pI, pJ, pK, VAR, VALUE)
Definition: param_mod.f:2
double precision, dimension(:,:), allocatable ro_s
Definition: fldvar_mod.f:45
integer, dimension(0:dim_m) nmax
Definition: physprop_mod.f:119
double precision, dimension(:), allocatable p_star
Definition: fldvar_mod.f:142
double precision mw_avg
Definition: physprop_mod.f:71
logical energy_eq
Definition: run_mod.f:100
character(len=line_length), dimension(line_count) err_msg
double precision, dimension(:), allocatable u_g
Definition: fldvar_mod.f:87
integer nscalar
Definition: scalars_mod.f:7
double precision, dimension(:,:), allocatable rop_s
Definition: fldvar_mod.f:51
integer smax
Definition: physprop_mod.f:22
integer, dimension(:), allocatable flag
Definition: geometry_mod.f:99
subroutine report_error(ABORT, pI, pJ, pK, VAR, LC1, LC2)
logical granular_energy
Definition: run_mod.f:112
double precision, dimension(:), allocatable ro_g
Definition: fldvar_mod.f:32
double precision, dimension(:), allocatable rop_g
Definition: fldvar_mod.f:38
double precision, dimension(:), allocatable x
Definition: geometry_mod.f:129
double precision, parameter zero
Definition: param1_mod.f:27
subroutine flush_err_msg(DEBUG, HEADER, FOOTER, ABORT, LOG, CALL_TREE)