Go to the documentation of this file.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 "faceSource.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <OpenFOAM/cyclicPolyPatch.H>
00029 #include <OpenFOAM/emptyPolyPatch.H>
00030 #include <OpenFOAM/processorPolyPatch.H>
00031 #include <sampling/sampledSurface.H>
00032
00033
00034
00035 namespace Foam
00036 {
00037 namespace fieldValues
00038 {
00039 defineTypeNameAndDebug(faceSource, 0);
00040 }
00041
00042 template<>
00043 const char* NamedEnum<fieldValues::faceSource::sourceType, 3>::
00044 names[] = {"faceZone", "patch", "sampledSurface"};
00045
00046 const NamedEnum<fieldValues::faceSource::sourceType, 3>
00047 fieldValues::faceSource::sourceTypeNames_;
00048
00049 template<>
00050 const char* NamedEnum<fieldValues::faceSource::operationType, 7>::
00051 names[] =
00052 {
00053 "none", "sum", "areaAverage",
00054 "areaIntegrate", "weightedAverage", "min", "max"
00055 };
00056
00057 const NamedEnum<fieldValues::faceSource::operationType, 7>
00058 fieldValues::faceSource::operationTypeNames_;
00059
00060 }
00061
00062
00063
00064
00065 void Foam::fieldValues::faceSource::setFaceZoneFaces()
00066 {
00067 label zoneId = mesh().faceZones().findZoneID(sourceName_);
00068
00069 if (zoneId < 0)
00070 {
00071 FatalErrorIn("faceSource::faceSource::setFaceZoneFaces()")
00072 << type() << " " << name_ << ": "
00073 << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
00074 << " Unknown face zone name: " << sourceName_
00075 << ". Valid face zones are: " << mesh().faceZones().names()
00076 << nl << exit(FatalError);
00077 }
00078
00079 const faceZone& fZone = mesh().faceZones()[zoneId];
00080
00081 faceId_.setSize(fZone.size());
00082 facePatchId_.setSize(fZone.size());
00083 faceSign_.setSize(fZone.size());
00084
00085 label count = 0;
00086 forAll(fZone, i)
00087 {
00088 label faceI = fZone[i];
00089
00090 label faceId = -1;
00091 label facePatchId = -1;
00092 if (mesh().isInternalFace(faceI))
00093 {
00094 faceId = faceI;
00095 facePatchId = -1;
00096 }
00097 else
00098 {
00099 facePatchId = mesh().boundaryMesh().whichPatch(faceI);
00100 const polyPatch& pp = mesh().boundaryMesh()[facePatchId];
00101 if (isA<processorPolyPatch>(pp))
00102 {
00103 if (refCast<const processorPolyPatch>(pp).owner())
00104 {
00105 faceId = pp.whichFace(faceI);
00106 }
00107 else
00108 {
00109 faceId = -1;
00110 }
00111 }
00112 else if (isA<cyclicPolyPatch>(pp))
00113 {
00114 label patchFaceI = faceI - pp.start();
00115 if (patchFaceI < pp.size()/2)
00116 {
00117 faceId = patchFaceI;
00118 }
00119 else
00120 {
00121 faceId = -1;
00122 }
00123 }
00124 else if (!isA<emptyPolyPatch>(pp))
00125 {
00126 faceId = faceI - pp.start();
00127 }
00128 else
00129 {
00130 faceId = -1;
00131 facePatchId = -1;
00132 }
00133 }
00134
00135 if (faceId >= 0)
00136 {
00137 if (fZone.flipMap()[i])
00138 {
00139 faceSign_[count] = -1;
00140 }
00141 else
00142 {
00143 faceSign_[count] = 1;
00144 }
00145 faceId_[count] = faceId;
00146 facePatchId_[count] = facePatchId;
00147 count++;
00148 }
00149 }
00150
00151 faceId_.setSize(count);
00152 facePatchId_.setSize(count);
00153 faceSign_.setSize(count);
00154 nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
00155
00156 if (debug)
00157 {
00158 Pout<< "Original face zone size = " << fZone.size()
00159 << ", new size = " << count << endl;
00160 }
00161 }
00162
00163
00164 void Foam::fieldValues::faceSource::setPatchFaces()
00165 {
00166 label patchId = mesh().boundaryMesh().findPatchID(sourceName_);
00167
00168 if (patchId < 0)
00169 {
00170 FatalErrorIn("faceSource::constructFaceAddressing()")
00171 << type() << " " << name_ << ": "
00172 << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
00173 << " Unknown patch name: " << sourceName_
00174 << ". Valid patch names are: "
00175 << mesh().boundaryMesh().names() << nl
00176 << exit(FatalError);
00177 }
00178
00179 const polyPatch& pp = mesh().boundaryMesh()[patchId];
00180
00181 label nFaces = pp.size();
00182 if (isA<cyclicPolyPatch>(pp))
00183 {
00184 nFaces /= 2;
00185 }
00186 else if (isA<emptyPolyPatch>(pp))
00187 {
00188 nFaces = 0;
00189 }
00190
00191 faceId_.setSize(nFaces);
00192 facePatchId_.setSize(nFaces);
00193 faceSign_.setSize(nFaces);
00194 nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
00195
00196 forAll(faceId_, faceI)
00197 {
00198 faceId_[faceI] = faceI;
00199 facePatchId_[faceI] = patchId;
00200 faceSign_[faceI] = 1;
00201 }
00202 }
00203
00204
00205 void Foam::fieldValues::faceSource::sampledSurfaceFaces(const dictionary& dict)
00206 {
00207 surfacePtr_ = sampledSurface::New
00208 (
00209 name_,
00210 mesh(),
00211 dict.subDict("sampledSurfaceDict")
00212 );
00213 surfacePtr_().update();
00214 nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
00215 }
00216
00217
00218
00219
00220 void Foam::fieldValues::faceSource::initialise(const dictionary& dict)
00221 {
00222 switch (source_)
00223 {
00224 case stFaceZone:
00225 {
00226 setFaceZoneFaces();
00227 break;
00228 }
00229 case stPatch:
00230 {
00231 setPatchFaces();
00232 break;
00233 }
00234 case stSampledSurface:
00235 {
00236 sampledSurfaceFaces(dict);
00237 break;
00238 }
00239 default:
00240 {
00241 FatalErrorIn("faceSource::initialise()")
00242 << type() << " " << name_ << ": "
00243 << sourceTypeNames_[source_] << "(" << sourceName_ << "):"
00244 << nl << " Unknown source type. Valid source types are:"
00245 << sourceTypeNames_ << nl << exit(FatalError);
00246 }
00247 }
00248
00249 scalar totalArea;
00250
00251 if (surfacePtr_.valid())
00252 {
00253 surfacePtr_().update();
00254 totalArea = gSum(surfacePtr_().magSf());
00255 }
00256 else
00257 {
00258 totalArea = gSum(filterField(mesh().magSf(), false));
00259 }
00260
00261 Info<< type() << " " << name_ << ":" << nl
00262 << " total faces = " << nFaces_
00263 << nl
00264 << " total area = " << totalArea
00265 << nl;
00266
00267 if (operation_ == opWeightedAverage)
00268 {
00269 dict.lookup("weightField") >> weightFieldName_;
00270 if
00271 (
00272 obr().foundObject<volScalarField>(weightFieldName_)
00273 || obr().foundObject<surfaceScalarField>(weightFieldName_)
00274 )
00275 {
00276 Info<< " weight field = " << weightFieldName_;
00277 }
00278 else
00279 {
00280 FatalErrorIn("faceSource::initialise()")
00281 << type() << " " << name_ << ": "
00282 << sourceTypeNames_[source_] << "(" << sourceName_ << "):"
00283 << nl << " Weight field " << weightFieldName_
00284 << " must be either a " << volScalarField::typeName << " or "
00285 << surfaceScalarField::typeName << nl << exit(FatalError);
00286 }
00287 }
00288
00289 Info<< nl << endl;
00290 }
00291
00292
00293 void Foam::fieldValues::faceSource::writeFileHeader()
00294 {
00295 if (outputFilePtr_.valid())
00296 {
00297 outputFilePtr_()
00298 << "# Source : " << sourceTypeNames_[source_] << " "
00299 << sourceName_ << nl << "# Faces : " << nFaces_ << nl
00300 << "# Time" << tab << "sum(magSf)";
00301
00302 forAll(fields_, i)
00303 {
00304 outputFilePtr_()
00305 << tab << operationTypeNames_[operation_]
00306 << "(" << fields_[i] << ")";
00307 }
00308
00309 outputFilePtr_() << endl;
00310 }
00311 }
00312
00313
00314
00315
00316 Foam::fieldValues::faceSource::faceSource
00317 (
00318 const word& name,
00319 const objectRegistry& obr,
00320 const dictionary& dict,
00321 const bool loadFromFiles
00322 )
00323 :
00324 fieldValue(name, obr, dict, loadFromFiles),
00325 source_(sourceTypeNames_.read(dict.lookup("source"))),
00326 operation_(operationTypeNames_.read(dict.lookup("operation"))),
00327 weightFieldName_("undefinedWeightedFieldName"),
00328 nFaces_(0),
00329 faceId_(),
00330 facePatchId_(),
00331 faceSign_()
00332 {
00333 read(dict);
00334 }
00335
00336
00337
00338
00339 Foam::fieldValues::faceSource::~faceSource()
00340 {}
00341
00342
00343
00344
00345 void Foam::fieldValues::faceSource::read(const dictionary& dict)
00346 {
00347 fieldValue::read(dict);
00348
00349 if (active_)
00350 {
00351 initialise(dict);
00352 }
00353 }
00354
00355
00356 void Foam::fieldValues::faceSource::write()
00357 {
00358 fieldValue::write();
00359
00360 if (active_)
00361 {
00362 scalar totalArea;
00363
00364 if (surfacePtr_.valid())
00365 {
00366 surfacePtr_().update();
00367 totalArea = gSum(surfacePtr_().magSf());
00368 }
00369 else
00370 {
00371 totalArea = gSum(filterField(mesh().magSf(), false));
00372 }
00373
00374 if (Pstream::master())
00375 {
00376 outputFilePtr_() << obr_.time().value() << tab << totalArea;
00377 }
00378
00379 forAll(fields_, i)
00380 {
00381 writeValues<scalar>(fields_[i]);
00382 writeValues<vector>(fields_[i]);
00383 writeValues<sphericalTensor>(fields_[i]);
00384 writeValues<symmTensor>(fields_[i]);
00385 writeValues<tensor>(fields_[i]);
00386 }
00387
00388 if (Pstream::master())
00389 {
00390 outputFilePtr_()<< endl;
00391 }
00392
00393 if (log_)
00394 {
00395 Info<< endl;
00396 }
00397 }
00398 }
00399
00400
00401