get_active_region_info Subroutine

private subroutine get_active_region_info(self, unique_ids, ierr)

Type Bound

type_vtk

Arguments

Type IntentOptional Attributes Name
class(type_vtk), intent(in) :: self

VTK data

integer(kind=int32), intent(inout), allocatable :: unique_ids(:)
integer(kind=int32), intent(out) :: ierr

Calls

proc~~get_active_region_info~~CallsGraph proc~get_active_region_info type_vtk%get_active_region_info interface~unique unique proc~get_active_region_info->interface~unique proc~unique_int16 unique_int16 interface~unique->proc~unique_int16 proc~unique_int32 unique_int32 interface~unique->proc~unique_int32 proc~unique_int64 unique_int64 interface~unique->proc~unique_int64 proc~unique_int8 unique_int8 interface~unique->proc~unique_int8 interface~allocate_array allocate_array proc~unique_int16->interface~allocate_array sort sort proc~unique_int16->sort proc~unique_int32->interface~allocate_array proc~unique_int32->sort proc~unique_int64->interface~allocate_array proc~unique_int64->sort proc~unique_int8->interface~allocate_array proc~unique_int8->sort proc~allocate_rank1_int16 allocate_rank1_int16 interface~allocate_array->proc~allocate_rank1_int16 proc~allocate_rank1_int32 allocate_rank1_int32 interface~allocate_array->proc~allocate_rank1_int32 proc~allocate_rank1_int64 allocate_rank1_int64 interface~allocate_array->proc~allocate_rank1_int64 proc~allocate_rank1_int8 allocate_rank1_int8 interface~allocate_array->proc~allocate_rank1_int8 proc~allocate_rank1_logical1 allocate_rank1_logical1 interface~allocate_array->proc~allocate_rank1_logical1 proc~allocate_rank1_logical4 allocate_rank1_logical4 interface~allocate_array->proc~allocate_rank1_logical4 proc~allocate_rank1_logical8 allocate_rank1_logical8 interface~allocate_array->proc~allocate_rank1_logical8 proc~allocate_rank1_real128 allocate_rank1_real128 interface~allocate_array->proc~allocate_rank1_real128 proc~allocate_rank1_real32 allocate_rank1_real32 interface~allocate_array->proc~allocate_rank1_real32 proc~allocate_rank1_real64 allocate_rank1_real64 interface~allocate_array->proc~allocate_rank1_real64 proc~allocate_rank2_int16 allocate_rank2_int16 interface~allocate_array->proc~allocate_rank2_int16 proc~allocate_rank2_int32 allocate_rank2_int32 interface~allocate_array->proc~allocate_rank2_int32 proc~allocate_rank2_int64 allocate_rank2_int64 interface~allocate_array->proc~allocate_rank2_int64 proc~allocate_rank2_int8 allocate_rank2_int8 interface~allocate_array->proc~allocate_rank2_int8 proc~allocate_rank2_logical1 allocate_rank2_logical1 interface~allocate_array->proc~allocate_rank2_logical1 proc~allocate_rank2_logical4 allocate_rank2_logical4 interface~allocate_array->proc~allocate_rank2_logical4 proc~allocate_rank2_logical8 allocate_rank2_logical8 interface~allocate_array->proc~allocate_rank2_logical8 proc~allocate_rank2_real128 allocate_rank2_real128 interface~allocate_array->proc~allocate_rank2_real128 proc~allocate_rank2_real32 allocate_rank2_real32 interface~allocate_array->proc~allocate_rank2_real32 proc~allocate_rank2_real64 allocate_rank2_real64 interface~allocate_array->proc~allocate_rank2_real64 proc~error_message error_message proc~allocate_rank1_int16->proc~error_message proc~allocate_rank1_int32->proc~error_message proc~allocate_rank1_int64->proc~error_message proc~allocate_rank1_int8->proc~error_message proc~allocate_rank1_logical1->proc~error_message proc~allocate_rank1_logical4->proc~error_message proc~allocate_rank1_logical8->proc~error_message proc~allocate_rank1_real128->proc~error_message proc~allocate_rank1_real32->proc~error_message proc~allocate_rank1_real64->proc~error_message proc~allocate_rank2_int16->proc~error_message proc~allocate_rank2_int32->proc~error_message proc~allocate_rank2_int64->proc~error_message proc~allocate_rank2_int8->proc~error_message proc~allocate_rank2_logical1->proc~error_message proc~allocate_rank2_logical4->proc~error_message proc~allocate_rank2_logical8->proc~error_message proc~allocate_rank2_real128->proc~error_message proc~allocate_rank2_real32->proc~error_message proc~allocate_rank2_real64->proc~error_message log_error log_error proc~error_message->log_error

