Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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
38 use checkpoint, only : chkp_t
39 use num_types, only : rp, i8
41 use field, only : field_t
42 use space, only : gll
43 use dofmap, only : dofmap_t
44 use krylov, only : ksp_t, krylov_solver_factory, ksp_max_iter
45 use coefs, only : coef_t
46 use dirichlet, only : dirichlet_t
47 use jacobi, only : jacobi_t
48 use sx_jacobi, only : sx_jacobi_t
50 use hsmg, only : hsmg_t
51 use phmg, only : phmg_t
52 use precon, only : pc_t, precon_factory, precon_destroy
53 use fluid_stats, only : fluid_stats_t
54 use bc, only : bc_t
55 use bc_list, only : bc_list_t
56 use mesh, only : mesh_t
57 use math, only : glsum
58 use operators, only : cfl, rotate_cyc
60 use registry, only : neko_registry
63 use json_module, only : json_file
67 use utils, only : neko_error
71 use time_state, only : time_state_t
72 implicit none
73 private
74
79 class(ksp_t), allocatable :: ksp_vel
80 class(ksp_t), allocatable :: ksp_prs
81 class(pc_t), allocatable :: pc_vel
82 class(pc_t), allocatable :: pc_prs
83 integer :: vel_projection_dim
84 integer :: pr_projection_dim
85 integer :: vel_projection_activ_step
86 integer :: pr_projection_activ_step
87 logical :: pr_projection_reorthogonalize_basis
88 logical :: strict_convergence
89 logical :: allow_stabilization
91 type(field_t), pointer :: u_e => null()
92 type(field_t), pointer :: v_e => null()
93 type(field_t), pointer :: w_e => null()
94
96 logical :: forced_flow_rate = .false.
97
99 character(len=:), allocatable :: nut_field_name
100
101 ! The total viscosity field
102 type(field_t), pointer :: mu_tot => null()
103
105 integer(kind=i8) :: glb_n_points
107 integer(kind=i8) :: glb_unique_points
108 contains
110 procedure, pass(this) :: init_base => fluid_scheme_init_base
111 procedure, pass(this) :: scheme_free => fluid_scheme_free
113 procedure, pass(this) :: validate => fluid_scheme_validate
115 procedure, pass(this) :: bc_apply_vel => fluid_scheme_bc_apply_vel
117 procedure, pass(this) :: bc_apply_prs => fluid_scheme_bc_apply_prs
119 procedure, pass(this) :: compute_cfl => fluid_compute_cfl
121 procedure, pass(this) :: set_material_properties => &
123
125 procedure, pass(this) :: update_material_properties => &
128 procedure, nopass :: solver_factory => fluid_scheme_solver_factory
130 procedure, pass(this) :: precon_factory_ => fluid_scheme_precon_factory
132
133 interface
134
135 module subroutine fluid_scheme_factory(object, type_name)
136 class(fluid_scheme_base_t), intent(inout), allocatable :: object
137 character(len=*) :: type_name
138 end subroutine fluid_scheme_factory
139 end interface
140
141 public :: fluid_scheme_incompressible_t, fluid_scheme_factory
142
143contains
144
146 subroutine fluid_scheme_init_base(this, msh, lx, params, scheme, user, &
147 kspv_init)
148 implicit none
149 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
150 type(mesh_t), target, intent(inout) :: msh
151 integer, intent(in) :: lx
152 character(len=*), intent(in) :: scheme
153 type(json_file), target, intent(inout) :: params
154 type(user_t), target, intent(in) :: user
155 logical, intent(in) :: kspv_init
156 type(dirichlet_t) :: bdry_mask
157 character(len=LOG_SIZE) :: log_buf
158 real(kind=rp), allocatable :: real_vec(:)
159 real(kind=rp) :: real_val, kappa, b, z0
160 logical :: logical_val
161 integer :: integer_val, ierr
162 type(json_file) :: wm_json
163 character(len=:), allocatable :: string_val1, string_val2
164 type(json_file) :: json_subdict
165
166 !
167 ! SEM simulation fundamentals
168 !
169
170 this%msh => msh
171
172 if (msh%gdim .eq. 2) then
173 call this%Xh%init(gll, lx, lx)
174 else
175 call this%Xh%init(gll, lx, lx, lx)
176 end if
177
178 call this%dm_Xh%init(msh, this%Xh)
179
180 call this%gs_Xh%init(this%dm_Xh)
181
182 call this%c_Xh%init(this%gs_Xh)
183
184 ! Assign Dofmap to scratch registry
185 call neko_scratch_registry%set_dofmap(this%dm_Xh)
186
187 ! Assign a name
188 call json_get_or_default(params, 'case.fluid.name', this%name, "fluid")
189
190 !
191 ! First section of fluid log
192 !
193
194 call neko_log%section('Fluid')
195 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
196 call neko_log%message(log_buf)
197 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
198 call neko_log%message(log_buf)
199
200 ! Assign velocity fields
201 call neko_registry%add_field(this%dm_Xh, 'u')
202 call neko_registry%add_field(this%dm_Xh, 'v')
203 call neko_registry%add_field(this%dm_Xh, 'w')
204 this%u => neko_registry%get_field('u')
205 this%v => neko_registry%get_field('v')
206 this%w => neko_registry%get_field('w')
207
208 !
209 ! Material properties
210 !
211 call this%set_material_properties(params, user)
212
213 ! Projection spaces
214 call json_get_or_lookup_or_default(params, &
215 'case.fluid.velocity_solver.projection_space_size', &
216 this%vel_projection_dim, 0)
217 call json_get_or_lookup_or_default(params, &
218 'case.fluid.pressure_solver.projection_space_size', &
219 this%pr_projection_dim, 0)
220 call json_get_or_lookup_or_default(params, &
221 'case.fluid.velocity_solver.projection_hold_steps', &
222 this%vel_projection_activ_step, 5)
223 call json_get_or_lookup_or_default(params, &
224 'case.fluid.pressure_solver.projection_hold_steps', &
225 this%pr_projection_activ_step, 5)
226 call json_get_or_default(params, &
227 'case.fluid.pressure_solver.projection_reorthogonalize_basis', &
228 this%pr_projection_reorthogonalize_basis, .false.)
229
230 call json_get_or_default(params, 'case.fluid.freeze', this%freeze, .false.)
231
232 if (params%valid_path("case.fluid.flow_rate_force")) then
233 this%forced_flow_rate = .true.
234 end if
235
236
237 if (lx .lt. 10) then
238 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
239 else if (lx .ge. 10) then
240 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
241 else
242 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
243 end if
244 call neko_log%message(log_buf)
245 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
246 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
247
248 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
249 call neko_log%message(log_buf)
250 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
251 call neko_log%message(log_buf)
252
253
254 call json_get(params, 'case.numerics.dealias', logical_val)
255 write(log_buf, '(A, L1)') 'Dealias : ', logical_val
256 call neko_log%message(log_buf)
257
258
259 call json_get_or_default(params, 'case.output_boundary', logical_val, &
260 .false.)
261 write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
262 call neko_log%message(log_buf)
263
264 call json_get_or_default(params, "case.fluid.full_stress_formulation", &
265 logical_val, .false.)
266 write(log_buf, '(A, L1)') 'Full stress: ', logical_val
267 call neko_log%message(log_buf)
268
269
270 !
271 ! Setup right-hand side fields.
272 !
273 allocate(this%f_x)
274 allocate(this%f_y)
275 allocate(this%f_z)
276 call this%f_x%init(this%dm_Xh, fld_name = "fluid_rhs_x")
277 call this%f_y%init(this%dm_Xh, fld_name = "fluid_rhs_y")
278 call this%f_z%init(this%dm_Xh, fld_name = "fluid_rhs_z")
279
280 ! Initialize velocity solver
281 if (kspv_init) then
282 call neko_log%section("Velocity solver")
283 call json_get_or_lookup_or_default(params, &
284 'case.fluid.velocity_solver.max_iterations', &
285 integer_val, ksp_max_iter)
286 call json_get(params, 'case.fluid.velocity_solver.type', string_val1)
287 call json_get(params, 'case.fluid.velocity_solver.preconditioner.type', &
288 string_val2)
289 call json_get(params, &
290 'case.fluid.velocity_solver.preconditioner', json_subdict)
291 call json_get_or_lookup(params, &
292 'case.fluid.velocity_solver.absolute_tolerance', &
293 real_val)
294 call json_get_or_default(params, &
295 'case.fluid.velocity_solver.monitor', &
296 logical_val, .false.)
297
298 call neko_log%message('Type : ('// trim(string_val1) // &
299 ', ' // trim(string_val2) // ')')
300
301 write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
302 call neko_log%message(log_buf)
303 call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
304 string_val1, integer_val, real_val, logical_val)
305 call this%precon_factory_(this%pc_vel, this%ksp_vel, &
306 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, &
307 string_val2, json_subdict)
308 call neko_log%end_section()
309 end if
310
311 ! Strict convergence for the velocity solver
312 call json_get_or_default(params, 'case.fluid.strict_convergence', &
313 this%strict_convergence, .false.)
314 ! Allow stabilization period where we do not warn about non-convergence
315 call json_get_or_default(params, 'case.fluid.allow_stabilization', &
316 this%allow_stabilization, .false.)
317
318
319 !! Initialize time-lag fields
320 call this%ulag%init(this%u, 2)
321 call this%vlag%init(this%v, 2)
322 call this%wlag%init(this%w, 2)
323
324 call neko_registry%add_field(this%dm_Xh, 'u_e')
325 call neko_registry%add_field(this%dm_Xh, 'v_e')
326 call neko_registry%add_field(this%dm_Xh, 'w_e')
327 this%u_e => neko_registry%get_field('u_e')
328 this%v_e => neko_registry%get_field('v_e')
329 this%w_e => neko_registry%get_field('w_e')
330
331 ! Initialize the source term
332 call neko_log%section('Fluid Source term')
333 call this%source_term%init(this%f_x, this%f_y, this%f_z, this%c_Xh, user, &
334 this%name)
335 call this%source_term%add(params, 'case.fluid.source_terms')
336 call neko_log%end_section()
337
338 end subroutine fluid_scheme_init_base
339
340 subroutine fluid_scheme_free(this)
341 class(fluid_scheme_incompressible_t), intent(inout) :: this
342 class(bc_t), pointer :: bc
343 integer :: i
344
345
346 call this%Xh%free()
347
348 if (allocated(this%ksp_vel)) then
349 call this%ksp_vel%free()
350 deallocate(this%ksp_vel)
351 end if
352
353 if (allocated(this%ksp_prs)) then
354 call this%ksp_prs%free()
355 deallocate(this%ksp_prs)
356 end if
357
358 if (allocated(this%pc_vel)) then
359 call precon_destroy(this%pc_vel)
360 deallocate(this%pc_vel)
361 end if
362
363 if (allocated(this%pc_prs)) then
364 call precon_destroy(this%pc_prs)
365 deallocate(this%pc_prs)
366 end if
367
368 do i = 1, this%bcs_vel%size()
369 bc => this%bcs_vel%get(i)
370 call bc%free()
371 end do
372 call this%bcs_vel%free()
373
374 do i = 1, this%bcs_prs%size()
375 bc => this%bcs_prs%get(i)
376 call bc%free()
377 end do
378 call this%bcs_prs%free()
379
380 call this%source_term%free()
381
382 call this%gs_Xh%free()
383
384 call this%c_Xh%free()
385
386 nullify(this%u)
387 nullify(this%v)
388 nullify(this%w)
389 nullify(this%p)
390
391 nullify(this%u_e)
392 nullify(this%v_e)
393 nullify(this%w_e)
394
395 call this%ulag%free()
396 call this%vlag%free()
397 call this%wlag%free()
398
399
400 if (associated(this%f_x)) then
401 call this%f_x%free()
402 deallocate(this%f_x)
403 end if
404
405 if (associated(this%f_y)) then
406 call this%f_y%free()
407 deallocate(this%f_y)
408 end if
409
410 if (associated(this%f_z)) then
411 call this%f_z%free()
412 deallocate(this%f_z)
413 end if
414
415 nullify(this%f_x)
416 nullify(this%f_y)
417 nullify(this%f_z)
418 nullify(this%rho)
419 nullify(this%mu)
420 nullify(this%mu_tot)
421
422 call this%dm_Xh%free()
423 call this%Xh%free()
424 nullify(this%msh)
425
426 end subroutine fluid_scheme_free
427
430 subroutine fluid_scheme_validate(this)
431 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
432 ! Variables for retrieving json parameters
433 logical :: logical_val
434
435 if ( (.not. associated(this%u)) .or. &
436 (.not. associated(this%v)) .or. &
437 (.not. associated(this%w)) .or. &
438 (.not. associated(this%p))) then
439 call neko_error('Fields are not registered')
440 end if
441
442 if ( (.not. allocated(this%u%x)) .or. &
443 (.not. allocated(this%v%x)) .or. &
444 (.not. allocated(this%w%x)) .or. &
445 (.not. allocated(this%p%x))) then
446 call neko_error('Fields are not allocated')
447 end if
448
449 if (.not. allocated(this%ksp_vel)) then
450 call neko_error('No Krylov solver for velocity defined')
451 end if
452
453 if (.not. allocated(this%ksp_prs)) then
454 call neko_error('No Krylov solver for pressure defined')
455 end if
456
457 end subroutine fluid_scheme_validate
458
463 subroutine fluid_scheme_bc_apply_vel(this, time, strong)
464 class(fluid_scheme_incompressible_t), intent(inout) :: this
465 type(time_state_t), intent(in) :: time
466 logical, intent(in) :: strong
467 integer :: i
468 class(bc_t), pointer :: b
469 b => null()
470
471 call this%bcs_vel%apply_vector(&
472 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), time, strong)
473
474 call rotate_cyc(this%u%x, this%v%x, this%w%x, 1, this%c_Xh)
475 call this%gs_Xh%op(this%u, gs_op_min, glb_cmd_event)
477 call this%gs_Xh%op(this%v, gs_op_min, glb_cmd_event)
479 call this%gs_Xh%op(this%w, gs_op_min, glb_cmd_event)
481 call rotate_cyc(this%u%x, this%v%x, this%w%x, 0, this%c_Xh)
482
483
484 call this%bcs_vel%apply_vector(&
485 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), time, strong)
486
487 call rotate_cyc(this%u%x, this%v%x, this%w%x, 1, this%c_Xh)
488 call this%gs_Xh%op(this%u, gs_op_max, glb_cmd_event)
490 call this%gs_Xh%op(this%v, gs_op_max, glb_cmd_event)
492 call this%gs_Xh%op(this%w, gs_op_max, glb_cmd_event)
494 call rotate_cyc(this%u%x, this%v%x, this%w%x, 0, this%c_Xh)
495
496 do i = 1, this%bcs_vel%size()
497 b => this%bcs_vel%get(i)
498 b%updated = .false.
499 end do
500 nullify(b)
501
502 end subroutine fluid_scheme_bc_apply_vel
503
506 subroutine fluid_scheme_bc_apply_prs(this, time)
507 class(fluid_scheme_incompressible_t), intent(inout) :: this
508 type(time_state_t), intent(in) :: time
509
510 integer :: i
511 class(bc_t), pointer :: b
512 b => null()
513
514 call this%bcs_prs%apply(this%p, time)
515 call this%gs_Xh%op(this%p, gs_op_min, glb_cmd_event)
517
518 call this%bcs_prs%apply(this%p, time)
519 call this%gs_Xh%op(this%p, gs_op_max, glb_cmd_event)
521
522 do i = 1, this%bcs_prs%size()
523 b => this%bcs_prs%get(i)
524 b%updated = .false.
525 end do
526 nullify(b)
527
528 end subroutine fluid_scheme_bc_apply_prs
529
532 subroutine fluid_scheme_solver_factory(ksp, n, solver, &
533 max_iter, abstol, monitor)
534 class(ksp_t), allocatable, target, intent(inout) :: ksp
535 integer, intent(in), value :: n
536 character(len=*), intent(in) :: solver
537 integer, intent(in) :: max_iter
538 real(kind=rp), intent(in) :: abstol
539 logical, intent(in) :: monitor
540
541 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
542 monitor = monitor)
543
544 end subroutine fluid_scheme_solver_factory
545
547 subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, &
548 pctype, pcparams)
549 class(fluid_scheme_incompressible_t), intent(inout) :: this
550 class(pc_t), allocatable, target, intent(inout) :: pc
551 class(ksp_t), target, intent(inout) :: ksp
552 type(coef_t), target, intent(in) :: coef
553 type(dofmap_t), target, intent(in) :: dof
554 type(gs_t), target, intent(inout) :: gs
555 type(bc_list_t), target, intent(inout) :: bclst
556 character(len=*) :: pctype
557 type(json_file), intent(inout) :: pcparams
558
559 call precon_factory(pc, pctype)
560
561 select type (pcp => pc)
562 type is (jacobi_t)
563 call pcp%init(coef, dof, gs)
564 type is (sx_jacobi_t)
565 call pcp%init(coef, dof, gs)
566 type is (device_jacobi_t)
567 call pcp%init(coef, dof, gs)
568 type is (hsmg_t)
569 call pcp%init(coef, bclst, pcparams)
570 type is (phmg_t)
571 call pcp%init(coef, bclst, pcparams)
572 end select
573
574 call ksp%set_pc(pc)
575
576 end subroutine fluid_scheme_precon_factory
577
579 function fluid_compute_cfl(this, dt) result(c)
580 class(fluid_scheme_incompressible_t), intent(in) :: this
581 real(kind=rp), intent(in) :: dt
582 real(kind=rp) :: c
583
584 c = cfl(dt, this%u%x, this%v%x, this%w%x, &
585 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
586
587 end function fluid_compute_cfl
588
589
594 subroutine fluid_scheme_update_material_properties(this, time)
595 class(fluid_scheme_incompressible_t), intent(inout) :: this
596 type(time_state_t), intent(in) :: time
597 type(field_t), pointer :: nut
598
599 call this%user_material_properties(this%name, this%material_properties, &
600 time)
601
602 if (len(trim(this%nut_field_name)) > 0) then
603 nut => neko_registry%get_field(this%nut_field_name)
604 ! Copy material property
605 call field_copy(this%mu_tot, this%mu)
606 ! Add turbulent contribution
607 call field_addcol3(this%mu_tot, nut, this%rho)
608 end if
609
610 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
611 ! host array data to pass constant density and viscosity
612 ! to some routines, we need to make sure that the host
613 ! values are also filled
614 if (neko_bcknd_device .eq. 1) then
615 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
616 device_to_host, sync = .false.)
617 end if
618 end subroutine fluid_scheme_update_material_properties
619
623 subroutine fluid_scheme_set_material_properties(this, params, user)
624 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
625 type(json_file), intent(inout) :: params
626 type(user_t), target, intent(in) :: user
627 character(len=LOG_SIZE) :: log_buf
628 ! A local pointer that is needed to make Intel happy
629 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
630 logical :: nondimensional
631 real(kind=rp) :: dummy_lambda, dummy_cp
632 real(kind=rp) :: const_mu, const_rho
633 type(time_state_t) :: dummy_time_state
634
635
636 dummy_mp_ptr => dummy_user_material_properties
637
638 call neko_registry%add_field(this%dm_Xh, this%name // "_mu")
639 call neko_registry%add_field(this%dm_Xh, this%name // "_mu_tot")
640 call neko_registry%add_field(this%dm_Xh, this%name // "_rho")
641 this%mu => neko_registry%get_field(this%name // "_mu")
642 this%mu_tot => neko_registry%get_field(this%name // "_mu_tot")
643 this%rho => neko_registry%get_field(this%name // "_rho")
644
645 call this%material_properties%init(2)
646 call this%material_properties%assign(1, this%rho)
647 call this%material_properties%assign(2, this%mu)
648
649 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
650
651 write(log_buf, '(A)') 'Material properties must be set in the user' // &
652 ' file!'
653 call neko_log%message(log_buf)
654 this%user_material_properties => user%material_properties
655
656 call user%material_properties(this%name, this%material_properties, &
657 dummy_time_state)
658
659 else
660 this%user_material_properties => dummy_user_material_properties
661 ! Incorrect user input
662 if (params%valid_path('case.fluid.Re') .and. &
663 (params%valid_path('case.fluid.mu') .or. &
664 params%valid_path('case.fluid.rho'))) then
665 call neko_error("To set the material properties for the fluid, " // &
666 "either provide Re OR mu and rho in the case file.")
667
668 else if (params%valid_path('case.fluid.Re')) then
669 ! Non-dimensional case
670 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
671 & input.'
672 call neko_log%message(log_buf, lvl = neko_log_verbose)
673 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
674 & 1/Re.'
675 call neko_log%message(log_buf, lvl = neko_log_verbose)
676
677 ! Read Re into mu for further manipulation.
678 call json_get_or_lookup(params, 'case.fluid.Re', const_mu)
679 write(log_buf, '(A)') 'Read non-dimensional material properties'
680 call neko_log%message(log_buf)
681 write(log_buf, '(A,ES13.6)') 'Re :', const_mu
682 call neko_log%message(log_buf)
683
684 ! Set rho to 1 since the setup is non-dimensional.
685 const_rho = 1.0_rp
686 ! Invert the Re to get viscosity.
687 const_mu = 1.0_rp/const_mu
688 else
689 ! Dimensional case
690 call json_get_or_lookup(params, 'case.fluid.mu', const_mu)
691 call json_get_or_lookup(params, 'case.fluid.rho', const_rho)
692 end if
693 end if
694
695 ! We need to fill the fields based on the parsed const values
696 ! if the user routine is not used.
697 if (associated(user%material_properties, dummy_mp_ptr)) then
698 ! Fill mu and rho field with the physical value
699 call field_cfill(this%mu, const_mu)
700 call field_cfill(this%mu_tot, const_mu)
701 call field_cfill(this%rho, const_rho)
702
703
704 write(log_buf, '(A,ES13.6)') 'rho :', const_rho
705 call neko_log%message(log_buf)
706 write(log_buf, '(A,ES13.6)') 'mu :', const_mu
707 call neko_log%message(log_buf)
708 end if
709
710 ! Copy over material property to the total one
711 call field_copy(this%mu_tot, this%mu)
712
713 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
714 ! host array data to pass constant density and viscosity
715 ! to some routines, we need to make sure that the host
716 ! values are also filled
717 if (neko_bcknd_device .eq. 1) then
718 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
719 device_to_host, sync = .false.)
720 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
721 device_to_host, sync = .false.)
722 call device_memcpy(this%mu_tot%x, this%mu_tot%x_d, this%mu%size(), &
723 device_to_host, sync = .false.)
724 end if
725 end subroutine fluid_scheme_set_material_properties
726
Copy data between host and device (or device and device)
Definition device.F90:71
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.
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
Jacobi preconditioner accelerator backend.
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_event_sync(event)
Synchronize an event.
Definition device.F90:1468
integer, parameter, public device_to_host
Definition device.F90:47
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:62
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_addcol3(a, b, c, n)
Returns .
subroutine, public field_copy(a, b, n)
Copy a vector .
Defines a field.
Definition field.f90:34
subroutine fluid_scheme_set_material_properties(this, params, user)
Sets rho and mu.
subroutine fluid_scheme_update_material_properties(this, time)
Call user material properties routine and update the values of mu if necessary.
subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, pctype, pcparams)
Initialize a Krylov preconditioner.
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_bc_apply_vel(this, time, strong)
Apply all boundary conditions defined for velocity Here we perform additional gs operations to take c...
subroutine fluid_scheme_solver_factory(ksp, n, solver, max_iter, abstol, monitor)
Initialize a linear solver.
subroutine fluid_scheme_bc_apply_prs(this, time)
Apply all boundary conditions defined for pressure.
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.
Krylov preconditioner.
Definition pc_hsmg.f90:61
Jacobi preconditioner.
Definition pc_jacobi.f90:34
Utilities for retrieving parameters from the case files.
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:85
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:499
Defines a mesh.
Definition mesh.f90:34
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 solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
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:49
Defines a container for all statistics.
Jacobi preconditioner SX-Aurora backend.
Module with things related to the simulation time.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
subroutine, public dummy_user_material_properties(scheme_name, properties, time)
Utilities.
Definition utils.f90:35
Base type for a boundary condition.
Definition bc.f90:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:61
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
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:73
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
Defines a jacobi preconditioner for SX-Aurora.
A struct that contains all info about the time, expand as needed.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...