Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
fluid_scheme_incompressible.f90
Go to the documentation of this file.
1! Copyright (c) 2020-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!
36 use gather_scatter, only : gs_t, gs_op_min, gs_op_max
39 use checkpoint, only : chkp_t
40 use mean_flow, only : mean_flow_t
41 use num_types, only : rp, i8
42 use comm
44 use field, only : field_t
45 use space, only : space_t, gll
46 use dofmap, only : dofmap_t
48 use krylov, only : ksp_t, krylov_solver_factory, krylov_solver_destroy, &
50 use coefs, only: coef_t
52 use dirichlet, only : dirichlet_t
55 use jacobi, only : jacobi_t
56 use sx_jacobi, only : sx_jacobi_t
58 use hsmg, only : hsmg_t
59 use phmg, only : phmg_t
60 use precon, only : pc_t, precon_factory, precon_destroy
61 use fluid_stats, only : fluid_stats_t
62 use bc, only : bc_t
63 use bc_list, only : bc_list_t
65 use math, only : cfill, add2s2, glsum
67 use operators, only : cfl
72 use json_module, only : json_file, json_core, json_value
76 use utils, only : neko_error, neko_warning
84 implicit none
85 private
86
91 class(ksp_t), allocatable :: ksp_vel
92 class(ksp_t), allocatable :: ksp_prs
93 class(pc_t), allocatable :: pc_vel
94 class(pc_t), allocatable :: pc_prs
95 integer :: vel_projection_dim
96 integer :: pr_projection_dim
97 integer :: vel_projection_activ_step
98 integer :: pr_projection_activ_step
99 logical :: strict_convergence
101 logical :: if_gradient_jump_penalty
102 type(gradient_jump_penalty_t) :: gradient_jump_penalty_u
103 type(gradient_jump_penalty_t) :: gradient_jump_penalty_v
104 type(gradient_jump_penalty_t) :: gradient_jump_penalty_w
106 type(field_t), pointer :: u_e => null()
107 type(field_t), pointer :: v_e => null()
108 type(field_t), pointer :: w_e => null()
109
110 type(mean_flow_t) :: mean
112 type(mean_sqr_flow_t) :: mean_sqr
113 logical :: forced_flow_rate = .false.
114
116 character(len=:), allocatable :: nut_field_name
117
119 integer(kind=i8) :: glb_n_points
121 integer(kind=i8) :: glb_unique_points
122 type(scratch_registry_t) :: scratch
123 contains
125 procedure, pass(this) :: init_base => fluid_scheme_init_base
126 procedure, pass(this) :: scheme_free => fluid_scheme_free
128 procedure, pass(this) :: validate => fluid_scheme_validate
130 procedure, pass(this) :: bc_apply_vel => fluid_scheme_bc_apply_vel
132 procedure, pass(this) :: bc_apply_prs => fluid_scheme_bc_apply_prs
134 procedure, pass(this) :: compute_cfl => fluid_compute_cfl
136 procedure, pass(this) :: set_material_properties => &
138
140 procedure, pass(this) :: update_material_properties => &
143 procedure, nopass :: solver_factory => fluid_scheme_solver_factory
145 procedure, pass(this) :: precon_factory_ => fluid_scheme_precon_factory
147
148 interface
149
150 module subroutine fluid_scheme_factory(object, type_name)
151 class(fluid_scheme_base_t), intent(inout), allocatable :: object
152 character(len=*) :: type_name
153 end subroutine fluid_scheme_factory
154 end interface
155
156 public :: fluid_scheme_incompressible_t, fluid_scheme_factory
157
158contains
159
161 subroutine fluid_scheme_init_base(this, msh, lx, params, scheme, user, &
162 kspv_init)
163 implicit none
164 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
165 type(mesh_t), target, intent(inout) :: msh
166 integer, intent(in) :: lx
167 character(len=*), intent(in) :: scheme
168 type(json_file), target, intent(inout) :: params
169 type(user_t), target, intent(in) :: user
170 logical, intent(in) :: kspv_init
171 type(dirichlet_t) :: bdry_mask
172 character(len=LOG_SIZE) :: log_buf
173 real(kind=rp), allocatable :: real_vec(:)
174 real(kind=rp) :: real_val, kappa, b, z0
175 logical :: logical_val
176 integer :: integer_val, ierr
177 type(json_file) :: wm_json
178 character(len=:), allocatable :: string_val1, string_val2
179 real(kind=rp) :: gjp_param_a, gjp_param_b
180
181 !
182 ! SEM simulation fundamentals
183 !
184
185 this%msh => msh
186
187 if (msh%gdim .eq. 2) then
188 call this%Xh%init(gll, lx, lx)
189 else
190 call this%Xh%init(gll, lx, lx, lx)
191 end if
192
193 call this%dm_Xh%init(msh, this%Xh)
194
195 call this%gs_Xh%init(this%dm_Xh)
196
197 call this%c_Xh%init(this%gs_Xh)
198
199 ! Local scratch registry
200 this%scratch = scratch_registry_t(this%dm_Xh, 10, 2)
201
202 !
203 ! First section of fluid log
204 !
205
206 call neko_log%section('Fluid')
207 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
208 call neko_log%message(log_buf)
209
210 !
211 ! Material properties
212 !
213 call this%set_material_properties(params, user)
214
215 !
216 ! Turbulence modelling and variable material properties
217 !
218 if (params%valid_path('case.fluid.nut_field')) then
219 call json_get(params, 'case.fluid.nut_field', this%nut_field_name)
220 this%variable_material_properties = .true.
221 else
222 this%nut_field_name = ""
223 end if
224
225 ! Fill mu and rho field with the physical value
226
227 call this%mu_field%init(this%dm_Xh, "mu")
228 call this%rho_field%init(this%dm_Xh, "rho")
229 call field_cfill(this%mu_field, this%mu, this%mu_field%size())
230 call field_cfill(this%rho_field, this%rho, this%mu_field%size())
231
232 ! Since mu, rho is a field, and the none-stress simulation fetches
233 ! data from the host arrays, we need to mirror the constant
234 ! material properties on the host
235 if (neko_bcknd_device .eq. 1) then
236 call cfill(this%mu_field%x, this%mu, this%mu_field%size())
237 call cfill(this%rho_field%x, this%rho, this%rho_field%size())
238 end if
239
240
241 ! Projection spaces
242 call json_get_or_default(params, &
243 'case.fluid.velocity_solver.projection_space_size', &
244 this%vel_projection_dim, 20)
245 call json_get_or_default(params, &
246 'case.fluid.pressure_solver.projection_space_size', &
247 this%pr_projection_dim, 20)
248 call json_get_or_default(params, &
249 'case.fluid.velocity_solver.projection_hold_steps', &
250 this%vel_projection_activ_step, 5)
251 call json_get_or_default(params, &
252 'case.fluid.pressure_solver.projection_hold_steps', &
253 this%pr_projection_activ_step, 5)
254
255
256 call json_get_or_default(params, 'case.fluid.freeze', this%freeze, .false.)
257
258 if (params%valid_path("case.fluid.flow_rate_force")) then
259 this%forced_flow_rate = .true.
260 end if
261
262
263 if (lx .lt. 10) then
264 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
265 else if (lx .ge. 10) then
266 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
267 else
268 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
269 end if
270 call neko_log%message(log_buf)
271 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
272 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
273
274 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
275 call neko_log%message(log_buf)
276 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
277 call neko_log%message(log_buf)
278
279 write(log_buf, '(A,ES13.6)') 'rho :', this%rho
280 call neko_log%message(log_buf)
281 write(log_buf, '(A,ES13.6)') 'mu :', this%mu
282 call neko_log%message(log_buf)
283
284 call json_get(params, 'case.numerics.dealias', logical_val)
285 write(log_buf, '(A, L1)') 'Dealias : ', logical_val
286 call neko_log%message(log_buf)
287
288 write(log_buf, '(A, L1)') 'LES : ', this%variable_material_properties
289 call neko_log%message(log_buf)
290
291 call json_get_or_default(params, 'case.output_boundary', logical_val, &
292 .false.)
293 write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
294 call neko_log%message(log_buf)
295
296 !
297 ! Setup right-hand side fields.
298 !
299 allocate(this%f_x)
300 allocate(this%f_y)
301 allocate(this%f_z)
302 call this%f_x%init(this%dm_Xh, fld_name = "fluid_rhs_x")
303 call this%f_y%init(this%dm_Xh, fld_name = "fluid_rhs_y")
304 call this%f_z%init(this%dm_Xh, fld_name = "fluid_rhs_z")
305
306 ! Initialize the source term
307 call this%source_term%init(this%f_x, this%f_y, this%f_z, this%c_Xh, user)
308 call this%source_term%add(params, 'case.fluid.source_terms')
309
310 ! Initialize velocity solver
311 if (kspv_init) then
312 call neko_log%section("Velocity solver")
313 call json_get_or_default(params, &
314 'case.fluid.velocity_solver.max_iterations', &
315 integer_val, ksp_max_iter)
316 call json_get(params, 'case.fluid.velocity_solver.type', string_val1)
317 call json_get(params, 'case.fluid.velocity_solver.preconditioner', &
318 string_val2)
319 call json_get(params, 'case.fluid.velocity_solver.absolute_tolerance', &
320 real_val)
321 call json_get_or_default(params, &
322 'case.fluid.velocity_solver.monitor', &
323 logical_val, .false.)
324
325 call neko_log%message('Type : ('// trim(string_val1) // &
326 ', ' // trim(string_val2) // ')')
327
328 write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
329 call neko_log%message(log_buf)
330 call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
331 string_val1, integer_val, real_val, logical_val)
332 call this%precon_factory_(this%pc_vel, this%ksp_vel, &
333 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, string_val2)
334 call neko_log%end_section()
335 end if
336
337 ! Strict convergence for the velocity solver
338 call json_get_or_default(params, 'case.fluid.strict_convergence', &
339 this%strict_convergence, .false.)
340
341 ! Assign velocity fields
342 call neko_field_registry%add_field(this%dm_Xh, 'u')
343 call neko_field_registry%add_field(this%dm_Xh, 'v')
344 call neko_field_registry%add_field(this%dm_Xh, 'w')
345 this%u => neko_field_registry%get_field('u')
346 this%v => neko_field_registry%get_field('v')
347 this%w => neko_field_registry%get_field('w')
348
349 !! Initialize time-lag fields
350 call this%ulag%init(this%u, 2)
351 call this%vlag%init(this%v, 2)
352 call this%wlag%init(this%w, 2)
353
354 ! Initiate gradient jump penalty
355 call json_get_or_default(params, &
356 'case.fluid.gradient_jump_penalty.enabled',&
357 this%if_gradient_jump_penalty, .false.)
358
359 if (this%if_gradient_jump_penalty .eqv. .true.) then
360 if ((this%dm_Xh%xh%lx - 1) .eq. 1) then
361 call json_get_or_default(params, &
362 'case.fluid.gradient_jump_penalty.tau',&
363 gjp_param_a, 0.02_rp)
364 gjp_param_b = 0.0_rp
365 else
366 call json_get_or_default(params, &
367 'case.fluid.gradient_jump_penalty.scaling_factor',&
368 gjp_param_a, 0.8_rp)
369 call json_get_or_default(params, &
370 'case.fluid.gradient_jump_penalty.scaling_exponent',&
371 gjp_param_b, 4.0_rp)
372 end if
373 call this%gradient_jump_penalty_u%init(params, this%dm_Xh, this%c_Xh, &
374 gjp_param_a, gjp_param_b)
375 call this%gradient_jump_penalty_v%init(params, this%dm_Xh, this%c_Xh, &
376 gjp_param_a, gjp_param_b)
377 call this%gradient_jump_penalty_w%init(params, this%dm_Xh, this%c_Xh, &
378 gjp_param_a, gjp_param_b)
379 end if
380
381 call neko_field_registry%add_field(this%dm_Xh, 'u_e')
382 call neko_field_registry%add_field(this%dm_Xh, 'v_e')
383 call neko_field_registry%add_field(this%dm_Xh, 'w_e')
384 this%u_e => neko_field_registry%get_field('u_e')
385 this%v_e => neko_field_registry%get_field('v_e')
386 this%w_e => neko_field_registry%get_field('w_e')
387
388 call neko_log%end_section()
389
390 end subroutine fluid_scheme_init_base
391
392 subroutine fluid_scheme_free(this)
393 class(fluid_scheme_incompressible_t), intent(inout) :: this
394
395 !
396 ! Free everything related to field_dirichlet BCs
397 !
398
399 call this%Xh%free()
400
401 if (allocated(this%ksp_vel)) then
402 call krylov_solver_destroy(this%ksp_vel)
403 deallocate(this%ksp_vel)
404 end if
405
406 if (allocated(this%ksp_prs)) then
407 call krylov_solver_destroy(this%ksp_prs)
408 deallocate(this%ksp_prs)
409 end if
410
411 if (allocated(this%pc_vel)) then
412 call precon_destroy(this%pc_vel)
413 deallocate(this%pc_vel)
414 end if
415
416 if (allocated(this%pc_prs)) then
417 call precon_destroy(this%pc_prs)
418 deallocate(this%pc_prs)
419 end if
420
421 call this%source_term%free()
422
423 call this%gs_Xh%free()
424
425 call this%c_Xh%free()
426
427 call this%scratch%free()
428
429 nullify(this%u)
430 nullify(this%v)
431 nullify(this%w)
432 nullify(this%p)
433
434 if (this%variable_material_properties) then
435 nullify(this%u_e)
436 nullify(this%v_e)
437 nullify(this%w_e)
438 end if
439
440 call this%ulag%free()
441 call this%vlag%free()
442 call this%wlag%free()
443
444
445 if (associated(this%f_x)) then
446 call this%f_x%free()
447 end if
448
449 if (associated(this%f_y)) then
450 call this%f_y%free()
451 end if
452
453 if (associated(this%f_z)) then
454 call this%f_z%free()
455 end if
456
457 nullify(this%f_x)
458 nullify(this%f_y)
459 nullify(this%f_z)
460
461 call this%rho_field%free()
462 call this%mu_field%free()
463
464 ! Free gradient jump penalty
465 if (this%if_gradient_jump_penalty .eqv. .true.) then
466 call this%gradient_jump_penalty_u%free()
467 call this%gradient_jump_penalty_v%free()
468 call this%gradient_jump_penalty_w%free()
469 end if
470
471 end subroutine fluid_scheme_free
472
475 subroutine fluid_scheme_validate(this)
476 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
477 ! Variables for retrieving json parameters
478 logical :: logical_val
479
480 if ( (.not. associated(this%u)) .or. &
481 (.not. associated(this%v)) .or. &
482 (.not. associated(this%w)) .or. &
483 (.not. associated(this%p))) then
484 call neko_error('Fields are not registered')
485 end if
486
487 if ( (.not. allocated(this%u%x)) .or. &
488 (.not. allocated(this%v%x)) .or. &
489 (.not. allocated(this%w%x)) .or. &
490 (.not. allocated(this%p%x))) then
491 call neko_error('Fields are not allocated')
492 end if
493
494 if (.not. allocated(this%ksp_vel)) then
495 call neko_error('No Krylov solver for velocity defined')
496 end if
497
498 if (.not. allocated(this%ksp_prs)) then
499 call neko_error('No Krylov solver for pressure defined')
500 end if
501
502 !
503 ! Setup checkpoint structure (if everything is fine)
504 !
505 call this%chkp%init(this%u, this%v, this%w, this%p)
506
507 end subroutine fluid_scheme_validate
508
513 subroutine fluid_scheme_bc_apply_vel(this, t, tstep, strong)
514 class(fluid_scheme_incompressible_t), intent(inout) :: this
515 real(kind=rp), intent(in) :: t
516 integer, intent(in) :: tstep
517 logical, intent(in) :: strong
518
519 call this%bcs_vel%apply_vector(&
520 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), t, tstep, strong)
521 call this%gs_Xh%op(this%u, gs_op_min, glb_cmd_event)
522 call this%gs_Xh%op(this%v, gs_op_min, glb_cmd_event)
523 call this%gs_Xh%op(this%w, gs_op_min, glb_cmd_event)
524 if (neko_bcknd_device .eq. 1) then
526 end if
527
528 call this%bcs_vel%apply_vector(&
529 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), t, tstep, strong)
530 call this%gs_Xh%op(this%u, gs_op_max, glb_cmd_event)
531 call this%gs_Xh%op(this%v, gs_op_max, glb_cmd_event)
532 call this%gs_Xh%op(this%w, gs_op_max, glb_cmd_event)
533 if (neko_bcknd_device .eq. 1) then
535 end if
536
537 end subroutine fluid_scheme_bc_apply_vel
538
541 subroutine fluid_scheme_bc_apply_prs(this, t, tstep)
542 class(fluid_scheme_incompressible_t), intent(inout) :: this
543 real(kind=rp), intent(in) :: t
544 integer, intent(in) :: tstep
545
546 call this%bcs_prs%apply(this%p, t, tstep)
547 call this%gs_Xh%op(this%p,gs_op_min, glb_cmd_event)
548 if (neko_bcknd_device .eq. 1) then
550 end if
551
552 call this%bcs_prs%apply(this%p, t, tstep)
553 call this%gs_Xh%op(this%p,gs_op_max, glb_cmd_event)
554 if (neko_bcknd_device .eq. 1) then
556 end if
557
558
559 end subroutine fluid_scheme_bc_apply_prs
560
563 subroutine fluid_scheme_solver_factory(ksp, n, solver, &
564 max_iter, abstol, monitor)
565 class(ksp_t), allocatable, target, intent(inout) :: ksp
566 integer, intent(in), value :: n
567 character(len=*), intent(in) :: solver
568 integer, intent(in) :: max_iter
569 real(kind=rp), intent(in) :: abstol
570 logical, intent(in) :: monitor
571
572 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
573 monitor = monitor)
574
575 end subroutine fluid_scheme_solver_factory
576
578 subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, &
579 pctype)
580 class(fluid_scheme_incompressible_t), intent(inout) :: this
581 class(pc_t), allocatable, target, intent(inout) :: pc
582 class(ksp_t), target, intent(inout) :: ksp
583 type(coef_t), target, intent(in) :: coef
584 type(dofmap_t), target, intent(in) :: dof
585 type(gs_t), target, intent(inout) :: gs
586 type(bc_list_t), target, intent(inout) :: bclst
587 character(len=*) :: pctype
588
589 call precon_factory(pc, pctype)
590
591 select type (pcp => pc)
592 type is (jacobi_t)
593 call pcp%init(coef, dof, gs)
594 type is (sx_jacobi_t)
595 call pcp%init(coef, dof, gs)
596 type is (device_jacobi_t)
597 call pcp%init(coef, dof, gs)
598 type is (hsmg_t)
599 if (len_trim(pctype) .gt. 4) then
600 if (index(pctype, '+') .eq. 5) then
601 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst, &
602 trim(pctype(6:)))
603 else
604 call neko_error('Unknown coarse grid solver')
605 end if
606 else
607 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
608 end if
609 type is (phmg_t)
610 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
611 end select
612
613 call ksp%set_pc(pc)
614
615 end subroutine fluid_scheme_precon_factory
616
618 function fluid_compute_cfl(this, dt) result(c)
619 class(fluid_scheme_incompressible_t), intent(in) :: this
620 real(kind=rp), intent(in) :: dt
621 real(kind=rp) :: c
622
623 c = cfl(dt, this%u%x, this%v%x, this%w%x, &
624 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
625
626 end function fluid_compute_cfl
627
628
630 subroutine fluid_scheme_update_material_properties(this)
631 class(fluid_scheme_incompressible_t), intent(inout) :: this
632 type(field_t), pointer :: nut
633 integer :: n
634
635 this%mu_field = this%mu
636 if (this%variable_material_properties) then
637 nut => neko_field_registry%get_field(this%nut_field_name)
638 n = nut%size()
639 call field_add2s2(this%mu_field, nut, this%rho, n)
640 end if
641 end subroutine fluid_scheme_update_material_properties
642
646 subroutine fluid_scheme_set_material_properties(this, params, user)
647 class(fluid_scheme_incompressible_t), intent(inout) :: this
648 type(json_file), intent(inout) :: params
649 type(user_t), target, intent(in) :: user
650 character(len=LOG_SIZE) :: log_buf
651 ! A local pointer that is needed to make Intel happy
652 procedure(user_material_properties), pointer :: dummy_mp_ptr
653 logical :: nondimensional
654 real(kind=rp) :: dummy_lambda, dummy_cp
655
656 dummy_mp_ptr => dummy_user_material_properties
657
658 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
659
660 write(log_buf, '(A)') "Material properties must be set in the user&
661 & file!"
662 call neko_log%message(log_buf)
663 call user%material_properties(0.0_rp, 0, this%rho, this%mu, &
664 dummy_cp, dummy_lambda, params)
665 else
666 ! Incorrect user input
667 if (params%valid_path('case.fluid.Re') .and. &
668 (params%valid_path('case.fluid.mu') .or. &
669 params%valid_path('case.fluid.rho'))) then
670 call neko_error("To set the material properties for the fluid,&
671 & either provide Re OR mu and rho in the case file.")
672
673 ! Non-dimensional case
674 else if (params%valid_path('case.fluid.Re')) then
675
676 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
677 & input.'
678 call neko_log%message(log_buf, lvl = neko_log_verbose)
679 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
680 & 1/Re.'
681 call neko_log%message(log_buf, lvl = neko_log_verbose)
682
683 ! Read Re into mu for further manipulation.
684 call json_get(params, 'case.fluid.Re', this%mu)
685 write(log_buf, '(A)') 'Read non-dimensional material properties'
686 call neko_log%message(log_buf)
687 write(log_buf, '(A,ES13.6)') 'Re :', this%mu
688 call neko_log%message(log_buf)
689
690 ! Set rho to 1 since the setup is non-dimensional.
691 this%rho = 1.0_rp
692 ! Invert the Re to get viscosity.
693 this%mu = 1.0_rp/this%mu
694 ! Dimensional case
695 else
696 call json_get(params, 'case.fluid.mu', this%mu)
697 call json_get(params, 'case.fluid.rho', this%rho)
698 end if
699
700 end if
701 end subroutine fluid_scheme_set_material_properties
702
Abstract interface to sets rho and mu.
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.
Abstract interface for setting material properties.
Abstract interface defining a user defined inflow condition (pointwise)
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
Jacobi preconditioner accelerator backend.
subroutine, public device_add2s2(a_d, b_d, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_cfill(a_d, c, n)
Set all elements to a constant c .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_event_sync(event)
Synchronize an event.
Definition device.F90:1244
subroutine, public device_stream_wait_event(stream, event, flags)
Synchronize a device stream with an event.
Definition device.F90:1138
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:50
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:56
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines inflow dirichlet conditions.
Defines user dirichlet condition for a scalar field.
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
subroutine, public field_add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Stores a series fields.
Defines a field.
Definition field.f90:34
subroutine fluid_scheme_set_material_properties(this, params, user)
Sets rho and mu.
subroutine fluid_scheme_validate(this)
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
real(kind=rp) function fluid_compute_cfl(this, dt)
Compute CFL.
subroutine fluid_scheme_solver_factory(ksp, n, solver, max_iter, abstol, monitor)
Initialize a linear solver.
subroutine fluid_scheme_bc_apply_prs(this, t, tstep)
Apply all boundary conditions defined for pressure.
subroutine fluid_scheme_bc_apply_vel(this, t, tstep, strong)
Apply all boundary conditions defined for velocity Here we perform additional gs operations to take c...
subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, pctype)
Initialize a Krylov preconditioner.
subroutine fluid_scheme_update_material_properties(this)
Update the values of mu_field if necessary.
subroutine fluid_scheme_init_base(this, msh, lx, params, scheme, user, kspv_init)
Initialise a fluid scheme.
Implements the fluid_source_term_t type.
Computes various statistics for the fluid fields. We use the Reynolds decomposition for a field u = ...
Gather-scatter.
Implements gradient_jump_penalty_t.
Krylov preconditioner.
Definition pc_hsmg.f90:61
Jacobi preconditioner.
Definition pc_jacobi.f90:34
Utilities for retrieving parameters from the case files.
subroutine, public json_extract_object(json, name, object)
Extract object as a separate JSON dictionary.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
integer, parameter, public ksp_max_iter
Maximum number of iters.
Definition krylov.f90:51
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_verbose
Verbose log level.
Definition log.f90:71
type(log_t), public neko_log
Global log stream.
Definition log.f90:65
integer, parameter, public log_size
Definition log.f90:42
Definition math.f90:60
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:359
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:347
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:672
Defines a mean flow field.
Definition mean_flow.f90:34
Defines a mean squared flow field.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:56
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:58
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
real(kind=rp) function, public cfl(dt, u, v, w, xh, coef, nelv, gdim)
Hybrid ph-multigrid preconditioner.
Definition phmg.f90:34
Krylov preconditioner.
Definition precon.f90:34
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
Defines a shear stress boundary condition for a vector field. Maintainer: Timofey Mukha.
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:48
Defines a container for all statistics.
Jacobi preconditioner SX-Aurora backend.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
subroutine, public dummy_user_material_properties(t, tstep, rho, mu, cp, lambda, params)
Defines inflow dirichlet conditions.
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:266
Defines a zero-valued Dirichlet boundary condition.
Base type for a boundary condition.
Definition bc.f90:57
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:47
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:47
User defined dirichlet condition, for which the user can work with an entire field....
Extension of the user defined dirichlet condition field_dirichlet
Base type of all fluid formulations.
Wrapper contaning and executing the fluid source terms.
Defines a jacobi preconditioner.
Definition pc_jacobi.f90:45
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:68
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
A shear stress boundary condition.
The function space for the SEM solution fields.
Definition space.f90:62
Defines a jacobi preconditioner for SX-Aurora.
User defined dirichlet condition for velocity.
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...