sampledIsoSurface.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 "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  {
59  Info<< "sampledIsoSurface::getIsoFields() : lookup volField "
60  << 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  {
71  Info<< "sampledIsoSurface::getIsoFields() : checking "
72  << isoField_ << " for same time " << fvm.time().timeName()
73  << endl;
74  }
75 
76  if
77  (
78  storedVolFieldPtr_.empty()
79  || (fvm.time().timeName() != storedVolFieldPtr_().instance())
80  )
81  {
82  if (debug)
83  {
84  Info<< "sampledIsoSurface::getIsoFields() : reading volField "
85  << isoField_ << " from time " << fvm.time().timeName()
86  << endl;
87  }
88 
89  IOobject vfHeader
90  (
91  isoField_,
92  fvm.time().timeName(),
93  fvm,
96  false
97  );
98 
99  if (vfHeader.headerOk())
100  {
101  storedVolFieldPtr_.reset
102  (
103  new volScalarField
104  (
105  vfHeader,
106  fvm
107  )
108  );
109  volFieldPtr_ = storedVolFieldPtr_.operator->();
110  }
111  else
112  {
113  FatalErrorIn("sampledIsoSurface::getIsoFields()")
114  << "Cannot find isosurface field " << isoField_
115  << " in database or directory " << vfHeader.path()
116  << exit(FatalError);
117  }
118  }
119  }
120 
121 
122  // Get pointField
123  // ~~~~~~~~~~~~~~
124 
125  if (!subMeshPtr_.valid())
126  {
127  word pointFldName = "volPointInterpolate(" + isoField_ + ')';
128 
129  if (fvm.foundObject<pointScalarField>(pointFldName))
130  {
131  if (debug)
132  {
133  Info<< "sampledIsoSurface::getIsoFields() : lookup pointField "
134  << pointFldName << endl;
135  }
136  pointFieldPtr_ = &fvm.lookupObject<pointScalarField>(pointFldName);
137  }
138  else
139  {
140  // Not in registry. Interpolate.
141 
142  if (debug)
143  {
144  Info<< "sampledIsoSurface::getIsoFields() : "
145  << "checking pointField " << pointFldName
146  << " for same time " << fvm.time().timeName()
147  << endl;
148  }
149 
150  if
151  (
152  storedPointFieldPtr_.empty()
153  || (fvm.time().timeName() != storedPointFieldPtr_().instance())
154  )
155  {
156  if (debug)
157  {
158  Info<< "sampledIsoSurface::getIsoFields() :"
159  << " interpolating volField " << volFieldPtr_->name()
160  << " to get pointField " << pointFldName << endl;
161  }
162 
163  storedPointFieldPtr_.reset
164  (
166  .interpolate(*volFieldPtr_).ptr()
167  );
168  storedPointFieldPtr_->checkOut();
169  pointFieldPtr_ = storedPointFieldPtr_.operator->();
170  }
171  }
172 
173 
174  // If averaging redo the volField. Can only be done now since needs the
175  // point field.
176  if (average_)
177  {
178  storedVolFieldPtr_.reset
179  (
180  pointAverage(*pointFieldPtr_).ptr()
181  );
182  volFieldPtr_ = storedVolFieldPtr_.operator->();
183  }
184 
185 
186  if (debug)
187  {
188  Info<< "sampledIsoSurface::getIsoFields() : volField "
189  << volFieldPtr_->name() << " min:" << min(*volFieldPtr_).value()
190  << " max:" << max(*volFieldPtr_).value() << endl;
191  Info<< "sampledIsoSurface::getIsoFields() : pointField "
192  << pointFieldPtr_->name()
193  << " min:" << gMin(pointFieldPtr_->internalField())
194  << " max:" << gMax(pointFieldPtr_->internalField()) << endl;
195  }
196  }
197  else
198  {
199  // Get subMesh variants
200  const fvMesh& subFvm = subMeshPtr_().subMesh();
201 
202  // Either lookup on the submesh or subset the whole-mesh volField
203 
204  if (subFvm.foundObject<volScalarField>(isoField_))
205  {
206  if (debug)
207  {
208  Info<< "sampledIsoSurface::getIsoFields() :"
209  << " submesh lookup volField "
210  << isoField_ << endl;
211  }
212  storedVolSubFieldPtr_.clear();
213  volSubFieldPtr_ = &subFvm.lookupObject<volScalarField>(isoField_);
214  }
215  else
216  {
217  if (debug)
218  {
219  Info<< "sampledIsoSurface::getIsoFields() : "
220  << "subsetting volField " << isoField_ << endl;
221  }
222  storedVolSubFieldPtr_.reset
223  (
224  subMeshPtr_().interpolate
225  (
226  *volFieldPtr_
227  ).ptr()
228  );
229  storedVolSubFieldPtr_->checkOut();
230  volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
231  }
232 
233 
234  // Pointfield on submesh
235 
236  word pointFldName =
237  "volPointInterpolate("
238  + volSubFieldPtr_->name()
239  + ')';
240 
241  if (subFvm.foundObject<pointScalarField>(pointFldName))
242  {
243  if (debug)
244  {
245  Info<< "sampledIsoSurface::getIsoFields() :"
246  << " submesh lookup pointField " << pointFldName << endl;
247  }
248  storedPointSubFieldPtr_.clear();
249  pointSubFieldPtr_ = &subFvm.lookupObject<pointScalarField>
250  (
251  pointFldName
252  );
253  }
254  else
255  {
256  if (debug)
257  {
258  Info<< "sampledIsoSurface::getIsoFields() :"
259  << " interpolating submesh volField "
260  << volSubFieldPtr_->name()
261  << " to get submesh pointField " << pointFldName << endl;
262  }
263  storedPointSubFieldPtr_.reset
264  (
266  (
267  subFvm
268  ).interpolate(*volSubFieldPtr_).ptr()
269  );
270  storedPointSubFieldPtr_->checkOut();
271  pointSubFieldPtr_ = storedPointSubFieldPtr_.operator->();
272  }
273 
274 
275  // If averaging redo the volField. Can only be done now since needs the
276  // point field.
277  if (average_)
278  {
279  storedVolSubFieldPtr_.reset
280  (
281  pointAverage(*pointSubFieldPtr_).ptr()
282  );
283  volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
284  }
285 
286 
287  if (debug)
288  {
289  Info<< "sampledIsoSurface::getIsoFields() : volSubField "
290  << volSubFieldPtr_->name()
291  << " min:" << min(*volSubFieldPtr_).value()
292  << " max:" << max(*volSubFieldPtr_).value() << endl;
293  Info<< "sampledIsoSurface::getIsoFields() : pointSubField "
294  << pointSubFieldPtr_->name()
295  << " min:" << gMin(pointSubFieldPtr_->internalField())
296  << " max:" << gMax(pointSubFieldPtr_->internalField()) << endl;
297  }
298  }
299 }
300 
301 
302 bool Foam::sampledIsoSurface::updateGeometry() const
303 {
304  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
305 
306  // No update needed
307  if (fvm.time().timeIndex() == prevTimeIndex_)
308  {
309  return false;
310  }
311 
312  // Get any subMesh
313  if (zoneID_.index() != -1 && !subMeshPtr_.valid())
314  {
315  const polyBoundaryMesh& patches = mesh().boundaryMesh();
316 
317  // Patch to put exposed internal faces into
318  const label exposedPatchI = patches.findPatchID(exposedPatchName_);
319 
320  if (debug)
321  {
322  Info<< "Allocating subset of size "
323  << mesh().cellZones()[zoneID_.index()].size()
324  << " with exposed faces into patch "
325  << patches[exposedPatchI].name() << endl;
326  }
327 
328  subMeshPtr_.reset
329  (
330  new fvMeshSubset(fvm)
331  );
332  subMeshPtr_().setLargeCellSubset
333  (
334  labelHashSet(mesh().cellZones()[zoneID_.index()]),
335  exposedPatchI
336  );
337  }
338 
339 
340  prevTimeIndex_ = fvm.time().timeIndex();
341  getIsoFields();
342 
343  // Clear any stored topo
344  surfPtr_.clear();
345  facesPtr_.clear();
346 
347  // Clear derived data
348  clearGeom();
349 
350  if (subMeshPtr_.valid())
351  {
352  surfPtr_.reset
353  (
354  new isoSurface
355  (
356  *volSubFieldPtr_,
357  *pointSubFieldPtr_,
358  isoVal_,
359  regularise_,
360  mergeTol_
361  )
362  );
363  }
364  else
365  {
366  surfPtr_.reset
367  (
368  new isoSurface
369  (
370  *volFieldPtr_,
371  *pointFieldPtr_,
372  isoVal_,
373  regularise_,
374  mergeTol_
375  )
376  );
377  }
378 
379 
380  if (debug)
381  {
382  Pout<< "sampledIsoSurface::updateGeometry() : constructed iso:"
383  << nl
384  << " regularise : " << regularise_ << nl
385  << " average : " << average_ << nl
386  << " isoField : " << isoField_ << nl
387  << " isoValue : " << isoVal_ << nl;
388  if (subMeshPtr_.valid())
389  {
390  Pout<< " zone size : " << subMeshPtr_().subMesh().nCells()
391  << nl;
392  }
393  Pout<< " points : " << points().size() << nl
394  << " tris : " << surface().size() << nl
395  << " cut cells : " << surface().meshCells().size()
396  << endl;
397  }
398 
399  return true;
400 }
401 
402 
403 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
404 
406 (
407  const word& name,
408  const polyMesh& mesh,
409  const dictionary& dict
410 )
411 :
412  sampledSurface(name, mesh, dict),
413  isoField_(dict.lookup("isoField")),
414  isoVal_(readScalar(dict.lookup("isoValue"))),
415  mergeTol_(dict.lookupOrDefault("mergeTol", 1e-6)),
416  regularise_(dict.lookupOrDefault("regularise", true)),
417  average_(dict.lookupOrDefault("average", false)),
418  zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
419  exposedPatchName_(word::null),
420  surfPtr_(NULL),
421  facesPtr_(NULL),
422  prevTimeIndex_(-1),
423  storedVolFieldPtr_(NULL),
424  volFieldPtr_(NULL),
425  storedPointFieldPtr_(NULL),
426  pointFieldPtr_(NULL)
427 {
429  {
431  (
432  "sampledIsoSurface::sampledIsoSurface"
433  "(const word&, const polyMesh&, const dictionary&)",
434  dict
435  ) << "Non-interpolated iso surface not supported since triangles"
436  << " span across cells." << exit(FatalIOError);
437  }
438 
439  if (zoneID_.index() != -1)
440  {
441  dict.lookup("exposedPatchName") >> exposedPatchName_;
442 
443  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
444  {
446  (
447  "sampledIsoSurface::sampledIsoSurface"
448  "(const word&, const polyMesh&, const dictionary&)",
449  dict
450  ) << "Cannot find patch " << exposedPatchName_
451  << " in which to put exposed faces." << endl
452  << "Valid patches are " << mesh.boundaryMesh().names()
453  << exit(FatalIOError);
454  }
455 
456  if (debug && zoneID_.index() != -1)
457  {
458  Info<< "Restricting to cellZone " << zoneID_.name()
459  << " with exposed internal faces into patch "
460  << exposedPatchName_ << endl;
461  }
462  }
463 }
464 
465 
466 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
467 
469 {}
470 
471 
472 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
473 
475 {
476  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
477 
478  return fvm.time().timeIndex() != prevTimeIndex_;
479 }
480 
481 
483 {
484  surfPtr_.clear();
485  facesPtr_.clear();
486  subMeshPtr_.clear();
487 
488  // Clear derived data
489  clearGeom();
490 
491  // already marked as expired
492  if (prevTimeIndex_ == -1)
493  {
494  return false;
495  }
496 
497  // force update
498  prevTimeIndex_ = -1;
499  return true;
500 }
501 
502 
504 {
505  return updateGeometry();
506 }
507 
508 
510 (
511  const volScalarField& vField
512 ) const
513 {
514  return sampleField(vField);
515 }
516 
517 
519 (
520  const volVectorField& vField
521 ) const
522 {
523  return sampleField(vField);
524 }
525 
526 
528 (
529  const volSphericalTensorField& vField
530 ) const
531 {
532  return sampleField(vField);
533 }
534 
535 
537 (
538  const volSymmTensorField& vField
539 ) const
540 {
541  return sampleField(vField);
542 }
543 
544 
546 (
547  const volTensorField& vField
548 ) const
549 {
550  return sampleField(vField);
551 }
552 
553 
555 (
556  const interpolation<scalar>& interpolator
557 ) const
558 {
559  return interpolateField(interpolator);
560 }
561 
562 
564 (
565  const interpolation<vector>& interpolator
566 ) const
567 {
568  return interpolateField(interpolator);
569 }
570 
572 (
573  const interpolation<sphericalTensor>& interpolator
574 ) const
575 {
576  return interpolateField(interpolator);
577 }
578 
579 
581 (
582  const interpolation<symmTensor>& interpolator
583 ) const
584 {
585  return interpolateField(interpolator);
586 }
587 
588 
590 (
591  const interpolation<tensor>& interpolator
592 ) const
593 {
594  return interpolateField(interpolator);
595 }
596 
597 
599 {
600  os << "sampledIsoSurface: " << name() << " :"
601  << " field :" << isoField_
602  << " value :" << isoVal_;
603 }
604 
605 
606 // ************************************************************************* //
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
const pointField & points
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
An abstract class for surfaces with sampling.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const volPointInterpolation & New(const fvMesh &mesh)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
virtual bool update()
Update the surface as required.
A class for handling words, derived from string.
Definition: word.H:59
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
virtual ~sampledIsoSurface()
Destructor.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
sampledIsoSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
wordList names() const
Return a list of patch names.
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Macros for easy insertion into run-time selection tables.
Type gMin(const FieldField< Field, Type > &f)
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
bool interpolate() const
Interpolation requested for surface.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual bool needsUpdate() const
Does the surface need an update?
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
addNamedToRunTimeSelectionTable(GAMGProcAgglomeration, noneGAMGProcAgglomeration, GAMGAgglomeration, none)
label findPatchID(const word &patchName) const
Find patch index given a name.
virtual void print(Ostream &) const
Write.
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
Type gMax(const FieldField< Field, Type > &f)
A class for managing temporary objects.
Definition: PtrList.H:118
virtual bool expire()
Mark the surface as needing an update.
static const word null
An empty word.
Definition: word.H:77
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53