MONC
Functions/Subroutines | Variables
pstep_mod Module Reference

Stepping of the pressure field. Completes the time-stepping of the velocity fields by adding the pressure term (dp/dx_i). In addition, ensures that l_zu and l_zv satisfy the Galilean-transformed boundary condition. This does not do the flow field _p terms which are only needed for diagnostics, nore does it do field halo swapping which is again only needed for diagnostics. More...

Functions/Subroutines

type(component_descriptor_type) function, public pstep_get_descriptor ()
 Descriptor of this component for registration. 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...
 
subroutine initialisation_callback (current_state)
 Initialisation callback hook which will check the diverr component is enabled (as this allocates p) More...
 
subroutine finalisation_callback (current_state)
 
subroutine timestep_callback (current_state)
 Called each timestep, this will step the pressure field for the non halo columns. More...
 
subroutine step_pressure_field (current_state)
 Does the actual stepping of the pressure field. More...
 
subroutine perform_galilean_transformation (current_state, y_index, x_index)
 Performs Galilean transformation of flow current and z fields. More...
 
subroutine save_precomponent_tendencies (current_state, cxn, cyn, txn, tyn)
 Save the 3d tendencies coming into the component. More...
 
subroutine compute_component_tendencies (current_state, cxn, cyn, txn, tyn)
 Computation of component tendencies. More...
 
subroutine set_published_field_value (field_value, real_1d_field, real_2d_field, real_3d_field)
 Sets the published field value from the temporary diagnostic values held by this component. More...
 

Variables

real(kind=default_precision), dimension(:,:,:), allocatable tendp_3d_u
 
real(kind=default_precision), dimension(:,:,:), allocatable tendp_3d_v
 
real(kind=default_precision), dimension(:,:,:), allocatable tendp_3d_w
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_u
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_v
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_w
 
logical l_tendp_3d_u
 
logical l_tendp_3d_v
 
logical l_tendp_3d_w
 
logical l_tend_3d_u
 
logical l_tend_3d_v
 
logical l_tend_3d_w
 
real(kind=default_precision), dimension(:), allocatable tendp_pr_tot_u
 
real(kind=default_precision), dimension(:), allocatable tendp_pr_tot_v
 
real(kind=default_precision), dimension(:), allocatable tendp_pr_tot_w
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_u
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_v
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_w
 
logical l_tendp_pr_tot_u
 
logical l_tendp_pr_tot_v
 
logical l_tendp_pr_tot_w
 
logical l_tend_pr_tot_u
 
logical l_tend_pr_tot_v
 
logical l_tend_pr_tot_w
 
integer diagnostic_generation_frequency
 

Detailed Description

Stepping of the pressure field. Completes the time-stepping of the velocity fields by adding the pressure term (dp/dx_i). In addition, ensures that l_zu and l_zv satisfy the Galilean-transformed boundary condition. This does not do the flow field _p terms which are only needed for diagnostics, nore does it do field halo swapping which is again only needed for diagnostics.

Function/Subroutine Documentation

◆ compute_component_tendencies()

subroutine pstep_mod::compute_component_tendencies ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  cxn,
integer, intent(in)  cyn,
integer, intent(in)  txn,
integer, intent(in)  tyn 
)
private

Computation of component tendencies.

Parameters
current_stateCurrent model state
cxnThe current slice, x, index
cynThe current column, y, index.
txntarget_x_index
tyntarget_y_index

Definition at line 476 of file pstep.F90.

