Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
design_simplefield.f90
1! Copyright (c) 2024, 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
33! Implements the `simplefield_design_t` type.
34module simplefield_design
35 use num_types, only: rp, sp
36 use field, only: field_t
37 use json_module, only: json_file
38 use mapping, only: mapping_t
39 use pde_filter, only: pde_filter_t
41 use coefs, only: coef_t
42 use scratch_registry, only: neko_scratch_registry
43 use point_zone_registry, only: neko_point_zone_registry
44 use point_zone, only: point_zone_t
46 use neko_config, only: neko_bcknd_device
47 use device, only: device_memcpy, host_to_device
48 use device_math, only: device_copy
49 use design, only: design_t
50 use math, only: rzero
51 use simulation_m, only: simulation_t
52 use json_module, only: json_file
53 use json_utils, only: json_get
55 use vector, only: vector_t
56 use math, only: copy
57 use field_registry, only: neko_field_registry
58
59 use fld_file_output, only: fld_file_output_t
60
61 implicit none
62 private
63
65 type, extends(design_t), public :: simplefield_design_t
66 private
67
68 type(field_t) :: designfield
69 type(vector_t) :: x_coord
70 type(vector_t) :: y_coord
71 type(vector_t) :: z_coord
72
73 ! needed to write the design
74 type(fld_file_output_t) :: output
75
76
77 contains
78
79 ! ----------------------------------------------------------------------- !
80 ! Initializations
81
83 generic, public :: init => init_from_components
85 procedure, pass(this) :: init_from_components => &
86 design_simple_init_from_components
87
89 procedure, pass(this) :: add_mapping => design_simple_add_mapping
90
92 procedure, pass(this) :: get_values => design_simple_get_values
93
94 ! Overrides of base class deferred procedures
95 procedure, pass(this) :: design_get_x => design_simple_get_x
96 procedure, pass(this) :: design_get_y => design_simple_get_y
97 procedure, pass(this) :: design_get_z => design_simple_get_z
98
100 procedure, pass(this) :: update_design => design_simple_update_design
101
103 procedure, pass(this) :: map_forward => design_simple_map_forward
104
106 procedure, pass(this) :: map_backward => design_simple_map_backward
107
109 procedure, pass(this) :: write => design_simple_write
110
112 procedure, pass(this) :: free => design_simple_free
113
114 end type simplefield_design_t
115
116contains
117
118
119 subroutine design_simple_init_from_components(this, n, x, y, z, neko_field)
120 class(simplefield_design_t), intent(inout) :: this
121 integer, intent(in) :: n
122 type(vector_t), intent(in) :: x, y, z
123 type(field_t) :: neko_field
124
125 call this%init_base('simplefield_design', n)
126
127 this%x_coord = x
128 this%y_coord = y
129 this%z_coord = z
130 this%designfield = neko_field
131 call this%output%init(sp, 'design', 1)
132 call this%output%fields%assign_to_field(1, this%designfield)
133 end subroutine design_simple_init_from_components
134
136 subroutine design_simple_free(this)
137 class(simplefield_design_t), intent(inout) :: this
138
139 call this%free_base()
140 call this%x_coord%free()
141 call this%y_coord%free()
142 call this%z_coord%free()
143
144 call this%designfield%free()
145 end subroutine design_simple_free
146
148 subroutine design_simple_add_mapping(this, parameters, simulation)
149 class(simplefield_design_t), intent(inout) :: this
150 type(json_file), intent(inout) :: parameters
151 type(simulation_t), intent(inout) :: simulation
152
153 end subroutine design_simple_add_mapping
154
155
156 subroutine design_simple_map_forward(this)
157 class(simplefield_design_t), intent(inout) :: this
158
159
160 end subroutine design_simple_map_forward
161
162 subroutine design_simple_get_values(this, values)
163 class(simplefield_design_t), intent(in) :: this
164 type(vector_t), intent(inout) :: values
165 integer :: n
166
167 n = this%size()
168 if (neko_bcknd_device .eq. 1) then
169 call device_copy(values%x_d, this%designfield%x_d, n)
170 else
171 call copy(values%x, this%designfield%x, n)
172 end if
173 end subroutine design_simple_get_values
174
175 subroutine design_simple_get_x(this, x)
176 class(simplefield_design_t), intent(in) :: this
177 type(vector_t), intent(inout) :: x
178
179 x = this%x_coord
180 end subroutine design_simple_get_x
181
182 subroutine design_simple_get_y(this, y)
183 class(simplefield_design_t), intent(in) :: this
184 type(vector_t), intent(inout) :: y
185
186 y = this%y_coord
187 end subroutine design_simple_get_y
188
189 subroutine design_simple_get_z(this, z)
190 class(simplefield_design_t), intent(in) :: this
191 type(vector_t), intent(inout) :: z
192
193 z = this%z_coord
194 end subroutine design_simple_get_z
195
196
197
198 subroutine design_simple_update_design(this, values)
199 class(simplefield_design_t), intent(inout) :: this
200 type(vector_t), intent(inout) :: values
201 integer :: n
202
203 n = this%size()
204 if (neko_bcknd_device .eq. 1) then
205 call device_copy(this%designfield%x_d, values%x_d, n)
206 else
207 call copy(this%designfield%x, values%x, n)
208 end if
209 end subroutine design_simple_update_design
210
211 subroutine design_simple_map_backward(this, sensitivity)
212 class(simplefield_design_t), intent(inout) :: this
213 type(vector_t), intent(in) :: sensitivity
214
215 end subroutine design_simple_map_backward
216
217 subroutine design_simple_write(this, idx)
218 class(simplefield_design_t), intent(inout) :: this
219 integer, intent(in) :: idx
220 call this%output%sample(real(idx, kind=rp))
221 end subroutine design_simple_write
222
223end module simplefield_design
Implements the design_t.
Definition design.f90:34
Mappings to be applied to a scalar field.
Definition mapping.f90:35
Some common Masking operations we may need.
Definition mask_ops.f90:34
A PDE based filter.
A RAMP mapping of coefficients.
Implements the simple_brinkman_source_term_t type.
Implements the steady_problem_t type.
An abstract design type.
Definition design.f90:52
Base abstract class for mapping.
Definition mapping.f90:45
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .
A RAMP mapping of coefficients This is the standard RAMP described in https://doi....
A topology optimization design variable.