Neko  0.8.1
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
47  implicit none
48  private
49 
53  type, public, extends(source_term_t) :: brinkman_source_term_t
54  private
55 
57  type(field_t), pointer :: indicator => null()
59  type(field_t), pointer :: brinkman => null()
60  contains
62  procedure, public, pass(this) :: init => brinkman_source_term_init_from_json
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
90  use json_module, only: json_core, json_value
91  implicit none
92 
93  class(brinkman_source_term_t), intent(inout) :: this
94  type(json_file), intent(inout) :: json
95  type(field_list_t), intent(inout), target :: fields
96  type(coef_t), intent(inout) :: coef
97  real(kind=rp) :: start_time, end_time
98 
99  character(len=:), allocatable :: filter_type
100  real(kind=rp), dimension(:), allocatable :: brinkman_limits
101  real(kind=rp) :: brinkman_penalty
102 
103  type(json_value), pointer :: json_object_list
104  type(json_core) :: core
105 
106  character(len=:), allocatable :: object_type
107  type(json_file) :: object_settings
108  integer :: n_regions
109  integer :: i
110 
111  ! Mandatory fields for the general source term
112  call json_get_or_default(json, "start_time", start_time, 0.0_rp)
113  call json_get_or_default(json, "end_time", end_time, huge(0.0_rp))
114 
115  ! Read the options for the permeability field
116  call json_get(json, 'brinkman.limits', brinkman_limits)
117  call json_get(json, 'brinkman.penalty', brinkman_penalty)
118 
119  if (size(brinkman_limits) .ne. 2) then
120  call neko_error('brinkman_limits must be a 2 element array of reals')
121  end if
122 
123  call this%free()
124  call this%init_base(fields, coef, start_time, end_time)
125 
126  ! ------------------------------------------------------------------------ !
127  ! Allocate the permeability and indicator field
128 
129  if (neko_field_registry%field_exists('brinkman_indicator') &
130  .or. neko_field_registry%field_exists('brinkman')) then
131  call neko_error('Brinkman field already exists.')
132  end if
133 
134  call neko_field_registry%add_field(coef%dof, 'brinkman_indicator')
135  this%indicator => neko_field_registry%get_field_by_name('brinkman_indicator')
136 
137  call neko_field_registry%add_field(coef%dof, 'brinkman')
138  this%brinkman => neko_field_registry%get_field_by_name('brinkman')
139 
140  ! ------------------------------------------------------------------------ !
141  ! Select which constructor should be called
142 
143  call json%get('objects', json_object_list)
144  call json%info('objects', n_children=n_regions)
145  call json%get_core(core)
146 
147  do i=1, n_regions
148  call json_extract_item(core, json_object_list, i, object_settings)
149  call json_get_or_default(object_settings, 'type', object_type, 'none')
150 
151  select case (object_type)
152  case ('boundary_mesh')
153  call this%init_boundary_mesh(object_settings)
154  case ('point_zone')
155  call this%init_point_zone(object_settings)
156 
157  case ('none')
158  call object_settings%print()
159  call neko_error('Brinkman source term objects require a region type')
160  case default
161  call neko_error('Brinkman source term unknown region type')
162  end select
163 
164  end do
165 
166  ! Run filter on the full indicator field to smooth it out.
167  call json_get_or_default(json, 'filter.type', filter_type, 'none')
168 
169  select case (filter_type)
170  case ('none')
171  ! Do nothing
172  case default
173  call neko_error('Brinkman source term unknown filter type')
174  end select
175 
176  ! ------------------------------------------------------------------------ !
177  ! Compute the permeability field
178 
179  call permeability_field(this%brinkman, this%indicator, &
180  & brinkman_limits(1), brinkman_limits(2), brinkman_penalty)
181 
182  ! Copy the permeability field to the device
183  if (neko_bcknd_device .eq. 1) then
184  call device_memcpy(this%brinkman%x, this%brinkman%x_d, &
185  this%brinkman%dof%size(), host_to_device, .true.)
186  end if
187 
189 
191  subroutine brinkman_source_term_free(this)
192  class(brinkman_source_term_t), intent(inout) :: this
193 
194  this%brinkman => null()
195  call this%free_base()
196  end subroutine brinkman_source_term_free
197 
201  subroutine brinkman_source_term_compute(this, t, tstep)
203  implicit none
204 
205  class(brinkman_source_term_t), intent(inout) :: this
206  real(kind=rp), intent(in) :: t
207  integer, intent(in) :: tstep
208 
209  if (neko_bcknd_device .eq. 1) then
210  call brinkman_source_term_compute_device(this%fields, this%brinkman)
211  else
212  call brinkman_source_term_compute_cpu(this%fields, this%brinkman)
213  end if
214  end subroutine brinkman_source_term_compute
215 
216  ! ========================================================================== !
217  ! Private methods
218 
220  subroutine init_boundary_mesh(this, json)
221  use file, only: file_t
222  use tri_mesh, only: tri_mesh_t
227  use aabb
228  implicit none
229 
230  class(brinkman_source_term_t), intent(inout) :: this
231  type(json_file), intent(inout) :: json
232 
233  ! Options
234  character(len=:), allocatable :: mesh_file_name
235  character(len=:), allocatable :: distance_transform
236  character(len=:), allocatable :: filter_type
237  character(len=:), allocatable :: mesh_transform
238 
239  ! Read the options for the boundary mesh
240  type(file_t) :: mesh_file
241  type(tri_mesh_t) :: boundary_mesh
242  real(kind=rp) :: scalar_r
243  real(kind=dp) :: scalar_d
244 
245  ! Mesh transform options variables
246  real(kind=dp), dimension(:), allocatable :: box_min, box_max
247  logical :: keep_aspect_ratio
248  real(kind=dp), dimension(3) :: scaling
249  real(kind=dp), dimension(3) :: translation
250  type(field_t) :: temp_field
251  type(aabb_t) :: mesh_box, target_box
252  integer :: idx_p
253 
254  ! ------------------------------------------------------------------------ !
255  ! Read the options for the boundary mesh
256 
257  call json_get(json, 'name', mesh_file_name)
258 
259  ! Settings on how to filter the design field
260  call json_get(json, 'distance_transform.type', distance_transform)
261 
262  ! ------------------------------------------------------------------------ !
263  ! Load the immersed boundary mesh
264 
265  mesh_file = file_t(mesh_file_name)
266  call mesh_file%read(boundary_mesh)
267 
268  if (boundary_mesh%nelv .eq. 0) then
269  call neko_error('No elements in the boundary mesh')
270  end if
271 
272  ! ------------------------------------------------------------------------ !
273  ! Transform the mesh if specified.
274 
275  call json_get_or_default(json, 'mesh_transform.type', mesh_transform, 'none')
276 
277  select case (mesh_transform)
278  case ('none')
279  ! Do nothing
280  case ('bounding_box')
281  call json_get(json, 'mesh_transform.box_min', box_min)
282  call json_get(json, 'mesh_transform.box_max', box_max)
283  call json_get_or_default(json, 'mesh_transform.keep_aspect_ratio', &
284  keep_aspect_ratio, .true.)
285 
286  if (size(box_min) .ne. 3 .or. size(box_max) .ne. 3) then
287  call neko_error('Case file: mesh_transform. &
288  &box_min and box_max must be 3 element arrays of reals')
289  end if
290 
291  call target_box%init(box_min, box_max)
292 
293  mesh_box = get_aabb(boundary_mesh)
294 
295  scaling = target_box%get_diagonal() / mesh_box%get_diagonal()
296  if (keep_aspect_ratio) then
297  scaling = minval(scaling)
298  end if
299 
300  translation = - scaling * mesh_box%get_min() + target_box%get_min()
301 
302  do idx_p = 1, boundary_mesh%mpts
303  boundary_mesh%points(idx_p)%x = &
304  scaling * boundary_mesh%points(idx_p)%x + translation
305  end do
306 
307  case default
308  call neko_error('Unknown mesh transform')
309  end select
310 
311  ! ------------------------------------------------------------------------ !
312  ! Compute the permeability field
313 
314  ! Assign the signed distance field to all GLL points in the permeability
315  ! field. Initally we just run a brute force loop over all GLL points and
316  ! compute the signed distance function. This should be replaced with a
317  ! more efficient method, such as a tree search.
318 
319  call temp_field%init(this%indicator%dof)
320 
321  ! Select how to transform the distance field to a design field
322  select case (distance_transform)
323  case ('smooth_step')
324  call json_get(json, 'distance_transform.value', scalar_d)
325  scalar_r = real(scalar_d, kind=rp)
326 
327  call signed_distance_field(temp_field, boundary_mesh, scalar_d)
328  call smooth_step_field(temp_field, scalar_r, 0.0_rp)
329 
330  case ('step')
331 
332  call json_get(json, 'distance_transform.value', scalar_d)
333 
334  call signed_distance_field(temp_field, boundary_mesh, scalar_d)
335  call step_function_field(temp_field, scalar_r, 1.0_rp, 0.0_rp)
336 
337  case default
338  call neko_error('Unknown distance transform')
339  end select
340 
341  ! ------------------------------------------------------------------------ !
342  ! Run filter on the temporary indicator field to smooth it out.
343  call json_get_or_default(json, 'filter.type', filter_type, 'none')
344 
345  select case (filter_type)
346  case ('none')
347  ! Do nothing
348  case default
349  call neko_error('Unknown filter type')
350  end select
351 
352  ! Update the global indicator field by max operator
353  this%indicator%x = max(this%indicator%x, temp_field%x)
354 
355  end subroutine init_boundary_mesh
356 
358  subroutine init_point_zone(this, json)
362  use point_zone, only: point_zone_t
365  implicit none
366 
367  class(brinkman_source_term_t), intent(inout) :: this
368  type(json_file), intent(inout) :: json
369 
370  ! Options
371  character(len=:), allocatable :: zone_name
372  character(len=:), allocatable :: filter_type
373 
374  type(field_t) :: temp_field
375  class(point_zone_t), pointer :: my_point_zone
376  integer :: i
377 
378  ! ------------------------------------------------------------------------ !
379  ! Read the options for the point zone
380 
381  call json_get(json,'name', zone_name)
382  call json_get_or_default(json, 'filter.type', filter_type, 'none')
383 
384  ! Compute the indicator field
385 
386  call temp_field%init(this%indicator%dof)
387 
388  my_point_zone => neko_point_zone_registry%get_point_zone(zone_name)
389 
390  do i = 1, my_point_zone%size
391  temp_field%x(my_point_zone%mask(i), 1, 1, 1) = 1.0_rp
392  end do
393 
394  ! Run filter on the temporary indicator field to smooth it out.
395 
396  select case (filter_type)
397  case ('none')
398  ! Do nothing
399  case default
400  call neko_error('Unknown filter type')
401  end select
402 
403  ! Update the global indicator field by max operator
404  this%indicator%x = max(this%indicator%x, temp_field%x)
405 
406  end subroutine init_point_zone
407 
408 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 cpu kernel for the brinkman_source_term_t type.
subroutine, public brinkman_source_term_compute_cpu(fields, brinkman)
Computes the Brinkman source term on the cpu.
Implements the device kernel for the brinkman_source_term_t type.
subroutine, public brinkman_source_term_compute_device(fields, brinkman)
Computes the Brinkman source term on the device.
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
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
subroutine, public json_extract_item(core, array, i, item)
Extract ith item from a JSON array as a separate JSON object.
Definition: json_utils.f90:364
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_end_region
End the most recently started profiler region.
Definition: profiler.F90:95
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition: profiler.F90:76
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:54
field_list_t, To be able to group fields together
Definition: field_list.f90:7
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition: file.f90:53
Base abstract type for point zones.
Definition: point_zone.f90:47
Base abstract type for source terms.
Definition: source_term.f90:44
#define max(a, b)
Definition: tensor.cu:40