MOM6
mom_interface_heights::find_eta Interface Reference

Detailed Description

Calculates the heights of sruface or all interfaces from layer thicknesses.

Definition at line 21 of file MOM_interface_heights.F90.

Private functions

subroutine find_eta_2d (h, tv, G, GV, US, eta, eta_bt, halo_size, eta_to_m)
 Calculates the free surface height, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, the sea surface height may be adjusted for consistency with the corresponding time-average quantity from the barotropic calculation. More...
 
subroutine find_eta_3d (h, tv, G, GV, US, eta, eta_bt, halo_size, eta_to_m)
 Calculates the heights of all interfaces between layers, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, these height may be dilated for consistency with the corresponding time-average quantity from the barotropic calculation. More...
 

Functions and subroutines

◆ find_eta_2d()

subroutine mom_interface_heights::find_eta::find_eta_2d ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  eta,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in), optional  eta_bt,
integer, intent(in), optional  halo_size,
real, intent(in), optional  eta_to_m 
)
private

Calculates the free surface height, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, the sea surface height may be adjusted for consistency with the corresponding time-average quantity from the barotropic calculation.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]tvA structure pointing to various thermodynamic variables.
[out]etafree surface height relative to mean sea level (z=0) often [Z ~> m].
[in]eta_btoptional barotropic variable that gives the "correct" free surface height (Boussinesq) or total water column mass per unit area (non-Boussinesq) [H ~> m or kg m-2].
[in]halo_sizewidth of halo points on which to calculate eta.
[in]eta_to_mThe conversion factor from the units of eta to m; by default this is USZ_to_m.

Definition at line 145 of file MOM_interface_heights.F90.

145  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
146  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
147  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
148  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
149  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
150  !! thermodynamic variables.
151  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: eta !< free surface height relative to
152  !! mean sea level (z=0) often [Z ~> m].
153  real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: eta_bt !< optional barotropic
154  !! variable that gives the "correct" free surface height (Boussinesq) or total
155  !! water column mass per unit area (non-Boussinesq) [H ~> m or kg m-2].
156  integer, optional, intent(in) :: halo_size !< width of halo points on
157  !! which to calculate eta.
158  real, optional, intent(in) :: eta_to_m !< The conversion factor from
159  !! the units of eta to m; by default this is US%Z_to_m.
160  ! Local variables
161  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
162  p ! The pressure at interfaces [Pa].
163  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
164  dz_geo ! The change in geopotential height across a layer [m2 s-2].
165  real :: htot(SZI_(G)) ! The sum of all layers' thicknesses [H ~> m or kg m-2].
166  real :: I_gEarth
167  real :: Z_to_eta, H_to_eta, H_to_rho_eta ! Unit conversion factors with obvious names.
168  integer i, j, k, is, ie, js, je, nz, halo
169 
170  halo = 0 ; if (present(halo_size)) halo = max(0,halo_size)
171  is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
172  nz = g%ke
173 
174  z_to_eta = 1.0 ; if (present(eta_to_m)) z_to_eta = us%Z_to_m / eta_to_m
175  h_to_eta = gv%H_to_Z * z_to_eta
176  h_to_rho_eta = gv%H_to_RZ * z_to_eta
177  i_gearth = z_to_eta / (us%Z_to_m * gv%mks_g_Earth)
178 
179 !$OMP parallel default(shared) private(htot)
180 !$OMP do
181  do j=js,je ; do i=is,ie ; eta(i,j) = -z_to_eta*g%bathyT(i,j) ; enddo ; enddo
182 
183  if (gv%Boussinesq) then
184  if (present(eta_bt)) then
185 !$OMP do
186  do j=js,je ; do i=is,ie
187  eta(i,j) = h_to_eta*eta_bt(i,j)
188  enddo ; enddo
189  else
190 !$OMP do
191  do j=js,je ; do k=1,nz ; do i=is,ie
192  eta(i,j) = eta(i,j) + h(i,j,k)*h_to_eta
193  enddo ; enddo ; enddo
194  endif
195  else
196  if (associated(tv%eqn_of_state)) then
197 !$OMP do
198  do j=js,je
199  do i=is,ie ; p(i,j,1) = 0.0 ; enddo
200 
201  do k=1,nz ; do i=is,ie
202  p(i,j,k+1) = p(i,j,k) + gv%H_to_Pa*h(i,j,k)
203  enddo ; enddo
204  enddo
205 !$OMP do
206  do k = 1, nz
207  call int_specific_vol_dp(tv%T(:,:,k), tv%S(:,:,k), p(:,:,k), p(:,:,k+1), 0.0, &
208  g%HI, tv%eqn_of_state, dz_geo(:,:,k), halo_size=halo)
209  enddo
210 !$OMP do
211  do j=js,je ; do k=1,nz ; do i=is,ie
212  eta(i,j) = eta(i,j) + i_gearth * dz_geo(i,j,k)
213  enddo ; enddo ; enddo
214  else
215 !$OMP do
216  do j=js,je ; do k=1,nz ; do i=is,ie
217  eta(i,j) = eta(i,j) + h_to_rho_eta*h(i,j,k) / gv%Rlay(k)
218  enddo ; enddo ; enddo
219  endif
220  if (present(eta_bt)) then
221  ! Dilate the water column to agree with the time-averaged column
222  ! mass from the barotropic solution.
223 !$OMP do
224  do j=js,je
225  do i=is,ie ; htot(i) = gv%H_subroundoff ; enddo
226  do k=1,nz ; do i=is,ie ; htot(i) = htot(i) + h(i,j,k) ; enddo ; enddo
227  do i=is,ie
228  eta(i,j) = (eta_bt(i,j) / htot(i)) * (eta(i,j) + z_to_eta*g%bathyT(i,j)) - &
229  z_to_eta*g%bathyT(i,j)
230  enddo
231  enddo
232  endif
233  endif
234 !$OMP end parallel
235 

