Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
euler_bc_fctry.f90
Go to the documentation of this file.
1! Copyright (c) 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!
33!
34submodule(fluid_scheme_compressible_euler) euler_bc_fctry
35 use dirichlet, only : dirichlet_t
36 use inflow, only : inflow_t
38 use symmetry, only : symmetry_t
40 implicit none
41
42 ! List of all possible types created by the boundary condition factories
43 character(len=25) :: EULER_KNOWN_BCS(7) = [character(len=25) :: &
44 "velocity_value", &
45 "density_value", &
46 "pressure_value", &
47 "no_slip", &
48 "symmetry", &
49 "outflow", &
50 "normal_outflow"]
51
52contains
59 module subroutine density_bc_factory(object, scheme, json, coef, user)
60 class(bc_t), pointer, intent(inout) :: object
61 type(fluid_scheme_compressible_euler_t), intent(in) :: scheme
62 type(json_file), intent(inout) :: json
63 type(coef_t), intent(in) :: coef
64 type(user_t), intent(in) :: user
65 character(len=:), allocatable :: type
66 integer :: i, j, k
67 integer, allocatable :: zone_indices(:)
68
69 call json_get(json, "type", type)
70
71 select case (trim(type))
72 case ("density_value")
73 allocate(dirichlet_t::object)
74 case default
75 do i = 1, size(euler_known_bcs)
76 if (trim(type) .eq. trim(euler_known_bcs(i))) return
77 end do
78 call neko_type_error("compressible_euler boundary conditions", type, &
79 EULER_KNOWN_BCS)
80 end select
81
82 call json_get_or_lookup(json, "zone_indices", zone_indices)
83 call object%init(coef, json)
84
85 do i = 1, size(zone_indices)
86 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
87 end do
88 call object%finalize()
89 end subroutine density_bc_factory
90
97 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
98 class(bc_t), pointer, intent(inout) :: object
99 type(fluid_scheme_compressible_euler_t), intent(inout) :: scheme
100 type(json_file), intent(inout) :: json
101 type(coef_t), intent(in) :: coef
102 type(user_t), intent(in) :: user
103 character(len=:), allocatable :: type
104 integer :: i, j, k
105 integer, allocatable :: zone_indices(:)
106
107 call json_get(json, "type", type)
108
109 select case (trim(type))
110 case ("outflow", "normal_outflow")
111 allocate(zero_dirichlet_t::object)
112 case ("pressure_value")
113 allocate(dirichlet_t::object)
114 case default
115 do i = 1, size(euler_known_bcs)
116 if (trim(type) .eq. trim(euler_known_bcs(i))) return
117 end do
118 call neko_type_error("compressible_euler boundary conditions", type, &
119 EULER_KNOWN_BCS)
120 end select
121
122 call json_get_or_lookup(json, "zone_indices", zone_indices)
123 call object%init(coef, json)
124
125 do i = 1, size(zone_indices)
126 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
127 end do
128 call object%finalize()
129
130 ! All pressure bcs are currently strong, so for all of them we
131 ! mark with value 1 in the mesh
132 do i = 1, size(zone_indices)
133 do j = 1, scheme%msh%nelv
134 do k = 1, 2 * scheme%msh%gdim
135 if (scheme%msh%facet_type(k,j) .eq. -zone_indices(i)) then
136 scheme%msh%facet_type(k, j) = 1
137 end if
138 end do
139 end do
140 end do
141 end subroutine pressure_bc_factory
142
149 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
150 class(bc_t), pointer, intent(inout) :: object
151 type(fluid_scheme_compressible_euler_t), intent(in) :: scheme
152 type(json_file), intent(inout) :: json
153 type(coef_t), intent(in) :: coef
154 type(user_t), intent(in) :: user
155 character(len=:), allocatable :: type
156 integer :: i, j, k
157 integer, allocatable :: zone_indices(:)
158
159 call json_get(json, "type", type)
160
161 select case (trim(type))
162 case ("symmetry")
163 allocate(symmetry_t::object)
164 case ("no_slip")
165 allocate(zero_dirichlet_t::object)
166 case ("velocity_value")
167 allocate(inflow_t::object)
168 case default
169 do i = 1, size(euler_known_bcs)
170 if (trim(type) .eq. trim(euler_known_bcs(i))) return
171 end do
172 call neko_type_error("compressible_euler boundary conditions", type, &
173 EULER_KNOWN_BCS)
174 end select
175
176 call json_get_or_lookup(json, "zone_indices", zone_indices)
177 call object%init(coef, json)
178 do i = 1, size(zone_indices)
179 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
180 end do
181 call object%finalize()
182
183 end subroutine velocity_bc_factory
184
185end submodule euler_bc_fctry
Retrieves a parameter by name or throws an error.
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines inflow dirichlet conditions.
Definition inflow.f90:34
Utilities for retrieving parameters from the case files.
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition symmetry.f90:34
Defines a zero-valued Dirichlet boundary condition.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
Dirichlet condition for inlet (vector valued)
Definition inflow.f90:47
Mixed Dirichlet-Neumann symmetry plane condition.
Definition symmetry.f90:50
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...