Skip to content

Commit

Permalink
Merge work in fork _joseph8; first version of scribe diectated I/O. A…
Browse files Browse the repository at this point in the history
…t the moment this is implemented as a CPP option and the old I/O approach is invoked with OLDIO in make files
  • Loading branch information
josephzhang8 committed Aug 30, 2021
1 parent c2e760b commit 8efc374
Show file tree
Hide file tree
Showing 14 changed files with 320 additions and 53 deletions.
2 changes: 1 addition & 1 deletion cmake/README
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@ cd ../build; rm -rf *
cmake -C ../cmake/SCHISM.local.build -C ../cmake/SCHISM.local.<cluster> ../src/
(On Pleiades: cmake -DMPIVERSION=1 -C ../cmake/SCHISM.local.build -C ../cmake/SCHISM.local.pleiades ../src/)

make -j8 or make VERBOSE=1 >& tmp
make -j8 pschism or make VERBOSE=1 pschism >& tmp

(Executables are in bin/; libs in lib/)

Expand Down
3 changes: 3 additions & 0 deletions cmake/SCHISM.local.build
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,9 @@
# Algorithm choices
# TVD_LIM must be one of SB, VL, MM or OS for Superbee, Van Leer, Minmod, or Osher.")
set (TVD_LIM VL CACHE STRING "Flux limiter")
#Turn OLDIO off to use the new scribe based I/O
set (OLDIO OFF CACHE BOOLEAN "Old nc output (each rank dumps its own data)")

set (PREC_EVAP OFF CACHE BOOLEAN "Include precipitation and evaporation calculation")
##Older versions of GOTM (3.*) have issues with netcdf v4, so are not maintained
set (USE_GOTM OFF CACHE BOOLEAN "Use GOTM turbulence model. This just enables the build -- GOTM must still be selected in param.nml")
Expand Down
2 changes: 2 additions & 0 deletions cmake/SCHISM.local.build.example
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,8 @@
# Algorithm choices
# TVD_LIM must be one of SB, VL, MM or OS for Superbee, Van Leer, Minmod, or Osher.")
set (TVD_LIM VL CACHE STRING "Flux limiter")
#Turn OLDIO off to use the new scribe based I/O
set (OLDIO OFF CACHE BOOLEAN "Old nc output (each rank dumps its own data)")
set (PREC_EVAP OFF CACHE BOOLEAN "Include precipitation and evaporation calculation")
set (USE_GOTM OFF CACHE BOOLEAN "Use GOTM turbulence model. This just enables the build -- GOTM must still be selected in param.nml")
set (USE_HA OFF CACHE BOOLEAN "Enable harmonic analysis output modules") # Not implemented yet
Expand Down
4 changes: 4 additions & 0 deletions mk/include_modules
Original file line number Diff line number Diff line change
@@ -1,3 +1,7 @@
# Old nc output option (each rank dumps its own data)
# USE_OLDIO = yes
# EXEC := $(EXEC)_OLDIO

