43 use,
intrinsic :: iso_c_binding
48 real(kind=
rp),
allocatable :: x(:,:,:,:)
54 logical :: internal_dofmap = .false.
55 character(len=80) :: name =
""
56 type(c_ptr) :: x_d = c_null_ptr
59 procedure,
private, pass(this) :: init_external_dof => &
61 procedure,
private, pass(this) :: init_internal_dof => &
72 generic :: init => init_external_dof, init_internal_dof
74 generic ::
assignment(=) => assign_field, assign_scalar
78 generic :: add => add_field, add_scalar
90 class(
field_t),
intent(inout) :: this
91 type(
mesh_t),
target,
intent(in) :: msh
92 type(
space_t),
target,
intent(in) :: space
93 character(len=*),
optional :: fld_name
101 call this%dof%init(this%msh, this%Xh)
102 this%internal_dofmap = .true.
104 if (
present(fld_name))
then
105 call this%init_common(fld_name)
107 call this%init_common()
114 class(
field_t),
intent(inout) :: this
115 type(
dofmap_t),
target,
intent(in) :: dof
116 character(len=*),
optional :: fld_name
124 if (
present(fld_name))
then
125 call this%init_common(fld_name)
127 call this%init_common()
134 class(
field_t),
intent(inout) :: this
135 character(len=*),
optional :: fld_name
139 associate(lx => this%Xh%lx, ly => this%Xh%ly, &
140 lz => this%Xh%lz, nelv => this%msh%nelv)
142 if (.not.
allocated(this%x))
then
143 allocate(this%x(lx, ly, lz, nelv), stat = ierr)
147 if (
present(fld_name))
then
154 n = lx * ly * lz * nelv
157 real(c_rp) :: rp_dummy
158 integer(c_size_t) :: s
159 s = c_sizeof(rp_dummy) * n
169 class(
field_t),
intent(inout) :: this
172 if (
allocated(this%x))
then
176 if (this%internal_dofmap)
then
178 this%internal_dofmap = .false.
185 if (c_associated(this%x_d))
then
186 call device_free(this%x_d)
196 class(
field_t),
intent(inout) :: this
197 integer,
intent(in) :: memdir
198 logical,
intent(in) :: sync
200 if (neko_bcknd_device .eq. 1)
then
201 call device_memcpy(this%x, this%x_d, this%size(), memdir, sync)
211 class(
field_t),
intent(inout) :: this
214 if (
allocated(this%x))
then
215 if (.not.
associated(this%Xh, g%Xh))
then
222 if (len_trim(this%name) == 0)
then
226 if (.not. g%internal_dofmap)
then
229 if (this%internal_dofmap)
then
233 this%internal_dofmap = .true.
235 call this%dof%init(this%msh, this%Xh)
238 if (.not.
allocated(this%x))
then
240 allocate(this%x(this%Xh%lx, this%Xh%ly, this%Xh%lz, this%msh%nelv))
242 if (neko_bcknd_device .eq. 1)
then
243 call device_map(this%x, this%x_d, this%size())
248 if (neko_bcknd_device .eq. 1)
then
249 call device_copy(this%x_d, g%x_d, this%size())
251 call copy(this%x, g%x, this%dof%size())
258 class(
field_t),
intent(inout) :: this
259 real(kind=rp),
intent(in) :: a
260 integer :: i, j, k, l
262 if (neko_bcknd_device .eq. 1)
then
263 call device_cfill(this%x_d, a, this%size())
265 call cfill(this%x, a, this%size())
274 class(
field_t),
intent(inout) :: this
277 if (neko_bcknd_device .eq. 1)
then
278 call device_add2(this%x_d, g%x_d, this%size())
280 call add2(this%x, g%x, this%size())
289 class(
field_t),
intent(inout) :: this
290 real(kind=rp),
intent(in) :: a
292 if (neko_bcknd_device .eq. 1)
then
293 call device_cadd(this%x_d, a, this%size())
295 call cadd(this%x, a, this%size())
302 class(
field_t),
intent(in) :: this
305 size = this%dof%size()
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
Device abstraction, common interface for various accelerators.
subroutine, public device_free(x_d)
Deallocate memory on the device.
subroutine, public device_memset(x_d, v, s, sync, strm)
Set memory on the device to a value.
Defines a mapping of the degrees of freedom.
subroutine field_add_field(this, g)
Add .
pure integer function field_size(this)
Return the size of the field based on the underlying dofmap.
subroutine field_assign_field(this, g)
Assignment .
subroutine field_add_scalar(this, a)
Add .
subroutine field_copy_from(this, memdir, sync)
Easy way to copy between host and device.
subroutine field_assign_scalar(this, a)
Assignment .
subroutine field_init_common(this, fld_name)
Initialize a field this.
subroutine field_init_external_dof(this, dof, fld_name)
Initialize a field this on the mesh msh using an internal dofmap.
subroutine field_init_internal_dof(this, msh, space, fld_name)
Initialize a field this on the mesh msh using an internal dofmap.
subroutine field_free(this)
Deallocate a field f.
subroutine, public cadd(a, s, n)
Add a scalar to vector .
subroutine, public add2(a, b, n)
Vector addition .
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
subroutine, public copy(a, b, n)
Copy a vector .
integer, parameter neko_bcknd_device
integer, parameter, public c_rp
integer, parameter, public rp
Global precision used in computations.
Defines a function space.
field_ptr_t, To easily obtain a pointer to a field
The function space for the SEM solution fields.