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