Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_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(scalar_pnpn) scalar_pnpn_bc_fctry
36 use dirichlet, only : dirichlet_t
37 use neumann, only : neumann_t
38 use user_intf, only : user_t
39 use utils, only : neko_type_error
43 implicit none
44
45 ! List of all possible types created by the boundary condition factories
46 character(len=25) :: SCALAR_PNPN_KNOWN_BCS(5) = [character(len=25) :: &
47 "dirichlet", &
48 "user_dirichlet", &
49 "user_neumann", &
50 "neumann", &
51 "overset_interface"]
52
53contains
54
62 module subroutine bc_factory(object, scheme, json, coef, user)
63 class(bc_t), pointer, intent(inout) :: object
64 type(scalar_pnpn_t), intent(in) :: scheme
65 type(json_file), intent(inout) :: json
66 type(coef_t), intent(in) :: coef
67 type(user_t), intent(in) :: user
68 character(len=:), allocatable :: type
69 integer :: i
70 integer, allocatable :: zone_indices(:)
71 character(len=:), allocatable :: default_name
72 character(len=64) :: buf
73
74 if (associated(object)) then
75 call object%free()
76 nullify(object)
77 end if
78
79 call json_get(json, "type", type)
80
81 select case (trim(type))
82 case ("user_dirichlet")
83 allocate(field_dirichlet_t::object)
84 select type (obj => object)
85 type is (field_dirichlet_t)
86 obj%update => user%dirichlet_conditions
87 ! Add the name of the dummy field in the bc, matching the scalar
88 ! solved for.
89 call json%add("field_name", scheme%s%name)
90 end select
91 case ("dirichlet")
92 allocate(dirichlet_t::object)
93 case ("user_neumann")
94 allocate(field_neumann_t::object)
95 select type (obj => object)
96 type is (field_neumann_t)
97 obj%update => user%neumann_conditions
98 ! Add the name of the dummy field in the bc, matching the scalar
99 ! solved for.
100 call json%add("field_name", scheme%s%name)
101 end select
102 case ("neumann")
103 allocate(neumann_t::object)
104 case ("overset_interface")
105 allocate(overset_interface_t::object)
106 select type (obj => object)
107 type is (overset_interface_t)
108 call json%add("field_name", scheme%s%name)
109 end select
110 case default
111 call neko_type_error("scalar_pnpn boundary conditions", type, &
112 SCALAR_PNPN_KNOWN_BCS)
113 end select
114
115 call json_get(json, "zone_indices", zone_indices)
116 call object%init(coef, json)
117 do i = 1, size(zone_indices)
118 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
119 end do
120
121 write(buf,'("scalar_bc_",I0)') zone_indices(1)
122 default_name = trim(buf)
123 call json_get_or_default(json, "name", object%name, default_name)
124 object%zone_indices = zone_indices
125 call object%finalize()
126
127 end subroutine bc_factory
128
129
130end submodule scalar_pnpn_bc_fctry
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines user dirichlet condition for a scalar field.
Defines user neumann condition for a scalar field.
Defines a Neumann boundary condition.
Definition neumann.f90:34
Defines overset interface scalar boundary conditions.
Contains the scalar_pnpn_t type.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
subroutine, public neko_type_error(base_type, wrong_type, known_types)
Reports an error allocating a type for a particular base pointer class.
Definition utils.f90:359
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
User defined dirichlet condition, for which the user can work with an entire field....
User defined neumann condition, for which the user can work with an entire field. The type stores a s...
A Neumann boundary condition. Sets the flux of the field to the chosen values.
Definition neumann.f90:60
Overset interface BC for a scalar field.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...