BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / postProcessing / functionObjects / utilities / dsmcFields / dsmcFields.C
blobac8b9cca10ec129323237aa279ec3e575ddc2c5f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "dsmcFields.H"
27 #include "volFields.H"
28 #include "dictionary.H"
29 #include "dsmcCloud.H"
31 #include "constants.H"
33 using namespace Foam::constant;
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 defineTypeNameAndDebug(Foam::dsmcFields, 0);
40 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
42 Foam::dsmcFields::dsmcFields
44     const word& name,
45     const objectRegistry& obr,
46     const dictionary& dict,
47     const bool loadFromFiles
50     name_(name),
51     obr_(obr),
52     active_(true)
54     // Check if the available mesh is an fvMesh, otherwise deactivate
55     if (!isA<fvMesh>(obr_))
56     {
57         active_ = false;
58         WarningIn
59         (
60             "dsmcFields::dsmcFields"
61             "("
62                 "const word&, "
63                 "const objectRegistry&, "
64                 "const dictionary&, "
65                 "const bool"
66             ")"
67         )   << "No fvMesh available, deactivating." << nl
68             << endl;
69     }
71     read(dict);
75 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
77 Foam::dsmcFields::~dsmcFields()
81 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
83 void Foam::dsmcFields::read(const dictionary& dict)
85     if (active_)
86     {
88     }
92 void Foam::dsmcFields::execute()
94     // Do nothing - only valid on write
98 void Foam::dsmcFields::end()
100     // Do nothing - only valid on write
104 void Foam::dsmcFields::write()
106     if (active_)
107     {
108         word rhoNMeanName = "rhoNMean";
109         word rhoMMeanName = "rhoMMean";
110         word momentumMeanName = "momentumMean";
111         word linearKEMeanName = "linearKEMean";
112         word internalEMeanName = "internalEMean";
113         word iDofMeanName = "iDofMean";
114         word fDMeanName = "fDMean";
116         const volScalarField& rhoNMean = obr_.lookupObject<volScalarField>
117         (
118             rhoNMeanName
119         );
121         const volScalarField& rhoMMean = obr_.lookupObject<volScalarField>
122         (
123             rhoMMeanName
124         );
126         const volVectorField& momentumMean = obr_.lookupObject<volVectorField>
127         (
128             momentumMeanName
129         );
131         const volScalarField& linearKEMean = obr_.lookupObject<volScalarField>
132         (
133             linearKEMeanName
134         );
136         const volScalarField& internalEMean = obr_.lookupObject<volScalarField>
137         (
138             internalEMeanName
139         );
141         const volScalarField& iDofMean = obr_.lookupObject<volScalarField>
142         (
143             iDofMeanName
144         );
146         const volVectorField& fDMean = obr_.lookupObject<volVectorField>
147         (
148             fDMeanName
149         );
151         if (min(mag(rhoNMean)).value() > VSMALL)
152         {
153             Info<< "Calculating dsmcFields." << endl;
155             Info<< "    Calculating UMean field." << endl;
156             volVectorField UMean
157             (
158                 IOobject
159                 (
160                     "UMean",
161                     obr_.time().timeName(),
162                     obr_,
163                     IOobject::NO_READ
164                 ),
165                 momentumMean/rhoMMean
166             );
168             Info<< "    Calculating translationalT field." << endl;
169             volScalarField translationalT
170             (
171                 IOobject
172                 (
173                     "translationalT",
174                     obr_.time().timeName(),
175                     obr_,
176                     IOobject::NO_READ
177                 ),
179                 2.0/(3.0*physicoChemical::k.value()*rhoNMean)
180                *(linearKEMean - 0.5*rhoMMean*(UMean & UMean))
181             );
183             Info<< "    Calculating internalT field." << endl;
184             volScalarField internalT
185             (
186                 IOobject
187                 (
188                     "internalT",
189                     obr_.time().timeName(),
190                     obr_,
191                     IOobject::NO_READ
192                 ),
193                 (2.0/physicoChemical::k.value())*(internalEMean/iDofMean)
194             );
196             Info<< "    Calculating overallT field." << endl;
197             volScalarField overallT
198             (
199                 IOobject
200                 (
201                     "overallT",
202                     obr_.time().timeName(),
203                     obr_,
204                     IOobject::NO_READ
205                 ),
206                 2.0/(physicoChemical::k.value()*(3.0*rhoNMean + iDofMean))
207                *(linearKEMean - 0.5*rhoMMean*(UMean & UMean) + internalEMean)
208             );
210             Info<< "    Calculating pressure field." << endl;
211             volScalarField p
212             (
213                 IOobject
214                 (
215                     "p",
216                     obr_.time().timeName(),
217                     obr_,
218                     IOobject::NO_READ
219                 ),
220                 physicoChemical::k.value()*rhoNMean*translationalT
221             );
223             const fvMesh& mesh = fDMean.mesh();
225             forAll(mesh.boundaryMesh(), i)
226             {
227                 const polyPatch& patch = mesh.boundaryMesh()[i];
229                 if (isA<wallPolyPatch>(patch))
230                 {
231                     p.boundaryField()[i] =
232                         fDMean.boundaryField()[i]
233                       & (patch.faceAreas()/mag(patch.faceAreas()));
234                 }
235             }
237             Info<< "    mag(UMean) max/min : "
238                 << max(mag(UMean)).value() << " "
239                 << min(mag(UMean)).value() << endl;
241             Info<< "    translationalT max/min : "
242                 << max(translationalT).value() << " "
243                 << min(translationalT).value() << endl;
245             Info<< "    internalT max/min : "
246                 << max(internalT).value() << " "
247                 << min(internalT).value() << endl;
249             Info<< "    overallT max/min : "
250                 << max(overallT).value() << " "
251                 << min(overallT).value() << endl;
253             Info<< "    p max/min : "
254                 << max(p).value() << " "
255                 << min(p).value() << endl;
257             UMean.write();
259             translationalT.write();
261             internalT.write();
263             overallT.write();
265             p.write();
267             Info<< "dsmcFields written." << nl << endl;
268         }
269         else
270         {
271             Info<< "Small value (" << min(mag(rhoNMean))
272                 << ") found in rhoNMean field. "
273                 << "Not calculating dsmcFields to avoid division by zero."
274                 << endl;
275         }
276     }
280 // ************************************************************************* //