MONC
Functions/Subroutines | Variables
pdf_analysis_mod Module Reference

Calculates fields related to distributions of data on full-domain horizontal 2d slices. More...

Functions/Subroutines

type(component_descriptor_type) function, public pdf_analysis_get_descriptor ()
 Returns the component descriptor of pdf analysis module. More...
 
subroutine init_callback (current_state)
 Called on MONC initialisation, will allocate appropriate data structures. More...
 
subroutine timestep_callback (current_state)
 Will sort the values across the whole domain and calculate the value corresponding to the percentile threshold. More...
 
subroutine finalisation_callback (current_state)
 Frees up the temporary data for the tm_allp. More...
 
subroutine calculate_w_percentiles (current_state)
 Calculates the w percentiles over the whole domain and stores these in the w up/dwn percentile arrays of current_state. More...
 
subroutine mergesortmerge (A, NA, B, NB, C, NC)
 Combines with MergeSort sorting algorithm taken from: More...
 
recursive subroutine mergesort (A, N, T)
 Combines with MergeSortMerge sorting algorithm taken from: More...
 
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 start_x
 
integer end_x
 
integer start_y
 
integer end_y
 
integer xsize
 
integer ysize
 
real(kind=default_precision) uppercrit
 
real(kind=default_precision) dwnpercrit
 
real(kind=default_precision), dimension(:), allocatable tmp_all
 
integer, dimension(:), allocatable gpts_on_proc
 
integer, dimension(:), allocatable displacements
 
integer tpts
 
integer lpts
 
logical show_critical_w
 
integer diagnostic_generation_frequency
 

Detailed Description

Calculates fields related to distributions of data on full-domain horizontal 2d slices.

Function/Subroutine Documentation

◆ calculate_w_percentiles()

subroutine pdf_analysis_mod::calculate_w_percentiles ( type(model_state_type), intent(inout), target  current_state)
private

Calculates the w percentiles over the whole domain and stores these in the w up/dwn percentile arrays of current_state.

Parameters
current_stateThe current model state

initialize diagnostic thresholds

reset thresholds

Loop over levels

specify local data area

Reshape to 1-D array

Perform sort of data on local process

Gather 2d field to single process

Perform global operations

Sort the global data on single process

Determine threshold updraft and downdraft values

do some stdout diagnostic work

Inform all processes of calculated thresholds

Display some diagnostics, if requested

Definition at line 136 of file pdf_analysis.F90.