477  type(model_state_type), target, intent(inout) :: current_state
478  integer, intent(in) :: cxn, cyn, txn, tyn
479  real(kind=default_precision) :: dtmtmp
480 
481  dtmtmp=2.0_default_precision* &
482  merge(current_state%dtm, 0.5_default_precision*current_state%dtm, current_state%field_stepping == centred_stepping)
483 
484  ! Calculate change in tendency due to component
485  if (l_tendp_3d_u) then
486  tendp_3d_u(:,tyn,txn)=(current_state%zu%data(:,cyn,cxn) - tendp_3d_u(:,tyn,txn))/dtmtmp
487  endif
488  if (l_tendp_3d_v) then
489  tendp_3d_v(:,tyn,txn)=(current_state%zv%data(:,cyn,cxn) - tendp_3d_v(:,tyn,txn))/dtmtmp
490  endif
491  if (l_tendp_3d_w) then
492  tendp_3d_w(:,tyn,txn)=(current_state%zw%data(:,cyn,cxn) - tendp_3d_w(:,tyn,txn))/dtmtmp
493  endif
494  if (l_tend_3d_u) then
495  tend_3d_u(:,tyn,txn)=tend_3d_u(:,tyn,txn) + tendp_3d_u(:,tyn,txn)
496  endif
497  if (l_tend_3d_v) then
498  tend_3d_v(:,tyn,txn)=tend_3d_v(:,tyn,txn) + tendp_3d_v(:,tyn,txn)
499  endif
500  if (l_tend_3d_w) then
501  tend_3d_w(:,tyn,txn)=tend_3d_w(:,tyn,txn) + tendp_3d_w(:,tyn,txn)
502  endif
503 
504  ! Add local tendency fields to the profile total
505  if (l_tendp_pr_tot_u) then
506  tendp_pr_tot_u(:)=tendp_pr_tot_u(:) + tendp_3d_u(:,tyn,txn)
507  endif
508  if (l_tendp_pr_tot_v) then
509  tendp_pr_tot_v(:)=tendp_pr_tot_v(:) + tendp_3d_v(:,tyn,txn)
510  endif
511  if (l_tendp_pr_tot_w) then
512  tendp_pr_tot_w(:)=tendp_pr_tot_w(:) + tendp_3d_w(:,tyn,txn)
513  endif
514  if (l_tend_pr_tot_u) then
515  tend_pr_tot_u(:)=tend_pr_tot_u(:) + tend_3d_u(:,tyn,txn)
516  endif
517  if (l_tend_pr_tot_v) then
518  tend_pr_tot_v(:)=tend_pr_tot_v(:) + tend_3d_v(:,tyn,txn)
519  endif
520  if (l_tend_pr_tot_w) then
521  tend_pr_tot_w(:)=tend_pr_tot_w(:) + tend_3d_w(:,tyn,txn)
522  endif
523 
Here is the caller graph for this function:

◆ field_information_retrieval_callback()

subroutine pstep_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 75 of file pstep.F90.

76  type(model_state_type), target, intent(inout) :: current_state
77  character(len=*), intent(in) :: name
78  type(component_field_information_type), intent(out) :: field_information
79  integer :: strcomp
80 
81  ! Field information for 3d
82  strcomp=index(name, "_pstep_3d_local")
83  if (strcomp .ne. 0) then
84  field_information%field_type=component_array_field_type
85  field_information%number_dimensions=3
86  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
87  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
88  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
89  field_information%data_type=component_double_data_type
90 
91  if (name .eq. "tend_u_pstep_3d_local") then
92  field_information%enabled=l_tend_3d_u
93  else if (name .eq. "tend_v_pstep_3d_local") then
94  field_information%enabled=l_tend_3d_v
95  else if (name .eq. "tend_w_pstep_3d_local") then
96  field_information%enabled=l_tend_3d_w
97  else
98  field_information%enabled=.true.
99  end if
100 
101  end if !end 3d check
102 
103  ! Field information for profiles
104  strcomp=index(name, "_pstep_profile_total_local")
105  if (strcomp .ne. 0) then
106  field_information%field_type=component_array_field_type
107  field_information%number_dimensions=1
108  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
109  field_information%data_type=component_double_data_type
110 
111  if (name .eq. "tend_u_pstep_profile_total_local") then
112  field_information%enabled=l_tend_pr_tot_u
113  else if (name .eq. "tend_v_pstep_profile_total_local") then
114  field_information%enabled=l_tend_pr_tot_v
115  else if (name .eq. "tend_w_pstep_profile_total_local") then
116  field_information%enabled=l_tend_pr_tot_w
117  else
118  field_information%enabled=.true.
119  end if
120 
121  end if !end profile check
122 
123  ! Field information for 3d totals
124  strcomp=index(name, "_total_3d_local")
125  if (strcomp .ne. 0) then
126  field_information%field_type=component_array_field_type
127  field_information%number_dimensions=3
128  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
129  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
130  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
131  field_information%data_type=component_double_data_type
132 
133  if (name .eq. "tend_u_total_3d_local") then
134  field_information%enabled=l_tend_3d_u
135  else if (name .eq. "tend_v_total_3d_local") then
136  field_information%enabled=l_tend_3d_v
137  else if (name .eq. "tend_w_total_3d_local") then
138  field_information%enabled=l_tend_3d_w
139  else
140  field_information%enabled=.true.
141  end if
142 
143  end if !end 3d check
144 
145  ! Field information for total profiles
146  strcomp=index(name, "_total_profile_total_local")
147  if (strcomp .ne. 0) then
148  field_information%field_type=component_array_field_type
149  field_information%number_dimensions=1
150  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
151  field_information%data_type=component_double_data_type
152 
153  if (name .eq. "tend_u_total_profile_total_local") then
154  field_information%enabled=l_tend_pr_tot_u
155  else if (name .eq. "tend_v_total_profile_total_local") then
156  field_information%enabled=l_tend_pr_tot_v
157  else if (name .eq. "tend_w_total_profile_total_local") then
158  field_information%enabled=l_tend_pr_tot_w
159  else
160  field_information%enabled=.true.
161  end if
162 
163  end if !end profile check
164 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine pstep_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 172 of file pstep.F90.

