Skip to content

Commit

Permalink
Change paramter name and though up PR
Browse files Browse the repository at this point in the history
Change "MLD_PHA_FIXED" to "MLD_PHA_CALC" which is default false.
Touch up a few formatting issues.
Add initialization and comment for mld_pha.
  • Loading branch information
Theresa Morrison authored and Theresa Morrison committed Aug 19, 2024
1 parent 735cc9d commit 0a97d48
Show file tree
Hide file tree
Showing 2 changed files with 16 additions and 15 deletions.
11 changes: 5 additions & 6 deletions src/diagnostics/MOM_diagnose_MLD.F90
Original file line number Diff line number Diff line change
Expand Up @@ -30,8 +30,7 @@ module MOM_diagnose_mld
!> Diagnose a mixed layer depth (MLD) determined by a given density difference with the surface.
!> This routine is appropriate in MOM_diabatic_aux due to its position within the time stepping.
subroutine diagnoseMLDbyDensityDifference(id_MLD, h, tv, densityDiff, G, GV, US, diagPtr, &
ref_h_mld, id_ref_z, id_ref_rho, id_N2subML, id_MLDsq, dz_subML, &
MLD_out)
ref_h_mld, id_ref_z, id_ref_rho, id_N2subML, id_MLDsq, dz_subML, MLD_out)
type(ocean_grid_type), intent(in) :: G !< Grid type
type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
Expand All @@ -46,7 +45,7 @@ subroutine diagnoseMLDbyDensityDifference(id_MLD, h, tv, densityDiff, G, GV, US,
integer, intent(in) :: id_ref_z !< Handle (ID) of reference depth diagnostic
integer, intent(in) :: id_ref_rho !< Handle (ID) of reference density diagnostic
real, dimension(SZI_(G),SZJ_(G)), &
optional, intent(inout) :: MLD_out !< Send MLD to other routines [Z ~> m]
optional, intent(inout) :: MLD_out !< Send MLD to other routines [Z ~> m]
integer, optional, intent(in) :: id_N2subML !< Optional handle (ID) of subML stratification
integer, optional, intent(in) :: id_MLDsq !< Optional handle (ID) of squared MLD
real, optional, intent(in) :: dz_subML !< The distance over which to calculate N2subML
Expand Down Expand Up @@ -237,7 +236,7 @@ subroutine diagnoseMLDbyDensityDifference(id_MLD, h, tv, densityDiff, G, GV, US,
if ((id_ref_z > 0) .and. (pRef_MLD(is)/=0.)) call post_data(id_ref_z, z_ref_diag , diagPtr)
if (id_ref_rho > 0) call post_data(id_ref_rho, rhoSurf_2d , diagPtr)

if (present(MLD_out)) MLD_out(:,:)=MLD(:,:)
if (present(MLD_out)) MLD_out(:,:) = MLD(:,:)

end subroutine diagnoseMLDbyDensityDifference

Expand Down Expand Up @@ -276,7 +275,7 @@ subroutine diagnoseMLDbyEnergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr,
!! available thermodynamic fields.
type(diag_ctrl), pointer :: diagPtr !< Diagnostics structure
real, dimension(SZI_(G),SZJ_(G)), &
optional, intent(inout) :: MLD_out !< Send MLD to other routines [Z ~> m]
optional, intent(inout) :: MLD_out !< Send MLD to other routines [Z ~> m]

! Local variables
real, dimension(SZI_(G),SZJ_(G),3) :: MLD ! Diagnosed mixed layer depth [Z ~> m].
Expand Down Expand Up @@ -474,7 +473,7 @@ subroutine diagnoseMLDbyEnergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr,
if (id_MLD(2) > 0) call post_data(id_MLD(2), MLD(:,:,2), diagPtr)
if (id_MLD(3) > 0) call post_data(id_MLD(3), MLD(:,:,3), diagPtr)

if (present(MLD_out)) MLD_out(:,:)=MLD(:,:,1)
if (present(MLD_out)) MLD_out(:,:) = MLD(:,:,1)

end subroutine diagnoseMLDbyEnergy

Expand Down
20 changes: 11 additions & 9 deletions src/tracer/MOM_generic_tracer.F90
Original file line number Diff line number Diff line change
Expand Up @@ -84,12 +84,11 @@ module MOM_generic_tracer
logical :: tracers_may_reinit !< If true, tracers may go through the
!! initialization code if they are not found in the restart files.

!real, allocatable, dimension(:,:) :: mld_pha
logical :: mld_pha_fixed = .True. !< If true, use a fixed value for photoacclimation MLD
logical :: mld_pha_calc = .False. !< If true, use a fixed value for photoacclimation MLD
real :: mld_pha_val = 0.0 !< The value of fixed photoacclimation MLD
logical :: mld_pha_use_delta_rho = .False. !< If true, use a density diference to find the MLD
real :: mld_pha_href = 0.0 !< The reference depth for density difference based MLD
real :: mld_pha_drho = 0.03 !< The density thershold for a density difference based MLD
real :: mld_pha_drho = 0.03 !< The density thershold for a density difference based MLD
logical :: mld_pha_use_delta_eng = .False. !< If true, use an energy diference to find the MLD
real :: mld_pha_deng = 25.0 !< The energy threshold for an energy d ifference based MLD

Expand Down Expand Up @@ -433,10 +432,11 @@ subroutine initialize_MOM_generic_tracer(restart, day, G, GV, US, h, tv, param_f
enddo
!! end section to re-initialize generic tracers

call get_param(param_file, "MOM", "PHA_MLD_FIXED", CS%mld_pha_fixed, &
"If true, use a fixed value for the photoacclimation mixed layer depth within the "//&
"generic tracer update. This MLD is only used for photoacclimation.", default=.true.)
if (CS%mld_pha_fixed) then
call get_param(param_file, "MOM", "PHA_MLD_CALC", CS%mld_pha_calc, &
"If false, use a fixed value for the photoacclimation mixed layer depth within the "//&
"generic tracer update. This MLD is only used for photoacclimation. This variable should "//&
"be set to true if using COBALTv3 for the BGC.", default=.false.)
if (.not.CS%mld_pha_calc) then
call get_param(param_file, "MOM", "PHA_MLD_VAL", CS%mld_pha_val, &
"The depth of photoacclimation if fixed depth is used [m].", &
units='m', default=0.0, scale=US%m_to_Z)
Expand Down Expand Up @@ -544,7 +544,8 @@ subroutine MOM_generic_tracer_column_physics(h_old, h_new, ea, eb, fluxes, Hml,
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: rho_dzt ! Layer mass per unit area [kg m-2]
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: dzt ! Layer vertical extents [m]
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h_work ! A work array of thicknesses [H ~> m or kg m-2]
real, dimension(SZI_(G),SZJ_(G)) :: mld_pha !
real, dimension(SZI_(G),SZJ_(G)) :: mld_pha ! The mixed layer depth calculated for photoacclimation
! that is used in COBALTv3
integer :: i, j, k, isc, iec, jsc, jec, nk

isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; nk = GV%ke
Expand Down Expand Up @@ -610,7 +611,8 @@ subroutine MOM_generic_tracer_column_physics(h_old, h_new, ea, eb, fluxes, Hml,
enddo ; enddo
sosga = global_area_mean(surface_field, G, unscale=US%S_to_ppt)

if (CS%mld_pha_fixed) then
mld_pha(:,:) = 0.0
if (.not.CS%mld_pha_calc) then
mld_pha(:,:) = CS%mld_pha_val
else
if (CS%mld_pha_use_delta_rho) then
Expand Down

0 comments on commit 0a97d48

Please sign in to comment.