dsmcFields.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-2013 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 "dsmcFields.H"
27 #include "volFields.H"
28 #include "dictionary.H"
29 #include "dsmcCloud.H"
30 
31 #include "constants.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 defineTypeNameAndDebug(dsmcFields, 0);
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::dsmcFields::dsmcFields
46 (
47  const word& name,
48  const objectRegistry& obr,
49  const dictionary& dict,
50  const bool loadFromFiles
51 )
52 :
53  name_(name),
54  obr_(obr),
55  active_(true)
56 {
57  // Check if the available mesh is an fvMesh, otherwise deactivate
58  if (!isA<fvMesh>(obr_))
59  {
60  active_ = false;
61  WarningIn
62  (
63  "dsmcFields::dsmcFields"
64  "("
65  "const word&, "
66  "const objectRegistry&, "
67  "const dictionary&, "
68  "const bool"
69  ")"
70  ) << "No fvMesh available, deactivating " << name_ << nl
71  << endl;
72  }
73 
74  read(dict);
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
87 {
88  if (active_)
89  {
90 
91  }
92 }
93 
94 
96 {
97  // Do nothing - only valid on write
98 }
99 
100 
102 {
103  // Do nothing - only valid on write
104 }
105 
106 
108 {
109  // Do nothing - only valid on write
110 }
111 
112 
114 {
115  if (active_)
116  {
117  word rhoNMeanName = "rhoNMean";
118  word rhoMMeanName = "rhoMMean";
119  word momentumMeanName = "momentumMean";
120  word linearKEMeanName = "linearKEMean";
121  word internalEMeanName = "internalEMean";
122  word iDofMeanName = "iDofMean";
123  word fDMeanName = "fDMean";
124 
125  const volScalarField& rhoNMean = obr_.lookupObject<volScalarField>
126  (
127  rhoNMeanName
128  );
129 
130  const volScalarField& rhoMMean = obr_.lookupObject<volScalarField>
131  (
132  rhoMMeanName
133  );
134 
135  const volVectorField& momentumMean = obr_.lookupObject<volVectorField>
136  (
137  momentumMeanName
138  );
139 
140  const volScalarField& linearKEMean = obr_.lookupObject<volScalarField>
141  (
142  linearKEMeanName
143  );
144 
145  const volScalarField& internalEMean = obr_.lookupObject<volScalarField>
146  (
147  internalEMeanName
148  );
149 
150  const volScalarField& iDofMean = obr_.lookupObject<volScalarField>
151  (
152  iDofMeanName
153  );
154 
155  const volVectorField& fDMean = obr_.lookupObject<volVectorField>
156  (
157  fDMeanName
158  );
159 
160  if (min(mag(rhoNMean)).value() > VSMALL)
161  {
162  Info<< "Calculating dsmcFields." << endl;
163 
164  Info<< " Calculating UMean field." << endl;
166  (
167  IOobject
168  (
169  "UMean",
170  obr_.time().timeName(),
171  obr_,
173  ),
174  momentumMean/rhoMMean
175  );
176 
177  Info<< " Calculating translationalT field." << endl;
178  volScalarField translationalT
179  (
180  IOobject
181  (
182  "translationalT",
183  obr_.time().timeName(),
184  obr_,
186  ),
187 
188  2.0/(3.0*physicoChemical::k.value()*rhoNMean)
189  *(linearKEMean - 0.5*rhoMMean*(UMean & UMean))
190  );
191 
192  Info<< " Calculating internalT field." << endl;
193  volScalarField internalT
194  (
195  IOobject
196  (
197  "internalT",
198  obr_.time().timeName(),
199  obr_,
201  ),
202  (2.0/physicoChemical::k.value())*(internalEMean/iDofMean)
203  );
204 
205  Info<< " Calculating overallT field." << endl;
206  volScalarField overallT
207  (
208  IOobject
209  (
210  "overallT",
211  obr_.time().timeName(),
212  obr_,
214  ),
215  2.0/(physicoChemical::k.value()*(3.0*rhoNMean + iDofMean))
216  *(linearKEMean - 0.5*rhoMMean*(UMean & UMean) + internalEMean)
217  );
218 
219  Info<< " Calculating pressure field." << endl;
221  (
222  IOobject
223  (
224  "p",
225  obr_.time().timeName(),
226  obr_,
228  ),
229  physicoChemical::k.value()*rhoNMean*translationalT
230  );
231 
232  const fvMesh& mesh = fDMean.mesh();
233 
234  forAll(mesh.boundaryMesh(), i)
235  {
236  const polyPatch& patch = mesh.boundaryMesh()[i];
237 
238  if (isA<wallPolyPatch>(patch))
239  {
240  p.boundaryField()[i] =
241  fDMean.boundaryField()[i]
242  & (patch.faceAreas()/mag(patch.faceAreas()));
243  }
244  }
245 
246  Info<< " mag(UMean) max/min : "
247  << max(mag(UMean)).value() << " "
248  << min(mag(UMean)).value() << endl;
249 
250  Info<< " translationalT max/min : "
251  << max(translationalT).value() << " "
252  << min(translationalT).value() << endl;
253 
254  Info<< " internalT max/min : "
255  << max(internalT).value() << " "
256  << min(internalT).value() << endl;
257 
258  Info<< " overallT max/min : "
259  << max(overallT).value() << " "
260  << min(overallT).value() << endl;
261 
262  Info<< " p max/min : "
263  << max(p).value() << " "
264  << min(p).value() << endl;
265 
266  UMean.write();
267 
268  translationalT.write();
269 
270  internalT.write();
271 
272  overallT.write();
273 
274  p.write();
275 
276  Info<< "dsmcFields written." << nl << endl;
277  }
278  else
279  {
280  Info<< "Small value (" << min(mag(rhoNMean))
281  << ") found in rhoNMean field. "
282  << "Not calculating dsmcFields to avoid division by zero."
283  << endl;
284  }
285  }
286 }
287 
288 
289 // ************************************************************************* //
virtual void read(const dictionary &)
Read the dsmcFields data.
Definition: dsmcFields.C:86
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:290
dimensioned< scalar > mag(const dimensioned< Type > &)
Collection of constants.
A class for handling words, derived from string.
Definition: word.H:59
messageStream Info
virtual void execute()
Execute, currently does nothing.
Definition: dsmcFields.C:95
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void write()
Calculate the dsmcFields and write.
Definition: dsmcFields.C:113
Namespace for OpenFOAM.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual bool write() const
Write using setting from DB.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define WarningIn(functionName)
Report a warning using Foam::Warning.
volScalarField & p
Definition: createFields.H:51
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
virtual Ostream & write(const token &)=0
Write next token to stream.
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: dsmcFields.C:107
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: dsmcFields.C:101
Registry of regIOobjects.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
virtual ~dsmcFields()
Destructor.
Definition: dsmcFields.C:80
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volVectorField UMean(UMeanHeader, mesh)
const dimensionedScalar k
Boltzmann constant.
defineTypeNameAndDebug(combustionModel, 0)