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
61 type(field_t),
pointer :: chi => null()
63 type(field_t),
pointer :: u => null()
65 type(field_t),
pointer :: v => null()
67 type(field_t),
pointer :: w => null()
70 type(space_t),
pointer :: xh_gll
72 type(space_t),
pointer :: xh_gl
74 type(coef_t),
pointer :: c_xh_gl
76 type(interpolator_t),
pointer :: gll_to_gl
80 type(scratch_registry_t),
pointer :: scratch_gl
84 procedure, pass(this) :: init => &
87 procedure, pass(this) :: init_from_components => &
88 simple_brinkman_source_term_init_from_components
90 procedure, pass(this) :: free => simple_brinkman_source_term_free
92 procedure, pass(this) :: compute_ => simple_brinkman_source_term_compute
105 type(json_file),
intent(inout) :: json
106 type(field_list_t),
intent(in),
target :: fields
107 type(coef_t),
intent(in),
target :: coef
108 character(len=*),
intent(in) :: variable_name
127 subroutine simple_brinkman_source_term_init_from_components(this, &
128 f_x, f_y, f_z, chi, u, v, w, coef, c_Xh_GL, GLL_to_GL, dealias, &
131 type(field_t),
pointer,
intent(in) :: f_x, f_y, f_z
132 type(field_list_t) :: fields
133 type(coef_t),
intent(in) :: coef
134 type(coef_t),
intent(in),
target :: c_Xh_GL
135 type(interpolator_t),
intent(in),
target :: GLL_to_GL
136 logical,
intent(in) :: dealias
137 real(kind=rp) :: start_time
138 real(kind=rp) :: end_time
139 type(field_t),
intent(in),
target :: u, v, w
140 type(field_t),
intent(in),
target :: chi
141 type(scratch_registry_t),
intent(in),
target :: scratch_GL
146 end_time = 100000000.0_rp
153 call fields%assign(1, f_x)
154 call fields%assign(2, f_y)
155 call fields%assign(3, f_z)
157 call this%init_base(fields, coef, start_time, end_time)
166 this%dealias = dealias
167 this%c_Xh_GL => c_xh_gl
168 this%Xh_GL => this%c_Xh_GL%Xh
169 this%Xh_GLL => this%coef%Xh
170 this%GLL_to_GL => gll_to_gl
171 this%scratch_GL => scratch_gl
173 end subroutine simple_brinkman_source_term_init_from_components
176 subroutine simple_brinkman_source_term_free(this)
179 call this%free_base()
184 nullify(this%c_Xh_GL)
187 nullify(this%GLL_to_GL)
188 nullify(this%scratch_GL)
190 end subroutine simple_brinkman_source_term_free
195 subroutine simple_brinkman_source_term_compute(this, time)
197 type(time_state_t),
intent(in) :: time
198 type(field_t),
pointer :: fu, fv, fw
199 type(field_t),
pointer :: work
200 type(field_t),
pointer :: accumulate, fld_GL, chi_GL
201 integer :: temp_indices(1)
202 integer :: temp_indices_GL(3)
205 fu => this%fields%get_by_index(1)
206 fv => this%fields%get_by_index(2)
207 fw => this%fields%get_by_index(3)
209 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
211 if (this%dealias)
then
212 nel = this%coef%msh%nelv
213 n_gl = nel * this%Xh_GL%lxyz
214 call this%scratch_GL%request_field(accumulate, &
215 temp_indices_gl(1), .false.)
216 call this%scratch_GL%request_field(fld_gl, &
217 temp_indices_gl(2), .false.)
218 call this%scratch_GL%request_field(chi_gl, &
219 temp_indices_gl(3), .false.)
221 call this%GLL_to_GL%map(chi_gl%x, this%chi%x, nel, this%Xh_GL)
224 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
225 call field_col3(accumulate, chi_gl, fld_gl)
227 if (neko_bcknd_device .eq. 1)
then
228 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
229 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
230 call device_invcol2(work%x_d, this%coef%B_d, work%size())
232 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
233 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
234 call invcol2(work%x, this%coef%B, work%size())
236 call field_sub2(fu, work)
239 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
240 call field_col3(accumulate, chi_gl, fld_gl)
242 if (neko_bcknd_device .eq. 1)
then
243 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
244 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
245 call device_invcol2(work%x_d, this%coef%B_d, work%size())
247 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
248 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
249 call invcol2(work%x, this%coef%B, work%size())
251 call field_sub2(fv, work)
254 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
255 call field_col3(accumulate, chi_gl, fld_gl)
257 if (neko_bcknd_device .eq. 1)
then
258 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
259 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
260 call device_invcol2(work%x_d, this%coef%B_d, work%size())
262 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
263 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
264 call invcol2(work%x, this%coef%B, work%size())
266 call field_sub2(fw, work)
268 call this%scratch_GL%relinquish_field(temp_indices_gl)
271 call field_subcol3(fu, this%u, this%chi)
272 call field_subcol3(fv, this%v, this%chi)
273 call field_subcol3(fw, this%w, this%chi)
277 call neko_scratch_registry%relinquish_field(temp_indices)
279 end subroutine simple_brinkman_source_term_compute
Implements the simple_brinkman_source_term_t type.
subroutine simple_brinkman_source_term_init_from_json(this, json, fields, coef, variable_name)
The common constructor using a JSON object.
A simple Brinkman source term.