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: scratch_registry_t, neko_scratch_registry
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()
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
78 type(scratch_registry_t),
pointer :: scratch_gl
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
100 subroutine simple_brinkman_source_term_init_from_json(this, json, fields, &
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
113 end subroutine simple_brinkman_source_term_init_from_json
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, &
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
144 end_time = 100000000.0_rp
151 call fields%assign(1, f_x)
152 call fields%assign(2, f_y)
153 call fields%assign(3, f_z)
155 call this%init_base(fields, coef, start_time, end_time)
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
171 end subroutine simple_brinkman_source_term_init_from_components
174 subroutine simple_brinkman_source_term_free(this)
177 call this%free_base()
182 nullify(this%c_Xh_GL)
185 nullify(this%GLL_to_GL)
186 nullify(this%scratch_GL)
188 end subroutine simple_brinkman_source_term_free
193 subroutine simple_brinkman_source_term_compute(this, time)
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)
203 fu => this%fields%get_by_index(1)
204 fv => this%fields%get_by_index(2)
205 fw => this%fields%get_by_index(3)
207 call neko_scratch_registry%request_field(work, temp_indices(1))
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, temp_indices_gl(1))
213 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2))
214 call this%scratch_GL%request_field(chi_gl, temp_indices_gl(3))
216 call this%GLL_to_GL%map(chi_gl%x, this%chi%x, nel, this%Xh_GL)
219 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
220 call field_col3(accumulate, chi_gl, fld_gl)
222 if (neko_bcknd_device .eq. 1)
then
223 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
224 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
225 call device_invcol2(work%x_d, this%coef%B_d, work%size())
227 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
228 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
229 call invcol2(work%x, this%coef%B, work%size())
231 call field_sub2(fu, work)
234 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
235 call field_col3(accumulate, chi_gl, fld_gl)
237 if (neko_bcknd_device .eq. 1)
then
238 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
239 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
240 call device_invcol2(work%x_d, this%coef%B_d, work%size())
242 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
243 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
244 call invcol2(work%x, this%coef%B, work%size())
246 call field_sub2(fv, work)
249 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
250 call field_col3(accumulate, chi_gl, fld_gl)
252 if (neko_bcknd_device .eq. 1)
then
253 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
254 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
255 call device_invcol2(work%x_d, this%coef%B_d, work%size())
257 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
258 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
259 call invcol2(work%x, this%coef%B, work%size())
261 call field_sub2(fw, work)
263 call this%scratch_GL%relinquish_field(temp_indices_gl)
266 call field_subcol3(fu, this%u, this%chi)
267 call field_subcol3(fv, this%v, this%chi)
268 call field_subcol3(fw, this%w, this%chi)
272 call neko_scratch_registry%relinquish_field(temp_indices)
274 end subroutine simple_brinkman_source_term_compute
Implements the simple_brinkman_source_term_t type.
A simple Brinkman source term.