Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
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 implicit none
41
42 ! List of all possible types created by the factory routine
43 character(len=20) :: LES_KNOWN_TYPES(5) = [character(len=20) :: &
44 "vreman", &
45 "smagorinsky", &
46 "dymamic_smagorinsky", &
47 "sigma", &
48 "wale"]
49
50contains
58 module subroutine les_model_factory(object, type_name, fluid, json)
59 class(les_model_t), allocatable, intent(inout) :: object
60 character(len=*), intent(in) :: type_name
61 class(fluid_scheme_base_t), intent(inout) :: fluid
62 type(json_file), intent(inout) :: json
63 character(len=:), allocatable :: type_string
64
65 call les_model_allocator(object, type_name)
66 call object%init(fluid, json)
67 end subroutine les_model_factory
68
72 module subroutine les_model_allocator(object, type_name)
73 class(les_model_t), allocatable, intent(inout) :: object
74 character(len=*), intent(in) :: type_name
75 integer :: i
76
77 if (allocated(object)) deallocate(object)
78
79 select case (trim(type_name))
80 case ('vreman')
81 allocate(vreman_t::object)
82 case ('smagorinsky')
83 allocate(smagorinsky_t::object)
84 case ('dynamic_smagorinsky')
85 allocate(dynamic_smagorinsky_t::object)
86 case ('sigma')
87 allocate(sigma_t::object)
88 case ('wale')
89 allocate(wale_t::object)
90 case default
91 do i = 1, les_model_registry_size
92 if (trim(type_name) == trim(les_model_registry(i)%type_name)) then
93 call les_model_registry(i)%allocator(object)
94 return
95 end if
96 end do
97
98 call neko_type_error("LES model", type_name, les_known_types)
99 end select
100
101 end subroutine les_model_allocator
102
107 module subroutine register_les_model(type_name, allocator)
108 character(len=*), intent(in) :: type_name
109 procedure(les_model_allocate), pointer, intent(in) :: allocator
110 type(allocator_entry), allocatable :: temp(:)
111
112 ! Expand registry
113 if (les_model_registry_size == 0) then
114 allocate(les_model_registry(1))
115 else
116 allocate(temp(les_model_registry_size + 1))
117 temp(1:les_model_registry_size) = les_model_registry
118 call move_alloc(temp, les_model_registry)
119 end if
120
121 les_model_registry_size = les_model_registry_size + 1
122 les_model_registry(les_model_registry_size)%type_name = type_name
123 les_model_registry(les_model_registry_size)%allocator => allocator
124 end subroutine register_les_model
125
126end submodule les_model_fctry
Implements dynamic_smagorinsky_t.
Implements les_model_t.
Definition les_model.f90:35
Implements sigma_t.
Definition sigma.f90:35
Implements smagorinsky_t.
Implements vreman_t.
Definition vreman.f90:35
Implements wale_t.
Definition wale.f90:35
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:52
Implements the Wale LES model.
Definition wale.f90:52