Neko  0.8.99
A portable framework for high-order spectral element flow simulations
fluid_scheme.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
35  use gather_scatter, only : gs_t
36  use mean_sqr_flow, only : mean_sqr_flow_t
37  use neko_config, only : neko_bcknd_device
38  use checkpoint, only : chkp_t
39  use mean_flow, only : mean_flow_t
40  use num_types, only : rp
41  use comm
43  use field, only : field_t
44  use space, only : space_t, gll
45  use dofmap, only : dofmap_t
46  use krylov, only : ksp_t, krylov_solver_factory, krylov_solver_destroy, &
48  use coefs, only: coef_t
49  use wall, only : no_slip_wall_t
50  use inflow, only : inflow_t
52  use blasius, only : blasius_t
53  use dirichlet, only : dirichlet_t
54  use dong_outflow, only : dong_outflow_t
55  use symmetry, only : symmetry_t
58  use jacobi, only : jacobi_t
59  use sx_jacobi, only : sx_jacobi_t
60  use device_jacobi, only : device_jacobi_t
61  use hsmg, only : hsmg_t
62  use precon, only : pc_t, precon_factory, precon_destroy
63  use fluid_stats, only : fluid_stats_t
67  use math, only : cfill, add2s2
70  use operators, only : cfl
71  use logger, only : neko_log, log_size
74  use json_module, only : json_file
76  use user_intf, only : user_t
77  use utils, only : neko_error
79  use field_series, only : field_series_t
81  use field_math, only : field_cfill
82  implicit none
83  private
84 
86  type, abstract :: fluid_scheme_t
87  type(field_t), pointer :: u => null()
88  type(field_t), pointer :: v => null()
89  type(field_t), pointer :: w => null()
90  type(field_t), pointer :: p => null()
91  type(field_series_t) :: ulag, vlag, wlag
92  type(space_t) :: xh
93  type(dofmap_t) :: dm_xh
94  type(gs_t) :: gs_xh
95  type(coef_t) :: c_xh
99  type(field_t), pointer :: f_x => null()
101  type(field_t), pointer :: f_y => null()
103  type(field_t), pointer :: f_z => null()
104  class(ksp_t), allocatable :: ksp_vel
105  class(ksp_t), allocatable :: ksp_prs
106  class(pc_t), allocatable :: pc_vel
107  class(pc_t), allocatable :: pc_prs
108  integer :: vel_projection_dim
109  integer :: pr_projection_dim
110  integer :: vel_projection_activ_step
111  integer :: pr_projection_activ_step
112  type(no_slip_wall_t) :: bc_wall
113  class(bc_t), allocatable :: bc_inflow
114 
115  ! Attributes for field dirichlet BCs
116  type(field_dirichlet_vector_t) :: user_field_bc_vel
117  type(field_dirichlet_t) :: user_field_bc_prs
118  type(dirichlet_t) :: bc_prs
119  type(dong_outflow_t) :: bc_dong
120  type(symmetry_t) :: bc_sym
121  type(bc_list_t) :: bclst_vel
122  type(bc_list_t) :: bclst_prs
123  type(field_t) :: bdry
124  type(json_file), pointer :: params
125  type(mesh_t), pointer :: msh => null()
126  type(chkp_t) :: chkp
127  type(mean_flow_t) :: mean
129  type(mean_sqr_flow_t) :: mean_sqr
130  logical :: forced_flow_rate = .false.
131  logical :: freeze = .false.
133  real(kind=rp), pointer :: mu => null()
135  type(field_t) :: mu_field
137  character(len=:), allocatable :: nut_field_name
139  logical :: variable_material_properties = .false.
141  real(kind=rp), pointer :: rho => null()
143  type(field_t) :: rho_field
144  type(scratch_registry_t) :: scratch
146  character(len=NEKO_MSH_MAX_ZLBL_LEN), allocatable :: bc_labels(:)
147  contains
149  procedure, pass(this) :: fluid_scheme_init_all
150  procedure, pass(this) :: fluid_scheme_init_common
153  procedure, pass(this) :: scheme_free => fluid_scheme_free
155  procedure, pass(this) :: validate => fluid_scheme_validate
157  procedure, pass(this) :: bc_apply_vel => fluid_scheme_bc_apply_vel
159  procedure, pass(this) :: bc_apply_prs => fluid_scheme_bc_apply_prs
161  procedure, pass(this) :: set_usr_inflow => fluid_scheme_set_usr_inflow
163  procedure, pass(this) :: compute_cfl => fluid_compute_cfl
165  procedure(fluid_scheme_init_intrf), pass(this), deferred :: init
167  procedure(fluid_scheme_free_intrf), pass(this), deferred :: free
169  procedure(fluid_scheme_step_intrf), pass(this), deferred :: step
171  procedure(fluid_scheme_restart_intrf), pass(this), deferred :: restart
172  procedure, private, pass(this) :: set_bc_type_output => &
175  procedure, pass(this) :: update_material_properties => &
177  end type fluid_scheme_t
178 
180  abstract interface
181  subroutine fluid_scheme_init_intrf(this, msh, lx, params, user, &
182  material_properties, time_scheme)
183  import fluid_scheme_t
184  import json_file
185  import mesh_t
186  import user_t
187  import material_properties_t
189  class(fluid_scheme_t), target, intent(inout) :: this
190  type(mesh_t), target, intent(inout) :: msh
191  integer, intent(inout) :: lx
192  type(json_file), target, intent(inout) :: params
193  type(user_t), intent(in) :: user
194  type(material_properties_t), &
195  target, intent(inout) :: material_properties
196  type(time_scheme_controller_t), target, intent(in) :: time_scheme
197  end subroutine fluid_scheme_init_intrf
198  end interface
199 
201  abstract interface
202  subroutine fluid_scheme_free_intrf(this)
203  import fluid_scheme_t
204  class(fluid_scheme_t), intent(inout) :: this
205  end subroutine fluid_scheme_free_intrf
206  end interface
207 
209  abstract interface
210  subroutine fluid_scheme_step_intrf(this, t, tstep, dt, ext_bdf, &
211  dt_controller)
212  import fluid_scheme_t
215  import rp
216  class(fluid_scheme_t), target, intent(inout) :: this
217  real(kind=rp), intent(inout) :: t
218  integer, intent(inout) :: tstep
219  real(kind=rp), intent(in) :: dt
220  type(time_scheme_controller_t), intent(inout) :: ext_bdf
221  type(time_step_controller_t), intent(in) :: dt_controller
222  end subroutine fluid_scheme_step_intrf
223  end interface
224 
226  abstract interface
227  subroutine fluid_scheme_restart_intrf(this, dtlag, tlag)
228  import fluid_scheme_t
229  import rp
230  class(fluid_scheme_t), target, intent(inout) :: this
231  real(kind=rp) :: dtlag(10), tlag(10)
232 
233  end subroutine fluid_scheme_restart_intrf
234  end interface
235 
236  interface
237 
238  module subroutine fluid_scheme_factory(object, type_name)
239  class(fluid_scheme_t), intent(inout), allocatable :: object
240  character(len=*) :: type_name
241  end subroutine fluid_scheme_factory
242  end interface
243 
244  public :: fluid_scheme_t, fluid_scheme_factory
245 
246 contains
247 
249  subroutine fluid_scheme_init_common(this, msh, lx, params, scheme, user, &
250  material_properties, kspv_init)
251  implicit none
252  class(fluid_scheme_t), target, intent(inout) :: this
253  type(mesh_t), target, intent(inout) :: msh
254  integer, intent(inout) :: lx
255  character(len=*), intent(in) :: scheme
256  type(json_file), target, intent(inout) :: params
257  type(user_t), target, intent(in) :: user
258  type(material_properties_t), target, intent(inout) :: material_properties
259  logical, intent(in) :: kspv_init
260  type(dirichlet_t) :: bdry_mask
261  character(len=LOG_SIZE) :: log_buf
262  real(kind=rp), allocatable :: real_vec(:)
263  real(kind=rp) :: real_val
264  logical :: logical_val
265  integer :: integer_val, ierr
266  character(len=:), allocatable :: string_val1, string_val2
267 
268 
269  !
270  ! SEM simulation fundamentals
271  !
272 
273  this%msh => msh
274 
275  if (msh%gdim .eq. 2) then
276  call this%Xh%init(gll, lx, lx)
277  else
278  call this%Xh%init(gll, lx, lx, lx)
279  end if
280 
281  this%dm_Xh = dofmap_t(msh, this%Xh)
282 
283  call this%gs_Xh%init(this%dm_Xh)
284 
285  call this%c_Xh%init(this%gs_Xh)
286 
287  ! Local scratch registry
288  this%scratch = scratch_registry_t(this%dm_Xh, 10, 2)
289 
290  ! Case parameters
291  this%params => params
292 
293  !
294  ! Material properties
295  !
296 
297  this%rho => material_properties%rho
298  this%mu => material_properties%mu
299 
300  !
301  ! Turbulence modelling and variable material properties
302  !
303  if (params%valid_path('case.fluid.nut_field')) then
304  call json_get(params, 'case.fluid.nut_field', this%nut_field_name)
305  this%variable_material_properties = .true.
306  else
307  this%nut_field_name = ""
308  end if
309 
310  ! Fill mu and rho field with the physical value
311 
312  call this%mu_field%init(this%dm_Xh, "mu")
313  call this%rho_field%init(this%dm_Xh, "mu")
314  call field_cfill(this%mu_field, this%mu, this%mu_field%size())
315  call field_cfill(this%rho_field, this%rho, this%mu_field%size())
316 
317  ! Since mu, rho is a field, and the none-stress simulation fetches
318  ! data from the host arrays, we need to mirror the constant
319  ! material properties on the host
320  if (neko_bcknd_device .eq. 1) then
321  call cfill(this%mu_field%x, this%mu, this%mu_field%size())
322  call cfill(this%rho_field%x, this%rho, this%rho_field%size())
323  end if
324 
325 
326  ! Projection spaces
327  call json_get_or_default(params, &
328  'case.fluid.velocity_solver.projection_space_size', &
329  this%vel_projection_dim, 20)
330  call json_get_or_default(params, &
331  'case.fluid.pressure_solver.projection_space_size', &
332  this%pr_projection_dim, 20)
333  call json_get_or_default(params, &
334  'case.fluid.velocity_solver.projection_hold_steps', &
335  this%vel_projection_activ_step, 5)
336  call json_get_or_default(params, &
337  'case.fluid.pressure_solver.projection_hold_steps', &
338  this%pr_projection_activ_step, 5)
339 
340 
341  call json_get_or_default(params, 'case.fluid.freeze', this%freeze, .false.)
342 
343  if (params%valid_path("case.fluid.flow_rate_force")) then
344  this%forced_flow_rate = .true.
345  end if
346 
347 
348  !
349  ! First section of fluid log
350  !
351 
352  call neko_log%section('Fluid')
353  write(log_buf, '(A, A)') 'Type : ', trim(scheme)
354  call neko_log%message(log_buf)
355  if (lx .lt. 10) then
356  write(log_buf, '(A, I1)') 'Poly order : ', lx-1
357  else if (lx .ge. 10) then
358  write(log_buf, '(A, I2)') 'Poly order : ', lx-1
359  else
360  write(log_buf, '(A, I3)') 'Poly order : ', lx-1
361  end if
362  call neko_log%message(log_buf)
363  write(log_buf, '(A, I0)') 'DoFs : ', this%dm_Xh%size()
364  call neko_log%message(log_buf)
365 
366  write(log_buf, '(A,ES13.6)') 'rho :', this%rho
367  call neko_log%message(log_buf)
368  write(log_buf, '(A,ES13.6)') 'mu :', this%mu
369  call neko_log%message(log_buf)
370 
371  call json_get(params, 'case.numerics.dealias', logical_val)
372  write(log_buf, '(A, L1)') 'Dealias : ', logical_val
373  call neko_log%message(log_buf)
374 
375  call json_get_or_default(params, 'case.output_boundary', logical_val, &
376  .false.)
377  write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
378  call neko_log%message(log_buf)
379  if (logical_val) then
380  write(log_buf, '(A)') 'bdry keys: '
381  call neko_log%message(log_buf)
382  write(log_buf, '(A)') 'No-slip wall ''w'' = 1'
383  call neko_log%message(log_buf)
384  write(log_buf, '(A)') 'Inlet/velocity dirichlet ''v'' = 2'
385  call neko_log%message(log_buf)
386  write(log_buf, '(A)') 'Outlet ''o'' = 3'
387  call neko_log%message(log_buf)
388  write(log_buf, '(A)') 'Symmetry ''sym'' = 4'
389  call neko_log%message(log_buf)
390  write(log_buf, '(A)') 'Outlet-normal ''on'' = 5'
391  call neko_log%message(log_buf)
392  write(log_buf, '(A)') 'Periodic = 6'
393  call neko_log%message(log_buf)
394  write(log_buf, '(A)') 'Dirichlet on velocity components: '
395  call neko_log%message(log_buf)
396  write(log_buf, '(A)') ' ''d_vel_u, d_vel_v, d_vel_w'' = 7'
397  call neko_log%message(log_buf)
398  write(log_buf, '(A)') 'Pressure dirichlet ''d_pres'' = 8'
399  call neko_log%message(log_buf)
400  write(log_buf, '(A)') '''d_vel_(u,v,w)'' and ''d_pres'' = 8'
401  call neko_log%message(log_buf)
402  write(log_buf, '(A)') 'No boundary condition set = 0'
403  call neko_log%message(log_buf)
404  end if
405 
406 
407  !
408  ! Setup velocity boundary conditions
409  !
410  allocate(this%bc_labels(neko_msh_max_zlbls))
411  this%bc_labels = "not"
412 
413  if (params%valid_path('case.fluid.boundary_types')) then
414  call json_get(params, &
415  'case.fluid.boundary_types', &
416  this%bc_labels)
417  end if
418 
419  call bc_list_init(this%bclst_vel)
420 
421  call this%bc_sym%init_base(this%c_Xh)
422  call this%bc_sym%mark_zone(msh%sympln)
423  call this%bc_sym%mark_zones_from_list(msh%labeled_zones, &
424  'sym', this%bc_labels)
425  call this%bc_sym%finalize()
426  call this%bc_sym%init(this%c_Xh)
427  call bc_list_add(this%bclst_vel, this%bc_sym)
428 
429  !
430  ! Inflow
431  !
432  if (params%valid_path('case.fluid.inflow_condition')) then
433  call json_get(params, 'case.fluid.inflow_condition.type', string_val1)
434  if (trim(string_val1) .eq. "uniform") then
435  allocate(inflow_t::this%bc_inflow)
436  else if (trim(string_val1) .eq. "blasius") then
437  allocate(blasius_t::this%bc_inflow)
438  else if (trim(string_val1) .eq. "user") then
439  allocate(usr_inflow_t::this%bc_inflow)
440  else
441  call neko_error('Invalid inflow condition '//string_val1)
442  end if
443 
444  call this%bc_inflow%init_base(this%c_Xh)
445  call this%bc_inflow%mark_zone(msh%inlet)
446  call this%bc_inflow%mark_zones_from_list(msh%labeled_zones, &
447  'v', this%bc_labels)
448  call this%bc_inflow%finalize()
449  call bc_list_add(this%bclst_vel, this%bc_inflow)
450 
451  if (trim(string_val1) .eq. "uniform") then
452  call json_get(params, 'case.fluid.inflow_condition.value', real_vec)
453  select type (bc_if => this%bc_inflow)
454  type is (inflow_t)
455  call bc_if%set_inflow(real_vec)
456  end select
457  else if (trim(string_val1) .eq. "blasius") then
458  select type (bc_if => this%bc_inflow)
459  type is (blasius_t)
460  call json_get(params, 'case.fluid.blasius.delta', real_val)
461  call json_get(params, 'case.fluid.blasius.approximation', &
462  string_val2)
463  call json_get(params, 'case.fluid.blasius.freestream_velocity', &
464  real_vec)
465 
466  call bc_if%set_params(real_vec, real_val, string_val2)
467 
468  end select
469  else if (trim(string_val1) .eq. "user") then
470  end if
471  end if
472 
473  call this%bc_wall%init_base(this%c_Xh)
474  call this%bc_wall%mark_zone(msh%wall)
475  call this%bc_wall%mark_zones_from_list(msh%labeled_zones, &
476  'w', this%bc_labels)
477  call this%bc_wall%finalize()
478  call bc_list_add(this%bclst_vel, this%bc_wall)
479 
480  ! Setup field dirichlet bc for u-velocity
481  call this%user_field_bc_vel%bc_u%init_base(this%c_Xh)
482  call this%user_field_bc_vel%bc_u%mark_zones_from_list(msh%labeled_zones, &
483  'd_vel_u', this%bc_labels)
484  call this%user_field_bc_vel%bc_u%finalize()
485 
486  call mpi_allreduce(this%user_field_bc_vel%bc_u%msk(0), integer_val, 1, &
487  mpi_integer, mpi_sum, neko_comm, ierr)
488  if (integer_val .gt. 0) then
489  call this%user_field_bc_vel%bc_u%init_field('d_vel_u')
490  end if
491 
492  ! Setup field dirichlet bc for v-velocity
493  call this%user_field_bc_vel%bc_v%init_base(this%c_Xh)
494  call this%user_field_bc_vel%bc_v%mark_zones_from_list(msh%labeled_zones, &
495  'd_vel_v', this%bc_labels)
496  call this%user_field_bc_vel%bc_v%finalize()
497 
498  call mpi_allreduce(this%user_field_bc_vel%bc_v%msk(0), integer_val, 1, &
499  mpi_integer, mpi_sum, neko_comm, ierr)
500  if (integer_val .gt. 0) then
501  call this%user_field_bc_vel%bc_v%init_field('d_vel_v')
502  end if
503 
504  ! Setup field dirichlet bc for w-velocity
505  call this%user_field_bc_vel%bc_w%init_base(this%c_Xh)
506  call this%user_field_bc_vel%bc_w%mark_zones_from_list(msh%labeled_zones, &
507  'd_vel_w', this%bc_labels)
508  call this%user_field_bc_vel%bc_w%finalize()
509 
510  call mpi_allreduce(this%user_field_bc_vel%bc_w%msk(0), integer_val, 1, &
511  mpi_integer, mpi_sum, neko_comm, ierr)
512  if (integer_val .gt. 0) then
513  call this%user_field_bc_vel%bc_w%init_field('d_vel_w')
514  end if
515 
516  ! Setup our global field dirichlet bc
517  call this%user_field_bc_vel%init_base(this%c_Xh)
518  call this%user_field_bc_vel%mark_zones_from_list(msh%labeled_zones, &
519  'd_vel_u', this%bc_labels)
520  call this%user_field_bc_vel%mark_zones_from_list(msh%labeled_zones, &
521  'd_vel_v', this%bc_labels)
522  call this%user_field_bc_vel%mark_zones_from_list(msh%labeled_zones, &
523  'd_vel_w', this%bc_labels)
524  call this%user_field_bc_vel%finalize()
525 
526  ! Add the field bc to velocity bcs
527  call bc_list_add(this%bclst_vel, this%user_field_bc_vel)
528 
529  !
530  ! Associate our field dirichlet update to the user one.
531  !
532  this%user_field_bc_vel%update => user%user_dirichlet_update
533 
534  !
535  ! Initialize field list and bc list for user_dirichlet_update
536  !
537 
538  ! Note, some of these are potentially not initialized !
539  call this%user_field_bc_vel%field_list%init(4)
540  call this%user_field_bc_vel%field_list%assign_to_field(1, &
541  this%user_field_bc_vel%bc_u%field_bc)
542  call this%user_field_bc_vel%field_list%assign_to_field(2, &
543  this%user_field_bc_vel%bc_v%field_bc)
544  call this%user_field_bc_vel%field_list%assign_to_field(3, &
545  this%user_field_bc_vel%bc_w%field_bc)
546  call this%user_field_bc_vel%field_list%assign_to_field(4, &
547  this%user_field_bc_prs%field_bc)
548 
549  call bc_list_init(this%user_field_bc_vel%bc_list, size = 4)
550  ! Note, bc_list_add only adds if the bc is not empty
551  call bc_list_add(this%user_field_bc_vel%bc_list, &
552  this%user_field_bc_vel%bc_u)
553  call bc_list_add(this%user_field_bc_vel%bc_list, &
554  this%user_field_bc_vel%bc_v)
555  call bc_list_add(this%user_field_bc_vel%bc_list, &
556  this%user_field_bc_vel%bc_w)
557 
558  !
559  ! Check if we need to output boundary types to a separate field
560  call fluid_scheme_set_bc_type_output(this, params)
561 
562  !
563  ! Setup right-hand side fields.
564  !
565  allocate(this%f_x)
566  allocate(this%f_y)
567  allocate(this%f_z)
568  call this%f_x%init(this%dm_Xh, fld_name = "fluid_rhs_x")
569  call this%f_y%init(this%dm_Xh, fld_name = "fluid_rhs_y")
570  call this%f_z%init(this%dm_Xh, fld_name = "fluid_rhs_z")
571 
572  ! Initialize the source term
573  call this%source_term%init(this%f_x, this%f_y, this%f_z, this%c_Xh, user)
574  call this%source_term%add(params, 'case.fluid.source_term')
575 
576  ! If case.output_boundary is true, set the values for the bc types in the
577  ! output of the field.
578  call this%set_bc_type_output(params)
579 
580  ! Initialize velocity solver
581  if (kspv_init) then
582  call neko_log%section("Velocity solver")
583  call json_get_or_default(params, &
584  'case.fluid.velocity_solver.max_iterations', &
585  integer_val, ksp_max_iter)
586  call json_get(params, 'case.fluid.velocity_solver.type', string_val1)
587  call json_get(params, 'case.fluid.velocity_solver.preconditioner', &
588  string_val2)
589  call json_get(params, 'case.fluid.velocity_solver.absolute_tolerance', &
590  real_val)
591  call json_get_or_default(params, &
592  'case.fluid.velocity_solver.monitor', &
593  logical_val, .false.)
594 
595  call neko_log%message('Type : ('// trim(string_val1) // &
596  ', ' // trim(string_val2) // ')')
597 
598  write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
599  call neko_log%message(log_buf)
600  call fluid_scheme_solver_factory(this%ksp_vel, this%dm_Xh%size(), &
601  string_val1, integer_val, real_val, logical_val)
602  call fluid_scheme_precon_factory(this%pc_vel, this%ksp_vel, &
603  this%c_Xh, this%dm_Xh, this%gs_Xh, this%bclst_vel, string_val2)
604  call neko_log%end_section()
605  end if
606 
607  ! Assign velocity fields
608  call neko_field_registry%add_field(this%dm_Xh, 'u')
609  call neko_field_registry%add_field(this%dm_Xh, 'v')
610  call neko_field_registry%add_field(this%dm_Xh, 'w')
611  this%u => neko_field_registry%get_field('u')
612  this%v => neko_field_registry%get_field('v')
613  this%w => neko_field_registry%get_field('w')
614 
615  !! Initialize time-lag fields
616  call this%ulag%init(this%u, 2)
617  call this%vlag%init(this%v, 2)
618  call this%wlag%init(this%w, 2)
619 
620 
621  end subroutine fluid_scheme_init_common
622 
624  subroutine fluid_scheme_init_all(this, msh, lx, params, kspv_init, &
625  kspp_init, scheme, user, &
626  material_properties)
627  implicit none
628  class(fluid_scheme_t), target, intent(inout) :: this
629  type(mesh_t), target, intent(inout) :: msh
630  integer, intent(inout) :: lx
631  type(json_file), target, intent(inout) :: params
632  type(user_t), target, intent(in) :: user
633  type(material_properties_t), target, intent(inout) :: material_properties
634  logical :: kspv_init
635  logical :: kspp_init
636  character(len=*), intent(in) :: scheme
637  real(kind=rp) :: abs_tol
638  integer :: integer_val, ierr
639  logical :: logical_val
640  character(len=:), allocatable :: solver_type, precon_type
641  character(len=LOG_SIZE) :: log_buf
642 
643  call fluid_scheme_init_common(this, msh, lx, params, scheme, user, &
644  material_properties, kspv_init)
645 
646  call neko_field_registry%add_field(this%dm_Xh, 'p')
647  this%p => neko_field_registry%get_field('p')
648 
649  !
650  ! Setup pressure boundary conditions
651  !
652  call bc_list_init(this%bclst_prs)
653  call this%bc_prs%init_base(this%c_Xh)
654  call this%bc_prs%mark_zones_from_list(msh%labeled_zones, &
655  'o', this%bc_labels)
656  call this%bc_prs%mark_zones_from_list(msh%labeled_zones, &
657  'on', this%bc_labels)
658 
659  ! Field dirichlet pressure bc
660  call this%user_field_bc_prs%init_base(this%c_Xh)
661  call this%user_field_bc_prs%mark_zones_from_list(msh%labeled_zones, &
662  'd_pres', this%bc_labels)
663  call this%user_field_bc_prs%finalize()
664  call mpi_allreduce(this%user_field_bc_prs%msk(0), integer_val, 1, &
665  mpi_integer, mpi_sum, neko_comm, ierr)
666 
667  if (integer_val .gt. 0) call this%user_field_bc_prs%init_field('d_pres')
668  call bc_list_add(this%bclst_prs, this%user_field_bc_prs)
669  call bc_list_add(this%user_field_bc_vel%bc_list, this%user_field_bc_prs)
670 
671  if (msh%outlet%size .gt. 0) then
672  call this%bc_prs%mark_zone(msh%outlet)
673  end if
674  if (msh%outlet_normal%size .gt. 0) then
675  call this%bc_prs%mark_zone(msh%outlet_normal)
676  end if
677 
678  call this%bc_prs%finalize()
679  call this%bc_prs%set_g(0.0_rp)
680  call bc_list_add(this%bclst_prs, this%bc_prs)
681  call this%bc_dong%init_base(this%c_Xh)
682  call this%bc_dong%mark_zones_from_list(msh%labeled_zones, &
683  'o+dong', this%bc_labels)
684  call this%bc_dong%mark_zones_from_list(msh%labeled_zones, &
685  'on+dong', this%bc_labels)
686  call this%bc_dong%finalize()
687 
688 
689  call this%bc_dong%init(this%c_Xh, params)
690 
691  call bc_list_add(this%bclst_prs, this%bc_dong)
692 
693  ! Pressure solver
694  if (kspp_init) then
695  call neko_log%section("Pressure solver")
696 
697  call json_get_or_default(params, &
698  'case.fluid.pressure_solver.max_iterations', &
699  integer_val, ksp_max_iter)
700  call json_get(params, 'case.fluid.pressure_solver.type', solver_type)
701  call json_get(params, 'case.fluid.pressure_solver.preconditioner', &
702  precon_type)
703  call json_get(params, 'case.fluid.pressure_solver.absolute_tolerance', &
704  abs_tol)
705  call json_get_or_default(params, &
706  'case.fluid.pressure_solver.monitor', &
707  logical_val, .false.)
708  call neko_log%message('Type : ('// trim(solver_type) // &
709  ', ' // trim(precon_type) // ')')
710  write(log_buf, '(A,ES13.6)') 'Abs tol :', abs_tol
711  call neko_log%message(log_buf)
712 
713  call fluid_scheme_solver_factory(this%ksp_prs, this%dm_Xh%size(), &
714  solver_type, integer_val, abs_tol, logical_val)
715  call fluid_scheme_precon_factory(this%pc_prs, this%ksp_prs, &
716  this%c_Xh, this%dm_Xh, this%gs_Xh, this%bclst_prs, precon_type)
717 
718  call neko_log%end_section()
719 
720  end if
721 
722 
723  call neko_log%end_section()
724 
725  end subroutine fluid_scheme_init_all
726 
728  subroutine fluid_scheme_free(this)
729  class(fluid_scheme_t), intent(inout) :: this
730 
731  call this%bdry%free()
732 
733  if (allocated(this%bc_inflow)) then
734  call this%bc_inflow%free()
735  end if
736 
737  call this%bc_wall%free()
738  call this%bc_sym%free()
739 
740  !
741  ! Free everything related to field_dirichlet BCs
742  !
743  call this%user_field_bc_prs%field_bc%free()
744  call this%user_field_bc_prs%free()
745  call this%user_field_bc_vel%bc_u%field_bc%free()
746  call this%user_field_bc_vel%bc_v%field_bc%free()
747  call this%user_field_bc_vel%bc_w%field_bc%free()
748  call this%user_field_bc_vel%free()
749 
750  call this%Xh%free()
751 
752  if (allocated(this%ksp_vel)) then
753  call krylov_solver_destroy(this%ksp_vel)
754  deallocate(this%ksp_vel)
755  end if
756 
757  if (allocated(this%ksp_prs)) then
758  call krylov_solver_destroy(this%ksp_prs)
759  deallocate(this%ksp_prs)
760  end if
761 
762  if (allocated(this%pc_vel)) then
763  call precon_destroy(this%pc_vel)
764  deallocate(this%pc_vel)
765  end if
766 
767  if (allocated(this%pc_prs)) then
768  call precon_destroy(this%pc_prs)
769  deallocate(this%pc_prs)
770  end if
771 
772  if (allocated(this%bc_labels)) then
773  deallocate(this%bc_labels)
774  end if
775 
776  call this%source_term%free()
777 
778  call this%gs_Xh%free()
779 
780  call this%c_Xh%free()
781 
782  call bc_list_free(this%bclst_vel)
783 
784  call this%scratch%free()
785 
786  nullify(this%params)
787 
788  nullify(this%u)
789  nullify(this%v)
790  nullify(this%w)
791  nullify(this%p)
792 
793  call this%ulag%free()
794  call this%vlag%free()
795  call this%wlag%free()
796 
797 
798  if (associated(this%f_x)) then
799  call this%f_x%free()
800  end if
801 
802  if (associated(this%f_y)) then
803  call this%f_y%free()
804  end if
805 
806  if (associated(this%f_z)) then
807  call this%f_z%free()
808  end if
809 
810  nullify(this%f_x)
811  nullify(this%f_y)
812  nullify(this%f_z)
813 
814  call this%mu_field%free()
815 
816 
817  end subroutine fluid_scheme_free
818 
821  subroutine fluid_scheme_validate(this)
822  class(fluid_scheme_t), target, intent(inout) :: this
823  ! Variables for retrieving json parameters
824  logical :: logical_val
825 
826  if ( (.not. associated(this%u)) .or. &
827  (.not. associated(this%v)) .or. &
828  (.not. associated(this%w)) .or. &
829  (.not. associated(this%p))) then
830  call neko_error('Fields are not registered')
831  end if
832 
833  if ( (.not. allocated(this%u%x)) .or. &
834  (.not. allocated(this%v%x)) .or. &
835  (.not. allocated(this%w%x)) .or. &
836  (.not. allocated(this%p%x))) then
837  call neko_error('Fields are not allocated')
838  end if
839 
840  if (.not. allocated(this%ksp_vel)) then
841  call neko_error('No Krylov solver for velocity defined')
842  end if
843 
844  if (.not. allocated(this%ksp_prs)) then
845  call neko_error('No Krylov solver for pressure defined')
846  end if
847 
848  if (.not. associated(this%params)) then
849  call neko_error('No parameters defined')
850  end if
851 
852  if (allocated(this%bc_inflow)) then
853  select type (ip => this%bc_inflow)
854  type is (usr_inflow_t)
855  call ip%validate
856  end select
857  end if
858 
859  !
860  ! Setup checkpoint structure (if everything is fine)
861  !
862  call this%chkp%init(this%u, this%v, this%w, this%p)
863 
864  !
865  ! Setup mean flow fields if requested
866  !
867  if (this%params%valid_path('case.statistics')) then
868  call json_get_or_default(this%params, 'case.statistics.enabled', &
869  logical_val, .true.)
870  if (logical_val) then
871  call this%mean%init(this%u, this%v, this%w, this%p)
872  call this%stats%init(this%c_Xh, this%mean%u, &
873  this%mean%v, this%mean%w, this%mean%p)
874  end if
875  end if
876 
877  end subroutine fluid_scheme_validate
878 
883  subroutine fluid_scheme_bc_apply_vel(this, t, tstep)
884  class(fluid_scheme_t), intent(inout) :: this
885  real(kind=rp), intent(in) :: t
886  integer, intent(in) :: tstep
887 
888  call bc_list_apply_vector(this%bclst_vel, &
889  this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), t, tstep)
890 
891  end subroutine fluid_scheme_bc_apply_vel
892 
895  subroutine fluid_scheme_bc_apply_prs(this, t, tstep)
896  class(fluid_scheme_t), intent(inout) :: this
897  real(kind=rp), intent(in) :: t
898  integer, intent(in) :: tstep
899 
900  call bc_list_apply_scalar(this%bclst_prs, this%p%x, &
901  this%p%dof%size(), t, tstep)
902 
903  end subroutine fluid_scheme_bc_apply_prs
904 
907  subroutine fluid_scheme_solver_factory(ksp, n, solver, &
908  max_iter, abstol, monitor)
909  class(ksp_t), allocatable, target, intent(inout) :: ksp
910  integer, intent(in), value :: n
911  character(len=*), intent(in) :: solver
912  integer, intent(in) :: max_iter
913  real(kind=rp), intent(in) :: abstol
914  logical, intent(in) :: monitor
915 
916  call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
917  monitor = monitor)
918 
919  end subroutine fluid_scheme_solver_factory
920 
922  subroutine fluid_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
923  pctype)
924  class(pc_t), allocatable, target, intent(inout) :: pc
925  class(ksp_t), target, intent(inout) :: ksp
926  type(coef_t), target, intent(inout) :: coef
927  type(dofmap_t), target, intent(inout) :: dof
928  type(gs_t), target, intent(inout) :: gs
929  type(bc_list_t), target, intent(inout) :: bclst
930  character(len=*) :: pctype
931 
932  call precon_factory(pc, pctype)
933 
934  select type (pcp => pc)
935  type is (jacobi_t)
936  call pcp%init(coef, dof, gs)
937  type is (sx_jacobi_t)
938  call pcp%init(coef, dof, gs)
939  type is (device_jacobi_t)
940  call pcp%init(coef, dof, gs)
941  type is (hsmg_t)
942  if (len_trim(pctype) .gt. 4) then
943  if (index(pctype, '+') .eq. 5) then
944  call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst, &
945  trim(pctype(6:)))
946  else
947  call neko_error('Unknown coarse grid solver')
948  end if
949  else
950  call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
951  end if
952  end select
953 
954  call ksp%set_pc(pc)
955 
956  end subroutine fluid_scheme_precon_factory
957 
959  subroutine fluid_scheme_set_usr_inflow(this, usr_eval)
960  class(fluid_scheme_t), intent(inout) :: this
961  procedure(usr_inflow_eval) :: usr_eval
962 
963  select type (bc_if => this%bc_inflow)
964  type is (usr_inflow_t)
965  call bc_if%set_eval(usr_eval)
966  class default
967  call neko_error("Not a user defined inflow condition")
968  end select
969  end subroutine fluid_scheme_set_usr_inflow
970 
972  function fluid_compute_cfl(this, dt) result(c)
973  class(fluid_scheme_t), intent(in) :: this
974  real(kind=rp), intent(in) :: dt
975  real(kind=rp) :: c
976 
977  c = cfl(dt, this%u%x, this%v%x, this%w%x, &
978  this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
979 
980  end function fluid_compute_cfl
981 
984  subroutine fluid_scheme_set_bc_type_output(this, params)
985  class(fluid_scheme_t), target, intent(inout) :: this
986  type(json_file), intent(inout) :: params
987  type(dirichlet_t) :: bdry_mask
988  logical :: found
989 
990  !
991  ! Check if we need to output boundaries
992  !
993  call json_get_or_default(params, 'case.output_boundary', found, .false.)
994 
995  if (found) then
996  call this%bdry%init(this%dm_Xh, 'bdry')
997  this%bdry = 0.0_rp
998 
999  call bdry_mask%init_base(this%c_Xh)
1000  call bdry_mask%mark_zone(this%msh%wall)
1001  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1002  'w', this%bc_labels)
1003  call bdry_mask%finalize()
1004  call bdry_mask%set_g(1.0_rp)
1005  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1006  call bdry_mask%free()
1007 
1008  call bdry_mask%init_base(this%c_Xh)
1009  call bdry_mask%mark_zone(this%msh%inlet)
1010  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1011  'v', this%bc_labels)
1012  call bdry_mask%finalize()
1013  call bdry_mask%set_g(2.0_rp)
1014  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1015  call bdry_mask%free()
1016 
1017  call bdry_mask%init_base(this%c_Xh)
1018  call bdry_mask%mark_zone(this%msh%outlet)
1019  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1020  'o', this%bc_labels)
1021  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1022  'o+dong', this%bc_labels)
1023  call bdry_mask%finalize()
1024  call bdry_mask%set_g(3.0_rp)
1025  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1026  call bdry_mask%free()
1027 
1028  call bdry_mask%init_base(this%c_Xh)
1029  call bdry_mask%mark_zone(this%msh%sympln)
1030  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1031  'sym', this%bc_labels)
1032  call bdry_mask%finalize()
1033  call bdry_mask%set_g(4.0_rp)
1034  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1035  call bdry_mask%free()
1036 
1037  call bdry_mask%init_base(this%c_Xh)
1038  call bdry_mask%mark_zone(this%msh%outlet_normal)
1039  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1040  'on', this%bc_labels)
1041  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1042  'on+dong', this%bc_labels)
1043  call bdry_mask%finalize()
1044  call bdry_mask%set_g(5.0_rp)
1045  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1046  call bdry_mask%free()
1047 
1048  call bdry_mask%init_base(this%c_Xh)
1049  call bdry_mask%mark_zone(this%msh%periodic)
1050  call bdry_mask%finalize()
1051  call bdry_mask%set_g(6.0_rp)
1052  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1053  call bdry_mask%free()
1054 
1055  call bdry_mask%init_base(this%c_Xh)
1056  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1057  'd_vel_u', this%bc_labels)
1058  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1059  'd_vel_v', this%bc_labels)
1060  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1061  'd_vel_w', this%bc_labels)
1062  call bdry_mask%finalize()
1063  call bdry_mask%set_g(7.0_rp)
1064  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1065  call bdry_mask%free()
1066 
1067  call bdry_mask%init_base(this%c_Xh)
1068  call bdry_mask%mark_zones_from_list(this%msh%labeled_zones, &
1069  'd_pres', this%bc_labels)
1070  call bdry_mask%finalize()
1071  call bdry_mask%set_g(8.0_rp)
1072  call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
1073  call bdry_mask%free()
1074 
1075  end if
1076 
1077  end subroutine fluid_scheme_set_bc_type_output
1078 
1080  subroutine fluid_scheme_update_material_properties(this)
1081  class(fluid_scheme_t), intent(inout) :: this
1082  type(field_t), pointer :: nut
1083  integer :: n
1084 
1085  if (this%variable_material_properties) then
1086  nut => neko_field_registry%get_field(this%nut_field_name)
1087  n = nut%size()
1088 
1089  if (neko_bcknd_device .eq. 1) then
1090  call device_cfill(this%mu_field%x_d, this%mu, n)
1091  call device_add2s2(this%mu_field%x_d, nut%x_d, this%rho, n)
1092  else
1093  call cfill(this%mu_field%x, this%mu, n)
1094  call add2s2(this%mu_field%x, nut%x, this%rho, n)
1095  end if
1096  end if
1097 
1098  end subroutine fluid_scheme_update_material_properties
1099 
1100 end module fluid_scheme
Abstract interface to dealocate a fluid formulation.
Abstract interface to initialize a fluid formulation.
Abstract interface to restart a fluid scheme.
Abstract interface to compute a time-step.
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:53
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:44
Abstract interface defining a user defined inflow condition (pointwise)
Definition: usr_inflow.f90:80
Defines a boundary condition.
Definition: bc.f90:34
subroutine, public bc_list_add(bclst, bc)
Add a condition to a list of boundary conditions.
Definition: bc.f90:501
subroutine, public bc_list_apply_vector(bclst, x, y, z, n, t, tstep)
Apply a list of boundary conditions to a vector field.
Definition: bc.f90:583
subroutine, public bc_list_init(bclst, size)
Constructor for a list of boundary conditions.
Definition: bc.f90:460
subroutine, public bc_list_free(bclst)
Destructor for a list of boundary conditions.
Definition: bc.f90:487
subroutine, public bc_list_apply_scalar(bclst, x, n, t, tstep)
Apply a list of boundary conditions to a scalar field.
Definition: bc.f90:526
Defines a Blasius profile dirichlet condition.
Definition: blasius.f90:34
Defines a checkpoint.
Definition: checkpoint.f90:34
Coefficients.
Definition: coef.f90:34
Definition: comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
Jacobi preconditioner accelerator backend.
subroutine, public device_add2s2(a_d, b_d, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_cfill(a_d, c, n)
Set all elements to a constant c .
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Defines a dong outflow condition.
Defines inflow dirichlet conditions.
Defines inflow dirichlet conditions.
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
Definition: field_math.f90:186
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Stores a series fields.
Defines a field.
Definition: field.f90:34
Fluid formulations.
subroutine fluid_scheme_init_all(this, msh, lx, params, kspv_init, kspp_init, scheme, user, material_properties)
Initialize all components of the current scheme.
subroutine fluid_scheme_free(this)
Deallocate a fluid formulation.
subroutine fluid_scheme_update_material_properties(this)
Update the values of mu_field if necessary.
subroutine fluid_scheme_set_bc_type_output(this, params)
Set boundary types for the diagnostic output.
real(kind=rp) function fluid_compute_cfl(this, dt)
Compute CFL.
subroutine fluid_scheme_set_usr_inflow(this, usr_eval)
Initialize a user defined inflow condition.
subroutine fluid_scheme_bc_apply_vel(this, t, tstep)
Apply all boundary conditions defined for velocity Here we perform additional gs operations to take c...
subroutine fluid_scheme_validate(this)
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
subroutine fluid_scheme_init_common(this, msh, lx, params, scheme, user, material_properties, kspv_init)
Initialise a fluid scheme.
subroutine fluid_scheme_bc_apply_prs(this, t, tstep)
Apply all boundary conditions defined for pressure.
Implements the fluid_source_term_t type.
Computes various statistics for the fluid fields. We use the Reynolds decomposition for a field u = ...
Definition: fluid_stats.f90:36
Gather-scatter.
Krylov preconditioner.
Definition: pc_hsmg.f90:61
Defines inflow dirichlet conditions.
Definition: inflow.f90:34
Jacobi preconditioner.
Definition: pc_jacobi.f90:34
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
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
type(log_t), public neko_log
Global log stream.
Definition: log.f90:61
integer, parameter, public log_size
Definition: log.f90:40
Implements material_properties_t type.
Definition: math.f90:60
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition: math.f90:314
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition: math.f90:633
Defines a mean flow field.
Definition: mean_flow.f90:34
Defines a mean squared flow field.
Defines a mesh.
Definition: mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition: mesh.f90:56
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition: mesh.f90:58
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
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)
Definition: operators.f90:392
Krylov preconditioner.
Definition: precon.f90:34
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Definition: source_term.f90:34
Defines a function space.
Definition: space.f90:34
integer, parameter, public gll
Definition: space.f90:48
Defines a container for all statistics.
Definition: statistics.f90:34
Jacobi preconditioner SX-Aurora backend.
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition: symmetry.f90:34
Compound scheme for the advection and diffusion operators in a transport equation.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition: user_intf.f90:34
Defines inflow dirichlet conditions.
Definition: usr_inflow.f90:34
Utilities.
Definition: utils.f90:35
Defines wall boundary conditions.
Definition: wall.f90:34
A list of boundary conditions.
Definition: bc.f90:104
Base type for a boundary condition.
Definition: bc.f90:51
Blasius profile for inlet (vector valued)
Definition: blasius.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44
Dong outflow condition Follows "A Convective-like Energy-Stable Open Boundary Condition for Simulati...
User defined dirichlet condition, for which the user can work with an entire field....
Extension of the user defined dirichlet condition field_dirichlet
Base type of all fluid formulations.
Wrapper contaning and executing the fluid source terms.
Dirichlet condition for inlet (vector valued)
Definition: inflow.f90:43
Defines a jacobi preconditioner.
Definition: pc_jacobi.f90:45
Base abstract type for a canonical Krylov method, solving .
Definition: krylov.f90:66
Contains all the material properties necessary in the simulation.
Defines a canonical Krylov preconditioner.
Definition: precon.f90:40
The function space for the SEM solution fields.
Definition: space.f90:62
Defines a jacobi preconditioner for SX-Aurora.
Mixed Dirichlet-Neumann symmetry plane condition.
Definition: symmetry.f90:46
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
Provides a tool to set time step dt.
User defined dirichlet condition for inlet (vector valued)
Definition: usr_inflow.f90:46
No-slip Wall boundary condition.
Definition: wall.f90:43