Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_fluid_ic.f90
Go to the documentation of this file.
1! Copyright (c) 2021, 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!
35 use num_types, only : rp
36 use gather_scatter, only : gs_t, gs_op_add
37 use neko_config, only : neko_bcknd_device
38 use flow_profile, only : blasius_profile, blasius_linear, blasius_cubic, &
39 blasius_quadratic, blasius_quartic, blasius_sin
40 use device, only: device_memcpy, host_to_device
41 use field, only: field_t
42 use utils, only: neko_error
43 use coefs, only: coef_t
44 use math, only: col2, cfill, cfill_mask
45 use device_math, only: device_col2, device_cfill, device_cfill_mask
46 use user_intf, only: useric
47 use json_module, only: json_file
48 use json_utils, only: json_get
49 use point_zone, only: point_zone_t
50 use point_zone_registry, only: neko_point_zone_registry
51 implicit none
52 private
53
56 end interface set_adjoint_fluid_ic
57
58 public :: set_adjoint_fluid_ic
59
60contains
61
63 subroutine set_adjoint_fluid_ic_int(u_adj, v_adj, w_adj, p_adj, coef, gs, &
64 type, params)
65 type(field_t), intent(inout) :: u_adj
66 type(field_t), intent(inout) :: v_adj
67 type(field_t), intent(inout) :: w_adj
68 type(field_t), intent(inout) :: p_adj
69 type(coef_t), intent(in) :: coef
70 type(gs_t), intent(inout) :: gs
71 character(len=*) :: type
72 type(json_file), intent(inout) :: params
73 real(kind=rp) :: delta
74 real(kind=rp), allocatable :: uinf(:)
75 real(kind=rp), allocatable :: zone_value(:)
76 character(len=:), allocatable :: blasius_approximation
77 character(len=:), allocatable :: zone_name
78
79 character(len=:), allocatable :: type_
80 call json_get(params, 'type', type_)
81
82 if (trim(type_) .eq. 'uniform') then
83 call json_get(params, 'value', uinf)
84 call set_adjoint_fluid_ic_uniform(u_adj, v_adj, w_adj, uinf)
85 else if (trim(type_) .eq. 'blasius') then
86 call json_get(params, 'blasius.delta', delta)
87 call json_get(params, 'blasius.approximation', &
88 blasius_approximation)
89 call json_get(params, 'blasius.freestream_velocity', uinf)
90 call set_adjoint_fluid_ic_blasius(u_adj, v_adj, w_adj, delta, uinf, &
91 blasius_approximation)
92 else if (trim(type_) .eq. 'point_zone') then
93 call json_get(params, 'base_value', uinf)
94 call json_get(params, 'zone_name', &
95 zone_name)
96 call json_get(params, 'zone_value', &
97 zone_value)
98 call set_adjoint_fluid_ic_point_zone(u_adj, v_adj, w_adj, uinf, &
99 zone_name, zone_value)
100 else
101 call neko_error('Invalid initial condition')
102 end if
103
104 call set_adjoint_fluid_ic_common(u_adj, v_adj, w_adj, p_adj, coef, gs)
105
106 end subroutine set_adjoint_fluid_ic_int
107
109 subroutine set_adjoint_fluid_ic_usr(u_adj, v_adj, w_adj, p_adj, coef, gs, &
110 usr_ic, params)
111 type(field_t), intent(inout) :: u_adj
112 type(field_t), intent(inout) :: v_adj
113 type(field_t), intent(inout) :: w_adj
114 type(field_t), intent(inout) :: p_adj
115 type(coef_t), intent(in) :: coef
116 type(gs_t), intent(inout) :: gs
117 procedure(useric) :: usr_ic
118 type(json_file), intent(inout) :: params
119
120 call usr_ic(u_adj, v_adj, w_adj, p_adj, params)
121
122 call set_adjoint_fluid_ic_common(u_adj, v_adj, w_adj, p_adj, coef, gs)
123
124 end subroutine set_adjoint_fluid_ic_usr
125
126 subroutine set_adjoint_fluid_ic_common(u_adj, v_adj, w_adj, p_adj, coef, gs)
127 type(field_t), intent(inout) :: u_adj
128 type(field_t), intent(inout) :: v_adj
129 type(field_t), intent(inout) :: w_adj
130 type(field_t), intent(inout) :: p_adj
131 type(coef_t), intent(in) :: coef
132 type(gs_t), intent(inout) :: gs
133 integer :: n
134
135 n = u_adj%dof%size()
136
137 if (neko_bcknd_device .eq. 1) then
138 call device_memcpy(u_adj%x, u_adj%x_d, n, &
139 host_to_device, sync = .false.)
140 call device_memcpy(v_adj%x, v_adj%x_d, n, &
141 host_to_device, sync = .false.)
142 call device_memcpy(w_adj%x, w_adj%x_d, n, &
143 host_to_device, sync = .false.)
144 end if
145
146 ! Ensure continuity across elements for initial conditions
147 call gs%op(u_adj%x, u_adj%dof%size(), gs_op_add)
148 call gs%op(v_adj%x, v_adj%dof%size(), gs_op_add)
149 call gs%op(w_adj%x, w_adj%dof%size(), gs_op_add)
150
151 if (neko_bcknd_device .eq. 1) then
152 call device_col2(u_adj%x_d, coef%mult_d, u_adj%dof%size())
153 call device_col2(v_adj%x_d, coef%mult_d, v_adj%dof%size())
154 call device_col2(w_adj%x_d, coef%mult_d, w_adj%dof%size())
155 else
156 call col2(u_adj%x, coef%mult, u_adj%dof%size())
157 call col2(v_adj%x, coef%mult, v_adj%dof%size())
158 call col2(w_adj%x, coef%mult, w_adj%dof%size())
159 end if
160
161 end subroutine set_adjoint_fluid_ic_common
162
164 subroutine set_adjoint_fluid_ic_uniform(u_adj, v_adj, w_adj, uinf)
165 type(field_t), intent(inout) :: u_adj
166 type(field_t), intent(inout) :: v_adj
167 type(field_t), intent(inout) :: w_adj
168 real(kind=rp), intent(in) :: uinf(3)
169 integer :: n
170 u_adj = uinf(1)
171 v_adj = uinf(2)
172 w_adj = uinf(3)
173 n = u_adj%dof%size()
174
175 if (neko_bcknd_device .eq. 1) then
176 call device_cfill(u_adj%x_d, uinf(1), n)
177 call device_cfill(v_adj%x_d, uinf(2), n)
178 call device_cfill(w_adj%x_d, uinf(3), n)
179 else
180 call cfill(u_adj%x, uinf(1), n)
181 call cfill(v_adj%x, uinf(2), n)
182 call cfill(w_adj%x, uinf(3), n)
183 end if
184
185 end subroutine set_adjoint_fluid_ic_uniform
186
189 subroutine set_adjoint_fluid_ic_blasius(u_adj, v_adj, w_adj, delta, uinf, &
190 type)
191 type(field_t), intent(inout) :: u_adj
192 type(field_t), intent(inout) :: v_adj
193 type(field_t), intent(inout) :: w_adj
194 real(kind=rp), intent(in) :: delta
195 real(kind=rp), intent(in) :: uinf(3)
196 character(len=*), intent(in) :: type
197 procedure(blasius_profile), pointer :: bla => null()
198 integer :: i
199
200 select case (trim(type))
201 case ('linear')
202 bla => blasius_linear
203 case ('quadratic')
204 bla => blasius_quadratic
205 case ('cubic')
206 bla => blasius_cubic
207 case ('quartic')
208 bla => blasius_quartic
209 case ('sin')
210 bla => blasius_sin
211 case default
212 call neko_error('Invalid Blasius approximation')
213 end select
214
215 if ((uinf(1) .gt. 0.0_rp) .and. (uinf(2) .le. 0.0_rp) &
216 .and. (uinf(3) .le. 0.0_rp)) then
217 do i = 1, u_adj%dof%size()
218 u_adj%x(i,1,1,1) = bla(u_adj%dof%z(i,1,1,1), delta, uinf(1))
219 v_adj%x(i,1,1,1) = 0.0_rp
220 w_adj%x(i,1,1,1) = 0.0_rp
221 end do
222 else if ((uinf(1) .le. 0.0_rp) .and. (uinf(2) .gt. 0.0_rp) &
223 .and. (uinf(3) .le. 0.0_rp)) then
224 do i = 1, u_adj%dof%size()
225 u_adj%x(i,1,1,1) = 0.0_rp
226 v_adj%x(i,1,1,1) = bla(u_adj%dof%x(i,1,1,1), delta, uinf(2))
227 w_adj%x(i,1,1,1) = 0.0_rp
228 end do
229 else if ((uinf(1) .le. 0.0_rp) .and. (uinf(2) .le. 0.0_rp) &
230 .and. (uinf(3) .gt. 0.0_rp)) then
231 do i = 1, u_adj%dof%size()
232 u_adj%x(i,1,1,1) = 0.0_rp
233 v_adj%x(i,1,1,1) = 0.0_rp
234 w_adj%x(i,1,1,1) = bla(u_adj%dof%y(i,1,1,1), delta, uinf(3))
235 end do
236 end if
237
238 end subroutine set_adjoint_fluid_ic_blasius
239
249 subroutine set_adjoint_fluid_ic_point_zone(u_adj, v_adj, w_adj, base_value, &
250 zone_name, zone_value)
251 type(field_t), intent(inout) :: u_adj
252 type(field_t), intent(inout) :: v_adj
253 type(field_t), intent(inout) :: w_adj
254 real(kind=rp), intent(in), dimension(3) :: base_value
255 character(len=*), intent(in) :: zone_name
256 real(kind=rp), intent(in) :: zone_value(:)
257
258 ! Internal variables
259 class(point_zone_t), pointer :: zone
260 integer :: size
261
262 call set_adjoint_fluid_ic_uniform(u_adj, v_adj, w_adj, base_value)
263 size = u_adj%dof%size()
264
265 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
266
267 if (neko_bcknd_device .eq. 1) then
268 call device_cfill_mask(u_adj%x_d, zone_value(1), size, &
269 zone%mask_d, zone%size)
270 call device_cfill_mask(v_adj%x_d, zone_value(2), size, &
271 zone%mask_d, zone%size)
272 call device_cfill_mask(w_adj%x_d, zone_value(3), size, &
273 zone%mask_d, zone%size)
274 else
275 call cfill_mask(u_adj%x, zone_value(1), size, zone%mask, zone%size)
276 call cfill_mask(v_adj%x, zone_value(2), size, zone%mask, zone%size)
277 call cfill_mask(w_adj%x, zone_value(3), size, zone%mask, zone%size)
278
279 end if
281
282end module adjoint_fluid_ic
Initial flow condition.
subroutine set_adjoint_fluid_ic_int(u_adj, v_adj, w_adj, p_adj, coef, gs, type, params)
Set initial flow condition (builtin)
subroutine set_adjoint_fluid_ic_common(u_adj, v_adj, w_adj, p_adj, coef, gs)
subroutine set_adjoint_fluid_ic_uniform(u_adj, v_adj, w_adj, uinf)
Uniform initial condition.
subroutine set_adjoint_fluid_ic_usr(u_adj, v_adj, w_adj, p_adj, coef, gs, usr_ic, params)
Set intial flow condition (user defined)
subroutine set_adjoint_fluid_ic_blasius(u_adj, v_adj, w_adj, delta, uinf, type)
Set a Blasius profile as initial condition.
subroutine set_adjoint_fluid_ic_point_zone(u_adj, v_adj, w_adj, base_value, zone_name, zone_value)
Set the initial condition of the flow based on a point zone.