Neko  0.9.99
A portable framework for high-order spectral element flow simulations
pc_jacobi_device.F90
Go to the documentation of this file.
1 ! Copyright (c) 2021-2022, 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 precon
36  use coefs
37  use dofmap
38  use num_types
39  use device_math
40  use device
41  use gather_scatter
42  implicit none
43  private
44 
46  type, public, extends(pc_t) :: device_jacobi_t
47  real(kind=rp), allocatable :: d(:,:,:,:)
48  type(c_ptr) :: d_d = c_null_ptr
49  type(gs_t), pointer :: gs_h
50  type(dofmap_t), pointer :: dof
51  type(coef_t), pointer :: coef
52  contains
53  procedure, pass(this) :: init => device_jacobi_init
54  procedure, pass(this) :: free => device_jacobi_free
55  procedure, pass(this) :: solve => device_jacobi_solve
56  procedure, pass(this) :: update => device_jacobi_update
57  end type device_jacobi_t
58 
59  interface
60  subroutine hip_jacobi_update(d_d, dxt_d, dyt_d, dzt_d, &
61  G11_d, G22_d, G33_d, G12_d, G13_d, G23_d, nelv, lx) &
62  bind(c, name='hip_jacobi_update')
63  use, intrinsic :: iso_c_binding
64  type(c_ptr), value :: d_d, dxt_d, dyt_d, dzt_d
65  type(c_ptr), value :: G11_d, G22_d, G33_d, G12_d, G13_d, G23_d
66  integer(c_int) :: nelv, lx
67  end subroutine hip_jacobi_update
68  end interface
69 
70  interface
71  subroutine cuda_jacobi_update(d_d, dxt_d, dyt_d, dzt_d, &
72  G11_d, G22_d, G33_d, G12_d, G13_d, G23_d, nelv, lx) &
73  bind(c, name='cuda_jacobi_update')
74  use, intrinsic :: iso_c_binding
75  type(c_ptr), value :: d_d, dxt_d, dyt_d, dzt_d
76  type(c_ptr), value :: G11_d, G22_d, G33_d, G12_d, G13_d, G23_d
77  integer(c_int) :: nelv, lx
78  end subroutine cuda_jacobi_update
79  end interface
80 
81  interface
82  subroutine opencl_jacobi_update(d_d, dxt_d, dyt_d, dzt_d, &
83  G11_d, G22_d, G33_d, G12_d, G13_d, G23_d, nelv, lx) &
84  bind(c, name='opencl_jacobi_update')
85  use, intrinsic :: iso_c_binding
86  type(c_ptr), value :: d_d, dxt_d, dyt_d, dzt_d
87  type(c_ptr), value :: G11_d, G22_d, G33_d, G12_d, G13_d, G23_d
88  integer(c_int) :: nelv, lx
89  end subroutine opencl_jacobi_update
90  end interface
91 
92 contains
93 
94  subroutine device_jacobi_init(this, coef, dof, gs_h)
95  class(device_jacobi_t), intent(inout) :: this
96  type(coef_t), intent(inout), target :: coef
97  type(dofmap_t), intent(inout), target :: dof
98  type(gs_t), intent(inout), target :: gs_h
99 
100  call this%free()
101 
102  this%gs_h => gs_h
103  this%dof => dof
104  this%coef => coef
105 
106  allocate(this%d(dof%Xh%lx,dof%Xh%ly,dof%Xh%lz, dof%msh%nelv))
107 
108  call device_map(this%d, this%d_d, size(this%d))
109 
110  call device_jacobi_update(this)
111 
112  end subroutine device_jacobi_init
113 
114  subroutine device_jacobi_free(this)
115  class(device_jacobi_t), intent(inout) :: this
116 
117  if (c_associated(this%d_d)) then
118  call device_free(this%d_d)
119  this%d_d = c_null_ptr
120  end if
121 
122  if (allocated(this%d)) then
123  deallocate(this%d)
124  end if
125 
126  nullify(this%dof)
127  nullify(this%gs_h)
128  nullify(this%coef)
129  end subroutine device_jacobi_free
130 
133  subroutine device_jacobi_solve(this, z, r, n)
134  integer, intent(in) :: n
135  class(device_jacobi_t), intent(inout) :: this
136  real(kind=rp), dimension(n), intent(inout) :: z
137  real(kind=rp), dimension(n), intent(inout) :: r
138  type(c_ptr) :: z_d, r_d
139 
140  z_d = device_get_ptr(z)
141  r_d = device_get_ptr(r)
142 
143  call device_col3(z_d, r_d, this%d_d, n)
144 
145  end subroutine device_jacobi_solve
146 
147  subroutine device_jacobi_update(this)
148  class(device_jacobi_t), intent(inout) :: this
149  integer :: lz, ly, lx
150  associate(dof => this%dof, coef => this%coef, xh => this%dof%Xh, &
151  gs_h => this%gs_h, nelv => this%dof%msh%nelv)
152 
153  lx = xh%lx
154  ly = xh%ly
155  lz = xh%lz
156 
157 #ifdef HAVE_HIP
158  call hip_jacobi_update(this%d_d, xh%dxt_d, xh%dyt_d, xh%dzt_d, &
159  coef%G11_d, coef%G22_d, coef%G33_d, &
160  coef%G12_d, coef%G13_d, coef%G23_d, &
161  nelv, lx)
162 #elif HAVE_CUDA
163  call cuda_jacobi_update(this%d_d, xh%dxt_d, xh%dyt_d, xh%dzt_d, &
164  coef%G11_d, coef%G22_d, coef%G33_d, &
165  coef%G12_d, coef%G13_d, coef%G23_d, &
166  nelv, lx)
167 #elif HAVE_OPENCL
168  call opencl_jacobi_update(this%d_d, xh%dxt_d, xh%dyt_d, xh%dzt_d, &
169  coef%G11_d, coef%G22_d, coef%G33_d, &
170  coef%G12_d, coef%G13_d, coef%G23_d, &
171  nelv, lx)
172 #endif
173 
174  call device_col2(this%d_d, coef%h1_d, coef%dof%size())
175 
176  if (coef%ifh2) then
177  call device_addcol3(this%d_d, coef%h2_d, coef%B_d, coef%dof%size())
178  end if
179 
180  call gs_h%op(this%d, dof%size(), gs_op_add)
181 
182  call device_invcol1(this%d_d, dof%size())
183  end associate
184  end subroutine device_jacobi_update
185 
186 end module device_jacobi
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
Coefficients.
Definition: coef.f90:34
Jacobi preconditioner accelerator backend.
subroutine device_jacobi_free(this)
subroutine device_jacobi_init(this, coef, dof, gs_h)
subroutine device_jacobi_update(this)
subroutine device_jacobi_solve(this, z, r, n)
The jacobi preconditioner where .
subroutine, public device_addcol3(a_d, b_d, c_d, n)
Returns .
subroutine, public device_col2(a_d, b_d, n)
Vector multiplication .
subroutine, public device_invcol1(a_d, n)
Invert a vector .
subroutine, public device_col3(a_d, b_d, c_d, n)
Vector multiplication with 3 vectors .
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition: device.F90:185
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Gather-scatter.
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Krylov preconditioner.
Definition: precon.f90:34
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
Defines a jacobi preconditioner.
Defines a canonical Krylov preconditioner.
Definition: precon.f90:40