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)
 
subroutine, public set_zlike_params (CS, min_thickness)
 
subroutine, public build_zstar_column (CS, nz, depth, total_thickness, zInterface, z_rigid_top, eta_orig)
 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,
integer, intent(in)  nz,
real, intent(in)  depth,
real, intent(in)  total_thickness,
real, dimension(nz+1), intent(inout)  zInterface,
real, intent(in), optional  z_rigid_top,
real, intent(in), optional  eta_orig 
)

Builds a z* coordinate with a minimum thickness.

Parameters
[in]csCoordinate control structure
[in]nzNumber of levels
[in]depthDepth of ocean bottom (positive in m)
[in]total_thicknessColumn thickness (positive in m)
[in,out]zinterfaceAbsolute positions of interfaces
[in]z_rigid_topThe height of a rigid top (negative in m)
[in]eta_origThe actual original height of the top (m)

Definition at line 60 of file coord_zlike.F90.

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

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

◆ end_coord_zlike()

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

Definition at line 41 of file coord_zlike.F90.

Referenced by mom_regridding::end_regridding().

41  type(zlike_cs), pointer :: cs
42 
43  ! nothing to do
44  if (.not. associated(cs)) return
45  deallocate(cs%coordinateResolution)
46  deallocate(cs)
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

Definition at line 28 of file coord_zlike.F90.

References mom_error_handler::mom_error().

28  type(zlike_cs), pointer :: cs !< Unassociated pointer to hold the control structure
29  integer, intent(in) :: nk
30  real, dimension(:), intent(in) :: coordinateresolution
31 
32  if (associated(cs)) call mom_error(fatal, "init_coord_zlike: CS already associated!")
33  allocate(cs)
34  allocate(cs%coordinateResolution(nk))
35 
36  cs%nk = nk
37  cs%coordinateResolution = coordinateresolution
Here is the call 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 
)

Definition at line 50 of file coord_zlike.F90.

References mom_error_handler::mom_error().

50  type(zlike_cs), pointer :: cs
51  real, optional, intent(in) :: min_thickness
52 
53  if (.not. associated(cs)) call mom_error(fatal, "set_zlike_params: CS not associated")
54 
55  if (present(min_thickness)) cs%min_thickness = min_thickness
Here is the call graph for this function: