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 implicit none
60 private
61
65 type, public, extends(objective_t) :: scalar_mixing_objective_t
66 private
67
69 type(field_t), pointer :: phi
71 real(kind=rp) :: phi_ref
73 class(coef_t), pointer :: coef
75 real(kind=rp) :: domain_volume
77 character(len=:), allocatable :: scalar_name
78
79 ! -----------------------------------------------------------------
80 ! THESE SHOULD BE DELETED WHEN THE DESIGN UPDATE COMES IN.
81 ! These sensitivity should return zero, and the simulation
82 ! should compute the u u_adj component
83 type(field_t), pointer :: u, v, w, u_adj, v_adj, w_adj
84
85 contains
86
88 procedure, public, pass(this) :: init_json_sim => &
91 procedure, public, pass(this) :: init_from_attributes => &
92 scalar_mixing_init_attributes
94 procedure, public, pass(this) :: free => scalar_mixing_free
96 procedure, public, pass(this) :: update_value => &
97 scalar_mixing_update_value
99 procedure, public, pass(this) :: update_sensitivity => &
100 scalar_mixing_update_sensitivity
101
103
104contains
105
111 subroutine scalar_mixing_init_json_sim(this, json, design, simulation)
112 class(scalar_mixing_objective_t), intent(inout) :: this
113 type(json_file), intent(inout) :: json
114 class(design_t), intent(in) :: design
115 type(simulation_t), target, intent(inout) :: simulation
116 real(kind=rp) :: weight
117 real(kind=rp) :: phi_ref
118 character(len=:), allocatable :: name
119 character(len=:), allocatable :: mask_name
120 character(len=:), allocatable :: scalar_name
121
122 call json_get_or_default(json, "weight", weight, 1.0_rp)
123 call json_get_or_default(json, "mask_name", mask_name, "")
124 call json_get_or_default(json, "target_concentration", phi_ref, 0.5_rp)
125 call json_get_or_default(json, "name", name, "Scalar Mixing")
126 call json_get_or_default(json, "scalar_name", scalar_name, "s")
127
128 ! initialize
129 call this%init_from_attributes(design, simulation, weight, name, &
130 mask_name, phi_ref, scalar_name)
131 end subroutine scalar_mixing_init_json_sim
132
142 subroutine scalar_mixing_init_attributes(this, design, simulation, weight, &
143 name, mask_name, phi_ref, scalar_name)
144 class(scalar_mixing_objective_t), intent(inout) :: this
145 class(design_t), intent(in) :: design
146 type(simulation_t), target, intent(inout) :: simulation
147 real(kind=rp), intent(in) :: weight
148 real(kind=rp), intent(in) :: phi_ref
149 character(len=*), intent(in) :: mask_name
150 character(len=*), intent(in) :: name
151 character(len=*), intent(in) :: scalar_name
152
153 type(adjoint_mixing_scalar_source_term_t) :: adjoint_forcing
154 integer :: i_scalar, i_adjoint_scalar
155
156 ! Start by checking if the adjoint scalar has been initialized
157 if (.not.allocated(simulation%adjoint_case%adjoint_scalars)) then
158 call neko_error("adjoint passive scalar not initialized")
159 end if
160
161 ! Call the base initializer
162 call this%init_base(name, design%size(), weight, mask_name)
163
164 ! Associate the integration weights
165 this%coef => simulation%fluid%c_Xh
166
167 ! Compute the masked volume
168 if (this%has_mask) then
169 this%domain_volume = compute_masked_volume(this%mask, this%coef)
170 else
171 this%domain_volume = this%coef%volume
172 end if
173
174 this%scalar_name = trim(scalar_name)
175
176 ! figure out the index associated with the scalar and adjoint scalar.
177 call get_scalar_indicies(i_scalar, i_adjoint_scalar, simulation%scalars, &
178 simulation%adjoint_scalars, this%scalar_name)
179
182 associate(f_phi_adj => &
183 simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar)%f_Xh)
184
185 ! Associate json parameters
186 this%phi_ref = phi_ref
187
188 ! Associate forward passive scalar
189 this%phi => simulation%scalars%scalar_fields(i_scalar)%s
190
191 ! Initialize the scalar mixing adjoint source term
192 call adjoint_forcing%init_from_components(f_phi_adj, this%phi, &
193 this%get_weight(), this%phi_ref, this%mask, this%has_mask, this%coef)
194
195 end associate
196
197 ! append adjoint source term to the adjoint passive scalar equation
198 call simulation%adjoint_scalars%adjoint_scalar_fields(i_adjoint_scalar) &
199 %source_term%add_source_term(adjoint_forcing)
200
201 !--------------------------------------------------------------------------
202 ! THIS SHOULD BE REPLACED WHEN THE DESIGN UPDATE OCCURS
203 this%u => simulation%fluid%u
204 this%v => simulation%fluid%v
205 this%w => simulation%fluid%w
206 this%u_adj => simulation%adjoint_fluid%u_adj
207 this%v_adj => simulation%adjoint_fluid%v_adj
208 this%w_adj => simulation%adjoint_fluid%w_adj
209 end subroutine scalar_mixing_init_attributes
210
212 subroutine scalar_mixing_free(this)
213 class(scalar_mixing_objective_t), intent(inout) :: this
214 call this%free_base()
215
216 this%u => null()
217 this%v => null()
218 this%w => null()
219 this%phi => null()
220 this%u_adj => null()
221 this%v_adj => null()
222 this%w_adj => null()
223
224 end subroutine scalar_mixing_free
225
229 subroutine scalar_mixing_update_value(this, design)
230 class(scalar_mixing_objective_t), intent(inout) :: this
231 class(design_t), intent(in) :: design
232 type(field_t), pointer :: work
233 integer :: temp_indices(1), n
234
235
236 ! The objective being computed is
237 !
238 ! \f$1/2 \int_\Omega_{obj} (\phi - \phi_ref)^2 d\Omega\f$
239
240 ! get a working array
241 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
242 n = work%size()
243
244 ! \f$ \phi \f$
245 call field_copy(work, this%phi)
246 ! \f$ \phi - \phi_ref \f$
247 call field_cadd(work, -this%phi_ref)
248 ! \f$ (\phi - \phi_ref)^2 \f$
249 call field_col2(work, work)
250 ! \f$ mask to \Omega_{obj} \f$
251 if (this%has_mask) then
252 call mask_exterior_const(work, this%mask, 0.0_rp)
253 end if
254 ! integrate
255 if (neko_bcknd_device .eq. 1) then
256 this%value = device_glsc2(work%x_d, this%coef%B_d, n)
257 else
258 this%value = glsc2(work%x, this%coef%B, n)
259 end if
260
261 ! \f$ scale by 1/2 and 1/|\Omega_{obj}| \f$
262 this%value = 0.5_rp * this%value / this%domain_volume
263
264 ! relingush the scratch field
265 call neko_scratch_registry%relinquish_field(temp_indices)
266 end subroutine scalar_mixing_update_value
267
272 subroutine scalar_mixing_update_sensitivity(this, design)
273 class(scalar_mixing_objective_t), intent(inout) :: this
274 class(design_t), intent(in) :: design
275
276 end subroutine scalar_mixing_update_sensitivity
277
Implements the adjoint_mixing_scalar_source_term type.
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:432
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:54
The abstract objective type.
Definition objective.f90:52
An objective function corresponding to the mixing of a passive scalar .