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