Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_scheme.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2024, 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!
35 use gather_scatter, only : gs_t
38 use checkpoint, only : chkp_t
39 use mean_flow, only : mean_flow_t
40 use num_types, only : rp, i8
41 use comm
43 use field, only : field_t
44 use space, only : space_t, gll
45 use dofmap, only : dofmap_t
46 use krylov, only : ksp_t, krylov_solver_factory, krylov_solver_destroy, &
48 use coefs, only: coef_t
49 use wall, only : no_slip_wall_t
50 use inflow, only : inflow_t
52 use blasius, only : blasius_t
53 use dirichlet, only : dirichlet_t
55 use symmetry, only : symmetry_t
58 use jacobi, only : jacobi_t
59 use sx_jacobi, only : sx_jacobi_t
61 use hsmg, only : hsmg_t
62 use phmg, only : phmg_t
63 use precon, only : pc_t, precon_factory, precon_destroy
64 use fluid_stats, only : fluid_stats_t
65 use bc, only : bc_t
66 use bc_list, only : bc_list_t
68 use math, only : cfill, add2s2, glsum
71 use operators, only : cfl
75 use json_module, only : json_file
79 use utils, only : neko_error, neko_warning
82 use field_math, only : field_cfill
86 implicit none
87 private
88
90 type, abstract :: fluid_scheme_t
91 type(field_t), pointer :: u => null()
92 type(field_t), pointer :: v => null()
93 type(field_t), pointer :: w => null()
94 type(field_t), pointer :: p => null()
95 type(field_series_t) :: ulag, vlag, wlag
96 type(space_t) :: xh
97 type(dofmap_t) :: dm_xh
98 type(gs_t) :: gs_xh
99 type(coef_t) :: c_xh
103 type(field_t), pointer :: f_x => null()
105 type(field_t), pointer :: f_y => null()
107 type(field_t), pointer :: f_z => null()
108
109 ! Krylov solvers and settings
110 class(ksp_t), allocatable :: ksp_vel
111 class(ksp_t), allocatable :: ksp_prs
112 class(pc_t), allocatable :: pc_vel
113 class(pc_t), allocatable :: pc_prs
114 integer :: vel_projection_dim
115 integer :: pr_projection_dim
116 integer :: vel_projection_activ_step
117 integer :: pr_projection_activ_step
118 logical :: strict_convergence
119
120 type(no_slip_wall_t) :: bc_wall
121 class(bc_t), allocatable :: bc_inflow
122 type(wall_model_bc_t) :: bc_wallmodel
124 logical :: if_gradient_jump_penalty
125 type(gradient_jump_penalty_t) :: gradient_jump_penalty_u
126 type(gradient_jump_penalty_t) :: gradient_jump_penalty_v
127 type(gradient_jump_penalty_t) :: gradient_jump_penalty_w
128
129 ! Attributes for field dirichlet BCs
130 type(field_dirichlet_vector_t) :: user_field_bc_vel
131 type(field_dirichlet_t) :: user_field_bc_prs
132 type(dirichlet_t) :: bc_prs
133 type(dong_outflow_t) :: bc_dong
134 type(symmetry_t) :: bc_sym
135 type(shear_stress_t) :: bc_sh
136 type(bc_list_t) :: bclst_vel
137 type(bc_list_t) :: bclst_vel_neumann
138 type(bc_list_t) :: bclst_prs
139 type(field_t) :: bdry
140 type(json_file), pointer :: params
141 type(mesh_t), pointer :: msh => null()
142 type(chkp_t) :: chkp
143 type(mean_flow_t) :: mean
145 type(mean_sqr_flow_t) :: mean_sqr
146 logical :: forced_flow_rate = .false.
147 logical :: freeze = .false.
149 real(kind=rp) :: mu
151 type(field_t) :: mu_field
153 character(len=:), allocatable :: nut_field_name
155 logical :: variable_material_properties = .false.
157 real(kind=rp) :: rho
159 type(field_t) :: rho_field
161 integer(kind=i8) :: glb_n_points
163 integer(kind=i8) :: glb_unique_points
164 type(scratch_registry_t) :: scratch
166 character(len=NEKO_MSH_MAX_ZLBL_LEN), allocatable :: bc_labels(:)
167 contains
169 procedure, pass(this) :: fluid_scheme_init_all
170 procedure, pass(this) :: fluid_scheme_init_common
173 procedure, pass(this) :: scheme_free => fluid_scheme_free
175 procedure, pass(this) :: validate => fluid_scheme_validate
177 procedure, pass(this) :: bc_apply_vel => fluid_scheme_bc_apply_vel
179 procedure, pass(this) :: bc_apply_prs => fluid_scheme_bc_apply_prs
181 procedure, pass(this) :: set_usr_inflow => fluid_scheme_set_usr_inflow
183 procedure, pass(this) :: compute_cfl => fluid_compute_cfl
185 procedure, pass(this) :: set_material_properties => &
188 procedure(fluid_scheme_init_intrf), pass(this), deferred :: init
190 procedure(fluid_scheme_free_intrf), pass(this), deferred :: free
192 procedure(fluid_scheme_step_intrf), pass(this), deferred :: step
194 procedure(fluid_scheme_restart_intrf), pass(this), deferred :: restart
195 procedure, private, pass(this) :: set_bc_type_output => &
198 procedure, pass(this) :: update_material_properties => &
200 end type fluid_scheme_t
201
203 abstract interface
204 subroutine fluid_scheme_init_intrf(this, msh, lx, params, user, &
205 time_scheme)
206 import fluid_scheme_t
207 import json_file
208 import mesh_t
209 import user_t
211 class(fluid_scheme_t), target, intent(inout) :: this
212 type(mesh_t), target, intent(inout) :: msh
213 integer, intent(in) :: lx
214 type(json_file), target, intent(inout) :: params
215 type(user_t), target, intent(in) :: user
216 type(time_scheme_controller_t), target, intent(in) :: time_scheme
217 end subroutine fluid_scheme_init_intrf
218 end interface
219
221 abstract interface
222 subroutine fluid_scheme_free_intrf(this)
223 import fluid_scheme_t
224 class(fluid_scheme_t), intent(inout) :: this
225 end subroutine fluid_scheme_free_intrf
226 end interface
227
229 abstract interface
230 subroutine fluid_scheme_step_intrf(this, t, tstep, dt, ext_bdf, &
231 dt_controller)
232 import fluid_scheme_t
235 import rp
236 class(fluid_scheme_t), target, intent(inout) :: this
237 real(kind=rp), intent(in) :: t
238 integer, intent(in) :: tstep
239 real(kind=rp), intent(in) :: dt
240 type(time_scheme_controller_t), intent(in) :: ext_bdf
241 type(time_step_controller_t), intent(in) :: dt_controller
242 end subroutine fluid_scheme_step_intrf
243 end interface
244
246 abstract interface
247 subroutine fluid_scheme_restart_intrf(this, dtlag, tlag)
248 import fluid_scheme_t
249 import rp
250 class(fluid_scheme_t), target, intent(inout) :: this
251 real(kind=rp) :: dtlag(10), tlag(10)
252
253 end subroutine fluid_scheme_restart_intrf
254 end interface
255
256 interface
257
258 module subroutine fluid_scheme_factory(object, type_name)
259 class(fluid_scheme_t), intent(inout), allocatable :: object
260 character(len=*) :: type_name
261 end subroutine fluid_scheme_factory
262 end interface
263
264 public :: fluid_scheme_t, fluid_scheme_factory
265
266contains
267
269 subroutine fluid_scheme_init_common(this, msh, lx, params, scheme, user, &
270 kspv_init)
271 implicit none
272 class(fluid_scheme_t), target, intent(inout) :: this
273 type(mesh_t), target, intent(inout) :: msh
274 integer, intent(in) :: lx
275 character(len=*), intent(in) :: scheme
276 type(json_file), target, intent(inout) :: params
277 type(user_t), target, intent(in) :: user
278 logical, intent(in) :: kspv_init
279 type(dirichlet_t) :: bdry_mask
280 character(len=LOG_SIZE) :: log_buf
281 real(kind=rp), allocatable :: real_vec(:)
282 real(kind=rp) :: real_val, kappa, b, z0
283 logical :: logical_val
284 integer :: integer_val, ierr
285 type(json_file) :: wm_json
286 character(len=:), allocatable :: string_val1, string_val2
287
288 !
289 ! SEM simulation fundamentals
290 !
291
292 this%msh => msh
293
294 if (msh%gdim .eq. 2) then
295 call this%Xh%init(gll, lx, lx)
296 else
297 call this%Xh%init(gll, lx, lx, lx)
298 end if
299
300 call this%dm_Xh%init(msh, this%Xh)
301
302 call this%gs_Xh%init(this%dm_Xh)
303
304 call this%c_Xh%init(this%gs_Xh)
305
306 ! Local scratch registry
307 this%scratch = scratch_registry_t(this%dm_Xh, 10, 2)
308
309 ! Case parameters
310 this%params => params
311
312
313 !
314 ! First section of fluid log
315 !
316
317 call neko_log%section('Fluid')
318 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
319 call neko_log%message(log_buf)
320
321 !
322 ! Material properties
323 !
324 call this%set_material_properties(params, user)
325
326 !
327 ! Turbulence modelling and variable material properties
328 !
329 if (params%valid_path('case.fluid.nut_field')) then
330 call json_get(params, 'case.fluid.nut_field', this%nut_field_name)
331 this%variable_material_properties = .true.
332 else
333 this%nut_field_name = ""
334 end if
335
336 ! Fill mu and rho field with the physical value
337
338 call this%mu_field%init(this%dm_Xh, "mu")
339 call this%rho_field%init(this%dm_Xh, "mu")
340 call field_cfill(this%mu_field, this%mu, this%mu_field%size())
341 call field_cfill(this%rho_field, this%rho, this%mu_field%size())
342
343 ! Since mu, rho is a field, and the none-stress simulation fetches
344 ! data from the host arrays, we need to mirror the constant
345 ! material properties on the host
346 if (neko_bcknd_device .eq. 1) then
347 call cfill(this%mu_field%x, this%mu, this%mu_field%size())
348 call cfill(this%rho_field%x, this%rho, this%rho_field%size())
349 end if
350
351
352 ! Projection spaces
353 call json_get_or_default(params, &
354 'case.fluid.velocity_solver.projection_space_size', &
355 this%vel_projection_dim, 20)
356 call json_get_or_default(params, &
357 'case.fluid.pressure_solver.projection_space_size', &
358 this%pr_projection_dim, 20)
359 call json_get_or_default(params, &
360 'case.fluid.velocity_solver.projection_hold_steps', &
361 this%vel_projection_activ_step, 5)
362 call json_get_or_default(params, &
363 'case.fluid.pressure_solver.projection_hold_steps', &
364 this%pr_projection_activ_step, 5)
365
366
367 call json_get_or_default(params, 'case.fluid.freeze', this%freeze, .false.)
368
369 if (params%valid_path("case.fluid.flow_rate_force")) then
370 this%forced_flow_rate = .true.
371 end if
372
373
374 if (lx .lt. 10) then
375 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
376 else if (lx .ge. 10) then
377 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
378 else
379 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
380 end if
381 call neko_log%message(log_buf)
382 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
383 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
384
385 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
386 call neko_log%message(log_buf)
387 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
388 call neko_log%message(log_buf)
389
390 write(log_buf, '(A,ES13.6)') 'rho :', this%rho
391 call neko_log%message(log_buf)
392 write(log_buf, '(A,ES13.6)') 'mu :', this%mu
393 call neko_log%message(log_buf)
394
395 call json_get(params, 'case.numerics.dealias', logical_val)
396 write(log_buf, '(A, L1)') 'Dealias : ', logical_val
397 call neko_log%message(log_buf)
398
399 write(log_buf, '(A, L1)') 'LES : ', this%variable_material_properties
400 call neko_log%message(log_buf)
401
402 call json_get_or_default(params, 'case.output_boundary', logical_val, &
403 .false.)
404 write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
405 call neko_log%message(log_buf)
406 if (logical_val) then
407 write(log_buf, '(A)') 'bdry keys: '
408 call neko_log%message(log_buf)
409 write(log_buf, '(A)') 'No-slip wall ''w'' = 1'
410 call neko_log%message(log_buf)
411 write(log_buf, '(A)') 'Inlet/velocity dirichlet ''v'' = 2'
412 call neko_log%message(log_buf)
413 write(log_buf, '(A)') 'Outlet ''o'' = 3'
414 call neko_log%message(log_buf)
415 write(log_buf, '(A)') 'Symmetry ''sym'' = 4'
416 call neko_log%message(log_buf)
417 write(log_buf, '(A)') 'Outlet-normal ''on'' = 5'
418 call neko_log%message(log_buf)
419 write(log_buf, '(A)') 'Periodic = 6'
420 call neko_log%message(log_buf)
421 write(log_buf, '(A)') 'Dirichlet on velocity components: '
422 call neko_log%message(log_buf)
423 write(log_buf, '(A)') ' ''d_vel_u, d_vel_v, d_vel_w'' = 7'
424 call neko_log%message(log_buf)
425 write(log_buf, '(A)') 'Pressure dirichlet ''d_pres'' = 8'
426 call neko_log%message(log_buf)
427 write(log_buf, '(A)') '''d_vel_(u,v,w)'' and ''d_pres'' = 8'
428 call neko_log%message(log_buf)
429 write(log_buf, '(A)') 'Shear stress ''sh'' = 9'
430 call neko_log%message(log_buf)
431 write(log_buf, '(A)') 'Wall modelling ''wm'' = 10'
432 call neko_log%message(log_buf)
433 write(log_buf, '(A)') 'No boundary condition set = 0'
434 call neko_log%message(log_buf)
435 end if
436
437
438 !
439 ! Setup velocity boundary conditions
440 !
441 allocate(this%bc_labels(neko_msh_max_zlbls))
442 this%bc_labels = "not"
443
444 if (params%valid_path('case.fluid.boundary_types')) then
445 call json_get(params, &
446 'case.fluid.boundary_types', &
447 this%bc_labels)
448 end if
449
450 call this%bclst_vel%init()
451
452 call this%bc_sym%init_base(this%c_Xh)
453 call this%bc_sym%mark_zone(msh%sympln)
454 call this%bc_sym%mark_zones_from_list(msh%labeled_zones, &
455 'sym', this%bc_labels)
456 call this%bc_sym%finalize()
457 call this%bc_sym%init(this%c_Xh)
458 call this%bclst_vel%append(this%bc_sym)
459
460 ! Shear stress conditions
461 call this%bc_sh%init_base(this%c_Xh)
462 call this%bc_sh%mark_zones_from_list(msh%labeled_zones, &
463 'sh', this%bc_labels)
464 call this%bc_sh%finalize()
465 ! This marks the dirichlet and neumann conditions inside, and finalizes
466 ! them.
467 call this%bc_sh%init_shear_stress(this%c_Xh)
468
469 call this%bclst_vel%append(this%bc_sh%symmetry)
470
471 ! Read stress value, default to [0 0 0]
472 if (this%bc_sh%msk(0) .gt. 0) then
473 call params%get('case.fluid.shear_stress.value', real_vec, logical_val)
474 if (.not. logical_val .and. this%bc_sh%msk(0) .gt. 0) then
475 call neko_warning("No stress values provided for sh boundaries, &
476 & defaulting to 0. Use fluid.shear_stress.value to set.")
477 allocate(real_vec(3))
478 real_vec = 0.0_rp
479 else if (size(real_vec) .ne. 3) then
480 call neko_error ("The shear stress vector provided in &
481 &fluid.shear_stress.value should have 3 components.")
482 end if
483 call this%bc_sh%set_stress(real_vec(1), real_vec(2), real_vec(3))
484 end if
485
486 call this%bclst_vel_neumann%init()
487 call this%bclst_vel_neumann%append(this%bc_sh)
488
489 !
490 ! Inflow
491 !
492 if (params%valid_path('case.fluid.inflow_condition')) then
493 call json_get(params, 'case.fluid.inflow_condition.type', string_val1)
494 if (trim(string_val1) .eq. "uniform") then
495 allocate(inflow_t::this%bc_inflow)
496 else if (trim(string_val1) .eq. "blasius") then
497 allocate(blasius_t::this%bc_inflow)
498 else if (trim(string_val1) .eq. "user") then
499 allocate(usr_inflow_t::this%bc_inflow)
500 else
501 call neko_error('Invalid inflow condition '//string_val1)
502 end if
503
504 call this%bc_inflow%init_base(this%c_Xh)
505 call this%bc_inflow%mark_zone(msh%inlet)
506 call this%bc_inflow%mark_zones_from_list(msh%labeled_zones, &
507 'v', this%bc_labels)
508 call this%bc_inflow%finalize()
509 call this%bclst_vel%append(this%bc_inflow)
510
511 if (trim(string_val1) .eq. "uniform") then
512 call json_get(params, 'case.fluid.inflow_condition.value', real_vec)
513 select type (bc_if => this%bc_inflow)
514 type is (inflow_t)
515 call bc_if%set_inflow(real_vec)
516 end select
517 else if (trim(string_val1) .eq. "blasius") then
518 select type (bc_if => this%bc_inflow)
519 type is (blasius_t)
520 call json_get(params, 'case.fluid.blasius.delta', real_val)
521 call json_get(params, 'case.fluid.blasius.approximation', &
522 string_val2)
523 call json_get(params, 'case.fluid.blasius.freestream_velocity', &
524 real_vec)
525
526 call bc_if%set_params(real_vec, real_val, string_val2)
527
528 end select
529 else if (trim(string_val1) .eq. "user") then
530 end if
531 end if
532
533 !
534 ! Wall models
535 !
536
537 call this%bc_wallmodel%init_base(this%c_Xh)
538 call this%bc_wallmodel%mark_zones_from_list(msh%labeled_zones,&
539 'wm', this%bc_labels)
540 call this%bc_wallmodel%finalize()
541
542 if (this%bc_wallmodel%msk(0) .gt. 0) then
543 call json_extract_object(params, 'case.fluid.wall_modelling', wm_json)
544 call this%bc_wallmodel%init_wall_model_bc(wm_json, this%mu / this%rho)
545 else
546 call this%bc_wallmodel%shear_stress_t%init_shear_stress(this%c_Xh)
547 end if
548
549 call this%bclst_vel%append(this%bc_wallmodel%symmetry)
550 call this%bclst_vel_neumann%append(this%bc_wallmodel)
551
552 call this%bc_wall%init_base(this%c_Xh)
553 call this%bc_wall%mark_zone(msh%wall)
554 call this%bc_wall%mark_zones_from_list(msh%labeled_zones, &
555 'w', this%bc_labels)
556 call this%bc_wall%finalize()
557 call this%bclst_vel%append(this%bc_wall)
558
559 ! Setup field dirichlet bc for u-velocity
560 call this%user_field_bc_vel%bc_u%init_base(this%c_Xh)
561 call this%user_field_bc_vel%bc_u%mark_zones_from_list(msh%labeled_zones, &
562 'd_vel_u', this%bc_labels)
563 call this%user_field_bc_vel%bc_u%finalize()
564
565 call mpi_allreduce(this%user_field_bc_vel%bc_u%msk(0), integer_val, 1, &
566 mpi_integer, mpi_sum, neko_comm, ierr)
567 if (integer_val .gt. 0) then
568 call this%user_field_bc_vel%bc_u%init_field('d_vel_u')
569 end if
570
571 ! Setup field dirichlet bc for v-velocity
572 call this%user_field_bc_vel%bc_v%init_base(this%c_Xh)
573 call this%user_field_bc_vel%bc_v%mark_zones_from_list(msh%labeled_zones, &
574 'd_vel_v', this%bc_labels)
575 call this%user_field_bc_vel%bc_v%finalize()
576
577 call mpi_allreduce(this%user_field_bc_vel%bc_v%msk(0), integer_val, 1, &
578 mpi_integer, mpi_sum, neko_comm, ierr)
579 if (integer_val .gt. 0) then
580 call this%user_field_bc_vel%bc_v%init_field('d_vel_v')
581 end if
582
583 ! Setup field dirichlet bc for w-velocity
584 call this%user_field_bc_vel%bc_w%init_base(this%c_Xh)
585 call this%user_field_bc_vel%bc_w%mark_zones_from_list(msh%labeled_zones, &
586 'd_vel_w', this%bc_labels)
587 call this%user_field_bc_vel%bc_w%finalize()
588
589 call mpi_allreduce(this%user_field_bc_vel%bc_w%msk(0), integer_val, 1, &
590 mpi_integer, mpi_sum, neko_comm, ierr)
591 if (integer_val .gt. 0) then
592 call this%user_field_bc_vel%bc_w%init_field('d_vel_w')
593 end if
594
595 ! Setup our global field dirichlet bc
596 call this%user_field_bc_vel%init_base(this%c_Xh)
597 call this%user_field_bc_vel%mark_zones_from_list(msh%labeled_zones, &
598 'd_vel_u', this%bc_labels)
599 call this%user_field_bc_vel%mark_zones_from_list(msh%labeled_zones, &
600 'd_vel_v', this%bc_labels)
601 call this%user_field_bc_vel%mark_zones_from_list(msh%labeled_zones, &
602 'd_vel_w', this%bc_labels)
603 call this%user_field_bc_vel%finalize()
604
605 ! Add the field bc to velocity bcs
606 call this%bclst_vel%append(this%user_field_bc_vel)
607
608 !
609 ! Associate our field dirichlet update to the user one.
610 !
611 this%user_field_bc_vel%update => user%user_dirichlet_update
612
613 !
614 ! Initialize field list and bc list for user_dirichlet_update
615 !
616
617 ! Note, some of these are potentially not initialized !
618 call this%user_field_bc_vel%field_list%init(4)
619 call this%user_field_bc_vel%field_list%assign_to_field(1, &
620 this%user_field_bc_vel%bc_u%field_bc)
621 call this%user_field_bc_vel%field_list%assign_to_field(2, &
622 this%user_field_bc_vel%bc_v%field_bc)
623 call this%user_field_bc_vel%field_list%assign_to_field(3, &
624 this%user_field_bc_vel%bc_w%field_bc)
625 call this%user_field_bc_vel%field_list%assign_to_field(4, &
626 this%user_field_bc_prs%field_bc)
627
628 call this%user_field_bc_vel%bc_list%init(size = 4)
629 ! Note, bc_list_add only adds if the bc is not empty
630 call this%user_field_bc_vel%bc_list%append(this%user_field_bc_vel%bc_u)
631 call this%user_field_bc_vel%bc_list%append(this%user_field_bc_vel%bc_v)
632 call this%user_field_bc_vel%bc_list%append(this%user_field_bc_vel%bc_w)
633
634 !
635 ! Check if we need to output boundary types to a separate field
636 call fluid_scheme_set_bc_type_output(this, params)
637
638 !
639 ! Setup right-hand side fields.
640 !
641 allocate(this%f_x)
642 allocate(this%f_y)
643 allocate(this%f_z)
644 call this%f_x%init(this%dm_Xh, fld_name = "fluid_rhs_x")
645 call this%f_y%init(this%dm_Xh, fld_name = "fluid_rhs_y")
646 call this%f_z%init(this%dm_Xh, fld_name = "fluid_rhs_z")
647
648 ! Initialize the source term
649 call this%source_term%init(this%f_x, this%f_y, this%f_z, this%c_Xh, user)
650 call this%source_term%add(params, 'case.fluid.source_terms')
651
652 ! If case.output_boundary is true, set the values for the bc types in the
653 ! output of the field.
654 call this%set_bc_type_output(params)
655
656 ! Initialize velocity solver
657 if (kspv_init) then
658 call neko_log%section("Velocity solver")
659 call json_get_or_default(params, &
660 'case.fluid.velocity_solver.max_iterations', &
661 integer_val, ksp_max_iter)
662 call json_get(params, 'case.fluid.velocity_solver.type', string_val1)
663 call json_get(params, 'case.fluid.velocity_solver.preconditioner', &
664 string_val2)
665 call json_get(params, 'case.fluid.velocity_solver.absolute_tolerance', &
666 real_val)
667 call json_get_or_default(params, &
668 'case.fluid.velocity_solver.monitor', &
669 logical_val, .false.)
670
671 call neko_log%message('Type : ('// trim(string_val1) // &
672 ', ' // trim(string_val2) // ')')
673
674 write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
675 call neko_log%message(log_buf)
676 call fluid_scheme_solver_factory(this%ksp_vel, this%dm_Xh%size(), &
677 string_val1, integer_val, real_val, logical_val)
678 call fluid_scheme_precon_factory(this%pc_vel, this%ksp_vel, &
679 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bclst_vel, string_val2)
680 call neko_log%end_section()
681 end if
682
683 ! Strict convergence for the velocity solver
684 call json_get_or_default(params, 'case.fluid.strict_convergence', &
685 this%strict_convergence, .false.)
686
687 ! Assign velocity fields
688 call neko_field_registry%add_field(this%dm_Xh, 'u')
689 call neko_field_registry%add_field(this%dm_Xh, 'v')
690 call neko_field_registry%add_field(this%dm_Xh, 'w')
691 this%u => neko_field_registry%get_field('u')
692 this%v => neko_field_registry%get_field('v')
693 this%w => neko_field_registry%get_field('w')
694
695 !! Initialize time-lag fields
696 call this%ulag%init(this%u, 2)
697 call this%vlag%init(this%v, 2)
698 call this%wlag%init(this%w, 2)
699
700
701 end subroutine fluid_scheme_init_common
702
704 subroutine fluid_scheme_init_all(this, msh, lx, params, kspv_init, &
705 kspp_init, scheme, user)
706 implicit none
707 class(fluid_scheme_t), target, intent(inout) :: this
708 type(mesh_t), target, intent(inout) :: msh
709 integer, intent(in) :: lx
710 type(json_file), target, intent(inout) :: params
711 type(user_t), target, intent(in) :: user
712 logical :: kspv_init
713 logical :: kspp_init
714 character(len=*), intent(in) :: scheme
715 real(kind=rp) :: abs_tol
716 integer :: integer_val, ierr
717 logical :: logical_val
718 character(len=:), allocatable :: solver_type, precon_type
719 character(len=LOG_SIZE) :: log_buf
720 real(kind=rp) :: gjp_param_a, gjp_param_b
721
722 call fluid_scheme_init_common(this, msh, lx, params, scheme, user, &
723 kspv_init)
724
725 call neko_field_registry%add_field(this%dm_Xh, 'p')
726 this%p => neko_field_registry%get_field('p')
727
728 !
729 ! Setup pressure boundary conditions
730 !
731 call this%bclst_prs%init()
732 call this%bc_prs%init_base(this%c_Xh)
733 call this%bc_prs%mark_zones_from_list(msh%labeled_zones, &
734 'o', this%bc_labels)
735 call this%bc_prs%mark_zones_from_list(msh%labeled_zones, &
736 'on', this%bc_labels)
737
738 ! Field dirichlet pressure bc
739 call this%user_field_bc_prs%init_base(this%c_Xh)
740 call this%user_field_bc_prs%mark_zones_from_list(msh%labeled_zones, &
741 'd_pres', this%bc_labels)
742 call this%user_field_bc_prs%finalize()
743 call mpi_allreduce(this%user_field_bc_prs%msk(0), integer_val, 1, &
744 mpi_integer, mpi_sum, neko_comm, ierr)
745
746 if (integer_val .gt. 0) call this%user_field_bc_prs%init_field('d_pres')
747 call this%bclst_prs%append(this%user_field_bc_prs)
748 call this%user_field_bc_vel%bc_list%append(this%user_field_bc_prs)
749
750 if (msh%outlet%size .gt. 0) then
751 call this%bc_prs%mark_zone(msh%outlet)
752 end if
753 if (msh%outlet_normal%size .gt. 0) then
754 call this%bc_prs%mark_zone(msh%outlet_normal)
755 end if
756
757 call this%bc_prs%finalize()
758 call this%bc_prs%set_g(0.0_rp)
759 call this%bclst_prs%append(this%bc_prs)
760 call this%bc_dong%init_base(this%c_Xh)
761 call this%bc_dong%mark_zones_from_list(msh%labeled_zones, &
762 'o+dong', this%bc_labels)
763 call this%bc_dong%mark_zones_from_list(msh%labeled_zones, &
764 'on+dong', this%bc_labels)
765 call this%bc_dong%finalize()
766
767
768 call this%bc_dong%init(this%c_Xh, params)
769
770 call this%bclst_prs%append(this%bc_dong)
771
772 ! Pressure solver
773 if (kspp_init) then
774 call neko_log%section("Pressure solver")
775
776 call json_get_or_default(params, &
777 'case.fluid.pressure_solver.max_iterations', &
778 integer_val, ksp_max_iter)
779 call json_get(params, 'case.fluid.pressure_solver.type', solver_type)
780 call json_get(params, 'case.fluid.pressure_solver.preconditioner', &
781 precon_type)
782 call json_get(params, 'case.fluid.pressure_solver.absolute_tolerance', &
783 abs_tol)
784 call json_get_or_default(params, &
785 'case.fluid.pressure_solver.monitor', &
786 logical_val, .false.)
787 call neko_log%message('Type : ('// trim(solver_type) // &
788 ', ' // trim(precon_type) // ')')
789 write(log_buf, '(A,ES13.6)') 'Abs tol :', abs_tol
790 call neko_log%message(log_buf)
791
792 call fluid_scheme_solver_factory(this%ksp_prs, this%dm_Xh%size(), &
793 solver_type, integer_val, abs_tol, logical_val)
794 call fluid_scheme_precon_factory(this%pc_prs, this%ksp_prs, &
795 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bclst_prs, precon_type)
796
797 call neko_log%end_section()
798
799 end if
800
801 ! Initiate gradient jump penalty
802 call json_get_or_default(params, &
803 'case.fluid.gradient_jump_penalty.enabled',&
804 this%if_gradient_jump_penalty, .false.)
805
806 if (this%if_gradient_jump_penalty .eqv. .true.) then
807 if ((this%dm_Xh%xh%lx - 1) .eq. 1) then
808 call json_get_or_default(params, &
809 'case.fluid.gradient_jump_penalty.tau',&
810 gjp_param_a, 0.02_rp)
811 gjp_param_b = 0.0_rp
812 else
813 call json_get_or_default(params, &
814 'case.fluid.gradient_jump_penalty.scaling_factor',&
815 gjp_param_a, 0.8_rp)
816 call json_get_or_default(params, &
817 'case.fluid.gradient_jump_penalty.scaling_exponent',&
818 gjp_param_b, 4.0_rp)
819 end if
820 call this%gradient_jump_penalty_u%init(params, this%dm_Xh, this%c_Xh, &
821 gjp_param_a, gjp_param_b)
822 call this%gradient_jump_penalty_v%init(params, this%dm_Xh, this%c_Xh, &
823 gjp_param_a, gjp_param_b)
824 call this%gradient_jump_penalty_w%init(params, this%dm_Xh, this%c_Xh, &
825 gjp_param_a, gjp_param_b)
826 end if
827
828 call neko_log%end_section()
829
830 end subroutine fluid_scheme_init_all
831
833 subroutine fluid_scheme_free(this)
834 class(fluid_scheme_t), intent(inout) :: this
835
836 call this%bdry%free()
837
838 if (allocated(this%bc_inflow)) then
839 call this%bc_inflow%free()
840 end if
841
842 call this%bc_wall%free()
843 call this%bc_sym%free()
844 call this%bc_sh%free()
845
846 !
847 ! Free everything related to field_dirichlet BCs
848 !
849 call this%user_field_bc_prs%field_bc%free()
850 call this%user_field_bc_prs%free()
851 call this%user_field_bc_vel%bc_u%field_bc%free()
852 call this%user_field_bc_vel%bc_v%field_bc%free()
853 call this%user_field_bc_vel%bc_w%field_bc%free()
854 call this%user_field_bc_vel%free()
855
856 call this%Xh%free()
857
858 if (allocated(this%ksp_vel)) then
859 call krylov_solver_destroy(this%ksp_vel)
860 deallocate(this%ksp_vel)
861 end if
862
863 if (allocated(this%ksp_prs)) then
864 call krylov_solver_destroy(this%ksp_prs)
865 deallocate(this%ksp_prs)
866 end if
867
868 if (allocated(this%pc_vel)) then
869 call precon_destroy(this%pc_vel)
870 deallocate(this%pc_vel)
871 end if
872
873 if (allocated(this%pc_prs)) then
874 call precon_destroy(this%pc_prs)
875 deallocate(this%pc_prs)
876 end if
877
878 if (allocated(this%bc_labels)) then
879 deallocate(this%bc_labels)
880 end if
881
882 call this%source_term%free()
883
884 call this%gs_Xh%free()
885
886 call this%c_Xh%free()
887
888 call this%bclst_vel%free()
889
890 call this%scratch%free()
891
892 nullify(this%params)
893
894 nullify(this%u)
895 nullify(this%v)
896 nullify(this%w)
897 nullify(this%p)
898
899 call this%ulag%free()
900 call this%vlag%free()
901 call this%wlag%free()
902
903
904 if (associated(this%f_x)) then
905 call this%f_x%free()
906 end if
907
908 if (associated(this%f_y)) then
909 call this%f_y%free()
910 end if
911
912 if (associated(this%f_z)) then
913 call this%f_z%free()
914 end if
915
916 nullify(this%f_x)
917 nullify(this%f_y)
918 nullify(this%f_z)
919
920 call this%mu_field%free()
921
922 ! Free gradient jump penalty
923 if (this%if_gradient_jump_penalty .eqv. .true.) then
924 call this%gradient_jump_penalty_u%free()
925 call this%gradient_jump_penalty_v%free()
926 call this%gradient_jump_penalty_w%free()
927 end if
928
929
930 end subroutine fluid_scheme_free
931
934 subroutine fluid_scheme_validate(this)
935 class(fluid_scheme_t), target, intent(inout) :: this
936 ! Variables for retrieving json parameters
937 logical :: logical_val
938
939 if ( (.not. associated(this%u)) .or. &
940 (.not. associated(this%v)) .or. &
941 (.not. associated(this%w)) .or. &
942 (.not. associated(this%p))) then
943 call neko_error('Fields are not registered')
944 end if
945
946 if ( (.not. allocated(this%u%x)) .or. &
947 (.not. allocated(this%v%x)) .or. &
948 (.not. allocated(this%w%x)) .or. &
949 (.not. allocated(this%p%x))) then
950 call neko_error('Fields are not allocated')
951 end if
952
953 if (.not. allocated(this%ksp_vel)) then
954 call neko_error('No Krylov solver for velocity defined')
955 end if
956
957 if (.not. allocated(this%ksp_prs)) then
958 call neko_error('No Krylov solver for pressure defined')
959 end if
960
961 if (.not. associated(this%params)) then
962 call neko_error('No parameters defined')
963 end if
964
965 if (allocated(this%bc_inflow)) then
966 select type (ip => this%bc_inflow)
967 type is (usr_inflow_t)
968 call ip%validate
969 end select
970 end if
971
972 !
973 ! Setup checkpoint structure (if everything is fine)
974 !
975 call this%chkp%init(this%u, this%v, this%w, this%p)
976
977 end subroutine fluid_scheme_validate
978
983 subroutine fluid_scheme_bc_apply_vel(this, t, tstep)
984 class(fluid_scheme_t), intent(inout) :: this
985 real(kind=rp), intent(in) :: t
986 integer, intent(in) :: tstep
987
988 call this%bclst_vel%apply_vector( &
989 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), t, tstep)
990
991 end subroutine fluid_scheme_bc_apply_vel
992
995 subroutine fluid_scheme_bc_apply_prs(this, t, tstep)
996 class(fluid_scheme_t), intent(inout) :: this
997 real(kind=rp), intent(in) :: t
998 integer, intent(in) :: tstep
999
1000 call this%bclst_prs%apply_scalar(this%p%x, this%p%dof%size(), t, tstep)
1001
1002 end subroutine fluid_scheme_bc_apply_prs
1003
1006 subroutine fluid_scheme_solver_factory(ksp, n, solver, &
1007 max_iter, abstol, monitor)
1008 class(ksp_t), allocatable, target, intent(inout) :: ksp
1009 integer, intent(in), value :: n
1010 character(len=*), intent(in) :: solver
1011 integer, intent(in) :: max_iter
1012 real(kind=rp), intent(in) :: abstol
1013 logical, intent(in) :: monitor
1014
1015 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
1016 monitor = monitor)
1017
1018 end subroutine fluid_scheme_solver_factory
1019
1021 subroutine fluid_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
1022 pctype)
1023 class(pc_t), allocatable, target, intent(inout) :: pc
1024 class(ksp_t), target, intent(inout) :: ksp
1025 type(coef_t), target, intent(in) :: coef
1026 type(dofmap_t), target, intent(in) :: dof
1027 type(gs_t), target, intent(inout) :: gs
1028 type(bc_list_t), target, intent(inout) :: bclst
1029 character(len=*) :: pctype
1030
1031 call precon_factory(pc, pctype)
1032
1033 select type (pcp => pc)
1034 type is (jacobi_t)
1035 call pcp%init(coef, dof, gs)
1036 type is (sx_jacobi_t)
1037 call pcp%init(coef, dof, gs)
1038 type is (device_jacobi_t)
1039 call pcp%init(coef, dof, gs)
1040 type is (hsmg_t)
1041 if (len_trim(pctype) .gt. 4) then
1042 if (index(pctype, '+') .eq. 5) then
1043 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst, &
1044 trim(pctype(6:)))
1045 else
1046 call neko_error('Unknown coarse grid solver')
1047 end if
1048 else
1049 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
1050 end if
1051 type is (phmg_t)
1052 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
1053 end select
1054
1055 call ksp%set_pc(pc)
1056
1057 end subroutine fluid_scheme_precon_factory
1058
1060 subroutine fluid_scheme_set_usr_inflow(this, usr_eval)
1061 class(fluid_scheme_t), intent(inout) :: this
1062 procedure(usr_inflow_eval) :: usr_eval
1063
1064 select type (bc_if => this%bc_inflow)
1065 type is (usr_inflow_t)
1066 call bc_if%set_eval(usr_eval)
1067 class default
1068 call neko_error("Not a user defined inflow condition")
1069 end select
1070 end subroutine fluid_scheme_set_usr_inflow
1071
1073 function fluid_compute_cfl(this, dt) result(c)
1074 class(fluid_scheme_t), intent(in) :: this
1075 real(kind=rp), intent(in) :: dt
1076 real(kind=rp) :: c
1077
1078 c = cfl(dt, this%u%x, this%v%x, this%w%x, &
1079 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
1080
1081 end function fluid_compute_cfl
1082
1085 subroutine fluid_scheme_set_bc_type_output(this, params)
1086 class(fluid_scheme_t), target, intent(inout) :: this
1087 type(json_file), intent(inout) :: params
1088 type(dirichlet_t) :: bdry_mask
1089 logical :: found
1090
1091 !
1092 ! Check if we need to output boundaries
1093 !
1094 call json_get_or_default(params, 'case.output_boundary', found, .false.)
1095
1096 if (found) then
1097 call this%bdry%init(this%dm_Xh, 'bdry')
1098 this%bdry = 0.0_rp
1099
1100 call bdry_mask%init_base(this%c_Xh)
1101 call bdry_mask%mark_zone(this%msh%wall)
1102 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1103 'w', this%bc_labels)
1104 call bdry_mask%finalize()
1105 call bdry_mask%set_g(1.0_rp)
1106 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1107 call bdry_mask%free()
1108
1109 call bdry_mask%init_base(this%c_Xh)
1110 call bdry_mask%mark_zone(this%msh%inlet)
1111 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1112 'v', this%bc_labels)
1113 call bdry_mask%finalize()
1114 call bdry_mask%set_g(2.0_rp)
1115 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1116 call bdry_mask%free()
1117
1118 call bdry_mask%init_base(this%c_Xh)
1119 call bdry_mask%mark_zone(this%msh%outlet)
1120 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1121 'o', this%bc_labels)
1122 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1123 'o+dong', this%bc_labels)
1124 call bdry_mask%finalize()
1125 call bdry_mask%set_g(3.0_rp)
1126 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1127 call bdry_mask%free()
1128
1129 call bdry_mask%init_base(this%c_Xh)
1130 call bdry_mask%mark_zone(this%msh%sympln)
1131 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1132 'sym', this%bc_labels)
1133 call bdry_mask%finalize()
1134 call bdry_mask%set_g(4.0_rp)
1135 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1136 call bdry_mask%free()
1137
1138 call bdry_mask%init_base(this%c_Xh)
1139 call bdry_mask%mark_zone(this%msh%outlet_normal)
1140 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1141 'on', this%bc_labels)
1142 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1143 'on+dong', this%bc_labels)
1144 call bdry_mask%finalize()
1145 call bdry_mask%set_g(5.0_rp)
1146 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1147 call bdry_mask%free()
1148
1149 call bdry_mask%init_base(this%c_Xh)
1150 call bdry_mask%mark_zone(this%msh%periodic)
1151 call bdry_mask%finalize()
1152 call bdry_mask%set_g(6.0_rp)
1153 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1154 call bdry_mask%free()
1155
1156 call bdry_mask%init_base(this%c_Xh)
1157 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1158 'd_vel_u', this%bc_labels)
1159 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1160 'd_vel_v', this%bc_labels)
1161 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1162 'd_vel_w', this%bc_labels)
1163 call bdry_mask%finalize()
1164 call bdry_mask%set_g(7.0_rp)
1165 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1166 call bdry_mask%free()
1167
1168 call bdry_mask%init_base(this%c_Xh)
1169 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1170 'd_pres', this%bc_labels)
1171 call bdry_mask%finalize()
1172 call bdry_mask%set_g(8.0_rp)
1173 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1174 call bdry_mask%free()
1175
1176 call bdry_mask%init_base(this%c_Xh)
1177 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1178 'sh', this%bc_labels)
1179 call bdry_mask%finalize()
1180 call bdry_mask%set_g(9.0_rp)
1181 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1182 call bdry_mask%free()
1183
1184 call bdry_mask%init_base(this%c_Xh)
1185 call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1186 'wm', this%bc_labels)
1187 call bdry_mask%finalize()
1188 call bdry_mask%set_g(10.0_rp)
1189 call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1190 call bdry_mask%free()
1191
1192 end if
1193
1194 end subroutine fluid_scheme_set_bc_type_output
1195
1197 subroutine fluid_scheme_update_material_properties(this)
1198 class(fluid_scheme_t), intent(inout) :: this
1199 type(field_t), pointer :: nut
1200 integer :: n
1201
1202 if (this%variable_material_properties) then
1203 nut => neko_field_registry%get_field(this%nut_field_name)
1204 n = nut%size()
1205
1206 if (neko_bcknd_device .eq. 1) then
1207 call device_cfill(this%mu_field%x_d, this%mu, n)
1208 call device_add2s2(this%mu_field%x_d, nut%x_d, this%rho, n)
1209 else
1210 call cfill(this%mu_field%x, this%mu, n)
1211 call add2s2(this%mu_field%x, nut%x, this%rho, n)
1212 end if
1213 end if
1214
1215 end subroutine fluid_scheme_update_material_properties
1216
1220 subroutine fluid_scheme_set_material_properties(this, params, user)
1221 class(fluid_scheme_t), intent(inout) :: this
1222 type(json_file), intent(inout) :: params
1223 type(user_t), target, intent(in) :: user
1224 character(len=LOG_SIZE) :: log_buf
1225 ! A local pointer that is needed to make Intel happy
1226 procedure(user_material_properties), pointer :: dummy_mp_ptr
1227 logical :: nondimensional
1228 real(kind=rp) :: dummy_lambda, dummy_cp
1229
1230 dummy_mp_ptr => dummy_user_material_properties
1231
1232 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
1233
1234 write(log_buf, '(A)') "Material properties must be set in the user&
1235 & file!"
1236 call neko_log%message(log_buf)
1237 call user%material_properties(0.0_rp, 0, this%rho, this%mu, &
1238 dummy_cp, dummy_lambda, params)
1239 else
1240 ! Incorrect user input
1241 if (params%valid_path('case.fluid.Re') .and. &
1242 (params%valid_path('case.fluid.mu') .or. &
1243 params%valid_path('case.fluid.rho'))) then
1244 call neko_error("To set the material properties for the fluid,&
1245 & either provide Re OR mu and rho in the case file.")
1246
1247 ! Non-dimensional case
1248 else if (params%valid_path('case.fluid.Re')) then
1249
1250 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
1251 & input.'
1252 call neko_log%message(log_buf, lvl = neko_log_verbose)
1253 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
1254 & 1/Re.'
1255 call neko_log%message(log_buf, lvl = neko_log_verbose)
1256
1257 ! Read Re into mu for further manipulation.
1258 call json_get(params, 'case.fluid.Re', this%mu)
1259 write(log_buf, '(A)') 'Read non-dimensional material properties'
1260 call neko_log%message(log_buf)
1261 write(log_buf, '(A,ES13.6)') 'Re :', this%mu
1262 call neko_log%message(log_buf)
1263
1264 ! Set rho to 1 since the setup is non-dimensional.
1265 this%rho = 1.0_rp
1266 ! Invert the Re to get viscosity.
1267 this%mu = 1.0_rp/this%mu
1268 ! Dimensional case
1269 else
1270 call json_get(params, 'case.fluid.mu', this%mu)
1271 call json_get(params, 'case.fluid.rho', this%rho)
1272 end if
1273
1274 end if
1275 end subroutine fluid_scheme_set_material_properties
1276
1277end module fluid_scheme
Abstract interface to dealocate a fluid formulation.
Abstract interface to initialize a fluid formulation.
Abstract interface to restart a fluid scheme.
Abstract interface to compute a time-step.
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 Blasius profile dirichlet condition.
Definition blasius.f90:34
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:16
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 .
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a dong outflow condition.
Defines inflow dirichlet conditions.
Defines inflow dirichlet conditions.
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
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
Fluid formulations.
subroutine fluid_scheme_set_material_properties(this, params, user)
Sets rho and mu.
subroutine fluid_scheme_free(this)
Deallocate a fluid formulation.
subroutine fluid_scheme_update_material_properties(this)
Update the values of mu_field if necessary.
subroutine fluid_scheme_set_bc_type_output(this, params)
Set boundary types for the diagnostic output.
real(kind=rp) function fluid_compute_cfl(this, dt)
Compute CFL.
subroutine fluid_scheme_init_all(this, msh, lx, params, kspv_init, kspp_init, scheme, user)
Initialize all components of the current scheme.
subroutine fluid_scheme_set_usr_inflow(this, usr_eval)
Initialize a user defined inflow condition.
subroutine fluid_scheme_init_common(this, msh, lx, params, scheme, user, kspv_init)
Initialise a fluid scheme.
subroutine fluid_scheme_bc_apply_vel(this, t, tstep)
Apply all boundary conditions defined for velocity Here we perform additional gs operations to take c...
subroutine fluid_scheme_validate(this)
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
subroutine fluid_scheme_bc_apply_prs(this, t, tstep)
Apply all boundary conditions defined for pressure.
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
Defines inflow dirichlet conditions.
Definition inflow.f90:34
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...
type(scratch_registry_t) function init(dof, size, expansion_size)
Constructor, optionally taking initial registry and expansion size as argument.
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.
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition symmetry.f90:34
Compound scheme for the advection and diffusion operators in a transport equation.
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 the wall_model_bc_t type. Maintainer: Timofey Mukha.
Defines wall boundary conditions.
Definition wall.f90:34
Base type for a boundary condition.
Definition bc.f90:51
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:46
Blasius profile for inlet (vector valued)
Definition blasius.f90:48
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:44
Dong outflow condition Follows "A Convective-like Energy-Stable Open Boundary Condition for Simulati...
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.
Dirichlet condition for inlet (vector valued)
Definition inflow.f90:43
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.
Mixed Dirichlet-Neumann symmetry plane condition.
Definition symmetry.f90:46
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
User defined dirichlet condition for inlet (vector valued)
No-slip Wall boundary condition.
Definition wall.f90:43
A shear stress boundary condition, computing the stress values using a wall model.