Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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 aabb, only : aabb_t, get_aabb
36 use coefs, only: coef_t
39 use field, only: field_t
40 use field_list, only: field_list_t
41 use field_math, only: field_subcol3
44 use file, only: file_t
45 use json_module, only: json_file, json_core, json_value
47 use logger, only: neko_log, log_size
48 use math, only: pwmax, cfill_mask
49 use tri_mesh, only: tri_mesh_t
51 use num_types, only: rp, dp
52 use point_zone, only: point_zone_t
56 use source_term, only: source_term_t
57 use utils, only: neko_error
58 implicit none
59 private
60
64 type, public, extends(source_term_t) :: brinkman_source_term_t
65 private
66
68 type(field_t) :: indicator
70 type(field_t) :: brinkman
71 contains
73 procedure, public, pass(this) :: init => &
76 procedure, public, pass(this) :: free => brinkman_source_term_free
78 procedure, public, pass(this) :: compute_ => brinkman_source_term_compute
79
80 ! ----------------------------------------------------------------------- !
81 ! Private methods
82 procedure, pass(this) :: init_boundary_mesh
83 procedure, pass(this) :: init_point_zone
85
86contains
87
88 ! ========================================================================== !
89 ! Public methods
90
95 subroutine brinkman_source_term_init_from_json(this, json, fields, coef)
96 class(brinkman_source_term_t), intent(inout) :: this
97 type(json_file), intent(inout) :: json
98 type(field_list_t), intent(in), target :: fields
99 type(coef_t), intent(in), target :: coef
100 real(kind=rp) :: start_time, end_time
101
102 character(len=:), allocatable :: filter_type
103 real(kind=rp), dimension(:), allocatable :: brinkman_limits
104 real(kind=rp) :: brinkman_penalty
105
106 type(json_value), pointer :: json_object_list
107 type(json_core) :: core
108
109 character(len=:), allocatable :: object_type
110 type(json_file) :: object_settings
111 integer :: n_regions
112 integer :: i
113
114 ! Mandatory fields for the general source term
115 call json_get_or_default(json, "start_time", start_time, 0.0_rp)
116 call json_get_or_default(json, "end_time", end_time, huge(0.0_rp))
117
118 ! Read the options for the permeability field
119 call json_get(json, 'brinkman.limits', brinkman_limits)
120 call json_get(json, 'brinkman.penalty', brinkman_penalty)
121
122 if (size(brinkman_limits) .ne. 2) then
123 call neko_error('brinkman_limits must be a 2 element array of reals')
124 end if
125
126 call this%free()
127 call this%init_base(fields, coef, start_time, end_time)
128
129 ! ------------------------------------------------------------------------ !
130 ! Allocate the permeability and indicator field
131
132 if (neko_field_registry%field_exists('brinkman_indicator') &
133 .or. neko_field_registry%field_exists('brinkman')) then
134 call neko_error('Brinkman field already exists.')
135 end if
136
137 call this%indicator%init(coef%dof)
138 call this%brinkman%init(coef%dof)
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 this%brinkman = this%indicator
180 call permeability_field(this%brinkman, &
181 brinkman_limits(1), brinkman_limits(2), brinkman_penalty)
182
184
187 class(brinkman_source_term_t), intent(inout) :: this
188
189 call this%indicator%free()
190 call this%brinkman%free()
191 call this%free_base()
192
193 end subroutine brinkman_source_term_free
194
198 subroutine brinkman_source_term_compute(this, t, tstep)
199 class(brinkman_source_term_t), intent(inout) :: this
200 real(kind=rp), intent(in) :: t
201 integer, intent(in) :: tstep
202 type(field_t), pointer :: u, v, w, fu, fv, fw
203 integer :: n
204
205 n = this%fields%item_size(1)
206
207 u => neko_field_registry%get_field('u')
208 v => neko_field_registry%get_field('v')
209 w => neko_field_registry%get_field('w')
210
211 fu => this%fields%get(1)
212 fv => this%fields%get(2)
213 fw => this%fields%get(3)
214
215 call field_subcol3(fu, u, this%brinkman, n)
216 call field_subcol3(fv, v, this%brinkman, n)
217 call field_subcol3(fw, w, this%brinkman, n)
218
219 end subroutine brinkman_source_term_compute
220
221 ! ========================================================================== !
222 ! Private methods
223
225 subroutine init_boundary_mesh(this, json)
226 class(brinkman_source_term_t), intent(inout) :: this
227 type(json_file), intent(inout) :: json
228
229 ! Options
230 character(len=:), allocatable :: mesh_file_name
231 character(len=:), allocatable :: distance_transform
232 character(len=:), allocatable :: filter_type
233 character(len=:), allocatable :: mesh_transform
234
235 ! Read the options for the boundary mesh
236 type(file_t) :: mesh_file
237 type(tri_mesh_t) :: boundary_mesh
238 real(kind=rp) :: scalar_r
239 real(kind=dp) :: scalar_d
240
241 ! Mesh transform options variables
242 real(kind=dp), dimension(:), allocatable :: box_min, box_max
243 logical :: keep_aspect_ratio
244 real(kind=dp), dimension(3) :: scaling
245 real(kind=dp), dimension(3) :: translation
246 type(field_t) :: temp_field
247 type(aabb_t) :: mesh_box, target_box
248 integer :: idx_p
249 character(len=LOG_SIZE) :: log_msg
250
251 ! ------------------------------------------------------------------------ !
252 ! Read the options for the boundary mesh
253
254 call json_get(json, 'name', mesh_file_name)
255
256 ! Settings on how to filter the design field
257 call json_get(json, 'distance_transform.type', distance_transform)
258
259 ! ------------------------------------------------------------------------ !
260 ! Load the immersed boundary mesh
261
262 mesh_file = file_t(mesh_file_name)
263 call mesh_file%read(boundary_mesh)
264
265 if (boundary_mesh%nelv .eq. 0) then
266 call neko_error('No elements in the boundary mesh')
267 end if
268
269 ! ------------------------------------------------------------------------ !
270 ! Transform the mesh if specified.
271
272 call json_get_or_default(json, 'mesh_transform.type', &
273 mesh_transform, 'none')
274
275 select case (mesh_transform)
276 case ('none')
277 ! Do nothing
278 case ('bounding_box')
279 call json_get(json, 'mesh_transform.box_min', box_min)
280 call json_get(json, 'mesh_transform.box_max', box_max)
281 call json_get_or_default(json, 'mesh_transform.keep_aspect_ratio', &
282 keep_aspect_ratio, .true.)
283
284 if (size(box_min) .ne. 3 .or. size(box_max) .ne. 3) then
285 call neko_error('Case file: mesh_transform. &
286 &box_min and box_max must be 3 element arrays of reals')
287 end if
288
289 call target_box%init(box_min, box_max)
290
291 mesh_box = get_aabb(boundary_mesh)
292
293 scaling = target_box%get_diagonal() / mesh_box%get_diagonal()
294 if (keep_aspect_ratio) then
295 scaling = minval(scaling)
296 end if
297
298 translation = - scaling * mesh_box%get_min() + target_box%get_min()
299
300 do idx_p = 1, boundary_mesh%mpts
301 boundary_mesh%points(idx_p)%x = &
302 scaling * boundary_mesh%points(idx_p)%x + translation
303 end do
304
305 ! Report the transformation applied
306 write(log_msg, '(A)') "The following transformation was applied:"
307 call neko_log%message(log_msg)
308 write(log_msg, '(A, 3F12.6)') "Scaling: ", scaling
309 call neko_log%message(log_msg)
310 write(log_msg, '(A, 3F12.6)') "Translation: ", translation
311 call neko_log%message(log_msg)
312
313 case default
314 call neko_error('Unknown mesh transform')
315 end select
316
317 ! ------------------------------------------------------------------------ !
318 ! Compute the permeability field
319
320 ! Assign the signed distance field to all GLL points in the permeability
321 ! field. Initally we just run a brute force loop over all GLL points and
322 ! compute the signed distance function. This should be replaced with a
323 ! more efficient method, such as a tree search.
324
325 call temp_field%init(this%coef%dof)
326
327 ! Select how to transform the distance field to a design field
328 select case (distance_transform)
329 case ('smooth_step')
330 call json_get(json, 'distance_transform.value', scalar_d)
331 scalar_r = real(scalar_d, kind=rp)
332
333 call signed_distance_field(temp_field, boundary_mesh, scalar_d)
334 call smooth_step_field(temp_field, scalar_r, 0.0_rp)
335
336 case ('step')
337
338 call json_get(json, 'distance_transform.value', scalar_d)
339 scalar_r = real(scalar_d, kind=rp)
340
341 call signed_distance_field(temp_field, boundary_mesh, scalar_d)
342 call step_function_field(temp_field, scalar_r, 1.0_rp, 0.0_rp)
343
344 case default
345 call neko_error('Unknown distance transform')
346 end select
347
348 ! ------------------------------------------------------------------------ !
349 ! Run filter on the temporary indicator field to smooth it out.
350 call json_get_or_default(json, 'filter.type', filter_type, 'none')
351
352 select case (filter_type)
353 case ('none')
354 ! Do nothing
355 case default
356 call neko_error('Unknown filter type')
357 end select
358
359 ! Update the global indicator field by max operator
360 if (neko_bcknd_device .eq. 1) then
361 call device_pwmax(this%indicator%x_d, temp_field%x_d, &
362 this%indicator%size())
363 else
364 this%indicator%x = max(this%indicator%x, temp_field%x)
365 end if
366
367 end subroutine init_boundary_mesh
368
370 subroutine init_point_zone(this, json)
371 class(brinkman_source_term_t), intent(inout) :: this
372 type(json_file), intent(inout) :: json
373
374 ! Options
375 character(len=:), allocatable :: zone_name
376 character(len=:), allocatable :: filter_type
377
378 type(field_t) :: temp_field
379 class(point_zone_t), pointer :: my_point_zone
380 integer :: i
381
382 ! ------------------------------------------------------------------------ !
383 ! Read the options for the point zone
384
385 call json_get(json, 'name', zone_name)
386 call json_get_or_default(json, 'filter.type', filter_type, 'none')
387
388 ! Compute the indicator field
389 call temp_field%init(this%coef%dof)
390
391 my_point_zone => neko_point_zone_registry%get_point_zone(zone_name)
392
393 if (neko_bcknd_device .eq. 1) then
394 call device_cfill_mask(temp_field%x_d, 1.0_rp, temp_field%size(), &
395 my_point_zone%mask_d, my_point_zone%size)
396 else
397 call cfill_mask(temp_field%x, 1.0_rp, temp_field%size(), &
398 my_point_zone%mask, my_point_zone%size)
399 end if
400
401 ! Run filter on the temporary indicator field to smooth it out.
402
403 select case (filter_type)
404 case ('none')
405 ! Do nothing
406 case default
407 call neko_error('Unknown filter type')
408 end select
409
410 ! Update the global indicator field by max operator
411 if (neko_bcknd_device .eq. 1) then
412 call device_pwmax(this%indicator%x_d, temp_field%x_d, &
413 this%indicator%size())
414 else
415 this%indicator%x = max(this%indicator%x, temp_field%x)
416 end if
417
418 end subroutine init_point_zone
419
420end module brinkman_source_term
double real
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...
Retrieves a parameter by name or throws an error.
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
subroutine, public device_cfill_mask(a_d, c, size, mask_d, mask_size)
Fill a constant to a masked vector. .
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 .
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 permeability_field(f, k_0, k_1, q)
Apply a permeability function to a field.
Definition filters.f90:90
subroutine, public smooth_step_field(f, edge0, edge1)
Apply a smooth step function to a field.
Definition filters.f90:70
subroutine, public step_function_field(f, x0, value0, value1)
Apply a step function to a field.
Definition filters.f90:109
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:65
integer, parameter, public log_size
Definition log.f90:42
Definition math.f90:60
subroutine, public cfill_mask(a, c, size, mask, mask_size)
Fill a constant to a masked vector. .
Definition math.f90:296
Build configurations.
integer, parameter neko_bcknd_device
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, public signed_distance_field(field_data, object, max_distance)
Signed distance field.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
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
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.
Base abstract type for source terms.
#define max(a, b)
Definition tensor.cu:40