TABLE OF CONTENTS


ABINIT/m_haydock_io [ Modules ]

[ Top ] [ Modules ]

NAME

 m_haydock_io

FUNCTION

  This module provides routines to read the Haydock file

COPYRIGHT

  Copyright (C) 2013-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

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 MODULE m_haydock_io
23 
24  use defs_basis
25  use m_abicore
26  use m_errors
27 
28  use m_io_tools,       only : open_file
29 
30  implicit none
31 
32  private

m_haydock_io/close_haydock [ Functions ]

[ Top ] [ m_haydock_io ] [ Functions ]

NAME

 close_haydock

FUNCTION

  Closes the haydock file and free the memory related to the file descriptor

INPUTS

  haydock_file = haydock file descriptor

SOURCE

405 subroutine close_haydock(haydock_file)
406 
407 !Arguments ------------------------------------
408  type(haydock_type),intent(inout) :: haydock_file
409 
410 ! *************************************************************************
411 
412  close(haydock_file%unt)
413 
414  ABI_SFREE(haydock_file%qpoints)
415  ABI_SFREE(haydock_file%niter)
416 
417 end subroutine close_haydock

m_haydock_io/haydock_type [ Types ]

[ Top ] [ m_haydock_io ] [ Types ]

NAME

  haydock_type

FUNCTION

  The structure defining the content of the haydock file

SOURCE

46  type,public :: haydock_type
47 
48    integer :: version
49    ! Version of the file
50 
51    integer :: hsize
52    ! Size of the hamiltonian
53 
54    integer :: use_coupling
55    ! 1 if coupling is used, 0 if not
56 
57    integer :: op
58    ! Type of the file
59 
60    integer :: nq
61    ! Number of q-points in the file
62 
63    integer :: unt
64    ! Fortran unit number.
65 
66    real(dp) :: broad
67    ! Broadening factor
68 
69    complex(dpc) :: factor
70 
71    real(dp),allocatable :: qpoints(:,:)
72    ! qpoints(3,nq)
73 
74    integer,allocatable :: niter(:)
75    ! niter(nq)
76 
77  end type haydock_type
78 
79  public :: open_haydock             ! Init the Haydock file
80  public :: read_haydock             ! Reads the data associated with 1 q-point
81  public :: read_dim_haydock         ! Read dimensions
82  public :: write_dim_haydock        ! Write dimensions.
83  public :: skip_dim_haydock         ! Skip the section with dimensions.
84  public :: write_haydock            ! Writes data related to a q-point inside the file
85  public :: close_haydock            ! Close the Haydock file and release memory

m_haydock_io/open_haydock [ Functions ]

[ Top ] [ m_haydock_io ] [ Functions ]

NAME

 open_haydock

FUNCTION

 Open the file and initialize haydock file descriptor that will
 be used for later access

INPUTS

  filename = Name of the file to be opened

OUTPUT

  haydock_file = file descriptor for the haydock file

SOURCE

107 subroutine open_haydock(filename, haydock_file)
108 
109 !Arguments ------------------------------------
110 !scalars
111  character(len=*),intent(in) :: filename
112  type(haydock_type),intent(out) :: haydock_file
113 !arrays
114 
115 !Local variables ------------------------------
116 !scalars
117  character(len=500) :: msg
118 
119 !************************************************************************
120 
121  if (open_file(filename,msg,newunit=haydock_file%unt,form="unformatted") /= 0) then
122    ABI_ERROR(msg)
123  end if
124 
125  haydock_file%version = CUR_VERSION
126 
127 end subroutine open_haydock

m_haydock_io/read_dim_haydock [ Functions ]

[ Top ] [ m_haydock_io ] [ Functions ]

NAME

 read_dim_haydock

FUNCTION

  Reads the header dimensions of the haydock file and store them in
  the haydock file descriptor

 INPUT/OUTPUT
  haydock_file = haydock file descriptor

SOURCE

