From de59adf02748b001a1b06860c5bd050f3cfb2f39 Mon Sep 17 00:00:00 2001 From: Robert Hallberg Date: Wed, 20 Mar 2024 09:29:38 -0400 Subject: [PATCH] Make the US argument to MOM_domains_init optional This commit makes the unit_scale_type argument US to MOM_domains_init and gen_auto_mask_table optional and moves it to the end of the argument list, so that coupled or ice-ocean models using SIS2 will compile with the proposed updates to the main branch of MOM6 from dev/ncar. Because MOM6 and SIS2 use some common framework code but are managed in separate github repositories, we need to use optional argument to allow a single version of SIS2 to work across changes to MOM6 interfaces. Because the TOPO_CONFIG parameter as used in SIS2 has a default value, there is an alternative call to get_param for TOPO_CONFIG with a default when MOM_domains_init is called with a domain_name argument. Also added missing scale arguments to get_param calls for MINIMUM_DEPTH and MASKING_DEPTH. This commit also adds or corrects units in the comments describing 4 recently added or modified variables. All answers are bitwise identical in any cases that worked before (noting that some cases using SIS2 would not even compile). --- src/core/MOM.F90 | 12 ++--- src/framework/MOM_domains.F90 | 54 ++++++++++++------- src/ice_shelf/MOM_ice_shelf.F90 | 4 +- src/ocean_data_assim/MOM_oda_driver.F90 | 2 +- .../lateral/MOM_hor_visc.F90 | 2 +- 5 files changed, 45 insertions(+), 29 deletions(-) diff --git a/src/core/MOM.F90 b/src/core/MOM.F90 index b7f8bd3f66..a9c8c5cd9e 100644 --- a/src/core/MOM.F90 +++ b/src/core/MOM.F90 @@ -2538,13 +2538,13 @@ subroutine initialize_MOM(Time, Time_init, param_file, dirs, CS, & #endif G_in => CS%G_in #ifdef STATIC_MEMORY_ - call MOM_domains_init(G_in%domain, US, param_file, symmetric=symmetric, & - static_memory=.true., NIHALO=NIHALO_, NJHALO=NJHALO_, & - NIGLOBAL=NIGLOBAL_, NJGLOBAL=NJGLOBAL_, NIPROC=NIPROC_, & - NJPROC=NJPROC_) + call MOM_domains_init(G_in%domain, param_file, symmetric=symmetric, & + static_memory=.true., NIHALO=NIHALO_, NJHALO=NJHALO_, & + NIGLOBAL=NIGLOBAL_, NJGLOBAL=NJGLOBAL_, NIPROC=NIPROC_, & + NJPROC=NJPROC_, US=US) #else - call MOM_domains_init(G_in%domain, US, param_file, symmetric=symmetric, & - domain_name="MOM_in") + call MOM_domains_init(G_in%domain, param_file, symmetric=symmetric, & + domain_name="MOM_in", US=US) #endif ! Copy input grid (G_in) domain to active grid G diff --git a/src/framework/MOM_domains.F90 b/src/framework/MOM_domains.F90 index f2c3225025..d937ed7b0c 100644 --- a/src/framework/MOM_domains.F90 +++ b/src/framework/MOM_domains.F90 @@ -63,12 +63,11 @@ module MOM_domains !> MOM_domains_init initializes a MOM_domain_type variable, based on the information !! read in from a param_file_type, and optionally returns data describing various !! properties of the domain type. -subroutine MOM_domains_init(MOM_dom, US, param_file, symmetric, static_memory, & +subroutine MOM_domains_init(MOM_dom, param_file, symmetric, static_memory, & NIHALO, NJHALO, NIGLOBAL, NJGLOBAL, NIPROC, NJPROC, & - min_halo, domain_name, include_name, param_suffix) + min_halo, domain_name, include_name, param_suffix, US) type(MOM_domain_type), pointer :: MOM_dom !< A pointer to the MOM_domain_type !! being defined here. - type(unit_scale_type), pointer :: US !< A dimensional unit scaling type type(param_file_type), intent(in) :: param_file !< A structure to parse for !! run-time parameters logical, optional, intent(in) :: symmetric !< If present, this specifies @@ -99,6 +98,7 @@ subroutine MOM_domains_init(MOM_dom, US, param_file, symmetric, static_memory, & !! "MOM_memory.h" if missing. character(len=*), optional, intent(in) :: param_suffix !< A suffix to apply to !! layout-specific parameters. + type(unit_scale_type), optional, pointer :: US !< A dimensional unit scaling type ! Local variables integer, dimension(2) :: layout ! The number of logical processors in the i- and j- directions @@ -120,6 +120,7 @@ subroutine MOM_domains_init(MOM_dom, US, param_file, symmetric, static_memory, & logical :: auto_mask_table ! Runtime flag that turns on automatic mask table generator integer :: auto_io_layout_fac ! Used to compute IO layout when auto_mask_table is True. logical :: mask_table_exists ! True if there is a mask table file + logical :: is_MOM_domain ! True if this domain is being set for MOM, and not another component like SIS2. character(len=128) :: inputdir ! The directory in which to find the diag table character(len=200) :: mask_table ! The file name and later the full path to the diag table character(len=64) :: inc_nm ! The name of the memory include file @@ -288,7 +289,16 @@ subroutine MOM_domains_init(MOM_dom, US, param_file, symmetric, static_memory, & call get_param(param_file, mdl, "INPUTDIR", inputdir, do_not_log=.true., default=".") inputdir = slasher(inputdir) - call get_param(param_file, mdl, "TOPO_CONFIG", topo_config, do_not_log=.true., fail_if_missing=.true.) + is_MOM_domain = .true. + if (present(domain_name)) then + is_MOM_domain = (index(domain_name, "MOM") > 1) + endif + + if (is_MOM_domain) then + call get_param(param_file, mdl, "TOPO_CONFIG", topo_config, do_not_log=.true., fail_if_missing=.true.) + else ! SIS2 has a default value for TOPO_CONFIG. + call get_param(param_file, mdl, "TOPO_CONFIG", topo_config, default="file", do_not_log=.true.) + endif auto_mask_table = .false. if (.not. present(param_suffix) .and. .not. is_static .and. trim(topo_config) == 'file') then @@ -314,7 +324,7 @@ subroutine MOM_domains_init(MOM_dom, US, param_file, symmetric, static_memory, & call cpu_clock_begin(id_clock_auto_mask) if (is_root_PE()) then call gen_auto_mask_table(n_global, reentrant, tripolar_N, PEs_used, param_file, inputdir, & - auto_mask_table_fname, US, auto_layout) + auto_mask_table_fname, auto_layout, US) endif call broadcast(auto_layout, length=2) call cpu_clock_end(id_clock_auto_mask) @@ -462,17 +472,18 @@ subroutine MOM_define_layout(n_global, ndivs, layout) end subroutine MOM_define_layout !> Given a desired number of active npes, generate a layout and mask_table -subroutine gen_auto_mask_table(n_global, reentrant, tripolar_N, npes, param_file, inputdir, filename, US, layout) +subroutine gen_auto_mask_table(n_global, reentrant, tripolar_N, npes, param_file, inputdir, filename, layout, US) integer, dimension(2), intent(in) :: n_global !< The total number of gridpoints in 2 directions logical, dimension(2), intent(in) :: reentrant !< True if the x- and y- directions are periodic. - logical :: tripolar_N !< A flag indicating whether there is n. tripolar connectivity + logical, intent(in) :: tripolar_N !< A flag indicating whether there is n. tripolar connectivity integer, intent(in) :: npes !< The desired number of active PEs. type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters - character(len=128), intent(in) :: inputdir !< INPUTDIR parameter + character(len=128), intent(in) :: inputdir !< INPUTDIR parameter character(len=:), allocatable, intent(in) :: filename !< Mask table file path (to be auto-generated.) - type(unit_scale_type), pointer :: US !< A dimensional unit scaling type integer, dimension(2), intent(out) :: layout !< The generated layout of PEs (incl. masked blocks) - !local + type(unit_scale_type), optional, pointer :: US !< A dimensional unit scaling type + + ! Local variables real, dimension(n_global(1), n_global(2)) :: D ! Bathymetric depth (to be read in from TOPO_FILE) [Z ~> m] integer, dimension(:,:), allocatable :: mask ! Cell masks (based on D and MINIMUM_DEPTH) character(len=200) :: topo_filepath, topo_file ! Strings for file/path @@ -483,18 +494,23 @@ subroutine gen_auto_mask_table(n_global, reentrant, tripolar_N, npes, param_file real :: Dmask ! The depth for masking in the same units as D [Z ~> m] real :: min_depth ! The minimum ocean depth in the same units as D [Z ~> m] real :: mask_depth ! The depth shallower than which to mask a point as land. [Z ~> m] - real :: glob_ocn_frac ! ratio of ocean points to total number of points - real :: r_p ! aspect ratio for division count p. + real :: glob_ocn_frac ! ratio of ocean points to total number of points [nondim] + real :: r_p ! aspect ratio for division count p. [nondim] + real :: m_to_Z ! A conversion factor from m to height units [Z m-1 ~> 1] integer :: nx, ny ! global domain sizes integer, parameter :: ibuf=2, jbuf=2 - real, parameter :: r_extreme = 4.0 ! aspect ratio limit (>1) for a layout to be considered. + real, parameter :: r_extreme = 4.0 ! aspect ratio limit (>1) for a layout to be considered [nondim] integer :: num_masked_blocks integer, allocatable :: mask_table(:,:) + m_to_Z = 1.0 ; if (present(US)) m_to_Z = US%m_to_Z + ! Read in params necessary for auto-masking - call get_param(param_file, mdl, "MINIMUM_DEPTH", min_depth, do_not_log=.true., units="m", default=0.0) - call get_param(param_file, mdl, "MASKING_DEPTH", mask_depth, do_not_log=.true., units="m", default=-9999.0) - call get_param(param_file, mdl, "TOPO_CONFIG", topo_config, do_not_log=.true., fail_if_missing=.true.) + call get_param(param_file, mdl, "MINIMUM_DEPTH", min_depth, & + units="m", default=0.0, scale=m_to_Z, do_not_log=.true.) + call get_param(param_file, mdl, "MASKING_DEPTH", mask_depth, & + units="m", default=-9999.0, scale=m_to_Z, do_not_log=.true.) + call get_param(param_file, mdl, "TOPO_CONFIG", topo_config, default="file", do_not_log=.true.) call get_param(param_file, mdl, "TOPO_FILE", topo_file, do_not_log=.true., default="topog.nc") call get_param(param_file, mdl, "TOPO_VARNAME", topo_varname, do_not_log=.true., default="depth") topo_filepath = trim(inputdir)//trim(topo_file) @@ -514,15 +530,15 @@ subroutine gen_auto_mask_table(n_global, reentrant, tripolar_N, npes, param_file ny = n_global(2) ! Read in bathymetric depth. - D(:,:) = -9.0e30 * US%m_to_Z ! Initializing to a very large negative depth (tall mountains) everywhere. + D(:,:) = -9.0e30 * m_to_Z ! Initializing to a very large negative depth (tall mountains) everywhere. call read_field(topo_filepath, trim(topo_varname), D, start=(/1, 1/), nread=n_global, no_domain=.true., & - scale=US%m_to_Z) + scale=m_to_Z) allocate(mask(nx+2*ibuf, ny+2*jbuf), source=0) ! Determine cell masks Dmask = mask_depth - if (mask_depth == -9999.0) Dmask = min_depth + if (mask_depth == -9999.0*m_to_Z) Dmask = min_depth do i=1,nx ; do j=1,ny if (D(i,j) <= Dmask) then mask(i+ibuf,j+jbuf) = 0 diff --git a/src/ice_shelf/MOM_ice_shelf.F90 b/src/ice_shelf/MOM_ice_shelf.F90 index eab178280c..b9640502d2 100644 --- a/src/ice_shelf/MOM_ice_shelf.F90 +++ b/src/ice_shelf/MOM_ice_shelf.F90 @@ -1322,8 +1322,8 @@ subroutine initialize_ice_shelf(param_file, ocn_grid, Time, CS, diag, Time_init, ! Set up the ice-shelf domain and grid wd_halos(:)=0 allocate(CS%Grid) - call MOM_domains_init(CS%Grid%domain, CS%US, param_file, min_halo=wd_halos, symmetric=GRID_SYM_,& - domain_name='MOM_Ice_Shelf_in') + call MOM_domains_init(CS%Grid%domain, param_file, min_halo=wd_halos, symmetric=GRID_SYM_,& + domain_name='MOM_Ice_Shelf_in', US=CS%US) !allocate(CS%Grid_in%HI) !call hor_index_init(CS%Grid%Domain, CS%Grid%HI, param_file, & ! local_indexing=.not.global_indexing) diff --git a/src/ocean_data_assim/MOM_oda_driver.F90 b/src/ocean_data_assim/MOM_oda_driver.F90 index f45939d007..6e24b9faee 100644 --- a/src/ocean_data_assim/MOM_oda_driver.F90 +++ b/src/ocean_data_assim/MOM_oda_driver.F90 @@ -291,7 +291,7 @@ subroutine init_oda(Time, G, GV, US, diag_CS, CS) CS%G => G allocate(CS%Grid) ! params NIHALO_ODA, NJHALO_ODA set the DA halo size - call MOM_domains_init(CS%Grid%Domain,CS%US,PF,param_suffix='_ODA') + call MOM_domains_init(CS%Grid%Domain, PF, param_suffix='_ODA', US=CS%US) allocate(HI) call hor_index_init(CS%Grid%Domain, HI, PF) call verticalGridInit( PF, CS%GV, CS%US ) diff --git a/src/parameterizations/lateral/MOM_hor_visc.F90 b/src/parameterizations/lateral/MOM_hor_visc.F90 index 9b1d81348e..1b73a7ec12 100644 --- a/src/parameterizations/lateral/MOM_hor_visc.F90 +++ b/src/parameterizations/lateral/MOM_hor_visc.F90 @@ -300,7 +300,7 @@ subroutine horizontal_viscosity(u, v, h, diffu, diffv, MEKE, VarMix, G, GV, US, GME_effic_h, & ! The filtered efficiency of the GME terms at h points [nondim] m_leithy, & ! Kh=m_leithy*Ah in Leith+E parameterization [L-2 ~> m-2] Ah_sq, & ! The square of the biharmonic viscosity [L8 T-2 ~> m8 s-2] - htot ! The total thickness of all layers [Z ~> m] + htot ! The total thickness of all layers [H ~> m or kg m-2] real :: Del2vort_h ! Laplacian of vorticity at h-points [L-2 T-1 ~> m-2 s-1] real :: grad_vel_mag_bt_h ! Magnitude of the barotropic velocity gradient tensor squared at h-points [T-2 ~> s-2] real :: boundary_mask_h ! A mask that zeroes out cells with at least one land edge [nondim]