MONC
Functions/Subroutines | Variables
simplecloud_mod Module Reference

A very simple saturation adjustment scheme without any microphysics. More...

Functions/Subroutines

type(component_descriptor_type) function, public simplecloud_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 moisture fields. More...
 
subroutine timestep_callback (current_state)
 Called for each column per timestep this will apply a forcing term to the aerosol fields. More...
 

Variables

integer iqv
 
integer iql
 
integer k_cloudmax
 
real(kind=default_precision) max_height_cloud
 

Detailed Description

A very simple saturation adjustment scheme without any microphysics.

Function/Subroutine Documentation

◆ initialisation_callback()

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

The initialisation callback sets up the moisture fields.

Parameters
current_stateThe current model state

Definition at line 41 of file simplecloud.F90.

41  type(model_state_type), target, intent(inout) :: current_state
42 
43  integer :: k ! look counter
44 
45  if (is_component_enabled(current_state%options_database, "casim")) then
46  call log_master_log(log_error, "Casim and Simplecloud are enabled, this does not work yet. Please disable one")
47  end if
48 
49  iqv=get_q_index(standard_q_names%VAPOUR, 'simplecloud')
50  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'simplecloud')
51 
52  ! set buoyancy coefficient (value for vapour should be set
53  ! elsewhere for a moist model
54  if (.not. allocated(current_state%cq))then
55  allocate(current_state%cq(current_state%number_q_fields))
56  current_state%cq=0.0_default_precision
57  end if
58  current_state%cq(iql) = -1.0
59 
60  max_height_cloud=options_get_real(current_state%options_database, "max_height_cloud")
61  do k=2, current_state%local_grid%size(z_index)-1
62  if (current_state%global_grid%configuration%vertical%zn(k) > max_height_cloud) exit
63  end do
64  k_cloudmax=k-1
65 
Here is the caller graph for this function:

◆ simplecloud_get_descriptor()

type(component_descriptor_type) function, public simplecloud_mod::simplecloud_get_descriptor ( )

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

Returns
The termination check component descriptor

Definition at line 32 of file simplecloud.F90.

32  simplecloud_get_descriptor%name="simplecloud"
33  simplecloud_get_descriptor%version=0.1
34  simplecloud_get_descriptor%initialisation=>initialisation_callback
35  simplecloud_get_descriptor%timestep=>timestep_callback
Here is the call graph for this function:

◆ timestep_callback()

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

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 72 of file simplecloud.F90.

72  type(model_state_type), target, intent(inout) :: current_state
73 
74  real(DEFAULT_PRECISION) :: tdegk ! Temperature in Kelvin
75  real(DEFAULT_PRECISION) :: pmb ! Pressure in mb
76  real(DEFAULT_PRECISION) :: exner ! Exner pressure
77  real(DEFAULT_PRECISION) :: one_over_exner ! Reciprocal of Exner pressure
78  real(DEFAULT_PRECISION) :: qv,qc ! Shorthand for vapour and cloud mass mixing ratio
79  real(DEFAULT_PRECISION) :: qs ! Saturation mixing ratio
80  real(DEFAULT_PRECISION) :: dqsdt ! Rate of change of qs with temperature
81  real(DEFAULT_PRECISION) :: qsatfac ! Multiplicative factor
82  real(DEFAULT_PRECISION) :: dmass ! Mass transfer mixing ratio
83 
84  integer :: k ! Loop counter
85  integer :: icol, jcol ! Shorthand column indices
86 
87  real(DEFAULT_PRECISION) :: dtm ! Local timestep variable
88 
89  if (current_state%halo_column) return
90 
91  dtm = current_state%dtm*2.0
92  if (current_state%field_stepping == forward_stepping) dtm=current_state%dtm! Should this be revised to scalar_stepping
93 
94  icol=current_state%column_local_x
95  jcol=current_state%column_local_y
96 
97  do k=2,k_cloudmax
98 
99  exner = current_state%global_grid%configuration%vertical%rprefrcp(k)
100  one_over_exner = current_state%global_grid%configuration%vertical%prefrcp(k)
101  pmb = (current_state%global_grid%configuration%vertical%prefn(k)/100.)
102 
103  if (current_state%field_stepping == forward_stepping) then ! Should this be revised to scalar_stepping
104  qv = current_state%q(iqv)%data(k, jcol, icol) + current_state%sq(iqv)%data(k, jcol, icol)*dtm
105  qc = current_state%q(iql)%data(k, jcol, icol) + current_state%sq(iql)%data(k, jcol, icol)*dtm
106  tdegk = (current_state%th%data(k, jcol, icol) + current_state%sth%data(k, jcol, icol)*dtm &
107  + current_state%global_grid%configuration%vertical%thref(k))*exner
108  else
109  qv = current_state%zq(iqv)%data(k, jcol, icol) + current_state%sq(iqv)%data(k, jcol, icol)*dtm
110  qc = current_state%zq(iql)%data(k, jcol, icol) + current_state%sq(iql)%data(k, jcol, icol)*dtm
111  tdegk = (current_state%zth%data(k, jcol, icol) + current_state%sth%data(k, jcol, icol)*dtm &
112  + current_state%global_grid%configuration%vertical%thref(k))*exner
113  end if
114  ! Calculate the cloud/vapour increments
115 
116  qs = qsaturation(tdegk, pmb)
117 
118  if (qv > qs .or. qc >0.0)then
119  dqsdt = dqwsatdt(qs, tdegk)
120 
121  qsatfac = 1.0/(1.0 + rlvap_over_cp*dqsdt)
122 
123  dmass = max(-qc,(qv-qs)*qsatfac)/dtm
124 
125  current_state%sq(iqv)%data(k, jcol, icol) = current_state%sq(iqv)%data(k, jcol, icol) - dmass
126  current_state%sq(iql)%data(k, jcol, icol) = current_state%sq(iql)%data(k, jcol, icol) + dmass
127 
128  current_state%sth%data(k, jcol, icol) = current_state%sth%data(k, jcol, icol) &
129  + rlvap_over_cp*dmass*one_over_exner
130 
131  end if
132 
133  end do
134 
135  ! If there's any cloud above then evaporate it
136  do k=k_cloudmax+1, current_state%local_grid%size(z_index)
137  if (current_state%scalar_stepping == forward_stepping) then
138  qv = current_state%q(iqv)%data(k, jcol, icol) + current_state%sq(iqv)%data(k, jcol, icol)*dtm
139  qc = current_state%q(iql)%data(k, jcol, icol) + current_state%sq(iql)%data(k, jcol, icol)*dtm
140  tdegk = (current_state%th%data(k, jcol, icol) + current_state%sth%data(k, jcol, icol)*dtm &
141  + current_state%global_grid%configuration%vertical%thref(k))*exner
142  else
143  qv = current_state%zq(iqv)%data(k, jcol, icol) + current_state%sq(iqv)%data(k, jcol, icol)*dtm
144  qc = current_state%zq(iql)%data(k, jcol, icol) + current_state%sq(iql)%data(k, jcol, icol)*dtm
145  tdegk = (current_state%zth%data(k, jcol, icol) + current_state%sth%data(k, jcol, icol)*dtm &
146  + current_state%global_grid%configuration%vertical%thref(k))*exner
147  end if
148  if (qc >0.0)then
149  dmass = -qc/dtm
150 
151  current_state%sq(iqv)%data(k, jcol, icol) = current_state%sq(iqv)%data(k, jcol, icol) - dmass
152  current_state%sq(iql)%data(k, jcol, icol) = current_state%sq(iql)%data(k, jcol, icol) + dmass
153 
154  current_state%sth%data(k, jcol, icol) = current_state%sth%data(k, jcol, icol) &
155  + rlvap_over_cp*dmass*one_over_exner
156 
157  end if
158  end do
159 
Here is the caller graph for this function:

Variable Documentation

◆ iql

integer simplecloud_mod::iql
private

Definition at line 21 of file simplecloud.F90.

◆ iqv

integer simplecloud_mod::iqv
private

Definition at line 21 of file simplecloud.F90.

21  integer :: iqv, iql

◆ k_cloudmax

integer simplecloud_mod::k_cloudmax
private

Definition at line 23 of file simplecloud.F90.

23  integer :: k_cloudmax ! max k index for height

◆ max_height_cloud

real(kind=default_precision) simplecloud_mod::max_height_cloud
private

Definition at line 24 of file simplecloud.F90.

24  real(kind=DEFAULT_PRECISION) :: max_height_cloud