Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
simulation.f90
Go to the documentation of this file.
1
34!
36! Here, we simply march forward to steady state solutions
38 use case, only: case_t
39 use neko, only: neko_solve
40 use adjoint_case, only: adjoint_case_t
41 use fluid_scheme_incompressible, only: fluid_scheme_incompressible_t
42 use adjoint_fluid_scheme, only: adjoint_fluid_scheme_t
44 use scalar_pnpn, only: scalar_pnpn_t
47 use scalars, only: scalars_t
48 use scalar_scheme, only: scalar_scheme_t
49 use fluid_pnpn, only: fluid_pnpn_t
50 use time_step_controller, only: time_step_controller_t
51 use time_state, only: time_state_t
52 use fld_file_output, only: fld_file_output_t
53 use chkp_output, only: chkp_output_t
54 use simcomp_executor, only: neko_simcomps
55 use neko_ext, only: reset, reset_adjoint
56 use field, only: field_t
57 use registry, only: neko_registry
58 use field_math, only: field_rzero, field_copy
59 use checkpoint, only: chkp_t
60 use file, only: file_t
61 use utils, only: neko_warning, neko_error
62 use comm, only: pe_rank
63 use json_file_module, only: json_file
64 use json_utils, only: json_get, json_get_or_default
65 use num_types, only: rp, sp, dp
66 use logger, only: log_size, neko_log
67 use mpi_f08, only: mpi_wtime
68 use jobctrl, only: jobctrl_time_limit
69 use profiler, only: profiler_start, profiler_stop, &
70 profiler_start_region, profiler_end_region
73 use simulation, only: simulation_init, simulation_step, simulation_finalize, &
74 simulation_restart
75 use simulation_checkpoint, only: simulation_checkpoint_t
76 use runtime_stats, only: neko_rt_stats
77 use scratch_registry, only: neko_scratch_registry
78 use registry, only: neko_registry
79 implicit none
80 private
81
83
85 type(case_t), public :: neko_case
87 type(adjoint_case_t), public :: adjoint_case
89 class(fluid_scheme_incompressible_t), public, pointer :: fluid => null()
91 type(scalars_t), public, pointer :: scalars => null()
93 class(adjoint_fluid_scheme_t), public, pointer :: adjoint_fluid => null()
95 type(adjoint_scalars_t), public, pointer :: adjoint_scalars => null()
98 type(fld_file_output_t), public :: output_forward
101 type(fld_file_output_t), public :: output_adjoint
103 logical :: unsteady = .false.
104
105 logical :: have_scalar = .false.
106 integer :: n_timesteps = 0
107
108 ! ----------------------------------------------------------------------- !
109 ! Checkpoint system
110
112 type(simulation_checkpoint_t) :: checkpoint
113
114 contains
116 procedure, pass(this) :: init => simulation_initialize
118 procedure, pass(this) :: free => simulation_free
120 procedure, pass(this) :: run_forward => simulation_run_forward
122 procedure, pass(this) :: run_backward => simulation_run_backward
124 procedure, pass(this) :: reset => simulation_reset
126 procedure, pass(this) :: write => simulation_write
128 procedure, pass(this) :: write_forward => simulation_write_forward
130 procedure, pass(this) :: write_adjoint => simulation_write_adjoint
131
132 end type simulation_t
133 public :: simulation_t
134contains
135
137 subroutine simulation_initialize(this, parameters)
138 class(simulation_t), intent(inout), target :: this
139 type(json_file), intent(inout) :: parameters
140 type(json_file) :: checkpoint_params
141 integer :: i, n_scalars, unsteady_support
142 logical :: unsteady
143
144 ! initialize the primal Neko objects
145 call this%neko_case%init(parameters)
146 call neko_rt_stats%init(parameters)
147 call neko_simcomps%init(this%neko_case)
148
149 ! initialize the adjoint
150 call this%adjoint_case%init(this%neko_case)
151
152 ! Start the profiler
153 call profiler_start
154
155 select type (fluid => this%neko_case%fluid)
156 type is (fluid_pnpn_t)
157 this%fluid => fluid
158 end select
159
160 select type (adjoint_fluid => this%adjoint_case%fluid_adj)
161 type is (adjoint_fluid_pnpn_t)
162 this%adjoint_fluid => adjoint_fluid
163 end select
164
165 if (allocated(this%neko_case%scalars)) then
166 this%scalars => this%neko_case%scalars
167 end if
168
169 if (allocated(this%adjoint_case%adjoint_scalars)) then
170 this%adjoint_scalars => this%adjoint_case%adjoint_scalars
171 end if
172
173 ! init the sampler
174 !---------------------------------------------------------
175 ! Allocate the output type
176 n_scalars = 0
177 if (allocated(this%neko_case%scalars)) then
178 n_scalars = size(this%neko_case%scalars%scalar_fields)
179 end if
180 call this%output_forward%init(sp, 'forward_fields', 4 + n_scalars)
181
182 call this%output_forward%fields%assign(1, this%fluid%p)
183 call this%output_forward%fields%assign(2, this%fluid%u)
184 call this%output_forward%fields%assign(3, this%fluid%v)
185 call this%output_forward%fields%assign(4, this%fluid%w)
186
187 ! Assign all scalar fields
188 if (allocated(this%neko_case%scalars)) then
189 do i = 1, n_scalars
190 call this%output_forward%fields%assign(4 + i, &
191 this%scalars%scalar_fields(i)%s)
192 end do
193 end if
194
195 n_scalars = 0
196 if (allocated(this%adjoint_case%adjoint_scalars)) then
197 n_scalars = size(this%adjoint_case%adjoint_scalars%adjoint_scalar_fields)
198 end if
199 call this%output_adjoint%init(sp, 'adjoint_fields', 4 + n_scalars)
200 call this%output_adjoint%fields%assign(1, this%adjoint_fluid%p_adj)
201 call this%output_adjoint%fields%assign(2, this%adjoint_fluid%u_adj)
202 call this%output_adjoint%fields%assign(3, this%adjoint_fluid%v_adj)
203 call this%output_adjoint%fields%assign(4, this%adjoint_fluid%w_adj)
204
205 ! Assign all scalar fields
206 if (allocated(this%adjoint_case%adjoint_scalars)) then
207 do i = 1, n_scalars
208 call this%output_adjoint%fields%assign(4 + i, &
209 this%adjoint_scalars%adjoint_scalar_fields(i)%s_adj)
210 end do
211 end if
212
213 ! Check if the simulation is steady or unsteady
214 call json_get_or_default(parameters, "unsteady", unsteady, .false.)
215 this%unsteady = unsteady
216
217 ! Ensure there is a means to deal with unsteadiness
218 if (this%unsteady) then
219 unsteady_support = 0
220 if ("checkpoints" .in. parameters) then
221 unsteady_support = unsteady_support + 1
222 end if
223
224 if (unsteady_support .eq. 0) then
225 call neko_error("No support for unsteady simulation provided, \\ &
226 & \\ current options include enabling checkpoints.")
227 end if
228
229 if (unsteady_support .gt. 1) then
230 call neko_error("Too many supports for unsteady simulation \\ &
231 & \\ provided, please select one.")
232 end if
233 end if
234
235 if ("checkpoints" .in. parameters) then
236 call json_get(parameters, 'checkpoints', checkpoint_params)
237 call this%checkpoint%init(this%neko_case, checkpoint_params)
238 end if
239
240 end subroutine simulation_initialize
241
243 subroutine simulation_free(this)
244 class(simulation_t), intent(inout) :: this
245
246 ! Stop the profiler
247 call profiler_stop
248
249 ! Free the objects
250 call this%neko_case%free()
251 call this%adjoint_case%free()
252 call this%output_forward%free()
253 call this%output_adjoint%free()
254 call this%checkpoint%free()
255
256 ! Nullify pointers
257 nullify(this%fluid)
258 nullify(this%scalars)
259 nullify(this%adjoint_fluid)
260 nullify(this%adjoint_scalars)
261
262 ! Reset flags and counters
263 this%unsteady = .false.
264 this%have_scalar = .false.
265 this%n_timesteps = 0
266
267 ! Close global objects
268 call neko_simcomps%free()
269 call neko_rt_stats%free()
270
271 end subroutine simulation_free
272
274 subroutine simulation_run_forward(this)
275 class(simulation_t), intent(inout) :: this
276 type(time_step_controller_t) :: dt_controller
277 real(kind=dp) :: loop_start
278
279 call dt_controller%init(this%neko_case%params)
280
281 call simulation_init(this%neko_case, dt_controller)
282
283 call profiler_start_region("Forward simulation")
284 loop_start = mpi_wtime()
285 this%n_timesteps = 0
286 do while (this%neko_case%time%t .lt. this%neko_case%time%end_time)
287 this%n_timesteps = this%n_timesteps + 1
288
289 call simulation_step(this%neko_case, dt_controller, loop_start)
290
291 call this%checkpoint%save(this%neko_case)
292 end do
293 call profiler_end_region("Forward simulation")
294
295 call simulation_finalize(this%neko_case)
296
297 end subroutine simulation_run_forward
298
300 subroutine simulation_run_backward(this)
301 class(simulation_t), intent(inout) :: this
302 type(time_step_controller_t) :: dt_controller
303 real(kind=dp) :: loop_start
304 real(kind=rp) :: cfl
305 integer :: i
306
307 call dt_controller%init(this%neko_case%params)
308
309 call simulation_adjoint_init(this%adjoint_case, dt_controller)
310
311 call profiler_start_region("Adjoint simulation")
312 cfl = this%adjoint_case%fluid_adj%compute_cfl(this%adjoint_case%time%dt)
313 loop_start = mpi_wtime()
314 do i = this%n_timesteps, 1, -1
315 call this%checkpoint%restore(this%neko_case, i)
316
317 call simulation_adjoint_step(this%adjoint_case, dt_controller, cfl, &
318 loop_start)
319 end do
320 call profiler_end_region("Adjoint simulation")
321
322 call simulation_adjoint_finalize(this%adjoint_case)
323
324 end subroutine simulation_run_backward
325
327 subroutine simulation_reset(this)
328 class(simulation_t), intent(inout) :: this
329
330 call reset(this%neko_case)
331 call reset_adjoint(this%adjoint_case, this%neko_case)
332 call this%checkpoint%reset()
333
334 end subroutine simulation_reset
335
337 subroutine simulation_write(this, idx)
338 class(simulation_t), intent(inout) :: this
339 integer, intent(in) :: idx
340
341 call this%output_forward%sample(real(idx, kind=rp))
342 call this%output_adjoint%sample(real(idx, kind=rp))
343
344 end subroutine simulation_write
345
347 subroutine simulation_write_forward(this, idx)
348 class(simulation_t), intent(inout) :: this
349 integer, intent(in) :: idx
350
351 call this%output_forward%sample(real(idx, kind=rp))
352
353 end subroutine simulation_write_forward
354
356 subroutine simulation_write_adjoint(this, idx)
357 class(simulation_t), intent(inout) :: this
358 integer, intent(in) :: idx
359
360 call this%output_adjoint%sample(real(idx, kind=rp))
361
362 end subroutine simulation_write_adjoint
363
364end module simulation_m
Adjoint Pn/Pn formulation.
Contains the adjoint_scalar_pnpn_t type.
Contains the adjoint_scalars_t type that manages multiple scalar fields.
Contains extensions to the neko library required to run the topology optimization code.
Definition neko_ext.f90:42
subroutine, public reset(neko_case)
Reset the case data structure.
Definition neko_ext.f90:91
subroutine, public reset_adjoint(adjoint_case, neko_case)
Reset the adjoint case data structure.
Definition neko_ext.f90:243
Adjoint simulation driver.
subroutine, public simulation_adjoint_init(c, dt_controller)
Initialise a simulation_adjoint of a case.
subroutine, public simulation_adjoint_step(c, dt_controller, cfl, tstep_loop_start_time, final_time)
Compute a single time-step of an adjoint case.
subroutine, public simulation_adjoint_finalize(c)
Finalize a simulation of a case.
Implements the steady_problem_t type.
subroutine simulation_initialize(this, parameters)
Initialize the simulation.
Adjoint case type. Todo: This should Ideally be a subclass of case_t, however, that is not yet suppor...
Type to manage multiple adjoint scalar transport equations.