37 use json_module,
only : json_file
74 character(len=1024) :: bf_rgstry_pref
76 logical :: baseflow_set = .false.
79 logical :: baseflow_is_ic = .false.
82 type(
field_t),
pointer :: fringe => null()
84 real(kind=
rp) :: amplitudes(3)
86 character(len=1024) :: fringe_registry_name
88 logical :: dump_fields = .false.
90 character(NEKO_FNAME_LEN) :: dump_fname
92 logical :: check = .true.
97 procedure, pass(this) :: init_constant => &
100 procedure, pass(this) :: init_user => &
103 procedure, pass(this) :: init_field => &
106 procedure, pass(this) :: init_common => &
119 type(json_file),
intent(inout) :: json
121 type(
coef_t),
intent(in),
target :: coef
122 character(len=*),
intent(in) :: variable_name
124 real(kind=
rp),
allocatable :: amplitudes(:)
125 character(len=:),
allocatable :: baseflow_method
126 character(len=:),
allocatable :: read_str, fringe_registry_name, &
127 bf_registry_pref, dump_fname
128 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
129 logical :: interpolate, dump_fields
130 real(kind=
rp),
allocatable :: constant_value(:)
131 real(kind=
rp) :: start_time, end_time, tolerance
134 type(json_file) :: baseflow_subdict
143 fringe_registry_name,
"sponge_fringe")
146 if (
size(amplitudes) .ne. 3)
then
147 call neko_error(
"(SPONGE) Expected 3 elements for 'amplitudes'")
156 call json_get(baseflow_subdict,
"method", baseflow_method)
158 bf_registry_pref,
"sponge_bf")
160 select case (trim(baseflow_method))
172 call json_get(baseflow_subdict,
'file_name', read_str)
173 fname = trim(read_str)
178 mesh_fname = trim(read_str)
181 type(json_file) :: interp_subdict
185 call this%init_field(fields, coef, start_time, end_time, amplitudes, &
186 fringe_registry_name, bf_registry_pref, dump_fields, dump_fname, &
187 fname, interpolate, mesh_fname, interp_subdict)
194 if (
size(constant_value) .lt. 3)
then
195 call neko_error(
"(SPONGE) Expected 3 elements for 'value'")
198 call this%init_constant(fields, coef, start_time, end_time, &
199 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
200 dump_fname, constant_value)
205 call this%init_user(fields, coef, start_time, end_time, amplitudes, &
206 fringe_registry_name, bf_registry_pref, dump_fields, dump_fname)
209 call neko_error(
"(SPONGE) " // trim(baseflow_method) // &
210 " is not a valid method")
219 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
220 dump_fname, constant_values)
222 type(field_list_t),
intent(in),
target :: fields
223 type(coef_t),
intent(in),
target :: coef
224 real(kind=rp),
intent(in) :: start_time, end_time
225 real(kind=rp),
intent(in) :: amplitudes(:)
226 character(len=*),
intent(in) :: fringe_registry_name, dump_fname, &
228 logical,
intent(in) :: dump_fields
229 real(kind=rp),
intent(in) :: constant_values(:)
235 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
241 call neko_log%message(
"Initializing bf fields", &
242 lvl = neko_log_debug)
244 call neko_registry%add_field(this%u%dof, &
245 trim(bf_registry_pref) //
"_u")
246 call neko_registry%add_field(this%v%dof, &
247 trim(bf_registry_pref) //
"_v")
248 call neko_registry%add_field(this%w%dof, &
249 trim(bf_registry_pref) //
"_w")
251 this%u_bf => neko_registry%get_field(trim(bf_registry_pref) //
"_u")
252 this%v_bf => neko_registry%get_field(trim(bf_registry_pref) //
"_v")
253 this%w_bf => neko_registry%get_field(trim(bf_registry_pref) //
"_w")
258 this%u_bf = constant_values(1)
259 this%v_bf = constant_values(2)
260 this%w_bf = constant_values(3)
262 this%baseflow_set = .true.
268 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
269 dump_fname, file_name, interpolate, mesh_file_name, interp_subdict)
271 type(field_list_t),
intent(in),
target :: fields
272 type(coef_t),
intent(in),
target :: coef
273 real(kind=rp),
intent(in) :: start_time, end_time
274 real(kind=rp),
intent(in) :: amplitudes(:)
275 character(len=*),
intent(in) :: fringe_registry_name, dump_fname, &
277 logical,
intent(in) :: dump_fields
278 character(len=*),
intent(in) :: file_name
279 logical,
intent(in) :: interpolate
280 character(len=*),
intent(inout) :: mesh_file_name
281 type(json_file),
intent(inout) :: interp_subdict
290 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
296 call neko_log%message(
"Initializing bf fields", &
297 lvl = neko_log_debug)
299 call neko_registry%add_field(this%u%dof, &
300 trim(bf_registry_pref) //
"_u")
301 call neko_registry%add_field(this%v%dof, &
302 trim(bf_registry_pref) //
"_v")
303 call neko_registry%add_field(this%w%dof, &
304 trim(bf_registry_pref) //
"_w")
306 this%u_bf => neko_registry%get_field(trim(bf_registry_pref) //
"_u")
307 this%v_bf => neko_registry%get_field(trim(bf_registry_pref) //
"_v")
308 this%w_bf => neko_registry%get_field(trim(bf_registry_pref) //
"_w")
313 call import_fields(file_name, interp_subdict, mesh_file_name, &
314 u = this%u_bf, v = this%v_bf, w = this%w_bf, &
315 interpolate = interpolate)
317 this%baseflow_set = .true.
323 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
326 type(field_list_t),
intent(in),
target :: fields
327 type(coef_t),
intent(in),
target :: coef
328 real(kind=rp),
intent(in) :: start_time, end_time
329 real(kind=rp),
intent(in) :: amplitudes(:)
330 character(len=*),
intent(in) :: fringe_registry_name, dump_fname, &
332 logical,
intent(in) :: dump_fields
338 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
345 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
348 type(field_list_t),
intent(in),
target :: fields
349 type(coef_t),
intent(in),
target :: coef
350 real(kind=rp),
intent(in) :: start_time, end_time
351 real(kind=rp),
intent(in) :: amplitudes(:)
352 character(len=*),
intent(in) :: fringe_registry_name, dump_fname, &
354 logical,
intent(in) :: dump_fields
359 call this%init_base(fields, coef, start_time, end_time)
361 this%amplitudes(1) = amplitudes(1)
362 this%amplitudes(2) = amplitudes(2)
363 this%amplitudes(3) = amplitudes(3)
365 this%fringe_registry_name = trim(fringe_registry_name)
366 this%bf_rgstry_pref = trim(bf_registry_pref)
367 this%dump_fields = dump_fields
368 this%dump_fname = trim(dump_fname)
370 call neko_log%message(
"Initializing sponge", lvl = neko_log_debug)
372 call neko_log%message(
"Pointing at fields u,v,w", &
373 lvl = neko_log_debug)
374 this%u => neko_registry%get_field_by_name(
"u")
375 this%v => neko_registry%get_field_by_name(
"v")
376 this%w => neko_registry%get_field_by_name(
"w")
384 call this%free_base()
394 this%fringe_registry_name =
""
395 this%bf_rgstry_pref =
""
396 this%baseflow_set = .false.
414 type(time_state_t),
intent(in) :: time
416 integer :: i, nfields
417 type(field_t),
pointer :: fu, fv, fw
419 type(json_file) :: json_subdict
420 character(len=:),
allocatable :: string_val
422 character(len=1024) :: u_name, v_name, w_name
423 type(fld_file_output_t) :: fout
425 type(field_t),
pointer :: wk
432 u_name = trim(this%bf_rgstry_pref) //
"_u"
433 v_name = trim(this%bf_rgstry_pref) //
"_v"
434 w_name = trim(this%bf_rgstry_pref) //
"_w"
437 this%baseflow_set = neko_registry%field_exists(trim(u_name)) &
438 .and. neko_registry%field_exists(trim(v_name)) .and. &
439 neko_registry%field_exists(trim(w_name))
441 if (.not. this%baseflow_set)
then
442 call neko_error(
"SPONGE: No baseflow set (searching for " // &
443 trim(this%bf_rgstry_pref) //
"_u)")
447 if (.not. neko_registry%field_exists( &
448 trim(this%fringe_registry_name)))
then
449 call neko_error(
"SPONGE: No fringe field set (" // &
450 this%fringe_registry_name //
" not found)")
456 neko_registry%get_field(trim(this%fringe_registry_name))
460 this%u_bf => neko_registry%get_field(trim(u_name))
461 this%v_bf => neko_registry%get_field(trim(v_name))
462 this%w_bf => neko_registry%get_field(trim(w_name))
470 if (this%dump_fields)
then
471 call fout%init(sp, trim(this%dump_fname), 4)
472 call fout%fields%assign_to_ptr(1, this%fringe)
473 call fout%fields%assign_to_field(2, this%u_bf)
474 call fout%fields%assign_to_field(3, this%v_bf)
475 call fout%fields%assign_to_field(4, this%w_bf)
476 call fout%sample(time%t)
486 fu => this%fields%get(1)
487 fv => this%fields%get(2)
488 fw => this%fields%get(3)
490 call neko_scratch_registry%request_field(wk, tmp_index, .false.)
492 if (neko_bcknd_device .eq. 1)
then
494 call device_sub3(wk%x_d, this%u_bf%x_d, this%u%x_d, this%u%size())
496 call device_col2(wk%x_d, this%fringe%x_d, this%fringe%size())
498 call device_add2s2(fu%x_d, wk%x_d, this%amplitudes(1), &
501 call device_sub3(wk%x_d, this%v_bf%x_d, this%v%x_d, this%v%size())
502 call device_col2(wk%x_d, this%fringe%x_d, this%fringe%size())
503 call device_add2s2(fv%x_d, wk%x_d, this%amplitudes(2), &
506 call device_sub3(wk%x_d, this%w_bf%x_d, this%w%x_d, this%w%size())
507 call device_col2(wk%x_d, this%fringe%x_d, this%fringe%size())
508 call device_add2s2(fw%x_d, wk%x_d, this%amplitudes(3), &
511 call sub3(wk%x, this%u_bf%x, this%u%x, this%u%size())
512 call col2(wk%x, this%fringe%x, this%fringe%size())
513 call add2s2(fu%x, wk%x, this%amplitudes(1), fu%dof%size())
515 call sub3(wk%x, this%v_bf%x, this%v%x, this%v%size())
516 call col2(wk%x, this%fringe%x, this%fringe%size())
517 call add2s2(fv%x, wk%x, this%amplitudes(2), fv%dof%size())
519 call sub3(wk%x, this%w_bf%x, this%w%x, this%w%size())
520 call col2(wk%x, this%fringe%x, this%fringe%size())
521 call add2s2(fw%x, wk%x, this%amplitudes(3), fw%dof%size())
524 call neko_scratch_registry%relinquish_field(tmp_index)
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.
Defines a simulation case.
integer, public pe_rank
MPI rank.
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
Module for file I/O operations.
Implements fld_file_output_t.
Importation of fields from fld files.
Utilities for retrieving parameters from the case files.
subroutine, public json_get_subdict_or_empty(json, key, output)
Extract a sub-object from a json object and returns an empty object if the key is missing.
integer, parameter, public neko_log_debug
Debug log level.
type(log_t), public neko_log
Global log stream.
subroutine, public sub3(a, b, c, n)
Vector subtraction .
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
integer, parameter neko_bcknd_device
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.
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.
Contains the simcomp_executor_t type.
type(simcomp_executor_t), target, public neko_simcomps
Global variable for the simulation component driver.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Implements the sponge_t type.
subroutine sponge_init_from_json(this, json, fields, coef, variable_name)
Constructor from json.
subroutine sponge_init_user(this, fields, coef, start_time, end_time, amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, dump_fname)
Initialize a sponge with a baseflow set by the user.
subroutine sponge_free(this)
Destructor.
subroutine sponge_compute(this, time)
Compute the sponge field. The sponge is applied according to the following definition: f_x = amp_x * ...
subroutine sponge_init_constant(this, fields, coef, start_time, end_time, amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, dump_fname, constant_values)
Initialize a sponge with a constant baseflow.
subroutine sponge_init_field(this, fields, coef, start_time, end_time, amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, dump_fname, file_name, interpolate, mesh_file_name, interp_subdict)
Initialize a sponge with a baseflow imported from a field file.
subroutine sponge_init_common(this, fields, coef, start_time, end_time, amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, dump_fname)
Common constructor.
Module with things related to the simulation time.
integer, parameter, public neko_fname_len
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...
A simple output saving a list of fields to a .fld file.
Base abstract type for source terms.
A struct that contains all info about the time, expand as needed.