MONC
Functions/Subroutines
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 initialisation_callback (current_state)
 Initialisation callback hook which will check the diverr component is enabled (as this allocates p) More...
 
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...
 

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

◆ 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 34 of file pstep.F90.

34  type(model_state_type), target, intent(inout) :: current_state
35 
36  if (.not. is_component_enabled(current_state%options_database, "diverr")) then
37  call log_master_log(log_error, "The pstep component requires the diverr component to be enabled")
38  end if
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 108 of file pstep.F90.

108  type(model_state_type), target, intent(inout) :: current_state
109  integer, intent(in) :: y_index, x_index
110 
111  integer :: k
112 
113  do k=1,current_state%local_grid%size(z_index)
114 #ifdef U_ACTIVE
115  current_state%u%data(k, y_index, x_index)= current_state%u%data(k, y_index, x_index)-current_state%ugal
116  current_state%zu%data(k, y_index, x_index)= current_state%zu%data(k, y_index, x_index)-current_state%ugal
117 #endif
118 #ifdef V_ACTIVE
119  current_state%v%data(k, y_index, x_index)= current_state%v%data(k, y_index, x_index)-current_state%vgal
120  current_state%zv%data(k, y_index, x_index)= current_state%zv%data(k, y_index, x_index)-current_state%vgal
121 #endif
122  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 25 of file pstep.F90.

25  pstep_get_descriptor%name="pstep"
26  pstep_get_descriptor%version=0.1
27  pstep_get_descriptor%initialisation=>initialisation_callback
28  pstep_get_descriptor%timestep=>timestep_callback
Here is the call 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 54 of file pstep.F90.

54  type(model_state_type), target, intent(inout) :: current_state
55 
56  integer :: k, x_index, y_index
57  real(kind=DEFAULT_PRECISION) :: dtmtmp
58 
59  x_index=current_state%column_local_x
60  y_index=current_state%column_local_y
61 
62  dtmtmp=merge(current_state%dtm, 0.5_default_precision*current_state%dtm, current_state%field_stepping == centred_stepping)
63  do k=2,current_state%local_grid%size(z_index)
64 
65 #ifdef U_ACTIVE
66  current_state%zu%data(k, y_index, x_index)= current_state%zu%data(k, y_index, x_index)+ 2.0_default_precision*&
67  current_state%global_grid%configuration%horizontal%cx*dtmtmp*(current_state%p%data(k, y_index, x_index)-&
68  current_state%p%data(k, y_index, x_index+1))
69 #endif
70 #ifdef V_ACTIVE
71  current_state%zv%data(k, y_index, x_index)=&
72  current_state%zv%data(k, y_index, x_index)+2.0_default_precision*&
73  current_state%global_grid%configuration%horizontal%cy*dtmtmp*&
74  (current_state%p%data(k, y_index, x_index) - current_state%p%data(k, y_index+1, x_index))
75 #endif
76 #ifdef W_ACTIVE
77  if (k .lt. current_state%local_grid%size(z_index)) then
78  current_state%zw%data(k, y_index, x_index)=current_state%zw%data(k, y_index, x_index)+2.0_default_precision*&
79  current_state%global_grid%configuration%vertical%rdzn(k+1)*dtmtmp*(current_state%p%data(k, y_index, x_index)-&
80  current_state%p%data(k+1, y_index, x_index))
81  end if
82 #endif
83  end do
84  if (current_state%use_viscosity_and_diffusion .and. current_state%use_surface_boundary_conditions) then
85 #ifdef U_ACTIVE
86  current_state%zu%data(1, y_index, x_index)=-current_state%zu%data(2, y_index, x_index)-&
87  2.0_default_precision*current_state%ugal
88 #endif
89 #ifdef V_ACTIVE
90  current_state%zv%data(1, y_index, x_index)=-current_state%zv%data(2, y_index, x_index)-&
91  2.0_default_precision*current_state%vgal
92 #endif
93  else
94 #ifdef U_ACTIVE
95  current_state%zu%data(1, y_index, x_index)=current_state%zu%data(2, y_index, x_index)
96 #endif
97 #ifdef V_ACTIVE
98  current_state%zv%data(1, y_index, x_index)=current_state%zv%data(2, y_index, x_index)
99 #endif
100  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 44 of file pstep.F90.

44  type(model_state_type), target, intent(inout) :: current_state
45 
46  if (current_state%galilean_transformation) call perform_galilean_transformation(current_state, &
47  current_state%column_local_y, current_state%column_local_x)
48  if (.not. current_state%halo_column) call step_pressure_field(current_state)
Here is the call graph for this function:
Here is the caller graph for this function: