TABLE OF CONTENTS


ABINIT/outxfhist [ Functions ]

[ Top ] [ Functions ]

NAME

 outxfhist

FUNCTION

  read/write xfhist

COPYRIGHT

 Copyright (C) 2003-2018 ABINIT group (MB)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  option =
   1: write
   2: read only nxfh
   3: read xfhist
  response =
   0: GS wavefunctions
   1: RF wavefunctions
  natom = number of atoms in unit cell
  mxfh = last dimension of the xfhist array

OUTPUT

  ios = error code returned by read operations

SIDE EFFECTS

  nxfh = actual number of (x,f) history pairs, see xfhist array
  wff2 = structured info for wavefunctions
  xfhist(3,natom+4,2,ab_xfh%mxfh) = (x,f) history array, also including
   rprim and stress

PARENTS

      gstate

CHILDREN

      xderiveread,xderiverrecend,xderiverrecinit,xderivewrecend
      xderivewrecinit,xderivewrite

SOURCE

 45 #if defined HAVE_CONFIG_H
 46 #include "config.h"
 47 #endif
 48 
 49 #include "abi_common.h"
 50 
 51 subroutine outxfhist(ab_xfh,natom,option,wff2,ios)
 52 
 53  use defs_basis
 54  use m_profiling_abi
 55  use m_abimover
 56  use m_xmpi
 57  use m_wffile
 58  use m_errors
 59 #if defined HAVE_NETCDF
 60  use netcdf
 61 #endif
 62 
 63 !This section has been created automatically by the script Abilint (TD).
 64 !Do not modify the following lines by hand.
 65 #undef ABI_FUNC
 66 #define ABI_FUNC 'outxfhist'
 67 !End of the abilint section
 68 
 69  implicit none
 70 
 71 !Arguments ------------------------------------
 72  integer          ,intent(in)    :: natom,option
 73  integer          ,intent(out)   :: ios
 74  type(wffile_type),intent(inout)    :: wff2
 75  type(ab_xfh_type),intent(inout) :: ab_xfh
 76 
 77 !Local variables-------------------------------
 78  integer :: ierr,ixfh,ncid_hdr,spaceComm,xfdim2
 79  real(dp),allocatable :: xfhist_tmp(:)
 80  character(len=500) :: message
 81 !no_abirules
 82 #if defined HAVE_NETCDF
 83  integer :: ncerr
 84  integer :: nxfh_id, mxfh_id, xfdim2_id, dim2inout_id, dimr3_id,xfhist_id
 85  integer :: nxfh_tmp,mxfh_tmp,xfdim2_tmp,dim2inout_tmp
 86 #endif
 87 
 88 
 89 ! *************************************************************************
 90 
 91 !DEBUG
 92 !write(std_out,*)'outxfhist  : enter, option = ', option
 93 !ENDDEBUG
 94  ncid_hdr = wff2%unwff
 95  xfdim2 = natom+4
 96 
 97  ios = 0
 98 
 99 !### (Option=1) Write out content of all iterations
