PaSR.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) 2011-2020 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 "PaSR.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ReactionThermo>
32 (
33  const word& modelType,
34  const ReactionThermo& thermo,
36  const word& combustionProperties
37 )
38 :
39  laminar<ReactionThermo>(modelType, thermo, turb, combustionProperties),
40  Cmix_(this->coeffs().template lookup<scalar>("Cmix")),
41  kappa_
42  (
43  IOobject
44  (
45  thermo.phasePropertyName(typeName + ":kappa"),
46  this->mesh().time().timeName(),
47  this->mesh(),
48  IOobject::NO_READ,
49  IOobject::AUTO_WRITE
50  ),
51  this->mesh(),
53  )
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
59 template<class ReactionThermo>
61 {}
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
66 template<class ReactionThermo>
68 {
70 
71  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
72  const scalarField& epsilon = tepsilon();
73 
74  tmp<volScalarField> tmuEff(this->turbulence().muEff());
75  const scalarField& muEff = tmuEff();
76 
77  tmp<volScalarField> ttc(this->tc());
78  const scalarField& tc = ttc();
79 
80  tmp<volScalarField> trho(this->rho());
81  const scalarField& rho = trho();
82 
83  forAll(epsilon, i)
84  {
85  const scalar tk =
86  Cmix_*sqrt(max(muEff[i]/rho[i]/(epsilon[i] + small), 0));
87 
88  if (tk > small)
89  {
90  kappa_[i] = tc[i]/(tc[i] + tk);
91  }
92  else
93  {
94  kappa_[i] = 1.0;
95  }
96  }
97 }
98 
99 
100 template<class ReactionThermo>
103 {
104  return kappa_*laminar<ReactionThermo>::R(Y);
105 }
106 
107 
108 template<class ReactionThermo>
111 {
112  return volScalarField::New
113  (
114  this->thermo().phasePropertyName(typeName + ":Qdot"),
116  );
117 }
118 
119 
120 template<class ReactionThermo>
122 {
124  {
125  this->coeffs().lookup("Cmix") >> Cmix_;
126  return true;
127  }
128  else
129  {
130  return false;
131  }
132 }
133 
134 
135 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: PaSR.C:102
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
Definition: PaSR.C:110
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
PaSR(const word &modelType, const ReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: PaSR.C:32
virtual bool read()
Update properties from given dictionary.
Definition: PaSR.C:121
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > trho
rhoReactionThermo & thermo
Definition: createFields.H:28
Laminar combustion model.
Definition: laminar.H:51
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
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
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
virtual void correct()
Correct combustion rate.
Definition: PaSR.C:67
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define R(A, B, C, D, E, F, K, M)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
scalar epsilon
virtual ~PaSR()
Destructor.
Definition: PaSR.C:60
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Abstract base class for turbulence models (RAS, LES and laminar).