MOM6
mom_ale Module Reference

Detailed Description

This module contains the main regridding routines.

Regridding comprises two steps:

  1. Interpolation and creation of a new grid based on target interface densities (or any other criterion).
  2. Remapping of quantities between old grid and new grid.

Original module written by Laurent White, 2008.06.09

Data Types

type  ale_cs
 ALE control structure. More...
 

Functions/Subroutines

subroutine, public ale_init (param_file, GV, US, max_depth, CS)
 This routine is typically called (from initialize_MOM in file MOM.F90) before the main time integration loop to initialize the regridding stuff. We read the MOM_input file to register the values of different regridding/remapping parameters. More...
 
subroutine, public ale_register_diags (Time, G, GV, US, diag, CS)
 Initialize diagnostics for the ALE module. More...
 
subroutine, public adjustgridforintegrity (CS, G, GV, h)
 Crudely adjust (initial) grid for integrity. This routine is typically called (from initialize_MOM in file MOM.F90) before the main time integration loop to initialize the regridding stuff. We read the MOM_input file to register the values of different regridding/remapping parameters. More...
 
subroutine, public ale_end (CS)
 End of regridding (memory deallocation). This routine is typically called (from MOM_end in file MOM.F90) after the main time integration loop to deallocate the regridding stuff. More...
 
subroutine, public ale_main (G, GV, US, h, u, v, tv, Reg, CS, OBC, dt, frac_shelf_h)
 Takes care of (1) building a new grid and (2) remapping all variables between the old grid and the new grid. The creation of the new grid can be based on z coordinates, target interface densities, sigma coordinates or any arbitrary coordinate system. More...
 
subroutine, public ale_main_offline (G, GV, h, tv, Reg, CS, OBC, dt)
 Takes care of (1) building a new grid and (2) remapping all variables between the old grid and the new grid. The creation of the new grid can be based on z coordinates, target interface densities, sigma coordinates or any arbitrary coordinate system. More...
 
subroutine, public ale_offline_inputs (CS, G, GV, h, tv, Reg, uhtr, vhtr, Kd, debug, OBC)
 Regrid/remap stored fields used for offline tracer integrations. These input fields are assumed to have the same layer thicknesses at the end of the last offline interval (which should be a Zstar grid). This routine builds a grid on the runtime specified vertical coordinate. More...
 
subroutine, public ale_offline_tracer_final (G, GV, h, tv, h_target, Reg, CS, OBC)
 Remaps all tracers from h onto h_target. This is intended to be called when tracers are done offline. In the case where transports don't quite conserve, we still want to make sure that layer thicknesses offline do not drift too far away from the online model. More...
 
subroutine check_grid (G, GV, h, threshold)
 Check grid for negative thicknesses. More...
 
subroutine, public ale_build_grid (G, GV, regridCS, remapCS, h, tv, debug, frac_shelf_h)
 Generates new grid. More...
 
subroutine, public ale_regrid_accelerated (CS, G, GV, h, tv, n, u, v, OBC, Reg, dt, dzRegrid, initial)
 For a state-based coordinate, accelerate the process of regridding by repeatedly applying the grid calculation algorithm. More...
 
subroutine remap_all_state_vars (CS_remapping, CS_ALE, G, GV, h_old, h_new, Reg, OBC, dxInterface, u, v, debug, dt)
 This routine takes care of remapping all variable between the old and the new grids. When velocity components need to be remapped, thicknesses at velocity points are taken to be arithmetic averages of tracer thicknesses. This routine is called during initialization of the model at time=0, to remap initiali conditions to the model grid. It is also called during a time step to update the state. More...
 
subroutine, public ale_remap_scalar (CS, G, GV, nk_src, h_src, s_src, h_dst, s_dst, all_cells, old_remap, answers_2018)
 Remaps a single scalar between grids described by thicknesses h_src and h_dst. h_dst must be dimensioned as a model array with GVke layers while h_src can have an arbitrary number of layers specified by nk_src. More...
 
subroutine, public ts_plm_edge_values (CS, S_t, S_b, T_t, T_b, G, GV, tv, h, bdry_extrap)
 Calculate edge values (top and bottom of layer) for T and S consistent with a PLM reconstruction in the vertical direction. Boundary reconstructions are PCM unless bdry_extrap is true. More...
 
subroutine, public ale_plm_edge_values (CS, G, GV, h, Q, bdry_extrap, Q_t, Q_b)
 Calculate edge values (top and bottom of layer) 3d scalar array. Boundary reconstructions are PCM unless bdry_extrap is true. More...
 
subroutine, public ts_ppm_edge_values (CS, S_t, S_b, T_t, T_b, G, GV, tv, h, bdry_extrap)
 Calculate edge values (top and bottom of layer) for T and S consistent with a PPM reconstruction in the vertical direction. Boundary reconstructions are PCM unless bdry_extrap is true. More...
 
subroutine, public ale_initregridding (GV, US, max_depth, param_file, mdl, regridCS)
 Initializes regridding for the main ALE algorithm. More...
 
real function, dimension(cs%nk+1), public ale_getcoordinate (CS)
 Query the target coordinate interfaces positions. More...
 
character(len=20) function, public ale_getcoordinateunits (CS)
 Query the target coordinate units. More...
 
logical function, public ale_remap_init_conds (CS)
 Returns true if initial conditions should be regridded and remapped. More...
 
subroutine, public ale_update_regrid_weights (dt, CS)
 Updates the weights for time filtering the new grid generated in regridding. More...
 
subroutine, public ale_updateverticalgridtype (CS, GV)
 Update the vertical grid type with ALE information. This subroutine sets information in the verticalGrid_type to be consistent with the use of ALE mode. More...
 
subroutine, public ale_writecoordinatefile (CS, GV, directory)
 Write the vertical coordinate information into a file. This subroutine writes out a file containing any available data related to the vertical grid used by the MOM ocean model when in ALE mode. More...
 
subroutine, public ale_initthicknesstocoord (CS, G, GV, h)
 Set h to coordinate values for fixed coordinate systems. More...
 

Function/Subroutine Documentation

◆ adjustgridforintegrity()

subroutine, public mom_ale::adjustgridforintegrity ( type(ale_cs), pointer  CS,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  h 
)

Crudely adjust (initial) grid for integrity. This routine is typically called (from initialize_MOM in file MOM.F90) before the main time integration loop to initialize the regridding stuff. We read the MOM_input file to register the values of different regridding/remapping parameters.

Parameters
csRegridding parameters and options
[in]gOcean grid informations
[in]gvOcean vertical grid structure
[in,out]hCurrent 3D grid thickness that are to be adjusted [H ~> m or kg-2]

Definition at line 292 of file MOM_ALE.F90.

293  type(ALE_CS), pointer :: CS !< Regridding parameters and options
294  type(ocean_grid_type), intent(in) :: G !< Ocean grid informations
295  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
296  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(inout) :: h !< Current 3D grid thickness that
297  !! are to be adjusted [H ~> m or kg-2]
298  call inflate_vanished_layers_old( cs%regridCS, g, gv, h(:,:,:) )
299 

◆ ale_build_grid()

subroutine, public mom_ale::ale_build_grid ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(regridding_cs), intent(in)  regridCS,
type(remapping_cs), intent(in)  remapCS,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(inout)  h,
type(thermo_var_ptrs), intent(inout)  tv,
logical, intent(in), optional  debug,
real, dimension(:,:), optional, pointer  frac_shelf_h 
)

Generates new grid.

Parameters
[in]gOcean grid structure
[in]gvOcean vertical grid structure
[in]regridcsRegridding parameters and options
[in]remapcsRemapping parameters and options
[in,out]tvThermodynamical variable structure
[in,out]hCurrent 3D grid obtained after the last time step [H ~> m or kg-2]
[in]debugIf true, show the call tree
frac_shelf_hFractional ice shelf coverage

Definition at line 615 of file MOM_ALE.F90.

616  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
617  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
618  type(regridding_CS), intent(in) :: regridCS !< Regridding parameters and options
619  type(remapping_CS), intent(in) :: remapCS !< Remapping parameters and options
620  type(thermo_var_ptrs), intent(inout) :: tv !< Thermodynamical variable structure
621  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)), intent(inout) :: h !< Current 3D grid obtained after the
622  !! last time step [H ~> m or kg-2]
623  logical, optional, intent(in) :: debug !< If true, show the call tree
624  real, dimension(:,:), optional, pointer :: frac_shelf_h !< Fractional ice shelf coverage
625  ! Local variables
626  integer :: nk, i, j, k
627  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)+1) :: dzRegrid ! The change in grid interface positions
628  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)) :: h_new ! The new grid thicknesses
629  logical :: show_call_tree, use_ice_shelf
630 
631  show_call_tree = .false.
632  if (present(debug)) show_call_tree = debug
633  if (show_call_tree) call calltree_enter("ALE_build_grid(), MOM_ALE.F90")
634  use_ice_shelf = .false.
635  if (present(frac_shelf_h)) then
636  if (associated(frac_shelf_h)) use_ice_shelf = .true.
637  endif
638 
639  ! Build new grid. The new grid is stored in h_new. The old grid is h.
640  ! Both are needed for the subsequent remapping of variables.
641  if (use_ice_shelf) then
642  call regridding_main( remapcs, regridcs, g, gv, h, tv, h_new, dzregrid, frac_shelf_h )
643  else
644  call regridding_main( remapcs, regridcs, g, gv, h, tv, h_new, dzregrid )
645  endif
646 
647  ! Override old grid with new one. The new grid 'h_new' is built in
648  ! one of the 'build_...' routines above.
649 !$OMP parallel do default(none) shared(G,h,h_new)
650  do j = g%jsc,g%jec ; do i = g%isc,g%iec
651  if (g%mask2dT(i,j)>0.) h(i,j,:) = h_new(i,j,:)
652  enddo ; enddo
653 
654  if (show_call_tree) call calltree_leave("ALE_build_grid()")

◆ ale_end()

subroutine, public mom_ale::ale_end ( type(ale_cs), pointer  CS)

End of regridding (memory deallocation). This routine is typically called (from MOM_end in file MOM.F90) after the main time integration loop to deallocate the regridding stuff.

Parameters
csmodule control structure

Definition at line 306 of file MOM_ALE.F90.

307  type(ALE_CS), pointer :: CS !< module control structure
308 
309  ! Deallocate memory used for the regridding
310  call end_remapping( cs%remapCS )
311  call end_regridding( cs%regridCS )
312 
313  deallocate(cs)
314 

◆ ale_getcoordinate()

real function, dimension(cs%nk+1), public mom_ale::ale_getcoordinate ( type(ale_cs), pointer  CS)

Query the target coordinate interfaces positions.

Parameters
csmodule control structure

Definition at line 1206 of file MOM_ALE.F90.

1207  type(ALE_CS), pointer :: CS !< module control structure
1208 
1209  real, dimension(CS%nk+1) :: ALE_getCoordinate
1210  ale_getcoordinate(:) = getcoordinateinterfaces( cs%regridCS, undo_scaling=.true. )
1211 

◆ ale_getcoordinateunits()

character(len=20) function, public mom_ale::ale_getcoordinateunits ( type(ale_cs), pointer  CS)

Query the target coordinate units.

Parameters
csmodule control structure

Definition at line 1216 of file MOM_ALE.F90.

1217  type(ALE_CS), pointer :: CS !< module control structure
1218 
1219  character(len=20) :: ALE_getCoordinateUnits
1220 
1221  ale_getcoordinateunits = getcoordinateunits( cs%regridCS )
1222 

◆ ale_init()

subroutine, public mom_ale::ale_init ( type(param_file_type), intent(in)  param_file,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, intent(in)  max_depth,
type(ale_cs), pointer  CS 
)

This routine is typically called (from initialize_MOM in file MOM.F90) before the main time integration loop to initialize the regridding stuff. We read the MOM_input file to register the values of different regridding/remapping parameters.

Parameters
[in]param_fileParameter file
[in]gvOcean vertical grid structure
[in]usA dimensional unit scaling type
[in]max_depthThe maximum depth of the ocean [Z ~> m].
csModule control structure

Definition at line 139 of file MOM_ALE.F90.

140  type(param_file_type), intent(in) :: param_file !< Parameter file
141  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
142  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
143  real, intent(in) :: max_depth !< The maximum depth of the ocean [Z ~> m].
144  type(ALE_CS), pointer :: CS !< Module control structure
145 
146  ! Local variables
147  real, dimension(:), allocatable :: dz
148  character(len=40) :: mdl = "MOM_ALE" ! This module's name.
149  character(len=80) :: string ! Temporary strings
150  real :: filter_shallow_depth, filter_deep_depth
151  logical :: default_2018_answers
152  logical :: check_reconstruction
153  logical :: check_remapping
154  logical :: force_bounds_in_subcell
155  logical :: local_logical
156  logical :: remap_boundary_extrap
157 
158  if (associated(cs)) then
159  call mom_error(warning, "ALE_init called with an associated "// &
160  "control structure.")
161  return
162  endif
163  allocate(cs)
164 
165  cs%show_call_tree = calltree_showquery()
166  if (cs%show_call_tree) call calltree_enter("ALE_init(), MOM_ALE.F90")
167 
168  call get_param(param_file, mdl, "REMAP_UV_USING_OLD_ALG", cs%remap_uv_using_old_alg, &
169  "If true, uses the old remapping-via-a-delta-z method for "//&
170  "remapping u and v. If false, uses the new method that remaps "//&
171  "between grids described by an old and new thickness.", &
172  default=.false.)
173 
174  ! Initialize and configure regridding
175  call ale_initregridding(gv, us, max_depth, param_file, mdl, cs%regridCS)
176 
177  ! Initialize and configure remapping
178  call get_param(param_file, mdl, "REMAPPING_SCHEME", string, &
179  "This sets the reconstruction scheme used "//&
180  "for vertical remapping for all variables. "//&
181  "It can be one of the following schemes: "//&
182  trim(remappingschemesdoc), default=remappingdefaultscheme)
183  call get_param(param_file, mdl, "FATAL_CHECK_RECONSTRUCTIONS", check_reconstruction, &
184  "If true, cell-by-cell reconstructions are checked for "//&
185  "consistency and if non-monotonicity or an inconsistency is "//&
186  "detected then a FATAL error is issued.", default=.false.)
187  call get_param(param_file, mdl, "FATAL_CHECK_REMAPPING", check_remapping, &
188  "If true, the results of remapping are checked for "//&
189  "conservation and new extrema and if an inconsistency is "//&
190  "detected then a FATAL error is issued.", default=.false.)
191  call get_param(param_file, mdl, "REMAP_BOUND_INTERMEDIATE_VALUES", force_bounds_in_subcell, &
192  "If true, the values on the intermediate grid used for remapping "//&
193  "are forced to be bounded, which might not be the case due to "//&
194  "round off.", default=.false.)
195  call get_param(param_file, mdl, "REMAP_BOUNDARY_EXTRAP", remap_boundary_extrap, &
196  "If true, values at the interfaces of boundary cells are "//&
197  "extrapolated instead of piecewise constant", default=.false.)
198  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
199  "This sets the default value for the various _2018_ANSWERS parameters.", &
200  default=.false.)
201  call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", cs%answers_2018, &
202  "If true, use the order of arithmetic and expressions that recover the "//&
203  "answers from the end of 2018. Otherwise, use updated and more robust "//&
204  "forms of the same expressions.", default=default_2018_answers)
205  call initialize_remapping( cs%remapCS, string, &
206  boundary_extrapolation=remap_boundary_extrap, &
207  check_reconstruction=check_reconstruction, &
208  check_remapping=check_remapping, &
209  force_bounds_in_subcell=force_bounds_in_subcell, &
210  answers_2018=cs%answers_2018)
211 
212  call get_param(param_file, mdl, "REMAP_AFTER_INITIALIZATION", cs%remap_after_initialization, &
213  "If true, applies regridding and remapping immediately after "//&
214  "initialization so that the state is ALE consistent. This is a "//&
215  "legacy step and should not be needed if the initialization is "//&
216  "consistent with the coordinate mode.", default=.true.)
217 
218  call get_param(param_file, mdl, "REGRID_TIME_SCALE", cs%regrid_time_scale, &
219  "The time-scale used in blending between the current (old) grid "//&
220  "and the target (new) grid. A short time-scale favors the target "//&
221  "grid (0. or anything less than DT_THERM) has no memory of the old "//&
222  "grid. A very long time-scale makes the model more Lagrangian.", &
223  units="s", default=0., scale=us%s_to_T)
224  call get_param(param_file, mdl, "REGRID_FILTER_SHALLOW_DEPTH", filter_shallow_depth, &
225  "The depth above which no time-filtering is applied. Above this depth "//&
226  "final grid exactly matches the target (new) grid.", &
227  units="m", default=0., scale=gv%m_to_H)
228  call get_param(param_file, mdl, "REGRID_FILTER_DEEP_DEPTH", filter_deep_depth, &
229  "The depth below which full time-filtering is applied with time-scale "//&
230  "REGRID_TIME_SCALE. Between depths REGRID_FILTER_SHALLOW_DEPTH and "//&
231  "REGRID_FILTER_SHALLOW_DEPTH the filter weights adopt a cubic profile.", &
232  units="m", default=0., scale=gv%m_to_H)
233  call set_regrid_params(cs%regridCS, depth_of_time_filter_shallow=filter_shallow_depth, &
234  depth_of_time_filter_deep=filter_deep_depth)
235  call get_param(param_file, mdl, "REGRID_USE_OLD_DIRECTION", local_logical, &
236  "If true, the regridding ntegrates upwards from the bottom for "//&
237  "interface positions, much as the main model does. If false "//&
238  "regridding integrates downward, consistant with the remapping "//&
239  "code.", default=.true., do_not_log=.true.)
240  call set_regrid_params(cs%regridCS, integrate_downward_for_e=.not.local_logical)
241 
242  ! Keep a record of values for subsequent queries
243  cs%nk = gv%ke
244 
245  if (cs%show_call_tree) call calltree_leave("ALE_init()")

◆ ale_initregridding()

subroutine, public mom_ale::ale_initregridding ( type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, intent(in)  max_depth,
type(param_file_type), intent(in)  param_file,
character(len=*), intent(in)  mdl,
type(regridding_cs), intent(out)  regridCS 
)

Initializes regridding for the main ALE algorithm.

Parameters
[in]gvOcean vertical grid structure
[in]usA dimensional unit scaling type
[in]max_depthThe maximum depth of the ocean [Z ~> m].
[in]param_fileparameter file
[in]mdlName of calling module
[out]regridcsRegridding parameters and work arrays

Definition at line 1185 of file MOM_ALE.F90.

1186  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1187  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1188  real, intent(in) :: max_depth !< The maximum depth of the ocean [Z ~> m].
1189  type(param_file_type), intent(in) :: param_file !< parameter file
1190  character(len=*), intent(in) :: mdl !< Name of calling module
1191  type(regridding_CS), intent(out) :: regridCS !< Regridding parameters and work arrays
1192  ! Local variables
1193  character(len=30) :: coord_mode
1194 
1195  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_MODE", coord_mode, &
1196  "Coordinate mode for vertical regridding. "//&
1197  "Choose among the following possibilities: "//&
1198  trim(regriddingcoordinatemodedoc), &
1199  default=default_coordinate_mode, fail_if_missing=.true.)
1200 
1201  call initialize_regridding(regridcs, gv, us, max_depth, param_file, mdl, coord_mode, '', '')
1202 

◆ ale_initthicknesstocoord()

subroutine, public mom_ale::ale_initthicknesstocoord ( type(ale_cs), intent(inout)  CS,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(out)  h 
)

