46 use json_module,
only : json_file, json_core, json_value
83 procedure,
public, pass(this) :: init => &
110 type(json_file),
intent(inout) :: json
112 type(
coef_t),
intent(in),
target :: coef
113 character(len=*),
intent(in) :: variable_name
114 real(kind=
rp) :: start_time, end_time
116 character(len=:),
allocatable :: filter_type
117 real(kind=
rp) :: filter_radius
118 real(kind=
rp),
dimension(:),
allocatable :: brinkman_limits
119 real(kind=
rp) :: brinkman_penalty
121 type(json_value),
pointer :: json_object_list
122 type(json_core) :: core
124 character(len=:),
allocatable :: object_type
125 type(json_file) :: object_settings
136 call json_get(json,
'brinkman.limits', brinkman_limits)
137 call json_get(json,
'brinkman.penalty', brinkman_penalty)
139 if (
size(brinkman_limits) .ne. 2)
then
140 call neko_error(
'brinkman_limits must be a 2 element array of reals')
144 call this%init_base(fields, coef, start_time, end_time)
151 call neko_error(
'Brinkman field already exists.')
154 call this%indicator%init(coef%dof)
155 call this%brinkman%init(coef%dof)
160 call json%get(
'objects', json_object_list)
161 call json%info(
'objects', n_children = n_regions)
162 call json%get_core(core)
168 select case (object_type)
169 case (
'boundary_mesh')
170 call this%init_boundary_mesh(object_settings)
172 call this%init_point_zone(object_settings)
175 call object_settings%print()
176 call neko_error(
'Brinkman source term objects require a region type')
178 call neko_error(
'Brinkman source term unknown region type')
187 select case (filter_type)
190 call this%indicator_unfiltered%init(coef%dof)
196 call this%filter%init(json, coef)
199 call field_copy(this%indicator_unfiltered, this%indicator)
202 call this%filter%apply(this%indicator, this%indicator_unfiltered)
206 call output%fields%assign_to_field(1, this%indicator_unfiltered)
207 call output%fields%assign_to_field(2, this%indicator)
208 call output%fields%assign_to_field(3, this%brinkman)
213 call output%fields%assign_to_field(1, this%indicator)
214 call output%fields%assign_to_field(2, this%brinkman)
217 call neko_error(
'Brinkman source term unknown filter type')
223 this%brinkman = this%indicator
225 brinkman_limits(1), brinkman_limits(2), brinkman_penalty)
228 call output%sample(0.0_rp)
236 call this%indicator%free()
237 call this%brinkman%free()
238 call this%free_base()
247 type(
field_t),
pointer :: u, v, w, fu, fv, fw
250 n = this%fields%item_size(1)
256 fu => this%fields%get(1)
257 fv => this%fields%get(2)
258 fw => this%fields%get(3)
272 type(json_file),
intent(inout) :: json
275 character(len=:),
allocatable :: mesh_file_name
276 character(len=:),
allocatable :: distance_transform
277 character(len=:),
allocatable :: filter_type
278 character(len=:),
allocatable :: mesh_transform
283 real(kind=
rp) :: scalar_r
284 real(kind=
dp) :: scalar_d
287 real(kind=
dp),
dimension(:),
allocatable :: box_min, box_max
288 logical :: keep_aspect_ratio
289 real(kind=
dp),
dimension(3) :: scaling
290 real(kind=
dp),
dimension(3) :: translation
292 type(
aabb_t) :: mesh_box, target_box
294 character(len=LOG_SIZE) :: log_msg
299 call json_get(json,
'name', mesh_file_name)
302 call json_get(json,
'distance_transform.type', distance_transform)
307 call mesh_file%init(mesh_file_name)
308 call mesh_file%read(boundary_mesh)
310 if (boundary_mesh%nelv .eq. 0)
then
311 call neko_error(
'No elements in the boundary mesh')
318 mesh_transform,
'none')
320 select case (mesh_transform)
323 case (
'bounding_box')
324 call json_get(json,
'mesh_transform.box_min', box_min)
325 call json_get(json,
'mesh_transform.box_max', box_max)
327 keep_aspect_ratio, .true.)
329 if (
size(box_min) .ne. 3 .or.
size(box_max) .ne. 3)
then
331 &box_min and box_max must be 3 element arrays of reals')
334 call target_box%init(box_min, box_max)
338 scaling = target_box%get_diagonal() / mesh_box%get_diagonal()
339 if (keep_aspect_ratio)
then
340 scaling = minval(scaling)
343 translation = - scaling * mesh_box%get_min() + target_box%get_min()
345 do idx_p = 1, boundary_mesh%mpts
346 boundary_mesh%points(idx_p)%x = &
347 scaling * boundary_mesh%points(idx_p)%x + translation
351 write(log_msg,
'(A)')
"The following transformation was applied:"
353 write(log_msg,
'(A, 3F12.6)')
"Scaling: ", scaling
355 write(log_msg,
'(A, 3F12.6)')
"Translation: ", translation
370 call temp_field%init(this%coef%dof)
373 select case (distance_transform)
375 call json_get(json,
'distance_transform.value', scalar_d)
376 scalar_r =
real(scalar_d, kind=
rp)
383 call json_get(json,
'distance_transform.value', scalar_d)
384 scalar_r =
real(scalar_d, kind=
rp)
396 this%indicator%size())
398 this%indicator%x =
max(this%indicator%x, temp_field%x)
406 type(json_file),
intent(inout) :: json
409 character(len=:),
allocatable :: zone_name
418 call json_get(json,
'name', zone_name)
421 call temp_field%init(this%coef%dof)
427 zone%mask%get_d(), zone%size)
429 call cfill_mask(temp_field%x, 1.0_rp, temp_field%size(), &
430 zone%mask%get(), zone%size)
436 this%indicator%size())
438 this%indicator%x =
max(this%indicator%x, temp_field%x)
Copy data between host and device (or device and device)
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.
type(aabb_t) function, public get_aabb(object, padding)
Construct the aabb of a predefined object.
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_init_from_json(this, json, fields, coef, variable_name)
The common constructor using a JSON object.
subroutine brinkman_source_term_free(this)
Destructor.
subroutine brinkman_source_term_compute(this, time)
Computes the source term and adds the result to fields.
subroutine, public device_cfill_mask(a_d, c, n, mask_d, n_mask, strm)
Fill a constant to a masked vector. .
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
subroutine, public field_copy(a, b, n)
Copy a vector .
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.
Module for file I/O operations.
Filter to be applied to a scalar field.
Implements fld_file_output_t.
Utilities for retrieving parameters from the case files.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
A module containing mapping functions and subroutines. These functions are used to modify fields in a...
subroutine, public permeability_field(f, k_0, k_1, q)
Apply a permeability function to a field.
subroutine, public smooth_step_field(f, edge0, edge1)
Apply a smooth step function to a field.
subroutine, public step_function_field(f, x0, value0, value1)
Apply a step function to a field.
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
integer, parameter neko_bcknd_device
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
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.
Module with things related to the simulation time.
Defines a triangular surface mesh.
Axis Aligned Bounding Box (aabb) data structure.
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,...
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...
Base abstract class for filter.
A simple output saving a list of fields to a .fld file.
A PDE based filter mapping $\rho \mapsto \tilde{\rho}$, see Lazarov & O. Sigmund 2010,...
Base abstract type for point zones.
Base abstract type for source terms.
A struct that contains all info about the time, expand as needed.