tem_varSys_append_derVar Subroutine

public subroutine tem_varSys_append_derVar(me, varname, operType, nComponents, input_varname, input_varIndex, method_data, get_point, get_element, set_params, get_params, setup_indices, get_valOfIndex, pos, wasAdded)

Append a new (non-state) variable to the variable system.

Non-state variables might depend on other variables, and do not modify the internal state variable counter.

Arguments

Type IntentOptional Attributes Name
type(tem_varSys_type), intent(inout) :: me

Variable system to append the state variable to.

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

Variable to append to the state.

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

Operation type

integer, intent(in) :: nComponents

Number of components in this variable.

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

List of variable names, this variable depends on.

(Might be empty). The variable will only be appended, if all invars are already found in the variable system.

integer, intent(in), optional :: input_varIndex(:)

Component index to access from input variable nComponents

(Might be empty). The variable will only be appended, if all index are available in input_varname

type(c_ptr), intent(in) :: method_data

Data that is required by the methods to obtain the variable.

procedure(tem_varSys_proc_point), pointer :: get_point

Procedure which allows the retrieval of the variable at given points.

procedure(tem_varSys_proc_element), pointer :: get_element

Procedure which allows the retrieval of the variable in an element.

procedure(tem_varSys_proc_setParams), optional, pointer :: set_params

Procedure which allows to set parameter in method_data

procedure(tem_varSys_proc_getParams), optional, pointer :: get_params

Procedure which allows to get parameter in method_data

procedure(tem_varSys_proc_setupIndices), pointer :: setup_indices

Procedure to setup growing array of points, variable value in method_data and return index of points set

procedure(tem_varSys_proc_getValOfIndex), pointer :: get_valOfIndex

Procedure which allows to retrieval of the variable at point or val array index

integer, intent(out), optional :: pos

Position of the variable in the system.

logical, intent(out), optional :: wasAdded

Indicator, if the variable was actually added to the system.


Calls

proc~~tem_varsys_append_dervar~~CallsGraph proc~tem_varsys_append_dervar tem_varSys_append_derVar interface~append~29 append proc~tem_varsys_append_dervar->interface~append~29 interface~positionofval~5 positionofval proc~tem_varsys_append_dervar->interface~positionofval~5 proc~append_da_label append_da_label interface~append~29->proc~append_da_label proc~append_da_veclabel append_da_veclabel interface~append~29->proc~append_da_veclabel proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label interface~expand~27 expand proc~append_da_label->interface~expand~27 interface~sortedposofval~5 sortedposofval proc~append_da_label->interface~sortedposofval~5 proc~append_da_veclabel->interface~expand~27 proc~posofval_label->interface~sortedposofval~5 proc~expand_da_label expand_da_label interface~expand~27->proc~expand_da_label proc~sortposofval_label sortposofval_label interface~sortedposofval~5->proc~sortposofval_label

Called by

proc~~tem_varsys_append_dervar~~CalledByGraph proc~tem_varsys_append_dervar tem_varSys_append_derVar proc~tem_varsys_append_meshinfovar tem_varSys_append_meshInfoVar proc~tem_varsys_append_meshinfovar->proc~tem_varsys_append_dervar proc~tem_varsys_append_opervar tem_varSys_append_operVar proc~tem_varsys_append_opervar->proc~tem_varsys_append_dervar proc~tem_varsys_append_stfun_raw tem_varSys_append_stFun_raw proc~tem_varsys_append_stfun_raw->proc~tem_varsys_append_dervar proc~tem_varsys_append_stfunvar tem_varSys_append_stFunVar proc~tem_varsys_append_stfunvar->proc~tem_varsys_append_dervar interface~tem_varsys_append_stfun tem_varSys_append_stfun interface~tem_varsys_append_stfun->proc~tem_varsys_append_stfun_raw interface~tem_varsys_append_stfun->proc~tem_varsys_append_stfunvar proc~tem_varsys_append_luavar tem_varSys_append_luaVar proc~tem_varsys_append_luavar->proc~tem_varsys_append_opervar proc~tem_varsys_append_luavar->interface~tem_varsys_append_stfun proc~tem_variable_loadmapping_single tem_variable_loadMapping_single proc~tem_variable_loadmapping_single->interface~tem_varsys_append_stfun interface~tem_variable_loadmapping tem_variable_loadMapping interface~tem_variable_loadmapping->proc~tem_variable_loadmapping_single proc~tem_variable_loadmapping_vector tem_variable_loadMapping_vector proc~tem_variable_loadmapping_vector->proc~tem_variable_loadmapping_single

