37 use num_types,
only: rp
41 use json_module,
only: json_file
42 use json_utils,
only: json_get_or_default
43 use field,
only: field_t
44 use field_math,
only: field_copy, field_cadd, field_col2, field_cmult
45 use neko_config,
only: neko_bcknd_device
46 use math,
only: glsc2, copy
47 use device_math,
only: device_glsc2, device_copy
49 use coefs,
only: coef_t
50 use scratch_registry,
only: neko_scratch_registry
51 use utils,
only: neko_error
55 use field_math,
only: field_addcol3, field_col3
66 type(field_t),
pointer :: phi
68 real(kind=rp) :: phi_ref
70 class(coef_t),
pointer :: coef
72 real(kind=rp) :: domain_volume
78 type(field_t),
pointer :: u, v, w, u_adj, v_adj, w_adj
83 procedure,
public, pass(this) :: init_json_sim => &
84 scalar_mixing_init_json_sim
86 procedure,
public, pass(this) :: init_from_attributes => &
87 scalar_mixing_init_attributes
89 procedure,
public, pass(this) :: free => scalar_mixing_free
91 procedure,
public, pass(this) :: update_value => &
92 scalar_mixing_update_value
94 procedure,
public, pass(this) :: update_sensitivity => &
95 scalar_mixing_update_sensitivity
106 subroutine scalar_mixing_init_json_sim(this, json, design, simulation)
108 type(json_file),
intent(inout) :: json
109 class(
design_t),
intent(in) :: design
111 real(kind=rp) :: weight
112 real(kind=rp) :: phi_ref
113 character(len=:),
allocatable :: name
114 character(len=:),
allocatable :: mask_name
116 call json_get_or_default(json,
"weight", weight, 1.0_rp)
117 call json_get_or_default(json,
"mask_name", mask_name,
"")
118 call json_get_or_default(json,
"target_concentration", phi_ref, 0.5_rp)
119 call json_get_or_default(json,
"name", name,
"Scalar Mixing")
122 call this%init_from_attributes(
design, simulation, weight, name, &
124 end subroutine scalar_mixing_init_json_sim
134 subroutine scalar_mixing_init_attributes(this, design, simulation, weight, &
135 name, mask_name, phi_ref)
137 class(
design_t),
intent(in) :: design
139 real(kind=rp),
intent(in) :: weight
140 real(kind=rp),
intent(in) :: phi_ref
141 character(len=*),
intent(in) :: mask_name
142 character(len=*),
intent(in) :: name
146 if (.not.
allocated(simulation%adjoint_case%scalar_adj))
then
147 call neko_error(
"adjoint passive scalar not initialized")
151 call this%init_base(name,
design%size(), weight, mask_name)
154 this%coef => simulation%fluid%c_Xh
157 if (this%has_mask)
then
158 this%domain_volume = compute_masked_volume(this%mask, this%coef)
160 this%domain_volume = this%coef%volume
165 associate(f_phi_adj => simulation%adjoint_scalar%f_Xh)
168 this%phi_ref = phi_ref
171 this%phi => simulation%scalar%s
174 call adjoint_forcing%init_from_components(f_phi_adj, this%phi, &
175 this%weight, this%phi_ref, this%mask, this%has_mask, this%coef)
180 call simulation%adjoint_scalar%source_term%add_source_term(adjoint_forcing)
184 this%u => simulation%fluid%u
185 this%v => simulation%fluid%v
186 this%w => simulation%fluid%w
187 this%u_adj => simulation%adjoint_fluid%u_adj
188 this%v_adj => simulation%adjoint_fluid%v_adj
189 this%w_adj => simulation%adjoint_fluid%w_adj
190 end subroutine scalar_mixing_init_attributes
193 subroutine scalar_mixing_free(this)
195 call this%free_base()
205 end subroutine scalar_mixing_free
210 subroutine scalar_mixing_update_value(this, design)
212 class(design_t),
intent(in) :: design
213 type(field_t),
pointer :: work
214 integer :: temp_indices(1), n
222 call neko_scratch_registry%request_field(work, temp_indices(1))
226 call field_copy(work, this%phi)
228 call field_cadd(work, -this%phi_ref)
230 call field_col2(work, work)
232 if (this%has_mask)
then
233 call mask_exterior_const(work, this%mask, 0.0_rp)
236 if (neko_bcknd_device .eq. 1)
then
237 this%value = device_glsc2(work%x_d, this%coef%B_d, n)
239 this%value = glsc2(work%x, this%coef%B, n)
243 this%value = 0.5_rp * this%value / this%domain_volume
246 call neko_scratch_registry%relinquish_field(temp_indices)
247 end subroutine scalar_mixing_update_value
253 subroutine scalar_mixing_update_sensitivity(this, design)
255 class(design_t),
intent(in) :: design
259 type(field_t),
pointer :: work
260 integer :: temp_indices(1)
262 call neko_scratch_registry%request_field(work, temp_indices(1))
265 call field_col3(work, this%u, this%u_adj)
266 call field_addcol3(work, this%v, this%v_adj)
267 call field_addcol3(work, this%w, this%w_adj)
269 call field_cmult(work, -1.0_rp)
271 if (neko_bcknd_device .eq. 1)
then
272 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
274 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
277 call neko_scratch_registry%relinquish_field(temp_indices)
279 end subroutine scalar_mixing_update_sensitivity
Implements the adjoint_mixing_scalar_source_term type.
Some common Masking operations we may need.
Implements the objective_t type.
An objective function corresponding to the mixing of a passive scalar .
Implements the steady_problem_t type.
The abstract objective type.
An objective function corresponding to the mixing of a passive scalar .