Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

dimension rescaling fixes and CFC cleanup #195

Merged
merged 6 commits into from
Aug 27, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 2 additions & 1 deletion config_src/drivers/mct_cap/mom_surface_forcing_mct.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1245,7 +1245,8 @@ subroutine surface_forcing_init(Time, G, US, param_file, diag, CS, restore_salt,
"If true, use a 2-dimensional gustiness supplied from "//&
"an input file", default=.false.)
call get_param(param_file, mdl, "GUST_CONST", CS%gust_const, &
"The background gustiness in the winds.", units="Pa", default=0.0)
"The background gustiness in the winds.", units="Pa", default=0.0, &
scale=US%kg_m3_to_R*US%m_s_to_L_T**2*US%L_to_Z)
if (CS%read_gust_2d) then
call get_param(param_file, mdl, "GUST_2D_FILE", gust_file, &
"The file in which the wind gustiness is found in "//&
Expand Down
2 changes: 1 addition & 1 deletion config_src/drivers/nuopc_cap/mom_cap.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1113,7 +1113,7 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
k = k + 1 ! Increment position within gindex
if (mask(k) /= 0) then
mesh_areas(k) = dataPtr_mesh_areas(k)
model_areas(k) = ocean_grid%AreaT(i,j) / ocean_grid%Rad_Earth**2
model_areas(k) = ocean_grid%US%L_to_m**2 * ocean_grid%AreaT(i,j) / ocean_grid%Rad_Earth**2
mod2med_areacor(k) = model_areas(k) / mesh_areas(k)
med2mod_areacor(k) = mesh_areas(k) / model_areas(k)
end if
Expand Down
6 changes: 4 additions & 2 deletions config_src/drivers/nuopc_cap/mom_surface_forcing_nuopc.F90
Original file line number Diff line number Diff line change
Expand Up @@ -576,7 +576,7 @@ subroutine convert_IOB_to_fluxes(IOB, fluxes, index_bounds, Time, valid_time, G,

! CFCs
if (CS%use_CFC) then
call CFC_cap_fluxes(fluxes, sfc_state, G, CS%Rho0, Time, CS%id_cfc11_atm, CS%id_cfc11_atm)
call CFC_cap_fluxes(fluxes, sfc_state, G, US, CS%Rho0, Time, CS%id_cfc11_atm, CS%id_cfc11_atm)
endif

if (associated(IOB%salt_flux)) then
Expand Down Expand Up @@ -1396,7 +1396,9 @@ subroutine surface_forcing_init(Time, G, US, param_file, diag, CS, restore_salt,
"internal BC generation (TODO).", default=" ", do_not_log=.true.)
if ((len_trim(CS%CFC_BC_file) > 0) .and. (scan(CS%CFC_BC_file,'/') == 0)) then
! Add the directory if CFC_BC_file is not already a complete path.
CS%CFC_BC_file = trim(slasher(CS%inputdir))//trim(CS%CFC_BC_file)
CS%CFC_BC_file = trim(CS%inputdir) // trim(CS%CFC_BC_file)
endif
if (len_trim(CS%CFC_BC_file) > 0) then
call get_param(param_file, mdl, "CFC11_VARIABLE", CS%cfc11_var_name, &
"The name of the variable representing CFC-11 in "//&
"CFC_BC_FILE.", default="CFC_11", do_not_log=.true.)
Expand Down
3 changes: 2 additions & 1 deletion config_src/drivers/solo_driver/MESO_surface_forcing.F90
Original file line number Diff line number Diff line change
Expand Up @@ -243,7 +243,8 @@ subroutine MESO_surface_forcing_init(Time, G, US, param_file, diag, CS)
"parameters from vertical units of m to kg m-2.", &
units="kg m-3", default=1035.0, scale=US%kg_m3_to_R)
call get_param(param_file, mdl, "GUST_CONST", CS%gust_const, &
"The background gustiness in the winds.", units="Pa", default=0.0)
"The background gustiness in the winds.", units="Pa", default=0.0, &
scale=US%kg_m3_to_R*US%m_s_to_L_T**2*US%L_to_Z)

