Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_fluid_pnpn_bc_fctry.f90
Go to the documentation of this file.
1! Copyright (c) 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!
33!
35submodule(adjoint_fluid_pnpn) adjoint_fluid_pnpn_bc_fctry
36 use user_intf, only: user_t
37 use utils, only: neko_type_error
38 use field_dirichlet, only: field_dirichlet_t
39 use inflow, only: inflow_t
40 use usr_inflow, only: usr_inflow_t, usr_inflow_eval
41 use blasius, only: blasius_t
42 use dirichlet, only: dirichlet_t
43 use dong_outflow, only: dong_outflow_t
44 use symmetry, only: symmetry_t
45 use non_normal, only: non_normal_t
46 use field_dirichlet_vector, only: field_dirichlet_vector_t
47 implicit none
48
49 ! List of all possible types created by the boundary condition factories
50 character(len=25) :: ADJOINT_FLUID_PNPN_KNOWN_BCS(13) = &
51 [character(len=25) :: &
52 "symmetry", &
53 "velocity_value", &
54 "no_slip", &
55 "outflow", &
56 "normal_outflow", &
57 "outflow+dong", &
58 "normal_outflow+dong", &
59 "shear_stress", &
60 "user_velocity", &
61 "user_pressure", &
62 "blasius_profile", &
63 "user_velocity_pointwise", &
64 "wall_model"]
65
66contains
67
74 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
75 class(bc_t), pointer, intent(inout) :: object
76 type(adjoint_fluid_pnpn_t), intent(in) :: scheme
77 type(json_file), intent(inout) :: json
78 type(coef_t), intent(in) :: coef
79 type(user_t), intent(in) :: user
80 character(len=:), allocatable :: type
81 integer :: i, j, k
82 integer, allocatable :: zone_indices(:)
83
84 call json_get(json, "type", type)
85
86 select case (trim(type))
87 case ("outflow", "normal_outflow")
88 allocate(zero_dirichlet_t::object)
89
90 case ("outflow+dong", "normal_outflow+dong")
91 allocate(dong_outflow_t::object)
92
93 case ("user_pressure")
94 allocate(field_dirichlet_t::object)
95 select type (obj => object)
96 type is (field_dirichlet_t)
97 obj%update => user%user_dirichlet_update
98 call json%add("field_name", scheme%p_adj%name)
99 end select
100
101 case default
102 do i = 1, size(adjoint_fluid_pnpn_known_bcs)
103 if (trim(type) .eq. trim(adjoint_fluid_pnpn_known_bcs(i))) return
104 end do
105 call neko_type_error("adjoint_fluid_pnpn boundary conditions", type, &
106 ADJOINT_FLUID_PNPN_KNOWN_BCS)
107 end select
108
109 call json_get(json, "zone_indices", zone_indices)
110 call object%init(coef, json)
111
112 do i = 1, size(zone_indices)
113 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
114 end do
115 call object%finalize()
116
117 ! All pressure bcs are currently strong, so for all of them we
118 ! mark with value 1 in the mesh
119 do i = 1, size(zone_indices)
120 do j = 1, scheme%msh%nelv
121 do k = 1, 2 * scheme%msh%gdim
122 if (scheme%msh%facet_type(k,j) .eq. -zone_indices(i)) then
123 scheme%msh%facet_type(k, j) = 1
124 end if
125 end do
126 end do
127 end do
128 end subroutine pressure_bc_factory
129
136 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
137 class(bc_t), pointer, intent(inout) :: object
138 type(adjoint_fluid_pnpn_t), intent(in) :: scheme
139 type(json_file), intent(inout) :: json
140 type(coef_t), intent(in) :: coef
141 type(user_t), intent(in) :: user
142 character(len=:), allocatable :: type
143 integer :: i, j, k
144 integer, allocatable :: zone_indices(:)
145
146 call json_get(json, "type", type)
147
148 select case (trim(type))
149 case ("symmetry")
150 allocate(symmetry_t::object)
151 case ("velocity_value")
152 allocate(inflow_t::object)
153 case ("no_slip")
154 allocate(zero_dirichlet_t::object)
155 case ("normal_outflow", "normal_outflow+dong")
156 allocate(non_normal_t::object)
157 case ("blasius_profile")
158 allocate(blasius_t::object)
159 case ("shear_stress")
160 allocate(shear_stress_t::object)
161 case ("wall_model")
162 allocate(wall_model_bc_t::object)
163 ! Kind of hack, but maybe OK? The thing is, we need the nu for
164 ! initing the wall model, and forcing the user duplicate that there
165 ! would be a nightmare.
166 call json%add("nu", scheme%mu / scheme%rho)
167
168 case ("user_velocity")
169 allocate(field_dirichlet_vector_t::object)
170 select type (obj => object)
171 type is (field_dirichlet_vector_t)
172 obj%update => user%user_dirichlet_update
173 end select
174
175 ! case ("user_velocity_pointwise")
176 ! allocate(usr_inflow_t::object)
177 ! select type (obj => object)
178 ! type is (usr_inflow_t)
179 ! call obj%set_eval(user%adjoint_user_if)
180 ! call obj%validate()
181 ! end select
182
183 case default
184 do i = 1, size(adjoint_fluid_pnpn_known_bcs)
185 if (trim(type) .eq. trim(adjoint_fluid_pnpn_known_bcs(i))) return
186 end do
187 call neko_type_error("adjoint_fluid_pnpn boundary conditions", type, &
188 ADJOINT_FLUID_PNPN_KNOWN_BCS)
189 end select
190
191 call json_get(json, "zone_indices", zone_indices)
192 call object%init(coef, json)
193 do i = 1, size(zone_indices)
194 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
195 end do
196 call object%finalize()
197
198 ! Exclude these two because they are bcs for the residual, not velocity
199 if (trim(type) .ne. "normal_outflow" .and. &
200 trim(type) .ne. "normal_outflow+dong") then
201 do i = 1, size(zone_indices)
202 do j = 1, scheme%msh%nelv
203 do k = 1, 2 * scheme%msh%gdim
204 if (scheme%msh%facet_type(k,j) .eq. -zone_indices(i)) then
205 scheme%msh%facet_type(k, j) = 2
206 end if
207 end do
208 end do
209 end do
210 end if
211 end subroutine velocity_bc_factory
212
213end submodule adjoint_fluid_pnpn_bc_fctry
Adjoint Pn/Pn formulation.
User defined user region.
Definition user.f90:2