MFIX  2016-1
check_point_sources.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! !
3 ! SUBROUTINE: CHECK_POINT_SOURCES !
4 ! Author: J. Musser Date: 10-JUN-13 !
5 ! !
6 ! Purpose: Check point source specifications. !
7 ! !
8 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
9  SUBROUTINE check_point_sources
10 
11 ! Global Variables:
12 !---------------------------------------------------------------------//
13 ! Flag: PS geometry was detected.
14  use ps, only: ps_defined
15 
16 ! Global Parameters:
17 !---------------------------------------------------------------------//
18 ! Maximum number of PS.
19  use param, only: dimension_ps
20 
21 ! Use the error manager for posting error messages.
22 !---------------------------------------------------------------------//
23  use error_manager
24 
25  implicit none
26 
27 ! Local Variables:
28 !---------------------------------------------------------------------//
29 ! Loop counter for BCs
30  INTEGER :: PSV
31 !......................................................................!
32 
33 
34 ! Initialize the error manager.
35  CALL init_err_msg("CHECK_POINT_SOURCES")
36 
37 ! Determine which PSs are DEFINED
39 
40 ! Loop over all PS arrays.
41  DO psv = 1, dimension_ps
42 
43 ! Verify user input for defined defined PS.
44  IF(ps_defined(psv)) THEN
45  CALL get_ps(psv)
46  CALL check_ps_gas_phase(psv)
47  CALL check_ps_solids_phases(psv)
48  ELSE
49 ! Verify that no data was defined for unspecified PS.
50  CALL check_ps_overflow(psv)
51  ENDIF
52  ENDDO
53 
54 ! Clear the error manager.
55  CALL finl_err_msg
56 
57  RETURN
58  END SUBROUTINE check_point_sources
59 
60 
61 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
62 ! !
63 ! SUBROUTINE: CHECK_PS_GEOMETRY !
64 ! Author: J. Musser Date: 10-JUN-13 !
65 ! !
66 ! Purpose: Check point source specifications. !
67 ! !
68 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
69  SUBROUTINE check_ps_geometry
70 
71 
72 ! Global Variables:
73 !---------------------------------------------------------------------//
74 ! Flag: PS contains geometric data and/or specified type
75  use ps, only: ps_defined, point_source
76 ! User specified: PS geometry
77  use ps, only: ps_x_e, ps_x_w, ps_i_e, ps_i_w
78  use ps, only: ps_y_n, ps_y_s, ps_j_n, ps_j_s
79  use ps, only: ps_z_t, ps_z_b, ps_k_t, ps_k_b
80 ! User specified: System geometry
81  use geometry, only: no_i, xlength
82  use geometry, only: no_j, ylength
83  use geometry, only: no_k, zlength
84 
85 ! Global Parameters:
86 !---------------------------------------------------------------------//
87 ! The max number of BCs.
88  use param, only: dimension_ps
89 ! Parameter constants
90  use param1, only: zero, undefined, undefined_i
91 
92 ! Use the error manager for posting error messages.
93 !---------------------------------------------------------------------//
94  use error_manager
95 
96 
97  implicit none
98 
99 
100 ! Local Variables:
101 !---------------------------------------------------------------------//
102 ! PS loop counter.
103  INTEGER :: PSV
104 !......................................................................!
105 
106 ! Initialize the error manager.
107  CALL init_err_msg("CHECK_PS_GEOMETRY")
108 
109 ! Initialize the PS runtime flag.
110  point_source = .false.
111 
112 ! Determine which point source indices have values.
113  psv_lp: do psv = 1, dimension_ps
114 
115  IF (ps_x_w(psv) /= undefined) ps_defined(psv) = .true.
116  IF (ps_x_e(psv) /= undefined) ps_defined(psv) = .true.
117  IF (ps_y_s(psv) /= undefined) ps_defined(psv) = .true.
118  IF (ps_y_n(psv) /= undefined) ps_defined(psv) = .true.
119  IF (ps_z_b(psv) /= undefined) ps_defined(psv) = .true.
120  IF (ps_z_t(psv) /= undefined) ps_defined(psv) = .true.
121  IF (ps_i_w(psv) /= undefined_i) ps_defined(psv) = .true.
122  IF (ps_i_e(psv) /= undefined_i) ps_defined(psv) = .true.
123  IF (ps_j_s(psv) /= undefined_i) ps_defined(psv) = .true.
124  IF (ps_j_n(psv) /= undefined_i) ps_defined(psv) = .true.
125  IF (ps_k_b(psv) /= undefined_i) ps_defined(psv) = .true.
126  IF (ps_k_t(psv) /= undefined_i) ps_defined(psv) = .true.
127 
128 ! Skip consistency checks if nothing was defined.
129  IF (.NOT.ps_defined(psv)) cycle psv_lp
130 
131 ! Flag that one or more point sources has been detected.
132  point_source = .true.
133 
134  IF(ps_x_w(psv)==undefined .AND. ps_i_w(psv)==undefined_i) THEN
135  IF(no_i) THEN
136  ps_x_w(psv) = zero
137  ELSE
138  WRITE(err_msg,1101) psv, 'PS_X_w and PS_I_w '
139  CALL flush_err_msg(abort=.true.)
140  ENDIF
141  ENDIF
142  IF(ps_x_e(psv)==undefined .AND. ps_i_e(psv)==undefined_i) THEN
143  IF(no_i) THEN
144  ps_x_e(psv) = xlength
145  ELSE
146  WRITE(err_msg,1101) psv, 'PS_X_e and PS_I_e '
147  CALL flush_err_msg(abort=.true.)
148  ENDIF
149  ENDIF
150  IF(ps_y_s(psv)==undefined .AND. ps_j_s(psv)==undefined_i) THEN
151  IF(no_j) THEN
152  ps_y_s(psv) = zero
153  ELSE
154  WRITE(err_msg,1101) psv, 'PS_Y_s and PS_J_s '
155  CALL flush_err_msg(abort=.true.)
156  ENDIF
157  ENDIF
158  IF(ps_y_n(psv)==undefined .AND. ps_j_n(psv)==undefined_i) THEN
159  IF(no_j) THEN
160  ps_y_n(psv) = ylength
161  ELSE
162  WRITE(err_msg,1101) psv, 'PS_Y_n and PS_J_n '
163  CALL flush_err_msg(abort=.true.)
164  ENDIF
165  ENDIF
166  IF(ps_z_b(psv)==undefined .AND. ps_k_b(psv)==undefined_i) THEN
167  IF(no_k) THEN
168  ps_z_b(psv) = zero
169  ELSE
170  WRITE(err_msg,1101) psv, 'PS_Z_b and PS_K_b '
171  CALL flush_err_msg(abort=.true.)
172  ENDIF
173  ENDIF
174  IF(ps_z_t(psv)==undefined .AND. ps_k_t(psv)==undefined_i) THEN
175  IF(no_k) THEN
176  ps_z_t(psv) = zlength
177  ELSE
178  WRITE(err_msg,1101) psv, 'PS_Z_t and PS_K_t '
179  CALL flush_err_msg(abort=.true.)
180  ENDIF
181  ENDIF
182 
183  1101 FORMAT('Error 1101: Point source ',i3,' is ill-defined.',/a, &
184  ' are not specified.',/'Please correct the mfix.dat file.')
185 
186  ENDDO psv_lp
187 
188  CALL finl_err_msg
189 
190  RETURN
191  END SUBROUTINE check_ps_geometry
192 
193 
194 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
195 ! !
196 ! SUBROUTINE: CHECK_PS_GAS_PHASE !
197 ! Author: J. Musser Date: 10-JUN-13 !
198 ! !
199 ! Purpose: Check point source specifications. !
200 ! !
201 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
202  SUBROUTINE check_ps_gas_phase(PSV)
204 ! Global Variables:
205 !---------------------------------------------------------------------//
206 ! Gas phase mass flowrate for PS
207  use ps, only: ps_massflow_g
208 ! Gas phase velocity for PS
209  use ps, only: ps_u_g, ps_v_g, ps_w_g
210 ! Gas phase tempture and species mass fractions
211  use ps, only: ps_t_g, ps_x_g
212 ! Flag: Solve energy equations.
213  use run, only: energy_eq
214 ! Flag: Solve species equations.
215  use run, only: species_eq
216 ! Number of species.
217  use physprop, only: nmax
218 
219 ! Global Parameters:
220 !---------------------------------------------------------------------//
221 ! Parameter constants
222  use param1, only: zero, one, undefined
223 
224 ! Use the error manager for posting error messages.
225 !---------------------------------------------------------------------//
226  use error_manager
227 
228  use toleranc
229 
230  implicit none
231 
232 ! Dummy Arguments:
233 !---------------------------------------------------------------------//
234  INTEGER, INTENT(in) :: PSV
235 
236 ! Local Variables:
237 !---------------------------------------------------------------------//
238 ! Loop counter
239  INTEGER :: N
240 ! Sum of solids mass fractions.
241  DOUBLE PRECISION :: SUM
242 !......................................................................!
243 
244 ! Initialze the error manager.
245  CALL init_err_msg("CHECK_PS_GAS_PHASE")
246 
247 
248 ! Check mass flow and velocity
249  IF(ps_massflow_g(psv) == undefined) THEN
250  IF(ps_u_g(psv) /= undefined .OR. &
251  ps_v_g(psv) /= undefined .OR. &
252  ps_w_g(psv) /= undefined) THEN
253 
254  WRITE(err_msg,1100) psv, trim(ivar('PS_MASSFLOW_G',psv))
255  CALL flush_err_msg(abort=.true.)
256 
257  1100 FORMAT('Error 1100: Invalid specification for point source ',i3,&
258  '.',/a,' is undefined but velocity is given.',/'Please ', &
259  'correct the mfix.dat file.')
260 
261  ELSE
262  ps_massflow_g(psv) = zero
263  ps_u_g(psv) = zero
264  ps_v_g(psv) = zero
265  ps_w_g(psv) = zero
266  ENDIF
267 
268  ELSEIF(ps_massflow_g(psv) == zero) THEN
269  IF(ps_u_g(psv) /= zero .OR. &
270  ps_v_g(psv) /= zero .OR. &
271  ps_w_g(psv) /= zero) THEN
272 
273  WRITE(err_msg,1101) psv, trim(ivar('PS_MASSFLOW_G',psv))
274  CALL flush_err_msg(abort=.true.)
275  ENDIF
276 
277  1101 FORMAT('Error 1101: Invalid specification for point source ',i3,&
278  '.',/a,' is zero but velocity is given.',/'Please correct ', &
279  'the mfix.dat file.')
280 
281 ! Verify a physical mass flow
282  ELSEIF(ps_massflow_g(psv) < zero) THEN
283  WRITE(err_msg,1102) psv, trim(ivar('PS_MASSFLOW_G',psv))
284  CALL flush_err_msg(abort=.true.)
285 
286  1102 FORMAT('Error 1102: Invalid specifications for point source ',i3,&
287  '.',/a,' < 0.0. Point sources can only add mass to a system',/&
288  'Please correct the mfix.dat file.')
289 
290 
291 ! Mass flow is specified:
292  ELSE
293 
294 ! Velocity does not have to be defined (no momentum source). If the
295 ! components are UNDEFINED, zero them out.
296  IF(ps_u_g(psv) == undefined) ps_u_g(psv) = zero
297  IF(ps_v_g(psv) == undefined) ps_v_g(psv) = zero
298  IF(ps_w_g(psv) == undefined) ps_w_g(psv) = zero
299 
300 ! Sum together defiend gas phase species mass fractions.
301  sum = zero
302  DO n = 1, nmax(0)
303  IF(ps_x_g(psv,n) /= undefined) THEN
304  sum = sum + ps_x_g(psv,n)
305  ELSE
306  ps_x_g(psv,n) = zero
307  ENDIF
308  ENDDO
309 
310 ! Enforce that the species mass fractions must sum to one.
311  IF(.NOT.compare(one,sum)) THEN
312 
313  IF(species_eq(0)) THEN
314  WRITE(err_msg, 1110) psv
315  CALL flush_err_msg(abort=.true.)
316 
317  1110 FORMAT('Error 1110: PS_X_g(',i3,',:) do NOT sum to ONE and the ',&
318  'gas phase',/'species equations are solved. Please correct ', &
319  'the mfix.dat file.')
320 
321  ELSEIF(.NOT.compare(sum,zero)) THEN
322  WRITE(err_msg, 1111) psv
323  CALL flush_err_msg(abort=.true.)
324 
325  1111 FORMAT('Error 1111: PS_X_g(',i3,',:) do not sum to ONE or ZERO ',&
326  'and they',/'are not needed. Please correct the mfix.dat ', &
327  'the mfix.dat file.')
328 
329  ELSE
330  ps_x_g(psv,:) = zero
331  ps_x_g(psv,1) = one
332  ENDIF
333 
334  ENDIF
335 
336 ! Verify that a temperature is provided.
337  IF(energy_eq)THEN
338  IF(ps_t_g(psv) == undefined) THEN
339  WRITE(err_msg,1000) trim(ivar('PS_T_g',psv))
340  CALL flush_err_msg(abort=.true.)
341 
342 ! Verify that a given temperature is physical.
343  ELSEIF(ps_t_g(psv) <= zero) THEN
344  WRITE(err_msg,1001) psv, trim(ivar('PS_T_g',psv)), &
345  trim(ival(ps_t_g(psv)))
346  CALL flush_err_msg(abort=.true.)
347  ENDIF
348  ENDIF
349 
350  ENDIF
351 
352  CALL finl_err_msg
353 
354  RETURN
355 
356  1000 FORMAT('Error 1000: Required input not specified: ',a,/'Please ',&
357  'correct the mfix.dat file.')
358 
359  1001 FORMAT('Error 1001: Illegal or unknown input: ',a,' = ',a,/ &
360  'Please correct the mfix.dat file.')
361 
362  END SUBROUTINE check_ps_gas_phase
363 
364 
365 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
366 ! !
367 ! SUBROUTINE: CHECK_PS_SOLIDS_PHASES !
368 ! Author: J. Musser Date: 10-JUN-13 !
369 ! !
370 ! Purpose: Check point source specifications. !
371 ! !
372 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
373  SUBROUTINE check_ps_solids_phases(PSV)
375 ! Global Variables:
376 !---------------------------------------------------------------------//
377 ! Solids phase mass flowrate for PS
378  use ps, only: ps_massflow_s
379 ! Solids phase velocity for PS
380  use ps, only: ps_u_s, ps_v_s, ps_w_s
381 ! Solids phase tempture and species mass fractions
382  use ps, only: ps_t_s, ps_x_s
383 ! Flag: Solve energy equations.
384  use run, only: energy_eq
385 ! Flag: Solve species equations.
386  use run, only: species_eq
387 ! Type of each solids phase.
388  use run, only: solids_model
389 ! Number of (TFM) solids.
390  use physprop, only: smax
391 ! Number of discrete solids phases.
392  use discretelement, only: des_mmax
393 ! Number of solids species.
394  use physprop, only: smax, nmax
395 
396 ! Global Parameters:
397 !---------------------------------------------------------------------//
398 ! Parameter constants
399  use param1, only: zero, one, undefined
400 
401 ! Use the error manager for posting error messages.
402 !---------------------------------------------------------------------//
403  use error_manager
404 
405  use toleranc
406 
407  implicit none
408 
409 ! Dummy Arguments:
410 !---------------------------------------------------------------------//
411  INTEGER, INTENT(in) :: PSV
412 
413 ! Local Variables:
414 !---------------------------------------------------------------------//
415 ! Total number of solid phases
416  INTEGER :: MMAX_TOT
417 ! Loop counters
418  INTEGER :: M, N
419 ! Sum of solids mass fractions.
420  DOUBLE PRECISION :: SUM
421 !......................................................................!
422 
423 ! Initialize the error manager.
424  CALL init_err_msg("CHECK_PS_SOLIDS_PHASES")
425 
426 ! The total number of solids phases (all models).
427  mmax_tot = smax + des_mmax
428 
429  DO m=1, mmax_tot
430 
431 ! Check mass flow and velocity
432  IF(ps_massflow_s(psv,m) == undefined) THEN
433  IF(ps_u_s(psv,m) /= undefined .OR. &
434  ps_v_s(psv,m) /= undefined .OR. &
435  ps_w_s(psv,m) /= undefined) THEN
436 
437  WRITE(err_msg,1100)psv, trim(ivar('PS_MASSFLOW_S',psv,m))
438  CALL flush_err_msg(abort=.true.)
439 
440  1100 FORMAT('Error 1100: Invalid specification for point source ',i3,&
441  '.',/a,' is undefined but velocity is given.',/'Please ', &
442  'correct the mfix.dat file.')
443 
444  ELSE
445  ps_massflow_s(psv,m) = zero
446  ps_u_s(psv,m) = zero
447  ps_v_s(psv,m) = zero
448  ps_w_s(psv,m) = zero
449  ENDIF
450 
451  ELSEIF(ps_massflow_s(psv,m) == zero) THEN
452  IF(ps_u_s(psv,m) /= zero .OR. &
453  ps_v_s(psv,m) /= zero .OR. &
454  ps_w_s(psv,m) /= zero) THEN
455 
456  WRITE(err_msg,1101)psv, trim(ivar('PS_MASSFLOW_S',psv,m))
457  CALL flush_err_msg(abort=.true.)
458  ENDIF
459 
460  1101 FORMAT('Error 1100: Invalid specification for point source ',i3,&
461  '.',/a,' is zero but velocity is given.',/'Please correct ', &
462  'the mfix.dat file.')
463 
464  ELSEIF(ps_massflow_s(psv,m) < zero) THEN
465  WRITE(err_msg,1102) psv, trim(ivar('PS_MASSFLOW_S',psv,m))
466  CALL flush_err_msg(abort=.true.)
467 
468  1102 FORMAT('Error 1102: Invalid specifications for point source ',i3,&
469  '.',/a,' < 0.0. Point sources can only add mass to a system',/&
470  'Please correct the mfix.dat file.')
471 
472 
473 ! Mass flow is specified:
474  ELSE
475 
476 ! Currently, only TFM solids can be used with point sources. However,
477 ! the could be implemented for PIC solids as well.
478  SELECT CASE(solids_model(m))
479  CASE ('DEM','PIC')
480  WRITE(err_msg, 1110) psv, solids_model(m)
481  CALL flush_err_msg(abort=.true.)
482  CASE DEFAULT
483  END SELECT
484 
485  1110 FORMAT('Error 1110: Invalid specifications for point source ',i3,&
486  '.',/'Point sources are not supported for ',a,' solids.',/ &
487  'Please correct the mfix.dat file.')
488 
489 ! Velocity does not have to be defined (no momentum source). If the
490 ! components are UNDEFINED, zero them out.
491  IF(ps_u_s(psv,m) == undefined) ps_u_s(psv,m) = zero
492  IF(ps_v_s(psv,m) == undefined) ps_v_s(psv,m) = zero
493  IF(ps_w_s(psv,m) == undefined) ps_w_s(psv,m) = zero
494 
495 ! Sum together defiend gas phase species mass fractions.
496  sum = zero
497  DO n = 1, nmax(m)
498  IF(ps_x_s(psv,m,n) /= undefined) THEN
499  sum = sum + ps_x_s(psv,m,n)
500  ELSE
501  ps_x_s(psv,m,n) = zero
502  ENDIF
503  ENDDO
504 
505 ! Enforce that the species mass fractions must sum to one.
506  IF(.NOT.compare(one,sum)) THEN
507 
508  IF(species_eq(m)) THEN
509  WRITE(err_msg, 1120) psv,m
510  CALL flush_err_msg(abort=.true.)
511 
512  1120 FORMAT('Error 1120: PS_X_s(',i3,',',i2,',:) do NOT sum to ONE ', &
513  'and the solids phase',/'species equations are solved. ', &
514  'Please correct the mfix.dat file.')
515 
516  ELSEIF(.NOT.compare(sum,zero)) THEN
517  WRITE(err_msg, 1121) psv,m
518  CALL flush_err_msg(abort=.true.)
519 
520  1121 FORMAT('Error 1121: PS_X_s(',i3,',',i2,',:) do not sum to ONE ', &
521  'or ZERO and they',/'are not needed. Please correct the ', &
522  'mfix.dat the mfix.dat file.')
523 
524  ELSE
525  ps_x_s(psv,m,1) = one
526  ps_x_s(psv,m,2:) = zero
527  ENDIF
528 
529  ENDIF
530 
531 ! Verify that a temperature is provided.
532  IF(energy_eq)THEN
533  IF(ps_t_s(psv,m) == undefined) THEN
534  WRITE(err_msg,1000) trim(ivar('PS_T_s',psv,m))
535  CALL flush_err_msg(abort=.true.)
536 
537 ! Verify that a given temperature is physical.
538  ELSEIF(ps_t_s(psv,m) <= zero) THEN
539  WRITE(err_msg,1001) psv, trim(ivar('PS_T_s',psv,m)), &
540  trim(ival(ps_t_s(psv,m)))
541  CALL flush_err_msg(abort=.true.)
542  ENDIF
543  ENDIF
544  ENDIF
545  ENDDO
546 
547  CALL finl_err_msg
548 
549  RETURN
550 
551  1000 FORMAT('Error 1000: Required input not specified: ',a,/'Please ',&
552  'correct the mfix.dat file.')
553 
554  1001 FORMAT('Error 1001: Illegal or unknown input: ',a,' = ',a,/ &
555  'Please correct the mfix.dat file.')
556 
557  END SUBROUTINE check_ps_solids_phases
558 
559 
560 
561 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
562 ! !
563 ! SUBROUTINE: CHECK_PS_OVERFLOW !
564 ! Author: J. Musser Date: 10-JUN-13 !
565 ! !
566 ! Purpose: Check point source specifications. !
567 ! !
568 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
569  SUBROUTINE check_ps_overflow(PSV)
571 ! Global Variables:
572 !---------------------------------------------------------------------//
573 ! Gas phase mass flowrate for PS and velocities
574  use ps, only: ps_massflow_g, ps_u_g, ps_v_g, ps_w_g
575 ! Gas phase tempture and species mass fractions
576  use ps, only: ps_t_g, ps_x_g
577 ! Solids phase mass flowrate and velocity
578  use ps, only: ps_massflow_s, ps_u_s, ps_v_s, ps_w_s
579 ! Solids phase tempture and species mass fractions
580  use ps, only: ps_t_s, ps_x_s
581 
582 ! Global Parameters:
583 !---------------------------------------------------------------------//
584 ! Maximum input array sizes.
585  use param, only: dim_m, dim_n_g, dim_n_s
586 ! Parameter constants
587  use param1, only: undefined
588 
589 ! Use the error manager for posting error messages.
590 !---------------------------------------------------------------------//
591  use error_manager
592 
593  implicit none
594 
595 ! Dummy Arguments:
596 !---------------------------------------------------------------------//
597  INTEGER, INTENT(in) :: PSV
598 
599 ! Local Variables:
600 !---------------------------------------------------------------------//
601 ! Loop counters
602  INTEGER :: M, N
603 !......................................................................!
604 
605 
606 ! Initialize the error manager.
607  CALL init_err_msg("CHECK_PS_OVERFLOW")
608 
609 
610  IF(ps_massflow_g(psv) /= undefined) THEN
611  WRITE(err_msg,1010) trim(ivar('PS_MASSFLOW_G',psv))
612  CALL flush_err_msg(abort=.true.)
613  ELSEIF(ps_u_g(psv) /= undefined) THEN
614  WRITE(err_msg,1010) trim(ivar('PS_U_g',psv))
615  CALL flush_err_msg(abort=.true.)
616  ELSEIF(ps_v_g(psv) /= undefined) THEN
617  WRITE(err_msg,1010) trim(ivar('PS_V_g',psv))
618  CALL flush_err_msg(abort=.true.)
619  ELSEIF(ps_w_g(psv) /= undefined) THEN
620  WRITE(err_msg,1010) trim(ivar('PS_W_g',psv))
621  CALL flush_err_msg(abort=.true.)
622  ELSEIF(ps_t_g(psv) /= undefined) THEN
623  WRITE(err_msg,1010) trim(ivar('PS_T_g',psv))
624  CALL flush_err_msg(abort=.true.)
625  ENDIF
626  DO n = 1, dim_n_g
627  IF(ps_x_g(psv,n) /= undefined) THEN
628  WRITE(err_msg,1010) trim(ivar('PS_X_G',psv,n))
629  CALL flush_err_msg(abort=.true.)
630  ENDIF
631  ENDDO
632 
633  DO m=1, dim_m
634  IF(ps_massflow_s(psv,m) /= undefined) THEN
635  WRITE(err_msg,1010) trim(ivar('PS_MASSFLOW_S',psv,m))
636  CALL flush_err_msg(abort=.true.)
637  ELSEIF(ps_u_s(psv,m) /= undefined) THEN
638  WRITE(err_msg,1010) trim(ivar('PS_U_s',psv,m))
639  CALL flush_err_msg(abort=.true.)
640  ELSEIF(ps_v_s(psv,m) /= undefined) THEN
641  WRITE(err_msg,1010) trim(ivar('PS_V_s',psv,m))
642  CALL flush_err_msg(abort=.true.)
643  ELSEIF(ps_w_s(psv,m) /= undefined) THEN
644  WRITE(err_msg,1010) trim(ivar('PS_W_s',psv,m))
645  CALL flush_err_msg(abort=.true.)
646  ELSEIF(ps_t_s(psv,m) /= undefined) THEN
647  WRITE(err_msg,1010) trim(ivar('PS_T_s',psv,m))
648  CALL flush_err_msg(abort=.true.)
649  ENDIF
650  DO n = 1, dim_n_s
651  IF(ps_x_s(psv,m,n) /= undefined) THEN
652  WRITE(err_msg,1010) trim(ivar('PS_X_S',psv,m,n))
653  CALL flush_err_msg(abort=.true.)
654  ENDIF
655  ENDDO
656  ENDDO
657 
658  CALL finl_err_msg
659 
660  RETURN
661 
662  1010 FORMAT('Error 1010: ',a,' specified in an undefined PS region.')
663 
664  END SUBROUTINE check_ps_overflow
integer, dimension(dimension_ps) ps_i_w
Definition: ps_mod.f:40
integer, parameter dim_n_g
Definition: param_mod.f:69
double precision, dimension(dimension_ps) ps_v_g
Definition: ps_mod.f:52
double precision, dimension(dimension_ps) ps_t_g
Definition: ps_mod.f:62
subroutine check_ps_overflow(PSV)
double precision, dimension(dimension_ps, dim_m) ps_v_s
Definition: ps_mod.f:70
character(len=32) function ivar(VAR, i1, i2, i3)
subroutine finl_err_msg
logical function compare(V1, V2)
Definition: toleranc_mod.f:94
logical no_i
Definition: geometry_mod.f:20
subroutine check_ps_gas_phase(PSV)
double precision, parameter one
Definition: param1_mod.f:29
logical, dimension(0:dim_m) species_eq
Definition: run_mod.f:115
integer, parameter dim_m
Definition: param_mod.f:67
integer, dimension(dimension_ps) ps_j_n
Definition: ps_mod.f:43
character(len=3), dimension(dim_m) solids_model
Definition: run_mod.f:253
double precision, dimension(dimension_ps) ps_y_n
Definition: ps_mod.f:35
double precision, parameter undefined
Definition: param1_mod.f:18
double precision, dimension(dimension_ps) ps_z_b
Definition: ps_mod.f:36
double precision, dimension(dimension_ps) ps_massflow_g
Definition: ps_mod.f:48
logical, dimension(dimension_ps) ps_defined
Definition: ps_mod.f:29
double precision, dimension(dimension_ps) ps_w_g
Definition: ps_mod.f:53
double precision, dimension(dimension_ps, dim_m) ps_t_s
Definition: ps_mod.f:80
double precision, dimension(dimension_ps) ps_x_e
Definition: ps_mod.f:33
double precision, dimension(dimension_ps) ps_x_w
Definition: ps_mod.f:32
subroutine init_err_msg(CALLER)
integer, dimension(dimension_ps) ps_k_b
Definition: ps_mod.f:44
double precision, dimension(dimension_ps, dim_m) ps_u_s
Definition: ps_mod.f:69
subroutine get_ps(PSV)
Definition: get_ps.f:10
double precision, dimension(dimension_ps) ps_y_s
Definition: ps_mod.f:34
integer, dimension(dimension_ps) ps_k_t
Definition: ps_mod.f:45
double precision xlength
Definition: geometry_mod.f:33
Definition: run_mod.f:13
subroutine check_ps_geometry
double precision, dimension(dimension_ps, dim_n_g) ps_x_g
Definition: ps_mod.f:59
Definition: param_mod.f:2
integer, parameter dimension_ps
Definition: param_mod.f:65
logical no_k
Definition: geometry_mod.f:28
integer, dimension(0:dim_m) nmax
Definition: physprop_mod.f:119
logical no_j
Definition: geometry_mod.f:24
logical energy_eq
Definition: run_mod.f:100
integer, parameter undefined_i
Definition: param1_mod.f:19
double precision, dimension(dimension_ps, dim_m) ps_massflow_s
Definition: ps_mod.f:66
character(len=line_length), dimension(line_count) err_msg
integer, parameter dim_n_s
Definition: param_mod.f:71
double precision, dimension(dimension_ps) ps_z_t
Definition: ps_mod.f:37
double precision, dimension(dimension_ps, dim_m) ps_w_s
Definition: ps_mod.f:71
double precision, dimension(dimension_ps, dim_m, dim_n_s) ps_x_s
Definition: ps_mod.f:77
Definition: ps_mod.f:22
subroutine check_point_sources
integer smax
Definition: physprop_mod.f:22
double precision ylength
Definition: geometry_mod.f:35
logical point_source
Definition: ps_mod.f:27
double precision, dimension(dimension_ps) ps_u_g
Definition: ps_mod.f:51
integer, dimension(dimension_ps) ps_j_s
Definition: ps_mod.f:42
integer, dimension(dimension_ps) ps_i_e
Definition: ps_mod.f:41
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)
subroutine check_ps_solids_phases(PSV)