48 use json_module,
only : json_file, json_core, json_value
81 type(
field_t),
pointer :: unfiltered => null()
83 type(
field_t),
pointer :: indicator => null()
85 type(
field_t),
pointer :: brinkman => null()
91 procedure,
public, pass(this) :: init => &
122 type(json_file),
intent(inout) :: json
124 type(
coef_t),
intent(in),
target :: coef
125 character(len=*),
intent(in) :: variable_name
126 real(kind=
rp) :: start_time, end_time
127 character(len=LOG_SIZE) :: error_msg
129 character(len=:),
allocatable :: filter_type
130 real(kind=
rp),
dimension(:),
allocatable :: brinkman_limits
131 real(kind=
rp) :: brinkman_penalty
133 character(len=:),
allocatable :: object_type
134 type(json_file) :: object_settings, filter_settings
135 integer :: n_objects, i
139 logical :: output_enable
140 character(len=:),
allocatable :: output_path, output_format, &
141 output_precision, fname, suffix
147 call neko_log%section(
'Brinkman source term')
163 if (
size(brinkman_limits) .ne. 2)
then
164 call neko_error(
'brinkman_limits must be a 2 element array of reals')
167 if (output_path(len_trim(output_path):len_trim(output_path)) .ne.
'/')
then
168 output_path = trim(output_path) //
'/'
172 call this%init_base(fields, coef, start_time, end_time)
177 call neko_registry%add_field(coef%dof,
'brinkman_indicator', .true.)
178 call neko_registry%add_field(coef%dof,
'brinkman_permeability', .true.)
180 this%indicator =>
neko_registry%get_field(
'brinkman_indicator')
181 this%brinkman =>
neko_registry%get_field(
'brinkman_permeability')
186 call json%info(
'objects', n_children = n_objects)
190 call json_get(object_settings,
'type', object_type)
192 select case (object_type)
193 case (
'boundary_mesh')
194 call this%add_boundary_mesh(object_settings)
196 call this%add_point_zone(object_settings)
198 call this%add_field(object_settings)
200 call this%add_file(object_settings)
203 write(error_msg,
'(A,I0,A,A)')
'Brinkman object: ', i, &
204 ' unknown type: ', trim(object_type)
214 select case (filter_type)
218 call neko_registry%add_field(coef%dof,
'brinkman_unfiltered', .true.)
219 this%unfiltered =>
neko_registry%get_field(
'brinkman_unfiltered')
220 call field_copy(this%unfiltered, this%indicator)
224 call json_get(json,
'filter', filter_settings)
225 call this%filter%init(filter_settings, coef)
226 call this%filter%apply(this%indicator, this%unfiltered)
231 call neko_error(
'Brinkman source term unknown filter type')
237 this%brinkman = this%indicator
239 brinkman_limits(1), brinkman_limits(2), brinkman_penalty)
244 if (output_enable)
then
245 fname = trim(output_path) //
'brinkman'
246 select case (trim(output_format))
254 suffix =
'.' // trim(output_format)
256 fname = trim(fname) // trim(suffix)
258 select case (trim(output_precision))
259 case (
'sp',
'single')
261 case (
'dp',
'double')
267 call output%init(fname, precision = precision)
269 call output_fields%init(2)
270 call output_fields%assign_to_field(1, this%indicator)
271 call output_fields%assign_to_field(2, this%brinkman)
273 call neko_log%message(
'Brinkman output')
274 call neko_log%message(
' 1: Indicator')
275 call neko_log%message(
' 2: Permeability')
276 if (
associated(this%unfiltered))
then
277 call output_fields%append(this%unfiltered)
278 call neko_log%message(
' 3: Unfiltered Indicator')
281 call output%write(output_fields)
284 call output_fields%free()
294 if (
associated(this%indicator))
nullify(this%indicator)
295 if (
associated(this%unfiltered))
nullify(this%unfiltered)
296 if (
associated(this%brinkman))
nullify(this%brinkman)
298 if (
allocated(this%filter))
then
299 call this%filter%free()
300 deallocate(this%filter)
302 call this%free_base()
311 type(
field_t),
pointer :: u, v, w, fu, fv, fw
314 n = this%fields%item_size(1)
320 fu => this%fields%get(1)
321 fv => this%fields%get(2)
322 fw => this%fields%get(3)
336 type(json_file),
intent(inout) :: json
339 character(len=:),
allocatable :: mesh_file_name
340 character(len=:),
allocatable :: distance_transform
341 character(len=:),
allocatable :: filter_type
342 character(len=:),
allocatable :: mesh_transform
347 real(kind=
rp) :: scalar_r
348 real(kind=
dp) :: scalar_d
349 logical :: cache, cache_exist
350 character(len=:),
allocatable :: cache_filename
351 type(
file_t) :: cache_file
359 real(kind=
dp),
dimension(:),
allocatable :: box_min, box_max
360 logical :: keep_aspect_ratio
361 real(kind=
dp),
dimension(3) :: scaling
362 real(kind=
dp),
dimension(3) :: translation
363 type(
field_t),
pointer :: temp_field
364 type(
aabb_t) :: mesh_box, target_box
365 integer :: idx_p, temp_idx
366 character(len=LOG_SIZE) :: log_msg
371 call json_get(json,
'name', mesh_file_name)
375 call json_get(json,
'distance_transform.type', distance_transform)
380 call json_get(json,
'cache_file', cache_filename)
382 inquire(
file=trim(cache_filename) //
"0.nek5000", exist=cache_exist)
383 write(log_msg,
'(A)')
"Checking for Brinkman source term cache."
386 if (cache_exist)
then
387 write(log_msg,
'(A)')
"Loading Brinkman source term from cache."
392 call cache_data%init()
393 call cache_file%init(cache_filename //
"0.fld")
394 call cache_file%set_counter(0)
395 call cache_file%read(cache_data)
396 call cache_data%import_fields(p = temp_field)
403 call cache_data%free()
404 call cache_file%free()
412 call mesh_file%init(mesh_file_name)
413 call mesh_file%read(boundary_mesh)
415 if (boundary_mesh%nelv .eq. 0)
then
416 call neko_error(
'No elements in the boundary mesh')
423 mesh_transform,
'none')
425 select case (mesh_transform)
428 case (
'bounding_box')
432 keep_aspect_ratio, .true.)
434 if (
size(box_min) .ne. 3 .or.
size(box_max) .ne. 3)
then
436 &box_min and box_max must be 3 element arrays of reals')
439 call target_box%init(box_min, box_max)
443 scaling = target_box%get_diagonal() / mesh_box%get_diagonal()
444 if (keep_aspect_ratio)
then
445 scaling = minval(scaling)
448 translation = - scaling * mesh_box%get_min() + target_box%get_min()
450 do idx_p = 1, boundary_mesh%mpts
451 boundary_mesh%points(idx_p)%x = &
452 scaling * boundary_mesh%points(idx_p)%x + translation
456 write(log_msg,
'(A)')
"The following transformation was applied:"
458 write(log_msg,
'(A, 3F12.6)')
"Scaling: ", scaling
460 write(log_msg,
'(A, 3F12.6)')
"Translation: ", translation
478 select case (distance_transform)
481 scalar_r =
real(scalar_d, kind=
rp)
489 scalar_r =
real(scalar_d, kind=
rp)
500 write(log_msg,
'(A)')
"Writing Brinkman source term to cache."
502 call cache_output%init(
dp, cache_filename, 1)
503 call cache_output%fields%assign_to_field(1, temp_field)
504 call cache_output%sample(0.0_rp)
517 type(json_file),
intent(inout) :: json
520 character(len=:),
allocatable :: zone_name
522 type(
field_t),
pointer :: temp_field
524 integer :: i, temp_idx
529 call json_get(json,
'name', zone_name)
538 zone%mask%get_d(), zone%size)
540 call cfill_mask(temp_field%x, 1.0_rp, temp_field%size(), &
541 zone%mask%get(), zone%size)
553 type(json_file),
intent(inout) :: json
555 type(
field_t),
pointer :: temp_field
556 character(len=:),
allocatable :: file_name, field_name, tmp_str
557 character(len=80) :: suffix
558 integer :: file_idx, temp_idx
560 call json_get(json,
'file_name', file_name)
562 'brinkman_indicator')
567 call file%init(file_name)
568 call file%set_counter(file_idx)
571 select case (trim(suffix))
579 call file%read(fld_data)
580 select case (field_name(1:1))
582 call fld_data%import_fields(p = temp_field)
584 call fld_data%import_fields(u = temp_field)
586 call fld_data%import_fields(v = temp_field)
588 call fld_data%import_fields(w = temp_field)
590 call fld_data%import_fields(t = temp_field)
593 if (len_trim(field_name) .eq. 3)
then
594 read(field_name(2:3),
'(I2)') idx(1)
595 else if (len_trim(field_name) .eq. 2)
then
596 read(field_name(2:2),
'(I1)') idx(1)
598 call neko_error(
'For fields with prefix s, the field name ' // &
599 'must be in the format sXX, where XX is the index of ' // &
600 'the field in the fld file')
603 call fld_fields%init(1)
604 call fld_fields%assign(1, temp_field)
606 call fld_data%import_fields(s_target_list = fld_fields, &
609 call fld_fields%free()
611 call neko_error(
'Unknown field prefix in field name: ' // &
621 tmp_str = trim(temp_field%name)
622 temp_field%name = trim(field_name)
624 call file%read(temp_field)
626 temp_field%name = trim(tmp_str)
629 call neko_error(
"Brinkman cannot read file: " // trim(file_name))
643 type(json_file),
intent(inout) :: json
644 character(len=:),
allocatable :: field_name
645 type(field_t),
pointer :: temp_field
648 call json_get(json,
'name', field_name)
649 temp_field => neko_registry%get_field(field_name)
652 call field_pwmax2(this%indicator, temp_field)
654 if (
associated(temp_field))
nullify(temp_field)
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 add_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 add_point_zone(this, json)
Initializes the source term from a point zone.
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 add_field(this, json)
Initializes the source term from a field.
subroutine add_file(this, json)
Initializes the source term from a file.
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
integer, parameter, public device_to_host
subroutine, public field_pwmax2(a, b, n)
Point-wise max operation .
subroutine, public field_copy(a, b, n)
Copy a vector .
subroutine, public field_subcol3(a, b, c, n)
Returns .
Implements field_output_t.
Module for file I/O operations.
Filter to be applied to a scalar field.
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
Implements fld_file_output_t.
Implements global_interpolation given a dofmap.
Routines to interpolate between different spaces.
Utilities for retrieving parameters from the case files.
integer, parameter, public neko_log_debug
Debug log level.
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 copy(a, b, n)
Copy a vector .
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.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
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 function space.
integer, parameter, public gll
Module with things related to the simulation time.
Defines a triangular surface mesh.
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
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 simple output saving a list of fields to a file.
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.
Implements global interpolation for arbitrary points in the domain.
Interpolation between two space::space_t.
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .
Base abstract type for point zones.
Base abstract type for source terms.
The function space for the SEM solution fields.
A struct that contains all info about the time, expand as needed.