Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
cg_coupled.f90
Go to the documentation of this file.
1! Copyright (c) 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!
34module cg_cpld
35 use num_types, only: rp
37 use precon, only : pc_t
38 use ax_product, only : ax_t
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, glsc2, abscmp
44 use utils, only : neko_error
45 use operators, only : rotate_cyc
46 implicit none
47 private
48
50 type, public, extends(ksp_t) :: cg_cpld_t
51 real(kind=rp), allocatable :: w1(:)
52 real(kind=rp), allocatable :: w2(:)
53 real(kind=rp), allocatable :: w3(:)
54 real(kind=rp), allocatable :: r1(:)
55 real(kind=rp), allocatable :: r2(:)
56 real(kind=rp), allocatable :: r3(:)
57 real(kind=rp), allocatable :: p1(:)
58 real(kind=rp), allocatable :: p2(:)
59 real(kind=rp), allocatable :: p3(:)
60 real(kind=rp), allocatable :: z1(:)
61 real(kind=rp), allocatable :: z2(:)
62 real(kind=rp), allocatable :: z3(:)
63 real(kind=rp), allocatable :: tmp(:)
64 contains
65 procedure, pass(this) :: init => cg_cpld_init
66 procedure, pass(this) :: free => cg_cpld_free
67 procedure, pass(this) :: solve => cg_cpld_nop
68 procedure, pass(this) :: solve_coupled => cg_cpld_solve
69 end type cg_cpld_t
70
71contains
72
74 subroutine cg_cpld_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
75 class(cg_cpld_t), target, intent(inout) :: this
76 integer, intent(in) :: max_iter
77 class(pc_t), optional, intent(in), target :: M
78 integer, intent(in) :: n
79 real(kind=rp), optional, intent(in) :: rel_tol
80 real(kind=rp), optional, intent(in) :: abs_tol
81 logical, optional, intent(in) :: monitor
82
83 call this%free()
84
85 allocate(this%w1(n))
86 allocate(this%w2(n))
87 allocate(this%w3(n))
88 allocate(this%r1(n))
89 allocate(this%r2(n))
90 allocate(this%r3(n))
91 allocate(this%p1(n))
92 allocate(this%p2(n))
93 allocate(this%p3(n))
94 allocate(this%z1(n))
95 allocate(this%z2(n))
96 allocate(this%z3(n))
97 allocate(this%tmp(n))
98
99 if (present(m)) then
100 this%M => m
101 end if
102
103 if (present(rel_tol) .and. present(abs_tol) .and. present(monitor)) then
104 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
105 else if (present(rel_tol) .and. present(abs_tol)) then
106 call this%ksp_init(max_iter, rel_tol, abs_tol)
107 else if (present(monitor) .and. present(abs_tol)) then
108 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
109 else if (present(rel_tol) .and. present(monitor)) then
110 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
111 else if (present(rel_tol)) then
112 call this%ksp_init(max_iter, rel_tol = rel_tol)
113 else if (present(abs_tol)) then
114 call this%ksp_init(max_iter, abs_tol = abs_tol)
115 else if (present(monitor)) then
116 call this%ksp_init(max_iter, monitor = monitor)
117 else
118 call this%ksp_init(max_iter)
119 end if
120
121 end subroutine cg_cpld_init
122
124 subroutine cg_cpld_free(this)
125 class(cg_cpld_t), intent(inout) :: this
126
127 call this%ksp_free()
128
129 if (allocated(this%w1)) then
130 deallocate(this%w1)
131 end if
132
133 if (allocated(this%w2)) then
134 deallocate(this%w2)
135 end if
136
137 if (allocated(this%w3)) then
138 deallocate(this%w3)
139 end if
140
141 if (allocated(this%r1)) then
142 deallocate(this%r1)
143 end if
144
145 if (allocated(this%r2)) then
146 deallocate(this%r2)
147 end if
148
149 if (allocated(this%r3)) then
150 deallocate(this%r3)
151 end if
152
153 if (allocated(this%p1)) then
154 deallocate(this%p1)
155 end if
156
157 if (allocated(this%p2)) then
158 deallocate(this%p2)
159 end if
160
161 if (allocated(this%p3)) then
162 deallocate(this%p3)
163 end if
164
165 if (allocated(this%z1)) then
166 deallocate(this%z1)
167 end if
168
169 if (allocated(this%z2)) then
170 deallocate(this%z2)
171 end if
172
173 if (allocated(this%z3)) then
174 deallocate(this%z3)
175 end if
176
177 if (allocated(this%tmp)) then
178 deallocate(this%tmp)
179 end if
180
181 nullify(this%M)
182
183 end subroutine cg_cpld_free
184
185 function cg_cpld_nop(this, Ax, x, f, n, coef, blst, gs_h, niter) &
186 result(ksp_results)
187 class(cg_cpld_t), intent(inout) :: this
188 class(ax_t), intent(in) :: ax
189 type(field_t), intent(inout) :: x
190 integer, intent(in) :: n
191 real(kind=rp), dimension(n), intent(in) :: f
192 type(coef_t), intent(inout) :: coef
193 type(bc_list_t), intent(inout) :: blst
194 type(gs_t), intent(inout) :: gs_h
195 type(ksp_monitor_t) :: ksp_results
196 integer, optional, intent(in) :: niter
197
198 ! Throw and error
199 call neko_error('The cpldcg solver is only defined for coupled solves')
200
201 ksp_results%res_final = 0.0
202 ksp_results%iter = 0
203 end function cg_cpld_nop
204
206 function cg_cpld_solve(this, Ax, x, y, z, fx, fy, fz, &
207 n, coef, blstx, blsty, blstz, gs_h, niter) result(ksp_results)
208 class(cg_cpld_t), intent(inout) :: this
209 class(ax_t), intent(in) :: ax
210 type(field_t), intent(inout) :: x
211 type(field_t), intent(inout) :: y
212 type(field_t), intent(inout) :: z
213 integer, intent(in) :: n
214 real(kind=rp), dimension(n), intent(in) :: fx
215 real(kind=rp), dimension(n), intent(in) :: fy
216 real(kind=rp), dimension(n), intent(in) :: fz
217 type(coef_t), intent(inout) :: coef
218 type(bc_list_t), intent(inout) :: blstx
219 type(bc_list_t), intent(inout) :: blsty
220 type(bc_list_t), intent(inout) :: blstz
221 type(gs_t), intent(inout) :: gs_h
222 type(ksp_monitor_t), dimension(3) :: ksp_results
223 integer, optional, intent(in) :: niter
224 integer :: i, iter, max_iter
225 real(kind=rp) :: rnorm, rtr, rtr0, rtz2, rtz1
226 real(kind=rp) :: beta, pap, alpha, alphm, norm_fac
227
228 if (present(niter)) then
229 max_iter = niter
230 else
231 max_iter = this%max_iter
232 end if
233 norm_fac = 1.0_rp / coef%volume
234
235 associate(p1 => this%p1, p2 => this%p2, p3 => this%p3, z1 => this%z1, &
236 z2 => this%z2, z3 => this%z3, r1 => this%r1, r2 => this%r2, &
237 r3 => this%r3, tmp => this%tmp, w1 => this%w1, w2 => this%w2, &
238 w3 => this%w3)
239
240 rtz1 = 1.0_rp
241 do concurrent(i = 1:n)
242 x%x(i,1,1,1) = 0.0_rp
243 y%x(i,1,1,1) = 0.0_rp
244 z%x(i,1,1,1) = 0.0_rp
245 p1(i) = 0.0_rp
246 p2(i) = 0.0_rp
247 p3(i) = 0.0_rp
248 z1(i) = 0.0_rp
249 z2(i) = 0.0_rp
250 z3(i) = 0.0_rp
251 r1(i) = fx(i)
252 r2(i) = fy(i)
253 r3(i) = fz(i)
254 tmp(i) = r1(i)**2 + r2(i)**2 + r3(i)**2
255 end do
256
257 rtr = glsc3(tmp, coef%mult, coef%binv, n)
258 rnorm = sqrt(rtr)*norm_fac
259 ksp_results%res_start = rnorm
260 ksp_results%res_final = rnorm
261 ksp_results%iter = 0
262 if (abscmp(rnorm, 0.0_rp)) then
263 ksp_results%converged = .true.
264 return
265 end if
266
267 call this%monitor_start('cpldCG')
268 do iter = 1, max_iter
269 call this%M%solve(z1, this%r1, n)
270 call this%M%solve(z2, this%r2, n)
271 call this%M%solve(z3, this%r3, n)
272 rtz2 = rtz1
273
274 do concurrent(i = 1:n)
275 this%tmp(i) = z1(i) * r1(i) &
276 + z2(i) * r2(i) &
277 + z3(i) * r3(i)
278 end do
279
280 rtz1 = glsc2(tmp, coef%mult, n)
281
282 beta = rtz1 / rtz2
283 if (iter .eq. 1) beta = 0.0_rp
284 do concurrent(i = 1:n)
285 p1(i) = p1(i) * beta + z1(i)
286 p2(i) = p2(i) * beta + z2(i)
287 p3(i) = p3(i) * beta + z3(i)
288 end do
289
290 call ax%compute_vector(w1, w2, w3, p1, p2, p3, coef, x%msh, x%Xh)
291
292 call rotate_cyc(w1, w2, w3, 1, coef)
293 call gs_h%op(w1, n, gs_op_add)
294 call gs_h%op(w2, n, gs_op_add)
295 call gs_h%op(w3, n, gs_op_add)
296 call rotate_cyc(w1, w2, w3, 0, coef)
297
298 call blstx%apply_scalar(w1, n)
299 call blsty%apply_scalar(w2, n)
300 call blstz%apply_scalar(w3, n)
301
302 do concurrent(i = 1:n)
303 tmp(i) = w1(i) * p1(i) &
304 + w2(i) * p2(i) &
305 + w3(i) * p3(i)
306 end do
307
308 pap = glsc2(tmp, coef%mult, n)
309
310 alpha = rtz1 / pap
311 alphm = -alpha
312 do concurrent(i = 1:n)
313 x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * p1(i)
314 y%x(i,1,1,1) = y%x(i,1,1,1) + alpha * p2(i)
315 z%x(i,1,1,1) = z%x(i,1,1,1) + alpha * p3(i)
316 r1(i) = r1(i) + alphm * w1(i)
317 r2(i) = r2(i) + alphm * w2(i)
318 r3(i) = r3(i) + alphm * w3(i)
319 tmp(i) = r1(i)**2 + r2(i)**2 + r3(i)**2
320 end do
321
322 rtr = glsc3(tmp, coef%mult, coef%binv, n)
323 if (iter .eq. 1) rtr0 = rtr
324 rnorm = sqrt(rtr) * norm_fac
325 call this%monitor_iter(iter, rnorm)
326 if (rnorm .lt. this%abs_tol) then
327 exit
328 end if
329 end do
330 end associate
331 call this%monitor_stop()
332 ksp_results%res_final = rnorm
333 ksp_results%iter = iter
334 ksp_results%converged = this%is_converged(iter, rnorm)
335 end function cg_cpld_solve
336
337end module cg_cpld
__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
Defines a coupled Conjugate Gradient methods.
type(ksp_monitor_t) function, dimension(3) cg_cpld_solve(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Coupled PCG solve.
type(ksp_monitor_t) function cg_cpld_nop(this, ax, x, f, n, coef, blst, gs_h, niter)
subroutine cg_cpld_free(this)
Deallocate a coupled PCG solver.
subroutine cg_cpld_init(this, n, max_iter, m, rel_tol, abs_tol, monitor)
Initialise a coupled PCG solver.
Coefficients.
Definition coef.f90:34
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:1067
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
Krylov preconditioner.
Definition precon.f90:34
Utilities.
Definition utils.f90:35
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
Coupled preconditioned conjugate gradient method.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
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
Defines a canonical Krylov preconditioner.
Definition precon.f90:40