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)
 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 pressure_gradient_plm (CS, S_t, S_b, T_t, T_b, G, GV, tv, h, bdry_extrap)
 Use plm reconstruction for pressure gradient (determine edge values) By using a PLM (limited piecewise linear method) reconstruction, this routine determines the edge values for the salinity and temperature within each layer. These edge values are returned and are used to compute the pressure gradient (by computing the densities). More...
 
subroutine, public pressure_gradient_ppm (CS, S_t, S_b, T_t, T_b, G, GV, tv, h, bdry_extrap)
 Use ppm reconstruction for pressure gradient (determine edge values) By using a PPM (limited piecewise linear method) reconstruction, this routine determines the edge values for the salinity and temperature within each layer. These edge values are returned and are used to compute the pressure gradient (by computing the densities). 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 295 of file MOM_ALE.F90.

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

◆ 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 618 of file MOM_ALE.F90.

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

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Here is the call graph for this function:

◆ 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 309 of file MOM_ALE.F90.

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

◆ 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 1187 of file MOM_ALE.F90.

1187  type(ALE_CS), pointer :: CS !< module control structure
1188 
1189  real, dimension(CS%nk+1) :: ALE_getCoordinate
1190  ale_getcoordinate(:) = getcoordinateinterfaces( cs%regridCS, undo_scaling=.true. )
1191 

◆ 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 1197 of file MOM_ALE.F90.

1197  type(ALE_CS), pointer :: CS !< module control structure
1198 
1199  character(len=20) :: ALE_getCoordinateUnits
1200 
1201  ale_getcoordinateunits = getcoordinateunits( cs%regridCS )
1202 

◆ 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 141 of file MOM_ALE.F90.

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

References ale_initregridding(), mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_oda_driver_mod::init_oda().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 1166 of file MOM_ALE.F90.

1166  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1167  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1168  real, intent(in) :: max_depth !< The maximum depth of the ocean [Z ~> m].
1169  type(param_file_type), intent(in) :: param_file !< parameter file
1170  character(len=*), intent(in) :: mdl !< Name of calling module
1171  type(regridding_CS), intent(out) :: regridCS !< Regridding parameters and work arrays
1172  ! Local variables
1173  character(len=30) :: coord_mode
1174 
1175  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_MODE", coord_mode, &
1176  "Coordinate mode for vertical regridding. "//&
1177  "Choose among the following possibilities: "//&
1178  trim(regriddingcoordinatemodedoc), &
1179  default=default_coordinate_mode, fail_if_missing=.true.)
1180 
1181  call initialize_regridding(regridcs, gv, us, max_depth, param_file, mdl, coord_mode, '', '')
1182 

Referenced by ale_init().

Here is the caller graph for this function:

◆ 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 1286 of file MOM_ALE.F90.

1286  type(ALE_CS), intent(inout) :: CS !< module control structure
1287  type(ocean_grid_type), intent(in) :: G !< module grid structure
1288  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1289  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(out) :: h !< layer thickness [H ~> m or kg m-2]
1290 
1291  ! Local variables
1292  integer :: i, j, k
1293 
1294  do j = g%jsd,g%jed ; do i = g%isd,g%ied
1295  h(i,j,:) = gv%Z_to_H * getstaticthickness( cs%regridCS, 0., g%bathyT(i,j) )
1296  enddo ; enddo
1297 

References mom_regridding::getstaticthickness().

Referenced by mom_oda_driver_mod::init_oda().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 324 of file MOM_ALE.F90.

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

References ale_update_regrid_weights(), mom_error_handler::calltree_enter(), mom_error_handler::calltree_leave(), mom_error_handler::calltree_waypoint(), check_grid(), mom_diag_mediator::diag_update_remap_grids(), and remap_all_state_vars().

Here is the call graph for this function:

◆ 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 411 of file MOM_ALE.F90.

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

References ale_update_regrid_weights(), mom_error_handler::calltree_enter(), mom_error_handler::calltree_leave(), mom_error_handler::calltree_waypoint(), check_grid(), and remap_all_state_vars().

Referenced by mom_offline_main::offline_advection_ale().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 465 of file MOM_ALE.F90.

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

References ale_remap_scalar(), mom_error_handler::calltree_leave(), mom_error_handler::calltree_waypoint(), mom_debugging::check_column_integrals(), check_grid(), mom_diag_vkernels::interpolate_column(), mom_tracer_registry::mom_tracer_chkinv(), mom_diag_vkernels::reintegrate_column(), and remap_all_state_vars().

Referenced by mom_offline_main::update_offline_fields().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 546 of file MOM_ALE.F90.

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

