Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
simple_brinkman_source_term.f90
1! Copyright (c) 2023, 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!
34! a term in the form $\chi \mathbf{u}$
36 use num_types, only: rp
37 use field_list, only: field_list_t
38 use json_module, only: json_file
39 use json_utils, only: json_get, json_get_or_default
40 use source_term, only: source_term_t
41 use coefs, only: coef_t
42 use neko_config, only: neko_bcknd_device
43 use time_state, only: time_state_t
44 use utils, only: neko_error
45 use field, only: field_t
46 use field_math, only: field_subcol3, field_sub2, field_col3
47 use interpolation, only: interpolator_t
48 use space, only: space_t, gl
49 use math, only: col2, invcol2
50 use device_math, only: device_col2, device_invcol2
51 use scratch_registry, only: neko_scratch_registry, scratch_registry_t
52 implicit none
53 private
54
56 ! We have a source term of the form $\chi \mathbf{u}$
57 type, public, extends(source_term_t) :: simple_brinkman_source_term_t
59 type(field_t), pointer :: chi => null()
61 type(field_t), pointer :: u => null()
63 type(field_t), pointer :: v => null()
65 type(field_t), pointer :: w => null()
66 ! --- for over-integration
68 type(space_t), pointer :: xh_gll
70 type(space_t), pointer :: xh_gl
72 type(coef_t), pointer :: c_xh_gl
74 type(interpolator_t), pointer :: gll_to_gl
76 logical :: dealias
78 type(scratch_registry_t), pointer :: scratch_gl
79
80 contains
82 procedure, pass(this) :: init => &
83 simple_brinkman_source_term_init_from_json
85 procedure, pass(this) :: init_from_components => &
86 simple_brinkman_source_term_init_from_components
88 procedure, pass(this) :: free => simple_brinkman_source_term_free
90 procedure, pass(this) :: compute_ => simple_brinkman_source_term_compute
92
93contains
100 subroutine simple_brinkman_source_term_init_from_json(this, json, fields, &
101 coef, variable_name)
102 class(simple_brinkman_source_term_t), intent(inout) :: this
103 type(json_file), intent(inout) :: json
104 type(field_list_t), intent(in), target :: fields
105 type(coef_t), intent(in), target :: coef
106 character(len=*), intent(in) :: variable_name
107
108
109 ! we shouldn't be initializing this from JSON
110 ! maybe throw an error?
111
112
113 end subroutine simple_brinkman_source_term_init_from_json
114
125 subroutine simple_brinkman_source_term_init_from_components(this, &
126 f_x, f_y, f_z, chi, u, v, w, coef, c_Xh_GL, GLL_to_GL, dealias, &
127 scratch_GL)
128 class(simple_brinkman_source_term_t), intent(inout) :: this
129 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
130 type(field_list_t) :: fields
131 type(coef_t), intent(in) :: coef
132 type(coef_t), intent(in), target :: c_Xh_GL
133 type(interpolator_t), intent(in), target :: GLL_to_GL
134 logical, intent(in) :: dealias
135 real(kind=rp) :: start_time
136 real(kind=rp) :: end_time
137 type(field_t), intent(in), target :: u, v, w
138 type(field_t), intent(in), target :: chi
139 type(scratch_registry_t), intent(in), target :: scratch_GL
140
141 ! I wish you didn't need a start time and end time...
142 ! but I'm just going to set a super big number...
143 start_time = 0.0_rp
144 end_time = 100000000.0_rp
145
146 call this%free()
147
148 ! this is copying the fluid source term init
149 ! We package the fields for the source term to operate on in a field list.
150 call fields%init(3)
151 call fields%assign(1, f_x)
152 call fields%assign(2, f_y)
153 call fields%assign(3, f_z)
154
155 call this%init_base(fields, coef, start_time, end_time)
156
157 ! point everything in the correct places
158 this%u => u
159 this%v => v
160 this%w => w
161 ! and get chi out of the design
162 this%chi => chi
163 ! for over integration
164 this%dealias = dealias
165 this%c_Xh_GL => c_xh_gl
166 this%Xh_GL => this%c_Xh_GL%Xh
167 this%Xh_GLL => this%coef%Xh
168 this%GLL_to_GL => gll_to_gl
169 this%scratch_GL => scratch_gl
170
171 end subroutine simple_brinkman_source_term_init_from_components
172
174 subroutine simple_brinkman_source_term_free(this)
175 class(simple_brinkman_source_term_t), intent(inout) :: this
176
177 call this%free_base()
178 nullify(this%u)
179 nullify(this%v)
180 nullify(this%w)
181 nullify(this%chi)
182 nullify(this%c_Xh_GL)
183 nullify(this%Xh_GL)
184 nullify(this%Xh_GLL)
185 nullify(this%GLL_to_GL)
186 nullify(this%scratch_GL)
187
188 end subroutine simple_brinkman_source_term_free
189
193 subroutine simple_brinkman_source_term_compute(this, time)
194 class(simple_brinkman_source_term_t), intent(inout) :: this
195 type(time_state_t), intent(in) :: time
196 type(field_t), pointer :: fu, fv, fw
197 type(field_t), pointer :: work
198 type(field_t), pointer :: accumulate, fld_GL, chi_GL
199 integer :: temp_indices(1)
200 integer :: temp_indices_GL(3)
201 integer :: n_GL, nel
202
203 fu => this%fields%get_by_index(1)
204 fv => this%fields%get_by_index(2)
205 fw => this%fields%get_by_index(3)
206
207 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
208
209 if (this%dealias) then
210 nel = this%coef%msh%nelv
211 n_gl = nel * this%Xh_GL%lxyz
212 call this%scratch_GL%request_field(accumulate, &
213 temp_indices_gl(1), .false.)
214 call this%scratch_GL%request_field(fld_gl, &
215 temp_indices_gl(2), .false.)
216 call this%scratch_GL%request_field(chi_gl, &
217 temp_indices_gl(3), .false.)
218
219 call this%GLL_to_GL%map(chi_gl%x, this%chi%x, nel, this%Xh_GL)
220
221 ! u
222 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
223 call field_col3(accumulate, chi_gl, fld_gl)
224 ! Evaluate term on GL and preempt the GLL premultiplication
225 if (neko_bcknd_device .eq. 1) then
226 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
227 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
228 call device_invcol2(work%x_d, this%coef%B_d, work%size())
229 else
230 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
231 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
232 call invcol2(work%x, this%coef%B, work%size())
233 end if
234 call field_sub2(fu, work)
235
236 ! v
237 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
238 call field_col3(accumulate, chi_gl, fld_gl)
239 ! Evaluate term on GL and preempt the GLL premultiplication
240 if (neko_bcknd_device .eq. 1) then
241 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
242 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
243 call device_invcol2(work%x_d, this%coef%B_d, work%size())
244 else
245 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
246 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
247 call invcol2(work%x, this%coef%B, work%size())
248 end if
249 call field_sub2(fv, work)
250
251 ! w
252 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
253 call field_col3(accumulate, chi_gl, fld_gl)
254 ! Evaluate term on GL and preempt the GLL premultiplication
255 if (neko_bcknd_device .eq. 1) then
256 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
257 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
258 call device_invcol2(work%x_d, this%coef%B_d, work%size())
259 else
260 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
261 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
262 call invcol2(work%x, this%coef%B, work%size())
263 end if
264 call field_sub2(fw, work)
265
266 call this%scratch_GL%relinquish_field(temp_indices_gl)
267 else
268
269 call field_subcol3(fu, this%u, this%chi)
270 call field_subcol3(fv, this%v, this%chi)
271 call field_subcol3(fw, this%w, this%chi)
272
273 end if
274
275 call neko_scratch_registry%relinquish_field(temp_indices)
276
277 end subroutine simple_brinkman_source_term_compute
278
Implements the simple_brinkman_source_term_t type.