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