Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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
44 use utils, only : neko_error
45 implicit none
46 private
47
50 type, public, extends(simulation_component_t) :: les_simcomp_t
52 class(les_model_t), allocatable :: les_model
54 type(field_writer_t) :: writer
55 contains
57 procedure, pass(this) :: init => les_simcomp_init_from_json
59 procedure, pass(this) :: free => les_simcomp_free
61 procedure, pass(this) :: compute_ => les_simcomp_compute
62 end type les_simcomp_t
63
64contains
65
67 subroutine les_simcomp_init_from_json(this, json, case)
68 class(les_simcomp_t), intent(inout) :: this
69 type(json_file), intent(inout) :: json
70 class(case_t), intent(inout), target :: case
71 character(len=:), allocatable :: name
72 character(len=:), allocatable :: nut_field
73 character(len=20) :: fields(2)
74
75 call this%free()
76
77 ! Check for whether eddy viscosity is enabled in fluid_scheme
78 if (case%fluid%variable_material_properties .eqv. .false.) then
79 call neko_error("Eddy viscosity is not acting &
80 &on the equations. &
81 &Please set up a nut_field option &
82 &in the fluid solver")
83 end if
84
85 ! Add fields keyword to the json so that the field_writer picks it up.
86 ! Will also add fields to the registry if missing.
87 call json_get_or_default(json, "nut_field", nut_field, "nut")
88 fields(1) = "les_delta"
89 fields(2) = nut_field
90
91 call json%add("fields", fields)
92
93 call this%init_base(json, case)
94 call this%writer%init(json, case)
95
96 call json_get(json, "model", name)
97
98 call les_model_factory(this%les_model, name, case%fluid%dm_Xh,&
99 case%fluid%c_Xh, json)
100
101 end subroutine les_simcomp_init_from_json
102
104 subroutine les_simcomp_free(this)
105 class(les_simcomp_t), intent(inout) :: this
106 call this%free_base()
107 call this%writer%free()
108
109 if (allocated(this%les_model)) then
110 call this%les_model%free()
111 deallocate(this%les_model)
112 end if
113 end subroutine les_simcomp_free
114
118 subroutine les_simcomp_compute(this, t, tstep)
119 class(les_simcomp_t), intent(inout) :: this
120 real(kind=rp), intent(in) :: t
121 integer, intent(in) :: tstep
122
123 call this%les_model%compute(t, tstep)
124 end subroutine les_simcomp_compute
125
126end module les_simcomp
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
Implements the field_writer_t type.
Utilities for retrieving parameters from the case files.
Implements les_model_t.
Definition les_model.f90:35
Implements the les_simcomp_t type.
subroutine les_simcomp_init_from_json(this, json, case)
Constructor from json.
subroutine les_simcomp_compute(this, t, tstep)
Compute the les_simcomp field.
subroutine les_simcomp_free(this)
Destructor.
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...
subroutine compute_(this, t, tstep)
Dummy compute function.
Utilities.
Definition utils.f90:35
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.
Base abstract class for simulation components.