Called by

proc~~get_active_region_info~~CalledByGraph proc~get_active_region_info type_vtk%get_active_region_info proc~initialize type_material_manager%initialize proc~initialize->proc~get_active_region_info proc~type_proereties_manager_initialize type_proereties_manager%type_proereties_manager_initialize proc~type_proereties_manager_initialize->proc~initialize

Source Code

    subroutine get_active_region_info(self, unique_ids, ierr)
        ! --- 引数 ---
        implicit none
        class(Type_VTK), intent(in) :: self !! VTK data
        integer(int32), allocatable, intent(inout) :: unique_ids(:)
        integer(int32), intent(out) :: ierr

        ! --- ローカル変数 ---
        integer(int32) :: max_dim
        integer(int32), allocatable :: collected_ids(:)
        integer(int32) :: i_cell, count
        logical(4) :: is_max_dim_element

        max_dim = 0
        ierr = 0

        ! --- ステップ1: メッシュ内の最大次元を判定 ---
        do i_cell = 1, self%num_total_cells
            select case (self%CELLS(i_cell)%cell_type)
            case (VTK_TETRA, VTK_HEXAHEDRON, &
                  VTK_WEDGE, VTK_PYRAMID, &
                  VTK_QUADRATIC_TETRA, VTK_QUADRATIC_HEXAHEDRON)
                max_dim = 3
                exit ! 3Dが見つかったら、それ以上探す必要はない
            case (VTK_TRIANGLE, VTK_PIXEL, &
                  VTK_QUAD, VTK_QUADRATIC_TRIANGLE, &
                  VTK_QUADRATIC_QUAD)
                max_dim = max(max_dim, 2)
            case (VTK_LINE, VTK_QUADRATIC_EDGE)
                max_dim = max(max_dim, 1)
            end select
        end do

        if (max_dim == 0) then
            ierr = -1
            return ! アクティブな要素がない
        end if

        ! --- ステップ2: 最大次元を持つ要素から、すべてのCellEntityIdを収集 ---
        allocate (collected_ids(self%num_total_cells))
        count = 0
        do i_cell = 1, self%num_total_cells
            is_max_dim_element = .false.
            select case (self%CELLS(i_cell)%cell_type)
            case (VTK_TETRA, VTK_HEXAHEDRON, &
                  VTK_WEDGE, VTK_PYRAMID, &
                  VTK_QUADRATIC_TETRA, VTK_QUADRATIC_HEXAHEDRON)
                if (max_dim == 3) is_max_dim_element = .true.
            case (VTK_TRIANGLE, VTK_PIXEL, &
                  VTK_QUAD, VTK_QUADRATIC_TRIANGLE, &
                  VTK_QUADRATIC_QUAD)
                if (max_dim == 2) is_max_dim_element = .true.
            case (VTK_LINE, VTK_QUADRATIC_EDGE)
                if (max_dim == 1) is_max_dim_element = .true.
            end select

            if (is_max_dim_element) then
                count = count + 1
                collected_ids(count) = self%CELLS(i_cell)%cell_entity_id
            end if
        end do

        ! --- ステップ3: 収集したIDリストから、ユニークなものだけを抽出 ---
        ! (これはFortranの標準的なユニーク化のアルゴリズム)
        if (count > 0) then
            call unique(collected_ids(1:count), unique_ids)
        else
            allocate (unique_ids(0))
        end if

    end subroutine get_active_region_info