All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
transport.H
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) 2011-2021 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 Class
25  Foam::XiModels::transport
26 
27 Description
28  Simple transport model for Xi based on Gulders correlation
29  with a linear correction function to give a plausible profile for Xi.
30  See report TR/HGW/10 for details on the Weller two equations model.
31  See \link XiModel.H \endlink for more details on flame wrinkling modelling.
32 
33 SourceFiles
34  transport.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef transport_H
39 #define transport_H
40 
41 #include "XiModel.H"
42 #include "XiEqModel.H"
43 #include "XiGModel.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace XiModels
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class transport Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class transport
57 :
58  public XiModel
59 {
60  // Private Data
61 
62  scalar XiShapeCoef;
63 
64  autoPtr<XiEqModel> XiEqModel_;
65  autoPtr<XiGModel> XiGModel_;
66 
67 
68 public:
69 
70  //- Runtime type information
71  TypeName("transport");
72 
73 
74  // Constructors
75 
76  //- Construct from components
77  transport
78  (
79  const dictionary& XiProperties,
82  const volScalarField& Su,
83  const volScalarField& rho,
84  const volScalarField& b,
85  const surfaceScalarField& phi
86  );
87 
88  //- Disallow default bitwise copy construction
89  transport(const transport&);
90 
91 
92  //- Destructor
93  virtual ~transport();
94 
95 
96  // Member Functions
97 
98  //- Return the flame diffusivity
99  virtual tmp<volScalarField> Db() const;
100 
101  //- Add Xi to the multivariateSurfaceInterpolationScheme table
102  virtual void addXi
103  (
105  )
106  {
107  fields.add(Xi_);
108  }
109 
110  //- Correct the flame-wrinkling Xi
111  virtual void correct()
112  {
114  }
115 
116  //- Correct the flame-wrinkling Xi using the given convection scheme
118 
119  //- Update properties from given dictionary
120  virtual bool read(const dictionary& XiProperties);
121 
122  //- Write fields of the XiEq model
123  virtual void writeFields()
124  {
125  XiEqModel_().writeFields();
126  }
127 
128 
129  // Member Operators
130 
131  //- Disallow default bitwise assignment
132  void operator=(const transport&) = delete;
133 };
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace XiModels
139 } // End namespace Foam
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 #endif
144 
145 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for de...
Definition: XiModel.H:108
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
volScalarField Xi_
Flame wrinkling field.
Definition: XiModel.H:125
virtual void addXi(multivariateSurfaceInterpolationScheme< scalar >::fieldTable &fields)
Add Xi to the multivariateSurfaceInterpolationScheme table.
Definition: transport.H:102
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Abstract base class for convection schemes.
void operator=(const transport &)=delete
Disallow default bitwise assignment.
phi
Definition: correctPhi.H:3
Simple transport model for Xi based on Gulders correlation with a linear correction function to give ...
Definition: transport.H:55
Base-class for combustion fluid thermodynamic properties based on compressibility.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
TypeName("transport")
Runtime type information.
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
virtual void correct()
Correct the flame-wrinkling Xi.
Definition: transport.H:110
virtual bool read(const dictionary &XiProperties)
Update properties from given dictionary.
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.schemes().div("div(phi,Yi_h)")))
virtual void writeFields()
Write fields of the XiEq model.
Definition: transport.H:122
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
transport(const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi)
Construct from components.
A class for managing temporary objects.
Definition: PtrList.H:53
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
void add(const GeometricField< Type, fvPatchField, volMesh > &f)
virtual ~transport()
Destructor.
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Namespace for OpenFOAM.