MONC
Functions/Subroutines | Variables
swapsmooth_mod Module Reference

Does the swapping and smoothing which is called for each column as part of the pressure-terms group of components Note that this does not currently implement smoothing for mean profiles. More...

Functions/Subroutines

type(component_descriptor_type) function, public swapsmooth_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine initialisation_callback (current_state)
 Initialises the swap and smooth component. More...
 
subroutine timestep_callback (current_state)
 Called for each non halo timestep column and will perform swapping and smoothing as required on that column. More...
 
subroutine swap_and_smooth_classic (current_state, old_smoother)
 Classic swap and smooth based upon the old or no smoothing. More...
 
subroutine swap_and_smooth_robert_filter (current_state)
 Swap and smooth with a Robert filter. More...
 
subroutine classic_for_average_profiles (current_state, old_smoother)
 Does swapping and smoothing (using classic algorithm) for the average profiles (the bars) More...
 
subroutine robert_filter_for_average_profiles (current_state)
 Does swapping and smoothing (using robert filter algorithm) for the average profiles (the bars) More...
 

Variables

logical mean_profiles_active =.false.
 Whether or not mean profiles need smoothing. More...
 

Detailed Description

Does the swapping and smoothing which is called for each column as part of the pressure-terms group of components Note that this does not currently implement smoothing for mean profiles.

Function/Subroutine Documentation

◆ classic_for_average_profiles()

subroutine swapsmooth_mod::classic_for_average_profiles ( type(model_state_type), intent(inout)  current_state,
logical, intent(in)  old_smoother 
)
private

Does swapping and smoothing (using classic algorithm) for the average profiles (the bars)

Parameters
current_stateThe current model state

Definition at line 180 of file swapsmooth.F90.

181  type(model_state_type), intent(inout) :: current_state
182  logical, intent(in) :: old_smoother
183 
184  integer :: k, n
185  real(kind=default_precision) :: c1, c2
186 
187  if (old_smoother) then
188  c1 = 1.0_default_precision - current_state%tsmth
189  c2 = 2.0_default_precision * current_state%tsmth - 1.0_default_precision
190  else
191  c1 = 1.0_default_precision
192  c2 = -1.0_default_precision
193  end if
194 
195  do k=1,current_state%global_grid%size(z_index)
196 #ifdef U_ACTIVE
197  current_state%global_grid%configuration%vertical%olzubar(k)=current_state%global_grid%configuration%vertical%olubar(k) +&
198  current_state%global_grid%configuration%vertical%olzubar(k)
199  current_state%global_grid%configuration%vertical%olubar(k)=current_state%global_grid%configuration%vertical%olzubar(k) *&
200  c1 + current_state%global_grid%configuration%vertical%olubar(k) * c2
201  current_state%global_grid%configuration%vertical%olzubar(k)=current_state%global_grid%configuration%vertical%olzubar(k) -&
202  current_state%global_grid%configuration%vertical%olubar(k)
203 #endif
204 #ifdef V_ACTIVE
205  current_state%global_grid%configuration%vertical%olzvbar(k)=current_state%global_grid%configuration%vertical%olvbar(k) +&
206  current_state%global_grid%configuration%vertical%olzvbar(k)
207  current_state%global_grid%configuration%vertical%olvbar(k)=current_state%global_grid%configuration%vertical%olzvbar(k) *&
208  c1 + current_state%global_grid%configuration%vertical%olvbar(k) * c2
209  current_state%global_grid%configuration%vertical%olzvbar(k)=current_state%global_grid%configuration%vertical%olzvbar(k) -&
210  current_state%global_grid%configuration%vertical%olvbar(k)
211 #endif
212  if (current_state%th%active) then
213  current_state%global_grid%configuration%vertical%olzthbar(k)=current_state%global_grid%configuration%vertical%olthbar(k)+&
214  current_state%global_grid%configuration%vertical%olzthbar(k)
215  current_state%global_grid%configuration%vertical%olthbar(k)=current_state%global_grid%configuration%vertical%olzthbar(k)*&
216  c1 + current_state%global_grid%configuration%vertical%olthbar(k) * c2
217  current_state%global_grid%configuration%vertical%olzthbar(k)=&
218  current_state%global_grid%configuration%vertical%olzthbar(k)-&
219  current_state%global_grid%configuration%vertical%olthbar(k)
220  end if
221  if (current_state%number_q_fields .gt. 0) then
222  do n=1, current_state%number_q_fields
223  current_state%global_grid%configuration%vertical%olzqbar(k,n)=&
224  current_state%global_grid%configuration%vertical%olqbar(k,n)+&
225  current_state%global_grid%configuration%vertical%olzqbar(k,n)
226  current_state%global_grid%configuration%vertical%olqbar(k,n)=&
227  current_state%global_grid%configuration%vertical%olzqbar(k,n)*&
228  c1 + current_state%global_grid%configuration%vertical%olqbar(k,n) * c2
229  current_state%global_grid%configuration%vertical%olzqbar(k,n)=&
230  current_state%global_grid%configuration%vertical%olzqbar(k,n)-&
231  current_state%global_grid%configuration%vertical%olqbar(k,n)
232  end do
233  end if
234  end do
Here is the caller graph for this function:

