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