MONC
Functions/Subroutines | Variables
meanprofiles_mod Module Reference

Calculates the mean profiles of prognostic variables which are then used in smoothing and other areas. More...

Functions/Subroutines

type(component_descriptor_type) function, public meanprofiles_get_descriptor ()
 Returns the component descriptor of the mean profiles module. More...
 
subroutine init_callback (current_state)
 Called on MONC initialisation, will allocate appropriate data structures. More...
 
subroutine timestep_callback (current_state)
 Will recalculate the mean profiles of each prognostic when called (for the entire local domain) More...
 
subroutine finalisation_callback (current_state)
 Frees up the temporary data for the bars. More...
 
subroutine calculate_mean_profiles (current_state)
 Calculates the global mean profiles and stores these in the ol bar arrays. More...
 
subroutine calculate_sum_profiles (current_state)
 Calculates the sum profiles for the bars for each level globally. More...
 

Variables

integer start_x
 
integer end_x
 
integer start_y
 
integer end_y
 
integer bar_fields
 
real(kind=default_precision) rnhpts
 
real(kind=default_precision), dimension(:,:), allocatable bartmp
 

Detailed Description

Calculates the mean profiles of prognostic variables which are then used in smoothing and other areas.

Function/Subroutine Documentation

◆ calculate_mean_profiles()

subroutine meanprofiles_mod::calculate_mean_profiles ( type(model_state_type), intent(inout), target  current_state)
private

Calculates the global mean profiles and stores these in the ol bar arrays.

Parameters
current_stateThe current model state

Definition at line 105 of file meanprofiles.F90.

105  type(model_state_type), target, intent(inout) :: current_state
106 
107  integer :: bar_index, i
108 
109  call calculate_sum_profiles(current_state)
110 
111  bar_index=1
112 #ifdef U_ACTIVE
113  current_state%global_grid%configuration%vertical%olubar(:)=bartmp(:, bar_index)*rnhpts
114  current_state%global_grid%configuration%vertical%olzubar(:)=bartmp(:, bar_index+1)*rnhpts
115  bar_index=bar_index+2
116 #endif
117 #ifdef V_ACTIVE
118  current_state%global_grid%configuration%vertical%olvbar(:)=bartmp(:, bar_index)*rnhpts
119  current_state%global_grid%configuration%vertical%olzvbar(:)=bartmp(:, bar_index+1)*rnhpts
120  bar_index=bar_index+2
121 #endif
122  if (current_state%th%active) then
123  current_state%global_grid%configuration%vertical%olthbar(:)=bartmp(:, bar_index)*rnhpts
124  current_state%global_grid%configuration%vertical%olzthbar(:)=bartmp(:, bar_index+1)*rnhpts
125  bar_index=bar_index+2
126  end if
127  do i=1,current_state%number_q_fields
128  if (current_state%q(i)%active) then
129  current_state%global_grid%configuration%vertical%olqbar(:, i)=bartmp(:, bar_index)*rnhpts
130  current_state%global_grid%configuration%vertical%olzqbar(:, i)=bartmp(:, bar_index+1)*rnhpts
131  bar_index=bar_index+2
132  end if
133  end do
134 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_sum_profiles()

subroutine meanprofiles_mod::calculate_sum_profiles ( type(model_state_type), intent(inout)  current_state)
private

Calculates the sum profiles for the bars for each level globally.

Parameters
current_stateThe current model state_mod

Definition at line 140 of file meanprofiles.F90.

140  type(model_state_type), intent(inout) :: current_state
141 
142  integer :: k, n, bar_index, ierr
143 
144  do k=current_state%local_grid%local_domain_start_index(z_index), current_state%local_grid%local_domain_end_index(z_index)
145  bar_index=1
146 #ifdef U_ACTIVE
147  bartmp(k, bar_index)=sum(current_state%u%data(k, start_y:end_y, start_x:end_x))
148  bartmp(k, bar_index+1)=sum(current_state%zu%data(k, start_y:end_y, start_x:end_x))
149  bar_index=bar_index+2
150 #endif
151 #ifdef V_ACTIVE
152  bartmp(k, bar_index)=sum(current_state%v%data(k, start_y:end_y, start_x:end_x))
153  bartmp(k, bar_index+1)=sum(current_state%zv%data(k, start_y:end_y, start_x:end_x))
154  bar_index=bar_index+2
155 #endif
156  if (current_state%th%active) then
157  bartmp(k, bar_index)=sum(current_state%th%data(k, start_y:end_y, start_x:end_x))
158  bartmp(k, bar_index+1)=sum(current_state%zth%data(k, start_y:end_y, start_x:end_x))
159  bar_index=bar_index+2
160  end if
161  do n=1,current_state%number_q_fields
162  if (current_state%q(n)%active) then
163  bartmp(k, bar_index)=sum(current_state%q(n)%data(k, start_y:end_y, start_x:end_x))
164  bartmp(k, bar_index+1)=sum(current_state%zq(n)%data(k, start_y:end_y, start_x:end_x))
165  bar_index=bar_index+2
166  end if
167  end do
168  end do
169 
170  call mpi_allreduce(mpi_in_place, bartmp, bar_fields*current_state%local_grid%size(z_index), precision_type, mpi_sum, &
171  current_state%parallel%monc_communicator, ierr)
Here is the caller graph for this function:

