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

sensible+latent heatfluxes using linear bulk formula #371

Merged
merged 21 commits into from
Sep 23, 2021
Merged
Show file tree
Hide file tree
Changes from 16 commits
Commits
Show all changes
21 commits
Select commit Hold shift + click to select a range
cd863c4
'heatflux_linear' flag: sensible+latent heatfluxes using traditional …
mhrib Sep 14, 2021
aa871ef
ERROR correction: heatflux_linear is logical
mhrib Sep 14, 2021
d602c2c
ERROR/syntax correction.
mhrib Sep 14, 2021
2da0b98
Add option atmbndy='mixed' boundary layer condition
mhrib Sep 15, 2021
dd63009
For 'atmbndy': change 'default' (obsolete) to 'similarity'. Same physics
mhrib Sep 15, 2021
bd91a7d
For 'atmbndy': change 'default' (obsolete) to 'similarity'. Same physics
mhrib Sep 15, 2021
b785a57
New options for 'atmbndy'
mhrib Sep 15, 2021
ffb374d
"atmbndy" options: 'similarity', 'constant', 'mixed'
mhrib Sep 16, 2021
4cfa3ab
Abort if "atmbndy" option is unknown + set default=similarity for bac…
mhrib Sep 16, 2021
644d324
Update doc/source/science_guide/sg_boundary_forcing.rst
mhrib Sep 16, 2021
ef391d5
atmbndy option 'ccsm' not an option. Previous given as 'default'
mhrib Sep 16, 2021
2c10349
namelist parameter plus string inside double quotes
mhrib Sep 16, 2021
eb58e3f
Merge branch 'Icepack' of https://github.com/mhrib/Icepack into Icepack
mhrib Sep 16, 2021
a44c806
Re-interduce file, after accidentally git rm command
mhrib Sep 16, 2021
56ccc77
Remove "ccsm" as an option to atmbndy
mhrib Sep 16, 2021
4807cc4
atmbndy: 'constant' or 'mixed'. Anything else end in 'similarity' beh…
mhrib Sep 17, 2021
78660ed
Error correction. Revert atmbndy to constant as before
mhrib Sep 20, 2021
1ec0b3f
Put constants intp icepack_parameters
mhrib Sep 20, 2021
bee1627
Introduce p0012, p0015 from icepack_parameters
mhrib Sep 20, 2021
2640405
Changed p0012/p0015 to specific constants senscoef/latncoef
mhrib Sep 23, 2021
81b80ba
senscoef/latncoef introduced
mhrib Sep 23, 2021
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
14 changes: 11 additions & 3 deletions columnphysics/icepack_atmo.F90
Original file line number Diff line number Diff line change
Expand Up @@ -359,8 +359,16 @@ subroutine atmo_boundary_layer (sfctype, &
! as in Jordan et al (JGR, 1999)
!------------------------------------------------------------

shcoef = rhoa * ustar * cp * rh + c1
lhcoef = rhoa * ustar * Lheat * re
if (trim(atmbndy) == 'mixed') then
!- Use constant coefficients for sensible and latent heat fluxes
! similar to atmo_boundary_const but using vmag instead of wind
shcoef = (1.20e-3_dbl_kind)*cp_air*rhoa*vmag
lhcoef = (1.50e-3_dbl_kind)*Lheat *rhoa*vmag
Copy link
Contributor

@proteanplanet proteanplanet Sep 17, 2021

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

If possible, please remove literal constants from the code (here 1.20e-3_dbl_kind and 1.50e-3_dbl_kind) and place them in the Icepack constants file. I appreciate this is already done in atmo_boundary_const, but if we leave two sets of identical literal constants in the code, the code is susceptible to errors creeping in if further changes are made by another developer to just one set of constants.

Copy link
Contributor Author

@mhrib mhrib Sep 20, 2021

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Good point w.r.t. secure code for future errors by placing the identical parameters in one place.
Will introduce p0012 and p0015.

else ! 'similarity'
!- Monin-Obukhov similarity theory for boundary layer
shcoef = rhoa * ustar * cp * rh + c1
lhcoef = rhoa * ustar * Lheat * re
endif

!------------------------------------------------------------
! Compute diagnostics: 2m ref T, Q, U
Expand Down Expand Up @@ -952,7 +960,7 @@ subroutine icepack_atm_boundary(sfctype, &
delt, delq, &
lhcoef, shcoef )
if (icepack_warnings_aborted(subname)) return
else ! default
else
call atmo_boundary_layer (sfctype, &
calc_strair, formdrag, &
Tsf, potT, &
Expand Down
18 changes: 9 additions & 9 deletions columnphysics/icepack_parameters.F90
Original file line number Diff line number Diff line change
Expand Up @@ -250,7 +250,7 @@ module icepack_parameters
TTTocn = 5107.4_dbl_kind ! for qsat over ocn

character (len=char_len), public :: &
atmbndy = 'default' ! atmo boundary method, 'default' ('ccsm3') or 'constant'
atmbndy = 'similarity' ! atmo boundary method, 'similarity', 'constant' or 'mixed'

logical (kind=log_kind), public :: &
calc_strair = .true. , & ! if true, calculate wind stress
Expand Down Expand Up @@ -646,12 +646,12 @@ subroutine icepack_init_parameters( &
TTTocn_in ! for qsat over ocn

character (len=*), intent(in), optional :: &
atmbndy_in ! atmo boundary method, 'default' ('ccsm3') or 'constant'
atmbndy_in ! atmo boundary method, 'similarity', 'constant' or 'mixed'

logical (kind=log_kind), intent(in), optional :: &
calc_strair_in, & ! if true, calculate wind stress components
formdrag_in, & ! if true, calculate form drag
highfreq_in ! if true, use high frequency coupling
calc_strair_in, & ! if true, calculate wind stress components
formdrag_in, & ! if true, calculate form drag
highfreq_in ! if true, use high frequency coupling

integer (kind=int_kind), intent(in), optional :: &
natmiter_in ! number of iterations for boundary layer calculations
Expand Down Expand Up @@ -1326,12 +1326,12 @@ subroutine icepack_query_parameters( &
TTTocn_out ! for qsat over ocn

character (len=*), intent(out), optional :: &
atmbndy_out ! atmo boundary method, 'default' ('ccsm3') or 'constant'
atmbndy_out ! atmo boundary method, 'similarity', 'constant' or 'mixed'

logical (kind=log_kind), intent(out), optional :: &
calc_strair_out, & ! if true, calculate wind stress components
formdrag_out, & ! if true, calculate form drag
highfreq_out ! if true, use high frequency coupling
calc_strair_out, & ! if true, calculate wind stress components
formdrag_out, & ! if true, calculate form drag
highfreq_out ! if true, use high frequency coupling

integer (kind=int_kind), intent(out), optional :: &
natmiter_out ! number of iterations for boundary layer calculations
Expand Down
4 changes: 2 additions & 2 deletions configuration/driver/icedrv_init.F90
Original file line number Diff line number Diff line change
Expand Up @@ -563,8 +563,8 @@ subroutine input_data
if (formdrag) then
if (trim(atmbndy) == 'constant') then
write (nu_diag,*) 'WARNING: atmbndy = constant not allowed with formdrag'
write (nu_diag,*) 'WARNING: Setting atmbndy = default'
atmbndy = 'default'
write (nu_diag,*) 'WARNING: Setting atmbndy = similarity'
atmbndy = 'similarity'
endif

if (.not. calc_strair) then
Expand Down
2 changes: 1 addition & 1 deletion configuration/scripts/icepack_in
Original file line number Diff line number Diff line change
Expand Up @@ -97,7 +97,7 @@

&forcing_nml
formdrag = .false.
atmbndy = 'default'
atmbndy = 'similarity'
calc_strair = .true.
calc_Tsfc = .true.
highfreq = .false.
Expand Down
2 changes: 1 addition & 1 deletion configuration/scripts/options/set_nml.dyn
Original file line number Diff line number Diff line change
Expand Up @@ -3,5 +3,5 @@ kstrength = 0
krdg_partic = 0
krdg_redist = 0
formdrag = .false.
atmbndy = 'constant'
atmbndy = 'similarity'
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This changes the behavior of at least one test in the test suites. Is that intended?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This was an error from my site. Will revert atmbndy to constant

highfreq = .false.
2 changes: 1 addition & 1 deletion doc/source/icepack_index.rst
Original file line number Diff line number Diff line change
Expand Up @@ -52,7 +52,7 @@ either Celsius or Kelvin units).
"atmiter_conv", ":math:`\bullet` convergence criteria for ustar", "0.0"
"atm_data_format", ":math:`\bullet` format of atmospheric forcing files", ""
"atm_data_type", ":math:`\bullet` type of atmospheric forcing", ""
"atmbndy", ":math:`\bullet` atmo boundary layer parameterization (‘default’ or ‘constant’)", ""
"atmbndy", ":math:`\bullet` atmo boundary layer parameterization (‘similarity’,‘constant’ or 'mixed')", ""
"awtidf", "weighting factor for near-ir, diffuse albedo", "0.36218"
"awtidr", "weighting factor for near-ir, direct albedo", "0.00182"
"awtvdf", "weighting factor for visible, diffuse albedo", "0.63282"
Expand Down
7 changes: 6 additions & 1 deletion doc/source/science_guide/sg_boundary_forcing.rst
Original file line number Diff line number Diff line change
Expand Up @@ -120,7 +120,12 @@ parameterizations. Rain and all melted snow end up in the ocean.
Wind stress and transfer coefficients for the
turbulent heat fluxes are computed in subroutine
*atmo\_boundary\_layer* following :cite:`Kauffman02`, with additions and changes as detailed in Appendix A of :cite:`Roberts15` for high frequency coupling (namelist variable ``highfreq``).
The resulting equations are provided here.
The resulting equations are provided here for the default boundary layer
scheme, which is based on Monin-Obukhov theory (``atmbndy = ‘stability’``).
Alternatively, ``atmbndy = ‘constant’`` provides constant coefficients for
wind stress, sensible heat and latent heat calculations (computed in subroutine
*atmo\_boundary\_const*); ``atmbndy = ‘mixed’`` uses the stability based
calculation for wind stress and constant coefficients for sensible and latent heat fluxes.

The wind stress and turbulent heat flux calculation accounts for both
stable and unstable atmosphere–ice boundary layers. We first define the
Expand Down
4 changes: 2 additions & 2 deletions doc/source/user_guide/interfaces.include
Original file line number Diff line number Diff line change
Expand Up @@ -944,7 +944,7 @@ icepack_init_parameters
TTTocn_in ! for qsat over ocn

character (len=char_len), intent(in), optional :: &
atmbndy_in ! atmo boundary method, 'default' ('ccsm3') or 'constant'
atmbndy_in ! atmo boundary method, 'similarity', 'constant' or 'mixed'

logical (kind=log_kind), intent(in), optional :: &
calc_strair_in, & ! if true, calculate wind stress components
Expand Down Expand Up @@ -1306,7 +1306,7 @@ icepack_query_parameters
TTTocn_out ! for qsat over ocn

character (len=char_len), intent(out), optional :: &
atmbndy_out ! atmo boundary method, 'default' ('ccsm3') or 'constant'
atmbndy_out ! atmo boundary method, 'similarity', 'constant' or 'mixed'

logical (kind=log_kind), intent(out), optional :: &
calc_strair_out, & ! if true, calculate wind stress components
Expand Down
6 changes: 4 additions & 2 deletions doc/source/user_guide/ug_case_settings.rst
Original file line number Diff line number Diff line change
Expand Up @@ -314,8 +314,10 @@ forcing_nml
:widths: 15, 15, 30, 15

"", "", "", ""
"``atmbndy``", "``constant``", "bulk transfer coefficients", "``default``"
"", "``default``", "stability-based boundary layer", ""
"``atmbndy``", "string", "bulk transfer coefficients", "``similarity``"
"", "``similarity``", "stability-based boundary layer", ""
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Maybe here we could add something like

Suggested change
"", "``similarity``", "stability-based boundary layer", ""
"", "``similarity``", "stability-based boundary layer (previously: `default`)", ""

to help users that could be confused by their existing namelist having a now undocumented default value ?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

From CICE this is automatically changed from 'default' to 'similarity' together with a WARNING message printed to the log/stdout.

"", "``constant``", "constant-based boundary layer", ""
"", "``mixed``", "stability-based, but constant for sensible+latent heatfluxes", ""
"``atmiter_conv``", "real", "convergence criteria for ustar", "0.0"
"``atm_data_file``", "string", "file containing atmospheric data", "' '"
"``atm_data_format``", "``bin``", "read direct access binary forcing files", "``bin``"
Expand Down