Neko  0.9.99
A portable framework for high-order spectral element flow simulations
pc_hsmg.f90
Go to the documentation of this file.
1 ! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2 !
3 ! The UChicago Argonne, LLC as Operator of Argonne National
4 ! Laboratory holds copyright in the Software. The copyright holder
5 ! reserves all rights except those expressly granted to licensees,
6 ! and U.S. Government license rights.
7 !
8 ! Redistribution and use in source and binary forms, with or without
9 ! modification, are permitted provided that the following conditions
10 ! are met:
11 !
12 ! 1. Redistributions of source code must retain the above copyright
13 ! notice, this list of conditions and the disclaimer below.
14 !
15 ! 2. Redistributions in binary form must reproduce the above copyright
16 ! notice, this list of conditions and the disclaimer (as noted below)
17 ! in the documentation and/or other materials provided with the
18 ! distribution.
19 !
20 ! 3. Neither the name of ANL nor the names of its contributors
21 ! may be used to endorse or promote products derived from this software
22 ! without specific prior written permission.
23 !
24 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
28 ! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
29 ! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
30 ! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
31 ! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
32 ! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
33 ! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
34 ! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
35 ! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 !
37 ! Additional BSD Notice
38 ! ---------------------
39 ! 1. This notice is required to be provided under our contract with
40 ! the U.S. Department of Energy (DOE). This work was produced at
41 ! Argonne National Laboratory under Contract
42 ! No. DE-AC02-06CH11357 with the DOE.
43 !
44 ! 2. Neither the United States Government nor UCHICAGO ARGONNE,
45 ! LLC nor any of their employees, makes any warranty,
46 ! express or implied, or assumes any liability or responsibility for the
47 ! accuracy, completeness, or usefulness of any information, apparatus,
48 ! product, or process disclosed, or represents that its use would not
49 ! infringe privately-owned rights.
50 !
51 ! 3. Also, reference herein to any specific commercial products, process,
52 ! or services by trade name, trademark, manufacturer or otherwise does
53 ! not necessarily constitute or imply its endorsement, recommendation,
54 ! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
55 ! The views and opinions of authors expressed
56 ! herein do not necessarily state or reflect those of the United States
57 ! Government or UCHICAGO ARGONNE, LLC, and shall
58 ! not be used for advertising or product endorsement purposes.
59 !
61 module hsmg
62  use neko_config, only : neko_bcknd_device
63  use num_types, only : rp
64  use math, only : copy, col2, add2
65  use utils, only : neko_error
66  use precon, only : pc_t, precon_factory, precon_destroy
67  use ax_product, only : ax_t, ax_helm_factory
68  use gather_scatter, only : gs_t, gs_op_add
69  use interpolation, only : interpolator_t
72  use dirichlet, only : dirichlet_t
73  use schwarz, only : schwarz_t
74  use jacobi, only : jacobi_t
75  use sx_jacobi, only : sx_jacobi_t
76  use device_jacobi, only : device_jacobi_t
77  use device
80  use space, only : space_t, gll
81  use dofmap, only : dofmap_t
82  use field, only : field_t
83  use coefs, only : coef_t
84  use mesh, only : mesh_t
85  use krylov, only : ksp_t, ksp_monitor_t, ksp_max_iter, &
86  krylov_solver_factory, krylov_solver_destroy
87  !$ use omp_lib
88  implicit none
89  private
90 
91  !Struct to arrange our multigridlevels
92  type, private :: multigrid_t
93  type(dofmap_t), pointer :: dof
94  type(gs_t), pointer :: gs_h
95  type(space_t), pointer :: xh
96  type(coef_t), pointer :: coef
97  type(bc_list_t), pointer :: bclst
98  type(schwarz_t), pointer :: schwarz
99  type(field_t), pointer :: e
100  end type multigrid_t
101 
102  type, public, extends(pc_t) :: hsmg_t
103  type(mesh_t), pointer :: msh
104  integer :: nlvls
105  type(multigrid_t), allocatable :: grids(:)
106  type(gs_t) :: gs_crs, gs_mg
107  type(space_t) :: xh_crs, xh_mg
108  type(dofmap_t) :: dm_crs, dm_mg
109  type(coef_t) :: c_crs, c_mg
110  type(dirichlet_t) :: bc_crs, bc_mg, bc_reg
111  type(bc_list_t) :: bclst_crs, bclst_mg, bclst_reg
112  type(schwarz_t) :: schwarz, schwarz_mg, schwarz_crs
113  type(field_t) :: e, e_mg, e_crs
114  type(field_t) :: wf
115  class(ksp_t), allocatable :: crs_solver
116  integer :: niter = 10
117  class(pc_t), allocatable :: pc_crs
118  class(ax_t), allocatable :: ax
119  real(kind=rp), allocatable :: r(:)
120  type(interpolator_t) :: interp_fine_mid
121  type(interpolator_t) :: interp_mid_crs
122  real(kind=rp), allocatable :: w(:)
123  type(c_ptr) :: w_d = c_null_ptr
124  type(c_ptr) :: r_d = c_null_ptr
125  type(c_ptr) :: hsmg_event
126  type(c_ptr) :: gs_event
127  contains
128  procedure, pass(this) :: init => hsmg_init
129  procedure, pass(this) :: free => hsmg_free
130  procedure, pass(this) :: solve => hsmg_solve
131  procedure, pass(this) :: update => hsmg_set_h
132  end type hsmg_t
133 
134 contains
135 
137  subroutine hsmg_init(this, msh, Xh, coef, dof, gs_h, bclst, crs_pctype)
138  class(hsmg_t), intent(inout), target :: this
139  type(mesh_t), intent(inout), target :: msh
140  type(space_t), intent(inout), target :: Xh
141  type(coef_t), intent(inout), target :: coef
142  type(dofmap_t), intent(inout), target :: dof
143  type(gs_t), intent(inout), target :: gs_h
144  type(bc_list_t), intent(inout), target :: bclst
145  character(len=*), optional :: crs_pctype
146  integer :: n, i
147  integer :: lx_crs, lx_mid
148 
149  call this%free()
150  this%nlvls = 3
151  lx_crs = 2
152  if (xh%lx .lt. 5) then
153  lx_mid = max(xh%lx-1,3)
154 
155  if(xh%lx .le. 2) then
156  call neko_error('Polynomial order < 2 not supported for hsmg precon')
157  end if
158 
159  else
160  lx_mid = 4
161  end if
162  this%msh => msh
163  allocate(this%grids(this%nlvls))
164  allocate(this%w(dof%size()))
165  allocate(this%r(dof%size()))
166 
167 
168  ! Compute all elements as if they are deformed
169  call msh%all_deformed()
170 
171  n = dof%size()
172  call this%e%init(dof, 'work array')
173  call this%wf%init(dof, 'work 2')
174 
175  call this%Xh_crs%init(gll, lx_crs, lx_crs, lx_crs)
176  call this%dm_crs%init(msh, this%Xh_crs)
177  call this%gs_crs%init(this%dm_crs)
178  call this%e_crs%init(this%dm_crs, 'work crs')
179  call this%c_crs%init(this%gs_crs)
180 
181  call this%Xh_mg%init(gll, lx_mid, lx_mid, lx_mid)
182  call this%dm_mg%init(msh, this%Xh_mg)
183  call this%gs_mg%init(this%dm_mg)
184  call this%e_mg%init(this%dm_mg, 'work midl')
185  call this%c_mg%init(this%gs_mg)
186 
187  ! Create backend specific Ax operator
188  call ax_helm_factory(this%ax, full_formulation = .false.)
189 
190  ! Create a backend specific preconditioner
191  call precon_factory(this%pc_crs, 'jacobi')
192 
193  ! Create a backend specific krylov solver
194  if (present(crs_pctype)) then
195  call krylov_solver_factory(this%crs_solver, &
196  this%dm_crs%size(), trim(crs_pctype), ksp_max_iter, m = this%pc_crs)
197  else
198  call krylov_solver_factory(this%crs_solver, &
199  this%dm_crs%size(), 'cg', ksp_max_iter, m = this%pc_crs)
200  end if
201 
202  call this%bc_crs%init_base(this%c_crs)
203  call this%bc_mg%init_base(this%c_mg)
204  call this%bc_reg%init_base(coef)
205  if (bclst%n .gt. 0) then
206  do i = 1, bclst%n
207  call this%bc_reg%mark_facets(bclst%bc(i)%bcp%marked_facet)
208  call this%bc_crs%mark_facets(bclst%bc(i)%bcp%marked_facet)
209  call this%bc_mg%mark_facets(bclst%bc(i)%bcp%marked_facet)
210  end do
211  end if
212  call this%bc_reg%finalize()
213  call this%bc_reg%set_g(real(0d0, rp))
214  call bc_list_init(this%bclst_reg)
215  call bc_list_add(this%bclst_reg, this%bc_reg)
216 
217  call this%bc_crs%finalize()
218  call this%bc_crs%set_g(real(0d0, rp))
219  call bc_list_init(this%bclst_crs)
220  call bc_list_add(this%bclst_crs, this%bc_crs)
221 
222 
223  call this%bc_mg%finalize()
224  call this%bc_mg%set_g(0.0_rp)
225  call bc_list_init(this%bclst_mg)
226  call bc_list_add(this%bclst_mg, this%bc_mg)
227 
228  call this%schwarz%init(xh, dof, gs_h, this%bclst_reg, msh)
229  call this%schwarz_mg%init(this%Xh_mg, this%dm_mg, this%gs_mg,&
230  this%bclst_mg, msh)
231 
232  call this%interp_fine_mid%init(xh, this%Xh_mg)
233  call this%interp_mid_crs%init(this%Xh_mg, this%Xh_crs)
234 
235  call hsmg_fill_grid(dof, gs_h, xh, coef, this%bclst_reg, this%schwarz, &
236  this%e, this%grids, 3)
237  call hsmg_fill_grid(this%dm_mg, this%gs_mg, this%Xh_mg, this%c_mg, &
238  this%bclst_mg, this%schwarz_mg, this%e_mg, &
239  this%grids, 2)
240  call hsmg_fill_grid(this%dm_crs, this%gs_crs, this%Xh_crs, &
241  this%c_crs, this%bclst_crs, this%schwarz_crs, &
242  this%e_crs, this%grids, 1)
243 
244  call hsmg_set_h(this)
245  if (neko_bcknd_device .eq. 1) then
246  call device_map(this%w, this%w_d, n)
247  call device_map(this%r, this%r_d, n)
248  end if
249 
250  select type(pc => this%pc_crs)
251  type is (jacobi_t)
252  call pc%init(this%c_crs, this%dm_crs, this%gs_crs)
253  type is (sx_jacobi_t)
254  call pc%init(this%c_crs, this%dm_crs, this%gs_crs)
255  type is (device_jacobi_t)
256  call pc%init(this%c_crs, this%dm_crs, this%gs_crs)
257  end select
258 
259  call device_event_create(this%hsmg_event, 2)
260  call device_event_create(this%gs_event, 2)
261  end subroutine hsmg_init
262 
263  subroutine hsmg_set_h(this)
264  class(hsmg_t), intent(inout) :: this
265 ! integer :: i
266  !Yeah I dont really know what to do here. For incompressible flow not much happens
267  this%grids(1)%coef%ifh2 = .false.
268  call copy(this%grids(1)%coef%h1, this%grids(3)%coef%h1, &
269  this%grids(1)%dof%size())
270  if (neko_bcknd_device .eq. 1) then
271  call device_copy(this%grids(1)%coef%h1_d, this%grids(3)%coef%h1_d, &
272  this%grids(1)%dof%size())
273  end if
274  end subroutine hsmg_set_h
275 
276 
277  subroutine hsmg_fill_grid(dof, gs_h, Xh, coef, bclst, schwarz, e, grids, l)
278  type(dofmap_t), target, intent(in) :: dof
279  type(gs_t), target, intent(in) :: gs_h
280  type(space_t), target, intent(in) :: Xh
281  type(coef_t), target, intent(in) :: coef
282  type(bc_list_t), target, intent(in) :: bclst
283  type(schwarz_t), target, intent(in) :: schwarz
284  type(field_t), target, intent(in) :: e
285  integer, intent(in) :: l
286  type(multigrid_t), intent(inout), dimension(l) :: grids
287 
288 
289  grids(l)%dof => dof
290  grids(l)%gs_h => gs_h
291  grids(l)%Xh => xh
292  grids(l)%coef => coef
293  grids(l)%bclst => bclst
294  grids(l)%schwarz => schwarz
295  grids(l)%e => e
296 
297  end subroutine hsmg_fill_grid
298 
299  subroutine hsmg_free(this)
300  class(hsmg_t), intent(inout) :: this
301 
302  if (allocated(this%ax)) then
303  deallocate(this%ax)
304  end if
305 
306  if (allocated(this%grids)) then
307  deallocate(this%grids)
308  end if
309 
310  if (allocated(this%w)) then
311  deallocate(this%w)
312  end if
313 
314  if (allocated(this%r)) then
315  deallocate(this%r)
316  end if
317 
318  call this%schwarz%free()
319  call this%schwarz_mg%free()
320 
321  call this%c_crs%free()
322  call this%c_mg%free()
323  call this%e%free()
324  call this%e_mg%free()
325  call this%e_crs%free()
326 
327  call this%gs_crs%free()
328  call this%gs_mg%free()
329  call this%interp_mid_crs%free()
330  call this%interp_fine_mid%free()
331 
332  if (allocated(this%crs_solver)) then
333  call krylov_solver_destroy(this%crs_solver)
334  deallocate(this%crs_solver)
335  end if
336 
337  if (allocated(this%pc_crs)) then
338  call precon_destroy(this%pc_crs)
339  end if
340 
341  end subroutine hsmg_free
342 
344  subroutine hsmg_solve(this, z, r, n)
345  integer, intent(in) :: n
346  class(hsmg_t), intent(inout) :: this
347  real(kind=rp), dimension(n), intent(inout) :: z
348  real(kind=rp), dimension(n), intent(inout) :: r
349  type(c_ptr) :: z_d, r_d
350  type(ksp_monitor_t) :: crs_info
351  integer :: thrdid, nthrds
352 
353  call profiler_start_region('HSMG_solve', 8)
354  if (neko_bcknd_device .eq. 1) then
355  z_d = device_get_ptr(z)
356  r_d = device_get_ptr(r)
357  !We should not work with the input
358  call device_copy(this%r_d, r_d, n)
359  call bc_list_apply_scalar(this%bclst_reg, r, n)
360 
361  !OVERLAPPING Schwarz exchange and solve
362  !! DOWNWARD Leg of V-cycle, we are pretty hardcoded here but w/e
363  call device_col2(this%r_d, this%grids(3)%coef%mult_d, &
364  this%grids(3)%dof%size())
365  !Restrict to middle level
366  call this%interp_fine_mid%map(this%e%x, this%r, &
367  this%msh%nelv, this%grids(2)%Xh)
368  call this%grids(2)%gs_h%op(this%e%x, &
369  this%grids(2)%dof%size(), gs_op_add, this%gs_event)
370  call device_event_sync(this%gs_event)
371  call device_copy(this%r_d, r_d, n)
372  call bc_list_apply_scalar(this%bclst_reg, r, n)
373  call device_copy(this%w_d, this%e%x_d, this%grids(2)%dof%size())
374  call bc_list_apply_scalar(this%bclst_mg, this%w, &
375  this%grids(2)%dof%size())
376  !OVERLAPPING Schwarz exchange and solve
377  call device_col2(this%w_d, this%grids(2)%coef%mult_d, &
378  this%grids(2)%dof%size())
379  !restrict residual to crs
380  call this%interp_mid_crs%map(this%wf%x, this%w, this%msh%nelv, &
381  this%grids(1)%Xh)
382  !Crs solve
383  call device_copy(this%w_d, this%e%x_d, this%grids(2)%dof%size())
384  call bc_list_apply_scalar(this%bclst_mg, this%w, &
385  this%grids(2)%dof%size())
386 
387  !$omp parallel private(thrdid, nthrds)
388 
389  thrdid = 0
390  nthrds = 1
391  !$ thrdid = omp_get_thread_num()
392  !$ nthrds = omp_get_num_threads()
393 
394  if (thrdid .eq. 0) then
395  call profiler_start_region('HSMG_schwarz', 9)
396  call this%grids(3)%schwarz%compute(z, this%r)
397  call this%grids(2)%schwarz%compute(this%grids(2)%e%x, this%w)
398  call profiler_end_region('HSMG_schwarz', 9)
399  end if
400  if (nthrds .eq. 1 .or. thrdid .eq. 1) then
401  call profiler_start_region('HSMG_coarse_grid', 10)
402  call this%grids(1)%gs_h%op(this%wf%x, &
403  this%grids(1)%dof%size(), gs_op_add, this%gs_event)
404  call device_event_sync(this%gs_event)
405  call bc_list_apply_scalar(this%grids(1)%bclst, this%wf%x, &
406  this%grids(1)%dof%size())
407  call profiler_start_region('HSMG_coarse_solve', 11)
408  crs_info = this%crs_solver%solve(this%Ax, this%grids(1)%e, &
409  this%wf%x, &
410  this%grids(1)%dof%size(), &
411  this%grids(1)%coef, &
412  this%grids(1)%bclst, &
413  this%grids(1)%gs_h, this%niter)
414  call profiler_end_region('HSMG_coarse_solve', 11)
415  call bc_list_apply_scalar(this%grids(1)%bclst, this%grids(1)%e%x,&
416  this%grids(1)%dof%size())
417  call profiler_end_region('HSMG_coarse_grid', 10)
418  end if
419  !$omp end parallel
420 
421  call this%interp_mid_crs%map(this%w, this%grids(1)%e%x, &
422  this%msh%nelv, this%grids(2)%Xh)
423  call device_add2(this%grids(2)%e%x_d, this%w_d, this%grids(2)%dof%size())
424 
425  call this%interp_fine_mid%map(this%w, this%grids(2)%e%x, &
426  this%msh%nelv, this%grids(3)%Xh)
427  call device_add2(z_d, this%w_d, this%grids(3)%dof%size())
428  call this%grids(3)%gs_h%op(z, this%grids(3)%dof%size(), &
429  gs_op_add, this%gs_event)
430  call device_event_sync(this%gs_event)
431  call device_col2(z_d, this%grids(3)%coef%mult_d, &
432  this%grids(3)%dof%size())
433  else
434  !We should not work with the input
435  call copy(this%r, r, n)
436 
437  !OVERLAPPING Schwarz exchange and solve
438  call this%grids(3)%schwarz%compute(z, this%r)
439  ! DOWNWARD Leg of V-cycle, we are pretty hardcoded here but w/e
440  call col2(this%r, this%grids(3)%coef%mult, &
441  this%grids(3)%dof%size())
442  !Restrict to middle level
443  call this%interp_fine_mid%map(this%w, this%r, &
444  this%msh%nelv, this%grids(2)%Xh)
445  call this%grids(2)%gs_h%op(this%w, this%grids(2)%dof%size(), gs_op_add)
446  !OVERLAPPING Schwarz exchange and solve
447  call this%grids(2)%schwarz%compute(this%grids(2)%e%x, this%w)
448  call col2(this%w, this%grids(2)%coef%mult, this%grids(2)%dof%size())
449  !restrict residual to crs
450  call this%interp_mid_crs%map(this%r, this%w, &
451  this%msh%nelv, this%grids(1)%Xh)
452  !Crs solve
453 
454  call this%grids(1)%gs_h%op(this%r, this%grids(1)%dof%size(), gs_op_add)
455  call bc_list_apply_scalar(this%grids(1)%bclst, this%r, &
456  this%grids(1)%dof%size())
457  call profiler_start_region('HSMG_coarse-solve', 11)
458  crs_info = this%crs_solver%solve(this%Ax, this%grids(1)%e, this%r, &
459  this%grids(1)%dof%size(), &
460  this%grids(1)%coef, &
461  this%grids(1)%bclst, &
462  this%grids(1)%gs_h, this%niter)
463  call profiler_end_region('HSMG_coarse-solve', 11)
464  call bc_list_apply_scalar(this%grids(1)%bclst, this%grids(1)%e%x,&
465  this%grids(1)%dof%size())
466 
467 
468  call this%interp_mid_crs%map(this%w, this%grids(1)%e%x, &
469  this%msh%nelv, this%grids(2)%Xh)
470  call add2(this%grids(2)%e%x, this%w, this%grids(2)%dof%size())
471 
472  call this%interp_fine_mid%map(this%w, this%grids(2)%e%x, &
473  this%msh%nelv, this%grids(3)%Xh)
474  call add2(z, this%w, this%grids(3)%dof%size())
475  call this%grids(3)%gs_h%op(z, this%grids(3)%dof%size(), gs_op_add)
476  call col2(z, this%grids(3)%coef%mult, this%grids(3)%dof%size())
477 
478  end if
479  call profiler_end_region('HSMG_solve', 8)
480  end subroutine hsmg_solve
481 end module hsmg
double real
Definition: device_config.h:12
Return the device pointer for an associated Fortran array.
Definition: device.F90:81
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
Defines a Matrix-vector product.
Definition: ax.f90:34
Defines a boundary condition.
Definition: bc.f90:34
subroutine, public bc_list_add(bclst, bc)
Add a condition to a list of boundary conditions.
Definition: bc.f90:505
subroutine, public bc_list_init(bclst, size)
Constructor for a list of boundary conditions.
Definition: bc.f90:464
subroutine, public bc_list_apply_scalar(bclst, x, n, t, tstep)
Apply a list of boundary conditions to a scalar field.
Definition: bc.f90:530
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_copy(a_d, b_d, n)
Copy a vector .
Definition: device_math.F90:76
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
subroutine, public device_event_sync(event)
Synchronize an event.
Definition: device.F90:1229
subroutine, public device_event_create(event, flags)
Create a device event queue.
Definition: device.F90:1164
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.
Krylov preconditioner.
Definition: pc_hsmg.f90:61
subroutine hsmg_solve(this, z, r, n)
The h1mg preconditioner from Nek5000.
Definition: pc_hsmg.f90:345
subroutine hsmg_init(this, msh, Xh, coef, dof, gs_h, bclst, crs_pctype)
Definition: pc_hsmg.f90:138
subroutine hsmg_set_h(this)
Definition: pc_hsmg.f90:264
subroutine hsmg_free(this)
Definition: pc_hsmg.f90:300
subroutine hsmg_fill_grid(dof, gs_h, Xh, coef, bclst, schwarz, e, grids, l)
Definition: pc_hsmg.f90:278
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:587
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:729
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:239
Defines a mesh.
Definition: mesh.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
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
Overlapping schwarz solves.
Definition: schwarz.f90:61
Defines a function space.
Definition: space.f90:34
integer, parameter, public gll
Definition: space.f90:48
Jacobi preconditioner SX-Aurora backend.
Utilities.
Definition: utils.f90:35
Base type for a matrix-vector product providing .
Definition: ax.f90:43
A list of boundary conditions.
Definition: bc.f90:104
Base type for a boundary condition.
Definition: bc.f90:51
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: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:66
Defines a canonical Krylov preconditioner.
Definition: precon.f90:40
The function space for the SEM solution fields.
Definition: space.f90:62
Defines a jacobi preconditioner for SX-Aurora.
#define max(a, b)
Definition: tensor.cu:40