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, adjoint_init, adjoint_free
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 implicit none
78 private
79
81
83 type(case_t), public :: neko_case
85 type(adjoint_case_t), public :: adjoint_case
87 class(fluid_scheme_incompressible_t), public, pointer :: fluid => null()
89 type(scalars_t), public, pointer :: scalars => null()
91 class(adjoint_fluid_scheme_t), public, pointer :: adjoint_fluid => null()
93 type(adjoint_scalars_t), public, pointer :: adjoint_scalars => null()
96 type(fld_file_output_t), public :: output_forward
99 type(fld_file_output_t), public :: output_adjoint
101 logical :: unsteady = .false.
102
103 logical :: have_scalar = .false.
104 integer :: n_timesteps = 0
105
106 ! ----------------------------------------------------------------------- !
107 ! Checkpoint system
108
110 type(simulation_checkpoint_t) :: checkpoint
111
112 contains
114 procedure, pass(this) :: init => simulation_initialize
116 procedure, pass(this) :: free => simulation_free
118 procedure, pass(this) :: run_forward => simulation_run_forward
120 procedure, pass(this) :: run_backward => simulation_run_backward
122 procedure, pass(this) :: reset => simulation_reset
124 procedure, pass(this) :: write => simulation_write
126 procedure, pass(this) :: write_forward => simulation_write_forward
128 procedure, pass(this) :: write_adjoint => simulation_write_adjoint
129
130 end type simulation_t
131 public :: simulation_t
132contains
133
135 subroutine simulation_initialize(this, parameters)
136 class(simulation_t), intent(inout), target :: this
137 type(json_file), intent(inout) :: parameters
138 type(json_file) :: checkpoint_params
139 integer :: i, n_scalars, unsteady_support
140 logical :: unsteady
141
142 ! initialize the primal Neko objects
143 call this%neko_case%init(parameters)
144 call neko_rt_stats%init(parameters)
145 call neko_simcomps%init(this%neko_case)
146
147 ! initialize the adjoint
148 call adjoint_init(this%adjoint_case, this%neko_case)
149
150 ! Start the profiler
151 call profiler_start
152
153 select type (fluid => this%neko_case%fluid)
154 type is (fluid_pnpn_t)
155 this%fluid => fluid
156 end select
157
158 select type (adjoint_fluid => this%adjoint_case%fluid_adj)
159 type is (adjoint_fluid_pnpn_t)
160 this%adjoint_fluid => adjoint_fluid
161 end select
162
163 if (allocated(this%neko_case%scalars)) then
164 this%scalars => this%neko_case%scalars
165 end if
166
167 if (allocated(this%adjoint_case%adjoint_scalars)) then
168 this%adjoint_scalars => this%adjoint_case%adjoint_scalars
169 end if
170
171 ! init the sampler
172 !---------------------------------------------------------
173 ! Allocate the output type
174 n_scalars = 0
175 if (allocated(this%neko_case%scalars)) then
176 n_scalars = size(this%neko_case%scalars%scalar_fields)
177 end if
178 call this%output_forward%init(sp, 'forward_fields', 4 + n_scalars)
179
180 call this%output_forward%fields%assign(1, this%fluid%p)
181 call this%output_forward%fields%assign(2, this%fluid%u)
182 call this%output_forward%fields%assign(3, this%fluid%v)
183 call this%output_forward%fields%assign(4, this%fluid%w)
184
185 ! Assign all scalar fields
186 if (allocated(this%neko_case%scalars)) then
187 do i = 1, n_scalars
188 call this%output_forward%fields%assign(4 + i, &
189 this%scalars%scalar_fields(i)%s)
190 end do
191 end if
192
193 n_scalars = 0
194 if (allocated(this%adjoint_case%adjoint_scalars)) then
195 n_scalars = size(this%adjoint_case%adjoint_scalars%adjoint_scalar_fields)
196 end if
197 call this%output_adjoint%init(sp, 'adjoint_fields', 4 + n_scalars)
198 call this%output_adjoint%fields%assign(1, this%adjoint_fluid%p_adj)
199 call this%output_adjoint%fields%assign(2, this%adjoint_fluid%u_adj)
200 call this%output_adjoint%fields%assign(3, this%adjoint_fluid%v_adj)
201 call this%output_adjoint%fields%assign(4, this%adjoint_fluid%w_adj)
202
203 ! Assign all scalar fields
204 if (allocated(this%adjoint_case%adjoint_scalars)) then
205 do i = 1, n_scalars
206 call this%output_adjoint%fields%assign(4 + i, &
207 this%adjoint_scalars%adjoint_scalar_fields(i)%s_adj)
208 end do
209 end if
210
211 ! Check if the simulation is steady or unsteady
212 call json_get_or_default(parameters, "unsteady", unsteady, .false.)
213 this%unsteady = unsteady
214
215 ! Ensure there is a means to deal with unsteadiness
216 if (this%unsteady) then
217 unsteady_support = 0
218 if ("checkpoints" .in. parameters) then
219 unsteady_support = unsteady_support + 1
220 end if
221
222 if (unsteady_support .eq. 0) then
223 call neko_error("No support for unsteady simulation provided, \\ &
224 & \\ current options include enabling checkpoints.")
225 end if
226
227 if (unsteady_support .gt. 1) then
228 call neko_error("Too many supports for unsteady simulation \\ &
229 & \\ provided, please select one.")
230 end if
231 end if
232
233 if ("checkpoints" .in. parameters) then
234 call json_get(parameters, 'checkpoints', checkpoint_params)
235 call this%checkpoint%init(this%neko_case, checkpoint_params)
236 end if
237
238 end subroutine simulation_initialize
239
241 subroutine simulation_free(this)
242 class(simulation_t), intent(inout) :: this
243
244 ! Stop the profiler
245 call profiler_stop
246
247 call this%checkpoint%free()
248 call adjoint_free(this%adjoint_case)
249
250 call neko_simcomps%free()
251 call neko_rt_stats%free()
252 call this%neko_case%free()
253
254 end subroutine simulation_free
255
257 subroutine simulation_run_forward(this)
258 class(simulation_t), intent(inout) :: this
259 type(time_step_controller_t) :: dt_controller
260 real(kind=dp) :: loop_start
261
262 call dt_controller%init(this%neko_case%params)
263
264 call simulation_init(this%neko_case, dt_controller)
265
266 call profiler_start_region("Forward simulation")
267 loop_start = mpi_wtime()
268 this%n_timesteps = 0
269 do while (this%neko_case%time%t .lt. this%neko_case%time%end_time)
270 this%n_timesteps = this%n_timesteps + 1
271
272 call simulation_step(this%neko_case, dt_controller, loop_start)
273
274 call this%checkpoint%save(this%neko_case)
275 end do
276 call profiler_end_region("Forward simulation")
277
278 call simulation_finalize(this%neko_case)
279
280 end subroutine simulation_run_forward
281
283 subroutine simulation_run_backward(this)
284 class(simulation_t), intent(inout) :: this
285 type(time_step_controller_t) :: dt_controller
286 real(kind=dp) :: loop_start
287 real(kind=rp) :: cfl
288 integer :: i
289
290 call dt_controller%init(this%neko_case%params)
291
292 call simulation_adjoint_init(this%adjoint_case, dt_controller)
293
294 call profiler_start_region("Adjoint simulation")
295 cfl = this%adjoint_case%fluid_adj%compute_cfl(this%adjoint_case%time%dt)
296 loop_start = mpi_wtime()
297 do i = this%n_timesteps, 1, -1
298 call this%checkpoint%restore(this%neko_case, i)
299
300 call simulation_adjoint_step(this%adjoint_case, dt_controller, cfl, &
301 loop_start)
302 end do
303 call profiler_end_region("Adjoint simulation")
304
305 call simulation_adjoint_finalize(this%adjoint_case)
306
307 end subroutine simulation_run_backward
308
310 subroutine simulation_reset(this)
311 class(simulation_t), intent(inout) :: this
312
313 call reset(this%neko_case)
314 call reset_adjoint(this%adjoint_case, this%neko_case)
315 call this%checkpoint%reset()
316
317 end subroutine simulation_reset
318
320 subroutine simulation_write(this, idx)
321 class(simulation_t), intent(inout) :: this
322 integer, intent(in) :: idx
323
324 call this%output_forward%sample(real(idx, kind=rp))
325 call this%output_adjoint%sample(real(idx, kind=rp))
326
327 end subroutine simulation_write
328
330 subroutine simulation_write_forward(this, idx)
331 class(simulation_t), intent(inout) :: this
332 integer, intent(in) :: idx
333
334 call this%output_forward%sample(real(idx, kind=rp))
335
336 end subroutine simulation_write_forward
337
339 subroutine simulation_write_adjoint(this, idx)
340 class(simulation_t), intent(inout) :: this
341 integer, intent(in) :: idx
342
343 call this%output_adjoint%sample(real(idx, kind=rp))
344
345 end subroutine simulation_write_adjoint
346
347end 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.