References mom_error_handler::calltree_enter(), mom_error_handler::calltree_leave(), mom_error_handler::calltree_waypoint(), check_grid(), and remap_all_state_vars().

Referenced by mom::step_offline().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 252 of file MOM_ALE.F90.

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

References mom_verticalgrid::get_thickness_units().

Here is the call graph for this function:

◆ 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 662 of file MOM_ALE.F90.

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

References ale_update_regrid_weights(), mom_domains::do_group_pass(), and remap_all_state_vars().

Referenced by mom_state_initialization::mom_initialize_state().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 1208 of file MOM_ALE.F90.

1208  type(ALE_CS), pointer :: CS !< module control structure
1209 
1210  ale_remap_init_conds = .false.
1211  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 
)

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".

Definition at line 946 of file MOM_ALE.F90.

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

References mom_remapping::dzfromh1h2().

Referenced by ale_offline_inputs(), and mom_tracer_initialization_from_z::mom_initialize_tracer_from_z().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 1216 of file MOM_ALE.F90.

1216  real, intent(in) :: dt !< Time-step used between ALE calls [T ~> s]
1217  type(ALE_CS), pointer :: CS !< ALE control structure
1218  ! Local variables
1219  real :: w ! An implicit weighting estimate.
1220 
1221  if (associated(cs)) then
1222  w = 0.0
1223  if (cs%regrid_time_scale > 0.0) then
1224  w = cs%regrid_time_scale / (cs%regrid_time_scale + dt)
1225  endif
1226  call set_regrid_params(cs%regridCS, old_grid_weight=w)
1227  endif
1228 

Referenced by ale_main(), ale_main_offline(), and ale_regrid_accelerated().

Here is the caller graph for this function:

◆ 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 1235 of file MOM_ALE.F90.

1235  type(ALE_CS), pointer :: CS !< ALE control structure
1236  type(verticalGrid_type), pointer :: GV !< vertical grid information
1237 
1238  integer :: nk
1239 
1240  nk = gv%ke
1241  gv%sInterface(1:nk+1) = getcoordinateinterfaces( cs%regridCS, undo_scaling=.true. )
1242  gv%sLayer(1:nk) = 0.5*( gv%sInterface(1:nk) + gv%sInterface(2:nk+1) )
1243  gv%zAxisUnits = getcoordinateunits( cs%regridCS )
1244  gv%zAxisLongName = getcoordinateshortname( cs%regridCS )
1245  gv%direction = -1 ! Because of ferret in z* mode. Need method to set
1246  ! as function of coordinae mode.
1247 

Referenced by mom_oda_driver_mod::init_oda().

Here is the caller graph for this function:

◆ 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 1255 of file MOM_ALE.F90.

1255  type(ALE_CS), pointer :: CS !< module control structure
1256  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
1257  character(len=*), intent(in) :: directory !< directory for writing grid info
1258 
1259  character(len=240) :: filepath
1260  type(vardesc) :: vars(2)
1261  type(fieldtype) :: fields(2)
1262  integer :: unit
1263  real :: ds(GV%ke), dsi(GV%ke+1)
1264 
1265  filepath = trim(directory) // trim("Vertical_coordinate")
1266  ds(:) = getcoordinateresolution( cs%regridCS, undo_scaling=.true. )
1267  dsi(1) = 0.5*ds(1)
1268  dsi(2:gv%ke) = 0.5*( ds(1:gv%ke-1) + ds(2:gv%ke) )
1269  dsi(gv%ke+1) = 0.5*ds(gv%ke)
1270 
1271  vars(1) = var_desc('ds', getcoordinateunits( cs%regridCS ), &
1272  'Layer Coordinate Thickness','1','L','1')
1273  vars(2) = var_desc('ds_interface', getcoordinateunits( cs%regridCS ), &
1274  'Layer Center Coordinate Separation','1','i','1')
1275 
1276  call create_file(unit, trim(filepath), vars, 2, fields, single_file, gv=gv)
1277  call write_field(unit, fields(1), ds)
1278  call write_field(unit, fields(2), dsi)
1279  call close_file(unit)
1280 

References mom_io::create_file().

Here is the call graph for this function:

◆ 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 591 of file MOM_ALE.F90.

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

Referenced by ale_main(), ale_main_offline(), ale_offline_inputs(), and ale_offline_tracer_final().

Here is the caller graph for this function:

◆ pressure_gradient_plm()

subroutine, public mom_ale::pressure_gradient_plm ( 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 
)

