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
54 ! delete after the simulation computes u u_adj
55 use field_math, only: field_addcol3, field_col3
56 implicit none
57 private
58
62 type, public, extends(objective_t) :: scalar_mixing_objective_t
63 private
64
66 type(field_t), pointer :: phi
68 real(kind=rp) :: phi_ref
70 class(coef_t), pointer :: coef
72 real(kind=rp) :: domain_volume
73
74 ! -----------------------------------------------------------------
75 ! THESE SHOULD BE DELETED WHEN THE DESIGN UPDATE COMES IN.
76 ! These sensitivity should return zero, and the simulation
77 ! should compute the u u_adj component
78 type(field_t), pointer :: u, v, w, u_adj, v_adj, w_adj
79
80 contains
81
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
96
98
99contains
100
106 subroutine scalar_mixing_init_json_sim(this, json, design, simulation)
107 class(scalar_mixing_objective_t), intent(inout) :: this
108 type(json_file), intent(inout) :: json
109 class(design_t), intent(in) :: design
110 type(simulation_t), target, intent(inout) :: simulation
111 real(kind=rp) :: weight
112 real(kind=rp) :: phi_ref
113 character(len=:), allocatable :: name
114 character(len=:), allocatable :: mask_name
115
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")
120
121 ! initialize
122 call this%init_from_attributes(design, simulation, weight, name, &
123 mask_name, phi_ref)
124 end subroutine scalar_mixing_init_json_sim
125
134 subroutine scalar_mixing_init_attributes(this, design, simulation, weight, &
135 name, mask_name, phi_ref)
136 class(scalar_mixing_objective_t), intent(inout) :: this
137 class(design_t), intent(in) :: design
138 type(simulation_t), target, intent(inout) :: simulation
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
143 type(adjoint_mixing_scalar_source_term_t) :: adjoint_forcing
144
145 ! Start by checking if the adjoint scalar has been initialized
146 if (.not.allocated(simulation%adjoint_case%scalar_adj)) then
147 call neko_error("adjoint passive scalar not initialized")
148 end if
149
150 ! Call the base initializer
151 call this%init_base(name, design%size(), weight, mask_name)
152
153 ! Associate the integration weights
154 this%coef => simulation%fluid%c_Xh
155
156 ! Compute the masked volume
157 if (this%has_mask) then
158 this%domain_volume = compute_masked_volume(this%mask, this%coef)
159 else
160 this%domain_volume = this%coef%volume
161 end if
162
165 associate(f_phi_adj => simulation%adjoint_scalar%f_Xh)
166
167 ! Associate json parameters
168 this%phi_ref = phi_ref
169
170 ! Associate forward passive scalar
171 this%phi => simulation%scalar%s
172
173 ! Initialize the scalar mixing adjoint source term
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)
176
177 end associate
178
179 ! append adjoint source term to the adjoint passive scalar equation
180 call simulation%adjoint_scalar%source_term%add_source_term(adjoint_forcing)
181
182 !--------------------------------------------------------------------------
183 ! THIS SHOULD BE REPLACED WHEN THE DESIGN UPDATE OCCURS
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
191
193 subroutine scalar_mixing_free(this)
194 class(scalar_mixing_objective_t), intent(inout) :: this
195 call this%free_base()
196
197 this%u => null()
198 this%v => null()
199 this%w => null()
200 this%phi => null()
201 this%u_adj => null()
202 this%v_adj => null()
203 this%w_adj => null()
204
205 end subroutine scalar_mixing_free
206
210 subroutine scalar_mixing_update_value(this, design)
211 class(scalar_mixing_objective_t), intent(inout) :: this
212 class(design_t), intent(in) :: design
213 type(field_t), pointer :: work
214 integer :: temp_indices(1), n
215
216
217 ! The objective being computed is
218 !
219 ! \f$1/2 \int_\Omega_{obj} (\phi - \phi_ref)^2 d\Omega\f$
220
221 ! get a working array
222 call neko_scratch_registry%request_field(work, temp_indices(1))
223 n = work%size()
224
225 ! \f$ \phi \f$
226 call field_copy(work, this%phi)
227 ! \f$ \phi - \phi_ref \f$
228 call field_cadd(work, -this%phi_ref)
229 ! \f$ (\phi - \phi_ref)^2 \f$
230 call field_col2(work, work)
231 ! \f$ mask to \Omega_{obj} \f$
232 if (this%has_mask) then
233 call mask_exterior_const(work, this%mask, 0.0_rp)
234 end if
235 ! integrate
236 if (neko_bcknd_device .eq. 1) then
237 this%value = device_glsc2(work%x_d, this%coef%B_d, n)
238 else
239 this%value = glsc2(work%x, this%coef%B, n)
240 end if
241
242 ! \f$ scale by 1/2 and 1/|\Omega_{obj}| \f$
243 this%value = 0.5_rp * this%value / this%domain_volume
244
245 ! relingush the scratch field
246 call neko_scratch_registry%relinquish_field(temp_indices)
247 end subroutine scalar_mixing_update_value
248
253 subroutine scalar_mixing_update_sensitivity(this, design)
254 class(scalar_mixing_objective_t), intent(inout) :: this
255 class(design_t), intent(in) :: design
256
257 !--------------------------------------------------------------------------
258 ! THIS SHOULD BE REPLACED WITH A ZERO CONTRIBUTION AFTER THE DESIGN UPDATE
259 type(field_t), pointer :: work
260 integer :: temp_indices(1)
261
262 call neko_scratch_registry%request_field(work, temp_indices(1))
263
264 ! here it should just be an inner product between the forward and adjoint
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)
268 ! but negative
269 call field_cmult(work, -1.0_rp)
270
271 if (neko_bcknd_device .eq. 1) then
272 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
273 else
274 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
275 end if
276
277 call neko_scratch_registry%relinquish_field(temp_indices)
278
279 end subroutine scalar_mixing_update_sensitivity
280
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
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 .