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
56 use field_math,
only: field_addcol3, field_col3
67 type(field_t),
pointer :: phi
69 real(kind=rp) :: phi_ref
71 class(coef_t),
pointer :: coef
73 real(kind=rp) :: domain_volume
75 character(len=:),
allocatable :: scalar_name
81 type(field_t),
pointer :: u, v, w, u_adj, v_adj, w_adj
86 procedure,
public, pass(this) :: init_json_sim => &
87 scalar_mixing_init_json_sim
89 procedure,
public, pass(this) :: init_from_attributes => &
90 scalar_mixing_init_attributes
92 procedure,
public, pass(this) :: free => scalar_mixing_free
94 procedure,
public, pass(this) :: update_value => &
95 scalar_mixing_update_value
97 procedure,
public, pass(this) :: update_sensitivity => &
98 scalar_mixing_update_sensitivity
109 subroutine scalar_mixing_init_json_sim(this, json, design, simulation)
111 type(json_file),
intent(inout) :: json
112 class(
design_t),
intent(in) :: design
114 real(kind=rp) :: weight
115 real(kind=rp) :: phi_ref
116 character(len=:),
allocatable :: name
117 character(len=:),
allocatable :: mask_name
118 character(len=:),
allocatable :: scalar_name
120 call json_get_or_default(json,
"weight", weight, 1.0_rp)
121 call json_get_or_default(json,
"mask_name", mask_name,
"")
122 call json_get_or_default(json,
"target_concentration", phi_ref, 0.5_rp)
123 call json_get_or_default(json,
"name", name,
"Scalar Mixing")
124 call json_get_or_default(json,
"scalar_name", scalar_name,
"s")
127 call this%init_from_attributes(
design, simulation, weight, name, &
128 mask_name, phi_ref, scalar_name)
129 end subroutine scalar_mixing_init_json_sim
140 subroutine scalar_mixing_init_attributes(this, design, simulation, weight, &
141 name, mask_name, phi_ref, scalar_name)
143 class(
design_t),
intent(in) :: design
145 real(kind=rp),
intent(in) :: weight
146 real(kind=rp),
intent(in) :: phi_ref
147 character(len=*),
intent(in) :: mask_name
148 character(len=*),
intent(in) :: name
149 character(len=*),
intent(in) :: scalar_name
152 integer :: i_scalar, i_adjoint_scalar
155 if (.not.
allocated(simulation%adjoint_case%adjoint_scalars))
then
156 call neko_error(
"adjoint passive scalar not initialized")
160 call this%init_base(name,
design%size(), weight, mask_name)
163 this%coef => simulation%fluid%c_Xh
166 if (this%has_mask)
then
167 this%domain_volume = compute_masked_volume(this%mask, this%coef)
169 this%domain_volume = this%coef%volume
172 this%scalar_name = trim(scalar_name)
176 simulation%adjoint_scalars, this%scalar_name)
180 associate(f_phi_adj => &
181 simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar)%f_Xh)
184 this%phi_ref = phi_ref
187 this%phi => simulation%scalars%scalar_fields(i_scalar)%s
190 call adjoint_forcing%init_from_components(f_phi_adj, this%phi, &
191 this%weight, this%phi_ref, this%mask, this%has_mask, this%coef)
196 call simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar) &
197 %source_term%add_source_term(adjoint_forcing)
201 this%u => simulation%fluid%u
202 this%v => simulation%fluid%v
203 this%w => simulation%fluid%w
204 this%u_adj => simulation%adjoint_fluid%u_adj
205 this%v_adj => simulation%adjoint_fluid%v_adj
206 this%w_adj => simulation%adjoint_fluid%w_adj
207 end subroutine scalar_mixing_init_attributes
210 subroutine scalar_mixing_free(this)
212 call this%free_base()
222 end subroutine scalar_mixing_free
227 subroutine scalar_mixing_update_value(this, design)
229 class(design_t),
intent(in) :: design
230 type(field_t),
pointer :: work
231 integer :: temp_indices(1), n
239 call neko_scratch_registry%request_field(work, temp_indices(1))
243 call field_copy(work, this%phi)
245 call field_cadd(work, -this%phi_ref)
247 call field_col2(work, work)
249 if (this%has_mask)
then
250 call mask_exterior_const(work, this%mask, 0.0_rp)
253 if (neko_bcknd_device .eq. 1)
then
254 this%value = device_glsc2(work%x_d, this%coef%B_d, n)
256 this%value = glsc2(work%x, this%coef%B, n)
260 this%value = 0.5_rp * this%value / this%domain_volume
263 call neko_scratch_registry%relinquish_field(temp_indices)
264 end subroutine scalar_mixing_update_value
270 subroutine scalar_mixing_update_sensitivity(this, design)
272 class(design_t),
intent(in) :: design
276 type(field_t),
pointer :: work
277 integer :: temp_indices(1)
279 call neko_scratch_registry%request_field(work, temp_indices(1))
282 call field_col3(work, this%u, this%u_adj)
283 call field_addcol3(work, this%v, this%v_adj)
284 call field_addcol3(work, this%w, this%w_adj)
286 call field_cmult(work, -1.0_rp)
288 if (neko_bcknd_device .eq. 1)
then
289 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
291 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
294 call neko_scratch_registry%relinquish_field(temp_indices)
296 end subroutine scalar_mixing_update_sensitivity
Implements the adjoint_mixing_scalar_source_term type.
Some common Masking operations we may need.
Contains extensions to the neko library required to run the topology optimization code.
subroutine, public get_scalar_indicies(i_primal, i_adjoint, scalars, adjoint_scalars, primal_name)
get scalar indices
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 .