Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
case.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2023, 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
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 comm
52 use checkpoint, only: chkp_t
54 use logger, only : neko_log, neko_log_quiet
56 use user_intf, only : user_t
57 use scalar_pnpn, only : scalar_pnpn_t
58 use time_state, only : time_state_t
59 use json_module, only : json_file
63 implicit none
64 private
65 type, public :: case_t
66 type(mesh_t) :: msh
67 type(json_file) :: params
68 character(len=:), allocatable :: output_directory
70 type(fluid_output_t) :: f_out
71 type(time_state_t) :: time
72 type(chkp_output_t) :: chkp_out
73 type(chkp_t) :: chkp
74 type(user_t) :: usr
75 class(fluid_scheme_base_t), allocatable :: fluid
76 type(scalar_pnpn_t), allocatable :: scalar
77 end type case_t
78
79 interface case_init
81 end interface case_init
82
83 public :: case_init, case_free
84
85contains
86
88 subroutine case_init_from_file(this, case_file)
89 type(case_t), target, intent(inout) :: this
90 character(len=*), intent(in) :: case_file
91 integer :: ierr, integer_val
92 character(len=:), allocatable :: json_buffer
93 logical :: exist
94
95 ! Check if the file exists
96 inquire(file = trim(case_file), exist = exist)
97 if (.not. exist) then
98 call neko_error('The case file '//trim(case_file)//' does not exist.')
99 end if
100
101 call neko_log%section('Case')
102 call neko_log%message('Reading case file ' // trim(case_file), &
104
105 if (pe_rank .eq. 0) then
106 call this%params%load_file(filename = trim(case_file))
107 call this%params%print_to_string(json_buffer)
108 integer_val = len(json_buffer)
109 end if
110
111 call mpi_bcast(integer_val, 1, mpi_integer, 0, neko_comm, ierr)
112 if (pe_rank .ne. 0) allocate(character(len = integer_val) :: json_buffer)
113 call mpi_bcast(json_buffer, integer_val, mpi_character, 0, neko_comm, ierr)
114 call this%params%load_from_string(json_buffer)
115
116 deallocate(json_buffer)
117
118 call case_init_common(this)
119
120 end subroutine case_init_from_file
121
123 subroutine case_init_from_json(this, case_json)
124 type(case_t), target, intent(inout) :: this
125 type(json_file), intent(in) :: case_json
126
127 call neko_log%section('Case')
128 call neko_log%message('Creating case from JSON object', neko_log_quiet)
129
130 this%params = case_json
131
132 call case_init_common(this)
133
134 end subroutine case_init_from_json
135
137 subroutine case_init_common(this)
138 type(case_t), target, intent(inout) :: this
139 integer :: lx = 0
140 logical :: scalar = .false.
141 type(file_t) :: msh_file, bdry_file, part_file
142 type(mesh_fld_t) :: msh_part, parts
143 logical :: found, logical_val
144 integer :: integer_val
145 real(kind=rp) :: real_val
146 character(len = :), allocatable :: string_val, name
147 integer :: output_dir_len
148 integer :: precision
149 type(json_file) :: scalar_params, numerics_params
150 type(json_file) :: json_subdict
151
152 !
153 ! Setup user defined functions
154 !
155 call this%usr%init()
156
157 ! Run user startup routine
158 call this%usr%user_startup(this%params)
159
160 !
161 ! Load mesh
162 !
163 call json_get_or_default(this%params, 'case.mesh_file', string_val, &
164 'no mesh')
165 if (trim(string_val) .eq. 'no mesh') then
166 call neko_error('The mesh_file keyword could not be found in the .' // &
167 'case file. Often caused by incorrectly formatted json.')
168 end if
169 msh_file = file_t(string_val)
170
171 call msh_file%read(this%msh)
172
173 !
174 ! Load Balancing
175 !
176 call json_get_or_default(this%params, 'case.load_balancing', logical_val,&
177 .false.)
178
179 if (pe_size .gt. 1 .and. logical_val) then
180 call neko_log%section('Load Balancing')
181 call parmetis_partmeshkway(this%msh, parts)
182 call redist_mesh(this%msh, parts)
183
184 ! store the balanced mesh (for e.g. restarts)
185 string_val = trim(string_val(1:scan(trim(string_val), &
186 '.', back = .true.) - 1)) // '_lb.nmsh'
187 msh_file = file_t(string_val)
188 call msh_file%write(this%msh)
189
190 call neko_log%end_section()
191 end if
192
193 !
194 ! Time step
195 !
196 call json_get_or_default(this%params, 'case.variable_timestep', &
197 logical_val, .false.)
198 if (.not. logical_val) then
199 call json_get(this%params, 'case.timestep', this%time%dt)
200 else
201 ! randomly set an initial dt to get cfl when dt is variable
202 this%time%dt = 1.0_rp
203 end if
204
205 !
206 ! End time
207 !
208 call json_get(this%params, 'case.end_time', this%time%end_time)
209
210 !
211 ! Initialize point_zones registry
212 !
213 call neko_point_zone_registry%init(this%params, this%msh)
214
215 ! Run user mesh motion routine
216 call this%usr%user_mesh_setup(this%msh)
217
218 call json_extract_object(this%params, 'case.numerics', numerics_params)
219
220 !
221 ! Setup fluid scheme
222 !
223 call json_get(this%params, 'case.fluid.scheme', string_val)
224 call fluid_scheme_base_factory(this%fluid, trim(string_val))
225
226 call json_get(this%params, 'case.numerics.polynomial_order', lx)
227 lx = lx + 1 ! add 1 to get number of gll points
228 ! Set time lags in chkp
229 this%chkp%tlag => this%time%tlag
230 this%chkp%dtlag => this%time%dtlag
231 call this%fluid%init(this%msh, lx, this%params, this%usr, this%chkp)
232
233
234 !
235 ! Setup scratch registry
236 !
237 neko_scratch_registry = scratch_registry_t(this%fluid%dm_Xh, 10, 10)
238
239 !
240 ! Setup scalar scheme
241 !
242 ! @todo no scalar factory for now, probably not needed
243 if (this%params%valid_path('case.scalar')) then
244 call json_get_or_default(this%params, 'case.scalar.enabled', scalar,&
245 .true.)
246 end if
247
248 if (scalar) then
249 allocate(this%scalar)
250 call json_extract_object(this%params, 'case.scalar', scalar_params)
251 call this%scalar%init(this%msh, this%fluid%c_Xh, this%fluid%gs_Xh, &
252 scalar_params, numerics_params, this%usr, this%chkp, this%fluid%ulag, &
253 this%fluid%vlag, this%fluid%wlag, this%fluid%ext_bdf, &
254 this%fluid%rho)
255
256 end if
257
258 !
259 ! Setup initial conditions
260 !
261 call json_get(this%params, 'case.fluid.initial_condition.type', &
262 string_val)
263 call json_extract_object(this%params, 'case.fluid.initial_condition', &
264 json_subdict)
265
266 call neko_log%section("Fluid initial condition ")
267
268 if (trim(string_val) .ne. 'user') then
269 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, &
270 this%fluid%p, this%fluid%c_Xh, this%fluid%gs_Xh, string_val, &
271 json_subdict)
272 else
273 call json_get(this%params, 'case.fluid.scheme', string_val)
274 if (trim(string_val) .eq. 'compressible') then
275 call set_flow_ic(this%fluid%rho_field, &
276 this%fluid%u, this%fluid%v, this%fluid%w, this%fluid%p, &
277 this%fluid%c_Xh, this%fluid%gs_Xh, this%usr%fluid_compressible_user_ic, &
278 this%params)
279 else
280 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, this%fluid%p,&
281 this%fluid%c_Xh, this%fluid%gs_Xh, this%usr%fluid_user_ic, &
282 this%params)
283 end if
284 end if
285
286 call neko_log%end_section()
287
288 if (scalar) then
289
290 call json_get(this%params, 'case.scalar.initial_condition.type', &
291 string_val)
292 call json_extract_object(this%params, 'case.scalar.initial_condition', &
293 json_subdict)
294
295 call neko_log%section("Scalar initial condition ")
296
297 if (trim(string_val) .ne. 'user') then
298 call set_scalar_ic(this%scalar%s, &
299 this%scalar%c_Xh, this%scalar%gs_Xh, string_val, json_subdict)
300 else
301 call set_scalar_ic(this%scalar%s, &
302 this%scalar%c_Xh, this%scalar%gs_Xh, this%usr%scalar_user_ic, &
303 this%params)
304 end if
305
306 call neko_log%end_section()
307
308 end if
309
310 ! Add initial conditions to BDF scheme (if present)
311 select type (f => this%fluid)
312 type is (fluid_pnpn_t)
313 call f%ulag%set(f%u)
314 call f%vlag%set(f%v)
315 call f%wlag%set(f%w)
316 end select
317
318 !
319 ! Validate that the case is properly setup for time-stepping
320 !
321 call this%fluid%validate
322
323 if (scalar) then
324 call this%scalar%slag%set(this%scalar%s)
325 call this%scalar%validate
326 end if
327
328 !
329 ! Get and process output directory
330 !
331 call json_get_or_default(this%params, 'case.output_directory',&
332 this%output_directory, '')
333
334 output_dir_len = len(trim(this%output_directory))
335 if (output_dir_len .gt. 0) then
336 if (this%output_directory(output_dir_len:output_dir_len) .ne. "/") then
337 this%output_directory = trim(this%output_directory)//"/"
338 if (pe_rank .eq. 0) then
339 call execute_command_line('mkdir -p '//this%output_directory)
340 end if
341 end if
342 end if
343
344 !
345 ! Save mesh partitions (if requested)
346 !
347 call json_get_or_default(this%params, 'case.output_partitions',&
348 logical_val, .false.)
349 if (logical_val) then
350 call mesh_field_init(msh_part, this%msh, 'MPI_Rank')
351 msh_part%data = pe_rank
352 part_file = file_t(trim(this%output_directory)//'partitions.vtk')
353 call part_file%write(msh_part)
354 call mesh_field_free(msh_part)
355 end if
356
357 !
358 ! Setup output precision of the field files
359 !
360 call json_get_or_default(this%params, 'case.output_precision', string_val,&
361 'single')
362
363 if (trim(string_val) .eq. 'double') then
364 precision = dp
365 else
366 precision = sp
367 end if
368
369 !
370 ! Setup output_controller
371 !
372 call this%output_controller%init(this%time%end_time)
373 call json_get_or_default(this%params, 'case.fluid.output_filename', &
374 name, "field")
375 if (scalar) then
376 call this%f_out%init(precision, this%fluid, this%scalar, name = name, &
377 path = trim(this%output_directory))
378 else
379 call this%f_out%init(precision, this%fluid, name = name, &
380 path = trim(this%output_directory))
381 end if
382
383 call json_get_or_default(this%params, 'case.fluid.output_control',&
384 string_val, 'org')
385
386 if (trim(string_val) .eq. 'org') then
387 ! yes, it should be real_val below for type compatibility
388 call json_get(this%params, 'case.nsamples', real_val)
389 call this%output_controller%add(this%f_out, real_val, 'nsamples')
390 else if (trim(string_val) .eq. 'never') then
391 ! Fix a dummy 0.0 output_value
392 call json_get_or_default(this%params, 'case.fluid.output_value', &
393 real_val, 0.0_rp)
394 call this%output_controller%add(this%f_out, 0.0_rp, string_val)
395 else
396 call json_get(this%params, 'case.fluid.output_value', real_val)
397 call this%output_controller%add(this%f_out, real_val, string_val)
398 end if
399
400 !
401 ! Save checkpoints (if nothing specified, default to saving at end of sim)
402 !
403 call json_get_or_default(this%params, 'case.output_checkpoints',&
404 logical_val, .true.)
405 if (logical_val) then
406 call json_get_or_default(this%params, 'case.checkpoint_filename', &
407 name, "fluid")
408 call json_get_or_default(this%params, 'case.checkpoint_format', &
409 string_val, "chkp")
410 this%chkp_out = chkp_output_t(this%chkp, name = name,&
411 path = this%output_directory, fmt = trim(string_val))
412 call json_get_or_default(this%params, 'case.checkpoint_control', &
413 string_val, "simulationtime")
414 call json_get_or_default(this%params, 'case.checkpoint_value', real_val,&
415 1e10_rp)
416 call this%output_controller%add(this%chkp_out, real_val, string_val, &
417 neko_eps)
418 end if
419
420 !
421 ! Setup joblimit
422 !
423 if (this%params%valid_path('case.job_timelimit')) then
424 call json_get(this%params, 'case.job_timelimit', string_val)
425 call jobctrl_set_time_limit(string_val)
426 end if
427
428 call neko_log%end_section()
429
430 end subroutine case_init_common
431
433 subroutine case_free(this)
434 type(case_t), intent(inout) :: this
435
436 if (allocated(this%fluid)) then
437 call this%fluid%free()
438 deallocate(this%fluid)
439 end if
440
441 if (allocated(this%scalar)) then
442 call this%scalar%free()
443 deallocate(this%scalar)
444 end if
445
446 call this%msh%free()
447
448 call this%f_out%free()
449
450 call this%output_controller%free()
451
452 end subroutine case_free
453
454end 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:89
subroutine, public case_free(this)
Deallocate a case.
Definition case.f90:434
subroutine case_init_from_json(this, case_json)
Initialize a case from a JSON object describing a case.
Definition case.f90:124
subroutine case_init_common(this)
Initialize a case from its (loaded) params object.
Definition case.f90:138
Defines a checkpoint.
Defines an output for a checkpoint.
Definition comm.F90:1
integer pe_rank
MPI rank.
Definition comm.F90:51
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:38
integer pe_size
MPI size of communicator.
Definition comm.F90:54
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.
subroutine, public json_extract_object(json, name, object)
Extract object as a separate JSON dictionary.
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_quiet
Always logged.
Definition log.f90:67
type(log_t), public neko_log
Global log stream.
Definition log.f90:65
Definition math.f90:60
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
Defines a mesh field.
subroutine, public mesh_field_free(fld)
subroutine, public mesh_field_init(fld, msh, fld_name)
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:111
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:56
Scalar initial condition.
Definition scalar_ic.f90:34
Contains the scalar_pnpn_t type.
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:54
Base type of all fluid formulations.
Centralized controller for a list of outputs.
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.