Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
phmg.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 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
48 use jacobi, only : jacobi_t
49 use ax_product, only : ax_t, ax_helm_factory
52 use math, only : copy, col2, add2
53 use krylov, only : ksp_t, ksp_monitor_t, ksp_max_iter, &
54 krylov_solver_factory, krylov_solver_destroy
55 implicit none
56 private
57
58
59 type, private :: phmg_lvl_t
60 integer :: lvl = -1
61 type(space_t), pointer :: xh
62 type(dofmap_t), pointer :: dm_xh
63 type(gs_t), pointer :: gs_h
64 type(cheby_t) :: cheby
66 type(coef_t), pointer :: coef
67 type(bc_list_t) :: bclst
68 type(dirichlet_t) :: bc
69 type(field_t) :: r, w, z
70 end type phmg_lvl_t
71
72 type, public :: phmg_hrchy_t
73 type(phmg_lvl_t), allocatable :: lvl(:)
74 end type phmg_hrchy_t
75
76
77 type, public, extends(pc_t) :: phmg_t
78 type(tamg_solver_t) :: amg_solver
79 integer :: nlvls
80 type(phmg_hrchy_t) :: phmg_hrchy
81 class(ax_t), allocatable :: ax
82 type(interpolator_t), allocatable :: intrp(:)
83 type(mesh_t), pointer :: msh
84 contains
85 procedure, pass(this) :: init => phmg_init
86 procedure, pass(this) :: free => phmg_free
87 procedure, pass(this) :: solve => phmg_solve
88 procedure, pass(this) :: update => phmg_update
89 end type phmg_t
90
91contains
92
93 subroutine phmg_init(this, msh, Xh, coef, dof, gs_h, bclst)
94 class(phmg_t), intent(inout), target :: this
95 type(mesh_t), intent(inout), target :: msh
96 type(space_t), intent(inout), target :: Xh
97 type(coef_t), intent(in), target :: coef
98 type(dofmap_t), intent(in), target :: dof
99 type(gs_t), intent(inout), target :: gs_h
100 type(bc_list_t), intent(inout), target :: bclst
101 integer :: lx_crs, lx_mid
102 integer, allocatable :: lx_lvls(:)
103 integer :: n, i, j
104
105 this%msh => msh
106
107 this%nlvls = xh%lx - 1
108
109 allocate(lx_lvls(0:this%nlvls - 1))
110 lx_lvls(0) = xh%lx
111 do i = 1, this%nlvls -1
112 lx_lvls(i) = xh%lx - i
113 end do
114
115 allocate(this%phmg_hrchy%lvl(0:this%nlvls - 1))
116
117 this%phmg_hrchy%lvl(0)%lvl = 0
118 this%phmg_hrchy%lvl(0)%Xh => xh
119 this%phmg_hrchy%lvl(0)%coef => coef
120 this%phmg_hrchy%lvl(0)%dm_Xh => dof
121 this%phmg_hrchy%lvl(0)%gs_h => gs_h
122
123 do i = 1, this%nlvls - 1
124 allocate(this%phmg_hrchy%lvl(i)%Xh)
125 allocate(this%phmg_hrchy%lvl(i)%dm_Xh)
126 allocate(this%phmg_hrchy%lvl(i)%gs_h)
127 allocate(this%phmg_hrchy%lvl(i)%coef)
128
129 this%phmg_hrchy%lvl(i)%lvl = i
130 call this%phmg_hrchy%lvl(i)%Xh%init(gll, lx_lvls(i), lx_lvls(i), lx_lvls(i))
131 call this%phmg_hrchy%lvl(i)%dm_Xh%init(msh, this%phmg_hrchy%lvl(i)%Xh)
132 call this%phmg_hrchy%lvl(i)%gs_h%init(this%phmg_hrchy%lvl(i)%dm_Xh)
133 call this%phmg_hrchy%lvl(i)%coef%init(this%phmg_hrchy%lvl(i)%gs_h)
134 end do
135
136 do i = 0, this%nlvls - 1
137 call this%phmg_hrchy%lvl(i)%r%init(this%phmg_hrchy%lvl(i)%dm_Xh)
138 call this%phmg_hrchy%lvl(i)%w%init(this%phmg_hrchy%lvl(i)%dm_Xh)
139 call this%phmg_hrchy%lvl(i)%z%init(this%phmg_hrchy%lvl(i)%dm_Xh)
140
141 call this%phmg_hrchy%lvl(i)%cheby%init(this%phmg_hrchy%lvl(i)%dm_Xh%size(), ksp_max_iter)
142
143 call this%phmg_hrchy%lvl(i)%jacobi%init(this%phmg_hrchy%lvl(i)%coef, &
144 this%phmg_hrchy%lvl(i)%dm_Xh, &
145 this%phmg_hrchy%lvl(i)%gs_h)
146
147 this%phmg_hrchy%lvl(i)%coef%ifh2 = coef%ifh2
148 call copy(this%phmg_hrchy%lvl(i)%coef%h1, coef%h1, &
149 this%phmg_hrchy%lvl(i)%dm_Xh%size())
150
151 call this%phmg_hrchy%lvl(i)%bc%init_base(this%phmg_hrchy%lvl(i)%coef)
152 if (bclst%size .gt. 0 ) then
153 do j = 1, bclst%size
154 call this%phmg_hrchy%lvl(i)%bc%mark_facets(&
155 bclst%items(j)%ptr%marked_facet)
156 end do
157 end if
158 call this%phmg_hrchy%lvl(i)%bc%finalize()
159 call this%phmg_hrchy%lvl(i)%bc%set_g(0.0_rp)
160 call this%phmg_hrchy%lvl(i)%bclst%init()
161 call this%phmg_hrchy%lvl(i)%bclst%append(this%phmg_hrchy%lvl(i)%bc)
162 end do
163
164 ! Create backend specific Ax operator
165 call ax_helm_factory(this%ax, full_formulation = .false.)
166
167 ! Interpolator Fine + mg levels
168 allocate(this%intrp(this%nlvls - 1))
169 do i = 1, this%nlvls -1
170 call this%intrp(i)%init(this%phmg_hrchy%lvl(i-1)%Xh, this%phmg_hrchy%lvl(i)%Xh)
171 end do
172
173 call this%amg_solver%init(this%ax, this%phmg_hrchy%lvl(this%nlvls -1)%Xh, &
174 this%phmg_hrchy%lvl(this%nlvls -1)%coef, this%msh, &
175 this%phmg_hrchy%lvl(this%nlvls-1)%gs_h, 4, &
176 this%phmg_hrchy%lvl(this%nlvls -1)%bclst, 1)
177
178 end subroutine phmg_init
179
180 subroutine phmg_free(this)
181 class(phmg_t), intent(inout) :: this
182
183
184 end subroutine phmg_free
185
186 subroutine phmg_solve(this, z, r, n)
187 class(phmg_t), intent(inout) :: this
188 integer, intent(in) :: n
189 real(kind=rp), dimension(n), intent(inout) :: z
190 real(kind=rp), dimension(n), intent(inout) :: r
191 type(ksp_monitor_t) :: ksp_results
192
193
194 associate( mglvl => this%phmg_hrchy%lvl)
195 !We should not work with the input
196 call copy(mglvl(0)%r%x, r, n)
197
198 mglvl(0)%z%x = 0.0_rp
199 mglvl(0)%w%x = 0.0_rp
200
201 call phmg_mg_cycle(mglvl(0)%z, mglvl(0)%r, mglvl(0)%w, 0, this%nlvls -1, &
202 mglvl, this%intrp, this%msh, this%Ax, this%amg_solver)
203
204 call copy(z, mglvl(0)%z%x, n)
205
206 end associate
207
208 end subroutine phmg_solve
209
210 subroutine phmg_update(this)
211 class(phmg_t), intent(inout) :: this
212 end subroutine phmg_update
213
214
215 recursive subroutine phmg_mg_cycle(z, r, w, lvl, clvl, &
216 mg, intrp, msh, Ax, amg_solver)
217 type(ksp_monitor_t) :: ksp_results
218 integer :: lvl, clvl
219 type(phmg_lvl_t) :: mg(0:clvl)
220 type(interpolator_t) :: intrp(1:clvl)
221 type(tamg_solver_t), intent(inout) :: amg_solver
222 class(ax_t), intent(inout) :: ax
223 type(mesh_t), intent(inout) :: msh
224 type(field_t) :: z, r, w
225 integer :: i
226
227
228 ksp_results = mg(lvl)%cheby%solve(ax, z, &
229 r%x, mg(lvl)%dm_Xh%size(), &
230 mg(lvl)%coef, mg(lvl)%bclst, &
231 mg(lvl)%gs_h, niter = 15)
232
233 call ax%compute(w%x, z%x, mg(lvl)%coef, msh, mg(lvl)%Xh)
234
235 w%x = r%x - w%x
236
237 call intrp(lvl+1)%map(mg(lvl+1)%r%x, w%x, msh%nelv, mg(lvl+1)%Xh)
238
239 call mg(lvl+1)%gs_h%op(mg(lvl+1)%r%x, mg(lvl+1)%dm_Xh%size(), gs_op_add)
240
241 call col2(mg(lvl+1)%r%x, mg(lvl+1)%coef%mult, mg(lvl+1)%dm_Xh%size())
242
243
244 call mg(lvl+1)%bclst%apply_scalar( &
245 mg(lvl+1)%r%x, &
246 mg(lvl+1)%dm_Xh%size())
247
248 mg(lvl+1)%z%x = 0.0_rp
249 if (lvl+1 .eq. clvl) then
250
251 call amg_solver%solve(mg(lvl+1)%z%x, &
252 mg(lvl+1)%r%x, &
253 mg(lvl+1)%dm_Xh%size())
254
255 call mg(lvl+1)%bclst%apply_scalar( &
256 mg(lvl+1)%z%x,&
257 mg(lvl+1)%dm_Xh%size())
258 else
259 call phmg_mg_cycle(mg(lvl+1)%z, mg(lvl+1)%r, mg(lvl+1)%w, lvl+1, &
260 clvl, mg, intrp, msh, ax, amg_solver)
261 end if
262
263 call intrp(lvl+1)%map(w%x, mg(lvl+1)%z%x, msh%nelv, mg(lvl)%Xh)
264
265 call mg(lvl)%gs_h%op(w%x, mg(lvl)%dm_Xh%size(), gs_op_add)
266
267 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
268
269 call mg(lvl)%bclst%apply_scalar(w%x, mg(lvl)%dm_Xh%size())
270
271 z%x = z%x + w%x
272
273 ksp_results = mg(lvl)%cheby%solve(ax, z, &
274 r%x, mg(lvl)%dm_Xh%size(), &
275 mg(lvl)%coef, mg(lvl)%bclst, &
276 mg(lvl)%gs_h, niter = 15)
277
278
279
280 end subroutine phmg_mg_cycle
281
282end module phmg
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.
Definition cheby.f90:34
Coefficients.
Definition coef.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
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
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Hybrid ph-multigrid preconditioner.
Definition phmg.f90:34
subroutine phmg_solve(this, z, r, n)
Definition phmg.f90:187
subroutine phmg_init(this, msh, xh, coef, dof, gs_h, bclst)
Definition phmg.f90:94
subroutine phmg_update(this)
Definition phmg.f90:211
subroutine phmg_free(this)
Definition phmg.f90:181
recursive subroutine phmg_mg_cycle(z, r, w, lvl, clvl, mg, intrp, msh, ax, amg_solver)
Definition phmg.f90:217
Krylov preconditioner.
Definition precon.f90:34
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:51
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:46
Defines a Chebyshev preconditioner.
Definition cheby.f90:52
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:44
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.