MONC
Functions/Subroutines | Variables
scalar_diagnostics_mod Module Reference

Functions/Subroutines

type(component_descriptor_type) function, public scalar_diagnostics_get_descriptor ()
 Provides the component descriptor for the core to register. More...
 
subroutine initialisation_callback (current_state)
 
subroutine timestep_callback (current_state)
 
subroutine field_information_retrieval_callback (current_state, name, field_information)
 Field information retrieval callback, this returns information for a specific components published field. More...
 
subroutine field_value_retrieval_callback (current_state, name, field_value)
 Field value retrieval callback, this returns the value of a specific published field. More...
 

Variables

integer total_points
 
integer ncl_col
 
integer iqv =0
 
integer iql =0
 
integer iqr =0
 
integer iqi =0
 
integer iqs =0
 
integer iqg =0
 
real(kind=default_precision), dimension(:), allocatable dz_rhon_fac
 
real(kind=default_precision), dimension(:), allocatable ww_prime_res
 
real(kind=default_precision), dimension(:), allocatable uu_prime_res
 
real(kind=default_precision), dimension(:), allocatable vv_prime_res
 
real(kind=default_precision), dimension(:), allocatable cloud_content
 
real(kind=default_precision) qlcrit
 
real(kind=default_precision), dimension(:,:), allocatable vwp
 
real(kind=default_precision), dimension(:,:), allocatable lwp
 
real(kind=default_precision), dimension(:,:), allocatable wmax
 
real(kind=default_precision), dimension(:,:), allocatable wmin
 
real(kind=default_precision), dimension(:,:), allocatable qlmax
 
real(kind=default_precision), dimension(:,:), allocatable hqlmax
 
real(kind=default_precision), dimension(:,:), allocatable cltop
 
real(kind=default_precision), dimension(:,:), allocatable clbas
 
real(kind=default_precision), dimension(:,:), allocatable senhf
 
real(kind=default_precision), dimension(:,:), allocatable lathf
 
real(kind=default_precision), dimension(:,:), allocatable rwp
 
real(kind=default_precision), dimension(:,:), allocatable iwp
 
real(kind=default_precision), dimension(:,:), allocatable swp
 
real(kind=default_precision), dimension(:,:), allocatable gwp
 
real(kind=default_precision), dimension(:,:), allocatable tot_iwp
 
real(kind=default_precision), dimension(:,:), allocatable reske
 

Function/Subroutine Documentation

◆ field_information_retrieval_callback()

subroutine scalar_diagnostics_mod::field_information_retrieval_callback ( type(model_state_type), intent(inout), target  current_state,
character(len=*), intent(in)  name,
type(component_field_information_type), intent(out)  field_information 
)
private

Field information retrieval callback, this returns information for a specific components published field.

Parameters
current_stateCurrent model state
nameThe name of the field to retrieve information for
field_informationPopulated with information about the field

Definition at line 322 of file scalar_diagnostics.F90.

