Skip to content

Commit

Permalink
Add REMAPPING_SCHEME for OBC, ODA and SPONGE (#751)
Browse files Browse the repository at this point in the history
Most modules that use remapping have their own parameter e.g.
Z_INIT_REMAPPING_SCHEME.  However, OBC, ODA and SPONGE currently
sliently use the primary ALE REMAPPING_SCHEME parameter.

Added logged parameters OBC_REMAPPING_SCHEME, ODA_REMAPPING_SCHEME
and SPONGE_REMAPPING_SCHEME to allow customization.  All take
REMAPPING_SCHEME as their default to maintain compatibility.

For ODA, also added ODA_REMAPPING_USE_OM4_SUBCELLS that takes
REMAPPING_USE_OM4_SUBCELLS as its default to maintain compatibility,
and ODA_BOUNDARY_EXTRAP that takes BOUNDARY_EXTRAPOLATION as its
default to maintain compatibility.  Note that BOUNDARY_EXTRAPOLATION
is a bug, it should have been REMAP_BOUNDARY_EXTRAP, but the former
remains the default for compatibility.

For SPONGE, also added SPONGE_BOUNDARY_EXTRAP that takes
BOUNDARY_EXTRAPOLATION as its default to maintain compatibility.
Note that BOUNDARY_EXTRAPOLATION is a bug, it should have been
REMAP_BOUNDARY_EXTRAP, but the former remains the default for
compatibility.

Answers are unchanged, but there are new logged parameters.
  • Loading branch information
awallcraft authored Nov 8, 2024
1 parent 13cc946 commit 31a4d8b
Show file tree
Hide file tree
Showing 4 changed files with 38 additions and 24 deletions.
9 changes: 5 additions & 4 deletions src/core/MOM_open_boundary.F90
Original file line number Diff line number Diff line change
Expand Up @@ -448,9 +448,8 @@ subroutine open_boundary_config(G, US, param_file, OBC)
character(len=200) :: config1 ! String for OBC_USER_CONFIG
real :: Lscale_in, Lscale_out ! parameters controlling tracer values at the boundaries [L ~> m]
integer :: default_answer_date ! The default setting for the various ANSWER_DATE flags.
logical :: check_reconstruction, check_remapping, force_bounds_in_subcell
logical :: check_remapping, force_bounds_in_subcell
logical :: om4_remap_via_sub_cells ! If true, use the OM4 remapping algorithm
character(len=64) :: remappingScheme
! This include declares and sets the variable "version".
# include "version_variable.h"

Expand Down Expand Up @@ -676,10 +675,12 @@ subroutine open_boundary_config(G, US, param_file, OBC)
enddo

call get_param(param_file, mdl, "REMAPPING_SCHEME", OBC%remappingScheme, &
default=remappingDefaultScheme, do_not_log=.true.)
call get_param(param_file, mdl, "OBC_REMAPPING_SCHEME", OBC%remappingScheme, &
"This sets the reconstruction scheme used "//&
"for vertical remapping for all variables. "//&
"for OBC vertical remapping for all variables. "//&
"It can be one of the following schemes: \n"//&
trim(remappingSchemesDoc), default=remappingDefaultScheme,do_not_log=.true.)
trim(remappingSchemesDoc), default=OBC%remappingScheme)
call get_param(param_file, mdl, "FATAL_CHECK_RECONSTRUCTIONS", OBC%check_reconstruction, &
"If true, cell-by-cell reconstructions are checked for "//&
"consistency and if non-monotonicity or an inconsistency is "//&
Expand Down
7 changes: 5 additions & 2 deletions src/ocean_data_assim/MOM_oda_driver.F90
Original file line number Diff line number Diff line change
Expand Up @@ -245,7 +245,7 @@ subroutine init_oda(Time, G, GV, US, diag_CS, CS)
call get_param(PF, mdl, "INPUTDIR", inputdir)
call get_param(PF, mdl, "ODA_REMAPPING_SCHEME", remap_scheme, &
"This sets the reconstruction scheme used "//&
"for vertical remapping for all variables. "//&
"for vertical remapping for all ODA variables. "//&
"It can be one of the following schemes: "//&
trim(remappingSchemesDoc), default="PPM_H4")
call get_param(PF, mdl, "DEFAULT_ANSWER_DATE", default_answer_date, &
Expand Down Expand Up @@ -324,7 +324,10 @@ subroutine init_oda(Time, G, GV, US, diag_CS, CS)
default="ZSTAR", fail_if_missing=.false.)
call get_param(PF, mdl, "REMAPPING_USE_OM4_SUBCELLS", om4_remap_via_sub_cells, &
do_not_log=.true., default=.true.)

call get_param(PF, mdl, "ODA_REMAPPING_USE_OM4_SUBCELLS", om4_remap_via_sub_cells, &
"If true, use the OM4 remapping-via-subcells algorithm for ODA. "//&
"See REMAPPING_USE_OM4_SUBCELLS for more details. "//&
"We recommend setting this option to false.", default=om4_remap_via_sub_cells)
call initialize_regridding(CS%regridCS, CS%GV, CS%US, dG%max_depth,PF,'oda_driver',coord_mode,'','')

h_neglect = set_h_neglect(GV, CS%answer_date, h_neglect_edge)
Expand Down
20 changes: 14 additions & 6 deletions src/ocean_data_assim/MOM_oda_incupd.F90
Original file line number Diff line number Diff line change
Expand Up @@ -25,6 +25,7 @@ module MOM_oda_incupd
use MOM_grid, only : ocean_grid_type
use MOM_io, only : vardesc, var_desc
use MOM_remapping, only : remapping_cs, remapping_core_h, initialize_remapping
use MOM_remapping, only : remappingSchemesDoc
use MOM_restart, only : register_restart_field, register_restart_pair, MOM_restart_CS
use MOM_restart, only : restart_init, save_restart, query_initialized
use MOM_spatial_means, only : global_i_mean
Expand Down Expand Up @@ -184,22 +185,29 @@ subroutine initialize_oda_incupd( G, GV, US, param_file, CS, data_h, nz_data, re
"use U,V increments.", &
default=.true.)
call get_param(param_file, mdl, "REMAPPING_SCHEME", remapScheme, &
"This sets the reconstruction scheme used "//&
" for vertical remapping for all variables.", &
default="PLM", do_not_log=.true.)
call get_param(param_file, mdl, "ODA_REMAPPING_SCHEME", remapScheme, &
"This sets the reconstruction scheme used "//&
"for vertical remapping for all ODA variables. "//&
"It can be one of the following schemes: "//&
trim(remappingSchemesDoc), default=remapScheme)

!The default should be REMAP_BOUNDARY_EXTRAP
call get_param(param_file, mdl, "BOUNDARY_EXTRAPOLATION", bndExtrapolation, &
"When defined, a proper high-order reconstruction "//&
"scheme is used within boundary cells rather "//&
"than PCM. E.g., if PPM is used for remapping, a "//&
"PPM reconstruction will also be used within boundary cells.", &
default=.false., do_not_log=.true.)
call get_param(param_file, mdl, "ODA_BOUNDARY_EXTRAP", bndExtrapolation, &
"If true, values at the interfaces of boundary cells are "//&
"extrapolated instead of piecewise constant", default=bndExtrapolation)
call get_param(param_file, mdl, "ODA_INCUPD_DATA_ONGRID", CS%incupdDataOngrid, &
"When defined, the incoming oda_incupd data are "//&
"assumed to be on the model horizontal grid " , &
default=.true.)
call get_param(param_file, mdl, "REMAPPING_USE_OM4_SUBCELLS", om4_remap_via_sub_cells, &
do_not_log=.true., default=.true.)
call get_param(param_file, mdl, "ODA_REMAPPING_USE_OM4_SUBCELLS", om4_remap_via_sub_cells, &
"If true, use the OM4 remapping-via-subcells algorithm for ODA. "//&
"See REMAPPING_USE_OM4_SUBCELLS for more details. "//&
"We recommend setting this option to false.", default=om4_remap_via_sub_cells)

CS%nz = GV%ke

Expand Down
26 changes: 14 additions & 12 deletions src/parameterizations/vertical/MOM_ALE_sponge.F90
Original file line number Diff line number Diff line change
Expand Up @@ -215,16 +215,17 @@ subroutine initialize_ALE_sponge_fixed(Iresttime, G, GV, param_file, CS, data_h,
default=.false.)

call get_param(param_file, mdl, "REMAPPING_SCHEME", remapScheme, &
"This sets the reconstruction scheme used "//&
" for vertical remapping for all variables.", &
default="PLM", do_not_log=.true.)
call get_param(param_file, mdl, "SPONGE_REMAPPING_SCHEME", remapScheme, &
"This sets the reconstruction scheme used "//&
"for vertical remapping for all SPONGE variables.", default=remapScheme)

!This default should be from REMAP_BOUNDARY_EXTRAP
call get_param(param_file, mdl, "BOUNDARY_EXTRAPOLATION", bndExtrapolation, &
"When defined, a proper high-order reconstruction "//&
"scheme is used within boundary cells rather "//&
"than PCM. E.g., if PPM is used for remapping, a "//&
"PPM reconstruction will also be used within boundary cells.", &
default=.false., do_not_log=.true.)
call get_param(param_file, mdl, "SPONGE_BOUNDARY_EXTRAP", bndExtrapolation, &
"If true, values at the interfaces of SPONGE boundary cells are "//&
"extrapolated instead of piecewise constant", default=bndExtrapolation)
call get_param(param_file, mdl, "DEFAULT_ANSWER_DATE", default_answer_date, &
"This sets the default value for the various _ANSWER_DATE parameters.", &
default=99991231)
Expand Down Expand Up @@ -495,15 +496,16 @@ subroutine initialize_ALE_sponge_varying(Iresttime, G, GV, US, param_file, CS, I
"Apply sponges in u and v, in addition to tracers.", &
default=.false.)
call get_param(param_file, mdl, "REMAPPING_SCHEME", remapScheme, &
"This sets the reconstruction scheme used "//&
" for vertical remapping for all variables.", &
default="PLM", do_not_log=.true.)
call get_param(param_file, mdl, "SPONGE_REMAPPING_SCHEME", remapScheme, &
"This sets the reconstruction scheme used "//&
"for vertical remapping for all SPONGE variables.", default=remapScheme)
!This default should be from REMAP_BOUNDARY_EXTRAP
call get_param(param_file, mdl, "BOUNDARY_EXTRAPOLATION", bndExtrapolation, &
"When defined, a proper high-order reconstruction "//&
"scheme is used within boundary cells rather "//&
"than PCM. E.g., if PPM is used for remapping, a "//&
"PPM reconstruction will also be used within boundary cells.", &
default=.false., do_not_log=.true.)
call get_param(param_file, mdl, "SPONGE_BOUNDARY_EXTRAP", bndExtrapolation, &
"If true, values at the interfaces of SPONGE boundary cells are "//&
"extrapolated instead of piecewise constant", default=bndExtrapolation)
call get_param(param_file, mdl, "VARYING_SPONGE_MASK_THICKNESS", CS%varying_input_dz_mask, &
"An input file thickness below which the target values with "//&
"time-varying sponges are replaced by the value above.", &
Expand Down

0 comments on commit 31a4d8b

Please sign in to comment.