dgefa Subroutine

public subroutine dgefa(A, Lda, N, Ipvt, Info)

dgefa(3f) - [M_odepack::matrix] Factor a matrix using Gaussian elimination.

  subroutine dgefa(A,Lda,N,Ipvt,Info)
  real(kind=dp),intent(inout) :: A(Lda,*)
  integer,intent(in)          :: Lda
  integer,intent(in)          :: N
  integer,intent(out)         :: Ipvt(*)
  integer,intent(out)         :: Info

DESCRIPTION

DGEFA factors a double precision matrix by Gaussian elimination.

DGEFA is usually called by DGECO, but it can be called directly with a saving in time if RCOND is not needed. (Time for DGECO) = (1 + 9/N)*(Time for DGEFA) .

OPTIONS

A

the matrix to be factored.

LDA

the leading dimension of the array A .

N

the order of the matrix A .

RETURNS

A

an upper triangular matrix and the multipliers which were used to obtain it. The factorization can be written A = L*U where

L

is a product of permutation and unit lower triangular matrices and U is upper triangular.

IPVT

an integer vector of pivot indices.

INFO

: = 0 normal value. = K if U(K,K) .EQ. 0.0 . This is not an error condition for this subroutine, but it does indicate that DGESL or DGEDI will divide by zero if called. Use RCOND in DGECO for a reliable indication of singularity.

REFERENCES

J. J. Dongarra, J. R. Bunch, C. B. Moler, and G. W. Stewart, LINPACK Users’ Guide, SIAM, 1979.

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout) :: A(Lda,*)
integer, intent(in) :: Lda
integer, intent(in) :: N
integer, intent(out) :: Ipvt(*)
integer, intent(out) :: Info

Calls

proc~~dgefa~2~~CallsGraph proc~dgefa~2 M_odepack::dgefa proc~daxpy~2 M_odepack::daxpy proc~dgefa~2->proc~daxpy~2 proc~dscal M_odepack::dscal proc~dgefa~2->proc~dscal

Called by

proc~~dgefa~2~~CalledByGraph proc~dgefa~2 M_odepack::dgefa proc~daigbt~2 M_odepack::daigbt proc~daigbt~2->proc~dgefa~2 proc~dainvg~2 M_odepack::dainvg proc~dainvg~2->proc~dgefa~2 proc~dpjibt~2 M_odepack::dpjibt proc~dpjibt~2->proc~dgefa~2 proc~dprepji~2 M_odepack::dprepji proc~dprepji~2->proc~dgefa~2 proc~dprepj~2 M_odepack::dprepj proc~dprepj~2->proc~dgefa~2 proc~dprja~2 M_odepack::dprja proc~dprja~2->proc~dgefa~2 proc~dlsodi~2 M_odepack::dlsodi proc~dlsodi~2->proc~dainvg~2 proc~dlsoibt~2 M_odepack::dlsoibt proc~dlsoibt~2->proc~daigbt~2