MOM6
|
This module contains the main regridding routines.
Regridding comprises two steps:
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... | |
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.
cs | Regridding parameters and options | |
[in] | g | Ocean grid informations |
[in] | gv | Ocean vertical grid structure |
[in,out] | h | Current 3D grid thickness that are to be adjusted [H ~> m or kg-2] |
Definition at line 295 of file MOM_ALE.F90.
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.
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | regridcs | Regridding parameters and options |
[in] | remapcs | Remapping parameters and options |
[in,out] | tv | Thermodynamical variable structure |
[in,out] | h | Current 3D grid obtained after the last time step [H ~> m or kg-2] |
[in] | debug | If true, show the call tree |
frac_shelf_h | Fractional ice shelf coverage |
Definition at line 618 of file MOM_ALE.F90.
References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().
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.
cs | module control structure |
Definition at line 309 of file MOM_ALE.F90.
real function, dimension(cs%nk+1), public mom_ale::ale_getcoordinate | ( | type(ale_cs), pointer | CS | ) |
Query the target coordinate interfaces positions.
cs | module control structure |
Definition at line 1187 of file MOM_ALE.F90.
Query the target coordinate units.
cs | module control structure |
Definition at line 1197 of file MOM_ALE.F90.
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.
[in] | param_file | Parameter file |
[in] | gv | Ocean vertical grid structure |
[in] | us | A dimensional unit scaling type |
[in] | max_depth | The maximum depth of the ocean [Z ~> m]. |
cs | Module control structure |
Definition at line 141 of file MOM_ALE.F90.
References ale_initregridding(), mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().
Referenced by mom_oda_driver_mod::init_oda().
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.
[in] | gv | Ocean vertical grid structure |
[in] | us | A dimensional unit scaling type |
[in] | max_depth | The maximum depth of the ocean [Z ~> m]. |
[in] | param_file | parameter file |
[in] | mdl | Name of calling module |
[out] | regridcs | Regridding parameters and work arrays |
Definition at line 1166 of file MOM_ALE.F90.
Referenced by ale_init().
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.
[in,out] | cs | module control structure |
[in] | g | module grid structure |
[in] | gv | Ocean vertical grid structure |
[out] | h | layer thickness [H ~> m or kg m-2] |
Definition at line 1286 of file MOM_ALE.F90.
References mom_regridding::getstaticthickness().
Referenced by mom_oda_driver_mod::init_oda().
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.
[in] | g | Ocean grid informations |
[in] | gv | Ocean vertical grid structure |
[in] | us | A dimensional unit scaling type |
[in,out] | h | Current 3D grid obtained after the last time step [H ~> m or kg m-2] |
[in,out] | u | Zonal velocity field [L T-1 ~> m s-1] |
[in,out] | v | Meridional velocity field [L T-1 ~> m s-1] |
[in,out] | tv | Thermodynamic variable structure |
reg | Tracer registry structure | |
cs | Regridding parameters and options | |
obc | Open boundary structure | |
[in] | dt | Time step between calls to ALE_main [T ~> s] |
frac_shelf_h | Fractional ice shelf coverage |
Definition at line 324 of file MOM_ALE.F90.
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().
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.
[in] | g | Ocean grid informations |
[in] | gv | Ocean vertical grid structure |
[in,out] | h | Current 3D grid obtained after the last time step [H ~> m or kg-2] |
[in,out] | tv | Thermodynamic variable structure |
reg | Tracer registry structure | |
cs | Regridding parameters and options | |
obc | Open boundary structure | |
[in] | dt | Time step between calls to ALE_main [T ~> s] |
Definition at line 411 of file MOM_ALE.F90.
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().
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.
cs | Regridding parameters and options | |
[in] | g | Ocean grid informations |
[in] | gv | Ocean vertical grid structure |
[in,out] | h | Layer thicknesses |
[in,out] | tv | Thermodynamic variable structure |
reg | Tracer registry structure | |
[in,out] | uhtr | Zonal mass fluxes |
[in,out] | vhtr | Meridional mass fluxes |
[in,out] | kd | Input diffusivites |
[in] | debug | If true, then turn checksums |
obc | Open boundary structure |
Definition at line 465 of file MOM_ALE.F90.
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().
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.
[in] | g | Ocean grid informations |
[in] | gv | Ocean vertical grid structure |
[in,out] | h | Current 3D grid obtained after the last time step [H ~> m or kg-2] |
[in,out] | tv | Thermodynamic variable structure |
[in,out] | h_target | Current 3D grid obtained after last time step [H ~> m or kg-2] |
reg | Tracer registry structure | |
cs | Regridding parameters and options | |
obc | Open boundary structure |
Definition at line 546 of file MOM_ALE.F90.
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().
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.
[in] | time | Time structure |
[in] | g | Grid structure |
[in] | us | A dimensional unit scaling type |
[in] | gv | Ocean vertical grid structure |
[in] | diag | Diagnostics control structure |
cs | Module control structure |
Definition at line 252 of file MOM_ALE.F90.
References mom_verticalgrid::get_thickness_units().
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.
cs | ALE control structure | |
[in,out] | g | Ocean grid |
[in] | gv | Vertical grid |
[in,out] | h | Original thicknesses [H ~> m or kg-2] |
[in,out] | tv | Thermo vars (T/S/EOS) |
[in] | n | Number of times to regrid |
[in,out] | u | Zonal velocity [L T-1 ~> m s-1] |
[in,out] | v | Meridional velocity [L T-1 ~> m s-1] |
obc | Open boundary structure | |
reg | Tracer registry to remap onto new grid | |
[in] | dt | Model timestep to provide a timescale for regridding [T ~> s] |
[in,out] | dzregrid | Final change in interface positions |
[in] | initial | Whether we're being called from an initialization routine (and expect diagnostics to work) |
Definition at line 662 of file MOM_ALE.F90.
References ale_update_regrid_weights(), mom_domains::do_group_pass(), and remap_all_state_vars().
Referenced by mom_state_initialization::mom_initialize_state().
logical function, public mom_ale::ale_remap_init_conds | ( | type(ale_cs), pointer | CS | ) |
Returns true if initial conditions should be regridded and remapped.
cs | module control structure |
Definition at line 1208 of file MOM_ALE.F90.
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.
[in] | cs | Remapping control structure |
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | nk_src | Number of levels on source grid |
[in] | h_src | Level thickness of source grid [H ~> m or kg-2] |
[in] | s_src | Scalar on source grid |
[in] | h_dst | Level thickness of destination grid [H ~> m or kg-2] |
[in,out] | s_dst | Scalar on destination grid |
[in] | all_cells | If false, only reconstruct for non-vanished cells. Use all vanished layers otherwise (default). |
[in] | old_remap | If true, use the old "remapping_core_w" method, otherwise use "remapping_core_h". |
Definition at line 946 of file MOM_ALE.F90.
References mom_remapping::dzfromh1h2().
Referenced by ale_offline_inputs(), and mom_tracer_initialization_from_z::mom_initialize_tracer_from_z().
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.
[in] | dt | Time-step used between ALE calls [T ~> s] |
cs | ALE control structure |
Definition at line 1216 of file MOM_ALE.F90.
Referenced by ale_main(), ale_main_offline(), and ale_regrid_accelerated().
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.
cs | ALE control structure |
gv | vertical grid information |
Definition at line 1235 of file MOM_ALE.F90.
Referenced by mom_oda_driver_mod::init_oda().
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.
cs | module control structure | |
[in] | gv | ocean vertical grid structure |
[in] | directory | directory for writing grid info |
Definition at line 1255 of file MOM_ALE.F90.
References mom_io::create_file().
|
private |
Check grid for negative thicknesses.
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h | Current 3D grid obtained after the last time step [H ~> m or kg m-2] |
[in] | threshold | Value below which to flag issues, [H ~> m or kg m-2] |
Definition at line 591 of file MOM_ALE.F90.
Referenced by ale_main(), ale_main_offline(), ale_offline_inputs(), and ale_offline_tracer_final().
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).
[in] | g | ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in,out] | cs | module control structure |
[in,out] | s_t | Salinity at the top edge of each layer |
[in,out] | s_b | Salinity at the bottom edge of each layer |
[in,out] | t_t | Temperature at the top edge of each layer |
[in,out] | t_b | Temperature at the bottom edge of each layer |
[in] | tv | thermodynamics structure |
[in] | h | layer thickness [H ~> m or kg m-2] |
[in] | bdry_extrap | If true, use high-order boundary extrapolation within boundary cells |
Definition at line 1012 of file MOM_ALE.F90.
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().
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).
[in] | g | ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in,out] | cs | module control structure |
[in,out] | s_t | Salinity at the top edge of each layer |
[in,out] | s_b | Salinity at the bottom edge of each layer |
[in,out] | t_t | Temperature at the top edge of each layer |
[in,out] | t_b | Temperature at the bottom edge of each layer |
[in] | tv | thermodynamics structure |
[in] | h | layer thicknesses [H ~> m or kg m-2] |
[in] | bdry_extrap | If true, use high-order boundary extrapolation within boundary cells |
Definition at line 1086 of file MOM_ALE.F90.
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().
|
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.
[in] | cs_remapping | Remapping control structure |
[in] | cs_ale | ALE control structure |
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h_old | Thickness of source grid [H ~> m or kg-2] |
[in] | h_new | Thickness of destination grid [H ~> m or kg-2] |
reg | Tracer registry structure | |
obc | Open boundary structure | |
[in] | dxinterface | Change in interface position |
[in,out] | u | Zonal velocity [L T-1 ~> m s-1] |
[in,out] | v | Meridional velocity [L T-1 ~> m s-1] |
[in] | debug | If true, show the call tree |
[in] | dt | time step for diagnostics [T ~> s] |
Definition at line 748 of file MOM_ALE.F90.
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().