Set h to coordinate values for fixed coordinate systems.

Parameters
[in,out]csmodule control structure
[in]gmodule grid structure
[in]gvOcean vertical grid structure
[out]hlayer thickness [H ~> m or kg m-2]

Definition at line 1305 of file MOM_ALE.F90.

1306  type(ALE_CS), intent(inout) :: CS !< module control structure
1307  type(ocean_grid_type), intent(in) :: G !< module grid structure
1308  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1309  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(out) :: h !< layer thickness [H ~> m or kg m-2]
1310 
1311  ! Local variables
1312  integer :: i, j, k
1313 
1314  do j = g%jsd,g%jed ; do i = g%isd,g%ied
1315  h(i,j,:) = gv%Z_to_H * getstaticthickness( cs%regridCS, 0., g%bathyT(i,j) )
1316  enddo ; enddo
1317 

◆ ale_main()

subroutine, public mom_ale::ale_main ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  h,
real, dimension(szib_(g),szj_(g),szk_(gv)), intent(inout)  u,
real, dimension(szi_(g),szjb_(g),szk_(gv)), intent(inout)  v,
type(thermo_var_ptrs), intent(inout)  tv,
type(tracer_registry_type), pointer  Reg,
type(ale_cs), pointer  CS,
type(ocean_obc_type), pointer  OBC,
real, intent(in), optional  dt,
real, dimension(:,:), optional, pointer  frac_shelf_h 
)

Takes care of (1) building a new grid and (2) remapping all variables between the old grid and the new grid. The creation of the new grid can be based on z coordinates, target interface densities, sigma coordinates or any arbitrary coordinate system.

Parameters
[in]gOcean grid informations
[in]gvOcean vertical grid structure
[in]usA dimensional unit scaling type
[in,out]hCurrent 3D grid obtained after the last time step [H ~> m or kg m-2]
[in,out]uZonal velocity field [L T-1 ~> m s-1]
[in,out]vMeridional velocity field [L T-1 ~> m s-1]
[in,out]tvThermodynamic variable structure
regTracer registry structure
csRegridding parameters and options
obcOpen boundary structure
[in]dtTime step between calls to ALE_main [T ~> s]
frac_shelf_hFractional ice shelf coverage

Definition at line 321 of file MOM_ALE.F90.

322  type(ocean_grid_type), intent(in) :: G !< Ocean grid informations
323  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
324  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
325  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(inout) :: h !< Current 3D grid obtained after the
326  !! last time step [H ~> m or kg m-2]
327  real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(inout) :: u !< Zonal velocity field [L T-1 ~> m s-1]
328  real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(inout) :: v !< Meridional velocity field [L T-1 ~> m s-1]
329  type(thermo_var_ptrs), intent(inout) :: tv !< Thermodynamic variable structure
330  type(tracer_registry_type), pointer :: Reg !< Tracer registry structure
331  type(ALE_CS), pointer :: CS !< Regridding parameters and options
332  type(ocean_OBC_type), pointer :: OBC !< Open boundary structure
333  real, optional, intent(in) :: dt !< Time step between calls to ALE_main [T ~> s]
334  real, dimension(:,:), optional, pointer :: frac_shelf_h !< Fractional ice shelf coverage
335  ! Local variables
336  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: dzRegrid ! The change in grid interface positions
337  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: eta_preale
338  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h_new ! New 3D grid obtained after last time step [H ~> m or kg-2]
339  integer :: nk, i, j, k, isc, iec, jsc, jec
340  logical :: ice_shelf
341 
342  nk = gv%ke; isc = g%isc; iec = g%iec; jsc = g%jsc; jec = g%jec
343 
344  ice_shelf = .false.
345  if (present(frac_shelf_h)) then
346  if (associated(frac_shelf_h)) ice_shelf = .true.
347  endif
348 
349  if (cs%show_call_tree) call calltree_enter("ALE_main(), MOM_ALE.F90")
350 
351  ! These diagnostics of the state before ALE is applied are mostly used for debugging.
352  if (cs%id_u_preale > 0) call post_data(cs%id_u_preale, u, cs%diag)
353  if (cs%id_v_preale > 0) call post_data(cs%id_v_preale, v, cs%diag)
354  if (cs%id_h_preale > 0) call post_data(cs%id_h_preale, h, cs%diag)
355  if (cs%id_T_preale > 0) call post_data(cs%id_T_preale, tv%T, cs%diag)
356  if (cs%id_S_preale > 0) call post_data(cs%id_S_preale, tv%S, cs%diag)
357  if (cs%id_e_preale > 0) then
358  call find_eta(h, tv, g, gv, us, eta_preale)
359  call post_data(cs%id_e_preale, eta_preale, cs%diag)
360  endif
361 
362  if (present(dt)) then
363  call ale_update_regrid_weights( dt, cs )
364  endif
365  dzregrid(:,:,:) = 0.0
366 
367  ! Build new grid. The new grid is stored in h_new. The old grid is h.
368  ! Both are needed for the subsequent remapping of variables.
369  if (ice_shelf) then
370  call regridding_main( cs%remapCS, cs%regridCS, g, gv, h, tv, h_new, dzregrid, frac_shelf_h)
371  else
372  call regridding_main( cs%remapCS, cs%regridCS, g, gv, h, tv, h_new, dzregrid)
373  endif
374 
375  call check_grid( g, gv, h, 0. )
376 
377  if (cs%show_call_tree) call calltree_waypoint("new grid generated (ALE_main)")
378 
379  ! The presence of dt is used for expediency to distinguish whether ALE_main is being called during init
380  ! or in the main loop. Tendency diagnostics in remap_all_state_vars also rely on this logic.
381  if (present(dt)) then
382  call diag_update_remap_grids(cs%diag)
383  endif
384  ! Remap all variables from old grid h onto new grid h_new
385  call remap_all_state_vars( cs%remapCS, cs, g, gv, h, h_new, reg, obc, -dzregrid, &
386  u, v, cs%show_call_tree, dt )
387 
388  if (cs%show_call_tree) call calltree_waypoint("state remapped (ALE_main)")
389 
390  ! Override old grid with new one. The new grid 'h_new' is built in
391  ! one of the 'build_...' routines above.
392  !$OMP parallel do default(shared)
393  do k = 1,nk ; do j = jsc-1,jec+1 ; do i = isc-1,iec+1
394  h(i,j,k) = h_new(i,j,k)
395  enddo ; enddo ; enddo
396 
397  if (cs%show_call_tree) call calltree_leave("ALE_main()")
398 
399  if (cs%id_dzRegrid>0 .and. present(dt)) call post_data(cs%id_dzRegrid, dzregrid, cs%diag)
400 
401 

◆ ale_main_offline()

subroutine, public mom_ale::ale_main_offline ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  h,
type(thermo_var_ptrs), intent(inout)  tv,
type(tracer_registry_type), pointer  Reg,
type(ale_cs), pointer  CS,
type(ocean_obc_type), pointer  OBC,
real, intent(in), optional  dt 
)

Takes care of (1) building a new grid and (2) remapping all variables between the old grid and the new grid. The creation of the new grid can be based on z coordinates, target interface densities, sigma coordinates or any arbitrary coordinate system.

Parameters
[in]gOcean grid informations
[in]gvOcean vertical grid structure
[in,out]hCurrent 3D grid obtained after the last time step [H ~> m or kg-2]
[in,out]tvThermodynamic variable structure
regTracer registry structure
csRegridding parameters and options
obcOpen boundary structure
[in]dtTime step between calls to ALE_main [T ~> s]

Definition at line 408 of file MOM_ALE.F90.

409  type(ocean_grid_type), intent(in) :: G !< Ocean grid informations
410  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
411  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(inout) :: h !< Current 3D grid obtained after the
412  !! last time step [H ~> m or kg-2]
413  type(thermo_var_ptrs), intent(inout) :: tv !< Thermodynamic variable structure
414  type(tracer_registry_type), pointer :: Reg !< Tracer registry structure
415  type(ALE_CS), pointer :: CS !< Regridding parameters and options
416  type(ocean_OBC_type), pointer :: OBC !< Open boundary structure
417  real, optional, intent(in) :: dt !< Time step between calls to ALE_main [T ~> s]
418  ! Local variables
419  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)+1) :: dzRegrid ! The change in grid interface positions
420  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h_new ! New 3D grid obtained after last time step [H ~> m or kg-2]
421  integer :: nk, i, j, k, isc, iec, jsc, jec
422 
423  nk = gv%ke; isc = g%isc; iec = g%iec; jsc = g%jsc; jec = g%jec
424 
425  if (cs%show_call_tree) call calltree_enter("ALE_main_offline(), MOM_ALE.F90")
426 
427  if (present(dt)) then
428  call ale_update_regrid_weights( dt, cs )
429  endif
430  dzregrid(:,:,:) = 0.0
431 
432  ! Build new grid. The new grid is stored in h_new. The old grid is h.
433  ! Both are needed for the subsequent remapping of variables.
434  call regridding_main( cs%remapCS, cs%regridCS, g, gv, h, tv, h_new, dzregrid )
435 
436  call check_grid( g, gv, h, 0. )
437 
438  if (cs%show_call_tree) call calltree_waypoint("new grid generated (ALE_main)")
439 
440  ! Remap all variables from old grid h onto new grid h_new
441 
442  call remap_all_state_vars(cs%remapCS, cs, g, gv, h, h_new, reg, obc, &
443  debug=cs%show_call_tree, dt=dt )
444 
445  if (cs%show_call_tree) call calltree_waypoint("state remapped (ALE_main)")
446 
447  ! Override old grid with new one. The new grid 'h_new' is built in
448  ! one of the 'build_...' routines above.
449  !$OMP parallel do default(shared)
450  do k = 1,nk ; do j = jsc-1,jec+1 ; do i = isc-1,iec+1
451  h(i,j,k) = h_new(i,j,k)
452  enddo ; enddo ; enddo
453 
454  if (cs%show_call_tree) call calltree_leave("ALE_main()")
455  if (cs%id_dzRegrid>0 .and. present(dt)) call post_data(cs%id_dzRegrid, dzregrid, cs%diag)
456 

