Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_scalar_convection_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 scratch_registry, only: neko_scratch_registry
40 use json_module, only: json_file
41 use time_state, only: time_state_t
42 use source_term, only: source_term_t
43 use coefs, only: coef_t
44 use field_math, only: field_subcol3
45 use operators, only: grad
46 implicit none
47 private
48
49 ! I don't know how to name this term, but when you have a passive
50 ! scalar you get an extra term in the adjoint velocity equation, which comes
51 ! from the convective term in the passive scalar equation.
52 ! Maybe this should be called just `adjoint_scalar_convection`?
53 ! but it does come in a source term...
54
55 ! In any case,
56 ! it's a source term acting on the adjoint velocity equations, of the form:
57 ! \f$\nabla s s_adj\f$
58 type, public, extends(source_term_t) :: &
61 type(field_t), pointer :: s_adj => null()
63 type(field_t), pointer :: s => null()
64 contains
66 procedure, pass(this) :: init => &
67 adjoint_scalar_convection_source_term_init_from_json
69 procedure, pass(this) :: init_from_components => &
70 adjoint_scalar_convection_source_term_init_from_components
72 procedure, pass(this) :: free => adjoint_scalar_convection_source_term_free
74 procedure, pass(this) :: compute_ => &
75 adjoint_scalar_convection_source_term_compute
77
78contains
85 subroutine adjoint_scalar_convection_source_term_init_from_json(this, &
86 json, fields, coef, variable_name)
87 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
88 type(json_file), intent(inout) :: json
89 type(field_list_t), intent(in), target :: fields
90 type(coef_t), intent(in), target :: coef
91 character(len=*), intent(in) :: variable_name
92
93 ! this is a bit weird... because I don't think this should come from the
94 ! JSON.
95 ! Maybe we should think of all these source terms as only "appendable"
96 !
97 ! Because we'll never have the whole case here, so we'll never be able
98 ! init from components anyway...
99
100
101 end subroutine adjoint_scalar_convection_source_term_init_from_json
102
103
104 subroutine adjoint_scalar_convection_source_term_init_from_components(this,&
105 f_x, f_y, f_z, s, s_adj, coef)
106 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
107 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
108 type(field_list_t) :: fields
109 type(coef_t) :: coef
110 real(kind=rp) :: start_time
111 real(kind=rp) :: end_time
112 type(field_t), intent(in), target :: s, s_adj
113
114 ! I wish you didn't need a start time and end time...
115 ! but I'm just going to set a super big number...
116 start_time = 0.0_rp
117 end_time = 100000000.0_rp
118
119 call this%free()
120
121 ! this is copying the fluid source term init
122 ! We package the fields for the source term to operate on in a field list.
123 call fields%init(3)
124 call fields%assign(1, f_x)
125 call fields%assign(2, f_y)
126 call fields%assign(3, f_z)
127
128 call this%init_base(fields, coef, start_time, end_time)
129
130 ! point everything in the correct places
131 this%s_adj => s_adj
132 this%s => s
133
134 end subroutine adjoint_scalar_convection_source_term_init_from_components
135
137 subroutine adjoint_scalar_convection_source_term_free(this)
138 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
139
140 call this%free_base()
141 end subroutine adjoint_scalar_convection_source_term_free
142
146 subroutine adjoint_scalar_convection_source_term_compute(this, time)
147 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
148 type(time_state_t), intent(in) :: time
149 type(field_t), pointer :: fu, fv, fw
150 integer :: temp_indices(3)
151 type(field_t), pointer :: dsdx, dsdy, dsdz
152
153
154 call neko_scratch_registry%request_field(dsdx, temp_indices(1))
155 call neko_scratch_registry%request_field(dsdy, temp_indices(2))
156 call neko_scratch_registry%request_field(dsdz, temp_indices(3))
157
158 fu => this%fields%get(1)
159 fv => this%fields%get(2)
160 fw => this%fields%get(3)
161
162 ! we basically just need the term
163 ! \f$\nabla s s_adj\f$
164 ! TODO
165 ! In principle, this should have to option to be evaluated on the dealiased
166 ! mesh!
167
168 ! So should the Brinkman term actually....
169
170 call grad(dsdx%x, dsdy%x, dsdz%x, this%s%x, this%coef)
171
172 ! TODO
173 ! So in principal, the derivatives could have kinks now.
174 ! I don't think a gsop will remedy this (or even whether it's a good idea)
175 ! But I want to leave this todo as a reminder.
176
177 call field_subcol3(fu, this%s_adj, dsdx)
178 call field_subcol3(fv, this%s_adj, dsdy)
179 call field_subcol3(fw, this%s_adj, dsdz)
180
181 ! free the scratch
182 call neko_scratch_registry%relinquish_field(temp_indices)
183 end subroutine adjoint_scalar_convection_source_term_compute
184
Implements the adjoint_scalar_convection_source_term type.