Loading [MathJax]/jax/output/HTML-CSS/config.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, ksp_max_iter
49 use coefs, only: coef_t
51 use dirichlet, only : dirichlet_t
52 use jacobi, only : jacobi_t
53 use sx_jacobi, only : sx_jacobi_t
55 use hsmg, only : hsmg_t
56 use phmg, only : phmg_t
57 use precon, only : pc_t, precon_factory, precon_destroy
58 use fluid_stats, only : fluid_stats_t
59 use bc, only : bc_t
60 use bc_list, only : bc_list_t
62 use math, only : cfill, add2s2, glsum
64 use operators, only : cfl
69 use json_module, only : json_file, json_core, json_value
73 use utils, only : neko_error, neko_warning
80 implicit none
81 private
82
87 class(ksp_t), allocatable :: ksp_vel
88 class(ksp_t), allocatable :: ksp_prs
89 class(pc_t), allocatable :: pc_vel
90 class(pc_t), allocatable :: pc_prs
91 integer :: vel_projection_dim
92 integer :: pr_projection_dim
93 integer :: vel_projection_activ_step
94 integer :: pr_projection_activ_step
95 logical :: strict_convergence
97 logical :: if_gradient_jump_penalty
98 type(gradient_jump_penalty_t) :: gradient_jump_penalty_u
99 type(gradient_jump_penalty_t) :: gradient_jump_penalty_v
100 type(gradient_jump_penalty_t) :: gradient_jump_penalty_w
102 type(field_t), pointer :: u_e => null()
103 type(field_t), pointer :: v_e => null()
104 type(field_t), pointer :: w_e => null()
105
106 type(mean_flow_t) :: mean
108 type(mean_sqr_flow_t) :: mean_sqr
109 logical :: forced_flow_rate = .false.
110
112 character(len=:), allocatable :: nut_field_name
113
115 integer(kind=i8) :: glb_n_points
117 integer(kind=i8) :: glb_unique_points
118 type(scratch_registry_t) :: scratch
119 contains
121 procedure, pass(this) :: init_base => fluid_scheme_init_base
122 procedure, pass(this) :: scheme_free => fluid_scheme_free
124 procedure, pass(this) :: validate => fluid_scheme_validate
126 procedure, pass(this) :: bc_apply_vel => fluid_scheme_bc_apply_vel
128 procedure, pass(this) :: bc_apply_prs => fluid_scheme_bc_apply_prs
130 procedure, pass(this) :: compute_cfl => fluid_compute_cfl
132 procedure, pass(this) :: set_material_properties => &
134
136 procedure, pass(this) :: update_material_properties => &
139 procedure, nopass :: solver_factory => fluid_scheme_solver_factory
141 procedure, pass(this) :: precon_factory_ => fluid_scheme_precon_factory
143
144 interface
145
146 module subroutine fluid_scheme_factory(object, type_name)
147 class(fluid_scheme_base_t), intent(inout), allocatable :: object
148 character(len=*) :: type_name
149 end subroutine fluid_scheme_factory
150 end interface
151
152 public :: fluid_scheme_incompressible_t, fluid_scheme_factory
153
154contains
155
157 subroutine fluid_scheme_init_base(this, msh, lx, params, scheme, user, &
158 kspv_init)
159 implicit none
160 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
161 type(mesh_t), target, intent(inout) :: msh
162 integer, intent(in) :: lx
163 character(len=*), intent(in) :: scheme
164 type(json_file), target, intent(inout) :: params
165 type(user_t), target, intent(in) :: user
166 logical, intent(in) :: kspv_init
167 type(dirichlet_t) :: bdry_mask
168 character(len=LOG_SIZE) :: log_buf
169 real(kind=rp), allocatable :: real_vec(:)
170 real(kind=rp) :: real_val, kappa, b, z0
171 logical :: logical_val
172 integer :: integer_val, ierr
173 type(json_file) :: wm_json
174 character(len=:), allocatable :: string_val1, string_val2
175 real(kind=rp) :: gjp_param_a, gjp_param_b
176
177 !
178 ! SEM simulation fundamentals
179 !
180
181 this%msh => msh
182
183 if (msh%gdim .eq. 2) then
184 call this%Xh%init(gll, lx, lx)
185 else
186 call this%Xh%init(gll, lx, lx, lx)
187 end if
188
189 call this%dm_Xh%init(msh, this%Xh)
190
191 call this%gs_Xh%init(this%dm_Xh)
192
193 call this%c_Xh%init(this%gs_Xh)
194
195 ! Local scratch registry
196 this%scratch = scratch_registry_t(this%dm_Xh, 10, 2)
197
198 ! Assign a name
199 call json_get_or_default(params, 'case.fluid.name', this%name, "fluid")
200
201 !
202 ! First section of fluid log
203 !
204
205 call neko_log%section('Fluid')
206 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
207 call neko_log%message(log_buf)
208 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
209 call neko_log%message(log_buf)
210
211 !
212 ! Material properties
213 !
214 call this%set_material_properties(params, user)
215
216 !
217 ! Turbulence modelling and variable material properties
218 !
219 if (params%valid_path('case.fluid.nut_field')) then
220 call json_get(params, 'case.fluid.nut_field', this%nut_field_name)
221 this%variable_material_properties = .true.
222 else
223 this%nut_field_name = ""
224 end if
225
226 ! Fill mu and rho field with the physical value
227
228 call this%mu_field%init(this%dm_Xh, "mu")
229 call this%rho_field%init(this%dm_Xh, "rho")
230 call field_cfill(this%mu_field, this%mu, this%mu_field%size())
231 call field_cfill(this%rho_field, this%rho, this%mu_field%size())
232
233 ! Since mu, rho is a field, and the none-stress simulation fetches
234 ! data from the host arrays, we need to mirror the constant
235 ! material properties on the host
236 if (neko_bcknd_device .eq. 1) then
237 call cfill(this%mu_field%x, this%mu, this%mu_field%size())
238 call cfill(this%rho_field%x, this%rho, this%rho_field%size())
239 end if
240
241
242 ! Projection spaces
243 call json_get_or_default(params, &
244 'case.fluid.velocity_solver.projection_space_size', &
245 this%vel_projection_dim, 0)
246 call json_get_or_default(params, &
247 'case.fluid.pressure_solver.projection_space_size', &
248 this%pr_projection_dim, 0)
249 call json_get_or_default(params, &
250 'case.fluid.velocity_solver.projection_hold_steps', &
251 this%vel_projection_activ_step, 5)
252 call json_get_or_default(params, &
253 'case.fluid.pressure_solver.projection_hold_steps', &
254 this%pr_projection_activ_step, 5)
255
256
257 call json_get_or_default(params, 'case.fluid.freeze', this%freeze, .false.)
258
259 if (params%valid_path("case.fluid.flow_rate_force")) then
260 this%forced_flow_rate = .true.
261 end if
262
263
264 if (lx .lt. 10) then
265 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
266 else if (lx .ge. 10) then
267 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
268 else
269 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
270 end if
271 call neko_log%message(log_buf)
272 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
273 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
274
275 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
276 call neko_log%message(log_buf)
277 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
278 call neko_log%message(log_buf)
279
280 write(log_buf, '(A,ES13.6)') 'rho :', this%rho
281 call neko_log%message(log_buf)
282 write(log_buf, '(A,ES13.6)') 'mu :', this%mu
283 call neko_log%message(log_buf)
284
285 call json_get(params, 'case.numerics.dealias', logical_val)
286 write(log_buf, '(A, L1)') 'Dealias : ', logical_val
287 call neko_log%message(log_buf)
288
289 write(log_buf, '(A, L1)') 'LES : ', this%variable_material_properties
290 call neko_log%message(log_buf)
291
292 call json_get_or_default(params, 'case.output_boundary', logical_val, &
293 .false.)
294 write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
295 call neko_log%message(log_buf)
296
297 !
298 ! Setup right-hand side fields.
299 !
300 allocate(this%f_x)
301 allocate(this%f_y)
302 allocate(this%f_z)
303 call this%f_x%init(this%dm_Xh, fld_name = "fluid_rhs_x")
304 call this%f_y%init(this%dm_Xh, fld_name = "fluid_rhs_y")
305 call this%f_z%init(this%dm_Xh, fld_name = "fluid_rhs_z")
306
307 ! Initialize the source term
308 call this%source_term%init(this%f_x, this%f_y, this%f_z, this%c_Xh, user)
309 call this%source_term%add(params, 'case.fluid.source_terms')
310
311 ! Initialize velocity solver
312 if (kspv_init) then
313 call neko_log%section("Velocity solver")
314 call json_get_or_default(params, &
315 'case.fluid.velocity_solver.max_iterations', &
316 integer_val, ksp_max_iter)
317 call json_get(params, 'case.fluid.velocity_solver.type', string_val1)
318 call json_get(params, 'case.fluid.velocity_solver.preconditioner', &
319 string_val2)
320 call json_get(params, 'case.fluid.velocity_solver.absolute_tolerance', &
321 real_val)
322 call json_get_or_default(params, &
323 'case.fluid.velocity_solver.monitor', &
324 logical_val, .false.)
325
326 call neko_log%message('Type : ('// trim(string_val1) // &
327 ', ' // trim(string_val2) // ')')
328
329 write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
330 call neko_log%message(log_buf)
331 call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
332 string_val1, integer_val, real_val, logical_val)
333 call this%precon_factory_(this%pc_vel, this%ksp_vel, &
334 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, string_val2)
335 call neko_log%end_section()
336 end if
337
338 ! Strict convergence for the velocity solver
339 call json_get_or_default(params, 'case.fluid.strict_convergence', &
340 this%strict_convergence, .false.)
341
342 ! Assign velocity fields
343 call neko_field_registry%add_field(this%dm_Xh, 'u')
344 call neko_field_registry%add_field(this%dm_Xh, 'v')
345 call neko_field_registry%add_field(this%dm_Xh, 'w')
346 this%u => neko_field_registry%get_field('u')
347 this%v => neko_field_registry%get_field('v')
348 this%w => neko_field_registry%get_field('w')
349
350 !! Initialize time-lag fields
351 call this%ulag%init(this%u, 2)
352 call this%vlag%init(this%v, 2)
353 call this%wlag%init(this%w, 2)
354
355 ! Initiate gradient jump penalty
356 call json_get_or_default(params, &
357 'case.fluid.gradient_jump_penalty.enabled',&
358 this%if_gradient_jump_penalty, .false.)
359
360 if (this%if_gradient_jump_penalty .eqv. .true.) then
361 if ((this%dm_Xh%xh%lx - 1) .eq. 1) then
362 call json_get_or_default(params, &
363 'case.fluid.gradient_jump_penalty.tau',&
364 gjp_param_a, 0.02_rp)
365 gjp_param_b = 0.0_rp
366 else
367 call json_get_or_default(params, &
368 'case.fluid.gradient_jump_penalty.scaling_factor',&
369 gjp_param_a, 0.8_rp)
370 call json_get_or_default(params, &
371 'case.fluid.gradient_jump_penalty.scaling_exponent',&
372 gjp_param_b, 4.0_rp)
373 end if
374 call this%gradient_jump_penalty_u%init(params, this%dm_Xh, this%c_Xh, &
375 gjp_param_a, gjp_param_b)
376 call this%gradient_jump_penalty_v%init(params, this%dm_Xh, this%c_Xh, &
377 gjp_param_a, gjp_param_b)
378 call this%gradient_jump_penalty_w%init(params, this%dm_Xh, this%c_Xh, &
379 gjp_param_a, gjp_param_b)
380 end if
381
382 call neko_field_registry%add_field(this%dm_Xh, 'u_e')
383 call neko_field_registry%add_field(this%dm_Xh, 'v_e')
384 call neko_field_registry%add_field(this%dm_Xh, 'w_e')
385 this%u_e => neko_field_registry%get_field('u_e')
386 this%v_e => neko_field_registry%get_field('v_e')
387 this%w_e => neko_field_registry%get_field('w_e')
388
389 call neko_log%end_section()
390
391 end subroutine fluid_scheme_init_base
392
393 subroutine fluid_scheme_free(this)
394 class(fluid_scheme_incompressible_t), intent(inout) :: this
395
396 call this%Xh%free()
397
398 if (allocated(this%ksp_vel)) then
399 call this%ksp_vel%free()
400 deallocate(this%ksp_vel)
401 end if
402
403 if (allocated(this%ksp_prs)) then
404 call this%ksp_prs%free()
405 deallocate(this%ksp_prs)
406 end if
407
408 if (allocated(this%pc_vel)) then
409 call precon_destroy(this%pc_vel)
410 deallocate(this%pc_vel)
411 end if
412
413 if (allocated(this%pc_prs)) then
414 call precon_destroy(this%pc_prs)
415 deallocate(this%pc_prs)
416 end if
417
418 call this%source_term%free()
419
420 call this%gs_Xh%free()
421
422 call this%c_Xh%free()
423
424 call this%scratch%free()
425
426 nullify(this%u)
427 nullify(this%v)
428 nullify(this%w)
429 nullify(this%p)
430
431 if (this%variable_material_properties) then
432 nullify(this%u_e)
433 nullify(this%v_e)
434 nullify(this%w_e)
435 end if
436
437 call this%ulag%free()
438 call this%vlag%free()
439 call this%wlag%free()
440
441
442 if (associated(this%f_x)) then
443 call this%f_x%free()
444 end if
445
446 if (associated(this%f_y)) then
447 call this%f_y%free()
448 end if
449
450 if (associated(this%f_z)) then
451 call this%f_z%free()
452 end if
453
454 nullify(this%f_x)
455 nullify(this%f_y)
456 nullify(this%f_z)
457
458 call this%rho_field%free()
459 call this%mu_field%free()
460
461 ! Free gradient jump penalty
462 if (this%if_gradient_jump_penalty .eqv. .true.) then
463 call this%gradient_jump_penalty_u%free()
464 call this%gradient_jump_penalty_v%free()
465 call this%gradient_jump_penalty_w%free()
466 end if
467
468 end subroutine fluid_scheme_free
469
472 subroutine fluid_scheme_validate(this)
473 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
474 ! Variables for retrieving json parameters
475 logical :: logical_val
476
477 if ( (.not. associated(this%u)) .or. &
478 (.not. associated(this%v)) .or. &
479 (.not. associated(this%w)) .or. &
480 (.not. associated(this%p))) then
481 call neko_error('Fields are not registered')
482 end if
483
484 if ( (.not. allocated(this%u%x)) .or. &
485 (.not. allocated(this%v%x)) .or. &
486 (.not. allocated(this%w%x)) .or. &
487 (.not. allocated(this%p%x))) then
488 call neko_error('Fields are not allocated')
489 end if
490
491 if (.not. allocated(this%ksp_vel)) then
492 call neko_error('No Krylov solver for velocity defined')
493 end if
494
495 if (.not. allocated(this%ksp_prs)) then
496 call neko_error('No Krylov solver for pressure defined')
497 end if
498
499 end subroutine fluid_scheme_validate
500
505 subroutine fluid_scheme_bc_apply_vel(this, t, tstep, strong)
506 class(fluid_scheme_incompressible_t), intent(inout) :: this
507 real(kind=rp), intent(in) :: t
508 integer, intent(in) :: tstep
509 logical, intent(in) :: strong
510
511 integer :: i
512 class(bc_t), pointer :: b
513 b => null()
514
515 call this%bcs_vel%apply_vector(&
516 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), t, tstep, strong)
517 call this%gs_Xh%op(this%u, gs_op_min, glb_cmd_event)
519 call this%gs_Xh%op(this%v, gs_op_min, glb_cmd_event)
521 call this%gs_Xh%op(this%w, gs_op_min, glb_cmd_event)
523
524
525 call this%bcs_vel%apply_vector(&
526 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), t, tstep, strong)
527 call this%gs_Xh%op(this%u, gs_op_max, glb_cmd_event)
529 call this%gs_Xh%op(this%v, gs_op_max, glb_cmd_event)
531 call this%gs_Xh%op(this%w, gs_op_max, glb_cmd_event)
533
534 do i = 1, this%bcs_vel%size()
535 b => this%bcs_vel%get(i)
536 b%updated = .false.
537 end do
538 nullify(b)
539
540 end subroutine fluid_scheme_bc_apply_vel
541
544 subroutine fluid_scheme_bc_apply_prs(this, t, tstep)
545 class(fluid_scheme_incompressible_t), intent(inout) :: this
546 real(kind=rp), intent(in) :: t
547 integer, intent(in) :: tstep
548
549 integer :: i
550 class(bc_t), pointer :: b
551 b => null()
552
553 call this%bcs_prs%apply(this%p, t, tstep)
554 call this%gs_Xh%op(this%p,gs_op_min, glb_cmd_event)
556
557 call this%bcs_prs%apply(this%p, t, tstep)
558 call this%gs_Xh%op(this%p,gs_op_max, glb_cmd_event)
560
561 do i = 1, this%bcs_prs%size()
562 b => this%bcs_prs%get(i)
563 b%updated = .false.
564 end do
565 nullify(b)
566
567 end subroutine fluid_scheme_bc_apply_prs
568
571 subroutine fluid_scheme_solver_factory(ksp, n, solver, &
572 max_iter, abstol, monitor)
573 class(ksp_t), allocatable, target, intent(inout) :: ksp
574 integer, intent(in), value :: n
575 character(len=*), intent(in) :: solver
576 integer, intent(in) :: max_iter
577 real(kind=rp), intent(in) :: abstol
578 logical, intent(in) :: monitor
579
580 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
581 monitor = monitor)
582
583 end subroutine fluid_scheme_solver_factory
584
586 subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, &
587 pctype)
588 class(fluid_scheme_incompressible_t), intent(inout) :: this
589 class(pc_t), allocatable, target, intent(inout) :: pc
590 class(ksp_t), target, intent(inout) :: ksp
591 type(coef_t), target, intent(in) :: coef
592 type(dofmap_t), target, intent(in) :: dof
593 type(gs_t), target, intent(inout) :: gs
594 type(bc_list_t), target, intent(inout) :: bclst
595 character(len=*) :: pctype
596
597 call precon_factory(pc, pctype)
598
599 select type (pcp => pc)
600 type is (jacobi_t)
601 call pcp%init(coef, dof, gs)
602 type is (sx_jacobi_t)
603 call pcp%init(coef, dof, gs)
604 type is (device_jacobi_t)
605 call pcp%init(coef, dof, gs)
606 type is (hsmg_t)
607 if (len_trim(pctype) .gt. 4) then
608 if (index(pctype, '+') .eq. 5) then
609 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst, &
610 trim(pctype(6:)))
611 else
612 call neko_error('Unknown coarse grid solver')
613 end if
614 else
615 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
616 end if
617 type is (phmg_t)
618 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
619 end select
620
621 call ksp%set_pc(pc)
622
623 end subroutine fluid_scheme_precon_factory
624
626 function fluid_compute_cfl(this, dt) result(c)
627 class(fluid_scheme_incompressible_t), intent(in) :: this
628 real(kind=rp), intent(in) :: dt
629 real(kind=rp) :: c
630
631 c = cfl(dt, this%u%x, this%v%x, this%w%x, &
632 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
633
634 end function fluid_compute_cfl
635
636
638 subroutine fluid_scheme_update_material_properties(this)
639 class(fluid_scheme_incompressible_t), intent(inout) :: this
640 type(field_t), pointer :: nut
641 integer :: n
642
643 this%mu_field = this%mu
644 if (this%variable_material_properties) then
645 nut => neko_field_registry%get_field(this%nut_field_name)
646 n = nut%size()
647 call field_add2s2(this%mu_field, nut, this%rho, n)
648 end if
649 end subroutine fluid_scheme_update_material_properties
650
654 subroutine fluid_scheme_set_material_properties(this, params, user)
655 class(fluid_scheme_incompressible_t), intent(inout) :: this
656 type(json_file), intent(inout) :: params
657 type(user_t), target, intent(in) :: user
658 character(len=LOG_SIZE) :: log_buf
659 ! A local pointer that is needed to make Intel happy
660 procedure(user_material_properties), pointer :: dummy_mp_ptr
661 logical :: nondimensional
662 real(kind=rp) :: dummy_lambda, dummy_cp
663
664 dummy_mp_ptr => dummy_user_material_properties
665
666 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
667
668 write(log_buf, '(A)') "Material properties must be set in the user&
669 & file!"
670 call neko_log%message(log_buf)
671 call user%material_properties(0.0_rp, 0, this%rho, this%mu, &
672 dummy_cp, dummy_lambda, params)
673 else
674 ! Incorrect user input
675 if (params%valid_path('case.fluid.Re') .and. &
676 (params%valid_path('case.fluid.mu') .or. &
677 params%valid_path('case.fluid.rho'))) then
678 call neko_error("To set the material properties for the fluid," // &
679 " either provide Re OR mu and rho in the case file.")
680
681 ! Non-dimensional case
682 else if (params%valid_path('case.fluid.Re')) then
683
684 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
685 & input.'
686 call neko_log%message(log_buf, lvl = neko_log_verbose)
687 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
688 & 1/Re.'
689 call neko_log%message(log_buf, lvl = neko_log_verbose)
690
691 ! Read Re into mu for further manipulation.
692 call json_get(params, 'case.fluid.Re', this%mu)
693 write(log_buf, '(A)') 'Read non-dimensional material properties'
694 call neko_log%message(log_buf)
695 write(log_buf, '(A,ES13.6)') 'Re :', this%mu
696 call neko_log%message(log_buf)
697
698 ! Set rho to 1 since the setup is non-dimensional.
699 this%rho = 1.0_rp
700 ! Invert the Re to get viscosity.
701 this%mu = 1.0_rp/this%mu
702 ! Dimensional case
703 else
704 call json_get(params, 'case.fluid.mu', this%mu)
705 call json_get(params, 'case.fluid.rho', this%rho)
706 end if
707
708 end if
709 end subroutine fluid_scheme_set_material_properties
710
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
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
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
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.
A type collecting all the overridable user routines.
User defined dirichlet condition for velocity.
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...