This subroutine load data for radial sponge layer Example:
spatial = {
--supported options: 'spongelayer_radial','sponge_radial_2d',
-- 'viscous_spongelayer_radial',
-- 'viscous_spongelayer_radial_2d'
predefined = 'viscous_spongelayer_radial',
origin = {0.0,0.0,0.0},
radius = 1.0, -- Sponge start
thickness = 0.3,
damp_profile = 'linear', --'exponential', 'polynomial_n5', 'polynomial_n6'
damp_factor = 0.5,
damp_exponent = 1.0,
target_state = {
Default: density, velocityX, velocityY, velocityZ and pressure
viscosity = 1e-3
}
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_spongeLayer_radial_type), | intent(out) | :: | me |
Radial spongeLayer data type |
||
type(flu_State) | :: | conf |
lua state type |
|||
integer, | intent(in) | :: | thandle |
aotus parent handle |
||
integer, | intent(in) | :: | nDim |
number of Dimension for nonViscous sponges |
||
integer, | intent(in) | :: | nComp |
Number of component of St-Fun variable under which this spatial function is defined |
||
character(len=*), | intent(in), | optional | :: | stateName |
Load stateName from target_state table |
subroutine tem_load_spongeLayer_radial(me, conf, thandle, nDim, nComp, &
& stateName)
! --------------------------------------------------------------------------
!> Radial spongeLayer data type
type(tem_spongeLayer_radial_type), intent(out) :: me
!> lua state type
type(flu_State) :: conf
!> aotus parent handle
integer, intent(in) :: thandle
!> number of Dimension for nonViscous sponges
integer, intent(in) :: nDim
!> Number of component of St-Fun variable under which this spatial function
!! is defined
integer, intent(in) :: nComp
!> Load stateName from target_state table
character(len=*), intent(in), optional :: stateName
! --------------------------------------------------------------------------
integer :: iError
integer :: vError(3), errfatal(3)
! --------------------------------------------------------------------------
errfatal = aotErr_Fatal
! Sponge origin
call aot_get_val( L = conf, &
& thandle = thandle, &
& key = 'origin', &
& val = me%origin, &
& ErrCode = vError )
if (any(btest( vError, errFatal )) ) then
write(logUnit(1),*) 'ERROR reading the sponge origin, ' &
& // 'origin is not well defined. ' &
& // 'It should have 3 entries for each coordinate.'
call tem_abort()
end if
write(logUnit(1),*) ' * Origin =', me%origin
! Sponge inner radius
call aot_get_val( L = conf, &
& thandle = thandle, &
& key = 'radius', &
& val = me%radius, &
& ErrCode = iError )
if (btest(iError, aotErr_Fatal)) then
write(logUnit(1),*) 'ERROR reading the sponge inner radius. '
call tem_abort()
end if
write(logUnit(1),*) ' * Inner radius =', me%radius
! Load base information required for sponge layer definition like
! damp_factor, damp_exponent and target_state
call load_spongeLayer( conf = conf, &
& thandle = thandle, &
& me = me%spongeLayer_base_type, &
& nDim = nDim, &
& nComp = nComp, &
& stateName = stateName )
end subroutine tem_load_spongeLayer_radial