MONC
Data Types | Functions/Subroutines | Variables
iterativesolver_mod Module Reference

This is the iterative pressure solver and uses a Jacobi preconditioned BiCGStab which we implement here. More...

Data Types

type  matrix_type
 A helper type to abstract the concrete details of the matrix. More...
 

Functions/Subroutines

type(component_descriptor_type) function, public iterativesolver_get_descriptor ()
 Descriptor of the iterative solver component used by the registry. More...
 
subroutine initialisation_callback (current_state)
 Initialisation callback hook which will set up the halo swapping state and allocate some data. More...
 
subroutine timestep_callback (current_state)
 Timestep callback, this ignores all but the last column where it calls the solver. More...
 
subroutine finalisation_callback (current_state)
 Called as MONC is shutting down and frees the halo swap state and deallocates local data. More...
 
subroutine bicgstab (current_state, A, x, b, i_strt, i_end, j_strt, j_end, k_end)
 Performs the BiCGStab KS method. More...
 
subroutine cg_solver (current_state, A, x, b, i_strt, i_end, j_strt, j_end, k_end)
 Performs the preconditioned conjugate gradient method. More...
 
subroutine precond (current_state, A, s, r, preits)
 Jacobi preconditioner. More...
 
subroutine calc_ax (current_state, A, x, Ax)
 Calculates A * x. More...
 
real(kind=default_precision) function inner_prod (current_state, x, y, i_strt, i_end, j_strt, j_end, k_end)
 Returns the global inner product of two vectors, ignoring the halo cells. More...
 
real(kind=default_precision) function, dimension(3) inner_prod_three_way (current_state, t, s, i_strt, i_end, j_strt, j_end, k_end)
 Returns the global inner product of a pair of vectors, ignoring the halo cells for three separate pairs. This call is for optimisation to bunch up the comms in a BiCGStab solver per iteration. More...
 
subroutine set_matrix_for_poisson (grid_configuration, A, z_size)
 Sets the values of the provided matrix to solve the poisson equation. More...
 
subroutine deduce_global_divmax (current_state)
 Determines the global divmax which is written into the current state. More...
 
subroutine copy_p_to_halo_buffer (current_state, neighbour_description, dim, source_index, pid_location, current_page, source_data)
 Copies the p field data to halo buffers for a specific process in a dimension and halo cell. More...
 
subroutine copy_calc_ax_to_halo_buffer (current_state, neighbour_description, dim, source_index, pid_location, current_page, source_data)
 Copies the source field data to halo buffers for a specific process in a dimension and halo cell - for the calc_Ax halo swaps. More...
 
subroutine copy_halo_buffer_to_p (current_state, neighbour_description, dim, target_index, neighbour_location, current_page, source_data)
 Copies the halo buffer to halo location for the p field. More...
 
subroutine copy_halo_buffer_to_calc_ax (current_state, neighbour_description, dim, target_index, neighbour_location, current_page, source_data)
 Copies the halo buffer to halo location for the source field as required in the calc_Ax procedure. More...
 
subroutine perform_local_data_copy_for_p (current_state, halo_depth, involve_corners, source_data)
 Does local data copying for P variable halo swap. More...
 
subroutine perform_local_data_copy_for_calc_ax (current_state, halo_depth, involve_corners, source_data)
 Does a local data copy for halo swapping cells with wrap around (to maintain periodic boundary condition) More...
 
type(matrix_type) function create_problem_matrix (z_size)
 Creates a problem matrix, allocating the required data based upon the column size. More...
 
subroutine complete_psrce_calculation (current_state, y_halo_size, x_halo_size)
 Completes the psrce calculation by waiting on all outstanding psrce communications to complete and then combine the received values with the P field for U and V. More...
 

Variables

real(kind=default_precision) tolerance
 
real(kind=default_precision) relaxation
 Solving tollerance. More...
 
integer max_iterations
 Maximum number of BiCGStab iterations. More...
 
integer preconditioner_iterations
 Number of preconditioner iterations to perform per call. More...
 
logical symm_prob
 
real(kind=default_precision), parameter tiny = 1.0e-16
 Minimum residual - if we go below this then something has gone wrong. More...
 
type(halo_communication_type), save halo_swap_state
 The halo swap state as initialised by that module. More...
 
real(kind=default_precision), dimension(:,:,:), allocatable psource
 
real(kind=default_precision), dimension(:,:,:), allocatable prev_p
 Passed to BiCGStab as the RHS. More...
 
logical first_run =.true.
 
type(matrix_typea
 

Detailed Description

This is the iterative pressure solver and uses a Jacobi preconditioned BiCGStab which we implement here.

Function/Subroutine Documentation

◆ bicgstab()

subroutine iterativesolver_mod::bicgstab ( type(model_state_type), intent(inout), target  current_state,
type(matrix_type), intent(inout)  A,
real(kind=default_precision), dimension(:,:,:), intent(inout)  x,
real(kind=default_precision), dimension(:,:,:), intent(in)  b,
integer, intent(in)  i_strt,
integer, intent(in)  i_end,
integer, intent(in)  j_strt,
integer, intent(in)  j_end,
integer, intent(in)  k_end 
)
private

Performs the BiCGStab KS method.

Parameters
current_stateThe current model state
AThe matrix
xThe solution
bThe RHS

Definition at line 136 of file iterativesolver.F90.

136  type(model_state_type), target, intent(inout) :: current_state
137  type(matrix_type), intent(inout) :: a
138  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(inout) :: x
139  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(in) :: b
140  integer, intent(in) :: i_strt, i_end, j_strt, j_end, k_end
141 
142  integer :: it, i, j, k
143  real(kind=DEFAULT_PRECISION) :: sc_err, alf, omg, nrm, my_rho, bet, tt, ts, ss, err, init_err, inner_prod_results(3)
144  real(kind=DEFAULT_PRECISION), dimension(current_state%local_grid%size(Z_INDEX) + & current_state%local_grid%halo_size(Z_INDEX) * 2, & current_state%local_grid%size(Y_INDEX) + current_state%local_grid%halo_size(Y_INDEX) * 2, & current_state%local_grid%size(X_INDEX) + current_state%local_grid%halo_size(X_INDEX) * 2) :: ax, r, cr, pp, v, t, s, cs
145 
146  ! Calculate scale factor for error
147  sc_err = sqrt(inner_prod(current_state, b, b, i_strt, i_end, j_strt, j_end, k_end))
148  sc_err = max(sc_err, 0.0001_default_precision)
149 
150  ! Calculate initial residual
151  call calc_ax(current_state, a, x, ax)
152 
153  do i = i_strt, i_end
154  do j = j_strt, j_end
155  do k = 2, k_end
156  r(k,j,i) = b(k,j,i) - ax(k,j,i)
157  cr(k,j,i) = r(k,j,i)
158  end do
159  end do
160  end do
161 
162  my_rho = inner_prod(current_state, r, r, i_strt, i_end, j_strt, j_end, k_end)
163  err = sqrt(my_rho)/sc_err
164  init_err = err
165 
166  alf = 1.0_default_precision
167  omg = 1.0_default_precision
168  nrm = 1.0_default_precision
169 
170  if (err .ge. tolerance) then
171  do it=1, max_iterations
172  if (it > 1) my_rho = inner_prod(current_state, r, cr, i_strt, i_end, j_strt, j_end, k_end)
173  bet = (my_rho/nrm) * (alf/omg)
174  if (it == 1) then
175  call precond(current_state, a, pp, r, preconditioner_iterations)
176  else
177  do i = i_strt, i_end
178  do j = j_strt, j_end
179  do k = 2, k_end
180  t(k,j,i) = r(k,j,i) - bet*omg*v(k,j,i)
181  end do
182  end do
183  end do
184  call precond(current_state, a, s, t, preconditioner_iterations)
185  do i = i_strt, i_end
186  do j = j_strt, j_end
187  do k = 2, k_end
188  pp(k,j,i) = s(k,j,i) + bet*pp(k,j,i)
189  end do
190  end do
191  end do
192  end if
193  call calc_ax(current_state, a, pp, v)
194  nrm = inner_prod(current_state, cr, v, i_strt, i_end, j_strt, j_end, k_end)
195  alf = my_rho / nrm
196 
197  do i = i_strt, i_end
198  do j = j_strt, j_end
199  do k = 2, k_end
200  s(k,j,i) = r(k,j,i) - alf*v(k,j,i)
201  end do
202  end do
203  end do
204 
205  call precond(current_state, a, cs, s, preconditioner_iterations)
206  call calc_ax(current_state, a, cs, t)
207 
208  inner_prod_results=inner_prod_three_way(current_state, t, s, i_strt, i_end, j_strt, j_end, k_end)
209  tt = inner_prod_results(1)
210  ts = inner_prod_results(2)
211  ss = inner_prod_results(3)
212  omg = ts/tt
213  x = x + alf*pp + omg*cs
214  do i = i_strt, i_end
215  do j = j_strt, j_end
216  do k = 2, k_end
217  r(k,j,i) = s(k,j,i) - omg*t(k,j,i)
218  end do
219  end do
220  end do
221  nrm = my_rho
222 
223  if (abs(omg) < tiny) then
224  call log_log(log_warn, "Convergence problem, omega="//conv_to_string(omg))
225  endif
226 
227  err = sqrt(ss - 2*omg*ts + omg**2 *tt)/sc_err
228  if (err < tolerance) exit
229  end do
230  end if
231 
232  if (err > tolerance) then
233  call log_log(log_warn, "Convergence failed, RNorm="//conv_to_string(err, exponent=.true.))
234  else if (current_state%parallel%my_rank==0 .and. log_get_logging_level() .eq. log_debug) then
235  call log_log(log_debug, "Converged in "//trim(conv_to_string(it))//" iterations with RNorm="//&
236  trim(conv_to_string(err, 5, .true.))//" initial norm="//trim(conv_to_string(init_err, 5, .true.)))
237  end if
238 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_ax()

subroutine iterativesolver_mod::calc_ax ( type(model_state_type), intent(inout), target  current_state,
type(matrix_type), intent(in)  A,
real(kind=default_precision), dimension(:,:,:), intent(inout), target  x,
real(kind=default_precision), dimension(:,:,:), intent(inout), target  Ax 
)
private

Calculates A * x.

Parameters
current_stateThe current model state
AThe matrix
xVector to multiply with
AxResult of A*x

Definition at line 379 of file iterativesolver.F90.

379  type(model_state_type), target, intent(inout) :: current_state
380  type(matrix_type), intent(in) :: a
381  real(kind=DEFAULT_PRECISION), dimension(:,:,:), target, intent(inout) :: x, ax
382 
383  integer :: i, k, j, n, istart, iend, jstart, jend
384  type(field_data_wrapper_type) :: source_data
385 
386  source_data%data=>x
387 
388  call initiate_nonblocking_halo_swap(current_state, halo_swap_state, &
389  copy_calc_ax_to_halo_buffer, source_data=(/source_data/))
390 
391  ax(1,:,:) = 0.0_default_precision
392  if (symm_prob) then
393  do n=1, 5
394  if (n==1) then
395  istart=current_state%local_grid%local_domain_start_index(x_index)+1
396  iend=current_state%local_grid%local_domain_end_index(x_index)-1
397  jstart=current_state%local_grid%local_domain_start_index(y_index)+1
398  jend=current_state%local_grid%local_domain_end_index(y_index)-1
399  else if (n==2) then
400  istart=current_state%local_grid%local_domain_start_index(x_index)
401  iend=current_state%local_grid%local_domain_start_index(x_index)
402  else if (n==3) then
403  istart=current_state%local_grid%local_domain_end_index(x_index)
404  iend=current_state%local_grid%local_domain_end_index(x_index)
405  else if (n==4) then
406  jstart=current_state%local_grid%local_domain_start_index(y_index)
407  jend=current_state%local_grid%local_domain_start_index(y_index)
408  istart=current_state%local_grid%local_domain_start_index(x_index)
409  iend=current_state%local_grid%local_domain_end_index(x_index)
410  else if (n==5) then
411  jstart=current_state%local_grid%local_domain_end_index(y_index)
412  jend=current_state%local_grid%local_domain_end_index(y_index)
413  end if
414  do i=istart, iend
415  do j=jstart, jend
416  k=2
417  ax(k,j,i)=a%vol(k)*(a%n*(x(k,j,i+1)+x(k,j,i-1))+a%e*(x(k,j+1,i)+x(k,j-1,i)))+ a%u(k)*x(k+1,j,i)+a%p(k)*x(k,j,i)
418  do k=3,current_state%local_grid%size(z_index)-1
419  ax(k,j,i)=a%vol(k)*(a%n*(x(k,j,i+1)+x(k,j,i-1))+a%e*(x(k,j+1,i)+x(k,j-1,i)))+&
420  a%u(k)*x(k+1,j,i)+a%d(k)*x(k-1,j,i)+a%p(k)*x(k,j,i)
421  end do
422  k=current_state%local_grid%size(z_index)
423  ax(k,j,i) = a%vol(k)*(a%n*(x(k,j,i+1)+x(k,j,i-1))+a%e*(x(k,j+1,i)+x(k,j-1,i)))+ a%d(k)*x(k-1,j,i)+a%p(k)*x(k,j,i)
424  end do
425  end do
426  if (n==1) then
427  call complete_nonblocking_halo_swap(current_state, halo_swap_state, perform_local_data_copy_for_calc_ax, &
428  copy_halo_buffer_to_calc_ax, source_data=(/source_data/))
429  end if
430  end do
431  else
432  do n=1, 5
433  if (n==1) then
434  istart=current_state%local_grid%local_domain_start_index(x_index)+1
435  iend=current_state%local_grid%local_domain_end_index(x_index)-1
436  jstart=current_state%local_grid%local_domain_start_index(y_index)+1
437  jend=current_state%local_grid%local_domain_end_index(y_index)-1
438  else if (n==2) then
439  istart=current_state%local_grid%local_domain_start_index(x_index)
440  iend=current_state%local_grid%local_domain_start_index(x_index)
441  else if (n==3) then
442  istart=current_state%local_grid%local_domain_end_index(x_index)
443  iend=current_state%local_grid%local_domain_end_index(x_index)
444  else if (n==4) then
445  jstart=current_state%local_grid%local_domain_start_index(y_index)
446  jend=current_state%local_grid%local_domain_start_index(y_index)
447  istart=current_state%local_grid%local_domain_start_index(x_index)
448  iend=current_state%local_grid%local_domain_end_index(x_index)
449  else if (n==5) then
450  jstart=current_state%local_grid%local_domain_end_index(y_index)
451  jend=current_state%local_grid%local_domain_end_index(y_index)
452  end if
453  do i=istart, iend
454  do j=jstart, jend
455  k=2
456  ax(k,j,i)=a%n*(x(k,j,i+1)+x(k,j,i-1))+a%e*(x(k,j+1,i)+x(k,j-1,i))+ a%u(k)*x(k+1,j,i)+a%p(k)*x(k,j,i)
457  do k=3,current_state%local_grid%size(z_index)-1
458  ax(k,j,i)=a%n*(x(k,j,i+1)+x(k,j,i-1))+a%e*(x(k,j+1,i)+x(k,j-1,i))+&
459  a%u(k)*x(k+1,j,i)+a%d(k)*x(k-1,j,i)+a%p(k)*x(k,j,i)
460  end do
461  k=current_state%local_grid%size(z_index)
462  ax(k,j,i) = a%n*(x(k,j,i+1)+x(k,j,i-1))+a%e*(x(k,j+1,i)+x(k,j-1,i))+ a%d(k)*x(k-1,j,i)+a%p(k)*x(k,j,i)
463  end do
464  end do
465  if (n==1) then
466  call complete_nonblocking_halo_swap(current_state, halo_swap_state, perform_local_data_copy_for_calc_ax, &
467  copy_halo_buffer_to_calc_ax, source_data=(/source_data/))
468  end if
469  end do
470  endif
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cg_solver()

subroutine iterativesolver_mod::cg_solver ( type(model_state_type), intent(inout), target  current_state,
type(matrix_type), intent(inout)  A,
real(kind=default_precision), dimension(:,:,:), intent(inout)  x,
real(kind=default_precision), dimension(:,:,:), intent(in)  b,
integer, intent(in)  i_strt,
integer, intent(in)  i_end,
integer, intent(in)  j_strt,
integer, intent(in)  j_end,
integer, intent(in)  k_end 
)
private

Performs the preconditioned conjugate gradient method.

Parameters
current_stateThe current model state
AThe matrix
xThe solution
bThe RHS

Definition at line 249 of file iterativesolver.F90.

249  type(model_state_type), target, intent(inout) :: current_state
250  type(matrix_type), intent(inout) :: a
251  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(inout) :: x
252  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(in) :: b
253  integer, intent(in) :: i_strt, i_end, j_strt, j_end, k_end
254 
255  integer :: it, k, i, j
256  real(kind=DEFAULT_PRECISION) :: sc_err, alf, bet, err, init_err, rho
257  real(kind=DEFAULT_PRECISION), dimension(current_state%local_grid%size(Z_INDEX) + & current_state%local_grid%halo_size(Z_INDEX) * 2, & current_state%local_grid%size(Y_INDEX) + current_state%local_grid%halo_size(Y_INDEX) * 2, & current_state%local_grid%size(X_INDEX) + current_state%local_grid%halo_size(X_INDEX) * 2) :: ax, r, z, p
258 
259  ! first rescale RHS for symmetry (this could be done when p_source is calculated
260  do i=current_state%local_grid%local_domain_start_index(x_index), current_state%local_grid%local_domain_end_index(x_index)
261  do j=current_state%local_grid%local_domain_start_index(y_index), current_state%local_grid%local_domain_end_index(y_index)
262 
263  r(1,j,i) = 0.0_default_precision
264  do k=2,current_state%local_grid%size(z_index)
265  r(k,j,i) = b(k,j,i) * a%vol(k)
266  end do
267  end do
268  end do
269 
270  ! Calculate scale factor for error
271 
272  call calc_ax(current_state, a, x, ax)
273 
274  sc_err = sqrt(inner_prod(current_state, r, r, i_strt, i_end, j_strt, j_end, k_end))
275  sc_err = max(sc_err, 0.0001_default_precision)
276  r = r - ax
277  init_err = sqrt(inner_prod(current_state, r, r, i_strt, i_end, j_strt, j_end, k_end))/sc_err
278 
279  do it=1, max_iterations
280  if( it == 1 ) then
281  call precond(current_state, a, p, r, preconditioner_iterations)
282  rho = inner_prod(current_state, p, r, i_strt, i_end, j_strt, j_end, k_end)
283  alf = rho
284  else
285  call precond(current_state, a, z, r, preconditioner_iterations)
286  alf = inner_prod(current_state, z, r, i_strt, i_end, j_strt, j_end, k_end)
287  bet = alf/rho
288  rho = alf
289  p = z + bet*p
290  end if
291 
292  call calc_ax(current_state, a, p, ax)
293  alf = alf/inner_prod(current_state, p, ax, i_strt, i_end, j_strt, j_end, k_end)
294  x = x + alf*p
295  r = r - alf*ax
296 
297  err = sqrt(inner_prod(current_state, r, r, i_strt, i_end, j_strt, j_end, k_end))/sc_err
298  if (err < tolerance) exit
299  end do
300 
301  if( current_state%parallel%my_rank == 0 ) print*,it, err, init_err
302 
303  if (err > tolerance) then
304  call log_log(log_warn, "Convergence failed, RNorm="//conv_to_string(err, exponent=.true.))
305  else if (current_state%parallel%my_rank==0 .and. log_get_logging_level() .eq. log_debug) then
306  call log_log(log_debug, "Converged in "//trim(conv_to_string(it))//" iterations with RNorm="//&
307  trim(conv_to_string(err, 5, .true.))//" initial norm="//trim(conv_to_string(init_err, 5, .true.)))
308  end if
309 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ complete_psrce_calculation()

subroutine iterativesolver_mod::complete_psrce_calculation ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  y_halo_size,
integer, intent(in)  x_halo_size 
)
private

Completes the psrce calculation by waiting on all outstanding psrce communications to complete and then combine the received values with the P field for U and V.

Parameters
current_stateThe current model state
y_halo_sizeThe halo size in the Y dimension
x_halo_sizeThe halo size in the X dimension

Definition at line 736 of file iterativesolver.F90.

736  type(model_state_type), target, intent(inout) :: current_state
737  integer, intent(in) :: y_halo_size, x_halo_size
738 
739  integer :: ierr, combined_handles(2), i, j, k
740 
741  combined_handles(1)=current_state%psrce_x_hs_recv_request
742  combined_handles(2)=current_state%psrce_y_hs_recv_request
743  call mpi_waitall(2, combined_handles, mpi_statuses_ignore, ierr)
744 
745  do j=current_state%local_grid%local_domain_start_index(y_index), current_state%local_grid%local_domain_end_index(y_index)
746  do k=2,current_state%local_grid%size(z_index)
747 #ifdef U_ACTIVE
748  current_state%p%data(k,j,x_halo_size+1)=current_state%p%data(k,j,x_halo_size+1)-&
749  current_state%psrce_recv_buffer_x(k-1,j-x_halo_size)
750 #endif
751 #ifdef V_ACTIVE
752  if (j .gt. y_halo_size+1) current_state%p%data(k, j, x_halo_size+1)=current_state%p%data(k, j, x_halo_size+1)-&
753  current_state%global_grid%configuration%horizontal%cy * current_state%sv%data(k, j-1, x_halo_size+1)
754 #endif
755  end do
756  end do
757 
758 #ifdef V_ACTIVE
759  do i=current_state%local_grid%local_domain_start_index(x_index), current_state%local_grid%local_domain_end_index(x_index)
760  do k=2,current_state%local_grid%size(z_index)
761  current_state%p%data(k,y_halo_size+1,i)=current_state%p%data(k,y_halo_size+1,i)-&
762  current_state%psrce_recv_buffer_y(k-1,i-y_halo_size)
763  end do
764  end do
765 #endif
766 
767  combined_handles(1)=current_state%psrce_x_hs_send_request
768  combined_handles(2)=current_state%psrce_y_hs_send_request
769  call mpi_waitall(2, combined_handles, mpi_statuses_ignore, ierr)
Here is the caller graph for this function:

◆ copy_calc_ax_to_halo_buffer()

subroutine iterativesolver_mod::copy_calc_ax_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the source field data to halo buffers for a specific process in a dimension and halo cell - for the calc_Ax halo swaps.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
dimDimension to copy from
source_indexThe source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from

Definition at line 625 of file iterativesolver.F90.

625  type(model_state_type), intent(inout) :: current_state
626  integer, intent(in) :: dim, pid_location, source_index
627  integer, intent(inout) :: current_page(:)
628  type(neighbour_description_type), intent(inout) :: neighbour_description
629  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
630 
631  type(field_data_wrapper_type) :: selected_source
632 
633  selected_source=source_data(1)
634 
635  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, selected_source%data, &
636  dim, source_index, current_page(pid_location))
637 
638  current_page(pid_location)=current_page(pid_location)+1
Here is the caller graph for this function:

◆ copy_halo_buffer_to_calc_ax()

subroutine iterativesolver_mod::copy_halo_buffer_to_calc_ax ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  target_index,
integer, intent(in)  neighbour_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the halo buffer to halo location for the source field as required in the calc_Ax procedure.

Parameters
current_stateThe current model state
neighbour_descriptionThe halo swapping description of the neighbour we are accessing the buffer of
dimThe dimension we receive for
target_indexThe target index for the dimension we are receiving for
neighbour_locationThe location in the local neighbour data stores of this neighbour
current_pageThe current, next, halo swap page to read from (all previous have been read and copied already)
source_dataOptional source data which is written into

Definition at line 673 of file iterativesolver.F90.

673  type(model_state_type), intent(inout) :: current_state
674  integer, intent(in) :: dim, target_index, neighbour_location
675  integer, intent(inout) :: current_page(:)
676  type(neighbour_description_type), intent(inout) :: neighbour_description
677  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
678 
679  type(field_data_wrapper_type) :: selected_source
680 
681  selected_source=source_data(1)
682 
683  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, selected_source%data, &
684  dim, target_index, current_page(neighbour_location))
685 
686  current_page(neighbour_location)=current_page(neighbour_location)+1
Here is the caller graph for this function:

◆ copy_halo_buffer_to_p()

subroutine iterativesolver_mod::copy_halo_buffer_to_p ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  target_index,
integer, intent(in)  neighbour_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the halo buffer to halo location for the p field.

Parameters
current_stateThe current model state
neighbour_descriptionThe halo swapping description of the neighbour we are accessing the buffer of
dimThe dimension we receive for
target_indexThe target index for the dimension we are receiving for
neighbour_locationThe location in the local neighbour data stores of this neighbour
current_pageThe current, next, halo swap page to read from (all previous have been read and copied already)
source_dataOptional source data which is written into

Definition at line 651 of file iterativesolver.F90.

651  type(model_state_type), intent(inout) :: current_state
652  integer, intent(in) :: dim, target_index, neighbour_location
653  integer, intent(inout) :: current_page(:)
654  type(neighbour_description_type), intent(inout) :: neighbour_description
655  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
656 
657  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, current_state%p%data, &
658  dim, target_index, current_page(neighbour_location))
659 
660  current_page(neighbour_location)=current_page(neighbour_location)+1
Here is the caller graph for this function:

◆ copy_p_to_halo_buffer()

subroutine iterativesolver_mod::copy_p_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the p field data to halo buffers for a specific process in a dimension and halo cell.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
dimDimension to copy from
source_indexThe source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from

Definition at line 603 of file iterativesolver.F90.

603  type(model_state_type), intent(inout) :: current_state
604  integer, intent(in) :: dim, pid_location, source_index
605  integer, intent(inout) :: current_page(:)
606  type(neighbour_description_type), intent(inout) :: neighbour_description
607  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
608 
609  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, current_state%p%data, &
610  dim, source_index, current_page(pid_location))
611 
612  current_page(pid_location)=current_page(pid_location)+1
Here is the caller graph for this function:

◆ create_problem_matrix()

type(matrix_type) function iterativesolver_mod::create_problem_matrix ( integer, intent(in)  z_size)
private

Creates a problem matrix, allocating the required data based upon the column size.

Parameters
z_sizeNumber of elements in the vertical column
Returns
The allocated matrix ready to be used

Definition at line 723 of file iterativesolver.F90.

723  integer, intent(in) :: z_size
724  type(matrix_type) :: create_problem_matrix
725 
726  allocate(create_problem_matrix%u(z_size), create_problem_matrix%d(z_size), create_problem_matrix%p(z_size), &
727  create_problem_matrix%lu_u(z_size), create_problem_matrix%lu_d(z_size), create_problem_matrix%vol(z_size))
Here is the caller graph for this function:

◆ deduce_global_divmax()

subroutine iterativesolver_mod::deduce_global_divmax ( type(model_state_type), intent(inout), target  current_state)
private

Determines the global divmax which is written into the current state.

Parameters
current_stateThe current model state

Definition at line 585 of file iterativesolver.F90.

585  type(model_state_type), target, intent(inout) :: current_state
586 
587  integer :: ierr
588 
589  call mpi_allreduce(current_state%local_divmax, current_state%global_divmax, 1, precision_type, mpi_max, &
590  current_state%parallel%monc_communicator, ierr)
Here is the caller graph for this function:

◆ finalisation_callback()

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

Called as MONC is shutting down and frees the halo swap state and deallocates local data.

Parameters
current_stateThe current model state

Definition at line 124 of file iterativesolver.F90.

124  type(model_state_type), target, intent(inout) :: current_state
125 
126  call finalise_halo_communication(halo_swap_state)
127  deallocate(psource, prev_p, a%u, a%d, a%p, a%lu_u, a%lu_d, a%vol)
Here is the caller graph for this function:

◆ initialisation_callback()

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

Initialisation callback hook which will set up the halo swapping state and allocate some data.

Parameters
current_stateThe current model state

Definition at line 58 of file iterativesolver.F90.

58  type(model_state_type), target, intent(inout) :: current_state
59 
60  if (.not. is_component_enabled(current_state%options_database, "diverr")) then
61  call log_master_log(log_error, "The iterative solver component requires the diverr component to be enabled")
62  end if
63 
64  tolerance=options_get_real(current_state%options_database, "tolerance")
65  max_iterations=options_get_integer(current_state%options_database, "max_iterations")
66  preconditioner_iterations=options_get_integer(current_state%options_database, "preconditioner_iterations")
67  symm_prob=options_get_logical(current_state%options_database, "symm_prob")
68 
69  call init_halo_communication(current_state, get_single_field_per_halo_cell, halo_swap_state, 1, .false.)
70 
71  allocate(psource(current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2, &
72  current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2, &
73  current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2),&
74  prev_p(current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2, &
75  current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2, &
76  current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2))
77 
78  a=create_problem_matrix(current_state%local_grid%size(z_index))
79  call set_matrix_for_poisson(current_state%global_grid%configuration, a, current_state%local_grid%size(z_index))
Here is the call graph for this function:
Here is the caller graph for this function:

◆ inner_prod()

real(kind=default_precision) function iterativesolver_mod::inner_prod ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), dimension(:,:,:), intent(in)  x,
real(kind=default_precision), dimension(:,:,:), intent(in)  y,
integer, intent(in)  i_strt,
integer, intent(in)  i_end,
integer, intent(in)  j_strt,
integer, intent(in)  j_end,
integer, intent(in)  k_end 
)
private

Returns the global inner product of two vectors, ignoring the halo cells.

Parameters
current_stateThe current model state
xFirst vector y Second vector
Returns
Global inner product of the two input vectors

Definition at line 479 of file iterativesolver.F90.

479  type(model_state_type), target, intent(inout) :: current_state
480  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(in) :: x, y
481  integer, intent(in) :: i_strt, i_end, j_strt, j_end, k_end
482 
483  real(kind=DEFAULT_PRECISION) :: local_sum, global_sum
484  integer :: ierr, i, j, k
485 
486  local_sum=0.0_default_precision
487 
488  do i=i_strt, i_end
489  do j=j_strt, j_end
490  do k=2, k_end
491  local_sum=local_sum+x(k,j,i)*y(k,j,i)
492  end do
493  end do
494  end do
495 
496  call mpi_allreduce(local_sum, global_sum, 1, precision_type, mpi_sum, current_state%parallel%monc_communicator, ierr)
497  inner_prod=global_sum
Here is the caller graph for this function:

◆ inner_prod_three_way()

real(kind=default_precision) function, dimension(3) iterativesolver_mod::inner_prod_three_way ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), dimension(:,:,:), intent(in)  t,
real(kind=default_precision), dimension(:,:,:), intent(in)  s,
integer, intent(in)  i_strt,
integer, intent(in)  i_end,
integer, intent(in)  j_strt,
integer, intent(in)  j_end,
integer, intent(in)  k_end 
)
private

Returns the global inner product of a pair of vectors, ignoring the halo cells for three separate pairs. This call is for optimisation to bunch up the comms in a BiCGStab solver per iteration.

Parameters
current_stateThe current model state
xFirst vector y Second vector
Returns
Global inner product of the two input vectors

Definition at line 507 of file iterativesolver.F90.

507  type(model_state_type), target, intent(inout) :: current_state
508  integer, intent(in) :: i_strt, i_end, j_strt, j_end, k_end
509  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(in) :: t, s
510  real(kind=DEFAULT_PRECISION), dimension(3) :: inner_prod_three_way
511 
512  real(kind=DEFAULT_PRECISION), dimension(3) :: local_sum, global_sum
513  integer :: ierr, i, j, k
514 
515  local_sum(1)=0.0_default_precision
516  local_sum(2)=0.0_default_precision
517  local_sum(3)=0.0_default_precision
518 
519  do i=i_strt, i_end
520  do j=j_strt, j_end
521  do k=2, k_end
522  local_sum(1)=local_sum(1)+t(k,j,i)*t(k,j,i)
523  local_sum(2)=local_sum(2)+t(k,j,i)*s(k,j,i)
524  local_sum(3)=local_sum(3)+s(k,j,i)*s(k,j,i)
525  end do
526  end do
527  end do
528 
529  call mpi_allreduce(local_sum, global_sum, 3, precision_type, mpi_sum, current_state%parallel%monc_communicator, ierr)
530  inner_prod_three_way=global_sum
Here is the caller graph for this function:

◆ iterativesolver_get_descriptor()

type(component_descriptor_type) function, public iterativesolver_mod::iterativesolver_get_descriptor ( )

Descriptor of the iterative solver component used by the registry.

Returns
The iterative solver component descriptor

Definition at line 48 of file iterativesolver.F90.

48  iterativesolver_get_descriptor%name="iterativesolver"
49  iterativesolver_get_descriptor%version=0.1
50  iterativesolver_get_descriptor%initialisation=>initialisation_callback
51  iterativesolver_get_descriptor%timestep=>timestep_callback
52  iterativesolver_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ perform_local_data_copy_for_calc_ax()

subroutine iterativesolver_mod::perform_local_data_copy_for_calc_ax ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  halo_depth,
logical, intent(in)  involve_corners,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Does a local data copy for halo swapping cells with wrap around (to maintain periodic boundary condition)

Parameters
current_stateThe current model state
source_dataOptional source data which is written into

Definition at line 706 of file iterativesolver.F90.

706  type(model_state_type), intent(inout) :: current_state
707  integer, intent(in) :: halo_depth
708  logical, intent(in) :: involve_corners
709  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
710 
711  type(field_data_wrapper_type) :: selected_source
712 
713  selected_source=source_data(1)
714 
715  call perform_local_data_copy_for_field(selected_source%data, current_state%local_grid, &
716  current_state%parallel%my_rank, halo_depth, involve_corners)
Here is the caller graph for this function:

◆ perform_local_data_copy_for_p()

subroutine iterativesolver_mod::perform_local_data_copy_for_p ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  halo_depth,
logical, intent(in)  involve_corners,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Does local data copying for P variable halo swap.

Parameters
current_stateThe current model state_mod
source_dataOptional source data which is written into

Definition at line 693 of file iterativesolver.F90.

693  type(model_state_type), intent(inout) :: current_state
694  integer, intent(in) :: halo_depth
695  logical, intent(in) :: involve_corners
696  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
697 
698  call perform_local_data_copy_for_field(current_state%p%data, current_state%local_grid, &
699  current_state%parallel%my_rank, halo_depth, involve_corners)
Here is the caller graph for this function:

◆ precond()

subroutine iterativesolver_mod::precond ( type(model_state_type), intent(inout), target  current_state,
type(matrix_type), intent(inout)  A,
real(kind=default_precision), dimension(:,:,:), intent(inout)  s,
real(kind=default_precision), dimension(:,:,:), intent(in)  r,
integer, intent(in)  preits 
)
private

Jacobi preconditioner.

Parameters
current_stateThe current model state
AThe matrix
sWritten into as result of preconditioning
rInput values to preconditioner
preitsNumber of iterations of the preconditioner to perform per call

Definition at line 321 of file iterativesolver.F90.