◆ ale_offline_inputs()

subroutine, public mom_ale::ale_offline_inputs ( type(ale_cs), pointer  CS,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  h,
type(thermo_var_ptrs), intent(inout)  tv,
type(tracer_registry_type), pointer  Reg,
real, dimension(szib_(g),szj_(g),szk_(gv)), intent(inout)  uhtr,
real, dimension(szi_(g),szjb_(g),szk_(gv)), intent(inout)  vhtr,
real, dimension(szi_(g),szj_(g),szk_(gv)+1), intent(inout)  Kd,
logical, intent(in)  debug,
type(ocean_obc_type), pointer  OBC 
)

Regrid/remap stored fields used for offline tracer integrations. These input fields are assumed to have the same layer thicknesses at the end of the last offline interval (which should be a Zstar grid). This routine builds a grid on the runtime specified vertical coordinate.

Parameters
csRegridding parameters and options
[in]gOcean grid informations
[in]gvOcean vertical grid structure
[in,out]hLayer thicknesses
[in,out]tvThermodynamic variable structure
regTracer registry structure
[in,out]uhtrZonal mass fluxes
[in,out]vhtrMeridional mass fluxes
[in,out]kdInput diffusivites
[in]debugIf true, then turn checksums
obcOpen boundary structure

Definition at line 462 of file MOM_ALE.F90.

463  type(ALE_CS), pointer :: CS !< Regridding parameters and options
464  type(ocean_grid_type), intent(in ) :: G !< Ocean grid informations
465  type(verticalGrid_type), intent(in ) :: GV !< Ocean vertical grid structure
466  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(inout) :: h !< Layer thicknesses
467  type(thermo_var_ptrs), intent(inout) :: tv !< Thermodynamic variable structure
468  type(tracer_registry_type), pointer :: Reg !< Tracer registry structure
469  real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(inout) :: uhtr !< Zonal mass fluxes
470  real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(inout) :: vhtr !< Meridional mass fluxes
471  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), intent(inout) :: Kd !< Input diffusivites
472  logical, intent(in ) :: debug !< If true, then turn checksums
473  type(ocean_OBC_type), pointer :: OBC !< Open boundary structure
474  ! Local variables
475  integer :: nk, i, j, k, isc, iec, jsc, jec
476  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)) :: h_new ! Layer thicknesses after regridding
477  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)+1) :: dzRegrid ! The change in grid interface positions
478  real, dimension(SZK_(GV)) :: h_src
479  real, dimension(SZK_(GV)) :: h_dest, uh_dest
480  real, dimension(SZK_(GV)) :: temp_vec
481 
482  nk = gv%ke; isc = g%isc; iec = g%iec; jsc = g%jsc; jec = g%jec
483  dzregrid(:,:,:) = 0.0
484  h_new(:,:,:) = 0.0
485 
486  if (debug) call mom_tracer_chkinv("Before ALE_offline_inputs", g, h, reg%Tr, reg%ntr)
487 
488  ! Build new grid from the Zstar state onto the requested vertical coordinate. The new grid is stored
489  ! in h_new. The old grid is h. Both are needed for the subsequent remapping of variables. Convective
490  ! adjustment right now is not used because it is unclear what to do with vanished layers
491  call regridding_main( cs%remapCS, cs%regridCS, g, gv, h, tv, h_new, dzregrid, conv_adjust = .false. )
492  call check_grid( g, gv, h_new, 0. )
493  if (cs%show_call_tree) call calltree_waypoint("new grid generated (ALE_offline_inputs)")
494 
495  ! Remap all variables from old grid h onto new grid h_new
496  call remap_all_state_vars( cs%remapCS, cs, g, gv, h, h_new, reg, obc, debug=cs%show_call_tree )
497  if (cs%show_call_tree) call calltree_waypoint("state remapped (ALE_inputs)")
498 
499  ! Reintegrate mass transports from Zstar to the offline vertical coordinate
500  do j=jsc,jec ; do i=g%iscB,g%iecB
501  if (g%mask2dCu(i,j)>0.) then
502  h_src(:) = 0.5 * (h(i,j,:) + h(i+1,j,:))
503  h_dest(:) = 0.5 * (h_new(i,j,:) + h_new(i+1,j,:))
504  call reintegrate_column(nk, h_src, uhtr(i,j,:), nk, h_dest, 0., temp_vec)
505  uhtr(i,j,:) = temp_vec
506  endif
507  enddo ; enddo
508  do j=g%jscB,g%jecB ; do i=isc,iec
509  if (g%mask2dCv(i,j)>0.) then
510  h_src(:) = 0.5 * (h(i,j,:) + h(i,j+1,:))
511  h_dest(:) = 0.5 * (h_new(i,j,:) + h_new(i,j+1,:))
512  call reintegrate_column(nk, h_src, vhtr(i,j,:), nk, h_dest, 0., temp_vec)
513  vhtr(i,j,:) = temp_vec
514  endif
515  enddo ; enddo
516 
517  do j = jsc,jec ; do i=isc,iec
518  if (g%mask2dT(i,j)>0.) then
519  if (check_column_integrals(nk, h_src, nk, h_dest)) then
520  call mom_error(fatal, "ALE_offline_inputs: Kd interpolation columns do not match")
521  endif
522  call interpolate_column(nk, h(i,j,:), kd(i,j,:), nk, h_new(i,j,:), 0., kd(i,j,:))
523  endif
524  enddo ; enddo
525 
526  call ale_remap_scalar(cs%remapCS, g, gv, nk, h, tv%T, h_new, tv%T, answers_2018=cs%answers_2018)
527  call ale_remap_scalar(cs%remapCS, g, gv, nk, h, tv%S, h_new, tv%S, answers_2018=cs%answers_2018)
528 
529  if (debug) call mom_tracer_chkinv("After ALE_offline_inputs", g, h_new, reg%Tr, reg%ntr)
530 
531  ! Copy over the new layer thicknesses
532  do k = 1,nk ; do j = jsc-1,jec+1 ; do i = isc-1,iec+1
533  h(i,j,k) = h_new(i,j,k)
534  enddo ; enddo ; enddo
535 
536  if (cs%show_call_tree) call calltree_leave("ALE_offline_inputs()")

◆ ale_offline_tracer_final()

subroutine, public mom_ale::ale_offline_tracer_final ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(inout)  h,
type(thermo_var_ptrs), intent(inout)  tv,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(inout)  h_target,
type(tracer_registry_type), pointer  Reg,
type(ale_cs), pointer  CS,
type(ocean_obc_type), pointer  OBC 
)

Remaps all tracers from h onto h_target. This is intended to be called when tracers are done offline. In the case where transports don't quite conserve, we still want to make sure that layer thicknesses offline do not drift too far away from the online model.

Parameters
[in]gOcean grid informations
[in]gvOcean vertical grid structure
[in,out]hCurrent 3D grid obtained after the last time step [H ~> m or kg-2]
[in,out]tvThermodynamic variable structure
[in,out]h_targetCurrent 3D grid obtained after last time step [H ~> m or kg-2]
regTracer registry structure
csRegridding parameters and options
obcOpen boundary structure

Definition at line 543 of file MOM_ALE.F90.

544  type(ocean_grid_type), intent(in) :: G !< Ocean grid informations
545  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
546  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(inout) :: h !< Current 3D grid obtained after the
547  !! last time step [H ~> m or kg-2]
548  type(thermo_var_ptrs), intent(inout) :: tv !< Thermodynamic variable structure
549  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(inout) :: h_target !< Current 3D grid obtained after
550  !! last time step [H ~> m or kg-2]
551  type(tracer_registry_type), pointer :: Reg !< Tracer registry structure
552  type(ALE_CS), pointer :: CS !< Regridding parameters and options
553  type(ocean_OBC_type), pointer :: OBC !< Open boundary structure
554  ! Local variables
555 
556  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)+1) :: dzRegrid !< The change in grid interface positions
557  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)) :: h_new !< Regridded target thicknesses
558  integer :: nk, i, j, k, isc, iec, jsc, jec
559 
560  nk = gv%ke; isc = g%isc; iec = g%iec; jsc = g%jsc; jec = g%jec
561 
562  if (cs%show_call_tree) call calltree_enter("ALE_offline_tracer_final(), MOM_ALE.F90")
563  ! Need to make sure that h_target is consistent with the current offline ALE confiuration
564  call regridding_main( cs%remapCS, cs%regridCS, g, gv, h_target, tv, h_new, dzregrid )
565  call check_grid( g, gv, h_target, 0. )
566 
567 
568  if (cs%show_call_tree) call calltree_waypoint("Source and target grids checked (ALE_offline_tracer_final)")
569 
570  ! Remap all variables from old grid h onto new grid h_new
571 
572  call remap_all_state_vars( cs%remapCS, cs, g, gv, h, h_new, reg, obc, debug=cs%show_call_tree )
573 
574  if (cs%show_call_tree) call calltree_waypoint("state remapped (ALE_offline_tracer_final)")
575 
576  ! Override old grid with new one. The new grid 'h_new' is built in
577  ! one of the 'build_...' routines above.
578  !$OMP parallel do default(shared)
579  do k = 1,nk
580  do j = jsc-1,jec+1 ; do i = isc-1,iec+1
581  h(i,j,k) = h_new(i,j,k)
582  enddo ; enddo
583  enddo
584  if (cs%show_call_tree) call calltree_leave("ALE_offline_tracer_final()")

◆ ale_plm_edge_values()

subroutine, public mom_ale::ale_plm_edge_values ( type(ale_cs), intent(in)  CS,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  h,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  Q,
logical, intent(in)  bdry_extrap,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  Q_t,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  Q_b 
)

