domain_manager.F90 Source File


This file depends on

sourcefile~~domain_manager.f90~~EfferentGraph sourcefile~domain_manager.f90 domain_manager.F90 sourcefile~adjacency.f90 adjacency.F90 sourcefile~domain_manager.f90->sourcefile~adjacency.f90 sourcefile~core.f90 core.F90 sourcefile~domain_manager.f90->sourcefile~core.f90 sourcefile~element.f90 element.F90 sourcefile~domain_manager.f90->sourcefile~element.f90 sourcefile~element_factory.f90 element_factory.F90 sourcefile~domain_manager.f90->sourcefile~element_factory.f90 sourcefile~input.f90 input.F90 sourcefile~domain_manager.f90->sourcefile~input.f90 sourcefile~multicoloring.f90 multicoloring.F90 sourcefile~domain_manager.f90->sourcefile~multicoloring.f90 sourcefile~reordering.f90 reordering.F90 sourcefile~domain_manager.f90->sourcefile~reordering.f90 sourcefile~side.f90 side.F90 sourcefile~domain_manager.f90->sourcefile~side.f90 sourcefile~side_factory.f90 side_factory.F90 sourcefile~domain_manager.f90->sourcefile~side_factory.f90 sourcefile~adjacency_element.f90 adjacency_element.F90 sourcefile~adjacency.f90->sourcefile~adjacency_element.f90 sourcefile~adjacency_node.f90 adjacency_node.F90 sourcefile~adjacency.f90->sourcefile~adjacency_node.f90 sourcefile~allocate.f90 allocate.F90 sourcefile~core.f90->sourcefile~allocate.f90 sourcefile~check_range.f90 check_range.F90 sourcefile~core.f90->sourcefile~check_range.f90 sourcefile~deallocate.f90 deallocate.F90 sourcefile~core.f90->sourcefile~deallocate.f90 sourcefile~error.f90 error.F90 sourcefile~core.f90->sourcefile~error.f90 sourcefile~fortran_utils.f90 fortran_utils.F90 sourcefile~core.f90->sourcefile~fortran_utils.f90 sourcefile~string_utils.f90 string_utils.F90 sourcefile~core.f90->sourcefile~string_utils.f90 sourcefile~types.f90 types.F90 sourcefile~core.f90->sourcefile~types.f90 sourcefile~unique.f90 unique.F90 sourcefile~core.f90->sourcefile~unique.f90 sourcefile~vtk.f90 vtk.F90 sourcefile~core.f90->sourcefile~vtk.f90 sourcefile~vtk_constants.f90 vtk_constants.F90 sourcefile~core.f90->sourcefile~vtk_constants.f90 sourcefile~element.f90->sourcefile~core.f90 sourcefile~element.f90->sourcefile~input.f90 sourcefile~element_factory.f90->sourcefile~core.f90 sourcefile~element_factory.f90->sourcefile~element.f90 sourcefile~element_factory.f90->sourcefile~input.f90 sourcefile~input_interface.f90 input_interface.F90 sourcefile~input.f90->sourcefile~input_interface.f90 sourcefile~multicoloring.f90->sourcefile~core.f90 sourcefile~multicoloring.f90->sourcefile~adjacency_element.f90 sourcefile~reordering.f90->sourcefile~core.f90 sourcefile~reordering.f90->sourcefile~element.f90 sourcefile~reordering.f90->sourcefile~adjacency_node.f90 sourcefile~side.f90->sourcefile~core.f90 sourcefile~side.f90->sourcefile~input.f90 sourcefile~side_factory.f90->sourcefile~core.f90 sourcefile~side_factory.f90->sourcefile~input.f90 sourcefile~side_factory.f90->sourcefile~side.f90 sourcefile~adjacency_element.f90->sourcefile~core.f90 sourcefile~adjacency_element.f90->sourcefile~element.f90 sourcefile~adjacency_node.f90->sourcefile~core.f90 sourcefile~allocate.f90->sourcefile~error.f90 sourcefile~deallocate.f90->sourcefile~error.f90 sourcefile~memory_stats_wrapper.f90 memory_stats_wrapper.F90 sourcefile~fortran_utils.f90->sourcefile~memory_stats_wrapper.f90 sourcefile~signal_flag_wrapper.f90 signal_flag_wrapper.F90 sourcefile~fortran_utils.f90->sourcefile~signal_flag_wrapper.f90 sourcefile~system_info_wrapper.f90 system_info_wrapper.F90 sourcefile~fortran_utils.f90->sourcefile~system_info_wrapper.f90 sourcefile~input_interface.f90->sourcefile~core.f90 sourcefile~project_settings.f90 project_settings.F90 sourcefile~input_interface.f90->sourcefile~project_settings.f90 sourcefile~string_utils.f90->sourcefile~allocate.f90 sourcefile~array.f90 array.F90 sourcefile~types.f90->sourcefile~array.f90 sourcefile~gauss.f90 gauss.F90 sourcefile~types.f90->sourcefile~gauss.f90 sourcefile~pointer.f90 pointer.F90 sourcefile~types.f90->sourcefile~pointer.f90 sourcefile~variable.f90 variable.F90 sourcefile~types.f90->sourcefile~variable.f90 sourcefile~vector.f90 vector.F90 sourcefile~types.f90->sourcefile~vector.f90 sourcefile~unique.f90->sourcefile~allocate.f90 sourcefile~vtk.f90->sourcefile~allocate.f90 sourcefile~vtk.f90->sourcefile~deallocate.f90 sourcefile~vtk.f90->sourcefile~unique.f90 sourcefile~vtk.f90->sourcefile~vtk_constants.f90 sourcefile~vtk.f90->sourcefile~array.f90 sourcefile~vtk_wrapper.f90 vtk_wrapper.F90 sourcefile~vtk.f90->sourcefile~vtk_wrapper.f90 sourcefile~vtu_wrapper.f90 vtu_wrapper.F90 sourcefile~vtk.f90->sourcefile~vtu_wrapper.f90 sourcefile~array.f90->sourcefile~allocate.f90 sourcefile~array.f90->sourcefile~deallocate.f90 sourcefile~c_utils.f90 c_utils.F90 sourcefile~memory_stats_wrapper.f90->sourcefile~c_utils.f90 sourcefile~project_settings.f90->sourcefile~core.f90 sourcefile~signal_flag.f90 signal_flag.F90 sourcefile~signal_flag_wrapper.f90->sourcefile~signal_flag.f90 sourcefile~system_info_wrapper.f90->sourcefile~c_utils.f90 sourcefile~variable.f90->sourcefile~allocate.f90 sourcefile~c_utils.f90->sourcefile~signal_flag.f90 sourcefile~memory_stats.f90 memory_stats.F90 sourcefile~c_utils.f90->sourcefile~memory_stats.f90 sourcefile~system_info.f90 system_info.F90 sourcefile~c_utils.f90->sourcefile~system_info.f90