173  type(model_state_type), target, intent(inout) :: current_state
174  character(len=*), intent(in) :: name
175  type(component_field_value_type), intent(out) :: field_value
176  ! 3d Tendency Fields
177 
178  if (name .eq. "tend_u_pstep_3d_local" .and. allocated(tendp_3d_u)) then
179  call set_published_field_value(field_value, real_3d_field=tendp_3d_u)
180  else if (name .eq. "tend_v_pstep_3d_local" .and. allocated(tendp_3d_v)) then
181  call set_published_field_value(field_value, real_3d_field=tendp_3d_v)
182  else if (name .eq. "tend_w_pstep_3d_local" .and. allocated(tendp_3d_w)) then
183  call set_published_field_value(field_value, real_3d_field=tendp_3d_w)
184 
185  ! Profile Tendency Fields
186  else if (name .eq. "tend_u_pstep_profile_total_local" .and. allocated(tendp_pr_tot_u)) then
187  call set_published_field_value(field_value, real_1d_field=tendp_pr_tot_u)
188  else if (name .eq. "tend_v_pstep_profile_total_local" .and. allocated(tendp_pr_tot_v)) then
189  call set_published_field_value(field_value, real_1d_field=tendp_pr_tot_v)
190  else if (name .eq. "tend_w_pstep_profile_total_local" .and. allocated(tendp_pr_tot_w)) then
191  call set_published_field_value(field_value, real_1d_field=tendp_pr_tot_w)
192 
193  ! 3d Tendency Fields
194  else if (name .eq. "tend_u_total_3d_local" .and. allocated(tend_3d_u)) then
195  call set_published_field_value(field_value, real_3d_field=tend_3d_u)
196  else if (name .eq. "tend_v_total_3d_local" .and. allocated(tend_3d_v)) then
197  call set_published_field_value(field_value, real_3d_field=tend_3d_v)
198  else if (name .eq. "tend_w_total_3d_local" .and. allocated(tend_3d_w)) then
199  call set_published_field_value(field_value, real_3d_field=tend_3d_w)
200 
201  ! Profile Tendency Fields
202  else if (name .eq. "tend_u_total_profile_total_local" .and. allocated(tend_pr_tot_u)) then
203  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_u)
204  else if (name .eq. "tend_v_total_profile_total_local" .and. allocated(tend_pr_tot_v)) then
205  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_v)
206  else if (name .eq. "tend_w_total_profile_total_local" .and. allocated(tend_pr_tot_w)) then
207  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_w)
208  end if
209 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

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

Definition at line 294 of file pstep.F90.

