Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
sigma.f90
Go to the documentation of this file.
1! Copyright (c) 2024, 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 sigma
36 use num_types, only : rp
37 use field, only : field_t
38 use les_model, only : les_model_t
39 use dofmap , only : dofmap_t
41 use json_module, only : json_file
42 use utils, only : neko_error
46 use coefs, only : coef_t
47 use logger, only : log_size, neko_log
48 implicit none
49 private
50
53 type, public, extends(les_model_t) :: sigma_t
55 real(kind=rp) :: c
56 contains
58 procedure, pass(this) :: init => sigma_init
60 procedure, pass(this) :: init_from_components => sigma_init_from_components
62 procedure, pass(this) :: free => sigma_free
64 procedure, pass(this) :: compute => sigma_compute
65 end type sigma_t
66
67contains
72 subroutine sigma_init(this, dofmap, coef, json)
73 class(sigma_t), intent(inout) :: this
74 type(dofmap_t), intent(in) :: dofmap
75 type(coef_t), intent(in) :: coef
76 type(json_file), intent(inout) :: json
77 character(len=:), allocatable :: nut_name
78 real(kind=rp) :: c
79 character(len=:), allocatable :: delta_type
80 character(len=LOG_SIZE) :: log_buf
81
82 call json_get_or_default(json, "nut_field", nut_name, "nut")
83 call json_get_or_default(json, "delta_type", delta_type, "pointwise")
84 ! Based on C = 1.35 as default values
85 call json_get_or_default(json, "c", c, 1.35_rp)
86
87 call neko_log%section('LES model')
88 write(log_buf, '(A)') 'Model : Sigma'
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 call neko_log%end_section()
95
96 call sigma_init_from_components(this, dofmap, coef, c, nut_name, delta_type)
97 end subroutine sigma_init
98
104 subroutine sigma_init_from_components(this, dofmap, coef, c, nut_name, &
105 delta_type)
106 class(sigma_t), intent(inout) :: this
107 type(dofmap_t), intent(in) :: dofmap
108 type(coef_t), intent(in) :: coef
109 real(kind=rp) :: c
110 character(len=*), intent(in) :: nut_name
111 character(len=*), intent(in) :: delta_type
112
113 call this%free()
114
115 call this%init_base(dofmap, coef, nut_name, delta_type)
116
117 this%c = c
118
119 end subroutine sigma_init_from_components
120
122 subroutine sigma_free(this)
123 class(sigma_t), intent(inout) :: this
124
125 call this%free_base()
126 end subroutine sigma_free
127
131 subroutine sigma_compute(this, t, tstep)
132 class(sigma_t), intent(inout) :: this
133 real(kind=rp), intent(in) :: t
134 integer, intent(in) :: tstep
135
136 if (neko_bcknd_device .eq. 1) then
137 call sigma_compute_device(t, tstep, this%coef, this%nut, this%delta, &
138 this%c)
139 else
140 call sigma_compute_cpu(t, tstep, this%coef, this%nut, this%delta, &
141 this%c)
142 end if
143
144 end subroutine sigma_compute
145
146end module sigma
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.
Coefficients.
Definition coef.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
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 sigma_t type. Following Nicoud et al. "Using singular values to bui...
Definition sigma_cpu.f90:38
subroutine, public sigma_compute_cpu(t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the CPU.
Definition sigma_cpu.f90:63
Implements the device kernel for the sigma_t type.
subroutine, public sigma_compute_device(t, tstep, coef, nut, delta, c)
Compute eddy viscosity on device.
Implements sigma_t.
Definition sigma.f90:35
subroutine sigma_compute(this, t, tstep)
Compute eddy viscosity.
Definition sigma.f90:132
subroutine sigma_init_from_components(this, dofmap, coef, c, nut_name, delta_type)
Constructor from components.
Definition sigma.f90:106
subroutine sigma_init(this, dofmap, coef, json)
Constructor.
Definition sigma.f90:73
subroutine sigma_free(this)
Destructor for the les_model_t (base) class.
Definition sigma.f90:123
Utilities.
Definition utils.f90:35
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:51
Implements the Sigma LES model.
Definition sigma.f90:53