tem_load_bc_state Subroutine

public subroutine tem_load_bc_state(bc, state_name, nComp, style, conf, bc_handle, varDict, varSys, solverData_evalElem, ErrCode)

Function to load spatial, constant or temporal boundary conditions to the boundary state type (or combinations of them). or Valid definitions: Example given for state variable velocityX in variable table

 variable = {
   { name = 'bc_pressure',
     ncomponents = 1,
     var_type = 'st_fun',
     st_fun = 0.01
   },
   { name = 'bc_velocity',
     ncomponents = 3,
     var_type = 'st_fun',
     st_fun = {0.01,0.0,0.0}
   },
   { name = 'bc_velocity_2',
     ncomponents = 3,
     var_type = 'st_fun',
     st_fun = {
                predefined = 'combined',
                temporal = {predefined="linear", min_factor = 0.0,
                             max_factor = 1.0, from_time = 0, to_time = 1000},
                spatial = {predefined='parabol', shape = {kind = 'canoND',
                           object = { origin={-2.0,0.0,0.0},vec={0.0,1.0,0.0}}},
                           amplitude = {0.0,1.0,0.0}
                          }
                 }
       }
 }
 boundary_condition = {
  {
    label = 'from_seeder',
    kind = 'bc_kind',
    style = 'dirichlet',
    pressure = 'bc_pressure',
    velocity = 'bc_velocity' --'bc_velocity_2'
  }
 }

Arguments

Type IntentOptional Attributes Name
type(tem_bc_state_type), intent(inout) :: bc

The boundary to fill

character(len=*), intent(in) :: state_name

The state variable to set with this boundary condition

integer, intent(in), optional :: nComp

Number of Components in this boundary variable.

character(len=*), intent(in), optional :: style

Style of this boundary condition dirichlet = set value itself neumann = set derivative of value

type(flu_State), intent(in) :: conf
integer, intent(in) :: bc_handle

Handle to the table describing the boundary

type(grw_stringkeyvaluepairarray_type), intent(inout) :: varDict

The dictionary that contains the mapping between expected variables and the actual variables defined by the user.

type(tem_varSys_type), intent(inout) :: varSys
type(tem_varSys_solverData_evalElem_type), intent(in), optional :: solverData_evalElem

A routine that allows the construction of an element representation from a point values.

integer, intent(out), optional :: ErrCode

Error code


Calls

proc~~tem_load_bc_state~~CallsGraph proc~tem_load_bc_state tem_load_bc_state aot_exists aot_exists proc~tem_load_bc_state->aot_exists aot_get_val aot_get_val proc~tem_load_bc_state->aot_get_val aot_table_close aot_table_close proc~tem_load_bc_state->aot_table_close aot_table_open aot_table_open proc~tem_load_bc_state->aot_table_open interface~append~4 append proc~tem_load_bc_state->interface~append~4 interface~tem_variable_loadmapping tem_variable_loadMapping proc~tem_load_bc_state->interface~tem_variable_loadmapping proc~tem_abort tem_abort proc~tem_load_bc_state->proc~tem_abort proc~append_ga_stringkeyvaluepair append_ga_stringkeyvaluepair interface~append~4->proc~append_ga_stringkeyvaluepair proc~append_ga_stringkeyvaluepair_vec append_ga_stringkeyvaluepair_vec interface~append~4->proc~append_ga_stringkeyvaluepair_vec proc~tem_variable_loadmapping_single tem_variable_loadMapping_single interface~tem_variable_loadmapping->proc~tem_variable_loadmapping_single proc~tem_variable_loadmapping_vector tem_variable_loadMapping_vector interface~tem_variable_loadmapping->proc~tem_variable_loadmapping_vector mpi_abort mpi_abort proc~tem_abort->mpi_abort interface~expand~3 expand proc~append_ga_stringkeyvaluepair->interface~expand~3 proc~append_ga_stringkeyvaluepair_vec->interface~expand~3 proc~tem_variable_loadmapping_single->aot_exists proc~tem_variable_loadmapping_single->aot_get_val proc~tem_variable_loadmapping_single->interface~append~4 interface~tem_load_spacetime tem_load_spacetime proc~tem_variable_loadmapping_single->interface~tem_load_spacetime interface~tem_varsys_append_stfun tem_varSys_append_stfun proc~tem_variable_loadmapping_single->interface~tem_varsys_append_stfun proc~tem_spacetime_hash_id tem_spacetime_hash_id proc~tem_variable_loadmapping_single->proc~tem_spacetime_hash_id proc~tem_variable_loadmapping_vector->aot_table_close proc~tem_variable_loadmapping_vector->aot_table_open proc~tem_variable_loadmapping_vector->proc~tem_variable_loadmapping_single proc~expand_ga_stringkeyvaluepair expand_ga_stringkeyvaluepair interface~expand~3->proc~expand_ga_stringkeyvaluepair proc~tem_load_spacetime_single tem_load_spacetime_single interface~tem_load_spacetime->proc~tem_load_spacetime_single proc~tem_load_spacetime_table tem_load_spacetime_table interface~tem_load_spacetime->proc~tem_load_spacetime_table proc~tem_varsys_append_stfun_raw tem_varSys_append_stFun_raw interface~tem_varsys_append_stfun->proc~tem_varsys_append_stfun_raw proc~tem_varsys_append_stfunvar tem_varSys_append_stFunVar interface~tem_varsys_append_stfun->proc~tem_varsys_append_stfunvar aot_fun_close aot_fun_close proc~tem_spacetime_hash_id->aot_fun_close aot_fun_id aot_fun_id proc~tem_spacetime_hash_id->aot_fun_id aot_fun_open aot_fun_open proc~tem_spacetime_hash_id->aot_fun_open

