Skip to content

Commit 1344e7c

Browse files
theresa-corderoTheresa MorrisonTheresa Morrison
authored
Add MLD_out to diagnose MLD (#832)
* Add MLD_out to diagnose MLD call Add an optional argument to pass the MLD out of the diagnose MLD routine. * Change MLD_out to be only computaional domain Includes three changes: - Make MLD_out optional - initalize MLD_out to 0. everywhere - copy MLD to MLD_out only in the computational domain. * Remove unnedded MLD=0 initalization Remove an extra initalization of MLD in diagnoseMLDbyEnergy --------- Co-authored-by: Theresa Morrison <Theresa.Morrison@gaea68.ncrc.gov> Co-authored-by: Theresa Morrison <Theresa.Morrison@gaea63.ncrc.gov>
1 parent 51b1515 commit 1344e7c

File tree

1 file changed

+17
-4
lines changed

1 file changed

+17
-4
lines changed

src/diagnostics/MOM_diagnose_MLD.F90

Lines changed: 17 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -30,7 +30,8 @@ module MOM_diagnose_mld
3030
!> Diagnose a mixed layer depth (MLD) determined by a given density difference with the surface.
3131
!> This routine is appropriate in MOM_diabatic_aux due to its position within the time stepping.
3232
subroutine diagnoseMLDbyDensityDifference(id_MLD, h, tv, densityDiff, G, GV, US, diagPtr, &
33-
ref_h_mld, id_ref_z, id_ref_rho, id_N2subML, id_MLDsq, dz_subML)
33+
ref_h_mld, id_ref_z, id_ref_rho, id_N2subML, id_MLDsq, &
34+
dz_subML, MLD_out)
3435
type(ocean_grid_type), intent(in) :: G !< Grid type
3536
type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
3637
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
@@ -48,6 +49,8 @@ subroutine diagnoseMLDbyDensityDifference(id_MLD, h, tv, densityDiff, G, GV, US,
4849
integer, optional, intent(in) :: id_MLDsq !< Optional handle (ID) of squared MLD
4950
real, optional, intent(in) :: dz_subML !< The distance over which to calculate N2subML
5051
!! or 50 m if missing [Z ~> m]
52+
real, dimension(SZI_(G),SZJ_(G)), &
53+
optional, intent(out) :: MLD_out !< Send MLD to other routines [Z ~> m]
5154

5255
! Local variables
5356
real, dimension(SZI_(G)) :: deltaRhoAtKm1, deltaRhoAtK ! Density differences [R ~> kg m-3].
@@ -234,11 +237,16 @@ subroutine diagnoseMLDbyDensityDifference(id_MLD, h, tv, densityDiff, G, GV, US,
234237
if ((id_ref_z > 0) .and. (pRef_MLD(is)/=0.)) call post_data(id_ref_z, z_ref_diag , diagPtr)
235238
if (id_ref_rho > 0) call post_data(id_ref_rho, rhoSurf_2d , diagPtr)
236239

240+
if (present(MLD_out)) then
241+
MLD_out(:,:) = 0.0
242+
MLD_out(is:ie,js:je) = MLD(is:ie,js:je)
243+
endif
244+
237245
end subroutine diagnoseMLDbyDensityDifference
238246

239247
!> Diagnose a mixed layer depth (MLD) determined by the depth a given energy value would mix.
240248
!> This routine is appropriate in MOM_diabatic_aux due to its position within the time stepping.
241-
subroutine diagnoseMLDbyEnergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr)
249+
subroutine diagnoseMLDbyEnergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr, MLD_out)
242250
! Author: Brandon Reichl
243251
! Date: October 2, 2020
244252
! //
@@ -270,6 +278,8 @@ subroutine diagnoseMLDbyEnergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr)
270278
type(thermo_var_ptrs), intent(in) :: tv !< Structure containing pointers to any
271279
!! available thermodynamic fields.
272280
type(diag_ctrl), pointer :: diagPtr !< Diagnostics structure
281+
real, dimension(SZI_(G),SZJ_(G)), &
282+
optional, intent(out) :: MLD_out !< Send MLD to other routines [Z ~> m]
273283

274284
! Local variables
275285
real, dimension(SZI_(G),SZJ_(G),3) :: MLD ! Diagnosed mixed layer depth [Z ~> m].
@@ -325,8 +335,6 @@ subroutine diagnoseMLDbyEnergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr)
325335
PE_threshold(iM) = Mixing_Energy(iM) / GV%g_Earth_Z_T2
326336
enddo
327337

328-
MLD(:,:,:) = 0.0
329-
330338
EOSdom(:) = EOS_domain(G%HI)
331339

332340
do j=js,je
@@ -467,6 +475,11 @@ subroutine diagnoseMLDbyEnergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr)
467475
if (id_MLD(2) > 0) call post_data(id_MLD(2), MLD(:,:,2), diagPtr)
468476
if (id_MLD(3) > 0) call post_data(id_MLD(3), MLD(:,:,3), diagPtr)
469477

478+
if (present(MLD_out)) then
479+
MLD_out(:,:) = 0.0
480+
MLD_out(is:ie,js:je) = MLD(is:ie,js:je,1)
481+
endif
482+
470483
end subroutine diagnoseMLDbyEnergy
471484

472485
!> \namespace mom_diagnose_mld

0 commit comments

Comments
 (0)