Neko  0.8.99
A portable framework for high-order spectral element flow simulations
les_simcomp.f90
Go to the documentation of this file.
1 ! Copyright (c) 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 !
33 !
35 
37  use num_types, only : rp
38  use json_module, only : json_file
40  use case, only : case_t
41  use les_model, only : les_model_t, les_model_factory
43  use field_writer, only : field_writer_t
44  implicit none
45  private
46 
49  type, public, extends(simulation_component_t) :: les_simcomp_t
51  class(les_model_t), allocatable :: les_model
53  type(field_writer_t) :: writer
54  contains
56  procedure, pass(this) :: init => les_simcomp_init_from_json
58  procedure, pass(this) :: free => les_simcomp_free
60  procedure, pass(this) :: compute_ => les_simcomp_compute
61  end type les_simcomp_t
62 
63 contains
64 
66  subroutine les_simcomp_init_from_json(this, json, case)
67  class(les_simcomp_t), intent(inout) :: this
68  type(json_file), intent(inout) :: json
69  class(case_t), intent(inout), target :: case
70  character(len=:), allocatable :: name
71  character(len=:), allocatable :: nut_field
72  character(len=20) :: fields(2)
73 
74  call this%free()
75 
76  ! Add fields keyword to the json so that the field_writer picks it up.
77  ! Will also add fields to the registry if missing.
78  call json_get_or_default(json, "nut_field", nut_field, "nut")
79  fields(1) = "les_delta"
80  fields(2) = nut_field
81 
82  call json%add("fields", fields)
83 
84  call this%init_base(json, case)
85  call this%writer%init(json, case)
86 
87  call json_get(json, "model", name)
88 
89  call les_model_factory(this%les_model, name, case%fluid%dm_Xh,&
90  case%fluid%c_Xh, json)
91 
92  end subroutine les_simcomp_init_from_json
93 
95  subroutine les_simcomp_free(this)
96  class(les_simcomp_t), intent(inout) :: this
97  call this%free_base()
98  call this%writer%free()
99 
100  if (allocated(this%les_model)) then
101  call this%les_model%free()
102  deallocate(this%les_model)
103  end if
104  end subroutine les_simcomp_free
105 
109  subroutine les_simcomp_compute(this, t, tstep)
110  class(les_simcomp_t), intent(inout) :: this
111  real(kind=rp), intent(in) :: t
112  integer, intent(in) :: tstep
113 
114  call this%les_model%compute(t, tstep)
115  end subroutine les_simcomp_compute
116 
117 end module les_simcomp
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
Implements the field_writer_t type.
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Implements les_model_t.
Definition: les_model.f90:35
Implements the les_simcomp_t type.
Definition: les_simcomp.f90:36
subroutine les_simcomp_init_from_json(this, json, case)
Constructor from json.
Definition: les_simcomp.f90:67
subroutine les_simcomp_compute(this, t, tstep)
Compute the les_simcomp field.
subroutine les_simcomp_free(this)
Destructor.
Definition: les_simcomp.f90:96
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
A simulation component that writes a 3d field to a file.
Base abstract type for LES models based on the Boussinesq approximation.
Definition: les_model.f90:51
A simulation component that drives the computation of the SGS viscosity.
Definition: les_simcomp.f90:49
Base abstract class for simulation components.