◆ initialisation_callback()

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

Initialises the swap and smooth component.

Parameters
current_stateThe current model state

Definition at line 31 of file swapsmooth.F90.

32  type(model_state_type), target, intent(inout) :: current_state
33 
34 #ifdef U_ACTIVE
35  mean_profiles_active=allocated(current_state%global_grid%configuration%vertical%olubar)
36  return
37 #endif
38 #ifdef V_ACTIVE
39  mean_profiles_active=allocated(current_state%global_grid%configuration%vertical%olvbar)
40  return
41 #endif
42  if (current_state%th%active) then
43  mean_profiles_active=allocated(current_state%global_grid%configuration%vertical%olthbar)
44  return
45  end if
46  if (current_state%number_q_fields .gt. 0) then
47  mean_profiles_active=allocated(current_state%global_grid%configuration%vertical%olqbar)
48  return
49  end if
Here is the caller graph for this function:

◆ robert_filter_for_average_profiles()

subroutine swapsmooth_mod::robert_filter_for_average_profiles ( type(model_state_type), intent(inout)  current_state)
private

Does swapping and smoothing (using robert filter algorithm) for the average profiles (the bars)

Parameters
current_stateThe current model state

Definition at line 239 of file swapsmooth.F90.

240  type(model_state_type), intent(inout) :: current_state
241 
242  integer :: k, n
243  real(kind=default_precision) :: c1, c2, existing_value
244 
245  c1 = 1.0_default_precision - 2.0_default_precision*current_state%tsmth
246  c2 = current_state%tsmth
247 
248  do k=1,current_state%global_grid%size(z_index)
249 #ifdef U_ACTIVE
250  existing_value=current_state%global_grid%configuration%vertical%olubar(k)
251  current_state%global_grid%configuration%vertical%olubar(k)=current_state%global_grid%configuration%vertical%olzubar(k)
252  current_state%global_grid%configuration%vertical%olzubar(k)=c1*existing_value+c2*&
253  (current_state%global_grid%configuration%vertical%olubar(k) + &
254  current_state%global_grid%configuration%vertical%savolubar(k))
255 #endif
256 #ifdef V_ACTIVE
257  existing_value=current_state%global_grid%configuration%vertical%olvbar(k)
258  current_state%global_grid%configuration%vertical%olvbar(k)=current_state%global_grid%configuration%vertical%olzvbar(k)
259  current_state%global_grid%configuration%vertical%olzvbar(k)=c1*existing_value+c2*&
260  (current_state%global_grid%configuration%vertical%olvbar(k) + &
261  current_state%global_grid%configuration%vertical%savolvbar(k))
262 #endif
263  if (current_state%th%active) then
264  existing_value=current_state%global_grid%configuration%vertical%olzthbar(k)
265  current_state%global_grid%configuration%vertical%olzthbar(k)=&
266  current_state%global_grid%configuration%vertical%olthbar(k) + current_state%tsmth * existing_value
267  current_state%global_grid%configuration%vertical%olthbar(k)=existing_value
268  end if
269  if (current_state%number_q_fields .gt. 0) then
270  do n=1, current_state%number_q_fields
271  existing_value=current_state%global_grid%configuration%vertical%olzqbar(k,n)
272  current_state%global_grid%configuration%vertical%olzqbar(k,n)=&
273  current_state%global_grid%configuration%vertical%olqbar(k,n) + current_state%tsmth * existing_value
274  current_state%global_grid%configuration%vertical%olqbar(k,n)=existing_value
275  end do
276  end if
277  end do
Here is the caller graph for this function:

