Neko 1.99.1
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 end type case_t
83
84 interface case_init
86 end interface case_init
87
88 public :: case_init, case_free
89
90contains
91
93 subroutine case_init_from_file(this, case_file)
94 type(case_t), target, intent(inout) :: this
95 character(len=*), intent(in) :: case_file
96 integer :: ierr, integer_val
97 character(len=:), allocatable :: json_buffer
98 logical :: exist
99
100 ! Check if the file exists
101 inquire(file = trim(case_file), exist = exist)
102 if (.not. exist) then
103 call neko_error('The case file '//trim(case_file)//' does not exist.')
104 end if
105
106 call neko_log%section('Case')
107 call neko_log%message('Reading case file ' // trim(case_file), &
109
110 if (pe_rank .eq. 0) then
111 call this%params%load_file(filename = trim(case_file))
112 call this%params%print_to_string(json_buffer)
113 integer_val = len(json_buffer)
114 end if
115
116 call mpi_bcast(integer_val, 1, mpi_integer, 0, neko_comm, ierr)
117 if (pe_rank .ne. 0) allocate(character(len = integer_val) :: json_buffer)
118 call mpi_bcast(json_buffer, integer_val, mpi_character, 0, neko_comm, ierr)
119 call this%params%load_from_string(json_buffer)
120
121 deallocate(json_buffer)
122
123 call case_init_common(this)
124
125 end subroutine case_init_from_file
126
128 subroutine case_init_from_json(this, case_json)
129 type(case_t), target, intent(inout) :: this
130 type(json_file), intent(in) :: case_json
131
132 call neko_log%section('Case')
133 call neko_log%message('Creating case from JSON object', neko_log_quiet)
134
135 this%params = case_json
136
137 call case_init_common(this)
138
139 end subroutine case_init_from_json
140
142 subroutine case_init_common(this)
143 type(case_t), target, intent(inout) :: this
144 integer :: lx = 0
145 logical :: scalar = .false.
146 type(file_t) :: msh_file, bdry_file, part_file
147 type(mesh_fld_t) :: msh_part, parts
148 logical :: found, logical_val
149 logical :: temperature_found = .false.
150 integer :: integer_val
151 real(kind=rp) :: real_val
152 character(len = :), allocatable :: string_val, name, file_format
153 integer :: output_dir_len
154 integer :: precision, layout
155 type(json_file) :: scalar_params, numerics_params
156 type(json_file) :: json_subdict
157 integer :: n_scalars, i
158
159 !
160 ! Setup user defined functions
161 !
162 call this%user%init()
163
164 ! Run user startup routine
165 call this%user%startup(this%params)
166
167 ! Check if default value fill-in is allowed
168 if (this%params%valid_path('case.no_defaults')) then
169 call json_get(this%params, 'case.no_defaults', json_no_defaults)
170 end if
171
172 !
173 ! Load mesh
174 !
175 call json_get_or_default(this%params, 'case.mesh_file', string_val, &
176 'no mesh')
177 if (trim(string_val) .eq. 'no mesh') then
178 call neko_error('The mesh_file keyword could not be found in the .' // &
179 'case file. Often caused by incorrectly formatted json.')
180 end if
181 call msh_file%init(string_val)
182
183 call msh_file%read(this%msh)
184
185 !
186 ! Load Balancing
187 !
188 call json_get_or_default(this%params, 'case.load_balancing', logical_val,&
189 .false.)
190
191 if (pe_size .gt. 1 .and. logical_val) then
192 call neko_log%section('Load Balancing')
193 call parmetis_partmeshkway(this%msh, parts)
194 call redist_mesh(this%msh, parts)
195
196 ! store the balanced mesh (for e.g. restarts)
197 string_val = trim(string_val(1:scan(trim(string_val), &
198 '.', back = .true.) - 1)) // '_lb.nmsh'
199 call msh_file%init(string_val)
200 call msh_file%write(this%msh)
201
202 call neko_log%end_section()
203 end if
204
205 !
206 ! Time control
207 !
208 call json_get(this%params, 'case.time', json_subdict)
209 call this%time%init(json_subdict)
210
211 !
212 ! Initialize point_zones registry
213 !
214 call neko_point_zone_registry%init(this%params, this%msh)
215
216 ! Run user mesh motion routine
217 call this%user%mesh_setup(this%msh, this%time)
218
219 call json_get(this%params, 'case.numerics', numerics_params)
220
221 !
222 ! Setup fluid scheme
223 !
224 call json_get(this%params, 'case.fluid.scheme', string_val)
225 call fluid_scheme_base_factory(this%fluid, trim(string_val))
226
227 call json_get(this%params, 'case.numerics.polynomial_order', lx)
228 lx = lx + 1 ! add 1 to get number of gll points
229 ! Set time lags in chkp
230 this%chkp%tlag => this%time%tlag
231 this%chkp%dtlag => this%time%dtlag
232 call this%fluid%init(this%msh, lx, this%params, this%user, this%chkp)
233
234
235 !
236 ! Setup scratch registry
237 !
238 call neko_scratch_registry%init(this%fluid%dm_Xh, 10, 10)
239
240 !
241 ! Setup scalar scheme
242 !
243 ! @todo no scalar factory for now, probably not needed
244 scalar = .false.
245 n_scalars = 0
246 if (this%params%valid_path('case.scalar')) then
247 call json_get_or_default(this%params, 'case.scalar.enabled', scalar, &
248 .true.)
249 n_scalars = 1
250 else if (this%params%valid_path('case.scalars')) then
251 call this%params%info('case.scalars', n_children = n_scalars)
252 if (n_scalars > 0) then
253 scalar = .true.
254 end if
255 end if
256
257 if (scalar) then
258 allocate(this%scalars)
259 if (this%params%valid_path('case.scalar')) then
260 ! For backward compatibility
261 call json_get(this%params, 'case.scalar', scalar_params)
262 call this%scalars%init(this%msh, this%fluid%c_Xh, this%fluid%gs_Xh, &
263 scalar_params, numerics_params, this%user, this%chkp, this%fluid%ulag, &
264 this%fluid%vlag, this%fluid%wlag, this%fluid%ext_bdf, &
265 this%fluid%rho)
266 else
267 ! Multiple scalars
268 call json_get(this%params, 'case.scalars', json_subdict)
269 call this%scalars%init(n_scalars, this%msh, this%fluid%c_Xh, this%fluid%gs_Xh, &
270 json_subdict, numerics_params, this%user, this%chkp, this%fluid%ulag, &
271 this%fluid%vlag, this%fluid%wlag, this%fluid%ext_bdf, &
272 this%fluid%rho)
273 end if
274 end if
275
276 !
277 ! Setup initial conditions
278 !
279 call json_get(this%params, 'case.fluid.initial_condition.type', &
280 string_val)
281 call json_get(this%params, 'case.fluid.initial_condition', &
282 json_subdict)
283
284 call neko_log%section("Fluid initial condition ")
285
286 if (this%params%valid_path('case.restart_file')) then
287 call neko_log%message("Restart file specified, " // &
288 "initial conditions ignored")
289 else if (trim(string_val) .ne. 'user') then
290 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, &
291 this%fluid%p, this%fluid%c_Xh, this%fluid%gs_Xh, string_val, &
292 json_subdict)
293 else
294 call json_get(this%params, 'case.fluid.scheme', string_val)
295 if (trim(string_val) .eq. 'compressible') then
296 call set_flow_ic(this%fluid%rho, &
297 this%fluid%u, this%fluid%v, this%fluid%w, this%fluid%p, &
298 this%fluid%c_Xh, this%fluid%gs_Xh, &
299 this%user%initial_conditions, this%fluid%name)
300 else
301 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, &
302 this%fluid%p, this%fluid%c_Xh, this%fluid%gs_Xh, &
303 this%user%initial_conditions, this%fluid%name)
304 end if
305 end if
306
307 call neko_log%end_section()
308
309 if (scalar) then
310 call neko_log%section("Scalar initial condition ")
311
312 if (this%params%valid_path('case.restart_file')) then
313 call neko_log%message("Restart file specified, " // &
314 "initial conditions ignored")
315 else if (this%params%valid_path('case.scalar')) then
316 ! For backward compatibility with single scalar
317 call json_get(this%params, 'case.scalar.initial_condition.type', &
318 string_val)
319 call json_get(this%params, &
320 'case.scalar.initial_condition', json_subdict)
321
322 if (trim(string_val) .ne. 'user') then
323 if (trim(this%scalars%scalar_fields(1)%name) .eq. 'temperature') then
324 call set_scalar_ic(this%scalars%scalar_fields(1)%s, &
325 this%scalars%scalar_fields(1)%c_Xh, &
326 this%scalars%scalar_fields(1)%gs_Xh, &
327 string_val, json_subdict, 0)
328 else
329 call set_scalar_ic(this%scalars%scalar_fields(1)%s, &
330 this%scalars%scalar_fields(1)%c_Xh, &
331 this%scalars%scalar_fields(1)%gs_Xh, &
332 string_val, json_subdict, 1)
333 end if
334 else
335 call set_scalar_ic(this%scalars%scalar_fields(1)%name, &
336 this%scalars%scalar_fields(1)%s, &
337 this%scalars%scalar_fields(1)%c_Xh, &
338 this%scalars%scalar_fields(1)%gs_Xh, &
339 this%user%initial_conditions)
340 end if
341
342 else
343 ! Handle multiple scalars
344 do i = 1, n_scalars
345 call json_extract_item(this%params, 'case.scalars', i, &
346 scalar_params)
347 call json_get(scalar_params, 'initial_condition.type', string_val)
348 call json_get(scalar_params, 'initial_condition', &
349 json_subdict)
350
351 if (trim(string_val) .ne. 'user') then
352 if (trim(this%scalars%scalar_fields(i)%name) .eq. 'temperature') then
353 call set_scalar_ic(this%scalars%scalar_fields(i)%s, &
354 this%scalars%scalar_fields(i)%c_Xh, &
355 this%scalars%scalar_fields(i)%gs_Xh, &
356 string_val, json_subdict, 0)
357 temperature_found = .true.
358 else
359 if (temperature_found) then
360 ! if temperature is found, other scalars start from index 1
361 call set_scalar_ic(this%scalars%scalar_fields(i)%s, &
362 this%scalars%scalar_fields(i)%c_Xh, &
363 this%scalars%scalar_fields(i)%gs_Xh, &
364 string_val, json_subdict, i - 1)
365 else
366 ! if temperature is not found, other scalars start from index 0
367 call set_scalar_ic(this%scalars%scalar_fields(i)%s, &
368 this%scalars%scalar_fields(i)%c_Xh, &
369 this%scalars%scalar_fields(i)%gs_Xh, &
370 string_val, json_subdict, i)
371 end if
372 end if
373 else
374 call set_scalar_ic(this%scalars%scalar_fields(i)%name,&
375 this%scalars%scalar_fields(i)%s, &
376 this%scalars%scalar_fields(i)%c_Xh, &
377 this%scalars%scalar_fields(i)%gs_Xh, &
378 this%user%initial_conditions)
379 end if
380 end do
381 end if
382
383 call neko_log%end_section()
384 end if
385
386 ! Add initial conditions to BDF scheme (if present)
387 select type (f => this%fluid)
388 type is (fluid_pnpn_t)
389 call f%ulag%set(f%u)
390 call f%vlag%set(f%v)
391 call f%wlag%set(f%w)
392 end select
393
394 !
395 ! Validate that the case is properly setup for time-stepping
396 !
397 call this%fluid%validate
398
399 if (scalar) then
400 call this%scalars%validate()
401 end if
402
403 !
404 ! Get and process output directory
405 !
406 call json_get_or_default(this%params, 'case.output_directory',&
407 this%output_directory, '')
408
409 output_dir_len = len(trim(this%output_directory))
410 if (output_dir_len .gt. 0) then
411 if (this%output_directory(output_dir_len:output_dir_len) .ne. "/") then
412 this%output_directory = trim(this%output_directory)//"/"
413 if (pe_rank .eq. 0) then
414 call execute_command_line('mkdir -p '//this%output_directory)
415 end if
416 end if
417 end if
418
419 !
420 ! Save mesh partitions (if requested)
421 !
422 call json_get_or_default(this%params, 'case.output_partitions',&
423 logical_val, .false.)
424 if (logical_val) then
425 call msh_part%init(this%msh, 'MPI_Rank')
426 msh_part%data = pe_rank
427 call part_file%init(trim(this%output_directory)//'partitions.vtk')
428 call part_file%write(msh_part)
429 call msh_part%free()
430 end if
431
432 !
433 ! Setup output precision of the field files
434 !
435 call json_get_or_default(this%params, 'case.output_precision', string_val,&
436 'single')
437
438 if (trim(string_val) .eq. 'double') then
439 precision = dp
440 else
441 precision = sp
442 end if
443
444 !
445 ! Setup output layout of the field bp file
446 !
447 call json_get_or_default(this%params, 'case.output_layout', layout, 1)
448
449 !
450 ! Setup output_controller
451 !
452 call json_get_or_default(this%params, 'case.fluid.output_filename', &
453 name, "field")
454 call json_get_or_default(this%params, 'case.fluid.output_format', &
455 file_format, 'fld')
456 call this%output_controller%init(this%time%end_time)
457 if (scalar) then
458 call this%f_out%init(precision, this%fluid, this%scalars, name = name, &
459 path = trim(this%output_directory), &
460 fmt = trim(file_format), layout = layout)
461 else
462 call this%f_out%init(precision, this%fluid, name = name, &
463 path = trim(this%output_directory), &
464 fmt = trim(file_format), layout = layout)
465 end if
466
467 call json_get_or_default(this%params, 'case.fluid.output_control',&
468 string_val, 'org')
469
470 if (trim(string_val) .eq. 'org') then
471 ! yes, it should be real_val below for type compatibility
472 call json_get(this%params, 'case.nsamples', real_val)
473 call this%output_controller%add(this%f_out, real_val, 'nsamples')
474 else if (trim(string_val) .eq. 'never') then
475 ! Fix a dummy 0.0 output_value
476 call json_get_or_default(this%params, 'case.fluid.output_value', &
477 real_val, 0.0_rp)
478 call this%output_controller%add(this%f_out, 0.0_rp, string_val)
479 else
480 call json_get(this%params, 'case.fluid.output_value', real_val)
481 call this%output_controller%add(this%f_out, real_val, string_val)
482 end if
483
484 !
485 ! Save checkpoints (if nothing specified, default to saving at end of sim)
486 !
487 call json_get_or_default(this%params, 'case.output_checkpoints',&
488 logical_val, .true.)
489 if (logical_val) then
490 call json_get_or_default(this%params, 'case.checkpoint_filename', &
491 name, "fluid")
492 call json_get_or_default(this%params, 'case.checkpoint_format', &
493 string_val, "chkp")
494 call this%chkp_out%init(this%chkp, name = name,&
495 path = this%output_directory, fmt = trim(string_val))
496 call json_get_or_default(this%params, 'case.checkpoint_control', &
497 string_val, "simulationtime")
498 call json_get_or_default(this%params, 'case.checkpoint_value', &
499 real_val, 1e10_rp)
500 call this%output_controller%add(this%chkp_out, real_val, string_val, &
501 neko_eps)
502 end if
503
504 !
505 ! Setup joblimit
506 !
507 if (this%params%valid_path('case.job_timelimit')) then
508 call json_get(this%params, 'case.job_timelimit', string_val)
509 call jobctrl_set_time_limit(string_val)
510 end if
511
512 call neko_log%end_section()
513
514 end subroutine case_init_common
515
517 subroutine case_free(this)
518 type(case_t), intent(inout) :: this
519
520 if (allocated(this%fluid)) then
521 call this%fluid%free()
522 deallocate(this%fluid)
523 end if
524
525 if (allocated(this%scalars)) then
526 call this%scalars%free()
527 deallocate(this%scalars)
528 end if
529
530 call this%msh%free()
531
532 call this%f_out%free()
533
534 call this%output_controller%free()
535
536 end subroutine case_free
537
538end 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:94
subroutine, public case_free(this)
Deallocate a case.
Definition case.f90:518
subroutine case_init_from_json(this, case_json)
Initialize a case from a JSON object describing a case.
Definition case.f90:129
subroutine case_init_common(this)
Initialize a case from its (loaded) params object.
Definition case.f90:143
Defines a checkpoint.
Defines an output for a checkpoint.
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:58
integer, public pe_rank
MPI rank.
Definition comm.F90:55
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
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:72
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
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 fields This can be used when you have a funct...
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:62
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...