TABLE OF CONTENTS


ABINIT/m_io_redirect [ Modules ]

[ Top ] [ Modules ]

NAME

  m_io_redirect

FUNCTION

  management of output and log files when parallelisation on cells is activated

COPYRIGHT

  Copyright (C) 2001-2024 ABINIT group (FJ,MT)
  This file is distributed under the terms of the
  GNU General Public License, see ~ABINIT/Infos/copyright
  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_io_redirect
23 
24  use defs_basis
25  use m_abicore
26  use m_errors
27 
28  use m_xmpi,        only : xmpi_comm_rank, xmpi_barrier
29  use m_libpaw_tools,only : libpaw_write_comm_set
30  use m_io_tools,    only : delete_file, flush_unit, open_file, get_unit
31  use m_fstrings,    only : int2char4
32 
33  implicit none
34 
35  private
36 
37  character(len=fnlen),allocatable,save :: fillog(:),filout(:)
38 
39 !Procedures ------------------------------------
40  public :: localfilnam    ! Define new appended name
41  public :: localwrfile    ! write local files
42  public :: localrdfile    ! read local files
43  public :: localredirect  ! redirect communicators for local files
44 
45 contains

m_io_redirect/localfilnam [ Functions ]

[ Top ] [ m_io_redirect ] [ Functions ]

NAME

 localfilnam

FUNCTION

 Define new appended name

INPUTS

   commspace= communicator between cells
   commworld= communicator of the world
   filnam(5)=character strings giving file names
   nam= name of the extension of the files
   nfil= number of files to be named

SOURCE

 64  subroutine localfilnam(commspace,commspace1,commworld,filnam,nam,nfil)
 65 
 66 !Arguments ------------------------------------
 67 !scalars
 68  integer,intent(in) :: commspace,commspace1,commworld,nfil
 69  character(len=4) :: nam
 70 !arrays
 71  character(len=fnlen),intent(in) :: filnam(:)
 72 
 73 !Local variables-------------------------------
 74 !scalars
 75  integer :: ii,ios,me,me_loc
 76  character(len=10) :: appen,tag
 77 
 78 ! *********************************************************************
 79  
 80  if (nfil<=1) return
 81 
 82  ABI_MALLOC(filout,(nfil))
 83  ABI_MALLOC(fillog,(nfil))
 84  me=xmpi_comm_rank(commworld)
 85  me_loc=xmpi_comm_rank(commspace1)
 86  call int2char4(me,tag)
 87  ABI_CHECK((tag(1:1)/='#'),'Bug: string length too short!')
 88  do ii=1,nfil
 89    call appdig(ii,'',appen)
 90    filout(ii)=trim(filnam(2))//nam//trim(appen)
 91    if (me_loc==0) then
 92      fillog(ii)=trim(filnam(5))//nam//trim(appen)//"_LOG"
 93    else
 94      fillog(ii)=trim(filnam(5))//nam//trim(appen)//"_LOG_"//trim(tag)
 95    end if
 96  end do
 97 
 98  if (me==0) then
 99    do ii=1,nfil
100      call delete_file(filout(ii),ios)
101    end do
102  end if
103 
104  call xmpi_barrier(commspace)
105 
106  end subroutine localfilnam

m_io_redirect/localrdfile [ Functions ]

[ Top ] [ m_io_redirect ] [ Functions ]

NAME

   localrdfile

FUNCTION

   read local (for each cell) output and log files
   to gather in one output or log file

INPUTS

   commspace= communicator between cells
   commworld= communicator of the world
   compute_all= flag to activate the reading on all cells
   [dyn(nfil)]= --optional-- indicate if this cell is to be taken into account
   nfil= number of files to be named
   paral= flag to activate parallelisation
   prtvol= flag to activate printing

SOURCE

