Browse Source

Fixed bug with redistributing radiation fluxes.

Multiplication by area fraction lead to a too weak radiation flux.
Since flux is a density-like quantity there must not be a multiplication
with fraction of area.

The too cold bottom model seems to be fixed by this.
1.00.00
Sven Karsten 2 years ago
parent
commit
5811994ca4
  1. 7
      src/flux_calculator_calculate.F90
  2. 10
      src/flux_lib/radiation/distribute_radiation_flux.F90

7
src/flux_calculator_calculate.F90

@ -97,7 +97,7 @@ MODULE flux_calculator_calculate
method = methods(my_bottom_model,i)
IF (trim(method) /= 'none') THEN
IF (trim(method)=='zero') THEN
local_field(i,1)%var(idx_MEVA)%field(:) = 0.0
local_field(i,1)%var(idx_HLAT)%field(:) = 0.0
ELSEIF (trim(method)=='water') THEN
DO j=1,grid_size(1)
CALL flux_heat_latent_water(local_field(i,1)%var(idx_HLAT)%field(j), &
@ -129,7 +129,7 @@ MODULE flux_calculator_calculate
method = methods(my_bottom_model,i)
IF (trim(method) /= 'none') THEN
IF (trim(method)=='zero') THEN
local_field(i,1)%var(idx_MEVA)%field(:) = 0.0
local_field(i,1)%var(idx_HSEN)%field(:) = 0.0
ELSEIF (trim(method)=='CCLM') THEN
DO j=1,grid_size(1)
CALL flux_heat_sensible_cclm(local_field(i,1)%var(idx_HSEN)%field(j), &
@ -262,8 +262,7 @@ MODULE flux_calculator_calculate
CALL distribute_radiation_flux(local_field(i,1)%var(idx_RSDR)%field(j), &
local_field(0,1)%var(idx_RSDD)%field(j), &
local_field(0,1)%var(idx_ALBA)%field(j), &
local_field(i,1)%var(idx_ALBE)%field(j), &
local_field(i,1)%var(idx_FARE)%field(j))
local_field(i,1)%var(idx_ALBE)%field(j))
ENDDO
ENDDO
!CALL average_across_surface_types(1,idx_RSDD,num_surface_types,grid_size,local_field)

10
src/flux_lib/radiation/distribute_radiation_flux.F90

@ -13,16 +13,16 @@ module distribute_radiation_flux_mod
flux_radiation_surface_type, & ! RESULT, surface-dependent radiation flux sent to bottom
flux_radiation_averaged, & ! surface-independent radiation flux from atmosphere
albedo_averaged, & ! surface-independent albedo from atmosphere
albedo_surface_type, & ! surface-dependent albedo received from ocean
fraction_surface_type & ! fraction of surface in exchange-grid cell
albedo_surface_type & ! surface-dependent albedo received from ocean
)
real(prec), intent(out) :: flux_radiation_surface_type ! RESULT
real(prec), intent(in) :: flux_radiation_averaged, albedo_averaged, albedo_surface_type, fraction_surface_type
real(prec), intent(in) :: flux_radiation_averaged, albedo_averaged, albedo_surface_type
flux_radiation_surface_type = fraction_surface_type * flux_radiation_averaged * (1.0 - albedo_surface_type) / (1.0 - albedo_averaged)
! apply surface-type-dependent albedo and get rid of averaged albedo
flux_radiation_surface_type = flux_radiation_averaged * (1.0 - albedo_surface_type) / (1.0 - albedo_averaged)
end subroutine distribute_radiation_flux
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

Loading…
Cancel
Save