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