tem_opVar_reduction_transient_init Subroutine

public subroutine tem_opVar_reduction_transient_init(varSys, tree, redTransVarMap, nDofs, time)

Initialize time reduction operation variable Loop over all variable in varSys and allocate redTrans%val for reduction_transient operation variable with nElems

Arguments

Type IntentOptional Attributes Name
type(tem_varSys_type), intent(in) :: varSys

Global variable system

type(treelmesh_type), intent(in) :: tree

treelmesh_type

type(tem_varMap_type), intent(out) :: redTransVarMap

position of time reduction variable in varSys

integer, intent(in), optional :: nDofs

Solver nDegrees of freedom

type(tem_time_type), intent(in) :: time

Current time


Calls

proc~~tem_opvar_reduction_transient_init~~CallsGraph proc~tem_opvar_reduction_transient_init tem_opVar_reduction_transient_init interface~append~24 append proc~tem_opvar_reduction_transient_init->interface~append~24 interface~init~24 init proc~tem_opvar_reduction_transient_init->interface~init~24 proc~tem_create_varmap tem_create_varMap proc~tem_opvar_reduction_transient_init->proc~tem_create_varmap proc~tem_reduction_transient_init tem_reduction_transient_init proc~tem_opvar_reduction_transient_init->proc~tem_reduction_transient_init proc~append_arrayga2d_real append_arrayga2d_real interface~append~24->proc~append_arrayga2d_real proc~append_singlega2d_real append_singlega2d_real interface~append~24->proc~append_singlega2d_real proc~init_ga2d_real init_ga2d_real interface~init~24->proc~init_ga2d_real interface~append~4 append proc~tem_create_varmap->interface~append~4 interface~init~29 init proc~tem_create_varmap->interface~init~29 interface~positionofval~5 positionofval proc~tem_create_varmap->interface~positionofval~5 interface~truncate~23 truncate proc~tem_create_varmap->interface~truncate~23 proc~tem_abort tem_abort proc~tem_reduction_transient_init->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~init_da_label init_da_label interface~init~29->proc~init_da_label proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~truncate_da_label truncate_da_label interface~truncate~23->proc~truncate_da_label interface~expand~22 expand proc~append_arrayga2d_real->interface~expand~22 proc~append_singlega2d_real->interface~expand~22 mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~expand_ga2d_real expand_ga2d_real interface~expand~22->proc~expand_ga2d_real interface~expand~3 expand proc~append_ga_stringkeyvaluepair->interface~expand~3 proc~append_ga_stringkeyvaluepair_vec->interface~expand~3 interface~sortedposofval~5 sortedposofval proc~posofval_label->interface~sortedposofval~5

Source Code

  subroutine tem_opVar_reduction_transient_init(varSys, tree, redTransVarMap,&
    &                                           nDofs, time)
    ! -------------------------------------------------------------------------
    !> Global variable system
    type(tem_varSys_type), intent(in) :: varSys
    !> treelmesh_type
    type(treelmesh_type), intent(in) :: tree
    !> position of time reduction variable in varSys
    type(tem_varMap_type), intent(out) :: redTransVarMap
    !> Solver nDegrees of freedom
    integer, intent(in), optional :: nDofs
    !> Current time
    type(tem_time_type), intent(in) :: time
    ! -------------------------------------------------------------------------
    type(tem_varSys_op_data_type), pointer :: opData
    integer :: iVar, iElem, varPos, nDofs_loc, posDepVar, nCompMax, idxMax
    integer :: nRedVars
    type(grw_labelArray_type) :: redTransVarName
    integer :: elemPos(tree%nElems)
    real(kind=rk), allocatable :: input_varRes(:)
    ! -------------------------------------------------------------------------
    if (present(nDofs)) then
      nDofs_loc = nDofs
    else
      nDofs_loc = 1
    end if

    ! Gather list of variable names which has reduction_transient operation
    call init(redTransVarName)
    do iVar = 1, varSys%varName%nVals
      if (trim(varSys%method%val(iVar)%operType) == 'reduction_transient') then
        call append(me=redTransVarName, val=varSys%varName%val(iVar))
      end if
    end do

    ! create varMap to store position of reduction_transient variable in varSys
    call tem_create_varMap(                                                 &
      &             varName = redTransVarName%val(1:redTransVarName%nVals), &
      &             varSys  = varSys,                                       &
      &             varMap  = redTransVarMap                                )

    elemPos(1:tree%nElems) = (/ (iElem, iElem=1, tree%nElems) /)
    nRedVars =  redTransVarMap%varPos%nVals
    nCompMax = maxval(varSys%method%val(redTransVarMap%varPos%val(1:nRedVars)) &
      &                            %nComponents)
    allocate(input_varRes(tree%nElems*nCompMax*nDofs_loc))

    ! Initialize time reduction
    do iVar = 1, redTransVarMap%varPos%nVals
      varPos = redTransVarmap%varPos%val(iVar)
      call C_F_POINTER(varSys%method%val(varPos)%method_data, opData)

      call tem_reduction_transient_init(                            &
        &      me          = opData%redTrans,                       &
        &      nElems      = tree%nElems,                           &
        &      nComponents = varSys%method%val(varPos)%nComponents, &
        &      nDofs       = nDofs_loc                              )
      ! Fill last
      posDepVar = varSys%method%val(varPos)%input_varPos(1)
      call varSys%method%val(posDepVar)%get_element(                        &
        &                                   varSys  = varSys,               &
        &                                   elemPos = elemPos,              &
        &                                   time    = time,                 &
        &                                   tree    = tree,                 &
        &                                   nElems  = tree%nElems,          &
        &                                   nDofs   = nDofs_loc,            &
        &                                   res     = input_varRes(:)       )

      idxMax = opData%redTrans%nEntries
      opData%redTrans%val(:, opData%redTrans%last) = input_varRes(1:idxMax)
    end do

  end subroutine tem_opVar_reduction_transient_init