Files dependent on this one

sourcefile~~domain_manager.f90~~AfferentGraph sourcefile~domain_manager.f90 domain_manager.F90 sourcefile~domain.f90 domain.F90 sourcefile~domain.f90->sourcefile~domain_manager.f90 sourcefile~boundary_interface.f90 boundary_interface.F90 sourcefile~boundary_interface.f90->sourcefile~domain.f90 sourcefile~matrix.f90 matrix.F90 sourcefile~boundary_interface.f90->sourcefile~matrix.f90 sourcefile~boundary_manager.f90 boundary_manager.F90 sourcefile~boundary_manager.f90->sourcefile~domain.f90 sourcefile~boundary_manager.f90->sourcefile~boundary_interface.f90 sourcefile~boundary_manager.f90->sourcefile~matrix.f90 sourcefile~ftdss.f90 ftdss.F90 sourcefile~ftdss.f90->sourcefile~domain.f90 sourcefile~thermal_interface.f90 thermal_interface.F90 sourcefile~ftdss.f90->sourcefile~thermal_interface.f90 sourcefile~boundary.f90 boundary.F90 sourcefile~ftdss.f90->sourcefile~boundary.f90 sourcefile~initial.f90 initial.F90 sourcefile~ftdss.f90->sourcefile~initial.f90 sourcefile~output.f90 output.F90 sourcefile~ftdss.f90->sourcefile~output.f90 sourcefile~initial_interface.f90 initial_interface.F90 sourcefile~initial_interface.f90->sourcefile~domain.f90 sourcefile~initial_interface.f90->sourcefile~boundary.f90 sourcefile~initial_manager.f90 initial_manager.F90 sourcefile~initial_manager.f90->sourcefile~domain.f90 sourcefile~initial_manager.f90->sourcefile~initial_interface.f90 sourcefile~matrix_base.f90 matrix_base.F90 sourcefile~matrix_base.f90->sourcefile~domain.f90 sourcefile~matrix_coo.f90 matrix_coo.F90 sourcefile~matrix_coo.f90->sourcefile~domain.f90 sourcefile~matrix_coo.f90->sourcefile~matrix_base.f90 sourcefile~matrix_crs.f90 matrix_crs.F90 sourcefile~matrix_crs.f90->sourcefile~domain.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_base.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_coo.f90 sourcefile~matrix_dense.f90 matrix_dense.F90 sourcefile~matrix_dense.f90->sourcefile~domain.f90 sourcefile~matrix_dense.f90->sourcefile~matrix_base.f90 sourcefile~output_interface.f90 output_interface.F90 sourcefile~output_interface.f90->sourcefile~domain.f90 sourcefile~output_interface.f90->sourcefile~matrix.f90 sourcefile~thermal_interface.f90->sourcefile~domain.f90 sourcefile~thermal_interface.f90->sourcefile~boundary.f90 sourcefile~thermal_interface.f90->sourcefile~matrix.f90 sourcefile~solver.f90 solver.F90 sourcefile~thermal_interface.f90->sourcefile~solver.f90 sourcefile~boundary.f90->sourcefile~boundary_interface.f90 sourcefile~boundary.f90->sourcefile~boundary_manager.f90 sourcefile~boundary_adiabatic.f90 boundary_adiabatic.F90 sourcefile~boundary_adiabatic.f90->sourcefile~boundary_interface.f90 sourcefile~boundary_base.f90 boundary_base.F90 sourcefile~boundary_base.f90->sourcefile~boundary_interface.f90 sourcefile~boundary_dirichlet.f90 boundary_dirichlet.F90 sourcefile~boundary_dirichlet.f90->sourcefile~boundary_interface.f90 sourcefile~initial.f90->sourcefile~initial_interface.f90 sourcefile~initial.f90->sourcefile~initial_manager.f90 sourcefile~initial_laplace.f90 initial_laplace.F90 sourcefile~initial_laplace.f90->sourcefile~initial_interface.f90 sourcefile~initial_uniform.f90 initial_uniform.F90 sourcefile~initial_uniform.f90->sourcefile~initial_interface.f90 sourcefile~matrix.f90->sourcefile~matrix_base.f90 sourcefile~matrix.f90->sourcefile~matrix_coo.f90 sourcefile~matrix.f90->sourcefile~matrix_crs.f90 sourcefile~matrix.f90->sourcefile~matrix_dense.f90 sourcefile~output.f90->sourcefile~output_interface.f90 sourcefile~output_base.f90 output_base.F90 sourcefile~output_base.f90->sourcefile~output_interface.f90 sourcefile~output_observation.f90 output_observation.F90 sourcefile~output_observation.f90->sourcefile~output_interface.f90 sourcefile~output_overall_base.f90 output_overall_base.F90 sourcefile~output_overall_base.f90->sourcefile~output_interface.f90 sourcefile~output_overall_vtk.f90 output_overall_vtk.F90 sourcefile~output_overall_vtk.f90->sourcefile~output_interface.f90 sourcefile~output_overall_vtu.f90 output_overall_vtu.F90 sourcefile~output_overall_vtu.f90->sourcefile~output_interface.f90 sourcefile~output_system_logger.f90 output_system_logger.F90 sourcefile~output_system_logger.f90->sourcefile~output_interface.f90 sourcefile~thermal.f90 thermal.F90 sourcefile~thermal.f90->sourcefile~thermal_interface.f90 sourcefile~thermal_3phase.f90 thermal_3phase.F90 sourcefile~thermal_3phase.f90->sourcefile~thermal_interface.f90 sourcefile~solver_factory.f90 solver_factory.F90 sourcefile~solver_factory.f90->sourcefile~matrix.f90 sourcefile~solver.f90->sourcefile~solver_factory.f90