Use plm reconstruction for pressure gradient (determine edge values) By using a PLM (limited piecewise linear method) reconstruction, this routine determines the edge values for the salinity and temperature within each layer. These edge values are returned and are used to compute the pressure gradient (by computing the densities).

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 1012 of file MOM_ALE.F90.

1012  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
1013  type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure
1014  type(ALE_CS), intent(inout) :: CS !< module control structure
1015  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1016  intent(inout) :: S_t !< Salinity at the top edge of each layer
1017  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1018  intent(inout) :: S_b !< Salinity at the bottom edge of each layer
1019  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1020  intent(inout) :: T_t !< Temperature at the top edge of each layer
1021  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1022  intent(inout) :: T_b !< Temperature at the bottom edge of each layer
1023  type(thermo_var_ptrs), intent(in) :: tv !< thermodynamics structure
1024  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1025  intent(in) :: h !< layer thickness [H ~> m or kg m-2]
1026  logical, intent(in) :: bdry_extrap !< If true, use high-order boundary
1027  !! extrapolation within boundary cells
1028 
1029  ! Local variables
1030  integer :: i, j, k
1031  real :: hTmp(GV%ke)
1032  real :: tmp(GV%ke)
1033  real, dimension(CS%nk,2) :: ppol_E !Edge value of polynomial
1034  real, dimension(CS%nk,2) :: ppol_coefs !Coefficients of polynomial
1035  real :: h_neglect
1036 
1037  !### Replace this with GV%H_subroundoff
1038  if (gv%Boussinesq) then
1039  h_neglect = gv%m_to_H*1.0e-30
1040  else
1041  h_neglect = gv%kg_m2_to_H*1.0e-30
1042  endif
1043 
1044  ! Determine reconstruction within each column
1045  !$OMP parallel do default(shared) private(hTmp,ppol_E,ppol_coefs,tmp)
1046  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
1047  ! Build current grid
1048  htmp(:) = h(i,j,:)
1049  tmp(:) = tv%S(i,j,:)
1050  ! Reconstruct salinity profile
1051  ppol_e(:,:) = 0.0
1052  ppol_coefs(:,:) = 0.0
1053  call plm_reconstruction( gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect )
1054  if (bdry_extrap) &
1055  call plm_boundary_extrapolation( gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect )
1056 
1057  do k = 1,gv%ke
1058  s_t(i,j,k) = ppol_e(k,1)
1059  s_b(i,j,k) = ppol_e(k,2)
1060  enddo
1061 
1062  ! Reconstruct temperature profile
1063  ppol_e(:,:) = 0.0
1064  ppol_coefs(:,:) = 0.0
1065  tmp(:) = tv%T(i,j,:)
1066  call plm_reconstruction( gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect )
1067  if (bdry_extrap) &
1068  call plm_boundary_extrapolation( gv%ke, htmp, tmp, ppol_e, ppol_coefs, h_neglect )
1069 
1070  do k = 1,gv%ke
1071  t_t(i,j,k) = ppol_e(k,1)
1072  t_b(i,j,k) = ppol_e(k,2)
1073  enddo
1074 
1075  enddo ; enddo
1076 

References plm_functions::plm_boundary_extrapolation(), and plm_functions::plm_reconstruction().

Referenced by mom_pressureforce_afv::pressureforce_afv_bouss(), mom_pressureforce_afv::pressureforce_afv_nonbouss(), and mom_pressureforce_blk_afv::pressureforce_blk_afv_bouss().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ pressure_gradient_ppm()

subroutine, public mom_ale::pressure_gradient_ppm ( 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 
)

Use ppm reconstruction for pressure gradient (determine edge values) By using a PPM (limited piecewise linear method) reconstruction, this routine determines the edge values for the salinity and temperature within each layer. These edge values are returned and are used to compute the pressure gradient (by computing the densities).

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 1086 of file MOM_ALE.F90.

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

References regrid_edge_values::edge_values_implicit_h4(), ppm_functions::ppm_boundary_extrapolation(), and ppm_functions::ppm_reconstruction().

Referenced by mom_pressureforce_afv::pressureforce_afv_bouss(), mom_pressureforce_afv::pressureforce_afv_nonbouss(), and mom_pressureforce_blk_afv::pressureforce_blk_afv_bouss().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 748 of file MOM_ALE.F90.

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

References mom_error_handler::calltree_enter(), mom_error_handler::calltree_leave(), mom_error_handler::calltree_waypoint(), and mom_open_boundary::obc_direction_n.

Referenced by ale_main(), ale_main_offline(), ale_offline_inputs(), ale_offline_tracer_final(), and ale_regrid_accelerated().

Here is the call graph for this function:
Here is the caller graph for this function: