TABLE OF CONTENTS


ABINIT/diag_occ [ Functions ]

[ Top ] [ Functions ]

NAME

 diag_occ

FUNCTION

 Use for DMFT in KGB parallelisation. Diagonalise the occupation matrix
 and return diagonalised occupations and associated eigen vectors sorted
 with descending occupation

INPUTS

   occ_nd_cpx(nband, nband) = matrix of non diagonal occupations for DMFT
   nband = number of band to be processed

OUTPUT

   occ_diag(2, nband) = diagonal occupation in the new band space
   cwavef_rot(2, npw, nband) = fourier coefficient of wave functions for all bands rotated in the band space

SOURCE

 59 !! TODO /!\ No parallel computing yet !
 60 !! TODO add the possibility of using ScaLAPACK to do computation in parallel
 61 
 62 subroutine diag_occ(occ_nd_cpx, nband, occ_diag)
 63 
 64 !Arguments ------------------------------------
 65 !scalars
 66   integer,intent(in) :: nband
 67 !! type(MPI_type),intent(inout) :: mpi_enreg
 68 !! type(dataset_type),intent(in) :: dtset
 69 !! type(paw_dmft_type), intent(in)  :: paw_dmft
 70 !no_abirules
 71   complex(kind=dpc), intent(inout) :: occ_nd_cpx(nband, nband)
 72   real(kind=dp), intent(out) :: occ_diag(nband)
 73 
 74 !Local variables-------------------------------
 75 
 76 !scalars
 77   integer :: info, lwork
 78   character(len=500) :: message
 79 
 80 !arrays
 81   real(kind=dp) :: rwork(3*nband-1)
 82   complex(kind=dpc), allocatable :: work(:)
 83 
 84 ! *************************************************************************
 85 
 86   DBG_ENTER("COLL")
 87 
 88 !! Initialisation
 89   rwork = zero
 90   info = 0
 91 
 92 !! use the opposite to have zheev orders the eigenvalues descending (once the
 93 !! opposite have been taken again)
 94   occ_nd_cpx = -occ_nd_cpx
 95 
 96 !! Get diagonal occupations and associeted base
 97 
 98 ! Compute the optimal working array size
 99   ABI_MALLOC(work,(1))
