Neko 1.99.1
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, 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
38 use lambda2, only : lambda2_t
39 use probes, only : probes_t
40 use les_simcomp, only : les_simcomp_t
43 use curl_simcomp, only : curl_t
50 implicit none
51
52 ! List of all possible types created by the factory routine
53 character(len=20) :: SIMCOMPS_KNOWN_TYPES(12) = [character(len=20) :: &
54 "lambda2", &
55 "probes", &
56 "les_model", &
57 "field_writer", &
58 "fluid_stats", &
59 "grad", &
60 "div", &
61 "curl", &
62 "derivative", &
63 "weak_grad", &
64 "force_torque", &
65 "spectral_error"]
66
67contains
68
73 module subroutine simulation_component_factory(object, json, case)
74 class(simulation_component_t), allocatable, intent(inout) :: object
75 type(json_file), intent(inout) :: json
76 class(case_t), intent(inout), target :: case
77 character(len=:), allocatable :: type_name
78 character(len=:), allocatable :: type_string
79 logical :: is_user
80
81 ! Check if this is a user-defined component
82 call json_get_or_default(json, "is_user", is_user, .false.)
83 if (is_user) return
84
85 ! Get the type name
86 call json_get(json, "type", type_name)
87
88 ! Allocate
89 call simulation_component_allocator(object, type_name)
90
91 ! Initialize
92 call object%init(json, case)
93
94 end subroutine simulation_component_factory
95
99 module subroutine simulation_component_allocator(object, type_name)
100 class(simulation_component_t), allocatable, intent(inout) :: object
101 character(len=*), intent(in):: type_name
102 integer :: i
103
104 select case (trim(type_name))
105 case ("lambda2")
106 allocate(lambda2_t::object)
107 case ("probes")
108 allocate(probes_t::object)
109 case ("les_model")
110 allocate(les_simcomp_t::object)
111 case ("field_writer")
112 allocate(field_writer_t::object)
113 case ("weak_grad")
114 allocate(weak_gradient_t::object)
115 case ("grad")
116 allocate(gradient_t::object)
117 case ("derivative")
118 allocate(derivative_t::object)
119 case ("curl")
120 allocate(curl_t::object)
121 case ("div")
122 allocate(divergence_t::object)
123 case ("force_torque")
124 allocate(force_torque_t::object)
125 case ("fluid_stats")
126 allocate(fluid_stats_simcomp_t::object)
127 case ("spectral_error")
128 allocate(spectral_error_t::object)
129 case default
130 do i = 1, simcomp_registry_size
131 if (trim(type_name) == &
132 trim(simcomp_registry(i)%type_name)) then
133 call simcomp_registry(i)%allocator(object)
134 return
135 end if
136 end do
137 call neko_type_error("simulation component", trim(type_name), &
138 simcomps_known_types)
139 end select
140
141 end subroutine simulation_component_allocator
142
148 module subroutine register_simulation_component(type_name, allocator)
149 character(len=*), intent(in) :: type_name
150 procedure(simulation_component_allocate), pointer, intent(in) :: allocator
151 type(allocator_entry), allocatable :: temp(:)
152 integer :: i
153
154 do i = 1, size(simcomps_known_types)
155 if (trim(type_name) .eq. trim(simcomps_known_types(i))) then
156 call neko_type_registration_error("simulation component", type_name, &
157 .true.)
158 end if
159 end do
160
161 do i = 1, simcomp_registry_size
162 if (trim(type_name) .eq. &
163 trim(simcomp_registry(i)%type_name)) then
164 call neko_type_registration_error("simulation component", type_name, &
165 .false.)
166 end if
167 end do
168
169 ! Expand registry
170 if (simcomp_registry_size == 0) then
171 allocate(simcomp_registry(1))
172 else
173 allocate(temp(simcomp_registry_size + 1))
174 temp(1:simcomp_registry_size) = simcomp_registry
175 call move_alloc(temp, simcomp_registry)
176 end if
177
178 simcomp_registry_size = simcomp_registry_size + 1
179 simcomp_registry(simcomp_registry_size)%type_name = type_name
180 simcomp_registry(simcomp_registry_size)%allocator => allocator
181 end subroutine register_simulation_component
182
183end 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_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
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
Implements type spectral_error_t.
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:294
subroutine, public neko_type_registration_error(base_type, wrong_type, known)
Definition utils.f90:266
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:251
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 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.
Provides tools to calculate the spectral error indicator.
A simulation component that computes the weak gradient of a field. Wraps the opgradient operator.