Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
bc.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2024, 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!
34module bc
35 use neko_config
36 use num_types
37 use device
38 use dofmap, only : dofmap_t
39 use coefs, only : coef_t
40 use space, only : space_t
42 use facet_zone, only : facet_zone_t
43 use stack, only : stack_i4t2_t
44 use tuple, only : tuple_i4_t
46 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
47 implicit none
48 private
49
51 type, public, abstract :: bc_t
53 integer, allocatable :: msk(:)
55 integer, allocatable :: facet(:)
57 type(dofmap_t), pointer :: dof
59 type(coef_t), pointer :: coef
61 type(mesh_t), pointer :: msh
63 type(space_t), pointer :: xh
65 type(stack_i4t2_t) :: marked_facet
67 type(c_ptr) :: msk_d = c_null_ptr
69 type(c_ptr) :: facet_d = c_null_ptr
70 contains
72 procedure, pass(this) :: init_base => bc_init_base
74 procedure, pass(this) :: free_base => bc_free_base
76 procedure, pass(this) :: mark_facet => bc_mark_facet
78 procedure, pass(this) :: mark_facets => bc_mark_facets
80 procedure, pass(this) :: mark_zones_from_list => bc_mark_zones_from_list
82 procedure, pass(this) :: mark_zone => bc_mark_zone
85 procedure, pass(this) :: finalize => bc_finalize
87 procedure(bc_apply_scalar), pass(this), deferred :: apply_scalar
89 procedure(bc_apply_vector), pass(this), deferred :: apply_vector
91 procedure(bc_apply_scalar_dev), pass(this), deferred :: apply_scalar_dev
93 procedure(bc_apply_vector_dev), pass(this), deferred :: apply_vector_dev
95 procedure(bc_destructor), pass(this), deferred :: free
96 end type bc_t
97
99 type, public :: bc_ptr_t
100 class(bc_t), pointer :: ptr
101 end type bc_ptr_t
102
104 type, public :: bc_alloc_t
105 class(bc_t), allocatable :: obj
106 end type bc_alloc_t
107
108 abstract interface
109
114 subroutine bc_apply_scalar(this, x, n, t, tstep)
115 import :: bc_t
116 import :: rp
117 class(bc_t), intent(inout) :: this
118 integer, intent(in) :: n
119 real(kind=rp), intent(inout), dimension(n) :: x
120 real(kind=rp), intent(in), optional :: t
121 integer, intent(in), optional :: tstep
122 end subroutine bc_apply_scalar
123 end interface
124
125 abstract interface
126
133 subroutine bc_apply_vector(this, x, y, z, n, t, tstep)
134 import :: bc_t
135 import :: rp
136 class(bc_t), intent(inout) :: this
137 integer, intent(in) :: n
138 real(kind=rp), intent(inout), dimension(n) :: x
139 real(kind=rp), intent(inout), dimension(n) :: y
140 real(kind=rp), intent(inout), dimension(n) :: z
141 real(kind=rp), intent(in), optional :: t
142 integer, intent(in), optional :: tstep
143 end subroutine bc_apply_vector
144 end interface
145
146 abstract interface
147
148 subroutine bc_destructor(this)
149 import :: bc_t
150 class(bc_t), intent(inout), target :: this
151 end subroutine bc_destructor
152 end interface
153
154 abstract interface
155
157 subroutine bc_apply_scalar_dev(this, x_d, t, tstep)
158 import :: c_ptr
159 import :: bc_t
160 import :: rp
161 class(bc_t), intent(inout), target :: this
162 type(c_ptr) :: x_d
163 real(kind=rp), intent(in), optional :: t
164 integer, intent(in), optional :: tstep
165 end subroutine bc_apply_scalar_dev
166 end interface
167
168 abstract interface
169
173 subroutine bc_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
174 import :: c_ptr
175 import :: bc_t
176 import :: rp
177 class(bc_t), intent(inout), target :: this
178 type(c_ptr) :: x_d
179 type(c_ptr) :: y_d
180 type(c_ptr) :: z_d
181 real(kind=rp), intent(in), optional :: t
182 integer, intent(in), optional :: tstep
183 end subroutine bc_apply_vector_dev
184 end interface
185
186contains
187
190 subroutine bc_init_base(this, coef)
191 class(bc_t), intent(inout) :: this
192 type(coef_t), target, intent(in) :: coef
193
194 call this%free_base
195
196 this%dof => coef%dof
197 this%coef => coef
198 this%Xh => this%dof%Xh
199 this%msh => this%dof%msh
200
201 call this%marked_facet%init()
202
203 end subroutine bc_init_base
204
206 subroutine bc_free_base(this)
207 class(bc_t), intent(inout) :: this
208
209 call this%marked_facet%free()
210
211 nullify(this%Xh)
212 nullify(this%msh)
213 nullify(this%dof)
214 nullify(this%coef)
215
216 if (allocated(this%msk)) then
217 deallocate(this%msk)
218 end if
219
220 if (allocated(this%facet)) then
221 deallocate(this%facet)
222 end if
223
224 if (c_associated(this%msk_d)) then
225 call device_free(this%msk_d)
226 this%msk_d = c_null_ptr
227 end if
228
229 if (c_associated(this%facet_d)) then
230 call device_free(this%facet_d)
231 this%facet_d = c_null_ptr
232 end if
233
234 end subroutine bc_free_base
235
239 subroutine bc_mark_facet(this, facet, el)
240 class(bc_t), intent(inout) :: this
241 integer, intent(in) :: facet
242 integer, intent(in) :: el
243 type(tuple_i4_t) :: t
244
245 t%x = (/facet, el/)
246 call this%marked_facet%push(t)
247
248 end subroutine bc_mark_facet
249
252 subroutine bc_mark_facets(this, facet_list)
253 class(bc_t), intent(inout) :: this
254 type(stack_i4t2_t), intent(inout) :: facet_list
255 type(tuple_i4_t), pointer :: fp(:)
256 integer :: i
257
258 fp => facet_list%array()
259 do i = 1, facet_list%size()
260 call this%marked_facet%push(fp(i))
261 end do
262
263 end subroutine bc_mark_facets
264
267 subroutine bc_mark_zone(this, bc_zone)
268 class(bc_t), intent(inout) :: this
269 class(facet_zone_t), intent(in) :: bc_zone
270 integer :: i
271 do i = 1, bc_zone%size
272 call this%marked_facet%push(bc_zone%facet_el(i))
273 end do
274 end subroutine bc_mark_zone
275
282 subroutine bc_mark_zones_from_list(this, bc_zones, bc_key, bc_labels)
283 class(bc_t), intent(inout) :: this
284 class(facet_zone_t), intent(in) :: bc_zones(:)
285 character(len=*) :: bc_key
286 character(len=100), allocatable :: split_key(:)
287 character(len=NEKO_MSH_MAX_ZLBL_LEN) :: bc_labels(:)
288 integer :: i, j, k, l, msh_bc_type
289
290 msh_bc_type = 0
291 if(trim(bc_key) .eq. 'o' .or. trim(bc_key) .eq. 'on' &
292 .or. trim(bc_key) .eq. 'o+dong' .or. trim(bc_key) .eq. 'on+dong') then
293 msh_bc_type = 1
294 else if(trim(bc_key) .eq. 'd_pres') then
295 msh_bc_type = 1
296 else if(trim(bc_key) .eq. 'w') then
297 msh_bc_type = 2
298 else if(trim(bc_key) .eq. 'v') then
299 msh_bc_type = 2
300 else if(trim(bc_key) .eq. 'd_vel_u') then
301 msh_bc_type = 2
302 else if(trim(bc_key) .eq. 'd_vel_v') then
303 msh_bc_type = 2
304 else if(trim(bc_key) .eq. 'd_vel_w') then
305 msh_bc_type = 2
306 else if(trim(bc_key) .eq. 'sym') then
307 msh_bc_type = 2
308 else if(trim(bc_key) .eq. 'sh') then
309 msh_bc_type = 2
310 else if(trim(bc_key) .eq. 'wm') then
311 msh_bc_type = 2
312 end if
313
314 do i = 1, size(bc_labels)
315 !Check if several bcs are defined for this zone
316 !bcs are seperated by /, but we could use something else
317 if (index(trim(bc_labels(i)), '/') .eq. 0) then
318 if (trim(bc_key) .eq. trim(bc_labels(i))) then
319 call bc_mark_zone(this, bc_zones(i))
320 ! Loop across all faces in the mesh
321 do j = 1,this%msh%nelv
322 do k = 1, 2 * this%msh%gdim
323 if (this%msh%facet_type(k,j) .eq. -i) then
324 this%msh%facet_type(k,j) = msh_bc_type
325 end if
326 end do
327 end do
328 end if
329 else
330 split_key = split_string(trim(bc_labels(i)),'/')
331 do l = 1, size(split_key)
332 if (trim(split_key(l)) .eq. trim(bc_key)) then
333 call bc_mark_zone(this, bc_zones(i))
334 ! Loop across all faces in the mesh
335 do j = 1,this%msh%nelv
336 do k = 1, 2 * this%msh%gdim
337 if (this%msh%facet_type(k,j) .eq. -i) then
338 this%msh%facet_type(k,j) = msh_bc_type
339 end if
340 end do
341 end do
342 end if
343 end do
344 end if
345 end do
346 end subroutine bc_mark_zones_from_list
347
348
352 subroutine bc_finalize(this)
353 class(bc_t), target, intent(inout) :: this
354 type(tuple_i4_t), pointer :: bfp(:)
355 type(tuple_i4_t) :: bc_facet
356 integer :: facet_size, facet, el
357 integer :: i, j, k, l, msk_c
358 integer :: lx, ly, lz, n
359
360 lx = this%Xh%lx
361 ly = this%Xh%ly
362 lz = this%Xh%lz
363
365
366 ! Note we assume that lx = ly = lz
367 facet_size = lx**2
368 allocate(this%msk(0:facet_size * this%marked_facet%size()))
369 allocate(this%facet(0:facet_size * this%marked_facet%size()))
370
371 msk_c = 0
372 bfp => this%marked_facet%array()
373
374 ! Loop through each (facet, element) id tuple
375 ! Then loop over all the nodes of the face and compute their linear index
376 ! This index goes into this%msk, whereas the corresponding face id goes into
377 ! this%facet
378 do i = 1, this%marked_facet%size()
379 bc_facet = bfp(i)
380 facet = bc_facet%x(1)
381 el = bc_facet%x(2)
382 select case (facet)
383 case (1)
384 do l = 1, lz
385 do k = 1, ly
386 msk_c = msk_c + 1
387 this%msk(msk_c) = linear_index(1,k,l,el,lx,ly,lz)
388 this%facet(msk_c) = 1
389 end do
390 end do
391 case (2)
392 do l = 1, lz
393 do k = 1, ly
394 msk_c = msk_c + 1
395 this%msk(msk_c) = linear_index(lx,k,l,el,lx,ly,lz)
396 this%facet(msk_c) = 2
397 end do
398 end do
399 case(3)
400 do l = 1, lz
401 do j = 1, lx
402 msk_c = msk_c + 1
403 this%msk(msk_c) = linear_index(j,1,l,el,lx,ly,lz)
404 this%facet(msk_c) = 3
405 end do
406 end do
407 case(4)
408 do l = 1, lz
409 do j = 1, lx
410 msk_c = msk_c + 1
411 this%msk(msk_c) = linear_index(j,ly,l,el,lx,ly,lz)
412 this%facet(msk_c) = 4
413 end do
414 end do
415 case(5)
416 do k = 1, ly
417 do j = 1, lx
418 msk_c = msk_c + 1
419 this%msk(msk_c) = linear_index(j,k,1,el,lx,ly,lz)
420 this%facet(msk_c) = 5
421 end do
422 end do
423 case(6)
424 do k = 1, ly
425 do j = 1, lx
426 msk_c = msk_c + 1
427 this%msk(msk_c) = linear_index(j,k,lz,el,lx,ly,lz)
428 this%facet(msk_c) = 6
429 end do
430 end do
431 end select
432 end do
433
434 this%msk(0) = msk_c
435 this%facet(0) = msk_c
436
437 if (neko_bcknd_device .eq. 1) then
438 n = facet_size * this%marked_facet%size() + 1
439 call device_map(this%msk, this%msk_d, n)
440 call device_map(this%facet, this%facet_d, n)
441
442 call device_memcpy(this%msk, this%msk_d, n, &
443 host_to_device, sync=.false.)
444 call device_memcpy(this%facet, this%facet_d, n, &
445 host_to_device, sync=.false.)
446 end if
447
448 end subroutine bc_finalize
449
450end module bc
Apply the boundary condition to a scalar field on the device.
Definition bc.f90:157
Apply the boundary condition to a scalar field.
Definition bc.f90:114
Apply the boundary condition to a vector field on the device.
Definition bc.f90:173
Apply the boundary condition to a vector field.
Definition bc.f90:133
Destructor.
Definition bc.f90:148
Map a Fortran array to a device (allocate and associate)
Definition device.F90:57
Copy data between host and device (or device and device)
Definition device.F90:51
Defines a boundary condition.
Definition bc.f90:34
subroutine bc_mark_zone(this, bc_zone)
Mark all facets from a zone.
Definition bc.f90:268
subroutine bc_finalize(this)
Finalize the construction of the bc by populting the msk and facet arrays.
Definition bc.f90:353
subroutine bc_free_base(this)
Destructor for the base type, bc_t.
Definition bc.f90:207
subroutine bc_init_base(this, coef)
Constructor.
Definition bc.f90:191
subroutine bc_mark_zones_from_list(this, bc_zones, bc_key, bc_labels)
Mark all facets from a list of zones, also marks type of bc in the mesh. The facet_type in mesh is be...
Definition bc.f90:283
subroutine bc_mark_facet(this, facet, el)
Mark facet on element el as part of the boundary condition.
Definition bc.f90:240
subroutine bc_mark_facets(this, facet_list)
Mark all facets from a (facet, el) tuple list.
Definition bc.f90:253
Coefficients.
Definition coef.f90:34
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
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
Defines a zone as a subset of facets in a mesh.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:56
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:58
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a function space.
Definition space.f90:34
Implements a dynamic stack ADT.
Definition stack.f90:35
Implements a n-tuple.
Definition tuple.f90:34
Utilities.
Definition utils.f90:35
character(len=100) function, dimension(:), allocatable, public split_string(string, delimiter)
Split a string based on delimiter (tokenizer) OBS: very hacky, this should really be improved,...
Definition utils.f90:136
pure integer function, public linear_index(i, j, k, l, lx, ly, lz)
Compute the address of a (i,j,k,l) array with sizes (1:lx, 1:ly, 1:lz, :)
Definition utils.f90:174
Wrapper around an allocatable bc.
Definition bc.f90:104
Pointer to boundary condtiion.
Definition bc.f90:99
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
The function space for the SEM solution fields.
Definition space.f90:62
Integer 2-tuple based stack.
Definition stack.f90:84
Integer based 2-tuple.
Definition tuple.f90:51