295  type(model_state_type), target, intent(inout) :: current_state
296 
297  if (allocated(tendp_3d_u)) deallocate(tendp_3d_u)
298  if (allocated(tendp_3d_v)) deallocate(tendp_3d_v)
299  if (allocated(tendp_3d_w)) deallocate(tendp_3d_w)
300  if (allocated(tend_3d_u)) deallocate(tend_3d_u)
301  if (allocated(tend_3d_v)) deallocate(tend_3d_v)
302  if (allocated(tend_3d_w)) deallocate(tend_3d_w)
303 
304  if (allocated(tendp_pr_tot_u)) deallocate(tendp_pr_tot_u)
305  if (allocated(tendp_pr_tot_v)) deallocate(tendp_pr_tot_v)
306  if (allocated(tendp_pr_tot_w)) deallocate(tendp_pr_tot_w)
307  if (allocated(tend_pr_tot_u)) deallocate(tend_pr_tot_u)
308  if (allocated(tend_pr_tot_v)) deallocate(tend_pr_tot_v)
309  if (allocated(tend_pr_tot_w)) deallocate(tend_pr_tot_w)
310 
Here is the caller graph for this function:

◆ initialisation_callback()

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

Initialisation callback hook which will check the diverr component is enabled (as this allocates p)

Parameters
current_stateThe current model state

Definition at line 214 of file pstep.F90.

215  type(model_state_type), target, intent(inout) :: current_state
216 
217  if (.not. is_component_enabled(current_state%options_database, "diverr")) then
218  call log_master_log(log_error, "The pstep component requires the diverr component to be enabled")
219  end if
220 
221  ! Tendency Logicals
222  l_tendp_pr_tot_u = current_state%u%active
223  l_tendp_pr_tot_v = current_state%v%active
224  l_tendp_pr_tot_w = current_state%w%active
225  l_tend_pr_tot_u = current_state%u%active
226  l_tend_pr_tot_v = current_state%v%active
227  l_tend_pr_tot_w = current_state%w%active
228 
229  l_tendp_3d_u = current_state%u%active .or. l_tendp_pr_tot_u
230  l_tendp_3d_v = current_state%v%active .or. l_tendp_pr_tot_v
231  l_tendp_3d_w = current_state%w%active .or. l_tendp_pr_tot_w
232  l_tend_3d_u = current_state%u%active .or. l_tend_pr_tot_u
233  l_tend_3d_v = current_state%v%active .or. l_tend_pr_tot_v
234  l_tend_3d_w = current_state%w%active .or. l_tend_pr_tot_w
235 
236  ! Allocate 3d tendency fields upon availability
237  if (l_tendp_3d_u) then
238  allocate( tendp_3d_u(current_state%local_grid%size(z_index), &
239  current_state%local_grid%size(y_index), &
240  current_state%local_grid%size(x_index) ) )
241  endif
242  if (l_tendp_3d_v) then
243  allocate( tendp_3d_v(current_state%local_grid%size(z_index), &
244  current_state%local_grid%size(y_index), &
245  current_state%local_grid%size(x_index) ) )
246  endif
247  if (l_tendp_3d_w) then
248  allocate( tendp_3d_w(current_state%local_grid%size(z_index), &
249  current_state%local_grid%size(y_index), &
250  current_state%local_grid%size(x_index) ) )
251  endif
252  if (l_tend_3d_u) then
253  allocate( tend_3d_u(current_state%local_grid%size(z_index), &
254  current_state%local_grid%size(y_index), &
255  current_state%local_grid%size(x_index) ) )
256  endif
257  if (l_tend_3d_v) then
258  allocate( tend_3d_v(current_state%local_grid%size(z_index), &
259  current_state%local_grid%size(y_index), &
260  current_state%local_grid%size(x_index) ) )
261  endif
262  if (l_tend_3d_w) then
263  allocate( tend_3d_w(current_state%local_grid%size(z_index), &
264  current_state%local_grid%size(y_index), &
265  current_state%local_grid%size(x_index) ) )
266  endif
267 
268  ! Allocate profile tendency fields upon availability
269  if (l_tendp_pr_tot_u) then
270  allocate( tendp_pr_tot_u(current_state%local_grid%size(z_index)) )
271  endif
272  if (l_tendp_pr_tot_v) then
273  allocate( tendp_pr_tot_v(current_state%local_grid%size(z_index)) )
274  endif
275  if (l_tendp_pr_tot_w) then
276  allocate( tendp_pr_tot_w(current_state%local_grid%size(z_index)) )
277  endif
278  if (l_tend_pr_tot_u) then
279  allocate( tend_pr_tot_u(current_state%local_grid%size(z_index)) )
280  endif
281  if (l_tend_pr_tot_v) then
282  allocate( tend_pr_tot_v(current_state%local_grid%size(z_index)) )
283  endif
284  if (l_tend_pr_tot_w) then
285  allocate( tend_pr_tot_w(current_state%local_grid%size(z_index)) )
286  endif
287 
288  ! Save the sampling_frequency to force diagnostic calculation on select time steps
289  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
290 
Here is the caller graph for this function:

