lammps/lib/linalg/dpotrf2.f

238 lines
6.2 KiB
Fortran

*> \brief \b DPOTRF2
*
* =========== DOCUMENTATION ===========
*
* Online html documentation available at
* http://www.netlib.org/lapack/explore-html/
*
* Definition:
* ===========
*
* RECURSIVE SUBROUTINE DPOTRF2( UPLO, N, A, LDA, INFO )
*
* .. Scalar Arguments ..
* CHARACTER UPLO
* INTEGER INFO, LDA, N
* ..
* .. Array Arguments ..
* REAL A( LDA, * )
* ..
*
*
*> \par Purpose:
* =============
*>
*> \verbatim
*>
*> DPOTRF2 computes the Cholesky factorization of a real symmetric
*> positive definite matrix A using the recursive algorithm.
*>
*> The factorization has the form
*> A = U**T * U, if UPLO = 'U', or
*> A = L * L**T, if UPLO = 'L',
*> where U is an upper triangular matrix and L is lower triangular.
*>
*> This is the recursive version of the algorithm. It divides
*> the matrix into four submatrices:
*>
*> [ A11 | A12 ] where A11 is n1 by n1 and A22 is n2 by n2
*> A = [ -----|----- ] with n1 = n/2
*> [ A21 | A22 ] n2 = n-n1
*>
*> The subroutine calls itself to factor A11. Update and scale A21
*> or A12, update A22 then calls itself to factor A22.
*>
*> \endverbatim
*
* Arguments:
* ==========
*
*> \param[in] UPLO
*> \verbatim
*> UPLO is CHARACTER*1
*> = 'U': Upper triangle of A is stored;
*> = 'L': Lower triangle of A is stored.
*> \endverbatim
*>
*> \param[in] N
*> \verbatim
*> N is INTEGER
*> The order of the matrix A. N >= 0.
*> \endverbatim
*>
*> \param[in,out] A
*> \verbatim
*> A is DOUBLE PRECISION array, dimension (LDA,N)
*> On entry, the symmetric matrix A. If UPLO = 'U', the leading
*> N-by-N upper triangular part of A contains the upper
*> triangular part of the matrix A, and the strictly lower
*> triangular part of A is not referenced. If UPLO = 'L', the
*> leading N-by-N lower triangular part of A contains the lower
*> triangular part of the matrix A, and the strictly upper
*> triangular part of A is not referenced.
*>
*> On exit, if INFO = 0, the factor U or L from the Cholesky
*> factorization A = U**T*U or A = L*L**T.
*> \endverbatim
*>
*> \param[in] LDA
*> \verbatim
*> LDA is INTEGER
*> The leading dimension of the array A. LDA >= max(1,N).
*> \endverbatim
*>
*> \param[out] INFO
*> \verbatim
*> INFO is INTEGER
*> = 0: successful exit
*> < 0: if INFO = -i, the i-th argument had an illegal value
*> > 0: if INFO = i, the leading minor of order i is not
*> positive definite, and the factorization could not be
*> completed.
*> \endverbatim
*
* Authors:
* ========
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Univ. of Colorado Denver
*> \author NAG Ltd.
*
*> \date December 2016
*
*> \ingroup doublePOcomputational
*
* =====================================================================
RECURSIVE SUBROUTINE DPOTRF2( UPLO, N, A, LDA, INFO )
*
* -- LAPACK computational routine (version 3.7.0) --
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
* December 2016
*
* .. Scalar Arguments ..
CHARACTER UPLO
INTEGER INFO, LDA, N
* ..
* .. Array Arguments ..
DOUBLE PRECISION A( LDA, * )
* ..
*
* =====================================================================
*
* .. Parameters ..
DOUBLE PRECISION ONE, ZERO
PARAMETER ( ONE = 1.0D+0, ZERO = 0.0D+0 )
* ..
* .. Local Scalars ..
LOGICAL UPPER
INTEGER N1, N2, IINFO
* ..
* .. External Functions ..
LOGICAL LSAME, DISNAN
EXTERNAL LSAME, DISNAN
* ..
* .. External Subroutines ..
EXTERNAL DSYRK, DTRSM, XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, SQRT
* ..
* .. Executable Statements ..
*
* Test the input parameters
*
INFO = 0
UPPER = LSAME( UPLO, 'U' )
IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
INFO = -1
ELSE IF( N.LT.0 ) THEN
INFO = -2
ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
INFO = -4
END IF
IF( INFO.NE.0 ) THEN
CALL XERBLA( 'DPOTRF2', -INFO )
RETURN
END IF
*
* Quick return if possible
*
IF( N.EQ.0 )
$ RETURN
*
* N=1 case
*
IF( N.EQ.1 ) THEN
*
* Test for non-positive-definiteness
*
IF( A( 1, 1 ).LE.ZERO.OR.DISNAN( A( 1, 1 ) ) ) THEN
INFO = 1
RETURN
END IF
*
* Factor
*
A( 1, 1 ) = SQRT( A( 1, 1 ) )
*
* Use recursive code
*
ELSE
N1 = N/2
N2 = N-N1
*
* Factor A11
*
CALL DPOTRF2( UPLO, N1, A( 1, 1 ), LDA, IINFO )
IF ( IINFO.NE.0 ) THEN
INFO = IINFO
RETURN
END IF
*
* Compute the Cholesky factorization A = U**T*U
*
IF( UPPER ) THEN
*
* Update and scale A12
*
CALL DTRSM( 'L', 'U', 'T', 'N', N1, N2, ONE,
$ A( 1, 1 ), LDA, A( 1, N1+1 ), LDA )
*
* Update and factor A22
*
CALL DSYRK( UPLO, 'T', N2, N1, -ONE, A( 1, N1+1 ), LDA,
$ ONE, A( N1+1, N1+1 ), LDA )
CALL DPOTRF2( UPLO, N2, A( N1+1, N1+1 ), LDA, IINFO )
IF ( IINFO.NE.0 ) THEN
INFO = IINFO + N1
RETURN
END IF
*
* Compute the Cholesky factorization A = L*L**T
*
ELSE
*
* Update and scale A21
*
CALL DTRSM( 'R', 'L', 'T', 'N', N2, N1, ONE,
$ A( 1, 1 ), LDA, A( N1+1, 1 ), LDA )
*
* Update and factor A22
*
CALL DSYRK( UPLO, 'N', N2, N1, -ONE, A( N1+1, 1 ), LDA,
$ ONE, A( N1+1, N1+1 ), LDA )
CALL DPOTRF2( UPLO, N2, A( N1+1, N1+1 ), LDA, IINFO )
IF ( IINFO.NE.0 ) THEN
INFO = IINFO + N1
RETURN
END IF
END IF
END IF
RETURN
*
* End of DPOTRF2
*
END