Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
phmg.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!
34module phmg
35 use num_types, only : rp
36 use precon, only : pc_t
37 use gather_scatter, only : gs_t, gs_op_add
38 use space, only : space_t, gll
39 use dofmap, only : dofmap_t
40 use field, only : field_t
41 use coefs, only : coef_t
42 use mesh, only : mesh_t
43 use bc, only : bc_t
44 use bc_list, only : bc_list_t
45 use dirichlet, only : dirichlet_t
46 use utils, only : neko_error
47 use cheby, only : cheby_t
49 use jacobi, only : jacobi_t
51 use ax_product, only : ax_t, ax_helm_factory
54 use math, only : copy, col2, add2
55 use device, only : device_get_ptr
60 use krylov, only : ksp_t, ksp_monitor_t, ksp_max_iter, &
61 krylov_solver_factory, krylov_solver_destroy
62 use logger, only : neko_log, log_size
63 use, intrinsic :: iso_c_binding
64 implicit none
65 private
66
67
68 type, private :: phmg_lvl_t
69 integer :: lvl = -1
70 type(space_t), pointer :: xh
71 type(dofmap_t), pointer :: dm_xh
72 type(gs_t), pointer :: gs_h
73 type(cheby_t) :: cheby
77 type(coef_t), pointer :: coef
78 type(bc_list_t) :: bclst
79 type(dirichlet_t) :: bc
80 type(field_t) :: r, w, z
81 end type phmg_lvl_t
82
83 type, public :: phmg_hrchy_t
84 type(phmg_lvl_t), allocatable :: lvl(:)
85 end type phmg_hrchy_t
86
87
88 type, public, extends(pc_t) :: phmg_t
89 type(tamg_solver_t) :: amg_solver
90 integer :: nlvls
91 type(phmg_hrchy_t) :: phmg_hrchy
92 class(ax_t), allocatable :: ax
93 type(interpolator_t), allocatable :: intrp(:)
94 type(mesh_t), pointer :: msh
95 contains
96 procedure, pass(this) :: init => phmg_init
97 procedure, pass(this) :: free => phmg_free
98 procedure, pass(this) :: solve => phmg_solve
99 procedure, pass(this) :: update => phmg_update
100 end type phmg_t
101
102contains
103
104 subroutine phmg_init(this, msh, Xh, coef, dof, gs_h, bclst)
105 class(phmg_t), intent(inout), target :: this
106 type(mesh_t), intent(inout), target :: msh
107 type(space_t), intent(inout), target :: Xh
108 type(coef_t), intent(in), target :: coef
109 type(dofmap_t), intent(in), target :: dof
110 type(gs_t), intent(inout), target :: gs_h
111 type(bc_list_t), intent(inout), target :: bclst
112 integer :: lx_crs, lx_mid
113 integer, allocatable :: lx_lvls(:)
114 integer :: n, i, j
115 class(bc_t), pointer :: bc_j
116 logical :: use_jacobi, use_cheby
117 use_jacobi = .false.
118 use_cheby = .true.
119
120 this%msh => msh
121
122 this%nlvls = xh%lx - 1
123
124 allocate(lx_lvls(0:this%nlvls - 1))
125 do i = 1, this%nlvls -1
126 lx_lvls(i) = xh%lx - i
127 end do
128
129 allocate(this%phmg_hrchy%lvl(0:this%nlvls - 1))
130
131 this%phmg_hrchy%lvl(0)%lvl = 0
132 this%phmg_hrchy%lvl(0)%Xh => xh
133 this%phmg_hrchy%lvl(0)%coef => coef
134 this%phmg_hrchy%lvl(0)%dm_Xh => dof
135 this%phmg_hrchy%lvl(0)%gs_h => gs_h
136
137 do i = 1, this%nlvls - 1
138 allocate(this%phmg_hrchy%lvl(i)%Xh)
139 allocate(this%phmg_hrchy%lvl(i)%dm_Xh)
140 allocate(this%phmg_hrchy%lvl(i)%gs_h)
141 allocate(this%phmg_hrchy%lvl(i)%coef)
142
143 this%phmg_hrchy%lvl(i)%lvl = i
144 call this%phmg_hrchy%lvl(i)%Xh%init(gll, lx_lvls(i), lx_lvls(i), &
145 lx_lvls(i))
146 call this%phmg_hrchy%lvl(i)%dm_Xh%init(msh, this%phmg_hrchy%lvl(i)%Xh)
147 call this%phmg_hrchy%lvl(i)%gs_h%init(this%phmg_hrchy%lvl(i)%dm_Xh)
148 call this%phmg_hrchy%lvl(i)%coef%init(this%phmg_hrchy%lvl(i)%gs_h)
149 end do
150
151 do i = 0, this%nlvls - 1
152 call this%phmg_hrchy%lvl(i)%r%init(this%phmg_hrchy%lvl(i)%dm_Xh)
153 call this%phmg_hrchy%lvl(i)%w%init(this%phmg_hrchy%lvl(i)%dm_Xh)
154 call this%phmg_hrchy%lvl(i)%z%init(this%phmg_hrchy%lvl(i)%dm_Xh)
155
156 if (use_cheby) then
157 if (neko_bcknd_device .eq. 1) then
158 call this%phmg_hrchy%lvl(i)%cheby_device%init(this%phmg_hrchy%lvl(i)%dm_Xh%size(), ksp_max_iter)
159 else
160 call this%phmg_hrchy%lvl(i)%cheby%init(this%phmg_hrchy%lvl(i)%dm_Xh%size(), ksp_max_iter)
161 end if
162 end if
163
164 if (use_jacobi) then
165 if (neko_bcknd_device .eq. 1) then
166 call this%phmg_hrchy%lvl(i)%device_jacobi%init(this%phmg_hrchy%lvl(i)%coef, &
167 this%phmg_hrchy%lvl(i)%dm_Xh, &
168 this%phmg_hrchy%lvl(i)%gs_h)
169 else
170 call this%phmg_hrchy%lvl(i)%jacobi%init(this%phmg_hrchy%lvl(i)%coef, &
171 this%phmg_hrchy%lvl(i)%dm_Xh, &
172 this%phmg_hrchy%lvl(i)%gs_h)
173 end if
174 end if
175
176 this%phmg_hrchy%lvl(i)%coef%ifh2 = coef%ifh2
177 call copy(this%phmg_hrchy%lvl(i)%coef%h1, coef%h1, &
178 this%phmg_hrchy%lvl(i)%dm_Xh%size())
179
180 call this%phmg_hrchy%lvl(i)%bc%init_base(this%phmg_hrchy%lvl(i)%coef)
181 if (bclst%size() .gt. 0 ) then
182 do j = 1, bclst%size()
183 bc_j => bclst%get(j)
184 call this%phmg_hrchy%lvl(i)%bc%mark_facets(bc_j%marked_facet)
185 end do
186 end if
187 call this%phmg_hrchy%lvl(i)%bc%finalize()
188 call this%phmg_hrchy%lvl(i)%bc%set_g(0.0_rp)
189 call this%phmg_hrchy%lvl(i)%bclst%init()
190 call this%phmg_hrchy%lvl(i)%bclst%append(this%phmg_hrchy%lvl(i)%bc)
191 end do
192
193 ! Create backend specific Ax operator
194 call ax_helm_factory(this%ax, full_formulation = .false.)
195
196 ! Interpolator Fine + mg levels
197 allocate(this%intrp(this%nlvls - 1))
198 do i = 1, this%nlvls -1
199 call this%intrp(i)%init(this%phmg_hrchy%lvl(i-1)%Xh, &
200 this%phmg_hrchy%lvl(i)%Xh)
201 end do
202
203 call this%amg_solver%init(this%ax, this%phmg_hrchy%lvl(this%nlvls -1)%Xh, &
204 this%phmg_hrchy%lvl(this%nlvls -1)%coef, this%msh, &
205 this%phmg_hrchy%lvl(this%nlvls-1)%gs_h, 4, &
206 this%phmg_hrchy%lvl(this%nlvls -1)%bclst, 1)
207
208 end subroutine phmg_init
209
210 subroutine phmg_free(this)
211 class(phmg_t), intent(inout) :: this
212 end subroutine phmg_free
213
214 subroutine phmg_solve(this, z, r, n)
215 class(phmg_t), intent(inout) :: this
216 integer, intent(in) :: n
217 real(kind=rp), dimension(n), intent(inout) :: z
218 real(kind=rp), dimension(n), intent(inout) :: r
219 type(c_ptr) :: z_d, r_d
220 type(ksp_monitor_t) :: ksp_results
221
222
223 associate( mglvl => this%phmg_hrchy%lvl)
224 if (neko_bcknd_device .eq. 1) then
225 z_d = device_get_ptr(z)
226 r_d = device_get_ptr(r)
227 !We should not work with the input
228 call device_copy(mglvl(0)%r%x_d, r_d, n)
229 call device_rzero(mglvl(0)%z%x_d, n)
230 call device_rzero(mglvl(0)%w%x_d, n)
231 call phmg_mg_cycle(mglvl(0)%z, mglvl(0)%r, mglvl(0)%w, 0, this%nlvls -1, &
232 mglvl, this%intrp, this%msh, this%Ax, this%amg_solver)
233 call device_copy(z_d, mglvl(0)%z%x_d, n)
234 else
235 !We should not work with the input
236 call copy(mglvl(0)%r%x, r, n)
237
238 mglvl(0)%z%x = 0.0_rp
239 mglvl(0)%w%x = 0.0_rp
240
241 call phmg_mg_cycle(mglvl(0)%z, mglvl(0)%r, mglvl(0)%w, 0, this%nlvls -1, &
242 mglvl, this%intrp, this%msh, this%Ax, this%amg_solver)
243
244 call copy(z, mglvl(0)%z%x, n)
245 end if
246 end associate
247
248 end subroutine phmg_solve
249
250 subroutine phmg_update(this)
251 class(phmg_t), intent(inout) :: this
252 end subroutine phmg_update
253
254
255 recursive subroutine phmg_mg_cycle(z, r, w, lvl, clvl, &
256 mg, intrp, msh, Ax, amg_solver)
257 type(ksp_monitor_t) :: ksp_results
258 integer :: lvl, clvl
259 type(phmg_lvl_t) :: mg(0:clvl)
260 type(interpolator_t) :: intrp(1:clvl)
261 type(tamg_solver_t), intent(inout) :: amg_solver
262 class(ax_t), intent(inout) :: ax
263 type(mesh_t), intent(inout) :: msh
264 type(field_t) :: z, r, w
265 integer :: i
266 real(kind=rp) :: val
267
268 call profiler_start_region('PHMG_cycle', 8)
272 call profiler_start_region('PHMG_PreSmooth', 9)
273 !!!!!!--------------------------------------------------------------------------------------------------------------
274 !call phmg_jacobi_smoother(z, r, w, mg(lvl), msh, Ax, mg(lvl)%dm_Xh%size(), lvl)
275 !!!!!!--------------------------------------------------------------------------------------------------------------
276 if (neko_bcknd_device .eq. 1) then
277 ksp_results = mg(lvl)%cheby_device%solve(ax, z, &
278 r%x, mg(lvl)%dm_Xh%size(), &
279 mg(lvl)%coef, mg(lvl)%bclst, &
280 mg(lvl)%gs_h, niter = 11)
281 else
282 ksp_results = mg(lvl)%cheby%solve(ax, z, &
283 r%x, mg(lvl)%dm_Xh%size(), &
284 mg(lvl)%coef, mg(lvl)%bclst, &
285 mg(lvl)%gs_h, niter = 15)
286 end if
287 call profiler_end_region('PHMG_PreSmooth', 9)
288
292 call ax%compute(w%x, z%x, mg(lvl)%coef, msh, mg(lvl)%Xh)
293 call mg(lvl)%gs_h%op(w%x, mg(lvl)%dm_Xh%size(), gs_op_add)
294 !call mg(lvl)%bclst%apply_scalar(w%x, mg(lvl)%dm_Xh%size())
295
296 if (neko_bcknd_device .eq. 1) then
297 call device_sub3(w%x_d, r%x_d, w%x_d, mg(lvl)%dm_Xh%size())
298 else
299 w%x = r%x - w%x
300 end if
301
305 if (neko_bcknd_device .eq. 1) then
306 call device_col2(w%x_d, mg(lvl)%coef%mult_d, mg(lvl)%dm_Xh%size())
307 else
308 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
309 end if
310
311 call profiler_start_region('PHMG_map_to_coarse', 9)
312 call intrp(lvl+1)%map(mg(lvl+1)%r%x, w%x, msh%nelv, mg(lvl+1)%Xh)
313 call profiler_end_region('PHMG_map_to_coarse', 9)
314
315 call mg(lvl+1)%gs_h%op(mg(lvl+1)%r%x, mg(lvl+1)%dm_Xh%size(), gs_op_add)
316
317 call mg(lvl+1)%bclst%apply_scalar( &
318 mg(lvl+1)%r%x, &
319 mg(lvl+1)%dm_Xh%size())
323 if (neko_bcknd_device .eq. 1) then
324 call device_rzero(mg(lvl+1)%z%x_d, mg(lvl+1)%dm_Xh%size())
325 else
326 mg(lvl+1)%z%x = 0.0_rp
327 end if
328 if (lvl+1 .eq. clvl) then
329 call profiler_start_region('PHMG_tAMG_coarse_grid', 10)
330
331 if (neko_bcknd_device .eq. 1) then
332 call amg_solver%device_solve(mg(lvl+1)%z%x, &
333 mg(lvl+1)%r%x, &
334 mg(lvl+1)%z%x_d, &
335 mg(lvl+1)%r%x_d, &
336 mg(lvl+1)%dm_Xh%size())
337 else
338 call amg_solver%solve(mg(lvl+1)%z%x, &
339 mg(lvl+1)%r%x, &
340 mg(lvl+1)%dm_Xh%size())
341 end if
342 call profiler_end_region('PHMG_tAMG_coarse_grid', 10)
343
344 call mg(lvl+1)%bclst%apply_scalar( &
345 mg(lvl+1)%z%x,&
346 mg(lvl+1)%dm_Xh%size())
347 else
348 call phmg_mg_cycle(mg(lvl+1)%z, mg(lvl+1)%r, mg(lvl+1)%w, lvl+1, &
349 clvl, mg, intrp, msh, ax, amg_solver)
350 end if
351
355 call profiler_start_region('PHMG_map_to_fine', 9)
356 call intrp(lvl+1)%map(w%x, mg(lvl+1)%z%x, msh%nelv, mg(lvl)%Xh)
357 call profiler_end_region('PHMG_map_to_fine', 9)
358
359 call mg(lvl)%gs_h%op(w%x, mg(lvl)%dm_Xh%size(), gs_op_add)
360
361 if (neko_bcknd_device .eq. 1) then
362 call device_col2(w%x_d, mg(lvl)%coef%mult_d, mg(lvl)%dm_Xh%size())
363 else
364 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
365 end if
366
367 call mg(lvl)%bclst%apply_scalar(w%x, mg(lvl)%dm_Xh%size())
368
372 if (neko_bcknd_device .eq. 1) then
373 call device_add2(z%x_d, w%x_d, mg(lvl)%dm_Xh%size())
374 else
375 z%x = z%x + w%x
376 end if
377
381 call profiler_start_region('PHMG_PostSmooth', 9)
382 !!!!!!--------------------------------------------------------------------------------------------------------------
383 !call phmg_jacobi_smoother(z, r, w, mg(lvl), msh, Ax, mg(lvl)%dm_Xh%size(), lvl)
384 !!!!!!--------------------------------------------------------------------------------------------------------------
385 if (neko_bcknd_device .eq. 1) then
386 ksp_results = mg(lvl)%cheby_device%solve(ax, z, &
387 r%x, mg(lvl)%dm_Xh%size(), &
388 mg(lvl)%coef, mg(lvl)%bclst, &
389 mg(lvl)%gs_h, niter = 11)
390 else
391 ksp_results = mg(lvl)%cheby%solve(ax, z, &
392 r%x, mg(lvl)%dm_Xh%size(), &
393 mg(lvl)%coef, mg(lvl)%bclst, &
394 mg(lvl)%gs_h, niter = 15)
395 end if
396 call profiler_end_region('PHMG_PostSmooth', 9)
397
398 call profiler_end_region('PHMG_cycle', 8)
399 end subroutine phmg_mg_cycle
400
401 subroutine phmg_jacobi_smoother(z, r, w, mg, msh, Ax, n, lvl)
402 type(phmg_lvl_t) :: mg
403 class(ax_t), intent(inout) :: Ax
404 type(mesh_t), intent(inout) :: msh
405 type(field_t), intent(inout) :: z, r, w
406 integer, intent(in) :: n, lvl
407 integer :: i, iblk, ni, niblk
408
409 ni = 1
410 niblk = 3
411
412 do i = 1, ni
413 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
414 call mg%gs_h%op(w%x, n, gs_op_add)
415 call mg%bclst%apply_scalar(w%x, n)
416 call device_sub3(w%x_d, r%x_d, w%x_d, n)
417
418 call mg%device_jacobi%solve(w%x, w%x, n)
419
420 call device_add2s2(z%x_d, w%x_d, 0.8_rp, n)
421
422 do iblk = 1, niblk
423 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
424 call device_invcol2(w%x_d, mg%coef%mult_d, n)
425 call device_sub3(w%x_d, r%x_d, w%x_d, n)
426 call mg%device_jacobi%solve(w%x, w%x, n)
427 call device_add2s2(z%x_d, w%x_d, 0.8_rp, n)
428 end do
429 end do
430 end subroutine phmg_jacobi_smoother
431
432
433 subroutine phmg_resid_monitor(z, r, w, mg, msh, Ax, lvl, typ)
434 integer :: lvl, typ
435 type(phmg_lvl_t) :: mg
436 class(ax_t), intent(inout) :: Ax
437 type(mesh_t), intent(inout) :: msh
438 type(field_t) :: z, r, w
439 real(kind=rp) :: val
440 character(len=LOG_SIZE) :: log_buf
441 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
442 call mg%gs_h%op(w%x, mg%dm_Xh%size(), gs_op_add)
443 call mg%bclst%apply_scalar(w%x, mg%dm_Xh%size())
444 call device_sub3(w%x_d, r%x_d, w%x_d, mg%dm_Xh%size())
445 val = device_glsc2(w%x_d, w%x_d, mg%dm_Xh%size())
446 if (typ .eq. 1) then
447 write(log_buf, '(A15,I4,F12.6)') 'PRESMOO - PRE', lvl, val
448 else if (typ .eq. 2) then
449 write(log_buf, '(A15,I4,F12.6)') 'PRESMOO -POST', lvl, val
450 else if (typ .eq. 3) then
451 write(log_buf, '(A15,I4,F12.6)') 'POSTSMOO- PRE', lvl, val
452 else if (typ .eq. 4) then
453 write(log_buf, '(A15,I4,F12.6)') 'POSTSMOO-POST', lvl, val
454 else if (typ .eq. 5) then
455 write(log_buf, '(A15,I4,F12.6)') 'TAMG - PRE', lvl, val
456 else if (typ .eq. 6) then
457 write(log_buf, '(A15,I4,F12.6)') 'TAMG -POST', lvl, val
458 else
459 write(log_buf, '(A15,I4,F12.6)') 'RESID', lvl, val
460 end if
461 call neko_log%message(log_buf)
462 end subroutine phmg_resid_monitor
463
464end module phmg
Return the device pointer for an associated Fortran array.
Definition device.F90:95
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Chebyshev preconditioner.
Chebyshev preconditioner.
Definition cheby.f90:34
Coefficients.
Definition coef.f90:34
Jacobi preconditioner accelerator backend.
subroutine, public device_add2(a_d, b_d, n)
Vector addition .
subroutine, public device_col2(a_d, b_d, n)
Vector multiplication .
subroutine, public device_rzero(a_d, n)
Zero a real vector.
subroutine, public device_add2s2(a_d, b_d, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
real(kind=rp) function, public device_glsc2(a_d, b_d, n)
Weighted inner product .
subroutine, public device_sub3(a_d, b_d, c_d, n)
Vector subtraction .
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
subroutine, public device_invcol2(a_d, b_d, n)
Vector division .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Gather-scatter.
Routines to interpolate between different spaces.
Jacobi preconditioner.
Definition pc_jacobi.f90:34
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
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:65
integer, parameter, public log_size
Definition log.f90:42
Definition math.f90:60
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:586
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:728
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:238
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Hybrid ph-multigrid preconditioner.
Definition phmg.f90:34
subroutine phmg_jacobi_smoother(z, r, w, mg, msh, ax, n, lvl)
Definition phmg.f90:402
subroutine phmg_solve(this, z, r, n)
Definition phmg.f90:215
subroutine phmg_init(this, msh, xh, coef, dof, gs_h, bclst)
Definition phmg.f90:105
subroutine phmg_update(this)
Definition phmg.f90:251
subroutine phmg_free(this)
Definition phmg.f90:211
subroutine phmg_resid_monitor(z, r, w, mg, msh, ax, lvl, typ)
Definition phmg.f90:434
recursive subroutine phmg_mg_cycle(z, r, w, lvl, clvl, mg, intrp, msh, ax, amg_solver)
Definition phmg.f90:257
Krylov preconditioner.
Definition precon.f90:34
Profiling interface.
Definition profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:109
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:48
Implements multigrid using the TreeAMG hierarchy structure. USE:
Utilities.
Definition utils.f90:35
Base type for a matrix-vector product providing .
Definition ax.f90:43
Base type for a boundary condition.
Definition bc.f90:57
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:47
Defines a Chebyshev preconditioner.
Definition cheby.f90:52
Defines a Chebyshev preconditioner.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:47
Interpolation between two space::space_t.
Defines a jacobi preconditioner.
Definition pc_jacobi.f90:45
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
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
The function space for the SEM solution fields.
Definition space.f90:62
Type for the TreeAMG solver.