compute_jacobian_dense Subroutine

private subroutine compute_jacobian_dense(me, x, jac)

just a wrapper for compute_jacobian, that returns a dense (m x n) matrix.

Note

This one will include the constant elements if the linear pattern is available.

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

the jacobian matrix


Calls

proc~~compute_jacobian_dense~~CallsGraph proc~compute_jacobian_dense numerical_differentiation_module::numdiff_type%compute_jacobian_dense proc~compute_jacobian numerical_differentiation_module::numdiff_type%compute_jacobian proc~compute_jacobian_dense->proc~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

Source Code

    subroutine compute_jacobian_dense(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 !! the jacobian matrix

    real(wp),dimension(:),allocatable :: jac_vec  !! sparse jacobian representation
    integer :: i !! counter

    if (me%exception_raised) return ! check for exceptions

    ! size output matrix:
    allocate(jac(me%m,me%n))

    ! convert to dense form:
    jac = zero

    ! compute sparse form of jacobian:
    call me%compute_jacobian(x,jac_vec)
    if (me%exception_raised) return ! check for exceptions

    if (allocated(jac_vec)) then
        ! add the nonlinear elements:
        do i = 1, me%sparsity%num_nonzero_elements
            jac(me%sparsity%irow(i),me%sparsity%icol(i)) = jac_vec(i)
        end do
        deallocate(jac_vec)
    end if

    ! add the constant elements if necessary:
    do i = 1, me%sparsity%num_nonzero_linear_elements
        jac(me%sparsity%linear_irow(i),me%sparsity%linear_icol(i)) = me%sparsity%linear_vals(i)
    end do

    end subroutine compute_jacobian_dense