137  type(model_state_type), target, intent(inout) :: current_state
138  real(kind=default_precision), dimension(lpts) :: tmp_var
139 
140  integer :: i, j, k, num_neg, num_pos, dd_thresh_pos, ud_thresh_pos
141  integer :: max_up_k, min_dwn_k
142  real(kind=default_precision), dimension((lpts+1)/2) :: t
143  real(kind=default_precision), dimension((tpts+1)/2) :: tall
144  real(kind=default_precision) :: max_up, min_dwn, &
145  max_up_th, min_dwn_th
146  integer :: ierr
147  real(kind=default_precision), dimension(ysize,xsize) :: l2d ! local 2d data
148 
150  max_up_th = 0.0_default_precision
151  min_dwn_th = 0.0_default_precision
152  max_up = max_up_th
153  min_dwn = min_dwn_th
154  max_up_k = 0
155  min_dwn_k = 0
156 
158  current_state%global_grid%configuration%vertical%w_dwn(:) = 0.0_default_precision
159  current_state%global_grid%configuration%vertical%w_up(:) = 0.0_default_precision
160 
162  do k = 2, current_state%local_grid%size(z_index)
163 
165  l2d(:,:)=0.5_default_precision*( current_state%w%data(k , start_y:end_y, start_x:end_x) &
166  + current_state%w%data(k-1, start_y:end_y, start_x:end_x) )
167 
169  tmp_var=pack(l2d,.true.)
170 
172  call mergesort(tmp_var,lpts,t)
173 
175  call mpi_gatherv(tmp_var, lpts, precision_type, tmp_all, gpts_on_proc, displacements, precision_type, &
176  0, current_state%parallel%monc_communicator, ierr )
177 
179  if (current_state%parallel%my_rank == 0) then
180 
182  call mergesort(tmp_all,tpts,tall)
183 
185  num_neg = count(tmp_all < 0.0_default_precision)
186  num_pos = count(tmp_all > 0.0_default_precision)
187 
188  dd_thresh_pos = int(num_neg * dwnpercrit)
189  ud_thresh_pos = tpts - int(num_pos * uppercrit) + 1
190 
191  if ( dd_thresh_pos == 0 ) dd_thresh_pos = 1
192  if ( ud_thresh_pos == 0 .or. num_pos == 0 ) ud_thresh_pos = tpts
193 
194  current_state%global_grid%configuration%vertical%w_dwn(k) = tmp_all(dd_thresh_pos)
195  current_state%global_grid%configuration%vertical%w_up(k) = tmp_all(ud_thresh_pos)
196 
198  if (show_critical_w) then
199  if ( tmp_all(dd_thresh_pos) < min_dwn_th ) then
200  min_dwn_th = tmp_all(dd_thresh_pos)
201  min_dwn = tmp_all(1) ! sorted array
202  min_dwn_k = k
203  end if
204  if ( tmp_all(ud_thresh_pos) > max_up_th ) then
205  max_up_th = tmp_all(ud_thresh_pos)
206  max_up = tmp_all(tpts) ! sorted array
207  max_up_k = k
208  end if
209  end if ! show_critical_w
210 
211  end if ! global operations section
212 
213  end do ! loop over k
214 
216  call mpi_bcast(current_state%global_grid%configuration%vertical%w_dwn(:), current_state%local_grid%size(z_index), &
217  precision_type, 0, current_state%parallel%monc_communicator, ierr)
218  call mpi_bcast(current_state%global_grid%configuration%vertical%w_up(:), current_state%local_grid%size(z_index), &
219  precision_type, 0, current_state%parallel%monc_communicator, ierr)
220 
221 
223  if (show_critical_w) then
224  call log_master_log(log_info, 'Time: '//trim(conv_to_string(current_state%time))//' s')
225  call log_master_log(log_info, 'Maximum updraft threshold: '&
226  //trim(conv_to_string(max_up_th))//' found at level '//trim(conv_to_string(max_up_k)) )
227  call log_master_log(log_info, 'Maximum updraft: '&
228  //trim(conv_to_string(max_up))//' at level '//trim(conv_to_string(max_up_k)) )
229  call log_master_log(log_info, 'Minimum downdraft threshold: '&
230  //trim(conv_to_string(min_dwn_th))//' found at level '//trim(conv_to_string(min_dwn_k)) )
231  call log_master_log(log_info, 'Minimum downdraft: '&
232  //trim(conv_to_string(min_dwn))//' at level '//trim(conv_to_string(min_dwn_k)) )
233  end if ! show_critical_w
234 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ field_information_retrieval_callback()

subroutine pdf_analysis_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 314 of file pdf_analysis.F90.

315  type(model_state_type), target, intent(inout) :: current_state
316  character(len=*), intent(in) :: name
317  type(component_field_information_type), intent(out) :: field_information
318 
319  field_information%field_type=component_array_field_type
320  field_information%number_dimensions=1
321  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
322  field_information%data_type=component_double_data_type
323  field_information%enabled=.true.
324 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine pdf_analysis_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 333 of file pdf_analysis.F90.

334  type(model_state_type), target, intent(inout) :: current_state
335  character(len=*), intent(in) :: name
336  type(component_field_value_type), intent(out) :: field_value
337 
338  if (name .eq. "critical_updraft_local") then
339  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)), &
340  source=current_state%global_grid%configuration%vertical%w_up(:))
341  else if (name .eq. "critical_downdraft_local") then
342  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)), &
343  source=current_state%global_grid%configuration%vertical%w_dwn(:))
344  end if
Here is the caller graph for this function:

◆ finalisation_callback()

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

Frees up the temporary data for the tm_allp.

Parameters
current_stateThe current model state

Definition at line 127 of file pdf_analysis.F90.

128  type(model_state_type), target, intent(inout) :: current_state
129 
130  if (allocated(tmp_all)) deallocate(tmp_all)
Here is the caller graph for this function:

◆ init_callback()

subroutine pdf_analysis_mod::init_callback ( type(model_state_type), intent(inout), target  current_state)
private

Called on MONC initialisation, will allocate appropriate data structures.

Parameters
current_stateThe current model state

Allocate space for the global 2d field only on a single process

Allocate and collect horizontal local sizes, send to all proceses

Allocate and initialize displacement values

Allocate critial fields in current_state if a cold start

Definition at line 59 of file pdf_analysis.F90.