◆ swap_and_smooth_classic()

subroutine swapsmooth_mod::swap_and_smooth_classic ( type(model_state_type), intent(inout)  current_state,
logical, intent(in)  old_smoother 
)
private

Classic swap and smooth based upon the old or no smoothing.

Parameters
current_stateThe current model state_mod
old_smootherWhether to use the old smoother or not (not means no smoothing)

Definition at line 78 of file swapsmooth.F90.

79  type(model_state_type), intent(inout) :: current_state
80  logical, intent(in) :: old_smoother
81 
82  integer :: y_index, x_index, k, n
83  real(kind=default_precision) :: c1, c2, existing_value
84 
85  if (old_smoother) then
86  c1 = 1.0_default_precision - current_state%tsmth
87  c2 = 2.0_default_precision * current_state%tsmth - 1.0_default_precision
88  else
89  c1 = 1.0_default_precision
90  c2 = -1.0_default_precision
91  end if
92 
93  x_index=current_state%column_local_x
94  y_index=current_state%column_local_y
95 
96  do k=1,current_state%global_grid%size(z_index)
97 #ifdef U_ACTIVE
98  existing_value = current_state%u%data(k,y_index,x_index) + current_state%zu%data(k,y_index,x_index)
99  current_state%u%data(k,y_index,x_index)=existing_value * c1 + current_state%u%data(k,y_index,x_index) * c2
100  current_state%zu%data(k,y_index,x_index)=existing_value - current_state%u%data(k,y_index,x_index)
101 #endif
102 #ifdef V_ACTIVE
103  existing_value = current_state%v%data(k,y_index,x_index) + current_state%zv%data(k,y_index,x_index)
104  current_state%v%data(k,y_index,x_index)=existing_value * c1 + current_state%v%data(k,y_index,x_index) * c2
105  current_state%zv%data(k,y_index,x_index)=existing_value - current_state%v%data(k,y_index,x_index)
106 #endif
107 #ifdef W_ACTIVE
108  existing_value = current_state%w%data(k,y_index,x_index) + current_state%zw%data(k,y_index,x_index)
109  current_state%w%data(k,y_index,x_index)=existing_value * c1 + current_state%w%data(k,y_index,x_index) * c2
110  current_state%zw%data(k,y_index,x_index)=existing_value - current_state%w%data(k,y_index,x_index)
111 #endif
112  if (current_state%th%active) then
113  existing_value = current_state%th%data(k,y_index,x_index) + current_state%zth%data(k,y_index,x_index)
114  current_state%th%data(k,y_index,x_index)=existing_value * c1 + current_state%th%data(k,y_index,x_index) * c2
115  current_state%zth%data(k,y_index,x_index)=existing_value - current_state%th%data(k,y_index,x_index)
116  end if
117  do n=1,current_state%number_q_fields
118  if (current_state%q(n)%active) then
119  existing_value = current_state%q(n)%data(k,y_index,x_index) + current_state%zq(n)%data(k,y_index,x_index)
120  current_state%q(n)%data(k,y_index,x_index)=existing_value * c1 + current_state%q(n)%data(k,y_index,x_index) * c2
121  current_state%zq(n)%data(k,y_index,x_index)=existing_value - current_state%q(n)%data(k,y_index,x_index)
122  end if
123  end do
124  end do
Here is the caller graph for this function:

◆ swap_and_smooth_robert_filter()

subroutine swapsmooth_mod::swap_and_smooth_robert_filter ( type(model_state_type), intent(inout)  current_state)
private

Swap and smooth with a Robert filter.

Parameters
current_stateThe current model state

Definition at line 129 of file swapsmooth.F90.

