Neko 0.9.99
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-2024, 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
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 procedure, pass(this) :: solve_coupled => sx_pipecg_solve_coupled
64 end type sx_pipecg_t
65
66contains
67
69 subroutine sx_pipecg_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
70 class(sx_pipecg_t), intent(inout) :: this
71 class(pc_t), optional, intent(in), target :: M
72 integer, intent(in) :: n
73 integer, intent(in) :: max_iter
74 real(kind=rp), optional, intent(in) :: rel_tol
75 real(kind=rp), optional, intent(in) :: abs_tol
76 logical, optional, intent(in) :: monitor
77
78 call this%free()
79
80 allocate(this%p(n))
81 allocate(this%q(n))
82 allocate(this%r(n))
83 allocate(this%s(n))
84 allocate(this%u(n))
85 allocate(this%w(n))
86 allocate(this%z(n))
87 allocate(this%mi(n))
88 allocate(this%ni(n))
89 if (present(m)) then
90 this%M => m
91 end if
92
93 if (present(rel_tol) .and. present(abs_tol) .and. present(monitor)) then
94 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
95 else if (present(rel_tol) .and. present(abs_tol)) then
96 call this%ksp_init(max_iter, rel_tol, abs_tol)
97 else if (present(monitor) .and. present(abs_tol)) then
98 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
99 else if (present(rel_tol) .and. present(monitor)) then
100 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
101 else if (present(rel_tol)) then
102 call this%ksp_init(max_iter, rel_tol = rel_tol)
103 else if (present(abs_tol)) then
104 call this%ksp_init(max_iter, abs_tol = abs_tol)
105 else if (present(monitor)) then
106 call this%ksp_init(max_iter, monitor = monitor)
107 else
108 call this%ksp_init(max_iter)
109 end if
110
111 end subroutine sx_pipecg_init
112
114 subroutine sx_pipecg_free(this)
115 class(sx_pipecg_t), intent(inout) :: this
116
117 call this%ksp_free()
118
119 if (allocated(this%p)) then
120 deallocate(this%p)
121 end if
122 if (allocated(this%q)) then
123 deallocate(this%q)
124 end if
125 if (allocated(this%r)) then
126 deallocate(this%r)
127 end if
128 if (allocated(this%s)) then
129 deallocate(this%s)
130 end if
131 if (allocated(this%u)) then
132 deallocate(this%u)
133 end if
134 if (allocated(this%w)) then
135 deallocate(this%w)
136 end if
137 if (allocated(this%z)) then
138 deallocate(this%z)
139 end if
140 if (allocated(this%mi)) then
141 deallocate(this%mi)
142 end if
143 if (allocated(this%ni)) then
144 deallocate(this%ni)
145 end if
146
147 nullify(this%M)
148
149
150 end subroutine sx_pipecg_free
151
153 function sx_pipecg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) result(ksp_results)
154 class(sx_pipecg_t), intent(inout) :: this
155 class(ax_t), intent(in) :: ax
156 type(field_t), intent(inout) :: x
157 integer, intent(in) :: n
158 real(kind=rp), dimension(n), intent(in) :: f
159 type(coef_t), intent(inout) :: coef
160 type(bc_list_t), intent(inout) :: blst
161 type(gs_t), intent(inout) :: gs_h
162 type(ksp_monitor_t) :: ksp_results
163 integer, optional, intent(in) :: niter
164 integer :: iter, max_iter, i, ierr
165 real(kind=rp) :: rnorm, rtr, reduction(3), norm_fac
166 real(kind=rp) :: alpha, beta, gamma1, gamma2, delta
167 real(kind=rp) :: tmp1, tmp2, tmp3
168 type(mpi_request) :: request
169 type(mpi_status) :: status
170
171 if (present(niter)) then
172 max_iter = niter
173 else
174 max_iter = this%max_iter
175 end if
176 norm_fac = 1.0_rp / sqrt(coef%volume)
177
178 do i = 1, n
179 x%x(i,1,1,1) = 0.0_rp
180 this%z(i) = 0.0_rp
181 this%q(i) = 0.0_rp
182 this%p(i) = 0.0_rp
183 this%s(i) = 0.0_rp
184 this%r(i) = f(i)
185 end do
186
187 call this%M%solve(this%u, this%r, n)
188 call ax%compute(this%w, this%u, coef, x%msh, x%Xh)
189 call gs_h%op(this%w, n, gs_op_add)
190 call blst%apply_scalar(this%w, n)
191
192 rtr = glsc3(this%r, coef%mult, this%r, n)
193 rnorm = sqrt(rtr)*norm_fac
194 ksp_results%res_start = rnorm
195 ksp_results%res_final = rnorm
196 ksp_results%iter = 0
197 if(abscmp(rnorm, 0.0_rp)) return
198
199 gamma1 = 0.0_rp
200
201 call this%monitor_start('PipeCG')
202 do iter = 1, max_iter
203
204 tmp1 = 0.0_rp
205 tmp2 = 0.0_rp
206 tmp3 = 0.0_rp
207 do i = 1, n
208 tmp1 = tmp1 + this%r(i) * coef%mult(i,1,1,1) * this%u(i)
209 tmp2 = tmp2 + this%w(i) * coef%mult(i,1,1,1) * this%u(i)
210 tmp3 = tmp3 + this%r(i) * coef%mult(i,1,1,1) * this%r(i)
211 end do
212 reduction(1) = tmp1
213 reduction(2) = tmp2
214 reduction(3) = tmp3
215
216 call mpi_iallreduce(mpi_in_place, reduction, 3, &
217 mpi_real_precision, mpi_sum, neko_comm, request, ierr)
218
219 call this%M%solve(this%mi, this%w, n)
220 call ax%compute(this%ni, this%mi, coef, x%msh, x%Xh)
221 call gs_h%op(this%ni, n, gs_op_add)
222 call blst%apply(this%ni, n)
223
224 call mpi_wait(request, status, ierr)
225 gamma2 = gamma1
226 gamma1 = reduction(1)
227 delta = reduction(2)
228 rtr = reduction(3)
229
230 rnorm = sqrt(rtr)*norm_fac
231 call this%monitor_iter(iter, rnorm)
232 if (rnorm .lt. this%abs_tol) then
233 exit
234 end if
235
236 if (iter .gt. 1) then
237 beta = gamma1 / gamma2
238 alpha = gamma1 / (delta - (beta * gamma1/alpha))
239 else
240 beta = 0.0_rp
241 alpha = gamma1/delta
242 end if
243
244 do i = 1, n
245 this%z(i) = beta * this%z(i) + this%ni(i)
246 this%q(i) = beta * this%q(i) + this%mi(i)
247 this%s(i) = beta * this%s(i) + this%w(i)
248 this%p(i) = beta * this%p(i) + this%u(i)
249 end do
250
251 do i = 1, n
252 x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * this%p(i)
253 this%r(i) = this%r(i) - alpha * this%s(i)
254 this%u(i) = this%u(i) - alpha * this%q(i)
255 this%w(i) = this%w(i) - alpha * this%z(i)
256 end do
257
258 end do
259 call this%monitor_stop()
260 ksp_results%res_final = rnorm
261 ksp_results%iter = iter
262
263 end function sx_pipecg_solve
264
266 function sx_pipecg_solve_coupled(this, Ax, x, y, z, fx, fy, fz, &
267 n, coef, blstx, blsty, blstz, gs_h, niter) result(ksp_results)
268 class(sx_pipecg_t), intent(inout) :: this
269 class(ax_t), intent(in) :: ax
270 type(field_t), intent(inout) :: x
271 type(field_t), intent(inout) :: y
272 type(field_t), intent(inout) :: z
273 integer, intent(in) :: n
274 real(kind=rp), dimension(n), intent(in) :: fx
275 real(kind=rp), dimension(n), intent(in) :: fy
276 real(kind=rp), dimension(n), intent(in) :: fz
277 type(coef_t), intent(inout) :: coef
278 type(bc_list_t), intent(inout) :: blstx
279 type(bc_list_t), intent(inout) :: blsty
280 type(bc_list_t), intent(inout) :: blstz
281 type(gs_t), intent(inout) :: gs_h
282 type(ksp_monitor_t), dimension(3) :: ksp_results
283 integer, optional, intent(in) :: niter
284
285 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
286 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
287 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
288
289 end function sx_pipecg_solve_coupled
290
291end module pipecg_sx
292
293
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_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:23
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:894
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:70
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:46
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:68
Pipelined preconditioned conjugate gradient method for SX-Aurora.
Definition pipecg_sx.f90:49
Defines a canonical Krylov preconditioner.
Definition precon.f90:40