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
50 character(len=25) :: ADJOINT_FLUID_PNPN_KNOWN_BCS(13) = &
51 [character(len=25) :: &
58 "normal_outflow+dong", &
63 "user_velocity_pointwise", &
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
82 integer,
allocatable :: zone_indices(:)
84 call json_get(json,
"type", type)
86 select case (trim(type))
87 case (
"outflow",
"normal_outflow")
88 allocate(zero_dirichlet_t::object)
90 case (
"outflow+dong",
"normal_outflow+dong")
91 allocate(dong_outflow_t::object)
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)
102 do i = 1,
size(adjoint_fluid_pnpn_known_bcs)
103 if (trim(type) .eq. trim(adjoint_fluid_pnpn_known_bcs(i)))
return
105 call neko_type_error(
"adjoint_fluid_pnpn boundary conditions",
type, &
106 ADJOINT_FLUID_PNPN_KNOWN_BCS)
109 call json_get(json,
"zone_indices", zone_indices)
110 call object%init(coef, json)
112 do i = 1,
size(zone_indices)
113 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
115 call object%finalize()
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
128 end subroutine pressure_bc_factory
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
144 integer,
allocatable :: zone_indices(:)
146 call json_get(json,
"type", type)
148 select case (trim(type))
150 allocate(symmetry_t::object)
151 case (
"velocity_value")
152 allocate(inflow_t::object)
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)
162 allocate(wall_model_bc_t::object)
166 call json%add(
"nu", scheme%mu / scheme%rho)
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
184 do i = 1,
size(adjoint_fluid_pnpn_known_bcs)
185 if (trim(type) .eq. trim(adjoint_fluid_pnpn_known_bcs(i)))
return
187 call neko_type_error(
"adjoint_fluid_pnpn boundary conditions",
type, &
188 ADJOINT_FLUID_PNPN_KNOWN_BCS)
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)))
196 call object%finalize()
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
211 end subroutine velocity_bc_factory
213end submodule adjoint_fluid_pnpn_bc_fctry
Adjoint Pn/Pn formulation.
User defined user region.