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
42 implicit none
43
44 ! List of all possible types created by the boundary condition factories
45 character(len=25) :: SCALAR_PNPN_KNOWN_BCS(4) = [character(len=25) :: &
46 "dirichlet", &
47 "user_dirichlet", &
48 "user_neumann", &
49 "neumann"]
50
51contains
52
60 module subroutine bc_factory(object, scheme, json, coef, user)
61 class(bc_t), pointer, intent(inout) :: object
62 type(scalar_pnpn_t), intent(in) :: scheme
63 type(json_file), intent(inout) :: json
64 type(coef_t), intent(in) :: coef
65 type(user_t), intent(in) :: user
66 character(len=:), allocatable :: type
67 integer :: i
68 integer, allocatable :: zone_indices(:)
69 character(len=:), allocatable :: default_name
70 character(len=64) :: buf
71
72 if (associated(object)) then
73 call object%free()
74 nullify(object)
75 end if
76
77 call json_get(json, "type", type)
78
79 select case (trim(type))
80 case ("user_dirichlet")
81 allocate(field_dirichlet_t::object)
82 select type (obj => object)
83 type is (field_dirichlet_t)
84 obj%update => user%dirichlet_conditions
85 ! Add the name of the dummy field in the bc, matching the scalar
86 ! solved for.
87 call json%add("field_name", scheme%s%name)
88 end select
89 case ("dirichlet")
90 allocate(dirichlet_t::object)
91 case ("user_neumann")
92 allocate(field_neumann_t::object)
93 select type (obj => object)
94 type is (field_neumann_t)
95 obj%update => user%neumann_conditions
96 ! Add the name of the dummy field in the bc, matching the scalar
97 ! solved for.
98 call json%add("field_name", scheme%s%name)
99 end select
100 case ("neumann")
101 allocate(neumann_t::object)
102 case default
103 call neko_type_error("scalar_pnpn boundary conditions", type, &
104 SCALAR_PNPN_KNOWN_BCS)
105 end select
106
107 call json_get(json, "zone_indices", zone_indices)
108 call object%init(coef, json)
109 do i = 1, size(zone_indices)
110 call object%mark_zone(coef%msh%labeled_zones(zone_indices(i)))
111 end do
112
113 write(buf,'("scalar_bc_",I0)') zone_indices(1)
114 default_name = trim(buf)
115 call json_get_or_default(json, "name", object%name, default_name)
116 object%zone_indices = zone_indices
117 call object%finalize()
118
119 end subroutine bc_factory
120
121
122end 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
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:313
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
A type collecting all the overridable user routines and flag to suppress type injection from custom m...