Neko 1.99.1
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
39 implicit none
40
41 ! List of all possible types created by the boundary condition factories
42 character(len=25) :: EULER_KNOWN_BCS(7) = [character(len=25) :: &
43 "velocity_value", &
44 "density_value", &
45 "pressure_value", &
46 "no_slip", &
47 "symmetry", &
48 "outflow", &
49 "normal_outflow"]
50
51contains
58 module subroutine density_bc_factory(object, scheme, json, coef, user)
59 class(bc_t), pointer, intent(inout) :: object
60 type(fluid_scheme_compressible_euler_t), intent(in) :: scheme
61 type(json_file), intent(inout) :: json
62 type(coef_t), intent(in) :: coef
63 type(user_t), intent(in) :: user
64 character(len=:), allocatable :: type
65 integer :: i, j, k
66 integer, allocatable :: zone_indices(:)
67
68 call json_get(json, "type", type)
69
70 select case (trim(type))
71 case ("density_value")
72 allocate(dirichlet_t::object)
73 case default
74 do i = 1, size(euler_known_bcs)
75 if (trim(type) .eq. trim(euler_known_bcs(i))) return
76 end do
77 call neko_type_error("compressible_euler boundary conditions", type, &
78 EULER_KNOWN_BCS)
79 end select
80
81 call json_get(json, "zone_indices", zone_indices)
82 call object%init(coef, json)
83
84 do i = 1, size(zone_indices)
85 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
86 end do
87 call object%finalize()
88 end subroutine density_bc_factory
89
96 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
97 class(bc_t), pointer, intent(inout) :: object
98 type(fluid_scheme_compressible_euler_t), intent(inout) :: scheme
99 type(json_file), intent(inout) :: json
100 type(coef_t), intent(in) :: coef
101 type(user_t), intent(in) :: user
102 character(len=:), allocatable :: type
103 integer :: i, j, k
104 integer, allocatable :: zone_indices(:)
105
106 call json_get(json, "type", type)
107
108 select case (trim(type))
109 case ("outflow", "normal_outflow")
110 allocate(zero_dirichlet_t::object)
111 case ("pressure_value")
112 allocate(dirichlet_t::object)
113 case default
114 do i = 1, size(euler_known_bcs)
115 if (trim(type) .eq. trim(euler_known_bcs(i))) return
116 end do
117 call neko_type_error("compressible_euler boundary conditions", type, &
118 EULER_KNOWN_BCS)
119 end select
120
121 call json_get(json, "zone_indices", zone_indices)
122 call object%init(coef, json)
123
124 do i = 1, size(zone_indices)
125 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
126 end do
127 call object%finalize()
128
129 ! All pressure bcs are currently strong, so for all of them we
130 ! mark with value 1 in the mesh
131 do i = 1, size(zone_indices)
132 do j = 1, scheme%msh%nelv
133 do k = 1, 2 * scheme%msh%gdim
134 if (scheme%msh%facet_type(k,j) .eq. -zone_indices(i)) then
135 scheme%msh%facet_type(k, j) = 1
136 end if
137 end do
138 end do
139 end do
140 end subroutine pressure_bc_factory
141
148 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
149 class(bc_t), pointer, intent(inout) :: object
150 type(fluid_scheme_compressible_euler_t), intent(in) :: scheme
151 type(json_file), intent(inout) :: json
152 type(coef_t), intent(in) :: coef
153 type(user_t), intent(in) :: user
154 character(len=:), allocatable :: type
155 integer :: i, j, k
156 integer, allocatable :: zone_indices(:)
157
158 call json_get(json, "type", type)
159
160 select case (trim(type))
161 case ("symmetry")
162 allocate(symmetry_t::object)
163 case ("no_slip")
164 allocate(zero_dirichlet_t::object)
165 case ("velocity_value")
166 allocate(inflow_t::object)
167 case default
168 do i = 1, size(euler_known_bcs)
169 if (trim(type) .eq. trim(euler_known_bcs(i))) return
170 end do
171 call neko_type_error("compressible_euler boundary conditions", type, &
172 EULER_KNOWN_BCS)
173 end select
174
175 call json_get(json, "zone_indices", zone_indices)
176 call object%init(coef, json)
177 do i = 1, size(zone_indices)
178 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
179 end do
180 call object%finalize()
181
182 end subroutine velocity_bc_factory
183
184end submodule euler_bc_fctry
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines inflow dirichlet conditions.
Definition inflow.f90:34
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:48
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...