Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ReadTbcSection.F90
Go to the documentation of this file.
1 ! *********************************************************************
2 ! * Rocstar Simulation Suite *
3 ! * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4 ! * *
5 ! * Illinois Rocstar LLC *
6 ! * Champaign, IL *
7 ! * www.illinoisrocstar.com *
8 ! * sales@illinoisrocstar.com *
9 ! * *
10 ! * License: See LICENSE file in top level of distribution package or *
11 ! * http://opensource.org/licenses/NCSA *
12 ! *********************************************************************
13 ! *********************************************************************
14 ! * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15 ! * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16 ! * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17 ! * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18 ! * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19 ! * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20 ! * Arising FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21 ! * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22 ! *********************************************************************
23 ! ******************************************************************************
24 !
25 ! Purpose: Read in user input related to time-dependent boundary conditions.
26 !
27 ! Description: None.
28 !
29 ! Input:
30 ! pRegion Pointer to region type
31 ! tbcType tbc type
32 !
33 ! Output: None.
34 !
35 ! Notes:
36 ! 1. Sample input sections for TBCs are given in the files updateTbc*.F90
37 !
38 ! ******************************************************************************
39 !
40 ! $Id: RFLU_ReadTbcSection.F90,v 1.12 2008/12/06 08:44:12 mtcampbe Exp $
41 !
42 ! Copyright: (c) 2002-2005 by the University of Illinois
43 !
44 ! ******************************************************************************
45 
46 SUBROUTINE rflu_readtbcsection(pRegion,tbcType)
47 
48  USE moddatatypes
49  USE modbndpatch
50  USE moddatastruct, ONLY: t_region
51  USE modglobal, ONLY: t_global
52  USE modgrid, ONLY: t_grid
53  USE moderror
54  USE modparameters
55 
57 
58  IMPLICIT NONE
59 
60 ! ******************************************************************************
61 ! Definitions and declarations
62 ! ******************************************************************************
63 
64 ! ==============================================================================
65 ! Arguments
66 ! ==============================================================================
67 
68  INTEGER, INTENT(IN) :: tbctype
69  TYPE(t_region), POINTER :: pregion
70 
71 ! ==============================================================================
72 ! Locals
73 ! ==============================================================================
74 
75  INTEGER, PARAMETER :: nvals_max = 10
76 
77  LOGICAL :: defined(nvals_max)
78  CHARACTER(10) :: keys(nvals_max),bcvar,pwisestr
79  CHARACTER(32) :: bctitle
80  CHARACTER(256) :: fname, line
81  CHARACTER(CHRLEN) :: bcname,rcsidentstring
82  INTEGER :: bccode,distrib,errorflag,ftype,ipatch,ireg,n,nbvals,ndata, &
83  njumps,nparams,nsvals,nswitches,nvals,prbeg,prend,var
84  INTEGER :: switches(nvals_max)
85  REAL(RFREAL) :: vals(nvals_max)
86  REAL(RFREAL), DIMENSION(:), ALLOCATABLE :: holdvals,morevals
87  TYPE(t_patch), POINTER :: ppatch
88  TYPE(t_bcvalues), POINTER :: bc
89  TYPE(t_tbcvalues), POINTER :: tbc
90  TYPE(t_global), POINTER :: global
91 
92 ! ******************************************************************************
93 ! Start
94 ! ******************************************************************************
95 
96  rcsidentstring = '$RCSfile: RFLU_ReadTbcSection.F90,v $ $Revision: 1.12 $'
97 
98  global => pregion%global
99 
100  CALL registerfunction( global,'RFLU_ReadTbcSection',&
101  'RFLU_ReadTbcSection.F90' )
102 
103  IF ( global%flowType == flow_steady ) THEN
104  CALL errorstop( global,err_steady,__line__ )
105  END IF ! global%flowType
106 
107 ! ******************************************************************************
108 ! Read BC type and variable name
109 ! ******************************************************************************
110 
111  DO ! find first non-blank line of section
112  READ(if_input,'(A256)',err=10,end=10) line
113  IF ( len_trim(line) > 0 ) THEN
114  EXIT
115  END IF ! LEN_TRIM
116  END DO ! <empty>
117 
118  READ(line,*) bctitle,bcvar
119 
120 ! ******************************************************************************
121 ! Find the integers describing module type and variable name
122 ! ******************************************************************************
123 
124  bccode = -1
125  ftype = -1
126 
127  SELECT CASE( trim(bctitle) )
128 
129 ! ==============================================================================
130 ! Mixture BCs: ftype = FTYPE_MIXT
131 ! ==============================================================================
132 
133 ! ------------------------------------------------------------------------------
134 ! Slip wall
135 ! ------------------------------------------------------------------------------
136 
137  CASE ('SLIPW')
138  ftype = ftype_mixt
139  bccode = bc_slipwall
140  CALL errorstop(global,err_bc_varname,__line__,'SLIPW')
141 
142 ! ------------------------------------------------------------------------------
143 ! No-slip wall
144 ! ------------------------------------------------------------------------------
145 
146  CASE ('NOSLIP')
147  ftype = ftype_mixt
148  bccode = bc_noslipwall
149 
150  SELECT CASE(trim(bcvar))
151  CASE ('TWALL')
152  var = bcdat_noslip_twall
153  CASE default
154  CALL errorstop(global,err_bc_varname,__line__,'NOSLIP')
155  END SELECT
156 
157 ! ------------------------------------------------------------------------------
158 ! Inflow (based on total quantities and flow angles)
159 ! ------------------------------------------------------------------------------
160 
161  CASE ('INFLOW_TOTANG')
162  ftype = ftype_mixt
163  bccode = bc_inflow_totang
164 
165  SELECT CASE(trim(bcvar))
166  CASE ('PTOT')
167  var = bcdat_inflow_ptot
168  CASE ('TTOT')
169  var = bcdat_inflow_ttot
170  CASE ('BETAH')
171  var = bcdat_inflow_betah
172  CASE ('BETAV')
173  var = bcdat_inflow_betav
174  CASE ('MACH')
175  var = bcdat_inflow_mach
176  CASE default
177  CALL errorstop(global,err_bc_varname,__line__,'INFLOW_TOTANG')
178  END SELECT
179 
180 ! ------------------------------------------------------------------------------
181 ! Inflow (based on velocity and temperature)
182 ! ------------------------------------------------------------------------------
183 
184  CASE ('INFLOW_VELTEMP')
185  ftype = ftype_mixt
186  bccode = bc_inflow_veltemp
187 
188  SELECT CASE(trim(bcvar))
189  CASE ('VELX')
190  var = bcdat_inflow_u
191  CASE ('VELY')
192  var = bcdat_inflow_v
193  CASE ('VELZ')
194  var = bcdat_inflow_w
195  CASE ('TEMP')
196  var = bcdat_inflow_t
197  CASE ('PRESS')
198  var = bcdat_inflow_p
199  CASE default
200  CALL errorstop(global,err_bc_varname,__line__,'INFLOW_VELTEMP')
201  END SELECT
202 
203 ! ------------------------------------------------------------------------------
204 ! Outflow
205 ! ------------------------------------------------------------------------------
206 
207  CASE ('OUTFLOW')
208  ftype = ftype_mixt
209  bccode = bc_outflow
210 
211  SELECT CASE(trim(bcvar))
212  CASE ('PRESS')
213  var = bcdat_outflow_press
214  CASE default
215  CALL errorstop(global,err_bc_varname,__line__,'OUTFLOW')
216  END SELECT
217 
218 ! ------------------------------------------------------------------------------
219 ! Farfield
220 ! ------------------------------------------------------------------------------
221 
222  CASE ('FARF')
223  ftype = ftype_mixt
224  bccode = bc_farfield
225 
226  SELECT CASE(trim(bcvar))
227  CASE ('MACH')
228  var = bcdat_farf_mach
229  CASE ('ATTACK')
230  var = bcdat_farf_attack
231  CASE ('SLIP')
232  var = bcdat_farf_slip
233  CASE ('PRESS')
234  var = bcdat_farf_press
235  CASE ('TEMP')
236  var = bcdat_farf_temp
237  CASE default
238  CALL errorstop(global,err_bc_varname,__line__,'FARF')
239  END SELECT
240 
241 ! ------------------------------------------------------------------------------
242 ! Injection
243 ! ------------------------------------------------------------------------------
244 
245  CASE ('INJECT')
246  ftype = ftype_mixt
247  bccode = bc_injection
248 
249  SELECT CASE(trim(bcvar))
250  CASE ('MFRATE')
251  var = bcdat_inject_mfrate
252  CASE ('TEMP')
253  var = bcdat_inject_temp
254  CASE default
255  CALL errorstop(global,err_bc_varname,__line__,'INJECT')
256  END SELECT
257 
258 #ifdef TURB
259 ! ==============================================================================
260 ! Turbulence BCs: ftype = FTYPE_TURB
261 ! ==============================================================================
262 
263 ! TO DO
264 ! END TO DO
265 #endif
266 
267 #ifdef SPEC
268 ! ==============================================================================
269 ! Species BCs: ftype = FTYPE_SPEC
270 ! ==============================================================================
271 
272 ! ------------------------------------------------------------------------------
273 ! Inflow
274 ! ------------------------------------------------------------------------------
275 
276  CASE ('SPEC_INFLOW')
277  ftype = ftype_spec
278  bccode = bc_inflow
279 
280  IF ( bcvar(1:4) == 'SPEC' ) THEN
281  IF ( len_trim(bcvar) == 4 ) THEN
282  global%warnCounter = global%warnCounter + 1
283 
284  WRITE(stdout,'(A,1X,A)') solver_name, &
285  '*** WARNING *** Improper specification in SPEC_INFLOW TBC.'
286  WRITE(stdout,'(A,1X,A)') solver_name, &
287  ' Assuming SPEC to mean SPEC1.'
288 
289  var = 1
290  ELSE
291  READ(bcvar(5:),*) var
292  END IF ! LEN_TRIM
293  ELSE
294  CALL errorstop(global,err_bc_varname,__line__,'SPEC_INFLOW')
295  END IF ! bcvar
296 
297 ! ------------------------------------------------------------------------------
298 ! Farfield
299 ! ------------------------------------------------------------------------------
300 
301  CASE ('SPEC_FARF')
302  ftype = ftype_spec
303  bccode = bc_farfield
304 
305  IF ( bcvar(1:4) == 'SPEC' ) THEN
306  IF ( len_trim(bcvar) == 4 ) THEN
307  global%warnCounter = global%warnCounter + 1
308 
309  WRITE(stdout,'(A,1X,A)') solver_name, &
310  '*** WARNING *** Improper specification in SPEC_FARF TBC.'
311  WRITE(stdout,'(A,1X,A)') solver_name, &
312  ' Assuming SPEC to mean SPEC1.'
313 
314  var = 1
315  ELSE
316  READ(bcvar(5:),*) var
317  END IF ! LEN_TRIM
318  ELSE
319  CALL errorstop(global,err_bc_varname,__line__,'SPEC_FARF')
320  END IF ! bcvar
321 
322 ! ------------------------------------------------------------------------------
323 ! Injection
324 ! ------------------------------------------------------------------------------
325 
326  CASE ('SPEC_INJECT')
327  ftype = ftype_spec
328  bccode = bc_injection
329 
330  IF ( bcvar(1:4) == 'SPEC' ) THEN
331  IF ( len_trim(bcvar) == 4 ) THEN
332  global%warnCounter = global%warnCounter + 1
333 
334  WRITE(stdout,'(A,1X,A)') solver_name, &
335  '*** WARNING *** Improper specification in SPEC_INJECT TBC.'
336  WRITE(stdout,'(A,1X,A)') solver_name, &
337  ' Assuming SPEC to mean SPEC1.'
338 
339  var = 1
340  ELSE
341  READ(bcvar(5:),*) var
342  END IF ! LEN_TRIM
343  ELSE
344  CALL errorstop(global,err_bc_varname,__line__,'SPEC_INJECT')
345  END IF ! bcvar
346 #endif
347 
348 #ifdef RADI
349 ! ==============================================================================
350 ! Radiation BCs: ftype = FTYPE_RADI
351 ! ==============================================================================
352 
353 ! TO DO
354 ! END TO DO
355 #endif
356 
357  CASE default
358  CALL errorstop(global,err_unknown_bc,__line__)
359  END SELECT ! TRIM(bctitle)
360 
361 ! ******************************************************************************
362 ! Specify keywords and search for them
363 ! ******************************************************************************
364 
365  keys(tbcdat_ontime) = 'ONTIME'
366  keys(tbcdat_offtime) = 'OFFTIME'
367  keys(tbcdat_amp) = 'AMP'
368 
369  SELECT CASE ( tbctype )
370  CASE ( tbc_sinusoidal )
371  keys(tbcdat_freq) = 'FREQ'
372  keys(tbcdat_phase) = 'PHASE'
373  nvals = 5
374  CASE ( tbc_stochastic )
375  keys(tbcdat_timecor) = 'TIMECOR'
376  keys(tbcdat_shape) = 'SHAPE'
377  keys(tbcdat_mincut) = 'MINCUT'
378  keys(tbcdat_maxcut) = 'MAXCUT'
379  nvals = 7
380  CASE ( tbc_whitenoise )
381  keys(4) = 'SUBSTEP' ! becomes a switch
382  nvals = 4
383  CASE ( tbc_piecewise )
384  keys(3) = 'ORDER' ! becomes a switch (overwrites 'AMP')
385  keys(4) = 'NJUMPS' ! becomes a switch
386  nvals = 4
387  CASE default
388  CALL errorstop(global,err_reached_default,__line__)
389  END SELECT ! tbcType
390 
391  CALL readpatchsection(global,if_input,nvals,keys,vals,prbeg,prend, &
392  distrib,fname,bcname,defined)
393 
394  IF (.NOT. defined(tbcdat_ontime)) vals(tbcdat_ontime) = -1.e20_rfreal
395  IF (vals(tbcdat_ontime) < 0.0_rfreal) vals(tbcdat_ontime) = -1.e20_rfreal
396 
397  IF (.NOT. defined(tbcdat_offtime)) vals(tbcdat_offtime) = 1.e20_rfreal
398  IF (vals(tbcdat_offtime) < 0.0_rfreal) vals(tbcdat_offtime) = 1.e20_rfreal
399 
400  IF (tbctype /= tbc_piecewise) THEN
401  IF (.NOT. defined(tbcdat_amp)) &
402  CALL errorstop(global,err_bcval_missing,__line__)
403 
404  IF (vals(tbcdat_amp) <= 0.0_rfreal) &
405  CALL errorstop(global,err_val_bcval,__line__)
406  END IF ! tbcType
407 
408 ! ******************************************************************************
409 ! Set values of vals and switches
410 ! ******************************************************************************
411 
412  SELECT CASE ( tbctype )
413 
414 ! ==============================================================================
415 ! Sinusoidal variation
416 ! ==============================================================================
417 
418  CASE ( tbc_sinusoidal )
419  nswitches = 0
420  nparams = 5
421  nsvals = 1
422  nbvals = 0
423 
424  IF (.NOT. defined(tbcdat_freq)) &
425  CALL errorstop(global,err_bcval_missing,__line__)
426 
427  IF (.NOT. defined(tbcdat_phase)) vals(tbcdat_phase) = 0.0_rfreal
428 
429  vals(tbcdat_phase) = vals(tbcdat_phase)*global%deg2rad
430 
431 ! ==============================================================================
432 ! Stochastic variation
433 ! ==============================================================================
434 
435  CASE ( tbc_stochastic )
436  nswitches = 0
437  nparams = 7
438  nsvals = 0
439  nbvals = 3
440 
441  IF (.NOT. defined(tbcdat_timecor)) &
442  CALL errorstop(global,err_bcval_missing,__line__)
443 
444  IF (vals(tbcdat_timecor) <= 0.0_rfreal) &
445  CALL errorstop(global,err_val_bcval,__line__)
446 
447  IF (.NOT. defined(tbcdat_shape)) vals(tbcdat_shape) = 1.0_rfreal
448 
449  IF (vals(tbcdat_shape) < 0.1_rfreal .OR. vals(tbcdat_shape) > 10.0_rfreal) &
450  CALL errorstop(global,err_val_bcval,__line__)
451 
452  IF (.NOT. defined(tbcdat_mincut)) vals(tbcdat_mincut) = -1.0_rfreal
453 
454  IF (.NOT. defined(tbcdat_maxcut)) vals(tbcdat_maxcut) = -1.0_rfreal
455 
456 ! ==============================================================================
457 ! Whitenoise variation
458 ! ==============================================================================
459 
460  CASE ( tbc_whitenoise )
461  nswitches = 1
462  nparams = 3
463  nsvals = 0
464  nbvals = 1
465 
466  IF (.NOT. defined(4)) switches(tbcswi_substep) = tbcopt_step
467 
468  IF ( vals(4) < 0.1_rfreal ) THEN
469  switches(tbcswi_substep) = tbcopt_step
470  ELSE
471  switches(tbcswi_substep) = tbcopt_substep
472  END IF ! vals
473 
474 ! ==============================================================================
475 ! Piecewise variation
476 ! ==============================================================================
477 
478  CASE ( tbc_piecewise )
479  nswitches = 2
480  nsvals = 1
481  nbvals = 0
482 
483  IF (.NOT. defined(3)) switches(tbcswi_order) = tbcopt_constant
484 
485  IF (vals(3) < 0.1_rfreal) THEN
486  switches(tbcswi_order) = tbcopt_constant
487  ELSE
488  switches(tbcswi_order) = tbcopt_linear
489  END IF ! vals
490 
491  IF (.NOT. defined(4)) THEN
492  CALL errorstop( global,err_bcval_missing,__line__ )
493  END IF
494 
495  njumps = nint(vals(4))
496  switches(tbcswi_njumps) = njumps
497 
498  IF (njumps < 1 .OR. njumps > 1000000) &
499  CALL errorstop(global,err_val_bcval,__line__)
500 
501  SELECT CASE (switches(tbcswi_order))
502  CASE (tbcopt_constant)
503  nparams = tbcdat_dat0 + 2*njumps+1
504  CASE (tbcopt_linear)
505  nparams = tbcdat_dat0 + 2*njumps
506  END SELECT ! switches(TBCSWI_ORDER)
507 
508  ALLOCATE(morevals(nparams),stat=errorflag)
509  global%error = errorflag
510  IF ( global%error /= err_none ) THEN
511  CALL errorstop(global,err_allocate,__line__)
512  END IF ! global%error
513 
514  morevals(1:tbcdat_dat0) = vals(1:tbcdat_dat0)
515 
516  DO n=tbcdat_dat0+1,nparams
517  DO ! find first non-blank line
518  READ(if_input,'(A256)',err=10,end=10) line
519  IF (len_trim(line) > 0) EXIT
520  END DO ! <empty>
521 
522  IF (line(1:1) == '#') goto 10
523 
524  READ(line,*) pwisestr, morevals(n)
525  END DO ! n
526 
527 ! - check that the times entered are increasing
528 
529  DO n=tbcdat_dat0+4,nparams,2
530  IF (morevals(n) .le. morevals(n-2)) &
531  CALL errorstop(global,err_val_bcval,__line__)
532  END DO ! n
533 
534  DO ! find end of section
535  READ(if_input,'(A256)',err=10,end=10) line
536  IF (line(1:1) == '#') EXIT
537  END DO ! <empty>
538 
539 ! ==============================================================================
540 ! Default
541 ! ==============================================================================
542 
543  CASE default
544  CALL errorstop(global,err_reached_default,__line__)
545  END SELECT ! tbctype
546 
547 ! ******************************************************************************
548 ! Copy values to variables
549 ! ******************************************************************************
550 
551  DO ipatch = 1,pregion%grid%nPatches
552  ppatch => pregion%patches(ipatch)
553 
554 ! ==============================================================================
555 ! Check whether this global patch exists in this region
556 ! ==============================================================================
557 
558  IF ( (ppatch%bcType >= bccode .AND. &
559  ppatch%bcType <= bccode+bc_range) .AND. &
560  (ppatch%iPatchGlobal >= prbeg .AND. &
561  ppatch%iPatchGlobal <= prend) ) THEN
562  SELECT CASE (ftype)
563  CASE (ftype_mixt)
564  bc => ppatch%mixt
565  CASE (ftype_turb)
566  bc => ppatch%turb
567  CASE (ftype_spec)
568  bc => ppatch%spec
569  CASE (ftype_radi)
570  bc => ppatch%valRadi
571  CASE default
572  CALL errorstop(global,err_reached_default,__line__)
573  END SELECT ! ftype
574 
575  ndata = bc%nData
576  tbc => bc%tbcs(var)
577 
578  IF ( var > ndata ) THEN
579  CALL errorstop(global,err_bc_varname,__line__)
580  END IF ! var
581 
582  IF ( tbc%tbcType /= tbc_none ) THEN
583  CALL errorstop(global,err_patch_overspec,__line__,'TBC')
584  END IF ! tbc%tbcType
585 
586  tbc%tbcType = tbctype
587 
588  IF ( nswitches > 0 ) THEN
589  ALLOCATE(tbc%switches(nswitches),stat=errorflag)
590  global%error = errorflag
591  IF ( global%error /= err_none ) THEN
592  CALL errorstop(global,err_allocate,__line__)
593  END IF ! global%error
594 
595  tbc%switches = switches(1:nswitches)
596  END IF ! nswitches
597 
598  IF ( nparams > 0 ) THEN
599  ALLOCATE(tbc%params(nparams),stat=errorflag)
600  global%error = errorflag
601  IF ( global%error /= err_none ) THEN
602  CALL errorstop(global,err_allocate,__line__)
603  END IF ! global%error
604 
605  IF (tbctype == tbc_piecewise) THEN
606  tbc%params = morevals(1:nparams)
607  ELSE
608  tbc%params = vals(1:nparams)
609  END IF ! tbcType
610  END IF ! nparams
611 
612  IF (nsvals > 0) THEN
613  ALLOCATE(tbc%svals(nsvals),stat=errorflag)
614  global%error = errorflag
615  IF ( global%error /= err_none ) THEN
616  CALL errorstop(global,err_allocate,__line__)
617  END IF ! global%error
618 
619  tbc%svals = 0.0_rfreal
620  END IF ! nsvals
621 
622  IF (nbvals > 0) THEN
623  ALLOCATE(tbc%bvals(nbvals,ppatch%nBFaces),stat=errorflag)
624  global%error = errorflag
625  IF ( global%error /= err_none ) THEN
626  CALL errorstop(global,err_allocate,__line__)
627  END IF ! global%error
628 
629  tbc%bvals = 0.0_rfreal
630 
631  IF ( tbctype == tbc_stochastic ) THEN
632  tbc%bvals(tbcsto_val,:) = 0.5_rfreal * tbc%params(tbcdat_amp)**2
633  END IF ! tbcType
634  END IF ! nbvals
635 
636  distrib = bc%distrib ! original value of bc%distrib
637 
638  IF (tbctype == tbc_stochastic .OR. & ! tbc types forcing BCDAT_DISTRIB
639  tbctype == tbc_whitenoise) bc%distrib = bcdat_distrib
640 
641  IF ( bc%distrib == bcdat_distrib ) THEN
642  ALLOCATE(tbc%mean(ppatch%nBFaces),stat=errorflag)
643  ELSE
644  ALLOCATE(tbc%mean(0:1),stat=errorflag)
645  END IF ! bc%distrib
646 
647  global%error = errorflag
648  IF ( global%error /= err_none ) THEN
649  CALL errorstop(global,err_allocate,__line__)
650  END IF ! global%error
651 
652 ! ==============================================================================
653 ! If distrib needs to change for BC, allocate and fill expanded bc%vals
654 ! arrays
655 ! ==============================================================================
656 
657  IF ( bc%distrib == bcdat_distrib .AND. distrib == bcdat_constant ) THEN
658  ALLOCATE( holdvals(ndata),stat=errorflag )
659  global%error = errorflag
660  IF ( global%error /= err_none ) THEN
661  CALL errorstop(global,err_allocate,__line__)
662  END IF ! global%error
663 
664  holdvals = bc%vals(:,1)
665 
666  DEALLOCATE(bc%vals,stat=errorflag)
667  global%error = errorflag
668  IF ( global%error /= err_none ) THEN
669  CALL errorstop(global,err_deallocate,__line__)
670  END IF ! global%error
671 
672  ALLOCATE(bc%vals(ndata,ppatch%nBFaces),stat=errorflag)
673  global%error = errorflag
674  IF ( global%error /= err_none ) THEN
675  CALL errorstop(global,err_allocate,__line__)
676  END IF ! global%error
677 
678  DO n = 1, ndata
679  bc%vals(n,:) = holdvals(n)
680  END DO ! n
681 
682  DEALLOCATE (holdvals,stat=errorflag)
683  global%error = errorflag
684  IF ( global%error /= err_none ) THEN
685  CALL errorstop(global,err_deallocate,__line__)
686  END IF ! global%error
687 
688  END IF ! bc%distrib
689 
690  tbc%mean = bc%vals(var,:)
691 
692  IF ( global%error /= err_none ) THEN
693  CALL errorstop(global,err_allocate,__line__)
694  END IF ! global%error
695 
696  END IF ! patch in region, correct boundary type
697  END DO ! iPatch
698 
699  IF ( ALLOCATED(morevals) .EQV. .true. ) THEN
700  DEALLOCATE(morevals,stat=errorflag)
701  global%error = errorflag
702  IF ( global%error /= err_none ) THEN
703  CALL errorstop(global,err_deallocate,__line__)
704  END IF ! global%error
705  END IF ! ALLOCATED
706 
707  goto 999
708 
709 ! ******************************************************************************
710 ! Error handling
711 ! ******************************************************************************
712 
713 10 CONTINUE
714  CALL errorstop(global,err_file_read,__line__)
715 
716 999 CONTINUE
717 
718 ! ******************************************************************************
719 ! End
720 ! ******************************************************************************
721 
722  CALL deregisterfunction( global )
723 
724 END SUBROUTINE rflu_readtbcsection
725 
726 ! ******************************************************************************
727 !
728 ! RCS Revision history:
729 !
730 ! $Log: RFLU_ReadTbcSection.F90,v $
731 ! Revision 1.12 2008/12/06 08:44:12 mtcampbe
732 ! Updated license.
733 !
734 ! Revision 1.11 2008/11/19 22:17:26 mtcampbe
735 ! Added Illinois Open Source License/Copyright
736 !
737 ! Revision 1.10 2006/08/19 15:38:48 mparmar
738 ! Renamed patch variables
739 !
740 ! Revision 1.9 2005/04/27 02:09:00 haselbac
741 ! Added support for INFLOW_VELTEMP
742 !
743 ! Revision 1.8 2004/10/19 19:25:17 haselbac
744 ! Removed RFVFx on injecting boundaries, clean-up
745 !
746 ! Revision 1.7 2004/07/23 22:43:15 jferry
747 ! Integrated rocspecies into rocinteract
748 !
749 ! Revision 1.6 2004/02/26 21:01:54 haselbac
750 ! Removed stray TEMPORARY statement
751 !
752 ! Revision 1.5 2004/01/29 22:56:46 haselbac
753 ! Added species stuff, clean-up
754 !
755 ! Revision 1.4 2003/07/22 01:59:42 haselbac
756 ! Added global%warnCounter
757 !
758 ! Revision 1.3 2003/06/10 22:54:42 jferry
759 ! Added Piecewise TBC
760 !
761 ! Revision 1.2 2003/06/04 20:06:38 jferry
762 ! added capability for PEUL BCs with TBCs
763 !
764 ! Revision 1.1 2002/11/16 00:05:31 haselbac
765 ! Moved here from rocflu directory
766 !
767 ! Revision 1.4 2002/10/08 15:49:29 haselbac
768 ! {IO}STAT=global%error replaced by {IO}STAT=errorFlag - SGI problem
769 !
770 ! Revision 1.3 2002/09/28 14:33:24 jferry
771 ! added relative velocity to injection BC
772 !
773 ! Revision 1.2 2002/09/27 00:57:10 jblazek
774 ! Changed makefiles - no makelinks needed.
775 !
776 ! Revision 1.1 2002/09/25 18:31:57 jferry
777 ! Implemented TBCs for ROCFLU
778 !
779 ! ******************************************************************************
780 
781 
782 
783 
784 
785 
786 
CImg< T > & line(const unsigned int y0)
Get a line.
Definition: CImg.h:18421
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine readpatchsection(global, fileID, nvals, keys, vals, brbeg, brend, prbeg, prend, distrib, profType, fname, defined)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE USE ModDataTypes USE prend
const NT & n
subroutine rflu_readtbcsection(pRegion, tbcType)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE USE ModDataTypes USE prbeg
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE USE ModDataTypes USE nvals
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469