Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
cheby_device.F90
Go to the documentation of this file.
1! Copyright (c) 2024-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!
35 use krylov, only : ksp_t, ksp_monitor_t
36 use precon, only : pc_t
37 use ax_product, only : ax_t
38 use num_types, only : rp, c_rp
39 use field, only : field_t
40 use coefs, only : coef_t
41 use mesh, only : mesh_t
42 use space, only : space_t
43 use gather_scatter, only : gs_t, gs_op_add
44 use bc_list, only : bc_list_t
45 use schwarz, only : schwarz_t
49 use device
50 use, intrinsic :: iso_c_binding, only : c_ptr, c_int, &
51 c_null_ptr, c_associated
52 implicit none
53 private
54
56 type, public, extends(ksp_t) :: cheby_device_t
57 real(kind=rp), allocatable :: d(:)
58 real(kind=rp), allocatable :: w(:)
59 real(kind=rp), allocatable :: r(:)
60 type(c_ptr) :: d_d = c_null_ptr
61 type(c_ptr) :: w_d = c_null_ptr
62 type(c_ptr) :: r_d = c_null_ptr
63 type(c_ptr) :: gs_event = c_null_ptr
64 real(kind=rp) :: tha, dlt
65 integer :: power_its = 150
66 logical :: recompute_eigs = .true.
67 logical :: zero_initial_guess = .false.
68 type(schwarz_t), pointer :: schwarz => null()
69 contains
70 procedure, pass(this) :: init => cheby_device_init
71 procedure, pass(this) :: free => cheby_device_free
72 procedure, pass(this) :: solve => cheby_device_impl
73 procedure, pass(this) :: solve_coupled => cheby_device_solve_coupled
74 end type cheby_device_t
75
76#ifdef HAVE_HIP
77 interface
78 subroutine hip_cheby_device_part1(d_d, x_d, inv_tha, n, strm) &
79 bind(c, name = 'hip_cheby_part1')
80 use, intrinsic :: iso_c_binding
81 import c_rp
82 implicit none
83 type(c_ptr), value :: d_d, x_d, strm
84 real(c_rp) :: inv_tha
85 integer(c_int) :: n
86 end subroutine hip_cheby_device_part1
87 end interface
88
89 interface
90 subroutine hip_cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n, strm) &
91 bind(c, name = 'hip_cheby_part2')
92 use, intrinsic :: iso_c_binding
93 import c_rp
94 implicit none
95 type(c_ptr), value :: d_d, w_d, x_d, strm
96 real(c_rp) :: tmp1, tmp2
97 integer(c_int) :: n
98 end subroutine hip_cheby_device_part2
99 end interface
100#elif HAVE_CUDA
101 interface
102 subroutine cuda_cheby_device_part1(d_d, x_d, inv_tha, n, strm) &
103 bind(c, name = 'cuda_cheby_part1')
104 use, intrinsic :: iso_c_binding
105 import c_rp
106 implicit none
107 type(c_ptr), value :: d_d, x_d, strm
108 real(c_rp) :: inv_tha
109 integer(c_int) :: n
110 end subroutine cuda_cheby_device_part1
111 end interface
112
113 interface
114 subroutine cuda_cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n, strm) &
115 bind(c, name = 'cuda_cheby_part2')
116 use, intrinsic :: iso_c_binding
117 import c_rp
118 implicit none
119 type(c_ptr), value :: d_d, w_d, x_d, strm
120 real(c_rp) :: tmp1, tmp2
121 integer(c_int) :: n
122 end subroutine cuda_cheby_device_part2
123 end interface
124#elif HAVE_METAL
125 interface
126 subroutine metal_cheby_device_part1(d_d, x_d, inv_tha, n, strm) &
127 bind(c, name = 'metal_cheby_part1')
128 use, intrinsic :: iso_c_binding
129 import c_rp
130 implicit none
131 type(c_ptr), value :: d_d, x_d, strm
132 real(c_rp) :: inv_tha
133 integer(c_int) :: n
134 end subroutine metal_cheby_device_part1
135 end interface
136
137 interface
138 subroutine metal_cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n, strm) &
139 bind(c, name = 'metal_cheby_part2')
140 use, intrinsic :: iso_c_binding
141 import c_rp
142 implicit none
143 type(c_ptr), value :: d_d, w_d, x_d, strm
144 real(c_rp) :: tmp1, tmp2
145 integer(c_int) :: n
146 end subroutine metal_cheby_device_part2
147 end interface
148#endif
149
150contains
151 subroutine cheby_device_part1(d_d, x_d, inv_tha, n)
152 type(c_ptr) :: d_d, x_d
153 real(c_rp) :: inv_tha
154 integer(c_int) :: n
155#ifdef HAVE_HIP
156 call hip_cheby_device_part1(d_d, x_d, inv_tha, n, glb_cmd_queue)
157#elif HAVE_CUDA
158 call cuda_cheby_device_part1(d_d, x_d, inv_tha, n, glb_cmd_queue)
159#elif HAVE_METAL
160 call metal_cheby_device_part1(d_d, x_d, inv_tha, n, glb_cmd_queue)
161#else !Fallback to device_math for missing device kernels
162 call device_cmult( d_d, inv_tha, n)
163 call device_add2( x_d, d_d, n)
164#endif
165 end subroutine cheby_device_part1
166
167
168
169 subroutine cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n)
170 type(c_ptr) :: d_d, w_d, x_d
171 real(c_rp) :: tmp1, tmp2
172 integer(c_int) :: n
173#ifdef HAVE_HIP
174 call hip_cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n, glb_cmd_queue)
175#elif HAVE_CUDA
176 call cuda_cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n, glb_cmd_queue)
177#elif HAVE_METAL
178 call metal_cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n, glb_cmd_queue)
179#else !Fallback to device_math for missing device kernels
180 call device_cmult( d_d, tmp1, n)
181 call device_add2s2( d_d, w_d, tmp2, n)
182 call device_add2( x_d, d_d, n)
183#endif
184 end subroutine cheby_device_part2
185
187 subroutine cheby_device_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
188 class(cheby_device_t), intent(inout), target :: this
189 integer, intent(in) :: max_iter
190 class(pc_t), optional, intent(in), target :: M
191 integer, intent(in) :: n
192 real(kind=rp), optional, intent(in) :: rel_tol
193 real(kind=rp), optional, intent(in) :: abs_tol
194 logical, optional, intent(in) :: monitor
195
196 call this%free()
197 allocate(this%d(n))
198 allocate(this%w(n))
199 allocate(this%r(n))
200
201 call device_map(this%d, this%d_d, n)
202 call device_map(this%w, this%w_d, n)
203 call device_map(this%r, this%r_d, n)
204
205 if (present(m)) then
206 this%M => m
207 end if
208
209 if (present(rel_tol) .and. present(abs_tol) .and. present(monitor)) then
210 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
211 else if (present(rel_tol) .and. present(abs_tol)) then
212 call this%ksp_init(max_iter, rel_tol, abs_tol)
213 else if (present(monitor) .and. present(abs_tol)) then
214 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
215 else if (present(rel_tol) .and. present(monitor)) then
216 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
217 else if (present(rel_tol)) then
218 call this%ksp_init(max_iter, rel_tol = rel_tol)
219 else if (present(abs_tol)) then
220 call this%ksp_init(max_iter, abs_tol = abs_tol)
221 else if (present(monitor)) then
222 call this%ksp_init(max_iter, monitor = monitor)
223 else
224 call this%ksp_init(max_iter)
225 end if
226
227 call device_event_create(this%gs_event, 2)
228
229 end subroutine cheby_device_init
230
231 subroutine cheby_device_free(this)
232 class(cheby_device_t), intent(inout) :: this
233
234 call this%ksp_free()
235
236 if (allocated(this%d)) then
237 if (c_associated(this%d_d)) then
238 call device_unmap(this%d, this%d_d)
239 end if
240 deallocate(this%d)
241 end if
242
243 if (allocated(this%w)) then
244 if (c_associated(this%w_d)) then
245 call device_unmap(this%w, this%w_d)
246 end if
247 deallocate(this%w)
248 end if
249
250 if (allocated(this%r)) then
251 if (c_associated(this%r_d)) then
252 call device_unmap(this%r, this%r_d)
253 end if
254 deallocate(this%r)
255 end if
256
257 nullify(this%M)
258
259 if (c_associated(this%gs_event)) then
260 call device_event_destroy(this%gs_event)
261 end if
262
263 end subroutine cheby_device_free
264
265 subroutine cheby_device_power(this, Ax, x, n, coef, blst, gs_h)
266 class(cheby_device_t), intent(inout) :: this
267 class(ax_t), intent(in) :: Ax
268 type(field_t), intent(inout) :: x
269 integer, intent(in) :: n
270 type(coef_t), intent(inout) :: coef
271 type(bc_list_t), intent(inout) :: blst
272 type(gs_t), intent(inout) :: gs_h
273 real(kind=rp) :: lam, b, a, rn
274 real(kind=rp) :: boost = 1.1_rp
275 real(kind=rp) :: lam_factor = 30.0_rp
276 real(kind=rp) :: wtw, dtw, dtd
277 integer :: i
278
279 associate(w => this%w, w_d => this%w_d, d => this%d, d_d => this%d_d)
280
281 do i = 1, n
282 !TODO: replace with a better way to initialize power method
283 call random_number(rn)
284 d(i) = rn + 10.0_rp
285 end do
286 call device_memcpy(d, d_d, n, host_to_device, sync = .true.)
287
288 call gs_h%op(d, n, gs_op_add, this%gs_event)
289 call blst%apply(d, n)
290
291 !Power method to get lamba max
292 do i = 1, this%power_its
293 call ax%compute(w, d, coef, x%msh, x%Xh)
294 call gs_h%op(w, n, gs_op_add, this%gs_event)
295 call blst%apply(w, n)
296 if (associated(this%schwarz)) then
297 call this%schwarz%compute(this%r, w)
298 call device_copy(w_d, this%r_d, n)
299 else
300 call this%M%solve(this%r, w, n)
301 call device_copy(w_d, this%r_d, n)
302 end if
303
304 wtw = device_glsc3(w_d, coef%mult_d, w_d, n)
305 call device_cmult2(d_d, w_d, 1.0_rp/sqrt(wtw), n)
306 call blst%apply(d, n)
307 end do
308
309 call ax%compute(w, d, coef, x%msh, x%Xh)
310 call gs_h%op(w, n, gs_op_add, this%gs_event)
311 call blst%apply(w, n)
312 if (associated(this%schwarz)) then
313 call this%schwarz%compute(this%r, w)
314 call device_copy(w_d, this%r_d, n)
315 else
316 call this%M%solve(this%r, w, n)
317 call device_copy(w_d, this%r_d, n)
318 end if
319
320 dtw = device_glsc3(d_d, coef%mult_d, w_d, n)
321 dtd = device_glsc3(d_d, coef%mult_d, d_d, n)
322 lam = dtw / dtd
323 b = lam * boost
324 a = lam / lam_factor
325 this%tha = (b+a)/2.0_rp
326 this%dlt = (b-a)/2.0_rp
327
328 this%recompute_eigs = .false.
329 end associate
330 end subroutine cheby_device_power
331
333 function cheby_device_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) &
334 result(ksp_results)
335 class(cheby_device_t), intent(inout) :: this
336 class(ax_t), intent(in) :: ax
337 type(field_t), intent(inout) :: x
338 integer, intent(in) :: n
339 real(kind=rp), dimension(n), intent(in) :: f
340 type(coef_t), intent(inout) :: coef
341 type(bc_list_t), intent(inout) :: blst
342 type(gs_t), intent(inout) :: gs_h
343 type(ksp_monitor_t) :: ksp_results
344 integer, optional, intent(in) :: niter
345 integer :: iter, max_iter
346 real(kind=rp) :: a, b, rtr, rnorm, norm_fac
347 type(c_ptr) :: f_d
348
349 f_d = device_get_ptr(f)
350
351 if (this%recompute_eigs) then
352 call cheby_device_power(this, ax, x, n, coef, blst, gs_h)
353 end if
354
355 if (present(niter)) then
356 max_iter = niter
357 else
358 max_iter = this%max_iter
359 end if
360 norm_fac = 1.0_rp / sqrt(coef%volume)
361
362 associate( w => this%w, r => this%r, d => this%d, &
363 w_d => this%w_d, r_d => this%r_d, d_d => this%d_d)
364 ! calculate residual
365 call device_copy(r_d, f_d, n)
366 call ax%compute(w, x%x, coef, x%msh, x%Xh)
367 call gs_h%op(w, n, gs_op_add, this%gs_event)
368 call blst%apply(w, n)
369 call device_sub2(r_d, w_d, n)
370
371 rtr = device_glsc3(r_d, coef%mult_d, r_d, n)
372 rnorm = sqrt(rtr) * norm_fac
373 ksp_results%res_start = rnorm
374 ksp_results%res_final = rnorm
375 ksp_results%iter = 0
376
377 ! First iteration
378 call this%M%solve(w, r, n)
379 call device_copy(d_d, w_d, n)
380 a = 2.0_rp / this%tha
381 call device_add2s2(x%x_d, d_d, a, n)! x = x + a*d
382
383 ! Rest of the iterations
384 do iter = 2, max_iter
385 ! calculate residual
386 call device_copy(r_d, f_d, n)
387 call ax%compute(w, x%x, coef, x%msh, x%Xh)
388 call gs_h%op(w, n, gs_op_add, this%gs_event)
389 call blst%apply(w, n)
390 call device_sub2(r_d, w_d, n)
391
392 call this%M%solve(w, r, n)
393
394 if (iter .eq. 2) then
395 b = 0.5_rp * (this%dlt * a)**2
396 else
397 b = (this%dlt * a / 2.0_rp)**2
398 end if
399 a = 1.0_rp/(this%tha - b/a)
400 call device_add2s1(d_d, w_d, b, n)! d = w + b*d
401
402 call device_add2s2(x%x_d, d_d, a, n)! x = x + a*d
403 end do
404
405 ! calculate residual
406 call device_copy(r_d, f_d, n)
407 call ax%compute(w, x%x, coef, x%msh, x%Xh)
408 call gs_h%op(w, n, gs_op_add, this%gs_event)
409 call blst%apply(w, n)
410 call device_sub2(r_d, w_d, n)
411 rtr = device_glsc3(r_d, coef%mult_d, r_d, n)
412 rnorm = sqrt(rtr) * norm_fac
413
414
415 ksp_results%res_final = rnorm
416 ksp_results%iter = iter
417 ksp_results%converged = this%is_converged(iter, rnorm)
418 end associate
419 end function cheby_device_solve
420
422 function cheby_device_impl(this, Ax, x, f, n, coef, blst, gs_h, niter) &
423 result(ksp_results)
424 class(cheby_device_t), intent(inout) :: this
425 class(ax_t), intent(in) :: ax
426 type(field_t), intent(inout) :: x
427 integer, intent(in) :: n
428 real(kind=rp), dimension(n), intent(in) :: f
429 type(coef_t), intent(inout) :: coef
430 type(bc_list_t), intent(inout) :: blst
431 type(gs_t), intent(inout) :: gs_h
432 type(ksp_monitor_t) :: ksp_results
433 integer, optional, intent(in) :: niter
434 integer :: iter, max_iter
435 real(kind=rp) :: a, b, rtr, rnorm, norm_fac
436 real(kind=rp) :: rhok, rhokp1, sig1, tmp1, tmp2
437 type(c_ptr) :: f_d
438
439 f_d = device_get_ptr(f)
440
441 if (this%recompute_eigs) then
442 call cheby_device_power(this, ax, x, n, coef, blst, gs_h)
443 end if
444
445 if (present(niter)) then
446 max_iter = niter
447 else
448 max_iter = this%max_iter
449 end if
450 norm_fac = 1.0_rp / sqrt(coef%volume)
451
452 associate( w => this%w, r => this%r, d => this%d, &
453 w_d => this%w_d, r_d => this%r_d, d_d => this%d_d)
454 ! calculate residual
455 if (.not.this%zero_initial_guess) then
456 call ax%compute(w, x%x, coef, x%msh, x%Xh)
457 call gs_h%op(w, n, gs_op_add, this%gs_event)
458 call blst%apply(w, n)
459 call device_sub3(r_d, f_d, w_d, n)
460 else
461 call device_copy(r_d, f_d, n)
462 this%zero_initial_guess = .false.
463 end if
464
465 ! First iteration
466 if (associated(this%schwarz)) then
467 call this%schwarz%compute(d, r)
468 else
469 call this%M%solve(d, r, n)
470 end if
471
472 tmp1 = 1.0_rp / this%tha
473 call cheby_device_part1(d_d, x%x_d, tmp1, n)
474
475 sig1 = this%tha / this%dlt
476 rhok = 1.0_rp / sig1
477
478 ! Rest of the iterations
479 do iter = 2, max_iter
480 rhokp1 = 1.0_rp / (2.0_rp * sig1 - rhok)
481 tmp1 = rhokp1 * rhok
482 tmp2 = 2.0_rp * rhokp1 / this%dlt
483 rhok = rhokp1
484 ! calculate residual
485 call ax%compute(w, x%x, coef, x%msh, x%Xh)
486 call gs_h%op(w, n, gs_op_add, this%gs_event)
487 call blst%apply(w, n)
488 call device_sub3(r_d, f_d, w_d, n)
489
490 if (associated(this%schwarz)) then
491 call this%schwarz%compute(w, r)
492 else
493 call this%M%solve(w, r, n)
494 end if
495
496 call cheby_device_part2(d_d, w_d, x%x_d, tmp1, tmp2, n)
497
498 end do
499
500 end associate
501 end function cheby_device_impl
502
504 function cheby_device_solve_coupled(this, Ax, x, y, z, fx, fy, fz, &
505 n, coef, blstx, blsty, blstz, gs_h, niter) result(ksp_results)
506 class(cheby_device_t), intent(inout) :: this
507 class(ax_t), intent(in) :: ax
508 type(field_t), intent(inout) :: x
509 type(field_t), intent(inout) :: y
510 type(field_t), intent(inout) :: z
511 integer, intent(in) :: n
512 real(kind=rp), dimension(n), intent(in) :: fx
513 real(kind=rp), dimension(n), intent(in) :: fy
514 real(kind=rp), dimension(n), intent(in) :: fz
515 type(coef_t), intent(inout) :: coef
516 type(bc_list_t), intent(inout) :: blstx
517 type(bc_list_t), intent(inout) :: blsty
518 type(bc_list_t), intent(inout) :: blstz
519 type(gs_t), intent(inout) :: gs_h
520 type(ksp_monitor_t), dimension(3) :: ksp_results
521 integer, optional, intent(in) :: niter
522
523 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
524 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
525 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
526
527 end function cheby_device_solve_coupled
528
529end module cheby_device
__device__ T solve(const T u, const T y, const T guess, const T nu, const T kappa, const T B)
Map a Fortran array to a device (allocate and associate)
Definition device.F90:78
Copy data between host and device (or device and device)
Definition device.F90:72
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:84
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Chebyshev preconditioner.
subroutine cheby_device_init(this, n, max_iter, m, rel_tol, abs_tol, monitor)
Initialise a standard solver.
subroutine cheby_device_part2(d_d, w_d, x_d, tmp1, tmp2, n)
type(ksp_monitor_t) function, dimension(3) cheby_device_solve_coupled(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Standard Cheby_Deviceshev coupled solve.
type(ksp_monitor_t) function cheby_device_solve(this, ax, x, f, n, coef, blst, gs_h, niter)
A chebyshev preconditioner.
subroutine cheby_device_free(this)
subroutine cheby_device_part1(d_d, x_d, inv_tha, n)
subroutine cheby_device_power(this, ax, x, n, coef, blst, gs_h)
type(ksp_monitor_t) function cheby_device_impl(this, ax, x, f, n, coef, blst, gs_h, niter)
A chebyshev preconditioner.
Coefficients.
Definition coef.f90:34
subroutine, public device_add2s1(a_d, b_d, c1, n, strm)
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_sub2(a_d, b_d, n, strm)
Vector substraction .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
real(kind=rp) function, public device_glsc3(a_d, b_d, c_d, n, strm)
Weighted inner product .
subroutine, public device_cmult2(a_d, b_d, c, n, strm)
Multiplication by constant c .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:48
subroutine, public device_event_destroy(event)
Destroy a device event.
Definition device.F90:1550
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:52
subroutine, public device_event_create(event, flags)
Create a device event queue.
Definition device.F90:1516
Defines a field.
Definition field.f90:34
Gather-scatter.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public c_rp
Definition num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Krylov preconditioner.
Definition precon.f90:34
Overlapping schwarz solves.
Definition schwarz.f90:61
Defines a function space.
Definition space.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:49
Defines a Chebyshev preconditioner.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
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
The function space for the SEM solution fields.
Definition space.f90:63