Neko  0.9.0
A portable framework for high-order spectral element flow simulations
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 !
34 module 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
38  use fluid_output, only : fluid_output_t
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 stats, only : stats_t
47  use file, only : file_t
48  use utils, only : neko_error
49  use mesh, only : mesh_t
50  use comm
54  use user_intf, only : user_t
55  use scalar_pnpn, only : scalar_pnpn_t
56  use json_module, only : json_file, json_core, json_value
60  implicit none
61  private
62 
63  type, public :: case_t
64  type(mesh_t) :: msh
65  type(json_file) :: params
66  type(time_scheme_controller_t) :: ext_bdf
67  real(kind=rp), dimension(10) :: tlag
68  real(kind=rp), dimension(10) :: dtlag
69  real(kind=rp) :: dt
70  real(kind=rp) :: end_time
71  character(len=:), allocatable :: output_directory
73  type(fluid_output_t) :: f_out
74  type(chkp_output_t) :: f_chkp
75  type(user_t) :: usr
76  class(fluid_scheme_t), allocatable :: fluid
77  type(scalar_pnpn_t), allocatable :: scalar
78  end type case_t
79 
80  interface case_init
82  end interface case_init
83 
84  public :: case_init, case_free
85 
86 contains
87 
89  subroutine case_init_from_file(this, case_file)
90  type(case_t), target, intent(inout) :: this
91  character(len=*), intent(in) :: case_file
92  integer :: ierr, integer_val
93  character(len=:), allocatable :: json_buffer
94 
95  call neko_log%section('Case')
96  call neko_log%message('Reading case file ' // trim(case_file), &
98 
99  if (pe_rank .eq. 0) then
100  call this%params%load_file(filename = trim(case_file))
101  call this%params%print_to_string(json_buffer)
102  integer_val = len(json_buffer)
103  end if
104 
105  call mpi_bcast(integer_val, 1, mpi_integer, 0, neko_comm, ierr)
106  if (pe_rank .ne. 0) allocate(character(len = integer_val) :: json_buffer)
107  call mpi_bcast(json_buffer, integer_val, mpi_character, 0, neko_comm, ierr)
108  call this%params%load_from_string(json_buffer)
109 
110  deallocate(json_buffer)
111 
112  call case_init_common(this)
113 
114  end subroutine case_init_from_file
115 
117  subroutine case_init_from_json(this, case_json)
118  type(case_t), target, intent(inout) :: this
119  type(json_file), intent(in) :: case_json
120 
121  call neko_log%section('Case')
122  call neko_log%message('Creating case from JSON object', neko_log_quiet)
123 
124  this%params = case_json
125 
126  call case_init_common(this)
127 
128  end subroutine case_init_from_json
129 
131  subroutine case_init_common(this)
132  type(case_t), target, intent(inout) :: this
133  integer :: lx = 0
134  logical :: scalar = .false.
135  type(file_t) :: msh_file, bdry_file, part_file
136  type(mesh_fld_t) :: msh_part, parts
137  logical :: found, logical_val
138  integer :: integer_val
139  real(kind=rp) :: real_val
140  character(len = :), allocatable :: string_val
141  real(kind=rp) :: stats_start_time, stats_output_val
142  integer :: stats_sampling_interval
143  integer :: output_dir_len
144  integer :: precision
145 
146  !
147  ! Load mesh
148  !
149  call json_get_or_default(this%params, 'case.mesh_file', string_val, &
150  'no mesh')
151  if (trim(string_val) .eq. 'no mesh') then
152  call neko_error('The mesh_file keyword could not be found in the .' // &
153  'case file. Often caused by incorrectly formatted json.')
154  end if
155  msh_file = file_t(string_val)
156 
157  call msh_file%read(this%msh)
158 
159  !
160  ! Load Balancing
161  !
162  call json_get_or_default(this%params, 'case.load_balancing', logical_val,&
163  .false.)
164 
165  if (pe_size .gt. 1 .and. logical_val) then
166  call neko_log%section('Load Balancing')
167  call parmetis_partmeshkway(this%msh, parts)
168  call redist_mesh(this%msh, parts)
169  call neko_log%end_section()
170  end if
171 
172  !
173  ! Time step
174  !
175  call this%params%get('case.variable_timestep', logical_val, found)
176  if (.not. logical_val) then
177  call json_get(this%params, 'case.timestep', this%dt)
178  else
179  ! randomly set an initial dt to get cfl when dt is variable
180  this%dt = 1.0_rp
181  end if
182 
183  !
184  ! End time
185  !
186  call json_get(this%params, 'case.end_time', this%end_time)
187 
188  !
189  ! Initialize point_zones registry
190  !
191  call neko_point_zone_registry%init(this%params, this%msh)
192 
193  !
194  ! Setup user defined functions
195  !
196  call this%usr%init()
197  call this%usr%user_mesh_setup(this%msh)
198 
199  !
200  ! Set order of timestepper
201  !
202  call json_get(this%params, 'case.numerics.time_order', integer_val)
203  call this%ext_bdf%init(integer_val)
204 
205  !
206  ! Setup fluid scheme
207  !
208  call json_get(this%params, 'case.fluid.scheme', string_val)
209  call fluid_scheme_factory(this%fluid, trim(string_val))
210 
211  call json_get(this%params, 'case.numerics.polynomial_order', lx)
212  lx = lx + 1 ! add 1 to get number of gll points
213  this%fluid%chkp%tlag => this%tlag
214  this%fluid%chkp%dtlag => this%dtlag
215  call this%fluid%init(this%msh, lx, this%params, this%usr, this%ext_bdf)
216  select type (f => this%fluid)
217  type is (fluid_pnpn_t)
218  f%chkp%abx1 => f%abx1
219  f%chkp%abx2 => f%abx2
220  f%chkp%aby1 => f%aby1
221  f%chkp%aby2 => f%aby2
222  f%chkp%abz1 => f%abz1
223  f%chkp%abz2 => f%abz2
224  end select
225 
226 
227  !
228  ! Setup scratch registry
229  !
230  neko_scratch_registry = scratch_registry_t(this%fluid%dm_Xh, 10, 10)
231 
232  !
233  ! Setup scalar scheme
234  !
235  ! @todo no scalar factory for now, probably not needed
236  if (this%params%valid_path('case.scalar')) then
237  call json_get_or_default(this%params, 'case.scalar.enabled', scalar,&
238  .true.)
239  end if
240 
241  if (scalar) then
242  allocate(this%scalar)
243  this%scalar%chkp%tlag => this%tlag
244  this%scalar%chkp%dtlag => this%dtlag
245  call this%scalar%init(this%msh, this%fluid%c_Xh, this%fluid%gs_Xh, &
246  this%params, this%usr, this%fluid%ulag, this%fluid%vlag, &
247  this%fluid%wlag, this%ext_bdf, this%fluid%rho)
248 
249  call this%fluid%chkp%add_scalar(this%scalar%s)
250 
251  this%fluid%chkp%abs1 => this%scalar%abx1
252  this%fluid%chkp%abs2 => this%scalar%abx2
253  this%fluid%chkp%slag => this%scalar%slag
254  end if
255 
256  !
257  ! Setup user defined conditions
258  !
259  if (this%params%valid_path('case.fluid.inflow_condition')) then
260  call json_get(this%params, 'case.fluid.inflow_condition.type',&
261  string_val)
262  if (trim(string_val) .eq. 'user') then
263  call this%fluid%set_usr_inflow(this%usr%fluid_user_if)
264  end if
265  end if
266 
267  ! Setup user boundary conditions for the scalar.
268  if (scalar) then
269  call this%scalar%set_user_bc(this%usr%scalar_user_bc)
270  end if
271 
272  !
273  ! Setup initial conditions
274  !
275  call json_get(this%params, 'case.fluid.initial_condition.type',&
276  string_val)
277 
278 
279  call neko_log%section("Fluid initial condition ")
280 
281  if (trim(string_val) .ne. 'user') then
282  call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, &
283  this%fluid%p, this%fluid%c_Xh, this%fluid%gs_Xh, string_val, &
284  this%params)
285 
286  else
287  call set_flow_ic(this%fluid%u, this%fluid%v, this%fluid%w, this%fluid%p,&
288  this%fluid%c_Xh, this%fluid%gs_Xh, this%usr%fluid_user_ic, &
289  this%params)
290  end if
291 
292  call neko_log%end_section()
293 
294  if (scalar) then
295 
296  call json_get(this%params, 'case.scalar.initial_condition.type', &
297  string_val)
298 
299  call neko_log%section("Scalar initial condition ")
300 
301  if (trim(string_val) .ne. 'user') then
302  call set_scalar_ic(this%scalar%s, &
303  this%scalar%c_Xh, this%scalar%gs_Xh, string_val, this%params)
304  else
305  call set_scalar_ic(this%scalar%s, &
306  this%scalar%c_Xh, this%scalar%gs_Xh, this%usr%scalar_user_ic, &
307  this%params)
308  end if
309 
310  call neko_log%end_section()
311 
312  end if
313 
314  ! Add initial conditions to BDF scheme (if present)
315  select type (f => this%fluid)
316  type is (fluid_pnpn_t)
317  call f%ulag%set(f%u)
318  call f%vlag%set(f%v)
319  call f%wlag%set(f%w)
320  end select
321 
322  !
323  ! Validate that the case is properly setup for time-stepping
324  !
325  call this%fluid%validate
326 
327  if (scalar) then
328  call this%scalar%slag%set(this%scalar%s)
329  call this%scalar%validate
330  end if
331 
332  !
333  ! Get and process output directory
334  !
335  call json_get_or_default(this%params, 'case.output_directory',&
336  this%output_directory, '')
337 
338  output_dir_len = len(trim(this%output_directory))
339  if (output_dir_len .gt. 0) then
340  if (this%output_directory(output_dir_len:output_dir_len) .ne. "/") then
341  this%output_directory = trim(this%output_directory)//"/"
342  if (pe_rank .eq. 0) then
343  call execute_command_line('mkdir -p '//this%output_directory)
344  end if
345  end if
346  end if
347 
348  !
349  ! Save boundary markings for fluid (if requested)
350  !
351  call json_get_or_default(this%params, 'case.output_boundary',&
352  logical_val, .false.)
353  if (logical_val) then
354  bdry_file = file_t(trim(this%output_directory)//'bdry.fld')
355  call bdry_file%write(this%fluid%bdry)
356  end if
357 
358  !
359  ! Save mesh partitions (if requested)
360  !
361  call json_get_or_default(this%params, 'case.output_partitions',&
362  logical_val, .false.)
363  if (logical_val) then
364  call mesh_field_init(msh_part, this%msh, 'MPI_Rank')
365  msh_part%data = pe_rank
366  part_file = file_t(trim(this%output_directory)//'partitions.vtk')
367  call part_file%write(msh_part)
368  call mesh_field_free(msh_part)
369  end if
370 
371  !
372  ! Setup output precision of the field files
373  !
374  call json_get_or_default(this%params, 'case.output_precision', string_val,&
375  'single')
376 
377  if (trim(string_val) .eq. 'double') then
378  precision = dp
379  else
380  precision = sp
381  end if
382 
383  !
384  ! Setup output_controller
385  !
386  call this%output_controller%init(this%end_time)
387  if (scalar) then
388  this%f_out = fluid_output_t(precision, this%fluid, this%scalar, &
389  path = trim(this%output_directory))
390  else
391  this%f_out = fluid_output_t(precision, this%fluid, &
392  path = trim(this%output_directory))
393  end if
394 
395  call json_get_or_default(this%params, 'case.fluid.output_control',&
396  string_val, 'org')
397 
398  if (trim(string_val) .eq. 'org') then
399  ! yes, it should be real_val below for type compatibility
400  call json_get(this%params, 'case.nsamples', real_val)
401  call this%output_controller%add(this%f_out, real_val, 'nsamples')
402  else if (trim(string_val) .eq. 'never') then
403  ! Fix a dummy 0.0 output_value
404  call json_get_or_default(this%params, 'case.fluid.output_value', &
405  real_val, 0.0_rp)
406  call this%output_controller%add(this%f_out, 0.0_rp, string_val)
407  else
408  call json_get(this%params, 'case.fluid.output_value', real_val)
409  call this%output_controller%add(this%f_out, real_val, string_val)
410  end if
411 
412  !
413  ! Save checkpoints (if nothing specified, default to saving at end of sim)
414  !
415  call json_get_or_default(this%params, 'case.output_checkpoints',&
416  logical_val, .true.)
417  if (logical_val) then
418  call json_get_or_default(this%params, 'case.checkpoint_format', &
419  string_val, "chkp")
420  this%f_chkp = chkp_output_t(this%fluid%chkp, &
421  path = this%output_directory, fmt = trim(string_val))
422  call json_get_or_default(this%params, 'case.checkpoint_control', &
423  string_val, "simulationtime")
424  call json_get_or_default(this%params, 'case.checkpoint_value', real_val,&
425  1e10_rp)
426  call this%output_controller%add(this%f_chkp, real_val, string_val)
427  end if
428 
429  !
430  ! Setup joblimit
431  !
432  if (this%params%valid_path('case.job_timelimit')) then
433  call json_get(this%params, 'case.job_timelimit', string_val)
434  call jobctrl_set_time_limit(string_val)
435  end if
436 
437  call neko_log%end_section()
438 
439  end subroutine case_init_common
440 
442  subroutine case_free(this)
443  type(case_t), intent(inout) :: this
444 
445  if (allocated(this%fluid)) then
446  call this%fluid%free()
447  deallocate(this%fluid)
448  end if
449 
450  if (allocated(this%scalar)) then
451  call this%scalar%free()
452  deallocate(this%scalar)
453  end if
454 
455  call this%msh%free()
456 
457  call this%f_out%free()
458 
459  call this%output_controller%free()
460 
461  end subroutine case_free
462 
463 end module case
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:54
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:45
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:90
subroutine, public case_free(this)
Deallocate a case.
Definition: case.f90:443
subroutine case_init_from_json(this, case_json)
Initialize a case from a JSON object describing a case.
Definition: case.f90:118
subroutine case_init_common(this)
Initialize a case from its (loaded) params object.
Definition: case.f90:132
Defines an output for a checkpoint.
Definition: chkp_output.f90:34
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.
Definition: fluid_pnpn.f90:34
Fluid formulations.
Job control.
Definition: jobctrl.f90:34
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
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
integer, parameter, public log_size
Definition: log.f90:42
Defines a mesh field.
Definition: mesh_field.f90:35
subroutine, public mesh_field_free(fld)
Definition: mesh_field.f90:73
subroutine, public mesh_field_init(fld, msh, fld_name)
Definition: mesh_field.f90:52
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.
Definition: scalar_pnpn.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.
Defines a container for all statistics.
Definition: statistics.f90:34
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.
Statistics backend.
Definition: statistics.f90:48
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...