Neko 1.99.3
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-2026, 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
38 use field, only : field_t
39 use field_list, only : field_list_t
40 use math, only : cfill_mask, copy
43 use registry, only : neko_registry
47 use file, only : file_t
48 use json_module, only : json_file, json_core, json_value
52 use tri_mesh, only : tri_mesh_t
54 use num_types, only : rp, dp, sp
55 use point_zone, only : point_zone_t
59 use source_term, only : source_term_t
61 use filter, only : filter_t
62 use pde_filter, only : pde_filter_t
66 use time_state, only : time_state_t
67
70 use space, only: space_t, gll
71 implicit none
72 private
73
77 type, public, extends(source_term_t) :: brinkman_source_term_t
78 private
79
81 type(field_t), pointer :: unfiltered => null()
83 type(field_t), pointer :: indicator => null()
85 type(field_t), pointer :: brinkman => null()
87 class(filter_t), allocatable :: filter
88
89 contains
91 procedure, public, pass(this) :: init => &
94 procedure, public, pass(this) :: free => brinkman_source_term_free
96 procedure, public, pass(this) :: compute_ => brinkman_source_term_compute
97
98 ! ----------------------------------------------------------------------- !
99 ! Internal subroutines for adding a given object type to the Brinkman
100 ! source term
101
102 procedure, pass(this) :: add_boundary_mesh
103 procedure, pass(this) :: add_point_zone
104 procedure, pass(this) :: add_file
105 procedure, pass(this) :: add_field
106
108
109contains
110
111 ! ========================================================================== !
112 ! Public methods
113
119 subroutine brinkman_source_term_init_from_json(this, json, fields, coef, &
120 variable_name)
121 class(brinkman_source_term_t), intent(inout) :: this
122 type(json_file), intent(inout) :: json
123 type(field_list_t), intent(in), target :: fields
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
128
129 character(len=:), allocatable :: filter_type
130 real(kind=rp), dimension(:), allocatable :: brinkman_limits
131 real(kind=rp) :: brinkman_penalty
132
133 character(len=:), allocatable :: object_type
134 type(json_file) :: object_settings, filter_settings
135 integer :: n_objects, i
136
137 type(file_t) :: output
138 type(field_list_t) :: output_fields
139 logical :: output_enable
140 character(len=:), allocatable :: output_path, output_format, &
141 output_precision, fname, suffix
142 integer :: precision
143
144 ! ------------------------------------------------------------------------ !
145 ! Read the general options for the Brinkman source term
146
147 call neko_log%section('Brinkman source term')
148
149 ! Read the options for the permeability field
150 call json_get_or_lookup(json, 'limits', brinkman_limits)
151 call json_get_or_lookup(json, 'penalty', brinkman_penalty)
152
153 ! Mandatory fields for the general source term
154 call json_get_or_lookup_or_default(json, "start_time", start_time, 0.0_rp)
155 call json_get_or_lookup_or_default(json, "end_time", end_time, huge(0.0_rp))
156
157 ! Output settings
158 call json_get_or_default(json, 'output.enable', output_enable, .false.)
159 call json_get_or_default(json, 'output.precision', output_precision, 'sp')
160 call json_get_or_default(json, 'output.path', output_path, './')
161 call json_get_or_default(json, 'output.format', output_format, 'fld')
162
163 if (size(brinkman_limits) .ne. 2) then
164 call neko_error('brinkman_limits must be a 2 element array of reals')
165 end if
166
167 if (output_path(len_trim(output_path):len_trim(output_path)) .ne. '/') then
168 output_path = trim(output_path) // '/'
169 end if
170
171 call this%free()
172 call this%init_base(fields, coef, start_time, end_time)
173
174 ! ------------------------------------------------------------------------ !
175 ! Allocate the permeability and indicator field
176
177 call neko_registry%add_field(coef%dof, 'brinkman_indicator', .true.)
178 call neko_registry%add_field(coef%dof, 'brinkman_permeability', .true.)
179
180 this%indicator => neko_registry%get_field('brinkman_indicator')
181 this%brinkman => neko_registry%get_field('brinkman_permeability')
182
183 ! ------------------------------------------------------------------------ !
184 ! Select which constructor should be called
185
186 call json%info('objects', n_children = n_objects)
187
188 do i = 1, n_objects
189 call json_extract_item(json, "objects", i, object_settings)
190 call json_get(object_settings, 'type', object_type)
191
192 select case (object_type)
193 case ('boundary_mesh')
194 call this%add_boundary_mesh(object_settings)
195 case ('point_zone')
196 call this%add_point_zone(object_settings)
197 case ('field')
198 call this%add_field(object_settings)
199 case ('file')
200 call this%add_file(object_settings)
201
202 case default
203 write(error_msg, '(A,I0,A,A)') 'Brinkman object: ', i, &
204 ' unknown type: ', trim(object_type)
205 call neko_error(trim(error_msg))
206 end select
207
208 end do
209
210 ! ------------------------------------------------------------------------ !
211 ! Filter the indicator field
212
213 call json_get_or_default(json, 'filter.type', filter_type, 'none')
214 select case (filter_type)
215 case ('PDE')
216
217 ! Copy the current indicator to unfiltered (essentially a rename)
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)
221
222 ! Allocate a PDE filter
223 allocate(pde_filter_t::this%filter)
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)
227
228 case ('none')
229 ! Do nothing
230 case default
231 call neko_error('Brinkman source term unknown filter type')
232 end select
233
234 ! ------------------------------------------------------------------------ !
235 ! Compute the permeability field
236
237 this%brinkman = this%indicator
238 call permeability_field(this%brinkman, &
239 brinkman_limits(1), brinkman_limits(2), brinkman_penalty)
240
241 ! ------------------------------------------------------------------------ !
242 ! Set up output the brinkman fields if enabled
243
244 if (output_enable) then
245 fname = trim(output_path) // 'brinkman'
246 select case (trim(output_format))
247 case ('nek5000')
248 suffix = '.fld'
249 case ('adios2')
250 suffix = '.bp'
251 case ('hdf5')
252 suffix = '.h5'
253 case default
254 suffix = '.' // trim(output_format)
255 end select
256 fname = trim(fname) // trim(suffix)
257
258 select case (trim(output_precision))
259 case ('sp', 'single')
260 precision = sp
261 case ('dp', 'double')
262 precision = dp
263 case default
264 call neko_error('Unknown output precision')
265 end select
266
267 call output%init(fname, precision = precision)
268
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)
272
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')
279 end if
280
281 call output%write(output_fields)
282
283 call output%free()
284 call output_fields%free()
285 end if
286
287 call neko_log%end_section()
289
292 class(brinkman_source_term_t), intent(inout) :: this
293
294 if (associated(this%indicator)) nullify(this%indicator)
295 if (associated(this%unfiltered)) nullify(this%unfiltered)
296 if (associated(this%brinkman)) nullify(this%brinkman)
297
298 if (allocated(this%filter)) then
299 call this%filter%free()
300 deallocate(this%filter)
301 end if
302 call this%free_base()
303
304 end subroutine brinkman_source_term_free
305
308 subroutine brinkman_source_term_compute(this, time)
309 class(brinkman_source_term_t), intent(inout) :: this
310 type(time_state_t), intent(in) :: time
311 type(field_t), pointer :: u, v, w, fu, fv, fw
312 integer :: n
313
314 n = this%fields%item_size(1)
315
316 u => neko_registry%get_field('u')
317 v => neko_registry%get_field('v')
318 w => neko_registry%get_field('w')
319
320 fu => this%fields%get(1)
321 fv => this%fields%get(2)
322 fw => this%fields%get(3)
323
324 call field_subcol3(fu, u, this%brinkman, n)
325 call field_subcol3(fv, v, this%brinkman, n)
326 call field_subcol3(fw, w, this%brinkman, n)
327
328 end subroutine brinkman_source_term_compute
329
330 ! ========================================================================== !
331 ! Private methods
332
334 subroutine add_boundary_mesh(this, json)
335 class(brinkman_source_term_t), intent(inout) :: this
336 type(json_file), intent(inout) :: json
337
338 ! Options
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
343
344 ! Read the options for the boundary mesh
345 type(file_t) :: mesh_file
346 type(tri_mesh_t) :: boundary_mesh
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
352 type(fld_file_output_t) :: cache_output
353 type(fld_file_data_t) :: cache_data
354 type(global_interpolation_t) :: global_interp
355 type(space_t) :: prev_Xh
356 type(interpolator_t) :: space_interp
357
358 ! Mesh transform options variables
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
367
368 ! ------------------------------------------------------------------------ !
369 ! Read the options for the boundary mesh
370
371 call json_get(json, 'name', mesh_file_name)
372 call json_get_or_default(json, 'cache', cache, .false.)
373
374 ! Settings on how to filter the design field
375 call json_get(json, 'distance_transform.type', distance_transform)
376
377 ! ------------------------------------------------------------------------ !
378 ! Check if we can load from cache
379 if (cache) then
380 call json_get(json, 'cache_file', cache_filename)
381
382 inquire(file=trim(cache_filename) // "0.nek5000", exist=cache_exist)
383 write(log_msg, '(A)') "Checking for Brinkman source term cache."
384 call neko_log%message(log_msg, neko_log_debug)
385
386 if (cache_exist) then
387 write(log_msg, '(A)') "Loading Brinkman source term from cache."
388 call neko_log%message(log_msg, neko_log_debug)
389
390 call neko_scratch_registry%request_field(temp_field, temp_idx, .true.)
391
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)
397
398 ! Update the global indicator field by max operator
399 call field_pwmax2(this%indicator, temp_field)
400
401 ! Clean up
402 call neko_scratch_registry%relinquish(temp_idx)
403 call cache_data%free()
404 call cache_file%free()
405 return
406 end if
407 end if
408
409 ! ------------------------------------------------------------------------ !
410 ! Load the immersed boundary mesh
411
412 call mesh_file%init(mesh_file_name)
413 call mesh_file%read(boundary_mesh)
414
415 if (boundary_mesh%nelv .eq. 0) then
416 call neko_error('No elements in the boundary mesh')
417 end if
418
419 ! ------------------------------------------------------------------------ !
420 ! Transform the mesh if specified.
421
422 call json_get_or_default(json, 'mesh_transform.type', &
423 mesh_transform, 'none')
424
425 select case (mesh_transform)
426 case ('none')
427 ! Do nothing
428 case ('bounding_box')
429 call json_get_or_lookup(json, 'mesh_transform.box_min', box_min)
430 call json_get_or_lookup(json, 'mesh_transform.box_max', box_max)
431 call json_get_or_default(json, 'mesh_transform.keep_aspect_ratio', &
432 keep_aspect_ratio, .true.)
433
434 if (size(box_min) .ne. 3 .or. size(box_max) .ne. 3) then
435 call neko_error('Case file: mesh_transform. &
436 &box_min and box_max must be 3 element arrays of reals')
437 end if
438
439 call target_box%init(box_min, box_max)
440
441 mesh_box = get_aabb(boundary_mesh)
442
443 scaling = target_box%get_diagonal() / mesh_box%get_diagonal()
444 if (keep_aspect_ratio) then
445 scaling = minval(scaling)
446 end if
447
448 translation = - scaling * mesh_box%get_min() + target_box%get_min()
449
450 do idx_p = 1, boundary_mesh%mpts
451 boundary_mesh%points(idx_p)%x = &
452 scaling * boundary_mesh%points(idx_p)%x + translation
453 end do
454
455 ! Report the transformation applied
456 write(log_msg, '(A)') "The following transformation was applied:"
457 call neko_log%message(log_msg)
458 write(log_msg, '(A, 3F12.6)') "Scaling: ", scaling
459 call neko_log%message(log_msg)
460 write(log_msg, '(A, 3F12.6)') "Translation: ", translation
461 call neko_log%message(log_msg)
462
463 case default
464 call neko_error('Unknown mesh transform')
465 end select
466
467 ! ------------------------------------------------------------------------ !
468 ! Compute the permeability field
469
470 ! Assign the signed distance field to all GLL points in the permeability
471 ! field. Initially we just run a brute force loop over all GLL points and
472 ! compute the signed distance function. This should be replaced with a
473 ! more efficient method, such as a tree search.
474
475 call neko_scratch_registry%request_field(temp_field, temp_idx, .true.)
476
477 ! Select how to transform the distance field to a design field
478 select case (distance_transform)
479 case ('smooth_step')
480 call json_get_or_lookup(json, 'distance_transform.value', scalar_d)
481 scalar_r = real(scalar_d, kind=rp)
482
483 call signed_distance_field(temp_field, boundary_mesh, scalar_d)
484 call smooth_step_field(temp_field, scalar_r, 0.0_rp)
485
486 case ('step')
487
488 call json_get_or_lookup(json, 'distance_transform.value', scalar_d)
489 scalar_r = real(scalar_d, kind=rp)
490
491 call signed_distance_field(temp_field, boundary_mesh, scalar_d)
492 call step_function_field(temp_field, scalar_r, 1.0_rp, 0.0_rp)
493
494 case default
495 call neko_error('Unknown distance transform')
496 end select
497
498 ! Write the field to cache
499 if (cache) then
500 write(log_msg, '(A)') "Writing Brinkman source term to cache."
501 call neko_log%message(log_msg, neko_log_debug)
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)
505 end if
506
507 ! Update the global indicator field by max operator
508 call field_pwmax2(this%indicator, temp_field)
509
510 call neko_scratch_registry%relinquish(temp_idx)
511
512 end subroutine add_boundary_mesh
513
515 subroutine add_point_zone(this, json)
516 class(brinkman_source_term_t), intent(inout) :: this
517 type(json_file), intent(inout) :: json
518
519 ! Options
520 character(len=:), allocatable :: zone_name
521
522 type(field_t), pointer :: temp_field
523 class(point_zone_t), pointer :: zone
524 integer :: i, temp_idx
525
526 ! ------------------------------------------------------------------------ !
527 ! Read the options for the point zone
528
529 call json_get(json, 'name', zone_name)
530
531 ! Compute the indicator field
532 call neko_scratch_registry%request_field(temp_field, temp_idx, .true.)
533
534 zone => neko_point_zone_registry%get_point_zone(zone_name)
535
536 if (neko_bcknd_device .eq. 1) then
537 call device_cfill_mask(temp_field%x_d, 1.0_rp, temp_field%size(), &
538 zone%mask%get_d(), zone%size)
539 else
540 call cfill_mask(temp_field%x, 1.0_rp, temp_field%size(), &
541 zone%mask%get(), zone%size)
542 end if
543
544 ! Update the global indicator field by max operator
545 call field_pwmax2(this%indicator, temp_field)
546
547 call neko_scratch_registry%relinquish(temp_idx)
548 end subroutine add_point_zone
549
551 subroutine add_file(this, json)
552 class(brinkman_source_term_t), intent(inout) :: this
553 type(json_file), intent(inout) :: json
554 type(file_t) :: file
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
559
560 call json_get(json, 'file_name', file_name)
561 call json_get_or_default(json, 'field_name', field_name, &
562 'brinkman_indicator')
563 call json_get_or_default(json, 'file_index', file_idx, 0)
564
565 call neko_scratch_registry%request_field(temp_field, temp_idx, .true.)
566
567 call file%init(file_name)
568 call file%set_counter(file_idx)
569
570 call filename_suffix(file_name, suffix)
571 select case (trim(suffix))
572 case ('fld')
573 block
574 type(fld_file_data_t) :: fld_data
575 type(field_list_t) :: fld_fields
576 integer :: idx(1)
577
578 call fld_data%init()
579 call file%read(fld_data)
580 select case (field_name(1:1))
581 case ('p')
582 call fld_data%import_fields(p = temp_field)
583 case ('u')
584 call fld_data%import_fields(u = temp_field)
585 case ('v')
586 call fld_data%import_fields(v = temp_field)
587 case ('w')
588 call fld_data%import_fields(w = temp_field)
589 case ('t')
590 call fld_data%import_fields(t = temp_field)
591 case ('s')
592
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)
597 else
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')
601 end if
602
603 call fld_fields%init(1)
604 call fld_fields%assign(1, temp_field)
605
606 call fld_data%import_fields(s_target_list = fld_fields, &
607 s_index_list = idx)
608
609 call fld_fields%free()
610 case default
611 call neko_error('Unknown field prefix in field name: ' // &
612 trim(field_name))
613 end select
614
615 call fld_data%free()
616 end block
617
618 case ('vtkhdf')
619
620 ! VTKHDF will read the name of the field object.
621 tmp_str = trim(temp_field%name)
622 temp_field%name = trim(field_name)
623
624 call file%read(temp_field)
625
626 temp_field%name = trim(tmp_str)
627
628 case default
629 call neko_error("Brinkman cannot read file: " // trim(file_name))
630 end select
631
632 ! Update the global indicator field by max operator
633 call field_pwmax2(this%indicator, temp_field)
634
635 call neko_scratch_registry%relinquish(temp_idx)
636 call file%free()
637
638 end subroutine add_file
639
641 subroutine add_field(this, json)
642 class(brinkman_source_term_t), intent(inout) :: this
643 type(json_file), intent(inout) :: json
644 character(len=:), allocatable :: field_name
645 type(field_t), pointer :: temp_field
646
647 ! Read the field name
648 call json_get(json, 'name', field_name)
649 temp_field => neko_registry%get_field(field_name)
650
651 ! Update the global indicator field by max operator
652 call field_pwmax2(this%indicator, temp_field)
653
654 if (associated(temp_field)) nullify(temp_field)
655
656 end subroutine add_field
657
658end module brinkman_source_term
double real
Copy data between host and device (or device and device)
Definition device.F90:71
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 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.
Coefficients.
Definition coef.f90:34
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.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
integer, parameter, public device_to_host
Definition device.F90:47
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.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Filter to be applied to a scalar field.
Definition filter.f90:38
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.
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:89
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
A module containing mapping functions and subroutines. These functions are used to modify fields in a...
Definition mappings.f90:36
subroutine, public permeability_field(f, k_0, k_1, q)
Apply a permeability function to a field.
Definition mappings.f90:90
subroutine, public smooth_step_field(f, edge0, edge1)
Apply a smooth step function to a field.
Definition mappings.f90:70
subroutine, public step_function_field(f, x0, value0, value1)
Apply a step function to a field.
Definition mappings.f90:109
Definition math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:252
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:433
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines an output.
Definition output.f90:34
A PDE based filter.
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:107
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
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.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
Module with things related to the simulation time.
Defines a triangular surface mesh.
Definition tri_mesh.f90:35
Utilities.
Definition utils.f90:35
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Definition utils.f90:118
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:63
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...
Definition file.f90:56
Base abstract class for filter.
Definition filter.f90:47
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.
Definition space.f90:63
A struct that contains all info about the time, expand as needed.