Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_viscous_dissipation_source_term.f90
Go to the documentation of this file.
1
34!
36!
37!
38! If the objective function $\frac{\mu}{2} \int |\nabla u|^2$, the
39! corresponding adjoint forcing is
40! $ \mu \int \nabla v \cdot \nabla u $ in weak form.
42 use num_types, only: rp
43 use field_list, only: field_list_t
44 use json_module, only: json_file
45 use json_utils, only: json_get, json_get_or_default
46 use source_term, only: source_term_t
47 use coefs, only: coef_t
48 use neko_config, only: neko_bcknd_device
49 use utils, only: neko_error
50 use field, only: field_t
51 use field_math, only: field_subcol3, field_add2, field_add2s2, field_rzero
52 use json_module, only: json_file
53 use time_state, only: time_state_t
55 use simcomp_executor, only: neko_simcomps
56 use user_source_term, only: user_source_term_t
57 use num_types, only: rp
58 use field, only: field_t
59 use registry, only: neko_registry
60 use math, only: rzero, copy, chsign, cfill, invcol2
61 use device_math, only: device_copy, device_cmult, device_cfill, device_invcol2
62 use neko_config, only: neko_bcknd_device
63 use scratch_registry, only: neko_scratch_registry
65 use point_zone, only: point_zone_t
66 use ax_product, only : ax_t, ax_helm_factory
67
68 implicit none
69 private
71
73 ! $\int \nabla v \cdot \nabla u $
74 type, public, extends(source_term_t) :: &
77 type(field_t), pointer :: u => null()
79 type(field_t), pointer :: v => null()
81 type(field_t), pointer :: w => null()
83 real(kind=rp) :: obj_scale
85 class(point_zone_t), pointer :: mask => null()
87 logical :: if_mask
89 class(ax_t), allocatable :: ax
91 real(kind=rp) :: volume
92
93 contains
95 procedure, pass(this) :: init => &
96 adjoint_viscous_dissipation_source_term_init_from_json
98 procedure, pass(this) :: init_from_components => &
99 adjoint_viscous_dissipation_source_term_init_from_components
101 procedure, pass(this) :: free => &
102 adjoint_viscous_dissipation_source_term_free
104 procedure, pass(this) :: compute_ => &
105 adjoint_viscous_dissipation_source_term_compute
107
108contains
109
112 class(source_term_t), allocatable, intent(inout) :: obj
115
122 subroutine adjoint_viscous_dissipation_source_term_init_from_json(this, &
123 json, fields, coef, variable_name)
124 class(adjoint_viscous_dissipation_source_term_t), intent(inout) :: this
125 type(json_file), intent(inout) :: json
126 type(field_list_t), intent(in), target :: fields
127 type(coef_t), intent(in), target :: coef
128 character(len=*), intent(in) :: variable_name
129
130 ! we shouldn't be initializing this from JSON
131 ! maybe throw an error?
132
133
134 end subroutine adjoint_viscous_dissipation_source_term_init_from_json
135
147 subroutine adjoint_viscous_dissipation_source_term_init_from_components( &
148 this, f_x, f_y, f_z, u, v, w, obj_scale, mask, if_mask, coef, &
149 volume, start_time, end_time)
150 class(adjoint_viscous_dissipation_source_term_t), intent(inout) :: this
151 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
152 type(field_t), intent(in), target :: u, v, w
153 real(kind=rp), intent(in) :: obj_scale
154 class(point_zone_t), intent(in), target :: mask
155 logical, intent(in) :: if_mask
156 type(coef_t), intent(in) :: coef
157 real(kind=rp), intent(in) :: volume
158 real(kind=rp), intent(in) :: start_time
159 real(kind=rp), intent(in) :: end_time
160
161 type(field_list_t) :: fields
162
163 call this%free()
164
165 ! this is copying the fluid source term init
166 ! We package the fields for the source term to operate on in a field list.
167 call fields%init(3)
168 call fields%assign(1, f_x)
169 call fields%assign(2, f_y)
170 call fields%assign(3, f_z)
171
172 call this%init_base(fields, coef, start_time, end_time)
173 call fields%free()
174
175 ! point everything in the correct places
176 this%u => u
177 this%v => v
178 this%w => w
179
180 this%obj_scale = obj_scale
181 this%volume = volume
182
183 this%if_mask = if_mask
184 if (this%if_mask) then
185 this%mask => mask
186 end if
187
188 ! Initialize the ax_helm object
189 call ax_helm_factory(this%Ax, full_formulation = .false.)
190
191 end subroutine adjoint_viscous_dissipation_source_term_init_from_components
192
194 subroutine adjoint_viscous_dissipation_source_term_free(this)
195 class(adjoint_viscous_dissipation_source_term_t), intent(inout) :: this
196
197 call this%free_base()
198 nullify(this%u)
199 nullify(this%v)
200 nullify(this%w)
201 nullify(this%mask)
202 if (allocated(this%Ax)) then
203 deallocate(this%Ax)
204 end if
205
206 end subroutine adjoint_viscous_dissipation_source_term_free
207
211 subroutine adjoint_viscous_dissipation_source_term_compute(this, time)
212 class(adjoint_viscous_dissipation_source_term_t), intent(inout) :: this
213 type(time_state_t), intent(in) :: time
214 type(field_t), pointer :: u, v, w
215 type(field_t), pointer :: fu, fv, fw
216 type(field_t), pointer :: work
217 integer :: temp_indices(1)
218 integer n
219
220 fu => this%fields%get_by_index(1)
221 fv => this%fields%get_by_index(2)
222 fw => this%fields%get_by_index(3)
223
224 n = fu%size()
225
226 u => this%u
227 v => this%v
228 w => this%w
229
230 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
231
232 associate(coef => this%coef)
233
234 ! Note that axhelm computes h1 * lap u + h2 * u
235 ! we set h1 = 1 and h2 = 0 to compute the weak laplacian.
236 if (neko_bcknd_device .eq. 1) then
237 call device_cfill(coef%h1_d, 1.0_rp, n)
238 call device_cfill(coef%h2_d, 0.0_rp, n)
239 else
240 call cfill(coef%h1, 1.0_rp, n)
241 call cfill(coef%h2, 0.0_rp, n)
242 end if
243 coef%ifh2 = .false.
244
245 ! ------------------------------------------------------------------------
246 ! u
247
248 call this%Ax%compute(work%x, u%x, coef, coef%msh, coef%xh)
249
250 ! pre-divide out the mass matrix to counteract it's multiplication
251 if (neko_bcknd_device .eq. 1) then
252 call device_invcol2(work%x_d, coef%B_d, work%size())
253 else
254 call invcol2(work%x, coef%B, work%size())
255 end if
256
257 ! mask
258 if (this%if_mask) then
259 call mask_exterior_const(work, this%mask, 0.0_rp)
260 end if
261
262 ! add to RHS
263 call field_add2s2(fu, work, this%obj_scale / this%volume)
264
265 ! ------------------------------------------------------------------------
266 ! v
267
268 call this%Ax%compute(work%x, v%x, coef, coef%msh, coef%xh)
269
270 ! pre-divide out the mass matrix to counteract it's multiplication
271 if (neko_bcknd_device .eq. 1) then
272 call device_invcol2(work%x_d, coef%B_d, work%size())
273 else
274 call invcol2(work%x, coef%B, work%size())
275 end if
276
277 ! mask
278 if (this%if_mask) then
279 call mask_exterior_const(work, this%mask, 0.0_rp)
280 end if
281
282 ! add to RHS
283 call field_add2s2(fv, work, this%obj_scale / this%volume)
284
285 ! ------------------------------------------------------------------------
286 ! w
287
288 call this%Ax%compute(work%x, w%x, coef, coef%msh, coef%xh)
289
290 ! pre-divide out the mass matrix to counteract it's multiplication
291 if (neko_bcknd_device .eq. 1) then
292 call device_invcol2(work%x_d, coef%B_d, work%size())
293 else
294 call invcol2(work%x, coef%B, work%size())
295 end if
296
297 ! mask
298 if (this%if_mask) then
299 call mask_exterior_const(work, this%mask, 0.0_rp)
300 end if
301
302 ! add to RHS
303 call field_add2s2(fw, work, this%obj_scale / this%volume)
304
305 call neko_scratch_registry%relinquish_field(temp_indices)
306
307 end associate
308
309 end subroutine adjoint_viscous_dissipation_source_term_compute
310
Implements the adjoint_viscous_dissipation_source_term_t type.
subroutine, public adjoint_viscous_dissipation_source_term_allocate(obj)
Allocator for the adjoint viscous dissipation source term.
Some common Masking operations we may need.
Definition mask_ops.f90:36
Implements the steady_simcomp_t type.
The steady_simcomp_t type is a simulation component that terminates a simulation when the normed diff...