boundary_manager.F90 Source File


This file depends on

sourcefile~~boundary_manager.f90~~EfferentGraph sourcefile~boundary_manager.f90 boundary_manager.F90 sourcefile~boundary_interface.f90 boundary_interface.F90 sourcefile~boundary_manager.f90->sourcefile~boundary_interface.f90 sourcefile~core.f90 core.F90 sourcefile~boundary_manager.f90->sourcefile~core.f90 sourcefile~domain.f90 domain.F90 sourcefile~boundary_manager.f90->sourcefile~domain.f90 sourcefile~input.f90 input.F90 sourcefile~boundary_manager.f90->sourcefile~input.f90 sourcefile~matrix.f90 matrix.F90 sourcefile~boundary_manager.f90->sourcefile~matrix.f90 sourcefile~boundary_interface.f90->sourcefile~core.f90 sourcefile~boundary_interface.f90->sourcefile~domain.f90 sourcefile~boundary_interface.f90->sourcefile~input.f90 sourcefile~boundary_interface.f90->sourcefile~matrix.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~adjacency.f90 adjacency.F90 sourcefile~domain.f90->sourcefile~adjacency.f90 sourcefile~domain_manager.f90 domain_manager.F90 sourcefile~domain.f90->sourcefile~domain_manager.f90 sourcefile~element.f90 element.F90 sourcefile~domain.f90->sourcefile~element.f90 sourcefile~element_factory.f90 element_factory.F90 sourcefile~domain.f90->sourcefile~element_factory.f90 sourcefile~multicoloring.f90 multicoloring.F90 sourcefile~domain.f90->sourcefile~multicoloring.f90 sourcefile~reordering.f90 reordering.F90 sourcefile~domain.f90->sourcefile~reordering.f90 sourcefile~side.f90 side.F90 sourcefile~domain.f90->sourcefile~side.f90 sourcefile~side_factory.f90 side_factory.F90 sourcefile~domain.f90->sourcefile~side_factory.f90 sourcefile~input_interface.f90 input_interface.F90 sourcefile~input.f90->sourcefile~input_interface.f90 sourcefile~matrix_base.f90 matrix_base.F90 sourcefile~matrix.f90->sourcefile~matrix_base.f90 sourcefile~matrix_coo.f90 matrix_coo.F90 sourcefile~matrix.f90->sourcefile~matrix_coo.f90 sourcefile~matrix_crs.f90 matrix_crs.F90 sourcefile~matrix.f90->sourcefile~matrix_crs.f90 sourcefile~matrix_dense.f90 matrix_dense.F90 sourcefile~matrix.f90->sourcefile~matrix_dense.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->sourcefile~error.f90 sourcefile~deallocate.f90->sourcefile~error.f90 sourcefile~domain_manager.f90->sourcefile~core.f90 sourcefile~domain_manager.f90->sourcefile~input.f90 sourcefile~domain_manager.f90->sourcefile~adjacency.f90 sourcefile~domain_manager.f90->sourcefile~element.f90 sourcefile~domain_manager.f90->sourcefile~element_factory.f90 sourcefile~domain_manager.f90->sourcefile~multicoloring.f90 sourcefile~domain_manager.f90->sourcefile~reordering.f90 sourcefile~domain_manager.f90->sourcefile~side.f90 sourcefile~domain_manager.f90->sourcefile~side_factory.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~input.f90 sourcefile~element_factory.f90->sourcefile~element.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~matrix_base.f90->sourcefile~domain.f90 sourcefile~matrix_coo.f90->sourcefile~core.f90 sourcefile~matrix_coo.f90->sourcefile~domain.f90 sourcefile~matrix_coo.f90->sourcefile~matrix_base.f90 sourcefile~matrix_crs.f90->sourcefile~core.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->sourcefile~core.f90 sourcefile~matrix_dense.f90->sourcefile~domain.f90 sourcefile~matrix_dense.f90->sourcefile~matrix_base.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~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~adjacency_element.f90->sourcefile~core.f90 sourcefile~adjacency_element.f90->sourcefile~element.f90 sourcefile~adjacency_node.f90->sourcefile~core.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~~boundary_manager.f90~~AfferentGraph sourcefile~boundary_manager.f90 boundary_manager.F90 sourcefile~boundary.f90 boundary.F90 sourcefile~boundary.f90->sourcefile~boundary_manager.f90 sourcefile~ftdss.f90 ftdss.F90 sourcefile~ftdss.f90->sourcefile~boundary.f90 sourcefile~thermal_interface.f90 thermal_interface.F90 sourcefile~ftdss.f90->sourcefile~thermal_interface.f90 sourcefile~initial.f90 initial.F90 sourcefile~ftdss.f90->sourcefile~initial.f90 sourcefile~initial_interface.f90 initial_interface.F90 sourcefile~initial_interface.f90->sourcefile~boundary.f90 sourcefile~thermal_interface.f90->sourcefile~boundary.f90 sourcefile~initial.f90->sourcefile~initial_interface.f90 sourcefile~initial_manager.f90 initial_manager.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_manager.f90->sourcefile~initial_interface.f90 sourcefile~initial_uniform.f90 initial_uniform.F90 sourcefile~initial_uniform.f90->sourcefile~initial_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

Source Code

module conditions_boundary_manager
    use :: iso_fortran_env
    use :: module_core, only:allocate_array, deallocate_array
    use :: module_domain, only:type_domain
    use :: module_matrix, only:type_crs
    use :: module_input, only:type_input
    use :: conditions_boundary, only:abst_bc_thermal, type_bc_thermal_adiabatic, type_bc_thermal_dirichlet
    implicit none
    private

    public :: holder_bcs
    public :: type_bc

    type :: holder_bcs
        class(abst_bc_thermal), allocatable :: t
    end type holder_bcs

    type :: type_bc
        integer(int32) :: num_boundaries
        type(holder_bcs), allocatable :: bc(:)
    contains
        procedure :: initialize => initialize_type_bc
        procedure :: apply_crs => apply_type_bc_crs
        ! generic :: apply => apply_crs
    end type type_bc

