Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
source_term_fctry.f90
Go to the documentation of this file.
1! Copyright (c) 2023-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(source_term) source_term_fctry
43 use json_utils, only : json_get
45 implicit none
46
47 ! List of all possible types created by the factory routine
48 character(len=25) :: SOURCE_KNOWN_TYPES(7) = [character(len=25) :: &
49 "constant", &
50 "boussinesq", &
51 "coriolis", &
52 "centrifugal", &
53 "gradient_jump_penalty", &
54 "brinkman", &
55 "sponge"]
56
57contains
58
63 module subroutine source_term_factory(object, json, fields, coef, &
64 variable_name)
65 class(source_term_t), allocatable, intent(inout) :: object
66 type(json_file), intent(inout) :: json
67 type(field_list_t), intent(inout) :: fields
68 type(coef_t), intent(inout) :: coef
69 character(len=*), intent(in) :: variable_name
70 character(len=:), allocatable :: type_name
71 character(len=:), allocatable :: type_string
72
73 call json_get(json, "type", type_name)
74
75 ! Allocate
76 call source_term_allocator(object, type_name)
77
78 ! Initialize
79 call object%init(json, fields, coef, variable_name)
80
81 end subroutine source_term_factory
82
86 module subroutine source_term_allocator(object, type_name)
87 class(source_term_t), allocatable, intent(inout) :: object
88 character(len=:), allocatable, intent(in) :: type_name
89 integer :: i
90
91 if (allocated(object)) then
92 call object%free()
93 deallocate(object)
94 end if
95
96 select case (trim(type_name))
97 case ("constant")
98 allocate(const_source_term_t::object)
99 case ("boussinesq")
100 allocate(boussinesq_source_term_t::object)
101 case ("coriolis")
102 allocate(coriolis_source_term_t::object)
103 case ("centrifugal")
104 allocate(centrifugal_source_term_t::object)
105 case ("brinkman")
106 allocate(brinkman_source_term_t::object)
107 case ("sponge")
108 allocate(sponge_source_term_t::object)
109 case ("gradient_jump_penalty")
110 allocate(gradient_jump_penalty_t::object)
111 case default
112 do i = 1, source_term_registry_size
113 if (trim(type_name) .eq. trim(source_term_registry(i)%type_name)) then
114 call source_term_registry(i)%allocator(object)
115 return
116 end if
117 end do
118
119 call neko_type_error("source term", type_name, source_known_types)
120 end select
121 end subroutine source_term_allocator
122
128 module subroutine register_source_term(type_name, allocator)
129 character(len=*), intent(in) :: type_name
130 procedure(source_term_allocate), pointer, intent(in) :: allocator
131 type(allocator_entry), allocatable :: temp(:)
132 integer :: i
133
134 do i = 1, size(source_known_types)
135 if (trim(type_name) .eq. trim(source_known_types(i))) then
136 call neko_type_registration_error("source term", type_name, .true.)
137 end if
138 end do
139
140 do i = 1, source_term_registry_size
141 if (trim(type_name) .eq. trim(source_term_registry(i)%type_name)) then
142 call neko_type_registration_error("source term", type_name, .false.)
143 end if
144 end do
145
146 ! Expand registry
147 if (source_term_registry_size .eq. 0) then
148 allocate(source_term_registry(1))
149 else
150 allocate(temp(source_term_registry_size + 1))
151 temp(1:source_term_registry_size) = source_term_registry
152 call move_alloc(temp, source_term_registry)
153 end if
154
155 source_term_registry_size = source_term_registry_size + 1
156 source_term_registry(source_term_registry_size)%type_name = type_name
157 source_term_registry(source_term_registry_size)%allocator => allocator
158 end subroutine register_source_term
159
160end submodule source_term_fctry
Retrieves a parameter by name or throws an error.
Implements the boussinesq_source_term_t type.
Implements the brinkman_source_term_t type.
Implements the centrifugal_source_term_t type. Maintainer: Adam Peplinski.
Implements the const_source_term_t type.
Implements the coriolis_source_term_t type. Maintainer: Timofey Mukha.
Implements gradient_jump_penalty_t.
Utilities for retrieving parameters from the case files.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Implements the sponge_t type.
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
Bouyancy source term accroding to the Boussinesq approximation.
A Brinkman source term. The region and strength are controlled by assigning regions types and brinkma...
This source term adds the centrifugal force.
A constant source term. The strength is specified with the values keyword, which should be an array,...
This source term adds the Coriolis force.