00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "fieldAverage.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <OpenFOAM/Time.H>
00029
00030 #include <fieldFunctionObjects/fieldAverageItem.H>
00031
00032
00033
00034 defineTypeNameAndDebug(Foam::fieldAverage, 0);
00035
00036 const Foam::word Foam::fieldAverage::EXT_MEAN = "Mean";
00037 const Foam::word Foam::fieldAverage::EXT_PRIME2MEAN = "Prime2Mean";
00038
00039
00040
00041
00042 void Foam::fieldAverage::resetFields(wordList& names)
00043 {
00044 forAll(names, fieldI)
00045 {
00046 if (names[fieldI].size())
00047 {
00048 obr_.checkOut(*obr_[names[fieldI]]);
00049 }
00050 }
00051
00052 names.clear();
00053 names.setSize(faItems_.size());
00054 }
00055
00056
00057 void Foam::fieldAverage::initialize()
00058 {
00059 resetFields(meanScalarFields_);
00060 resetFields(meanVectorFields_);
00061 resetFields(meanSphericalTensorFields_);
00062 resetFields(meanSymmTensorFields_);
00063 resetFields(meanTensorFields_);
00064
00065 resetFields(prime2MeanScalarFields_);
00066 resetFields(prime2MeanSymmTensorFields_);
00067
00068 totalIter_.clear();
00069 totalIter_.setSize(faItems_.size(), 1);
00070
00071 totalTime_.clear();
00072 totalTime_.setSize(faItems_.size(), obr_.time().deltaT().value());
00073
00074
00075
00076 forAll(faItems_, fieldI)
00077 {
00078 const word& fieldName = faItems_[fieldI].fieldName();
00079 if (obr_.foundObject<volScalarField>(fieldName))
00080 {
00081 addMeanField<scalar>(fieldI, meanScalarFields_);
00082 }
00083 else if (obr_.foundObject<volVectorField>(fieldName))
00084 {
00085 addMeanField<vector>(fieldI, meanVectorFields_);
00086 }
00087 else if (obr_.foundObject<volSphericalTensorField>(fieldName))
00088 {
00089 addMeanField<sphericalTensor>(fieldI, meanSphericalTensorFields_);
00090 }
00091 else if (obr_.foundObject<volSymmTensorField>(fieldName))
00092 {
00093 addMeanField<symmTensor>(fieldI, meanSymmTensorFields_);
00094 }
00095 else if (obr_.foundObject<volTensorField>(fieldName))
00096 {
00097 addMeanField<tensor>(fieldI, meanTensorFields_);
00098 }
00099 else
00100 {
00101 FatalErrorIn("Foam::fieldAverage::initialize()")
00102 << "Requested field " << faItems_[fieldI].fieldName()
00103 << " does not exist in the database" << nl
00104 << exit(FatalError);
00105 }
00106 }
00107
00108
00109 forAll(faItems_, fieldI)
00110 {
00111 if (faItems_[fieldI].prime2Mean())
00112 {
00113 const word& fieldName = faItems_[fieldI].fieldName();
00114 if (!faItems_[fieldI].mean())
00115 {
00116 FatalErrorIn("Foam::fieldAverage::initialize()")
00117 << "To calculate the prime-squared average, the "
00118 << "mean average must also be selected for field "
00119 << fieldName << nl << exit(FatalError);
00120 }
00121
00122 if (obr_.foundObject<volScalarField>(fieldName))
00123 {
00124 addPrime2MeanField<scalar, scalar>
00125 (
00126 fieldI,
00127 meanScalarFields_,
00128 prime2MeanScalarFields_
00129 );
00130 }
00131 else if (obr_.foundObject<volVectorField>(fieldName))
00132 {
00133 addPrime2MeanField<vector, symmTensor>
00134 (
00135 fieldI,
00136 meanVectorFields_,
00137 prime2MeanSymmTensorFields_
00138 );
00139 }
00140 else
00141 {
00142 FatalErrorIn("Foam::fieldAverage::initialize()")
00143 << "prime2Mean average can only be applied to "
00144 << "volScalarFields and volVectorFields"
00145 << nl << " Field: " << fieldName << nl
00146 << exit(FatalError);
00147 }
00148 }
00149 }
00150 }
00151
00152
00153 void Foam::fieldAverage::calcAverages()
00154 {
00155 const label currentTimeIndex =
00156 static_cast<const fvMesh&>(obr_).time().timeIndex();
00157
00158 if (prevTimeIndex_ == currentTimeIndex)
00159 {
00160 return;
00161 }
00162 else
00163 {
00164 prevTimeIndex_ = currentTimeIndex;
00165 }
00166
00167
00168 Info<< "Calculating averages" << nl << endl;
00169 forAll(faItems_, fieldI)
00170 {
00171 totalIter_[fieldI]++;
00172 totalTime_[fieldI] += obr_.time().deltaT().value();
00173 }
00174
00175 addMeanSqrToPrime2Mean<scalar, scalar>
00176 (
00177 meanScalarFields_,
00178 prime2MeanScalarFields_
00179 );
00180 addMeanSqrToPrime2Mean<vector, symmTensor>
00181 (
00182 meanVectorFields_,
00183 prime2MeanSymmTensorFields_
00184 );
00185
00186 calculateMeanFields<scalar>(meanScalarFields_);
00187 calculateMeanFields<vector>(meanVectorFields_);
00188 calculateMeanFields<sphericalTensor>(meanSphericalTensorFields_);
00189 calculateMeanFields<symmTensor>(meanSymmTensorFields_);
00190 calculateMeanFields<tensor>(meanTensorFields_);
00191
00192 calculatePrime2MeanFields<scalar, scalar>
00193 (
00194 meanScalarFields_,
00195 prime2MeanScalarFields_
00196 );
00197 calculatePrime2MeanFields<vector, symmTensor>
00198 (
00199 meanVectorFields_,
00200 prime2MeanSymmTensorFields_
00201 );
00202 }
00203
00204
00205 void Foam::fieldAverage::writeAverages() const
00206 {
00207 writeFieldList<scalar>(meanScalarFields_);
00208 writeFieldList<vector>(meanVectorFields_);
00209 writeFieldList<sphericalTensor>(meanSphericalTensorFields_);
00210 writeFieldList<symmTensor>(meanSymmTensorFields_);
00211 writeFieldList<tensor>(meanTensorFields_);
00212
00213 writeFieldList<scalar>(prime2MeanScalarFields_);
00214 writeFieldList<symmTensor>(prime2MeanSymmTensorFields_);
00215 }
00216
00217
00218 void Foam::fieldAverage::writeAveragingProperties() const
00219 {
00220 IOdictionary propsDict
00221 (
00222 IOobject
00223 (
00224 "fieldAveragingProperties",
00225 obr_.time().timeName(),
00226 "uniform",
00227 obr_,
00228 IOobject::NO_READ,
00229 IOobject::NO_WRITE,
00230 false
00231 )
00232 );
00233
00234 forAll(faItems_, fieldI)
00235 {
00236 const word& fieldName = faItems_[fieldI].fieldName();
00237 propsDict.add(fieldName, dictionary());
00238 propsDict.subDict(fieldName).add("totalIter", totalIter_[fieldI]);
00239 propsDict.subDict(fieldName).add("totalTime", totalTime_[fieldI]);
00240 }
00241
00242 propsDict.regIOobject::write();
00243 }
00244
00245
00246 void Foam::fieldAverage::readAveragingProperties()
00247 {
00248 if (cleanRestart_)
00249 {
00250 Info<< "fieldAverage: starting averaging at time "
00251 << obr_.time().timeName() << nl << endl;
00252 }
00253 else
00254 {
00255 IOobject propsDictHeader
00256 (
00257 "fieldAveragingProperties",
00258 obr_.time().timeName(),
00259 "uniform",
00260 obr_,
00261 IOobject::MUST_READ,
00262 IOobject::NO_WRITE,
00263 false
00264 );
00265
00266 if (!propsDictHeader.headerOk())
00267 {
00268 Info<< "fieldAverage: starting averaging at time "
00269 << obr_.time().timeName() << nl << endl;
00270 return;
00271 }
00272
00273 IOdictionary propsDict(propsDictHeader);
00274
00275 Info<< "fieldAverage: restarting averaging for fields:" << endl;
00276 forAll(faItems_, fieldI)
00277 {
00278 const word& fieldName = faItems_[fieldI].fieldName();
00279 if (propsDict.found(fieldName))
00280 {
00281 dictionary fieldDict(propsDict.subDict(fieldName));
00282
00283 totalIter_[fieldI] = readLabel(fieldDict.lookup("totalIter"));
00284 totalTime_[fieldI] = readScalar(fieldDict.lookup("totalTime"));
00285 Info<< " " << fieldName
00286 << " iters = " << totalIter_[fieldI]
00287 << " time = " << totalTime_[fieldI] << endl;
00288 }
00289 }
00290 Info<< endl;
00291 }
00292 }
00293
00294
00295
00296
00297 Foam::fieldAverage::fieldAverage
00298 (
00299 const word& name,
00300 const objectRegistry& obr,
00301 const dictionary& dict,
00302 const bool loadFromFiles
00303 )
00304 :
00305 name_(name),
00306 obr_(obr),
00307 active_(true),
00308 prevTimeIndex_(-1),
00309 cleanRestart_(false),
00310 resetOnOutput_(false),
00311 faItems_(),
00312 meanScalarFields_(),
00313 meanVectorFields_(),
00314 meanSphericalTensorFields_(),
00315 meanSymmTensorFields_(),
00316 meanTensorFields_(),
00317 prime2MeanScalarFields_(),
00318 prime2MeanSymmTensorFields_(),
00319 totalIter_(),
00320 totalTime_()
00321 {
00322
00323 if (isA<fvMesh>(obr_))
00324 {
00325 read(dict);
00326 }
00327 else
00328 {
00329 active_ = false;
00330 WarningIn
00331 (
00332 "fieldAverage::fieldAverage\n"
00333 "(\n"
00334 "const word&,\n"
00335 "const objectRegistry&,\n"
00336 "const dictionary&,\n"
00337 "const bool\n"
00338 ")"
00339 ) << "No fvMesh available, deactivating."
00340 << nl << endl;
00341 }
00342 }
00343
00344
00345
00346
00347 Foam::fieldAverage::~fieldAverage()
00348 {}
00349
00350
00351
00352
00353 void Foam::fieldAverage::read(const dictionary& dict)
00354 {
00355 if (active_)
00356 {
00357 dict.readIfPresent("cleanRestart", cleanRestart_);
00358 dict.readIfPresent("resetOnOutput", resetOnOutput_);
00359 dict.lookup("fields") >> faItems_;
00360
00361 initialize();
00362 readAveragingProperties();
00363
00364
00365 prevTimeIndex_ = -1;
00366 }
00367 }
00368
00369
00370 void Foam::fieldAverage::execute()
00371 {
00372 if (active_)
00373 {
00374 calcAverages();
00375 }
00376 }
00377
00378
00379 void Foam::fieldAverage::end()
00380 {
00381 }
00382
00383
00384 void Foam::fieldAverage::write()
00385 {
00386 if (active_)
00387 {
00388 calcAverages();
00389 writeAverages();
00390 writeAveragingProperties();
00391
00392 if (resetOnOutput_)
00393 {
00394 Info<< "fieldAverage: restarting averaging at time "
00395 << obr_.time().timeName() << nl << endl;
00396
00397 initialize();
00398
00399
00400 prevTimeIndex_ = -1;
00401 }
00402 }
00403 }
00404
00405
00406 void Foam::fieldAverage::updateMesh(const mapPolyMesh&)
00407 {
00408
00409 }
00410
00411
00412 void Foam::fieldAverage::movePoints(const pointField&)
00413 {
00414
00415 }
00416
00417
00418