Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
turbulenceModels
compressible
RAS
derivedFvPatchFields
wallFunctions
kqRWallFunctions
kqRWallFunction
kqRWallFunctionFvPatchField.C
Go to the documentation of this file.
1
/*---------------------------------------------------------------------------*\
2
========= |
3
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4
\\ / O peration |
5
\\ / A nd | Copyright (C) 2008-2010 OpenCFD Ltd.
6
\\/ M anipulation |
7
-------------------------------------------------------------------------------
8
License
9
This file is part of OpenFOAM.
10
11
OpenFOAM is free software: you can redistribute it and/or modify it
12
under the terms of the GNU General Public License as published by
13
the Free Software Foundation, either version 3 of the License, or
14
(at your option) any later version.
15
16
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19
for more details.
20
21
You should have received a copy of the GNU General Public License
22
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23
24
\*---------------------------------------------------------------------------*/
25
26
#include "
kqRWallFunctionFvPatchField.H
"
27
#include <
finiteVolume/fvPatchFieldMapper.H
>
28
#include <
OpenFOAM/addToRunTimeSelectionTable.H
>
29
#include <
finiteVolume/wallFvPatch.H
>
30
31
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33
namespace
Foam
34
{
35
namespace
compressible
36
{
37
namespace
RASModels
38
{
39
40
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
42
template
<
class
Type>
43
void
kqRWallFunctionFvPatchField<Type>::checkType()
44
{
45
if
(!isA<wallFvPatch>(this->patch()))
46
{
47
FatalErrorIn
(
"kqRWallFunctionFvPatchField::checkType()"
)
48
<<
"Invalid wall function specification"
<<
nl
49
<<
" Patch type for patch "
<< this->patch().name()
50
<<
" must be wall"
<<
nl
51
<<
" Current patch type is "
<< this->patch().type() <<
nl
<<
endl
52
<<
abort
(
FatalError
);
53
}
54
}
55
56
57
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58
59
template
<
class
Type>
60
kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
61
(
62
const
fvPatch
&
p
,
63
const
DimensionedField<Type, volMesh>
& iF
64
)
65
:
66
zeroGradientFvPatchField<Type>
(
p
, iF)
67
{
68
checkType();
69
}
70
71
72
template
<
class
Type>
73
kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
74
(
75
const
kqRWallFunctionFvPatchField
& ptf,
76
const
fvPatch
& p,
77
const
DimensionedField<Type, volMesh>
& iF,
78
const
fvPatchFieldMapper
& mapper
79
)
80
:
81
zeroGradientFvPatchField<Type>
(ptf,
p
, iF, mapper)
82
{
83
checkType();
84
}
85
86
87
template
<
class
Type>
88
kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
89
(
90
const
fvPatch
& p,
91
const
DimensionedField<Type, volMesh>
& iF,
92
const
dictionary
& dict
93
)
94
:
95
zeroGradientFvPatchField<Type>
(
p
, iF, dict)
96
{
97
checkType();
98
}
99
100
101
template
<
class
Type>
102
kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
103
(
104
const
kqRWallFunctionFvPatchField
& tkqrwfpf
105
)
106
:
107
zeroGradientFvPatchField<Type>
(tkqrwfpf)
108
{
109
checkType();
110
}
111
112
113
template
<
class
Type>
114
kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
115
(
116
const
kqRWallFunctionFvPatchField
& tkqrwfpf,
117
const
DimensionedField<Type, volMesh>
& iF
118
)
119
:
120
zeroGradientFvPatchField<Type>
(tkqrwfpf, iF)
121
{
122
checkType();
123
}
124
125
126
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127
128
template
<
class
Type>
129
void
kqRWallFunctionFvPatchField<Type>::evaluate
130
(
131
const
Pstream::commsTypes
commsType
132
)
133
{
134
zeroGradientFvPatchField<Type>::evaluate
(commsType);
135
}
136
137
138
template
<
class
Type>
139
void
kqRWallFunctionFvPatchField<Type>::write
(
Ostream
& os)
const
140
{
141
zeroGradientFvPatchField<Type>::write
(os);
142
this->writeEntry(
"value"
, os);
143
}
144
145
146
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147
148
}
// End namespace RASModels
149
}
// End namespace compressible
150
}
// End namespace Foam
151
152
// ************************ vim: set sw=4 sts=4 et: ************************ //