The Elk Code
zminv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine zminv(n,a)
7 use modomp
8 implicit none
9 ! arguments
10 integer, intent(in) :: n
11 complex(8), intent(inout) :: a(n,n)
12 ! local variables
13 integer info,nthd,nts
14 ! automatic arrays
15 integer ipiv(n)
16 complex(8) work(n)
17 ! enable MKL parallelism
18 call holdthd(maxthdmkl,nthd)
19 nts=mkl_set_num_threads_local(nthd)
20 call zgetrf(n,n,a,n,ipiv,info)
21 if (info /= 0) then
22  write(*,*)
23  write(*,'("Error(zminv): unable to invert matrix")')
24  write(*,'(" ZGETRF returned INFO = ",I8)') info
25  write(*,*)
26  stop
27 end if
28 call zgetri(n,a,n,ipiv,work,n,info)
29 if (info /= 0) then
30  write(*,*)
31  write(*,'("Error(zminv): unable to invert matrix")')
32  write(*,'(" ZGETRI returned INFO = ",I8)') info
33  write(*,*)
34  stop
35 end if
36 nts=mkl_set_num_threads_local(0)
37 call freethd(nthd)
38 end subroutine
39 
Definition: modomp.f90:6
integer maxthdmkl
Definition: modomp.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine zminv(n, a)
Definition: zminv.f90:7