191  subroutine localrdfile(commspace,commworld,compute_all,nfil,paral,prtvol,dyn)
192 
193 !Arguments ------------------------------------
194 !scalars
195  integer,intent(in) :: commspace,commworld,nfil,paral,prtvol
196  logical,intent(in) :: compute_all
197 !arrays
198  integer,intent(in),target,optional :: dyn(nfil)
199 
200 !Local variables ------------------------------
201 !scalars
202  integer :: ii,ios,me,temp_unit
203  logical :: eof
204  character(len=500) :: msg
205  character(len=fnlen) :: line
206 !arrays
207  integer,pointer :: my_dyn(:)
208 
209 ! *********************************************************************
210 
211  if (nfil<=1) return
212  
213  me=xmpi_comm_rank(commworld)
214  call xmpi_barrier(commspace) !waiting for all files to be written and close
215 
216  if (prtvol<=0.and.(paral==1)) then
217    if (me==0) then
218      if (present(dyn)) then
219        my_dyn => dyn
220      else
221        ABI_MALLOC(my_dyn,(nfil))
222        my_dyn(:)=1
223      end if
224      do ii=1,nfil
225        if (open_file(filout(ii),msg,newunit=temp_unit,status='old') == 0) then
226          if (my_dyn(ii)==1.or.compute_all) then
227            eof=.false.
228            do while (.not.eof)
229              read(temp_unit,fmt='(a)',err=111,end=111,iostat=ios) line
230              if (ios==0) write(ab_out,'(a)') trim(line)
231              goto 112
232              111              eof=.true.
233              112              continue
234            end do
235          end if
236          close(unit=temp_unit,status='delete')
237        end if
238      end do
239      call flush_unit(ab_out)
240      if (.not.present(dyn)) then
241        ABI_FREE(my_dyn)
242      end if
243    end if
244    call xmpi_barrier(commspace)
245  end if
246 
247  if (paral==1.and.do_write_log) then
248    if (me==0) then
249      do ii=1,nfil
250        if (open_file(fillog(ii),msg,newunit=temp_unit,status='old') == 0) then
251          eof=.false.
252          do while (.not.eof)
253            read(temp_unit,fmt='(a)',err=113,end=113,iostat=ios) line
254            if (ios==0) write(std_out,'(a)') trim(line)
255            goto 114
256            113            eof=.true.
257            114            continue
258          end do
259          close(unit=temp_unit,status='delete')
260        end if
261      end do
262    end if
263  end if
264 
265  call xmpi_barrier(commspace)
266 
267  ABI_FREE(filout)
268  ABI_FREE(fillog)
269 
270  end subroutine localrdfile

m_io_redirect/localredirect [ Functions ]

[ Top ] [ m_io_redirect ] [ Functions ]

NAME

   localredirect

FUNCTION

   Close output units ; restore defaults

INPUTS

   commspace= communicator between cells
   commworld= communicator of the world
   nfil= number of files to be named
   paral= flag to activate parallelisation
   prtvol= flag to activate printing  

SOURCE

289  subroutine localredirect(commspace,commworld,nfil,paral,prtvol)
290 
291 !Arguments ------------------------------------
292 !scalars
293  integer, intent(in) :: commspace,commworld,nfil,paral,prtvol
294 
295 !Local variables ------------------------------
296  integer :: me
297 
298 ! *********************************************************************
299 
300  if (nfil<=1) return
301 
302  me=xmpi_comm_rank(commspace)
303  if (prtvol>0) then
304    if (do_write_log.and.me==0) close(unit=ab_out)
305  else if (paral==1) then
306    if (me==0) close(unit=ab_out)
307  end if
308  if (paral==1.and.me==0.and.do_write_log) close(unit=std_out)
309  call abi_io_redirect(new_ab_out=ab_out_default,new_std_out=std_out_default,new_io_comm=commworld)
310  call libpaw_write_comm_set(commworld)
311 
312  end subroutine localredirect

m_io_redirect/localwrfile [ Functions ]

[ Top ] [ m_io_redirect ] [ Functions ]

NAME

   localwrfile

FUNCTION

   Write local (for each cell) output and log files

INPUTS

   commspace= communicator inside cells
   ii=number of the iteration on the cells
   nfil= number of files to be named
   paral= flag to activate parallelisation
   prtvol= flag to activate printing

SOURCE

125  subroutine localwrfile(commspace,ii,nfil,paral,prtvol)
126 
127 !Arguments ------------------------------------
128  integer, intent(in) :: commspace,ii,nfil,paral,prtvol
129 
130 !Local variables ------------------------------
131  integer :: me
132  character(len=500) :: msg
133 
134 ! *********************************************************************
135 
136  if (nfil<=1) return
137 
138  me=xmpi_comm_rank(commspace)
139  if (paral==1) then
140    call abi_io_redirect(new_io_comm=commspace)
141    call libpaw_write_comm_set(commspace)
142  end if
143  if (prtvol>0) then
144    if (do_write_log) then
145      call abi_io_redirect(new_ab_out=get_unit())
146      if (me==0) then 
147        if (open_file(NULL_FILE,msg,unit=ab_out,status='unknown') /= 0) then
148          ABI_ERROR(msg)
149        end if
150      end if
151    else
152      call abi_io_redirect(new_ab_out=std_out)
153    end if
154  else if (paral==1) then
155    call abi_io_redirect(new_ab_out=get_unit())
156    if (me==0) then 
157      if (open_file(filout(ii),msg,unit=ab_out,status='unknown') /= 0) then
158        ABI_ERROR(msg)
159      end if
160    end if
161  end if
162  if (paral==1.and.me==0.and.do_write_log) then
163    call abi_io_redirect(new_std_out=get_unit())
164    if (open_file(fillog(ii),msg,unit=std_out,status='unknown') /= 0) then
165      ABI_ERROR(msg)
166    end if
167  end if
168 
169  end subroutine localwrfile