39 use num_types,
only: rp
43 use json_module,
only: json_file
44 use json_utils,
only: json_get_or_default
45 use field,
only: field_t
46 use field_math,
only: field_copy, field_cadd, field_col2, field_cmult
47 use neko_config,
only: neko_bcknd_device
48 use math,
only: glsc2, copy
49 use device_math,
only: device_glsc2, device_copy
51 use coefs,
only: coef_t
52 use scratch_registry,
only: neko_scratch_registry
53 use utils,
only: neko_error
58 use field_math,
only: field_addcol3, field_col3
70 type(field_t),
pointer :: phi
72 real(kind=rp) :: phi_ref
74 class(coef_t),
pointer :: coef
76 real(kind=rp) :: domain_volume
78 character(len=:),
allocatable :: scalar_name
84 type(field_t),
pointer :: u, v, w, u_adj, v_adj, w_adj
89 procedure,
public, pass(this) :: init_json_sim => &
92 procedure,
public, pass(this) :: init_from_attributes => &
93 scalar_mixing_init_attributes
95 procedure,
public, pass(this) :: free => scalar_mixing_free
97 procedure,
public, pass(this) :: update_value => &
98 scalar_mixing_update_value
100 procedure,
public, pass(this) :: update_sensitivity => &
101 scalar_mixing_update_sensitivity
114 type(json_file),
intent(inout) :: json
115 class(
design_t),
intent(in) :: design
117 real(kind=rp) :: weight
118 real(kind=rp) :: phi_ref
119 real(kind=rp) :: start_time, end_time
120 character(len=:),
allocatable :: name
121 character(len=:),
allocatable :: mask_name
122 character(len=:),
allocatable :: scalar_name
124 call nekotop_continuation%json_get_or_register(json,
'weight', &
125 this%weight, weight, 1.0_rp)
126 call json_get_or_default(json,
"mask_name", mask_name,
"")
127 call json_get_or_default(json,
"target_concentration", phi_ref, 0.5_rp)
128 call json_get_or_default(json,
"name", name,
"Scalar Mixing")
129 call json_get_or_default(json,
"scalar_name", scalar_name,
"s")
130 call json_get_or_default(json,
"start_time", start_time, 0.0_rp)
131 call json_get_or_default(json,
"end_time", end_time, huge(0.0_rp))
134 call this%init_from_attributes(
design, simulation, weight, name, &
135 mask_name, phi_ref, scalar_name, start_time, end_time)
149 subroutine scalar_mixing_init_attributes(this, design, simulation, weight, &
150 name, mask_name, phi_ref, scalar_name, start_time, end_time)
152 class(
design_t),
intent(in) :: design
154 real(kind=rp),
intent(in) :: weight
155 character(len=*),
intent(in) :: name
156 character(len=*),
intent(in) :: mask_name
157 real(kind=rp),
intent(in) :: phi_ref
158 character(len=*),
intent(in) :: scalar_name
159 real(kind=rp),
intent(in) :: start_time
160 real(kind=rp),
intent(in) :: end_time
163 integer :: i_scalar, i_adjoint_scalar
166 if (.not.
allocated(simulation%adjoint_case%adjoint_scalars))
then
167 call neko_error(
"adjoint passive scalar not initialized")
171 call this%init_base(name,
design%size(), weight, mask_name, &
172 start_time, end_time)
175 this%coef => simulation%fluid%c_Xh
178 if (this%has_mask)
then
181 this%domain_volume = this%coef%volume
184 this%scalar_name = trim(scalar_name)
188 simulation%adjoint_scalars, this%scalar_name)
192 associate(f_phi_adj => simulation%adjoint_scalars% &
193 adjoint_scalar_fields(i_adjoint_scalar)%f_Xh)
196 this%phi_ref = phi_ref
199 this%phi => simulation%scalars%scalar_fields(i_scalar)%scalar%s
202 call adjoint_forcing%init_from_components(f_phi_adj, this%phi, &
203 this%get_weight(), this%phi_ref, this%mask, this%has_mask, &
204 this%coef, this%start_time, this%end_time)
209 call simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar) &
210 %source_term%add_source_term(adjoint_forcing)
214 this%u => simulation%fluid%u
215 this%v => simulation%fluid%v
216 this%w => simulation%fluid%w
217 this%u_adj => simulation%adjoint_fluid%u_adj
218 this%v_adj => simulation%adjoint_fluid%v_adj
219 this%w_adj => simulation%adjoint_fluid%w_adj
220 end subroutine scalar_mixing_init_attributes
223 subroutine scalar_mixing_free(this)
225 call this%free_base()
235 end subroutine scalar_mixing_free
240 subroutine scalar_mixing_update_value(this, design)
242 class(design_t),
intent(in) :: design
243 type(field_t),
pointer :: work
244 integer :: temp_indices(1), n
252 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
256 call field_copy(work, this%phi)
258 call field_cadd(work, -this%phi_ref)
260 call field_col2(work, work)
262 if (this%has_mask)
then
263 call mask_exterior_const(work, this%mask, 0.0_rp)
266 if (neko_bcknd_device .eq. 1)
then
267 this%value = device_glsc2(work%x_d, this%coef%B_d, n)
269 this%value = glsc2(work%x, this%coef%B, n)
273 this%value = 0.5_rp * this%value / this%domain_volume
276 call neko_scratch_registry%relinquish_field(temp_indices)
277 end subroutine scalar_mixing_update_value
283 subroutine scalar_mixing_update_sensitivity(this, design)
285 class(design_t),
intent(in) :: design
287 end subroutine scalar_mixing_update_sensitivity
Implements the adjoint_mixing_scalar_source_term type.
Continuation scheduler for the optimization loop.
Some common Masking operations we may need.
real(kind=rp) function, public compute_masked_volume(mask, coef)
Compute the volume of the domain contained within the mask.
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 .
subroutine scalar_mixing_init_json_sim(this, json, design, simulation)
The common constructor using a JSON object.
Implements the steady_problem_t type.
The abstract objective type.
An objective function corresponding to the mixing of a passive scalar .