Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
develop.f90
Go to the documentation of this file.
1
34
35! ============================================================================ !
36! Supporting functions
37! ============================================================================ !
38module develop
39 use case, only: case_t
40 use global_interpolation, only: global_interpolation_t
41 use json_utils, only: json_get_or_default
42 use logger, only: neko_log, log_size
43 use num_types, only: rp
44 use tuple, only: tuple4_i4_t
45 use utils, only: neko_error
46 use hex, only: hex_t
47 implicit none
48
49 integer, public, parameter :: facet_type_interior = 0
50 integer, public, parameter :: facet_type_outlet = 1
51 integer, public, parameter :: facet_type_inlet = 2
52 integer, public, parameter :: facet_type_wall = 3
53 integer, public, parameter :: facet_type_sympln = 4
54
55contains
56
61 pure function cross(a, b) result(c)
62 real(kind=rp), dimension(3), intent(in) :: a
63 real(kind=rp), dimension(3), intent(in) :: b
64 real(kind=rp), dimension(3) :: c
65
66 c(1) = a(2) * b(3) - a(3) * b(2)
67 c(2) = a(3) * b(1) - a(1) * b(3)
68 c(3) = a(1) * b(2) - a(2) * b(1)
69
70 end function cross
71
76 pure function dot(a, b) result(d)
77 real(kind=rp), dimension(3), intent(in) :: a(3)
78 real(kind=rp), dimension(3), intent(in) :: b(3)
79 real(kind=rp) :: d
80
81 d = a(1) * b(1) + a(2) * b(2) + a(3) * b(3)
82
83 end function dot
84
91 pure function area(nv, vertices) result(a)
92 integer, intent(in) :: nv
93 real(kind=rp), dimension(3, nv), intent(in) :: vertices
94 real(kind=rp) :: a
95
96 ! Internal variables
97 real(kind=rp), dimension(3) :: v1, e1, e2, tmp
98 real(kind=rp) :: cp(3)
99 integer :: v
100
101 ! Initialize the origin to the first vertex, which reduce the
102 ! computational cost by allowing us to skip the first and last cross
103 ! products.
104 v1 = vertices(:, 1)
105 cp = 0.0
106 do v = 2, nv - 1
107 e1 = (vertices(:, v) - v1)
108 e2 = (vertices(:, v + 1) - v1)
109 tmp = cross(e1, e2)
110 cp = cp + 0.5 * tmp
111 end do
112
113 a = sqrt(dot(cp, cp))
114
115 end function area
116
123 subroutine get_facets(C, out, facet_type)
124 type(case_t), intent(in) :: C
125 integer, dimension(:, :), allocatable, intent(out) :: out
126 integer, intent(in), optional :: facet_type
127
128 ! Local variables
129 integer :: e, f, facet_id
130 integer :: N_elem, N_facets
131 integer :: facet_type_
132
133 if (present(facet_type)) then
134 facet_type_ = facet_type
135 else
136 facet_type_ = 1
137 end if
138
139 ! Total number of elements
140 n_elem = c%msh%nelv
141
142 ! Count number of output facets (Facet_type is 1)
143 n_facets = 0
144 do e = 1, n_elem
145 do f = 1, 6
146 if (c%msh%facet_type(f, e) .eq. facet_type_) then
147 n_facets = n_facets + 1
148 end if
149 end do
150 end do
151
152 ! Allocate the outlet_elements array
153 allocate (out(2, n_facets))
154
155 ! Loop across the elements and store the outlet elements
156 facet_id = 1
157 do e = 1, n_elem
158 do f = 1, 6
159 if (c%msh%facet_type(f, e) .eq. facet_type_) then
160 out(1, facet_id) = e
161 out(2, facet_id) = f
162 facet_id = facet_id + 1
163 end if
164 end do
165 end do
166
167 end subroutine get_facets
168
173 subroutine get_facets_outlet(C, out)
174 type(case_t), intent(in) :: C
175 integer, dimension(:, :), allocatable, intent(out) :: out
176
177 call get_facets(c, out, facet_type_outlet)
178 end subroutine get_facets_outlet
179
187 pure function sum_weighted(a, w) result(s)
188 real(kind=rp), dimension(:), intent(in) :: a
189 real(kind=rp), dimension(:), intent(in), optional :: w
190 real(kind=rp) :: s
191
192 if (present(w)) then
193 s = sum(a * w)
194 else
195 s = sum(a)
196 end if
197
198 end function sum_weighted
199
206 pure function average_weighted(a, w) result(m)
207 real(kind=rp), dimension(:), intent(in) :: a
208 real(kind=rp), dimension(:), intent(in), optional :: w
209 real(kind=rp) :: m
210
211 if (present(w)) then
212 m = sum_weighted(a, w) / sum(w)
213 else
214 m = sum(a) / size(a)
215 end if
216
217 end function average_weighted
218
224 subroutine get_facet_nodes(C, element_id, facet_id, nodes)
225 type(case_t), intent(in) :: c
226 integer, intent(in) :: element_id
227 integer, intent(in) :: facet_id
228 real(kind=rp), dimension(:, :), allocatable, intent(out) :: nodes
229
230 ! Local variables
231 integer :: n, v
232 integer :: n_nodes
233 type(tuple4_i4_t) :: t_hex
234
235 ! Determine the number of nodes in the facet
236 n_nodes = 0
237 select type (ele => c%msh%elements(element_id)%e)
238 type is (hex_t)
239 n_nodes = 4
240 end select
241
242 if (n_nodes .eq. 0) then
243 call neko_error("Error: the facet shape is not supported.")
244 end if
245
246 ! Allocate the nodes array
247 if (allocated(nodes) .and. size(nodes, 2) .eq. n_nodes) then
248 nodes = 0.0_rp
249 else if (.not. allocated(nodes)) then
250 allocate (nodes(3, n_nodes))
251 else
252 call neko_error("Error: the nodes array is not the correct size.")
253 end if
254
255 ! Get the nodes
256 select type (ele => c%msh%elements(element_id)%e)
257 type is (hex_t)
258 call ele%facet_order(t_hex, facet_id)
259 do n = 1, n_nodes
260 v = t_hex%x(n)
261 nodes(:, n) = real(c%msh%points(v)%x, kind=rp)
262 end do
263 end select
264
265 end subroutine get_facet_nodes
266
267 subroutine estimate_temperature(neko_case)
268 type(case_t), intent(inout) :: neko_case
269 type(global_interpolation_t) :: interpolator
270
271 real(kind=rp) :: target_temperature
272 real(kind=rp) :: temperature_mean
273 real(kind=rp), dimension(:), allocatable :: temperature_local
274
275 integer, allocatable :: facet_list(:, :)
276
277 real(kind=rp), dimension(:, :), allocatable :: facet_nodes
278 real(kind=rp), dimension(:, :), allocatable :: facet_centers
279 real(kind=rp), dimension(:), allocatable :: facet_area
280
281 integer :: element_id, facet_id
282 integer :: N_facets
283 integer :: f, n
284
285 character(len=LOG_SIZE) :: log_buf
286
287 ! Read the case file for options
288 call json_get_or_default(neko_case%params, 'topopt.target_temperature', &
289 target_temperature, 0.5_rp)
290
291 ! Initialize the global interpolation
292 call interpolator%init(neko_case%scalars%scalar_fields(1)%dm_xh)
293
294 ! Get the list of outlet facets
295 call get_facets(neko_case, facet_list)
296
297 ! Allocate the facet_list array
298 n_facets = size(facet_list, 2)
299 allocate (facet_centers(3, n_facets))
300 allocate (facet_area(n_facets))
301 allocate (temperature_local(n_facets))
302
303 ! Loop across the outlet elements and print the scalar field
304 do f = 1, n_facets
305 element_id = facet_list(1, f)
306 facet_id = facet_list(2, f)
307
308 ! Get the nodes for the facet
309 call get_facet_nodes(neko_case, element_id, facet_id, facet_nodes)
310
311 ! Compute the center of the facet
312 facet_centers(:, f) = 0.0
313 do n = 1, size(facet_nodes, 2)
314 facet_centers(:, f) = facet_centers(:, f) + facet_nodes(:, n)
315 end do
316 facet_centers(:, f) = facet_centers(:, f) / size(facet_nodes, 2)
317
318 ! Compute the area of the facet
319 facet_area(f) = area(size(facet_nodes, 2), facet_nodes)
320 end do
321
322 ! Find the outlet temperature at the supplied list of points
323 call interpolator%find_points_xyz(facet_centers, n_facets)
324 call interpolator%evaluate(temperature_local, &
325 neko_case%scalars%scalar_fields(1)%s%x, on_host=.false.)
326
327 temperature_local = temperature_local - target_temperature
328 temperature_mean = average_weighted(temperature_local, facet_area)
329
330 write (log_buf, '(a,f15.7)') &
331 "Outlet area-weighted average temperature deviation: ", &
332 temperature_mean
333 call neko_log%message(log_buf)
334
335 call interpolator%free()
336
337 end subroutine estimate_temperature
338
339end module develop