TABLE OF CONTENTS


ABINIT/m_eprenorms [ Modules ]

[ Top ] [ Modules ]

NAME

 m_eprenorms

FUNCTION

 This module contains datatypes to compute the renormalization of electronic states due to
 eph coupling and temperature effects

NOTES

 This code is still under development and the API will change in the next versions.
 Contact gmatteo

COPYRIGHT

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

SOURCE

21 #if defined HAVE_CONFIG_H
22 #include "config.h"
23 #endif
24 
25 #include "abi_common.h"
26 
27 module m_eprenorms
28 
29  use defs_basis
30  use m_abicore
31  use m_errors
32  use m_xmpi
33 #ifdef HAVE_NETCDF
34  use netcdf
35 #endif
36  use m_nctk
37 
38  use defs_datatypes, only : ebands_t
39  use m_crystal,      only : crystal_t
40  use m_kpts,         only : listkk
41 
42  implicit none
43 
44  private

m_eprenorms/eprenorms_bcast [ Functions ]

[ Top ] [ m_eprenorms ] [ Functions ]

NAME

 eprenorms_bcast

FUNCTION

  MPI broadcast all the content of eprenorms_t

INPUTS

  master = Rank of master
  comm = MPI communicator

SIDE EFFECTS

  Epren<eprenorms_t> = Data broadcasted on every node from master

SOURCE

259 subroutine eprenorms_bcast(Epren,master,comm)
260 
261 !Arguments -----------------------------------
262 !scalars
263  integer,intent(in) :: master, comm
264  type(eprenorms_t),intent(inout) :: Epren
265 
266 !Local variables------------------------------
267 !scalars
268  integer :: ierr
269 
270 ! ************************************************************************
271 
272  if (xmpi_comm_size(comm) == 1) return
273 
274  call xmpi_bcast(Epren%nkpt, master, comm, ierr)
275  call xmpi_bcast(Epren%nsppol, master, comm, ierr)
276  call xmpi_bcast(Epren%mband, master, comm, ierr)
277  call xmpi_bcast(Epren%ntemp, master, comm, ierr)
278 
279  if (xmpi_comm_rank(comm) /= master) then
280   call eprenorms_init(Epren, Epren%nkpt, Epren%nsppol, Epren%mband, Epren%ntemp)
281  end if
282 
283  call xmpi_bcast(Epren%kpts, master, comm, ierr)
284  call xmpi_bcast(Epren%temps, master, comm, ierr)
285  call xmpi_bcast(Epren%eigens, master, comm, ierr)
286  call xmpi_bcast(Epren%occs, master, comm, ierr)
287  call xmpi_bcast(Epren%renorms, master, comm, ierr)
288  call xmpi_bcast(Epren%linewidth, master, comm, ierr)
289 
290 end subroutine eprenorms_bcast

m_eprenorms/eprenorms_free [ Functions ]

[ Top ] [ m_eprenorms ] [ Functions ]

NAME

 eprenorms_free

FUNCTION

 Deallocate all memory associated with eprenorms

INPUTS

 Epren<eprenorms_t>=The datatype to be freed

SOURCE

166 subroutine eprenorms_free(Epren)
167 
168 !Arguments -----------------------------------
169 !scalars
170  type(eprenorms_t),intent(inout) :: Epren
171 
172 !*********************************************************************
173 
174  ABI_SFREE(Epren%temps)
175  ABI_SFREE(Epren%kpts)
176  ABI_SFREE(Epren%eigens)
177  ABI_SFREE(Epren%occs)
178  ABI_SFREE(Epren%renorms)
179  ABI_SFREE(Epren%linewidth)
180 
181 end subroutine eprenorms_free

m_eprenorms/eprenorms_from_epnc [ Functions ]

[ Top ] [ m_eprenorms ] [ Functions ]

NAME

 eprenorms_from_epnc

FUNCTION

 Allocates and initializes the datatype from a _EP.nc file

INPUTS

 filename = name of the file to be read

SIDE EFFECTS

 Epren<eprenorms_t> = fields are initialized and filled with data from filename

SOURCE