call get_param(param_file, mdl, "RESTOREBUOY", CS%restorebuoy, &
"If true, the buoyancy fluxes drive the model back "//&
Expand Down
16 changes: 10 additions & 6 deletions src/core/MOM_forcing_type.F90
Original file line number Diff line number Diff line change
Expand Up @@ -185,8 +185,8 @@ module MOM_forcing_type

! CFC-related arrays needed in the MOM_CFC_cap module
real, pointer, dimension(:,:) :: &
cfc11_flux => NULL(), & !< flux of cfc_11 into the ocean [mol m-2 s-1].
cfc12_flux => NULL(), & !< flux of cfc_12 into the ocean [mol m-2 s-1].
cfc11_flux => NULL(), & !< flux of cfc_11 into the ocean [CU Z T-1 kg m-3 = mol Z T-1 m-3 ~> mol m-2 s-1].
cfc12_flux => NULL(), & !< flux of cfc_12 into the ocean [CU Z T-1 kg m-3 = mol Z T-1 m-3 ~> mol m-2 s-1].
ice_fraction => NULL(), & !< fraction of sea ice coverage at h-cells, from 0 to 1 [nondim].
u10_sqr => NULL() !< wind magnitude at 10 m squared [L2 T-2 ~> m2 s-2]

Expand Down Expand Up @@ -1101,9 +1101,13 @@ subroutine MOM_forcing_chksum(mesg, fluxes, G, US, haloshift)
if (associated(fluxes%p_surf)) &
call hchksum(fluxes%p_surf, mesg//" fluxes%p_surf", G%HI, haloshift=hshift , scale=US%RL2_T2_to_Pa)
if (associated(fluxes%u10_sqr)) &
call hchksum(fluxes%u10_sqr, mesg//" fluxes%u10_sqr", G%HI, haloshift=hshift , scale=US%Z_to_m**2*US%s_to_T**2)
call hchksum(fluxes%u10_sqr, mesg//" fluxes%u10_sqr", G%HI, haloshift=hshift , scale=US%L_to_m**2*US%s_to_T**2)
if (associated(fluxes%ice_fraction)) &
call hchksum(fluxes%ice_fraction, mesg//" fluxes%ice_fraction", G%HI, haloshift=hshift)
if (associated(fluxes%cfc11_flux)) &
call hchksum(fluxes%cfc11_flux, mesg//" fluxes%cfc11_flux", G%HI, haloshift=hshift, scale=US%Z_to_m*US%s_to_T)
if (associated(fluxes%cfc12_flux)) &
call hchksum(fluxes%cfc12_flux, mesg//" fluxes%cfc12_flux", G%HI, haloshift=hshift, scale=US%Z_to_m*US%s_to_T)
if (associated(fluxes%salt_flux)) &
call hchksum(fluxes%salt_flux, mesg//" fluxes%salt_flux", G%HI, haloshift=hshift, scale=US%RZ_T_to_kg_m2s)
if (associated(fluxes%TKE_tidal)) &
Expand All @@ -1117,7 +1121,7 @@ subroutine MOM_forcing_chksum(mesg, fluxes, G, US, haloshift)
call hchksum(fluxes%frunoff, mesg//" fluxes%frunoff", G%HI, haloshift=hshift, scale=US%RZ_T_to_kg_m2s)
if (associated(fluxes%heat_content_lrunoff)) &
call hchksum(fluxes%heat_content_lrunoff, mesg//" fluxes%heat_content_lrunoff", G%HI, &
haloshift=hshift, scale=US%RZ_T_to_kg_m2s)
haloshift=hshift, scale=US%QRZ_T_to_W_m2)
if (associated(fluxes%heat_content_frunoff)) &
call hchksum(fluxes%heat_content_frunoff, mesg//" fluxes%heat_content_frunoff", G%HI, &
haloshift=hshift, scale=US%QRZ_T_to_W_m2)
Expand Down Expand Up @@ -1314,14 +1318,14 @@ subroutine register_forcing_type_diags(Time, diag, US, use_temperature, handles,
if (use_cfcs) then
handles%id_cfc11 = register_diag_field('ocean_model', 'cfc11_flux', diag%axesT1, Time, &
'Gas exchange flux of CFC11 into the ocean ', 'mol m-2 s-1', &
conversion= US%m_to_Z**2*US%s_to_T,&
conversion= US%Z_to_m*US%s_to_T,&
cmor_field_name='fgcfc11', &
cmor_long_name='Surface Downward CFC11 Flux', &
cmor_standard_name='surface_downward_cfc11_flux')

handles%id_cfc12 = register_diag_field('ocean_model', 'cfc12_flux', diag%axesT1, Time, &
'Gas exchange flux of CFC12 into the ocean ', 'mol m-2 s-1', &
conversion= US%m_to_Z**2*US%s_to_T,&
conversion= US%Z_to_m*US%s_to_T,&
cmor_field_name='fgcfc12', &
cmor_long_name='Surface Downward CFC12 Flux', &
cmor_standard_name='surface_downward_cfc12_flux')
Expand Down
18 changes: 13 additions & 5 deletions src/parameterizations/lateral/MOM_MEKE.F90
Original file line number Diff line number Diff line change
Expand Up @@ -81,7 +81,7 @@ module MOM_MEKE
real :: MEKE_advection_factor !< A scaling in front of the advection of MEKE [nondim]
real :: MEKE_topographic_beta !< Weight for how much topographic beta is considered
!! when computing beta in Rhines scale [nondim]
real :: MEKE_restoring_rate !< Inverse of the timescale used to nudge MEKE toward its equilibrium value [s-1].
real :: MEKE_restoring_rate !< Inverse of the timescale used to nudge MEKE toward its equilibrium value [T-1 ~> s-1].

logical :: kh_flux_enabled !< If true, lateral diffusive MEKE flux is enabled.
logical :: initialize !< If True, invokes a steady state solver to calculate MEKE.
Expand Down Expand Up @@ -340,6 +340,10 @@ subroutine step_forward_MEKE(MEKE, h, SN_u, SN_v, visc, dt, G, GV, US, CS, hu, h
enddo ; enddo
endif

if (CS%debug) then
call hchksum(src, "MEKE src", G%HI, haloshift=0, scale=US%L_to_m**2*US%s_to_T**3)
endif

! Increase EKE by a full time-steps worth of source
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
Expand Down Expand Up @@ -534,6 +538,10 @@ subroutine step_forward_MEKE(MEKE, h, SN_u, SN_v, visc, dt, G, GV, US, CS, hu, h
endif
endif ! MEKE_KH>=0

if (CS%debug) then
call hchksum(MEKE%MEKE, "MEKE post-update MEKE", G%HI, haloshift=0, scale=US%L_T_to_m_s**2)
endif

! do j=js,je ; do i=is,ie
! MEKE%MEKE(i,j) = MAX(MEKE%MEKE(i,j),0.0)
! enddo ; enddo
Expand Down Expand Up @@ -688,7 +696,7 @@ subroutine MEKE_equilibrium(CS, MEKE, G, GV, US, SN_u, SN_v, drag_rate_visc, I_m
SN = min(SN_u(I,j), SN_u(I-1,j), SN_v(i,J), SN_v(i,J-1))

if (CS%MEKE_equilibrium_alt) then
MEKE%MEKE(i,j) = (CS%MEKE_GEOMETRIC_alpha * SN * US%Z_to_m*G%bathyT(i,j))**2 / cd2
MEKE%MEKE(i,j) = (CS%MEKE_GEOMETRIC_alpha * SN * US%Z_to_L*G%bathyT(i,j))**2 / cd2
else
FatH = 0.25*((G%CoriolisBu(I,J) + G%CoriolisBu(I-1,J-1)) + &
(G%CoriolisBu(I-1,J) + G%CoriolisBu(I,J-1))) ! Coriolis parameter at h points
Expand Down Expand Up @@ -824,7 +832,7 @@ subroutine MEKE_equilibrium_restoring(CS, G, US, SN_u, SN_v)
! SN = 0.25*max( (SN_u(I,j) + SN_u(I-1,j)) + (SN_v(i,J) + SN_v(i,J-1)), 0.)
! This avoids extremes values in equilibrium solution due to bad values in SN_u, SN_v
SN = min(SN_u(I,j), SN_u(I-1,j), SN_v(i,J), SN_v(i,J-1))
CS%equilibrium_value(i,j) = (CS%MEKE_GEOMETRIC_alpha * SN * US%Z_to_m*G%bathyT(i,j))**2 / cd2
CS%equilibrium_value(i,j) = (CS%MEKE_GEOMETRIC_alpha * SN * US%Z_to_L*G%bathyT(i,j))**2 / cd2
enddo ; enddo

if (CS%id_MEKE_equilibrium>0) call post_data(CS%id_MEKE_equilibrium, CS%equilibrium_value, CS%diag)
Expand Down Expand Up @@ -957,7 +965,7 @@ subroutine MEKE_lengthScales_0d(CS, US, area, beta, depth, Rd_dx, SN, EKE, & ! Z
Leady = 0.
endif
if (CS%use_min_lscale) then
LmixScale = 1.e7
LmixScale = 1.e7*US%m_to_L
if (CS%aDeform*Ldeform > 0.) LmixScale = min(LmixScale,CS%aDeform*Ldeform)
if (CS%aFrict *Lfrict > 0.) LmixScale = min(LmixScale,CS%aFrict *Lfrict)
if (CS%aRhines*Lrhines > 0.) LmixScale = min(LmixScale,CS%aRhines*Lrhines)
Expand Down Expand Up @@ -1069,7 +1077,7 @@ logical function MEKE_init(Time, G, US, param_file, diag, CS, MEKE, restart_CS)
if (CS%MEKE_equilibrium_restoring) then
call get_param(param_file, mdl, "MEKE_RESTORING_TIMESCALE", MEKE_restoring_timescale, &
"The timescale used to nudge MEKE toward its equilibrium value.", units="s", &
default=1e6, scale=US%T_to_s)
default=1e6, scale=US%s_to_T)
CS%MEKE_restoring_rate = 1.0 / MEKE_restoring_timescale
endif

Expand Down
11 changes: 6 additions & 5 deletions src/parameterizations/lateral/MOM_lateral_mixing_coeffs.F90
Original file line number Diff line number Diff line change
Expand Up @@ -506,10 +506,10 @@ subroutine calc_Visbeck_coeffs_old(h, slope_x, slope_y, N2_u, N2_v, G, GV, US, C
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1), intent(in) :: slope_x !< Zonal isoneutral slope
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1), intent(in) :: N2_u !< Buoyancy (Brunt-Vaisala) frequency
!! at u-points [T-2 ~> s-2]
!! at u-points [L2 Z-2 T-2 ~> s-2]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1), intent(in) :: slope_y !< Meridional isoneutral slope
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1), intent(in) :: N2_v !< Buoyancy (Brunt-Vaisala) frequency
!! at v-points [T-2 ~> s-2]
!! at v-points [L2 Z-2 T-2 ~> s-2]
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(VarMix_CS), pointer :: CS !< Variable mixing coefficients
type(ocean_OBC_type), optional, pointer :: OBC !< Open boundaries control structure.
Expand Down Expand Up @@ -654,9 +654,10 @@ subroutine calc_Visbeck_coeffs_old(h, slope_x, slope_y, N2_u, N2_v, G, GV, US, C
endif

if (CS%debug) then
call uvchksum("calc_Visbeck_coeffs_old slope_[xy]", slope_x, slope_y, G%HI, haloshift=1)
call uvchksum("calc_Visbeck_coeffs_old slope_[xy]", slope_x, slope_y, G%HI, &
scale=US%Z_to_L, haloshift=1)
call uvchksum("calc_Visbeck_coeffs_old N2_u, N2_v", N2_u, N2_v, G%HI, &
scale=US%s_to_T**2, scalar_pair=.true.)
scale=US%L_to_Z**2 * US%s_to_T**2, scalar_pair=.true.)
call uvchksum("calc_Visbeck_coeffs_old SN_[uv]", CS%SN_u, CS%SN_v, G%HI, &
scale=US%s_to_T, scalar_pair=.true.)
endif
Expand Down Expand Up @@ -1484,7 +1485,7 @@ subroutine VarMix_init(Time, G, GV, US, param_file, diag, CS)
allocate(CS%Depth_fn_v(isd:ied,JsdB:JedB)) ; CS%Depth_fn_v(:,:) = 0.0
call get_param(param_file, mdl, "DEPTH_SCALED_KHTH_H0", CS%depth_scaled_khth_h0, &
"The depth above which KHTH is scaled away.",&
units="m", default=1000.)
units="m", scale=US%m_to_Z, default=1000.)
call get_param(param_file, mdl, "DEPTH_SCALED_KHTH_EXP", CS%depth_scaled_khth_exp, &
"The exponent used in the depth dependent scaling function for KHTH.",&
units="nondim", default=3.0)
Expand Down
4 changes: 2 additions & 2 deletions src/parameterizations/lateral/MOM_thickness_diffuse.F90
Original file line number Diff line number Diff line change
Expand Up @@ -433,7 +433,7 @@ subroutine thickness_diffuse(h, uhtr, vhtr, tv, dt, G, GV, US, MEKE, VarMix, CDp
call hchksum(e, "thickness_diffuse_1 e", G%HI, haloshift=1, scale=US%Z_to_m)
if (use_stored_slopes) then
call uvchksum("VarMix%slope_[xy]", VarMix%slope_x, VarMix%slope_y, &
G%HI, haloshift=0)
G%HI, haloshift=0, scale=US%Z_to_L)
endif
if (associated(tv%eqn_of_state)) then
call hchksum(tv%T, "thickness_diffuse T", G%HI, haloshift=1)
Expand Down Expand Up @@ -505,7 +505,7 @@ subroutine thickness_diffuse(h, uhtr, vhtr, tv, dt, G, GV, US, MEKE, VarMix, CDp
enddo

do j=js,je ; do i=is,ie
MEKE%Kh_diff(i,j) = MEKE%Kh_diff(i,j) / MAX(1.0,G%bathyT(i,j))
MEKE%Kh_diff(i,j) = MEKE%Kh_diff(i,j) / MAX(GV%m_to_H,GV%Z_to_H*G%bathyT(i,j))
enddo ; enddo
endif

Expand Down
6 changes: 3 additions & 3 deletions src/parameterizations/vertical/MOM_CVMix_KPP.F90
Original file line number Diff line number Diff line change
Expand Up @@ -644,7 +644,7 @@ subroutine KPP_calculate(CS, G, GV, US, h, uStar, &
if (CS%debug) then
call hchksum(h, "KPP in: h",G%HI,haloshift=0, scale=GV%H_to_m)
call hchksum(uStar, "KPP in: uStar",G%HI,haloshift=0, scale=US%Z_to_m*US%s_to_T)
call hchksum(buoyFlux, "KPP in: buoyFlux",G%HI,haloshift=0)
call hchksum(buoyFlux, "KPP in: buoyFlux",G%HI,haloshift=0, scale=US%L_to_m**2*US%s_to_T**3)
call hchksum(Kt, "KPP in: Kt",G%HI,haloshift=0, scale=US%Z2_T_to_m2_s)
call hchksum(Ks, "KPP in: Ks",G%HI,haloshift=0, scale=US%Z2_T_to_m2_s)
endif
Expand Down Expand Up @@ -967,8 +967,8 @@ subroutine KPP_compute_BLD(CS, G, GV, US, h, Temp, Salt, u, v, tv, uStar, buoyFl
if (CS%debug) then
call hchksum(Salt, "KPP in: S",G%HI,haloshift=0)
call hchksum(Temp, "KPP in: T",G%HI,haloshift=0)
call hchksum(u, "KPP in: u",G%HI,haloshift=0)
call hchksum(v, "KPP in: v",G%HI,haloshift=0)
call hchksum(u, "KPP in: u",G%HI,haloshift=0,scale=US%L_T_to_m_s)
call hchksum(v, "KPP in: v",G%HI,haloshift=0,scale=US%L_T_to_m_s)
endif

call cpu_clock_begin(id_clock_KPP_compute_BLD)
Expand Down
4 changes: 2 additions & 2 deletions src/parameterizations/vertical/MOM_CVMix_conv.F90
Original file line number Diff line number Diff line change
Expand Up @@ -288,8 +288,8 @@ subroutine calculate_CVMix_conv(h, tv, G, GV, US, CS, hbl, Kd, Kv, Kd_aux)
! call hchksum(Kd_conv, "MOM_CVMix_conv: Kd_conv", G%HI, haloshift=0, scale=US%Z2_T_to_m2_s)
! if (CS%id_kv_conv > 0) &
! call hchksum(Kv_conv, "MOM_CVMix_conv: Kv_conv", G%HI, haloshift=0, scale=US%m2_s_to_Z2_T)
if (present(Kd)) call hchksum(Kv, "MOM_CVMix_conv: Kd", G%HI, haloshift=0, scale=US%m2_s_to_Z2_T)
if (present(Kv)) call hchksum(Kv, "MOM_CVMix_conv: Kv", G%HI, haloshift=0, scale=US%m2_s_to_Z2_T)
if (present(Kd)) call hchksum(Kd, "MOM_CVMix_conv: Kd", G%HI, haloshift=0, scale=US%Z2_T_to_m2_s)
if (present(Kv)) call hchksum(Kv, "MOM_CVMix_conv: Kv", G%HI, haloshift=0, scale=US%Z2_T_to_m2_s)
endif

! send diagnostics to post_data
Expand Down
Loading