MFIX  2016-1
check_bc_inflow.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! !
3 ! Subroutine: CHECK_BC_INFLOW !
4 ! Author: J.Musser Date: 01-Mar-14 !
5 ! !
6 ! Purpose: Provided a detailed error message on common inflow BC !
7 ! !
8 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
9  SUBROUTINE check_bc_inflow(M_TOT, SKIP, BCV)
10 
11 ! Modules
12 !---------------------------------------------------------------------//
13  use bc, only: bc_x_g, bc_x_s
14  use bc, only: bc_t_g, bc_t_s, bc_theta_m
15  use bc, only: bc_k_turb_g, bc_e_turb_g
16  use bc, only: bc_scalar
17  use param, only: dim_m
18  use param1, only: undefined, one, zero
19  use physprop, only: inert_species
20  use physprop, only: mu_g0
21  use physprop, only: mw_avg, nmax
22  use physprop, only: ro_g0
25  use scalars, only: nscalar
26  use toleranc, only: compare
27 
28  use error_manager
29  IMPLICIT NONE
30 
31 ! Dummy arguments
32 !---------------------------------------------------------------------//
33  INTEGER, INTENT(in) :: BCV
34  INTEGER, INTENT(in) :: M_TOT
35  LOGICAL, INTENT(in) :: SKIP(dim_m)
36 
37 ! Local variables
38 !---------------------------------------------------------------------//
39 ! loop/variable indices
40  INTEGER :: M, N
41  DOUBLE PRECISION SUM_X
42 ! Index of inert species
43  INTEGER :: INERT
44 !---------------------------------------------------------------------//
45 
46  CALL init_err_msg("CHECK_BC_INFLOW")
47 
48 ! Check temperature dependency.
49  IF((energy_eq .OR. ro_g0==undefined .OR.mu_g0==undefined) .AND. &
50  bc_t_g(bcv)==undefined) THEN
51  WRITE(err_msg, 1000) trim(ivar('BC_T_g',bcv))
52  CALL flush_err_msg(abort=.true.)
53  ENDIF
54 
55 ! Sum together defined gas phase species mass fractions.
56  sum_x = zero
57  DO n = 1, nmax(0)
58  IF(bc_x_g(bcv,n) /= undefined) THEN
59  sum_x = sum_x + bc_x_g(bcv,n)
60  ELSE
61  bc_x_g(bcv,n) = zero
62  ENDIF
63  ENDDO
64 
65 ! Enforce that the species mass fractions must sum to one.
66  IF(.NOT.compare(one,sum_x)) THEN
67 
68  IF(species_eq(0)) THEN
69  WRITE(err_msg, 1110) bcv
70  CALL flush_err_msg(abort=.true.)
71  1110 FORMAT('Error 1110: BC_X_g(',i3,',:) do NOT sum to ONE and the ',&
72  'gas phase',/'species equations are solved. Please correct ', &
73  'the mfix.dat file.')
74 
75  ELSEIF(ro_g0 == undefined .AND. mw_avg == undefined) THEN
76  WRITE(err_msg, 1111) bcv
77  CALL flush_err_msg(abort=.true.)
78  1111 FORMAT('Error 1111: BC_X_g(',i3,',:) do NOT sum to ONE and the ',&
79  'gas phase',/'is compressible and MW_AVG is UNDEFINED.',/ &
80  'Please correct the mfix.dat the mfix.dat file.')
81 
82  ELSEIF(.NOT.compare(sum_x,zero)) THEN
83  WRITE(err_msg, 1112) bcv
84  CALL flush_err_msg(abort=.true.)
85  1112 FORMAT('Error 1112: BC_X_g(',i3,',:) do not sum to ONE or ZERO ',&
86  'and they',/'are not needed. Please correct the mfix.dat ', &
87  'the mfix.dat file.')
88 
89  ELSE
90  bc_x_g(bcv,:) = zero
91  bc_x_g(bcv,1) = one
92  ENDIF
93  ENDIF
94 
95 
96 ! Verify that species mass fractions are defined for mass flow BCs when
97 ! using variable solids density. Needed to calculation RO_s
98  DO m = 1, m_tot
99 
100 ! If this phase is not present, clear out x_s for the BC and
101 ! cycle the solids loop. No need to continue checks.
102  IF(skip(m)) THEN
103  IF(species_eq(m))THEN
104  bc_x_s(bcv,m,:) = zero
105  bc_x_s(bcv,m,1) = one
106  ENDIF
107  cycle
108  ENDIF
109 
110 ! Sum together defined species mass fractions.
111  sum_x = zero
112  DO n = 1, nmax(m)
113  IF(bc_x_s(bcv,m,n) /= undefined) THEN
114  sum_x = sum_x + bc_x_s(bcv,m,n)
115  ELSE
116  bc_x_s(bcv,m,n) = zero
117  ENDIF
118  ENDDO
119 
120 ! Enforce that the species mass fractions must sum to one.
121  IF(.NOT.compare(one,sum_x)) THEN
122  IF(species_eq(m)) THEN
123  WRITE(err_msg, 1210) bcv, m
124  CALL flush_err_msg(abort=.true.)
125  1210 FORMAT('Error 1210: BC_X_s(',i3,',',i2,':) do NOT sum to ONE ', &
126  'and the solids phase',/'species equations are solved. ', &
127  'Please correct the mfix.dat file.')
128 
129  ELSEIF(solve_ros(m)) THEN
130  WRITE(err_msg, 1211) bcv, m
131  CALL flush_err_msg(abort=.true.)
132  1211 FORMAT('Error 1211: BC_X_s(',i3,',',i2,':) do NOT sum to ONE ', &
133  'and the solids phase',/'density is calculated. Please ', &
134  'correct the mfix.dat file.')
135 
136  ELSEIF(.NOT.compare(sum_x,zero)) THEN
137  WRITE(err_msg, 1212) bcv, m
138  CALL flush_err_msg(abort=.true.)
139  1212 FORMAT('Error 1212: BC_X_s(',i3,',',i2,':) do NOT sum to ONE ', &
140  'or ZERO and',/'they are not needed. Please correct the ', &
141  'mfix.dat file.')
142 
143  ELSE
144  bc_x_s(bcv,m,:) = zero
145  bc_x_s(bcv,m,1) = one
146  ENDIF
147  ENDIF
148 
149 ! Set the solids density for the BC region.
150  IF(solve_ros(m)) THEN
151 ! Verify that the species mass fraction for the inert material is not
152 ! zero in the IC region when the solids is present.
153  inert = inert_species(m)
154  IF(bc_x_s(bcv,m,inert) == zero) THEN
155  WRITE(err_msg,1213) m, bcv
156  CALL flush_err_msg(abort=.true.)
157  1213 FORMAT('Error 1213: No inert species for phase ',i2,' in BC ', &
158  'region',i3,'.',/'Unable to calculate solids phase density. ',&
159  'Please refer to the Readme',/' file for required variable ', &
160  'solids density model input parameters and',/' make the ', &
161  'necessary corrections to the mfix.dat file.')
162 
163  ENDIF
164  ENDIF
165  ENDDO
166 
167 
168  DO m = 1, m_tot
169 ! Check solids phase temperature dependency.
170  IF(energy_eq .AND. bc_t_s(bcv,m)==undefined) THEN
171  IF(skip(m)) THEN
172  bc_t_s(bcv,m) = bc_t_g(bcv)
173  ELSE
174  WRITE(err_msg, 1000) trim(ivar('BC_T_s',bcv,m))
175  CALL flush_err_msg(abort=.true.)
176  ENDIF
177  ENDIF
178 
179 ! Check granular energy dependency
180  IF(granular_energy) THEN
181  IF(bc_theta_m(bcv,m) == undefined) THEN
182  IF(skip(m) .OR. solids_model(m) /= 'TFM') THEN
183  bc_theta_m(bcv,m) = zero
184  ELSE
185  WRITE(err_msg,1000) trim(ivar('BC_Theta_m',bcv,m))
186  CALL flush_err_msg(abort=.true.)
187  ENDIF
188  ENDIF
189  ENDIF
190  ENDDO
191 
192 ! Check K-Epsilon BCs.
193  IF(k_epsilon) THEN
194  IF(bc_k_turb_g(bcv) == undefined) THEN
195  WRITE(err_msg, 1000) trim(ivar('BC_K_Turb_G',bcv))
196  CALL flush_err_msg(abort=.true.)
197 
198  ELSEIF(bc_e_turb_g(bcv) == undefined) THEN
199  WRITE(err_msg, 1000) trim(ivar('BC_E_Turb_G',bcv))
200  CALL flush_err_msg(abort=.true.)
201  ENDIF
202  ENDIF
203 
204 ! Check scalar equation BCs.
205  DO n = 1, nscalar
206  IF(bc_scalar(bcv,n) == undefined) THEN
207  WRITE(err_msg, 1001) trim(ivar('BC_Scalar',bcv,n))
208  CALL flush_err_msg(abort=.true.)
209  ENDIF
210  ENDDO
211 
212  CALL finl_err_msg
213 
214  RETURN
215 
216  1000 FORMAT('Error 1000: Required input not specified: ',a,/'Please ',&
217  'correct the mfix.dat file.')
218 
219  1001 FORMAT('Error 1001: Illegal or unknown input: ',a,' = ',a,/ &
220  'Please correct the mfix.dat file.')
221 
222  END SUBROUTINE check_bc_inflow
223 
224 
225 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
226 ! !
227 ! Subroutine: CHECK_BC_MASS_INFLOW !
228 ! Author: J.Musser Date: 01-Mar-14 !
229 ! !
230 ! Purpose: Provided a detailed error message on BC !
231 ! !
232 ! Comments: !
233 ! The velocities at the inflow face are fixed and the momentum !
234 ! equations are not solved in the inflow cells. Since the flow is !
235 ! into the domain all other scalars that are used need to be !
236 ! specified (e.g., mass fractions, void fraction, etc.,) !
237 ! !
238 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
239  SUBROUTINE check_bc_mass_inflow(M_TOT, SKIP, BCV)
241 ! Modules
242 !---------------------------------------------------------------------//
243  use bc, only: bc_ep_g, bc_p_g
244  use bc, only: bc_rop_s, bc_ep_s, bc_x_s
245  use bc, only: bc_massflow_g
246  use eos, only: eoss
247  use param, only: dim_m
248  use param1, only: undefined, one, zero
249  use physprop, only: ro_g0, ro_s0
250  use physprop, only: inert_species, x_s0
251  use run, only: solve_ros
252  use toleranc, only: compare
253  use error_manager
254  IMPLICIT NONE
255 
256 ! Dummy arguments
257 !---------------------------------------------------------------------//
258  INTEGER, INTENT(in) :: BCV
259  INTEGER, INTENT(in) :: M_TOT
260  LOGICAL, INTENT(in) :: SKIP(dim_m)
261 
262 ! Local variables
263 !---------------------------------------------------------------------//
264 ! loop/variable indices
265  INTEGER :: M
266  DOUBLE PRECISION :: SUM_EP
267 ! Solids phase density in BC region.
268  DOUBLE PRECISION :: BC_ROs(dim_m)
269 ! Index of inert species
270  INTEGER :: INERT
271 
272 !---------------------------------------------------------------------//
273 
274  CALL init_err_msg("CHECK_BC_MASS_INFLOW")
275 
276 ! Check gas phase volume fraction.
277  IF(bc_ep_g(bcv) == undefined) THEN
278  WRITE(err_msg, 1000) trim(ivar('BC_EP_g',bcv))
279  CALL flush_err_msg(abort=.true.)
280  ENDIF
281 
282 ! Verify compressible boundary condition variables.
283  IF(ro_g0 == undefined) THEN
284  IF(bc_p_g(bcv) == undefined) THEN
285  IF(bc_massflow_g(bcv) /= undefined .AND. &
286  bc_massflow_g(bcv) /= zero) THEN
287  WRITE(err_msg, 1100) trim(ivar('BC_P_g',bcv))
288  CALL flush_err_msg(abort=.true.)
289  ENDIF
290  1100 FORMAT('Error 1100: ',a,' must be specified for compressible ', &
291  'flows',/'when specifying BC_MASSFLOW_g to make the ', &
292  'conversion to velocity.',/'Please correct the mfix.dat file.')
293 
294  ELSEIF(bc_p_g(bcv) <= zero) THEN
295  WRITE(err_msg, 1101) bcv, trim(ival(bc_p_g(bcv)))
296  CALL flush_err_msg(abort=.true.)
297  ENDIF
298  1101 FORMAT('Error 1101: Pressure must be greater than zero for ', &
299  'compressible flow',/' >>> BC_P_g(',i3,') = ',a,/'Please ', &
300  'correct the mfix.dat file.')
301  ENDIF
302 
303 ! Calculate the solids volume fraction from the gas phase if there is
304 ! only one solids phase.
305  IF(m_tot == 1 .AND. bc_ep_s(bcv,1) == undefined) THEN
306  bc_ep_s(bcv,1) = one - bc_ep_g(bcv)
307  ENDIF
308 
309 ! Bulk density or solids volume fraction must be explicitly defined
310 ! if there are more than one solids phase.
311  IF(m_tot > 1 .AND. .NOT.compare(bc_ep_g(bcv),one)) THEN
312  DO m = 1, m_tot
313  IF(bc_rop_s(bcv,m) == undefined .AND. &
314  bc_ep_s(bcv,m) == undefined) THEN
315  WRITE(err_msg, 1200) m, bcv, 'BC_ROP_s and BC_EP_s'
316  CALL flush_err_msg(abort=.true.)
317  ENDIF
318  ENDDO
319  ENDIF
320  1200 FORMAT('Error 1200: Insufficient solids phase ',i2,' data ', &
321  'for BC',i3,'. ',/a,' not specified.',/'Please correct the ', &
322  'mfix.dat file.')
323 
324 ! Initialize the sum of the total volume fraction.
325  sum_ep = bc_ep_g(bcv)
326  DO m = 1, m_tot
327 
328 ! If this phase is not present, clear out EPs and ROPs for the BC and
329 ! cycle the solids loop. No need to continue checks.
330  IF(skip(m)) THEN
331  bc_ep_s(bcv,m) = zero
332  bc_rop_s(bcv,m) = zero
333  cycle
334  ENDIF
335 
336 ! Set the solids density for the BC region.
337  IF(.NOT.solve_ros(m)) THEN
338  bc_ros(m) = ro_s0(m)
339  ELSE
340 ! presence of non-zero inert species is checked by bc_inflow
341  inert = inert_species(m)
342  bc_ros(m) = eoss(ro_s0(m), x_s0(m,inert),&
343  bc_x_s(bcv,m,inert))
344  ENDIF
345 
346 ! If both input parameters are defined. Make sure they are equivalent.
347  IF(bc_rop_s(bcv,m) /= undefined .AND. &
348  bc_ep_s(bcv,m) /= undefined) THEN
349 
350  IF(.NOT.compare(bc_ep_s(bcv,m)*bc_ros(m), &
351  bc_rop_s(bcv,m))) THEN
352  WRITE(err_msg,1214) bcv
353  CALL flush_err_msg(abort=.true.)
354  ENDIF
355  1214 FORMAT('Error 1214: Illegal initial condition region : ',i3,/ &
356  'BC_EP_s and BC_ROP_s are inconsistent. Please correct the ',/&
357  'mfix.dat file.')
358 
359 ! Compute BC_EP_s from BC_ROP_s
360  ELSEIF(bc_ep_s(bcv,m) == undefined) THEN
361  bc_ep_s(bcv,m) = bc_rop_s(bcv,m) / bc_ros(m)
362 
363 ! Compute BC_ROP_s from BC_EP_s and BC_ROs
364  ELSEIF(bc_rop_s(bcv,m) == undefined) THEN
365  bc_rop_s(bcv,m) = bc_ep_s(bcv,m) * bc_ros(m)
366 
367  ENDIF
368 ! Add this phase to the total volume fraction.
369  sum_ep = sum_ep + bc_ep_s(bcv,m)
370  ENDDO
371 
372 ! Verify that the volume fractions sum to one.
373  IF(.NOT.compare(sum_ep,one)) THEN
374  WRITE(err_msg,1215) bcv, trim(ival(sum_ep))
375  CALL flush_err_msg(abort=.true.)
376  ENDIF
377  1215 FORMAT('Error 1215: Illegal boundary condition region: ',i3,'. ',&
378  'Sum of volume',/'fractions does NOT equal ONE. (SUM = ',a, &
379  ')',/'Please correct the mfix.dat file.')
380 
381  CALL finl_err_msg
382 
383  RETURN
384 
385  1000 FORMAT('Error 1000: Required input not specified: ',a,/'Please ',&
386  'correct the mfix.dat file.')
387 
388  1001 FORMAT('Error 1001: Illegal or unknown input: ',a,' = ',a,/ &
389  'Please correct the mfix.dat file.')
390 
391  END SUBROUTINE check_bc_mass_inflow
392 
393 
394 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
395 ! !
396 ! Subroutine: CHECK_BC_P_INFLOW !
397 ! Author: J.Musser Date: 01-Mar-14 !
398 ! !
399 ! Purpose: Provided detailed error message on bc !
400 ! !
401 ! Comments: !
402 ! Unlike the MI boundary, for the PI boundary the velocities at !
403 ! the inflow face are calculated by solving the momentum eqns !
404 ! and are not fixed. In this way, the PI is similar to the PO !
405 ! except that the flow is into the domain and hence all other !
406 ! scalars (e.g., mass fractions, void fraction, temperature, !
407 ! etc.,) at the inflow cells need to be specified. To satisfy !
408 ! the error routines at the start of the simulation, both the !
409 ! tangential and normal components at the inflow also need to !
410 ! be specified. The velocities values essentially serve as IC. !
411 ! !
412 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
413  SUBROUTINE check_bc_p_inflow(M_TOT, SKIP, BCV)
415 ! Modules
416 !---------------------------------------------------------------------//
417  use bc, only: bc_p_g, bc_rop_s
418  use bc, only: bc_u_g, bc_v_g, bc_w_g
419  use bc, only: bc_u_s, bc_v_s, bc_w_s
420  use geometry, only: no_i, no_j, no_k
421  use param, only: dim_m
422  use param1, only: undefined, zero
423  use physprop, only: ro_g0
424  use error_manager
425  IMPLICIT NONE
426 
427 ! Dummy arguments
428 !---------------------------------------------------------------------//
429  INTEGER, INTENT(in) :: BCV
430  INTEGER, INTENT(in) :: M_TOT
431  LOGICAL, INTENT(in) :: SKIP(dim_m)
432 
433 ! Local variables
434 !---------------------------------------------------------------------//
435  INTEGER :: M
436 
437 !---------------------------------------------------------------------//
438 
439  CALL init_err_msg("CHECK_BC_P_INFLOW")
440 
441 ! Remove checks on bc_ep_g/bc_rop_s; using routine check_bc_outflow
442 
443  IF (bc_p_g(bcv) == undefined) THEN
444  WRITE(err_msg,1000) 'BC_P_g', bcv
445  CALL flush_err_msg(abort=.true.)
446  1000 FORMAT('Error 1000: Required input not specified: ',a,/'Please ',&
447  'correct the mfix.dat file.')
448 
449  ELSEIF (bc_p_g(bcv)<=zero .AND. ro_g0==undefined) THEN
450  WRITE(err_msg, 1101) bcv, trim(ival(bc_p_g(bcv)))
451  CALL flush_err_msg(abort=.true.)
452  1101 FORMAT('Error 1101: Pressure must be greater than zero for ', &
453  'compressible flow',/' >>> BC_P_g(',i3,') = ',a,/'Please ', &
454  'correct the mfix.dat file.')
455  ENDIF
456 
457 ! Check that velocities are also specified. These are essentially used
458 ! as initial conditions for the boundary region. If they are not
459 ! specified then a default value is set here otherwise check_data_20
460 ! will complain and cause MFIX to exit.
461  IF(bc_u_g(bcv) == undefined) THEN
462  bc_u_g(bcv) = zero
463  IF(.NOT.no_i) WRITE(err_msg, 1300) trim(ivar('BC_U_g',bcv))
464  ENDIF
465 
466  IF(bc_v_g(bcv) == undefined) THEN
467  bc_v_g(bcv) = zero
468  IF(.NOT.no_j) WRITE(err_msg, 1300) trim(ivar('BC_V_g',bcv))
469  ENDIF
470 
471  IF(bc_w_g(bcv) == undefined) THEN
472  bc_w_g(bcv) = zero
473  IF(.NOT.no_k) WRITE(err_msg, 1300) trim(ivar('BC_W_g',bcv))
474  ENDIF
475 
476  DO m = 1, m_tot
477  IF (skip(m)) THEN
478  bc_u_s(bcv,m) = zero
479  bc_v_s(bcv,m) = zero
480  bc_w_s(bcv,m) = zero
481  ELSE
482  IF(bc_u_s(bcv,m) == undefined) THEN
483  bc_u_s(bcv,m) = zero
484  IF(bc_rop_s(bcv,m) /= zero .AND. .NOT.no_i) &
485  WRITE(err_msg, 1300) trim(ivar('BC_U_s',bcv,m))
486  ENDIF
487 
488  IF(bc_v_s(bcv,m) == undefined) THEN
489  bc_v_s(bcv,m) = zero
490  IF(bc_rop_s(bcv,m) /= zero .AND. .NOT.no_j) &
491  WRITE(err_msg, 1300) trim(ivar('BC_V_s',bcv,m))
492  ENDIF
493 
494  IF(bc_w_s(bcv,m) == undefined) THEN
495  bc_w_s(bcv,m) = zero
496  IF(bc_rop_s(bcv,m) /= zero .AND. .NOT.no_k) &
497  WRITE(err_msg, 1300) trim(ivar('BC_W_s',bcv,m))
498  ENDIF
499  ENDIF
500  ENDDO
501 
502  1300 FORMAT('Warning 1300: ',a,' was undefined. This variable was ', &
503  'set ',/ 'to zero to be used as the initial value in the BC ',&
504  'region.')
505 
506  CALL finl_err_msg
507 
508  RETURN
509  END SUBROUTINE check_bc_p_inflow
character(len=32) function ivar(VAR, i1, i2, i3)
double precision, dimension(dimension_bc) bc_t_g
Definition: bc_mod.f:97
subroutine finl_err_msg
logical function compare(V1, V2)
Definition: toleranc_mod.f:94
logical no_i
Definition: geometry_mod.f:20
double precision, parameter one
Definition: param1_mod.f:29
logical, dimension(0:dim_m) species_eq
Definition: run_mod.f:115
double precision, dimension(dimension_bc, dim_m) bc_w_s
Definition: bc_mod.f:129
double precision, dimension(dimension_bc, dim_m, dim_n_s) bc_x_s
Definition: bc_mod.f:254
double precision mu_g0
Definition: physprop_mod.f:62
integer, parameter dim_m
Definition: param_mod.f:67
logical, dimension(dim_m) solve_ros
Definition: run_mod.f:250
double precision, dimension(dim_m, dim_n_s) x_s0
Definition: physprop_mod.f:32
character(len=3), dimension(dim_m) solids_model
Definition: run_mod.f:253
double precision, parameter undefined
Definition: param1_mod.f:18
double precision, dimension(dimension_bc) bc_v_g
Definition: bc_mod.f:117
subroutine init_err_msg(CALLER)
double precision, dimension(dimension_bc, dim_m) bc_t_s
Definition: bc_mod.f:101
double precision ro_g0
Definition: physprop_mod.f:59
double precision function eoss(pBase, Xs0_INERT, Xs_INERT)
Definition: eos_mod.f:155
integer, dimension(dim_m) inert_species
Definition: physprop_mod.f:39
Definition: eos_mod.f:10
subroutine check_bc_mass_inflow(M_TOT, SKIP, BCV)
double precision, dimension(dimension_bc, dim_scalar) bc_scalar
Definition: bc_mod.f:384
subroutine check_bc_inflow(M_TOT, SKIP, BCV)
double precision, dimension(dimension_bc) bc_p_g
Definition: bc_mod.f:80
subroutine check_bc_p_inflow(M_TOT, SKIP, BCV)
Definition: run_mod.f:13
Definition: param_mod.f:2
logical no_k
Definition: geometry_mod.f:28
double precision, dimension(dimension_bc, dim_m) bc_v_s
Definition: bc_mod.f:121
double precision, dimension(dimension_bc) bc_massflow_g
Definition: bc_mod.f:201
integer, dimension(0:dim_m) nmax
Definition: physprop_mod.f:119
logical k_epsilon
Definition: run_mod.f:97
double precision, dimension(dimension_bc) bc_u_g
Definition: bc_mod.f:109
logical no_j
Definition: geometry_mod.f:24
double precision mw_avg
Definition: physprop_mod.f:71
logical energy_eq
Definition: run_mod.f:100
double precision, dimension(dimension_bc, dim_m) bc_u_s
Definition: bc_mod.f:113
character(len=line_length), dimension(line_count) err_msg
double precision, dimension(dimension_bc) bc_e_turb_g
Definition: bc_mod.f:404
integer nscalar
Definition: scalars_mod.f:7
double precision, dimension(dimension_bc, dim_m) bc_theta_m
Definition: bc_mod.f:105
double precision, dimension(dimension_bc, dim_n_g) bc_x_g
Definition: bc_mod.f:251
double precision, dimension(dimension_bc) bc_ep_g
Definition: bc_mod.f:77
double precision, dimension(dimension_bc) bc_k_turb_g
Definition: bc_mod.f:403
double precision, dimension(dim_m) ro_s0
Definition: physprop_mod.f:28
logical granular_energy
Definition: run_mod.f:112
double precision, dimension(dimension_bc) bc_w_g
Definition: bc_mod.f:125
double precision, dimension(dimension_bc, dim_m) bc_ep_s
Definition: bc_mod.f:93
double precision, parameter zero
Definition: param1_mod.f:27
subroutine flush_err_msg(DEBUG, HEADER, FOOTER, ABORT, LOG, CALL_TREE)
double precision, dimension(dimension_bc, dim_m) bc_rop_s
Definition: bc_mod.f:92
Definition: bc_mod.f:23