MONC
Functions/Subroutines | Variables
conditional_diagnostics_whole_mod Module Reference

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

Functions/Subroutines

type(component_descriptor_type) function, public conditional_diagnostics_whole_get_descriptor ()
 Provides registry information for the component. More...
 
subroutine initialisation_callback (current_state)
 Initialisation hook: currently doesn't need to do anything. More...
 
subroutine timestep_callback (current_state)
 The timestep hook will perform averaging of the conditional diagnostics. More...
 
subroutine finalisation_callback (current_state)
 Called on termination: currently doesn't need to do anything. More...
 

Variables

integer diagnostic_generation_frequency
 

Detailed Description

Conditionally averaged diagnostics, Part 2 of 2.

Function/Subroutine Documentation

◆ conditional_diagnostics_whole_get_descriptor()

type(component_descriptor_type) function, public conditional_diagnostics_whole_mod::conditional_diagnostics_whole_get_descriptor

Provides registry information for the component.

Returns
The component descriptor that describes this component

Definition at line 33 of file conditional_diagnostics_whole.F90.

34  conditional_diagnostics_whole_get_descriptor%name="conditional_diagnostics_whole"
35  conditional_diagnostics_whole_get_descriptor%version=0.1
36  conditional_diagnostics_whole_get_descriptor%initialisation=>initialisation_callback
37  conditional_diagnostics_whole_get_descriptor%timestep=>timestep_callback
38  conditional_diagnostics_whole_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ finalisation_callback()

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

Called on termination: currently doesn't need to do anything.

Parameters
current_stateThe current model state

Definition at line 112 of file conditional_diagnostics_whole.F90.

113  type(model_state_type), target, intent(inout) :: current_state
114 
Here is the caller graph for this function:

◆ initialisation_callback()

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

Initialisation hook: currently doesn't need to do anything.

Parameters
current_stateThe current model state

Definition at line 44 of file conditional_diagnostics_whole.F90.

45  type(model_state_type), target, intent(inout) :: current_state
46 
47  ! Save the sampling_frequency to force diagnostic calculation on select time steps
48  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
49 
Here is the caller graph for this function:

◆ timestep_callback()

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

The timestep hook will perform averaging of the conditional diagnostics.

Parameters
current_stateThe current model state

Decide if conditions are appropriate to proceed with calculations

Sum conditional diagnostics total array (horizontally), placing the result on process 0 Reduction call on process 0 requires special MPI_IN_PLACE handling

Average summed diagnostics over the domain by dividing the total diagnostic for each condition by the total number of points for the associated conditions. This is NOT done for the area diagnostic, identified by the requested_area position in the array. Note: CondDiags_tot(k, ncond, ndiag)

Convert the area total number of points for each condition to fraction of the horizontal domain.

Apply missing data mask to top/bottom

Since the xml handling of CondDiags_tot will perform a sum over processes, divide by the number of proceses.

Broadcast the process-fractional solution to all processes.

Definition at line 55 of file conditional_diagnostics_whole.F90.

56  type(model_state_type), target, intent(inout) :: current_state
57  integer :: k,cnc,dnc,ierr
58  real(kind=default_precision) :: temp
59 
61  if (.not. mod(current_state%timestep, diagnostic_generation_frequency) == 0) return
62 
65  if (current_state%parallel%my_rank == 0) then
66  call mpi_reduce(mpi_in_place , conddiags_tot, ncond*2*ndiag*current_state%local_grid%size(z_index), &
67  precision_type, mpi_sum, 0, current_state%parallel%monc_communicator, ierr)
68  else
69  call mpi_reduce(conddiags_tot, conddiags_tot, ncond*2*ndiag*current_state%local_grid%size(z_index), &
70  precision_type, mpi_sum, 0, current_state%parallel%monc_communicator, ierr)
71  end if
72 
77  do dnc = 1, ndiag
78  if (dnc /= requested_area) then
79  do cnc = 1, ncond
80  do k = 2, current_state%local_grid%size(z_index) - 1
81  temp = conddiags_tot(k,cnc,requested_area)
82  if (temp .gt. 0) then
83  conddiags_tot(k,cnc,dnc) = conddiags_tot(k,cnc,dnc) / temp
84  else
85  conddiags_tot(k,cnc,dnc) = rmdi
86  end if
87  end do ! k
88  end do ! cnc over ncond*2
89  end if ! check requested_area
90  end do ! dnc over ndiag
91 
93  conddiags_tot(:,:,requested_area) = conddiags_tot(:,:,requested_area) / gpts_total
94 
96  conddiags_tot(1,:,:) = rmdi
97  conddiags_tot(current_state%local_grid%size(z_index),:,:) = rmdi
98 
101  conddiags_tot = conddiags_tot / current_state%parallel%processes
102 
104  call mpi_bcast(conddiags_tot, ncond*2*ndiag*current_state%local_grid%size(z_index), &
105  precision_type, 0, current_state%parallel%monc_communicator, ierr)
106 
Here is the caller graph for this function:

Variable Documentation

◆ diagnostic_generation_frequency

integer conditional_diagnostics_whole_mod::diagnostic_generation_frequency
private

Definition at line 24 of file conditional_diagnostics_whole.F90.

24  integer :: diagnostic_generation_frequency
missing_data_mod::rmdi
real, parameter rmdi
Definition: missing_data_mod.F90:29
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