ThermalDiffusivity.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2018 OpenFOAM Foundation
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 "ThermalDiffusivity.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class BasicTurbulenceModel>
32 (
33  const word& type,
34  const alphaField& alpha,
35  const volScalarField& rho,
36  const volVectorField& U,
37  const surfaceScalarField& alphaRhoPhi,
38  const surfaceScalarField& phi,
39  const transportModel& transport,
40  const word& propertiesName
41 )
42 :
43  BasicTurbulenceModel
44  (
45  type,
46  alpha,
47  rho,
48  U,
49  alphaRhoPhi,
50  phi,
51  transport,
52  propertiesName
53  )
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
58 
59 template<class BasicTurbulenceModel>
62 (
63  const alphaField& alpha,
64  const volScalarField& rho,
65  const volVectorField& U,
66  const surfaceScalarField& alphaRhoPhi,
67  const surfaceScalarField& phi,
68  const transportModel& transport,
69  const word& propertiesName
70 )
71 {
73  (
74  static_cast<ThermalDiffusivity*>(
76  (
77  alpha,
78  rho,
79  U,
80  alphaRhoPhi,
81  phi,
82  transport,
83  propertiesName
84  ).ptr())
85  );
86 }
87 
88 
89 template<class BasicTurbulenceModel>
92 (
93  const volScalarField& rho,
94  const volVectorField& U,
95  const surfaceScalarField& phi,
96  const transportModel& transport,
97  const word& propertiesName
98 )
99 {
101  (
102  static_cast<ThermalDiffusivity*>(
104  (
105  rho,
106  U,
107  phi,
108  transport,
109  propertiesName
110  ).ptr())
111  );
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 template<class BasicTurbulenceModel>
120 {
121  return volScalarField::New
122  (
123  IOobject::groupName("alphat", this->alphaRhoPhi_.group()),
124  this->mesh_,
126  );
127 }
128 
129 
130 template<class BasicTurbulenceModel>
133 (
134  const label patchi
135 ) const
136 {
137  return tmp<scalarField>
138  (
139  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
140  );
141 }
142 
143 
144 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimViscosity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
ThermalDiffusivity(const word &type, const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &trasport, const word &propertiesName)
Construct.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
A class for handling words, derived from string.
Definition: word.H:59
volScalarField scalarField(fieldObject, mesh)
const dimensionSet dimDensity
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static autoPtr< ThermalDiffusivity > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &trasportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual tmp< volScalarField > alphat() const
Turbulent thermal diffusivity for enthalpy [kg/m/s].
A class for managing temporary objects.
Definition: PtrList.H:53