Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
scalar_mixing_objective_function.f90
Go to the documentation of this file.
1
34!
39 use num_types, only: rp
40 use objective, only: objective_t
41 use simulation_m, only: simulation_t
42 use design, only: design_t
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
57 ! delete after the simulation computes u u_adj
58 use field_math, only: field_addcol3, field_col3
59 use continuation_scheduler, only: nekotop_continuation
60 implicit none
61 private
62
66 type, public, extends(objective_t) :: scalar_mixing_objective_t
67 private
68
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
79
80 ! -----------------------------------------------------------------
81 ! THESE SHOULD BE DELETED WHEN THE DESIGN UPDATE COMES IN.
82 ! These sensitivity should return zero, and the simulation
83 ! should compute the u u_adj component
84 type(field_t), pointer :: u, v, w, u_adj, v_adj, w_adj
85
86 contains
87
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
102
104
105contains
106
112 subroutine scalar_mixing_init_json_sim(this, json, design, simulation)
113 class(scalar_mixing_objective_t), intent(inout) :: this
114 type(json_file), intent(inout) :: json
115 class(design_t), intent(in) :: design
116 type(simulation_t), target, intent(inout) :: simulation
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
123
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))
132
133 ! initialize
134 call this%init_from_attributes(design, simulation, weight, name, &
135 mask_name, phi_ref, scalar_name, start_time, end_time)
136 end subroutine scalar_mixing_init_json_sim
137
149 subroutine scalar_mixing_init_attributes(this, design, simulation, weight, &
150 name, mask_name, phi_ref, scalar_name, start_time, end_time)
151 class(scalar_mixing_objective_t), intent(inout) :: this
152 class(design_t), intent(in) :: design
153 type(simulation_t), target, intent(inout) :: simulation
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
161
162 type(adjoint_mixing_scalar_source_term_t) :: adjoint_forcing
163 integer :: i_scalar, i_adjoint_scalar
164
165 ! Start by checking if the adjoint scalar has been initialized
166 if (.not.allocated(simulation%adjoint_case%adjoint_scalars)) then
167 call neko_error("adjoint passive scalar not initialized")
168 end if
169
170 ! Call the base initializer
171 call this%init_base(name, design%size(), weight, mask_name, &
172 start_time, end_time)
173
174 ! Associate the integration weights
175 this%coef => simulation%fluid%c_Xh
176
177 ! Compute the masked volume
178 if (this%has_mask) then
179 this%domain_volume = compute_masked_volume(this%mask, this%coef)
180 else
181 this%domain_volume = this%coef%volume
182 end if
183
184 this%scalar_name = trim(scalar_name)
185
186 ! figure out the index associated with the scalar and adjoint scalar.
187 call get_scalar_indicies(i_scalar, i_adjoint_scalar, simulation%scalars, &
188 simulation%adjoint_scalars, this%scalar_name)
189
192 associate(f_phi_adj => simulation%adjoint_scalars% &
193 adjoint_scalar_fields(i_adjoint_scalar)%f_Xh)
194
195 ! Associate json parameters
196 this%phi_ref = phi_ref
197
198 ! Associate forward passive scalar
199 this%phi => simulation%scalars%scalar_fields(i_scalar)%scalar%s
200
201 ! Initialize the scalar mixing adjoint source term
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)
205
206 end associate
207
208 ! append adjoint source term to the adjoint passive scalar equation
209 call simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar) &
210 %source_term%add_source_term(adjoint_forcing)
211
212 !--------------------------------------------------------------------------
213 ! THIS SHOULD BE REPLACED WHEN THE DESIGN UPDATE OCCURS
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
221
223 subroutine scalar_mixing_free(this)
224 class(scalar_mixing_objective_t), intent(inout) :: this
225 call this%free_base()
226
227 this%u => null()
228 this%v => null()
229 this%w => null()
230 this%phi => null()
231 this%u_adj => null()
232 this%v_adj => null()
233 this%w_adj => null()
234
235 end subroutine scalar_mixing_free
236
240 subroutine scalar_mixing_update_value(this, design)
241 class(scalar_mixing_objective_t), intent(inout) :: this
242 class(design_t), intent(in) :: design
243 type(field_t), pointer :: work
244 integer :: temp_indices(1), n
245
246
247 ! The objective being computed is
248 !
249 ! \f$1/2 \int_\Omega_{obj} (\phi - \phi_ref)^2 d\Omega\f$
250
251 ! get a working array
252 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
253 n = work%size()
254
255 ! \f$ \phi \f$
256 call field_copy(work, this%phi)
257 ! \f$ \phi - \phi_ref \f$
258 call field_cadd(work, -this%phi_ref)
259 ! \f$ (\phi - \phi_ref)^2 \f$
260 call field_col2(work, work)
261 ! \f$ mask to \Omega_{obj} \f$
262 if (this%has_mask) then
263 call mask_exterior_const(work, this%mask, 0.0_rp)
264 end if
265 ! integrate
266 if (neko_bcknd_device .eq. 1) then
267 this%value = device_glsc2(work%x_d, this%coef%B_d, n)
268 else
269 this%value = glsc2(work%x, this%coef%B, n)
270 end if
271
272 ! \f$ scale by 1/2 and 1/|\Omega_{obj}| \f$
273 this%value = 0.5_rp * this%value / this%domain_volume
274
275 ! relingush the scratch field
276 call neko_scratch_registry%relinquish_field(temp_indices)
277 end subroutine scalar_mixing_update_value
278
283 subroutine scalar_mixing_update_sensitivity(this, design)
284 class(scalar_mixing_objective_t), intent(inout) :: this
285 class(design_t), intent(in) :: design
286
287 end subroutine scalar_mixing_update_sensitivity
288
Implements the adjoint_mixing_scalar_source_term type.
Continuation scheduler for the optimization loop.
Implements the design_t.
Definition design.f90:36
Some common Masking operations we may need.
Definition mask_ops.f90:36
real(kind=rp) function, public compute_masked_volume(mask, coef)
Compute the volume of the domain contained within the mask.
Definition mask_ops.f90:182
Contains extensions to the neko library required to run the topology optimization code.
Definition neko_ext.f90:42
subroutine, public get_scalar_indicies(i_primal, i_adjoint, scalars, adjoint_scalars, primal_name)
get scalar indices
Definition neko_ext.f90:435
Implements the objective_t type.
Definition objective.f90:36
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.
An abstract design type.
Definition design.f90:53
The abstract objective type.
Definition objective.f90:52
An objective function corresponding to the mixing of a passive scalar .