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