Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
vreman.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!
35module vreman
36 use num_types, only : rp
37 use les_model, only : les_model_t
38 use field, only : field_t
41 use json_module, only : json_file
46 use logger, only : log_size, neko_log
47 implicit none
48 private
49
52 type, public, extends(les_model_t) :: vreman_t
54 real(kind=rp) :: c
55 contains
57 procedure, pass(this) :: init => vreman_init
59 procedure, pass(this) :: init_from_components => &
62 procedure, pass(this) :: free => vreman_free
64 procedure, pass(this) :: compute => vreman_compute
65 end type vreman_t
66
67contains
71 subroutine vreman_init(this, fluid, json)
72 class(vreman_t), intent(inout) :: this
73 class(fluid_scheme_base_t), intent(inout), target :: fluid
74 type(json_file), intent(inout) :: json
75 character(len=:), allocatable :: nut_name
76 real(kind=rp) :: c
77 character(len=:), allocatable :: delta_type
78 logical :: if_ext
79 character(len=LOG_SIZE) :: log_buf
80
81 call json_get_or_default(json, "nut_field", nut_name, "nut")
82 call json_get_or_default(json, "delta_type", delta_type, "pointwise")
83 ! Based on the Smagorinsky Cs = 0.17.
84 call json_get_or_default(json, "c", c, 0.07_rp)
85 call json_get_or_default(json, "extrapolation", if_ext, .true.)
86
87 call neko_log%section('LES model')
88 write(log_buf, '(A)') 'Model : Vreman'
89 call neko_log%message(log_buf)
90 write(log_buf, '(A, A)') 'Delta evaluation : ', delta_type
91 call neko_log%message(log_buf)
92 write(log_buf, '(A, E15.7)') 'c : ', c
93 call neko_log%message(log_buf)
94 write(log_buf, '(A, L1)') 'extrapolation : ', if_ext
95 call neko_log%message(log_buf)
96 call neko_log%end_section()
97
98 call vreman_init_from_components(this, fluid, c, nut_name, &
99 delta_type, if_ext)
100 end subroutine vreman_init
101
108 subroutine vreman_init_from_components(this, fluid, c, nut_name, &
109 delta_type, if_ext)
110 class(vreman_t), intent(inout) :: this
111 class(fluid_scheme_base_t), intent(inout), target :: fluid
112 real(kind=rp) :: c
113 character(len=*), intent(in) :: nut_name
114 character(len=*), intent(in) :: delta_type
115 logical, intent(in) :: if_ext
116
117 call this%free()
118
119 call this%init_base(fluid, nut_name, delta_type, if_ext)
120 this%c = c
121
122 end subroutine vreman_init_from_components
123
125 subroutine vreman_free(this)
126 class(vreman_t), intent(inout) :: this
127
128 call this%free_base()
129 end subroutine vreman_free
130
134 subroutine vreman_compute(this, t, tstep)
135 class(vreman_t), intent(inout) :: this
136 real(kind=rp), intent(in) :: t
137 integer, intent(in) :: tstep
138
139 type(field_t), pointer :: u, v, w, u_e, v_e, w_e
140
141 if (this%if_ext .eqv. .true.) then
142 ! Extrapolate the velocity fields
143 associate(ulag => this%ulag, vlag => this%vlag, &
144 wlag => this%wlag, ext_bdf => this%ext_bdf)
145
146 u => neko_field_registry%get_field_by_name("u")
147 v => neko_field_registry%get_field_by_name("v")
148 w => neko_field_registry%get_field_by_name("w")
149 u_e => neko_field_registry%get_field_by_name("u_e")
150 v_e => neko_field_registry%get_field_by_name("v_e")
151 w_e => neko_field_registry%get_field_by_name("w_e")
152
153 call this%sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
154 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
155
156 end associate
157 end if
158
159 ! Compute the eddy viscosity field
160 if (neko_bcknd_device .eq. 1) then
161 call vreman_compute_device(this%if_ext, t, tstep, this%coef, &
162 this%nut, this%delta, this%c)
163 else
164 call vreman_compute_cpu(this%if_ext, t, tstep, this%coef, &
165 this%nut, this%delta, this%c)
166 end if
167
168 end subroutine vreman_compute
169
170end module vreman
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
Implements les_model_t.
Definition les_model.f90:35
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:65
integer, parameter, public log_size
Definition log.f90:42
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements the CPU kernel for the vreman_t type.
subroutine, public vreman_compute_cpu(if_ext, t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the CPU.
Implements the device kernel for the vreman_t type.
subroutine, public vreman_compute_device(if_ext, t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the device.
Implements vreman_t.
Definition vreman.f90:35
subroutine vreman_init_from_components(this, fluid, c, nut_name, delta_type, if_ext)
Constructor from components.
Definition vreman.f90:110
subroutine vreman_free(this)
Destructor for the les_model_t (base) class.
Definition vreman.f90:126
subroutine vreman_compute(this, t, tstep)
Compute eddy viscosity.
Definition vreman.f90:135
subroutine vreman_init(this, fluid, json)
Constructor.
Definition vreman.f90:72
Base type of all fluid formulations.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:65
Implements the Vreman LES model.
Definition vreman.f90:52