# Use PetSc as solver for FE eq
#'v1': v3.5; 'v2': 3.6, 3.7; 'v3': 3.10
# PETSC_VERSION = v2
Expand Down
1 change: 1 addition & 0 deletions src/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -58,6 +58,7 @@ message(STATUS "\n### Processing algorithmic switches and build options involvin
define_opt(DEBUG "Enable diagnostics (slow and not suitable for production simulations" OFF)
define_opt(USE_ANALYSIS "Triggers detailed calculations and output options" OFF)
define_opt(INCLUDE_TIMING "Add timing profile instrumentation to simulation" OFF)
define_opt(OLDIO "Old nc output option" OFF)
define_opt(PREC_EVAP "Use precipitation/evaporation option" OFF)
define_opt(USE_HA "Enable Harmonic Analysis output" OFF)
define_opt(USE_SIMPLE_WIND "Simplified wind data interface" OFF)
Expand Down
1 change: 1 addition & 0 deletions src/Core/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -30,6 +30,7 @@ add_library(core ${local_extra_code}
schism_msgp.F90
schism_assert.F90
schism_io.F90
scribe_io.F90
misc_modules.F90
hydraulic_structures.F90
gen_modules_clock.F90
Expand Down
14 changes: 9 additions & 5 deletions src/Core/schism_glbl.F90
Original file line number Diff line number Diff line change
Expand Up @@ -55,7 +55,7 @@ module schism_glbl
integer,save :: nthreads

!In/out dirs
character(len=1000) :: in_dir,out_dir
character(len=1000),save :: in_dir,out_dir
integer,save :: len_in_dir,len_out_dir

! For timing
Expand Down Expand Up @@ -145,10 +145,14 @@ module schism_glbl
character(len=12),save :: ifile_char
! character(len=48),save,dimension(mnout) :: outfile !,variable_nm,variable_dim
integer,save :: ihfskip,nrec,nspool,ifile,ifile_len, &
&noutput,it_main,iths_main,id_out_var(2000)
! integer,save,dimension(mnout) :: iof
! integer,save,allocatable :: ichan_ns(:),iof_ns(:)
real(rkind) :: time_stamp !simulation time in sec
&noutput,noutvars,it_main,iths_main,id_out_var(2000),ncount_2dnode, &
&ncount_3dnode,nsend_varout
! integer,save,dimension(mnout) :: iof
integer,save,allocatable :: scribe_cat(:),srqst7(:)
real(rkind),save :: time_stamp !simulation time in sec
!Send var buffers
real(4),save,allocatable :: varout_3dnode(:,:,:),varout_3delem(:,:,:),varout_3dside(:,:,:)
real(4),save,allocatable :: varout_2dnode(:,:,:),varout_2delem(:,:,:),varout_2dside(:,:,:)
character(len=48),save,allocatable :: outfile_ns(:) !,varnm_ns(:)
character(len=48),save :: a_48
character(len=16),save :: a_16
Expand Down
52 changes: 43 additions & 9 deletions src/Core/schism_msgp.F90
Original file line number Diff line number Diff line change
Expand Up @@ -49,12 +49,20 @@ module schism_msgp
! Public data
!-----------------------------------------------------------------------------

integer,public,save :: myrank ! Rank of MPI task (0 base)
integer,public,save :: nproc ! Number of MPI tasks
integer,public,save :: nscribes ! # of I/O scribes (>=0)
integer,public,save :: task_id ! 1: compute; 2: I/O
integer,public,save :: myrank ! Rank of MPI compute task (0 base)
integer,public,save :: nproc ! Number of MPI compute tasks
integer,public,save :: myrank_schism ! Rank of MPI task from union world
integer,public,save :: nproc_schism ! Number of MPI tasks in union world
integer,public,save :: myrank_scribe ! Rank of MPI task from scribe world
integer,public,save :: nproc_scribe ! Number of MPI tasks from scribe world
integer,public,save :: ierr ! Return flag for MPI calls
integer,public,save :: istatus(MPI_STATUS_SIZE) ! Return status for MPI calls

integer,public,save :: comm ! MPI Communicator for ELCIRC
integer,public,save :: comm ! MPI Communicator for compute only
integer,public,save :: comm_scribe ! MPI Communicator for scribes only
integer,public,save :: comm_schism ! Union Communicator for compute & scribes
integer,public,save :: itype = MPI_INTEGER ! MPI Integer Type
!#ifdef USE_SINGLE
! integer,public,save :: rtype = MPI_REAL4 ! MPI Real Type -- to match rkind
Expand Down Expand Up @@ -368,28 +376,54 @@ subroutine parallel_init(communicator)
implicit none
integer, optional :: communicator

integer :: comm2,nproc3,myrank3,nproc_compute


if (present(communicator)) then
! Expect external call to mpi_init,
! use communicator as provided by the interface
call mpi_comm_dup(communicator,comm,ierr)
call mpi_comm_dup(communicator,comm_schism,ierr)
else
! Initialize MPI
call mpi_init(ierr)
if(ierr/=MPI_SUCCESS) call parallel_abort(error=ierr)
! Duplicate communication space to isolate ELCIRC communication
call mpi_comm_dup(MPI_COMM_WORLD,comm,ierr)
call mpi_comm_dup(MPI_COMM_WORLD,comm_schism,ierr)
endif
if(ierr/=MPI_SUCCESS) call parallel_abort(error=ierr)

! Get number of processors
call mpi_comm_size(comm,nproc,ierr)
call mpi_comm_size(comm_schism,nproc_schism,ierr)
if(ierr/=MPI_SUCCESS) call parallel_abort(error=ierr)

! Get rank
call mpi_comm_rank(comm,myrank,ierr)
call mpi_comm_rank(comm_schism,myrank_schism,ierr)
if(ierr/=MPI_SUCCESS) call parallel_abort(error=ierr)

nproc_compute=nproc_schism-nscribes
if(myrank_schism<nproc_schism-nscribes) then !compute ranks
task_id=1
else !IO ranks
task_id=2
endif

!Use original rank as key to order the new ranks
CALL MPI_Comm_split(comm_schism,task_id,myrank_schism,comm2,ierr)
CALL MPI_Comm_size(comm2, nproc3, ierr)
CALL MPI_Comm_rank(comm2, myrank3,ierr)

if(task_id==1) then !compute
comm=comm2
nproc=nproc3
myrank=myrank3
! print*, 'Compute:',myrank_schism,nproc_schism,myrank,nproc,nproc_compute,nscribes,task_id
else !I/O scribes
comm_scribe=comm2
nproc_scribe=nproc3
myrank_scribe=myrank3
! print*, 'Scribes:',myrank_schism,nproc_schism,myrank3,nproc3,nproc_compute,nscribes,task_id
endif

end subroutine parallel_init


Expand Down Expand Up @@ -424,10 +458,10 @@ subroutine parallel_abort(string,error)
if(lopen) write(11,'(i4,2a)') myrank,': MPI ERROR: ',s
endif
do i=1,500; inquire(i,opened=lopen); if(lopen) close(i); enddo;
call mpi_abort(comm,error,ierror)
call mpi_abort(comm_schism,error,ierror)
else
do i=1,500; inquire(i,opened=lopen); if(lopen) close(i); enddo;
call mpi_abort(comm,0,ierror)
call mpi_abort(comm_schism,0,ierror)
endif

end subroutine parallel_abort
Expand Down
18 changes: 9 additions & 9 deletions src/Core/scribe_io.F90
Original file line number Diff line number Diff line change
Expand Up @@ -52,7 +52,7 @@ module scribe_io
integer,save,allocatable :: np(:),ne(:),ns(:),iplg(:,:),ielg(:,:),islg(:,:),kbp00(:), &
&i34(:),elnode(:,:),rrqst2(:)
real(rkind),save,allocatable :: xnd(:),ynd(:),dp(:)
real(rkind),save,allocatable :: var2dnode(:,:),var3dnode(:,:),var3dnode_gb(:,:)
real(4),save,allocatable :: var2dnode(:,:),var3dnode(:,:),var3dnode_gb(:,:)

public :: scribe_init
public :: scribe_step
Expand Down Expand Up @@ -245,12 +245,12 @@ subroutine scribe_init(indir,iths,ntime)
!Alloc global output arrays for step
! if(ncount_2dnode>0) then
! allocate(var2dnode(ncount_2dnode,np_max))
! var2dnode(ncount_2dnode,np_max)=0.d0 !touch mem
! var2dnode(ncount_2dnode,np_max)=0. !touch mem
! endif

! if(ncount_3dnode>0) then
allocate(var3dnode(nvrt,np_max),var3dnode_gb(nvrt,np_global))
var3dnode_gb(nvrt,np_global)=0.d0 !touch mem
var3dnode_gb(nvrt,np_global)=0. !touch mem
! endif

! Define output var/file names (gfort requires equal length in assignment)
Expand Down Expand Up @@ -280,7 +280,7 @@ subroutine scribe_step(it)
if(myrank_schism==irank) then
!OK to fill partial arrays as long as respect column major
do i=1,nproc_compute
call mpi_irecv(var3dnode,np(i)*nvrt,rtype,i-1,200+itotal,comm_schism,rrqst,ierr)
call mpi_irecv(var3dnode,np(i)*nvrt,MPI_REAL4,i-1,200+itotal,comm_schism,rrqst,ierr)
call mpi_wait(rrqst,MPI_STATUSES_IGNORE,ierr)
var3dnode_gb(1:nvrt,iplg(1:np(i),i))=var3dnode(1:nvrt,1:np(i))
enddo !i
Expand Down Expand Up @@ -310,7 +310,7 @@ subroutine nc_writeout3D(it,idim1,idim2,var3dnode_gb,vname)
implicit none

integer, intent(in) :: it,idim1,idim2
real(rkind), intent(in) :: var3dnode_gb(idim1,idim2)
real(4), intent(in) :: var3dnode_gb(idim1,idim2)
character(len=*), intent(in) :: vname

integer :: irec,iret
Expand Down Expand Up @@ -339,9 +339,9 @@ subroutine nc_writeout3D(it,idim1,idim2,var3dnode_gb,vname)
iret=nf90_put_att(ncid_schism_io,itime_id,'i23d',0) !set i23d flag

time_dims(1)=node_dim
iret=nf90_def_var(ncid_schism_io,'SCHISM_hgrid_node_x',NF90_FLOAT,time_dims,ix_id)
iret=nf90_def_var(ncid_schism_io,'SCHISM_hgrid_node_x',NF90_DOUBLE,time_dims,ix_id)
if(iret.ne.NF90_NOERR) call parallel_abort('nc_writeout3D: xnd')
iret=nf90_def_var(ncid_schism_io,'SCHISM_hgrid_node_y',NF90_FLOAT,time_dims,iy_id)
iret=nf90_def_var(ncid_schism_io,'SCHISM_hgrid_node_y',NF90_DOUBLE,time_dims,iy_id)
if(iret.ne.NF90_NOERR) call parallel_abort('nc_writeout3D: ynd')
iret=nf90_def_var(ncid_schism_io,'depth',NF90_FLOAT,time_dims,ih_id)
if(iret.ne.NF90_NOERR) call parallel_abort('nc_writeout3D: dp')
Expand All @@ -359,8 +359,8 @@ subroutine nc_writeout3D(it,idim1,idim2,var3dnode_gb,vname)
iret=nf90_enddef(ncid_schism_io)

!Write static info (x,y...)
iret=nf90_put_var(ncid_schism_io,ix_id,real(xnd),(/1/),(/np_global/))
iret=nf90_put_var(ncid_schism_io,iy_id,real(ynd),(/1/),(/np_global/))
iret=nf90_put_var(ncid_schism_io,ix_id,xnd,(/1/),(/np_global/))
iret=nf90_put_var(ncid_schism_io,iy_id,ynd,(/1/),(/np_global/))
iret=nf90_put_var(ncid_schism_io,ih_id,real(dp),(/1/),(/np_global/))
! iret=nf90_put_var(ncid_schism_io,i34_id,i34,(/1/),(/ne_global/))
iret=nf90_put_var(ncid_schism_io,elnode_id,elnode,(/1,1/),(/4,ne_global/))
Expand Down
96 changes: 86 additions & 10 deletions src/Driver/schism_driver.F90
Original file line number Diff line number Diff line change
Expand Up @@ -65,41 +65,117 @@
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!! OMP: search for 'new21' for TODO
!'new35' for scribe I/O todo
! List of changed files from official repo:
! Driver/schism_driver.F90
! Corie: schism_[glbl,msgp], scribe_io
! Hydro: schism_[init,step]; renamed single_netcdf
! src/Makefile
!===============================================================================
!===============================================================================
! SCHISM main driver
!===============================================================================
!===============================================================================
program schism_driver
use schism_msgp, only: parallel_init,parallel_finalize,parallel_abort
use schism_msgp, only: nscribes,parallel_init,parallel_finalize,parallel_abort
use schism_version
implicit none
character*8 cli
character(len=20) :: cli

if(COMMAND_ARGUMENT_COUNT()<1) then
print*, 'Must have at least 1 cmd argument (# of scribes)'
stop
endif

call get_command_argument(1,cli)
if (cli(1:2) == "-v")then
print*, ""
call print_version
stop
read(cli,*)nscribes
if(nscribes<0) then
print*, 'nscribes<0:',nscribes
stop
endif
#ifdef OLDIO
if(nscribes/=0) then
print*, 'nscribes must be 0 in this option:',nscribes
stop
endif
#endif

if(COMMAND_ARGUMENT_COUNT()>1) then
call get_command_argument(1,cli)
if (cli(1:2) == "-v")then
print*, ""
call print_version
stop
endif
endif

call parallel_init

!Deal with command args
!call get_command_args
call schism_main
call parallel_finalize

end program schism_driver

subroutine schism_main
use schism_msgp, only: myrank !! debug only
! use schism_msgp, only: myrank !! debug only
implicit none
integer :: it,iths,ntime
call schism_init(0,'./',iths,ntime)

#ifndef OLDIO
call schism_init0(iths,ntime)
do it=iths+1,ntime
call schism_step0(it)
enddo !it
call schism_finalize0
#else
call schism_init(0,'./',iths,ntime)
do it=iths+1,ntime
call schism_step(it)
enddo !it
call schism_finalize
#endif

end subroutine schism_main

subroutine schism_init0(iths,ntime)
use schism_msgp, only: task_id
use scribe_io
implicit none
integer, intent(out) :: iths,ntime

if(task_id==1) then !compute
call schism_init(0,'./',iths,ntime)
else !I/O scribes
call scribe_init('./',iths,ntime)
endif !task_id

end subroutine schism_init0

subroutine schism_step0(it)
use schism_msgp, only: task_id
use scribe_io
implicit none
integer, intent(in) :: it

if(task_id==1) then !compute
call schism_step(it)
else !I/O scribes
call scribe_step(it)
endif !task_id

end subroutine schism_step0

subroutine schism_finalize0
use schism_msgp, only: task_id,comm_schism
use scribe_io
implicit none
integer :: ierr

if(task_id==1) then !compute
call schism_finalize
else !I/O scribes
call scribe_finalize
endif !task_id
call mpi_barrier(comm_schism,ierr)

end subroutine schism_finalize0

0 comments on commit 8efc374

Please sign in to comment.