Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_scalar_pnpn_bc_fctry.f90
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!
38submodule(adjoint_scalar_pnpn) adjoint_scalar_pnpn_bc_fctry
39 use dirichlet, only : dirichlet_t
40 use neumann, only : neumann_t
41 use usr_scalar, only : usr_scalar_t
42 use user_intf, only : user_t
43 use usr_scalar, only : usr_scalar_t
44 use utils, only : neko_type_error
45 use field_dirichlet, only : field_dirichlet_t
46 implicit none
47
48 ! List of all possible types created by the boundary condition factories
49 character(len=25) :: ADJOINT_SCALAR_KNOWN_BCS(4) = [character(len=25) :: &
50 "dirichlet", &
51 "user_pointwise", &
52 "user", &
53 "neumann"]
54
55contains
56
64 module subroutine adjoint_bc_factory(object, scheme, json, coef, user)
65 class(bc_t), pointer, intent(inout) :: object
66 type(adjoint_scalar_pnpn_t), intent(in) :: scheme
67 type(json_file), intent(inout) :: json
68 type(coef_t), intent(in) :: coef
69 type(user_t), intent(in) :: user
70 character(len=:), allocatable :: type
71 integer :: i
72 integer, allocatable :: zone_indices(:)
73
74 call json_get(json, "type", type)
75
76 select case (trim(type))
77 case ("user_pointwise")
78 allocate(usr_scalar_t::object)
79 select type (obj => object)
80 type is (usr_scalar_t)
81 call obj%set_eval(user%scalar_user_bc)
82 end select
83 case ("user")
84 allocate(field_dirichlet_t::object)
85 select type (obj => object)
86 type is (field_dirichlet_t)
87 obj%update => user%user_dirichlet_update
88 ! Add the name of the dummy field in the bc, matching the scalar
89 ! solved for.
90 call json%add("field_name", scheme%s_adj%name)
91 end select
92 case ("dirichlet")
93 allocate(dirichlet_t::object)
94 case ("neumann")
95 allocate(neumann_t::object)
96 case default
97 call neko_type_error("scalar_pnpn boundary conditions", type, &
98 ADJOINT_SCALAR_KNOWN_BCS)
99 end select
100
101 call json_get(json, "zone_indices", zone_indices)
102 call object%init(coef, json)
103 do i = 1, size(zone_indices)
104 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
105 end do
106 call object%finalize()
107
108 end subroutine adjoint_bc_factory
109
110
111end submodule adjoint_scalar_pnpn_bc_fctry
Contains the adjoint_scalar_pnpn_t type.