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 character(len=:), allocatable :: default_name
69 character(len=64) :: buf
70
71 call json_get(json, "type", type)
72
73 select case (trim(type))
74 case ("density_value")
75 allocate(dirichlet_t::object)
76 case default
77 do i = 1, size(euler_known_bcs)
78 if (trim(type) .eq. trim(euler_known_bcs(i))) return
79 end do
80 call neko_type_error("compressible_euler boundary conditions", type, &
81 EULER_KNOWN_BCS)
82 end select
83
84 call json_get_or_lookup(json, "zone_indices", zone_indices)
85 call object%init(coef, json)
86
87 do i = 1, size(zone_indices)
88 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
89 end do
90
91 write(buf,'("density_bc_",I0)') zone_indices(1)
92 default_name = trim(buf)
93 call json_get_or_default(json, "name", object%name, default_name)
94 object%zone_indices = zone_indices
95 call object%finalize()
96 end subroutine density_bc_factory
97
104 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
105 class(bc_t), pointer, intent(inout) :: object
106 type(fluid_scheme_compressible_euler_t), intent(inout) :: scheme
107 type(json_file), intent(inout) :: json
108 type(coef_t), intent(in) :: coef
109 type(user_t), intent(in) :: user
110 character(len=:), allocatable :: type
111 integer :: i, j, k
112 integer, allocatable :: zone_indices(:)
113 character(len=:), allocatable :: default_name
114 character(len=64) :: buf
115
116 call json_get(json, "type", type)
117
118 select case (trim(type))
119 case ("outflow", "normal_outflow")
120 allocate(zero_dirichlet_t::object)
121 case ("pressure_value")
122 allocate(dirichlet_t::object)
123 case default
124 do i = 1, size(euler_known_bcs)
125 if (trim(type) .eq. trim(euler_known_bcs(i))) return
126 end do
127 call neko_type_error("compressible_euler boundary conditions", type, &
128 EULER_KNOWN_BCS)
129 end select
130
131 call json_get_or_lookup(json, "zone_indices", zone_indices)
132 call object%init(coef, json)
133
134 do i = 1, size(zone_indices)
135 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
136 end do
137
138 write(buf,'("pressure_bc_",I0)') zone_indices(1)
139 default_name = trim(buf)
140 call json_get_or_default(json, "name", object%name, default_name)
141 object%zone_indices = zone_indices
142 call object%finalize()
143
144 ! All pressure bcs are currently strong, so for all of them we
145 ! mark with value 1 in the mesh
146 do i = 1, size(zone_indices)
147 do j = 1, scheme%msh%nelv
148 do k = 1, 2 * scheme%msh%gdim
149 if (scheme%msh%facet_type(k,j) .eq. -zone_indices(i)) then
150 scheme%msh%facet_type(k, j) = 1
151 end if
152 end do
153 end do
154 end do
155 end subroutine pressure_bc_factory
156
163 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
164 class(bc_t), pointer, intent(inout) :: object
165 type(fluid_scheme_compressible_euler_t), intent(in) :: scheme
166 type(json_file), intent(inout) :: json
167 type(coef_t), intent(in) :: coef
168 type(user_t), intent(in) :: user
169 character(len=:), allocatable :: type
170 integer :: i, j, k
171 integer, allocatable :: zone_indices(:)
172 character(len=:), allocatable :: default_name
173 character(len=64) :: buf
174
175 call json_get(json, "type", type)
176
177 select case (trim(type))
178 case ("symmetry")
179 allocate(symmetry_t::object)
180 case ("no_slip")
181 allocate(zero_dirichlet_t::object)
182 case ("velocity_value")
183 allocate(inflow_t::object)
184 case default
185 do i = 1, size(euler_known_bcs)
186 if (trim(type) .eq. trim(euler_known_bcs(i))) return
187 end do
188 call neko_type_error("compressible_euler boundary conditions", type, &
189 EULER_KNOWN_BCS)
190 end select
191
192 call json_get_or_lookup(json, "zone_indices", zone_indices)
193 call object%init(coef, json)
194 do i = 1, size(zone_indices)
195 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
196 end do
197
198 write(buf,'("velocity_bc_",I0)') zone_indices(1)
199 default_name = trim(buf)
200 call json_get_or_default(json, "name", object%name, default_name)
201 object%zone_indices = zone_indices
202 call object%finalize()
203
204 end subroutine velocity_bc_factory
205
206end submodule euler_bc_fctry
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.
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...