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 => &
71 generic :: init => init_external_dof, init_internal_dof
73 generic ::
assignment(=) => assign_field, assign_scalar
77 generic :: add => add_field, add_scalar
89 class(
field_t),
intent(inout) :: this
90 type(
mesh_t),
target,
intent(in) :: msh
91 type(
space_t),
target,
intent(in) :: space
92 character(len=*),
optional :: fld_name
100 call this%dof%init(this%msh, this%Xh)
101 this%internal_dofmap = .true.
103 if (
present(fld_name))
then
104 call this%init_common(fld_name)
106 call this%init_common()
113 class(
field_t),
intent(inout) :: this
114 type(
dofmap_t),
target,
intent(in) :: dof
115 character(len=*),
optional :: fld_name
123 if (
present(fld_name))
then
124 call this%init_common(fld_name)
126 call this%init_common()
133 class(
field_t),
intent(inout) :: this
134 character(len=*),
optional :: fld_name
138 associate(lx => this%Xh%lx, ly => this%Xh%ly, &
139 lz => this%Xh%lz, nelv => this%msh%nelv)
141 if (.not.
allocated(this%x))
then
142 allocate(this%x(lx, ly, lz, nelv), stat = ierr)
146 if (
present(fld_name))
then
153 n = lx * ly * lz * nelv
162 class(
field_t),
intent(inout) :: this
164 if (
allocated(this%x))
then
168 if (this%internal_dofmap)
then
170 this%internal_dofmap = .false.
177 if (c_associated(this%x_d))
then
178 call device_free(this%x_d)
187 class(
field_t),
intent(inout) :: this
190 if (
allocated(this%x))
then
191 if (this%Xh .ne. g%Xh)
then
205 if (.not.
allocated(this%x))
then
207 allocate(this%x(this%Xh%lx, this%Xh%ly, this%Xh%lz, this%msh%nelv))
209 if (neko_bcknd_device .eq. 1)
then
210 call device_map(this%x, this%x_d, this%size())
215 if (neko_bcknd_device .eq. 1)
then
216 call device_copy(this%x_d, g%x_d, this%size())
218 call copy(this%x, g%x, this%dof%size())
225 class(
field_t),
intent(inout) :: this
226 real(kind=rp),
intent(in) :: a
227 integer :: i, j, k, l
229 if (neko_bcknd_device .eq. 1)
then
230 call device_cfill(this%x_d, a, this%size())
232 do i = 1, this%msh%nelv
236 this%x(j, k, l, i) = a
249 class(
field_t),
intent(inout) :: this
252 if (neko_bcknd_device .eq. 1)
then
253 call device_add2(this%x_d, g%x_d, this%size())
255 call add2(this%x, g%x, this%size())
264 class(
field_t),
intent(inout) :: this
265 real(kind=rp),
intent(in) :: a
267 if (neko_bcknd_device .eq. 1)
then
268 call device_cadd(this%x_d, a, this%size())
270 call cadd(this%x, a, this%size())
277 class(
field_t),
intent(in) :: this
280 size = this%dof%size()
Map a Fortran array to a device (allocate and associate)
Device abstraction, common interface for various accelerators.
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_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 copy(a, b, n)
Copy a vector .
integer, parameter neko_bcknd_device
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.