This subroutine load standard temporal function variables from LUA file.
Default values: \li min_factor = 0.0 \li max_factor = 1.0 \li from_time = 0 \li to_time = tmax/2
Valid definition: \li linear function
temporal = {predefined="linear", min_factor = 0.0, max_factor = 1.0,
from_time = 0, to_time = 1000}
Example: Transient inlet velocity which starts from 0 to 1000 with maximum value 0.08 is shown below for linear and smooth with definition and image.
boundary_condition =
{
{ label = 'inlet',
kind = 'inlet_ubb',
velocityX = { spatial = 1.0
,temporal = {predefined="linear", min_factor = 0.0,
max_factor = 1.0, from_time = 0, to_time = 1000}}
-- ,temporal = {predefined="smooth", min_factor = 0.0,
-- max_factor = 1.0, from_time = 0, to_time = 1000}}
velocityY = 0.0
velocityZ = 0.0
}
}
\image html transient.png
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_linear_type), | intent(inout) | :: | me |
temporal predefined fun type |
||
type(flu_State) | :: | conf |
lua state type |
|||
integer, | intent(in) | :: | thandle |
aotus parent handle |
subroutine load_temporal_linear( me, conf, thandle ) ! --------------------------------------------------------------------------- !> temporal predefined fun type type(tem_linear_type),intent(inout) :: me !> lua state type type(flu_State) :: conf !> aotus parent handle integer, intent(in) :: thandle ! --------------------------------------------------------------------------- ! local variables integer :: iError ! --------------------------------------------------------------------------- !min_factor call aot_get_val( L = conf, & & thandle = thandle, & & key = 'min_factor', & & val = me%min_factor, & & ErrCode = iError, & & default = 0.0_rk ) if (btest(iError, aoterr_Fatal)) then write(logUnit(1),*)'FATAL Error occured, while retrieving min_factor :' if (btest(iError, aoterr_NonExistent)) & & write(logUnit(1),*)'Variable not existent!' if (btest(iError, aoterr_WrongType)) & & write(logUnit(1),*)'Variable has wrong type!' write(logUnit(1),*)'STOPPING' call tem_abort() end if !max_factor call aot_get_val( L = conf, & & thandle = thandle, & & key = 'max_factor', & & val = me%max_factor, & & ErrCode = iError, & & default = 1.0_rk ) if (btest(iError, aoterr_Fatal)) then write(logUnit(1),*)'FATAL Error occured, while retrieving max_factor :' if (btest(iError, aoterr_NonExistent)) & & write(logUnit(1),*)'Variable not existent!' if (btest(iError, aoterr_WrongType)) & & write(logUnit(1),*)'Variable has wrong type!' write(logUnit(1),*)'STOPPING' call tem_abort() end if !from_time call aot_get_val( L = conf, & & thandle = thandle, & & key = 'from_time', & & val = me%from_time, & & ErrCode = iError, & & default = 0.0_rk ) if (btest(iError, aoterr_Fatal)) then write(logUnit(1),*)'FATAL Error occured, while retrieving from_time :' if (btest(iError, aoterr_NonExistent)) & & write(logUnit(1),*)'Variable not existent!' if (btest(iError, aoterr_WrongType)) & & write(logUnit(1),*)'Variable has wrong type!' write(logUnit(1),*)'STOPPING' call tem_abort() end if ! to_time call aot_get_val( L = conf, & & thandle = thandle, & & key = 'to_time', & & val = me%to_time, & & ErrCode = iError, & & default = 1.0_rk ) if (btest(iError, aoterr_Fatal)) then write(logUnit(1),*)'FATAL Error occured, while retrieving to_time :' if (btest(iError, aoterr_NonExistent)) & & write(logUnit(1),*)'Variable not existent!' if (btest(iError, aoterr_WrongType)) & & write(logUnit(1),*)'Variable has wrong type!' write(logUnit(1),*)'STOPPING' call tem_abort() end if write(logUnit(3), "(A, F10.5)") " min factor: ", me%min_factor write(logUnit(3), "(A, F10.5)") " max factor: ", me%max_factor write(logUnit(3), "(A, E10.5)") " from time: ", me%from_time write(logUnit(3), "(A, E10.5)") " to time: ", me%to_time end subroutine load_temporal_linear