MONC
Functions/Subroutines | Variables
pressuresource_mod Module Reference

Calculates the gradient of the source flow fields (SU, SV, SW.) This is based upon the P field values already set for the divergence error in the diverr component. Note that some communication is required, as a process needs to know the x and y -1 values, which are held on the -1 neighbour in that dimension. These are computed here and sent to neighbour +1, the recvs from a process to neighbour -1 are also registered here and the handles are waited upon the the results combined into P in the solver. If the neighbour is local in any dimension then it is just a local memory copy. More...

Functions/Subroutines

type(component_descriptor_type) function, public pressuresource_get_descriptor ()
 Descriptor of this component for registration. More...
 
subroutine initialisation_callback (current_state)
 On initialisation this will allocate the buffer areas required and set communication handles to null. More...
 
subroutine timestep_callback (current_state)
 The timestep callback will update the values of P for each column. More...
 
subroutine finalisation_callback (current_state)
 Frees up the allocated buffers (if such were allocated) More...
 
subroutine calculate_psrce (current_state)
 Combines the source fields with the pressure values. For U and V, if this is on the low boundary then delay dealing with the -1 values as they will be communicated. Equally, for those fields if this is the high boundary then compute the values and send them on to the neighbouring process. More...
 
subroutine send_neighbouring_pressure_data (current_state)
 Sends the computed source pressure data terms to the p+1 process. More...
 
subroutine register_neighbouring_pressure_data_recv (current_state)
 Registers the receive requests for each neighbouring process if that is not local, is recieves from p-1. More...
 

Variables

real(kind=default_precision), dimension(:,:), allocatable send_buffer_x
 
real(kind=default_precision), dimension(:,:), allocatable send_buffer_y
 

Detailed Description

Calculates the gradient of the source flow fields (SU, SV, SW.) This is based upon the P field values already set for the divergence error in the diverr component. Note that some communication is required, as a process needs to know the x and y -1 values, which are held on the -1 neighbour in that dimension. These are computed here and sent to neighbour +1, the recvs from a process to neighbour -1 are also registered here and the handles are waited upon the the results combined into P in the solver. If the neighbour is local in any dimension then it is just a local memory copy.

Function/Subroutine Documentation

◆ calculate_psrce()

subroutine pressuresource_mod::calculate_psrce ( type(model_state_type), intent(inout), target  current_state)
private

Combines the source fields with the pressure values. For U and V, if this is on the low boundary then delay dealing with the -1 values as they will be communicated. Equally, for those fields if this is the high boundary then compute the values and send them on to the neighbouring process.

Parameters
current_stateThe current model state

Definition at line 84 of file pressuresource.F90.

85  type(model_state_type), target, intent(inout) :: current_state
86 
87  integer :: k, local_y, local_x, corrected_y, corrected_x
88  logical :: last_x, last_y
89 
90  local_y=current_state%column_local_y
91  local_x=current_state%column_local_x
92  last_y = local_y == current_state%local_grid%local_domain_end_index(y_index)
93  last_x = local_x == current_state%local_grid%local_domain_end_index(x_index)
94  if (last_x .or. last_y) then
95  corrected_x=local_x-current_state%local_grid%halo_size(x_index)
96  corrected_y=local_y-current_state%local_grid%halo_size(y_index)
97  end if
98  do k=2,current_state%local_grid%size(z_index)
99 #ifdef W_ACTIVE
100  current_state%p%data(k, local_y, local_x)=current_state%p%data(k, local_y, local_x)+&
101  4.0_default_precision*(current_state%global_grid%configuration%vertical%tzc2(k)*&
102  current_state%sw%data(k, local_y, local_x)-&
103  current_state%global_grid%configuration%vertical%tzc1(k)*current_state%sw%data(k-1, local_y, local_x))
104 #endif
105 #ifdef U_ACTIVE
106  current_state%p%data(k, local_y, local_x)=current_state%p%data(k, local_y, local_x)+&
107  current_state%global_grid%configuration%horizontal%cx * current_state%su%data(k, local_y, local_x)
108 #endif
109 #ifdef V_ACTIVE
110  current_state%p%data(k, local_y, local_x)=current_state%p%data(k, local_y, local_x)+&
111  current_state%global_grid%configuration%horizontal%cy * current_state%sv%data(k, local_y, local_x)
112 #endif
113 #ifdef U_ACTIVE
114  if (local_x .gt. 3) then
115  current_state%p%data(k, local_y, local_x)=current_state%p%data(k, local_y, local_x)-&
116  current_state%global_grid%configuration%horizontal%cx * current_state%su%data(k, local_y, local_x-1)
117  end if
118 
119  if (last_x) then
120  send_buffer_x(k-1, corrected_y)=&
121  current_state%global_grid%configuration%horizontal%cx * current_state%su%data(k, local_y, local_x)
122  end if
123 #endif
124 #ifdef V_ACTIVE
125  if (local_y .gt. 3 .and. local_x .gt. 3) then
126  current_state%p%data(k, local_y, local_x)=current_state%p%data(k, local_y, local_x)-&
127  current_state%global_grid%configuration%horizontal%cy * current_state%sv%data(k, local_y-1, local_x)
128  end if
129  if (last_y) then
130  send_buffer_y(k-1, corrected_x)=&
131  current_state%global_grid%configuration%horizontal%cy * current_state%sv%data(k, local_y, local_x)
132  end if
133 #endif
134  end do
Here is the caller graph for this function:

◆ finalisation_callback()

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

Frees up the allocated buffers (if such were allocated)

Parameters
current_stateThe current model state

Definition at line 71 of file pressuresource.F90.

72  type(model_state_type), target, intent(inout) :: current_state
73 
74  if (allocated(send_buffer_x)) deallocate(send_buffer_x)
75  if (allocated(current_state%psrce_recv_buffer_x)) deallocate(current_state%psrce_recv_buffer_x)
76  if (allocated(send_buffer_y)) deallocate(send_buffer_y)
77  if (allocated(current_state%psrce_recv_buffer_y)) deallocate(current_state%psrce_recv_buffer_y)
Here is the caller graph for this function:

◆ initialisation_callback()

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

On initialisation this will allocate the buffer areas required and set communication handles to null.

Parameters
current_stateThe current model state

Definition at line 38 of file pressuresource.F90.

39  type(model_state_type), target, intent(inout) :: current_state
40 
41  if (.not. is_component_enabled(current_state%options_database, "diverr")) then
42  call log_master_log(log_error, "The pressure source component requires the diverr component to be enabled")
43  end if
44 
45 #ifdef U_ACTIVE
46  allocate(send_buffer_x(current_state%local_grid%size(z_index)-1, current_state%local_grid%size(y_index)), &
47  current_state%psrce_recv_buffer_x(current_state%local_grid%size(z_index)-1, current_state%local_grid%size(y_index)))
48 #endif
49 #ifdef V_ACTIVE
50  allocate(send_buffer_y(current_state%local_grid%size(z_index)-1, current_state%local_grid%size(x_index)), &
51  current_state%psrce_recv_buffer_y(current_state%local_grid%size(z_index)-1, current_state%local_grid%size(x_index)))
52 #endif
53  current_state%psrce_x_hs_send_request=mpi_request_null
54  current_state%psrce_y_hs_send_request=mpi_request_null
55  current_state%psrce_x_hs_recv_request=mpi_request_null
56  current_state%psrce_y_hs_recv_request=mpi_request_null
Here is the caller graph for this function:

◆ pressuresource_get_descriptor()

type(component_descriptor_type) function, public pressuresource_mod::pressuresource_get_descriptor

Descriptor of this component for registration.

Returns
The pressure source component descriptor

Definition at line 28 of file pressuresource.F90.

29  pressuresource_get_descriptor%name="psrce"
30  pressuresource_get_descriptor%version=0.1
31  pressuresource_get_descriptor%initialisation=>initialisation_callback
32  pressuresource_get_descriptor%timestep=>timestep_callback
33  pressuresource_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ register_neighbouring_pressure_data_recv()

subroutine pressuresource_mod::register_neighbouring_pressure_data_recv ( type(model_state_type), intent(inout), target  current_state)
private

Registers the receive requests for each neighbouring process if that is not local, is recieves from p-1.

Parameters
current_stateThe current model state

Definition at line 164 of file pressuresource.F90.

165  type(model_state_type), target, intent(inout) :: current_state
166 
167  integer :: ierr
168 
169 #ifdef U_ACTIVE
170  if (current_state%local_grid%neighbours(x_index,2) .ne. current_state%parallel%my_rank) then
171  call mpi_irecv(current_state%psrce_recv_buffer_x, size(current_state%psrce_recv_buffer_x), precision_type, &
172  current_state%local_grid%neighbours(x_index,2), 10, current_state%parallel%neighbour_comm, &
173  current_state%psrce_x_hs_recv_request, ierr)
174  end if
175 #endif
176 #ifdef V_ACTIVE
177  if (current_state%local_grid%neighbours(y_index,2) .ne. current_state%parallel%my_rank) then
178  call mpi_irecv(current_state%psrce_recv_buffer_y, size(current_state%psrce_recv_buffer_y), precision_type, &
179  current_state%local_grid%neighbours(y_index,2), 10, current_state%parallel%neighbour_comm, &
180  current_state%psrce_y_hs_recv_request, ierr)
181  end if
182 #endif
Here is the caller graph for this function:

◆ send_neighbouring_pressure_data()

subroutine pressuresource_mod::send_neighbouring_pressure_data ( type(model_state_type), intent(inout), target  current_state)
private

Sends the computed source pressure data terms to the p+1 process.

Parameters
current_stateThe current model state

Definition at line 139 of file pressuresource.F90.

140  type(model_state_type), target, intent(inout) :: current_state
141 
142  integer :: ierr
143 
144 #ifdef U_ACTIVE
145  if (current_state%local_grid%neighbours(x_index,3) .eq. current_state%parallel%my_rank) then
146  current_state%psrce_recv_buffer_x=send_buffer_x
147  else
148  call mpi_isend(send_buffer_x, size(send_buffer_x), precision_type, current_state%local_grid%neighbours(x_index,3), &
149  10, current_state%parallel%neighbour_comm, current_state%psrce_x_hs_send_request, ierr)
150  end if
151 #endif
152 #ifdef V_ACTIVE
153  if (current_state%local_grid%neighbours(y_index,3) .eq. current_state%parallel%my_rank) then
154  current_state%psrce_recv_buffer_y=send_buffer_y
155  else
156  call mpi_isend(send_buffer_y, size(send_buffer_y), precision_type, current_state%local_grid%neighbours(y_index,3), &
157  10, current_state%parallel%neighbour_comm, current_state%psrce_y_hs_send_request, ierr)
158  end if
159 #endif
Here is the caller graph for this function:

◆ timestep_callback()

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

The timestep callback will update the values of P for each column.

Parameters
current_stateThe current model state

Definition at line 61 of file pressuresource.F90.

62  type(model_state_type), target, intent(inout) :: current_state
63 
64  if (current_state%first_timestep_column) call register_neighbouring_pressure_data_recv(current_state)
65  if (.not. current_state%halo_column) call calculate_psrce(current_state)
66  if (current_state%last_timestep_column) call send_neighbouring_pressure_data(current_state)
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ send_buffer_x

real(kind=default_precision), dimension(:,:), allocatable pressuresource_mod::send_buffer_x
private

Definition at line 21 of file pressuresource.F90.

21  real(kind=default_precision), dimension(:,:), allocatable :: send_buffer_x, send_buffer_y

◆ send_buffer_y

real(kind=default_precision), dimension(:,:), allocatable pressuresource_mod::send_buffer_y
private

Definition at line 21 of file pressuresource.F90.

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