321  type(model_state_type), target, intent(inout) :: current_state
322  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(in) :: r
323  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(inout) :: s
324  integer, intent(in) :: preits
325  type(matrix_type), intent(inout) :: a
326 
327  real(kind=DEFAULT_PRECISION), dimension(current_state%local_grid%size(Z_INDEX) + & current_state%local_grid%halo_size(Z_INDEX) * 2, current_state%local_grid%size(Y_INDEX) + & current_state%local_grid%halo_size(Y_INDEX) * 2, current_state%local_grid%size(X_INDEX) + & current_state%local_grid%halo_size(X_INDEX) * 2) :: t
328  real(kind=DEFAULT_PRECISION), dimension(current_state%local_grid%size(Z_INDEX)) :: s_k
329  integer :: it, i, j, k
330 
331  if (preits .lt. 0) then
332  s=r
333  return
334  end if
335 
336  do i=current_state%local_grid%local_domain_start_index(x_index), current_state%local_grid%local_domain_end_index(x_index)
337  do j=current_state%local_grid%local_domain_start_index(y_index), current_state%local_grid%local_domain_end_index(y_index)
338  s(1,j,i) = 0.0_default_precision
339  k=2
340  s(k,j,i)=r(k,j,i)*a%LU_d(k)
341  do k=3,current_state%local_grid%size(z_index)
342  s(k,j,i)=(r(k,j,i) - a%d(k)*s(k-1,j,i))*a%lu_d(k)
343  end do
344  do k=current_state%local_grid%size(z_index)-1, 2, -1
345  s(k,j,i)=s(k,j,i) - a%lu_u(k)*s(k+1,j,i)
346  end do
347  end do
348  end do
349 
350  do it=1, preits
351  call calc_ax(current_state, a, s, t)
352  do i=current_state%local_grid%local_domain_start_index(x_index), current_state%local_grid%local_domain_end_index(x_index)
353  do j=current_state%local_grid%local_domain_start_index(y_index), current_state%local_grid%local_domain_end_index(y_index)
354  k=2
355  s_k(k)=(r(k,j,i) - t(k,j,i))*a%lu_d(k)
356  do k=3,current_state%local_grid%size(z_index)
357  s_k(k)=(r(k,j,i) - t(k,j,i) - a%d(k)*s_k(k-1))*a%lu_d(k)
358  end do
359  k=current_state%local_grid%size(z_index)
360  s(k,j,i)=s(k,j,i)+s_k(k)
361  do k=current_state%local_grid%size(z_index)-1, 2, -1
362  s_k(k)=s_k(k) - a%lu_u(k)*s_k(k+1)
363  s(k,j,i)=s(k,j,i) + relaxation*s_k(k)
364  end do
365  end do
366  end do
367  end do
368 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_matrix_for_poisson()

subroutine iterativesolver_mod::set_matrix_for_poisson ( type(grid_configuration_type), intent(inout)  grid_configuration,
type(matrix_type), intent(inout)  A,
integer, intent(in)  z_size 
)
private

Sets the values of the provided matrix to solve the poisson equation.

Parameters
grid_configurationConfiguration of the vertical and horizontal grids
AThe matrix that the values are written into
z_sizeNumber of elements in a column

Definition at line 538 of file iterativesolver.F90.

