Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_scalars.f90
1! Copyright (c) 2022-2025, 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
36 use num_types, only: rp
39 use mesh, only: mesh_t
40 use space, only: space_t
41 use gather_scatter, only: gs_t
42 use time_scheme_controller, only: time_scheme_controller_t
43 use time_step_controller, only: time_step_controller_t
44 use json_module, only: json_file
45 use json_utils, only: json_get, json_get_or_default, json_extract_object, &
46 json_extract_item
47 use field, only: field_t
48 use field_series, only: field_series_t
49 use field_registry, only: neko_field_registry
50 use checkpoint, only: chkp_t
51 use krylov, only: ksp_t
52 use logger, only: neko_log, log_size, neko_log_verbose
53 use user_intf, only: user_t
54 use utils, only: neko_error
55 use coefs, only : coef_t
56 use time_state, only : time_state_t
57 implicit none
58 private
59
61 type, public :: adjoint_scalars_t
63 class(adjoint_scalar_scheme_t), allocatable :: adjoint_scalar_fields(:)
65 class(ksp_t), allocatable :: shared_ksp
66 contains
68 generic :: init => adjoint_scalars_init, adjoint_scalars_init_single
69 procedure, private :: adjoint_scalars_init
70 procedure, private :: adjoint_scalars_init_single
72 procedure :: step => adjoint_scalars_step
74 procedure :: restart => adjoint_scalars_restart
76 procedure :: validate => adjoint_scalars_validate
78 procedure :: free => adjoint_scalars_free
79 end type adjoint_scalars_t
80
81contains
82
101 subroutine adjoint_scalars_init(this, n_adjoint_scalars, n_primal_scalars, &
102 msh, coef, gs, params_adjoint, params_primal, &
103 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
104 class(adjoint_scalars_t), intent(inout) :: this
105 integer, intent(in) :: n_adjoint_scalars, n_primal_scalars
106 type(mesh_t), target, intent(in) :: msh
107 type(coef_t), target, intent(in) :: coef
108 type(gs_t), target, intent(inout) :: gs
109 type(json_file), target, intent(inout) :: params_adjoint, params_primal
110 type(json_file), target, intent(inout) :: numerics_params
111 type(user_t), target, intent(in) :: user
112 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
113 type(time_scheme_controller_t), target, intent(in) :: time_scheme
114 TYPE(field_t), TARGET, INTENT(IN) :: rho
115 type(chkp_t), target, intent(inout) :: chkp
116 type(json_file) :: json_subdict_adjoint, json_subdict_primal
117 integer :: i, j
118 character(len=:), allocatable :: field_name
119 character(len=256), allocatable :: field_names(:)
120 character(len=:), allocatable :: primal_name, trial_primal_name
121 logical :: found_primal
122
123 if (n_adjoint_scalars .gt. n_primal_scalars) then
124 call neko_error("More adjoint scalars than forward scalars")
125 ! Note. This assumes every adjoint scalar must have a corresponding
126 ! primal scalar, but not the other way around.
127 end if
128
129 ! Allocate the scalar fields
130 ! If there are more adjoint_scalar_scheme_t types, add a factory function
131 ! here
132 allocate( &
133 adjoint_scalar_pnpn_t::this%adjoint_scalar_fields(n_adjoint_scalars))
134
135 allocate(field_names(n_adjoint_scalars))
136
137 ! For multiple adjoint_scalars, collect and validate field names
138 if (n_adjoint_scalars > 1) then
139 do i = 1, n_adjoint_scalars
140 ! Extract element i from the "adjoint_scalars" array
141 call json_extract_item(params_adjoint, "", i, json_subdict_adjoint)
142
143 ! Try to get name from JSON, generate one if not found or empty
144 if (json_subdict_adjoint%valid_path('name')) then
145 call json_get(json_subdict_adjoint, 'name', field_name)
146 else
147 field_name = ''
148 end if
149
150 ! If name is empty or not provided, generate a default one
151 if (len_trim(field_name) == 0) then
152 write(field_name, '(A,I0,A)') 's_', i, '_adj'
153 end if
154
155 field_names(i) = trim(field_name)
156
157 ! If there's a duplicate, append a number until unique
158 j = 1
159 do while (j < i)
160 if (trim(field_names(i)) == trim(field_names(j))) then
161 write(field_name, '(A,I0)') trim(field_names(i))//'_', j
162 field_names(i) = trim(field_name)
163 j = 1 ! Start over to check if new name is unique
164 else
165 j = j + 1
166 end if
167 end do
168 end do
169 end if
170
171 do i = 1, n_adjoint_scalars
172 call json_extract_item(params_adjoint, "", i, json_subdict_adjoint)
173
174 ! Use the processed field names for multiple adjoint_scalars
175 if (n_adjoint_scalars > 1) then
176 call json_subdict_adjoint%add('name', trim(field_names(i)))
177 else
178 call json_subdict_adjoint%add('name', 's_adj')
179 end if
180
181 ! Before initializing, find the corresponding primal scalar, and find
182 ! its subdict
183 ! note: the default primal name is s.
184 ! So we search for it if not specified
185 call json_get_or_default(params_adjoint, 'primal_name', primal_name, 's')
186
187 found_primal = .false.
188 do j = 1, n_primal_scalars
189 ! Extract element i from the "adjoint_scalars" array
190 call json_extract_item(params_primal, "", j, json_subdict_primal)
191
192 ! Try to get name from JSON
193 if (json_subdict_primal%valid_path('name')) then
194 call json_get(json_subdict_primal, 'name', trial_primal_name)
195 if (trim(trial_primal_name) == trim(primal_name)) then
196 found_primal = .true.
197 exit
198 end if
199 end if
200 end do
201
202 if (.not. found_primal) then
203 call neko_error('Could not find corresponding primal scalar name')
204 else
205 ! now we modify the json in case a name wasn't specified
206 call json_subdict_adjoint%add('primal_name', primal_name)
207 end if
208
209 call this%adjoint_scalar_fields(i)%init(msh, coef, gs, &
210 json_subdict_adjoint, json_subdict_primal, numerics_params, &
211 user, chkp, ulag, vlag, wlag, time_scheme, rho)
212 end do
213 end subroutine adjoint_scalars_init
214
231 subroutine adjoint_scalars_init_single(this, msh, coef, gs, params_adjoint, &
232 params_primal, numerics_params, user, chkp, ulag, vlag, wlag, &
233 time_scheme, rho)
234 class(adjoint_scalars_t), intent(inout) :: this
235 type(mesh_t), target, intent(in) :: msh
236 type(coef_t), target, intent(in) :: coef
237 type(gs_t), target, intent(inout) :: gs
238 type(json_file), target, intent(inout) :: params_adjoint
239 type(json_file), target, intent(inout) :: params_primal
240 type(json_file), target, intent(inout) :: numerics_params
241 type(user_t), target, intent(in) :: user
242 type(chkp_t), target, intent(inout) :: chkp
243 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
244 type(time_scheme_controller_t), target, intent(in) :: time_scheme
245 TYPE(field_t), TARGET, INTENT(IN) :: rho
246
247 ! Allocate a single scalar field
248 allocate(adjoint_scalar_pnpn_t::this%adjoint_scalar_fields(1))
249
250 ! The default for primal scalar name will be s, so we force these two
251 call params_adjoint%add('name', 's_adj')
252 call params_adjoint%add('primal_name', 's')
253
254 ! TODO, there may be a corner case with multiple primal and a single adjoint
255 ! we'll need to catch that before entering this subroutine.
256
257 ! Initialize it directly with the params
258 call this%adjoint_scalar_fields(1)%init(msh, coef, gs, params_adjoint, &
259 params_primal, numerics_params, user, chkp, ulag, vlag, wlag, &
260 time_scheme, rho)
261 end subroutine adjoint_scalars_init_single
262
264 subroutine adjoint_scalars_step(this, time, ext_bdf, dt_controller)
265 class(adjoint_scalars_t), intent(inout) :: this
266 type(time_state_t), intent(in) :: time
267 type(time_scheme_controller_t), intent(inout) :: ext_bdf
268 type(time_step_controller_t), intent(inout) :: dt_controller
269 integer :: i
270
271 ! Iterate through all scalar fields
272 do i = 1, size(this%adjoint_scalar_fields)
273 call this%adjoint_scalar_fields(i)%step(time, ext_bdf, dt_controller)
274 end do
275 end subroutine adjoint_scalars_step
276
278 subroutine adjoint_scalars_restart(this, chkp)
279 class(adjoint_scalars_t), intent(inout) :: this
280 type(chkp_t), intent(inout) :: chkp
281 integer :: i
282 ! Iterate through all scalar fields
283 do i = 1, size(this%adjoint_scalar_fields)
284 call this%adjoint_scalar_fields(i)%restart(chkp)
285 end do
286 end subroutine adjoint_scalars_restart
287
289 subroutine adjoint_scalars_validate(this)
290 class(adjoint_scalars_t), intent(inout) :: this
291 integer :: i
292 ! Iterate through all scalar fields
293 do i = 1, size(this%adjoint_scalar_fields)
294 call this%adjoint_scalar_fields(i)%s_adj_lag%set(&
295 this%adjoint_scalar_fields(i)%s_adj)
296 call this%adjoint_scalar_fields(i)%validate()
297 end do
298 end subroutine adjoint_scalars_validate
299
301 subroutine adjoint_scalars_free(this)
302 class(adjoint_scalars_t), intent(inout) :: this
303 integer :: i
304
305 ! Iterate through all scalar fields
306 if (allocated(this%adjoint_scalar_fields)) then
307 do i = 1, size(this%adjoint_scalar_fields)
308 call this%adjoint_scalar_fields(i)%free()
309 end do
310 deallocate(this%adjoint_scalar_fields)
311 end if
312 end subroutine adjoint_scalars_free
313
314end module adjoint_scalars
Contains the adjoint_scalar_pnpn_t type.
Contains the adjoint_scalar_scheme_t type.
Contains the adjoint_scalars_t type that manages multiple scalar fields.
Base type for a scalar advection-diffusion solver.
Type to manage multiple adjoint scalar transport equations.