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 164 of file lwrad_exponential.F90.

165 
166  type(model_state_type), target, intent(inout) :: current_state
167  deallocate(lwrad_flux_top,lwrad_flux_base, lwrad_flux, density_factor, radiation_factor, sth_lw)
168 
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 57 of file lwrad_exponential.F90.

58  type(model_state_type), target, intent(inout) :: current_state
59 
60  integer :: kkp, k ! look counter
61 
62 
63  k_top = current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
64  x_local = current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
65  y_local = current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(x_index) * 2
66 
67  if (is_component_enabled(current_state%options_database, "socrates_couple")) then
68  call log_master_log &
69  (log_error, "Socrates and lwrad_exponential both enabled, switch off one in config - STOP")
70  endif
71 
72  allocate(lwrad_flux_top(k_top))
73  allocate(lwrad_flux_base(k_top))
74  allocate(lwrad_flux(k_top))
75  allocate(qc_col(k_top))
76  allocate(density_factor(k_top))
77  allocate(radiation_factor(k_top))
78  ! NOTE: this may be the wrong declaration as sth_lw may need to declared on the whole domain
79  allocate(sth_lw(k_top,y_local,x_local))
80 
81  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'lwrad_exponential')
82 
83  longwave_exp_decay=options_get_real(current_state%options_database, "longwave_exp_decay")
84  cltop_longwave_flux=options_get_real(current_state%options_database, "cltop_longwave_flux")
85  clbase_longwave_flux=options_get_real(current_state%options_database, "clbase_longwave_flux")
86 
87  density_factor(:) = current_state%global_grid%configuration%vertical%rhon(:)* &
88  current_state%global_grid%configuration%vertical%dz(:)
89 
90  radiation_factor(2:k_top) = 1.0/(density_factor(2:k_top)*cp)
91  radiation_factor(1) = radiation_factor(2)
92 
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 46 of file lwrad_exponential.F90.

47  lwrad_exponential_get_descriptor%name="lwrad_exponential"
48  lwrad_exponential_get_descriptor%version=0.1
49  lwrad_exponential_get_descriptor%initialisation=>initialisation_callback
50  lwrad_exponential_get_descriptor%timestep=>timestep_callback
51  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 98 of file lwrad_exponential.F90.

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

◆ cltop_longwave_flux

real(default_precision) lwrad_exponential_mod::cltop_longwave_flux

Definition at line 38 of file lwrad_exponential.F90.

◆ density_factor

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

Definition at line 31 of file lwrad_exponential.F90.

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

◆ iql

integer lwrad_exponential_mod::iql

Definition at line 21 of file lwrad_exponential.F90.

21  integer :: iql

◆ k_top

integer lwrad_exponential_mod::k_top

Definition at line 23 of file lwrad_exponential.F90.

23  integer :: k_top, x_local, y_local

◆ longwave_exp_decay

real(default_precision) lwrad_exponential_mod::longwave_exp_decay

Definition at line 38 of file lwrad_exponential.F90.

38  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 27 of file lwrad_exponential.F90.

◆ lwrad_flux_base

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

Definition at line 27 of file lwrad_exponential.F90.

◆ lwrad_flux_top

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

Definition at line 27 of file lwrad_exponential.F90.

27  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 29 of file lwrad_exponential.F90.

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

◆ radiation_factor

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

Definition at line 31 of file lwrad_exponential.F90.

◆ sth_lw

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

Definition at line 35 of file lwrad_exponential.F90.

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

◆ x_local

integer lwrad_exponential_mod::x_local

Definition at line 23 of file lwrad_exponential.F90.

◆ y_local

integer lwrad_exponential_mod::y_local

Definition at line 23 of file lwrad_exponential.F90.

logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
logging_mod::log_master_log
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
optionsdatabase_mod::options_get_real
real(kind=default_precision) function, public options_get_real(options_database, key, index)
Retrieves a real value from the database that matches the provided key.
Definition: optionsdatabase.F90:91