100 !#####################################################################
101  if ( option == 1 ) then
102 
103 !  Write the (x,f) history
104    if (wff2%iomode == IO_MODE_FORTRAN) then
105      write(unit=wff2%unwff)ab_xfh%nxfh
106      do ixfh=1,ab_xfh%nxfh
107        write(unit=wff2%unwff)ab_xfh%xfhist(:,:,:,ixfh)
108      end do
109 
110    else if (wff2%iomode == IO_MODE_FORTRAN_MASTER) then
111 !    FIXME: should copy the xfhist to other processors, and check that we are on the master to read in this case
112 !    if node is master
113      write(message, "(A,A,A,A)") ch10, " outxfhist: ERROR -", ch10, &
114 &     'iomode == -1 (localrdwf ) has not been coded yet for xfhist rereading.'
115      MSG_ERROR(message)
116 
117      write(unit=wff2%unwff)ab_xfh%nxfh
118      do ixfh=1,ab_xfh%nxfh
119        write(unit=wff2%unwff)ab_xfh%xfhist(:,:,:,ixfh)
120      end do
121 
122 !    insert mpi broadcast here
123 
124    else if(wff2%iomode==IO_MODE_MPI)then
125      ABI_ALLOCATE(xfhist_tmp,(3*(natom+4)*2))
126      spaceComm=xmpi_comm_self
127      call xderiveWRecInit(wff2,ierr)
128      call xderiveWrite(wff2,ab_xfh%nxfh,ierr)
129      call xderiveWRecEnd(wff2,ierr)
130      do ixfh=1,ab_xfh%nxfh
131        xfhist_tmp(:)=reshape(ab_xfh%xfhist(:,:,:,ixfh),(/3*(natom+4)*2/))
132        call xderiveWRecInit(wff2,ierr)
133        call xderiveWrite(wff2,xfhist_tmp,3*(natom+4)*2,spaceComm,ierr)
134        call xderiveWRecEnd(wff2,ierr)
135      end do
136      ABI_DEALLOCATE(xfhist_tmp)
137 
138 #if defined HAVE_NETCDF
139    else if (wff2%iomode == IO_MODE_NETCDF) then
140 !    check if nxfh and xfhist are defined
141      ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="nxfh",dimid=nxfh_id)
142 
143      if (ncerr /= NF90_NOERR) then
144 !      need to define everything
145        ncerr = nf90_redef (ncid=ncid_hdr)
146        NCF_CHECK_MSG(ncerr," outxfhist : going to define mode ")
147 
148        ncerr = nf90_def_dim(ncid=ncid_hdr,name="dim2inout",len=2,dimid=dim2inout_id)
149        NCF_CHECK_MSG(ncerr," outxfhist : define dim2inout")
150        ncerr = nf90_def_dim(ncid=ncid_hdr,name="mxfh",len=ab_xfh%mxfh,dimid=mxfh_id)
151        NCF_CHECK_MSG(ncerr," outxfhist : define mxfh")
152        ncerr = nf90_def_dim(ncid=ncid_hdr,name="nxfh",len=ab_xfh%nxfh,dimid=nxfh_id)
153        NCF_CHECK_MSG(ncerr," outxfhist : define nxfh")
154        ncerr = nf90_def_dim(ncid=ncid_hdr,name="xfdim2",len=xfdim2,dimid=xfdim2_id)
155        NCF_CHECK_MSG(ncerr," outxfhist : define xfdim2")
156 
157        ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="dimr3",dimid=dimr3_id)
158        NCF_CHECK_MSG(ncerr," outxfhist : inquire dimr3")
159 
160 !      ab_xfh%xfhist(3,natom+4,2,ab_xfh%mxfh)
161        ncerr = nf90_def_var(ncid=ncid_hdr,name="xfhist",xtype=NF90_DOUBLE,&
162 &       dimids=(/dimr3_id,xfdim2_id,dim2inout_id,mxfh_id/),varid=xfhist_id)
163        NCF_CHECK_MSG(ncerr," outxfhist : define xfhist")
164 
165 !      End define mode and go to data mode
166        ncerr = nf90_enddef(ncid=ncid_hdr)
167        NCF_CHECK_MSG(ncerr," outxfhist : enddef call ")
168      else
169 !      check that the dimensions are correct
170        ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="nxfh",dimid=nxfh_id)
171        NCF_CHECK_MSG(ncerr," outxfhist : inquire nxfh")
172        ncerr = nf90_Inquire_Dimension(ncid=ncid_hdr,dimid=nxfh_id,&
173 &       len=nxfh_tmp)
174        NCF_CHECK_MSG(ncerr,"  outxfhist : get nxfh")
175        ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="xfdim2",dimid=xfdim2_id)
176        NCF_CHECK_MSG(ncerr," outxfhist : inquire xfdim2")
177        ncerr = nf90_Inquire_Dimension(ncid=ncid_hdr,dimid=xfdim2_id,&
178 &       len=xfdim2_tmp)
179        NCF_CHECK_MSG(ncerr,"  outxfhist : get xfdim2")
180        ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="mxfh",dimid=mxfh_id)
181        NCF_CHECK_MSG(ncerr," outxfhist : inquire mxfh")
182        ncerr = nf90_Inquire_Dimension(ncid=ncid_hdr,dimid=mxfh_id,&
183 &       len=mxfh_tmp)
184        NCF_CHECK_MSG(ncerr,"  outxfhist : get mxfh")
185        ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="dim2inout",dimid=dim2inout_id)
186        NCF_CHECK_MSG(ncerr," outxfhist : inquire dim2inout")
187        ncerr = nf90_Inquire_Dimension(ncid=ncid_hdr,dimid=dim2inout_id,&
188 &       len=dim2inout_tmp)
189        NCF_CHECK_MSG(ncerr,"  outxfhist : get dim2inout")
190 
191        ncerr = nf90_inq_varid(ncid=ncid_hdr,name="xfhist",varid=xfhist_id)
192        NCF_CHECK_MSG(ncerr," outxfhist : inquire xfhist")
193 
194        if (mxfh_tmp /= ab_xfh%mxfh .or. dim2inout_tmp /= 2 .or. xfdim2_tmp /= xfdim2) then
195          write (message,"(A)") 'outxfhist : ERROR xfhist has bad dimensions in NetCDF file. Can not re-write it.'
196          MSG_ERROR(message)
197        end if
198 
199      end if
200 
201 !    Now fill the data
202      ncerr = nf90_put_var(ncid=ncid_hdr,varid=xfhist_id,values=ab_xfh%xfhist)
203      NCF_CHECK_MSG(ncerr," outxfhist : fill xfhist")
204 
205 !    end NETCDF definition ifdef
206 #endif
207    end if  ! end iomode if
208 
209 !  ### (Option=2) Read in number of iterations
210 !  #####################################################################
211  else if ( option == 2 ) then
212 
213    if (wff2%iomode == IO_MODE_FORTRAN) then
214      read(unit=wff2%unwff,iostat=ios)ab_xfh%nxfh
215 
216    else if (wff2%iomode == IO_MODE_FORTRAN_MASTER) then
217 !    FIXME: should copy the xfhist to other processors, and check that we are on the master to read in this case
218 !    if node is master
219      write(message, "(A,A,A,A)") ch10, " outxfhist: ERROR -", ch10, &
220 &     'iomode == -1 (localrdwf ) has not been coded yet for xfhist rereading.'
221      MSG_ERROR(message)
222 
223      read(unit=wff2%unwff,iostat=ios)ab_xfh%nxfh
224 
225    else if (wff2%iomode == IO_MODE_MPI) then
226      call xderiveRRecInit(wff2,ierr)
227      call xderiveRead(wff2,ab_xfh%nxfh,ierr)
228      call xderiveRRecEnd(wff2,ierr)
229 
230 #if defined HAVE_NETCDF
231    else if (wff2%iomode == IO_MODE_NETCDF) then
232      ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="nxfh",dimid=nxfh_id)
233      NCF_CHECK_MSG(ncerr," outxfhist : inquire nxfh")
234      ncerr = nf90_Inquire_Dimension(ncid=ncid_hdr,dimid=nxfh_id,&
235 &     len=ab_xfh%nxfh)
236      NCF_CHECK_MSG(ncerr,"  outxfhist : get nxfh")
237 #endif
238    end if
239 
240 !  ### (Option=3) Read in iteration content
241 !  #####################################################################
242  else if ( option == 3 ) then
243    if (wff2%iomode == IO_MODE_FORTRAN) then
244      do ixfh=1,ab_xfh%nxfhr
245        read(unit=wff2%unwff,iostat=ios)ab_xfh%xfhist(:,:,:,ixfh)
246      end do
247    else if (wff2%iomode == IO_MODE_FORTRAN_MASTER) then
248 !    FIXME: should copy the xfhist to other processors, and check that we are on the master to read in this case
249 !    if node is master
250      write(message, "(A,A,A,A)") ch10, " outxfhist: ERROR -", ch10, &
251 &     'iomode == -1 (localrdwf ) has not been coded yet for xfhist rereading.'
252      MSG_ERROR(message)
253 
254      do ixfh=1,ab_xfh%nxfhr
255        read(unit=wff2%unwff,iostat=ios)ab_xfh%xfhist(:,:,:,ixfh)
256      end do
257 
258    else if (wff2%iomode == IO_MODE_MPI) then
259      ABI_ALLOCATE(xfhist_tmp,(3*(natom+4)*2))
260      spaceComm=xmpi_comm_self
261      do ixfh=1,ab_xfh%nxfhr
262        call xderiveRRecInit(wff2,ierr)
263        call xderiveRead(wff2,xfhist_tmp,3*(natom+4)*2,spaceComm,ierr)
264        call xderiveRRecEnd(wff2,ierr)
265        xfhist_tmp(:)=xfhist_tmp(:)
266      end do
267      ABI_DEALLOCATE(xfhist_tmp)
268    end if
269 
270 !  FIXME: should this be inside the if not mpi as above for options 1 and 2?
271 !  it is placed here because the netcdf read is a single operation
272 #if defined HAVE_NETCDF
273    if (wff2%iomode == IO_MODE_NETCDF) then
274      ncerr = nf90_inq_dimid(ncid=ncid_hdr,name="nxfh",dimid=nxfh_id)
275      NCF_CHECK_MSG(ncerr," outxfhist : inquire nxfh")
276      ncerr = nf90_Inquire_Dimension(ncid=ncid_hdr,dimid=nxfh_id,&
277 &     len=ab_xfh%nxfhr)
278      NCF_CHECK_MSG(ncerr,"  outxfhist : get nxfh")
279 
280      ncerr = nf90_inq_varid(ncid=ncid_hdr,varid=xfhist_id,name="xfhist")
281      NCF_CHECK_MSG(ncerr," outxfhist : inquire xfhist")
282      ncerr = nf90_get_var(ncid=ncid_hdr,varid=xfhist_id,values=ab_xfh%xfhist,&
283 &     start=(/1,1,1,1/),count=(/3,natom+4,2,ab_xfh%nxfhr/))
284      NCF_CHECK_MSG(ncerr," outxfhist : read xfhist")
285    end if
286 #endif
287 
288  else
289 !  write(std_out,*)' outxfhist : option ', option , ' not available '
290    write(message, "(A,A,A,A,I3,A)") ch10, "outxfhist: ERROR -", ch10, &
291 &   "option ", option, " not available."
292    MSG_ERROR(message)
293  end if
294 
295 end subroutine outxfhist