145 subroutine read_dim_haydock(haydock_file)
146 
147 !Arguments ------------------------------------
148 !scalars
149  type(haydock_type),intent(inout) :: haydock_file
150 
151 !Local variables ------------------------------
152 !scalars
153  integer :: niter_file
154  integer :: iq, it
155  character(len=500) :: msg
156 !arrays
157  real(dp) :: q_file(3)
158 
159 !************************************************************************
160 
161  read(haydock_file%unt) haydock_file%version
162 
163  msg = "The haydock file has been produced by an uncompatible version of abinit"
164  ABI_CHECK(haydock_file%version == CUR_VERSION, msg)
165 
166  read(haydock_file%unt) haydock_file%hsize,haydock_file%use_coupling,&
167 &   haydock_file%op,haydock_file%nq,haydock_file%broad
168 
169  ABI_MALLOC(haydock_file%qpoints,(3,haydock_file%nq))
170  ABI_MALLOC(haydock_file%niter,(haydock_file%nq))
171 
172  do iq = 1,haydock_file%nq
173    read(haydock_file%unt) q_file(:)
174    read(haydock_file%unt) niter_file
175 
176    haydock_file%qpoints(:,iq) = q_file(:)
177    haydock_file%niter(iq) = niter_file
178    ! Skip data for this q.
179    do it=1,niter_file
180      read(haydock_file%unt) ! it,aa(it),bb(it)
181    end do
182    read(haydock_file%unt) ! phi_nm1
183    read(haydock_file%unt) ! phi_n
184    read(haydock_file%unt) ! factor
185  end do
186 
187 end subroutine read_dim_haydock

m_haydock_io/read_haydock [ Functions ]

[ Top ] [ m_haydock_io ] [ Functions ]

NAME

 read_haydock

FUNCTION

  Reads the data related to one q-point inside the file accessed by file descriptor

INPUTS

  haydock_file = haydock file descriptor
  q = q-point to be searched inside the file

OUTPUT

  aa = coefficients "a" of the lanczos chain
  bb = coefficients "b" of the lanczos chain
  phi_n = last vector of the lanczos chain
  phi_nm1 = penultimate vector of the lanczos chain
  niter = number of iterations done
  factor = pre-factor used to obtain the green function

NOTES

  niter = 0 if the q-point has not been found

SOURCE

274 subroutine read_haydock(haydock_file, q, aa, bb, phi_n, phi_nm1, niter, factor)
275 
276 !Arguments ------------------------------------
277 !scalars
278  integer,intent(out) :: niter
279  complex(dpc),intent(out) :: factor
280  type(haydock_type),intent(in) :: haydock_file
281 !arrays
282  real(dp),intent(in) :: q(3)
283  real(dp),allocatable,intent(out) :: bb(:)
284  complex(dpc),allocatable,intent(out) :: aa(:),phi_n(:),phi_nm1(:)
285 
286 !Local variables ------------------------------
287 !scalars
288  integer :: iq, it, inn
289  integer :: niter_file
290  logical :: found_q
291  character(len=500) :: msg
292 !arrays
293  real(dp) :: q_file(3)
294 
295 ! *************************************************************************
296 
297  rewind(haydock_file%unt)
298  call skip_dim_haydock(haydock_file)
299 
300  do iq = 1,haydock_file%nq
301    read(haydock_file%unt) q_file(:)
302    read(haydock_file%unt) niter_file
303 
304    if ( ALL(ABS(q_file - q) < tol6) ) then
305       found_q = .TRUE.; EXIT
306    else
307       ! Skip data for this q.
308       do it=1,niter_file
309         read(haydock_file%unt) ! it,aa(it),bb(it)
310       end do
311       read(haydock_file%unt) ! phi_nm1
312       read(haydock_file%unt) ! phi_n
313       read(haydock_file%unt) ! factor
314    end if
315  end do
316 
317  if(found_q) then
318    niter = niter_file
319    ABI_MALLOC(aa,(niter))
320    ABI_MALLOC(bb,(niter))
321    do inn=1,niter
322      read(haydock_file%unt)it,aa(inn),bb(inn)
323      if (inn/=it) then
324        write(msg,'(2(a,i0))')" Found it_file: ",it," while it should be: ",inn
325        ABI_ERROR(msg)
326      end if
327    end do
328    ABI_MALLOC(phi_nm1,(haydock_file%hsize))
329    ABI_MALLOC(phi_n,(haydock_file%hsize))
330    read(haydock_file%unt)phi_nm1
331    read(haydock_file%unt)phi_n
332    read(haydock_file%unt)factor
333  else
334    niter = 0
335  end if
336 
337 end subroutine read_haydock

m_haydock_io/skip_dim_haydock [ Functions ]

[ Top ] [ m_haydock_io ] [ Functions ]

NAME

 skip_dim_haydock

FUNCTION

  Skip the part of the file reading basic dimensions contained in the header

INPUTS

  haydock_file = haydock file descriptor

OUTPUT

SOURCE

234 subroutine skip_dim_haydock(haydock_file)
235 
236 !Arguments ------------------------------------
237 !scalars
238  type(haydock_type),intent(in) :: haydock_file
239 
240 ! *************************************************************************
241 
242  read(haydock_file%unt)
243  read(haydock_file%unt)
244 
245 end subroutine skip_dim_haydock

m_haydock_io/write_dim_haydock [ Functions ]

[ Top ] [ m_haydock_io ] [ Functions ]

NAME

 write_dim_haydock

FUNCTION

  Writes the basic dimensions stored inside the haydock descriptor inside the file

INPUTS

  haydock_file = haydock file descriptor

SOURCE

204 subroutine write_dim_haydock(haydock_file)
205 
206 !Arguments ------------------------------------
207  type(haydock_type),intent(in) :: haydock_file
208 
209 ! *************************************************************************
210 
211  write(haydock_file%unt) haydock_file%version
212  write(haydock_file%unt) haydock_file%hsize,haydock_file%use_coupling, &
213 &  haydock_file%op,haydock_file%nq,haydock_file%broad
214 
215 end subroutine write_dim_haydock

m_haydock_io/write_haydock [ Functions ]

[ Top ] [ m_haydock_io ] [ Functions ]

NAME

 write_haydock

FUNCTION

  Writes data related to a q-point inside the file

INPUTS

  haydock_file = haydock file descriptor
  q = q-point to be searched inside the file
  aa = coefficients "a" of the lanczos chain
  bb = coefficients "b" of the lanczos chain
  phi_n = last vector of the lanczos chain
  phi_nm1 = penultimate vector of the lanczos chain
  niter = number of iterations done
  factor = pre-factor used to obtain the green function

SOURCE

361 subroutine write_haydock(haydock_file, hsize, q, aa, bb, phi_n, phi_nm1, niter, factor)
362 
363 !Arguments ------------------------------------
364 !scalars
365  integer,intent(in) :: niter,hsize
366  complex(dpc),intent(in) :: factor
367  type(haydock_type),intent(in) :: haydock_file
368 !arrays
369  real(dp),intent(in) :: q(3)
370  real(dp),intent(in) :: bb(niter)
371  complex(dpc),intent(in) :: aa(niter),phi_n(hsize),phi_nm1(hsize)
372 
373 !Local variables -----------------------------
374 !scalars
375  integer :: it
376 
377 ! *************************************************************************
378 
379  write(haydock_file%unt) q
380  write(haydock_file%unt) niter  ! NB if the previous loop completed inn=niter_max+1
381  do it=1,niter        ! if we exited then inn is not incremented by one.
382    write(haydock_file%unt)it,aa(it),bb(it)
383  end do
384  write(haydock_file%unt) phi_nm1
385  write(haydock_file%unt) phi_n
386  write(haydock_file%unt) factor
387 
388 end subroutine write_haydock