Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
simulation_component_fctry.f90
Go to the documentation of this file.
1! Copyright (c) 2024-2026, 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(simulation_component) simulation_component_fctry
41 use user_stats, only : user_stats_t
42 use lambda2, only : lambda2_t
43 use probes, only : probes_t
44 use les_simcomp, only : les_simcomp_t
47 use curl_simcomp, only : curl_t
54 implicit none
55
56 ! List of all possible types created by the factory routine
57 character(len=20) :: SIMCOMPS_KNOWN_TYPES(16) = [character(len=20) :: &
58 "lambda2", &
59 "probes", &
60 "les_model", &
61 "field_writer", &
62 "fluid_stats", &
63 "fluid_sgs_stats", &
64 "scalar_stats", &
65 "scalar_sgs_stats", &
66 "grad", &
67 "div", &
68 "curl", &
69 "derivative", &
70 "weak_grad", &
71 "force_torque", &
72 "user_stats", &
73 "spectral_error"]
74
75contains
76
81 module subroutine simulation_component_factory(object, json, case)
82 class(simulation_component_t), allocatable, intent(inout) :: object
83 type(json_file), intent(inout) :: json
84 class(case_t), intent(inout), target :: case
85 character(len=:), allocatable :: type_name
86 character(len=:), allocatable :: type_string
87 logical :: is_user
88
89 ! Check if this is a user-defined component
90 call json_get_or_default(json, "is_user", is_user, .false.)
91 if (is_user) return
92
93 ! Get the type name
94 call json_get(json, "type", type_name)
95
96 ! Allocate
97 call simulation_component_allocator(object, type_name)
98
99 ! Initialize
100 call object%init(json, case)
101
102 end subroutine simulation_component_factory
103
107 module subroutine simulation_component_allocator(object, type_name)
108 class(simulation_component_t), allocatable, intent(inout) :: object
109 character(len=*), intent(in):: type_name
110 integer :: i
111
112 if (allocated(object)) then
113 call object%free()
114 deallocate(object)
115 end if
116
117 select case (trim(type_name))
118 case ("lambda2")
119 allocate(lambda2_t::object)
120 case ("probes")
121 allocate(probes_t::object)
122 case ("les_model")
123 allocate(les_simcomp_t::object)
124 case ("field_writer")
125 allocate(field_writer_t::object)
126 case ("weak_grad")
127 allocate(weak_gradient_t::object)
128 case ("grad")
129 allocate(gradient_t::object)
130 case ("derivative")
131 allocate(derivative_t::object)
132 case ("curl")
133 allocate(curl_t::object)
134 case ("div")
135 allocate(divergence_t::object)
136 case ("force_torque")
137 allocate(force_torque_t::object)
138 case ("fluid_stats")
139 allocate(fluid_stats_simcomp_t::object)
140 case ("fluid_sgs_stats")
141 allocate(fluid_sgs_stats_simcomp_t::object)
142 case ("scalar_stats")
143 allocate(scalar_stats_simcomp_t::object)
144 case ("scalar_sgs_stats")
145 allocate(scalar_sgs_stats_simcomp_t::object)
146 case ("user_stats")
147 allocate(user_stats_t::object)
148 case ("spectral_error")
149 allocate(spectral_error_t::object)
150 case default
151 do i = 1, simcomp_registry_size
152 if (trim(type_name) == &
153 trim(simcomp_registry(i)%type_name)) then
154 call simcomp_registry(i)%allocator(object)
155 return
156 end if
157 end do
158 call neko_type_error("simulation component", trim(type_name), &
159 simcomps_known_types)
160 end select
161
162 end subroutine simulation_component_allocator
163
169 module subroutine register_simulation_component(type_name, allocator)
170 character(len=*), intent(in) :: type_name
171 procedure(simulation_component_allocate), pointer, intent(in) :: allocator
172 type(allocator_entry), allocatable :: temp(:)
173 integer :: i
174
175 do i = 1, size(simcomps_known_types)
176 if (trim(type_name) .eq. trim(simcomps_known_types(i))) then
177 call neko_type_registration_error("simulation component", type_name, &
178 .true.)
179 end if
180 end do
181
182 do i = 1, simcomp_registry_size
183 if (trim(type_name) .eq. &
184 trim(simcomp_registry(i)%type_name)) then
185 call neko_type_registration_error("simulation component", type_name, &
186 .false.)
187 end if
188 end do
189
190 ! Expand registry
191 if (simcomp_registry_size == 0) then
192 allocate(simcomp_registry(1))
193 else
194 allocate(temp(simcomp_registry_size + 1))
195 temp(1:simcomp_registry_size) = simcomp_registry
196 call move_alloc(temp, simcomp_registry)
197 end if
198
199 simcomp_registry_size = simcomp_registry_size + 1
200 simcomp_registry(simcomp_registry_size)%type_name = type_name
201 simcomp_registry(simcomp_registry_size)%allocator => allocator
202 end subroutine register_simulation_component
203
204end submodule simulation_component_fctry
Defines a simulation case.
Definition case.f90:34
Implements the curl_t type.
Implements the derivative_t type.
Implements the divergence_t type.
Implements the field_writer_t type.
Implements the fluid_sgs_stats_simcomp_t type.
Implements the fluid_stats_simcomp_t type.
Implements the force_torque_t type.
Implements the gradient_t type.
A simulation component that computes lambda2 The values are stored in the field registry under the na...
Definition lambda2.f90:37
Implements the les_simcomp_t type.
Implements probes.
Definition probes.F90:37
Implements the scalar_sgs_stats_simcomp_t type.
Implements the scalar_stats_simcomp_t type.
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
Implements type spectral_error_t.
Implements the user_stats_t type.
Utilities.
Definition utils.f90:35
character(:) function, allocatable, public concat_string_array(array, sep, prepend)
Concatenate an array of strings into one string with array items separated by spaces.
Definition utils.f90:356
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 the weak_gradient_t type.
A simulation component that computes the curl of a vector field. Added to the field registry as curl_...
A simulation component that computes a derivative of a field. Wraps the duxyz operator.
A simulation component that computes the divergence of a vector field. Added to the field registry as...
A simulation component that writes a 3d field to a file.
A simulation component that computes the subgrid-scale contributions to the Reynolds stresses in LES.
A simulation component that computes the velocity and pressure statistics up to 4th order....
A simulation component that computes the force and torque on a given boundary zone.
A simulation component that computes the gradient of a field. Wraps the gradient operator.
A simulation component that drives the computation of the SGS viscosity.
A simulation component that computes the subgrid-scale contributions to the Reynolds stresses in LES.
A simulation component that computes the scalar statistics for the skewness, kurtosis,...
Provides tools to calculate the spectral error indicator.
A simulation component that computes the averages of fields in the registry.
A simulation component that computes the weak gradient of a field. Wraps the opgradient operator.