compute_jacobian Subroutine

private subroutine compute_jacobian(me, x, jac)

Compute the Jacobian.

Note

The output jac only includes the elements of the nonlinear Jacobian. If the constant elements are being handled separately (if the linear pattern is available), then those elements can be obtained by calling get_sparsity_pattern if required.

Type Bound

numdiff_type

Arguments

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

vector of variables (size n)

real(kind=wp), intent(out), dimension(:), allocatable :: jac

sparse jacobian vector


Calls

proc~~compute_jacobian~~CallsGraph proc~compute_jacobian numerical_differentiation_module::numdiff_type%compute_jacobian proc~compute_perturbation_vector numerical_differentiation_module::numdiff_type%compute_perturbation_vector proc~compute_jacobian->proc~compute_perturbation_vector proc~raise_exception numerical_differentiation_module::numdiff_type%raise_exception proc~compute_jacobian->proc~raise_exception proc~compute_perturb_vector numerical_differentiation_module::numdiff_type%compute_perturb_vector proc~compute_perturbation_vector->proc~compute_perturb_vector proc~compute_perturb_vector->proc~raise_exception

Called by

proc~~compute_jacobian~~CalledByGraph proc~compute_jacobian numerical_differentiation_module::numdiff_type%compute_jacobian proc~compute_jacobian_dense numerical_differentiation_module::numdiff_type%compute_jacobian_dense proc~compute_jacobian_dense->proc~compute_jacobian proc~compute_jacobian_times_vector numerical_differentiation_module::numdiff_type%compute_jacobian_times_vector proc~compute_jacobian_times_vector->proc~compute_jacobian

Source Code

    subroutine compute_jacobian(me,x,jac)

    implicit none

    class(numdiff_type),intent(inout)             :: me
    real(wp),dimension(:),intent(in)              :: x    !! vector of variables (size `n`)
    real(wp),dimension(:),allocatable,intent(out) :: jac  !! sparse jacobian vector

    real(wp),dimension(me%n) :: dx  !! absolute perturbation (>0) for each variable

    if (me%exception_raised) return ! check for exceptions

    ! if we don't have a sparsity pattern yet then compute it:
    ! [also computes the indices vector]
    if (.not. me%sparsity%sparsity_computed) call me%compute_sparsity(x)
    if (me%sparsity%num_nonzero_elements==0) return

    ! size the jacobian vector:
    allocate(jac(me%sparsity%num_nonzero_elements))

    ! compute dx vector:
    if (.not. me%use_diff) then
        ! only need this for the finite difference methods (not diff)
        call me%compute_perturbation_vector(x,dx)
        if (me%exception_raised) return ! check for exceptions
    end if

    ! compute the jacobian:
    if (associated(me%jacobian_function)) then
        call me%jacobian_function(x,dx,jac)
    else
        call me%raise_exception(22,'compute_jacobian',&
                                   'jacobian_function has not been associated.')
        return
    end if

    end subroutine compute_jacobian