Source Code

module domain_manager
    use, intrinsic :: iso_fortran_env, only: int32
    use :: stdlib_logger
    use :: module_core, only:type_dp_3d, allocate_array, deallocate_array
    use :: module_input, only:type_input
    use :: domain_element
    use :: domain_side, only:holder_sides
    use :: domain_element_factory, only:create_element
    use :: domain_side_factory, only:create_side
    use :: domain_adjacency, only:type_node_adjacency, type_crs_adjacency_element
    use :: domain_multicoloring, only:type_coloring, type_colored_info
    use :: domain_reordering, only:type_reordering
    implicit none
    private

    public :: type_domain

    type :: type_domain
        integer(int32), private :: num_sides
        integer(int32), private :: num_elements
        integer(int32), private :: num_volumes
        integer(int32), private :: num_nodes
        integer(int32), private :: num_materials
        type(holder_elements), allocatable :: elements(:)
        type(holder_sides), allocatable :: sides(:)

        type(type_coloring) :: colors
        type(type_reordering) :: reordering

        integer(int32), private :: computaion_dimension
        ! ...
    contains
        procedure, pass(self) :: initialize => initialize_type_domain
        procedure, pass(self) :: apply_reordering

        procedure, pass(self) :: get_num_elements
        procedure, pass(self) :: get_num_sides
        procedure, pass(self) :: get_num_nodes
        procedure, pass(self) :: get_num_materials
        procedure, pass(self) :: get_computation_dimension
    end type type_domain

