Neko  0.9.0
A portable framework for high-order spectral element flow simulations
field.f90
Go to the documentation of this file.
1 ! Copyright (c) 2018-2023, 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 !
34 module field
35  use neko_config, only : neko_bcknd_device
36  use device_math
37  use num_types, only : rp
38  use device
39  use math, only : add2, copy, cadd
40  use mesh, only : mesh_t
41  use space, only : space_t, operator(.ne.)
42  use dofmap, only : dofmap_t
43  implicit none
44  private
45 
46  type, public :: field_t
47  real(kind=rp), allocatable :: x(:,:,:,:)
48 
49  type(space_t), pointer :: xh
50  type(mesh_t), pointer :: msh
51  type(dofmap_t), pointer :: dof
52 
53  logical :: internal_dofmap = .false.
54  character(len=80) :: name
55  type(c_ptr) :: x_d = c_null_ptr
56  contains
57  procedure, private, pass(this) :: init_common => field_init_common
58  procedure, private, pass(this) :: init_external_dof => &
60  procedure, private, pass(this) :: init_internal_dof => &
62  procedure, private, pass(this) :: assign_field => field_assign_field
63  procedure, private, pass(this) :: assign_scalar => field_assign_scalar
64  procedure, private, pass(this) :: add_field => field_add_field
65  procedure, private, pass(this) :: add_scalar => field_add_scalar
66  procedure, pass(this) :: free => field_free
68  procedure, pass(this) :: size => field_size
70  generic :: init => init_external_dof, init_internal_dof
72  generic :: assignment(=) => assign_field, assign_scalar
76  generic :: add => add_field, add_scalar
77  end type field_t
78 
80  type, public :: field_ptr_t
81  type(field_t), pointer :: ptr => null()
82  end type field_ptr_t
83 
84 contains
85 
87  subroutine field_init_internal_dof(this, msh, space, fld_name)
88  class(field_t), intent(inout) :: this
89  type(mesh_t), target, intent(in) :: msh
90  type(space_t), target, intent(in) :: space
91  character(len=*), optional :: fld_name
92 
93  call this%free()
94 
95  this%Xh => space
96  this%msh => msh
97 
98  allocate(this%dof)
99  call this%dof%init(this%msh, this%Xh)
100  this%internal_dofmap = .true.
101 
102  if (present(fld_name)) then
103  call this%init_common(fld_name)
104  else
105  call this%init_common()
106  end if
107 
108  end subroutine field_init_internal_dof
109 
111  subroutine field_init_external_dof(this, dof, fld_name)
112  class(field_t), intent(inout) :: this
113  type(dofmap_t), target, intent(in) :: dof
114  character(len=*), optional :: fld_name
115 
116  call this%free()
117 
118  this%dof => dof
119  this%Xh => dof%Xh
120  this%msh => dof%msh
121 
122  if (present(fld_name)) then
123  call this%init_common(fld_name)
124  else
125  call this%init_common()
126  end if
127 
128  end subroutine field_init_external_dof
129 
131  subroutine field_init_common(this, fld_name)
132  class(field_t), intent(inout) :: this
133  character(len=*), optional :: fld_name
134  integer :: ierr
135  integer :: n
136 
137  associate(lx => this%Xh%lx, ly => this%Xh%ly, &
138  lz => this%Xh%lz, nelv => this%msh%nelv)
139 
140  if (.not. allocated(this%x)) then
141  allocate(this%x(lx, ly, lz, nelv), stat = ierr)
142  this%x = 0d0
143  end if
144 
145  if (present(fld_name)) then
146  this%name = fld_name
147  else
148  this%name = "Field"
149  end if
150 
151  if (neko_bcknd_device .eq. 1) then
152  n = lx * ly * lz * nelv
153  call device_map(this%x, this%x_d, n)
154  end if
155  end associate
156 
157  end subroutine field_init_common
158 
160  subroutine field_free(this)
161  class(field_t), intent(inout) :: this
162 
163  if (allocated(this%x)) then
164  deallocate(this%x)
165  end if
166 
167  if (this%internal_dofmap) then
168  deallocate(this%dof)
169  this%internal_dofmap = .false.
170  end if
171 
172  nullify(this%msh)
173  nullify(this%Xh)
174  nullify(this%dof)
175 
176  if (c_associated(this%x_d)) then
177  call device_free(this%x_d)
178  end if
179 
180  end subroutine field_free
181 
185  subroutine field_assign_field(this, g)
186  class(field_t), intent(inout) :: this
187  type(field_t), intent(in) :: g
188 
189  if (allocated(this%x)) then
190  if (this%Xh .ne. g%Xh) then
191  call this%free()
192  end if
193  end if
194 
195  this%Xh => g%Xh
196  this%msh => g%msh
197  this%dof => g%dof
198 
199 
200  this%Xh%lx = g%Xh%lx
201  this%Xh%ly = g%Xh%ly
202  this%Xh%lz = g%Xh%lz
203 
204  if (.not. allocated(this%x)) then
205 
206  allocate(this%x(this%Xh%lx, this%Xh%ly, this%Xh%lz, this%msh%nelv))
207 
208  if (neko_bcknd_device .eq. 1) then
209  call device_map(this%x, this%x_d, this%size())
210  end if
211 
212  end if
213 
214  if (neko_bcknd_device .eq. 1) then
215  call device_copy(this%x_d, g%x_d, this%size())
216  else
217  call copy(this%x, g%x, this%dof%size())
218  end if
219 
220  end subroutine field_assign_field
221 
223  subroutine field_assign_scalar(this, a)
224  class(field_t), intent(inout) :: this
225  real(kind=rp), intent(in) :: a
226  integer :: i, j, k, l
227 
228  if (neko_bcknd_device .eq. 1) then
229  call device_cfill(this%x_d, a, this%size())
230  else
231  do i = 1, this%msh%nelv
232  do l = 1, this%Xh%lz
233  do k = 1, this%Xh%ly
234  do j = 1, this%Xh%lx
235  this%x(j, k, l, i) = a
236  end do
237  end do
238  end do
239  end do
240  end if
241 
242  end subroutine field_assign_scalar
243 
247  subroutine field_add_field(this, g)
248  class(field_t), intent(inout) :: this
249  type(field_t), intent(in) :: g
250 
251  if (neko_bcknd_device .eq. 1) then
252  call device_add2(this%x_d, g%x_d, this%size())
253  else
254  call add2(this%x, g%x, this%size())
255  end if
256 
257  end subroutine field_add_field
258 
259 
262  subroutine field_add_scalar(this, a)
263  class(field_t), intent(inout) :: this
264  real(kind=rp), intent(in) :: a
265 
266  if (neko_bcknd_device .eq. 1) then
267  call device_cadd(this%x_d, a, this%size())
268  else
269  call cadd(this%x, a, this%size())
270  end if
271 
272  end subroutine field_add_scalar
273 
275  pure function field_size(this) result(size)
276  class(field_t), intent(in) :: this
277  integer :: size
278 
279  size = this%dof%size()
280  end function field_size
281 
282 end module field
283 
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Defines a field.
Definition: field.f90:34
subroutine field_add_field(this, g)
Add .
Definition: field.f90:248
pure integer function field_size(this)
Return the size of the field based on the underlying dofmap.
Definition: field.f90:276
subroutine field_assign_field(this, g)
Assignment .
Definition: field.f90:186
subroutine field_add_scalar(this, a)
Add .
Definition: field.f90:263
subroutine field_assign_scalar(this, a)
Assignment .
Definition: field.f90:224
subroutine field_init_common(this, fld_name)
Initialize a field this.
Definition: field.f90:132
subroutine field_init_external_dof(this, dof, fld_name)
Initialize a field this on the mesh msh using an internal dofmap.
Definition: field.f90:112
subroutine field_init_internal_dof(this, msh, space, fld_name)
Initialize a field this on the mesh msh using an internal dofmap.
Definition: field.f90:88
subroutine field_free(this)
Deallocate a field f.
Definition: field.f90:161
Definition: math.f90:60
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition: math.f90:323
subroutine, public add2(a, b, n)
Vector addition .
Definition: math.f90:587
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
Defines a function space.
Definition: space.f90:34
field_ptr_t, To easily obtain a pointer to a field
Definition: field.f90:80
The function space for the SEM solution fields.
Definition: space.f90:62