201 subroutine eprenorms_from_epnc(Epren,filename)
202 
203 !Arguments -----------------------------------
204 !scalars
205  character(len=fnlen),intent(in) :: filename
206  type(eprenorms_t),intent(inout) :: Epren
207 
208 !Local variables------------------------------
209  integer :: nkpt, mband, nsppol, ntemp
210 #ifdef HAVE_NETCDF
211  integer :: ncid
212 #endif
213 
214 ! ************************************************************************
215 
216 #ifdef HAVE_NETCDF
217  NCF_CHECK(nctk_open_read(ncid, filename, xmpi_comm_self))
218  NCF_CHECK(nctk_set_datamode(ncid))
219 
220  NCF_CHECK(nctk_get_dim(ncid, "number_of_kpoints", nkpt))
221  NCF_CHECK(nctk_get_dim(ncid, "number_of_spins",nsppol))
222  NCF_CHECK(nctk_get_dim(ncid, "max_number_of_states",mband))
223  NCF_CHECK(nctk_get_dim(ncid, "number_of_temperature",ntemp))
224 
225  call eprenorms_init(Epren, nkpt, nsppol, mband, ntemp)
226 
227  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid,"reduced_coordinates_of_kpoints"), Epren%kpts))
228  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid,"temperature"), Epren%temps))
229  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid,"eigenvalues"), Epren%eigens))
230  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid,"occupations"), Epren%occs))
231  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid,"zero_point_motion"), Epren%renorms))
232  ! TODO: This should be changed. What is stored is a linewidth, not a lifetime,
233  ! we postone the change so as to not break compatibility
234  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid,"lifetime"), Epren%linewidth))
235 
236 #endif
237 
238 end subroutine eprenorms_from_epnc

m_eprenorms/eprenorms_init [ Functions ]

[ Top ] [ m_eprenorms ] [ Functions ]

NAME

 eprenorms_init

FUNCTION

  Initializes an eprenorms_t datatype

INPUTS

OUTPUT

  Epren<eprenorms_t>=Datatype gathering electron-phonon renormalizations

SOURCE

123 subroutine eprenorms_init(Epren,nkpt,nsppol,mband,ntemp)
124 
125 !Arugments -----------------------------------
126 !scalars
127  integer,intent(in) :: nkpt, nsppol, mband, ntemp
128  type(eprenorms_t) :: Epren
129 !arrays
130 
131 !*************************************************************************
132 
133  DBG_ENTER("COLL")
134 
135  Epren%nkpt = nkpt
136  Epren%nsppol = nsppol
137  Epren%mband = mband
138  Epren%ntemp = ntemp
139 
140  ABI_MALLOC(Epren%temps,(Epren%ntemp))
141  ABI_MALLOC(Epren%kpts,(3,Epren%nkpt))
142  ABI_MALLOC(Epren%eigens,(Epren%mband,Epren%nkpt,Epren%nsppol))
143  ABI_MALLOC(Epren%occs,(Epren%mband,Epren%nkpt,Epren%nsppol))
144  ABI_MALLOC(Epren%renorms,(2,Epren%mband,Epren%nkpt,Epren%nsppol,Epren%ntemp))
145  ABI_MALLOC(Epren%linewidth,(2,Epren%mband,Epren%nkpt,Epren%nsppol,Epren%ntemp))
146 
147  DBG_EXIT("COLL")
148 
149 end subroutine eprenorms_init

m_eprenorms/eprenorms_t [ Types ]

[ Top ] [ m_eprenorms ] [ Types ]

NAME

 eprenorms_t

FUNCTION

 Datatype gathering data for electron-phonon renormalization of the band structure

SOURCE

 56  type,public :: eprenorms_t
 57 
 58   !scalars
 59   integer :: nkpt
 60   ! Number of kpoints
 61 
 62   integer :: nsppol
 63   ! Number of spin channels
 64 
 65   integer :: mband
 66   ! Maximum number of bands
 67 
 68   integer :: ntemp
 69   ! Number of temperatures
 70 
 71   !arrays
 72   real(dp), allocatable :: kpts(:,:)
 73   ! kpt(3,nkpt)
 74   ! Kpoints
 75 
 76   real(dp), allocatable :: temps(:)
 77   ! temps(ntemp)
 78   ! Temperatures
 79 
 80   real(dp), allocatable :: eigens(:,:,:)
 81   ! eigens(mband,nkpt,nsppol)
 82   ! Kohn-Sham eigenvalues
 83 
 84   real(dp), allocatable :: occs(:,:,:)
 85   ! occ(mband,nkpt,nsppol)
 86   ! Occupation numbers
 87 
 88   real(dp), allocatable :: renorms(:,:,:,:,:)
 89   ! renorms(2,mband,nkpt,nsppol,ntemp)
 90   ! Renormalization of the eigenvalues for each temperature
 91 
 92   real(dp), allocatable :: linewidth(:,:,:,:,:)
 93   ! linewidth(2,mband,nkpt,nsppol,ntemp)
 94   ! Electron-phonon induced linewidth of the eigens
 95 
 96  end type eprenorms_t
 97 
 98  public :: eprenorms_init
 99  public :: eprenorms_free
