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
41 use utils, only : neko_error
44 use time_state, only : time_state_t
45 use math, only : sub3, col2, add2s2
46 use logger, only : neko_log, neko_log_debug
48 use source_term, only : source_term_t
49 use case, only : case_t
52 use field_list, only : field_list_t
53 use coefs, only : coef_t
54 use utils, only : neko_fname_len
55 use file, only : file_t
57 use comm, only : pe_rank
59 implicit none
60 private
61
62 type, public, extends(source_term_t) :: sponge_source_term_t
63
65 type(field_t), pointer :: u => null()
66 type(field_t), pointer :: v => null()
67 type(field_t), pointer :: w => null()
69 type(field_t), pointer :: u_bf => null()
70 type(field_t), pointer :: v_bf => null()
71 type(field_t), pointer :: w_bf => null()
73 character(len=1024) :: bf_rgstry_pref
75 logical :: baseflow_set = .false.
78 logical :: baseflow_is_ic = .false.
81 type(field_t), pointer :: fringe => null()
83 real(kind=rp) :: amplitudes(3)
85 character(len=1024) :: fringe_registry_name
87 logical :: dump_fields = .false.
89 character(NEKO_FNAME_LEN) :: dump_fname
91 logical :: check = .true.
92 contains
94 procedure, pass(this) :: init => sponge_init_from_json
96 procedure, pass(this) :: init_constant => &
99 procedure, pass(this) :: init_user => &
102 procedure, pass(this) :: init_field => &
105 procedure, pass(this) :: init_common => &
108 procedure, pass(this) :: free => sponge_free
110 procedure, pass(this) :: compute_ => sponge_compute
111 end type sponge_source_term_t
112
113contains
114
116 subroutine sponge_init_from_json(this, json, fields, coef, variable_name)
117 class(sponge_source_term_t), intent(inout) :: this
118 type(json_file), intent(inout) :: json
119 type(field_list_t), intent(in), target :: fields
120 type(coef_t), intent(in), target :: coef
121 character(len=*), intent(in) :: variable_name
122
123 real(kind=rp), allocatable :: amplitudes(:)
124 character(len=:), allocatable :: baseflow_method
125 character(len=:), allocatable :: read_str, fringe_registry_name, &
126 bf_registry_pref, dump_fname
127 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
128 logical :: interpolate, dump_fields
129 real(kind=rp), allocatable :: constant_value(:)
130 real(kind=rp) :: start_time, end_time, tolerance
131 integer :: nzones
132
133 type(json_file) :: baseflow_subdict
134 integer :: i, izone
135
136 call neko_log%section("SPONGE SOURCE TERM", lvl = neko_log_debug)
137
138 call json_get_or_default(json, "dump_fields", dump_fields, .false.)
139 call json_get_or_default(json, "dump_file_name", dump_fname, &
140 "spng_fields")
141 call json_get_or_default(json, "fringe_registry_name", &
142 fringe_registry_name, "sponge_fringe")
143
144 call json_get_or_lookup(json, "amplitudes", amplitudes)
145 if (size(amplitudes) .ne. 3) then
146 call neko_error("(SPONGE) Expected 3 elements for 'amplitudes'")
147 end if
148
149 call json_get_or_default(json, "start_time", start_time, 0.0_rp)
150 call json_get_or_default(json, "end_time", end_time, huge(0.0_rp))
151
152
153 call json_get(json, 'baseflow', &
154 baseflow_subdict)
155 call json_get(baseflow_subdict, "method", baseflow_method)
156 call json_get_or_default(json, "baseflow_registry_prefix", &
157 bf_registry_pref, "sponge_bf")
158
159 select case (trim(baseflow_method))
160
161 ! Import the base flow from an fld file. The same parameters as field
162 ! initial condition apply since we use the same subroutine to read and
163 ! potentially interpolate from a field file. It would be nice to have
164 ! some kind of `gfldr` function to make this nicer.
165 case ("field")
166
167 ! The lines below are just copy pasted from set_flow_ic_int, because we
168 ! want to make sure to separate the JSON constructor
169 ! sponge_init_from_json and the constructor from attributes
170 ! sponge_init_field.
171 call json_get(baseflow_subdict, 'file_name', read_str)
172 fname = trim(read_str)
173 call json_get_or_default(baseflow_subdict, 'interpolate', interpolate, &
174 .false.)
175 call json_get_or_default(baseflow_subdict, 'tolerance', tolerance, &
176 0.000001_rp)
177 call json_get_or_default(baseflow_subdict, 'mesh_file_name', read_str, &
178 "none")
179 mesh_fname = trim(read_str)
180
181 call this%init_field(fields, coef, start_time, end_time, amplitudes, &
182 fringe_registry_name, bf_registry_pref, dump_fields, dump_fname, &
183 fname, interpolate, tolerance, mesh_fname)
184
185 ! Constant base flow
186 case ("constant")
187
188 call json_get(baseflow_subdict, "value", constant_value)
189 if (size(constant_value) .lt. 3) then
190 call neko_error("(SPONGE) Expected 3 elements for 'value'")
191 end if
192
193 call this%init_constant(fields, coef, start_time, end_time, &
194 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
195 dump_fname, constant_value)
196
197 ! Let the user set the base flow.
198 case ("user")
199
200 call this%init_user(fields, coef, start_time, end_time, amplitudes, &
201 fringe_registry_name, bf_registry_pref, dump_fields, dump_fname)
202
203 case default
204 call neko_error("(SPONGE) " // trim(baseflow_method) // &
205 " is not a valid method")
206 end select
207
208 call neko_log%end_section(lvl = neko_log_debug)
209
210 end subroutine sponge_init_from_json
211
213 subroutine sponge_init_constant(this, fields, coef, start_time, end_time, &
214 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
215 dump_fname, constant_values)
216 class(sponge_source_term_t), intent(inout) :: this
217 type(field_list_t), intent(in), target :: fields
218 type(coef_t), intent(in), target :: coef
219 real(kind=rp), intent(in) :: start_time, end_time
220 real(kind=rp), intent(in) :: amplitudes(:)
221 character(len=*), intent(in) :: fringe_registry_name, dump_fname, &
222 bf_registry_pref
223 logical, intent(in) :: dump_fields
224 real(kind=rp), intent(in) :: constant_values(:)
225
226 !
227 ! Common constructor
228 !
229 call sponge_init_common(this, fields, coef, start_time, end_time, &
230 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
231 dump_fname)
232
233 !
234 ! Create the base flow fields in the registry
235 !
236 call neko_log%message("Initializing bf fields", &
237 lvl = neko_log_debug)
238
239 call neko_registry%add_field(this%u%dof, &
240 trim(bf_registry_pref) // "_u")
241 call neko_registry%add_field(this%v%dof, &
242 trim(bf_registry_pref) // "_v")
243 call neko_registry%add_field(this%w%dof, &
244 trim(bf_registry_pref) // "_w")
245
246 this%u_bf => neko_registry%get_field(trim(bf_registry_pref) // "_u")
247 this%v_bf => neko_registry%get_field(trim(bf_registry_pref) // "_v")
248 this%w_bf => neko_registry%get_field(trim(bf_registry_pref) // "_w")
249
250 !
251 ! Assign constant values
252 !
253 this%u_bf = constant_values(1)
254 this%v_bf = constant_values(2)
255 this%w_bf = constant_values(3)
256
257 this%baseflow_set = .true.
258
259 end subroutine sponge_init_constant
260
262 subroutine sponge_init_field(this, fields, coef, start_time, end_time, &
263 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
264 dump_fname, file_name, interpolate, tolerance, mesh_file_name)
265 class(sponge_source_term_t), intent(inout) :: this
266 type(field_list_t), intent(in), target :: fields
267 type(coef_t), intent(in), target :: coef
268 real(kind=rp), intent(in) :: start_time, end_time
269 real(kind=rp), intent(in) :: amplitudes(:)
270 character(len=*), intent(in) :: fringe_registry_name, dump_fname, &
271 bf_registry_pref
272 logical, intent(in) :: dump_fields
273 character(len=*), intent(in) :: file_name
274 logical, intent(in) :: interpolate
275 real(kind=rp), intent(in) :: tolerance
276 character(len=*), intent(inout) :: mesh_file_name
277
278 type(field_t) :: wk
279 integer :: tmp_index
280
281 !
282 ! Common constructor
283 !
284 call sponge_init_common(this, fields, coef, start_time, end_time, &
285 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
286 dump_fname)
287
288 !
289 ! Create the base flow fields in the registry
290 !
291 call neko_log%message("Initializing bf fields", &
292 lvl = neko_log_debug)
293
294 call neko_registry%add_field(this%u%dof, &
295 trim(bf_registry_pref) // "_u")
296 call neko_registry%add_field(this%v%dof, &
297 trim(bf_registry_pref) // "_v")
298 call neko_registry%add_field(this%w%dof, &
299 trim(bf_registry_pref) // "_w")
300
301 this%u_bf => neko_registry%get_field(trim(bf_registry_pref) // "_u")
302 this%v_bf => neko_registry%get_field(trim(bf_registry_pref) // "_v")
303 this%w_bf => neko_registry%get_field(trim(bf_registry_pref) // "_w")
304
305 !
306 ! Import the u,v,w baseflows from fld
307 !
308 call import_fields(file_name, mesh_file_name, &
309 u = this%u_bf, v = this%v_bf, w = this%w_bf, &
310 interpolate = interpolate, tolerance = tolerance)
311
312 this%baseflow_set = .true.
313
314 end subroutine sponge_init_field
315
317 subroutine sponge_init_user(this, fields, coef, start_time, end_time, &
318 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
319 dump_fname)
320 class(sponge_source_term_t), intent(inout) :: this
321 type(field_list_t), intent(in), target :: fields
322 type(coef_t), intent(in), target :: coef
323 real(kind=rp), intent(in) :: start_time, end_time
324 real(kind=rp), intent(in) :: amplitudes(:)
325 character(len=*), intent(in) :: fringe_registry_name, dump_fname, &
326 bf_registry_pref
327 logical, intent(in) :: dump_fields
328
329 !
330 ! Common constructor
331 !
332 call sponge_init_common(this, fields, coef, start_time, end_time, &
333 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
334 dump_fname)
335
336 end subroutine sponge_init_user
337
339 subroutine sponge_init_common(this, fields, coef, start_time, end_time, &
340 amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, &
341 dump_fname)
342 class(sponge_source_term_t), intent(inout) :: this
343 type(field_list_t), intent(in), target :: fields
344 type(coef_t), intent(in), target :: coef
345 real(kind=rp), intent(in) :: start_time, end_time
346 real(kind=rp), intent(in) :: amplitudes(:)
347 character(len=*), intent(in) :: fringe_registry_name, dump_fname, &
348 bf_registry_pref
349 logical, intent(in) :: dump_fields
350
351 integer :: i
352
353 call this%free()
354 call this%init_base(fields, coef, start_time, end_time)
355
356 this%amplitudes(1) = amplitudes(1)
357 this%amplitudes(2) = amplitudes(2)
358 this%amplitudes(3) = amplitudes(3)
359
360 this%fringe_registry_name = trim(fringe_registry_name)
361 this%bf_rgstry_pref = trim(bf_registry_pref)
362 this%dump_fields = dump_fields
363 this%dump_fname = trim(dump_fname)
364
365 call neko_log%message("Initializing sponge", lvl = neko_log_debug)
366
367 call neko_log%message("Pointing at fields u,v,w", &
368 lvl = neko_log_debug)
369 this%u => neko_registry%get_field_by_name("u")
370 this%v => neko_registry%get_field_by_name("v")
371 this%w => neko_registry%get_field_by_name("w")
372
373 end subroutine sponge_init_common
374
376 subroutine sponge_free(this)
377 class(sponge_source_term_t), intent(inout) :: this
378
379 call this%free_base()
380
381 nullify(this%u)
382 nullify(this%v)
383 nullify(this%w)
384 nullify(this%u_bf)
385 nullify(this%v_bf)
386 nullify(this%w_bf)
387 nullify(this%fringe)
388
389 this%fringe_registry_name = ""
390 this%bf_rgstry_pref = ""
391 this%baseflow_set = .false.
392
393 end subroutine sponge_free
394
407 subroutine sponge_compute(this, time)
408 class(sponge_source_term_t), intent(inout) :: this
409 type(time_state_t), intent(in) :: time
410
411 integer :: i, nfields
412 type(field_t), pointer :: fu, fv, fw
413
414 type(json_file) :: json_subdict
415 character(len=:), allocatable :: string_val
416
417 character(len=1024) :: u_name, v_name, w_name
418 type(fld_file_output_t) :: fout
419 integer :: tmp_index
420 type(field_t), pointer :: wk
421
422 !
423 ! Do some checks at the first timestep
424 !
425 if (this%check) then
426
427 u_name = trim(this%bf_rgstry_pref) // "_u"
428 v_name = trim(this%bf_rgstry_pref) // "_v"
429 w_name = trim(this%bf_rgstry_pref) // "_w"
430
431 ! Check if all the base flow fields exist in the registry
432 this%baseflow_set = neko_registry%field_exists(trim(u_name)) &
433 .and. neko_registry%field_exists(trim(v_name)) .and. &
434 neko_registry%field_exists(trim(w_name))
435
436 if (.not. this%baseflow_set) then
437 call neko_error("SPONGE: No baseflow set (searching for " // &
438 trim(this%bf_rgstry_pref) // "_u)")
439 end if
440
441 ! Check if the user has added the fringe field in the registry
442 if (.not. neko_registry%field_exists( &
443 trim(this%fringe_registry_name))) then
444 call neko_error("SPONGE: No fringe field set (" // &
445 this%fringe_registry_name // " not found)")
446 end if
447
448 ! This will throw an error if the user hasn't added 'sponge_fringe'
449 ! to the registry.
450 this%fringe => &
451 neko_registry%get_field(trim(this%fringe_registry_name))
452
453 ! This will throw an error if the user hasn't added the base flow fields
454 ! to the registry.
455 this%u_bf => neko_registry%get_field(trim(u_name))
456 this%v_bf => neko_registry%get_field(trim(v_name))
457 this%w_bf => neko_registry%get_field(trim(w_name))
458
459 ! Reset the flag
460 this%check = .false.
461
462 !
463 ! Dump the fringe and/or baseflow fields for visualization
464 !
465 if (this%dump_fields) then
466 call fout%init(sp, trim(this%dump_fname), 4)
467 call fout%fields%assign_to_ptr(1, this%fringe)
468 call fout%fields%assign_to_field(2, this%u_bf)
469 call fout%fields%assign_to_field(3, this%v_bf)
470 call fout%fields%assign_to_field(4, this%w_bf)
471 call fout%sample(time%t)
472 end if
473
474 end if
475
476
477 !
478 ! Start computation of source term
479 !
480
481 fu => this%fields%get(1)
482 fv => this%fields%get(2)
483 fw => this%fields%get(3)
484
485 call neko_scratch_registry%request_field(wk, tmp_index, .false.)
486
487 if (neko_bcknd_device .eq. 1) then
488 ! wk = u_bf - u
489 call device_sub3(wk%x_d, this%u_bf%x_d, this%u%x_d, this%u%size())
490 ! wk = fringe * wk = fringe * (u_bf - u)
491 call device_col2(wk%x_d, this%fringe%x_d, this%fringe%size())
492 ! fu = fu + amplitude(1)*wk = fu + amplitude(1)*fringe*(u_bf - u)
493 call device_add2s2(fu%x_d, wk%x_d, this%amplitudes(1), &
494 fu%dof%size())
495
496 call device_sub3(wk%x_d, this%v_bf%x_d, this%v%x_d, this%v%size())
497 call device_col2(wk%x_d, this%fringe%x_d, this%fringe%size())
498 call device_add2s2(fv%x_d, wk%x_d, this%amplitudes(2), &
499 fv%dof%size())
500
501 call device_sub3(wk%x_d, this%w_bf%x_d, this%w%x_d, this%w%size())
502 call device_col2(wk%x_d, this%fringe%x_d, this%fringe%size())
503 call device_add2s2(fw%x_d, wk%x_d, this%amplitudes(3), &
504 fw%dof%size())
505 else
506 call sub3(wk%x, this%u_bf%x, this%u%x, this%u%size())
507 call col2(wk%x, this%fringe%x, this%fringe%size())
508 call add2s2(fu%x, wk%x, this%amplitudes(1), fu%dof%size())
509
510 call sub3(wk%x, this%v_bf%x, this%v%x, this%v%size())
511 call col2(wk%x, this%fringe%x, this%fringe%size())
512 call add2s2(fv%x, wk%x, this%amplitudes(2), fv%dof%size())
513
514 call sub3(wk%x, this%w_bf%x, this%w%x, this%w%size())
515 call col2(wk%x, this%fringe%x, this%fringe%size())
516 call add2s2(fw%x, wk%x, this%amplitudes(3), fw%dof%size())
517 end if
518
519 call neko_scratch_registry%relinquish_field(tmp_index)
520
521 end subroutine sponge_compute
522
523end 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.
subroutine, public import_fields(fname, mesh_fname, u, v, w, p, t, s_target_list, s_index_list, interpolate, tolerance)
Imports fields from an fld file, potentially with interpolation.
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
Definition math.f90:60
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:781
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:812
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_field(this, fields, coef, start_time, end_time, amplitudes, fringe_registry_name, bf_registry_pref, dump_fields, dump_fname, file_name, interpolate, tolerance, mesh_file_name)
Initialize a sponge with a baseflow imported from a field file.
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_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:58
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:55
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.