◆ find_eta_3d()

subroutine mom_interface_heights::find_eta::find_eta_3d ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke+1), intent(out)  eta,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in), optional  eta_bt,
integer, intent(in), optional  halo_size,
real, intent(in), optional  eta_to_m 
)
private

Calculates the heights of all interfaces between layers, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, these height may be dilated for consistency with the corresponding time-average quantity from the barotropic calculation.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]tvA structure pointing to various thermodynamic variables.
[out]etalayer interface heights [Z ~> m] or 1/eta_to_m m).
[in]eta_btoptional barotropic variable that gives the "correct" free surface height (Boussinesq) or total water column mass per unit area (non-Boussinesq). This is used to dilate the layer. thicknesses when calculating interfaceheights [H ~> m or kg m-2].
[in]halo_sizewidth of halo points on which to calculate eta.
[in]eta_to_mThe conversion factor from the units of eta to m; by default this is USZ_to_m.

Definition at line 32 of file MOM_interface_heights.F90.

32  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
33  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
34  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
35  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
36  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
37  !! thermodynamic variables.
38  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(out) :: eta !< layer interface heights
39  !! [Z ~> m] or 1/eta_to_m m).
40  real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: eta_bt !< optional barotropic
41  !! variable that gives the "correct" free surface height (Boussinesq) or total water
42  !! column mass per unit area (non-Boussinesq). This is used to dilate the layer.
43  !! thicknesses when calculating interfaceheights [H ~> m or kg m-2].
44  integer, optional, intent(in) :: halo_size !< width of halo points on
45  !! which to calculate eta.
46  real, optional, intent(in) :: eta_to_m !< The conversion factor from
47  !! the units of eta to m; by default this is US%Z_to_m.
48 
49  ! Local variables
50  real :: p(SZI_(G),SZJ_(G),SZK_(G)+1)
51  real :: dz_geo(SZI_(G),SZJ_(G),SZK_(G)) ! The change in geopotential height
52  ! across a layer [m2 s-2].
53  real :: dilate(SZI_(G)) ! non-dimensional dilation factor
54  real :: htot(SZI_(G)) ! total thickness H
55  real :: I_gEarth
56  real :: Z_to_eta, H_to_eta, H_to_rho_eta ! Unit conversion factors with obvious names.
57  integer i, j, k, isv, iev, jsv, jev, nz, halo
58 
59  halo = 0 ; if (present(halo_size)) halo = max(0,halo_size)
60 
61  isv = g%isc-halo ; iev = g%iec+halo ; jsv = g%jsc-halo ; jev = g%jec+halo
62  nz = g%ke
63 
64  if ((isv<g%isd) .or. (iev>g%ied) .or. (jsv<g%jsd) .or. (jev>g%jed)) &
65  call mom_error(fatal,"find_eta called with an overly large halo_size.")
66 
67  z_to_eta = 1.0 ; if (present(eta_to_m)) z_to_eta = us%Z_to_m / eta_to_m
68  h_to_eta = gv%H_to_Z * z_to_eta
69  h_to_rho_eta = gv%H_to_RZ * z_to_eta
70  i_gearth = z_to_eta / (us%Z_to_m * gv%mks_g_Earth)
71 
72 !$OMP parallel default(shared) private(dilate,htot)
73 !$OMP do
74  do j=jsv,jev ; do i=isv,iev ; eta(i,j,nz+1) = -z_to_eta*g%bathyT(i,j) ; enddo ; enddo
75 
76  if (gv%Boussinesq) then
77 !$OMP do
78  do j=jsv,jev ; do k=nz,1,-1; do i=isv,iev
79  eta(i,j,k) = eta(i,j,k+1) + h(i,j,k)*h_to_eta
80  enddo ; enddo ; enddo
81  if (present(eta_bt)) then
82  ! Dilate the water column to agree with the free surface height
83  ! that is used for the dynamics.
84 !$OMP do
85  do j=jsv,jev
86  do i=isv,iev
87  dilate(i) = (eta_bt(i,j)*h_to_eta + z_to_eta*g%bathyT(i,j)) / &
88  (eta(i,j,1) + z_to_eta*g%bathyT(i,j))
89  enddo
90  do k=1,nz ; do i=isv,iev
91  eta(i,j,k) = dilate(i) * (eta(i,j,k) + z_to_eta*g%bathyT(i,j)) - z_to_eta*g%bathyT(i,j)
92  enddo ; enddo
93  enddo
94  endif
95  else
96  if (associated(tv%eqn_of_state)) then
97 !$OMP do
98  do j=jsv,jev
99  ! ### THIS SHOULD BE P_SURF, IF AVAILABLE.
100  do i=isv,iev ; p(i,j,1) = 0.0 ; enddo
101  do k=1,nz ; do i=isv,iev
102  p(i,j,k+1) = p(i,j,k) + gv%H_to_Pa*h(i,j,k)
103  enddo ; enddo
104  enddo
105 !$OMP do
106  do k=1,nz
107  call int_specific_vol_dp(tv%T(:,:,k), tv%S(:,:,k), p(:,:,k), p(:,:,k+1), &
108  0.0, g%HI, tv%eqn_of_state, dz_geo(:,:,k), halo_size=halo)
109  enddo
110 !$OMP do
111  do j=jsv,jev
112  do k=nz,1,-1 ; do i=isv,iev
113  eta(i,j,k) = eta(i,j,k+1) + i_gearth * dz_geo(i,j,k)
114  enddo ; enddo
115  enddo
116  else
117 !$OMP do
118  do j=jsv,jev ; do k=nz,1,-1; do i=isv,iev
119  eta(i,j,k) = eta(i,j,k+1) + h_to_rho_eta*h(i,j,k) / gv%Rlay(k)
120  enddo ; enddo ; enddo
121  endif
122  if (present(eta_bt)) then
123  ! Dilate the water column to agree with the free surface height
124  ! from the time-averaged barotropic solution.
125 !$OMP do
126  do j=jsv,jev
127  do i=isv,iev ; htot(i) = gv%H_subroundoff ; enddo
128  do k=1,nz ; do i=isv,iev ; htot(i) = htot(i) + h(i,j,k) ; enddo ; enddo
129  do i=isv,iev ; dilate(i) = eta_bt(i,j) / htot(i) ; enddo
130  do k=1,nz ; do i=isv,iev
131  eta(i,j,k) = dilate(i) * (eta(i,j,k) + z_to_eta*g%bathyT(i,j)) - z_to_eta*g%bathyT(i,j)
132  enddo ; enddo
133  enddo
134  endif
135  endif
136 !$OMP end parallel
137 

The documentation for this interface was generated from the following file: