MONC
Functions/Subroutines | Variables
conditional_diagnostics_column_mod Module Reference

Conditionally averaged diagnostics, Part 1 of 2. More...

Functions/Subroutines

type(component_descriptor_type) function, public conditional_diagnostics_column_get_descriptor ()
 Provides the component descriptor for the core to register. More...
 
subroutine initialisation_callback (current_state)
 
subroutine finalisation_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 iqv
 
integer iql
 
integer iqr
 
integer iqi
 
integer iqs
 
integer iqg
 
real(kind=default_precision), dimension(:,:,:), allocatable, public conddiags_tot
 
real(kind=default_precision) qlcrit
 
real(kind=default_precision) qicrit
 
real(kind=default_precision) qpptcrit
 
real(kind=default_precision) vpptcrit
 
real(kind=default_precision) thvprcrit
 
real(kind=default_precision) wsdwncrit
 
real(kind=default_precision) wsupcrit
 
real(kind=default_precision) wupcrit
 
real(kind=default_precision) wdwncrit
 
integer, public ncond
 
integer, public ndiag
 
integer x_size
 
integer y_size
 
real(kind=default_precision), public gpts_total
 
logical thv_from_th_with_liqice
 
logical l_qi_qr_qs_qg
 
logical l_do_near
 
character(len=string_length), dimension(23) master_conditions_list
 
character(len=string_length), dimension(31) master_diagnostics_list
 
character(len=string_length), dimension(:), allocatable, public cond_request
 
character(len=string_length), dimension(:), allocatable, public cond_long
 
character(len=string_length), dimension(:), allocatable, public diag_request
 
character(len=string_length), dimension(:), allocatable, public diag_long
 
integer, dimension(:), allocatable diag_locations
 
integer, dimension(:), allocatable cond_locations
 
integer, public requested_area
 
integer diagnostic_generation_frequency
 

Detailed Description

Conditionally averaged diagnostics, Part 1 of 2.

Function/Subroutine Documentation

◆ conditional_diagnostics_column_get_descriptor()

type(component_descriptor_type) function, public conditional_diagnostics_column_mod::conditional_diagnostics_column_get_descriptor

Provides the component descriptor for the core to register.

Returns
The descriptor describing this component

Definition at line 68 of file conditional_diagnostics_column.F90.

69  conditional_diagnostics_column_get_descriptor%name="conditional_diagnostics_column"
70  conditional_diagnostics_column_get_descriptor%version=0.1
71  conditional_diagnostics_column_get_descriptor%initialisation=>initialisation_callback
72  conditional_diagnostics_column_get_descriptor%timestep=>timestep_callback
73  conditional_diagnostics_column_get_descriptor%finalisation=>finalisation_callback
74 
75  conditional_diagnostics_column_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
76  conditional_diagnostics_column_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
77  allocate(conditional_diagnostics_column_get_descriptor%published_fields(1))
78 
79  conditional_diagnostics_column_get_descriptor%published_fields(1)="CondDiags_tot"
Here is the call graph for this function:

◆ field_information_retrieval_callback()

subroutine conditional_diagnostics_column_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 768 of file conditional_diagnostics_column.F90.

769  type(model_state_type), target, intent(inout) :: current_state
770  character(len=*), intent(in) :: name
771  type(component_field_information_type), intent(out) :: field_information
772 
773 
774  field_information%field_type=component_array_field_type
775  field_information%number_dimensions=3
776  field_information%dimension_sizes(1) = current_state%local_grid%size(z_index)
777  field_information%dimension_sizes(2) = ncond * 2
778  field_information%dimension_sizes(3) = ndiag
779 
780  field_information%data_type = component_double_data_type
781 
782  if (name .eq. "CondDiags_tot") then
783  field_information%enabled = allocated(conddiags_tot)
784  end if
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine conditional_diagnostics_column_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 793 of file conditional_diagnostics_column.F90.

794  type(model_state_type), target, intent(inout) :: current_state
795  character(len=*), intent(in) :: name
796  type(component_field_value_type), intent(out) :: field_value
797 
798  if (name .eq. "CondDiags_tot") then
799  allocate(field_value%real_3d_array(current_state%local_grid%size(z_index), &
800  ncond * 2, &
801  ndiag))
802  field_value%real_3d_array(:,:,:) = conddiags_tot(:,:,:)
803  end if
Here is the caller graph for this function:

◆ finalisation_callback()

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

Definition at line 360 of file conditional_diagnostics_column.F90.

361  type(model_state_type), target, intent(inout) :: current_state
362 
363  if (allocated(conddiags_tot)) deallocate(conddiags_tot)
364  if (allocated(cond_request)) deallocate(cond_request)
365  if (allocated(cond_long)) deallocate(cond_long)
366  if (allocated(cond_locations)) deallocate(cond_locations)
367  if (allocated(diag_request)) deallocate(diag_request)
368  if (allocated(diag_long)) deallocate(diag_long)
369  if (allocated(diag_locations)) deallocate(diag_locations)
370 
Here is the caller graph for this function:

◆ initialisation_callback()

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

This is section contains the full list of available conditions and

Check for requested diagnostics.

Record which conditions and diagnostics have been requested

Definition at line 83 of file conditional_diagnostics_column.F90.

