Neko 0.9.99
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-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
37 use fluid_scheme, only : fluid_scheme_t, fluid_scheme_factory
39 use chkp_output, only : chkp_output_t
42 use redist, only : redist_mesh
44 use flow_ic, only : set_flow_ic
45 use scalar_ic, only : set_scalar_ic
46 use file, only : file_t
47 use utils, only : neko_error
48 use mesh, only : mesh_t
49 use comm
51 use logger, only : neko_log, neko_log_quiet
53 use user_intf, only : user_t
54 use scalar_pnpn, only : scalar_pnpn_t
55 use json_module, only : json_file
59 implicit none
60 private
61
62 type, public :: case_t
63 type(mesh_t) :: msh
64 type(json_file) :: params
65 type(time_scheme_controller_t) :: ext_bdf
66 real(kind=rp), dimension(10) :: tlag
67 real(kind=rp), dimension(10) :: dtlag
68 real(kind=rp) :: dt
69 real(kind=rp) :: end_time
70 character(len=:), allocatable :: output_directory
72 type(fluid_output_t) :: f_out
73 type(chkp_output_t) :: f_chkp
74 type(user_t) :: usr
75 class(fluid_scheme_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
147 integer :: output_dir_len
148 integer :: precision
149
150 !
151 ! Load mesh
152 !
153 call json_get_or_default(this%params, 'case.mesh_file', string_val, &
154 'no mesh')
155 if (trim(string_val) .eq. 'no mesh') then
156 call neko_error('The mesh_file keyword could not be found in the .' // &
157 'case file. Often caused by incorrectly formatted json.')
158 end if
159 msh_file = file_t(string_val)
160
161 call msh_file%read(this%msh)
162
163 !
164 ! Load Balancing
165 !
166 call json_get_or_default(this%params, 'case.load_balancing', logical_val,&
167 .false.)
168
169 if (pe_size .gt. 1 .and. logical_val) then
170 call neko_log%section('Load Balancing')
171 call parmetis_partmeshkway(this%msh, parts)
172 call redist_mesh(this%msh, parts)
173 call neko_log%end_section()
174 end if
175
176 !
177 ! Time step
178 !
179 call this%params%get('case.variable_timestep', logical_val, found)
180 if (.not. logical_val) then
181 call json_get(this%params, 'case.timestep', this%dt)
182 else
183 ! randomly set an initial dt to get cfl when dt is variable
184 this%dt = 1.0_rp
185 end if
186
187 !
188 ! End time
189 !
190 call json_get(this%params, 'case.end_time', this%end_time)
191
192 !
193 ! Initialize point_zones registry
194 !
195 call neko_point_zone_registry%init(this%params, this%msh)
196
197 !
198 ! Setup user defined functions
199 !
200 call this%usr%init()
201 call this%usr%user_mesh_setup(this%msh)
202
203 !
204 ! Set order of timestepper
205 !
206 call json_get(this%params, 'case.numerics.time_order', integer_val)
207 call this%ext_bdf%init(integer_val)
208
209 !
210 ! Setup fluid scheme
211 !
212 call json_get(this%params, 'case.fluid.scheme', string_val)
213 call fluid_scheme_factory(this%fluid, trim(string_val))
214
215 call json_get(this%params, 'case.numerics.polynomial_order', lx)
216 lx = lx + 1 ! add 1 to get number of gll points
217 this%fluid%chkp%tlag => this%tlag
218 this%fluid%chkp%dtlag => this%dtlag
219 call this%fluid%init(this%msh, lx, this%params, this%usr, this%ext_bdf)
220 select type (f => this%fluid)
221 type is (fluid_pnpn_t)
222 f%chkp%abx1 => f%abx1
223 f%chkp%abx2 => f%abx2
224 f%chkp%aby1 => f%aby1
225 f%chkp%aby2 => f%aby2
226 f%chkp%abz1 => f%abz1
227 f%chkp%abz2 => f%abz2
228 end select
229
230
231 !
232 ! Setup scratch registry
233 !
234 neko_scratch_registry = scratch_registry_t(this%fluid%dm_Xh, 10, 10)
235
236 !
237 ! Setup scalar scheme
238 !
239 ! @todo no scalar factory for now, probably not needed
240 if (this%params%valid_path('case.scalar')) then
241 call json_get_or_default(this%params, 'case.scalar.enabled', scalar,&
242 .true.)
243 end if
244
245 if (scalar) then
246 allocate(this%scalar)
247 this%scalar%chkp%tlag => this%tlag
248 this%scalar%chkp%dtlag => this%dtlag
249 call this%scalar%init(this%msh, this%fluid%c_Xh, this%fluid%gs_Xh, &
250 this%params, this%usr, this%fluid%ulag, this%fluid%vlag, &
251 this%fluid%wlag, this%ext_bdf, this%fluid%rho)
252
253 call this%fluid%chkp%add_scalar(this%scalar%s)
254
255 this%fluid%chkp%abs1 => this%scalar%abx1
256 this%fluid%chkp%abs2 => this%scalar%abx2
257 this%fluid%chkp%slag => this%scalar%slag
258 end if
259
260 !
261 ! Setup user defined conditions
262 !
263 if (this%params%valid_path('case.fluid.inflow_condition')) then
264 call json_get(this%params, 'case.fluid.inflow_condition.type',&
265 string_val)
266 if (trim(string_val) .eq. 'user') then
267 call this%fluid%set_usr_inflow(this%usr%fluid_user_if)
268 end if
269 end if
270
271 ! Setup user boundary conditions for the scalar.
272 if (scalar) then
273 call this%scalar%set_user_bc(this%usr%scalar_user_bc)
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
282 call neko_log%section("Fluid initial condition ")
283
284 if (trim(string_val) .ne. 'user') then
285 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, &
286 this%fluid%p, this%fluid%c_Xh, this%fluid%gs_Xh, string_val, &
287 this%params)
288
289 else
290 call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, this%fluid%p,&
291 this%fluid%c_Xh, this%fluid%gs_Xh, this%usr%fluid_user_ic, &
292 this%params)
293 end if
294
295 call neko_log%end_section()
296
297 if (scalar) then
298
299 call json_get(this%params, 'case.scalar.initial_condition.type', &
300 string_val)
301
302 call neko_log%section("Scalar initial condition ")
303
304 if (trim(string_val) .ne. 'user') then
305 call set_scalar_ic(this%scalar%s, &
306 this%scalar%c_Xh, this%scalar%gs_Xh, string_val, this%params)
307 else
308 call set_scalar_ic(this%scalar%s, &
309 this%scalar%c_Xh, this%scalar%gs_Xh, this%usr%scalar_user_ic, &
310 this%params)
311 end if
312
313 call neko_log%end_section()
314
315 end if
316
317 ! Add initial conditions to BDF scheme (if present)
318 select type (f => this%fluid)
319 type is (fluid_pnpn_t)
320 call f%ulag%set(f%u)
321 call f%vlag%set(f%v)
322 call f%wlag%set(f%w)
323 end select
324
325 !
326 ! Validate that the case is properly setup for time-stepping
327 !
328 call this%fluid%validate
329
330 if (scalar) then
331 call this%scalar%slag%set(this%scalar%s)
332 call this%scalar%validate
333 end if
334
335 !
336 ! Get and process output directory
337 !
338 call json_get_or_default(this%params, 'case.output_directory',&
339 this%output_directory, '')
340
341 output_dir_len = len(trim(this%output_directory))
342 if (output_dir_len .gt. 0) then
343 if (this%output_directory(output_dir_len:output_dir_len) .ne. "/") then
344 this%output_directory = trim(this%output_directory)//"/"
345 if (pe_rank .eq. 0) then
346 call execute_command_line('mkdir -p '//this%output_directory)
347 end if
348 end if
349 end if
350
351 !
352 ! Save boundary markings for fluid (if requested)
353 !
354 call json_get_or_default(this%params, 'case.output_boundary',&
355 logical_val, .false.)
356 if (logical_val) then
357 bdry_file = file_t(trim(this%output_directory)//'bdry.fld')
358 call bdry_file%write(this%fluid%bdry)
359 end if
360
361 !
362 ! Save mesh partitions (if requested)
363 !
364 call json_get_or_default(this%params, 'case.output_partitions',&
365 logical_val, .false.)
366 if (logical_val) then
367 call mesh_field_init(msh_part, this%msh, 'MPI_Rank')
368 msh_part%data = pe_rank
369 part_file = file_t(trim(this%output_directory)//'partitions.vtk')
370 call part_file%write(msh_part)
371 call mesh_field_free(msh_part)
372 end if
373
374 !
375 ! Setup output precision of the field files
376 !
377 call json_get_or_default(this%params, 'case.output_precision', string_val,&
378 'single')
379
380 if (trim(string_val) .eq. 'double') then
381 precision = dp
382 else
383 precision = sp
384 end if
385
386 !
387 ! Setup output_controller
388 !
389 call this%output_controller%init(this%end_time)
390 if (scalar) then
391 this%f_out = fluid_output_t(precision, this%fluid, this%scalar, &
392 path = trim(this%output_directory))
393 else
394 this%f_out = fluid_output_t(precision, this%fluid, &
395 path = trim(this%output_directory))
396 end if
397
398 call json_get_or_default(this%params, 'case.fluid.output_control',&
399 string_val, 'org')
400
401 if (trim(string_val) .eq. 'org') then
402 ! yes, it should be real_val below for type compatibility
403 call json_get(this%params, 'case.nsamples', real_val)
404 call this%output_controller%add(this%f_out, real_val, 'nsamples')
405 else if (trim(string_val) .eq. 'never') then
406 ! Fix a dummy 0.0 output_value
407 call json_get_or_default(this%params, 'case.fluid.output_value', &
408 real_val, 0.0_rp)
409 call this%output_controller%add(this%f_out, 0.0_rp, string_val)
410 else
411 call json_get(this%params, 'case.fluid.output_value', real_val)
412 call this%output_controller%add(this%f_out, real_val, string_val)
413 end if
414
415 !
416 ! Save checkpoints (if nothing specified, default to saving at end of sim)
417 !
418 call json_get_or_default(this%params, 'case.output_checkpoints',&
419 logical_val, .true.)
420 if (logical_val) then
421 call json_get_or_default(this%params, 'case.checkpoint_format', &
422 string_val, "chkp")
423 this%f_chkp = chkp_output_t(this%fluid%chkp, &
424 path = this%output_directory, fmt = trim(string_val))
425 call json_get_or_default(this%params, 'case.checkpoint_control', &
426 string_val, "simulationtime")
427 call json_get_or_default(this%params, 'case.checkpoint_value', real_val,&
428 1e10_rp)
429 call this%output_controller%add(this%f_chkp, real_val, string_val)
430 end if
431
432 !
433 ! Setup joblimit
434 !
435 if (this%params%valid_path('case.job_timelimit')) then
436 call json_get(this%params, 'case.job_timelimit', string_val)
437 call jobctrl_set_time_limit(string_val)
438 end if
439
440 call neko_log%end_section()
441
442 end subroutine case_init_common
443
445 subroutine case_free(this)
446 type(case_t), intent(inout) :: this
447
448 if (allocated(this%fluid)) then
449 call this%fluid%free()
450 deallocate(this%fluid)
451 end if
452
453 if (allocated(this%scalar)) then
454 call this%scalar%free()
455 deallocate(this%scalar)
456 end if
457
458 call this%msh%free()
459
460 call this%f_out%free()
461
462 call this%output_controller%free()
463
464 end subroutine case_free
465
466end 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:446
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 an output for a checkpoint.
Definition comm.F90:1
integer pe_rank
MPI rank.
Definition comm.F90:28
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:16
integer pe_size
MPI size of communicator.
Definition comm.F90:31
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.
Fluid formulations.
Job control.
Definition jobctrl.f90:34
Utilities for retrieving parameters from the case files.
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
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
Containts 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.
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 ...