Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
volume_constraint.f90
Go to the documentation of this file.
1
34!
41 use constraint, only: constraint_t
42
43 use design, only: design_t
44 use brinkman_design, only: brinkman_design_t
45 use simulation_m, only: simulation_t
47 use num_types, only: rp
48 use coefs, only: coef_t
49 use json_module, only: json_file
50 use json_utils, only: json_get, json_get_or_default
51 use field, only: field_t
52 use scratch_registry, only: neko_scratch_registry
53 use neko_config, only: neko_bcknd_device
55 use math, only: glsc2, copy, cmult
56 use device_math, only: device_glsc2, device_copy, device_cmult
57 use vector_math, only: vector_cmult
58 use math_ext, only: glsc2_mask
59 use field_math, only: field_rone, field_copy, field_cmult, field_cfill
60 use utils, only: neko_error
61 use vector, only: vector_t
62 use neko_ext, only: field_to_vector
63 implicit none
64 private
65
67 type, public, extends(constraint_t) :: volume_constraint_t
68 private
69
73 logical :: is_max
75 real(kind=rp) :: limit
77 real(kind=rp) :: volume_domain
79 class(coef_t), pointer :: c_xh => null()
83 logical :: if_mapping = .false.
84
85 contains
86
88 procedure, public, pass(this) :: init_json_sim => &
91 procedure, public, pass(this) :: init_from_components => &
92 volume_constraint_init_from_components
94 procedure, public, pass(this) :: free => volume_constraint_free
96 procedure, public, pass(this) :: update_value => &
97 volume_constraint_update_value
99 procedure, public, pass(this) :: update_sensitivity => &
100 volume_constraint_update_sensitivity
101
103 procedure, private, pass(this) :: compute_volume
104
105 end type volume_constraint_t
106
107contains
108
114 subroutine volume_constraint_init_json_sim(this, json, design, simulation)
115 class(volume_constraint_t), intent(inout) :: this
116 type(json_file), intent(inout) :: json
117 class(design_t), intent(in) :: design
118 type(simulation_t), target, intent(inout) :: simulation
119
120 character(len=:), allocatable :: mask_name
121 character(len=:), allocatable :: name
122 logical :: is_max
123 real(kind=rp) :: limit
124
125 call json_get_or_default(json, "name", name, "Volume constraint")
126 call json_get_or_default(json, "mask_name", mask_name, "")
127 call json_get_or_default(json, "is_max", is_max, .false.)
128 call json_get(json, "limit", limit)
129
130 call this%init_from_components(design, simulation, name, mask_name, &
131 is_max, limit)
132
133 ! check if we have a mapping
134 if ('mapping' .in. json) then
135 ! Initialize the mapper
136 call this%mapping%init_base(this%c_Xh)
137 call this%mapping%add(json, 'mapping')
138 this%if_mapping = .true.
139 ! recompute value after mapping
140 call this%update_value(design)
141 end if
142
144
153 subroutine volume_constraint_init_from_components(this, design, simulation, &
154 name, mask_name, is_max, limit)
155 class(volume_constraint_t), intent(inout) :: this
156 class(design_t), intent(in) :: design
157 type(simulation_t), target, intent(inout) :: simulation
158 character(len=*), intent(in) :: mask_name
159 character(len=*), intent(in) :: name
160 logical, intent(in) :: is_max
161 real(kind=rp), intent(in) :: limit
162
163 type(field_t), pointer :: work
164 integer :: temp_indices(1)
165
166 ! Initialize the base class
167 call this%init_base(name, design%size(), mask_name)
168
169 ! Store the attributes
170 this%is_max = is_max
171 this%limit = limit
172 this%c_Xh => simulation%neko_case%fluid%c_Xh
173
174 ! Now we can extract the mask/has_mask from the design
175 if (this%has_mask) then
176
177 ! calculate the volume of the optimization domain
178 call neko_scratch_registry%request(work, temp_indices(1), .false.)
179 call field_rone(work)
180 call mask_exterior_const(work, this%mask, 0.0_rp)
181
182 if (neko_bcknd_device .eq. 1) then
183 this%volume_domain = device_glsc2(work%x_d, this%c_xh%B_d, &
184 work%size())
185 else
186 this%volume_domain = glsc2_mask(work%x, this%c_Xh%B, &
187 design%size(), this%mask%mask%get(), this%mask%size)
188 end if
189
190 call neko_scratch_registry%relinquish(temp_indices)
191 else
192 this%volume_domain = this%c_Xh%volume
193 end if
194
195 ! ------------------------------------------------------------------------ !
196 ! Initialize the value of constraint
197
198 call this%update_value(design)
199
200 ! ------------------------------------------------------------------------ !
201 ! Initialize the sensitivity value
202
203 if (neko_bcknd_device .eq. 1) then
204 call device_copy(this%sensitivity%x_d, this%c_xh%B_d, design%size())
205 call device_cmult(this%sensitivity%x_d, -1.0_rp / this%volume_domain, &
206 design%size())
207
208 if (this%is_max) then
209 call device_cmult(this%sensitivity%x_d, -1.0_rp, design%size())
210 end if
211 else
212 call copy(this%sensitivity%x, this%c_Xh%B, design%size())
213 call cmult(this%sensitivity%x, -1.0_rp / this%volume_domain, &
214 design%size())
215
216 if (this%is_max) then
217 call vector_cmult(this%sensitivity, -1.0_rp)
218 end if
219 end if
220
221 if (this%has_mask) then
222 call mask_exterior_const(this%sensitivity, this%mask, 0.0_rp)
223 end if
224
225 end subroutine volume_constraint_init_from_components
226
228 subroutine volume_constraint_free(this)
229 class(volume_constraint_t), intent(inout) :: this
230
231 call this%free_base()
232 end subroutine volume_constraint_free
233
237 subroutine volume_constraint_update_value(this, design)
238 class(volume_constraint_t), intent(inout) :: this
239 class(design_t), intent(in) :: design
240 real(kind=rp) :: volume
241
242 volume = this%compute_volume(design)
243
244 ! Compute the distance to the target volume
245 this%value = this%limit - volume / this%volume_domain
246
247 ! Invert the sign if it is a maximum constraint
248 if (this%is_max) this%value = -this%value
249
250 end subroutine volume_constraint_update_value
251
255 subroutine volume_constraint_update_sensitivity(this, design)
256 class(volume_constraint_t), intent(inout) :: this
257 class(design_t), intent(in) :: design
259 type(field_t), pointer :: unmapped, mapped
260 integer :: temp_indices(2)
261
262 if (this%if_mapping) then
263 ! Recompute and map backward
264 call neko_scratch_registry%request(unmapped, temp_indices(1), &
265 .false.)
266 call neko_scratch_registry%request(mapped, temp_indices(2), &
267 .false.)
268 ! The mapping will handle the mass matrix
269 call field_cfill(unmapped, -1.0_rp / this%volume_domain)
270 if (this%is_max) then
271 call field_cmult(unmapped, -1.0_rp)
272 end if
273 ! mask
274 if (this%has_mask) then
275 call mask_exterior_const(unmapped, this%mask, 0.0_rp)
276 end if
277 ! map backwards
278 call this%mapping%apply_backward(mapped, unmapped)
279 call field_to_vector(this%sensitivity, mapped)
280
281 call neko_scratch_registry%relinquish(temp_indices)
282
283 else
284 ! Sensitivity is just a constant so it should not be updated
285 end if
286
287 end subroutine volume_constraint_update_sensitivity
288
289 ! ========================================================================== !
290 ! The actual volume computations for different types of designs
291
297 function compute_volume(this, design) result(volume)
298 class(volume_constraint_t), intent(inout) :: this
299 class(design_t), intent(in) :: design
300 real(kind=rp) :: volume
301
302 volume = 0.0_rp
303 select type (design)
304 type is (brinkman_design_t)
305 volume = volume_brinkman_design(this, design)
306
307 class default
308 call neko_error('Volume constraint only works with brinkman_design')
309 end select
310
311 end function compute_volume
312
316 function volume_brinkman_design(this, design) result(volume)
317 class(volume_constraint_t), intent(inout) :: this
318 type(brinkman_design_t), intent(in) :: design
319 real(kind=rp) :: volume
320 type(field_t), pointer :: work
321 type(vector_t), pointer :: values, unmapped_values
322 integer :: temp_indices, ind_value, ind_um_value
323
324 call neko_scratch_registry%request(values, ind_value, design%size(), &
325 .false.)
326
327 volume = 0.0_rp
328 if (this%if_mapping) then
329 call neko_scratch_registry%request(unmapped_values, ind_um_value, &
330 design%size(), .false.)
331 call design%get_values(unmapped_values)
332 call this%mapping%apply_forward(values, unmapped_values)
333 call neko_scratch_registry%relinquish(ind_um_value)
334 else
335 call design%get_values(values)
336 end if
337
338 if (this%has_mask) then
339
340 if (neko_bcknd_device .eq. 1) then
341 call neko_scratch_registry%request(work, temp_indices, .false.)
342 call device_copy(work%x_d, values%x_d, design%size())
343 call mask_exterior_const(work, this%mask, 0.0_rp)
344
345 volume = device_glsc2(work%x_d, this%c_xh%B_d, design%size())
346
347 call neko_scratch_registry%relinquish(temp_indices)
348 else
349 volume = glsc2_mask(values%x, this%c_Xh%B, design%size(), &
350 this%mask%mask%get(), this%mask%size)
351 end if
352
353 else
354
355 if (neko_bcknd_device .eq. 1) then
356 volume = device_glsc2(values%x_d, this%c_xh%B_d, design%size())
357 else
358 volume = glsc2(values%x, this%c_Xh%B, design%size())
359 end if
360
361 end if
362
363 call neko_scratch_registry%relinquish(ind_value)
364
365 end function volume_brinkman_design
366
367end module volume_constraint
Implements the constraint_t type.
Implements the design_t.
Definition design.f90:36
Implements the mapping_handler_t type.
Mappings to be applied to a scalar field.
Definition mapping.f90:36
Some common Masking operations we may need.
Definition mask_ops.f90:36
Contains extensions to the neko library required to run the topology optimization code.
Definition neko_ext.f90:42
subroutine, public field_to_vector(vector, field)
Field to vector.
Definition neko_ext.f90:406
Implements the steady_problem_t type.
Implements the volume_constraint_t type. Either .
subroutine volume_constraint_init_json_sim(this, json, design, simulation)
The common constructor using a JSON object.
A topology optimization design variable.
The abstract constraint type.
An abstract design type.
Definition design.f90:54
Abstract class for handling mapping_cascade.
A constraint on the volume of the design.