Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
objective_factory.f90
Go to the documentation of this file.
1! Copyright (c) 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
34submodule(objective) objective_factory_mod
35 use json_utils, only: json_get
36 use utils, only: neko_type_error
37
38 ! Import the objective function types
41
42 implicit none
43
45 character(len=25), parameter :: KNOWN_TYPES(2) = [ character(len=25) :: &
46 "minimum_dissipation", &
47 "lube_term"]
48
49contains
50
51 ! -------------------------------------------------------------------------- !
52 ! Factory function
53
60 module subroutine objective_factory(object, json, design, simulation)
61 class(objective_t), allocatable, intent(inout) :: object
62 type(json_file), intent(inout) :: json
63 class(design_t), intent(in) :: design
64 type(simulation_t), target, intent(inout) :: simulation
65 character(len=:), allocatable :: type
66
67 if (allocated(object)) then
68 call object%free()
69 deallocate(object)
70 end if
71
72 call json_get(json, "type", type)
73 select case (trim(type))
74 case ("minimum_dissipation")
75 allocate(minimum_dissipation_objective_t::object)
76 case ("lube_term")
77 allocate(lube_term_objective_t::object)
78
79 case default
80 call neko_type_error("Objective", type, KNOWN_TYPES)
81 end select
82
83 call object%init_json(json, design, simulation)
84 end subroutine objective_factory
85
86end submodule objective_factory_mod
Implements the design_t.
Definition design.f90:34
Implements the lube_term_objective_t type.
Implements the minimum_dissipation_objective_t type.
Implements the objective_t type.
Definition objective.f90:35
Implements the steady_problem_t type.
An objective function corresponding to minimum dissipation $ F = \int_\Omega |\nabla u|^2 d \Omega + ...
An objective function corresponding to minimum dissipation.