MONC
Functions/Subroutines | Variables
lwrad_exponential_mod Module Reference

Simple exponential scheme to calculate the longwave radiation associated with cloud. The scheme is based on the methods used in GASS intercomparison cases, e.g. DYCOMS, ISDAC. More...

Functions/Subroutines

type(component_descriptor_type) function, public lwrad_exponential_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine initialisation_callback (current_state)
 The initialisation callback sets up the prescribed longwave fluxes and the exponential decay factor. More...
 
subroutine timestep_callback (current_state)
 Called for each column per timestep this will apply a forcing term to the aerosol fields. More...
 
subroutine finalisation_callback (current_state)
 

Variables

integer iql
 
integer k_top
 
integer x_local
 
integer y_local
 
real(default_precision), dimension(:), allocatable lwrad_flux_top
 
real(default_precision), dimension(:), allocatable lwrad_flux_base
 
real(default_precision), dimension(:), allocatable lwrad_flux
 
real(default_precision), dimension(:), allocatable qc_col
 
real(default_precision), dimension(:), allocatable density_factor
 
real(default_precision), dimension(:), allocatable radiation_factor
 
real(default_precision), dimension(:,:,:), allocatable sth_lw
 
real(default_precision) longwave_exp_decay
 
real(default_precision) cltop_longwave_flux
 
real(default_precision) clbase_longwave_flux
 

Detailed Description

Simple exponential scheme to calculate the longwave radiation associated with cloud. The scheme is based on the methods used in GASS intercomparison cases, e.g. DYCOMS, ISDAC.

This scheme depends on exp_lw, lwtop_in, lwbase_in, which are set in the config file.

Function/Subroutine Documentation

◆ finalisation_callback()

subroutine lwrad_exponential_mod::finalisation_callback ( type(model_state_type), intent(inout), target  current_state)

Definition at line 156 of file lwrad_exponential.F90.

156 
157  type(model_state_type), target, intent(inout) :: current_state
158  deallocate(lwrad_flux_top,lwrad_flux_base, lwrad_flux, density_factor, radiation_factor, sth_lw)
159 
Here is the caller graph for this function:

◆ initialisation_callback()

subroutine lwrad_exponential_mod::initialisation_callback ( type(model_state_type), intent(inout), target  current_state)

The initialisation callback sets up the prescribed longwave fluxes and the exponential decay factor.

Parameters
current_stateThe current model state

Definition at line 56 of file lwrad_exponential.F90.

56  type(model_state_type), target, intent(inout) :: current_state
57 
58  integer :: kkp, k ! look counter
59  k_top = current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
60  x_local = current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
61  y_local = current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(x_index) * 2
62 
63  allocate(lwrad_flux_top(k_top))
64  allocate(lwrad_flux_base(k_top))
65  allocate(lwrad_flux(k_top))
66  allocate(qc_col(k_top))
67  allocate(density_factor(k_top))
68  allocate(radiation_factor(k_top))
69  ! NOTE: this may be the wrong declaration as sth_lw may need to declared on the whole domain
70  allocate(sth_lw(k_top,y_local,x_local))
71 
72  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'lwrad_exponential')
73 
74  longwave_exp_decay=options_get_real(current_state%options_database, "longwave_exp_decay")
75  cltop_longwave_flux=options_get_real(current_state%options_database, "cltop_longwave_flux")
76  clbase_longwave_flux=options_get_real(current_state%options_database, "clbase_longwave_flux")
77 
78  density_factor(:) = current_state%global_grid%configuration%vertical%rhon(:)* &
79  current_state%global_grid%configuration%vertical%dz(:)
80 
81  radiation_factor(2:k_top) = 1.0/(density_factor(2:k_top)*cp)
82  radiation_factor(1) = radiation_factor(2)
83 
Here is the caller graph for this function:

◆ lwrad_exponential_get_descriptor()

type(component_descriptor_type) function, public lwrad_exponential_mod::lwrad_exponential_get_descriptor ( )

Provides the descriptor back to the caller and is used in component registration.

Returns
The termination check component descriptor

Definition at line 45 of file lwrad_exponential.F90.

45  lwrad_exponential_get_descriptor%name="lwrad_exponential"
46  lwrad_exponential_get_descriptor%version=0.1
47  lwrad_exponential_get_descriptor%initialisation=>initialisation_callback
48  lwrad_exponential_get_descriptor%timestep=>timestep_callback
49  lwrad_exponential_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ timestep_callback()

subroutine lwrad_exponential_mod::timestep_callback ( type(model_state_type), intent(inout), target  current_state)

Called for each column per timestep this will apply a forcing term to the aerosol fields.

Parameters
current_stateThe current model state

Definition at line 90 of file lwrad_exponential.F90.

90  type(model_state_type), target, intent(inout) :: current_state
91 
92  integer :: k ! Loop counter
93  integer :: icol, jcol ! Shorthand column indices
94 
95  real(DEFAULT_PRECISION) :: dtm ! Local timestep variable
96 
97  if (current_state%halo_column) return
98 
99  dtm = current_state%dtm*2.0
100  if (current_state%field_stepping == forward_stepping) dtm=current_state%dtm! Should this be revised to scalar_stepping
101 
102  icol=current_state%column_local_x
103  jcol=current_state%column_local_y
104 
105  ! set the column liquid water content
106  if (current_state%field_stepping == forward_stepping) then ! Should this be revised to scalar_stepping
107  qc_col(:) = current_state%q(iql)%data(:, jcol, icol) + current_state%sq(iql)%data(:, jcol, icol)*dtm
108  else
109  qc_col(:)= current_state%zq(iql)%data(:, jcol, icol) + current_state%sq(iql)%data(:, jcol, icol)*dtm
110 
111  end if
112 
113  ! initialise the flux top and base to 0.0
114  lwrad_flux_top(:) = 0.0
115  lwrad_flux_base(:) = 0.0
116 
117  ! First workout cloud top cooling
118  lwrad_flux_top(k_top)=cltop_longwave_flux !units Wm-2
119 
120  do k = k_top-1, 1, -1
121  if (qc_col(k+1) > 1.e-10) then
122  lwrad_flux_top(k) = lwrad_flux_top(k+1)* &
123  exp(-qc_col(k+1)*density_factor(k+1)*longwave_exp_decay)
124  else
125  lwrad_flux_top(k) = lwrad_flux_top(k+1)
126  endif
127  enddo
128 
129  ! Second workout the cloud base warming
130  do k = 2, k_top
131  if (qc_col(k) > 1.e-10) then
132  lwrad_flux_base(k) = lwrad_flux_base(k-1)* &
133  exp(-qc_col(k)*density_factor(k)*longwave_exp_decay)
134  else
135  lwrad_flux_base(k) = lwrad_flux_base(k-1)
136  endif
137  enddo
138 
139  ! Next, sum the fluxes
140  do k = 1, k_top
141  lwrad_flux(k) = lwrad_flux_base(k) + lwrad_flux_top(k)
142  enddo
143 
144  ! workout radiative heating rate, and stored on the processors local grid - is this correct??
145  ! or should the declaration be on the global grid?
146  do k = 2, k_top
147  sth_lw(k, jcol, icol) = -(lwrad_flux(k) - lwrad_flux(k-1))*radiation_factor(k)
148  enddo
149 
150  ! update the current_state sth
151  current_state%sth%data(:, jcol, icol) = current_state%sth%data(:, jcol, icol) + sth_lw(:, jcol, icol)
152 
Here is the caller graph for this function:

Variable Documentation

◆ clbase_longwave_flux

real(default_precision) lwrad_exponential_mod::clbase_longwave_flux

Definition at line 36 of file lwrad_exponential.F90.

◆ cltop_longwave_flux

real(default_precision) lwrad_exponential_mod::cltop_longwave_flux

Definition at line 36 of file lwrad_exponential.F90.

◆ density_factor

real(default_precision), dimension(:), allocatable lwrad_exponential_mod::density_factor

Definition at line 29 of file lwrad_exponential.F90.

29  real(DEFAULT_PRECISION), allocatable :: density_factor(:), radiation_factor(:)

◆ iql

integer lwrad_exponential_mod::iql

Definition at line 19 of file lwrad_exponential.F90.

19  integer :: iql

◆ k_top

integer lwrad_exponential_mod::k_top

Definition at line 21 of file lwrad_exponential.F90.

21  integer :: k_top, x_local, y_local

◆ longwave_exp_decay

real(default_precision) lwrad_exponential_mod::longwave_exp_decay

Definition at line 36 of file lwrad_exponential.F90.

36  real(DEFAULT_PRECISION) :: longwave_exp_decay, cltop_longwave_flux, clbase_longwave_flux

◆ lwrad_flux

real(default_precision), dimension(:), allocatable lwrad_exponential_mod::lwrad_flux

Definition at line 25 of file lwrad_exponential.F90.

◆ lwrad_flux_base

real(default_precision), dimension(:), allocatable lwrad_exponential_mod::lwrad_flux_base

Definition at line 25 of file lwrad_exponential.F90.

◆ lwrad_flux_top

real(default_precision), dimension(:), allocatable lwrad_exponential_mod::lwrad_flux_top

Definition at line 25 of file lwrad_exponential.F90.

25  real(DEFAULT_PRECISION), allocatable :: lwrad_flux_top(:), lwrad_flux_base(:), lwrad_flux(:)

◆ qc_col

real(default_precision), dimension(:), allocatable lwrad_exponential_mod::qc_col

Definition at line 27 of file lwrad_exponential.F90.

27  real(DEFAULT_PRECISION), allocatable :: qc_col(:)

◆ radiation_factor

real(default_precision), dimension(:), allocatable lwrad_exponential_mod::radiation_factor

Definition at line 29 of file lwrad_exponential.F90.

◆ sth_lw

real(default_precision), dimension(:,:,:), allocatable lwrad_exponential_mod::sth_lw

Definition at line 33 of file lwrad_exponential.F90.

33  real(DEFAULT_PRECISION), allocatable :: sth_lw(:,:,:)

◆ x_local

integer lwrad_exponential_mod::x_local

Definition at line 21 of file lwrad_exponential.F90.

◆ y_local

integer lwrad_exponential_mod::y_local

Definition at line 21 of file lwrad_exponential.F90.