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
93 subroutine adjoint_mixing_scalar_source_term_init_from_json(this, &
94 json, fields, coef)
95 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
96 type(json_file), intent(inout) :: json
97 type(field_list_t), intent(in), target :: fields
98 type(coef_t), intent(in), target :: coef
99
100
101 end subroutine adjoint_mixing_scalar_source_term_init_from_json
102
103
104 subroutine adjoint_mixing_scalar_source_term_init_from_components(this, &
105 f_s, s, obj_scale, phi_ref, mask, if_mask, coef)
106 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
107 type(field_t), pointer, intent(in) :: f_s
108 real(kind=rp), intent(in) :: phi_ref
109 type(field_t), intent(in), target :: s
110 class(point_zone_t), intent(in), target :: mask
111 type(field_list_t) :: fields
112 type(coef_t) :: coef
113 real(kind=rp) :: start_time
114 real(kind=rp) :: end_time
115 real(kind=rp) :: obj_scale
116
117 logical :: if_mask
118
119 start_time = 0.0_rp
120 end_time = 100000000.0_rp
121
122
123 call this%free()
124
125 ! this is copying the fluid source term init
126 ! We package the fields for the source term to operate on in a field list.
127 call fields%init(1)
128 call fields%assign(1, f_s)
129
130 call this%init_base(fields, coef, start_time, end_time)
131
132 ! point everything in the correct places
133 this%s => s
134 this%obj_scale = obj_scale
135 this%phi_ref = phi_ref
136 this%if_mask = if_mask
137
138 if (this%if_mask) then
139 this%mask => mask
140 this%mask_volume = compute_masked_volume(this%mask, coef)
141 else
142 this%mask_volume = coef%volume
143 end if
144
145 end subroutine adjoint_mixing_scalar_source_term_init_from_components
146
148 subroutine adjoint_mixing_scalar_source_term_free(this)
149 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
150
151 call this%free_base()
152 end subroutine adjoint_mixing_scalar_source_term_free
153
158 subroutine adjoint_mixing_scalar_source_term_compute(this, t, tstep)
159 class(adjoint_mixing_scalar_source_term_t), intent(inout) :: this
160 real(kind=rp), intent(in) :: t
161 integer, intent(in) :: tstep
162 type(field_t), pointer :: fs
163 type(field_t), pointer :: work
164 integer :: temp_indices(1)
165
166
167 fs => this%fields%get(1)
168
169 call neko_scratch_registry%request_field(work, temp_indices(1))
170 ! \phi
171 call field_copy(work, this%s)
172 ! \phi - \phi_ref
173 call field_cadd(work, -this%phi_ref)
174 ! mask
175 if (this%if_mask) then
176 call mask_exterior_const(work, this%mask, 0.0_rp)
177 end if
178
179 ! append to RHS with scaling and mask volume
180 call field_add2s2(fs, work, this%obj_scale / this%mask_volume)
181 call neko_scratch_registry%relinquish_field(temp_indices)
182
183
184 end subroutine adjoint_mixing_scalar_source_term_compute
185
Implements the adjoint_mixing_scalar_source_term type.
Some common Masking operations we may need.
Definition mask_ops.f90:34