TABLE OF CONTENTS


ABINIT/m_wvl_rwwf [ Modules ]

[ Top ] [ Modules ]

NAME

  m_wvl_rwwf

FUNCTION

COPYRIGHT

  Copyright (C) 1998-2024 ABINIT group (DC)
  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_wvl_rwwf
23 
24  use defs_basis
25 
26  use defs_wvltypes
27  use m_wffile
28  use m_errors
29  use m_abicore
30  use m_hdr
31  use m_xmpi
32  use m_dtset
33 
34  use defs_abitypes,  only : MPI_type
35  use m_geometry,     only : xred2xcart
36 
37  implicit none
38 
39  private

ABINIT/wvl_read [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_read

FUNCTION

 Simple wrapper around the read disk methods of BigDFT for wavefunctions.

INPUTS

  dtset <type(dataset_type)>=input variables.
  mpi_enreg=information about MPI parallelization
  option= -2 for write with BigDFT format,
          -1 for reading wavelets coefficients with BigDFT format,
          2 for write,
          1 for read.
  rprimd(3,3)=dimensional primitive translations in real space (bohr)
  wff <type(wffile_type)>=struct info for wavefunction
  wfs <type(wvl_wf_type)>=wavefunctions information for wavelets.
  xred(3,natom)=reduced dimensionless atomic coordinates

OUTPUT

SOURCE

 72 subroutine wvl_read(dtset, hdr0, hdr, mpi_enreg, option, rprimd, wff, wfs, wvl, xred)
 73 
 74 #if defined HAVE_BIGDFT
 75   use BigDFT_API, only: readonewave, reformatonewave, readmywaves, &
 76 &                       WF_FORMAT_NONE
 77 #endif
 78   implicit none
 79 
 80 !Arguments -------------------------------
 81   !scalars
 82   integer, intent(in)                       :: option
 83   type(dataset_type), intent(in)            :: dtset
 84   type(hdr_type), intent(in)                :: hdr0
 85   type(hdr_type), intent(in)                :: hdr
 86   type(MPI_type), intent(in)                :: mpi_enreg
 87   type(wffile_type),intent(in)              :: wff
 88   type(wvl_wf_type), intent(inout)          :: wfs
 89   type(wvl_internal_type), intent(in)       :: wvl
 90   !arrays
 91   real(dp), intent(in)                      :: rprimd(3, 3)
 92   real(dp), intent(in)                      :: xred(3, dtset%natom)
 93 
 94 !Local variables-------------------------------
 95 #if defined HAVE_BIGDFT
 96   character(len = 500)  :: message
 97   integer               :: iBand, bandSize
 98   integer               :: comm,me
 99   real(dp), allocatable :: xcart(:,:), psifscf(:,:,:)
100   real(dp), allocatable :: xcart_old(:,:)
101 #endif
102 
103 ! *********************************************************************
104 
105 #if defined HAVE_BIGDFT
106 
107  if (abs(option) /= 1) then
108    write(message,'(a,a,a,i0,a)')&
109 &   '  Option argument is wrong,', ch10, &
110 &   '  awaited values are -1 or  1 but option = ', option, '.'
111    ABI_BUG(message)
112  end if
113 
114  comm=mpi_enreg%comm_wvl
115  me=xmpi_comm_rank(comm)
116 !Store xcart for each atom
117  ABI_MALLOC(xcart,(3, dtset%natom))
118  ABI_MALLOC(xcart_old,(3, dtset%natom))
119  call xred2xcart(dtset%natom, rprimd, xcart, xred)
120 
121  write(message,'(2a)') ch10,' wvl_read:  read wavefunctions from file.'
122  call wrtout(std_out,message,'COLL')
123 
124  if (option > 0) then
125    bandSize = wfs%ks%lzd%Glr%wfd%nvctr_c + 7 * wfs%ks%lzd%Glr%wfd%nvctr_f
126 !  Read in the ABINIT way.
127    if (wff%iomode == IO_MODE_FORTRAN .or. (wff%iomode == IO_MODE_FORTRAN_MASTER .and. wff%master==wff%me)) then
128      ABI_MALLOC(psifscf,(wvl%Glr%d%n1i, wvl%Glr%d%n2i, wvl%Glr%d%n3i))
129      do iBand = 1, dtset%mband * dtset%nsppol, 1
130        call readonewave(wff%unwff, .false., iBand, me, &
131 &       wvl%Glr%d%n1, wvl%Glr%d%n2, wvl%Glr%d%n3, &
132 &       wvl%h(1), wvl%h(2), wvl%h(3), wvl%atoms, &
133 &       wfs%ks%lzd%Glr%wfd, xcart_old, xcart, &
134 &       wfs%ks%psi(bandSize * (iBand - me * wfs%ks%orbs%norbp - 1) + 1: &
135 &       bandSize * (iBand - me * wfs%ks%orbs%norbp - 1) + bandSize), &
136 &       wfs%ks%orbs%eval(iBand), psifscf)
137      end do
138      ABI_FREE(psifscf)
139 
140    else
141      write(message,'(4a,i0,a)') ch10,&
142 &     '  wff%iomode argument is wrong,', ch10, &
143 &     '  awaited values are -1, 0 (or 3 if netcdf/etsf_io is available) but value = ', wff%iomode, '.'
144      ABI_BUG(message)
145    end if
146  else
147    call readmywaves(me, "wavefunctions", WF_FORMAT_NONE, wfs%ks%orbs, &
148 &   wvl%Glr%d%n1, wvl%Glr%d%n2, wvl%Glr%d%n3, &
149 &   wvl%h(1), wvl%h(2), wvl%h(3), wvl%atoms, &
150 &   xcart_old, xcart, wfs%ks%lzd%Glr%wfd, wfs%ks%psi)
151  end if
152 
153  ABI_FREE(xcart)
154  ABI_FREE(xcart_old)
155 #else
156  BIGDFT_NOTENABLED_ERROR()
157  if (.false.) write(std_out,*) option,dtset%nstep,hdr0%ecut,hdr%ecut,mpi_enreg%nproc,wff%me,&
158 & wfs%ks,wvl%h(1),rprimd(1,1),xred(1,1)
159 #endif
160 
161 end subroutine wvl_read

ABINIT/wvl_write [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_write

FUNCTION

 Simple wrapper around the write disk methods of BigDFT for wavefunctions.

INPUTS

  dtset <type(dataset_type)>=input variables.
  mpi_enreg=information about MPI parallelization
  option= -2 for write with BigDFT format,
          -1 for reading wavelets coefficients with BigDFT format,
          2 for write,
          1 for read.
  rprimd(3,3)=dimensional primitive translations in real space (bohr)
  wff <type(wffile_type)>=struct info for wavefunction
  wfs <type(wvl_wf_type)>=wavefunctions information for wavelets.

OUTPUT

SIDE EFFECTS

  xred(3,natom)=reduced dimensionless atomic coordinates

SOURCE

189 subroutine wvl_write(dtset, eigen, mpi_enreg, option, rprimd, wff, wfs, wvl, xred)
190 
191 #if defined HAVE_BIGDFT
192   use BigDFT_API, only : writeonewave,writemywaves,WF_FORMAT_NONE
193 #endif
194   implicit none
195 
196 !Arguments -------------------------------
197   !scalars
198   integer, intent(in)                       :: option
199   type(dataset_type), intent(in)            :: dtset
200   type(MPI_type), intent(in)                :: mpi_enreg
201   type(wffile_type),intent(in)              :: wff
202   type(wvl_wf_type), intent(in)             :: wfs
203   type(wvl_internal_type), intent(in)       :: wvl
204   !arrays
205   real(dp), intent(in), target              :: eigen(dtset%mband)
206   real(dp), intent(in)                      :: rprimd(3, 3)
207   real(dp), intent(in)                      :: xred(3, dtset%natom)
208 
209 !Local variables-------------------------------
210 #if defined HAVE_BIGDFT
211   character(len = 500)  :: message
212   integer               :: comm,me
213   integer               :: iorb
214   integer               :: iseg, nseg, ipsi, npsi
215   real(dp), allocatable :: xcart(:,:)
216 #endif
217 
218 ! *********************************************************************
219 
220 #if defined HAVE_BIGDFT
221 
222  if (abs(option) /= 2) then
223    write(message,'(a,a,a,i0,a)')&
224 &   '  Option argument is wrong,', ch10, &
225 &   '  awaited values are -2 or  2 but option = ', option, '.'
226    ABI_BUG(message)
227  end if
228 
229  comm=mpi_enreg%comm_wvl
230  me=xmpi_comm_rank(comm)
231 !Store xcart for each atom
232  ABI_MALLOC(xcart,(3, dtset%natom))
233  call xred2xcart(dtset%natom, rprimd, xcart, xred)
234 
235  write(message, '(a,a,a,a)' ) ch10,&
236 & ' wvl_write:  Write wavefunctions to file.'
237  call wrtout(std_out,message,'COLL')
238 
239  if (option > 0) then
240 !  Write in the ABINIT way.
241    if (wff%iomode == IO_MODE_FORTRAN .or. (wff%iomode == IO_MODE_FORTRAN_MASTER .and. wff%master==wff%me)) then
242      iseg = wfs%ks%lzd%Glr%wfd%nseg_c
243      nseg = wfs%ks%lzd%Glr%wfd%nseg_c + wfs%ks%lzd%Glr%wfd%nseg_f
244      ipsi = wfs%ks%lzd%Glr%wfd%nvctr_c
245      npsi = wfs%ks%lzd%Glr%wfd%nvctr_c + 7 * wfs%ks%lzd%Glr%wfd%nvctr_f
246      do iorb = 1, dtset%mband
247        call writeonewave(wff%unwff, .false., iorb, wvl%Glr%d%n1, &
248 &       wvl%Glr%d%n2, wvl%Glr%d%n3, &
249 &       wvl%h(1), wvl%h(2), wvl%h(3), dtset%natom, &
250 &       xcart, wfs%ks%lzd%Glr%wfd%nseg_c, wfs%ks%lzd%Glr%wfd%nvctr_c, &
251 &       wfs%ks%lzd%Glr%wfd%keygloc(:,1:iseg), wfs%ks%lzd%Glr%wfd%keyvloc(1:iseg), wfs%ks%lzd%Glr%wfd%nseg_f, &
252 &       wfs%ks%lzd%Glr%wfd%nvctr_f, wfs%ks%lzd%Glr%wfd%keygloc(:, iseg + 1:nseg), &
253 &       wfs%ks%lzd%Glr%wfd%keyvloc(iseg + 1:nseg), &
254 &       wfs%ks%psi(npsi * (iorb - me * wfs%ks%orbs%norbp - 1) + 1: &
255 &       npsi * (iorb - me * wfs%ks%orbs%norbp - 1) + ipsi), &
256 &       wfs%ks%psi(npsi * (iorb - me * wfs%ks%orbs%norbp - 1) + ipsi + 1: &
257 &       npsi * (iorb - me * wfs%ks%orbs%norbp - 1) + npsi), &
258 &       wfs%ks%orbs%eval(iorb))
259      end do
260    else
261      write(message,'(3a,i0,a)')&
262 &     '  wff%iomode argument is wrong,', ch10, &
263 &     '  awaited values are -1, 0 (or 3 if netcdf/etsf_io is available) but value = ', wff%iomode, '.'
264      ABI_BUG(message)
265    end if
266  else
267 !  Write in the BigDFT way.
268    call  writemywaves(me, "wavefunctions", WF_FORMAT_NONE, wfs%ks%orbs, &
269 &   wvl%Glr%d%n1, wvl%Glr%d%n2, wvl%Glr%d%n3, &
270 &   wvl%h(1), wvl%h(2), wvl%h(3),wvl%atoms, &
271 &   xcart, wfs%ks%lzd%Glr%wfd, wfs%ks%psi)
272  end if
273 
274  ABI_FREE(xcart)
275 
276 #else
277  BIGDFT_NOTENABLED_ERROR()
278  if (.false.) write(std_out,*) option,dtset%nstep,mpi_enreg%nproc,wff%me,&
279 & wfs%ks,wvl%h(1),eigen(1),rprimd(1,1),xred(1,1)
280 #endif
281 
282 end subroutine wvl_write