538  type(grid_configuration_type), intent(inout) :: grid_configuration
539  type(matrix_type), intent(inout) :: a
540  integer, intent(in) :: z_size
541 
542  integer :: k
543  real(kind=DEFAULT_PRECISION) :: d_sc, concat_scalars
544 
545  a%n=grid_configuration%horizontal%cx*grid_configuration%horizontal%cx
546  a%s=a%n
547  a%e=grid_configuration%horizontal%cy*grid_configuration%horizontal%cy
548  a%w=a%e
549  concat_scalars=a%n+a%s+a%e+a%w
550  do k=2, z_size
551  if (symm_prob) then
552  a%vol(k)=grid_configuration%vertical%dz(k)
553  d_sc=1.0/grid_configuration%vertical%rhon(k)
554  else
555  d_sc=grid_configuration%vertical%rdz(k) / grid_configuration%vertical%rhon(k)
556  a%vol(k)=1.0
557  endif
558 
559  if (k==z_size) then
560  a%u(k)=0.0_default_precision
561  else
562  a%u(k)=grid_configuration%vertical%rho(k)*grid_configuration%vertical%rdzn(k+1)
563  end if
564  if (k==2) then
565  a%d(k)=0.0_default_precision
566  else
567  a%d(k)=grid_configuration%vertical%rho(k-1)*grid_configuration%vertical%rdzn(k)
568  end if
569  a%p(k) = d_sc * (-(a%u(k) + a%d(k))) - concat_scalars * a%vol(k)
570  a%u(k)=d_sc * a%u(k)
571  a%d(k)=d_sc * a%d(k)
572  end do
573  k=2
574  a%lu_d(k)=1.0_default_precision/a%p(k)
575  a%lu_u(k)=a%lu_d(k)*a%u(k)
576  do k=3, z_size
577  a%lu_d(k)=1.0_default_precision/(a%p(k) - a%d(k)*a%lu_u(k-1))
578  a%lu_u(k)=a%u(k)*a%lu_d(k)
579  end do
Here is the caller graph for this function:

◆ timestep_callback()

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

Timestep callback, this ignores all but the last column where it calls the solver.

Parameters
current_stateThe current model state

Definition at line 85 of file iterativesolver.F90.

85  type(model_state_type), target, intent(inout) :: current_state
86 
87  integer :: i_strt, i_end, j_strt, j_end, k_end
88 
89  i_strt = current_state%local_grid%local_domain_start_index(x_index)
90  i_end = current_state%local_grid%local_domain_end_index(x_index)
91  j_strt = current_state%local_grid%local_domain_start_index(y_index)
92  j_end = current_state%local_grid%local_domain_end_index(y_index)
93  k_end = current_state%local_grid%size(z_index)
94 
95  call complete_psrce_calculation(current_state, current_state%local_grid%halo_size(y_index), &
96  current_state%local_grid%halo_size(x_index))
97 
98  call initiate_nonblocking_halo_swap(current_state, halo_swap_state, copy_p_to_halo_buffer)
99  call deduce_global_divmax(current_state)
100  call complete_nonblocking_halo_swap(current_state, halo_swap_state, perform_local_data_copy_for_p, copy_halo_buffer_to_p)
101 
102  psource=current_state%p%data
103  if (first_run) then
104  ! If first timestep then initial guess is zero
105  current_state%p%data=0.0_default_precision
106  first_run=.false.
107  else
108  ! Initial guess is set to previous timesteps p
109  current_state%p%data=prev_p
110  end if
111 
112  if (symm_prob) then
113  call cg_solver(current_state, a, current_state%p%data, psource, i_strt, i_end, j_strt, j_end, k_end)
114  else
115  call bicgstab(current_state, a, current_state%p%data, psource, i_strt, i_end, j_strt, j_end, k_end)
116  end if
117 
118  prev_p=current_state%p%data
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ a

type(matrix_type) iterativesolver_mod::a
private

Definition at line 40 of file iterativesolver.F90.

40  type(matrix_type) :: a

◆ first_run

logical iterativesolver_mod::first_run =.true.
private

Definition at line 39 of file iterativesolver.F90.

39  logical :: first_run=.true.

◆ halo_swap_state

type(halo_communication_type), save iterativesolver_mod::halo_swap_state
private

The halo swap state as initialised by that module.

Definition at line 37 of file iterativesolver.F90.

37  type(halo_communication_type), save :: halo_swap_state

◆ max_iterations

integer iterativesolver_mod::max_iterations
private

Maximum number of BiCGStab iterations.

Definition at line 31 of file iterativesolver.F90.

31  integer :: max_iterations, & !< Maximum number of BiCGStab iterations
32  preconditioner_iterations

◆ preconditioner_iterations

integer iterativesolver_mod::preconditioner_iterations
private

Number of preconditioner iterations to perform per call.

Definition at line 31 of file iterativesolver.F90.

◆ prev_p

real(kind=default_precision), dimension(:,:,:), allocatable iterativesolver_mod::prev_p
private

Passed to BiCGStab as the RHS.

Definition at line 38 of file iterativesolver.F90.

◆ psource

real(kind=default_precision), dimension(:,:,:), allocatable iterativesolver_mod::psource
private

Definition at line 38 of file iterativesolver.F90.

38  real(kind=DEFAULT_PRECISION), dimension(:,:,:), allocatable :: psource, prev_p

◆ relaxation

real(kind=default_precision) iterativesolver_mod::relaxation
private

Solving tollerance.

Definition at line 30 of file iterativesolver.F90.

◆ symm_prob

logical iterativesolver_mod::symm_prob
private

Definition at line 33 of file iterativesolver.F90.

33  logical :: symm_prob

◆ tiny

real(kind=default_precision), parameter iterativesolver_mod::tiny = 1.0e-16
private

Minimum residual - if we go below this then something has gone wrong.

Definition at line 35 of file iterativesolver.F90.

35  real(kind=DEFAULT_PRECISION), parameter :: tiny = 1.0e-16

◆ tolerance

real(kind=default_precision) iterativesolver_mod::tolerance
private

Definition at line 30 of file iterativesolver.F90.

30  real(kind=DEFAULT_PRECISION) :: tolerance, relaxation