Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
pipecg_sx.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2025, 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!
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_list, only : bc_list_t
43 use math, only : glsc3, abscmp
45 use mpi_f08, only : mpi_iallreduce, mpi_in_place, mpi_sum, mpi_wait, &
46 mpi_request, mpi_status
47 implicit none
48 private
49
51 type, public, extends(ksp_t) :: sx_pipecg_t
52 real(kind=rp), allocatable :: p(:)
53 real(kind=rp), allocatable :: q(:)
54 real(kind=rp), allocatable :: r(:)
55 real(kind=rp), allocatable :: s(:)
56 real(kind=rp), allocatable :: u(:)
57 real(kind=rp), allocatable :: w(:)
58 real(kind=rp), allocatable :: z(:)
59 real(kind=rp), allocatable :: mi(:)
60 real(kind=rp), allocatable :: ni(:)
61 contains
62 procedure, pass(this) :: init => sx_pipecg_init
63 procedure, pass(this) :: free => sx_pipecg_free
64 procedure, pass(this) :: solve => sx_pipecg_solve
65 procedure, pass(this) :: solve_coupled => sx_pipecg_solve_coupled
66 end type sx_pipecg_t
67
68contains
69
71 subroutine sx_pipecg_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
72 class(sx_pipecg_t), target, intent(inout) :: this
73 class(pc_t), optional, intent(in), target :: M
74 integer, intent(in) :: n
75 integer, intent(in) :: max_iter
76 real(kind=rp), optional, intent(in) :: rel_tol
77 real(kind=rp), optional, intent(in) :: abs_tol
78 logical, optional, intent(in) :: monitor
79
80 call this%free()
81
82 allocate(this%p(n))
83 allocate(this%q(n))
84 allocate(this%r(n))
85 allocate(this%s(n))
86 allocate(this%u(n))
87 allocate(this%w(n))
88 allocate(this%z(n))
89 allocate(this%mi(n))
90 allocate(this%ni(n))
91 if (present(m)) then
92 this%M => m
93 end if
94
95 if (present(rel_tol) .and. present(abs_tol) .and. present(monitor)) then
96 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
97 else if (present(rel_tol) .and. present(abs_tol)) then
98 call this%ksp_init(max_iter, rel_tol, abs_tol)
99 else if (present(monitor) .and. present(abs_tol)) then
100 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
101 else if (present(rel_tol) .and. present(monitor)) then
102 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
103 else if (present(rel_tol)) then
104 call this%ksp_init(max_iter, rel_tol = rel_tol)
105 else if (present(abs_tol)) then
106 call this%ksp_init(max_iter, abs_tol = abs_tol)
107 else if (present(monitor)) then
108 call this%ksp_init(max_iter, monitor = monitor)
109 else
110 call this%ksp_init(max_iter)
111 end if
112
113 end subroutine sx_pipecg_init
114
116 subroutine sx_pipecg_free(this)
117 class(sx_pipecg_t), intent(inout) :: this
118
119 call this%ksp_free()
120
121 if (allocated(this%p)) then
122 deallocate(this%p)
123 end if
124 if (allocated(this%q)) then
125 deallocate(this%q)
126 end if
127 if (allocated(this%r)) then
128 deallocate(this%r)
129 end if
130 if (allocated(this%s)) then
131 deallocate(this%s)
132 end if
133 if (allocated(this%u)) then
134 deallocate(this%u)
135 end if
136 if (allocated(this%w)) then
137 deallocate(this%w)
138 end if
139 if (allocated(this%z)) then
140 deallocate(this%z)
141 end if
142 if (allocated(this%mi)) then
143 deallocate(this%mi)
144 end if
145 if (allocated(this%ni)) then
146 deallocate(this%ni)
147 end if
148
149 nullify(this%M)
150
151
152 end subroutine sx_pipecg_free
153
155 function sx_pipecg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) result(ksp_results)
156 class(sx_pipecg_t), intent(inout) :: this
157 class(ax_t), intent(in) :: ax
158 type(field_t), intent(inout) :: x
159 integer, intent(in) :: n
160 real(kind=rp), dimension(n), intent(in) :: f
161 type(coef_t), intent(inout) :: coef
162 type(bc_list_t), intent(inout) :: blst
163 type(gs_t), intent(inout) :: gs_h
164 type(ksp_monitor_t) :: ksp_results
165 integer, optional, intent(in) :: niter
166 integer :: iter, max_iter, i, ierr
167 real(kind=rp) :: rnorm, rtr, reduction(3), norm_fac
168 real(kind=rp) :: alpha, beta, gamma1, gamma2, delta
169 real(kind=rp) :: tmp1, tmp2, tmp3
170 type(mpi_request) :: request
171 type(mpi_status) :: status
172
173 if (present(niter)) then
174 max_iter = niter
175 else
176 max_iter = this%max_iter
177 end if
178 norm_fac = 1.0_rp / sqrt(coef%volume)
179
180 do i = 1, n
181 x%x(i,1,1,1) = 0.0_rp
182 this%z(i) = 0.0_rp
183 this%q(i) = 0.0_rp
184 this%p(i) = 0.0_rp
185 this%s(i) = 0.0_rp
186 this%r(i) = f(i)
187 end do
188
189 call this%M%solve(this%u, this%r, n)
190 call ax%compute(this%w, this%u, coef, x%msh, x%Xh)
191 call gs_h%op(this%w, n, gs_op_add)
192 call blst%apply_scalar(this%w, n)
193
194 rtr = glsc3(this%r, coef%mult, this%r, n)
195 rnorm = sqrt(rtr)*norm_fac
196 ksp_results%res_start = rnorm
197 ksp_results%res_final = rnorm
198 ksp_results%iter = 0
199 if(abscmp(rnorm, 0.0_rp)) then
200 ksp_results%converged = .true.
201 return
202 end if
203
204 gamma1 = 0.0_rp
205
206 call this%monitor_start('PipeCG')
207 do iter = 1, max_iter
208
209 tmp1 = 0.0_rp
210 tmp2 = 0.0_rp
211 tmp3 = 0.0_rp
212 do i = 1, n
213 tmp1 = tmp1 + this%r(i) * coef%mult(i,1,1,1) * this%u(i)
214 tmp2 = tmp2 + this%w(i) * coef%mult(i,1,1,1) * this%u(i)
215 tmp3 = tmp3 + this%r(i) * coef%mult(i,1,1,1) * this%r(i)
216 end do
217 reduction(1) = tmp1
218 reduction(2) = tmp2
219 reduction(3) = tmp3
220
221 call mpi_iallreduce(mpi_in_place, reduction, 3, &
222 mpi_real_precision, mpi_sum, neko_comm, request, ierr)
223
224 call this%M%solve(this%mi, this%w, n)
225 call ax%compute(this%ni, this%mi, coef, x%msh, x%Xh)
226 call gs_h%op(this%ni, n, gs_op_add)
227 call blst%apply(this%ni, n)
228
229 call mpi_wait(request, status, ierr)
230 gamma2 = gamma1
231 gamma1 = reduction(1)
232 delta = reduction(2)
233 rtr = reduction(3)
234
235 rnorm = sqrt(rtr)*norm_fac
236 call this%monitor_iter(iter, rnorm)
237 if (rnorm .lt. this%abs_tol) then
238 exit
239 end if
240
241 if (iter .gt. 1) then
242 beta = gamma1 / gamma2
243 alpha = gamma1 / (delta - (beta * gamma1/alpha))
244 else
245 beta = 0.0_rp
246 alpha = gamma1/delta
247 end if
248
249 do i = 1, n
250 this%z(i) = beta * this%z(i) + this%ni(i)
251 this%q(i) = beta * this%q(i) + this%mi(i)
252 this%s(i) = beta * this%s(i) + this%w(i)
253 this%p(i) = beta * this%p(i) + this%u(i)
254 end do
255
256 do i = 1, n
257 x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * this%p(i)
258 this%r(i) = this%r(i) - alpha * this%s(i)
259 this%u(i) = this%u(i) - alpha * this%q(i)
260 this%w(i) = this%w(i) - alpha * this%z(i)
261 end do
262
263 end do
264 call this%monitor_stop()
265 ksp_results%res_final = rnorm
266 ksp_results%iter = iter
267
268 end function sx_pipecg_solve
269
271 function sx_pipecg_solve_coupled(this, Ax, x, y, z, fx, fy, fz, &
272 n, coef, blstx, blsty, blstz, gs_h, niter) result(ksp_results)
273 class(sx_pipecg_t), intent(inout) :: this
274 class(ax_t), intent(in) :: ax
275 type(field_t), intent(inout) :: x
276 type(field_t), intent(inout) :: y
277 type(field_t), intent(inout) :: z
278 integer, intent(in) :: n
279 real(kind=rp), dimension(n), intent(in) :: fx
280 real(kind=rp), dimension(n), intent(in) :: fy
281 real(kind=rp), dimension(n), intent(in) :: fz
282 type(coef_t), intent(inout) :: coef
283 type(bc_list_t), intent(inout) :: blstx
284 type(bc_list_t), intent(inout) :: blsty
285 type(bc_list_t), intent(inout) :: blstz
286 type(gs_t), intent(inout) :: gs_h
287 type(ksp_monitor_t), dimension(3) :: ksp_results
288 integer, optional, intent(in) :: niter
289
290 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
291 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
292 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
293
294 end function sx_pipecg_solve_coupled
295
296end module pipecg_sx
297
298
__device__ T solve(const T u, const T y, const T guess, const T nu, const T kappa, const T B)
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:50
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
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:51
Definition math.f90:60
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1026
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
type(ksp_monitor_t) function, dimension(3) sx_pipecg_solve_coupled(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Pipelined PCG coupled solve.
subroutine sx_pipecg_free(this)
Deallocate a pipelined PCG solver.
subroutine sx_pipecg_init(this, n, max_iter, m, rel_tol, abs_tol, monitor)
Initialise a pipelined PCG solver.
Definition pipecg_sx.f90:72
type(ksp_monitor_t) function sx_pipecg_solve(this, ax, x, f, n, coef, blst, gs_h, niter)
Pipelined PCG solve.
Krylov preconditioner.
Definition precon.f90:34
Base type for a matrix-vector product providing .
Definition ax.f90:43
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:73
Pipelined preconditioned conjugate gradient method for SX-Aurora.
Definition pipecg_sx.f90:51
Defines a canonical Krylov preconditioner.
Definition precon.f90:40