60  type(model_state_type), target, intent(inout) :: current_state
61  integer :: ierr, inc
62 
63  tpts = current_state%global_grid%size(x_index)*current_state%global_grid%size(y_index)
64 
65  start_x = current_state%local_grid%local_domain_start_index(x_index)
66  end_x = current_state%local_grid%local_domain_end_index(x_index)
67  start_y = current_state%local_grid%local_domain_start_index(y_index)
68  end_y = current_state%local_grid%local_domain_end_index(y_index)
69 
70  xsize = end_x - start_x + 1
71  ysize = end_y - start_y + 1
72 
73  lpts = xsize*ysize
74 
75  uppercrit = options_get_real(current_state%options_database, "uppercrit")
76  dwnpercrit = options_get_real(current_state%options_database, "dwnpercrit")
77 
78  show_critical_w = options_get_logical(current_state%options_database, "show_critical_w")
79 
81 ! if (current_state%parallel%my_rank == 0)
82  allocate(tmp_all(tpts))
83 ! else
84 ! allocate(tmp_all(1))
85 ! end if
86 
88  allocate(gpts_on_proc(current_state%parallel%processes))
89  call mpi_allgather(lpts, 1, mpi_int, gpts_on_proc, 1, mpi_int, current_state%parallel%monc_communicator, ierr)
90 
92  allocate(displacements(current_state%parallel%processes))
93  displacements(1) = 0
94  do inc = 2, current_state%parallel%processes
95  displacements(inc) = displacements(inc-1) + gpts_on_proc(inc-1)
96  end do ! loop over processes
97 
99  if (.not. current_state%continuation_run) then
100  allocate(current_state%global_grid%configuration%vertical%w_dwn(current_state%local_grid%size(z_index)),&
101  current_state%global_grid%configuration%vertical%w_up(current_state%local_grid%size(z_index)))
102  end if
103 
104  ! Save the sampling_frequency to force diagnostic calculation on select time steps
105  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
106 
Here is the caller graph for this function:

◆ mergesort()

recursive subroutine pdf_analysis_mod::mergesort ( real(kind=default_precision), dimension(n), intent(inout)  A,
integer, intent(in)  N,
real(kind=default_precision), dimension((n+1)/2), intent(out)  T 
)
private

Combines with MergeSortMerge sorting algorithm taken from:

Definition at line 277 of file pdf_analysis.F90.

278 
279  integer, intent(in) :: N
280  real(kind=default_precision), dimension(N), intent(in out) :: a
281  real(kind=default_precision), dimension((N+1)/2), intent (out) :: t
282 
283  integer :: NA,NB
284  real(kind=default_precision) :: v
285 
286  if (n < 2) return
287  if (n == 2) then
288  if (a(1) > a(2)) then
289  v = a(1)
290  a(1) = a(2)
291  a(2) = v
292  endif
293  return
294  endif
295  na=(n+1)/2
296  nb=n-na
297 
298  call mergesort(a,na,t)
299  call mergesort(a(na+1),nb,t)
300 
301  if (a(na) > a(na+1)) then
302  t(1:na)=a(1:na)
303  call mergesortmerge(t,na,a(na+1),nb,a,n)
304  endif
305  return
306 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mergesortmerge()

subroutine pdf_analysis_mod::mergesortmerge ( real(kind=default_precision), dimension(na), intent(inout)  A,
integer, intent(in)  NA,
real(kind=default_precision), dimension(nb), intent(in)  B,
integer, intent(in)  NB,
real(kind=default_precision), dimension(nc), intent(inout)  C,
integer, intent(in)  NC 
)
private

Combines with MergeSort sorting algorithm taken from:

Definition at line 242 of file pdf_analysis.F90.

243 
244  integer, intent(in) :: NA,NB,NC ! Normal usage: NA+NB = NC
245  real(kind=default_precision), intent(in out) :: a(na) ! B overlays C(NA+1:NC)
246  real(kind=default_precision), intent(in) :: b(nb)
247  real(kind=default_precision), intent(in out) :: c(nc)
248 
249  integer :: I,J,K
250 
251  i = 1; j = 1; k = 1;
252  do while(i <= na .and. j <= nb)
253  if (a(i) <= b(j)) then
254  c(k) = a(i)
255  i = i+1
256  else
257  c(k) = b(j)
258  j = j+1
259  endif
260  k = k + 1
261  enddo
262  do while (i <= na)
263  c(k) = a(i)
264  i = i + 1
265  k = k + 1
266  enddo
267  return
268 
Here is the caller graph for this function:

◆ pdf_analysis_get_descriptor()

type(component_descriptor_type) function, public pdf_analysis_mod::pdf_analysis_get_descriptor

Returns the component descriptor of pdf analysis module.

