Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
wall_model_fctry.f90
Go to the documentation of this file.
1! Copyright (c) 2021-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
34submodule(wall_model) wall_model_fctry
35 use spalding, only : spalding_t
37 use most, only : most_t
38 use richardson, only : richardson_t
39 use utils, only : neko_type_error
41 implicit none
42
43 ! List of all possible types created by the factory routine
44 character(len=20) :: WALLM_KNOWN_TYPES(4) = [character(len=20) :: &
45 "spalding", &
46 "rough_log_law", &
47 "most", &
48 "richardson"]
49
50contains
51
60 module subroutine wall_model_factory(object, scheme_name, coef, msk, facet, &
61 json)
62 class(wall_model_t), allocatable, intent(inout) :: object
63 character(len=*), intent(in) :: scheme_name
64 type(coef_t), intent(in) :: coef
65 integer, intent(in) :: msk(:)
66 integer, intent(in) :: facet(:)
67 type(json_file), intent(inout) :: json
68 character(len=:), allocatable :: type_name
69 character(len=:), allocatable :: type_string
70 integer :: h_index
71
72 call json_get(json, "model", type_name)
73 call json_get(json, "h_index", h_index)
74
75 call wall_model_allocator(object, type_name)
76
77 ! Initialize
78 call object%init(scheme_name, coef, msk, facet, h_index, json)
79
80 end subroutine wall_model_factory
81
85 module subroutine wall_model_allocator(object, type_name)
86 class(wall_model_t), allocatable, intent(inout) :: object
87 character(len=:), allocatable, intent(in) :: type_name
88 integer :: i
89
90 if (allocated(object)) then
91 call object%free()
92 deallocate(object)
93 end if
94
95 select case (trim(type_name) )
96 case ("spalding")
97 allocate(spalding_t::object)
98 case ("rough_log_law")
99 allocate(rough_log_law_t::object)
100 case ("most")
101 allocate(most_t::object)
102 case ("richardson")
103 allocate(richardson_t::object)
104 case default
105 do i = 1, wall_model_registry_size
106 if (trim(type_name) .eq. trim(wall_model_registry(i)%type_name)) then
107 call wall_model_registry(i)%allocator(object)
108 return
109 end if
110 end do
111 call neko_type_error("wall model", trim(type_name), wallm_known_types)
112 end select
113 end subroutine wall_model_allocator
114
120 module subroutine register_wall_model(type_name, allocator)
121 character(len=*), intent(in) :: type_name
122 procedure(wall_model_allocate), pointer, intent(in) :: allocator
123 type(allocator_entry), allocatable :: temp(:)
124 integer :: i
125
126 do i = 1, size(wallm_known_types)
127 if (trim(type_name) .eq. trim(wallm_known_types(i))) then
128 call neko_type_registration_error("wall model", type_name, .true.)
129 end if
130 end do
131
132 do i = 1, wall_model_registry_size
133 if (trim(type_name) .eq. trim(wall_model_registry(i)%type_name)) then
134 call neko_type_registration_error("wall model", type_name, .false.)
135 end if
136 end do
137
138 ! Expand registry
139 if (wall_model_registry_size .eq. 0) then
140 allocate(wall_model_registry(1))
141 else
142 allocate(temp(wall_model_registry_size + 1))
143 temp(1:wall_model_registry_size) = wall_model_registry
144 call move_alloc(temp, wall_model_registry)
145 end if
146
147 wall_model_registry_size = wall_model_registry_size + 1
148 wall_model_registry(wall_model_registry_size)%type_name = type_name
149 wall_model_registry(wall_model_registry_size)%allocator => allocator
150 end subroutine register_wall_model
151
152end submodule wall_model_fctry
Implements most_t.
Definition most.f90:35
Implements richardson_t.
Implements rough_log_law_t.
Implements spalding_t.
Definition spalding.f90:35
Utilities.
Definition utils.f90:35
subroutine, public neko_type_registration_error(base_type, wrong_type, known)
Definition utils.f90:374
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
Implements wall_model_t.
Wall model based on the Monin-Obukhov Similarity Theory for atmospheric boundary layer flows....
Definition most.f90:63
Wall model similar to the Monin-Obukhov Similarity Theory for atmospheric boundary layer flows,...
Wall model based on the log-law for a rough wall. The formula defining the law is ....
Wall model based on Spalding's law of the wall. Reference: http://dx.doi.org/10.1115/1....
Definition spalding.f90:58