1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "fieldAverage.H"
28 #include "volFields.H"
29 #include "dictionary.H"
34 #include "fieldAverageItem.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(fieldAverage, 0);
43 const Foam::word Foam::fieldAverage::EXT_MEAN = "Mean";
44 const Foam::word Foam::fieldAverage::EXT_PRIME2MEAN = "Prime2Mean";
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 void Foam::fieldAverage::checkoutFields(const wordList& fieldNames) const
53 if (fieldNames[i] != word::null)
55 obr_.checkOut(*obr_[fieldNames[i]]);
61 void Foam::fieldAverage::resetLists(const label nItems)
63 checkoutFields(meanScalarFields_);
64 meanScalarFields_.clear();
65 meanScalarFields_.setSize(nItems);
67 checkoutFields(meanVectorFields_);
68 meanVectorFields_.clear();
69 meanVectorFields_.setSize(nItems);
71 checkoutFields(meanSphericalTensorFields_);
72 meanSphericalTensorFields_.clear();
73 meanSphericalTensorFields_.setSize(nItems);
75 checkoutFields(meanSymmTensorFields_);
76 meanSymmTensorFields_.clear();
77 meanSymmTensorFields_.setSize(nItems);
79 checkoutFields(meanTensorFields_);
80 meanTensorFields_.clear();
81 meanTensorFields_.setSize(nItems);
83 checkoutFields(prime2MeanScalarFields_);
84 prime2MeanScalarFields_.clear();
85 prime2MeanScalarFields_.setSize(nItems);
87 checkoutFields(prime2MeanSymmTensorFields_);
88 prime2MeanSymmTensorFields_.clear();
89 prime2MeanSymmTensorFields_.setSize(nItems);
92 totalIter_.setSize(nItems, 1);
95 totalTime_.setSize(nItems, obr_.time().deltaT().value());
99 void Foam::fieldAverage::initialise()
101 // Add mean fields to the field lists
104 const word& fieldName = faItems_[i].fieldName();
105 if (obr_.foundObject<volScalarField>(fieldName))
107 addMeanField<scalar>(i, meanScalarFields_);
109 else if (obr_.foundObject<volVectorField>(fieldName))
111 addMeanField<vector>(i, meanVectorFields_);
113 else if (obr_.foundObject<volSphericalTensorField>(fieldName))
115 addMeanField<sphericalTensor>(i, meanSphericalTensorFields_);
117 else if (obr_.foundObject<volSymmTensorField>(fieldName))
119 addMeanField<symmTensor>(i, meanSymmTensorFields_);
121 else if (obr_.foundObject<volTensorField>(fieldName))
123 addMeanField<tensor>(i, meanTensorFields_);
127 FatalErrorIn("Foam::fieldAverage::initialise()")
128 << "Requested field " << faItems_[i].fieldName()
129 << " does not exist in the database" << nl
134 // Add prime-squared mean fields to the field lists
137 if (faItems_[i].prime2Mean())
139 const word& fieldName = faItems_[i].fieldName();
140 if (!faItems_[i].mean())
142 FatalErrorIn("Foam::fieldAverage::initialise()")
143 << "To calculate the prime-squared average, the "
144 << "mean average must also be selected for field "
145 << fieldName << nl << exit(FatalError);
148 if (obr_.foundObject<volScalarField>(fieldName))
150 addPrime2MeanField<scalar, scalar>
154 prime2MeanScalarFields_
157 else if (obr_.foundObject<volVectorField>(fieldName))
159 addPrime2MeanField<vector, symmTensor>
163 prime2MeanSymmTensorFields_
168 FatalErrorIn("Foam::fieldAverage::initialise()")
169 << "prime2Mean average can only be applied to "
170 << "volScalarFields and volVectorFields"
171 << nl << " Field: " << fieldName << nl
179 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
181 Foam::fieldAverage::fieldAverage
184 const objectRegistry& obr,
185 const dictionary& dict,
186 const bool loadFromFiles
192 faItems_(dict.lookup("fields")),
193 meanScalarFields_(faItems_.size()),
194 meanVectorFields_(faItems_.size()),
195 meanSphericalTensorFields_(faItems_.size()),
196 meanSymmTensorFields_(faItems_.size()),
197 meanTensorFields_(faItems_.size()),
198 prime2MeanScalarFields_(faItems_.size()),
199 prime2MeanSymmTensorFields_(faItems_.size()),
200 totalIter_(faItems_.size(), 1),
201 totalTime_(faItems_.size(), obr_.time().deltaT().value())
203 // Check if the available mesh is an fvMesh otherise deactivate
204 if (!isA<fvMesh>(obr_))
209 "fieldAverage::fieldAverage\n"
212 "const objectRegistry&,\n"
213 "const dictionary&,\n"
216 ) << "No fvMesh available, deactivating."
224 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
226 Foam::fieldAverage::~fieldAverage()
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
232 void Foam::fieldAverage::read(const dictionary& dict)
237 faItems_ = List<fieldAverageItem>(dict.lookup("fields"));
239 resetLists(faItems_.size());
243 readAveragingProperties();
248 void Foam::fieldAverage::write()
254 if (obr_.time().outputTime())
258 writeAveragingProperties();
264 void Foam::fieldAverage::calcAverages()
266 Info<< "Calculating averages" << nl << endl;
270 totalTime_[i] += obr_.time().deltaT().value();
273 addMeanSqrToPrime2Mean<scalar, scalar>
276 prime2MeanScalarFields_
278 addMeanSqrToPrime2Mean<vector, symmTensor>
281 prime2MeanSymmTensorFields_
284 calculateMeanFields<scalar>(meanScalarFields_);
285 calculateMeanFields<vector>(meanVectorFields_);
286 calculateMeanFields<sphericalTensor>(meanSphericalTensorFields_);
287 calculateMeanFields<symmTensor>(meanSymmTensorFields_);
288 calculateMeanFields<tensor>(meanTensorFields_);
290 calculatePrime2MeanFields<scalar, scalar>
293 prime2MeanScalarFields_
295 calculatePrime2MeanFields<vector, symmTensor>
298 prime2MeanSymmTensorFields_
303 void Foam::fieldAverage::writeAverages() const
305 writeFieldList<scalar>(meanScalarFields_);
306 writeFieldList<vector>(meanVectorFields_);
307 writeFieldList<sphericalTensor>(meanSphericalTensorFields_);
308 writeFieldList<symmTensor>(meanSymmTensorFields_);
309 writeFieldList<tensor>(meanTensorFields_);
311 writeFieldList<scalar>(prime2MeanScalarFields_);
312 writeFieldList<symmTensor>(prime2MeanSymmTensorFields_);
316 void Foam::fieldAverage::writeAveragingProperties() const
318 IOdictionary propsDict
322 "fieldAveragingProperties",
323 obr_.time().timeName(),
334 const word& fieldName = faItems_[i].fieldName();
335 propsDict.add(fieldName, dictionary());
336 propsDict.subDict(fieldName).add("totalIter", totalIter_[i]);
337 propsDict.subDict(fieldName).add("totalTime", totalTime_[i]);
340 propsDict.regIOobject::write();
344 void Foam::fieldAverage::readAveragingProperties()
348 "fieldAveragingProperties",
349 obr_.time().timeName(),
359 IOdictionary propsDict(io);
363 const word& fieldName = faItems_[i].fieldName();
364 if (propsDict.found(fieldName))
366 const dictionary& fieldDict = propsDict.subDict(fieldName);
368 totalIter_[i] = readLabel(fieldDict.lookup("totalIter"));
369 totalTime_[i] = readScalar(fieldDict.lookup("totalTime"));
376 void Foam::fieldAverage::updateMesh(const mapPolyMesh&)
382 void Foam::fieldAverage::movePoints(const pointField&)
388 // ************************************************************************* //