323  type(model_state_type), target, intent(inout) :: current_state
324  character(len=*), intent(in) :: name
325  type(component_field_information_type), intent(out) :: field_information
326 
327  field_information%field_type=component_array_field_type
328  field_information%number_dimensions=2
329  field_information%dimension_sizes(1)=current_state%local_grid%size(y_index)
330  field_information%dimension_sizes(2)=current_state%local_grid%size(x_index)
331  field_information%data_type=component_double_data_type
332 
333  if (name .eq. "senhf") then
334  field_information%enabled=current_state%use_surface_boundary_conditions .and. current_state%th%active
335  else if (name .eq. "lathf") then
336  field_information%enabled=current_state%use_surface_boundary_conditions .and. &
337  current_state%water_vapour_mixing_ratio_index .gt. 0 .and. &
338  current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index
339  else if (name .eq. "qlmax".or. name .eq. "cltop" .or. name .eq. "clbas") then
340  field_information%enabled=.not. current_state%passive_q .and. current_state%liquid_water_mixing_ratio_index .gt. 0 &
341  .and. current_state%number_q_fields .ge. current_state%liquid_water_mixing_ratio_index
342  else if (name .eq. "vwp" .or. name .eq. "lwp") then
343  field_information%enabled=current_state%number_q_fields .gt. 0 .and. current_state%water_vapour_mixing_ratio_index .gt. 0 &
344  .and. current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index
345  else if (name .eq. "rwp" ) then
346  field_information%enabled= current_state%rain_water_mixing_ratio_index .gt. 0
347  else if (name .eq. "iwp" .or. name .eq. 'tot_iwp') then
348  field_information%enabled= current_state%ice_water_mixing_ratio_index .gt. 0
349  else if (name .eq. "swp" ) then
350  field_information%enabled= current_state%snow_water_mixing_ratio_index .gt. 0
351  else if (name .eq. "gwp" ) then
352  field_information%enabled= current_state%graupel_water_mixing_ratio_index .gt. 0
353  else
354  field_information%enabled=.true.
355  end if
356 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine scalar_diagnostics_mod::field_value_retrieval_callback ( type(model_state_type), intent(inout), target  current_state,
character(len=*), intent(in)  name,
type(component_field_value_type), intent(out)  field_value 
)
private

Field value retrieval callback, this returns the value of a specific published field.

Parameters
current_stateCurrent model state
nameThe name of the field to retrieve the value for
field_valuePopulated with the value of the field

Definition at line 363 of file scalar_diagnostics.F90.

364  type(model_state_type), target, intent(inout) :: current_state
365  character(len=*), intent(in) :: name
366  type(component_field_value_type), intent(out) :: field_value
367 
368  integer :: i
369 
370  if (name .eq. "wmax") then
371  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
372  current_state%local_grid%size(x_index)))
373  field_value%real_2d_array(:,:)=wmax(:,:)
374  else if (name .eq. "wmin") then
375  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
376  current_state%local_grid%size(x_index)))
377  field_value%real_2d_array(:,:)=wmin(:,:)
378  else if (name .eq. 'reske') then
379  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
380  current_state%local_grid%size(x_index)))
381  field_value%real_2d_array(:,:)=reske(:,:)
382  else if (name .eq. "qlmax") then
383  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
384  current_state%local_grid%size(x_index)))
385  field_value%real_2d_array(:,:)=qlmax(:,:)
386  !else if (name .eq. "hqlmax") then
387  ! allocate(field_value%real_2d_array(current_state%local_grid%size(Y_INDEX), &
388  ! current_state%local_grid%size(X_INDEX)))
389  ! field_value%real_2d_array(:,:)=hqlmax(:,:)
390  else if (name .eq. "cltop") then
391  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
392  current_state%local_grid%size(x_index)))
393  field_value%real_2d_array(:,:)=cltop(:,:)
394  else if (name .eq. "clbas") then
395  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
396  current_state%local_grid%size(x_index)))
397  field_value%real_2d_array(:,:)=clbas(:,:)
398  else if (name .eq. "vwp") then
399  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
400  current_state%local_grid%size(x_index)))
401  field_value%real_2d_array(:,:)=vwp(:,:)
402  else if (name .eq. "lwp") then
403  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
404  current_state%local_grid%size(x_index)))
405  field_value%real_2d_array(:,:)=lwp(:,:)
406  else if (name .eq. "rwp") then
407  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
408  current_state%local_grid%size(x_index)))
409  field_value%real_2d_array(:,:)=rwp(:,:)
410  else if (name .eq. "iwp") then
411  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
412  current_state%local_grid%size(x_index)))
413  field_value%real_2d_array(:,:)=iwp(:,:)
414  else if (name .eq. "swp") then
415  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
416  current_state%local_grid%size(x_index)))
417  field_value%real_2d_array(:,:)=swp(:,:)
418  else if (name .eq. "gwp") then
419  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
420  current_state%local_grid%size(x_index)))
421  field_value%real_2d_array(:,:)=gwp(:,:)
422  else if (name .eq. "tot_iwp") then
423  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
424  current_state%local_grid%size(x_index)))
425  field_value%real_2d_array(:,:)=tot_iwp(:,:)
426  else if (name .eq. "senhf") then
427  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
428  current_state%local_grid%size(x_index)))
429  field_value%real_2d_array(:,:)=senhf(:,:)
430  else if (name .eq. "lathf") then
431  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
432  current_state%local_grid%size(x_index)))
433  field_value%real_2d_array(:,:)=lathf(:,:)
434  end if
435 
Here is the caller graph for this function:

