Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
les_model_fctry.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2025, 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!
33submodule(les_model) les_model_fctry
34 use vreman, only : vreman_t
35 use smagorinsky, only : smagorinsky_t
37 use sigma, only : sigma_t
39 use wale, only : wale_t
40 use deardorff, only : deardorff_t
42 implicit none
43
44 ! List of all possible types created by the factory routine
45 character(len=20) :: LES_KNOWN_TYPES(6) = [character(len=20) :: &
46 "vreman", &
47 "smagorinsky", &
48 "dymamic_smagorinsky", &
49 "sigma", &
50 "wale", &
51 "deardorff"]
52
53contains
61 module subroutine les_model_factory(object, type_name, fluid, json)
62 class(les_model_t), allocatable, intent(inout) :: object
63 character(len=*), intent(in) :: type_name
64 class(fluid_scheme_base_t), intent(inout) :: fluid
65 type(json_file), intent(inout) :: json
66 character(len=:), allocatable :: type_string
67
68 call les_model_allocator(object, type_name)
69 call object%init(fluid, json)
70 end subroutine les_model_factory
71
75 module subroutine les_model_allocator(object, type_name)
76 class(les_model_t), allocatable, intent(inout) :: object
77 character(len=*), intent(in) :: type_name
78 integer :: i
79
80 if (allocated(object)) then
81 call object%free()
82 deallocate(object)
83 end if
84
85 select case (trim(type_name))
86 case ('vreman')
87 allocate(vreman_t::object)
88 case ('smagorinsky')
89 allocate(smagorinsky_t::object)
90 case ('dynamic_smagorinsky')
91 allocate(dynamic_smagorinsky_t::object)
92 case ('sigma')
93 allocate(sigma_t::object)
94 case ('wale')
95 allocate(wale_t::object)
96 case ('deardorff')
97 allocate(deardorff_t::object)
98 case default
99 do i = 1, les_model_registry_size
100 if (trim(type_name) == trim(les_model_registry(i)%type_name)) then
101 call les_model_registry(i)%allocator(object)
102 return
103 end if
104 end do
105
106 call neko_type_error("LES model", type_name, les_known_types)
107 end select
108
109 end subroutine les_model_allocator
110
116 module subroutine register_les_model(type_name, allocator)
117 character(len=*), intent(in) :: type_name
118 procedure(les_model_allocate), pointer, intent(in) :: allocator
119 type(allocator_entry), allocatable :: temp(:)
120 integer :: i
121
122 do i = 1, size(les_known_types)
123 if (trim(type_name) .eq. trim(les_known_types(i))) then
124 call neko_type_registration_error("LES model", type_name, .true.)
125 end if
126 end do
127
128 do i = 1, les_model_registry_size
129 if (trim(type_name) .eq. trim(les_model_registry(i)%type_name)) then
130 call neko_type_registration_error("LES model", type_name, .false.)
131 end if
132 end do
133
134 ! Expand registry
135 if (les_model_registry_size == 0) then
136 allocate(les_model_registry(1))
137 else
138 allocate(temp(les_model_registry_size + 1))
139 temp(1:les_model_registry_size) = les_model_registry
140 call move_alloc(temp, les_model_registry)
141 end if
142
143 les_model_registry_size = les_model_registry_size + 1
144 les_model_registry(les_model_registry_size)%type_name = type_name
145 les_model_registry(les_model_registry_size)%allocator => allocator
146 end subroutine register_les_model
147
148end submodule les_model_fctry
Implements deardorff_t.
Definition deardorff.f90:35
Implements dynamic_smagorinsky_t.
Implements les_model_t.
Definition les_model.f90:35
Implements sigma_t.
Definition sigma.f90:35
Implements smagorinsky_t.
Utilities.
Definition utils.f90:35
subroutine, public neko_type_registration_error(base_type, wrong_type, known)
Definition utils.f90:328
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
Implements vreman_t.
Definition vreman.f90:35
Implements wale_t.
Definition wale.f90:35
Implements the deardorff LES model.
Definition deardorff.f90:53
Implements the dynamic Smagorinsky LES model.
Base type of all fluid formulations.
Implements the Sigma LES model.
Definition sigma.f90:52
Implements the smagorinsky LES model.
Implements the Vreman LES model.
Definition vreman.f90:54
Implements the Wale LES model.
Definition wale.f90:53