139 subroutine hsmg_init(this, msh, Xh, coef, dof, gs_h, bclst, crs_pctype)
140 class(
hsmg_t),
intent(inout),
target :: this
141 type(
mesh_t),
intent(inout),
target :: msh
142 type(
space_t),
intent(inout),
target :: Xh
143 type(
coef_t),
intent(in),
target :: coef
144 type(
dofmap_t),
intent(in),
target :: dof
145 type(
gs_t),
intent(inout),
target :: gs_h
146 type(
bc_list_t),
intent(inout),
target :: bclst
147 character(len=*),
optional :: crs_pctype
149 integer :: lx_crs, lx_mid
154 if (xh%lx .lt. 5)
then
155 lx_mid =
max(xh%lx-1,3)
157 if(xh%lx .le. 2)
then
158 call neko_error(
'Polynomial order < 2 not supported for hsmg precon')
165 allocate(this%grids(this%nlvls))
166 allocate(this%w(dof%size()))
167 allocate(this%r(dof%size()))
171 call msh%all_deformed()
174 call this%e%init(dof,
'work array')
175 call this%wf%init(dof,
'work 2')
177 call this%Xh_crs%init(
gll, lx_crs, lx_crs, lx_crs)
178 call this%dm_crs%init(msh, this%Xh_crs)
179 call this%gs_crs%init(this%dm_crs)
180 call this%e_crs%init(this%dm_crs,
'work crs')
181 call this%c_crs%init(this%gs_crs)
183 call this%Xh_mg%init(
gll, lx_mid, lx_mid, lx_mid)
184 call this%dm_mg%init(msh, this%Xh_mg)
185 call this%gs_mg%init(this%dm_mg)
186 call this%e_mg%init(this%dm_mg,
'work midl')
187 call this%c_mg%init(this%gs_mg)
190 call ax_helm_factory(this%ax, full_formulation = .false.)
193 call precon_factory(this%pc_crs,
'jacobi')
195 call this%bc_crs%init_base(this%c_crs)
196 call this%bc_mg%init_base(this%c_mg)
197 call this%bc_reg%init_base(coef)
198 if (bclst%size .gt. 0)
then
200 call this%bc_reg%mark_facets(bclst%items(i)%ptr%marked_facet)
201 call this%bc_crs%mark_facets(bclst%items(i)%ptr%marked_facet)
202 call this%bc_mg%mark_facets(bclst%items(i)%ptr%marked_facet)
205 call this%bc_reg%finalize()
206 call this%bc_reg%set_g(
real(0d0,
rp))
207 call this%bclst_reg%init()
208 call this%bclst_reg%append(this%bc_reg)
210 call this%bc_crs%finalize()
211 call this%bc_crs%set_g(
real(0d0,
rp))
212 call this%bclst_crs%init()
213 call this%bclst_crs%append(this%bc_crs)
215 call this%bc_mg%finalize()
216 call this%bc_mg%set_g(0.0_rp)
217 call this%bclst_mg%init()
218 call this%bclst_mg%append(this%bc_mg)
220 call this%schwarz%init(xh, dof, gs_h, this%bclst_reg, msh)
221 call this%schwarz_mg%init(this%Xh_mg, this%dm_mg, this%gs_mg,&
224 call this%interp_fine_mid%init(xh, this%Xh_mg)
225 call this%interp_mid_crs%init(this%Xh_mg, this%Xh_crs)
227 call hsmg_fill_grid(dof, gs_h, xh, coef, this%bclst_reg, this%schwarz, &
228 this%e, this%grids, 3)
229 call hsmg_fill_grid(this%dm_mg, this%gs_mg, this%Xh_mg, this%c_mg, &
230 this%bclst_mg, this%schwarz_mg, this%e_mg, &
233 this%c_crs, this%bclst_crs, this%schwarz_crs, &
234 this%e_crs, this%grids, 1)
242 select type(pc => this%pc_crs)
244 call pc%init(this%c_crs, this%dm_crs, this%gs_crs)
246 call pc%init(this%c_crs, this%dm_crs, this%gs_crs)
248 call pc%init(this%c_crs, this%dm_crs, this%gs_crs)
255 if (
present(crs_pctype))
then
256 if (trim(crs_pctype) .eq.
'tamg')
then
258 call neko_error(
'Tree-amg only supported for CPU')
261 allocate(this%amg_solver)
263 call this%amg_solver%init(this%ax, this%grids(1)%e%Xh, &
264 this%grids(1)%coef, this%msh, this%grids(1)%gs_h, 4, &
265 this%grids(1)%bclst, 1)
267 call krylov_solver_factory(this%crs_solver, &
268 this%dm_crs%size(), trim(crs_pctype),
ksp_max_iter, m = this%pc_crs)
271 call krylov_solver_factory(this%crs_solver, &
272 this%dm_crs%size(),
'cg',
ksp_max_iter, m = this%pc_crs)
362 integer,
intent(in) :: n
363 class(
hsmg_t),
intent(inout) :: this
364 real(kind=
rp),
dimension(n),
intent(inout) :: z
365 real(kind=
rp),
dimension(n),
intent(inout) :: r
366 type(c_ptr) :: z_d, r_d
368 integer :: thrdid, nthrds
376 call this%bclst_reg%apply_scalar(r, n)
380 call device_col2(this%r_d, this%grids(3)%coef%mult_d, &
381 this%grids(3)%dof%size())
383 call this%interp_fine_mid%map(this%e%x, this%r, &
384 this%msh%nelv, this%grids(2)%Xh)
385 call this%grids(2)%gs_h%op(this%e%x, &
386 this%grids(2)%dof%size(), gs_op_add, this%gs_event)
389 call this%bclst_reg%apply_scalar(r, n)
390 call device_copy(this%w_d, this%e%x_d, this%grids(2)%dof%size())
391 call this%bclst_mg%apply_scalar(this%w, this%grids(2)%dof%size())
393 call device_col2(this%w_d, this%grids(2)%coef%mult_d, &
394 this%grids(2)%dof%size())
396 call this%interp_mid_crs%map(this%wf%x, this%w, this%msh%nelv, &
399 call device_copy(this%w_d, this%e%x_d, this%grids(2)%dof%size())
400 call this%bclst_mg%apply_scalar(this%w, this%grids(2)%dof%size())
409 if (thrdid .eq. 0)
then
411 call this%grids(3)%schwarz%compute(z, this%r)
412 call this%grids(2)%schwarz%compute(this%grids(2)%e%x, this%w)
415 if (nthrds .eq. 1 .or. thrdid .eq. 1)
then
417 call this%grids(1)%gs_h%op(this%wf%x, &
418 this%grids(1)%dof%size(), gs_op_add, this%gs_event)
420 call this%grids(1)%bclst%apply_scalar(this%wf%x, &
421 this%grids(1)%dof%size())
423 crs_info = this%crs_solver%solve(this%Ax, this%grids(1)%e, &
425 this%grids(1)%dof%size(), &
426 this%grids(1)%coef, &
427 this%grids(1)%bclst, &
428 this%grids(1)%gs_h, this%niter)
430 call this%grids(1)%bclst%apply_scalar(this%grids(1)%e%x,&
431 this%grids(1)%dof%size())
436 call this%interp_mid_crs%map(this%w, this%grids(1)%e%x, &
437 this%msh%nelv, this%grids(2)%Xh)
438 call device_add2(this%grids(2)%e%x_d, this%w_d, this%grids(2)%dof%size())
440 call this%interp_fine_mid%map(this%w, this%grids(2)%e%x, &
441 this%msh%nelv, this%grids(3)%Xh)
442 call device_add2(z_d, this%w_d, this%grids(3)%dof%size())
443 call this%grids(3)%gs_h%op(z, this%grids(3)%dof%size(), &
444 gs_op_add, this%gs_event)
447 this%grids(3)%dof%size())
450 call copy(this%r, r, n)
453 call this%grids(3)%schwarz%compute(z, this%r)
455 call col2(this%r, this%grids(3)%coef%mult, &
456 this%grids(3)%dof%size())
458 call this%interp_fine_mid%map(this%w, this%r, &
459 this%msh%nelv, this%grids(2)%Xh)
460 call this%grids(2)%gs_h%op(this%w, this%grids(2)%dof%size(), gs_op_add)
462 call this%grids(2)%schwarz%compute(this%grids(2)%e%x, this%w)
463 call col2(this%w, this%grids(2)%coef%mult, this%grids(2)%dof%size())
465 call this%interp_mid_crs%map(this%r, this%w, &
466 this%msh%nelv, this%grids(1)%Xh)
469 call this%grids(1)%gs_h%op(this%r, this%grids(1)%dof%size(), gs_op_add)
470 call this%grids(1)%bclst%apply(this%r, this%grids(1)%dof%size())
473 if (
allocated(this%amg_solver))
then
474 call this%amg_solver%solve(this%grids(1)%e%x, this%r, this%grids(1)%dof%size())
476 crs_info = this%crs_solver%solve(this%Ax, this%grids(1)%e, this%r, &
477 this%grids(1)%dof%size(), &
478 this%grids(1)%coef, &
479 this%grids(1)%bclst, &
480 this%grids(1)%gs_h, this%niter)
484 call this%grids(1)%bclst%apply_scalar(this%grids(1)%e%x, &
485 this%grids(1)%dof%size())
488 call this%interp_mid_crs%map(this%w, this%grids(1)%e%x, &
489 this%msh%nelv, this%grids(2)%Xh)
490 call add2(this%grids(2)%e%x, this%w, this%grids(2)%dof%size())
492 call this%interp_fine_mid%map(this%w, this%grids(2)%e%x, &
493 this%msh%nelv, this%grids(3)%Xh)
494 call add2(z, this%w, this%grids(3)%dof%size())
495 call this%grids(3)%gs_h%op(z, this%grids(3)%dof%size(), gs_op_add)
496 call col2(z, this%grids(3)%coef%mult, this%grids(3)%dof%size())