◆ initialisation_callback()

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

Definition at line 64 of file scalar_diagnostics.F90.

65  type(model_state_type), target, intent(inout) :: current_state
66 
67  integer :: k
68  integer :: y_size_local, x_size_local
69 
70  ! qlcrit declared in the config file
71  qlcrit=options_get_real(current_state%options_database, "qlcrit")
72 
73  y_size_local = current_state%local_grid%size(y_index)
74  x_size_local = current_state%local_grid%size(x_index)
75 
76  ! allocate scalar diagnostics as 2-D fields (horizontal slices) that
77  ! are published to the ioserver so that the ioserver can do the manipulation
78  ! to obtain the scalar field
79  allocate(wmax(y_size_local, x_size_local), wmin(y_size_local, x_size_local), &
80  reske(y_size_local, x_size_local), &
81  senhf(y_size_local, x_size_local),lathf(y_size_local, x_size_local))
82  ! allocate the 1d velocity arrays for the kinetic energy calc
83  allocate(ww_prime_res(current_state%local_grid%size(z_index)), &
84  uu_prime_res(current_state%local_grid%size(z_index)), &
85  vv_prime_res(current_state%local_grid%size(z_index)))
86  ww_prime_res(:) = 0.0
87  uu_prime_res(:) = 0.0
88  vv_prime_res(:) = 0.0
89 
90  if (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) then
91  iqv = current_state%water_vapour_mixing_ratio_index
92  iql = current_state%liquid_water_mixing_ratio_index
93  allocate(vwp(y_size_local, x_size_local), lwp(y_size_local, x_size_local), &
94  qlmax(y_size_local, x_size_local), hqlmax(y_size_local, x_size_local), &
95  cltop(y_size_local, x_size_local), clbas(y_size_local, x_size_local))
96  allocate(cloud_content(current_state%local_grid%size(z_index)))
97  ! allocate other hydrometeors. Allocation dependent on index being set in
98  ! appropriate microphysics scheme (see casim component from example)
99  if (current_state%rain_water_mixing_ratio_index > 0) then
100  iqr = current_state%rain_water_mixing_ratio_index
101  allocate(rwp(y_size_local, x_size_local))
102  endif
103  if (current_state%ice_water_mixing_ratio_index > 0) then
104  iqi = current_state%ice_water_mixing_ratio_index
105  allocate(iwp(y_size_local, x_size_local))
106  allocate(tot_iwp(y_size_local, x_size_local))
107  endif
108  if (current_state%snow_water_mixing_ratio_index > 0) then
109  iqs = current_state%snow_water_mixing_ratio_index
110  allocate(swp(y_size_local, x_size_local))
111  endif
112  if (current_state%graupel_water_mixing_ratio_index > 0) then
113  iqg = current_state%graupel_water_mixing_ratio_index
114  allocate(gwp(y_size_local, x_size_local))
115  endif
116  endif
117 
118  allocate(dz_rhon_fac(current_state%local_grid%size(z_index)))
119  do k=2, current_state%local_grid%size(z_index)
120  ! used in the water path calculation
121  dz_rhon_fac(k)=current_state%global_grid%configuration%vertical%dz(k)*&
122  current_state%global_grid%configuration%vertical%rhon(k)
123  end do
124 
Here is the caller graph for this function:

◆ scalar_diagnostics_get_descriptor()

type(component_descriptor_type) function, public scalar_diagnostics_mod::scalar_diagnostics_get_descriptor

Provides the component descriptor for the core to register.

Returns
The descriptor describing this component

Definition at line 34 of file scalar_diagnostics.F90.

35  scalar_diagnostics_get_descriptor%name="scalar_diagnostics"
36  scalar_diagnostics_get_descriptor%version=0.1
37 
38  scalar_diagnostics_get_descriptor%initialisation=>initialisation_callback
39  scalar_diagnostics_get_descriptor%timestep=>timestep_callback
40 
41  scalar_diagnostics_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
42  scalar_diagnostics_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
43  allocate(scalar_diagnostics_get_descriptor%published_fields(16))
44 
45  scalar_diagnostics_get_descriptor%published_fields(1)="vwp"
46  scalar_diagnostics_get_descriptor%published_fields(2)="lwp"
47  scalar_diagnostics_get_descriptor%published_fields(3)="qlmax"
48  scalar_diagnostics_get_descriptor%published_fields(4)="hqlmax"
49  scalar_diagnostics_get_descriptor%published_fields(5)="cltop"
50  scalar_diagnostics_get_descriptor%published_fields(6)="clbas"
51  scalar_diagnostics_get_descriptor%published_fields(7)="wmax"
52  scalar_diagnostics_get_descriptor%published_fields(8)="wmin"
53  scalar_diagnostics_get_descriptor%published_fields(9)="senhf"
54  scalar_diagnostics_get_descriptor%published_fields(10)="lathf"
55  scalar_diagnostics_get_descriptor%published_fields(11)="rwp"
56  scalar_diagnostics_get_descriptor%published_fields(12)="iwp"
57  scalar_diagnostics_get_descriptor%published_fields(13)="swp"
58  scalar_diagnostics_get_descriptor%published_fields(14)="gwp"
59  scalar_diagnostics_get_descriptor%published_fields(15)="tot_iwp"
60  scalar_diagnostics_get_descriptor%published_fields(16)="reske"
61 
Here is the call graph for this function:

◆ timestep_callback()

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

Include ice if present.

Definition at line 127 of file scalar_diagnostics.F90.

128  type(model_state_type), target, intent(inout) :: current_state
129 
130  integer :: k, i
131  integer :: current_y_index, current_x_index, target_x_index, target_y_index
132 
133  current_y_index=current_state%column_local_y
134  current_x_index=current_state%column_local_x
135  target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
136  target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
137 
138  if (current_state%first_timestep_column) then
139  ! maximum vertical velocity for each column
140  wmax(:,:)=0.0
141  ! minimum vertical velocity for each column
142  wmin(:,:)=0.0
143  ! resolved ke
144  reske(:,:) = 0.0
145 
146  if (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) then
147  ! maximum liquid water content in a column
148  qlmax(:,:)=0.0
149  ! the height of the maximum liquid water content in a column
150  hqlmax(:,:)=0.0
151  ! cloud top height where liqud water content is greater than qlcrit
152  cltop(:,:)=0.0
153  ! minimum cloud base where liquid water content is greater than qlcrit
154  clbas(:,:)=0.0
155  ! water vapour path for each column
156  vwp(:,:)=0.0
157  ! liquid water path for each column
158  lwp(:,:)=0.0
159  ! rain water path for each column
160  if (current_state%rain_water_mixing_ratio_index > 0) rwp(:,:)=0.0
161  ! ice water path for each column
162  if (current_state%ice_water_mixing_ratio_index > 0) then
163  iwp(:,:)=0.0
164  ! total ice water path (iwp + swp + gwp) for each column
165  tot_iwp(:,:)=0.0
166  endif
167  ! snow water path for each column
168  if (current_state%snow_water_mixing_ratio_index > 0) swp(:,:)=0.0
169  ! graupel water path for each column
170  if (current_state%graupel_water_mixing_ratio_index > 0) gwp(:,:)=0.0
171  endif
172 
173  ! surface sensible heat flux
174  senhf(:,:)=0.0
175  ! surface latent heat flux
176  lathf(:,:)=0.0
177  end if
178 
179  if (.not. current_state%halo_column) then
180  ! maximum and minimum vertical velocity in each column
181  wmax(target_y_index, target_x_index)=maxval(current_state%w%data(:, current_y_index, current_x_index))
182  wmin(target_y_index, target_x_index)=minval(current_state%w%data(:, current_y_index, current_x_index))
183 
184  ! work out the column resolved ww, uu, vv
185  ww_prime_res(:) = &
186  (current_state%w%data(:,current_state%column_local_y,current_state%column_local_x)**2.)
187  if (allocated(current_state%global_grid%configuration%vertical%olubar)) then
188  uu_prime_res(:) = &
189  ((current_state%u%data(:,current_state%column_local_y,current_state%column_local_x) &
190  - (current_state%global_grid%configuration%vertical%olubar(:) - current_state%ugal))**2.)
191  else
192  uu_prime_res(:) = 0.0
193  endif
194  if (allocated(current_state%global_grid%configuration%vertical%olvbar)) then
195  vv_prime_res(:) = &
196  ((current_state%v%data(:,current_state%column_local_y,current_state%column_local_x) &
197  - (current_state%global_grid%configuration%vertical%olvbar(:) - current_state%vgal))**2.)
198  else
199  vv_prime_res(:) = 0.0
200  endif
201  ! use column resolved ww, uu, vv to derive total resolved KE for column
202  do k=2, current_state%local_grid%size(z_index)-1
203  reske(target_y_index, target_x_index) = reske(target_y_index, target_x_index) + &
204  (uu_prime_res(k) + uu_prime_res(k+1) &
205  + vv_prime_res(k) + vv_prime_res(k+1) &
206  + 2_default_precision * ww_prime_res(k))*0.25_default_precision* &
207  current_state%global_grid%configuration%vertical%dzn(k+1)* &
208  current_state%global_grid%configuration%vertical%rho(k)
209  enddo
210  ! divide reske by altitude to make it column mean (as in the LEM)
211  reske(target_y_index, target_x_index) = reske(target_y_index, target_x_index)/ &
212  current_state%global_grid%configuration%vertical%z(current_state%local_grid%size(z_index))
213  ! subke is derive in the subgrid_profile_diagnostics component
214 
215  if (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) then
216  !
217  ! calculate the lwc maximum and height of lwc max for each column
218  !
219  if (current_state%liquid_water_mixing_ratio_index .gt. 0 .and. &
220  current_state%number_q_fields .ge. current_state%liquid_water_mixing_ratio_index) then
221  qlmax(target_y_index, target_x_index) = &
222  maxval(current_state%q(current_state%liquid_water_mixing_ratio_index)%data &
223  (:,current_y_index, current_x_index))
224  !hqlmax(current_y_index, current_x_index) = &
225  ! current_state%global_grid%configuration%vertical% &
226  ! zn(maxloc(current_state%q(current_state%liquid_water_mixing_ratio_index)%data &
227  ! (:,current_y_index, current_x_index)))
228 
229 
230  ! calculate the cloud top maximum and minimum for each column
231  !
232  cloud_content(:) = current_state%q(current_state%liquid_water_mixing_ratio_index)%data( &
233  :,current_y_index, current_x_index)
235  if (current_state%ice_water_mixing_ratio_index .gt. 0) &
236  cloud_content(:) = cloud_content(:) + current_state%q(current_state%ice_water_mixing_ratio_index)%data( &
237  :,current_y_index, current_x_index)
238  do k = 2, current_state%local_grid%size(z_index)
239  if (cloud_content(k) .gt. qlcrit) then
240  cltop(target_y_index, target_x_index) = &
241  current_state%global_grid%configuration%vertical%zn(k)
242  endif
243  if (cloud_content(current_state%local_grid%size(z_index)+1-k) .gt. qlcrit) then
244  clbas(target_y_index, target_x_index) = &
245  current_state%global_grid%configuration%vertical%zn(current_state%local_grid%size(z_index)+1-k)
246  end if
247  enddo ! k loop over height
248  endif
249  !
250  ! calculate the vapour and liquid water path
251  !
252  if (current_state%water_vapour_mixing_ratio_index .gt. 0 .and. &
253  current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index) then
254  do k = 2, current_state%local_grid%size(z_index)
255  vwp(target_y_index, target_x_index)=vwp(target_y_index, target_x_index) &
256  +dz_rhon_fac(k)*current_state%q(current_state%water_vapour_mixing_ratio_index)%data(k, &
257  current_y_index, current_x_index)
258  lwp(target_y_index, target_x_index)=lwp(target_y_index, target_x_index) &
259  +dz_rhon_fac(k)*current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
260  current_y_index, current_x_index)
261  enddo
262  if (current_state%rain_water_mixing_ratio_index > 0) then
263  do k = 2, current_state%local_grid%size(z_index)
264  rwp(target_y_index, target_x_index)=rwp(target_y_index, target_x_index) &
265  +dz_rhon_fac(k)*current_state%q(iqr)%data(k,current_y_index, current_x_index)
266  enddo
267  endif
268  if (current_state%ice_water_mixing_ratio_index > 0) then
269  do k = 2, current_state%local_grid%size(z_index)
270  iwp(target_y_index, target_x_index)=iwp(target_y_index, target_x_index) &
271  +dz_rhon_fac(k)*current_state%q(iqi)%data(k,current_y_index, current_x_index)
272  enddo
273  tot_iwp(target_y_index, target_x_index)=tot_iwp(target_y_index, target_x_index)+ &
274  iwp(target_y_index, target_x_index)
275  endif
276  if (current_state%snow_water_mixing_ratio_index > 0) then
277  do k = 2, current_state%local_grid%size(z_index)
278  swp(target_y_index, target_x_index)=swp(target_y_index, target_x_index) &
279  +dz_rhon_fac(k)*current_state%q(iqs)%data(k,current_y_index, current_x_index)
280  enddo
281  tot_iwp(target_y_index, target_x_index)=tot_iwp(target_y_index, target_x_index)+ &
282  swp(target_y_index, target_x_index)
283  endif
284  if (current_state%graupel_water_mixing_ratio_index > 0) then
285  do k = 2, current_state%local_grid%size(z_index)
286  gwp(target_y_index, target_x_index)=gwp(target_y_index, target_x_index) &
287  +dz_rhon_fac(k)*current_state%q(iqg)%data(k,current_y_index, current_x_index)
288  enddo
289  tot_iwp(target_y_index, target_x_index)=tot_iwp(target_y_index, target_x_index)+ &
290  gwp(target_y_index, target_x_index)
291  endif
292  end if
293  end if
294 
295  ! surface flux diagnostics
296  if (current_state%use_surface_boundary_conditions) then
297  if (current_state%water_vapour_mixing_ratio_index .gt. 0 .and. &
298  current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index) then
299  lathf(target_y_index, target_x_index)= (current_state%diff_coefficient%data(1, current_y_index, current_x_index) * &
300  current_state%global_grid%configuration%vertical%rdzn(2) * &
301  (current_state%q(current_state%water_vapour_mixing_ratio_index)%data(1,current_y_index,current_x_index) - &
302  current_state%q(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index))) &
303  * rlvap * current_state%global_grid%configuration%vertical%rhon(1)
304  endif
305 
306  if (current_state%th%active) then
307  senhf(target_y_index, target_x_index)=(current_state%diff_coefficient%data(1, current_y_index, current_x_index) &
308  * current_state%global_grid%configuration%vertical%rdzn(2) &
309  * (current_state%th%data(1, current_y_index, current_x_index) &
310  - current_state%th%data(2, current_y_index, current_x_index) &
311  - current_state%global_grid%configuration%vertical%dthref(1))) &
312  * current_state%global_grid%configuration%vertical%rhon(1)*cp
313  endif
314  endif
315  end if
Here is the caller graph for this function:

Variable Documentation

◆ clbas

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::clbas
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ cloud_content

real(kind=default_precision), dimension(:), allocatable scalar_diagnostics_mod::cloud_content
private

Definition at line 21 of file scalar_diagnostics.F90.

◆ cltop

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::cltop
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ dz_rhon_fac

real(kind=default_precision), dimension(:), allocatable scalar_diagnostics_mod::dz_rhon_fac
private

Definition at line 20 of file scalar_diagnostics.F90.

20  real(kind=default_precision), dimension(:), allocatable :: dz_rhon_fac

◆ gwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::gwp
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ hqlmax

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::hqlmax
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ iqg

integer scalar_diagnostics_mod::iqg =0
private

Definition at line 18 of file scalar_diagnostics.F90.

◆ iqi

integer scalar_diagnostics_mod::iqi =0
private

Definition at line 18 of file scalar_diagnostics.F90.

◆ iql

integer scalar_diagnostics_mod::iql =0
private

Definition at line 18 of file scalar_diagnostics.F90.

◆ iqr

integer scalar_diagnostics_mod::iqr =0
private

Definition at line 18 of file scalar_diagnostics.F90.

◆ iqs

integer scalar_diagnostics_mod::iqs =0
private

Definition at line 18 of file scalar_diagnostics.F90.