100  public :: eprenorms_from_epnc
101  public :: eprenorms_bcast
102 
103  public :: renorm_bst

m_eprenorms/renorm_bst [ Functions ]

[ Top ] [ m_eprenorms ] [ Functions ]

NAME

 renorm_bst

FUNCTION

  Renormalize the band structure Bst from data contained Epren

INPUTS

  Epren<eprenorms_t> = datatype containing the elphon renormalization
  itemp = index of the temperature you want to use
  do_lifetime = .true. if we want to use imaginary eigenvalues (lifetime field)

SIDE EFFECTS

  Bst<bands_t> : eigens are changed according to epren
                 linewidth is allocated and filled with data if do_linewidth

SOURCE

313 subroutine renorm_bst(Epren,Bst,Cryst,itemp,do_lifetime,do_check)
314 
315 !Arguments -----------------------------------
316 !scalars
317  integer :: itemp
318  logical,intent(in) :: do_lifetime
319  logical,optional,intent(in) :: do_check
320  type(eprenorms_t),intent(in) :: Epren
321  type(ebands_t),intent(inout) :: Bst
322  type(crystal_t),intent(in) :: Cryst
323 
324 !Local variables------------------------------
325 !scalars
326  integer :: isppol,ikpt,comm
327  integer :: nband1, nband_tmp
328  integer :: timrev, sppoldbl
329  integer :: ik_eph
330  real(dp) :: dksqmax
331  logical :: check
332 !arrays
333  integer,allocatable :: bs2eph(:,:)
334 
335 ! ************************************************************************
336 
337  ABI_CHECK(Bst%nsppol == Epren%nsppol, "Nsppol should be the same")
338 
339  comm = xmpi_comm_self
340 
341  if(do_lifetime) then
342    ABI_MALLOC(Bst%linewidth,(1,Bst%mband,Bst%nkpt,Bst%nsppol))
343  end if
344 
345  check = .TRUE.
346  if(present(do_check)) then
347    check = do_check
348  end if
349 
350  sppoldbl = 1 !; if (any(Cryst%symafm == -1) .and. Epren%nsppol == 1) nsppoldbl=2
351  ABI_MALLOC(bs2eph, (BSt%nkpt*sppoldbl, 6))
352  timrev = 1
353  call listkk(dksqmax, Cryst%gmet, bs2eph, Epren%kpts, BSt%kptns, Epren%nkpt, Bst%nkpt, Cryst%nsym, &
354 &   sppoldbl, Cryst%symafm, Cryst%symrel, timrev, comm, use_symrec=.False.)
355 
356  do isppol=1,Bst%nsppol
357    do ikpt=1,Bst%nkpt
358      nband1 = Bst%nband(ikpt+(isppol-1)*Bst%nkpt)
359      nband_tmp=MIN(nband1,Epren%mband)
360 
361      ik_eph = bs2eph(ikpt,1)
362 
363      !FIXME change check
364      if (check) then
365        if (ANY(ABS(Bst%eig(1:MIN(10,nband_tmp),ikpt,isppol) - Epren%eigens(1:MIN(10,nband_tmp),ik_eph,isppol)) > tol3)) then
366          !write(stdout,*) "eig : ",BSt%eig(1:MIN(10,nband_tmp),ikpt,isppol)
367          !write(stdout,*) "eigens : ",Epren%eigens(1:MIN(10,nband_tmp),ikpt,isppol)
368          ABI_ERROR("Error in eigenvalues, check the _EP.nc file with respect to your input file !")
369        end if
370        if (ANY(ABS(Bst%occ(1:MIN(10,nband_tmp),ikpt,isppol) - Epren%occs(1:MIN(10,nband_tmp),ik_eph,isppol)) > tol3)) then
371          ABI_ERROR("Error in occupations, check the _EP.nc file with respect to your input file !")
372        end if
373      end if
374 
375      ! Upgrade energies
376      Bst%eig(1:nband_tmp,ikpt,isppol) = BSt%eig(1:nband_tmp,ikpt,isppol) + Epren%renorms(1,1:nband_tmp,ik_eph,isppol,itemp)
377 
378      if (do_lifetime) then
379        Bst%linewidth(1,1:nband_tmp,ikpt,isppol) = Epren%linewidth(1,1:nband_tmp,ik_eph,isppol,itemp)
380      end if
381    end do
382  end do
383 
384  ABI_FREE(bs2eph)
385 
386 end subroutine renorm_bst