sampledIsoSurface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 "sampledIsoSurface.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace sampledSurfaces
34 {
35  defineTypeNameAndDebug(isoSurface, 0);
36  addToRunTimeSelectionTable(sampledSurface, isoSurface, word);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 bool Foam::sampledSurfaces::isoSurface::updateGeometry() const
44 {
45  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
46 
47  // no update needed
48  if (fvm.time().timeIndex() == prevTimeIndex_)
49  {
50  return false;
51  }
52 
53  prevTimeIndex_ = fvm.time().timeIndex();
54 
55  // Clear derived data
57 
58  // Optionally read volScalarField
59  autoPtr<volScalarField> readFieldPtr_;
60 
61  // 1. see if field in database
62  // 2. see if field can be read
63  const volScalarField* cellFldPtr = nullptr;
64  if (fvm.foundObject<volScalarField>(isoField_))
65  {
66  if (debug)
67  {
68  InfoInFunction << "Lookup " << isoField_ << endl;
69  }
70 
71  cellFldPtr = &fvm.lookupObject<volScalarField>(isoField_);
72  }
73  else
74  {
75  // Bit of a hack. Read field and store.
76 
77  if (debug)
78  {
80  << "Reading " << isoField_
81  << " from time " <<fvm.time().timeName()
82  << endl;
83  }
84 
85  readFieldPtr_.reset
86  (
87  new volScalarField
88  (
89  IOobject
90  (
91  isoField_,
92  fvm.time().timeName(),
93  fvm,
96  false
97  ),
98  fvm
99  )
100  );
101 
102  cellFldPtr = readFieldPtr_.operator->();
103  }
104  const volScalarField& cellFld = *cellFldPtr;
105 
106  tmp<pointScalarField> pointFld
107  (
109  );
110 
111  PtrList<Foam::isoSurface> isos(isoVals_.size());
112  forAll(isos, isoi)
113  {
114  isos.set
115  (
116  isoi,
117  new Foam::isoSurface
118  (
119  fvm,
120  cellFld.primitiveField(),
121  pointFld().primitiveField(),
122  isoVals_[isoi],
123  filter_
124  )
125  );
126  }
127 
128  if (isos.size() == 1)
129  {
130  // Straight transfer
131  const_cast<isoSurface&>
132  (
133  *this
135  meshCells_ = isos[0].meshCells();
136  }
137  else
138  {
139  label nFaces = 0;
140  label nPoints = 0;
141  forAll(isos, isoi)
142  {
143  nFaces += isos[isoi].size();
144  nPoints += isos[isoi].points().size();
145  }
146 
147  faceList allFaces(nFaces);
148  labelList allMeshCells(nFaces);
149  pointField allPoints(nPoints);
150 
151  nFaces = 0;
152  nPoints = 0;
153  forAll(isos, isoi)
154  {
155  Foam::isoSurface& iso = isos[isoi];
156 
157  SubList<face> subAll(allFaces, iso.size(), nFaces);
158  subAll = iso;
159  // Offset point indices
160  if (nPoints > 0)
161  {
162  forAll(subAll, i)
163  {
164  face& f = subAll[i];
165  forAll(f, fp)
166  {
167  f[fp] += nPoints;
168  }
169  }
170  }
171  SubList<label>(allMeshCells, iso.size(), nFaces) = iso.meshCells();
172  nFaces += iso.size();
173 
174  const pointField& pts = iso.points();
175  SubList<point>(allPoints, pts.size(), nPoints) = pts;
176  nPoints += pts.size();
177 
178  // Clear out asap
179  isos.set(isoi, nullptr);
180  }
181 
182  if (nFaces != allFaces.size() || nPoints != allPoints.size())
183  {
184  FatalErrorInFunction << "nFaces:" << nFaces
185  << " nPoints:" << nPoints << exit(FatalError);
186  }
187 
188 
189  surfZoneList allZones(1);
190  allZones[0] = surfZone
191  (
192  "allFaces",
193  allFaces.size(), // size
194  0, // start
195  0 // index
196  );
197 
198  // Transfer
199  const_cast<isoSurface&>
200  (
201  *this
203  (
204  move(allPoints),
205  move(allFaces),
206  move(allZones)
207  );
208  meshCells_.transfer(allMeshCells);
209  }
210  if (debug)
211  {
212  Pout<< "sampledSurfaces::isoSurface::updateGeometry() : "
213  "constructed iso:"
214  << nl
215  << " filtering : "
216  << Foam::isoSurface::filterTypeNames_[filter_] << nl
217  << " isoField : " << isoField_ << nl;
218  if (isoVals_.size() == 1)
219  {
220  Pout<< " isoValue : " << isoVals_[0] << nl;
221  }
222  else
223  {
224  Pout<< " isoValues : " << isoVals_ << nl;
225  }
226  Pout<< " points : " << points().size() << nl
227  << " faces : " << faces().size() << nl
228  << " cut cells : " << meshCells_.size() << endl;
229  }
230 
231  return true;
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
236 
238 (
239  const word& name,
240  const polyMesh& mesh,
241  const dictionary& dict
242 )
243 :
244  sampledSurface(name, mesh, dict),
245  isoField_(dict.lookup("isoField")),
246  isoVals_
247  (
248  dict.found("isoValues")
249  ? scalarField(dict.lookup("isoValues"))
250  : scalarField(1, dict.lookup<scalar>("isoValue"))
251  ),
252  filter_
253  (
254  dict.found("filtering")
255  ? Foam::isoSurface::filterTypeNames_.read(dict.lookup("filtering"))
257  ),
258  prevTimeIndex_(-1),
259  meshCells_(0)
260 {}
261 
262 
263 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
264 
266 {}
267 
268 
269 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
270 
272 {
273  return wordList(isoField_);
274 }
275 
276 
278 {
279  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
280 
281  return fvm.time().timeIndex() != prevTimeIndex_;
282 }
283 
284 
286 {
287  // Clear derived data
290 
291  // already marked as expired
292  if (prevTimeIndex_ == -1)
293  {
294  return false;
295  }
296 
297  // force update
298  prevTimeIndex_ = -1;
299  return true;
300 }
301 
302 
304 {
305  return updateGeometry();
306 }
307 
308 
311 (
312  const volScalarField& vField
313 ) const
314 {
315  return sampleField(vField);
316 }
317 
318 
321 (
322  const volVectorField& vField
323 ) const
324 {
325  return sampleField(vField);
326 }
327 
328 
331 (
332  const volSphericalTensorField& vField
333 ) const
334 {
335  return sampleField(vField);
336 }
337 
338 
341 (
342  const volSymmTensorField& vField
343 ) const
344 {
345  return sampleField(vField);
346 }
347 
348 
351 (
352  const volTensorField& vField
353 ) const
354 {
355  return sampleField(vField);
356 }
357 
358 
361 (
362  const interpolation<scalar>& interpolator
363 ) const
364 {
365  return interpolateField(interpolator);
366 }
367 
368 
371 (
372  const interpolation<vector>& interpolator
373 ) const
374 {
375  return interpolateField(interpolator);
376 }
377 
380 (
381  const interpolation<sphericalTensor>& interpolator
382 ) const
383 {
384  return interpolateField(interpolator);
385 }
386 
387 
390 (
391  const interpolation<symmTensor>& interpolator
392 ) const
393 {
394  return interpolateField(interpolator);
395 }
396 
397 
400 (
401  const interpolation<tensor>& interpolator
402 ) const
403 {
404  return interpolateField(interpolator);
405 }
406 
407 
409 {
410  os << "isoSurface: " << name() << " :"
411  << " field:" << isoField_;
412  if (isoVals_.size() == 1)
413  {
414  os << " value:" << isoVals_[0];
415  }
416  else
417  {
418  os << " values:" << isoVals_;
419  }
420  os << " faces:" << faces().size()
421  << " points:" << points().size();
422 }
423 
424 
425 // ************************************************************************* //
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual void print(Ostream &) const
Write.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
An abstract class for surfaces with sampling.
defineTypeNameAndDebug(distanceSurface, 0)
static const NamedEnum< filterType, 4 > filterTypeNames_
Definition: isoSurface.H:92
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
bool interpolate() const
Interpolation requested for surface.
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual bool needsUpdate() const
Does the surface need an update?
void transfer(MeshedSurface< face > &)
Transfer the contents of the argument and annul the argument.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
isoSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
fvMesh & mesh
const labelList & meshCells() const
For every face original cell in mesh.
Definition: isoSurface.H:255
virtual bool expire()
Mark the surface as needing an update.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
virtual void clearGeom() const
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Field< PointType > & points() const
Return reference to global points.
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual wordList fields() const
Return the list of fields required.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void reset(pointField &&points, List< face > &&faces, surfZoneList &&zones)
Reset primitive data (points, faces and zones)
static const char nl
Definition: Ostream.H:260
labelList f(nPoints)
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< surfZone > surfZoneList
Definition: surfZoneList.H:45
List< word > wordList
A List of words.
Definition: fileName.H:54
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
label size() const
The surface size is the number of faces.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
virtual bool update()
Update the surface as required.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
Marching tet iso surface algorithm with filtering to remove unnecessary topology. ...
Definition: isoSurface.H:78
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
#define InfoInFunction
Report an information message using Foam::Info.