14 implicit none ;
private 16 #include <MOM_memory.h> 22 module procedure find_eta_2d, find_eta_3d
31 subroutine find_eta_3d(h, tv, G, GV, US, eta, eta_bt, halo_size, eta_to_m)
32 type(ocean_grid_type),
intent(in) :: G
33 type(verticalGrid_type),
intent(in) :: GV
34 type(unit_scale_type),
intent(in) :: US
35 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
36 type(thermo_var_ptrs),
intent(in) :: tv
38 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(out) :: eta
40 real,
dimension(SZI_(G),SZJ_(G)),
optional,
intent(in) :: eta_bt
44 integer,
optional,
intent(in) :: halo_size
46 real,
optional,
intent(in) :: eta_to_m
50 real :: p(SZI_(G),SZJ_(G),SZK_(G)+1)
51 real :: dz_geo(SZI_(G),SZJ_(G),SZK_(G))
53 real :: dilate(SZI_(G))
57 real :: Z_to_eta, H_to_eta, H_to_rho_eta
58 integer i, j, k, isv, iev, jsv, jev, nz, halo
60 halo = 0 ;
if (
present(halo_size)) halo = max(0,halo_size)
62 isv = g%isc-halo ; iev = g%iec+halo ; jsv = g%jsc-halo ; jev = g%jec+halo
65 if ((isv<g%isd) .or. (iev>g%ied) .or. (jsv<g%jsd) .or. (jev>g%jed)) &
66 call mom_error(fatal,
"find_eta called with an overly large halo_size.")
68 z_to_eta = 1.0 ;
if (
present(eta_to_m)) z_to_eta = us%Z_to_m / eta_to_m
69 h_to_eta = gv%H_to_Z * z_to_eta
70 h_to_rho_eta = gv%H_to_RZ * z_to_eta
71 i_gearth = z_to_eta / gv%g_Earth
75 do j=jsv,jev ;
do i=isv,iev ; eta(i,j,nz+1) = -z_to_eta*g%bathyT(i,j) ;
enddo ;
enddo 77 if (gv%Boussinesq)
then 79 do j=jsv,jev ;
do k=nz,1,-1;
do i=isv,iev
80 eta(i,j,k) = eta(i,j,k+1) + h(i,j,k)*h_to_eta
82 if (
present(eta_bt))
then 88 dilate(i) = (eta_bt(i,j)*h_to_eta + z_to_eta*g%bathyT(i,j)) / &
89 (eta(i,j,1) + z_to_eta*g%bathyT(i,j))
91 do k=1,nz ;
do i=isv,iev
92 eta(i,j,k) = dilate(i) * (eta(i,j,k) + z_to_eta*g%bathyT(i,j)) - z_to_eta*g%bathyT(i,j)
97 if (
associated(tv%eqn_of_state))
then 100 if (
associated(tv%p_surf))
then 101 do i=isv,iev ; p(i,j,1) = tv%p_surf(i,j) ;
enddo 103 do i=isv,iev ; p(i,j,1) = 0.0 ;
enddo 105 do k=1,nz ;
do i=isv,iev
106 p(i,j,k+1) = p(i,j,k) + gv%g_Earth*gv%H_to_RZ*h(i,j,k)
111 call int_specific_vol_dp(tv%T(:,:,k), tv%S(:,:,k), p(:,:,k), p(:,:,k+1), &
112 0.0, g%HI, tv%eqn_of_state, us, dz_geo(:,:,k), halo_size=halo)
116 do k=nz,1,-1 ;
do i=isv,iev
117 eta(i,j,k) = eta(i,j,k+1) + i_gearth * dz_geo(i,j,k)
122 do j=jsv,jev ;
do k=nz,1,-1;
do i=isv,iev
123 eta(i,j,k) = eta(i,j,k+1) + h_to_rho_eta*h(i,j,k) / gv%Rlay(k)
124 enddo ;
enddo ;
enddo 126 if (
present(eta_bt))
then 131 do i=isv,iev ; htot(i) = gv%H_subroundoff ;
enddo 132 do k=1,nz ;
do i=isv,iev ; htot(i) = htot(i) + h(i,j,k) ;
enddo ;
enddo 133 do i=isv,iev ; dilate(i) = eta_bt(i,j) / htot(i) ;
enddo 134 do k=1,nz ;
do i=isv,iev
135 eta(i,j,k) = dilate(i) * (eta(i,j,k) + z_to_eta*g%bathyT(i,j)) - z_to_eta*g%bathyT(i,j)
142 end subroutine find_eta_3d
148 subroutine find_eta_2d(h, tv, G, GV, US, eta, eta_bt, halo_size, eta_to_m)
149 type(ocean_grid_type),
intent(in) :: G
150 type(verticalGrid_type),
intent(in) :: GV
151 type(unit_scale_type),
intent(in) :: US
152 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
153 type(thermo_var_ptrs),
intent(in) :: tv
155 real,
dimension(SZI_(G),SZJ_(G)),
intent(out) :: eta
157 real,
dimension(SZI_(G),SZJ_(G)),
optional,
intent(in) :: eta_bt
160 integer,
optional,
intent(in) :: halo_size
162 real,
optional,
intent(in) :: eta_to_m
165 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
167 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
169 real :: htot(SZI_(G))
172 real :: Z_to_eta, H_to_eta, H_to_rho_eta
173 integer i, j, k, is, ie, js, je, nz, halo
175 halo = 0 ;
if (
present(halo_size)) halo = max(0,halo_size)
176 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
179 z_to_eta = 1.0 ;
if (
present(eta_to_m)) z_to_eta = us%Z_to_m / eta_to_m
180 h_to_eta = gv%H_to_Z * z_to_eta
181 h_to_rho_eta = gv%H_to_RZ * z_to_eta
182 i_gearth = z_to_eta / gv%g_Earth
186 do j=js,je ;
do i=is,ie ; eta(i,j) = -z_to_eta*g%bathyT(i,j) ;
enddo ;
enddo 188 if (gv%Boussinesq)
then 189 if (
present(eta_bt))
then 191 do j=js,je ;
do i=is,ie
192 eta(i,j) = h_to_eta*eta_bt(i,j)
196 do j=js,je ;
do k=1,nz ;
do i=is,ie
197 eta(i,j) = eta(i,j) + h(i,j,k)*h_to_eta
198 enddo ;
enddo ;
enddo 201 if (
associated(tv%eqn_of_state))
then 204 if (
associated(tv%p_surf))
then 205 do i=is,ie ; p(i,j,1) = tv%p_surf(i,j) ;
enddo 207 do i=is,ie ; p(i,j,1) = 0.0 ;
enddo 210 do k=1,nz ;
do i=is,ie
211 p(i,j,k+1) = p(i,j,k) + gv%g_Earth*gv%H_to_RZ*h(i,j,k)
216 call int_specific_vol_dp(tv%T(:,:,k), tv%S(:,:,k), p(:,:,k), p(:,:,k+1), 0.0, &
217 g%HI, tv%eqn_of_state, us, dz_geo(:,:,k), halo_size=halo)
220 do j=js,je ;
do k=1,nz ;
do i=is,ie
221 eta(i,j) = eta(i,j) + i_gearth * dz_geo(i,j,k)
222 enddo ;
enddo ;
enddo 225 do j=js,je ;
do k=1,nz ;
do i=is,ie
226 eta(i,j) = eta(i,j) + h_to_rho_eta*h(i,j,k) / gv%Rlay(k)
227 enddo ;
enddo ;
enddo 229 if (
present(eta_bt))
then 234 do i=is,ie ; htot(i) = gv%H_subroundoff ;
enddo 235 do k=1,nz ;
do i=is,ie ; htot(i) = htot(i) + h(i,j,k) ;
enddo ;
enddo 237 eta(i,j) = (eta_bt(i,j) / htot(i)) * (eta(i,j) + z_to_eta*g%bathyT(i,j)) - &
238 z_to_eta*g%bathyT(i,j)
245 end subroutine find_eta_2d
Calculates the heights of sruface or all interfaces from layer thicknesses.
Functions for calculating interface heights, including free surface height.
Ocean grid type. See mom_grid for details.
Provides the ocean grid type.
The MOM6 facility to parse input files for runtime parameters.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Provides integrals of density.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.