1!=============================================================================
2!
3!  Utilities:
4!
5!    epsmat_hdf5_upgrade  Originally by FHJ    Last Modified 08/2016 (FHJ)
6!
7!     This utility upgrades an epsmat.h5 file from BerkeleyGW-1.1-beta
8!     (version 1) to the version used in BerkeleyGW-1.2.0 (version 3).
9!
10!=============================================================================
11
12#include "f_defs.h"
13
14program epsmat_hdf5_upgrade
15
16  use global_m
17  use hdf5
18  use hdf5_io_m
19  use write_matrix_m
20  use epswrite_hdf5_m
21  use wfn_rho_vxc_io_m
22  use wfn_io_hdf5_m
23  use input_utils_m
24
25  implicit none
26
27  character(len=128) :: fname_eps, fname_wfn, tmpstr
28  integer :: nspin, ng
29  integer :: version, intinfo(9)
30  integer :: iw, narg, error
31  integer, allocatable :: gvecs(:,:), isorts(:,:,:)
32  logical, allocatable :: qpt_done(:)
33  real(DP), allocatable :: freqs_tmp(:,:), freqgrid(:)
34  logical :: skip_wfn, mf_header_exists, nspin_exists, exists
35  integer(HID_T) :: file_id
36  type(mf_header_t) :: mf
37  type(polarizability) :: pol
38
39  narg = iargc()
40  if (narg/=4) then
41    write(0,*)
42    write(0,*) 'Usage: epsmat_hdf5_upgrade.flavor.x epsmat.h5 WFN icutv [itruncval]'
43    write(0,*) '==================================================================='
44    write(0,*)
45    write(0,*) 'This utility will upgrade an epsmat.h5 file from version 1'
46    write(0,*) '(generated prior to BerkeleyGW 1.2.0)'
47    write(0,*)
48    write(0,*) 'Required arguments:'
49    write(0,*) '-------------------'
50    write(0,*) '  - epsmat.h5: epsmat HDF5 file in previous (beta) version'
51    write(0,*) '  - WFN: WFN file from which we`ll extract the header'
52    write(0,*) '      You may also specify a single dash (-) as an argument, but only'
53    write(0,*) '      if you generated your epsmat.h5 with BerkeleyGW 1.1-beta2, and not'
54    write(0,*) '      BerkeleyGW 1.1-beta!'
55    write(0,*) '  - icutv: truncation flag used in the calculation. Options are:'
56    write(0,*) '      0 -> no truncation'
57    write(0,*) '      2 -> spherical_truncation'
58    write(0,*) '      4 -> cell_wire_truncation'
59    write(0,*) '      5 -> cell_box_truncation'
60    write(0,*) '      6 -> cell_slab_truncation'
61    write(0,*) '  - itruncval: truncation radius, only needed if icutv==2.'
62    write(0,*)
63    write(0,*) 'Warning:'
64    write(0,*) '--------'
65    write(0,*) '  The upgrade process happens in place and modifies the source file.'
66    write(0,*) '  The process is instantaneous, but because there is always a chance that'
67    write(0,*) '  something will go wrong, BACK UP YOUR EPSMAT.H5 FILE FIRST!'
68    write(0,*)
69    stop
70  endif
71  call getarg(1, fname_eps)
72  call getarg(2, fname_wfn)
73  call getarg(3, tmpstr)
74  read(tmpstr,*) pol%icutv
75  call getarg(4, tmpstr)
76  read(tmpstr,*) pol%truncval(1)
77
78  write(6,*)
79  write(6,*) 'Epsmat HDF5 Upgrade Utility'
80  write(6,*) '==========================='
81  write(6,*)
82  call h5open_f(error)
83
84  skip_wfn = TRUNC(fname_wfn)=="-"
85  ! FHJ: Read WFN file, get the header
86  if (skip_wfn) then
87    write(6,*) 'We will not read any WFN file.'
88  else
89    write(6,*) 'Reading WFN file '//TRUNC(fname_wfn)//'.'
90      call open_file(unit=25,file=TRUNC(fname_wfn),form='unformatted',status='old')
91      call read_mf_header(25, mf, sheader='WFN', iflavor=SCALARSIZE, &
92        warn=.false., dont_warn_kgrid=.true.)
93      SAFE_ALLOCATE(mf%gvec%components, (3, mf%gvec%ng))
94      call read_binary_gvectors(25, mf%gvec%ng, mf%gvec%ng, mf%gvec%components)
95      call close_file(25)
96  endif
97
98  write(6,*) 'Reading epsmat file '//TRUNC(fname_eps)//'.'
99  call open_file(99, trim(fname_eps), status='old')
100  call close_file(99)
101  call h5fopen_f(trim(fname_eps), H5F_ACC_RDONLY_F, file_id, error)
102  write(6,*) 'Performing consistency checks.'
103  version = -1
104  call hdf5_read_int(file_id, 'versionnumber', version, error)
105  if (version/=1.or.error/=0) then
106    write(0,*)
107    write(0,*) 'ERROR: source file ',trim(fname_eps),' has an incorrect version tag.'
108    write(0,*) '  expecting: ', 1
109    write(0,*) '  read: ', version
110    write(0,*) '  error flag: ', error
111    write(0,*)
112    call die('Invalid file format in '//trim(fname_eps)//'.')
113  endif
114
115  !FHJ: mf_header is a recent addition, so it might not be present in older HDF5 file
116  call h5lexists_f(file_id, 'mf_header', mf_header_exists, error)
117  if (error/=0) call die('Could not search for mf_header.')
118  if (mf_header_exists) then
119    write(6,*) 'File '//TRUNC(fname_eps)//' contains mean-field information.'
120  else
121    write(6,*) 'File '//TRUNC(fname_eps)//' does not contain mean-field information.'
122    if (skip_wfn) &
123      call die('Could not find mf_header in '//trim(fname_eps)//'. You`ll need to specify the WFN argument.')
124  endif
125  !FHJ: nspin is a recent addition, so it might not be present in older HDF5 file
126  call h5lexists_f(file_id, 'nspin', nspin_exists, error)
127  if (error/=0) call die('Could not search for nspin.')
128  if (nspin_exists) then
129    call hdf5_read_int(file_id, 'nspin', nspin, error)
130  else
131    if (skip_wfn) &
132      call die('Could not find nspin in '//trim(fname_eps)//'. You`ll need to specify the WFN argument.')
133  endif
134
135  !FHJ: from the old epsmat.h5.spec file:
136  !Dataset: intinfo
137  !Rank: 1
138  !Dims(1): 9
139  !Value(1): 0 for REAL and 1 for CPLX
140  !Value(2): # of qpoints
141  !Value(3): freq_dep
142  !Value(4): # of frequencies
143  !Value(5): # of gvectors
144  !Value(6): Max of nmtx(q) over q
145  !Value(7): kgrid(1)
146  !Value(8): kgrid(2)
147  !Value(9): kgrid(3)
148  call hdf5_read_double(file_id, 'dblinfo', pol%ecuts, error)
149  call hdf5_read_int_array(file_id, 'intinfo', (/9/), intinfo, error)
150  pol%nq = intinfo(2)
151  pol%freq_dep = intinfo(3)
152  pol%nfreq = intinfo(4)
153  ng = intinfo(5)
154  if (pol%nfreq<1) pol%nfreq = 1
155  if (intinfo(1)+1/=SCALARSIZE) call die('Wrong flavor.')
156  if (.not.skip_wfn) then
157    ! FHJ: Only perform this consistency check if we are using an external WFN file
158    if (nspin_exists) then
159      if (nspin/=mf%kp%nspin) &
160        call die('nspin mismatch between '//trim(fname_wfn)//' and '//trim(fname_eps)//'.')
161    else
162      nspin = mf%kp%nspin
163    endif
164    write(6,*) 'Checking consistency between G-spaces from '//&
165      TRUNC(fname_eps)//' and '//TRUNC(fname_wfn)//'.'
166    if (ng/=mf%gvec%ng) &
167      call die('ng mismatch between '//trim(fname_wfn)//' and '//trim(fname_eps)//'.')
168    SAFE_ALLOCATE(gvecs, (3,ng))
169    call hdf5_read_int_array(file_id, 'gvecs', (/3,ng/), gvecs, error)
170    if (any(mf%gvec%components/=gvecs)) &
171      call die('gvectors mismatch between '//trim(fname_wfn)//' and '//trim(fname_eps)//'.')
172    SAFE_DEALLOCATE(gvecs)
173  endif
174  call h5fclose_f(file_id, error)
175
176  if (.not.mf_header_exists) then
177    ! FHJ: Only write mf_header if it wasn`t there in the first place.
178    write(6,*) 'Writing mean-field header.'
179    call setup_hdf5_mf_file(trim(fname_eps), create_file=.false.)
180    !mf%sheader = 'WFN'
181    !mf%iflavor = SCALARSIZE
182    call write_hdf5_mf_header(trim(fname_eps), mf)
183    call write_hdf5_gvectors(trim(fname_eps), ng, mf%gvec%components)
184  endif
185
186  call h5fopen_f(trim(fname_eps), H5F_ACC_RDWR_F, file_id, error)
187  if (mf_header_exists) then
188    call h5lexists_f(file_id, 'mf_header/flavor', exists, error)
189    ! FHJ: Some old mf_headers didn`t have the flavor/versionnumber.
190    if (.not.exists.or.error/=0) call hdf5_write_int(file_id, 'mf_header/flavor', SCALARSIZE, error)
191    ! FHJ: version hard coded to 1, we don`t want to auto-update it for now.
192    call h5lexists_f(file_id, 'mf_header/versionnumber', exists, error)
193    if (.not.exists.or.error/=0) call hdf5_write_int(file_id, 'mf_header/versionnumber', VER_WFN_HDF5, error)
194  else
195    ! FHJ: Re-read mf-header
196    call read_hdf5_mf_header(trim(fname_eps), mf, sheader='WFN', iflavor=SCALARSIZE, &
197      warn=.false., dont_warn_kgrid=.true.)
198  endif
199
200
201
202  write(6,*) 'Setting up new file layout.'
203  call hdf5_create_group(file_id, 'eps_header', error)
204  call hdf5_create_group(file_id, 'eps_header/params', error)
205  call hdf5_create_group(file_id, 'eps_header/qpoints', error)
206  call hdf5_create_group(file_id, 'eps_header/freqs', error)
207  call hdf5_create_group(file_id, 'eps_header/gspace', error)
208  call hdf5_create_group(file_id, 'mats', error)
209  call hdf5_write_int(file_id, 'eps_header/versionnumber', VER_EPS_HDF5, error)
210  call hdf5_write_int(file_id, 'eps_header/flavor', SCALARSIZE, error)
211
212  ! FHJ: some defaults..
213  if (index(fname_eps,'epsmat')/=0 .or. index(fname_eps,'eps0mat')/=0) then
214    pol%matrix_type = 0
215  elseif (index(fname_eps,'chimat')/=0 .or. index(fname_eps,'chi0mat')/=0) then
216    pol%matrix_type = 2
217  else
218    write(0,*)
219    write(0,*) 'WARNING: could not determine if we have a chimat or epsmat file.'
220    write(0,*) 'Assuming we are dealing with an epsmat file.'
221    write(0,*)
222    pol%matrix_type = 0
223  endif
224  call eps_setup_sizes(pol, SCALARSIZE, nspin)
225  pol%efermi = 0d0
226  pol%intraband_flag = 0
227  pol%intraband_overlap_min = 0.9d0
228  pol%subsample = .false.
229
230  ! FHJ: General datasets
231  write(6,*) 'Rewriting general datasets.'
232  call hdf5_write_int(file_id, 'eps_header/params/matrix_type', pol%matrix_type, error)
233  call hdf5_write_logical(file_id, 'eps_header/params/has_advanced', pol%has_advanced, error)
234  call hdf5_write_int(file_id, 'eps_header/params/nmatrix', pol%nmatrix, error)
235  call hdf5_write_int(file_id, 'eps_header/params/matrix_flavor', pol%matrix_flavor, error)
236  call hdf5_write_int(file_id, 'eps_header/params/icutv', pol%icutv, error)
237  call hdf5_write_double(file_id, 'eps_header/params/ecuts', pol%ecuts, error)
238  call hdf5_write_int(file_id, 'eps_header/params/nband', -1, error)
239  call hdf5_write_double(file_id, 'eps_header/params/efermi', pol%efermi/ryd, error)
240  call hdf5_write_int(file_id, 'eps_header/params/intraband_flag', pol%intraband_flag, error)
241  call hdf5_write_double(file_id, 'eps_header/params/intraband_overlap_min', pol%intraband_overlap_min, error)
242  call hdf5_write_logical(file_id, 'eps_header/params/subsample', pol%subsample, error)
243
244  ! FHJ: Q-points-related datasets
245  write(6,*) 'Rewriting q-points-related datasets.'
246  call hdf5_write_int(file_id, 'eps_header/qpoints/nq', pol%nq, error)
247  call move_dset('qpoints', 'eps_header/qpoints/qpts')
248  call hdf5_write_int_array(file_id, 'eps_header/qpoints/qgrid', (/3/), intinfo(7:9), error)
249  call h5lexists_f(file_id, 'qpt_done', exists, error)
250  if (exists) then
251    call move_dset('qpt_done', 'eps_header/qpoints/qpt_done')
252  else
253    SAFE_ALLOCATE(qpt_done, (pol%nq))
254    qpt_done(:) = .true.
255    call hdf5_write_logical_array(file_id, 'eps_header/qpoints/qpt_done', (/pol%nq/), qpt_done, error)
256    SAFE_DEALLOCATE(qpt_done)
257  endif
258
259  ! FHJ: Frequency-related datasets
260  write(6,*) 'Rewriting frequency-related datasets.'
261  call hdf5_write_int(file_id, 'eps_header/freqs/freq_dep', pol%freq_dep, error)
262  call hdf5_write_int(file_id, 'eps_header/freqs/nfreq', pol%nfreq, error)
263  if (pol%freq_dep==0) then
264    call hdf5_write_double_array(file_id, 'eps_header/freqs/freqs', (/2,pol%nfreq/), (/0d0,0d0/), error)
265  else
266    SAFE_ALLOCATE(freqgrid, (pol%nfreq))
267    SAFE_ALLOCATE(freqs_tmp, (2,pol%nfreq))
268    call hdf5_read_double_array(file_id, 'freqs', (/pol%nfreq/), freqgrid, error)
269    call hdf5_read_double_array(file_id, 'freqbrds', (/2,pol%nfreq/), freqs_tmp, error)
270    do iw=1,pol%nfreq
271      freqs_tmp(1,iw) = freqs_tmp(1,iw) + freqgrid(iw)
272    enddo
273    call hdf5_write_double_array(file_id, 'eps_header/freqs/freqs', (/2, pol%nfreq/), freqs_tmp, error)
274    SAFE_DEALLOCATE(freqgrid)
275    SAFE_DEALLOCATE(freqs_tmp)
276  endif
277
278  ! FHJ: G-vectors-related datasets
279  write(6,*) 'Rewriting G-vectors-related datasets.'
280  call move_dset('nmtx-of-q', 'eps_header/gspace/nmtx')
281  call hdf5_write_int(file_id, 'eps_header/gspace/nmtx_max',  intinfo(6), error)
282  call move_dset('q-gvec-ekin', 'eps_header/gspace/ekin')
283  SAFE_ALLOCATE(isorts, (ng,2,pol%nq))
284  call hdf5_read_int_array(file_id, 'q-gvec-index', shape(isorts), isorts, error)
285  call hdf5_write_int_array(file_id, 'eps_header/gspace/gind_eps2rho', &
286    (/ng,pol%nq/), isorts(:,1,:), error)
287  call hdf5_write_int_array(file_id, 'eps_header/gspace/gind_rho2eps', &
288    (/ng,pol%nq/), isorts(:,2,:), error)
289  SAFE_DEALLOCATE(isorts)
290  ! FHJ: Some old mf_headers don`t contain the gvec%components
291  call h5lexists_f(file_id, 'mf_header/gspace/components', exists, error)
292  if (.not.exists.or.error/=0) call move_dset('gvecs', 'mf_header/gspace/components')
293
294  ! FHJ: Move large matrix datasets
295  write(6,*) 'Moving matrices.'
296  call move_dset('matrix', 'mats/matrix')
297  call move_dset('matrix-diagonal', 'mats/matrix-diagonal')
298
299  write(6,*) 'Cleaning up old structures.'
300  call h5ldelete_f(file_id, 'q-gvec-index', error)
301  call h5ldelete_f(file_id, 'dblinfo', error)
302  call h5ldelete_f(file_id, 'intinfo', error)
303  call delete_if_exists('freqs') ! FHJ: Not present in freq_dep==0
304  call delete_if_exists('freqbrds') ! Same
305  call delete_if_exists('gvecs')
306  call delete_if_exists('info')
307  call delete_if_exists('nspin')
308  call delete_if_exists('versionnumber')
309
310  call h5fclose_f(file_id, error)
311  write(6,*) 'All done!'
312  write(6,*)
313  call h5close_f(error)
314
315contains
316
317  subroutine move_dset(src, dest)
318    character(len=*), intent(in) :: src, dest
319
320    PUSH_SUB(move_dset)
321
322    call h5lcreate_hard_f(file_id, src, file_id, dest, error)
323    call h5ldelete_f(file_id, src, error)
324
325    POP_SUB(move_dset)
326
327  end subroutine move_dset
328
329  subroutine delete_if_exists(dset)
330    character(len=*), intent(in) :: dset
331
332    PUSH_SUB(delete_if_exists)
333
334    exists = .false.
335    call h5lexists_f(file_id, dset, exists, error)
336    if (exists) call h5ldelete_f(file_id, dset, error)
337
338    POP_SUB(delete_if_exists)
339
340  end subroutine delete_if_exists
341
342
343end program epsmat_hdf5_upgrade
344