Source Code

  subroutine tem_varSys_append_derVar(me, varname, operType, nComponents,    &
    &                                 input_varname, input_varindex,         &
    &                                 method_data, get_point, get_element,   &
    &                                 set_params, get_params, setup_indices, &
    &                                 get_valOfIndex, pos, wasAdded          )
    ! --------------------------------------------------------------------------!
    !> Variable system to append the state variable to.
    type(tem_varSys_type), intent(inout) :: me

    !> Variable to append to the state.
    character(len=*), intent(in) :: varname

    !> Operation type
    character(len=*), optional, intent(in) :: operType

    !> Number of components in this variable.
    integer, intent(in) :: nComponents

    !> List of variable names, this variable depends on.
    !!
    !! (Might be empty). The variable will only be appended, if all invars are
    !! already found in the variable system.
    character(len=*), optional, intent(in) :: input_varname(:)

    !> Component index to access from input variable nComponents
    !!
    !! (Might be empty). The variable will only be appended, if all index are
    !! available in input_varname
    integer, optional, intent(in) :: input_varIndex(:)

    !> Data that is required by the methods to obtain the variable.
    type(c_ptr), intent(in) :: method_data

    !> Procedure which allows the retrieval of the variable at given points.
    procedure(tem_varSys_proc_point), pointer :: get_point

    !> Procedure which allows the retrieval of the variable in an element.
    procedure(tem_varSys_proc_element), pointer :: get_element

    !> Procedure which allows to set parameter in method_data
    procedure(tem_varSys_proc_setParams), optional, pointer :: set_params

    !> Procedure which allows to get parameter in method_data
    procedure(tem_varSys_proc_getParams), optional, pointer :: get_params

    !> Procedure to setup growing array of points, variable value
    !! in method_data and return index of points set
    procedure(tem_varSys_proc_setupIndices), pointer :: setup_indices

    !> Procedure which allows to retrieval of the variable at point or val
    !! array index
    procedure(tem_varSys_proc_getValOfIndex), pointer :: get_valOfIndex

    !> Position of the variable in the system.
    integer, optional, intent(out)   :: pos

    !> Indicator, if the variable was actually added to the system.
    logical,  optional, intent(out)   :: wasAdded
    ! --------------------------------------------------------------------------!
    logical :: newVar, nonexisting_input
    type(tem_varSys_op_type) :: var_access
    integer :: iIn, nInputs, varPos, in_nComp, nVarIndex
    integer, allocatable :: inpos(:), input_varIndex_loc(:)
    ! --------------------------------------------------------------------------!

    if (present(input_varname)) then
      nInputs = size(input_varname)
    else
      nInputs = 0
    end if

    if (present(input_varindex)) then
      nVarIndex = size(input_varindex)
    else
      nVarIndex = 0
    end if

    allocate(inpos(nInputs))
    allocate(input_varIndex_loc(nVarIndex))

    nonexisting_input = .false.

    ! Check whether all input variables are present. If not, the variable can't
    ! be added, because variables it depends on are not present.
    do iIn=1,nInputs
      inpos(iIn) = PositionOfVal(me%varname, input_varname(iIn))
      nonexisting_input = (inpos(iIn) <= 0)
      if (nonexisting_input) then
        write(logUnit(llerror),*) 'Variable ' // trim(varname) &
          & // ' cant be added due to missing input variable '   &
          & // trim(input_varname(iIn))
        EXIT
      end if

      ! allocate and set input_varIndex only if nVarIndex > 0
      ! input_varIndex is used only for one nInputs
      if ( nVarIndex > 0 .and. nInputs == 1 ) then
        ! Check whether the component index requested from input_varname
        ! is available
        in_nComp = me%method%val(inpos(iIn))%nComponents
        if ( (maxval(input_varIndex) > in_nComp) .or. &
          &  (nVarIndex > in_nComp) ) then
          write(logUnit(llerror),*) 'WARNING: Component index requested ' &
            & // 'is not extractable from variable: '  &
            & // trim(input_varname(iIn))
          write(logUnit(llerror),*) ' Input variable name: ' &
            & // trim(input_varname(iIn))
          write(logUnit(llerror),*) ' nComponents of input_varname: ', &
            & in_nComp
          write(logUnit(llerror),*) ' Requested component index of ' &
            & // 'input var: ', input_varIndex
          EXIT
        end if
        ! Copy input_varIndex
        input_varIndex_loc = input_varIndex
      end if

    end do

    if (nonexisting_input) then

      varPos = 0
      newVar = .false.

    else

      ! All inputs are found in the variable system, thus this derived
      ! quantity can be added now.
      call append( me       = me%varname, &
        &          val      = varname,    &
        &          pos      = varPos,     &
        &          wasAdded = newVar      )

      if (newVar) then
        write(logUnit(lldebug),*) 'Variable ' // trim(varname) &
          & // ' was newly added'
        var_access%mypos = varPos
        if (present(operType)) then
          var_access%operType = trim(operType)
        else
          var_access%operType = 'derived'
        end if
        var_access%nInputs = nInputs
        var_access%nComponents = nComponents
        allocate(var_access%state_varPos(0))
        allocate(var_access%auxField_varPos(0))
        allocate(var_access%input_varPos(nInputs))
        var_access%input_varPos = inpos
        allocate(var_access%input_varIndex(nVarIndex))
        var_access%input_varIndex = input_varIndex_loc

        var_access%method_data = method_data
        var_access%get_point => get_point
        var_access%get_element => get_element
        if (present(set_params)) then
          var_access%set_params => set_params
        else
          var_access%set_params => tem_varSys_setParams_dummy
        end if
        if (present(get_params)) then
          var_access%get_params => get_params
        else
          var_access%get_params => tem_varSys_getParams_dummy
        end if
        var_access%setup_indices => setup_indices
        var_access%get_valOfIndex => get_valOfIndex
        call append( me  = me%method,   &
          &          val = var_access )
      else
        write(logUnit(lldebug),*) 'Variable ' // trim(varname) &
          & // ' is already present'
      end if

    end if

    if (present(pos)) pos = varPos
    if (present(wasAdded)) wasAdded = newVar

  end subroutine tem_varSys_append_derVar