130  type(model_state_type), intent(inout) :: current_state
131 
132  integer :: y_index, x_index, k, n
133  real(kind=default_precision) :: c1, c2, existing_value
134 
135  x_index=current_state%column_local_x
136  y_index=current_state%column_local_y
137 
138  c1 = 1.0_default_precision - 2.0_default_precision*current_state%tsmth
139  c2 = current_state%tsmth
140 
141  do k=1,current_state%global_grid%size(z_index)
142 #ifdef U_ACTIVE
143  existing_value = current_state%u%data(k,y_index,x_index)
144  current_state%u%data(k,y_index,x_index)=current_state%zu%data(k,y_index,x_index)
145  current_state%zu%data(k,y_index,x_index)=c1*existing_value+c2*(current_state%u%data(k, y_index, x_index)+&
146  current_state%savu%data(k,y_index,x_index) -current_state%ugal)
147 #endif
148 #ifdef V_ACTIVE
149  existing_value = current_state%v%data(k,y_index,x_index)
150  current_state%v%data(k,y_index,x_index)=current_state%zv%data(k,y_index,x_index)
151  current_state%zv%data(k,y_index,x_index)=c1*existing_value+c2*(current_state%v%data(k, y_index, x_index)+&
152  current_state%savv%data(k,y_index,x_index)-current_state%vgal)
153 #endif
154 #ifdef W_ACTIVE
155  existing_value = current_state%w%data(k,y_index,x_index)
156  current_state%w%data(k,y_index,x_index)=current_state%zw%data(k,y_index,x_index)
157  current_state%zw%data(k,y_index,x_index)=c1*existing_value+c2*(current_state%w%data(k, y_index, x_index)+&
158  current_state%savw%data(k,y_index,x_index))
159 #endif
160  if (current_state%th%active) then
161  ! Uses the partial smooth of theta from stepfields
162  existing_value = current_state%zth%data(k,y_index,x_index)
163  current_state%zth%data(k,y_index,x_index)=current_state%th%data(k,y_index,x_index) + current_state%tsmth * existing_value
164  current_state%th%data(k,y_index,x_index)=existing_value
165  end if
166  do n=1, current_state%number_q_fields
167  if (current_state%q(n)%active) then
168  ! Uses the partial smooth of q from stepfields
169  existing_value = current_state%zq(n)%data(k,y_index,x_index)
170  current_state%zq(n)%data(k,y_index,x_index)=current_state%q(n)%data(k,y_index,x_index)+&
171  current_state%tsmth * existing_value
172  current_state%q(n)%data(k,y_index,x_index)=existing_value
173  end if
174  end do
175  end do
Here is the caller graph for this function:

◆ swapsmooth_get_descriptor()

type(component_descriptor_type) function, public swapsmooth_mod::swapsmooth_get_descriptor

Provides the descriptor back to the caller and is used in component registration.

Returns
The SwapSmooth component descriptor

Definition at line 22 of file swapsmooth.F90.

23  swapsmooth_get_descriptor%name="swap_smooth"
24  swapsmooth_get_descriptor%version=0.1
25  swapsmooth_get_descriptor%initialisation=>initialisation_callback
26  swapsmooth_get_descriptor%timestep=>timestep_callback
Here is the call graph for this function:

◆ timestep_callback()

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

Called for each non halo timestep column and will perform swapping and smoothing as required on that column.

Parameters
current_stateThe current model state_mod

Definition at line 54 of file swapsmooth.F90.

55  type(model_state_type), target, intent(inout) :: current_state
56 
57  if (.not. current_state%halo_column) then
58  if (current_state%field_stepping == forward_stepping) then
59  call swap_and_smooth_classic(current_state, .false.)
60  else
61  ! Centred stepping
62  call swap_and_smooth_robert_filter(current_state)
63  end if
64  end if
65 
66  if (mean_profiles_active .and. current_state%last_timestep_column) then
67  if (current_state%field_stepping == forward_stepping) then
68  call classic_for_average_profiles(current_state, .false.)
69  else
70  call robert_filter_for_average_profiles(current_state)
71  end if
72  end if
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ mean_profiles_active

logical swapsmooth_mod::mean_profiles_active =.false.
private

Whether or not mean profiles need smoothing.

Definition at line 14 of file swapsmooth.F90.

14  logical :: mean_profiles_active=.false.
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
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
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