Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
scalar_mixing_objective_function.f90
1! Copyright (c) 2023, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
37 use num_types, only: rp
38 use objective, only: objective_t
39 use simulation_m, only: simulation_t
40 use design, only: design_t
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
48 use mask_ops, only: mask_exterior_const, compute_masked_volume
49 use coefs, only: coef_t
50 use scratch_registry, only: neko_scratch_registry
51 use utils, only: neko_error
55 ! delete after the simulation computes u u_adj
56 use field_math, only: field_addcol3, field_col3
57 implicit none
58 private
59
63 type, public, extends(objective_t) :: scalar_mixing_objective_t
64 private
65
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
76
77 ! -----------------------------------------------------------------
78 ! THESE SHOULD BE DELETED WHEN THE DESIGN UPDATE COMES IN.
79 ! These sensitivity should return zero, and the simulation
80 ! should compute the u u_adj component
81 type(field_t), pointer :: u, v, w, u_adj, v_adj, w_adj
82
83 contains
84
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
99
101
102contains
103
109 subroutine scalar_mixing_init_json_sim(this, json, design, simulation)
110 class(scalar_mixing_objective_t), intent(inout) :: this
111 type(json_file), intent(inout) :: json
112 class(design_t), intent(in) :: design
113 type(simulation_t), target, intent(inout) :: simulation
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
119
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")
125
126 ! initialize
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
130
140 subroutine scalar_mixing_init_attributes(this, design, simulation, weight, &
141 name, mask_name, phi_ref, scalar_name)
142 class(scalar_mixing_objective_t), intent(inout) :: this
143 class(design_t), intent(in) :: design
144 type(simulation_t), target, intent(inout) :: simulation
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
150
151 type(adjoint_mixing_scalar_source_term_t) :: adjoint_forcing
152 integer :: i_scalar, i_adjoint_scalar
153
154 ! Start by checking if the adjoint scalar has been initialized
155 if (.not.allocated(simulation%adjoint_case%adjoint_scalars)) then
156 call neko_error("adjoint passive scalar not initialized")
157 end if
158
159 ! Call the base initializer
160 call this%init_base(name, design%size(), weight, mask_name)
161
162 ! Associate the integration weights
163 this%coef => simulation%fluid%c_Xh
164
165 ! Compute the masked volume
166 if (this%has_mask) then
167 this%domain_volume = compute_masked_volume(this%mask, this%coef)
168 else
169 this%domain_volume = this%coef%volume
170 end if
171
172 this%scalar_name = trim(scalar_name)
173
174 ! figure out the index associated with the scalar and adjoint scalar.
175 call get_scalar_indicies(i_scalar, i_adjoint_scalar, simulation%scalars, &
176 simulation%adjoint_scalars, this%scalar_name)
177
180 associate(f_phi_adj => &
181 simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar)%f_Xh)
182
183 ! Associate json parameters
184 this%phi_ref = phi_ref
185
186 ! Associate forward passive scalar
187 this%phi => simulation%scalars%scalar_fields(i_scalar)%s
188
189 ! Initialize the scalar mixing adjoint source term
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)
192
193 end associate
194
195 ! append adjoint source term to the adjoint passive scalar equation
196 call simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar) &
197 %source_term%add_source_term(adjoint_forcing)
198
199 !--------------------------------------------------------------------------
200 ! THIS SHOULD BE REPLACED WHEN THE DESIGN UPDATE OCCURS
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
208
210 subroutine scalar_mixing_free(this)
211 class(scalar_mixing_objective_t), intent(inout) :: this
212 call this%free_base()
213
214 this%u => null()
215 this%v => null()
216 this%w => null()
217 this%phi => null()
218 this%u_adj => null()
219 this%v_adj => null()
220 this%w_adj => null()
221
222 end subroutine scalar_mixing_free
223
227 subroutine scalar_mixing_update_value(this, design)
228 class(scalar_mixing_objective_t), intent(inout) :: this
229 class(design_t), intent(in) :: design
230 type(field_t), pointer :: work
231 integer :: temp_indices(1), n
232
233
234 ! The objective being computed is
235 !
236 ! \f$1/2 \int_\Omega_{obj} (\phi - \phi_ref)^2 d\Omega\f$
237
238 ! get a working array
239 call neko_scratch_registry%request_field(work, temp_indices(1))
240 n = work%size()
241
242 ! \f$ \phi \f$
243 call field_copy(work, this%phi)
244 ! \f$ \phi - \phi_ref \f$
245 call field_cadd(work, -this%phi_ref)
246 ! \f$ (\phi - \phi_ref)^2 \f$
247 call field_col2(work, work)
248 ! \f$ mask to \Omega_{obj} \f$
249 if (this%has_mask) then
250 call mask_exterior_const(work, this%mask, 0.0_rp)
251 end if
252 ! integrate
253 if (neko_bcknd_device .eq. 1) then
254 this%value = device_glsc2(work%x_d, this%coef%B_d, n)
255 else
256 this%value = glsc2(work%x, this%coef%B, n)
257 end if
258
259 ! \f$ scale by 1/2 and 1/|\Omega_{obj}| \f$
260 this%value = 0.5_rp * this%value / this%domain_volume
261
262 ! relingush the scratch field
263 call neko_scratch_registry%relinquish_field(temp_indices)
264 end subroutine scalar_mixing_update_value
265
270 subroutine scalar_mixing_update_sensitivity(this, design)
271 class(scalar_mixing_objective_t), intent(inout) :: this
272 class(design_t), intent(in) :: design
273
274 !--------------------------------------------------------------------------
275 ! THIS SHOULD BE REPLACED WITH A ZERO CONTRIBUTION AFTER THE DESIGN UPDATE
276 type(field_t), pointer :: work
277 integer :: temp_indices(1)
278
279 call neko_scratch_registry%request_field(work, temp_indices(1))
280
281 ! here it should just be an inner product between the forward and adjoint
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)
285 ! but negative
286 call field_cmult(work, -1.0_rp)
287
288 if (neko_bcknd_device .eq. 1) then
289 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
290 else
291 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
292 end if
293
294 call neko_scratch_registry%relinquish_field(temp_indices)
295
296 end subroutine scalar_mixing_update_sensitivity
297
Implements the adjoint_mixing_scalar_source_term type.
Implements the design_t.
Definition design.f90:34
Some common Masking operations we may need.
Definition mask_ops.f90:34
Contains extensions to the neko library required to run the topology optimization code.
Definition neko_ext.f90:9
subroutine, public get_scalar_indicies(i_primal, i_adjoint, scalars, adjoint_scalars, primal_name)
get scalar indices
Definition neko_ext.f90:284
Implements the objective_t type.
Definition objective.f90:35
An objective function corresponding to the mixing of a passive scalar .
Implements the steady_problem_t type.
An abstract design type.
Definition design.f90:50
The abstract objective type.
Definition objective.f90:51
An objective function corresponding to the mixing of a passive scalar .