MOM6
coord_zlike Module Reference

Detailed Description

Regrid columns for a z-like coordinate (z-star, z-level)

Data Types

type  zlike_cs
 Control structure containing required parameters for a z-like coordinate. More...
 

Functions/Subroutines

subroutine, public init_coord_zlike (CS, nk, coordinateResolution)
 Initialise a zlike_CS with pointers to parameters. More...
 
subroutine, public end_coord_zlike (CS)
 Deallocates the zlike control structure. More...
 
subroutine, public set_zlike_params (CS, min_thickness)
 Set parameters in the zlike structure. More...
 
subroutine, public build_zstar_column (CS, depth, total_thickness, zInterface, z_rigid_top, eta_orig, zScale)
 Builds a z* coordinate with a minimum thickness. More...
 

Function/Subroutine Documentation

◆ build_zstar_column()

subroutine, public coord_zlike::build_zstar_column ( type(zlike_cs), intent(in)  CS,
real, intent(in)  depth,
real, intent(in)  total_thickness,
real, dimension(cs%nk+1), intent(inout)  zInterface,
real, intent(in), optional  z_rigid_top,
real, intent(in), optional  eta_orig,
real, intent(in), optional  zScale 
)

Builds a z* coordinate with a minimum thickness.

Parameters
[in]csCoordinate control structure
[in]depthDepth of ocean bottom (positive in the output units)
[in]total_thicknessColumn thickness (positive in the same units as depth)
[in,out]zinterfaceAbsolute positions of interfaces
[in]z_rigid_topThe height of a rigid top (negative in the same units as depth)
[in]eta_origThe actual original height of the top in the same units as depth
[in]zscaleScaling factor from the target coordinate resolution in Z to desired units for zInterface, perhaps Z_to_H

Definition at line 65 of file coord_zlike.F90.

65  type(zlike_CS), intent(in) :: CS !< Coordinate control structure
66  real, intent(in) :: depth !< Depth of ocean bottom (positive in the output units)
67  real, intent(in) :: total_thickness !< Column thickness (positive in the same units as depth)
68  real, dimension(CS%nk+1), intent(inout) :: zInterface !< Absolute positions of interfaces
69  real, optional, intent(in) :: z_rigid_top !< The height of a rigid top (negative in the
70  !! same units as depth)
71  real, optional, intent(in) :: eta_orig !< The actual original height of the top in the
72  !! same units as depth
73  real, optional, intent(in) :: zScale !< Scaling factor from the target coordinate resolution
74  !! in Z to desired units for zInterface, perhaps Z_to_H
75  ! Local variables
76  real :: eta, stretching, dh, min_thickness, z0_top, z_star, z_scale
77  integer :: k
78  logical :: new_zstar_def
79 
80  z_scale = 1.0 ; if (present(zscale)) z_scale = zscale
81 
82  new_zstar_def = .false.
83  min_thickness = min( cs%min_thickness, total_thickness/real(cs%nk) )
84  z0_top = 0.
85  if (present(z_rigid_top)) then
86  z0_top = z_rigid_top
87  new_zstar_def = .true.
88  endif
89 
90  ! Position of free-surface (or the rigid top, for which eta ~ z0_top)
91  eta = total_thickness - depth
92  if (present(eta_orig)) eta = eta_orig
93 
94  ! Conventional z* coordinate:
95  ! z* = (z-eta) / stretching where stretching = (H+eta)/H
96  ! z = eta + stretching * z*
97  ! The above gives z*(z=eta) = 0, z*(z=-H) = -H.
98  ! With a rigid top boundary at eta = z0_top then
99  ! z* = z0 + (z-eta) / stretching where stretching = (H+eta)/(H+z0)
100  ! z = eta + stretching * (z*-z0) * stretching
101  stretching = total_thickness / ( depth + z0_top )
102 
103  if (new_zstar_def) then
104  ! z_star is the notional z* coordinate in absence of upper/lower topography
105  z_star = 0. ! z*=0 at the free-surface
106  zinterface(1) = eta ! The actual position of the top of the column
107  do k = 2,cs%nk
108  z_star = z_star - cs%coordinateResolution(k-1)*z_scale
109  ! This ensures that z is below a rigid upper surface (ice shelf bottom)
110  zinterface(k) = min( eta + stretching * ( z_star - z0_top ), z0_top )
111  ! This ensures that the layer in inflated
112  zinterface(k) = min( zinterface(k), zinterface(k-1) - min_thickness )
113  ! This ensures that z is above or at the topography
114  zinterface(k) = max( zinterface(k), -depth + real(cs%nk+1-k) * min_thickness )
115  enddo
116  zinterface(cs%nk+1) = -depth
117 
118  else
119  ! Integrate down from the top for a notional new grid, ignoring topography
120  ! The starting position is offset by z0_top which, if z0_top<0, will place
121  ! interfaces above the rigid boundary.
122  zinterface(1) = eta
123  do k = 1,cs%nk
124  dh = stretching * cs%coordinateResolution(k)*z_scale ! Notional grid spacing
125  zinterface(k+1) = zinterface(k) - dh
126  enddo
127 
128  ! Integrating up from the bottom adjusting interface position to accommodate
129  ! inflating layers without disturbing the interface above
130  zinterface(cs%nk+1) = -depth
131  do k = cs%nk,1,-1
132  if ( zinterface(k) < (zinterface(k+1) + min_thickness) ) then
133  zinterface(k) = zinterface(k+1) + min_thickness
134  endif
135  enddo
136  endif
137 