◆ iqv

integer scalar_diagnostics_mod::iqv =0
private

Definition at line 18 of file scalar_diagnostics.F90.

◆ iwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::iwp
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ lathf

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::lathf
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ lwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::lwp
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ ncl_col

integer scalar_diagnostics_mod::ncl_col
private

Definition at line 18 of file scalar_diagnostics.F90.

◆ qlcrit

real(kind=default_precision) scalar_diagnostics_mod::qlcrit
private

Definition at line 23 of file scalar_diagnostics.F90.

23  real(kind=default_precision) :: qlcrit

◆ qlmax

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::qlmax
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ reske

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::reske
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ rwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::rwp
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ senhf

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::senhf
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ swp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::swp
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ tot_iwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::tot_iwp
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ total_points

integer scalar_diagnostics_mod::total_points
private

Definition at line 18 of file scalar_diagnostics.F90.

18  integer :: total_points, ncl_col, iqv=0, iql=0, iqr=0, iqi=0, iqs=0, &
19  iqg=0

◆ uu_prime_res

real(kind=default_precision), dimension(:), allocatable scalar_diagnostics_mod::uu_prime_res
private

Definition at line 21 of file scalar_diagnostics.F90.

◆ vv_prime_res

real(kind=default_precision), dimension(:), allocatable scalar_diagnostics_mod::vv_prime_res
private

Definition at line 21 of file scalar_diagnostics.F90.

◆ vwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::vwp
private

Definition at line 24 of file scalar_diagnostics.F90.

24  real(kind=default_precision), dimension(:,:), allocatable :: vwp, lwp, wmax, wmin, &
25  qlmax, hqlmax, cltop, clbas, senhf, lathf, rwp, iwp, swp, gwp, tot_iwp, &
26  reske

◆ wmax

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::wmax
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ wmin

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::wmin
private

Definition at line 24 of file scalar_diagnostics.F90.

◆ ww_prime_res

real(kind=default_precision), dimension(:), allocatable scalar_diagnostics_mod::ww_prime_res
private

Definition at line 21 of file scalar_diagnostics.F90.

21  real(kind=default_precision), dimension(:), allocatable :: ww_prime_res, uu_prime_res, &
22  vv_prime_res, cloud_content
grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
science_constants_mod::cp
real(kind=default_precision), public cp
Definition: scienceconstants.F90:13
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
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
datadefn_mod::default_precision
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17