◆ perform_galilean_transformation()

subroutine pstep_mod::perform_galilean_transformation ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  y_index,
integer, intent(in)  x_index 
)
private

Performs Galilean transformation of flow current and z fields.

Parameters
current_stateThe current model state
y_indexLocal y index which we are working with
x_indexLocal x index which we are working with

Definition at line 418 of file pstep.F90.

419  type(model_state_type), target, intent(inout) :: current_state
420  integer, intent(in) :: y_index, x_index
421 
422  integer :: k
423 
424  do k=1,current_state%local_grid%size(z_index)
425 #ifdef U_ACTIVE
426  current_state%u%data(k, y_index, x_index)= current_state%u%data(k, y_index, x_index)-current_state%ugal
427  current_state%zu%data(k, y_index, x_index)= current_state%zu%data(k, y_index, x_index)-current_state%ugal
428 #endif
429 #ifdef V_ACTIVE
430  current_state%v%data(k, y_index, x_index)= current_state%v%data(k, y_index, x_index)-current_state%vgal
431  current_state%zv%data(k, y_index, x_index)= current_state%zv%data(k, y_index, x_index)-current_state%vgal
432 #endif
433  end do
Here is the caller graph for this function:

◆ pstep_get_descriptor()

type(component_descriptor_type) function, public pstep_mod::pstep_get_descriptor

Descriptor of this component for registration.

Returns
The pstep component descriptor

Definition at line 44 of file pstep.F90.

45  pstep_get_descriptor%name="pstep"
46  pstep_get_descriptor%version=0.1
47  pstep_get_descriptor%initialisation=>initialisation_callback
48  pstep_get_descriptor%timestep=>timestep_callback
49  pstep_get_descriptor%finalisation=>finalisation_callback
50 
51  pstep_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
52  pstep_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
53  allocate(pstep_get_descriptor%published_fields(6+6))
54 
55  pstep_get_descriptor%published_fields(1)= "tend_u_pstep_3d_local"
56  pstep_get_descriptor%published_fields(2)= "tend_v_pstep_3d_local"
57  pstep_get_descriptor%published_fields(3)= "tend_w_pstep_3d_local"
58  pstep_get_descriptor%published_fields(4)= "tend_u_total_3d_local"
59  pstep_get_descriptor%published_fields(5)= "tend_v_total_3d_local"
60  pstep_get_descriptor%published_fields(6)= "tend_w_total_3d_local"
61 
62  pstep_get_descriptor%published_fields(6+1)= "tend_u_pstep_profile_total_local"
63  pstep_get_descriptor%published_fields(6+2)= "tend_v_pstep_profile_total_local"
64  pstep_get_descriptor%published_fields(6+3)= "tend_w_pstep_profile_total_local"
65  pstep_get_descriptor%published_fields(6+4)= "tend_u_total_profile_total_local"
66  pstep_get_descriptor%published_fields(6+5)= "tend_v_total_profile_total_local"
67  pstep_get_descriptor%published_fields(6+6)= "tend_w_total_profile_total_local"
68 
Here is the call graph for this function:

◆ save_precomponent_tendencies()