Calculate edge values (top and bottom of layer) 3d scalar array. Boundary reconstructions are PCM unless bdry_extrap is true.

Parameters
[in]csmodule control structure
[in]gocean grid structure
[in]gvOcean vertical grid structure
[in]hlayer thickness [H ~> m or kg m-2]
[in]q3d scalar array
[in]bdry_extrapIf true, use high-order boundary extrapolation within boundary cells
[in,out]q_tScalar at the top edge of each layer
[in,out]q_bScalar at the bottom edge of each layer

Definition at line 1038 of file MOM_ALE.F90.

1039  type(ALE_CS), intent(in) :: CS !< module control structure
1040  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
1041  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1042  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1043  intent(in) :: h !< layer thickness [H ~> m or kg m-2]
1044  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1045  intent(in) :: Q !< 3d scalar array
1046  logical, intent(in) :: bdry_extrap !< If true, use high-order boundary
1047  !! extrapolation within boundary cells
1048  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1049  intent(inout) :: Q_t !< Scalar at the top edge of each layer
1050  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1051  intent(inout) :: Q_b !< Scalar at the bottom edge of each layer
1052  ! Local variables
1053  integer :: i, j, k
1054  real :: slp(GV%ke)
1055  real :: mslp
1056  real :: h_neglect
1057 
1058  if (.not.cs%answers_2018) then
1059  h_neglect = gv%H_subroundoff
1060  elseif (gv%Boussinesq) then
1061  h_neglect = gv%m_to_H*1.0e-30
1062  else
1063  h_neglect = gv%kg_m2_to_H*1.0e-30
1064  endif
1065 
1066  !$OMP parallel do default(shared) private(slp,mslp)
1067  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
1068  slp(1) = 0.
1069  do k = 2, gv%ke-1
1070  slp(k) = plm_slope_wa(h(i,j,k-1), h(i,j,k), h(i,j,k+1), h_neglect, q(i,j,k-1), q(i,j,k), q(i,j,k+1))
1071  enddo
1072  slp(gv%ke) = 0.
1073 
1074  do k = 2, gv%ke-1
1075  mslp = plm_monotonized_slope(q(i,j,k-1), q(i,j,k), q(i,j,k+1), slp(k-1), slp(k), slp(k+1))
1076  q_t(i,j,k) = q(i,j,k) - 0.5 * mslp
1077  q_b(i,j,k) = q(i,j,k) + 0.5 * mslp
1078  enddo
1079  if (bdry_extrap) then
1080  mslp = - plm_extrapolate_slope(h(i,j,2), h(i,j,1), h_neglect, q(i,j,2), q(i,j,1))
1081  q_t(i,j,1) = q(i,j,1) - 0.5 * mslp
1082  q_b(i,j,1) = q(i,j,1) + 0.5 * mslp
1083  mslp = plm_extrapolate_slope(h(i,j,gv%ke-1), h(i,j,gv%ke), h_neglect, q(i,j,gv%ke-1), q(i,j,gv%ke))
1084  q_t(i,j,gv%ke) = q(i,j,gv%ke) - 0.5 * mslp
1085  q_b(i,j,gv%ke) = q(i,j,gv%ke) + 0.5 * mslp
1086  else
1087  q_t(i,j,1) = q(i,j,1)
1088  q_b(i,j,1) = q(i,j,1)
1089  q_t(i,j,gv%ke) = q(i,j,gv%ke)
1090  q_b(i,j,gv%ke) = q(i,j,gv%ke)
1091  endif
1092 
1093  enddo ; enddo
1094 

◆ ale_register_diags()

subroutine, public mom_ale::ale_register_diags ( type(time_type), intent(in), target  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(diag_ctrl), intent(in), target  diag,
type(ale_cs), pointer  CS 
)

Initialize diagnostics for the ALE module.

Parameters
[in]timeTime structure
[in]gGrid structure
[in]usA dimensional unit scaling type
[in]gvOcean vertical grid structure
[in]diagDiagnostics control structure
csModule control structure

Definition at line 249 of file MOM_ALE.F90.

250  type(time_type),target, intent(in) :: Time !< Time structure
251  type(ocean_grid_type), intent(in) :: G !< Grid structure
252  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
253  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
254  type(diag_ctrl), target, intent(in) :: diag !< Diagnostics control structure
255  type(ALE_CS), pointer :: CS !< Module control structure
256 
257  cs%diag => diag
258 
259  ! These diagnostics of the state variables before ALE are useful for
260  ! debugging the ALE code.
261  cs%id_u_preale = register_diag_field('ocean_model', 'u_preale', diag%axesCuL, time, &
262  'Zonal velocity before remapping', 'm s-1', conversion=us%L_T_to_m_s)
263  cs%id_v_preale = register_diag_field('ocean_model', 'v_preale', diag%axesCvL, time, &
264  'Meridional velocity before remapping', 'm s-1', conversion=us%L_T_to_m_s)
265  cs%id_h_preale = register_diag_field('ocean_model', 'h_preale', diag%axesTL, time, &
266  'Layer Thickness before remapping', get_thickness_units(gv), &
267  conversion=gv%H_to_MKS, v_extensive=.true.)
268  cs%id_T_preale = register_diag_field('ocean_model', 'T_preale', diag%axesTL, time, &
269  'Temperature before remapping', 'degC')
270  cs%id_S_preale = register_diag_field('ocean_model', 'S_preale', diag%axesTL, time, &
271  'Salinity before remapping', 'PSU')
272  cs%id_e_preale = register_diag_field('ocean_model', 'e_preale', diag%axesTi, time, &
273  'Interface Heights before remapping', 'm', conversion=us%Z_to_m)
274 
275  cs%id_dzRegrid = register_diag_field('ocean_model','dzRegrid',diag%axesTi,time, &
276  'Change in interface height due to ALE regridding', 'm', &
277  conversion=gv%H_to_m)
278  cs%id_vert_remap_h = register_diag_field('ocean_model', 'vert_remap_h', &
279  diag%axestl, time, 'layer thicknesses after ALE regridding and remapping', 'm', &
280  conversion=gv%H_to_m, v_extensive=.true.)
281  cs%id_vert_remap_h_tendency = register_diag_field('ocean_model','vert_remap_h_tendency',diag%axestl,time, &
282  'Layer thicknesses tendency due to ALE regridding and remapping', 'm', &
283  conversion=gv%H_to_m*us%s_to_T, v_extensive = .true.)
284 

◆ ale_regrid_accelerated()

subroutine, public mom_ale::ale_regrid_accelerated ( type(ale_cs), pointer  CS,
type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(inout)  h,
type(thermo_var_ptrs), intent(inout)  tv,
integer, intent(in)  n,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, gv %ke), intent(inout)  u,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, gv %ke), intent(inout)  v,
type(ocean_obc_type), pointer  OBC,
type(tracer_registry_type), optional, pointer  Reg,
real, intent(in), optional  dt,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke+1), intent(inout), optional  dzRegrid,
logical, intent(in), optional  initial 
)

For a state-based coordinate, accelerate the process of regridding by repeatedly applying the grid calculation algorithm.

Parameters
csALE control structure
[in,out]gOcean grid
[in]gvVertical grid
[in,out]hOriginal thicknesses [H ~> m or kg-2]
[in,out]tvThermo vars (T/S/EOS)
[in]nNumber of times to regrid
[in,out]uZonal velocity [L T-1 ~> m s-1]
[in,out]vMeridional velocity [L T-1 ~> m s-1]
obcOpen boundary structure
regTracer registry to remap onto new grid
[in]dtModel timestep to provide a timescale for regridding [T ~> s]
[in,out]dzregridFinal change in interface positions
[in]initialWhether we're being called from an initialization routine (and expect diagnostics to work)

Definition at line 659 of file MOM_ALE.F90.

