Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.1
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 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 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
72 subroutine vreman_init(this, dofmap, coef, json)
73 class(vreman_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
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
86 call vreman_init_from_components(this, dofmap, coef, c, nut_name, &
87 delta_type)
88 end subroutine vreman_init
89
95 subroutine vreman_init_from_components(this, dofmap, coef, c, nut_name, &
96 delta_type)
97 class(vreman_t), intent(inout) :: this
98 type(dofmap_t), intent(in) :: dofmap
99 type(coef_t), intent(in) :: coef
100 real(kind=rp) :: c
101 character(len=*), intent(in) :: nut_name
102 character(len=*), intent(in) :: delta_type
103
104 call this%free()
105
106 call this%init_base(dofmap, coef, nut_name, delta_type)
107 this%c = c
108
109 end subroutine vreman_init_from_components
110
112 subroutine vreman_free(this)
113 class(vreman_t), intent(inout) :: this
114
115 call this%free_base()
116 end subroutine vreman_free
117
121 subroutine vreman_compute(this, t, tstep)
122 class(vreman_t), intent(inout) :: this
123 real(kind=rp), intent(in) :: t
124 integer, intent(in) :: tstep
125
126 if (neko_bcknd_device .eq. 1) then
127 call vreman_compute_device(t, tstep, this%coef, this%nut, this%delta,&
128 this%c)
129 else
130 call vreman_compute_cpu(t, tstep, this%coef, this%nut, this%delta,&
131 this%c)
132 end if
133
134 end subroutine vreman_compute
135
136end module vreman
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
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 vreman_t type.
subroutine, public vreman_compute_cpu(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(t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the device.
Implements vreman_t.
Definition vreman.f90:35
subroutine vreman_init(this, dofmap, coef, json)
Constructor.
Definition vreman.f90:73
subroutine vreman_free(this)
Destructor for the les_model_t (base) class.
Definition vreman.f90:113
subroutine vreman_compute(this, t, tstep)
Compute eddy viscosity.
Definition vreman.f90:122
subroutine vreman_init_from_components(this, dofmap, coef, c, nut_name, delta_type)
Constructor from components.
Definition vreman.f90:97
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 Vreman LES model.
Definition vreman.f90:52