iteration.F90 Source File


This file depends on

sourcefile~~iteration.f90~~EfferentGraph sourcefile~iteration.f90 iteration.F90 sourcefile~calculate.f90 calculate.F90 sourcefile~iteration.f90->sourcefile~calculate.f90 sourcefile~input.f90 input.F90 sourcefile~iteration.f90->sourcefile~input.f90 sourcefile~density_interface.f90 density_interface.F90 sourcefile~calculate.f90->sourcefile~density_interface.f90 sourcefile~gcc_interface.f90 gcc_interface.F90 sourcefile~calculate.f90->sourcefile~gcc_interface.f90 sourcefile~hcf_interface.f90 hcf_interface.F90 sourcefile~calculate.f90->sourcefile~hcf_interface.f90 sourcefile~heat_capacity_interface.f90 heat_capacity_interface.F90 sourcefile~calculate.f90->sourcefile~heat_capacity_interface.f90 sourcefile~linalg.f90 linalg.F90 sourcefile~calculate.f90->sourcefile~linalg.f90 sourcefile~specific_heat_interface.f90 specific_heat_interface.F90 sourcefile~calculate.f90->sourcefile~specific_heat_interface.f90 sourcefile~thermal_conductivity_interface.f90 thermal_conductivity_interface.F90 sourcefile~calculate.f90->sourcefile~thermal_conductivity_interface.f90 sourcefile~input_interface.f90 input_interface.F90 sourcefile~input.f90->sourcefile~input_interface.f90 sourcefile~density_interface.f90->sourcefile~input.f90 sourcefile~core.f90 core.F90 sourcefile~density_interface.f90->sourcefile~core.f90 sourcefile~gcc_interface.f90->sourcefile~input.f90 sourcefile~gcc_interface.f90->sourcefile~core.f90 sourcefile~hcf_interface.f90->sourcefile~input.f90 sourcefile~hcf_interface.f90->sourcefile~core.f90 sourcefile~heat_capacity_interface.f90->sourcefile~density_interface.f90 sourcefile~heat_capacity_interface.f90->sourcefile~input_interface.f90 sourcefile~heat_capacity_interface.f90->sourcefile~core.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_ops.f90 matrix_ops.F90 sourcefile~linalg.f90->sourcefile~matrix_ops.f90 sourcefile~matvec.f90 matvec.F90 sourcefile~linalg.f90->sourcefile~matvec.f90 sourcefile~vector_ops.f90 vector_ops.F90 sourcefile~linalg.f90->sourcefile~vector_ops.f90 sourcefile~specific_heat_interface.f90->sourcefile~input.f90 sourcefile~specific_heat_interface.f90->sourcefile~core.f90 sourcefile~thermal_conductivity_interface.f90->sourcefile~input.f90 sourcefile~thermal_conductivity_interface.f90->sourcefile~core.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~matrix_ops.f90->sourcefile~core.f90 sourcefile~matvec.f90->sourcefile~core.f90 sourcefile~project_settings.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~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~matrix.f90 matrix.F90 sourcefile~types.f90->sourcefile~matrix.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~types.f90 sourcefile~vtk.f90->sourcefile~unique.f90 sourcefile~vtk.f90->sourcefile~vtk_constants.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~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~matrix_interface.f90 matrix_interface.F90 sourcefile~matrix.f90->sourcefile~matrix_interface.f90 sourcefile~c_utils.f90 c_utils.F90 sourcefile~memory_stats_wrapper.f90->sourcefile~c_utils.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 sourcefile~matrix_coo.f90->sourcefile~allocate.f90 sourcefile~matrix_coo.f90->sourcefile~deallocate.f90 sourcefile~matrix_coo.f90->sourcefile~matrix_interface.f90 sourcefile~matrix_crs.f90->sourcefile~allocate.f90 sourcefile~matrix_crs.f90->sourcefile~deallocate.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_interface.f90 sourcefile~matrix_dense.f90->sourcefile~allocate.f90 sourcefile~matrix_dense.f90->sourcefile~deallocate.f90 sourcefile~matrix_dense.f90->sourcefile~matrix_interface.f90

Files dependent on this one