660  type(ALE_CS), pointer :: CS !< ALE control structure
661  type(ocean_grid_type), intent(inout) :: G !< Ocean grid
662  type(verticalGrid_type), intent(in) :: GV !< Vertical grid
663  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
664  intent(inout) :: h !< Original thicknesses [H ~> m or kg-2]
665  type(thermo_var_ptrs), intent(inout) :: tv !< Thermo vars (T/S/EOS)
666  integer, intent(in) :: n !< Number of times to regrid
667  real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
668  intent(inout) :: u !< Zonal velocity [L T-1 ~> m s-1]
669  real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
670  intent(inout) :: v !< Meridional velocity [L T-1 ~> m s-1]
671  type(ocean_OBC_type), pointer :: OBC !< Open boundary structure
672  type(tracer_registry_type), &
673  optional, pointer :: Reg !< Tracer registry to remap onto new grid
674  real, optional, intent(in) :: dt !< Model timestep to provide a timescale for regridding [T ~> s]
675  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), &
676  optional, intent(inout) :: dzRegrid !< Final change in interface positions
677  logical, optional, intent(in) :: initial !< Whether we're being called from an initialization
678  !! routine (and expect diagnostics to work)
679 
680  ! Local variables
681  integer :: i, j, k, nz
682  type(thermo_var_ptrs) :: tv_local ! local/intermediate temp/salt
683  type(group_pass_type) :: pass_T_S_h ! group pass if the coordinate has a stencil
684  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h_loc, h_orig ! A working copy of layer thicknesses
685  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), target :: T, S ! local temporary state
686  ! we have to keep track of the total dzInterface if for some reason
687  ! we're using the old remapping algorithm for u/v
688  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: dzInterface, dzIntTotal
689 
690  nz = gv%ke
691 
692  ! initial total interface displacement due to successive regridding
693  dzinttotal(:,:,:) = 0.
694 
695  call create_group_pass(pass_t_s_h, t, g%domain)
696  call create_group_pass(pass_t_s_h, s, g%domain)
697  call create_group_pass(pass_t_s_h, h_loc, g%domain)
698 
699  ! copy original temp/salt and set our local tv_pointers to them
700  tv_local = tv
701  t(:,:,:) = tv%T(:,:,:)
702  s(:,:,:) = tv%S(:,:,:)
703  tv_local%T => t
704  tv_local%S => s
705 
706  ! get local copy of thickness and save original state for remapping
707  h_loc(:,:,:) = h(:,:,:)
708  h_orig(:,:,:) = h(:,:,:)
709 
710  ! Apply timescale to regridding (for e.g. filtered_grid_motion)
711  if (present(dt)) &
712  call ale_update_regrid_weights(dt, cs)
713 
714  do k = 1, n
715  call do_group_pass(pass_t_s_h, g%domain)
716 
717  ! generate new grid
718  call regridding_main(cs%remapCS, cs%regridCS, g, gv, h_loc, tv_local, h, dzinterface)
719  dzinttotal(:,:,:) = dzinttotal(:,:,:) + dzinterface(:,:,:)
720 
721  ! remap from original grid onto new grid
722  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
723  call remapping_core_h(cs%remapCS, nz, h_orig(i,j,:), tv%S(i,j,:), nz, h(i,j,:), tv_local%S(i,j,:))
724  call remapping_core_h(cs%remapCS, nz, h_orig(i,j,:), tv%T(i,j,:), nz, h(i,j,:), tv_local%T(i,j,:))
725  enddo ; enddo
726 
727  ! starting grid for next iteration
728  h_loc(:,:,:) = h(:,:,:)
729  enddo
730 
731  ! remap all state variables (including those that weren't needed for regridding)
732  call remap_all_state_vars(cs%remapCS, cs, g, gv, h_orig, h, reg, obc, dzinttotal, u, v)
733 
734  ! save total dzregrid for diags if needed?
735  if (present(dzregrid)) dzregrid(:,:,:) = dzinttotal(:,:,:)

◆ ale_remap_init_conds()

logical function, public mom_ale::ale_remap_init_conds ( type(ale_cs), pointer  CS)

Returns true if initial conditions should be regridded and remapped.

Parameters
csmodule control structure

Definition at line 1227 of file MOM_ALE.F90.

1228  type(ALE_CS), pointer :: CS !< module control structure
1229 
1230  ale_remap_init_conds = .false.
1231  if (associated(cs)) ale_remap_init_conds = cs%remap_after_initialization

◆ ale_remap_scalar()

subroutine, public mom_ale::ale_remap_scalar ( type(remapping_cs), intent(in)  CS,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
integer, intent(in)  nk_src,
real, dimension(szi_(g),szj_(g),nk_src), intent(in)  h_src,
real, dimension(szi_(g),szj_(g),nk_src), intent(in)  s_src,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  h_dst,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  s_dst,
logical, intent(in), optional  all_cells,
logical, intent(in), optional  old_remap,
logical, intent(in), optional  answers_2018 
)

Remaps a single scalar between grids described by thicknesses h_src and h_dst. h_dst must be dimensioned as a model array with GVke layers while h_src can have an arbitrary number of layers specified by nk_src.

Parameters
[in]csRemapping control structure
[in]gOcean grid structure
[in]gvOcean vertical grid structure
[in]nk_srcNumber of levels on source grid
[in]h_srcLevel thickness of source grid [H ~> m or kg-2]
[in]s_srcScalar on source grid
[in]h_dstLevel thickness of destination grid [H ~> m or kg-2]
[in,out]s_dstScalar on destination grid
[in]all_cellsIf false, only reconstruct for non-vanished cells. Use all vanished layers otherwise (default).
[in]old_remapIf true, use the old "remapping_core_w" method, otherwise use "remapping_core_h".
[in]answers_2018If true, use the order of arithmetic and expressions that recover the answers for remapping from the end of 2018. Otherwise, use more robust forms of the same expressions.

Definition at line 944 of file MOM_ALE.F90.

945  type(remapping_CS), intent(in) :: CS !< Remapping control structure
946  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
947  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
948  integer, intent(in) :: nk_src !< Number of levels on source grid
949  real, dimension(SZI_(G),SZJ_(G),nk_src), intent(in) :: h_src !< Level thickness of source grid
950  !! [H ~> m or kg-2]
951  real, dimension(SZI_(G),SZJ_(G),nk_src), intent(in) :: s_src !< Scalar on source grid
952  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)),intent(in) :: h_dst !< Level thickness of destination grid
953  !! [H ~> m or kg-2]
954  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)),intent(inout) :: s_dst !< Scalar on destination grid
955  logical, optional, intent(in) :: all_cells !< If false, only reconstruct for
956  !! non-vanished cells. Use all vanished
957  !! layers otherwise (default).
958  logical, optional, intent(in) :: old_remap !< If true, use the old "remapping_core_w"
959  !! method, otherwise use "remapping_core_h".
960  logical, optional, intent(in) :: answers_2018 !< If true, use the order of arithmetic
961  !! and expressions that recover the answers for
962  !! remapping from the end of 2018. Otherwise,
963  !! use more robust forms of the same expressions.
964  ! Local variables
965  integer :: i, j, k, n_points
966  real :: dx(GV%ke+1)
967  real :: h_neglect, h_neglect_edge
968  logical :: ignore_vanished_layers, use_remapping_core_w, use_2018_remap
969 
970  ignore_vanished_layers = .false.
971  if (present(all_cells)) ignore_vanished_layers = .not. all_cells
972  use_remapping_core_w = .false.
973  if (present(old_remap)) use_remapping_core_w = old_remap
974  n_points = nk_src
975  use_2018_remap = .true. ; if (present(answers_2018)) use_2018_remap = answers_2018
976 
977  if (.not.use_2018_remap) then
978  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
979  elseif (gv%Boussinesq) then
980  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
981  else
982  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
983  endif
984 
985  !$OMP parallel do default(shared) firstprivate(n_points,dx)
986  do j = g%jsc,g%jec ; do i = g%isc,g%iec
987  if (g%mask2dT(i,j) > 0.) then
988  if (ignore_vanished_layers) then
989  n_points = 0
990  do k = 1, nk_src
991  if (h_src(i,j,k)>0.) n_points = n_points + 1
992  enddo
993  s_dst(i,j,:) = 0.
994  endif
995  if (use_remapping_core_w) then
996  call dzfromh1h2( n_points, h_src(i,j,1:n_points), gv%ke, h_dst(i,j,:), dx )
997  call remapping_core_w(cs, n_points, h_src(i,j,1:n_points), s_src(i,j,1:n_points), &
998  gv%ke, dx, s_dst(i,j,:), h_neglect, h_neglect_edge)
999  else
1000  call remapping_core_h(cs, n_points, h_src(i,j,1:n_points), s_src(i,j,1:n_points), &
1001  gv%ke, h_dst(i,j,:), s_dst(i,j,:), h_neglect, h_neglect_edge)
1002  endif
1003  else
1004  s_dst(i,j,:) = 0.
1005  endif
1006  enddo ; enddo
1007 

◆ ale_update_regrid_weights()

subroutine, public mom_ale::ale_update_regrid_weights ( real, intent(in)  dt,
type(ale_cs), pointer  CS 
)

Updates the weights for time filtering the new grid generated in regridding.

Parameters
[in]dtTime-step used between ALE calls [T ~> s]
csALE control structure

Definition at line 1235 of file MOM_ALE.F90.

1236  real, intent(in) :: dt !< Time-step used between ALE calls [T ~> s]
1237  type(ALE_CS), pointer :: CS !< ALE control structure
1238  ! Local variables
1239  real :: w ! An implicit weighting estimate.
1240 
1241  if (associated(cs)) then
1242  w = 0.0
1243  if (cs%regrid_time_scale > 0.0) then
1244  w = cs%regrid_time_scale / (cs%regrid_time_scale + dt)
1245  endif
1246  call set_regrid_params(cs%regridCS, old_grid_weight=w)
1247  endif
1248 

◆ ale_updateverticalgridtype()

subroutine, public mom_ale::ale_updateverticalgridtype ( type(ale_cs), pointer  CS,
type(verticalgrid_type), pointer  GV 
)

Update the vertical grid type with ALE information. This subroutine sets information in the verticalGrid_type to be consistent with the use of ALE mode.

Parameters
csALE control structure
gvvertical grid information

Definition at line 1254 of file MOM_ALE.F90.

1255  type(ALE_CS), pointer :: CS !< ALE control structure
1256  type(verticalGrid_type), pointer :: GV !< vertical grid information
1257 
1258  integer :: nk
1259 
1260  nk = gv%ke
1261  gv%sInterface(1:nk+1) = getcoordinateinterfaces( cs%regridCS, undo_scaling=.true. )
1262  gv%sLayer(1:nk) = 0.5*( gv%sInterface(1:nk) + gv%sInterface(2:nk+1) )
1263  gv%zAxisUnits = getcoordinateunits( cs%regridCS )
1264  gv%zAxisLongName = getcoordinateshortname( cs%regridCS )
1265  gv%direction = -1 ! Because of ferret in z* mode. Need method to set
1266  ! as function of coordinae mode.
1267 

◆ ale_writecoordinatefile()

subroutine, public mom_ale::ale_writecoordinatefile ( type(ale_cs), pointer  CS,
type(verticalgrid_type), intent(in)  GV,
character(len=*), intent(in)  directory 
)

Write the vertical coordinate information into a file. This subroutine writes out a file containing any available data related to the vertical grid used by the MOM ocean model when in ALE mode.

Parameters
csmodule control structure
[in]gvocean vertical grid structure
[in]directorydirectory for writing grid info

Definition at line 1274 of file MOM_ALE.F90.