subroutine pstep_mod::save_precomponent_tendencies ( type(model_state_type), intent(in), target  current_state,
integer, intent(in)  cxn,
integer, intent(in)  cyn,
integer, intent(in)  txn,
integer, intent(in)  tyn 
)
private

Save the 3d tendencies coming into the component.

Parameters
current_stateCurrent model state
cxnThe current slice, x, index
cynThe current column, y, index.
txntarget_x_index
tyntarget_y_index

Definition at line 443 of file pstep.F90.

444  type(model_state_type), target, intent(in) :: current_state
445  integer, intent(in) :: cxn, cyn, txn, tyn
446 
447  ! Save 3d tendency fields upon request (of 3d or profiles) and availability
448  if (l_tendp_3d_u) then
449  tendp_3d_u(:,tyn,txn)=current_state%zu%data(:,cyn,cxn)
450  endif
451  if (l_tendp_3d_v) then
452  tendp_3d_v(:,tyn,txn)=current_state%zv%data(:,cyn,cxn)
453  endif
454  if (l_tendp_3d_w) then
455  tendp_3d_w(:,tyn,txn)=current_state%zw%data(:,cyn,cxn)
456  endif
457  if (l_tend_3d_u) then
458  tend_3d_u(:,tyn,txn)=current_state%su%data(:,cyn,cxn)
459  endif
460  if (l_tend_3d_v) then
461  tend_3d_v(:,tyn,txn)=current_state%sv%data(:,cyn,cxn)
462  endif
463  if (l_tend_3d_w) then
464  tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn)
465  endif
466 
Here is the caller graph for this function:

◆ set_published_field_value()

subroutine pstep_mod::set_published_field_value ( type(component_field_value_type), intent(inout)  field_value,
real(kind=default_precision), dimension(:), optional  real_1d_field,
real(kind=default_precision), dimension(:,:), optional  real_2d_field,
real(kind=default_precision), dimension(:,:,:), optional  real_3d_field 
)
private

Sets the published field value from the temporary diagnostic values held by this component.

Parameters
field_valuePopulated with the value of the field
real_1d_fieldOptional one dimensional real of values to publish
real_2d_fieldOptional two dimensional real of values to publish

Definition at line 531 of file pstep.F90.

532  type(component_field_value_type), intent(inout) :: field_value
533  real(kind=default_precision), dimension(:), optional :: real_1d_field
534  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
535  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
536 
537  if (present(real_1d_field)) then
538  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
539  else if (present(real_2d_field)) then
540  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
541  else if (present(real_3d_field)) then
542  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
543  source=real_3d_field)
544  end if
Here is the caller graph for this function:

◆ step_pressure_field()

subroutine pstep_mod::step_pressure_field ( type(model_state_type), intent(inout), target  current_state)
private

Does the actual stepping of the pressure field.

Parameters
current_stateThe current model state

Definition at line 364 of file pstep.F90.