contains

    subroutine initialize_type_bc(self, input, domain)
        class(type_bc), intent(inout) :: self
        type(type_input), intent(in) :: input
        type(type_domain), intent(in) :: domain

        integer(int32), allocatable :: group_ids(:)
        integer(int32) :: i, i_material

        real(real64) :: time_conv

        self%num_boundaries = input%conditions%num_boundaries

        ! 2. 各グループのコンテナを確保
        if (allocated(self%bc)) deallocate (self%bc)
        allocate (self%bc(self%num_boundaries))

        select case (input%conditions%time_control%simulation_period%unit)
        case ("second")
            time_conv = 1.0d0
        case ("minute")
            time_conv = 60.0d0
        case ("hour")
            time_conv = 3600.0d0
        case ("day")
            time_conv = 86400.0d0
        case ("year")
            time_conv = 31557600.0d0
        end select

        ! 3. 各グループをループして、対応するBCオブジェクトを生成・セットアップ
        do i = 1, self%num_boundaries
            i_material = input%conditions%boundary_conditions(i)%id
            if (input%basic%analysis_controls%calculate_thermal) then
                select case (input%conditions%boundary_conditions(i)%thermal%type)
                case ("dirichlet")
                    allocate (type_bc_thermal_dirichlet :: self%bc(i)%t)
                    call self%bc(i)%t%initialize(input, domain, i, i_material, time_conv)
                case ("adiabatic")
                    allocate (type_bc_thermal_adiabatic :: self%bc(i)%t)
                    call self%bc(i)%t%initialize(input, domain, i, i_material, time_conv)
                end select
            end if
            if (input%basic%analysis_controls%calculate_hydraulic) then
                ! select case (input%conditions%boundary_conditions(i)%thermal%type)
                ! case ("dirichlet")
                !     allocate (type_bc_thermal_dirichlet :: self%bc(i)%t)
                !     call self%bc(i)%t%initialize(input, domain, i_material, time_conv)
                ! case ("adiabatic")
                !     allocate (type_bc_thermal_adiabatic :: self%bc(i)%t)
                !     call self%bc(i)%t%initialize(input, domain, i_material, time_conv)
                ! end select
            end if
        end do
    end subroutine initialize_type_bc

    subroutine apply_type_bc_crs(self, boundary_target, current_time, A, b, domain, mode)
        class(type_bc), intent(inout) :: self
        character(*), intent(in) :: boundary_target
        real(real64), intent(in) :: current_time
        type(type_crs), intent(inout), optional :: A
        real(real64), intent(inout) :: b(:)
        type(type_domain), intent(in) :: domain
        integer(int32), intent(in), optional :: mode

        integer(int32) :: i

        ! --------------------------------------------------------------------------
        ! 1st Pass: Apply all non-Dirichlet boundary conditions (e.g., Neumann, Adiabatic)
        ! --------------------------------------------------------------------------
        select case (trim(adjustl(boundary_target)))
        case ('thermal')
            do i = 1, self%num_boundaries
                if (allocated(self%bc(i)%t)) then
                    select type (bc => self%bc(i)%t)
                    type is (type_bc_thermal_dirichlet)
                        ! This is a Dirichlet BC, so we skip it in the first pass.
                        cycle
                    class default
                        ! This is any other type of BC, apply it now.
                        if (present(A)) then
                            ! Apply the BC using CRS matrix format.
                            if (present(mode)) then
                                call bc%apply_crs(current_time=current_time, &
                                                  A=A, &
                                                  b=b, &
                                                  domain=domain, &
                                                  mode=mode)
                            else
                                call bc%apply_crs(current_time=current_time, &
                                                  b=b, &
                                                  domain=domain)
                            end if
                        else
                            ! Apply the BC without matrix A.
                            if (present(mode)) then
                                call bc%apply_crs(current_time=current_time, &
                                                  b=b, &
                                                  domain=domain, &
                                                  mode=mode)
                            else
                                call bc%apply_crs(current_time=current_time, &
                                                  b=b, &
                                                  domain=domain)
                            end if
                        end if
                    end select
                end if
            end do

            ! --------------------------------------------------------------------------
            ! 2nd Pass: Apply ONLY the Dirichlet boundary conditions
            ! --------------------------------------------------------------------------
            do i = 1, self%num_boundaries
                if (allocated(self%bc(i)%t)) then
                    select type (bc => self%bc(i)%t)
                    type is (type_bc_thermal_dirichlet)
                        ! This is a Dirichlet BC, apply it in the final pass.
                        if (present(A)) then
                            ! Apply the BC using CRS matrix format.
                            if (present(mode)) then
                                call bc%apply_crs(current_time=current_time, &
                                                  A=A, &
                                                  b=b, &
                                                  domain=domain, &
                                                  mode=mode)
                            else
                                call bc%apply_crs(current_time=current_time, &
                                                  b=b, &
                                                  domain=domain)
                            end if
                        else
                            ! Apply the BC without matrix A.
                            if (present(mode)) then
                                call bc%apply_crs(current_time=current_time, &
                                                  b=b, &
                                                  domain=domain, &
                                                  mode=mode)
                            else
                                call bc%apply_crs(current_time=current_time, &
                                                  b=b, &
                                                  domain=domain)
                            end if
                        end if
                    class default
                        ! All other types were handled in the first pass, so do nothing.
                        cycle
                    end select
                end if
            end do
        end select

    end subroutine apply_type_bc_crs

end module conditions_boundary_manager