Neko  0.8.99
A portable framework for high-order spectral element flow simulations
brinkman_source_term.f90
Go to the documentation of this file.
1 ! Copyright (c) 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 !
35  use num_types, only: rp, dp
36  use field, only: field_t
37  use field_list, only: field_list_t
38  use json_module, only: json_file
41  use source_term, only: source_term_t
42  use coefs, only: coef_t
44  use utils, only: neko_error
45  use field_math, only: field_subcol3
46  implicit none
47  private
48 
52  type, public, extends(source_term_t) :: brinkman_source_term_t
53  private
54 
56  type(field_t) :: indicator
58  type(field_t) :: brinkman
59  contains
61  procedure, public, pass(this) :: init => &
64  procedure, public, pass(this) :: free => brinkman_source_term_free
66  procedure, public, pass(this) :: compute_ => brinkman_source_term_compute
67 
68  ! ----------------------------------------------------------------------- !
69  ! Private methods
70  procedure, pass(this) :: init_boundary_mesh
71  procedure, pass(this) :: init_point_zone
72  end type brinkman_source_term_t
73 
74 contains
75 
76  ! ========================================================================== !
77  ! Public methods
78 
83  subroutine brinkman_source_term_init_from_json(this, json, fields, coef)
84  use file, only: file_t
85  use tri_mesh, only: tri_mesh_t
91  use json_module, only: json_core, json_value
92  implicit none
93 
94  class(brinkman_source_term_t), intent(inout) :: this
95  type(json_file), intent(inout) :: json
96  type(field_list_t), intent(inout), target :: fields
97  type(coef_t), intent(inout), target :: coef
98  real(kind=rp) :: start_time, end_time
99 
100  character(len=:), allocatable :: filter_type
101  real(kind=rp), dimension(:), allocatable :: brinkman_limits
102  real(kind=rp) :: brinkman_penalty
103 
104  type(json_value), pointer :: json_object_list
105  type(json_core) :: core
106 
107  character(len=:), allocatable :: object_type
108  type(json_file) :: object_settings
109  integer :: n_regions
110  integer :: i
111 
112  ! Mandatory fields for the general source term
113  call json_get_or_default(json, "start_time", start_time, 0.0_rp)
114  call json_get_or_default(json, "end_time", end_time, huge(0.0_rp))
115 
116  ! Read the options for the permeability field
117  call json_get(json, 'brinkman.limits', brinkman_limits)
118  call json_get(json, 'brinkman.penalty', brinkman_penalty)
119 
120  if (size(brinkman_limits) .ne. 2) then
121  call neko_error('brinkman_limits must be a 2 element array of reals')
122  end if
123 
124  call this%free()
125  call this%init_base(fields, coef, start_time, end_time)
126 
127  ! ------------------------------------------------------------------------ !
128  ! Allocate the permeability and indicator field
129 
130  if (neko_field_registry%field_exists('brinkman_indicator') &
131  .or. neko_field_registry%field_exists('brinkman')) then
132  call neko_error('Brinkman field already exists.')
133  end if
134 
135  call this%indicator%init(coef%dof)
136  call this%brinkman%init(coef%dof)
137 
138  ! ------------------------------------------------------------------------ !
139  ! Select which constructor should be called
140 
141  call json%get('objects', json_object_list)
142  call json%info('objects', n_children = n_regions)
143  call json%get_core(core)
144 
145  do i = 1, n_regions
146  call json_extract_item(core, json_object_list, i, object_settings)
147  call json_get_or_default(object_settings, 'type', object_type, 'none')
148 
149  select case (object_type)
150  case ('boundary_mesh')
151  call this%init_boundary_mesh(object_settings)
152  case ('point_zone')
153  call this%init_point_zone(object_settings)
154 
155  case ('none')
156  call object_settings%print()
157  call neko_error('Brinkman source term objects require a region type')
158  case default
159  call neko_error('Brinkman source term unknown region type')
160  end select
161 
162  end do
163 
164  ! Run filter on the full indicator field to smooth it out.
165  call json_get_or_default(json, 'filter.type', filter_type, 'none')
166 
167  select case (filter_type)
168  case ('none')
169  ! Do nothing
170  case default
171  call neko_error('Brinkman source term unknown filter type')
172  end select
173 
174  ! ------------------------------------------------------------------------ !
175  ! Compute the permeability field
176 
177  call permeability_field(this%brinkman, this%indicator, &
178  & brinkman_limits(1), brinkman_limits(2), brinkman_penalty)
179 
180  ! Copy the permeability field to the device
181  if (neko_bcknd_device .eq. 1) then
182  call device_memcpy(this%brinkman%x, this%brinkman%x_d, &
183  this%brinkman%dof%size(), host_to_device, .true.)
184  end if
185 
187 
189  subroutine brinkman_source_term_free(this)
190  class(brinkman_source_term_t), intent(inout) :: this
191 
192  call this%indicator%free()
193  call this%brinkman%free()
194  call this%free_base()
195  end subroutine brinkman_source_term_free
196 
200  subroutine brinkman_source_term_compute(this, t, tstep)
201  class(brinkman_source_term_t), intent(inout) :: this
202  real(kind=rp), intent(in) :: t
203  integer, intent(in) :: tstep
204  type(field_t), pointer :: u, v, w, fu, fv, fw
205  integer :: n
206 
207  n = this%fields%item_size(1)
208 
209  u => neko_field_registry%get_field('u')
210  v => neko_field_registry%get_field('v')
211  w => neko_field_registry%get_field('w')
212 
213  fu => this%fields%get(1)
214  fv => this%fields%get(2)
215  fw => this%fields%get(3)
216 
217  call field_subcol3(fu, u, this%brinkman, n)
218  call field_subcol3(fv, v, this%brinkman, n)
219  call field_subcol3(fw, w, this%brinkman, n)
220 
221  end subroutine brinkman_source_term_compute
222 
223  ! ========================================================================== !
224  ! Private methods
225 
227  subroutine init_boundary_mesh(this, json)
228  use file, only: file_t
229  use tri_mesh, only: tri_mesh_t
235  use aabb, only : aabb_t, get_aabb
236  implicit none
237 
238  class(brinkman_source_term_t), intent(inout) :: this
239  type(json_file), intent(inout) :: json
240 
241  ! Options
242  character(len=:), allocatable :: mesh_file_name
243  character(len=:), allocatable :: distance_transform
244  character(len=:), allocatable :: filter_type
245  character(len=:), allocatable :: mesh_transform
246 
247  ! Read the options for the boundary mesh
248  type(file_t) :: mesh_file
249  type(tri_mesh_t) :: boundary_mesh
250  real(kind=rp) :: scalar_r
251  real(kind=dp) :: scalar_d
252 
253  ! Mesh transform options variables
254  real(kind=dp), dimension(:), allocatable :: box_min, box_max
255  logical :: keep_aspect_ratio
256  real(kind=dp), dimension(3) :: scaling
257  real(kind=dp), dimension(3) :: translation
258  type(field_t) :: temp_field
259  type(aabb_t) :: mesh_box, target_box
260  integer :: idx_p
261 
262  ! ------------------------------------------------------------------------ !
263  ! Read the options for the boundary mesh
264 
265  call json_get(json, 'name', mesh_file_name)
266 
267  ! Settings on how to filter the design field
268  call json_get(json, 'distance_transform.type', distance_transform)
269 
270  ! ------------------------------------------------------------------------ !
271  ! Load the immersed boundary mesh
272 
273  mesh_file = file_t(mesh_file_name)
274  call mesh_file%read(boundary_mesh)
275 
276  if (boundary_mesh%nelv .eq. 0) then
277  call neko_error('No elements in the boundary mesh')
278  end if
279 
280  ! ------------------------------------------------------------------------ !
281  ! Transform the mesh if specified.
282 
283  call json_get_or_default(json, 'mesh_transform.type', &
284  mesh_transform, 'none')
285 
286  select case (mesh_transform)
287  case ('none')
288  ! Do nothing
289  case ('bounding_box')
290  call json_get(json, 'mesh_transform.box_min', box_min)
291  call json_get(json, 'mesh_transform.box_max', box_max)
292  call json_get_or_default(json, 'mesh_transform.keep_aspect_ratio', &
293  keep_aspect_ratio, .true.)
294 
295  if (size(box_min) .ne. 3 .or. size(box_max) .ne. 3) then
296  call neko_error('Case file: mesh_transform. &
297  &box_min and box_max must be 3 element arrays of reals')
298  end if
299 
300  call target_box%init(box_min, box_max)
301 
302  mesh_box = get_aabb(boundary_mesh)
303 
304  scaling = target_box%get_diagonal() / mesh_box%get_diagonal()
305  if (keep_aspect_ratio) then
306  scaling = minval(scaling)
307  end if
308 
309  translation = - scaling * mesh_box%get_min() + target_box%get_min()
310 
311  do idx_p = 1, boundary_mesh%mpts
312  boundary_mesh%points(idx_p)%x = &
313  scaling * boundary_mesh%points(idx_p)%x + translation
314  end do
315 
316  case default
317  call neko_error('Unknown mesh transform')
318  end select
319 
320  ! ------------------------------------------------------------------------ !
321  ! Compute the permeability field
322 
323  ! Assign the signed distance field to all GLL points in the permeability
324  ! field. Initally we just run a brute force loop over all GLL points and
325  ! compute the signed distance function. This should be replaced with a
326  ! more efficient method, such as a tree search.
327 
328  call temp_field%init(this%indicator%dof)
329 
330  ! Select how to transform the distance field to a design field
331  select case (distance_transform)
332  case ('smooth_step')
333  call json_get(json, 'distance_transform.value', scalar_d)
334  scalar_r = real(scalar_d, kind=rp)
335 
336  call signed_distance_field(temp_field, boundary_mesh, scalar_d)
337  call smooth_step_field(temp_field, scalar_r, 0.0_rp)
338 
339  case ('step')
340 
341  call json_get(json, 'distance_transform.value', scalar_d)
342 
343  call signed_distance_field(temp_field, boundary_mesh, scalar_d)
344  call step_function_field(temp_field, scalar_r, 1.0_rp, 0.0_rp)
345 
346  case default
347  call neko_error('Unknown distance transform')
348  end select
349 
350  ! ------------------------------------------------------------------------ !
351  ! Run filter on the temporary indicator field to smooth it out.
352  call json_get_or_default(json, 'filter.type', filter_type, 'none')
353 
354  select case (filter_type)
355  case ('none')
356  ! Do nothing
357  case default
358  call neko_error('Unknown filter type')
359  end select
360 
361  ! Update the global indicator field by max operator
362  this%indicator%x = max(this%indicator%x, temp_field%x)
363 
364  end subroutine init_boundary_mesh
365 
367  subroutine init_point_zone(this, json)
372  use point_zone, only: point_zone_t
375  implicit none
376 
377  class(brinkman_source_term_t), intent(inout) :: this
378  type(json_file), intent(inout) :: json
379 
380  ! Options
381  character(len=:), allocatable :: zone_name
382  character(len=:), allocatable :: filter_type
383 
384  type(field_t) :: temp_field
385  class(point_zone_t), pointer :: my_point_zone
386  integer :: i
387 
388  ! ------------------------------------------------------------------------ !
389  ! Read the options for the point zone
390 
391  call json_get(json, 'name', zone_name)
392  call json_get_or_default(json, 'filter.type', filter_type, 'none')
393 
394  ! Compute the indicator field
395 
396  call temp_field%init(this%indicator%dof)
397 
398  my_point_zone => neko_point_zone_registry%get_point_zone(zone_name)
399 
400  do i = 1, my_point_zone%size
401  temp_field%x(my_point_zone%mask(i), 1, 1, 1) = 1.0_rp
402  end do
403 
404  ! Run filter on the temporary indicator field to smooth it out.
405 
406  select case (filter_type)
407  case ('none')
408  ! Do nothing
409  case default
410  call neko_error('Unknown filter type')
411  end select
412 
413  ! Update the global indicator field by max operator
414  this%indicator%x = max(this%indicator%x, temp_field%x)
415 
416  end subroutine init_point_zone
417 
418 end module brinkman_source_term
double real
Definition: device_config.h:12
Copy data between host and device (or device and device)
Definition: device.F90:51
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:53
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:44
Axis Aligned Bounding Box (aabb) implementation in Fortran.
Definition: aabb.f90:71
type(aabb_t) function, public get_aabb(object, padding)
Construct the aabb of a predefined object.
Definition: aabb.f90:181
Implements the brinkman_source_term_t type.
subroutine init_point_zone(this, json)
Initializes the source term from a point zone.
subroutine init_boundary_mesh(this, json)
Initializes the source term from a boundary mesh.
subroutine brinkman_source_term_compute(this, t, tstep)
Computes the source term and adds the result to fields.
subroutine brinkman_source_term_free(this)
Destructor.
subroutine brinkman_source_term_init_from_json(this, json, fields, coef)
The common constructor using a JSON object.
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 field_subcol3(a, b, c, n)
Returns .
Definition: field_math.f90:475
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines a field.
Definition: field.f90:34
Module for file I/O operations.
Definition: file.f90:34
A module containing filter functions and subroutines. These functions are used to modify fields in a ...
Definition: filters.f90:36
subroutine, public smooth_step_field(F, edge0, edge1)
Apply a smooth step function to a field.
Definition: filters.f90:62
subroutine, public step_function_field(F, x0, value0, value1)
Apply a step function to a field.
Definition: filters.f90:94
subroutine, public permeability_field(F_out, x, k_0, k_1, q)
Apply a permeability function to a field.
Definition: filters.f90:78
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public dp
Definition: num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Profiling interface.
Definition: profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition: profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition: profiler.F90:109
Module containing Signed Distance Functions.
subroutine signed_distance_field(field_data, object, max_distance)
Signed distance field.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Definition: source_term.f90:34
Defines a triangular surface mesh.
Definition: tri_mesh.f90:35
Utilities.
Definition: utils.f90:35
Axis Aligned Bounding Box (aabb) data structure.
Definition: aabb.f90:107
A Brinkman source term. The region and strength are controlled by assigning regions types and brinkma...
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
field_list_t, To be able to group fields together
Definition: field_list.f90:13
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition: file.f90:54
Base abstract type for point zones.
Definition: point_zone.f90:47
Base abstract type for source terms.
Definition: source_term.f90:43
#define max(a, b)
Definition: tensor.cu:40