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