Neko  0.8.1
A portable framework for high-order spectral element flow simulations
pipecg_sx.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module pipecg_sx
36  use precon, only : pc_t
37  use ax_product, only : ax_t
38  use num_types, only: rp
39  use field, only : field_t
40  use coefs, only : coef_t
41  use gather_scatter, only : gs_t, gs_op_add
42  use bc, only : bc_list_t, bc_list_apply
43  use math, only : glsc3, abscmp
44  use comm
45  implicit none
46  private
47 
49  type, public, extends(ksp_t) :: sx_pipecg_t
50  real(kind=rp), allocatable :: p(:)
51  real(kind=rp), allocatable :: q(:)
52  real(kind=rp), allocatable :: r(:)
53  real(kind=rp), allocatable :: s(:)
54  real(kind=rp), allocatable :: u(:)
55  real(kind=rp), allocatable :: w(:)
56  real(kind=rp), allocatable :: z(:)
57  real(kind=rp), allocatable :: mi(:)
58  real(kind=rp), allocatable :: ni(:)
59  contains
60  procedure, pass(this) :: init => sx_pipecg_init
61  procedure, pass(this) :: free => sx_pipecg_free
62  procedure, pass(this) :: solve => sx_pipecg_solve
63  end type sx_pipecg_t
64 
65 contains
66 
68  subroutine sx_pipecg_init(this, n, max_iter, M, rel_tol, abs_tol)
69  class(sx_pipecg_t), intent(inout) :: this
70  class(pc_t), optional, intent(inout), target :: M
71  integer, intent(in) :: n
72  integer, intent(in) :: max_iter
73  real(kind=rp), optional, intent(inout) :: rel_tol
74  real(kind=rp), optional, intent(inout) :: abs_tol
75 
76  call this%free()
77 
78  allocate(this%p(n))
79  allocate(this%q(n))
80  allocate(this%r(n))
81  allocate(this%s(n))
82  allocate(this%u(n))
83  allocate(this%w(n))
84  allocate(this%z(n))
85  allocate(this%mi(n))
86  allocate(this%ni(n))
87  if (present(m)) then
88  this%M => m
89  end if
90 
91  if (present(rel_tol) .and. present(abs_tol)) then
92  call this%ksp_init(max_iter, rel_tol, abs_tol)
93  else if (present(rel_tol)) then
94  call this%ksp_init(max_iter, rel_tol=rel_tol)
95  else if (present(abs_tol)) then
96  call this%ksp_init(max_iter, abs_tol=abs_tol)
97  else
98  call this%ksp_init(max_iter)
99  end if
100 
101  end subroutine sx_pipecg_init
102 
104  subroutine sx_pipecg_free(this)
105  class(sx_pipecg_t), intent(inout) :: this
106 
107  call this%ksp_free()
108 
109  if (allocated(this%p)) then
110  deallocate(this%p)
111  end if
112  if (allocated(this%q)) then
113  deallocate(this%q)
114  end if
115  if (allocated(this%r)) then
116  deallocate(this%r)
117  end if
118  if (allocated(this%s)) then
119  deallocate(this%s)
120  end if
121  if (allocated(this%u)) then
122  deallocate(this%u)
123  end if
124  if (allocated(this%w)) then
125  deallocate(this%w)
126  end if
127  if (allocated(this%z)) then
128  deallocate(this%z)
129  end if
130  if (allocated(this%mi)) then
131  deallocate(this%mi)
132  end if
133  if (allocated(this%ni)) then
134  deallocate(this%ni)
135  end if
136 
137  nullify(this%M)
138 
139 
140  end subroutine sx_pipecg_free
141 
143  function sx_pipecg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) result(ksp_results)
144  class(sx_pipecg_t), intent(inout) :: this
145  class(ax_t), intent(inout) :: ax
146  type(field_t), intent(inout) :: x
147  integer, intent(in) :: n
148  real(kind=rp), dimension(n), intent(inout) :: f
149  type(coef_t), intent(inout) :: coef
150  type(bc_list_t), intent(inout) :: blst
151  type(gs_t), intent(inout) :: gs_h
152  type(ksp_monitor_t) :: ksp_results
153  integer, optional, intent(in) :: niter
154  integer :: iter, max_iter, i, ierr
155  real(kind=rp) :: rnorm, rtr, reduction(3), norm_fac
156  real(kind=rp) :: alpha, beta, gamma1, gamma2, delta
157  real(kind=rp) :: tmp1, tmp2, tmp3
158  type(mpi_request) :: request
159  type(mpi_status) :: status
160 
161  if (present(niter)) then
162  max_iter = niter
163  else
164  max_iter = this%max_iter
165  end if
166  norm_fac = 1.0_rp / sqrt(coef%volume)
167 
168  do i = 1, n
169  x%x(i,1,1,1) = 0.0_rp
170  this%z(i) = 0.0_rp
171  this%q(i) = 0.0_rp
172  this%p(i) = 0.0_rp
173  this%s(i) = 0.0_rp
174  this%r(i) = f(i)
175  end do
176 
177  call this%M%solve(this%u, this%r, n)
178  call ax%compute(this%w, this%u, coef, x%msh, x%Xh)
179  call gs_h%op(this%w, n, gs_op_add)
180  call bc_list_apply(blst, this%w, n)
181 
182  rtr = glsc3(this%r, coef%mult, this%r, n)
183  rnorm = sqrt(rtr)*norm_fac
184  ksp_results%res_start = rnorm
185  ksp_results%res_final = rnorm
186  ksp_results%iter = 0
187  if(abscmp(rnorm, 0.0_rp)) return
188 
189  gamma1 = 0.0_rp
190 
191  do iter = 1, max_iter
192 
193  tmp1 = 0.0_rp
194  tmp2 = 0.0_rp
195  tmp3 = 0.0_rp
196  do i = 1, n
197  tmp1 = tmp1 + this%r(i) * coef%mult(i,1,1,1) * this%u(i)
198  tmp2 = tmp2 + this%w(i) * coef%mult(i,1,1,1) * this%u(i)
199  tmp3 = tmp3 + this%r(i) * coef%mult(i,1,1,1) * this%r(i)
200  end do
201  reduction(1) = tmp1
202  reduction(2) = tmp2
203  reduction(3) = tmp3
204 
205  call mpi_iallreduce(mpi_in_place, reduction, 3, &
206  mpi_real_precision, mpi_sum, neko_comm, request, ierr)
207 
208  call this%M%solve(this%mi, this%w, n)
209  call ax%compute(this%ni, this%mi, coef, x%msh, x%Xh)
210  call gs_h%op(this%ni, n, gs_op_add)
211  call bc_list_apply(blst, this%ni, n)
212 
213  call mpi_wait(request, status, ierr)
214  gamma2 = gamma1
215  gamma1 = reduction(1)
216  delta = reduction(2)
217  rtr = reduction(3)
218 
219  rnorm = sqrt(rtr)*norm_fac
220  if (rnorm .lt. this%abs_tol) then
221  exit
222  end if
223 
224  if (iter .gt. 1) then
225  beta = gamma1 / gamma2
226  alpha = gamma1 / (delta - (beta * gamma1/alpha))
227  else
228  beta = 0.0_rp
229  alpha = gamma1/delta
230  end if
231 
232  do i = 1, n
233  this%z(i) = beta * this%z(i) + this%ni(i)
234  this%q(i) = beta * this%q(i) + this%mi(i)
235  this%s(i) = beta * this%s(i) + this%w(i)
236  this%p(i) = beta * this%p(i) + this%u(i)
237  end do
238 
239  do i = 1, n
240  x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * this%p(i)
241  this%r(i) = this%r(i) - alpha * this%s(i)
242  this%u(i) = this%u(i) - alpha * this%q(i)
243  this%w(i) = this%w(i) - alpha * this%z(i)
244  end do
245 
246  end do
247 
248  ksp_results%res_final = rnorm
249  ksp_results%iter = iter
250 
251  end function sx_pipecg_solve
252 
253 end module pipecg_sx
254 
255 
Defines a Matrix-vector product.
Definition: ax.f90:34
Defines a boundary condition.
Definition: bc.f90:34
Coefficients.
Definition: coef.f90:34
Definition: comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
type(mpi_datatype) mpi_real_precision
MPI type for working precision of REAL types.
Definition: comm.F90:22
Defines a field.
Definition: field.f90:34
Gather-scatter.
Implements the base abstract type for Krylov solvers plus helper types.
Definition: krylov.f90:34
integer, parameter, public ksp_max_iter
Maximum number of iters.
Definition: krylov.f90:50
Definition: math.f90:60
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition: math.f90:810
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a pipelined Conjugate Gradient methods SX-Aurora backend.
Definition: pipecg_sx.f90:34
subroutine sx_pipecg_init(this, n, max_iter, M, rel_tol, abs_tol)
Initialise a pipelined PCG solver.
Definition: pipecg_sx.f90:69
subroutine sx_pipecg_free(this)
Deallocate a pipelined PCG solver.
Definition: pipecg_sx.f90:105
type(ksp_monitor_t) function sx_pipecg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter)
Pipelined PCG solve.
Definition: pipecg_sx.f90:144
Krylov preconditioner.
Definition: precon.f90:34
Base type for a matrix-vector product providing .
Definition: ax.f90:43
A list of boundary conditions.
Definition: bc.f90:102
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
Type for storing initial and final residuals in a Krylov solver.
Definition: krylov.f90:55
Base abstract type for a canonical Krylov method, solving .
Definition: krylov.f90:65
Pipelined preconditioned conjugate gradient method for SX-Aurora.
Definition: pipecg_sx.f90:49
Defines a canonical Krylov preconditioner.
Definition: precon.f90:40