Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_mixing_scalar_source_term.f90
1! Copyright (c) 2023, 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!
34! this is a such a dumb name
36 use num_types, only : rp
37 use field_list, only : field_list_t
38 use field, only: field_t
39 use field_registry, only: neko_field_registry
40 use scratch_registry, only: neko_scratch_registry
41 use json_module, only : json_file
42 use json_utils, only: json_get, json_get_or_default
43 use source_term, only : source_term_t
44 use coefs, only : coef_t
45 use field_math, only: field_add2s2, field_copy, field_cadd
46 use mask_ops, only: mask_exterior_const, compute_masked_volume
47 use point_zone, only: point_zone_t
48 implicit none
49 private
50
51 ! this will be the adjoint forcing from Casper's objective function
52 ! TODO
53 ! it actually isn't this. Infact, his is a surface integral on the outlet
54 ! Which means we get strange BC's in our adjoint problem.
55 ! This source term would be if we had a certain volume that we wanted more
56 ! mixed
57 ! the forcing is of the form:
58 ! \f$ \phi - \phi_ref \f$
59 ! ie, difference between it and the average.
60 type, public, extends(source_term_t) :: adjoint_mixing_scalar_source_term_t
62 type(field_t), pointer :: s => null()
64 real(kind=rp) :: obj_scale
66 real(kind=rp) :: phi_ref
68 class(point_zone_t), pointer :: mask => null()
70 logical :: if_mask
72 real(kind=rp) :: mask_volume
73 contains
75 procedure, pass(this) :: init => &
76 adjoint_mixing_scalar_source_term_init_from_json
78 procedure, pass(this) :: init_from_components => &
79 adjoint_mixing_scalar_source_term_init_from_components
81 procedure, pass(this) :: free => adjoint_mixing_scalar_source_term_free
83 procedure, pass(this) :: compute_ => &
84 adjoint_mixing_scalar_source_term_compute
86
87contains
94 subroutine adjoint_mixing_scalar_source_term_init_from_json(this, &
95 json, fields, coef, variable_name)
96 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
97 type(json_file), intent(inout) :: json
98 type(field_list_t), intent(in), target :: fields
99 type(coef_t), intent(in), target :: coef
100 character(len=*), intent(in) :: variable_name
101
102
103 end subroutine adjoint_mixing_scalar_source_term_init_from_json
104
105
106 subroutine adjoint_mixing_scalar_source_term_init_from_components(this, &
107 f_s, s, obj_scale, phi_ref, mask, if_mask, coef)
108 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
109 type(field_t), pointer, intent(in) :: f_s
110 real(kind=rp), intent(in) :: phi_ref
111 type(field_t), intent(in), target :: s
112 class(point_zone_t), intent(in), target :: mask
113 type(field_list_t) :: fields
114 type(coef_t) :: coef
115 real(kind=rp) :: start_time
116 real(kind=rp) :: end_time
117 real(kind=rp) :: obj_scale
118
119 logical :: if_mask
120
121 start_time = 0.0_rp
122 end_time = 100000000.0_rp
123
124
125 call this%free()
126
127 ! this is copying the fluid source term init
128 ! We package the fields for the source term to operate on in a field list.
129 call fields%init(1)
130 call fields%assign(1, f_s)
131
132 call this%init_base(fields, coef, start_time, end_time)
133
134 ! point everything in the correct places
135 this%s => s
136 this%obj_scale = obj_scale
137 this%phi_ref = phi_ref
138 this%if_mask = if_mask
139
140 if (this%if_mask) then
141 this%mask => mask
142 this%mask_volume = compute_masked_volume(this%mask, coef)
143 else
144 this%mask_volume = coef%volume
145 end if
146
147 end subroutine adjoint_mixing_scalar_source_term_init_from_components
148
150 subroutine adjoint_mixing_scalar_source_term_free(this)
151 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
152
153 call this%free_base()
154 end subroutine adjoint_mixing_scalar_source_term_free
155
160 subroutine adjoint_mixing_scalar_source_term_compute(this, t, tstep)
161 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
162 real(kind=rp), intent(in) :: t
163 integer, intent(in) :: tstep
164 type(field_t), pointer :: fs
165 type(field_t), pointer :: work
166 integer :: temp_indices(1)
167
168
169 fs => this%fields%get(1)
170
171 call neko_scratch_registry%request_field(work, temp_indices(1))
172 ! \phi
173 call field_copy(work, this%s)
174 ! \phi - \phi_ref
175 call field_cadd(work, -this%phi_ref)
176 ! mask
177 if (this%if_mask) then
178 call mask_exterior_const(work, this%mask, 0.0_rp)
179 end if
180
181 ! append to RHS with scaling and mask volume
182 call field_add2s2(fs, work, this%obj_scale / this%mask_volume)
183 call neko_scratch_registry%relinquish_field(temp_indices)
184
185
186 end subroutine adjoint_mixing_scalar_source_term_compute
187
Implements the adjoint_mixing_scalar_source_term type.
Some common Masking operations we may need.
Definition mask_ops.f90:34