sampledIsoSurfaceCell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledIsoSurfaceCell.H"
27 #include "dictionary.H"
28 #include "volFields.H"
29 #include "volPointInterpolation.H"
31 #include "fvMesh.H"
32 #include "isoSurfaceCell.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(sampledIsoSurfaceCell, 0);
40  (
41  sampledSurface,
42  sampledIsoSurfaceCell,
43  word,
44  isoSurfaceCell
45  );
46 }
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 {
52  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
53 
54  // no update needed
55  if (fvm.time().timeIndex() == prevTimeIndex_)
56  {
57  return false;
58  }
59 
60  prevTimeIndex_ = fvm.time().timeIndex();
61 
62  // Clear any stored topo
63  facesPtr_.clear();
64 
65  // Clear derived data
67 
68  // Optionally read volScalarField
69  autoPtr<volScalarField> readFieldPtr_;
70 
71  // 1. see if field in database
72  // 2. see if field can be read
73  const volScalarField* cellFldPtr = NULL;
75  {
76  if (debug)
77  {
78  Info<< "sampledIsoSurfaceCell::updateGeometry() : lookup "
79  << isoField_ << endl;
80  }
81 
82  cellFldPtr = &fvm.lookupObject<volScalarField>(isoField_);
83  }
84  else
85  {
86  // Bit of a hack. Read field and store.
87 
88  if (debug)
89  {
90  Info<< "sampledIsoSurfaceCell::updateGeometry() : reading "
91  << isoField_ << " from time " <<fvm.time().timeName()
92  << endl;
93  }
94 
95  readFieldPtr_.reset
96  (
97  new volScalarField
98  (
99  IOobject
100  (
101  isoField_,
102  fvm.time().timeName(),
103  fvm,
106  false
107  ),
108  fvm
109  )
110  );
111 
112  cellFldPtr = readFieldPtr_.operator->();
113  }
114  const volScalarField& cellFld = *cellFldPtr;
115 
116  tmp<pointScalarField> pointFld
117  (
119  );
120 
121  if (average_)
122  {
123  //- From point field and interpolated cell.
124  scalarField cellAvg(fvm.nCells(), scalar(0.0));
125  labelField nPointCells(fvm.nCells(), 0);
126  {
127  for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
128  {
129  const labelList& pCells = fvm.pointCells(pointI);
130 
131  forAll(pCells, i)
132  {
133  label cellI = pCells[i];
134 
135  cellAvg[cellI] += pointFld().internalField()[pointI];
136  nPointCells[cellI]++;
137  }
138  }
139  }
140  forAll(cellAvg, cellI)
141  {
142  cellAvg[cellI] /= nPointCells[cellI];
143  }
144 
145  const isoSurfaceCell iso
146  (
147  fvm,
148  cellAvg,
149  pointFld().internalField(),
150  isoVal_,
151  regularise_,
152  bounds_
153  );
154 
155  const_cast<sampledIsoSurfaceCell&>
156  (
157  *this
158  ).triSurface::operator=(iso);
159  meshCells_ = iso.meshCells();
160  }
161  else
162  {
163  //- Direct from cell field and point field. Gives bad continuity.
164  const isoSurfaceCell iso
165  (
166  fvm,
167  cellFld.internalField(),
168  pointFld().internalField(),
169  isoVal_,
170  regularise_,
171  bounds_
172  );
173 
174  const_cast<sampledIsoSurfaceCell&>
175  (
176  *this
177  ).triSurface::operator=(iso);
178  meshCells_ = iso.meshCells();
179  }
180 
181 
182  if (debug)
183  {
184  Pout<< "sampledIsoSurfaceCell::updateGeometry() : constructed iso:"
185  << nl
186  << " regularise : " << regularise_ << nl
187  << " average : " << average_ << nl
188  << " isoField : " << isoField_ << nl
189  << " isoValue : " << isoVal_ << nl
190  << " bounds : " << bounds_ << nl
191  << " points : " << points().size() << nl
192  << " tris : " << triSurface::size() << nl
193  << " cut cells : " << meshCells_.size() << endl;
194  }
195 
196  return true;
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
201 
203 (
204  const word& name,
205  const polyMesh& mesh,
206  const dictionary& dict
207 )
208 :
210  isoField_(dict.lookup("isoField")),
211  isoVal_(readScalar(dict.lookup("isoValue"))),
212  bounds_(dict.lookupOrDefault("bounds", boundBox::greatBox)),
213  regularise_(dict.lookupOrDefault("regularise", true)),
214  average_(dict.lookupOrDefault("average", true)),
215  zoneKey_(keyType::null),
216  facesPtr_(NULL),
217  prevTimeIndex_(-1),
218  meshCells_(0)
219 {
220 // dict.readIfPresent("zone", zoneKey_);
221 //
222 // if (debug && zoneKey_.size() && mesh.cellZones().findZoneID(zoneKey_) < 0)
223 // {
224 // Info<< "cellZone " << zoneKey_
225 // << " not found - using entire mesh" << endl;
226 // }
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
231 
233 {}
234 
235 
236 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
237 
239 {
240  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
241 
242  return fvm.time().timeIndex() != prevTimeIndex_;
243 }
244 
245 
247 {
248  facesPtr_.clear();
249 
250  // Clear derived data
252 
253  // already marked as expired
254  if (prevTimeIndex_ == -1)
255  {
256  return false;
257  }
258 
259  // force update
260  prevTimeIndex_ = -1;
261  return true;
262 }
263 
264 
266 {
267  return updateGeometry();
268 }
269 
270 
273 (
274  const volScalarField& vField
275 ) const
276 {
277  return sampleField(vField);
278 }
279 
280 
283 (
284  const volVectorField& vField
285 ) const
286 {
287  return sampleField(vField);
288 }
289 
290 
293 (
294  const volSphericalTensorField& vField
295 ) const
296 {
297  return sampleField(vField);
298 }
299 
300 
303 (
304  const volSymmTensorField& vField
305 ) const
306 {
307  return sampleField(vField);
308 }
309 
310 
313 (
314  const volTensorField& vField
315 ) const
316 {
317  return sampleField(vField);
318 }
319 
320 
323 (
324  const interpolation<scalar>& interpolator
325 ) const
326 {
327  return interpolateField(interpolator);
328 }
329 
330 
333 (
334  const interpolation<vector>& interpolator
335 ) const
336 {
337  return interpolateField(interpolator);
338 }
339 
342 (
343  const interpolation<sphericalTensor>& interpolator
344 ) const
345 {
346  return interpolateField(interpolator);
347 }
348 
349 
352 (
353  const interpolation<symmTensor>& interpolator
354 ) const
355 {
356  return interpolateField(interpolator);
357 }
358 
359 
362 (
363  const interpolation<tensor>& interpolator
364 ) const
365 {
366  return interpolateField(interpolator);
367 }
368 
369 
371 {
372  os << "sampledIsoSurfaceCell: " << name() << " :"
373  << " field:" << isoField_
374  << " value:" << isoVal_;
375  //<< " faces:" << faces().size() // possibly no geom yet
376  //<< " points:" << points().size();
377 }
378 
379 
380 // ************************************************************************* //
Foam::sampledIsoSurfaceCell::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledIsoSurfaceCell.C:246
volFields.H
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::sampledIsoSurfaceCell::print
virtual void print(Ostream &) const
Write.
Definition: sampledIsoSurfaceCell.C:370
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch, cyclicPeriodicAMI)
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::sampledIsoSurfaceCell::meshCells_
labelList meshCells_
For every triangle the original cell in mesh.
Definition: sampledIsoSurfaceCell.H:86
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::sampledIsoSurfaceCell::prevTimeIndex_
label prevTimeIndex_
Time at last call, also track it surface needs an update.
Definition: sampledIsoSurfaceCell.H:83
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
isoSurfaceCell.H
Foam::sampledIsoSurfaceCell::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledIsoSurfaceCell.C:238
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::New
static const volPointInterpolation & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Foam::isoSurfaceCell
A surface formed by the iso value. After "Polygonising A Scalar Field Using Tetrahedrons",...
Definition: isoSurfaceCell.H:64
Foam::sampledIsoSurfaceCell::updateGeometry
bool updateGeometry() const
Create iso surface (if time has changed)
Definition: sampledIsoSurfaceCell.C:50
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::sampledIsoSurfaceCell::facesPtr_
autoPtr< faceList > facesPtr_
Triangles converted to faceList.
Definition: sampledIsoSurfaceCell.H:77
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::isoSurfaceCell::meshCells
const labelList & meshCells() const
For every face original cell in mesh.
Definition: isoSurfaceCell.H:333
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::sampledIsoSurfaceCell::regularise_
const Switch regularise_
Whether to coarse.
Definition: sampledIsoSurfaceCell.H:68
Foam::sampledIsoSurfaceCell::points
virtual const pointField & points() const
Points of surface.
Definition: sampledIsoSurfaceCell.H:145
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::keyType::null
static const keyType null
An empty keyType.
Definition: keyType.H:75
Foam::interpolation< scalar >
Foam::sampledIsoSurfaceCell::update
virtual bool update()
Update the surface as required.
Definition: sampledIsoSurfaceCell.C:265
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::sampledIsoSurfaceCell::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledIsoSurfaceCell.C:273
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::sampledIsoSurfaceCell::isoVal_
const scalar isoVal_
Iso value.
Definition: sampledIsoSurfaceCell.H:62
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::List< labelledTri >::size
label size() const
Return the number of elements in the UList.
internalField
conserve internalField()+
volPointInterpolation.H
sampledIsoSurfaceCell.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
Foam::sampledIsoSurfaceCell::bounds_
const boundBox bounds_
Optional bounding box to trim triangles against.
Definition: sampledIsoSurfaceCell.H:65
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::boundBox::greatBox
static const boundBox greatBox
A very large boundBox: min/max == -/+ VGREAT.
Definition: boundBox.H:76
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sampledIsoSurfaceCell::~sampledIsoSurfaceCell
virtual ~sampledIsoSurfaceCell()
Destructor.
Definition: sampledIsoSurfaceCell.C:232
dictionary.H
Foam::primitiveMesh::pointCells
const labelListList & pointCells() const
Definition: primitiveMeshPointCells.C:108
Foam::sampledIsoSurfaceCell::isoField_
const word isoField_
Field to get isoSurface of.
Definition: sampledIsoSurfaceCell.H:59
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Definition: sampledSurface.C:42
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:244
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::TimeState::timeIndex
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation requested for surface.
Definition: sampledSurface.H:256
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::sampledIsoSurfaceCell::sampledIsoSurfaceCell
sampledIsoSurfaceCell(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledIsoSurfaceCell.C:203
Foam::sampledIsoSurfaceCell
A sampledSurface defined by a surface of iso value. Always triangulated. To be used in sampleSurfaces...
Definition: sampledIsoSurfaceCell.H:51
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::sampledIsoSurfaceCell::average_
const Switch average_
Whether to recalculate cell values as average of point values.
Definition: sampledIsoSurfaceCell.H:71