fdjs Subroutine

public subroutine fdjs(m, n, Col, Ind, Npntr, Ngrp, Numgrp, d, Fjacd, Fjac)

Given a consistent partition of the columns of an m by n jacobian matrix into groups, this subroutine computes approximations to those columns in a given group. the approximations are stored into either a column-oriented or a row-oriented pattern.

a partition is consistent if the columns in any group do not have a non-zero in the same row position.

approximations to the columns of the jacobian matrix in a given group can be obtained by specifying a difference parameter array d with d(jcol) non-zero if and only if jcol is a column in the group, and an approximation to jac*d where jac denotes the jacobian matrix of a mapping f.

d can be defined with the following segment of code.

  do jcol = 1, n
  d(jcol) = 0.0
  if (ngrp(jcol) == numgrp) d(jcol) = eta(jcol)
  end do

in the above code numgrp is the given group number, ngrp(jcol) is the group number of column jcol, and eta(jcol) is the difference parameter used to approximate column jcol of the jacobian matrix. suitable values for the array eta must be provided.

as mentioned above, an approximation to jac*d must also be provided. for example, the approximation

  f(x+d) - f(x)

corresponds to the forward difference formula at x.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: m

a positive integer input variable set to the number of rows of the jacobian matrix.

integer, intent(in) :: n

a positive integer input variable set to the number of columns of the jacobian matrix.

logical, intent(in) :: Col

a logical input variable. if col is set true, then the jacobian approximations are stored into a column-oriented pattern. if col is set false, then the jacobian approximations are stored into a row-oriented pattern.

integer, dimension(*) :: Ind

an integer input array which contains the row indices for the non-zeroes in the jacobian matrix if col is true, and contains the column indices for the non-zeroes in the jacobian matrix if col is false.

integer, dimension(*) :: Npntr

an integer input array which specifies the locations of the row indices in ind if col is true, and specifies the locations of the column indices in ind if col is false. if col is true, the indices for column j are ind(k), k = npntr(j),...,npntr(j+1)-1. if col is false, the indices for row i are ind(k), k = npntr(i),...,npntr(i+1)-1. Note that npntr(n+1)-1 if col is true, or npntr(m+1)-1 if col is false, is then the number of non-zero elements of the jacobian matrix.

integer, dimension(n) :: Ngrp

an integer input array of length n which specifies the partition of the columns of the jacobian matrix. column jcol belongs to group ngrp(jcol).

integer :: Numgrp

a positive integer input variable set to a group number in the partition. the columns of the jacobian matrix in this group are to be estimated on this call.

real(kind=wp), dimension(n) :: d

an input array of length n which contains the difference parameter vector for the estimate of the jacobian matrix columns in group numgrp.

real(kind=wp), dimension(m) :: Fjacd

an input array of length m which contains an approximation to the difference vector jac*d, where jac denotes the jacobian matrix.

real(kind=wp), dimension(*) :: Fjac

an output array of length nnz, where nnz is the number of its non-zero elements. at each call of fdjs, fjac is updated to include the non-zero elements of the jacobian matrix for those columns in group numgrp. fjac should not be altered between successive calls to fdjs.


Source Code

    subroutine fdjs(m,n,Col,Ind,Npntr,Ngrp,Numgrp,d,Fjacd,Fjac)

    implicit none

    integer,intent(in)    :: m        !! a positive integer input variable set to the number
                                      !! of rows of the jacobian matrix.
    integer,intent(in)    :: n        !! a positive integer input variable set to the number
                                      !! of columns of the jacobian matrix.
    integer               :: Numgrp   !! a positive integer input variable set to a group
                                      !! number in the partition. the columns of the jacobian
                                      !! matrix in this group are to be estimated on this call.
    integer,dimension(*)  :: Ind      !! an integer input array which contains the row
                                      !! indices for the non-zeroes in the jacobian matrix
                                      !! if `col` is true, and contains the column indices for
                                      !! the non-zeroes in the jacobian matrix if `col` is false.
    integer,dimension(*)  :: Npntr    !! an integer input array which specifies the
                                      !! locations of the row indices in `ind` if `col` is true, and
                                      !! specifies the locations of the column indices in `ind` if
                                      !! `col` is false. if `col` is true, the indices for column `j` are
                                      !!       `ind(k), k = npntr(j),...,npntr(j+1)-1`.
                                      !! if `col` is false, the indices for row `i` are
                                      !!       `ind(k), k = npntr(i),...,npntr(i+1)-1`.
                                      !! ***Note*** that `npntr(n+1)-1` if `col` is true, or `npntr(m+1)-1`
                                      !! if `col` is false, is then the number of non-zero elements
                                      !! of the jacobian matrix.
    integer,dimension(n)  :: Ngrp     !! an integer input array of length `n` which specifies
                                      !! the partition of the columns of the jacobian matrix.
                                      !! column `jcol` belongs to group `ngrp(jcol)`.
    real(wp),dimension(n) :: d        !! an input array of length `n` which contains the
                                      !! difference parameter vector for the estimate of
                                      !! the jacobian matrix columns in group `numgrp`.
    real(wp),dimension(m) :: Fjacd    !! an input array of length `m` which contains
                                      !! an approximation to the difference vector `jac*d`,
                                      !! where `jac` denotes the jacobian matrix.
    real(wp),dimension(*) :: Fjac     !! an output array of length `nnz`, where `nnz` is the
                                      !! number of its non-zero elements. at each call of `fdjs`,
                                      !! `fjac` is updated to include the non-zero elements of the
                                      !! jacobian matrix for those columns in group `numgrp`. `fjac`
                                      !! should not be altered between successive calls to `fdjs`.
    logical,intent(in)    :: Col      !! a logical input variable. if `col` is set true, then the
                                      !! jacobian approximations are stored into a column-oriented
                                      !! pattern. if `col` is set false, then the jacobian
                                      !! approximations are stored into a row-oriented pattern.

    integer :: ip , irow , jcol , jp

    ! compute estimates of jacobian matrix columns in group
    ! numgrp. the array fjacd must contain an approximation
    ! to jac*d, where jac denotes the jacobian matrix and d
    ! is a difference parameter vector with d(jcol) non-zero
    ! if and only if jcol is a column in group numgrp.

    if ( Col ) then  ! column orientation.

        do jcol = 1 , n
            if ( Ngrp(jcol)==Numgrp ) then
                do jp = Npntr(jcol) , Npntr(jcol+1) - 1
                    irow = Ind(jp)
                    Fjac(jp) = Fjacd(irow)/d(jcol)
                enddo
            endif
        enddo

    else  ! row orientation.

        do irow = 1 , m
            do ip = Npntr(irow) , Npntr(irow+1) - 1
                jcol = Ind(ip)
                if ( Ngrp(jcol)==Numgrp ) then
                    Fjac(ip) = Fjacd(irow)/d(jcol)
                    exit
                endif
            enddo
        enddo

    endif

    end subroutine fdjs