select_finite_diff_method_for_partition_group Subroutine

private subroutine select_finite_diff_method_for_partition_group(me, x, xlow, xhigh, dx, list_of_methods, fd, status_ok)

Select a finite diff method of a given class so that the perturbations of x will not violate the variable bounds for any variable in the group.

The x vector are only the variables in a group (not the full variable vector)

Type Bound

numdiff_type

Arguments

Type IntentOptional Attributes Name
class(numdiff_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

the variable values

real(kind=wp), intent(in), dimension(:) :: xlow

the variable lower bounds

real(kind=wp), intent(in), dimension(:) :: xhigh

the variable upper bounds

real(kind=wp), intent(in), dimension(:) :: dx

the perturbation values (>0)

type(meth_array), intent(in) :: list_of_methods

list of available methods to choose from

type(finite_diff_method), intent(out) :: fd

this method can be used

logical, intent(out) :: status_ok

true if it really doesn't violate the bounds (say, the bounds are very close to each other) if status_ok=False, then the first method in the given class is returned in fd.


Called by

proc~~select_finite_diff_method_for_partition_group~~CalledByGraph proc~select_finite_diff_method_for_partition_group numerical_differentiation_module::numdiff_type%select_finite_diff_method_for_partition_group proc~compute_jacobian_partitioned numerical_differentiation_module::compute_jacobian_partitioned proc~compute_jacobian_partitioned->proc~select_finite_diff_method_for_partition_group

Source Code

    subroutine select_finite_diff_method_for_partition_group(me,x,xlow,xhigh,dx,&
                                                             list_of_methods,fd,status_ok)

    implicit none

    class(numdiff_type),intent(inout)     :: me
    real(wp),dimension(:),intent(in)      :: x               !! the variable values
    real(wp),dimension(:),intent(in)      :: xlow            !! the variable lower bounds
    real(wp),dimension(:),intent(in)      :: xhigh           !! the variable upper bounds
    real(wp),dimension(:),intent(in)      :: dx              !! the perturbation values (>0)
    type(meth_array),intent(in)           :: list_of_methods !! list of available methods to choose from
    type(finite_diff_method),intent(out)  :: fd              !! this method can be used
    logical,intent(out)                   :: status_ok       !! true if it really doesn't violate the bounds
                                                             !! (say, the bounds are very close to each other)
                                                             !! if `status_ok=False`, then the first method in
                                                             !! the given class is returned in `fd`.

    integer  :: i   !! counter
    integer  :: j   !! counter
    real(wp),dimension(size(x)) :: xp  !! perturbed `x` values
    real(wp),dimension(size(x)) :: xs  !! if `x` is outside the bounds, this is the value
                                       !! on the nearest bound. otherwise equal to `x`.

    ! initialize:
    status_ok = .false.

    if (me%exception_raised) return ! check for exceptions

    ! make sure they are within the bounds
    xs = min(xhigh,max(xlow,x))

    ! try all the methods in the class:
    do i = 1, size(list_of_methods%meth)
        status_ok = .true. ! will be set to false if any
                           ! perturbation violates the bounds
        ! check each of the perturbations:
        do j = 1, size(list_of_methods%meth(i)%dx_factors)
            xp = xs + list_of_methods%meth(i)%dx_factors(j)*dx
            if (any(xp < xlow) .or. any(xp > xhigh)) then
                status_ok = .false.
                exit
            end if
        end do
        if (status_ok) then   ! this one is OK to use
            fd = list_of_methods%meth(i)
            exit
        end if
    end do

    if (.not. status_ok) then
        ! no method was found that doesn't violate the bounds,
        ! so just return the first one in the list.
        fd = list_of_methods%meth(1)
    end if

    end subroutine select_finite_diff_method_for_partition_group