Source Code

  subroutine tem_load_bc_state( bc, state_name, nComp, style, conf, bc_handle, &
    &                           varDict, varSys, solverData_evalElem, ErrCode        )
    !---------------------------------------------------------------------------
    !> The boundary to fill
    type(tem_bc_state_type), intent(inout) :: bc
    !> The state variable to set with this boundary condition
    character(len=*), intent(in) :: state_name
    !> Number of Components in this boundary variable.
    integer, intent(in), optional :: nComp
    !> Style of this boundary condition
    !! dirichlet = set value itself
    !! neumann = set derivative of value
    character(len=*), optional, intent(in) :: style
    type(flu_State),intent(in) :: conf !< Lua state
    !> Handle to the table describing the boundary
    integer, intent(in) :: bc_handle
    !> The dictionary that contains the mapping between expected variables
    !! and the actual variables defined by the user.
    type(grw_stringKeyValuePairArray_type), intent(inout) :: varDict
    type(tem_varSys_type), intent(inout) :: varSys
    !> A routine that allows the construction of an element representation
    !! from a point values.
    type(tem_varSys_solverData_evalElem_type), &
      &  optional, intent(in) :: solverData_evalElem
    !> Error code
    integer, optional, intent(out) :: ErrCode
    ! ---------------------------------------------------------------------------
    integer :: iError
    integer :: state_handle
    character(len=labelLen) :: def_style, varname
    logical :: varExist
    type(tem_stringKeyValuePair_type) :: kvp
    real(kind=rk) :: numtest
    ! ---------------------------------------------------------------------------
    ! Assume undefined boundary condition for this state.
    bc%isDefined = .false.
    bc%state_name = trim(state_name)

    if (present(style)) then
      def_style = trim(style)
    else
      def_style = 'dirichlet'
    end if

    if (present(nComp)) then
      bc%nComponents = nComp
    else
      bc%nComponents = 1
    end if

    write(logUnit(1),*)' Loading bc state for '// trim(state_name)

    ! If table is defined load spacetime function variable name from table
    ! else load from bc_handle
    call aot_table_open( L       = conf,         &
      &                  parent  = bc_handle,    &
      &                  thandle = state_handle, &
      &                  key     = state_name    )

    varExist = .false.
    ! If boundary variable is defined as a table, then load boundary
    ! style and check if varname exist.
    if (state_handle /= 0) then
      ! Found a table, now check for the style of this boundary condition
      ! Default to dirichlet.
      call aot_get_val( L       = conf,         &
        &               thandle = state_handle, &
        &               val     = bc%style,     &
        &               ErrCode = iError,       &
        &               key     = 'style',      &
        &               default = def_style     )

      ! try to load variable as string
      varExist = aot_exists( L       = conf,         &
        &                    thandle = state_handle, &
        &                    key     = 'varname'     )
    end if ! variable defined as table

    !If varname exist load refered variable name from 'varname'
    ! and append the value to dictionary
    if (varExist) then
      bc%isDefined = .true.

      ! First check, whether this variable definition is a number.
      ! (They also satisfy strings).
      ! We do not accept numbers as variable names, instead this
      ! will be read as constant stfun.
      call aot_get_val( L       = conf,         &
        &               thandle = state_handle, &
        &               key     = 'varname',    &
        &               val     = numtest,      &
        &               ErrCode = iError        )
      if (btest(iError, aoterr_WrongType)) then
        ! Not a number, try to interpret it as a string.
        call aot_get_val( L       = conf,         &
           &              thandle = state_handle, &
           &              key     = 'varname',    &
           &              val     = varname,      &
           &              ErrCode = iError        )
        if( iError == 0 ) then
          ! Found a string, use it to refer to a variable.
          write(logUnit(3),*) 'Corresponding variable for ' &
            & // trim(state_name) // ' found: ' // trim(varname)
          kvp%key = state_name
          kvp%value = varname
          call append( me = varDict, val = kvp )
        else
          write(logUnit(1),*) 'Error: Unable to load state name "'// &
            & trim(state_name)//'" for boundary "'
          write(logUnit(1),*) '"'//trim(state_name)//'" is defined as table with'
          write(logUnit(1),*) 'varname, but failed to load variable name '// &
            &                 'as string'
          call tem_abort()
        end if
      end if
    else
      bc%isDefined = .true.
      bc%style = def_style

      ! Load spacetime function variable name
      call tem_variable_loadMapping( expectedName        = trim(state_name),   &
        &                            conf                = conf,               &
        &                            thandle             = bc_handle,          &
        &                            varDict             = varDict,            &
        &                            varSys              = varSys,             &
        &                            nComp               = bc%nComponents,     &
        &                            solverData_evalElem = solverData_evalElem,&
        &                            ErrCode             = iError              )

      if (iError /= 0) then
        write(logUnit(1),*) 'Error: Unable to load state name "'// &
          & trim(state_name)//'" for boundary '
        call tem_abort()
      end if

    end if

    call aot_table_close(L=conf, thandle=state_handle)

    if (present(ErrCode)) ErrCode = iError

  end subroutine tem_load_bc_state