Returns
Component descriptor of pdf analysis

Definition at line 41 of file pdf_analysis.F90.

42  pdf_analysis_get_descriptor%name="pdf_analysis"
43  pdf_analysis_get_descriptor%version=0.1
44  pdf_analysis_get_descriptor%initialisation=>init_callback
45  pdf_analysis_get_descriptor%timestep=>timestep_callback
46  pdf_analysis_get_descriptor%finalisation=>finalisation_callback
47 
48  pdf_analysis_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
49  pdf_analysis_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
50  allocate(pdf_analysis_get_descriptor%published_fields(2))
51  pdf_analysis_get_descriptor%published_fields(1)="critical_updraft_local"
52  pdf_analysis_get_descriptor%published_fields(2)="critical_downdraft_local"
53 
Here is the call graph for this function:

◆ timestep_callback()

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

Will sort the values across the whole domain and calculate the value corresponding to the percentile threshold.

Parameters
current_stateThe current model state

Current forumulation only handles vertical velocity percentiles. Future enhancements may employ this component to perform additional operations that require access to full horizontal fields, such as pdf calculations.

Definition at line 112 of file pdf_analysis.F90.

113  type(model_state_type), target, intent(inout) :: current_state
114 
119 
120  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0) call calculate_w_percentiles(current_state)
121 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ diagnostic_generation_frequency

integer pdf_analysis_mod::diagnostic_generation_frequency
private

Definition at line 32 of file pdf_analysis.F90.

32  integer :: diagnostic_generation_frequency

◆ displacements

integer, dimension(:), allocatable pdf_analysis_mod::displacements
private

Definition at line 23 of file pdf_analysis.F90.

◆ dwnpercrit

real(kind=default_precision) pdf_analysis_mod::dwnpercrit
private

Definition at line 21 of file pdf_analysis.F90.

◆ end_x

integer pdf_analysis_mod::end_x
private

Definition at line 19 of file pdf_analysis.F90.

◆ end_y

integer pdf_analysis_mod::end_y
private

Definition at line 19 of file pdf_analysis.F90.

◆ gpts_on_proc

integer, dimension(:), allocatable pdf_analysis_mod::gpts_on_proc
private

Definition at line 23 of file pdf_analysis.F90.

23  integer, dimension(:), allocatable :: gpts_on_proc, & ! number of horizontal grid points on each process
24  displacements ! displacement for mpi_gatherv

◆ lpts

integer pdf_analysis_mod::lpts
private

Definition at line 28 of file pdf_analysis.F90.

28  integer :: lpts ! local number of horizontal grid points on

◆ show_critical_w

logical pdf_analysis_mod::show_critical_w
private

Definition at line 30 of file pdf_analysis.F90.

30  logical :: show_critical_w ! stdout diagnostic logical

◆ start_x

integer pdf_analysis_mod::start_x
private

Definition at line 19 of file pdf_analysis.F90.

19  integer :: start_x, end_x, start_y, end_y, xsize, ysize

◆ start_y

integer pdf_analysis_mod::start_y
private

Definition at line 19 of file pdf_analysis.F90.

◆ tmp_all

real(kind=default_precision), dimension(:), allocatable pdf_analysis_mod::tmp_all
private

Definition at line 22 of file pdf_analysis.F90.

22  real(kind=default_precision), dimension(:), allocatable :: tmp_all

◆ tpts

integer pdf_analysis_mod::tpts
private

Definition at line 27 of file pdf_analysis.F90.

27  integer :: tpts ! total number of horizontal grid points on full domain

◆ uppercrit

real(kind=default_precision) pdf_analysis_mod::uppercrit
private

Definition at line 21 of file pdf_analysis.F90.

21  real(kind=default_precision) :: uppercrit, dwnpercrit

◆ xsize

integer pdf_analysis_mod::xsize
private

Definition at line 19 of file pdf_analysis.F90.

◆ ysize

integer pdf_analysis_mod::ysize
private

Definition at line 19 of file pdf_analysis.F90.

optionsdatabase_mod::options_get_integer
integer function, public options_get_integer(options_database, key, index)
Retrieves an integer value from the database that matches the provided key.
Definition: optionsdatabase.F90:217
logging_mod::log_info
integer, parameter, public log_info
Log INFO, WARNING and ERROR messages.
Definition: logging.F90:13
optionsdatabase_mod::options_get_logical
logical function, public options_get_logical(options_database, key, index)
Retrieves a logical value from the database that matches the provided key.
Definition: optionsdatabase.F90:154
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
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