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