Referenced by mom_regridding::build_zstar_grid(), and mom_diag_remap::diag_remap_update().

Here is the caller graph for this function:

◆ end_coord_zlike()

subroutine, public coord_zlike::end_coord_zlike ( type(zlike_cs), pointer  CS)

Deallocates the zlike control structure.

Parameters
csCoordinate control structure

Definition at line 44 of file coord_zlike.F90.

44  type(zlike_CS), pointer :: CS !< Coordinate control structure
45 
46  ! Nothing to do
47  if (.not. associated(cs)) return
48  deallocate(cs%coordinateResolution)
49  deallocate(cs)

Referenced by mom_regridding::end_regridding().

Here is the caller graph for this function:

◆ init_coord_zlike()

subroutine, public coord_zlike::init_coord_zlike ( type(zlike_cs), pointer  CS,
integer, intent(in)  nk,
real, dimension(:), intent(in)  coordinateResolution 
)

Initialise a zlike_CS with pointers to parameters.

Parameters
csUnassociated pointer to hold the control structure
[in]nkNumber of levels in the grid
[in]coordinateresolutionTarget coordinate resolution [Z ~> m]

Definition at line 30 of file coord_zlike.F90.

30  type(zlike_CS), pointer :: CS !< Unassociated pointer to hold the control structure
31  integer, intent(in) :: nk !< Number of levels in the grid
32  real, dimension(:), intent(in) :: coordinateResolution !< Target coordinate resolution [Z ~> m]
33 
34  if (associated(cs)) call mom_error(fatal, "init_coord_zlike: CS already associated!")
35  allocate(cs)
36  allocate(cs%coordinateResolution(nk))
37 
38  cs%nk = nk
39  cs%coordinateResolution = coordinateresolution

References mom_error_handler::mom_error().

Referenced by mom_regridding::initcoord().

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

◆ set_zlike_params()

subroutine, public coord_zlike::set_zlike_params ( type(zlike_cs), pointer  CS,
real, intent(in), optional  min_thickness 
)

Set parameters in the zlike structure.

Parameters
csCoordinate control structure
[in]min_thicknessMinimum allowed thickness [H ~> m or kg m-2]

Definition at line 54 of file coord_zlike.F90.

54  type(zlike_CS), pointer :: CS !< Coordinate control structure
55  real, optional, intent(in) :: min_thickness !< Minimum allowed thickness [H ~> m or kg m-2]
56 
57  if (.not. associated(cs)) call mom_error(fatal, "set_zlike_params: CS not associated")
58 
59  if (present(min_thickness)) cs%min_thickness = min_thickness

References mom_error_handler::mom_error().

Referenced by mom_regridding::set_regrid_params().

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