sourcefile~~iteration.f90~~AfferentGraph sourcefile~iteration.f90 iteration.F90 sourcefile~control.f90 control.F90 sourcefile~control.f90->sourcefile~iteration.f90 sourcefile~ftdss.f90 ftdss.F90 sourcefile~ftdss.f90->sourcefile~control.f90 sourcefile~hydraulic.f90 hydraulic.F90 sourcefile~ftdss.f90->sourcefile~hydraulic.f90 sourcefile~output.f90 output.F90 sourcefile~ftdss.f90->sourcefile~output.f90 sourcefile~properties.f90 properties.F90 sourcefile~ftdss.f90->sourcefile~properties.f90 sourcefile~thermal.f90 thermal.F90 sourcefile~ftdss.f90->sourcefile~thermal.f90 sourcefile~hydraulic_assemble.f90 hydraulic_assemble.F90 sourcefile~hydraulic_assemble.f90->sourcefile~control.f90 sourcefile~hydraulic_assemble.f90->sourcefile~properties.f90 sourcefile~hydraulic_interface.f90 hydraulic_interface.F90 sourcefile~hydraulic_interface.f90->sourcefile~control.f90 sourcefile~hydraulic_interface.f90->sourcefile~hydraulic_assemble.f90 sourcefile~hydraulic_interface.f90->sourcefile~properties.f90 sourcefile~output_interface.f90 output_interface.F90 sourcefile~output_interface.f90->sourcefile~control.f90 sourcefile~output_interface.f90->sourcefile~properties.f90 sourcefile~properties_manager.f90 properties_manager.F90 sourcefile~properties_manager.f90->sourcefile~control.f90 sourcefile~thermal_interface.f90 thermal_interface.F90 sourcefile~thermal_interface.f90->sourcefile~control.f90 sourcefile~thermal_interface.f90->sourcefile~properties.f90 sourcefile~hydraulic.f90->sourcefile~hydraulic_interface.f90 sourcefile~hydraulic_crs.f90 hydraulic_crs.F90 sourcefile~hydraulic_crs.f90->sourcefile~hydraulic_interface.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~properties.f90->sourcefile~properties_manager.f90 sourcefile~thermal.f90->sourcefile~thermal_interface.f90 sourcefile~thermal_crs.f90 thermal_crs.F90 sourcefile~thermal_crs.f90->sourcefile~thermal_interface.f90

Source Code

module control_iteration
    use, intrinsic :: iso_fortran_env, only: int32, real64
    use :: module_input, only:type_input, type_convergence
    use :: module_calculate, only:norm_2, norm_inf
    implicit none
    private

    public :: type_iteration

    ! --- プライベート: イテレータ設定格納型 ---
    type, private :: type_iterator_config
        integer(int32) :: max_iterations
        integer(int32) :: update_frequency
        type(type_convergence) :: convergence
    end type type_iterator_config

    ! --- パブリック: イテレーション管理派生型 ---
    type :: type_iteration
        private
        integer(int32) :: iter = 0 ! 非線形反復回数
        integer(int32) :: step = 0 ! ステップカウンタ
        logical :: is_converged = .false.

        real(real64) :: init_res_norm_l2 = 0.0d0
        real(real64) :: init_res_norm_inf = 0.0d0
        real(real64) :: init_upd_norm_l2 = 0.0d0
        real(real64) :: init_upd_norm_inf = 0.0d0

        character(:), allocatable :: algorithm
        type(type_iterator_config) :: config
    contains
        procedure, pass(self), public :: initialize => initialize_type_iteration
        procedure, pass(self), public :: reset_step
        procedure, pass(self), public :: set_initial_norms
        procedure, pass(self), public :: reset_timestep
        procedure, pass(self), public :: check_convergence
        procedure, pass(self), public :: increment_iter
        procedure, pass(self), public :: increment_step
        procedure, pass(self), public :: should_continue => continue_loop
        procedure, pass(self), public :: get_iter
        procedure, pass(self), public :: get_step
        procedure, pass(self), public :: get_algorithm_name
        procedure, pass(self), public :: has_converged => get_status
    end type type_iteration

