MFIX  2016-1
check_initial_conditions.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! !
3 ! Subroutine: CHECK_INITIAL_CONDITIONS !
4 ! Author: P. Nicoletti Date: 02-DEC-91 !
5 ! Author: J.Musser Date: 01-MAR-14 !
6 ! !
7 ! Purpose: check the initial conditions input section !
8 ! - check geometry of any specified IC region !
9 ! - check specification of physical quantities !
10 ! !
11 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
12  SUBROUTINE check_initial_conditions
13 
14 ! Global Variables:
15 !---------------------------------------------------------------------//
16 ! Flag: IC geometry was detected.
17  use ic, only: ic_defined
18 ! Flag: DEM solids present.
19  use run, only: dem_solids
20 ! Flag: New run or a restart.
21  use run, only: run_type
22 ! Runtime flag specifying MPPIC solids
23  use run, only: pic_solids
24 ! Flag: Adjusting runtime parameters
25  use run, only: reinitializing
26 
27 ! Global Parameters:
28 !---------------------------------------------------------------------//
29 ! Maximum number of IC.
30  use param, only: dimension_ic
31 
32 ! Use the error manager for posting error messages.
33 !---------------------------------------------------------------------//
34  use error_manager
35 
36  IMPLICIT NONE
37 
38 ! Local Variables:
39 !---------------------------------------------------------------------//
40 ! Loop counter for ICs
41  INTEGER :: ICV
42 !......................................................................!
43 
44 ! Initialize the error manager.
45  CALL init_err_msg("CHECK_INITIAL_CONDITIONS")
46 
47 ! Determine which ICs are DEFINED
49 
50 ! Loop over all IC arrays.
51  DO icv=1, dimension_ic
52 
53 ! Verify user input for defined defined IC.
54  IF(ic_defined(icv)) THEN
55 ! Gas phase checks.
56  CALL check_ic_gas_phase(icv)
57 ! Generic solids phase checks.
58  CALL check_ic_solids_phases(icv)
59 
60 ! Verify that no data was defined for unspecified IC. ICs are only
61 ! defined for new runs, so these checks are restricted to new runs.
62  ELSEIF(run_type == 'NEW' .AND. .NOT.reinitializing) THEN
63  CALL check_ic_overflow(icv)
64  ENDIF
65  ENDDO
66 
67 
68 
69 ! Check the initial conditions for the DEM and MPPIC models as well
70  IF(dem_solids.OR.pic_solids) &
72  IF(dem_solids) CALL check_ic_dem
73  IF(pic_solids) CALL check_ic_mppic
74 
75 ! Finalize the error manager.
76  CALL finl_err_msg
77 
78  RETURN
79  END SUBROUTINE check_initial_conditions
80 
81 
82 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
83 ! !
84 ! Subroutine: CHECK_IC_GEOMETRY !
85 ! Author: J.Musser Date: 01-Mar-14 !
86 ! !
87 ! Purpose: Provided a detailed error message when the sum of volume !
88 ! !
89 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
90  SUBROUTINE check_ic_geometry
91 
92 
93 ! Global Variables:
94 !---------------------------------------------------------------------//
95 ! Flag: IC contains geometric data and/or specified type
96  use ic, only: ic_defined
97 ! Flag: IC type.
98  use ic, only: ic_type
99 ! User specified: IC geometry
100  use ic, only: ic_x_e, ic_x_w, ic_i_e, ic_i_w
101  use ic, only: ic_y_n, ic_y_s, ic_j_n, ic_j_s
102  use ic, only: ic_z_t, ic_z_b, ic_k_t, ic_k_b
103 ! User specified: System geometry
104  use geometry, only: no_i, imax, imin1, imax1, xlength, dx, xmin
105  use geometry, only: no_j, jmax, jmin1, jmax1, ylength, dy
106  use geometry, only: no_k, kmax, kmin1, kmax1, zlength, dz
107 ! Flag: New run or a restart.
108  use run, only: run_type
109  use run, only: reinitializing
110 
111 ! Global Parameters:
112 !---------------------------------------------------------------------//
113 ! The max number of ICs.
114  use param, only: dimension_ic
115 ! Parameter constants
116  use param1, only: zero, undefined, undefined_i
117 
118 ! Use the error manager for posting error messages.
119 !---------------------------------------------------------------------//
120  use error_manager
121 
122 
123  implicit none
124 
125 
126 ! Local Variables:
127 !---------------------------------------------------------------------//
128 ! Loop/variable indices
129  INTEGER :: ICV
130 ! Local spatial indices.
131  INTEGER :: I_w, I_e, J_s, J_n, K_b, K_t
132 !......................................................................!
133 
134 
135 ! Initialize the error manager.
136  CALL init_err_msg("CHECK_IC_GEOMETRY")
137 
138 ! Check geometry of any specified IC region
139  DO icv = 1, dimension_ic
140 
141  ic_defined(icv) = .false.
142  IF (ic_x_w(icv) /= undefined) ic_defined(icv) = .true.
143  IF (ic_x_e(icv) /= undefined) ic_defined(icv) = .true.
144  IF (ic_y_s(icv) /= undefined) ic_defined(icv) = .true.
145  IF (ic_y_n(icv) /= undefined) ic_defined(icv) = .true.
146  IF (ic_z_b(icv) /= undefined) ic_defined(icv) = .true.
147  IF (ic_z_t(icv) /= undefined) ic_defined(icv) = .true.
148  IF (ic_i_w(icv) /= undefined_i) ic_defined(icv) = .true.
149  IF (ic_i_e(icv) /= undefined_i) ic_defined(icv) = .true.
150  IF (ic_j_s(icv) /= undefined_i) ic_defined(icv) = .true.
151  IF (ic_j_n(icv) /= undefined_i) ic_defined(icv) = .true.
152  IF (ic_k_b(icv) /= undefined_i) ic_defined(icv) = .true.
153  IF (ic_k_t(icv) /= undefined_i) ic_defined(icv) = .true.
154 
155 ! An IC is defined for restart runs only if it is a 'PATCH'.
156  IF(run_type /= 'NEW' .AND. ic_type(icv) /= 'PATCH') &
157  ic_defined(icv) = .false.
158 
159 ! Ignore patched IC regions for new runs. It may be better to flag this as
160 ! and error to avoid user confusion.
161  IF(run_type == 'NEW' .AND. ic_type(icv) == 'PATCH') &
162  ic_defined(icv) = .false.
163 
164 ! Enable only PATCH IC regions when initializing.
165  IF(reinitializing) ic_defined(icv)=(ic_type(icv)=='PATCH')
166 
167  IF(.NOT.ic_defined(icv)) cycle
168 
169  IF (ic_x_w(icv)==undefined .AND. ic_i_w(icv)==undefined_i) THEN
170  IF (no_i) THEN
171  ic_x_w(icv) = zero
172  ELSE
173  WRITE(err_msg, 1100) icv, 'IC_X_w and IC_I_w'
174  CALL flush_err_msg(abort=.true.)
175  ENDIF
176  ENDIF
177 
178  IF (ic_x_e(icv)==undefined .AND. ic_i_e(icv)==undefined_i) THEN
179  IF (no_i) THEN
180  ic_x_e(icv) = xlength
181  ELSE
182  WRITE(err_msg, 1100) icv, 'IC_X_e and IC_I_e'
183  CALL flush_err_msg(abort=.true.)
184  ENDIF
185  ENDIF
186 
187  IF (ic_y_s(icv)==undefined .AND. ic_j_s(icv)==undefined_i) THEN
188  IF (no_j) THEN
189  ic_y_s(icv) = zero
190  ELSE
191  WRITE(err_msg, 1100) icv, 'IC_Y_s and IC_J_s'
192  CALL flush_err_msg(abort=.true.)
193  ENDIF
194  ENDIF
195 
196  IF (ic_y_n(icv)==undefined .AND. ic_j_n(icv)==undefined_i) THEN
197  IF (no_j) THEN
198  ic_y_n(icv) = ylength
199  ELSE
200  WRITE(err_msg, 1100) icv, 'IC_Y_n and IC_J_n'
201  CALL flush_err_msg(abort=.true.)
202  ENDIF
203  ENDIF
204 
205  IF (ic_z_b(icv)==undefined .AND. ic_k_b(icv)==undefined_i) THEN
206  IF (no_k) THEN
207  ic_z_b(icv) = zero
208  ELSE
209  WRITE(err_msg, 1100) icv, 'IC_Z_b and IC_K_b'
210  CALL flush_err_msg(abort=.true.)
211  ENDIF
212  ENDIF
213 
214  IF (ic_z_t(icv)==undefined .AND. ic_k_t(icv)==undefined_i) THEN
215  IF (no_k) THEN
216  ic_z_t(icv) = zlength
217  ELSE
218  WRITE(err_msg, 1100) icv, 'IC_Z_t and IC_K_t'
219  CALL flush_err_msg(abort=.true.)
220  ENDIF
221  ENDIF
222 
223  ENDDO ! end loop over (icv = 1,dimension_ic)
224 
225  1100 FORMAT('Error 1100: Initial condition region ',i3,' is ill-', &
226  'defined.',/' > ',a,' are not specified.',/'Please correct ', &
227  'the mfix.dat file.')
228 
229 
230  DO icv = 1, dimension_ic
231 
232 ! Skip this check if the IC region is not specified.
233  IF(.NOT.ic_defined(icv)) cycle
234 
235  IF (ic_x_w(icv)/=undefined .AND. ic_x_e(icv)/=undefined) THEN
236  IF (no_i) THEN
237  i_w = 1
238  i_e = 1
239  ELSE
240  CALL calc_cell (xmin, ic_x_w(icv), dx, imax, i_w)
241  i_w = i_w + 1
242  CALL calc_cell (xmin, ic_x_e(icv), dx, imax, i_e)
243  ENDIF
244  IF (ic_i_w(icv)/=undefined_i .OR. ic_i_e(icv)/=undefined_i) THEN
245  CALL location_check (ic_i_w(icv), i_w, icv, 'IC - west')
246  CALL location_check (ic_i_e(icv), i_e, icv, 'IC - east')
247  ELSE
248  ic_i_w(icv) = i_w
249  ic_i_e(icv) = i_e
250  ENDIF
251  ENDIF
252 
253 ! Report problems with calculated bounds.
254  IF(ic_i_w(icv) > ic_i_e(icv)) THEN
255  WRITE(err_msg, 1101) icv, 'IC_I_W > IC_I_E'
256  write(*,*)' dump:',ic_i_w(icv),ic_i_e(icv)
257  CALL flush_err_msg(abort=.true.)
258  ELSEIF(ic_i_w(icv) < imin1) THEN
259  WRITE(err_msg, 1101) icv, 'IC_I_W < IMIN1'
260  CALL flush_err_msg(abort=.true.)
261  ELSEIF(ic_i_w(icv) > imax1) THEN
262  WRITE(err_msg, 1101) icv, 'IC_I_W > IMAX1'
263  CALL flush_err_msg(abort=.true.)
264  ELSEIF(ic_i_e(icv) < imin1) THEN
265  WRITE(err_msg, 1101) icv, 'IC_I_E < IMIN1'
266  CALL flush_err_msg(abort=.true.)
267  ELSEIF(ic_i_e(icv) > imax1) THEN
268  WRITE(err_msg, 1101) icv, 'IC_Z_t and IC_K_t'
269  CALL flush_err_msg(abort=.true.)
270  ENDIF
271 
272  IF (ic_y_s(icv)/=undefined .AND. ic_y_n(icv)/=undefined) THEN
273  IF (no_j) THEN
274  j_s = 1
275  j_n = 1
276  ELSE
277  CALL calc_cell (zero, ic_y_s(icv), dy, jmax, j_s)
278  j_s = j_s + 1
279  CALL calc_cell (zero, ic_y_n(icv), dy, jmax, j_n)
280  ENDIF
281  IF (ic_j_s(icv)/=undefined_i .OR. ic_j_n(icv)/=undefined_i) THEN
282  CALL location_check (ic_j_s(icv), j_s, icv, 'IC - south')
283  CALL location_check (ic_j_n(icv), j_n, icv, 'IC - north')
284  ELSE
285  ic_j_s(icv) = j_s
286  ic_j_n(icv) = j_n
287  ENDIF
288  ENDIF
289 
290  IF(ic_j_s(icv) > ic_j_n(icv)) THEN
291  WRITE(err_msg, 1101) icv, 'IC_J_S > IC_J_N'
292  CALL flush_err_msg(abort=.true.)
293  ELSEIF(ic_j_s(icv)<jmin1) THEN
294  WRITE(err_msg, 1101) icv, 'IC_J_S < JMIN1'
295  CALL flush_err_msg(abort=.true.)
296  ELSEIF(ic_j_s(icv)>jmax1) THEN
297  WRITE(err_msg, 1101) icv, 'IC_J_S > JMAX1'
298  CALL flush_err_msg(abort=.true.)
299  ELSEIF(ic_j_n(icv)<jmin1) THEN
300  WRITE(err_msg, 1101) icv, 'IC_J_N < JMIN1'
301  CALL flush_err_msg(abort=.true.)
302  ELSEIF(ic_j_n(icv)>jmax1) THEN
303  WRITE(err_msg, 1101) icv, 'IC_J_N > JMAX1'
304  CALL flush_err_msg(abort=.true.)
305  ENDIF
306 
307 
308  IF (ic_z_b(icv)/=undefined .AND. ic_z_t(icv)/=undefined) THEN
309  IF (no_k) THEN
310  k_b = 1
311  k_t = 1
312  ELSE
313  CALL calc_cell (zero, ic_z_b(icv), dz, kmax, k_b)
314  k_b = k_b + 1
315  CALL calc_cell (zero, ic_z_t(icv), dz, kmax, k_t)
316  ENDIF
317  IF (ic_k_b(icv)/=undefined_i .OR. ic_k_t(icv)/=undefined_i) THEN
318  CALL location_check (ic_k_b(icv), k_b, icv, 'IC - bottom')
319  CALL location_check (ic_k_t(icv), k_t, icv, 'IC - top')
320  ELSE
321  ic_k_b(icv) = k_b
322  ic_k_t(icv) = k_t
323  ENDIF
324  ENDIF
325 
326  IF(ic_k_b(icv) > ic_k_t(icv)) THEN
327  WRITE(err_msg, 1101) icv, 'IC_K_B > IC_K_T'
328  CALL flush_err_msg(abort=.true.)
329  ELSEIF(ic_k_b(icv) < kmin1) THEN
330  WRITE(err_msg, 1101) icv, 'IC_K_B < KMIN1'
331  CALL flush_err_msg(abort=.true.)
332  ELSEIF(ic_k_b(icv) > kmax1) THEN
333  WRITE(err_msg, 1101) icv, 'IC_K_B > KMAX1'
334  CALL flush_err_msg(abort=.true.)
335  ELSEIF(ic_k_t(icv) < kmin1) THEN
336  WRITE(err_msg, 1101) icv, 'IC_K_T < KMIN1'
337  CALL flush_err_msg(abort=.true.)
338  ELSEIF(ic_k_t(icv) > kmax1) THEN
339  WRITE(err_msg, 1101) icv, 'IC_K_T > KMAX1'
340  CALL flush_err_msg(abort=.true.)
341  ENDIF
342 
343 
344  1101 FORMAT('Error 1101: Initial condition region ',i2,' is ill-', &
345  'defined.',/3x,a,/'Please correct the mfix.dat file.')
346 
347  ENDDO ! end loop over (icv=1,dimension_ic)
348 
349  CALL finl_err_msg
350 
351  RETURN
352  END SUBROUTINE check_ic_geometry
353 
354 
355 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
356 ! !
357 ! Subroutine: CHECK_IC_GAS_PHASE !
358 ! Author: J.Musser Date: 01-Mar-14 !
359 ! !
360 ! Purpose: Verify gas phase input variables in IC region. !
361 ! !
362 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
363  SUBROUTINE check_ic_gas_phase(ICV)
365 ! Global Variables:
366 !---------------------------------------------------------------------//
367  use constant, only: l_scale0
368 ! Gas phase volume fraction, pressure, temperature, species.
369  use ic, only: ic_ep_g, ic_p_g, ic_t_g, ic_x_g
370 ! Gas phase velocity components.
371  use ic, only: ic_u_g, ic_v_g, ic_w_g
372 ! Radiation model parameters.
373  use ic, only: ic_gama_rg, ic_t_rg
374 ! K-Epsilon model parameters.
375  use ic, only: ic_k_turb_g, ic_e_turb_g
376 ! L-scale model parameters
377  use ic, only: ic_l_scale
378 ! IC for user-defined scalar equation.
379  use ic, only: ic_scalar
380 ! IC Type: UNDEFINED or PATCH.
381  use ic, only: ic_type
382 
383 ! Flag. Solve Energy equations
384  use run, only: energy_eq
385 ! Flag. Solve Species equations
386  use run, only: species_eq
387 ! Flag: Solve K-Epsilon; K-E parameters
388  use run, only: k_epsilon
389 ! Specified constant gas density and viscosity.
390  use physprop, only: ro_g0, mu_g0
391 ! Specified average molecular weight
392  use physprop, only: mw_avg
393 ! Number of gas phase species
394  use physprop, only: nmax
395 ! Specified number of scalar equations.
396  use scalars, only: nscalar
397 ! Flag: Do not solve in specified direction.
398  use geometry, only: no_i, no_j, no_k
399 
400 ! Global Parameters:
401 !---------------------------------------------------------------------//
402 ! Parameter constants
403  use param1, only: zero, one, undefined
404 
405 ! Use the error manager for posting error messages.
406 !---------------------------------------------------------------------//
407  use error_manager
408 
409  use toleranc
410 
411  implicit none
412 
413 ! Dummy Arguments:
414 !---------------------------------------------------------------------//
415  INTEGER, INTENT(IN) :: ICV
416 
417 ! Local Variables:
418 !---------------------------------------------------------------------//
419 ! Loop index
420  INTEGER :: N
421 ! Sum of mass fraction.
422  DOUBLE PRECISION :: SUM
423 ! Flag for regular IC regions
424  LOGICAL :: BASIC_IC
425 !......................................................................!
426 
427 ! Initialize the error manager.
428  CALL init_err_msg("CHECK_IC_GAS_PHASE")
429 
430 ! Patch ICs skip various checks.
431  basic_ic = (ic_type(icv) /= 'PATCH')
432 
433 ! Check that gas phase velocity components are initialized.
434  IF(basic_ic) THEN
435  IF(ic_u_g(icv) == undefined) THEN
436  IF(no_i) THEN
437  ic_u_g(icv) = zero
438  ELSE
439  WRITE(err_msg, 1000) trim(ivar('IC_U_g',icv))
440  CALL flush_err_msg(abort=.true.)
441  ENDIF
442  ENDIF
443 
444  IF(ic_v_g(icv) == undefined) THEN
445  IF(no_j) THEN
446  ic_v_g(icv) = zero
447  ELSE
448  WRITE(err_msg, 1000) trim(ivar('IC_V_g',icv))
449  CALL flush_err_msg(abort=.true.)
450  ENDIF
451  ENDIF
452 
453  IF(ic_w_g(icv) == undefined) THEN
454  IF (no_k) THEN
455  ic_w_g(icv) = zero
456  ELSE
457  WRITE(err_msg, 1000) trim(ivar('IC_W_g',icv))
458  CALL flush_err_msg(abort=.true.)
459  ENDIF
460  ENDIF
461  ENDIF
462 
463 ! Check that gas phase void fraction is initialized. Patched ICs may
464 ! have an undefined volume fration. A second check is preformed on
465 ! the solids.
466  IF(ic_ep_g(icv) == undefined .AND. basic_ic) THEN
467  WRITE(err_msg, 1000) trim(ivar('IC_EP_g',icv))
468  CALL flush_err_msg(abort=.true.)
469  ENDIF
470 
471 ! Check that if the gas phase pressure is initialized and the gas is
472 ! compressible that the gas phase pressure is not zero or negative
473  IF(ic_p_g(icv) /= undefined) THEN
474  IF(ro_g0==undefined .AND. ic_p_g(icv)<=zero) THEN
475  WRITE(err_msg, 1100) trim(ivar('IC_P_g',icv)), &
476  ival(ic_p_g(icv))
477  CALL flush_err_msg(abort=.true.)
478  ENDIF
479  ENDIF
480 
481  1100 FORMAT('Error 1100: Pressure must be greater than 0.0 for ', &
482  'compressible flow',/'Illegal value: ',a,' = ',a,/'Please ', &
483  'correct the mfix.dat file.')
484 
485  IF(basic_ic) THEN
486  IF(energy_eq .OR. ro_g0==undefined .OR. mu_g0==undefined) THEN
487  IF(ic_t_g(icv)==undefined) THEN
488  WRITE(err_msg, 1000) trim(ivar('IC_T_g',icv))
489  CALL flush_err_msg(abort=.true.)
490  ENDIF
491  ENDIF
492  ENDIF
493 
494 ! Gas phase radiation values.
495  IF (energy_eq) THEN
496  IF (ic_gama_rg(icv) < zero) THEN
497  WRITE(err_msg, 1001) trim(ivar('IC_GAMA_Rg',icv)), &
498  ival(ic_gama_rg(icv))
499  CALL flush_err_msg(abort=.true.)
500  ELSEIF (ic_gama_rg(icv) > zero) THEN
501  IF (ic_t_rg(icv) == undefined) THEN
502  WRITE(err_msg, 1000) ivar('IC_T_Rg',icv)
503  CALL flush_err_msg(abort=.true.)
504  ENDIF
505  ENDIF
506  ENDIF
507 
508 
509 ! First: sum together defiend gas phase species mass fractions.
510  sum = zero
511  DO n = 1, nmax(0)
512  IF(ic_x_g(icv,n) /= undefined) THEN
513  sum = sum + ic_x_g(icv,n)
514  ELSEIF(basic_ic) THEN
515  ic_x_g(icv,n) = zero
516  ENDIF
517  ENDDO
518 
519 
520 ! Enforce that the species mass fractions must sum to one.
521  IF(.NOT.compare(one,sum)) THEN
522 
523 ! No error for ZERO sum PATCH IC data.
524  IF(.NOT.basic_ic .AND. compare(zero,sum))THEN
525 
526 ! Error for regular IC regions when solving species equations.
527  ELSEIF(species_eq(0)) THEN
528  WRITE(err_msg, 1110) icv
529  CALL flush_err_msg(abort=.true.)
530 
531  1110 FORMAT('Error 1110: IC_X_g(',i3,',:) do NOT sum to ONE and the ',&
532  'gas phase',/'species equations are solved. Please correct ', &
533  'the mfix.dat file.')
534 
535 ! Error for slightly compressible flows with MW_AVG specified.
536  ELSEIF(ro_g0 == undefined .AND. mw_avg == undefined) THEN
537  WRITE(err_msg, 1111) icv
538  CALL flush_err_msg(abort=.true.)
539 
540  1111 FORMAT('Error 1111: IC_X_g(',i3,',:) do NOT sum to ONE and the ',&
541  'gas phase',/'is compressible and MW_AVG is UNDEFINED.',/ &
542  'Please correct the mfix.dat file.')
543 
544  ELSEIF(.NOT.compare(sum,zero)) THEN
545  WRITE(err_msg, 1112) icv
546  CALL flush_err_msg(abort=.true.)
547 
548  1112 FORMAT('Error 1112: IC_X_g(',i3,',:) do not sum to ONE or ZERO ',&
549  'and they',/'are not needed. Please correct the mfix.dat ', &
550  'file.')
551 
552  ELSE
553  ic_x_g(icv,:) = zero
554  ic_x_g(icv,1) = one
555  ENDIF
556  ENDIF
557 
558 
559  DO n = 1, nscalar
560  IF(ic_scalar(icv,n) == undefined) ic_scalar(icv,n) = zero
561  ENDDO
562 
563 
564  IF(basic_ic) THEN
565  IF (l_scale0 == zero) THEN
566  IF (ic_l_scale(icv) /= undefined) THEN
567  WRITE(err_msg, 1113) icv
568  CALL flush_err_msg(abort=.true.)
569  1113 FORMAT('Error 1113: IC_L_SCALE(',i3,',:) is defined but L_SCALE0 ',&
570  'is equal',/,'to zero. A non-zero value must be specified to ',&
571  'activate this model.',/,'Please correct the mfix.dat file.')
572  ENDIF
573  ELSEIF (l_scale0 < zero) THEN
574  WRITE(err_msg, 1001) 'L_SCALE0', ival(l_scale0)
575  CALL flush_err_msg(abort=.true.)
576  ELSE ! l_scale0 is defined at not zero
577  IF (ic_l_scale(icv) < zero) THEN
578  WRITE(err_msg, 1001) ivar('IC_L_SCALE',icv), &
579  ival(ic_l_scale(icv))
580  CALL flush_err_msg(abort=.true.)
581  ENDIF
582  ENDIF
583  ENDIF
584 
585 
586  IF(k_epsilon .AND. basic_ic) THEN
587  IF (ic_k_turb_g(icv) == undefined) THEN
588  WRITE(err_msg, 1000) ivar('IC_K_Turb_G',icv)
589  CALL flush_err_msg(abort=.true.)
590  ENDIF
591  IF (ic_e_turb_g(icv) == undefined) THEN
592  WRITE(err_msg, 1000) ivar('IC_E_Turb_G',icv)
593  CALL flush_err_msg(abort=.true.)
594  ENDIF
595  ENDIF
596 
597  CALL finl_err_msg
598 
599  RETURN
600 
601  1000 FORMAT('Error 1000: Required input not specified: ',a,/'Please ',&
602  'correct the mfix.dat file.')
603 
604  1001 FORMAT('Error 1001: Illegal or unknown input: ',a,' = ',a,/ &
605  'Please correct the mfix.dat file.')
606 
607  END SUBROUTINE check_ic_gas_phase
608 
609 
610 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
611 ! !
612 ! Subroutine: CHECK_IC_SOLIDS_PHASES !
613 ! Author: P. Nicoletti Date: 02-DEC-91 !
614 ! Author: J.Musser Date: 01-MAR-14 !
615 ! !
616 ! Purpose: Verify solids phase(s) input variables in IC region. !
617 ! !
618 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
619  SUBROUTINE check_ic_solids_phases(ICV)
621 
622 ! Global Variables:
623 !---------------------------------------------------------------------//
624 ! Solids volume fraction, bulk density
625  use ic, only: ic_ep_s, ic_rop_s
626 ! Solids velocity components.
627  use ic, only: ic_u_s, ic_v_s, ic_w_s
628 ! Solids temperature, mass fractions, granular energy
629  use ic, only: ic_t_s, ic_x_s, ic_theta_m
630 ! Radiation model parameters.
631  use ic, only: ic_gama_rs, ic_t_rs
632 ! Gas phase volume fraction and temperature.
633  use ic, only: ic_ep_g, ic_t_g
634 ! IC Type: UNDEFINED or PATCH.
635  use ic, only: ic_type
636 
637 ! Flag. Solve Energy equations
638  use run, only: energy_eq
639 ! Flag. Solve Species equations
640  use run, only: species_eq
641 ! Flag. Solve Granular Energy PDE
642  use run, only: granular_energy
643 ! Flag. Solve variable solids density
644  use run, only: solve_ros
645 ! Baseline solids mass fraction, index of intert
646  use physprop, only: x_s0, inert_species
647 ! Specified constant solids density.
648  use physprop, only: ro_s0
649 ! Number of gas phase species
650  use physprop, only: nmax
651 ! Number of TFM solids phases.
652  use physprop, only: smax
653 ! Number of DEM or PIC solids.
654  use discretelement, only: des_mmax
655 ! Flag: Do not solve in specified direction.
656  use geometry, only: no_i, no_j, no_k
657 
658 ! Global Parameters:
659 !---------------------------------------------------------------------//
660 ! Maximum number of solids species.
661  use param, only: dim_m
662 ! Parameter constants
663  use param1, only: zero, one, undefined
664 ! function to calculate solids density.
665  USE eos, ONLY: eoss
666 
667 ! Use the error manager for posting error messages.
668 !---------------------------------------------------------------------//
669  use error_manager
670 
671  use toleranc
672 
673  IMPLICIT NONE
674 
675 ! Dummy Arguments:
676 !---------------------------------------------------------------------//
677 ! Index of IC region.
678  INTEGER, INTENT(IN) :: ICV
679 ! Loop/variable index
680  INTEGER :: M, N
681 ! Various sums.
682  DOUBLE PRECISION SUM, SUM_EP
683 ! Solids phase density in IC region.
684  DOUBLE PRECISION :: IC_ROs(1:dim_m)
685 ! Index of inert species
686  INTEGER :: INERT
687 ! Flag to skip checks on indexed solid phase.
688  LOGICAL :: SKIP(1:dim_m)
689 ! Total number of solids phases (TFM + DEM + MPPIC)
690  INTEGER :: MMAX_TOT
691 ! Flag for PATCH IC regions
692  LOGICAL :: BASIC_IC
693 !......................................................................!
694 
695 ! Initialize the error manager.
696  CALL init_err_msg("CHECK_IC_SOLIDS_PHASES")
697 
698 ! Patch ICs skip various checks.
699  basic_ic = (ic_type(icv) /= 'PATCH')
700 
701 ! The total number of solids phases (all models).
702  mmax_tot = smax + des_mmax
703 
704 ! Calculate EP_s from EP_g if there is only one solids phase.
705  IF(mmax_tot == 1 .AND. ic_ep_s(icv,1) == undefined) THEN
706  IF(ic_ep_g(icv) /= undefined) ic_ep_s(icv,1) = one-ic_ep_g(icv)
707  ENDIF
708 
709 ! Bulk density or solids volume fraction must be explicitly defined
710 ! if there are more than one solids phase.
711  IF(mmax_tot > 1 .AND. .NOT.compare(ic_ep_g(icv),one)) THEN
712 ! IC_EP_g may be undefined for PATCH IC regions.
713  IF(ic_ep_g(icv) /= undefined) THEN
714  DO m = 1, mmax_tot
715  IF(ic_rop_s(icv,m) == undefined .AND. &
716  ic_ep_s(icv,m) == undefined) THEN
717  WRITE(err_msg, 1400) m, icv, 'IC_ROP_s and IC_EP_s'
718  CALL flush_err_msg(abort=.true.)
719  ENDIF
720  ENDDO
721 
722 ! If IC_EP_G is undefined, then ROP_s and EP_s should be too.
723  ELSE
724  DO m = 1, mmax_tot
725  IF(ic_rop_s(icv,m) /= undefined .AND. &
726  ic_ep_s(icv,m) /= undefined) THEN
727  WRITE(err_msg, 1400) m, icv, 'IC_ROP_s and IC_EP_s'
728  CALL flush_err_msg(abort=.true.)
729  ENDIF
730  ENDDO
731  ENDIF
732  ENDIF
733 
734  1400 FORMAT('Error 1400: Insufficient solids phase ',i2,' ', &
735  'information for IC',/'region ',i3,'. ',a,' not specified.',/ &
736  'Please correct the mfix.dat file.')
737 
738 ! Determine which solids phases are present.
739  DO m = 1, mmax_tot
740  skip(m)=(ic_rop_s(icv,m)==undefined.OR.ic_rop_s(icv,m)==zero) &
741  .AND.(ic_ep_s(icv,m)==undefined .OR.ic_ep_s(icv,m)==zero)
742  ENDDO
743 
744  IF(mmax_tot == 1 .AND. ic_ep_g(icv)/=one) skip(1) = .false.
745 
746  DO m=1, mmax_tot
747 
748 ! check that solids phase m velocity components are initialized
749  IF(basic_ic) THEN
750  IF(ic_u_s(icv,m) == undefined) THEN
751  IF (skip(m) .OR. no_i) THEN
752  ic_u_s(icv,m) = zero
753  ELSE
754  WRITE(err_msg, 1000)trim(ivar('IC_U_s',icv,m))
755  CALL flush_err_msg(abort=.true.)
756  ENDIF
757  ENDIF
758 
759  IF(ic_v_s(icv,m) == undefined) THEN
760  IF(skip(m) .OR. no_j) THEN
761  ic_v_s(icv,m) = zero
762  ELSE
763  WRITE(err_msg, 1000)trim(ivar('IC_V_s',icv,m))
764  CALL flush_err_msg(abort=.true.)
765  ENDIF
766  ENDIF
767 
768  IF(ic_w_s(icv,m) == undefined) THEN
769  IF(skip(m) .OR. no_k) THEN
770  ic_w_s(icv,m) = zero
771  ELSE
772  WRITE(err_msg, 1000)trim(ivar('IC_W_s',icv,m))
773  CALL flush_err_msg(abort=.true.)
774  ENDIF
775  ENDIF
776 
777  IF(energy_eq .AND. ic_t_s(icv,m)==undefined) THEN
778  IF(skip(m)) THEN
779  ic_t_s(icv,m) = ic_t_g(icv)
780  ELSE
781  WRITE(err_msg, 1000)trim(ivar('IC_T_s',icv,m))
782  CALL flush_err_msg(abort=.true.)
783  ENDIF
784  ENDIF
785 
786  IF(granular_energy .AND. ic_theta_m(icv,m)==undefined) THEN
787  IF(skip(m)) THEN
788  ic_theta_m(icv,m) = zero
789  ELSE
790  WRITE(err_msg, 1000)trim(ivar('IC_Theta_M',icv,m))
791  CALL flush_err_msg(abort=.true.)
792  ENDIF
793  ENDIF
794 
795  IF(energy_eq) THEN
796  IF(ic_gama_rs(icv,m) < zero) THEN
797  WRITE(err_msg, 1001)trim(ivar('IC_GAMA_Rs',icv,m)), &
798  ival(ic_gama_rs(icv,m))
799  CALL flush_err_msg(abort=.true.)
800  ELSEIF (ic_gama_rs(icv,m) > zero) THEN
801  IF(ic_t_rs(icv,m) == undefined) THEN
802  WRITE(err_msg, 1001)trim(ivar('IC_T_Rs',icv,m))
803  CALL flush_err_msg(abort=.true.)
804  ENDIF
805  ENDIF
806  ENDIF
807  ENDIF
808 
809 
810 ! First: sum together defiend species mass fractions.
811  sum = zero
812  DO n = 1, nmax(m)
813  IF(ic_x_s(icv,m,n) /= undefined) THEN
814  sum = sum + ic_x_s(icv,m,n)
815  ELSEIF(basic_ic) THEN
816  ic_x_s(icv,m,n) = zero
817  ENDIF
818  ENDDO
819 
820 ! Enforce that the species mass fractions must sum to one.
821  IF(.NOT.compare(one,sum)) THEN
822 
823 ! Summing to zero is not an error for PATCH IC regions.
824  IF(.NOT.basic_ic .AND. compare(zero,sum)) THEN
825 
826  ELSEIF(species_eq(m) .AND. .NOT.skip(m)) THEN
827  WRITE(err_msg, 1402) icv, m
828  CALL flush_err_msg(abort=.true.)
829 
830  1402 FORMAT('Error 1402: IC_X_s(',i3,',',i2,',:) do NOT sum to ONE ', &
831  'and the solids phase',/'species equations are solved. ', &
832  'Please correct the mfix.dat file.')
833 
834  ELSEIF(solve_ros(m) .AND. .NOT.skip(m)) THEN
835  WRITE(err_msg, 1403) icv, m
836  CALL flush_err_msg(abort=.true.)
837 
838  1403 FORMAT('Error 1403: IC_X_s(',i3,',',i2,':) do NOT sum to ONE ', &
839  'and the solids phase',/'density is calculated. Please ', &
840  'correct the mfix.dat file.')
841 
842  ELSEIF(.NOT.compare(sum,zero)) THEN
843  WRITE(err_msg, 1404) icv
844  CALL flush_err_msg(abort=.true.)
845 
846  1404 FORMAT('Error 1404: IC_X_s(',i3,',',i2,':) do NOT sum to ONE ', &
847  'or ZERO and',/'they are not needed. Please correct the ', &
848  'mfix.dat file.')
849 
850  ELSE
851  ic_x_s(icv,m,:) = zero
852  ic_x_s(icv,m,1) = one
853  ENDIF
854  ENDIF
855 
856 ! Set the solids density for the IC region.
857  IF(skip(m)) THEN
858  ic_ros(m) = merge(ro_s0(m), ro_s0(m), solve_ros(m))
859 
860  ELSEIF(.NOT.solve_ros(m)) THEN
861  ic_ros(m) = ro_s0(m)
862 
863  ELSE
864 ! Verify that the species mass fraction for the inert material is not
865 ! zero in the IC region when the solids is present.
866  inert = inert_species(m)
867  IF(ic_x_s(icv,m,inert) == zero) THEN
868  WRITE(err_msg,1405) m, icv
869  CALL flush_err_msg(abort=.true.)
870  ENDIF
871 
872  1405 FORMAT('Error 1405: No inert species for phase ',i2,' in IC ', &
873  'region ',i3,'.',/'Unable to calculate solids phase density. ',&
874  'Please refer to the Readme',/' file for required variable ', &
875  'solids density model input parameters and',/' make the ', &
876  'necessary corrections to the data file.')
877 
878 ! Calculate the solids density.
879  ic_ros(m) = eoss(ro_s0(m), x_s0(m,inert), &
880  ic_x_s(icv,m,inert))
881  ENDIF
882 
883  ENDDO ! end loop over (m=1,smax)
884 
885 
886 ! Initialize the sum of the total volume fraction.
887  sum_ep = ic_ep_g(icv)
888 
889  DO m=1, mmax_tot
890 
891 ! Clear out both varaibles if this phase is skipped.
892  IF(basic_ic .AND. skip(m)) THEN
893  ic_ep_s(icv,m) = zero
894  ic_rop_s(icv,m) = zero
895 
896 ! Leave everything undefined for PATCH ICs that are not specifed.
897  ELSEIF(.NOT.basic_ic .AND. (ic_rop_s(icv,m) == undefined &
898  .AND. ic_ep_s(icv,m) == undefined)) THEN
899 
900 ! If both input parameters are defined. Make sure they are equivalent.
901  ELSEIF(ic_rop_s(icv,m) /= undefined .AND. &
902  ic_ep_s(icv,m) /= undefined) THEN
903 
904  IF(.NOT.compare(ic_ep_s(icv,m)*ic_ros(m), &
905  ic_rop_s(icv,m))) THEN
906 
907 ! BASIC_IC regions require that the IC_ROP_s and IC_EP_s specifications
908 ! match although it is unlikely that anyone would specify both.
909  IF(basic_ic) THEN
910 
911  WRITE(err_msg,1406) m, icv
912  CALL flush_err_msg(abort=.true.)
913 
914 1406 FORMAT('Error 1406: IC_EP_s and IC_ROP_s are inconsistent for ',&
915  'phase ',i2,/,'in IC region ', i3,'. Please correct the ',&
916  'mfix.dat file.')
917 
918 
919 ! PachedeIC regions defer to IC_EP_s if the values do not match. This
920 ! prevents a dead lock or the need to define both. This case is rather
921 ! common as a defined IC_EP_s is converted to IC_ROP_s. Therefore, if
922 ! a patch region is used more than once, these values may not match.
923  ELSE
924 
925  WRITE(err_msg,1407) trim(ivar('IC_ROP_s',icv,m)), &
926  trim(ival(ic_rop_s(icv,m))), trim(ivar('IC_EP_s',&
927  icv,m)), trim(ival(ic_ep_s(icv,m)))
928  CALL flush_err_msg()
929 
930 1407 FORMAT('Warning 1407: IC_EP_s and IC_ROP_s are inconsistent:', &
931  2(/3x,a,' = ',a),/'Deferring to IC_EP_s to overcome conflict.')
932 
933  ic_rop_s(icv,m) = ic_ep_s(icv,m)*ic_ros(m)
934 
935  ENDIF
936  ENDIF
937 
938 
939 ! Compute IC_EP_s from IC_ROP_s
940  ELSEIF(ic_ep_s(icv,m) == undefined)THEN
941  ic_ep_s(icv,m) = ic_rop_s(icv,m) / ic_ros(m)
942 
943 ! Compute IC_ROP_s from IC_EP_s and IC_ROs
944  ELSEIF(ic_rop_s(icv,m) == undefined) THEN
945  ic_rop_s(icv,m) = ic_ep_s(icv,m) * ic_ros(m)
946 ! This is a sanity check.
947  ELSE
948 
949  ENDIF
950 ! Add this phase to the total volume fraction.
951  sum_ep = sum_ep + ic_ep_s(icv,m)
952  ENDDO
953 
954 ! Verify that the volume fractions sum to one.
955  IF(basic_ic .AND. .NOT.compare(sum_ep,one)) THEN
956  WRITE(err_msg,1410) icv
957  CALL flush_err_msg(abort=.true.)
958  ENDIF
959 
960  1410 FORMAT('Error 1410: Illegal initial condition region : ',i3,/ &
961  'Sum of volume fractions does NOT equal ONE. Please correct',/&
962  'the mfix.dat file.')
963 
964 
965  CALL finl_err_msg
966 
967  RETURN
968 
969 
970  1000 FORMAT('Error 1000: Required input not specified: ',a,/'Please ',&
971  'correct the mfix.dat file.')
972 
973  1001 FORMAT('Error 1001: Illegal or unknown input: ',a,' = ',a,/ &
974  'Please correct the mfix.dat file.')
975 
976  1002 FORMAT('Error 1002: Illegal data in initial condition region ', &
977  i3,/a,' defined for unknown solids phase ',i2,'.',/ &
978  'Please correct the mfix.dat file.')
979 
980  END SUBROUTINE check_ic_solids_phases
981 
982 
983 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
984 ! !
985 ! Subroutine: CHECK_IC_OVERFLOW !
986 ! Author: J.Musser Date: 01-Mar-14 !
987 ! !
988 ! Purpose: Verify that no data was defined for unspecified IC. !
989 ! !
990 ! !
991 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
992 
993  SUBROUTINE check_ic_overflow(ICV)
995 ! Global Variables:
996 !---------------------------------------------------------------------//
997 ! Gas phase volume fraction, pressure, temperature, species.
998  use ic, only: ic_ep_g, ic_t_g, ic_x_g
999 ! Gas phase velocity components.
1000  use ic, only: ic_u_g, ic_v_g, ic_w_g
1001 ! Radiation model parameters.
1002  use ic, only: ic_t_rg
1003 ! K-Epsilon model parameters.
1004  use ic, only: ic_k_turb_g, ic_e_turb_g
1005 ! IC for user-defined scalar equation.
1006  use ic, only: ic_scalar
1007 ! Solids volume fraction, bulk density
1008  use ic, only: ic_rop_s
1009 ! Solids velocity components.
1010  use ic, only: ic_u_s, ic_v_s, ic_w_s
1011 ! Solids temperature, mass fractions, granular energy
1012  use ic, only: ic_t_s, ic_x_s
1013 ! Radiation model parameters.
1014  use ic, only: ic_t_rs
1015 ! IC Type: UNDEFINED or PATCH.
1016  use ic, only: ic_type
1017 
1018 ! Global Parameters:
1019 !---------------------------------------------------------------------//
1020 ! Parameter constant
1021  use param1, only: undefined
1022 ! Maximum number of solids phases and user-defined scalars.
1023  use param, only: dim_m, dim_scalar
1024 ! Maximum number of gas and solids phase species
1025  use param, only: dim_n_g, dim_n_s
1026 
1027 ! Use the error manager for posting error messages.
1028 !---------------------------------------------------------------------//
1029  use error_manager
1030 
1031  implicit none
1032 
1033 ! Dummy Arguments:
1034 !---------------------------------------------------------------------//
1035 ! IC region index.
1036  INTEGER, INTENT(IN) :: ICV
1037 
1038 ! Local Variables:
1039 !---------------------------------------------------------------------//
1040 ! Loop counters
1041  INTEGER :: M, N
1042 !......................................................................!
1043 
1044  IF (ic_type(icv) == 'PATCH') RETURN
1045 
1046 ! Initialize the error manager.
1047  CALL init_err_msg("CHECK_IC_OVERFLOW")
1048 
1049 ! GAS PHASE quantities
1050 ! -------------------------------------------->>>
1051  IF(ic_u_g(icv) /= undefined) THEN
1052  WRITE(err_msg, 1010) trim(ivar('IC_U_g',icv))
1053  CALL flush_err_msg(abort=.true.)
1054  ELSEIF(ic_v_g(icv) /= undefined) THEN
1055  WRITE(err_msg, 1010) trim(ivar('IC_V_g',icv))
1056  CALL flush_err_msg(abort=.true.)
1057  ELSEIF(ic_w_g(icv) /= undefined) THEN
1058  WRITE(err_msg, 1010) trim(ivar('IC_W_g',icv))
1059  CALL flush_err_msg(abort=.true.)
1060  ELSEIF(ic_ep_g(icv) /= undefined) THEN
1061  WRITE(err_msg, 1010) trim(ivar('IC_EP_g',icv))
1062  CALL flush_err_msg(abort=.true.)
1063  ELSEIF(ic_t_g(icv) /= undefined) THEN
1064  WRITE(err_msg, 1010) trim(ivar('IC_T_g',icv))
1065  CALL flush_err_msg(abort=.true.)
1066  ELSEIF(ic_t_rg(icv) /= undefined) THEN
1067  WRITE(err_msg, 1010) trim(ivar('IC_T_Rg',icv))
1068  CALL flush_err_msg(abort=.true.)
1069  ELSEIF(ic_k_turb_g(icv) /= undefined) THEN
1070  WRITE(err_msg, 1010) trim(ivar('IC_K_Turb_G',icv))
1071  CALL flush_err_msg(abort=.true.)
1072  ELSEIF(ic_e_turb_g(icv) /= undefined) THEN
1073  WRITE(err_msg, 1010) trim(ivar('IC_E_Turb_G',icv))
1074  CALL flush_err_msg(abort=.true.)
1075  ENDIF
1076  DO n = 1, dim_n_g
1077  IF(ic_x_g(icv,n) /= undefined) THEN
1078  WRITE(err_msg, 1010) trim(ivar('IC_X_g',icv))
1079  CALL flush_err_msg(abort=.true.)
1080  ENDIF
1081  ENDDO
1082  DO n = 1, dim_scalar
1083  IF(ic_scalar(icv,n) /= undefined) THEN
1084  WRITE(err_msg, 1010) trim(ivar('IC_Scalar',icv))
1085  CALL flush_err_msg(abort=.true.)
1086  ENDIF
1087  ENDDO
1088 ! --------------------------------------------<<<
1089 
1090 
1091 ! SOLIDS PHASE quantities
1092 ! -------------------------------------------->>>
1093  DO m=1, dim_m
1094  IF(ic_rop_s(icv,m) /= undefined) THEN
1095  WRITE(err_msg, 1010) trim(ivar('IC_ROP_s',icv,m))
1096  CALL flush_err_msg(abort=.true.)
1097  ELSEIF(ic_u_s(icv,m) /= undefined) THEN
1098  WRITE(err_msg, 1010) trim(ivar('IC_U_s',icv,m))
1099  CALL flush_err_msg(abort=.true.)
1100  ELSEIF(ic_v_s(icv,m) /= undefined) THEN
1101  WRITE(err_msg, 1010) trim(ivar('IC_V_s',icv,m))
1102  CALL flush_err_msg(abort=.true.)
1103  ELSEIF(ic_w_s(icv,m) /= undefined) THEN
1104  WRITE(err_msg, 1010) trim(ivar('IC_W_s',icv,m))
1105  CALL flush_err_msg(abort=.true.)
1106  ELSEIF(ic_t_s(icv,m) /= undefined) THEN
1107  WRITE(err_msg, 1010) trim(ivar('IC_T_s',icv,m))
1108  CALL flush_err_msg(abort=.true.)
1109  ELSEIF(ic_t_rs(icv,m) /= undefined) THEN
1110  WRITE(err_msg, 1010) trim(ivar('IC_T_Rs',icv,m))
1111  CALL flush_err_msg(abort=.true.)
1112  ENDIF
1113  DO n = 1, dim_n_s
1114  IF(ic_x_s(icv,m,n) /= undefined) THEN
1115  WRITE(err_msg, 1010) trim(ivar('IC_X_s',icv,m))
1116  CALL flush_err_msg(abort=.true.)
1117  ENDIF
1118  ENDDO
1119  ENDDO
1120 ! --------------------------------------------<<<
1121 
1122  CALL finl_err_msg
1123  RETURN
1124 
1125  1010 FORMAT('Error 1010: ',a,' specified in an undefined IC region')
1126 
1127  END SUBROUTINE check_ic_overflow
double precision l_scale0
Definition: constant_mod.f:177
double precision, dimension(dimension_ic) ic_e_turb_g
Definition: ic_mod.f:132
integer, parameter dim_n_g
Definition: param_mod.f:69
integer, parameter dimension_ic
Definition: param_mod.f:59
double precision, dimension(dimension_ic) ic_l_scale
Definition: ic_mod.f:71
subroutine calc_cell(RMIN, REACTOR_LOC, D_DIR, N_DIR, CELL_LOC)
Definition: calc_cell.f:14
logical dem_solids
Definition: run_mod.f:257
double precision, dimension(dimension_ic, dim_m) ic_rop_s
Definition: ic_mod.f:74
character(len=32) function ivar(VAR, i1, i2, i3)
double precision, dimension(dimension_ic) ic_t_g
Definition: ic_mod.f:80
subroutine finl_err_msg
subroutine check_ic_overflow(ICV)
logical function compare(V1, V2)
Definition: toleranc_mod.f:94
integer, dimension(dimension_ic) ic_j_s
Definition: ic_mod.f:47
logical no_i
Definition: geometry_mod.f:20
double precision, dimension(dimension_ic, dim_m) ic_theta_m
Definition: ic_mod.f:86
subroutine check_ic_mppic
double precision, parameter one
Definition: param1_mod.f:29
subroutine check_ic_common_discrete
double precision, dimension(dimension_ic, dim_scalar) ic_scalar
Definition: ic_mod.f:128
logical, dimension(0:dim_m) species_eq
Definition: run_mod.f:115
integer, dimension(dimension_ic) ic_j_n
Definition: ic_mod.f:50
integer, parameter dim_scalar
Definition: param_mod.f:85
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(0:dim_j) dy
Definition: geometry_mod.f:70
logical, dimension(dimension_ic) ic_defined
Definition: ic_mod.f:107
double precision, dimension(dimension_ic) ic_z_b
Definition: ic_mod.f:35
double precision, dimension(dim_m, dim_n_s) x_s0
Definition: physprop_mod.f:32
double precision, dimension(dimension_ic) ic_x_w
Definition: ic_mod.f:23
double precision, parameter undefined
Definition: param1_mod.f:18
character(len=16), dimension(dimension_ic) ic_type
Definition: ic_mod.f:59
double precision, dimension(dimension_ic, dim_m) ic_gama_rs
Definition: ic_mod.f:122
double precision, dimension(0:dim_k) dz
Definition: geometry_mod.f:72
subroutine check_ic_dem
Definition: check_ic_dem.f:12
subroutine check_ic_geometry
integer imax
Definition: geometry_mod.f:47
double precision, dimension(dimension_ic) ic_u_g
Definition: ic_mod.f:89
subroutine init_err_msg(CALLER)
Definition: ic_mod.f:9
double precision, dimension(dimension_ic) ic_z_t
Definition: ic_mod.f:38
double precision, dimension(dimension_ic) ic_k_turb_g
Definition: ic_mod.f:131
double precision, dimension(dimension_ic, dim_m) ic_w_s
Definition: ic_mod.f:104
integer kmax1
Definition: geometry_mod.f:58
integer, dimension(dimension_ic) ic_i_w
Definition: ic_mod.f:41
double precision ro_g0
Definition: physprop_mod.f:59
integer imax1
Definition: geometry_mod.f:54
double precision function eoss(pBase, Xs0_INERT, Xs_INERT)
Definition: eos_mod.f:155
integer, dimension(dimension_ic) ic_i_e
Definition: ic_mod.f:44
double precision, dimension(dimension_ic) ic_v_g
Definition: ic_mod.f:95
integer, dimension(dimension_ic) ic_k_b
Definition: ic_mod.f:53
integer, dimension(dim_m) inert_species
Definition: physprop_mod.f:39
Definition: eos_mod.f:10
character(len=16) run_type
Definition: run_mod.f:33
double precision, dimension(dimension_ic) ic_y_n
Definition: ic_mod.f:32
double precision, dimension(dimension_ic, dim_m, dim_n_s) ic_x_s
Definition: ic_mod.f:113
double precision, dimension(dimension_ic) ic_w_g
Definition: ic_mod.f:101
double precision, dimension(0:dim_i) dx
Definition: geometry_mod.f:68
double precision, dimension(dimension_ic) ic_gama_rg
Definition: ic_mod.f:116
double precision xlength
Definition: geometry_mod.f:33
integer jmax1
Definition: geometry_mod.f:56
Definition: run_mod.f:13
integer, dimension(dimension_ic) ic_k_t
Definition: ic_mod.f:56
Definition: param_mod.f:2
integer kmax
Definition: geometry_mod.f:51
double precision, dimension(dimension_ic, dim_m) ic_v_s
Definition: ic_mod.f:98
subroutine location_check(CELL_SPECIFIED, CELL_CALCULATED, COUNTER, MESSAGE)
logical no_k
Definition: geometry_mod.f:28
subroutine check_initial_conditions
integer, dimension(0:dim_m) nmax
Definition: physprop_mod.f:119
integer jmin1
Definition: geometry_mod.f:42
double precision, dimension(dimension_ic) ic_p_g
Definition: ic_mod.f:65
logical k_epsilon
Definition: run_mod.f:97
logical reinitializing
Definition: run_mod.f:208
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_ic) ic_x_e
Definition: ic_mod.f:26
subroutine check_ic_gas_phase(ICV)
integer, parameter undefined_i
Definition: param1_mod.f:19
character(len=line_length), dimension(line_count) err_msg
integer, parameter dim_n_s
Definition: param_mod.f:71
integer nscalar
Definition: scalars_mod.f:7
double precision, dimension(dimension_ic, dim_m) ic_u_s
Definition: ic_mod.f:92
double precision xmin
Definition: geometry_mod.f:75
subroutine check_ic_solids_phases(ICV)
integer jmax
Definition: geometry_mod.f:49
double precision, dimension(dimension_ic, dim_n_g) ic_x_g
Definition: ic_mod.f:110
double precision, dimension(dimension_ic) ic_ep_g
Definition: ic_mod.f:62
integer smax
Definition: physprop_mod.f:22
double precision, dimension(dim_m) ro_s0
Definition: physprop_mod.f:28
double precision ylength
Definition: geometry_mod.f:35
integer imin1
Definition: geometry_mod.f:40
logical granular_energy
Definition: run_mod.f:112
double precision, dimension(dimension_ic) ic_y_s
Definition: ic_mod.f:29
logical pic_solids
Definition: run_mod.f:258
double precision, dimension(dimension_ic) ic_t_rg
Definition: ic_mod.f:119
integer kmin1
Definition: geometry_mod.f:44
double precision, dimension(dimension_ic, dim_m) ic_t_s
Definition: ic_mod.f:83
double precision, dimension(dimension_ic, dim_m) ic_ep_s
Definition: ic_mod.f:77
double precision, parameter zero
Definition: param1_mod.f:27
double precision zlength
Definition: geometry_mod.f:37
subroutine flush_err_msg(DEBUG, HEADER, FOOTER, ABORT, LOG, CALL_TREE)
double precision, dimension(dimension_ic, dim_m) ic_t_rs
Definition: ic_mod.f:125