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! ============================================================================ !
2! Supporting functions
3! ============================================================================ !
4module develop
5 use case, only: case_t
6 use global_interpolation, only: global_interpolation_t
7 use json_utils, only: json_get_or_default
8 use logger, only: neko_log, log_size
9 use num_types, only: rp
10 use tuple, only: tuple4_i4_t
11 use utils, only: neko_error
12 use hex, only: hex_t
13 implicit none
14
15 integer, public, parameter :: facet_type_interior = 0
16 integer, public, parameter :: facet_type_outlet = 1
17 integer, public, parameter :: facet_type_inlet = 2
18 integer, public, parameter :: facet_type_wall = 3
19 integer, public, parameter :: facet_type_sympln = 4
20
21contains
22
27 pure function cross(a, b) result(c)
28 real(kind=rp), dimension(3), intent(in) :: a
29 real(kind=rp), dimension(3), intent(in) :: b
30 real(kind=rp), dimension(3) :: c
31
32 c(1) = a(2) * b(3) - a(3) * b(2)
33 c(2) = a(3) * b(1) - a(1) * b(3)
34 c(3) = a(1) * b(2) - a(2) * b(1)
35
36 end function cross
37
42 pure function dot(a, b) result(d)
43 real(kind=rp), dimension(3), intent(in) :: a(3)
44 real(kind=rp), dimension(3), intent(in) :: b(3)
45 real(kind=rp) :: d
46
47 d = a(1) * b(1) + a(2) * b(2) + a(3) * b(3)
48
49 end function dot
50
57 pure function area(nv, vertices) result(a)
58 integer, intent(in) :: nv
59 real(kind=rp), dimension(3, nv), intent(in) :: vertices
60 real(kind=rp) :: a
61
62 ! Internal variables
63 real(kind=rp), dimension(3) :: v1, e1, e2, tmp
64 real(kind=rp) :: cp(3)
65 integer :: v
66
67 ! Initialize the origin to the first vertex, which reduce the
68 ! computational cost by allowing us to skip the first and last cross
69 ! products.
70 v1 = vertices(:, 1)
71 cp = 0.0
72 do v = 2, nv - 1
73 e1 = (vertices(:, v) - v1)
74 e2 = (vertices(:, v + 1) - v1)
75 tmp = cross(e1, e2)
76 cp = cp + 0.5 * tmp
77 end do
78
79 a = sqrt(dot(cp, cp))
80
81 end function area
82
89 subroutine get_facets(C, out, facet_type)
90 type(case_t), intent(in) :: C
91 integer, dimension(:, :), allocatable, intent(out) :: out
92 integer, intent(in), optional :: facet_type
93
94 ! Local variables
95 integer :: e, f, facet_id
96 integer :: N_elem, N_facets
97 integer :: facet_type_
98
99 if (present(facet_type)) then
100 facet_type_ = facet_type
101 else
102 facet_type_ = 1
103 end if
104
105 ! Total number of elements
106 n_elem = c%msh%nelv
107
108 ! Count number of output facets (Facet_type is 1)
109 n_facets = 0
110 do e = 1, n_elem
111 do f = 1, 6
112 if (c%msh%facet_type(f, e) .eq. facet_type_) then
113 n_facets = n_facets + 1
114 end if
115 end do
116 end do
117
118 ! Allocate the outlet_elements array
119 allocate (out(2, n_facets))
120
121 ! Loop across the elements and store the outlet elements
122 facet_id = 1
123 do e = 1, n_elem
124 do f = 1, 6
125 if (c%msh%facet_type(f, e) .eq. facet_type_) then
126 out(1, facet_id) = e
127 out(2, facet_id) = f
128 facet_id = facet_id + 1
129 end if
130 end do
131 end do
132
133 end subroutine get_facets
134
139 subroutine get_facets_outlet(C, out)
140 type(case_t), intent(in) :: C
141 integer, dimension(:, :), allocatable, intent(out) :: out
142
143 call get_facets(c, out, facet_type_outlet)
144 end subroutine get_facets_outlet
145
153 pure function sum_weighted(a, w) result(s)
154 real(kind=rp), dimension(:), intent(in) :: a
155 real(kind=rp), dimension(:), intent(in), optional :: w
156 real(kind=rp) :: s
157
158 if (present(w)) then
159 s = sum(a * w)
160 else
161 s = sum(a)
162 end if
163
164 end function sum_weighted
165
172 pure function average_weighted(a, w) result(m)
173 real(kind=rp), dimension(:), intent(in) :: a
174 real(kind=rp), dimension(:), intent(in), optional :: w
175 real(kind=rp) :: m
176
177 if (present(w)) then
178 m = sum_weighted(a, w) / sum(w)
179 else
180 m = sum(a) / size(a)
181 end if
182
183 end function average_weighted
184
190 subroutine get_facet_nodes(C, element_id, facet_id, nodes)
191 type(case_t), intent(in) :: c
192 integer, intent(in) :: element_id
193 integer, intent(in) :: facet_id
194 real(kind=rp), dimension(:, :), allocatable, intent(out) :: nodes
195
196 ! Local variables
197 integer :: n, v
198 integer :: n_nodes
199 type(tuple4_i4_t) :: t_hex
200
201 ! Determine the number of nodes in the facet
202 n_nodes = 0
203 select type (ele => c%msh%elements(element_id)%e)
204 type is (hex_t)
205 n_nodes = 4
206 end select
207
208 if (n_nodes .eq. 0) then
209 call neko_error("Error: the facet shape is not supported.")
210 end if
211
212 ! Allocate the nodes array
213 if (allocated(nodes) .and. size(nodes, 2) .eq. n_nodes) then
214 nodes = 0.0_rp
215 else if (.not. allocated(nodes)) then
216 allocate (nodes(3, n_nodes))
217 else
218 call neko_error("Error: the nodes array is not the correct size.")
219 end if
220
221 ! Get the nodes
222 select type (ele => c%msh%elements(element_id)%e)
223 type is (hex_t)
224 call ele%facet_order(t_hex, facet_id)
225 do n = 1, n_nodes
226 v = t_hex%x(n)
227 nodes(:, n) = c%msh%points(v)%x
228 end do
229 end select
230
231 end subroutine get_facet_nodes
232
233 subroutine estimate_temperature(neko_case)
234 type(case_t), intent(inout) :: neko_case
235 type(global_interpolation_t) :: interpolator
236
237 real(kind=rp) :: target_temperature
238 real(kind=rp) :: temperature_mean
239 real(kind=rp), dimension(:), allocatable :: temperature_local
240
241 integer, allocatable :: facet_list(:, :)
242
243 real(kind=rp), dimension(:, :), allocatable :: facet_nodes
244 real(kind=rp), dimension(:, :), allocatable :: facet_centers
245 real(kind=rp), dimension(:), allocatable :: facet_area
246
247 integer :: element_id, facet_id
248 integer :: N_facets
249 integer :: f, n
250
251 character(len=LOG_SIZE) :: log_buf
252
253 ! Read the case file for options
254 call json_get_or_default(neko_case%params, 'topopt.target_temperature', &
255 target_temperature, 0.5_rp)
256
257 ! Initialize the global interpolation
258 call interpolator%init(neko_case%scalar%dm_xh)
259
260 ! Get the list of outlet facets
261 call get_facets(neko_case, facet_list)
262
263 ! Allocate the facet_list array
264 n_facets = size(facet_list, 2)
265 allocate (facet_centers(3, n_facets))
266 allocate (facet_area(n_facets))
267 allocate (temperature_local(n_facets))
268
269 ! Loop across the outlet elements and print the scalar field
270 do f = 1, n_facets
271 element_id = facet_list(1, f)
272 facet_id = facet_list(2, f)
273
274 ! Get the nodes for the facet
275 call get_facet_nodes(neko_case, element_id, facet_id, facet_nodes)
276
277 ! Compute the center of the facet
278 facet_centers(:, f) = 0.0
279 do n = 1, size(facet_nodes, 2)
280 facet_centers(:, f) = facet_centers(:, f) + facet_nodes(:, n)
281 end do
282 facet_centers(:, f) = facet_centers(:, f) / size(facet_nodes, 2)
283
284 ! Compute the area of the facet
285 facet_area(f) = area(size(facet_nodes, 2), facet_nodes)
286 end do
287
288 ! Find the outlet temperature at the supplied list of points
289 call interpolator%find_points_xyz(facet_centers, n_facets)
290 call interpolator%evaluate(temperature_local, neko_case%scalar%s%x)
291
292 temperature_local = temperature_local - target_temperature
293 temperature_mean = average_weighted(temperature_local, facet_area)
294
295 write (log_buf, '(a,f15.7)') &
296 "Outlet area-weighted average temperature deviation: ", &
297 temperature_mean
298 call neko_log%message(log_buf)
299
300 call interpolator%free()
301
302 end subroutine estimate_temperature
303
304end module develop
integer, parameter, public facet_type_outlet
Definition develop.f90:16
subroutine get_facets_outlet(c, out)
Construct a list of facets.
Definition develop.f90:140
subroutine get_facet_nodes(c, element_id, facet_id, nodes)
Construct a list of nodes for a facet.
Definition develop.f90:191
pure real(kind=rp) function dot(a, b)
Compute dot product of two vectors.
Definition develop.f90:43
integer, parameter, public facet_type_sympln
Definition develop.f90:19
integer, parameter, public facet_type_inlet
Definition develop.f90:17
pure real(kind=rp) function, dimension(3) cross(a, b)
Compute cross product of two vectors.
Definition develop.f90:28
subroutine estimate_temperature(neko_case)
Definition develop.f90:234
pure real(kind=rp) function area(nv, vertices)
Compute area of a polygon.
Definition develop.f90:58
subroutine get_facets(c, out, facet_type)
Construct a list of facets.
Definition develop.f90:90
integer, parameter, public facet_type_interior
Definition develop.f90:15
integer, parameter, public facet_type_wall
Definition develop.f90:18
pure real(kind=rp) function sum_weighted(a, w)
Compute the weighted sum of an array.
Definition develop.f90:154
pure real(kind=rp) function average_weighted(a, w)
Compute the weighted average of an array.
Definition develop.f90:173