84  type(model_state_type), target, intent(inout) :: current_state
85 
86  integer :: total_points, inc
87 
88  !----------------------------------------------------------------------
90  ! diagnostics.
91  !
92  ! If you add more, you will need to update some of these:
93  ! 1. master_conditions_list (and its declared size)
94  ! 2. master_diagnostics_list (and its declared size)
95  ! 3. cond_locations CASE checks
96  ! 4. diag_locations CASE checks
97  ! 5. local_diag value
98  ! 6. l_cond logical test
99  ! 7. add short name to global_config default list
100  !
101  ! TAKE CARE TO ENSURE that indices correspond directly between:
102  ! 1. master_conditions_list / cond_locations
103  ! 2. master_diagnostics_list / diag_locations / local_diag
104  !
105  !----------------------------------------------------------------------
106 
107  master_conditions_list(1) = "all points"
108  master_conditions_list(2) = "buoyant updraft"
109  master_conditions_list(3) = "buoyant cloudy updraft"
110  master_conditions_list(4) = "near buoyant cloudy updraft"
111  master_conditions_list(5) = "cloudy"
112  master_conditions_list(6) = "cloudy updraft"
113  master_conditions_list(7) = "cloudy downdraft"
114  master_conditions_list(8) = "strong updraft"
115  master_conditions_list(9) = "strong downdraft"
116  master_conditions_list(10) = "updraft"
117  master_conditions_list(11) = "downdraft"
118  master_conditions_list(12) = "liquid cloudy updraft"
119  master_conditions_list(13) = "liquid cloudy downdraft"
120  master_conditions_list(14) = "hydrometeors"
121  master_conditions_list(15) = "cloud liquid water"
122  master_conditions_list(16) = "cloud ice"
123  master_conditions_list(17) = "precipitating downdraft"
124  master_conditions_list(18) = "large precipitating downdraft"
125  master_conditions_list(19) = "precipitating strong downdraft"
126  master_conditions_list(20) = "saturated updraft"
127  master_conditions_list(21) = "buoyant saturated updraft"
128  master_conditions_list(22) = "cloud-free precipitating downdraft"
129  master_conditions_list(23) = "cloud-free large precipitating downdraft"
130 
131  master_diagnostics_list(1) = "Fractional area"
132  master_diagnostics_list(2) = "Vertical velocity"
133  master_diagnostics_list(3) = "Vertical velocity variance"
134  master_diagnostics_list(4) = "Potential temperature"
135  master_diagnostics_list(5) = "w * th"
136  master_diagnostics_list(6) = "Potential temperature anomalies"
137  master_diagnostics_list(7) = "Flux of potential temperature"
138  master_diagnostics_list(8) = "Virtual potential temperature anomalies"
139  master_diagnostics_list(9) = "Flux of virtual potential temperature"
140  master_diagnostics_list(10) = "Variance of potential temperature"
141  master_diagnostics_list(11) = "Diff_coeff * dth/dz"
142  master_diagnostics_list(12) = "w**3"
143  master_diagnostics_list(13) = "Relative humidity"
144  master_diagnostics_list(14) = "Zonal velocity"
145  master_diagnostics_list(15) = "Meridional velovity"
146  master_diagnostics_list(16) = "w * u"
147  master_diagnostics_list(17) = "w * v"
148  master_diagnostics_list(18) = "Vis_coeff * du/dz"
149  master_diagnostics_list(19) = "Vis_coeff * dv/dz"
150  master_diagnostics_list(20) = "Temperature"
151  master_diagnostics_list(21) = "th+lvap*qv/cp"
152  master_diagnostics_list(22) = "Anomalies of THL"
153  master_diagnostics_list(23) = "Variance of THL"
154  master_diagnostics_list(24) = "qv+ql+qi"
155  master_diagnostics_list(25) = "Anomalies of QVLI"
156  master_diagnostics_list(26) = "Variance of QVLI"
157  master_diagnostics_list(27) = "qr+qs+qg"
158  master_diagnostics_list(28) = "Anomalies of QRSG"
159  master_diagnostics_list(29) = "Variance of QRSG"
160  master_diagnostics_list(30) = "Flux of QVLI"
161  master_diagnostics_list(31) = "Flux of QRSG"
162 
163 
164  total_points=(current_state%local_grid%size(y_index) * current_state%local_grid%size(x_index))
165 
166  if (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) then
167  iqv =get_q_index(standard_q_names%VAPOUR, 'profile_diags')
168  iql =get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'profile_diags')
169  iqi =get_q_index(standard_q_names%ICE_MASS, 'profile_diags')
170  iqr =get_q_index(standard_q_names%RAIN_MASS, 'profile_diags')
171  iqs =get_q_index(standard_q_names%SNOW_MASS, 'profile_diags')
172  iqg =get_q_index(standard_q_names%GRAUPEL_MASS, 'profile_diags')
173 
174  qlcrit =options_get_real(current_state%options_database, "qlcrit")
175  qicrit =options_get_real(current_state%options_database, "qicrit")
176  qpptcrit =options_get_real(current_state%options_database, "qpptcrit")
177  vpptcrit =options_get_real(current_state%options_database, "vpptcrit")
178  endif
179 
180  !calculate qi, qr, qs and qg if and only if number_q_fields=11
181  l_qi_qr_qs_qg = (.not. current_state%passive_q .and. current_state%number_q_fields .ge. 11)
182 
183  x_size =options_get_integer(current_state%options_database, "x_size")
184  y_size =options_get_integer(current_state%options_database, "y_size")
185  thvprcrit =options_get_real(current_state%options_database, "thvprcrit")
186  wsdwncrit =options_get_real(current_state%options_database, "wSdwncrit")
187  wsupcrit =options_get_real(current_state%options_database, "wSupcrit")
188  wupcrit =options_get_real(current_state%options_database, "wupcrit")
189  wdwncrit =options_get_real(current_state%options_database, "wdwncrit")
190  ncond =options_get_array_size(current_state%options_database, "cond_request")
191  ndiag =options_get_array_size(current_state%options_database, "diag_request")
192  thv_from_th_with_liqice =options_get_logical(current_state%options_database, "thv_from_th_with_liqice")
193 
194  gpts_total = 1.0_default_precision * x_size * y_size
195 
197  if (ncond .lt. 1 .or. ndiag .lt. 1) call log_master_log(log_error, &
198  "When conditional_diagnostics_column are enabled, condition and diagnostic request lists must be provided.")
199 
200  if (current_state%th%active) then
201  allocate(conddiags_tot(current_state%local_grid%size(z_index), &
202  ncond * 2 , &
203  ndiag))
204  endif
205 
206 
207  l_do_near = .false.
208  requested_area = 0
210  allocate(cond_request(ncond))
211  allocate(cond_long(ncond))
212  allocate(cond_locations(ncond))
213  allocate(diag_request(ndiag))
214  allocate(diag_long(ndiag))
215  allocate(diag_locations(ndiag))
216 
217  do inc = 1, ncond
218  cond_request(inc) = trim(options_get_string(current_state%options_database, "cond_request", inc))
219 
220  SELECT CASE ( trim(cond_request(inc)) )
221  CASE DEFAULT
222  call log_master_log(log_error, &
223  "Condition '"//trim(cond_request(inc))//"' has not been set up.")
224  CASE ("ALL")
225  cond_locations(inc) = 1
226  CASE ("BYu")
227  cond_locations(inc) = 2
228  CASE ("BCu")
229  cond_locations(inc) = 3
230  CASE ("NrBCu")
231  cond_locations(inc) = 4
232  l_do_near = .true.
233  CASE ("AC")
234  cond_locations(inc) = 5
235  CASE ("ACu")
236  cond_locations(inc) = 6
237  CASE ("ACd")
238  cond_locations(inc) = 7
239  CASE ("WG1")
240  cond_locations(inc) = 8
241  CASE ("WL1")
242  cond_locations(inc) = 9
243  CASE ("ALu")
244  cond_locations(inc) = 10
245  CASE ("ALd")
246  cond_locations(inc) = 11
247  CASE ("CLu")
248  cond_locations(inc) = 12
249  CASE ("CLd")
250  cond_locations(inc) = 13
251  CASE ("AH")
252  cond_locations(inc) = 14
253  CASE ("AL")
254  cond_locations(inc) = 15
255  CASE ("AI")
256  cond_locations(inc) = 16
257  CASE ("PPd")
258  cond_locations(inc) = 17
259  CASE ("VPd")
260  cond_locations(inc) = 18
261  CASE ("PVd")
262  cond_locations(inc) = 19
263  CASE ("MO")
264  cond_locations(inc) = 20
265  CASE ("BM")
266  cond_locations(inc) = 21
267  CASE ("AA")
268  cond_locations(inc) = 22
269  CASE ("AV")
270  cond_locations(inc) = 23
271  END SELECT
272  end do ! loop over requested conditions
273 
274 
275  do inc = 1, ndiag
276  diag_request(inc) = trim(options_get_string(current_state%options_database, "diag_request", inc))
277  if ( trim(diag_request(inc)) .eq. "area" ) requested_area = inc
278 
279  SELECT CASE ( trim(diag_request(inc)) )
280  CASE DEFAULT
281  call log_master_log(log_error, &
282  "Diagnostic '"//trim(diag_request(inc))//"' has not been set up.")
283  CASE ("area")
284  diag_locations(inc) = 1
285  CASE ("W")
286  diag_locations(inc) = 2
287  CASE ("W2")
288  diag_locations(inc) = 3
289  CASE ("TH")
290  diag_locations(inc) = 4
291  CASE ("WTH")
292  diag_locations(inc) = 5
293  CASE ("THP")
294  diag_locations(inc) = 6
295  CASE ("WTHP")
296  diag_locations(inc) = 7
297  CASE ("THVP")
298  diag_locations(inc) = 8
299  CASE ("WTHVP")
300  diag_locations(inc) = 9
301  CASE ("THP2")
302  diag_locations(inc) = 10
303  CASE ("WTHSG")
304  diag_locations(inc) = 11
305  CASE ("W3")
306  diag_locations(inc) = 12
307  CASE ("RH")
308  diag_locations(inc) = 13
309  CASE ("U")
310  diag_locations(inc) = 14
311  CASE ("V")
312  diag_locations(inc) = 15
313  CASE ("WU")
314  diag_locations(inc) = 16
315  CASE ("WV")
316  diag_locations(inc) = 17
317  CASE ("WUSG")
318  diag_locations(inc) = 18
319  CASE ("WVSG")
320  diag_locations(inc) = 19
321  CASE ("TEMP")
322  diag_locations(inc) = 20
323  CASE ("THL")
324  diag_locations(inc) = 21
325  CASE ("THLP")
326  diag_locations(inc) = 22
327  CASE ("THLP2")
328  diag_locations(inc) = 23
329  CASE ("QVLI")
330  diag_locations(inc) = 24
331  CASE ("QVLIP")
332  diag_locations(inc) = 25
333  CASE ("QVLIP2")
334  diag_locations(inc) = 26
335  CASE ("QRSG")
336  diag_locations(inc) = 27
337  CASE ("QRSGP")
338  diag_locations(inc) = 28
339  CASE ("QRSGP2")
340  diag_locations(inc) = 29
341  CASE ("WQVLIP")
342  diag_locations(inc) = 30
343  CASE ("WQRSGP")
344  diag_locations(inc) = 31
345  END SELECT
346  end do ! loop over requested diagnostics
347 
348  if (requested_area == 0) call log_master_log(log_error, &
349  "The diagnostic 'area' must be provided to complete the conditional diagnostics process.")
350 
351  cond_long=master_conditions_list(cond_locations)
352  diag_long=master_diagnostics_list(diag_locations)
353 
354  ! Save the sampling_frequency to force diagnostic calculation on select time steps
355  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
356 
Here is the caller graph for this function:

◆ timestep_callback()

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

Reset sums on this process

Decide if conditions are appropriate to proceed with calculations

Begin the calculations Loop over levels

work out all potential terms at this location

end of potential term calculations

Store the diagnostics

Determine status of requested conditions

Near it, but not it.

Record the value for this diagnostic meeting this condition and its complement by adding to cumulative sum

Definition at line 374 of file conditional_diagnostics_column.F90.

375  type(model_state_type), target, intent(inout) :: current_state
376 
377  logical :: l_cond ! condition at local position (k,j,i)
378  real(kind=default_precision), dimension(size(master_diagnostics_list)) :: local_diag ! diagnostics at local position (k,j,i)
379 
380  real(kind=default_precision) :: qi, ql,qv,qli,qsat,qppt,qvl
381  real(kind=default_precision) :: qi_pr, ql_pr,qv_pr,qli_pr,qvli, qvli_pr, qppt_pr
382  real(kind=default_precision) :: thv_pr,thv,exner,tdegk,th_pr,th_pr2
383  real(kind=default_precision) :: tmp_th !=thref+th
384  real(kind=default_precision) :: pottemp !=thref+olzthbar
385  real(kind=default_precision) :: relhum,pmb,wth,wthsg, wthv_pr, wth_pr
386  real(kind=default_precision) :: w_zn ! vertical velocity at theta levels
387  real(kind=default_precision) :: tmp_u,tmp_v ! horizontal veleocities
388  real(kind=default_precision) :: w_zn2, w_zn3, wu, wv, wusg, wvsg, wqvli_pr,wqppt_pr
389  real(kind=default_precision) :: th_h,th_h_pr1,th_h_pr2
390  real(kind=default_precision) :: qvli_pr2, qppt_pr2
391 
392  real(kind=default_precision) :: qv_jip1,qv_jim1,qv_jp1i,qv_jm1i,qv_pr_jip1,qv_pr_jim1,qv_pr_jp1i,qv_pr_jm1i
393  real(kind=default_precision) :: qi_jip1,qi_jim1,qi_jp1i,qi_jm1i
394  real(kind=default_precision) :: ql_jip1,ql_jim1,ql_jp1i,ql_jm1i
395  real(kind=default_precision) :: qli_jip1,qli_jim1,qli_jp1i,qli_jm1i,qli_pr_jip1,qli_pr_jim1,qli_pr_jp1i,qli_pr_jm1i
396  real(kind=default_precision) :: w_zn_jip1,w_zn_jim1,w_zn_jp1i,w_zn_jm1i
397  real(kind=default_precision) :: thv_pr_jip1,thv_pr_jim1,thv_pr_jp1i,thv_pr_jm1i
398  real(kind=default_precision) :: tmp_th_jip1,tmp_th_jim1,tmp_th_jp1i,tmp_th_jm1i
399  real(kind=default_precision) :: th_pr_jip1,th_pr_jim1,th_pr_jp1i,th_pr_jm1i
400 
401  integer :: k, j, i
402  integer :: inc ! loop increment variable
403  integer :: local_y, locaL_x, target_x_index, target_y_index
404 
405  j=current_state%column_local_y
406  i=current_state%column_local_x
407 
409  if (current_state%first_timestep_column) then
410  conddiags_tot(:,:,:) = 0.0_default_precision
411  end if
412 
414  if (current_state%halo_column) return
415  if ( .not. (current_state%th%active .and. &
416  .not. current_state%passive_q .and. &
417  current_state%number_q_fields .gt. 0) ) return
418  if (.not. mod(current_state%timestep, diagnostic_generation_frequency) == 0) return
419 
422  ! Due to above/below level solution dependency, no data will exist at first/last levels,
423  ! but k-dimension size kept equal to size(Z_INDEX) for convenience
424  do k = 2, current_state%local_grid%size(z_index)-1
425 
427  qv = current_state%zq(iqv)%data(k,j,i)
428  ql = current_state%zq(iql)%data(k,j,i)
429  qvl = qv + ql
430  qv_pr = qv - current_state%global_grid%configuration%vertical%olzqbar(k,iqv)
431  ql_pr = ql - current_state%global_grid%configuration%vertical%olzqbar(k,iql)
432 
433  if (.not. l_qi_qr_qs_qg) then
434  qi = 0.0
435  qli = ql
436  qvli = qvl
437  qppt = 0.0
438  qi_pr = 0.0
439  qli_pr = ql_pr
440  qvli_pr = qv_pr + ql_pr
441  qppt_pr = 0.0
442  end if
443 
444  if (l_qi_qr_qs_qg) then
445  qi = current_state%zq(iqi)%data(k,j,i)
446  qli = ql + qi
447  qvli = qv + qli
448  qppt = current_state%zq(iqr)%data(k,j,i) + current_state%zq(iqs)%data(k,j,i)+current_state%zq(iqg)%data(k,j,i)
449  qi_pr = qli - current_state%global_grid%configuration%vertical%olzqbar(k,iqi)
450  qli_pr = ql_pr + qi_pr
451  qvli_pr = qv_pr + qli_pr
452  qppt_pr = current_state%zq(iqr)%data(k,j,i) - current_state%global_grid%configuration%vertical%olzqbar(k,iqr)+ &
453  current_state%zq(iqs)%data(k,j,i) - current_state%global_grid%configuration%vertical%olzqbar(k,iqs)+ &
454  current_state%zq(iqg)%data(k,j,i) - current_state%global_grid%configuration%vertical%olzqbar(k,iqg)
455  end if
456 
457  w_zn = 0.5 * (current_state%zw%data(k,j,i) + current_state%zw%data(k-1,j,i)) !w at theta levels
458  w_zn2 = w_zn * w_zn
459  w_zn3 = w_zn * w_zn * w_zn
460 
461  qvli_pr2 = qvli_pr * qvli_pr
462  qppt_pr2 = qppt_pr * qppt_pr
463  wqvli_pr = w_zn * qvli_pr
464  wqppt_pr = w_zn * qppt_pr
465 
466  tmp_u = 0.5 * (current_state%zu%data(k,j,i-1) + current_state%zu%data(k,j,i)) + current_state%ugal !u at theta points
467  tmp_v = 0.5 * (current_state%zv%data(k,j-1,i) + current_state%zv%data(k,j,i)) + current_state%vgal !v at theta points
468  wu = w_zn * tmp_u
469  wv = w_zn * tmp_v
470  tmp_th = current_state%global_grid%configuration%vertical%thref(k) + current_state%zth%data(k,j,i)
471  wth = w_zn * tmp_th
472  wthsg = -0.5 * current_state%diff_coefficient%data(k,j,i) * (current_state%zth%data(k+1,j,i) + &
473  current_state%global_grid%configuration%vertical%thref(k+1) - current_state%zth%data(k,j,i) - &
474  current_state%global_grid%configuration%vertical%thref(k)) * &
475  current_state%global_grid%configuration%vertical%rdzn(k+1) &
476  -0.5 * current_state%diff_coefficient%data(k-1,j,i) * (current_state%zth%data(k,j,i) + &
477  current_state%global_grid%configuration%vertical%thref(k) - current_state%zth%data(k-1,j,i) - &
478  current_state%global_grid%configuration%vertical%thref(k-1)) * &
479  current_state%global_grid%configuration%vertical%rdzn(k)
480 
481  wusg = -0.25 * current_state%vis_coefficient%data(k ,j,i-1) * &
482  (current_state%zu%data(k+1,j,i-1) - current_state%zu%data(k ,j,i-1)) * &
483  current_state%global_grid%configuration%vertical%rdzn(k+1) &
484  -0.25 * current_state%vis_coefficient%data(k ,j,i ) * &
485  (current_state%zu%data(k+1,j,i ) - current_state%zu%data(k ,j,i )) * &
486  current_state%global_grid%configuration%vertical%rdzn(k+1) &
487  -0.25 * current_state%vis_coefficient%data(k-1,j,i-1) * &
488  (current_state%zu%data(k ,j,i-1) - current_state%zu%data(k-1,j,i-1)) * &
489  current_state%global_grid%configuration%vertical%rdzn(k) &
490  -0.25 * current_state%vis_coefficient%data(k-1,j,i ) * &
491  (current_state%zu%data(k ,j,i ) - current_state%zu%data(k-1,j,i )) * &
492  current_state%global_grid%configuration%vertical%rdzn(k)
493 
494  wvsg = -0.25 * current_state%vis_coefficient%data(k ,j-1,i) * &
495  (current_state%zv%data(k+1,j-1,i) - current_state%zv%data(k ,j-1,i)) * &
496  current_state%global_grid%configuration%vertical%rdzn(k+1) &
497  -0.25 * current_state%vis_coefficient%data(k ,j ,i) * &
498  (current_state%zv%data(k+1,j ,i) - current_state%zv%data(k ,j ,i)) * &
499  current_state%global_grid%configuration%vertical%rdzn(k+1) &
500  -0.25 * current_state%vis_coefficient%data(k-1,j-1,i) * &
501  (current_state%zv%data(k ,j-1,i) - current_state%zv%data(k-1,j-1,i)) * &
502  current_state%global_grid%configuration%vertical%rdzn(k) &
503  -0.25 * current_state%vis_coefficient%data(k-1,j ,i) * &
504  (current_state%zv%data(k ,j ,i) - current_state%zv%data(k-1,j ,i)) * &
505  current_state%global_grid%configuration%vertical%rdzn(k)
506 
507  pottemp = current_state%global_grid%configuration%vertical%thref(k) + &
508  current_state%global_grid%configuration%vertical%olzthbar(k)
509 
510  th_pr = current_state%zth%data(k,j,i) - current_state%global_grid%configuration%vertical%olzthbar(k)
511  wth_pr = w_zn * th_pr
512  th_pr2 = th_pr * th_pr
513  exner = current_state%global_grid%configuration%vertical%rprefrcp(k)
514  tdegk = (current_state%global_grid%configuration%vertical%thref(k) + current_state%zth%data(k,j,i)) * exner
515  pmb = current_state%global_grid%configuration%vertical%prefn(k) / 100.
516  qsat = qsaturation(tdegk,pmb)
517  relhum = qv / qsat
518  th_h = tmp_th + rlvap * qv / cp
519  th_h_pr1 = th_pr + rlvap * qv_pr / cp
520  th_h_pr2 = th_h_pr1 * th_h_pr1
521 
522  !calculate values at nearby points
523  if (l_do_near) then
524  qv_jip1 = current_state%zq(iqv)%data(k,j,i+1) !qv at the position j and i+1
525  qv_jim1 = current_state%zq(iqv)%data(k,j,i-1) !qv at the position j and i-1
526  qv_jp1i = current_state%zq(iqv)%data(k,j+1,i) !qv at the position j+1 and i
527  qv_jm1i = current_state%zq(iqv)%data(k,j-1,i) !qv at the position j-1 and i
528  ql_jip1 = current_state%zq(iql)%data(k,j,i+1)
529  ql_jim1 = current_state%zq(iql)%data(k,j,i-1)
530  ql_jp1i = current_state%zq(iql)%data(k,j+1,i)
531  ql_jm1i = current_state%zq(iql)%data(k,j-1,i)
532  qv_pr_jip1 = current_state%zq(iqv)%data(k,j,i+1)-current_state%global_grid%configuration%vertical%olzqbar(k,iqv)
533  qv_pr_jim1 = current_state%zq(iqv)%data(k,j,i-1)-current_state%global_grid%configuration%vertical%olzqbar(k,iqv)
534  qv_pr_jp1i = current_state%zq(iqv)%data(k,j+1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iqv)
535  qv_pr_jm1i = current_state%zq(iqv)%data(k,j-1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iqv)
536 
537  if (.not. l_qi_qr_qs_qg) then
538  qi_jip1 = 0.0
539  qi_jim1 = 0.0
540  qi_jp1i = 0.0
541  qi_jm1i = 0.0
542  qli_jip1 = current_state%zq(iql)%data(k,j,i+1)+0.0
543  qli_jim1 = current_state%zq(iql)%data(k,j,i-1)+0.0
544  qli_jp1i = current_state%zq(iql)%data(k,j+1,i)+0.0
545  qli_jm1i = current_state%zq(iql)%data(k,j-1,i)+0.0
546  qli_pr_jip1 = current_state%zq(iql)%data(k,j,i+1)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+0.0
547  qli_pr_jim1 = current_state%zq(iql)%data(k,j,i-1)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+0.0
548  qli_pr_jp1i = current_state%zq(iql)%data(k,j+1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+0.0
549  qli_pr_jm1i = current_state%zq(iql)%data(k,j-1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+0.0
550  end if
551 
552  if (l_qi_qr_qs_qg) then
553  qi_jip1 = current_state%zq(iqi)%data(k,j,i+1)
554  qi_jim1 = current_state%zq(iqi)%data(k,j,i-1)
555  qi_jp1i = current_state%zq(iqi)%data(k,j+1,i)
556  qi_jm1i = current_state%zq(iqi)%data(k,j-1,i)
557  qli_jip1 = current_state%zq(iql)%data(k,j,i+1)+current_state%zq(iqi)%data(k,j,i+1)
558  qli_jim1 = current_state%zq(iql)%data(k,j,i-1)+current_state%zq(iqi)%data(k,j,i-1)
559  qli_jp1i = current_state%zq(iql)%data(k,j+1,i)+current_state%zq(iqi)%data(k,j+1,i)
560  qli_jm1i = current_state%zq(iql)%data(k,j-1,i)+current_state%zq(iqi)%data(k,j-1,i)
561  qli_pr_jip1 = current_state%zq(iql)%data(k,j,i+1)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+ &
562  current_state%zq(iqi)%data(k,j,i+1)-current_state%global_grid%configuration%vertical%olzqbar(k,iqi)
563  qli_pr_jim1 = current_state%zq(iql)%data(k,j,i-1)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+ &
564  current_state%zq(iqi)%data(k,j,i-1)-current_state%global_grid%configuration%vertical%olzqbar(k,iqi)
565  qli_pr_jp1i = current_state%zq(iql)%data(k,j+1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+ &
566  current_state%zq(iqi)%data(k,j+1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iqi)
567  qli_pr_jm1i = current_state%zq(iql)%data(k,j-1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iql)+ &
568  current_state%zq(iqi)%data(k,j-1,i)-current_state%global_grid%configuration%vertical%olzqbar(k,iqi)
569  end if
570 
571  tmp_th_jip1 = current_state%global_grid%configuration%vertical%thref(k) + current_state%zth%data(k,j,i+1)
572  tmp_th_jim1 = current_state%global_grid%configuration%vertical%thref(k) + current_state%zth%data(k,j,i-1)
573  tmp_th_jp1i = current_state%global_grid%configuration%vertical%thref(k) + current_state%zth%data(k,j+1,i)
574  tmp_th_jm1i = current_state%global_grid%configuration%vertical%thref(k) + current_state%zth%data(k,j-1,i)
575  th_pr_jip1 = current_state%zth%data(k,j,i+1)-current_state%global_grid%configuration%vertical%olzthbar(k)
576  th_pr_jim1 = current_state%zth%data(k,j,i-1)-current_state%global_grid%configuration%vertical%olzthbar(k)
577  th_pr_jp1i = current_state%zth%data(k,j-1,i)-current_state%global_grid%configuration%vertical%olzthbar(k)
578  th_pr_jm1i = current_state%zth%data(k,j+1,i)-current_state%global_grid%configuration%vertical%olzthbar(k)
579  w_zn_jip1 = 0.5*(current_state%zw%data(k,j,i+1)+current_state%zw%data(k-1,j,i+1))
580  w_zn_jim1 = 0.5*(current_state%zw%data(k,j,i-1)+current_state%zw%data(k-1,j,i-1))
581  w_zn_jp1i = 0.5*(current_state%zw%data(k,j+1,i)+current_state%zw%data(k-1,j+1,i))
582  w_zn_jm1i = 0.5*(current_state%zw%data(k,j-1,i)+current_state%zw%data(k-1,j-1,i))
583 
584  if (thv_from_th_with_liqice) then
585  thv = tmp_th*(1.0+0.61*qv-qli)
586  !the term th_pr*(0.61*qv-qli) is very small and can be neglected
587  thv_pr = th_pr+tmp_th*(0.61*qv_pr-qli_pr)+th_pr*(0.61*qv-qli)
588  wthv_pr = w_zn*thv_pr
589 
590  thv_pr_jip1 = th_pr_jip1+tmp_th_jip1*(0.61*qv_pr_jip1-qli_pr_jip1)+th_pr_jip1*(0.61*qv_jip1-qli_jip1)
591  thv_pr_jim1 = th_pr_jim1+tmp_th_jim1*(0.61*qv_pr_jim1-qli_pr_jim1)+th_pr_jim1*(0.61*qv_jim1-qli_jim1)
592  thv_pr_jp1i = th_pr_jp1i+tmp_th_jp1i*(0.61*qv_pr_jp1i-qli_pr_jp1i)+th_pr_jp1i*(0.61*qv_jp1i-qli_jp1i)
593  thv_pr_jm1i = th_pr_jm1i+tmp_th_jm1i*(0.61*qv_pr_jm1i-qli_pr_jm1i)+th_pr_jm1i*(0.61*qv_jm1i-qli_jm1i)
594  end if
595 
596  if (.not. thv_from_th_with_liqice) then
597  thv = tmp_th*(1.0+0.61*qv)
598  ! the 0.61*qv*th_pr is very small and can be neglected
599  thv_pr = th_pr + 0.61 * tmp_th * qv_pr + 0.61 * qv * th_pr
600  wthv_pr = w_zn * thv_pr
601  thv_pr_jip1 = th_pr_jip1+0.61*tmp_th_jip1*qv_pr_jip1+0.61*qv_jip1*th_pr_jip1
602  thv_pr_jim1 = th_pr_jim1+0.61*tmp_th_jim1*qv_pr_jim1+0.61*qv_jim1*th_pr_jim1
603  thv_pr_jp1i = th_pr_jp1i+0.61*tmp_th_jp1i*qv_pr_jp1i+0.61*qv_jp1i*th_pr_jp1i
604  thv_pr_jm1i = th_pr_jm1i+0.61*tmp_th_jm1i*qv_pr_jm1i+0.61*qv_jm1i*th_pr_jm1i
605  end if
606 
607  end if ! l_do_near
609 
610 
612  local_diag(1) = 1.0_default_precision
613  local_diag(2) = w_zn
614  local_diag(3) = w_zn2
615  local_diag(4) = tmp_th
616  local_diag(5) = wth
617  local_diag(6) = th_pr
618  local_diag(7) = wth_pr
619  local_diag(8) = thv_pr
620  local_diag(9) = wthv_pr
621  local_diag(10) = th_pr2
622  local_diag(11) = wthsg
623  local_diag(12) = w_zn3
624  local_diag(13) = relhum
625  local_diag(14) = tmp_u
626  local_diag(15) = tmp_v
627  local_diag(16) = wu
628  local_diag(17) = wv
629  local_diag(18) = wusg
630  local_diag(19) = wvsg
631  local_diag(20) = tdegk
632  local_diag(21) = th_h
633  local_diag(22) = th_h_pr1
634  local_diag(23) = th_h_pr2
635  local_diag(24) = qvli
636  local_diag(25) = qvli_pr
637  local_diag(26) = qvli_pr2
638  local_diag(27) = qppt
639  local_diag(28) = qppt_pr
640  local_diag(29) = qppt_pr2
641  local_diag(30) = wqvli_pr
642  local_diag(31) = wqppt_pr
643 
645  do inc = 1, ncond
646 
647  l_cond = .false. ! initialization
648 
649  !Condition processing on requested condition
650  SELECT CASE ( trim(cond_request(inc)) )
651 
652  CASE DEFAULT
653  call log_master_log(log_error, &
654  "Condition '"//trim(cond_request(inc))//"' has not been set up.")
655 
656  CASE ("ALL")
657  l_cond = .true.
658 
659  CASE ("BYu")
660  l_cond = ( (w_zn .gt. wupcrit) .and. (thv_pr .gt. thvprcrit) )
661 
662  CASE ("BCu")
663  l_cond = ( (w_zn .gt. wupcrit .and. thv_pr .gt. thvprcrit ) .and. &
664  ( ql .gt. qlcrit .or. qi .gt. qicrit ) )
665 
666  CASE ("NrBCu")
668  if ( .not. ( (w_zn .gt. wupcrit .and. thv_pr .gt. thvprcrit ) .and. &
669  ( ql .gt. qlcrit .or. qi .gt. qicrit ) ) ) then
670 
671  l_cond = ( ( (w_zn_jip1 .gt. wupcrit .and. thv_pr_jip1 .gt. thvprcrit) .and. &
672  (ql_jip1 .gt. qlcrit .or. qi_jip1 .gt. qicrit) ) &
673  .or. &
674  ( (w_zn_jim1 .gt. wupcrit .and. thv_pr_jim1 .gt. thvprcrit) .and. &
675  (ql_jim1 .gt. qlcrit .or. qi_jim1 .gt. qicrit) ) &
676  .or. &
677  ( (w_zn_jp1i .gt. wupcrit .and. thv_pr_jp1i .gt. thvprcrit) .and. &
678  (ql_jp1i .gt. qlcrit .or. qi_jp1i .gt. qicrit) ) &
679  .or. &
680  ( (w_zn_jm1i .gt. wupcrit .and. thv_pr_jm1i .gt. thvprcrit) .and. &
681  (ql_jm1i .gt. qlcrit .or. qi_jm1i .gt. qicrit) ) )
682 
683  end if
684 
685  CASE ("AC")
686  l_cond = ( (ql .gt. qlcrit) .or. (qi .gt. qicrit) )
687 
688  CASE ("ACu")
689  l_cond = ( (w_zn .gt. wupcrit) .and. &
690  ( (ql .gt. qlcrit) .or. (qi .gt. qicrit) ) )
691 
692  CASE ("ACd")
693  l_cond = ( (w_zn .lt. wdwncrit) .and. &
694  ( (ql .gt. qlcrit) .or. (qi .gt. qicrit) ) )
695 
696  CASE ("WG1")
697  l_cond = ( w_zn .gt. wsupcrit )
698 
699  CASE ("WL1")
700  l_cond = ( w_zn .lt. wsdwncrit )
701 
702  CASE ("ALu")
703  l_cond = ( w_zn .gt. wupcrit )
704 
705  CASE ("ALd")
706  l_cond = ( w_zn .lt. wdwncrit )
707 
708  CASE ("CLu")
709  l_cond = ( (w_zn .gt. wupcrit) .and. (ql .gt. qlcrit) )
710 
711  CASE ("CLd")
712  l_cond = ( (w_zn .lt. wdwncrit) .and. (ql .gt. qlcrit) )
713 
714  CASE ("AH")
715  l_cond = ( (ql .gt. qlcrit) .or. (qi .gt. qicrit) .or. (qppt .gt. qpptcrit) )
716 
717  CASE ("AL")
718  l_cond = ( ql .gt. qlcrit )
719 
720  CASE ("AI")
721  l_cond = ( qi .gt. qicrit )
722 
723  CASE ("PPd")
724  l_cond = ( (w_zn .lt. wdwncrit) .and. (qppt .gt. qpptcrit) )
725 
726  CASE ("VPd")
727  l_cond = ( (w_zn .lt. wdwncrit) .and. (qppt .gt. vpptcrit) )
728 
729  CASE ("PVd")
730  l_cond = ( (w_zn .lt. wsdwncrit) .and. (qppt .gt. qpptcrit) )
731 
732  CASE ("MO")
733  l_cond = ( (w_zn .gt. wupcrit) .and. (qvl .gt. qsat) )
734 
735  CASE ("BM")
736  l_cond = ( (w_zn .gt. wupcrit) .and. (thv_pr .gt. thvprcrit) .and. (qvl .gt. qsat) )
737 
738  CASE ("AA")
739  l_cond = ( (ql .lt. qlcrit) .and. (qi .lt. qicrit) .and. &
740  (w_zn .lt. wdwncrit) .and. (qppt .gt. qpptcrit) )
741 
742  CASE ("AV")
743  l_cond = ( (ql .lt. qlcrit) .and. (qi .lt. qicrit) .and. &
744  (w_zn .lt. wdwncrit) .and. (qppt .gt. vpptcrit) )
745 
746  END SELECT
747 
748 
750  if (l_cond) then
751  conddiags_tot(k,inc,:) = conddiags_tot(k,inc,:) + local_diag(diag_locations)
752  else
753  conddiags_tot(k,ncond+inc,:) = conddiags_tot(k,ncond+inc,:) + local_diag(diag_locations)
754  end if
755 
756  end do ! inc loop over requested conditions
757 
758  enddo ! k loop over vertical levels
759 
Here is the caller graph for this function:

Variable Documentation

◆ cond_locations

integer, dimension(:), allocatable conditional_diagnostics_column_mod::cond_locations
private

Definition at line 55 of file conditional_diagnostics_column.F90.

◆ cond_long

character(len=string_length), dimension(:), allocatable, public conditional_diagnostics_column_mod::cond_long

Definition at line 53 of file conditional_diagnostics_column.F90.

◆ cond_request

character(len=string_length), dimension(:), allocatable, public conditional_diagnostics_column_mod::cond_request

Definition at line 53 of file conditional_diagnostics_column.F90.

53  character(len=STRING_LENGTH), dimension(:), allocatable :: cond_request, cond_long

◆ conddiags_tot

real(kind=default_precision), dimension(:,:,:), allocatable, public conditional_diagnostics_column_mod::conddiags_tot

Definition at line 27 of file conditional_diagnostics_column.F90.

27  real(kind=default_precision), dimension(:,:,:), allocatable :: conddiags_tot

◆ diag_locations

integer, dimension(:), allocatable conditional_diagnostics_column_mod::diag_locations
private

Definition at line 55 of file conditional_diagnostics_column.F90.

55  integer, dimension(:), allocatable :: diag_locations, cond_locations

◆ diag_long

character(len=string_length), dimension(:), allocatable, public conditional_diagnostics_column_mod::diag_long

Definition at line 54 of file conditional_diagnostics_column.F90.

◆ diag_request

character(len=string_length), dimension(:), allocatable, public conditional_diagnostics_column_mod::diag_request

Definition at line 54 of file conditional_diagnostics_column.F90.

54  character(len=STRING_LENGTH), dimension(:), allocatable :: diag_request, diag_long

◆ diagnostic_generation_frequency

integer conditional_diagnostics_column_mod::diagnostic_generation_frequency
private

Definition at line 58 of file conditional_diagnostics_column.F90.

58  integer :: diagnostic_generation_frequency

◆ gpts_total

real(kind=default_precision), public conditional_diagnostics_column_mod::gpts_total

Definition at line 44 of file conditional_diagnostics_column.F90.

44  real(kind=default_precision) :: gpts_total

◆ iqg

integer conditional_diagnostics_column_mod::iqg
private

Definition at line 25 of file conditional_diagnostics_column.F90.

◆ iqi

integer conditional_diagnostics_column_mod::iqi
private

Definition at line 25 of file conditional_diagnostics_column.F90.

◆ iql

integer conditional_diagnostics_column_mod::iql
private

Definition at line 25 of file conditional_diagnostics_column.F90.

◆ iqr

integer conditional_diagnostics_column_mod::iqr
private

Definition at line 25 of file conditional_diagnostics_column.F90.

◆ iqs

integer conditional_diagnostics_column_mod::iqs
private

Definition at line 25 of file conditional_diagnostics_column.F90.

◆ iqv

integer conditional_diagnostics_column_mod::iqv
private

Definition at line 25 of file conditional_diagnostics_column.F90.

25  integer :: iqv, iql, iqr, iqi, iqs, iqg

◆ l_do_near

logical conditional_diagnostics_column_mod::l_do_near
private

Definition at line 48 of file conditional_diagnostics_column.F90.

48  logical :: l_do_near

◆ l_qi_qr_qs_qg

logical conditional_diagnostics_column_mod::l_qi_qr_qs_qg
private

Definition at line 47 of file conditional_diagnostics_column.F90.

47  logical :: l_qi_qr_qs_qg

◆ master_conditions_list

character(len=string_length), dimension(23) conditional_diagnostics_column_mod::master_conditions_list
private

Definition at line 50 of file conditional_diagnostics_column.F90.

50  character(len=STRING_LENGTH), dimension(23) :: master_conditions_list

◆ master_diagnostics_list

character(len=string_length), dimension(31) conditional_diagnostics_column_mod::master_diagnostics_list
private

Definition at line 51 of file conditional_diagnostics_column.F90.

51  character(len=STRING_LENGTH), dimension(31) :: master_diagnostics_list

◆ ncond

integer, public conditional_diagnostics_column_mod::ncond

Definition at line 39 of file conditional_diagnostics_column.F90.

39  integer :: ncond

◆ ndiag

integer, public conditional_diagnostics_column_mod::ndiag

Definition at line 40 of file conditional_diagnostics_column.F90.

40  integer :: ndiag

◆ qicrit

real(kind=default_precision) conditional_diagnostics_column_mod::qicrit
private

Definition at line 30 of file conditional_diagnostics_column.F90.

30  real(kind=default_precision) :: qicrit

◆ qlcrit

real(kind=default_precision) conditional_diagnostics_column_mod::qlcrit
private

Definition at line 29 of file conditional_diagnostics_column.F90.

29  real(kind=default_precision) :: qlcrit

◆ qpptcrit

real(kind=default_precision) conditional_diagnostics_column_mod::qpptcrit
private

Definition at line 31 of file conditional_diagnostics_column.F90.

31  real(kind=default_precision) :: qpptcrit

◆ requested_area

integer, public conditional_diagnostics_column_mod::requested_area

Definition at line 56 of file conditional_diagnostics_column.F90.

56  integer :: requested_area

◆ thv_from_th_with_liqice

logical conditional_diagnostics_column_mod::thv_from_th_with_liqice
private

Definition at line 46 of file conditional_diagnostics_column.F90.

46  logical :: thv_from_th_with_liqice

◆ thvprcrit

real(kind=default_precision) conditional_diagnostics_column_mod::thvprcrit
private

Definition at line 33 of file conditional_diagnostics_column.F90.

33  real(kind=default_precision) :: thvprcrit

◆ vpptcrit

real(kind=default_precision) conditional_diagnostics_column_mod::vpptcrit
private

Definition at line 32 of file conditional_diagnostics_column.F90.

32  real(kind=default_precision) :: vpptcrit

◆ wdwncrit

real(kind=default_precision) conditional_diagnostics_column_mod::wdwncrit
private

Definition at line 37 of file conditional_diagnostics_column.F90.

37  real(kind=default_precision) :: wdwncrit

◆ wsdwncrit

real(kind=default_precision) conditional_diagnostics_column_mod::wsdwncrit
private

Definition at line 34 of file conditional_diagnostics_column.F90.

34  real(kind=default_precision) :: wsdwncrit

◆ wsupcrit

real(kind=default_precision) conditional_diagnostics_column_mod::wsupcrit
private

Definition at line 35 of file conditional_diagnostics_column.F90.

35  real(kind=default_precision) :: wsupcrit

◆ wupcrit

real(kind=default_precision) conditional_diagnostics_column_mod::wupcrit
private

Definition at line 36 of file conditional_diagnostics_column.F90.

36  real(kind=default_precision) :: wupcrit

◆ x_size

integer conditional_diagnostics_column_mod::x_size
private

Definition at line 41 of file conditional_diagnostics_column.F90.

41  integer :: x_size

◆ y_size

integer conditional_diagnostics_column_mod::y_size
private

Definition at line 42 of file conditional_diagnostics_column.F90.

42  integer :: y_size
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
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