1275  type(ALE_CS), pointer :: CS !< module control structure
1276  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
1277  character(len=*), intent(in) :: directory !< directory for writing grid info
1278 
1279  character(len=240) :: filepath
1280  type(vardesc) :: vars(2)
1281  type(fieldtype) :: fields(2)
1282  integer :: unit
1283  real :: ds(GV%ke), dsi(GV%ke+1)
1284 
1285  filepath = trim(directory) // trim("Vertical_coordinate")
1286  ds(:) = getcoordinateresolution( cs%regridCS, undo_scaling=.true. )
1287  dsi(1) = 0.5*ds(1)
1288  dsi(2:gv%ke) = 0.5*( ds(1:gv%ke-1) + ds(2:gv%ke) )
1289  dsi(gv%ke+1) = 0.5*ds(gv%ke)
1290 
1291  vars(1) = var_desc('ds', getcoordinateunits( cs%regridCS ), &
1292  'Layer Coordinate Thickness','1','L','1')
1293  vars(2) = var_desc('ds_interface', getcoordinateunits( cs%regridCS ), &
1294  'Layer Center Coordinate Separation','1','i','1')
1295 
1296  call create_file(unit, trim(filepath), vars, 2, fields, single_file, gv=gv)
1297  call write_field(unit, fields(1), ds)
1298  call write_field(unit, fields(2), dsi)
1299  call close_file(unit)
1300 

◆ check_grid()

subroutine mom_ale::check_grid ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
real, intent(in)  threshold 
)
private

Check grid for negative thicknesses.

Parameters
[in]gOcean grid structure
[in]gvOcean vertical grid structure
[in]hCurrent 3D grid obtained after the last time step [H ~> m or kg m-2]
[in]thresholdValue below which to flag issues, [H ~> m or kg m-2]

Definition at line 588 of file MOM_ALE.F90.

589  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
590  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
591  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Current 3D grid obtained after the
592  !! last time step [H ~> m or kg m-2]
593  real, intent(in) :: threshold !< Value below which to flag issues,
594  !! [H ~> m or kg m-2]
595  ! Local variables
596  integer :: i, j
597 
598  do j = g%jsc,g%jec ; do i = g%isc,g%iec
599  if (g%mask2dT(i,j)>0.) then
600  if (minval(h(i,j,:)) < threshold) then
601  write(0,*) 'check_grid: i,j=',i,j,'h(i,j,:)=',h(i,j,:)
602  if (threshold <= 0.) then
603  call mom_error(fatal,"MOM_ALE, check_grid: negative thickness encountered.")
604  else
605  call mom_error(fatal,"MOM_ALE, check_grid: too tiny thickness encountered.")
606  endif
607  endif
608  endif
609  enddo ; enddo
610 
611 

◆ remap_all_state_vars()

subroutine mom_ale::remap_all_state_vars ( type(remapping_cs), intent(in)  CS_remapping,
type(ale_cs), intent(in)  CS_ALE,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  h_old,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  h_new,
type(tracer_registry_type), pointer  Reg,
type(ocean_obc_type), pointer  OBC,
real, dimension(szi_(g),szj_(g),szk_(gv)+1), intent(in), optional  dxInterface,
real, dimension(szib_(g),szj_(g),szk_(gv)), intent(inout), optional  u,
real, dimension(szi_(g),szjb_(g),szk_(gv)), intent(inout), optional  v,
logical, intent(in), optional  debug,
real, intent(in), optional  dt 
)
private

This routine takes care of remapping all variable between the old and the new grids. When velocity components need to be remapped, thicknesses at velocity points are taken to be arithmetic averages of tracer thicknesses. This routine is called during initialization of the model at time=0, to remap initiali conditions to the model grid. It is also called during a time step to update the state.

Parameters
[in]cs_remappingRemapping control structure
[in]cs_aleALE control structure
[in]gOcean grid structure
[in]gvOcean vertical grid structure
[in]h_oldThickness of source grid [H ~> m or kg-2]
[in]h_newThickness of destination grid [H ~> m or kg-2]
regTracer registry structure
obcOpen boundary structure
[in]dxinterfaceChange in interface position
[in,out]uZonal velocity [L T-1 ~> m s-1]
[in,out]vMeridional velocity [L T-1 ~> m s-1]
[in]debugIf true, show the call tree
[in]dttime step for diagnostics [T ~> s]

Definition at line 744 of file MOM_ALE.F90.

746  type(remapping_CS), intent(in) :: CS_remapping !< Remapping control structure
747  type(ALE_CS), intent(in) :: CS_ALE !< ALE control structure
748  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
749  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
750  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h_old !< Thickness of source grid
751  !! [H ~> m or kg-2]
752  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h_new !< Thickness of destination grid
753  !! [H ~> m or kg-2]
754  type(tracer_registry_type), pointer :: Reg !< Tracer registry structure
755  type(ocean_OBC_type), pointer :: OBC !< Open boundary structure
756  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), &
757  optional, intent(in) :: dxInterface !< Change in interface position
758  !! [H ~> m or kg-2]
759  real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
760  optional, intent(inout) :: u !< Zonal velocity [L T-1 ~> m s-1]
761  real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
762  optional, intent(inout) :: v !< Meridional velocity [L T-1 ~> m s-1]
763  logical, optional, intent(in) :: debug !< If true, show the call tree
764  real, optional, intent(in) :: dt !< time step for diagnostics [T ~> s]
765  ! Local variables
766  integer :: i, j, k, m
767  integer :: nz, ntr
768  real, dimension(GV%ke+1) :: dx
769  real, dimension(GV%ke) :: h1, u_column
770  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)) :: work_conc
771  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)) :: work_cont
772  real, dimension(SZI_(G), SZJ_(G)) :: work_2d
773  real :: Idt ! The inverse of the timestep [T-1 ~> s-1]
774  real :: ppt2mks
775  real, dimension(GV%ke) :: h2
776  real :: h_neglect, h_neglect_edge
777  logical :: show_call_tree
778  type(tracer_type), pointer :: Tr => null()
779 
780  show_call_tree = .false.
781  if (present(debug)) show_call_tree = debug
782  if (show_call_tree) call calltree_enter("remap_all_state_vars(), MOM_ALE.F90")
783 
784  ! If remap_uv_using_old_alg is .true. and u or v is requested, then we must have dxInterface. Otherwise,
785  ! u and v can be remapped without dxInterface
786  if ( .not. present(dxinterface) .and. (cs_ale%remap_uv_using_old_alg .and. (present(u) .or. present(v))) ) then
787  call mom_error(fatal, "remap_all_state_vars: dxInterface must be present if using old algorithm "// &
788  "and u/v are to be remapped")
789  endif
790 
791  if (.not.cs_ale%answers_2018) then
792  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
793  elseif (gv%Boussinesq) then
794  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
795  else
796  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
797  endif
798 
799  nz = gv%ke
800  ppt2mks = 0.001
801 
802  ntr = 0 ; if (associated(reg)) ntr = reg%ntr
803 
804  if (present(dt)) then
805  idt = 1.0/dt
806  work_conc(:,:,:) = 0.0
807  work_cont(:,:,:) = 0.0
808  endif
809 
810  ! Remap tracer
811  if (ntr>0) then
812  if (show_call_tree) call calltree_waypoint("remapping tracers (remap_all_state_vars)")
813  !$OMP parallel do default(shared) private(h1,h2,u_column,Tr)
814  do m=1,ntr ! For each tracer
815  tr => reg%Tr(m)
816  do j = g%jsc,g%jec ; do i = g%isc,g%iec ; if (g%mask2dT(i,j)>0.) then
817  ! Build the start and final grids
818  h1(:) = h_old(i,j,:)
819  h2(:) = h_new(i,j,:)
820  call remapping_core_h(cs_remapping, nz, h1, tr%t(i,j,:), nz, h2, &
821  u_column, h_neglect, h_neglect_edge)
822 
823  ! Intermediate steps for tendency of tracer concentration and tracer content.
824  if (present(dt)) then
825  if (tr%id_remap_conc > 0) then
826  do k=1,gv%ke
827  work_conc(i,j,k) = (u_column(k) - tr%t(i,j,k)) * idt
828  enddo
829  endif
830  if (tr%id_remap_cont > 0 .or. tr%id_remap_cont_2d > 0) then
831  do k=1,gv%ke
832  work_cont(i,j,k) = (u_column(k)*h2(k) - tr%t(i,j,k)*h1(k)) * idt
833  enddo
834  endif
835  endif
836  ! update tracer concentration
837  tr%t(i,j,:) = u_column(:)
838  endif ; enddo ; enddo
839 
840  ! tendency diagnostics.
841  if (present(dt)) then
842  if (tr%id_remap_conc > 0) then
843  call post_data(tr%id_remap_conc, work_conc, cs_ale%diag)
844  endif
845  if (tr%id_remap_cont > 0) then
846  call post_data(tr%id_remap_cont, work_cont, cs_ale%diag)
847  endif
848  if (tr%id_remap_cont_2d > 0) then
849  do j = g%jsc,g%jec ; do i = g%isc,g%iec
850  work_2d(i,j) = 0.0
851  do k = 1,gv%ke
852  work_2d(i,j) = work_2d(i,j) + work_cont(i,j,k)
853  enddo
854  enddo ; enddo
855  call post_data(tr%id_remap_cont_2d, work_2d, cs_ale%diag)
856  endif
857  endif
858  enddo ! m=1,ntr
859 
860  endif ! endif for ntr > 0
861 
862  if (show_call_tree) call calltree_waypoint("tracers remapped (remap_all_state_vars)")
863 
864  ! Remap u velocity component
865  if ( present(u) ) then
866  !$OMP parallel do default(shared) private(h1,h2,dx,u_column)
867  do j = g%jsc,g%jec ; do i = g%iscB,g%iecB ; if (g%mask2dCu(i,j)>0.) then
868  ! Build the start and final grids
869  h1(:) = 0.5 * ( h_old(i,j,:) + h_old(i+1,j,:) )
870  if (cs_ale%remap_uv_using_old_alg) then
871  dx(:) = 0.5 * ( dxinterface(i,j,:) + dxinterface(i+1,j,:) )
872  do k = 1, nz
873  h2(k) = max( 0., h1(k) + ( dx(k+1) - dx(k) ) )
874  enddo
875  else
876  h2(:) = 0.5 * ( h_new(i,j,:) + h_new(i+1,j,:) )
877  endif
878  if (associated(obc)) then
879  if (obc%segnum_u(i,j) .ne. 0) then
880  if (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_e) then
881  h1(:) = h_old(i,j,:)
882  h2(:) = h_new(i,j,:)
883  else ! (OBC%segment(n)%direction == OBC_DIRECTION_W)
884  h1(:) = h_old(i+1,j,:)
885  h2(:) = h_new(i+1,j,:)
886  endif
887  endif
888  endif
889  call remapping_core_h(cs_remapping, nz, h1, u(i,j,:), nz, h2, &
890  u_column, h_neglect, h_neglect_edge)
891  u(i,j,:) = u_column(:)
892  endif ; enddo ; enddo
893  endif
894 
895  if (show_call_tree) call calltree_waypoint("u remapped (remap_all_state_vars)")
896 
897  ! Remap v velocity component
898  if ( present(v) ) then
899  !$OMP parallel do default(shared) private(h1,h2,dx,u_column)
900  do j = g%jscB,g%jecB ; do i = g%isc,g%iec ; if (g%mask2dCv(i,j)>0.) then
901  ! Build the start and final grids
902  h1(:) = 0.5 * ( h_old(i,j,:) + h_old(i,j+1,:) )
903  if (cs_ale%remap_uv_using_old_alg) then
904  dx(:) = 0.5 * ( dxinterface(i,j,:) + dxinterface(i,j+1,:) )
905  do k = 1, nz
906  h2(k) = max( 0., h1(k) + ( dx(k+1) - dx(k) ) )
907  enddo
908  else
909  h2(:) = 0.5 * ( h_new(i,j,:) + h_new(i,j+1,:) )
910  endif
911  if (associated(obc)) then
912  if (obc%segnum_v(i,j) .ne. 0) then
913  if (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_n) then
914  h1(:) = h_old(i,j,:)
915  h2(:) = h_new(i,j,:)
916  else ! (OBC%segment(n)%direction == OBC_DIRECTION_S)
917  h1(:) = h_old(i,j+1,:)
918  h2(:) = h_new(i,j+1,:)
919  endif
920  endif
921  endif
922  call remapping_core_h(cs_remapping, nz, h1, v(i,j,:), nz, h2, &
923  u_column, h_neglect, h_neglect_edge)
924  v(i,j,:) = u_column(:)
925  endif ; enddo ; enddo
926  endif
927 
928  if (cs_ale%id_vert_remap_h > 0) call post_data(cs_ale%id_vert_remap_h, h_old, cs_ale%diag)
929  if ((cs_ale%id_vert_remap_h_tendency > 0) .and. present(dt)) then
930  do k = 1, nz ; do j = g%jsc,g%jec ; do i = g%isc,g%iec
931  work_cont(i,j,k) = (h_new(i,j,k) - h_old(i,j,k))*idt
932  enddo ; enddo ; enddo
933  call post_data(cs_ale%id_vert_remap_h_tendency, work_cont, cs_ale%diag)
934  endif
935  if (show_call_tree) call calltree_waypoint("v remapped (remap_all_state_vars)")
936  if (show_call_tree) call calltree_leave("remap_all_state_vars()")
937 

