Neko  0.8.1
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
36  use fluid_fctry
37  use fluid_output
38  use chkp_output
42  use mpi_f08
43  use mesh_field
44  use parmetis
45  use redist
46  use sampler
47  use flow_ic
48  use scalar_ic, only : set_scalar_ic
49  use stats
50  use file
51  use utils
52  use mesh
53  use comm
55  use logger
56  use jobctrl
57  use user_intf
58  use scalar_pnpn ! todo directly load the pnpn? can we have other
59  use json_module, only : json_file, json_core, json_value
64  implicit none
65 
66  type :: case_t
67  type(mesh_t) :: msh
68  type(json_file) :: params
69  type(time_scheme_controller_t) :: ext_bdf
70  real(kind=rp), dimension(10) :: tlag
71  real(kind=rp), dimension(10) :: dtlag
72  real(kind=rp) :: dt
73  real(kind=rp) :: end_time
74  type(sampler_t) :: s
75  type(fluid_output_t) :: f_out
76  type(fluid_stats_output_t) :: f_stats_output
77  type(chkp_output_t) :: f_chkp
78  type(mean_flow_output_t) :: f_mf
79  type(mean_sqr_flow_output_t) :: f_msqrf
80  type(stats_t) :: q
81  type(user_t) :: usr
82  class(fluid_scheme_t), allocatable :: fluid
83  type(scalar_pnpn_t), allocatable :: scalar
85  end type case_t
86 
87  interface case_init
89  end interface case_init
90 
92 
93 contains
94 
96  subroutine case_init_from_file(C, case_file)
97  type(case_t), target, intent(inout) :: C
98  character(len=*), intent(in) :: case_file
99  integer :: ierr, integer_val
100  character(len=:), allocatable :: json_buffer
101 
102  call neko_log%section('Case')
103  call neko_log%message('Reading case file ' // trim(case_file), &
105 
106  if (pe_rank .eq. 0) then
107  call c%params%load_file(filename=trim(case_file))
108  call c%params%print_to_string(json_buffer)
109  integer_val = len(json_buffer)
110  end if
111 
112  call mpi_bcast(integer_val, 1, mpi_integer, 0, neko_comm, ierr)
113  if (pe_rank .ne. 0) allocate(character(len=integer_val)::json_buffer)
114  call mpi_bcast(json_buffer, integer_val, mpi_character, 0, neko_comm, ierr)
115  call c%params%load_from_string(json_buffer)
116 
117  deallocate(json_buffer)
118 
119  call case_init_common(c)
120 
121  end subroutine case_init_from_file
122 
124  subroutine case_init_from_json(C, case_json)
125  type(case_t), target, intent(inout) :: C
126  type(json_file), intent(in) :: case_json
127 
128  call neko_log%section('Case')
129  call neko_log%message('Creating case from JSON object', neko_log_quiet)
130 
131  c%params = case_json
132 
133  call case_init_common(c)
134 
135  end subroutine case_init_from_json
136 
138  subroutine case_init_common(C)
139  type(case_t), target, intent(inout) :: c
140  character(len=:), allocatable :: output_directory
141  integer :: lx = 0
142  logical :: scalar = .false.
143  type(file_t) :: msh_file, bdry_file, part_file
144  type(mesh_fld_t) :: msh_part, parts
145  logical :: found, logical_val
146  integer :: integer_val
147  real(kind=rp) :: real_val
148  character(len=:), allocatable :: string_val
149  real(kind=rp) :: stats_start_time, stats_output_val
150  integer :: stats_sampling_interval
151  integer :: output_dir_len
152  integer :: n_simcomps
153  integer :: precision
154 
155  !
156  ! Load mesh
157  !
158  call json_get(c%params, 'case.mesh_file', string_val)
159  msh_file = file_t(string_val)
160 
161  call msh_file%read(c%msh)
162 
163  !
164  ! Load Balancing
165  !
166  call json_get_or_default(c%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(c%msh, parts)
172  call redist_mesh(c%msh, parts)
173  call neko_log%end_section()
174  end if
175 
176  !
177  ! Time step
178  !
179  call c%params%get('case.variable_timestep', logical_val, found)
180  if (.not. logical_val) then
181  call json_get(c%params, 'case.timestep', c%dt)
182  else
183  ! randomly set an initial dt to get cfl when dt is variable
184  c%dt = 1.0_rp
185  end if
186 
187  !
188  ! End time
189  !
190  call json_get(c%params, 'case.end_time', c%end_time)
191 
192  !
193  ! Initialize point_zones registry
194  !
195  call neko_point_zone_registry%init(c%params, c%msh)
196 
197  !
198  ! Setup user defined functions
199  !
200  call c%usr%init()
201  call c%usr%user_mesh_setup(c%msh)
202 
203 
204  !
205  ! Material properties
206  !
207  call c%material_properties%init(c%params, c%usr)
208 
209  !
210  ! Setup fluid scheme
211  !
212  call json_get(c%params, 'case.fluid.scheme', string_val)
213  call fluid_scheme_factory(c%fluid, trim(string_val))
214 
215  call json_get(c%params, 'case.numerics.polynomial_order', lx)
216  lx = lx + 1 ! add 1 to get poly order
217  call c%fluid%init(c%msh, lx, c%params, c%usr, c%material_properties)
218  c%fluid%chkp%tlag => c%tlag
219  c%fluid%chkp%dtlag => c%dtlag
220  select type(f => c%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(c%fluid%dm_Xh, 10, 10)
235 
236  !
237  ! Setup scalar scheme
238  !
239  ! @todo no scalar factory for now, probably not needed
240  if (c%params%valid_path('case.scalar')) then
241  call json_get_or_default(c%params, 'case.scalar.enabled', scalar,&
242  .true.)
243  end if
244 
245  if (scalar) then
246  allocate(c%scalar)
247  call c%scalar%init(c%msh, c%fluid%c_Xh, c%fluid%gs_Xh, c%params, c%usr,&
248  c%material_properties)
249  call c%fluid%chkp%add_scalar(c%scalar%s)
250  c%fluid%chkp%abs1 => c%scalar%abx1
251  c%fluid%chkp%abs2 => c%scalar%abx2
252  c%fluid%chkp%slag => c%scalar%slag
253  end if
254 
255  !
256  ! Setup user defined conditions
257  !
258  if (c%params%valid_path('case.fluid.inflow_condition')) then
259  call json_get(c%params, 'case.fluid.inflow_condition.type',&
260  string_val)
261  if (trim(string_val) .eq. 'user') then
262  call c%fluid%set_usr_inflow(c%usr%fluid_user_if)
263  end if
264  end if
265 
266  ! Setup user boundary conditions for the scalar.
267  if (scalar) then
268  call c%scalar%set_user_bc(c%usr%scalar_user_bc)
269  end if
270 
271  !
272  ! Setup initial conditions
273  !
274  call json_get(c%params, 'case.fluid.initial_condition.type',&
275  string_val)
276  if (trim(string_val) .ne. 'user') then
277  call set_flow_ic(c%fluid%u, c%fluid%v, c%fluid%w, c%fluid%p, &
278  c%fluid%c_Xh, c%fluid%gs_Xh, string_val, c%params)
279  else
280  call set_flow_ic(c%fluid%u, c%fluid%v, c%fluid%w, c%fluid%p, &
281  c%fluid%c_Xh, c%fluid%gs_Xh, c%usr%fluid_user_ic, c%params)
282  end if
283 
284  if (scalar) then
285  call json_get(c%params, 'case.scalar.initial_condition.type', string_val)
286  if (trim(string_val) .ne. 'user') then
287  call set_scalar_ic(c%scalar%s, &
288  c%scalar%c_Xh, c%scalar%gs_Xh, string_val, c%params)
289  else
290  call set_scalar_ic(c%scalar%s, &
291  c%scalar%c_Xh, c%scalar%gs_Xh, c%usr%scalar_user_ic, c%params)
292  end if
293  end if
294 
295  ! Add initial conditions to BDF scheme (if present)
296  select type(f => c%fluid)
297  type is(fluid_pnpn_t)
298  call f%ulag%set(f%u)
299  call f%vlag%set(f%v)
300  call f%wlag%set(f%w)
301  end select
302 
303  !
304  ! Validate that the case is properly setup for time-stepping
305  !
306  call c%fluid%validate
307 
308  if (scalar) then
309  call c%scalar%slag%set(c%scalar%s)
310  call c%scalar%validate
311  end if
312 
313  !
314  ! Set order of timestepper
315  !
316  call json_get(c%params, 'case.numerics.time_order', integer_val)
317  call c%ext_bdf%init(integer_val)
318 
319  !
320  ! Get and process output directory
321  !
322  call json_get_or_default(c%params, 'case.output_directory',&
323  output_directory, '')
324 
325  output_dir_len = len(trim(output_directory))
326  if (output_dir_len .gt. 0) then
327  if (output_directory(output_dir_len:output_dir_len) .ne. "/") then
328  output_directory = trim(output_directory)//"/"
329  if (pe_rank .eq. 0) then
330  call execute_command_line('mkdir -p '//output_directory)
331  end if
332  end if
333  end if
334 
335  !
336  ! Save boundary markings for fluid (if requested)
337  !
338  call json_get_or_default(c%params, 'case.output_boundary',&
339  logical_val, .false.)
340  if (logical_val) then
341  bdry_file = file_t(trim(output_directory)//'bdry.fld')
342  call bdry_file%write(c%fluid%bdry)
343  end if
344 
345  !
346  ! Save mesh partitions (if requested)
347  !
348  call json_get_or_default(c%params, 'case.output_partitions',&
349  logical_val, .false.)
350  if (logical_val) then
351  call mesh_field_init(msh_part, c%msh, 'MPI_Rank')
352  msh_part%data = pe_rank
353  part_file = file_t(trim(output_directory)//'partitions.vtk')
354  call part_file%write(msh_part)
355  call mesh_field_free(msh_part)
356  end if
357 
358  !
359  ! Setup output precision of the field files
360  !
361  call json_get_or_default(c%params, 'case.output_precision', string_val,&
362  'single')
363 
364  if (trim(string_val) .eq. 'double') then
365  precision = dp
366  else
367  precision = sp
368  end if
369 
370  !
371  ! Setup sampler
372  !
373  call c%s%init(c%end_time)
374  if (scalar) then
375  c%f_out = fluid_output_t(precision, c%fluid, c%scalar, &
376  path=trim(output_directory))
377  else
378  c%f_out = fluid_output_t(precision, c%fluid, &
379  path=trim(output_directory))
380  end if
381 
382  call json_get_or_default(c%params, 'case.fluid.output_control',&
383  string_val, 'org')
384 
385  if (trim(string_val) .eq. 'org') then
386  ! yes, it should be real_val below for type compatibility
387  call json_get(c%params, 'case.nsamples', real_val)
388  call c%s%add(c%f_out, real_val, 'nsamples')
389  else if (trim(string_val) .eq. 'never') then
390  ! Fix a dummy 0.0 output_value
391  call json_get_or_default(c%params, 'case.fluid.output_value', real_val, &
392  0.0_rp)
393  call c%s%add(c%f_out, 0.0_rp, string_val)
394  else
395  call json_get(c%params, 'case.fluid.output_value', real_val)
396  call c%s%add(c%f_out, real_val, string_val)
397  end if
398 
399  !
400  ! Save checkpoints (if nothing specified, default to saving at end of sim)
401  !
402  call json_get_or_default(c%params, 'case.output_checkpoints',&
403  logical_val, .true.)
404  if (logical_val) then
405  c%f_chkp = chkp_output_t(c%fluid%chkp, path=output_directory)
406  call json_get_or_default(c%params, 'case.checkpoint_control', &
407  string_val, "simulationtime")
408  call json_get_or_default(c%params, 'case.checkpoint_value', real_val,&
409  1e10_rp)
410  call c%s%add(c%f_chkp, real_val, string_val)
411  end if
412 
413  !
414  ! Setup statistics
415  !
416 
417  ! Always init, so that we can call eval in simulation.f90 with no if.
418  ! Note, don't use json_get_or_default here, because that will break the
419  ! valid_path if statement below (the path will become valid always).
420  call c%params%get('case.statistics.start_time', stats_start_time,&
421  found)
422  if (.not. found) stats_start_time = 0.0_rp
423 
424  call c%params%get('case.statistics.sampling_interval', &
425  stats_sampling_interval, found)
426  if (.not. found) stats_sampling_interval = 10
427 
428  call c%q%init(stats_start_time, stats_sampling_interval)
429 
430  found = c%params%valid_path('case.statistics')
431  if (found) then
432  call json_get_or_default(c%params, 'case.statistics.enabled',&
433  logical_val, .true.)
434  if (logical_val) then
435  call c%q%add(c%fluid%mean%u)
436  call c%q%add(c%fluid%mean%v)
437  call c%q%add(c%fluid%mean%w)
438  call c%q%add(c%fluid%mean%p)
439 
440  c%f_mf = mean_flow_output_t(c%fluid%mean, stats_start_time, &
441  path=output_directory)
442 
443  call json_get(c%params, 'case.statistics.output_control', &
444  string_val)
445  call json_get(c%params, 'case.statistics.output_value', &
446  stats_output_val)
447 
448  call c%s%add(c%f_mf, stats_output_val, string_val)
449  call c%q%add(c%fluid%stats)
450 
451  c%f_stats_output = fluid_stats_output_t(c%fluid%stats, &
452  stats_start_time, path=output_directory)
453  call c%s%add(c%f_stats_output, stats_output_val, string_val)
454  end if
455  end if
456 
457 ! if (C%params%stats_mean_sqr_flow) then
458 ! call C%q%add(C%fluid%mean_sqr%uu)
459 ! call C%q%add(C%fluid%mean_sqr%vv)
460 ! call C%q%add(C%fluid%mean_sqr%ww)
461 ! call C%q%add(C%fluid%mean_sqr%pp)
462 
463 ! if (C%params%output_mean_sqr_flow) then
464 ! C%f_msqrf = mean_sqr_flow_output_t(C%fluid%mean_sqr, &
465 ! C%params%stats_begin, &
466 ! path=output_directory)
467 ! call C%s%add(C%f_msqrf, C%params%stats_write_par, &
468 ! C%params%stats_write_control)
469 ! end if
470 ! end if
471 
472  !
473  ! Setup joblimit
474  !
475  if (c%params%valid_path('case.job_timelimit')) then
476  call json_get(c%params, 'case.job_timelimit', string_val)
477  call jobctrl_set_time_limit(string_val)
478  end if
479 
480  call neko_log%end_section()
481 
482  end subroutine case_init_common
483 
485  subroutine case_free(C)
486  type(case_t), intent(inout) :: C
487 
488  if (allocated(c%fluid)) then
489  call c%fluid%free()
490  deallocate(c%fluid)
491  end if
492 
493  if (allocated(c%scalar)) then
494  call c%scalar%free()
495  deallocate(c%scalar)
496  end if
497 
498  call c%msh%free()
499 
500  call c%s%free()
501 
502  call c%q%free()
503 
504  end subroutine case_free
505 
506 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:53
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:44
Defines a simulation case.
Definition: case.f90:34
subroutine case_free(C)
Deallocate a case.
Definition: case.f90:486
subroutine, private case_init_from_json(C, case_json)
Initialize a case from a JSON object describing a case.
Definition: case.f90:125
subroutine, private case_init_from_file(C, case_file)
Initialize a case from an input file case_file.
Definition: case.f90:97
subroutine, private case_init_common(C)
Initialize a case from its (loaded) params object.
Definition: case.f90:139
Defines an output for a checkpoint.
Definition: chkp_output.f90:34
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:26
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
integer pe_size
MPI size of communicator.
Definition: comm.F90:29
Module for file I/O operations.
Definition: file.f90:34
Initial flow condition.
Definition: flow_ic.f90:34
Factory for all fluid schemes.
Definition: fluid_fctry.f90:34
subroutine fluid_scheme_factory(fluid, fluid_scheme)
Initialise a fluid scheme.
Definition: fluid_fctry.f90:46
Defines an output for a fluid.
Defines an output for a mean flow field.
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:63
type(log_t), public neko_log
Global log stream.
Definition: log.f90:61
Implements material_properties_t type.
Defines an output for a mean flow field.
Defines an output for a mean squared flow field.
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
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
Defines a sampler.
Definition: sampler.f90:34
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:53
Contains all the material properties necessary in the simulation.
Statistics backend.
Definition: statistics.f90:47
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...