Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
case.f90
Go to the documentation of this file.
1! Copyright (c) 2020-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!
34module case
35 use num_types, only : rp, sp, dp
36 use fluid_pnpn, only : fluid_pnpn_t
38 use fluid_scheme_base, only: fluid_scheme_base_t, fluid_scheme_base_factory
40 use chkp_output, only : chkp_output_t
41 use mesh_field, only : mesh_fld_t
43 use redist, only : redist_mesh
45 use flow_ic, only : set_flow_ic
46 use scalar_ic, only : set_scalar_ic
47 use file, only : file_t
48 use utils, only : neko_error
49 use mesh, only : mesh_t
50 use math, only : neko_eps
51 use checkpoint, only: chkp_t
53 use logger, only : neko_log, neko_log_quiet
55 use user_intf, only : user_t
56 use scalar_pnpn, only : scalar_pnpn_t
58 use time_state, only : time_state_t
59 use json_module, only : json_file
63 use scalars, only : scalars_t
64 use comm, only : neko_comm, pe_rank, pe_size
65 use mpi_f08, only : mpi_bcast, mpi_character, mpi_integer
66
67 implicit none
68 private
69
70 type, public :: case_t
71 type(mesh_t) :: msh
72 type(json_file) :: params
73 character(len=:), allocatable :: output_directory
75 type(fluid_output_t) :: f_out
76 type(time_state_t) :: time
77 type(chkp_output_t) :: chkp_out
78 type(chkp_t) :: chkp
79 type(user_t) :: user
80 class(fluid_scheme_base_t), allocatable :: fluid
81 type(scalars_t), allocatable :: scalars
82 contains
83 procedure, private, pass(this) :: case_init_from_file
84 procedure, private, pass(this) :: case_init_from_json
85 procedure, pass(this) :: free => case_free
87 end type case_t
88
89contains
90
92 subroutine case_init_from_file(this, case_file)
93 class(case_t), target, intent(inout) :: this
94 character(len=*), intent(in) :: case_file
95 integer :: ierr, integer_val
96 character(len=:), allocatable :: json_buffer
97 logical :: exist
98
99 ! Check if the file exists
100 inquire(file = trim(case_file), exist = exist)
101 if (.not. exist) then
102 call neko_error('The case file '//trim(case_file)//' does not exist.')
103 end if
104
105 call neko_log%section('Case')
106 call neko_log%message('Reading case file ' // trim(case_file), &
108
109 if (pe_rank .eq. 0) then
110 call this%params%load_file(filename = trim(case_file))
111 call this%params%print_to_string(json_buffer)
112 integer_val = len(json_buffer)
113 end if
114
115 call mpi_bcast(integer_val, 1, mpi_integer, 0, neko_comm, ierr)
116 if (pe_rank .ne. 0) allocate(character(len = integer_val) :: json_buffer)
117 call mpi_bcast(json_buffer, integer_val, mpi_character, 0, neko_comm, ierr)
118 call this%params%load_from_string(json_buffer)
119
120 deallocate(json_buffer)
121
122 call case_init_common(this)
123
124 end subroutine case_init_from_file
125
127 subroutine case_init_from_json(this, case_json)
128 class(case_t), target, intent(inout) :: this
129 type(json_file), intent(in) :: case_json
130
131 call neko_log%section('Case')
132 call neko_log%message('Creating case from JSON object', neko_log_quiet)
133
134 this%params = case_json
135
136 call case_init_common(this)
137
138 end subroutine case_init_from_json
139
141 subroutine case_init_common(this)
142 type(case_t), target, intent(inout) :: this
143 integer :: lx = 0
144 logical :: scalar = .false.
145 type(file_t) :: msh_file, bdry_file, part_file
146 type(mesh_fld_t) :: msh_part, parts
147 logical :: found, logical_val
148 logical :: temperature_found = .false.
149 integer :: integer_val
150 real(kind=rp) :: real_val
151 character(len = :), allocatable :: string_val, name, file_format
152 integer :: output_dir_len
153 integer :: precision, layout
154 type(json_file) :: scalar_params, numerics_params
155 type(json_file) :: json_subdict
156 integer :: n_scalars, i
157
158 !
159 ! Setup user defined functions
160 !
161 call this%user%init()
162
163 ! Run user startup routine
164 call this%user%startup(this%params)
165
166 ! Check if default value fill-in is allowed
167 if (this%params%valid_path('case.no_defaults')) then
168 call json_get(this%params, 'case.no_defaults', json_no_defaults)
169 end if
170
171 !
172 ! Load mesh
173 !
174 call json_get_or_default(this%params, 'case.mesh_file', string_val, &
175 'no mesh')
176 if (trim(string_val) .eq. 'no mesh') then
177 call neko_error('The mesh_file keyword could not be found in the .' // &
178 'case file. Often caused by incorrectly formatted json.')
179 end if
180 call msh_file%init(string_val)
181 call msh_file%read(this%msh)
182
183 !
184 ! Load Balancing
185 !
186 call json_get_or_default(this%params, 'case.load_balancing', logical_val,&
187 .false.)
188
189 if (pe_size .gt. 1 .and. logical_val) then
190 call neko_log%section('Load Balancing')
191 call parmetis_partmeshkway(this%msh, parts)
192 call redist_mesh(this%msh, parts)
193
194 ! store the balanced mesh (for e.g. restarts)
195 string_val = trim(string_val(1:scan(trim(string_val), &
196 '.', back = .true.) - 1)) // '_lb.nmsh'
197 call msh_file%init(string_val)
198 call msh_file%write(this%msh)
199
200 call neko_log%end_section()
201 end if
202
203 ! Run user mesh motion routine
204 call this%user%mesh_setup(this%msh, this%time)
205
206 !
207 ! Time control
208 !
209 call json_get(this%params, 'case.time', json_subdict)
210 call this%time%init(json_subdict)
211
212 !
213 ! Initialize point_zones registry
214 !
215 call neko_point_zone_registry%init(this%params, this%msh)
216
217 !
218 ! Setup fluid scheme
219 !
220 call json_get(this%params, 'case.fluid.scheme', string_val)
221 call fluid_scheme_base_factory(this%fluid, trim(string_val))
222
223 call json_get(this%params, 'case.numerics.polynomial_order', lx)
224 lx = lx + 1 ! add 1 to get number of gll points
225 ! Set time lags in chkp
226 this%chkp%tlag => this%time%tlag
227 this%chkp%dtlag => this%time%dtlag
228 call this%fluid%init(this%msh, lx, this%params, this%user, this%chkp)
229
230
231 !
232 ! Setup scratch registry
233 !
234 call neko_scratch_registry%set_dofmap(this%fluid%dm_Xh)
235
236 !
237 ! Setup scalar scheme
238 !
239 ! @todo no scalar factory for now, probably not needed
240 scalar = .false.
241 n_scalars = 0
242 if (this%params%valid_path('case.scalar')) then
243 call json_get_or_default(this%params, 'case.scalar.enabled', scalar, &
244 .true.)
245 n_scalars = 1
246 else if (this%params%valid_path('case.scalars')) then
247 call this%params%info('case.scalars', n_children = n_scalars)
248 if (n_scalars > 0) then
249 scalar = .true.
250 end if
251 end if
252
253 if (scalar) then
254 allocate(this%scalars)
255 call json_get(this%params, 'case.numerics', numerics_params)
256 if (this%params%valid_path('case.scalar')) then
257 ! For backward compatibility
258 call json_get(this%params, 'case.scalar', scalar_params)
259 call this%scalars%init(this%msh, this%fluid%c_Xh, this%fluid%gs_Xh, &
260 scalar_params, numerics_params, this%user, this%chkp, this%fluid%ulag, &
261 this%fluid%vlag, this%fluid%wlag, this%fluid%ext_bdf, &
262 this%fluid%rho)
263 else
264 ! Multiple scalars
265 call json_get(this%params, 'case.scalars', json_subdict)
266 call this%scalars%init(n_scalars, this%msh, this%fluid%c_Xh, this%fluid%gs_Xh, &
267 json_subdict, numerics_params, this%user, this%chkp, this%fluid%ulag, &
268 this%fluid%vlag, this%fluid%wlag, this%fluid%ext_bdf, &
269 this%fluid%rho)
270 end if
271 end if
272
273 !
274 ! Setup initial conditions
275 !
276 call json_get(this%params, 'case.fluid.initial_condition.type', &
277 string_val)
278 call json_get(this%params, 'case.fluid.initial_condition', &
279 json_subdict)
280
281 call neko_log%section("Fluid initial condition ")
282
283 if (this%params%valid_path('case.restart_file')) then
284 call neko_log%message("Restart file specified, " // &
285 "initial conditions ignored")
286 else if (trim(string_val) .ne. 'user') then
287 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, &
288 this%fluid%p, this%fluid%c_Xh, this%fluid%gs_Xh, string_val, &
289 json_subdict)
290 else
291 call json_get(this%params, 'case.fluid.scheme', string_val)
292 if (trim(string_val) .eq. 'compressible') then
293 call set_flow_ic(this%fluid%rho, &
294 this%fluid%u, this%fluid%v, this%fluid%w, this%fluid%p, &
295 this%fluid%c_Xh, this%fluid%gs_Xh, &
296 this%user%initial_conditions, this%fluid%name)
297 else
298 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, &
299 this%fluid%p, this%fluid%c_Xh, this%fluid%gs_Xh, &
300 this%user%initial_conditions, this%fluid%name)
301 end if
302 end if
303
304 call neko_log%end_section()
305
306 if (scalar) then
307 call neko_log%section("Scalar initial condition ")
308
309 if (this%params%valid_path('case.restart_file')) then
310 call neko_log%message("Restart file specified, " // &
311 "initial conditions ignored")
312 else if (this%params%valid_path('case.scalar')) then
313 ! For backward compatibility with single scalar
314 call json_get(this%params, 'case.scalar.initial_condition.type', &
315 string_val)
316 call json_get(this%params, &
317 'case.scalar.initial_condition', json_subdict)
318
319 if (trim(string_val) .ne. 'user') then
320 if (trim(this%scalars%scalar_fields(1)%name) .eq. 'temperature') then
321 call set_scalar_ic(this%scalars%scalar_fields(1)%s, &
322 this%scalars%scalar_fields(1)%c_Xh, &
323 this%scalars%scalar_fields(1)%gs_Xh, &
324 string_val, json_subdict, 0)
325 else
326 call set_scalar_ic(this%scalars%scalar_fields(1)%s, &
327 this%scalars%scalar_fields(1)%c_Xh, &
328 this%scalars%scalar_fields(1)%gs_Xh, &
329 string_val, json_subdict, 1)
330 end if
331 else
332 call set_scalar_ic(this%scalars%scalar_fields(1)%name, &
333 this%scalars%scalar_fields(1)%s, &
334 this%scalars%scalar_fields(1)%c_Xh, &
335 this%scalars%scalar_fields(1)%gs_Xh, &
336 this%user%initial_conditions)
337 end if
338
339 else
340 ! Handle multiple scalars
341 do i = 1, n_scalars
342 call json_extract_item(this%params, 'case.scalars', i, &
343 scalar_params)
344 call json_get(scalar_params, 'initial_condition.type', string_val)
345 call json_get(scalar_params, 'initial_condition', &
346 json_subdict)
347
348 if (trim(string_val) .ne. 'user') then
349 if (trim(this%scalars%scalar_fields(i)%name) .eq. 'temperature') then
350 call set_scalar_ic(this%scalars%scalar_fields(i)%s, &
351 this%scalars%scalar_fields(i)%c_Xh, &
352 this%scalars%scalar_fields(i)%gs_Xh, &
353 string_val, json_subdict, 0)
354 temperature_found = .true.
355 else
356 if (temperature_found) then
357 ! if temperature is found, other scalars start from index 1
358 call set_scalar_ic(this%scalars%scalar_fields(i)%s, &
359 this%scalars%scalar_fields(i)%c_Xh, &
360 this%scalars%scalar_fields(i)%gs_Xh, &
361 string_val, json_subdict, i - 1)
362 else
363 ! if temperature is not found, other scalars start from index 0
364 call set_scalar_ic(this%scalars%scalar_fields(i)%s, &
365 this%scalars%scalar_fields(i)%c_Xh, &
366 this%scalars%scalar_fields(i)%gs_Xh, &
367 string_val, json_subdict, i)
368 end if
369 end if
370 else
371 call set_scalar_ic(this%scalars%scalar_fields(i)%name,&
372 this%scalars%scalar_fields(i)%s, &
373 this%scalars%scalar_fields(i)%c_Xh, &
374 this%scalars%scalar_fields(i)%gs_Xh, &
375 this%user%initial_conditions)
376 end if
377 end do
378 end if
379
380 call neko_log%end_section()
381 end if
382
383 ! Add initial conditions to BDF scheme (if present)
384 select type (f => this%fluid)
385 type is (fluid_pnpn_t)
386 call f%ulag%set(f%u)
387 call f%vlag%set(f%v)
388 call f%wlag%set(f%w)
389 end select
390
391 !
392 ! Validate that the case is properly setup for time-stepping
393 !
394 call this%fluid%validate
395
396 if (scalar) then
397 call this%scalars%validate()
398 end if
399
400 !
401 ! Get and process output directory
402 !
403 call json_get_or_default(this%params, 'case.output_directory',&
404 this%output_directory, '')
405
406 output_dir_len = len(trim(this%output_directory))
407 if (output_dir_len .gt. 0) then
408 if (this%output_directory(output_dir_len:output_dir_len) .ne. "/") then
409 this%output_directory = trim(this%output_directory)//"/"
410 if (pe_rank .eq. 0) then
411 call execute_command_line('mkdir -p '//this%output_directory)
412 end if
413 end if
414 end if
415
416 !
417 ! Save mesh partitions (if requested)
418 !
419 call json_get_or_default(this%params, 'case.output_partitions',&
420 logical_val, .false.)
421 if (logical_val) then
422 call msh_part%init(this%msh, 'MPI_Rank')
423 msh_part%data = pe_rank
424 call part_file%init(trim(this%output_directory)//'partitions.vtk')
425 call part_file%write(msh_part)
426 call msh_part%free()
427 end if
428
429 !
430 ! Setup output precision of the field files
431 !
432 call json_get_or_default(this%params, 'case.output_precision', string_val,&
433 'single')
434
435 if (trim(string_val) .eq. 'double') then
436 precision = dp
437 else
438 precision = sp
439 end if
440
441 !
442 ! Setup output layout of the field bp file
443 !
444 call json_get_or_default(this%params, 'case.output_layout', layout, 1)
445
446 !
447 ! Setup output_controller
448 !
449 call json_get_or_default(this%params, 'case.fluid.output_filename', &
450 name, "field")
451 call json_get_or_default(this%params, 'case.fluid.output_format', &
452 file_format, 'fld')
453 call this%output_controller%init(this%time%end_time)
454 if (scalar) then
455 call this%f_out%init(precision, this%fluid, this%scalars, name = name, &
456 path = trim(this%output_directory), &
457 fmt = trim(file_format), layout = layout)
458 else
459 call this%f_out%init(precision, this%fluid, name = name, &
460 path = trim(this%output_directory), &
461 fmt = trim(file_format), layout = layout)
462 end if
463
464 call json_get_or_default(this%params, 'case.fluid.output_control',&
465 string_val, 'org')
466
467 if (trim(string_val) .eq. 'org') then
468 ! yes, it should be real_val below for type compatibility
469 call json_get(this%params, 'case.nsamples', real_val)
470 call this%output_controller%add(this%f_out, real_val, 'nsamples')
471 else if (trim(string_val) .eq. 'never') then
472 ! Fix a dummy 0.0 output_value
473 call json_get_or_default(this%params, 'case.fluid.output_value', &
474 real_val, 0.0_rp)
475 call this%output_controller%add(this%f_out, 0.0_rp, string_val)
476 else
477 call json_get(this%params, 'case.fluid.output_value', real_val)
478 call this%output_controller%add(this%f_out, real_val, string_val)
479 end if
480
481 !
482 ! Save checkpoints (if nothing specified, default to saving at end of sim)
483 !
484 call json_get_or_default(this%params, 'case.output_checkpoints',&
485 logical_val, .true.)
486 if (logical_val) then
487 call json_get_or_default(this%params, 'case.checkpoint_filename', &
488 name, "fluid")
489 call json_get_or_default(this%params, 'case.checkpoint_format', &
490 string_val, "chkp")
491 call this%chkp_out%init(this%chkp, name = name,&
492 path = this%output_directory, fmt = trim(string_val))
493 call json_get_or_default(this%params, 'case.checkpoint_control', &
494 string_val, "simulationtime")
495 call json_get_or_default(this%params, 'case.checkpoint_value', &
496 real_val, 1e10_rp)
497 call this%output_controller%add(this%chkp_out, real_val, string_val, &
498 neko_eps)
499 end if
500
501 !
502 ! Setup joblimit
503 !
504 if (this%params%valid_path('case.job_timelimit')) then
505 call json_get(this%params, 'case.job_timelimit', string_val)
506 call jobctrl_set_time_limit(string_val)
507 end if
508
509 call neko_log%end_section()
510
511 end subroutine case_init_common
512
514 subroutine case_free(this)
515 class(case_t), intent(inout) :: this
516
517 if (allocated(this%fluid)) then
518 call this%fluid%free()
519 deallocate(this%fluid)
520 end if
521
522 if (allocated(this%scalars)) then
523 call this%scalars%free()
524 deallocate(this%scalars)
525 end if
526
527 call this%msh%free()
528
529 call this%f_out%free()
530
531 call this%output_controller%free()
532
533 end subroutine case_free
534
535end module case
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Defines a simulation case.
Definition case.f90:34
subroutine case_init_from_file(this, case_file)
Initialize a case from an input file case_file.
Definition case.f90:93
subroutine case_free(this)
Deallocate a case.
Definition case.f90:515
subroutine case_init_from_json(this, case_json)
Initialize a case from a JSON object describing a case.
Definition case.f90:128
subroutine case_init_common(this)
Initialize a case from its (loaded) params object.
Definition case.f90:142
Defines a checkpoint.
Defines an output for a checkpoint.
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Module for file I/O operations.
Definition file.f90:34
Initial flow condition.
Definition flow_ic.f90:34
Defines an output for a fluid.
Modular version of the Classic Nek5000 Pn/Pn formulation for fluids.
Job control.
Definition jobctrl.f90:34
Utilities for retrieving parameters from the case files.
logical, public json_no_defaults
If true, the json_get_or_default routines will not add missing parameters.
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_quiet
Always logged.
Definition log.f90:78
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
Definition math.f90:60
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
Defines a mesh field.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements output_controller_t
Interface to ParMETIS.
Definition parmetis.F90:34
subroutine, public parmetis_partmeshkway(msh, parts, weights, nprts)
Compute a k-way partitioning of a mesh msh.
Definition parmetis.F90:112
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Redistribution routines.
Definition redist.f90:34
subroutine, public redist_mesh(msh, parts)
Redistribute a mesh msh according to new partitions.
Definition redist.f90:59
Scalar initial condition.
Definition scalar_ic.f90:34
Contains the scalar_pnpn_t type.
Contains the scalar_scheme_t type.
Contains the scalars_t type that manages multiple scalar fields.
Definition scalars.f90:35
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Compound scheme for the advection and diffusion operators in a transport equation.
Module with things related to the simulation time.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:55
Base type of all fluid formulations.
Centralized controller for a list of outputs.
Base type for a scalar advection-diffusion solver.
Type to manage multiple scalar transport equations.
Definition scalars.f90:61
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
A struct that contains all info about the time, expand as needed.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...