◆ ts_plm_edge_values()

subroutine, public mom_ale::ts_plm_edge_values ( type(ale_cs), intent(inout)  CS,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  S_t,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  S_b,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  T_t,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  T_b,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  h,
logical, intent(in)  bdry_extrap 
)

Calculate edge values (top and bottom of layer) for T and S consistent with a PLM reconstruction in the vertical direction. Boundary reconstructions are PCM unless bdry_extrap is true.

Parameters
[in]gocean grid structure
[in]gvOcean vertical grid structure
[in,out]csmodule control structure
[in,out]s_tSalinity at the top edge of each layer
[in,out]s_bSalinity at the bottom edge of each layer
[in,out]t_tTemperature at the top edge of each layer
[in,out]t_bTemperature at the bottom edge of each layer
[in]tvthermodynamics structure
[in]hlayer thickness [H ~> m or kg m-2]
[in]bdry_extrapIf true, use high-order boundary extrapolation within boundary cells

Definition at line 1013 of file MOM_ALE.F90.

1014  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
1015  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1016  type(ALE_CS), intent(inout) :: CS !< module control structure
1017  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1018  intent(inout) :: S_t !< Salinity at the top edge of each layer
1019  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1020  intent(inout) :: S_b !< Salinity at the bottom edge of each layer
1021  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1022  intent(inout) :: T_t !< Temperature at the top edge of each layer
1023  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1024  intent(inout) :: T_b !< Temperature at the bottom edge of each layer
1025  type(thermo_var_ptrs), intent(in) :: tv !< thermodynamics structure
1026  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1027  intent(in) :: h !< layer thickness [H ~> m or kg m-2]
1028  logical, intent(in) :: bdry_extrap !< If true, use high-order boundary
1029  !! extrapolation within boundary cells
1030 
1031  call ale_plm_edge_values( cs, g, gv, h, tv%S, bdry_extrap, s_t, s_b )
1032  call ale_plm_edge_values( cs, g, gv, h, tv%T, bdry_extrap, t_t, t_b )
1033 

◆ ts_ppm_edge_values()

subroutine, public mom_ale::ts_ppm_edge_values ( type(ale_cs), intent(inout)  CS,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  S_t,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  S_b,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  T_t,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  T_b,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  h,
logical, intent(in)  bdry_extrap 
)

Calculate edge values (top and bottom of layer) for T and S consistent with a PPM reconstruction in the vertical direction. Boundary reconstructions are PCM unless bdry_extrap is true.

Parameters
[in]gocean grid structure
[in]gvOcean vertical grid structure
[in,out]csmodule control structure
[in,out]s_tSalinity at the top edge of each layer
[in,out]s_bSalinity at the bottom edge of each layer
[in,out]t_tTemperature at the top edge of each layer
[in,out]t_bTemperature at the bottom edge of each layer
[in]tvthermodynamics structure
[in]hlayer thicknesses [H ~> m or kg m-2]
[in]bdry_extrapIf true, use high-order boundary extrapolation within boundary cells

Definition at line 1099 of file MOM_ALE.F90.

1100  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
1101  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1102  type(ALE_CS), intent(inout) :: CS !< module control structure
1103  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1104  intent(inout) :: S_t !< Salinity at the top edge of each layer
1105  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1106  intent(inout) :: S_b !< Salinity at the bottom edge of each layer
1107  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1108  intent(inout) :: T_t !< Temperature at the top edge of each layer
1109  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1110  intent(inout) :: T_b !< Temperature at the bottom edge of each layer
1111  type(thermo_var_ptrs), intent(in) :: tv !< thermodynamics structure
1112  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1113  intent(in) :: h !< layer thicknesses [H ~> m or kg m-2]
1114  logical, intent(in) :: bdry_extrap !< If true, use high-order boundary
1115  !! extrapolation within boundary cells
1116 
1117  ! Local variables
1118  integer :: i, j, k
1119  real :: hTmp(GV%ke) ! A 1-d copy of h [H ~> m or kg m-2]
1120  real :: tmp(GV%ke) ! A 1-d copy of a column of temperature [degC] or salinity [ppt]
1121  real, dimension(CS%nk,2) :: &
1122  ppol_E ! Edge value of polynomial in [degC] or [ppt]
1123  real, dimension(CS%nk,3) :: &
1124  ppol_coefs ! Coefficients of polynomial, all in [degC] or [ppt]
1125  real :: h_neglect, h_neglect_edge ! Tiny thicknesses [H ~> m or kg m-2]
1126 
1127  if (.not.cs%answers_2018) then
1128  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
1129  elseif (gv%Boussinesq) then
1130  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
1131  else
1132  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
1133  endif
1134 
1135  ! Determine reconstruction within each column
1136  !$OMP parallel do default(shared) private(hTmp,tmp,ppol_E,ppol_coefs)
1137  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
1138 
1139  ! Build current grid
1140  htmp(:) = h(i,j,:)
1141  tmp(:) = tv%S(i,j,:)
1142 
1143  ! Reconstruct salinity profile
1144  ppol_e(:,:) = 0.0
1145  ppol_coefs(:,:) = 0.0
1146  call edge_values_implicit_h4( gv%ke, htmp, tmp, ppol_e, h_neglect=h_neglect_edge, &
1147  answers_2018=cs%answers_2018 )
1148  call ppm_reconstruction( gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect, &
1149  answers_2018=cs%answers_2018 )
1150  if (bdry_extrap) &
1151  call ppm_boundary_extrapolation( gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect )
1152 
1153  do k = 1,gv%ke
1154  s_t(i,j,k) = ppol_e(k,1)
1155  s_b(i,j,k) = ppol_e(k,2)
1156  enddo
1157 
1158  ! Reconstruct temperature profile
1159  ppol_e(:,:) = 0.0
1160  ppol_coefs(:,:) = 0.0
1161  tmp(:) = tv%T(i,j,:)
1162  if (cs%answers_2018) then
1163  call edge_values_implicit_h4( gv%ke, htmp, tmp, ppol_e, h_neglect=1.0e-10*gv%m_to_H, &
1164  answers_2018=cs%answers_2018 )
1165  else
1166  call edge_values_implicit_h4( gv%ke, htmp, tmp, ppol_e, h_neglect=gv%H_subroundoff, &
1167  answers_2018=cs%answers_2018 )
1168  endif
1169  call ppm_reconstruction( gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect, &
1170  answers_2018=cs%answers_2018 )
1171  if (bdry_extrap) &
1172  call ppm_boundary_extrapolation(gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect )
1173 
1174  do k = 1,gv%ke
1175  t_t(i,j,k) = ppol_e(k,1)
1176  t_b(i,j,k) = ppol_e(k,2)
1177  enddo
1178 
1179  enddo ; enddo
1180