compute_kernel_real64 Subroutine

private elemental subroutine compute_kernel_real64(agm, init_a, init_g)

Compute arithmetic–geometric mean (AGM) using the given arithmetic mean init_a and geometric mean init_g and store the complete iteration history.

This type-bound subroutine performs the AGM iteration starting from given initial values, storing each intermediate arithmetic and geometric mean in list_a and list_g respectively. The final iteration count is stored in n_iter.

Unlike the arithmetic_geometric_mean_kernel, this subroutine preserves the full convergence history for analysis.

Warning

  • This subroutine/interface assumes both inputs are positive.
  • No validation is performed on inputs.
  • Calling this subroutine overwrites any previous computation stored in the object.

Note

Convergence criterion
See is_not_converged

Type Bound

arithmetic_geometric_mean_real64_type

Arguments

Type IntentOptional Attributes Name
class(arithmetic_geometric_mean_real64_type), intent(inout) :: agm
real(kind=real64), intent(in) :: init_a

initial value: arithmetic mean

real(kind=real64), intent(in) :: init_g

initial value: geometric mean


Source Code

    elemental subroutine compute_kernel_real64(agm, init_a, init_g)
        !! Compute arithmetic–geometric mean (AGM)
        !! using the given arithmetic mean `init_a` and geometric mean `init_g`
        !! and store the complete iteration history.
        !!
        !! This type-bound subroutine performs the AGM iteration starting from
        !! given initial values, storing each intermediate arithmetic and geometric mean
        !! in `list_a` and `list_g` respectively.
        !! The final iteration count is stored in `n_iter`.
        !!
        !! Unlike the [[arithmetic_geometric_mean_kernel]],
        !! this subroutine preserves the full convergence history for analysis.
        !!
        !! @warning
        !! - This subroutine/interface assumes both inputs are positive.
        !! - No validation is performed on inputs.
        !! - Calling this subroutine overwrites any previous computation stored in the object.
        !! @endwarning
        !!
        !! @note
        !! **Convergence criterion**  
        !! See [[is_not_converged]]
        !! @endnote

        class(arithmetic_geometric_mean_real64_type), intent(inout) :: agm

        real(real64), intent(in) :: init_a
        !! initial value: arithmetic mean

        real(real64), intent(in) :: init_g
        !! initial value: geometric mean



        call initialize(agm)

        agm%list_a(0) = init_a
        agm%list_g(0) = init_g
        agm%n_iter_    = 0

        do

            associate(prev_iter => agm%n_iter_, next_iter => agm%n_iter_ + 1)

                associate( &!
                    prev_a => agm%list_a(prev_iter) , &!
                    prev_g => agm%list_g(prev_iter) , &!
                    next_a => agm%list_a(next_iter) , &!
                    next_g => agm%list_g(next_iter)   &!
                )

                    agm%n_iter_ = agm%n_iter_ + 1

                    call compute_step( &!
                        prev_a = prev_a , &!
                        prev_g = prev_g , &!
                        next_a = next_a , &!
                        next_g = next_g   &!
                    )

                    if ( is_not_converged(next_a, next_g) ) then

                        prev_iter = next_iter

                        cycle

                    else

                        prev_iter = next_iter

                        return

                    end if

                end associate

            end associate

        end do

    end subroutine compute_kernel_real64