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-2018 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"
27 #include "dictionary.H"
28 #include "volFields.H"
29 #include "volPointInterpolation.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(sampledIsoSurface, 0);
38  (
39  sampledSurface,
40  sampledIsoSurface,
41  word,
42  isoSurface
43  );
44 }
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::sampledIsoSurface::getIsoFields() const
49 {
50  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
51 
52  // Get volField
53  // ~~~~~~~~~~~~
54 
55  if (fvm.foundObject<volScalarField>(isoField_))
56  {
57  if (debug)
58  {
60  << "Lookup volField " << isoField_ << endl;
61  }
62  storedVolFieldPtr_.clear();
63  volFieldPtr_ = &fvm.lookupObject<volScalarField>(isoField_);
64  }
65  else
66  {
67  // Bit of a hack. Read field and store.
68 
69  if (debug)
70  {
72  << "Checking " << isoField_
73  << " for same time " << fvm.time().timeName()
74  << endl;
75  }
76 
77  if
78  (
79  storedVolFieldPtr_.empty()
80  || (fvm.time().timeName() != storedVolFieldPtr_().instance())
81  )
82  {
83  if (debug)
84  {
86  << "Reading volField " << isoField_
87  << " from time " << fvm.time().timeName()
88  << endl;
89  }
90 
91  IOobject vfHeader
92  (
93  isoField_,
94  fvm.time().timeName(),
95  fvm,
98  false
99  );
100 
101  if (vfHeader.typeHeaderOk<volScalarField>(true))
102  {
103  storedVolFieldPtr_.reset
104  (
105  new volScalarField
106  (
107  vfHeader,
108  fvm
109  )
110  );
111  volFieldPtr_ = storedVolFieldPtr_.operator->();
112  }
113  else
114  {
116  << "Cannot find isosurface field " << isoField_
117  << " in database or directory " << vfHeader.path()
118  << exit(FatalError);
119  }
120  }
121  }
122 
123 
124  // Get pointField
125  // ~~~~~~~~~~~~~~
126 
127  // In case of multiple iso values we don't want to calculate multiple e.g.
128  // "volPointInterpolate(p)" so register it and re-use it. This is the
129  // same as the 'cache' functionality from volPointInterpolate but
130  // unfortunately that one does not guarantee that the field pointer
131  // remain: e.g. some other functionObject might delete the cached version.
132  // (volPointInterpolation::interpolate with cache=false deletes any
133  // registered one or if mesh.changing())
134 
135  if (!subMeshPtr_.valid())
136  {
137  const word pointFldName =
138  "volPointInterpolate_"
139  + type()
140  + "("
141  + isoField_
142  + ')';
143 
144  if (fvm.foundObject<pointScalarField>(pointFldName))
145  {
146  if (debug)
147  {
149  << "lookup pointField " << pointFldName << endl;
150  }
151  const pointScalarField& pfld = fvm.lookupObject<pointScalarField>
152  (
153  pointFldName
154  );
155 
156  if (!pfld.upToDate(*volFieldPtr_))
157  {
158  if (debug)
159  {
161  << "updating pointField "
162  << pointFldName << endl;
163  }
164  // Update the interpolated value
166  (
167  *volFieldPtr_,
168  const_cast<pointScalarField&>(pfld)
169  );
170  }
171 
172  pointFieldPtr_ = &pfld;
173  }
174  else
175  {
176  // Not in registry. Interpolate.
177 
178  if (debug)
179  {
181  << "Checking pointField " << pointFldName
182  << " for same time " << fvm.time().timeName()
183  << endl;
184  }
185 
186  // Interpolate without cache. Note that we're registering it
187  // below so next time round it goes into the condition
188  // above.
189  tmp<pointScalarField> tpfld
190  (
192  (
193  *volFieldPtr_,
194  pointFldName,
195  false
196  )
197  );
198  pointFieldPtr_ = tpfld.ptr();
199  const_cast<pointScalarField*>(pointFieldPtr_)->store();
200  }
201 
202 
203  // If averaging redo the volField. Can only be done now since needs the
204  // point field.
205  if (average_)
206  {
207  storedVolFieldPtr_.reset
208  (
209  pointAverage(*pointFieldPtr_).ptr()
210  );
211  volFieldPtr_ = storedVolFieldPtr_.operator->();
212  }
213 
214 
215  if (debug)
216  {
218  << "volField " << volFieldPtr_->name()
219  << " min:" << min(*volFieldPtr_).value()
220  << " max:" << max(*volFieldPtr_).value() << endl;
222  << "pointField " << pointFieldPtr_->name()
223  << " min:" << gMin(pointFieldPtr_->primitiveField())
224  << " max:" << gMax(pointFieldPtr_->primitiveField()) << endl;
225  }
226  }
227  else
228  {
229  // Get subMesh variants
230  const fvMesh& subFvm = subMeshPtr_().subMesh();
231 
232  // Either lookup on the submesh or subset the whole-mesh volField
233 
234  if (subFvm.foundObject<volScalarField>(isoField_))
235  {
236  if (debug)
237  {
239  << "Sub-mesh lookup volField "
240  << isoField_ << endl;
241  }
242  storedVolSubFieldPtr_.clear();
243  volSubFieldPtr_ = &subFvm.lookupObject<volScalarField>(isoField_);
244  }
245  else
246  {
247  if (debug)
248  {
250  << "Sub-setting volField " << isoField_ << endl;
251  }
252  storedVolSubFieldPtr_.reset
253  (
254  subMeshPtr_().interpolate
255  (
256  *volFieldPtr_
257  ).ptr()
258  );
259  storedVolSubFieldPtr_->checkOut();
260  volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
261  }
262 
263 
264  // Pointfield on submesh
265 
266  word pointFldName =
267  "volPointInterpolate("
268  + volSubFieldPtr_->name()
269  + ')';
270 
271  if (subFvm.foundObject<pointScalarField>(pointFldName))
272  {
273  if (debug)
274  {
276  << "Sub-mesh lookup pointField " << pointFldName << endl;
277  }
278  storedPointSubFieldPtr_.clear();
279  pointSubFieldPtr_ = &subFvm.lookupObject<pointScalarField>
280  (
281  pointFldName
282  );
283  }
284  else
285  {
286  if (debug)
287  {
289  << "Interpolating submesh volField "
290  << volSubFieldPtr_->name()
291  << " to get submesh pointField " << pointFldName << endl;
292  }
293  storedPointSubFieldPtr_.reset
294  (
296  (
297  subFvm
298  ).interpolate(*volSubFieldPtr_).ptr()
299  );
300  storedPointSubFieldPtr_->checkOut();
301  pointSubFieldPtr_ = storedPointSubFieldPtr_.operator->();
302  }
303 
304 
305  // If averaging redo the volField. Can only be done now since needs the
306  // point field.
307  if (average_)
308  {
309  storedVolSubFieldPtr_.reset
310  (
311  pointAverage(*pointSubFieldPtr_).ptr()
312  );
313  volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
314  }
315 
316 
317  if (debug)
318  {
320  << "volSubField "
321  << volSubFieldPtr_->name()
322  << " min:" << min(*volSubFieldPtr_).value()
323  << " max:" << max(*volSubFieldPtr_).value() << endl;
325  << "pointSubField "
326  << pointSubFieldPtr_->name()
327  << " min:" << gMin(pointSubFieldPtr_->primitiveField())
328  << " max:" << gMax(pointSubFieldPtr_->primitiveField()) << endl;
329  }
330  }
331 }
332 
333 
334 bool Foam::sampledIsoSurface::updateGeometry() const
335 {
336  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
337 
338  // No update needed
339  if (fvm.time().timeIndex() == prevTimeIndex_)
340  {
341  return false;
342  }
343 
344  // Get any subMesh
345  if (zoneID_.index() != -1 && !subMeshPtr_.valid())
346  {
347  const polyBoundaryMesh& patches = mesh().boundaryMesh();
348 
349  // Patch to put exposed internal faces into
350  const label exposedPatchi = patches.findPatchID(exposedPatchName_);
351 
352  if (debug)
353  {
354  Info<< "Allocating subset of size "
355  << mesh().cellZones()[zoneID_.index()].size()
356  << " with exposed faces into patch "
357  << patches[exposedPatchi].name() << endl;
358  }
359 
360  subMeshPtr_.reset
361  (
362  new fvMeshSubset(fvm)
363  );
364  subMeshPtr_().setLargeCellSubset
365  (
366  labelHashSet(mesh().cellZones()[zoneID_.index()]),
367  exposedPatchi
368  );
369  }
370 
371 
372  prevTimeIndex_ = fvm.time().timeIndex();
373  getIsoFields();
374 
375  // Clear any stored topo
376  surfPtr_.clear();
377  facesPtr_.clear();
378 
379  // Clear derived data
380  clearGeom();
381 
382  if (subMeshPtr_.valid())
383  {
384  surfPtr_.reset
385  (
386  new isoSurface
387  (
388  *volSubFieldPtr_,
389  *pointSubFieldPtr_,
390  isoVal_,
391  regularise_,
392  mergeTol_
393  )
394  );
395  }
396  else
397  {
398  surfPtr_.reset
399  (
400  new isoSurface
401  (
402  *volFieldPtr_,
403  *pointFieldPtr_,
404  isoVal_,
405  regularise_,
406  mergeTol_
407  )
408  );
409  }
410 
411 
412  if (debug)
413  {
414  Pout<< "sampledIsoSurface::updateGeometry() : constructed iso:"
415  << nl
416  << " regularise : " << regularise_ << nl
417  << " average : " << average_ << nl
418  << " isoField : " << isoField_ << nl
419  << " isoValue : " << isoVal_ << nl;
420  if (subMeshPtr_.valid())
421  {
422  Pout<< " zone size : " << subMeshPtr_().subMesh().nCells()
423  << nl;
424  }
425  Pout<< " points : " << points().size() << nl
426  << " tris : " << surface().size() << nl
427  << " cut cells : " << surface().meshCells().size()
428  << endl;
429  }
430 
431  return true;
432 }
433 
434 
435 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
436 
438 (
439  const word& name,
440  const polyMesh& mesh,
441  const dictionary& dict
442 )
443 :
444  sampledSurface(name, mesh, dict),
445  isoField_(dict.lookup("isoField")),
446  isoVal_(readScalar(dict.lookup("isoValue"))),
447  mergeTol_(dict.lookupOrDefault("mergeTol", 1e-6)),
448  regularise_(dict.lookupOrDefault("regularise", true)),
449  average_(dict.lookupOrDefault("average", false)),
450  zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
451  exposedPatchName_(word::null),
452  surfPtr_(nullptr),
453  facesPtr_(nullptr),
454  prevTimeIndex_(-1),
455  storedVolFieldPtr_(nullptr),
456  volFieldPtr_(nullptr),
457  pointFieldPtr_(nullptr)
458 {
460  {
462  (
463  dict
464  ) << "Non-interpolated iso surface not supported since triangles"
465  << " span across cells." << exit(FatalIOError);
466  }
467 
468  if (zoneID_.index() != -1)
469  {
470  dict.lookup("exposedPatchName") >> exposedPatchName_;
471 
472  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
473  {
475  (
476  dict
477  ) << "Cannot find patch " << exposedPatchName_
478  << " in which to put exposed faces." << endl
479  << "Valid patches are " << mesh.boundaryMesh().names()
480  << exit(FatalIOError);
481  }
482 
483  if (debug && zoneID_.index() != -1)
484  {
485  Info<< "Restricting to cellZone " << zoneID_.name()
486  << " with exposed internal faces into patch "
487  << exposedPatchName_ << endl;
488  }
489  }
490 }
491 
492 
493 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
494 
496 {}
497 
498 
499 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
500 
502 {
503  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
504 
505  return fvm.time().timeIndex() != prevTimeIndex_;
506 }
507 
508 
510 {
511  surfPtr_.clear();
512  facesPtr_.clear();
513  subMeshPtr_.clear();
514 
515  // Clear derived data
516  clearGeom();
517 
518  // already marked as expired
519  if (prevTimeIndex_ == -1)
520  {
521  return false;
522  }
523 
524  // force update
525  prevTimeIndex_ = -1;
526  return true;
527 }
528 
529 
531 {
532  return updateGeometry();
533 }
534 
535 
537 (
538  const volScalarField& vField
539 ) const
540 {
541  return sampleField(vField);
542 }
543 
544 
546 (
547  const volVectorField& vField
548 ) const
549 {
550  return sampleField(vField);
551 }
552 
553 
555 (
556  const volSphericalTensorField& vField
557 ) const
558 {
559  return sampleField(vField);
560 }
561 
562 
564 (
565  const volSymmTensorField& vField
566 ) const
567 {
568  return sampleField(vField);
569 }
570 
571 
573 (
574  const volTensorField& vField
575 ) const
576 {
577  return sampleField(vField);
578 }
579 
580 
582 (
583  const interpolation<scalar>& interpolator
584 ) const
585 {
586  return interpolateField(interpolator);
587 }
588 
589 
591 (
592  const interpolation<vector>& interpolator
593 ) const
594 {
595  return interpolateField(interpolator);
596 }
597 
599 (
600  const interpolation<sphericalTensor>& interpolator
601 ) const
602 {
603  return interpolateField(interpolator);
604 }
605 
606 
608 (
609  const interpolation<symmTensor>& interpolator
610 ) const
611 {
612  return interpolateField(interpolator);
613 }
614 
615 
617 (
618  const interpolation<tensor>& interpolator
619 ) const
620 {
621  return interpolateField(interpolator);
622 }
623 
624 
626 {
627  os << "sampledIsoSurface: " << name() << " :"
628  << " field :" << isoField_
629  << " value :" << isoVal_;
630 }
631 
632 
633 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
An abstract class for surfaces with sampling.
Type gMin(const FieldField< Field, Type > &f)
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:163
bool interpolate() const
Interpolation requested for surface.
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
virtual ~sampledIsoSurface()
Destructor.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
label findPatchID(const word &patchName) const
Find patch index given a name.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Macros for easy insertion into run-time selection tables.
static const volPointInterpolation & New(const fvMesh &mesh)
Definition: MeshObject.C:44
virtual bool update()
Update the surface as required.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
virtual bool needsUpdate() const
Does the surface need an update?
dynamicFvMesh & mesh
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:472
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
wordList names() const
Return a list of patch names.
static const word null
An empty word.
Definition: word.H:77
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
virtual void print(Ostream &) const
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
addNamedToRunTimeSelectionTable(GAMGProcAgglomeration, noneGAMGProcAgglomeration, GAMGAgglomeration, none)
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
friend Ostream & operator(Ostream &, const GeometricField< Type, PatchField, GeoMesh > &)
messageStream Info
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
virtual bool expire()
Mark the surface as needing an update.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
sampledIsoSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
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:576
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.