TABLE OF CONTENTS


ABINIT/matrginv [ Functions ]

[ Top ] [ Functions ]

NAME

 matrginv

FUNCTION

 Invert a general matrix of real*8 elements.

COPYRIGHT

 Copyright (C) 2001-2017 ABINIT group (GMR)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

 lda=leading dimension of complex matrix a
 n=size of complex matrix a
 a=matrix of real elements

OUTPUT

 a=inverse of a input matrix

SIDE EFFECTS

 a(lda,n)= array of real elements, input, inverted at output

PARENTS

      calc_optical_mels,ddb_elast,ddb_piezo,get_tau_k,linear_optics_paw
      m_chi0,m_haydock,m_vcoul,matpointsym,mka2f_tr,mlwfovlp_ylmfar,setup_bse
      strainsym

CHILDREN

      dbgmdi,dbgmlu,dgeicd,dgetrf,dgetri

SOURCE

 36 #if defined HAVE_CONFIG_H
 37 #include "config.h"
 38 #endif
 39 
 40 #include "abi_common.h"
 41 
 42 subroutine matrginv(a,lda,n)
 43 
 44  use defs_basis
 45  use m_errors
 46  use m_profiling_abi
 47  use m_linalg_interfaces
 48 
 49 !This section has been created automatically by the script Abilint (TD).
 50 !Do not modify the following lines by hand.
 51 #undef ABI_FUNC
 52 #define ABI_FUNC 'matrginv'
 53 !End of the abilint section
 54 
 55  implicit none
 56 
 57 !Arguments ------------------------------------
 58 !scalars
 59  integer,intent(in) :: lda,n
 60 !arrays
 61  real(dp),intent(inout) :: a(lda,n)
 62 
 63 !Local variables-------------------------------
 64 !scalars
 65  integer :: ierr,nwork
 66 #if defined HAVE_LINALG_ESSL
 67  real(dp) :: rcond
 68 #endif
 69  character(len=500) :: message
 70 !arrays
 71  integer,allocatable :: ipvt(:)
 72 #if defined HAVE_LINALG_ESSL
 73  real(dp) :: det(2)
 74 #elif defined HAVE_LINALG_ASL
 75  real(dp) :: det(2)
 76 #endif
 77  real(dp),allocatable :: work(:)
 78 
 79 ! *************************************************************************
 80 
 81 #if defined HAVE_LINALG_ESSL
 82  nwork=200*n
 83 #else
 84  nwork=n
 85 #endif
 86 
 87  ABI_ALLOCATE(work,(nwork))
 88  ABI_ALLOCATE(ipvt,(n))
 89 
 90 
 91 #if defined HAVE_LINALG_ESSL
 92 
 93  call dgeicd(a,lda,n,0,rcond,det,work,nwork)
 94  if(abs(rcond)==zero) then
 95    write(message, '(10a)' ) ch10,&
 96 &   ' matrginv : BUG -',ch10,&
 97 &   '  The matrix that has been passed in argument of this subroutine',ch10,&
 98 &   '  is probably either singular or nearly singular.',ch10,&
 99 &   '  The ESSL routine dgeicd failed.',ch10,&
100 &   '  Action : Contact ABINIT group '
101    MSG_ERROR(message)
102  end if
103 
104 #elif defined HAVE_LINALG_ASL
105 
106  call dbgmlu(a,lda,n,ipvt,ierr)
107  if(ierr /= 0) then
108    write(message, '(10a)' ) ch10,&
109 &   ' matrginv : BUG -',ch10,&
110 &   '  The matrix that has been passed in argument of this subroutine',ch10,&
111 &   '  is probably either singular or nearly singular.',ch10,&
112 &   '  The ASL routine dbgmlu failed.',ch10,&
113 &   '  Action : Contact ABINIT group '
114    MSG_ERROR(message)
115  end if
116  call dbgmdi(a,lda,n,ipvt,det,-1,work,ierr)
117  if(ierr /= 0) then
118    write(message, '(10a)' ) ch10,&
119 &   ' matrginv : BUG -',ch10,&
120 &   '  The matrix that has been passed in argument of this subroutine',ch10,&
121 &   '  is probably either singular or nearly singular.',ch10,&
122 &   '  The ASL routine dbgmdi failed.',ch10,&
123 &   '  Action : Contact ABINIT group '
124    MSG_ERROR(message)
125  end if
126 
127 #else
128 
129  call dgetrf(n,n,a,lda,ipvt,ierr)
130  if(ierr /= 0) then
131    write(message, '(10a)' ) ch10,&
132 &   ' matrginv : BUG -',ch10,&
133 &   '  The matrix that has been passed in argument of this subroutine',ch10,&
134 &   '  is probably either singular or nearly singular.',ch10,&
135 &   '  The LAPACK routine dgetrf failed.',ch10,&
136 &   '  Action : Contact ABINIT group '
137    MSG_ERROR(message)
138  end if
139  call dgetri(n,a,lda,ipvt,work,n,ierr)
140  if(ierr /= 0) then
141    write(message, '(10a)' ) ch10,&
142 &   ' matrginv : BUG -',ch10,&
143 &   '  The matrix that has been passed in argument of this subroutine',ch10,&
144 &   '  is probably either singular or nearly singular.',ch10,&
145 &   '  The LAPACK routine dgetri failed.',ch10,&
146 &   '  Action : Contact ABINIT group '
147    MSG_ERROR(message)
148  end if
149 
150 #endif
151 
152  ABI_DEALLOCATE(work)
153  ABI_DEALLOCATE(ipvt)
154 
155 end subroutine matrginv