365  type(model_state_type), target, intent(inout) :: current_state
366 
367  integer :: k, x_index, y_index
368  real(kind=default_precision) :: dtmtmp
369 
370  x_index=current_state%column_local_x
371  y_index=current_state%column_local_y
372 
373  dtmtmp=merge(current_state%dtm, 0.5_default_precision*current_state%dtm, current_state%field_stepping == centred_stepping)
374  do k=2,current_state%local_grid%size(z_index)
375 
376 #ifdef U_ACTIVE
377  current_state%zu%data(k, y_index, x_index)= current_state%zu%data(k, y_index, x_index)+ 2.0_default_precision*&
378  current_state%global_grid%configuration%horizontal%cx*dtmtmp*(current_state%p%data(k, y_index, x_index)-&
379  current_state%p%data(k, y_index, x_index+1))
380 #endif
381 #ifdef V_ACTIVE
382  current_state%zv%data(k, y_index, x_index)=&
383  current_state%zv%data(k, y_index, x_index)+2.0_default_precision*&
384  current_state%global_grid%configuration%horizontal%cy*dtmtmp*&
385  (current_state%p%data(k, y_index, x_index) - current_state%p%data(k, y_index+1, x_index))
386 #endif
387 #ifdef W_ACTIVE
388  if (k .lt. current_state%local_grid%size(z_index)) then
389  current_state%zw%data(k, y_index, x_index)=current_state%zw%data(k, y_index, x_index)+2.0_default_precision*&
390  current_state%global_grid%configuration%vertical%rdzn(k+1)*dtmtmp*(current_state%p%data(k, y_index, x_index)-&
391  current_state%p%data(k+1, y_index, x_index))
392  end if
393 #endif
394  end do
395  if (current_state%use_viscosity_and_diffusion .and. current_state%use_surface_boundary_conditions) then
396 #ifdef U_ACTIVE
397  current_state%zu%data(1, y_index, x_index)=-current_state%zu%data(2, y_index, x_index)-&
398  2.0_default_precision*current_state%ugal
399 #endif
400 #ifdef V_ACTIVE
401  current_state%zv%data(1, y_index, x_index)=-current_state%zv%data(2, y_index, x_index)-&
402  2.0_default_precision*current_state%vgal
403 #endif
404  else
405 #ifdef U_ACTIVE
406  current_state%zu%data(1, y_index, x_index)=current_state%zu%data(2, y_index, x_index)
407 #endif
408 #ifdef V_ACTIVE
409  current_state%zv%data(1, y_index, x_index)=current_state%zv%data(2, y_index, x_index)
410 #endif
411  end if
Here is the caller graph for this function:

◆ timestep_callback()

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

Called each timestep, this will step the pressure field for the non halo columns.

Parameters
current_stateThe current model state

Definition at line 315 of file pstep.F90.

316  type(model_state_type), target, intent(inout) :: current_state
317  integer :: local_y, local_x, target_x_index, target_y_index
318 
319  local_y=current_state%column_local_y
320  local_x=current_state%column_local_x
321  target_y_index=local_y-current_state%local_grid%halo_size(y_index)
322  target_x_index=local_x-current_state%local_grid%halo_size(x_index)
323 
324  ! Zero profile tendency totals on first instance in the sum
325  if (current_state%first_timestep_column) then
326  if (l_tendp_pr_tot_u) then
327  tendp_pr_tot_u(:)= 0.0_default_precision
328  endif
329  if (l_tendp_pr_tot_v) then
330  tendp_pr_tot_v(:)= 0.0_default_precision
331  endif
332  if (l_tendp_pr_tot_w) then
333  tendp_pr_tot_w(:)= 0.0_default_precision
334  endif
335  if (l_tend_pr_tot_u) then
336  tend_pr_tot_u(:)= 0.0_default_precision
337  endif
338  if (l_tend_pr_tot_v) then
339  tend_pr_tot_v(:)= 0.0_default_precision
340  endif
341  if (l_tend_pr_tot_w) then
342  tend_pr_tot_w(:)= 0.0_default_precision
343  endif
344  endif ! zero totals
345 
346  if (current_state%galilean_transformation) call perform_galilean_transformation(current_state, &
347  current_state%column_local_y, current_state%column_local_x)
348 
349  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
350  call save_precomponent_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
351  end if
352 
353  if (.not. current_state%halo_column) call step_pressure_field(current_state)
354 
355  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
356  call compute_component_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
357  end if
358 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ diagnostic_generation_frequency

integer pstep_mod::diagnostic_generation_frequency
private

Definition at line 36 of file pstep.F90.

36  integer :: diagnostic_generation_frequency

◆ l_tend_3d_u

logical pstep_mod::l_tend_3d_u
private

Definition at line 30 of file pstep.F90.

◆ l_tend_3d_v

logical pstep_mod::l_tend_3d_v
private

Definition at line 30 of file pstep.F90.

◆ l_tend_3d_w

logical pstep_mod::l_tend_3d_w
private

Definition at line 30 of file pstep.F90.

◆ l_tend_pr_tot_u

logical pstep_mod::l_tend_pr_tot_u
private

Definition at line 34 of file pstep.F90.

◆ l_tend_pr_tot_v

logical pstep_mod::l_tend_pr_tot_v
private

