Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
volume_constraint.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!
34! $V = \frac{1}{|\Omega_O|}\int_{\Omega_O} \tilde{\rho} d\Omega$
35! Either
36! $V < V_\text{max}$
37! $V > V_\text{min}$
39 use constraint, only: constraint_t
40
41 use design, only: design_t
42 use brinkman_design, only: brinkman_design_t
43 use simulation_m, only: simulation_t
44
45 use num_types, only: rp
46 use coefs, only: coef_t
47 use json_module, only: json_file
48 use json_utils, only: json_get, json_get_or_default
49 use field, only: field_t
50 use field_registry, only: neko_field_registry
51 use scratch_registry, only: neko_scratch_registry
52 use neko_config, only: neko_bcknd_device
54 use math, only: glsc2
55 use device_math, only: device_glsc2
56 use math_ext, only: glsc2_mask
57 use field_math, only: field_rone, field_copy
58 use utils, only: neko_error
59 implicit none
60 private
61
63 type, public, extends(constraint_t) :: volume_constraint_t
64 private
65
69 logical :: is_max
71 real(kind=rp) :: limit
73 real(kind=rp) :: volume_domain
74
76 class(coef_t), pointer :: c_xh => null()
77
78 contains
79
81 procedure, public, pass(this) :: init_json_sim => &
82 volume_constraint_init_json_sim
84 procedure, public, pass(this) :: init_from_attributes => &
85 volume_constraint_init_attributes
87 procedure, public, pass(this) :: free => volume_constraint_free
89 procedure, public, pass(this) :: update_value => &
90 volume_constraint_update_value
92 procedure, public, pass(this) :: update_sensitivity => &
93 volume_constraint_update_sensitivity
94
96 procedure, private, pass(this) :: compute_volume
97
98 end type volume_constraint_t
99
100contains
101
107 subroutine volume_constraint_init_json_sim(this, json, design, simulation)
108 class(volume_constraint_t), intent(inout) :: this
109 type(json_file), intent(inout) :: json
110 class(design_t), intent(in) :: design
111 type(simulation_t), target, intent(inout) :: simulation
112
113 character(len=:), allocatable :: mask_name
114 character(len=:), allocatable :: name
115 logical :: is_max
116 real(kind=rp) :: limit
117
118 call json_get_or_default(json, "mask_name", mask_name, "")
119 call json_get_or_default(json, "name", name, "Volume constraint")
120 call json_get_or_default(json, "is_max", is_max, .false.)
121 call json_get(json, "limit", limit)
122
123 call this%init_from_attributes(design, simulation, name, mask_name, &
124 is_max, limit)
125 end subroutine volume_constraint_init_json_sim
126
135 subroutine volume_constraint_init_attributes(this, design, simulation, &
136 name, mask_name, is_max, limit)
137 class(volume_constraint_t), intent(inout) :: this
138 class(design_t), intent(in) :: design
139 type(simulation_t), target, intent(inout) :: simulation
140 character(len=*), intent(in) :: mask_name
141 character(len=*), intent(in) :: name
142 logical, intent(in) :: is_max
143 real(kind=rp), intent(in) :: limit
144
145 real(kind=rp) :: volume
146 type(field_t), pointer :: work
147 integer :: temp_indices(1)
148
149 ! Initialize the base class
150 call this%init_base(name, design%size(), mask_name)
151
152 ! Store the attributes
153 this%is_max = is_max
154 this%limit = limit
155 this%c_Xh => simulation%neko_case%fluid%c_Xh
156
157 ! Now we can extract the mask/has_mask from the design
158 if (this%has_mask) then
159
160 ! calculate the volume of the optimization domain
161 call neko_scratch_registry%request_field(work, temp_indices(1))
162 call field_rone(work)
163
164 if (neko_bcknd_device .eq. 1) then
165 call mask_exterior_const(work, this%mask, 0.0_rp)
166 this%volume_domain = device_glsc2(work%x_d, this%c_xh%B_d, &
167 work%size())
168 else
169 this%volume_domain = glsc2_mask(work%x, this%c_Xh%B, &
170 design%size(), this%mask%mask, this%mask%size)
171 end if
172
173 call neko_scratch_registry%relinquish_field(temp_indices)
174 else
175 this%volume_domain = this%c_Xh%volume
176 end if
177
178 ! ------------------------------------------------------------------------ !
179 ! Initialize the value of constraint
180
181 ! Compute the volume of the design
182 volume = this%compute_volume(design)
183
184 ! Compute the distance to the target volume
185 this%value = this%limit - volume / this%volume_domain
186
187 ! Invert the sign if it is a maximum constraint
188 if (this%is_max) this%value = - ( this%value )
189
190 ! ------------------------------------------------------------------------ !
191 ! Initialize the sensitivity value
192
193 this%sensitivity = 1.0_rp / this%volume_domain
194
195 ! Invert the sign if it is a maximum constraint
196 if (.not. this%is_max) this%sensitivity = (-1.0_rp) * this%sensitivity
197
198 if (this%has_mask) then
199 call mask_exterior_const(this%sensitivity, this%mask, 0.0_rp)
200 end if
201
202 end subroutine volume_constraint_init_attributes
203
205 subroutine volume_constraint_free(this)
206 class(volume_constraint_t), intent(inout) :: this
207
208 call this%free_base()
209 end subroutine volume_constraint_free
210
214 subroutine volume_constraint_update_value(this, design)
215 class(volume_constraint_t), intent(inout) :: this
216 class(design_t), intent(in) :: design
217 real(kind=rp) :: volume
218
219 volume = this%compute_volume(design)
220
221 ! Compute the distance to the target volume
222 this%value = this%limit - volume / this%volume_domain
223
224 ! Invert the sign if it is a maximum constraint
225 if (this%is_max) this%value = - ( this%value )
226
227 end subroutine volume_constraint_update_value
228
232 subroutine volume_constraint_update_sensitivity(this, design)
233 class(volume_constraint_t), intent(inout) :: this
234 class(design_t), intent(in) :: design
235
236 ! Sensitivity is just a constant so it should not be updated
237
238 end subroutine volume_constraint_update_sensitivity
239
240
241 ! ========================================================================== !
242 ! The actual volume computations for different types of designs
243
244
250 function compute_volume(this, design) result(volume)
251 class(volume_constraint_t), intent(inout) :: this
252 class(design_t), intent(in) :: design
253 real(kind=rp) :: volume
254
255 volume = 0.0_rp
256 select type (design)
257 type is (brinkman_design_t)
258 volume = volume_brinkman_design(this, design)
259
260 class default
261 call neko_error('Volume constraint only works with brinkman_design')
262 end select
263
264 end function compute_volume
265
269 function volume_brinkman_design(this, design) result(volume)
270 class(volume_constraint_t), intent(inout) :: this
271 type(brinkman_design_t), intent(in) :: design
272 real(kind=rp) :: volume
273 type(field_t), pointer :: work, design_indicator
274 integer :: temp_indices(1)
275
276 volume = 0.0_rp
277 design_indicator => neko_field_registry%get_field("design_indicator")
278
279 ! in the future we should be using the mapped design variable
280 ! corresponding to this constraint!!!
281 if (this%has_mask) then
282
283 if (neko_bcknd_device .eq. 1) then
284 call neko_scratch_registry%request_field(work, temp_indices(1))
285 call field_copy(work, design_indicator)
286 call mask_exterior_const(work, this%mask, 0.0_rp)
287 volume = device_glsc2(work%x_d, this%c_xh%B_d, design%size())
288 call neko_scratch_registry%relinquish_field(temp_indices)
289 else
290 volume = glsc2_mask(design_indicator%x, &
291 this%c_Xh%B, design%size(), this%mask%mask, this%mask%size)
292 end if
293
294 else
295
296 if (neko_bcknd_device .eq. 1) then
297 volume = device_glsc2(design_indicator%x_d, &
298 this%c_xh%B_d, design%size())
299 else
300 volume = glsc2(design_indicator%x, &
301 this%c_Xh%B, design%size())
302 end if
303
304 end if
305
306 end function volume_brinkman_design
307
308end module volume_constraint
Implements the constraint_t type.
Implements the design_t.
Definition design.f90:34
Some common Masking operations we may need.
Definition mask_ops.f90:34
Implements the steady_problem_t type.
Implements the volume_constraint_t type.
A topology optimization design variable.
The abstract constraint type.
An abstract design type.
Definition design.f90:50
A constraint on the volume of the design.