contains

    subroutine initialize_type_iteration(self, input)
        implicit none
        class(type_iteration), intent(out) :: self
        type(type_input), intent(in) :: input

        self%iter = 0
        self%step = 0
        self%is_converged = .false.

        self%algorithm = input%basic%solver_settings%nonlinear_solver%method
        select case (trim(self%algorithm))
        case ("newton", "modified_newton", "picard")
            self%config%max_iterations = input%basic%solver_settings%nonlinear_solver%max_iterations
            self%config%update_frequency = input%basic%solver_settings%nonlinear_solver%update_frequency
            self%config%convergence = input%basic%solver_settings%nonlinear_solver%convergence
        end select
    end subroutine initialize_type_iteration

    subroutine reset_step(self)
        implicit none
        class(type_iteration), intent(inout) :: self

        self%step = 0
        self%is_converged = .false.
    end subroutine reset_step

    subroutine set_initial_norms(self, res_vec, upd_vec)
        implicit none
        class(type_iteration), intent(inout) :: self
        real(real64), intent(in), optional :: res_vec(:), upd_vec(:)

        ! 初期ノルム値設定
        if (present(res_vec)) then
            self%init_res_norm_l2 = norm_2(res_vec)
            self%init_res_norm_inf = norm_inf(res_vec)
        end if
        if (present(upd_vec)) then
            self%init_upd_norm_l2 = norm_2(upd_vec)
            self%init_upd_norm_inf = norm_inf(upd_vec)
        end if
    end subroutine set_initial_norms

    subroutine reset_timestep(self)
        implicit none
        class(type_iteration), intent(inout) :: self

        self%step = 0
        self%iter = 0
        self%is_converged = .false.
        self%init_res_norm_l2 = 0.0d0
        self%init_res_norm_inf = 0.0d0
        self%init_upd_norm_l2 = 0.0d0
        self%init_upd_norm_inf = 0.0d0
    end subroutine reset_timestep

    subroutine check_convergence(self, res_vec, upd_vec)
        implicit none
        class(type_iteration), intent(inout) :: self
        real(real64), intent(in) :: res_vec(:), upd_vec(:)
        logical :: ok_res, ok_upd
        real(real64) :: abs_res, rel_res, abs_upd, rel_upd

        if (trim(self%algorithm) == "none") then
            self%is_converged = .true.
            return
        end if

        ! --- ノルム計算 ---
        abs_res = norm_2(res_vec)
        abs_upd = norm_2(upd_vec)
        if (self%init_res_norm_l2 > 1.0d-12) then
            rel_res = abs_res / self%init_res_norm_l2
        else
            rel_res = 0.0d0
        end if
        if (self%init_upd_norm_l2 > 1.0d-12) then
            rel_upd = abs_upd / self%init_upd_norm_l2
        else
            rel_upd = 0.0d0
        end if

        ! --- 残差基準 ---
        select case (trim(self%config%convergence%residual%criteria))
        case ("absolute")
            ok_res = abs_res < self%config%convergence%residual%absolute_tolerance
        case ("relative")
            ok_res = rel_res < self%config%convergence%residual%relative_tolerance
        case ("both")
            if (trim(self%config%convergence%residual%logic) == "and") then
                ok_res = (abs_res < self%config%convergence%residual%absolute_tolerance) .and. &
                         (rel_res < self%config%convergence%residual%relative_tolerance)
            else
                ok_res = (abs_res < self%config%convergence%residual%absolute_tolerance) .or. &
                         (rel_res < self%config%convergence%residual%relative_tolerance)
            end if
        end select

        ! --- 更新基準 ---
        select case (trim(self%config%convergence%update%criteria))
        case ("absolute")
            ok_upd = abs_upd < self%config%convergence%update%absolute_tolerance
        case ("relative")
            ok_upd = rel_upd < self%config%convergence%update%relative_tolerance
        case ("both")
            if (trim(self%config%convergence%update%logic) == "and") then
                ok_upd = (abs_upd < self%config%convergence%update%absolute_tolerance) .and. &
                         (rel_upd < self%config%convergence%update%relative_tolerance)
            else
                ok_upd = (abs_upd < self%config%convergence%update%absolute_tolerance) .or. &
                         (rel_upd < self%config%convergence%update%relative_tolerance)
            end if
        end select

        ! --- 総合判定 ---
        select case (trim(self%config%convergence%use_criteria))
        case ("residual")
            self%is_converged = ok_res
        case ("update")
            self%is_converged = ok_upd
        case ("both")
            if (trim(self%config%convergence%use_logic) == "and") then
                self%is_converged = ok_res .and. ok_upd
            else
                self%is_converged = ok_res .or. ok_upd
            end if
        end select
    end subroutine check_convergence

    subroutine increment_iter(self)
        implicit none
        class(type_iteration), intent(inout) :: self

        self%iter = self%iter + 1
    end subroutine increment_iter

    subroutine increment_step(self)
        implicit none
        class(type_iteration), intent(inout) :: self

        self%step = self%step + 1
    end subroutine increment_step

    function continue_loop(self) result(should_continue_loop)
        implicit none
        class(type_iteration), intent(in) :: self
        logical :: should_continue_loop

        if (self%is_converged) then
            should_continue_loop = .false.
        else
            if (self%step == 0) then
                should_continue_loop = .true.
            else
                should_continue_loop = (self%step < self%config%max_iterations)
            end if
        end if
    end function continue_loop

    pure function get_iter(self) result(iter)
        implicit none
        class(type_iteration), intent(in) :: self
        integer(int32) :: iter

        iter = self%iter
    end function get_iter

    pure function get_step(self) result(step)
        implicit none
        class(type_iteration), intent(in) :: self
        integer(int32) :: step

        step = self%step
    end function get_step

    pure function get_status(self) result(is_converged)
        implicit none
        class(type_iteration), intent(in) :: self
        logical :: is_converged

        is_converged = self%is_converged
    end function get_status

    pure function get_algorithm_name(self) result(algorithm_name)
        implicit none
        class(type_iteration), intent(in) :: self
        character(:), allocatable :: algorithm_name

        algorithm_name = self%algorithm
    end function get_algorithm_name

end module control_iteration