Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
wale.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 wale
36 use num_types, only : rp
37 use field, only : field_t
39 use les_model, only : les_model_t
41 use json_module, only : json_file
42 use utils, only : neko_error
44 use wale_cpu, only : wale_compute_cpu
47 use logger, only : log_size, neko_log
48 implicit none
49 private
50
53 type, public, extends(les_model_t) :: wale_t
55 real(kind=rp) :: c_w
56 contains
58 procedure, pass(this) :: init => wale_init
60 procedure, pass(this) :: init_from_components => &
63 procedure, pass(this) :: free => wale_free
65 procedure, pass(this) :: compute => wale_compute
66 end type wale_t
67
68contains
72 subroutine wale_init(this, fluid, json)
73 class(wale_t), intent(inout) :: this
74 class(fluid_scheme_base_t), intent(inout), target :: fluid
75 type(json_file), intent(inout) :: json
76 character(len=:), allocatable :: nut_name
77 real(kind=rp) :: c_w
78 character(len=:), allocatable :: delta_type
79 logical :: if_ext
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 call json_get_or_default(json, "c_w", c_w, 0.55_rp)
85 call json_get_or_default(json, "extrapolation", if_ext, .false.)
86
87 call neko_log%section('LES model')
88 write(log_buf, '(A)') 'Model : Wale'
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_w : ', c_w
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 wale_init_from_components(this, fluid, c_w, nut_name, &
99 delta_type, if_ext)
100
101 end subroutine wale_init
102
109 subroutine wale_init_from_components(this, fluid, c_w, &
110 nut_name, delta_type, if_ext)
111 class(wale_t), intent(inout) :: this
112 class(fluid_scheme_base_t), intent(inout), target :: fluid
113 real(kind=rp) :: c_w
114 character(len=*), intent(in) :: nut_name
115 character(len=*), intent(in) :: delta_type
116 logical, intent(in) :: if_ext
117
118 call this%free()
119
120 call this%init_base(fluid, nut_name, delta_type, if_ext)
121 this%c_w = c_w
122
123 end subroutine wale_init_from_components
124
126 subroutine wale_free(this)
127 class(wale_t), intent(inout) :: this
128
129 call this%free_base()
130 end subroutine wale_free
131
135 subroutine wale_compute(this, t, tstep)
136 class(wale_t), intent(inout) :: this
137 real(kind=rp), intent(in) :: t
138 integer, intent(in) :: tstep
139
140 type(field_t), pointer :: u, v, w, u_e, v_e, w_e
141
142 if (this%if_ext .eqv. .true.) then
143 ! Extrapolate the velocity fields
144 associate(ulag => this%ulag, vlag => this%vlag, &
145 wlag => this%wlag, ext_bdf => this%ext_bdf)
146
147 u => neko_field_registry%get_field_by_name("u")
148 v => neko_field_registry%get_field_by_name("v")
149 w => neko_field_registry%get_field_by_name("w")
150 u_e => neko_field_registry%get_field_by_name("u_e")
151 v_e => neko_field_registry%get_field_by_name("v_e")
152 w_e => neko_field_registry%get_field_by_name("w_e")
153
154 call this%sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
155 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
156
157 end associate
158 end if
159
160 ! Compute the eddy viscosity field
161 if (neko_bcknd_device .eq. 1) then
162 call wale_compute_device(this%if_ext, t, tstep, this%coef, &
163 this%nut, this%delta, this%c_w)
164 else
165 call wale_compute_cpu(this%if_ext, t, tstep, this%coef, &
166 this%nut, this%delta, this%c_w)
167 end if
168
169 end subroutine wale_compute
170
171end module wale
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:70
integer, parameter, public log_size
Definition log.f90:46
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35
Implements the CPU kernel for the wale_t type.
Definition wale_cpu.f90:34
subroutine, public wale_compute_cpu(if_ext, t, tstep, coef, nut, delta, c_w)
Compute eddy viscosity on the CPU.
Definition wale_cpu.f90:60
Implements the device kernel for the wale_t type.
subroutine, public wale_compute_device(if_ext, t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the device.
Implements wale_t.
Definition wale.f90:35
subroutine wale_init(this, fluid, json)
Constructor.
Definition wale.f90:73
subroutine wale_compute(this, t, tstep)
Compute eddy viscosity.
Definition wale.f90:136
subroutine wale_init_from_components(this, fluid, c_w, nut_name, delta_type, if_ext)
Constructor from components.
Definition wale.f90:111
subroutine wale_free(this)
Destructor for the les_model_t (base) class.
Definition wale.f90:127
Base type of all fluid formulations.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:64
Implements the Wale LES model.
Definition wale.f90:53