1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
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 "sampledSurfaces.H"
28 #include "volFields.H"
29 #include "dictionary.H"
32 #include "ListListOps.H"
33 #include "mergePoints.H"
34 #include "volPointInterpolation.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 //- Used to offset faces in Pstream::combineOffset
53 face result(x.size());
57 result[xI] = x[xI] + offset;
64 defineTypeNameAndDebug(sampledSurfaces, 0);
68 bool Foam::sampledSurfaces::verbose_(false);
69 Foam::scalar Foam::sampledSurfaces::mergeTol_(1e-10);
71 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
73 Foam::label Foam::sampledSurfaces::classifyFieldTypes()
77 scalarFields_.clear();
78 vectorFields_.clear();
79 sphericalTensorFields_.clear();
80 symmTensorFields_.clear();
81 tensorFields_.clear();
83 forAll(fieldNames_, fieldI)
85 const word& fieldName = fieldNames_[fieldI];
88 // check files for a particular time
94 mesh_.time().timeName(),
103 fieldType = io.headerClassName();
112 // check objectRegistry
113 objectRegistry::const_iterator iter = mesh_.find(fieldName);
115 if (iter != mesh_.objectRegistry::end())
117 fieldType = iter()->type();
126 if (fieldType == volScalarField::typeName)
128 scalarFields_.append(fieldName);
131 else if (fieldType == volVectorField::typeName)
133 vectorFields_.append(fieldName);
136 else if (fieldType == volSphericalTensorField::typeName)
138 sphericalTensorFields_.append(fieldName);
141 else if (fieldType == volSymmTensorField::typeName)
143 symmTensorFields_.append(fieldName);
146 else if (fieldType == volTensorField::typeName)
148 tensorFields_.append(fieldName);
158 void Foam::sampledSurfaces::writeGeometry() const
160 // Write to time directory under outputPath_
161 // skip surface without faces (eg, a failed cut-plane)
163 const fileName outputDir = outputPath_/mesh_.time().timeName();
167 const sampledSurface& s = operator[](surfI);
169 if (Pstream::parRun())
171 if (Pstream::master() && mergeList_[surfI].faces.size())
173 genericFormatter_->write
177 mergeList_[surfI].points,
178 mergeList_[surfI].faces
182 else if (s.faces().size())
184 genericFormatter_->write
196 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198 Foam::sampledSurfaces::sampledSurfaces
201 const objectRegistry& obr,
202 const dictionary& dict,
203 const bool loadFromFiles
206 PtrList<sampledSurface>(),
208 mesh_(refCast<const fvMesh>(obr)),
209 loadFromFiles_(loadFromFiles),
210 outputPath_(fileName::null),
212 interpolationScheme_(word::null),
213 writeFormat_(word::null),
215 genericFormatter_(NULL),
218 sphericalTensorFields_(),
222 if (Pstream::parRun())
224 outputPath_ = mesh_.time().path()/".."/name_;
228 outputPath_ = mesh_.time().path()/name_;
235 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
237 Foam::sampledSurfaces::~sampledSurfaces()
241 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
243 void Foam::sampledSurfaces::verbose(const bool verbosity)
245 verbose_ = verbosity;
249 void Foam::sampledSurfaces::execute()
251 // Do nothing - only valid on write
255 void Foam::sampledSurfaces::end()
257 // Do nothing - only valid on write
261 void Foam::sampledSurfaces::write()
265 // finalize surfaces, merge points etc.
268 const label nFields = classifyFieldTypes();
270 if (Pstream::master())
274 Pout<< "timeName = " << mesh_.time().timeName() << nl
275 << "scalarFields " << scalarFields_ << nl
276 << "vectorFields " << vectorFields_ << nl
277 << "sphTensorFields " << sphericalTensorFields_ << nl
278 << "symTensorFields " << symmTensorFields_ <<nl
279 << "tensorFields " << tensorFields_ <<nl;
281 Pout<< "Creating directory "
282 << outputPath_/mesh_.time().timeName() << nl << endl;
286 mkDir(outputPath_/mesh_.time().timeName());
289 // write geometry first if required, or when no fields would otherwise
291 if (nFields == 0 || genericFormatter_->separateFiles())
296 sampleAndWrite(scalarFields_);
297 sampleAndWrite(vectorFields_);
298 sampleAndWrite(sphericalTensorFields_);
299 sampleAndWrite(symmTensorFields_);
300 sampleAndWrite(tensorFields_);
305 void Foam::sampledSurfaces::read(const dictionary& dict)
307 fieldNames_ = wordList(dict.lookup("fields"));
309 const label nFields = fieldNames_.size();
311 scalarFields_.reset(nFields);
312 vectorFields_.reset(nFields);
313 sphericalTensorFields_.reset(nFields);
314 symmTensorFields_.reset(nFields);
315 tensorFields_.reset(nFields);
317 interpolationScheme_ = dict.lookupOrDefault<word>
319 "interpolationScheme",
322 writeFormat_ = dict.lookupOrDefault<word>
329 // define the generic (geometry) writer
330 genericFormatter_ = surfaceWriter<bool>::New(writeFormat_);
333 PtrList<sampledSurface> newList
335 dict.lookup("surfaces"),
336 sampledSurface::iNew(mesh_)
341 if (Pstream::parRun())
343 mergeList_.setSize(size());
346 // ensure all surfaces and merge information are expired
349 if (Pstream::master() && debug)
351 Pout<< "sample fields:" << fieldNames_ << nl
352 << "sample surfaces:" << nl << "(" << nl;
356 Pout<< " " << operator[](surfI) << endl;
363 void Foam::sampledSurfaces::updateMesh(const mapPolyMesh&)
369 void Foam::sampledSurfaces::movePoints(const pointField&)
375 void Foam::sampledSurfaces::readUpdate(const polyMesh::readUpdateState state)
377 if (state != polyMesh::UNCHANGED)
384 bool Foam::sampledSurfaces::needsUpdate() const
388 if (operator[](surfI).needsUpdate())
398 bool Foam::sampledSurfaces::expire()
400 bool justExpired = false;
404 if (operator[](surfI).expire())
409 // clear merge information
410 if (Pstream::parRun())
412 mergeList_[surfI].clear();
416 // reset interpolation
417 pointMesh::Delete(mesh_);
418 volPointInterpolation::Delete(mesh_);
420 // true if any surfaces just expired
425 bool Foam::sampledSurfaces::update()
427 bool updated = false;
434 // serial: quick and easy, no merging required
435 if (!Pstream::parRun())
439 if (operator[](surfI).update())
448 // dimension as fraction of mesh bounding box
449 scalar mergeDim = mergeTol_ * mesh_.globalData().bb().mag();
451 if (Pstream::master() && debug)
453 Pout<< nl << "Merging all points within "
454 << mergeDim << " meter" << endl;
459 sampledSurface& s = operator[](surfI);
471 // Collect points from all processors
472 List<pointField> gatheredPoints(Pstream::nProcs());
473 gatheredPoints[Pstream::myProcNo()] = s.points();
474 Pstream::gatherList(gatheredPoints);
476 if (Pstream::master())
478 mergeList_[surfI].points = ListListOps::combine<pointField>
481 accessOp<pointField>()
485 // Collect faces from all processors and renumber using sizes of
487 List<faceList> gatheredFaces(Pstream::nProcs());
488 gatheredFaces[Pstream::myProcNo()] = s.faces();
489 Pstream::gatherList(gatheredFaces);
491 if (Pstream::master())
493 mergeList_[surfI].faces = static_cast<const faceList&>
495 ListListOps::combineOffset<faceList>
498 ListListOps::subSizes
501 accessOp<pointField>()
503 accessOp<faceList>(),
509 pointField newPoints;
512 bool hasMerged = mergePoints
514 mergeList_[surfI].points,
523 // Store point mapping
524 mergeList_[surfI].pointsMap.transfer(oldToNew);
527 mergeList_[surfI].points.transfer(newPoints);
530 faceList& faces = mergeList_[surfI].faces;
534 inplaceRenumber(mergeList_[surfI].pointsMap, faces[faceI]);
537 if (Pstream::master() && debug)
539 Pout<< "For surface " << surfI << " merged from "
540 << mergeList_[surfI].pointsMap.size() << " points down to "
541 << mergeList_[surfI].points.size() << " points" << endl;
550 // ************************************************************************* //