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 if (associated(bc)) then
371 call bc%free()
372 deallocate(bc)
373 end if
374 end do
375 call this%bcs_vel%free()
376
377 do i = 1, this%bcs_prs%size()
378 bc => this%bcs_prs%get(i)
379 if (associated(bc)) then
380 call bc%free()
381 deallocate(bc)
382 end if
383 end do
384 call this%bcs_prs%free()
385
386 call this%source_term%free()
387
388 call this%gs_Xh%free()
389
390 call this%c_Xh%free()
391
392 nullify(this%u)
393 nullify(this%v)
394 nullify(this%w)
395 nullify(this%p)
396
397 nullify(this%u_e)
398 nullify(this%v_e)
399 nullify(this%w_e)
400
401 call this%ulag%free()
402 call this%vlag%free()
403 call this%wlag%free()
404
405
406 if (associated(this%f_x)) then
407 call this%f_x%free()
408 deallocate(this%f_x)
409 end if
410
411 if (associated(this%f_y)) then
412 call this%f_y%free()
413 deallocate(this%f_y)
414 end if
415
416 if (associated(this%f_z)) then
417 call this%f_z%free()
418 deallocate(this%f_z)
419 end if
420
421 nullify(this%f_x)
422 nullify(this%f_y)
423 nullify(this%f_z)
424 nullify(this%rho)
425 nullify(this%mu)
426 nullify(this%mu_tot)
427
428 call this%dm_Xh%free()
429 call this%Xh%free()
430 nullify(this%msh)
431
432 end subroutine fluid_scheme_free
433
436 subroutine fluid_scheme_validate(this)
437 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
438 ! Variables for retrieving json parameters
439 logical :: logical_val
440
441 if ( (.not. associated(this%u)) .or. &
442 (.not. associated(this%v)) .or. &
443 (.not. associated(this%w)) .or. &
444 (.not. associated(this%p))) then
445 call neko_error('Fields are not registered')
446 end if
447
448 if ( (.not. allocated(this%u%x)) .or. &
449 (.not. allocated(this%v%x)) .or. &
450 (.not. allocated(this%w%x)) .or. &
451 (.not. allocated(this%p%x))) then
452 call neko_error('Fields are not allocated')
453 end if
454
455 if (.not. allocated(this%ksp_vel)) then
456 call neko_error('No Krylov solver for velocity defined')
457 end if
458
459 if (.not. allocated(this%ksp_prs)) then
460 call neko_error('No Krylov solver for pressure defined')
461 end if
462
463 end subroutine fluid_scheme_validate
464
469 subroutine fluid_scheme_bc_apply_vel(this, time, strong)
470 class(fluid_scheme_incompressible_t), intent(inout) :: this
471 type(time_state_t), intent(in) :: time
472 logical, intent(in) :: strong
473 integer :: i
474 class(bc_t), pointer :: b
475 b => null()
476
477 call this%bcs_vel%apply_vector(&
478 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), time, strong)
479
480 call rotate_cyc(this%u, this%v, this%w, 1, this%c_Xh)
481 call this%gs_Xh%op(this%u, gs_op_min, glb_cmd_event)
483 call this%gs_Xh%op(this%v, gs_op_min, glb_cmd_event)
485 call this%gs_Xh%op(this%w, gs_op_min, glb_cmd_event)
487 call rotate_cyc(this%u, this%v, this%w, 0, this%c_Xh)
488
489
490 call this%bcs_vel%apply_vector(&
491 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), time, strong)
492
493 call rotate_cyc(this%u, this%v, this%w, 1, this%c_Xh)
494 call this%gs_Xh%op(this%u, gs_op_max, glb_cmd_event)
496 call this%gs_Xh%op(this%v, gs_op_max, glb_cmd_event)
498 call this%gs_Xh%op(this%w, gs_op_max, glb_cmd_event)
500 call rotate_cyc(this%u, this%v, this%w, 0, this%c_Xh)
501
502 do i = 1, this%bcs_vel%size()
503 b => this%bcs_vel%get(i)
504 b%updated = .false.
505 end do
506 nullify(b)
507
508 end subroutine fluid_scheme_bc_apply_vel
509
512 subroutine fluid_scheme_bc_apply_prs(this, time)
513 class(fluid_scheme_incompressible_t), intent(inout) :: this
514 type(time_state_t), intent(in) :: time
515
516 integer :: i
517 class(bc_t), pointer :: b
518 b => null()
519
520 call this%bcs_prs%apply(this%p, time)
521 call this%gs_Xh%op(this%p, gs_op_min, glb_cmd_event)
523
524 call this%bcs_prs%apply(this%p, time)
525 call this%gs_Xh%op(this%p, gs_op_max, glb_cmd_event)
527
528 do i = 1, this%bcs_prs%size()
529 b => this%bcs_prs%get(i)
530 b%updated = .false.
531 end do
532 nullify(b)
533
534 end subroutine fluid_scheme_bc_apply_prs
535
538 subroutine fluid_scheme_solver_factory(ksp, n, solver, &
539 max_iter, abstol, monitor)
540 class(ksp_t), allocatable, target, intent(inout) :: ksp
541 integer, intent(in), value :: n
542 character(len=*), intent(in) :: solver
543 integer, intent(in) :: max_iter
544 real(kind=rp), intent(in) :: abstol
545 logical, intent(in) :: monitor
546
547 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
548 monitor = monitor)
549
550 end subroutine fluid_scheme_solver_factory
551
553 subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, &
554 pctype, pcparams)
555 class(fluid_scheme_incompressible_t), intent(inout) :: this
556 class(pc_t), allocatable, target, intent(inout) :: pc
557 class(ksp_t), target, intent(inout) :: ksp
558 type(coef_t), target, intent(in) :: coef
559 type(dofmap_t), target, intent(in) :: dof
560 type(gs_t), target, intent(inout) :: gs
561 type(bc_list_t), target, intent(inout) :: bclst
562 character(len=*) :: pctype
563 type(json_file), intent(inout) :: pcparams
564
565 call precon_factory(pc, pctype)
566
567 select type (pcp => pc)
568 type is (jacobi_t)
569 call pcp%init(coef, dof, gs)
570 type is (sx_jacobi_t)
571 call pcp%init(coef, dof, gs)
572 type is (device_jacobi_t)
573 call pcp%init(coef, dof, gs)
574 type is (hsmg_t)
575 call pcp%init(coef, bclst, pcparams)
576 type is (phmg_t)
577 call pcp%init(coef, bclst, pcparams)
578 end select
579
580 call ksp%set_pc(pc)
581
582 end subroutine fluid_scheme_precon_factory
583
585 function fluid_compute_cfl(this, dt) result(c)
586 class(fluid_scheme_incompressible_t), intent(in) :: this
587 real(kind=rp), intent(in) :: dt
588 real(kind=rp) :: c
589
590 c = cfl(dt, this%u, this%v, this%w, &
591 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
592
593 end function fluid_compute_cfl
594
595
600 subroutine fluid_scheme_update_material_properties(this, time)
601 class(fluid_scheme_incompressible_t), intent(inout) :: this
602 type(time_state_t), intent(in) :: time
603 type(field_t), pointer :: nut
604
605 call this%user_material_properties(this%name, this%material_properties, &
606 time)
607
608 if (len(trim(this%nut_field_name)) > 0) then
609 nut => neko_registry%get_field(this%nut_field_name)
610 ! Copy material property
611 call field_copy(this%mu_tot, this%mu)
612 ! Add turbulent contribution
613 call field_addcol3(this%mu_tot, nut, this%rho)
614 end if
615
616 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
617 ! host array data to pass constant density and viscosity
618 ! to some routines, we need to make sure that the host
619 ! values are also filled
620 if (neko_bcknd_device .eq. 1) then
621 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
622 device_to_host, sync = .false.)
623 end if
624 end subroutine fluid_scheme_update_material_properties
625
629 subroutine fluid_scheme_set_material_properties(this, params, user)
630 class(fluid_scheme_incompressible_t), target, intent(inout) :: this
631 type(json_file), intent(inout) :: params
632 type(user_t), target, intent(in) :: user
633 character(len=LOG_SIZE) :: log_buf
634 ! A local pointer that is needed to make Intel happy
635 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
636 logical :: nondimensional
637 real(kind=rp) :: dummy_lambda, dummy_cp
638 real(kind=rp) :: const_mu, const_rho
639 type(time_state_t) :: dummy_time_state
640
641
642 dummy_mp_ptr => dummy_user_material_properties
643
644 call neko_registry%add_field(this%dm_Xh, this%name // "_mu")
645 call neko_registry%add_field(this%dm_Xh, this%name // "_mu_tot")
646 call neko_registry%add_field(this%dm_Xh, this%name // "_rho")
647 this%mu => neko_registry%get_field(this%name // "_mu")
648 this%mu_tot => neko_registry%get_field(this%name // "_mu_tot")
649 this%rho => neko_registry%get_field(this%name // "_rho")
650
651 call this%material_properties%init(2)
652 call this%material_properties%assign(1, this%rho)
653 call this%material_properties%assign(2, this%mu)
654
655 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
656
657 write(log_buf, '(A)') 'Material properties must be set in the user' // &
658 ' file!'
659 call neko_log%message(log_buf)
660 this%user_material_properties => user%material_properties
661
662 call user%material_properties(this%name, this%material_properties, &
663 dummy_time_state)
664
665 else
666 this%user_material_properties => dummy_user_material_properties
667 ! Incorrect user input
668 if (params%valid_path('case.fluid.Re') .and. &
669 (params%valid_path('case.fluid.mu') .or. &
670 params%valid_path('case.fluid.rho'))) then
671 call neko_error("To set the material properties for the fluid, " // &
672 "either provide Re OR mu and rho in the case file.")
673
674 else if (params%valid_path('case.fluid.Re')) then
675 ! Non-dimensional case
676 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
677 & input.'
678 call neko_log%message(log_buf, lvl = neko_log_verbose)
679 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
680 & 1/Re.'
681 call neko_log%message(log_buf, lvl = neko_log_verbose)
682
683 ! Read Re into mu for further manipulation.
684 call json_get_or_lookup(params, 'case.fluid.Re', const_mu)
685 write(log_buf, '(A)') 'Read non-dimensional material properties'
686 call neko_log%message(log_buf)
687 write(log_buf, '(A,ES13.6)') 'Re :', const_mu
688 call neko_log%message(log_buf)
689
690 ! Set rho to 1 since the setup is non-dimensional.
691 const_rho = 1.0_rp
692 ! Invert the Re to get viscosity.
693 const_mu = 1.0_rp/const_mu
694 else
695 ! Dimensional case
696 call json_get_or_lookup(params, 'case.fluid.mu', const_mu)
697 call json_get_or_lookup(params, 'case.fluid.rho', const_rho)
698 end if
699 end if
700
701 ! We need to fill the fields based on the parsed const values
702 ! if the user routine is not used.
703 if (associated(user%material_properties, dummy_mp_ptr)) then
704 ! Fill mu and rho field with the physical value
705 call field_cfill(this%mu, const_mu)
706 call field_cfill(this%mu_tot, const_mu)
707 call field_cfill(this%rho, const_rho)
708
709
710 write(log_buf, '(A,ES13.6)') 'rho :', const_rho
711 call neko_log%message(log_buf)
712 write(log_buf, '(A,ES13.6)') 'mu :', const_mu
713 call neko_log%message(log_buf)
714 end if
715
716 ! Copy over material property to the total one
717 call field_copy(this%mu_tot, this%mu)
718
719 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
720 ! host array data to pass constant density and viscosity
721 ! to some routines, we need to make sure that the host
722 ! values are also filled
723 if (neko_bcknd_device .eq. 1) then
724 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
725 device_to_host, sync = .false.)
726 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
727 device_to_host, sync = .false.)
728 call device_memcpy(this%mu_tot%x, this%mu_tot%x_d, this%mu%size(), &
729 device_to_host, sync = .false.)
730 end if
731 end subroutine fluid_scheme_set_material_properties
732
Copy data between host and device (or device and device)
Definition device.F90:72
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.
Compute CFL condition.
Definition operators.f90:97
Apply cyclic boundary condition to a vector field.
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:1594
integer, parameter, public device_to_host
Definition device.F90:48
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:63
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:86
type(log_t), public neko_log
Global log stream.
Definition log.f90:80
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:627
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
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:49
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
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...