Definition at line 34 of file pstep.F90.

◆ l_tend_pr_tot_w

logical pstep_mod::l_tend_pr_tot_w
private

Definition at line 34 of file pstep.F90.

◆ l_tendp_3d_u

logical pstep_mod::l_tendp_3d_u
private

Definition at line 30 of file pstep.F90.

30  logical :: l_tendp_3d_u, l_tendp_3d_v, l_tendp_3d_w, l_tend_3d_u, l_tend_3d_v, l_tend_3d_w

◆ l_tendp_3d_v

logical pstep_mod::l_tendp_3d_v
private

Definition at line 30 of file pstep.F90.

◆ l_tendp_3d_w

logical pstep_mod::l_tendp_3d_w
private

Definition at line 30 of file pstep.F90.

◆ l_tendp_pr_tot_u

logical pstep_mod::l_tendp_pr_tot_u
private

Definition at line 34 of file pstep.F90.

34  logical :: l_tendp_pr_tot_u, l_tendp_pr_tot_v, l_tendp_pr_tot_w, l_tend_pr_tot_u, l_tend_pr_tot_v, l_tend_pr_tot_w

◆ l_tendp_pr_tot_v

logical pstep_mod::l_tendp_pr_tot_v
private

Definition at line 34 of file pstep.F90.

◆ l_tendp_pr_tot_w

logical pstep_mod::l_tendp_pr_tot_w
private

Definition at line 34 of file pstep.F90.

◆ tend_3d_u

real(kind=default_precision), dimension(:,:,:), allocatable pstep_mod::tend_3d_u
private

Definition at line 28 of file pstep.F90.

◆ tend_3d_v

real(kind=default_precision), dimension(:,:,:), allocatable pstep_mod::tend_3d_v
private

Definition at line 28 of file pstep.F90.

◆ tend_3d_w

real(kind=default_precision), dimension(:,:,:), allocatable pstep_mod::tend_3d_w
private

Definition at line 28 of file pstep.F90.

◆ tend_pr_tot_u

real(kind=default_precision), dimension(:), allocatable pstep_mod::tend_pr_tot_u
private

Definition at line 32 of file pstep.F90.

◆ tend_pr_tot_v

real(kind=default_precision), dimension(:), allocatable pstep_mod::tend_pr_tot_v
private

Definition at line 32 of file pstep.F90.

◆ tend_pr_tot_w

real(kind=default_precision), dimension(:), allocatable pstep_mod::tend_pr_tot_w
private

Definition at line 32 of file pstep.F90.

◆ tendp_3d_u

real(kind=default_precision), dimension(:,:,:), allocatable pstep_mod::tendp_3d_u
private

Definition at line 28 of file pstep.F90.

28  real(kind=default_precision), dimension(:,:,:), allocatable :: &
29  tendp_3d_u, tendp_3d_v, tendp_3d_w, tend_3d_u, tend_3d_v, tend_3d_w

◆ tendp_3d_v

real(kind=default_precision), dimension(:,:,:), allocatable pstep_mod::tendp_3d_v
private

Definition at line 28 of file pstep.F90.

◆ tendp_3d_w

real(kind=default_precision), dimension(:,:,:), allocatable pstep_mod::tendp_3d_w
private

Definition at line 28 of file pstep.F90.

◆ tendp_pr_tot_u

real(kind=default_precision), dimension(:), allocatable pstep_mod::tendp_pr_tot_u
private

Definition at line 32 of file pstep.F90.

32  real(kind=default_precision), dimension(:), allocatable :: &
33  tendp_pr_tot_u, tendp_pr_tot_v, tendp_pr_tot_w, tend_pr_tot_u, tend_pr_tot_v, tend_pr_tot_w

◆ tendp_pr_tot_v

real(kind=default_precision), dimension(:), allocatable pstep_mod::tendp_pr_tot_v
private

Definition at line 32 of file pstep.F90.

◆ tendp_pr_tot_w

real(kind=default_precision), dimension(:), allocatable pstep_mod::tendp_pr_tot_w
private

Definition at line 32 of file pstep.F90.

logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
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_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