◆ finalisation_callback()

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

Frees up the temporary data for the bars.

Parameters
current_stateThe current model state

Definition at line 97 of file meanprofiles.F90.

97  type(model_state_type), target, intent(inout) :: current_state
98 
99  if (allocated(bartmp)) deallocate(bartmp)
Here is the caller graph for this function:

◆ init_callback()

subroutine meanprofiles_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

Definition at line 35 of file meanprofiles.F90.

35  type(model_state_type), target, intent(inout) :: current_state
36 
37  bar_fields=0
38 
39  rnhpts=1.0_default_precision/real(current_state%global_grid%size(x_index)*current_state%global_grid%size(y_index))
40 
41  start_x=current_state%local_grid%local_domain_start_index(x_index)
42  end_x=current_state%local_grid%local_domain_end_index(X_INDEX)
43  start_y=current_state%local_grid%local_domain_start_index(y_index)
44  end_y=current_state%local_grid%local_domain_end_index(Y_INDEX)
45 
46 #ifdef U_ACTIVE
47  if (.not. current_state%continuation_run) then
48  allocate(current_state%global_grid%configuration%vertical%olubar(current_state%local_grid%size(z_index)),&
49  current_state%global_grid%configuration%vertical%olzubar(current_state%local_grid%size(z_index)))
50  end if
51  allocate(current_state%global_grid%configuration%vertical%savolubar(current_state%local_grid%size(z_index)))
52  bar_fields=bar_fields+2
53 #endif
54 #ifdef V_ACTIVE
55  if (.not. current_state%continuation_run) then
56  allocate(current_state%global_grid%configuration%vertical%olvbar(current_state%local_grid%size(z_index)),&
57  current_state%global_grid%configuration%vertical%olzvbar(current_state%local_grid%size(z_index)))
58  end if
59  allocate(current_state%global_grid%configuration%vertical%savolvbar(current_state%local_grid%size(z_index)))
60  bar_fields=bar_fields+2
61 #endif
62  if (current_state%th%active) then
63  if (.not. current_state%continuation_run) then
64  allocate(current_state%global_grid%configuration%vertical%olthbar(current_state%local_grid%size(z_index)),&
65  current_state%global_grid%configuration%vertical%olzthbar(current_state%local_grid%size(z_index)))
66  end if
67  bar_fields=bar_fields+2
68  end if
69  if (current_state%number_q_fields .gt. 0) then
70  bar_fields=bar_fields+(current_state%number_q_fields*2)
71  if (.not. current_state%continuation_run) then
72  allocate(current_state%global_grid%configuration%vertical%olqbar(current_state%local_grid%size(z_index), &
73  current_state%number_q_fields), current_state%global_grid%configuration%vertical%olzqbar(&
74  current_state%local_grid%size(z_index), current_state%number_q_fields))
75  end if
76  end if
77  allocate(bartmp(current_state%local_grid%size(z_index), bar_fields))
78 
79  ! Do the initial calculation for the first timestep
80  if (.not. current_state%continuation_run) call calculate_mean_profiles(current_state)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ meanprofiles_get_descriptor()

type(component_descriptor_type) function, public meanprofiles_mod::meanprofiles_get_descriptor ( )

Returns the component descriptor of the mean profiles module.

Returns
Component descriptor of mean profiles

Definition at line 25 of file meanprofiles.F90.

25  meanprofiles_get_descriptor%name="mean_profiles"
26  meanprofiles_get_descriptor%version=0.1
27  meanprofiles_get_descriptor%initialisation=>init_callback
28  meanprofiles_get_descriptor%timestep=>timestep_callback
29  meanprofiles_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ timestep_callback()

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

Will recalculate the mean profiles of each prognostic when called (for the entire local domain)

Parameters
current_stateThe current model state

Definition at line 86 of file meanprofiles.F90.

86  type(model_state_type), target, intent(inout) :: current_state
87 
88  current_state%global_grid%configuration%vertical%savolubar=current_state%global_grid%configuration%vertical%olubar
89  current_state%global_grid%configuration%vertical%savolvbar=current_state%global_grid%configuration%vertical%olvbar
90 
91  call calculate_mean_profiles(current_state)
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ bar_fields

integer meanprofiles_mod::bar_fields
private

Definition at line 14 of file meanprofiles.F90.

◆ bartmp

real(kind=default_precision), dimension(:,:), allocatable meanprofiles_mod::bartmp
private

Definition at line 17 of file meanprofiles.F90.

17  real(kind=DEFAULT_PRECISION), dimension(:,:), allocatable :: bartmp

◆ end_x

integer meanprofiles_mod::end_x
private

Definition at line 14 of file meanprofiles.F90.

◆ end_y

integer meanprofiles_mod::end_y
private

Definition at line 14 of file meanprofiles.F90.

◆ rnhpts

real(kind=default_precision) meanprofiles_mod::rnhpts
private

Definition at line 16 of file meanprofiles.F90.

16  real(kind=DEFAULT_PRECISION) :: rnhpts

◆ start_x

integer meanprofiles_mod::start_x
private

Definition at line 14 of file meanprofiles.F90.

14  integer :: start_x, end_x, start_y, end_y, bar_fields

◆ start_y

integer meanprofiles_mod::start_y
private

Definition at line 14 of file meanprofiles.F90.