Skip to content

Commit

Permalink
Offline transport: add T,S to required inputs (overwrite these with
Browse files Browse the repository at this point in the history
hydro only)
  • Loading branch information
josephzhang8 committed Apr 21, 2021
1 parent ec86e0d commit d07d75d
Show file tree
Hide file tree
Showing 5 changed files with 51 additions and 14 deletions.
3 changes: 2 additions & 1 deletion cmake/SCHISM.local.femto
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,7 @@ set(CMAKE_Fortran_COMPILER ifort CACHE PATH "Path to serial Fortran compiler")
set(NetCDF_FORTRAN_DIR "$ENV{NETCDF_FORTRAN}" CACHE PATH "Path to NetCDF Fortran library")
set(NetCDF_C_DIR "$ENV{NETCDF}" CACHE PATH "Path to NetCDF C library")
###Compile flags
#set(CMAKE_Fortran_FLAGS_RELEASE "-O2 -axCORE-AVX2 -xSSE4.2" CACHE STRING "Fortran flags" FORCE)
set(CMAKE_Fortran_FLAGS_RELEASE "-O3 -xCORE-AVX2 -mcmodel=medium -fma -align array64byte -finline-functions" CACHE STRING "Fortran flags" FORCE)

#Hybrid (plz also update exec name above)
#set(CMAKE_Fortran_FLAGS_RELEASE "-O3 -xCORE-AVX2 -mcmodel=medium -fma -align array64byte -finline-functions -qopenmp" CACHE STRING "Fortran flags" FORCE)
2 changes: 1 addition & 1 deletion cmake/SCHISM.local.pleiades
Original file line number Diff line number Diff line change
Expand Up @@ -14,5 +14,5 @@ set(CMAKE_Fortran_FLAGS_RELEASE "-O3 -no-prec-div -mcmodel=medium -assume bytere
##set(CMAKE_EXE_LINKER_FLAGS "-O3 -no-prec-div -mcmodel=medium -assume byterecl -ipo -axCORE-AVX512 -xSSE4.2 -lmpi++ -lmpi -lstdc++ -lcpuset -lbitmask" CACHE STRING "Link flags" FORCE)
set(CMAKE_EXE_LINKER_FLAGS "-O3 -no-prec-div -mcmodel=medium -ipo -axCORE-AVX512 -xSSE4.2 -lmpi++ -lmpi -lstdc++ -lcpuset -lbitmask" CACHE STRING "Link flags" FORCE)

##Hybrid
##Hybrid (plz also update exec name above)
##set(CMAKE_EXE_LINKER_FLAGS "-O3 -no-prec-div -mcmodel=medium -ipo -axCORE-AVX512 -xSSE4.2 -lmpi++ -lmpi -lstdc++ -lcpuset -lbitmask -qopenmp" CACHE STRING "Link flags" FORCE)
4 changes: 2 additions & 2 deletions sample_inputs/param.nml
Original file line number Diff line number Diff line change
Expand Up @@ -47,10 +47,10 @@

!-----------------------------------------------------------------------
! Option to only solve tracer transport (and bypass most b-tropic solver)
! Usage: turn the flag on and turn on your tracer modules and make sure
! Usage: turn the flag on ('1') and turn on your tracer modules and make sure
! (1) your inputs are consistent with the original hydro-only run (with additional tracers
! of course); (2) hydro-only run results are in hydro_out/schout*.nc, which must
! have 'hvel_side' (not normal hvel), 'elev', 'diffusivity'; (3) dt above is
! have 'hvel_side' (not normal hvel), 'elev', 'diffusivity', 'temp_elem', 'salt_elem'; (3) dt above is
! multiple of _output_ step used in the original hydro-only run
! (as found in in hydro_out/schout*.nc); e.g. dt = 1 hour.
! Hotstart should work also, but you'd probably not use an aggressively large dt especially
Expand Down
54 changes: 45 additions & 9 deletions src/Hydro/schism_step.F90
Original file line number Diff line number Diff line change
Expand Up @@ -245,7 +245,7 @@ subroutine schism_step(it)
real(rkind),allocatable :: swild96(:,:,:),swild97(:,:,:) !used in ELAD (deallocate immediately afterwards)
real(rkind),allocatable :: swild95(:,:,:) !for analysis module
real(rkind),allocatable :: swild13(:)
real(4),allocatable :: swild11(:),swild12(:,:) !reading schout*
real(4),allocatable :: swild11(:),swild12(:,:),swild14(:,:,:),ts_offline(:,:,:) !reading schout*
real(rkind),allocatable :: hp_int(:,:,:),buf1(:,:),buf2(:,:),buf3(:),msource(:,:)
real(rkind),allocatable :: fluxes_tr(:,:),fluxes_tr_gb(:,:) !fluxes output between regions
logical :: ltmp,ltmp1(1),ltmp2(1)
Expand Down Expand Up @@ -349,6 +349,12 @@ subroutine schism_step(it)
endif
! End alloc.

! Offline transport
if(itransport_only/=0) then
allocate(ts_offline(2,nvrt,nea),stat=istat)
if(istat/=0) call parallel_abort('MAIN: failed to alloc. (73)')
endif

! do it=iths+1,ntime

#ifdef INCLUDE_TIMING
Expand Down Expand Up @@ -1735,7 +1741,7 @@ subroutine schism_step(it)
!new28: bypass solver for transport only option
if(itransport_only/=0) then
!=================================================================================
!Read in schout (saved hydro outputs), and update new soln: eta2, s[uv]2, dfh.
!Read in schout (saved hydro outputs), and update new soln: eta2, s[uv]2, dfh, tr_el(1:2,:,:).
!Other vars: zcor and dry flags are computed either from schism_init or from levels*() after
!transport solver; similarly for tr_nd*

Expand Down Expand Up @@ -1805,7 +1811,6 @@ subroutine schism_step(it)
!write(16,*)'done reading elev...'
j=nf90_inq_varid(ncid_schout, "diffusivity",mm)
if(j/=NF90_NOERR) call parallel_abort('STEP: nc dfh')
! do k=1,nvrt
!j=nf90_get_var(ncid_schout,mm,swild11(1:np_global),(/k,1,irec2/),(/1,np_global,1/))
j=nf90_get_var(ncid_schout,mm,swild12(:,1:np_global),(/1,1,irec2/),(/nvrt,np_global,1/))
if(j/=NF90_NOERR) call parallel_abort('STEP: nc get dfh')
Expand All @@ -1817,13 +1822,11 @@ subroutine schism_step(it)
dfh(:,ip)=swild12(:,i)
endif
enddo !i
! enddo !k

if(myrank==0) then
!write(16,*)'done reading dfh...'
j=nf90_inq_varid(ncid_schout, "hvel_side",mm)
if(j/=NF90_NOERR) call parallel_abort('STEP: nc hvel')
! do k=1,nvrt
j=nf90_get_var(ncid_schout,mm,swild12,(/1,1,1,irec2/),(/1,nvrt,ns_global,1/))
if(j/=NF90_NOERR) call parallel_abort('STEP: nc get hvel')
!write(16,*)'done reading su2'
Expand All @@ -1835,13 +1838,11 @@ subroutine schism_step(it)
su2(:,isd)=swild12(:,i)
endif
enddo !i
! enddo !k

! do k=1,nvrt
if(myrank==0) then
j=nf90_get_var(ncid_schout,mm,swild12,(/2,1,1,irec2/),(/1,nvrt,ns_global,1/))
if(j/=NF90_NOERR) call parallel_abort('STEP: nc get hvel2')
write(16,*)'finished reading schout...'
! write(16,*)'finished reading sv2...'
endif !myrank=0
call mpi_bcast(swild12,nvrt*ns_global,mpi_real,0,comm,istat)
do i=1,ns_global
Expand All @@ -1850,9 +1851,36 @@ subroutine schism_step(it)
sv2(:,isd)=swild12(:,i)
endif
enddo !i
! enddo !k
deallocate(swild12)

! T,S
allocate(swild14(nvrt,ne_global,2),stat=istat)
if(istat/=0) call parallel_abort('STEP: alloc swild14')
swild14(nvrt,ne_global,2)=0 !test mem

if(myrank==0) then
j=nf90_inq_varid(ncid_schout, "temp_elem",mm)
if(j/=NF90_NOERR) call parallel_abort('STEP: nc temp_elem')
j=nf90_inq_varid(ncid_schout, "salt_elem",jj)
if(j/=NF90_NOERR) call parallel_abort('STEP: nc salt_elem')
j=nf90_get_var(ncid_schout,mm,swild14(:,:,1),(/1,1,irec2/),(/nvrt,ne_global,1/))
if(j/=NF90_NOERR) call parallel_abort('STEP: nc get temp_elem')
j=nf90_get_var(ncid_schout,jj,swild14(:,:,2),(/1,1,irec2/),(/nvrt,ne_global,1/))
if(j/=NF90_NOERR) call parallel_abort('STEP: nc get salt_elem')
write(16,*)'done reading T,S and schout_'
endif !myrank=0
call mpi_bcast(swild14,2*nvrt*ne_global,mpi_real,0,comm,istat)
if(istat/=MPI_SUCCESS) call parallel_abort('STEP: mpi_bcast in reading T,S')
do i=1,ne_global
if(iegl(i)%rank==myrank) then
ie=iegl(i)%id !up to nea
ts_offline(1,:,ie)=swild14(:,i,1)
ts_offline(2,:,ie)=swild14(:,i,2)
endif
enddo !i

deallocate(swild14)

! Deal with junks
where(abs(dfh)>1.d3) dfh=1.d-6
where(abs(su2)>1.d2) su2=0.d0
Expand Down Expand Up @@ -7051,6 +7079,13 @@ subroutine schism_step(it)
!$OMP end single
#endif /*USE_AGE*/

! Overwrite T,S for offline transport option
if(itransport_only/=0) then
!$OMP workshare
tr_el(1:2,:,1:nea)=ts_offline(1:2,:,:)
!$OMP end workshare
endif !itransport_only/

! Convert to nodes and whole levels
!$OMP do
do i=1,nea
Expand Down Expand Up @@ -9361,6 +9396,7 @@ subroutine schism_step(it)
deallocate(hp_int,uth,vth,d2uv,dr_dxy,bcc)
if(allocated(rwild)) deallocate(rwild)
deallocate(swild9)
if(allocated(ts_offline)) deallocate(ts_offline)

#ifdef USE_NAPZD
deallocate(Bio_bdefp)
Expand Down
2 changes: 1 addition & 1 deletion src/Readme.beta_notes
Original file line number Diff line number Diff line change
Expand Up @@ -480,7 +480,7 @@ git versions:
(118) 28ff9d1 (Dec 16, 2020): added optional self-attraction loading tides. The option shares some constants
with tidal potential: freq names and cut-off depth;
(119) b1bcaa0 (April 20, 2021): removed most of 'goto'. Remaining ones: harm.F90, lap.F90, WWM
(120) 0cec024 (April 21, 2021): bug fix in misc_subs c/o of Fei (nodeA in final extrap stage may be interface node)
(120) 0cec024 (April 21, 2021): bug fix in misc_subs c/o of Fei (inunfl=1: nodeA in final extrap stage may be interface node)
---------------------------------------------------------------
Auto-test history:
R100 (Nov 2011); R168 (minus WWM); R224 (basic); R240 (full); R306 (basic); R370 (basic); R408 (basic); R430 (basic); R601 (basic); R730 (basic); R747 (all hydro); R1120 (hydro+ SF BayDelta); R1305: (all hydro; ICM); R1532 (branch/selfe_opt1): all; R1641: hydro; R2018: hydro; R2232 (all);
Expand Down

0 comments on commit d07d75d

Please sign in to comment.