100   work = czero
101   call zheev('V', 'U', nband, occ_nd_cpx, nband, occ_diag, work, -1, rwork, info)
102   lwork = int(work(1))
103   ABI_FREE(work)
104 
105 ! Compute the eigenvalues (occ_diag) and vectors
106   ABI_MALLOC(work,(lwork))
107   work = czero
108 
109   call zheev('V', 'U', nband, occ_nd_cpx, nband, occ_diag, work, lwork, rwork, info)
110 
111 !! obtain the true eigen values of occupation matrix in descending order
112   occ_diag = -occ_diag
113 
114   ABI_FREE(work)
115 
116   if (info > 0) then
117     message=""
118     write(message, "(a,i5)") " something wrong happened with the diagonalisation of &
119 &the occupation matrix (did't converge), info=",info
120     ABI_ERROR(message)
121   else if (info < 0) then
122     message=""
123     write(message, "(a,i5)") " something wrong happened with the diagonalisation of &
124 &the occupation matrix (bad input argument), info=",info
125     ABI_ERROR(message)
126   end if
127 
128   DBG_EXIT("COLL")
129 
130 end subroutine diag_occ

ABINIT/m_rot_cg [ Functions ]

[ Top ] [ Functions ]

NAME

 m_rot_cg

FUNCTION

  Rotate the cg coefficient with the rotation matrix obtained from the
  diagonalisation of the non-diagonal occupation matrix produced by DMFT.

INPUTS

OUTPUT

SOURCE

15 !! TODO /!\ No parallel computing yet !
16 
17 #if defined HAVE_CONFIG_H
18 #include "config.h"
19 #endif
20 
21 #include "abi_common.h"
22 
23 module m_rot_cg
24 
25   use defs_basis
26   use m_xmpi
27   use m_errors
28   use m_abicore
29 
30   implicit none
31 
32   private
33 
34   public :: rot_cg
35 
36   contains

ABINIT/rot_cg [ Functions ]

[ Top ] [ Functions ]

NAME

 rot_cg

FUNCTION

 Use for DMFT in KGB parallelisation. Diagonalise the occupation matrix
 and use the resulting base to represent the wave functions.

INPUTS

   occ_nd(2, nband, nband) = matrix of non diagonal occupations for DMFT
   cwavef(2, npw, nband) = fourier coefficient of wave functions for all bands
   npw = number of G vectors computed in this iteration
   nband = number of band to be processed
   blocksize = size of the block for th LO.. algorithm
               still have to be equal to nband
   nspinor = number of spinor components
   first_bandc = index of the first correlated band
   nbandc = number of bands correlated

OUTPUT

   occ_diag(2, nband) = diagonal occupation in the new band space

 SIDE EFFECT
   cwavef is rotated with the unitary matrix obtained from the diagonalisation
   of occupations (occ_nd)

SOURCE

161 !! TODO /!\ No parallel computing yet !
162 !! TODO add the possibility of using ScaLAPACK to do computation in parallel
163 !! TODO Make the computation of the new wf parallel
164 
165 subroutine rot_cg(occ_nd, cwavef, npw, nband, blocksize, nspinor, first_bandc, nbandc, occ_diag)
166 
167 !Arguments ------------------------------------
168 !scalars
169   integer,intent(in) :: npw, nband, blocksize, nspinor, first_bandc, nbandc
170 !! type(MPI_type),intent(inout) :: mpi_enreg
171 !! type(dataset_type),intent(in) :: dtset
172 !! type(paw_dmft_type), intent(in)  :: band_in
173 !no_abirules
174   real(kind=dp), intent(in) :: occ_nd(2, blocksize, blocksize)
175   real(kind=dp), intent(inout) :: cwavef(2, npw, blocksize, nspinor)
176   real(kind=dp), intent(out) :: occ_diag(blocksize)
177 
178 !Local variables-------------------------------
179 
180 !scalars
181   integer :: n, np, ig
182   character(len=500) :: message
183 
184 !arrays
185   real(kind=dp) :: occ_diag_red(nbandc)
186   complex(kind=dpc) :: occ_nd_cpx(nbandc, nbandc)
187   complex(kind=dpc) :: cwavef_rot_g(nbandc, nspinor)
188 
189 ! *************************************************************************
190 
191   DBG_ENTER("COLL")
192 
193   if(nband /= blocksize) then
194     message = " DMFT in KGB cannot be used with multiple blocks yet. Make sure that bandpp*npband = nband."
195     ABI_ERROR(message)
196   end if
197 
198 !! Initialisation
199 
200   do n=1,nbandc
201     do np=1,nbandc
202       occ_nd_cpx(n,np) = cmplx(occ_nd(1,n+first_bandc-1,np+first_bandc-1), occ_nd(2,n+first_bandc-1,np+first_bandc-1), kind=dpc)
203     end do
204   end do
205 
206 !! Get diagonal occupations and associeted base
207 
208   call diag_occ(occ_nd_cpx, nbandc, occ_diag_red)
209 
210   do n=1,nband
211     if (n < first_bandc .or. n >= first_bandc+nbandc) then
212       occ_diag(n) = occ_nd(1, n, n)
213     else
214       occ_diag(n) = occ_diag_red(n-first_bandc+1)
215     end if
216   end do
217 
218 !! Compute the corresponding wave functions if nothing wrong happened
219   ! $c^{rot}_{n,k}(g) =  \sum_{n'} [\bar{f_{n',n}} * c_{n',k}(g)]$
220   do ig=1,npw
221     cwavef_rot_g(:,:) = czero
222     do n=1,nbandc
223       do np=1,nbandc
224         cwavef_rot_g(n,:) = cwavef_rot_g(n,:) + occ_nd_cpx(np, n) * &
225 &                           cmplx(cwavef(1,ig,np+first_bandc-1,:), cwavef(2,ig,np+first_bandc-1,:), kind=dpc)
226       end do
227     end do
228     cwavef(1,ig,first_bandc:first_bandc+nbandc-1,:) = dreal(cwavef_rot_g)
229     cwavef(2,ig,first_bandc:first_bandc+nbandc-1,:) = dimag(cwavef_rot_g)
230   end do
231 
232   DBG_EXIT("COLL")
233 
234  end subroutine rot_cg
235 
236 end module m_rot_cg