Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
sponge_source_term.f90
Go to the documentation of this file.
1! Copyright (c) 2025, 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!
33!
36 use num_types, only : rp, dp, sp
37 use json_module, only : json_file
38 use registry, only : neko_registry
39 use field, only : field_t
42 use utils, only : neko_error
45 use time_state, only : time_state_t
46 use math, only : sub3, col2, add2s2
47 use logger, only : neko_log, neko_log_debug
49 use source_term, only : source_term_t
50 use case, only : case_t
53 use field_list, only : field_list_t
54 use coefs, only : coef_t
55 use utils, only : neko_fname_len
56 use file, only : file_t
58 use comm, only : pe_rank
60 implicit none
61 private
62
63 type, public, extends(source_term_t) :: sponge_source_term_t
64
66 type(field_t), pointer :: u => null()
67 type(field_t), pointer :: v => null()
68 type(field_t), pointer :: w => null()
70 type(field_t), pointer :: u_bf => null()
71 type(field_t), pointer :: v_bf => null()
72 type(field_t), pointer :: w_bf => null()
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.
93 contains
95 procedure, pass(this) :: init => sponge_init_from_json
97 procedure, pass(this) :: init_constant => &
100 procedure, pass(this) :: init_user => &
103 procedure, pass(this) :: init_field => &
106 procedure, pass(this) :: init_common => &
109 procedure, pass(this) :: free => sponge_free
111 procedure, pass(this) :: compute_ => sponge_compute
112 end type sponge_source_term_t
113
114contains
115
117 subroutine sponge_init_from_json(this, json, fields, coef, variable_name)
118 class(sponge_source_term_t), intent(inout) :: this
119 type(json_file), intent(inout) :: json
120 type(field_list_t), intent(in), target :: fields
121 type(coef_t), intent(in), target :: coef
122 character(len=*), intent(in) :: variable_name
123
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
132 integer :: nzones
133
134 type(json_file) :: baseflow_subdict
135 integer :: i, izone
136
137 call neko_log%section("SPONGE SOURCE TERM", lvl = neko_log_debug)
138
139 call json_get_or_default(json, "dump_fields", dump_fields, .false.)
140 call json_get_or_default(json, "dump_file_name", dump_fname, &
141 "spng_fields")
142 call json_get_or_default(json, "fringe_registry_name", &
143 fringe_registry_name, "sponge_fringe")
144
145 call json_get_or_lookup(json, "amplitudes", amplitudes)
146 if (size(amplitudes) .ne. 3) then
147 call neko_error("(SPONGE) Expected 3 elements for 'amplitudes'")
148 end if
149
150 call json_get_or_default(json, "start_time", start_time, 0.0_rp)
151 call json_get_or_default(json, "end_time", end_time, huge(0.0_rp))
152
153
154 call json_get(json, 'baseflow', &
155 baseflow_subdict)
156 call json_get(baseflow_subdict, "method", baseflow_method)
157 call json_get_or_default(json, "baseflow_registry_prefix", &
158 bf_registry_pref, "sponge_bf")
159
160 select case (trim(baseflow_method))
161
162 ! Import the base flow from an fld file. The same parameters as field
163 ! initial condition apply since we use the same subroutine to read and
164 ! potentially interpolate from a field file. It would be nice to have
165 ! some kind of `gfldr` function to make this nicer.
166 case ("field")
167
168 ! The lines below are just copy pasted from set_flow_ic_int, because we
169 ! want to make sure to separate the JSON constructor
170 ! sponge_init_from_json and the constructor from attributes
171 ! sponge_init_field.
172 call json_get(baseflow_subdict, 'file_name', read_str)
173 fname = trim(read_str)
174 call json_get_or_default(baseflow_subdict, 'interpolate', interpolate, &
175 .false.)
176 call json_get_or_default(baseflow_subdict, 'mesh_file_name', read_str, &
177 "none")
178 mesh_fname = trim(read_str)
179
180 block
181 type(json_file) :: interp_subdict
182 call json_get_subdict_or_empty(baseflow_subdict, "interpolation", &
183 interp_subdict)
184
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)
188 end block
189
190 ! Constant base flow
191 case ("constant")
192
193 call json_get_or_lookup(baseflow_subdict, "value", constant_value)
194 if (size(constant_value) .lt. 3) then
195 call neko_error("(SPONGE) Expected 3 elements for 'value'")
196 end if
197
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)
201
202 ! Let the user set the base flow.
203 case ("user")
204
205 call this%init_user(fields, coef, start_time, end_time, amplitudes, &
206 fringe_registry_name, bf_registry_pref, dump_fields, dump_fname)
207
208 case default
209 call neko_error("(SPONGE) " // trim(baseflow_method) // &
210 " is not a valid method")
211 end select
212
213 call neko_log%end_section(lvl = neko_log_debug)
214
215 end subroutine sponge_init_from_json
216
218 subroutine sponge_init_constant(this, fields, coef, start_time, end_time, &
219 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
220 dump_fname, constant_values)
221 class(sponge_source_term_t), intent(inout) :: this
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, &
227 bf_registry_pref
228 logical, intent(in) :: dump_fields
229 real(kind=rp), intent(in) :: constant_values(:)
230
231 !
232 ! Common constructor
233 !
234 call sponge_init_common(this, fields, coef, start_time, end_time, &
235 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
236 dump_fname)
237
238 !
239 ! Create the base flow fields in the registry
240 !
241 call neko_log%message("Initializing bf fields", &
242 lvl = neko_log_debug)
243
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")
250
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")
254
255 !
256 ! Assign constant values
257 !
258 this%u_bf = constant_values(1)
259 this%v_bf = constant_values(2)
260 this%w_bf = constant_values(3)
261
262 this%baseflow_set = .true.
263
264 end subroutine sponge_init_constant
265
267 subroutine sponge_init_field(this, fields, coef, start_time, end_time, &
268 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
269 dump_fname, file_name, interpolate, mesh_file_name, interp_subdict)
270 class(sponge_source_term_t), intent(inout) :: this
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, &
276 bf_registry_pref
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
282
283 type(field_t) :: wk
284 integer :: tmp_index
285
286 !
287 ! Common constructor
288 !
289 call sponge_init_common(this, fields, coef, start_time, end_time, &
290 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
291 dump_fname)
292
293 !
294 ! Create the base flow fields in the registry
295 !
296 call neko_log%message("Initializing bf fields", &
297 lvl = neko_log_debug)
298
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")
305
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")
309
310 !
311 ! Import the u,v,w baseflows from fld
312 !
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)
316
317 this%baseflow_set = .true.
318
319 end subroutine sponge_init_field
320
322 subroutine sponge_init_user(this, fields, coef, start_time, end_time, &
323 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
324 dump_fname)
325 class(sponge_source_term_t), intent(inout) :: this
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, &
331 bf_registry_pref
332 logical, intent(in) :: dump_fields
333
334 !
335 ! Common constructor
336 !
337 call sponge_init_common(this, fields, coef, start_time, end_time, &
338 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
339 dump_fname)
340
341 end subroutine sponge_init_user
342
344 subroutine sponge_init_common(this, fields, coef, start_time, end_time, &
345 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
346 dump_fname)
347 class(sponge_source_term_t), intent(inout) :: this
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, &
353 bf_registry_pref
354 logical, intent(in) :: dump_fields
355
356 integer :: i
357
358 call this%free()
359 call this%init_base(fields, coef, start_time, end_time)
360
361 this%amplitudes(1) = amplitudes(1)
362 this%amplitudes(2) = amplitudes(2)
363 this%amplitudes(3) = amplitudes(3)
364
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)
369
370 call neko_log%message("Initializing sponge", lvl = neko_log_debug)
371
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")
377
378 end subroutine sponge_init_common
379
381 subroutine sponge_free(this)
382 class(sponge_source_term_t), intent(inout) :: this
383
384 call this%free_base()
385
386 nullify(this%u)
387 nullify(this%v)
388 nullify(this%w)
389 nullify(this%u_bf)
390 nullify(this%v_bf)
391 nullify(this%w_bf)
392 nullify(this%fringe)
393
394 this%fringe_registry_name = ""
395 this%bf_rgstry_pref = ""
396 this%baseflow_set = .false.
397
398 end subroutine sponge_free
399
412 subroutine sponge_compute(this, time)
413 class(sponge_source_term_t), intent(inout) :: this
414 type(time_state_t), intent(in) :: time
415
416 integer :: i, nfields
417 type(field_t), pointer :: fu, fv, fw
418
419 type(json_file) :: json_subdict
420 character(len=:), allocatable :: string_val
421
422 character(len=1024) :: u_name, v_name, w_name
423 type(fld_file_output_t) :: fout
424 integer :: tmp_index
425 type(field_t), pointer :: wk
426
427 !
428 ! Do some checks at the first timestep
429 !
430 if (this%check) then
431
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"
435
436 ! Check if all the base flow fields exist in the registry
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))
440
441 if (.not. this%baseflow_set) then
442 call neko_error("SPONGE: No baseflow set (searching for " // &
443 trim(this%bf_rgstry_pref) // "_u)")
444 end if
445
446 ! Check if the user has added the fringe field in the registry
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)")
451 end if
452
453 ! This will throw an error if the user hasn't added 'sponge_fringe'
454 ! to the registry.
455 this%fringe => &
456 neko_registry%get_field(trim(this%fringe_registry_name))
457
458 ! This will throw an error if the user hasn't added the base flow fields
459 ! to the registry.
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))
463
464 ! Reset the flag
465 this%check = .false.
466
467 !
468 ! Dump the fringe and/or baseflow fields for visualization
469 !
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)
477 end if
478
479 end if
480
481
482 !
483 ! Start computation of source term
484 !
485
486 fu => this%fields%get(1)
487 fv => this%fields%get(2)
488 fw => this%fields%get(3)
489
490 call neko_scratch_registry%request_field(wk, tmp_index, .false.)
491
492 if (neko_bcknd_device .eq. 1) then
493 ! wk = u_bf - u
494 call device_sub3(wk%x_d, this%u_bf%x_d, this%u%x_d, this%u%size())
495 ! wk = fringe * wk = fringe * (u_bf - u)
496 call device_col2(wk%x_d, this%fringe%x_d, this%fringe%size())
497 ! fu = fu + amplitude(1)*wk = fu + amplitude(1)*fringe*(u_bf - u)
498 call device_add2s2(fu%x_d, wk%x_d, this%amplitudes(1), &
499 fu%dof%size())
500
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), &
504 fv%dof%size())
505
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), &
509 fw%dof%size())
510 else
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())
514
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())
518
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())
522 end if
523
524 call neko_scratch_registry%relinquish_field(tmp_index)
525
526 end subroutine sponge_compute
527
528end module sponge_source_term
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.
Defines a simulation case.
Definition case.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
integer, public pe_rank
MPI rank.
Definition comm.F90:56
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.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
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.
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
Definition math.f90:60
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:782
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:855
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:813
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 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.
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.
Utilities.
Definition utils.f90:35
integer, parameter, public neko_fname_len
Definition utils.f90:40
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:62
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:56
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.