tem_shape_findElemByBCLabels Subroutine

private subroutine tem_shape_findElemByBCLabels(bcLabels, bc_prop, foundAny, map2global, bcIDs, stencil)

This routine identify elements that belong to certain bounaries. Labels of required boundaries are given by bcLabels. bc_prop contains boudnary_ID of all local elements. Firstly, bcLabels are converted into bcIDs. Then all elements in bc_prop are looped over to check if it matches one of the required bcID. If match, its position is save in map2global. Number of elements found on each level is saved in countElems.

Arguments

Type IntentOptional AttributesName
character(len=labelLen), intent(in) :: bcLabels(:)

bcLabels

type(tem_BC_prop_type), intent(in) :: bc_prop

bc property

logical, intent(out) :: foundAny

if any element be identified

type(dyn_intarray_type), intent(inout) :: map2global

dynamic array. Elements positions in bc_prop%property%elemID

integer, intent(out), allocatable:: bcIDs(:)

id of boundary condition to be tracked

type(tem_stencilHeader_type), intent(in), optional :: stencil

stencil required to get useful links


Calls

proc~~tem_shape_findelembybclabels~~CallsGraph proc~tem_shape_findelembybclabels tem_shape_findElemByBCLabels interface~append~25 append proc~tem_shape_findelembybclabels->interface~append~25 proc~append_arrayga2d_real append_arrayga2d_real interface~append~25->proc~append_arrayga2d_real proc~append_singlega2d_real append_singlega2d_real interface~append~25->proc~append_singlega2d_real interface~expand~23 expand proc~append_arrayga2d_real->interface~expand~23 proc~append_singlega2d_real->interface~expand~23 proc~expand_ga2d_real expand_ga2d_real interface~expand~23->proc~expand_ga2d_real

Called by

proc~~tem_shape_findelembybclabels~~CalledByGraph proc~tem_shape_findelembybclabels tem_shape_findElemByBCLabels proc~tem_create_subtree_of tem_create_subTree_of proc~tem_create_subtree_of->proc~tem_shape_findelembybclabels proc~tem_create_subtree_of_st_funlist tem_create_subTree_of_st_funList proc~tem_create_subtree_of_st_funlist->proc~tem_create_subtree_of proc~tem_write_debugmesh tem_write_debugMesh proc~tem_write_debugmesh->proc~tem_create_subtree_of program~tem_varsys_test tem_varSys_test program~tem_varsys_test->proc~tem_create_subtree_of proc~tem_init_tracker_subtree tem_init_tracker_subTree proc~tem_init_tracker_subtree->proc~tem_create_subtree_of proc~tem_init_convergence tem_init_convergence proc~tem_init_convergence->proc~tem_create_subtree_of program~tem_variable_evaltype_test tem_variable_evaltype_test program~tem_variable_evaltype_test->proc~tem_create_subtree_of_st_funlist program~tem_varsys_stfunvar_test tem_varSys_stfunVar_test program~tem_varsys_stfunvar_test->proc~tem_create_subtree_of_st_funlist program~tem_variable_extract_test tem_variable_extract_test program~tem_variable_extract_test->proc~tem_create_subtree_of_st_funlist proc~check_variableoperations check_variableOperations proc~check_variableoperations->proc~tem_create_subtree_of_st_funlist program~tem_variable_combine_test tem_variable_combine_Test program~tem_variable_combine_test->proc~tem_create_subtree_of_st_funlist program~tem_varsys_opvar_test tem_varSys_opVar_test program~tem_varsys_opvar_test->proc~tem_create_subtree_of_st_funlist program~tem_logical_opertor_test tem_logical_opertor_test program~tem_logical_opertor_test->proc~check_variableoperations

Contents


Source Code

  subroutine tem_shape_findElemByBCLabels( bcLabels, bc_prop, &
    &                                      foundAny, map2global, bcIDs, stencil)
    ! ---------------------------------------------------------------------------
    !> bcLabels
    character(len=labelLen), intent(in) :: bcLabels(:)
    !> bc property
    type( tem_bc_prop_type ), intent(in) :: bc_prop
    !> if any element be identified
    logical,                 intent(out) :: foundAny
    !> dynamic array. Elements positions in bc_prop%property%elemID
    type(dyn_intArray_type), intent(inout) :: map2global
    !> id of boundary condition to be tracked
    integer, allocatable, intent(out) :: bcIDs(:)
    !> stencil required to get useful links
    type( tem_stencilHeader_type ), optional, intent(in) :: stencil
    ! ---------------------------------------------------------------------------
    integer :: dPos, iElem, iBC, nBCs, iBCtype, posInTree, QQN
    logical :: wasAdded, found
    integer, allocatable :: map(:)
    ! ---------------------------------------------------------------------------

    if( present( stencil ) ) then
      allocate( map( size(stencil%map) ) )
      map = stencil%map
      qqn = stencil%QQN
    else
      allocate( map(qQQQ) )
      map = qDir
      qqn = qQQQ
    end if

! write(*,*) 'inside tem_shape_findElemByBCLabels routine'

    foundAny = .false.
    ! convert bcLabels to bcIDs
    nBCs = size(bcLabels)
    allocate( bcIDs( nBCs ) )

! write(*,*) 'size of bcLabels: ', nBCs

    ! loop over all required bcLabels
    do iBC = 1, nBCs

      found = .false.

      ! loop over all BCs in the mesh
      do iBCtype = 1, bc_prop%nBCtypes
        if ( trim(bcLabels(iBC)) == bc_prop%bc_label( iBCtype ) ) then
          bcIDs( iBC ) = iBCtype
          found = .true.
! write(*,*) 'required bc label: '//trim(bcLabels(iBC))
! write(*,*) 'its bcID: ', bcIDs(iBC)
          exit
        end if
      end do

      if ( .not. found ) then
        write(logUnit(1),*) 'Required BC label: '//trim(bcLabels(iBC))
        write(logUnit(1),*) 'can not be found in given mesh!'
        stop
      end if
    end do ! iBC
! write(*,*) 'finished find bcIDs'
! write(*,*) 'nElems of bc_prop%property:', bc_prop%property%nElems

    ! Loop over all element with boundary property
    do iElem = 1, bc_prop%property%nElems

      do iBC = 1, nBCs
        if ( any(int(bc_prop%boundary_ID( map(:qqn), iElem )) == bcIDs(iBC) ) ) then

          posInTree = bc_prop%property%elemID(iElem)
          ! Append to treeID list (note that already existing ones are
          ! omitted)
          call append( me       = map2global,   &
            &          pos      = dPos,         &
            &          val      = posInTree,    &
            &          wasAdded = wasAdded )

          ! Count up if it was added
          if( wasAdded ) then
            ! tLevel = tem_levelOf( inTree%treeID(posInTree) )
            ! countElems( tLevel ) = countElems( tLevel ) + 1
            foundAny = .true.

! write(*,*) 'Found a element on level: ', tLevel
            exit ! continue with next element
          end if ! wasAdded

        end if ! boundary_ID == bcIDs
      end do

    end do ! iElem

  end subroutine tem_shape_findElemByBCLabels