PaSR.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) 2011-2015 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 #include "fvmSup.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const word& modelType,
35  const fvMesh& mesh,
36  const word& phaseName
37 )
38 :
39  laminar<Type>(modelType, mesh, phaseName),
40  Cmix_(readScalar(this->coeffs().lookup("Cmix"))),
41  turbulentReaction_(this->coeffs().lookup("turbulentReaction")),
42  kappa_
43  (
44  IOobject
45  (
46  IOobject::groupName("PaSR:kappa", phaseName),
47  mesh.time().timeName(),
48  mesh,
49  IOobject::NO_READ,
50  IOobject::AUTO_WRITE
51  ),
52  mesh,
53  dimensionedScalar("kappa", dimless, 0.0)
54  )
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
60 template<class Type>
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
67 template<class Type>
69 {
70  if (this->active())
71  {
73 
74  if (turbulentReaction_)
75  {
76  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
77  const volScalarField& epsilon = tepsilon();
78  tmp<volScalarField> tmuEff(this->turbulence().muEff());
79  const volScalarField& muEff = tmuEff();
80  tmp<volScalarField> ttc(this->tc());
81  const volScalarField& tc = ttc();
82  tmp<volScalarField> trho(this->rho());
83  const volScalarField& rho = trho();
84 
85  forAll(epsilon, i)
86  {
87  scalar tk =
88  Cmix_*sqrt(max(muEff[i]/rho[i]/(epsilon[i] + SMALL), 0));
89 
90  if (tk > SMALL)
91  {
92  kappa_[i] = tc[i]/(tc[i] + tk);
93  }
94  else
95  {
96  kappa_[i] = 1.0;
97  }
98  }
99  }
100  else
101  {
102  kappa_ = 1.0;
103  }
104  }
105 }
106 
107 
108 template<class Type>
111 {
112  return kappa_*laminar<Type>::R(Y);
113 }
114 
115 
116 template<class Type>
119 {
120  return tmp<volScalarField>
121  (
122  new volScalarField
123  (
124  IOobject::groupName("PaSR:dQ", this->phaseName_),
125  kappa_*laminar<Type>::dQ()
126  )
127  );
128 }
129 
130 
131 template<class Type>
134 {
135  return tmp<volScalarField>
136  (
137  new volScalarField
138  (
139  IOobject::groupName("PaSR:Sh", this->phaseName_),
140  kappa_*laminar<Type>::Sh()
141  )
142  );
143 }
144 
145 
146 template<class Type>
148 {
149  if (laminar<Type>::read())
150  {
151  this->coeffs().lookup("Cmix") >> Cmix_;
152  this->coeffs().lookup("turbulentReaction") >> turbulentReaction_;
153  return true;
154  }
155  else
156  {
157  return false;
158  }
159 }
160 
161 
162 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
#define readScalar
Definition: doubleScalar.C:38
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool read()
Update properties from given dictionary.
Definition: PaSR.C:147
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > trho
Partially stirred reactor combustion model. The model calculates a finite rate, based on both turbule...
Definition: PaSR.H:54
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: PaSR.C:110
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Laminar combustion model.
Definition: laminar.H:49
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil)-psi)*pSat, rhoMin);#1"/home/ubuntu/OpenFOAM-4.1/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:72
virtual tmp< volScalarField > Sh() const
Return source for enthalpy equation [kg/m/s3].
Definition: PaSR.C:133
virtual void correct()
Correct combustion rate.
Definition: PaSR.C:68
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
scalar epsilon
virtual ~PaSR()
Destructor.
Definition: PaSR.C:61
A class for managing temporary objects.
Definition: PtrList.H:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Calculate the matrix for implicit and explicit sources.
virtual tmp< volScalarField > dQ() const
Heat release rate calculated from fuel consumption rate matrix.
Definition: PaSR.C:118