contains
    subroutine initialize_type_domain(self, input, Coordinate, ierr)
        implicit none
        class(type_domain), intent(inout) :: self
        type(type_input), intent(in) :: input
        type(type_dp_3d), intent(inout), pointer :: Coordinate
        integer(int32), intent(inout) :: ierr

        type(type_crs_adjacency_element) :: element_adjacency
        type(type_node_adjacency) :: node_adjacency

        integer(int32) :: count_sides, count_elements, count_volumes
        integer(int32) :: iCell, iElem, iSide
        integer(int32) :: factory_ierr
        integer(int32) :: cell_dimension

        ! -----------------------------------------------------------------------------------
        ! 初期化処理
        ! -----------------------------------------------------------------------------------
        ierr = 0
        count_sides = 0
        count_elements = 0
        count_volumes = 0

        do iCell = 1, input%geometry%vtk%num_total_cells
            cell_dimension = input%geometry%vtk%cells(iCell)%get_dimension()
            select case (cell_dimension)
            case (1)
                count_sides = count_sides + 1
            case (2)
                count_elements = count_elements + 1
            case (3)
                count_volumes = count_volumes + 1
            end select
        end do

        self%num_elements = count_elements
        self%num_sides = count_sides
        self%num_nodes = input%geometry%vtk%num_points
        self%num_materials = input%basic%num_materials

        if (allocated(self%elements)) deallocate (self%elements)
        allocate (self%elements(self%num_elements))
        if (allocated(self%sides)) deallocate (self%sides)
        allocate (self%sides(self%num_sides))

        iElem = 1
        iSide = 1
        do iCell = 1, input%geometry%vtk%num_total_cells
            cell_dimension = input%geometry%vtk%cells(iCell)%get_dimension()
            select case (cell_dimension)
            case (1)
                call create_side(new_side=self%sides(iSide)%s, &
                                 id=iCell, &
                                 global_coordinate=Coordinate, &
                                 cell_info=input%geometry%vtk%cells(iCell), &
                                 integration=input%basic%geometry_settings, &
                                 ierr=factory_ierr)
                if (factory_ierr /= 0) then
                    ierr = -1
                    return
                end if
                iSide = iSide + 1
            case (2)
                call create_element(new_element=self%elements(iElem)%e, &
                                    id=iCell, &
                                    global_coordinate=Coordinate, &
                                    cell_info=input%geometry%vtk%cells(iCell), &
                                    integration=input%basic%geometry_settings, &
                                    ierr=factory_ierr)
                if (factory_ierr /= 0) then
                    ierr = -1
                    return
                end if
                iElem = iElem + 1
            case (3)
                !!TBI
            end select

        end do

        self%computaion_dimension = input%basic%simulation_settings%calculate_dimension

        !===============================================================
        ! 3. 隣接行列の構築
        !===============================================================
        call element_adjacency%initialize(self%elements)
        print *, "Step 3a: Element adjacency matrix created."

        !===============================================================
        ! 4. RCM並べ替えの実行
        !===============================================================
        call self%reordering%initialize(input%basic%solver_settings%reordering, self%elements)
        if (input%basic%solver_settings%reordering /= "none") then
            call self%apply_reordering()
            call global_logger%log_information(message="RCM reordering completed.")
        end if

        !===============================================================
        ! 5. グラフ彩色の実行
        !===============================================================
        call self%colors%initialize(input%basic%solver_settings%coloring, element_adjacency)
        call global_logger%log_information(message="Graph coloring completed using " &
                                           //trim(self%colors%algorithm_name)//" algorithm.")

        !===============================================================
        ! 6. 後片付け
        !===============================================================
        call global_logger%log_information(message="Initialization process completed successfully.")

    end subroutine initialize_type_domain

    function get_num_elements(self) result(num_elements)
        implicit none
        class(type_domain), intent(in) :: self
        integer(int32) :: num_elements

        num_elements = self%num_elements

    end function get_num_elements

    function get_num_sides(self) result(num_sides)
        implicit none
        class(type_domain), intent(in) :: self
        integer(int32) :: num_sides

        num_sides = self%num_sides

    end function get_num_sides

    function get_num_nodes(self) result(num_nodea)
        implicit none
        class(type_domain), intent(in) :: self
        integer(int32) :: num_nodea

        num_nodea = self%num_nodes

    end function get_num_nodes

    function get_num_materials(self) result(num_materials)
        implicit none
        class(type_domain), intent(in) :: self
        integer(int32) :: num_materials

        num_materials = self%num_materials

    end function get_num_materials

    function get_computation_dimension(self) result(computaion_dimension)
        implicit none
        class(type_domain), intent(in) :: self
        integer(int32) :: computaion_dimension

        computaion_dimension = self%computaion_dimension

    end function get_computation_dimension

    subroutine apply_reordering(self)
        implicit none
        class(type_domain), intent(inout) :: self

        integer(int32) :: iElem, iSide

        if (self%computaion_dimension >= 3) then
            !! TBI: Handle 3D reordering if necessary
        end if
        if (self%computaion_dimension >= 2) then
            do iElem = 1, self%num_elements
                call allocate_array(self%elements(iElem)%e%connectivity_reordered, self%elements(iElem)%e%get_num_nodes())
                call self%reordering%to_reordered(self%elements(iElem)%e%connectivity, &
                                                  self%elements(iElem)%e%connectivity_reordered)
                if (associated(self%elements(iElem)%e%interpolate)) then
                    nullify (self%elements(iElem)%e%interpolate)
                end if
                self%elements(iElem)%e%interpolate => interpolate_reordered

                if (associated(self%elements(iElem)%e%get_connectivity)) then
                    nullify (self%elements(iElem)%e%get_connectivity)
                end if
                self%elements(iElem)%e%get_connectivity => get_connectivity_reordered
            end do
        end if
        if (self%computaion_dimension >= 1) then
            do iSide = 1, self%num_sides
                call allocate_array(self%sides(iSide)%s%connectivity_reordered, self%sides(iSide)%s%get_num_nodes())
                call self%reordering%to_reordered(self%sides(iSide)%s%connectivity, &
                                                  self%sides(iSide)%s%connectivity_reordered)
            end do
        end if

    end subroutine apply_reordering

end module domain_manager