distanceSurface.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 "distanceSurface.H"
27 #include "dictionary.H"
28 #include "volFields.H"
29 #include "volPointInterpolation.H"
31 #include "fvMesh.H"
32 #include "volumeType.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(distanceSurface, 0);
39  addToRunTimeSelectionTable(sampledSurface, distanceSurface, word);
40 }
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 {
46  if (debug)
47  {
48  Pout<< "distanceSurface::createGeometry :updating geometry." << endl;
49  }
50 
51  // Clear any stored topologies
52  facesPtr_.clear();
53  isoSurfCellPtr_.clear();
54  isoSurfPtr_.clear();
55 
56  // Clear derived data
57  clearGeom();
58 
59  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
60 
61  // Distance to cell centres
62  // ~~~~~~~~~~~~~~~~~~~~~~~~
63 
64  cellDistancePtr_.reset
65  (
66  new volScalarField
67  (
68  IOobject
69  (
70  "cellDistance",
71  fvm.time().timeName(),
72  fvm.time(),
75  false
76  ),
77  fvm,
78  dimensionedScalar("zero", dimLength, 0)
79  )
80  );
81  volScalarField& cellDistance = cellDistancePtr_();
82 
83  // Internal field
84  {
85  const pointField& cc = fvm.C();
86  scalarField& fld = cellDistance.internalField();
87 
88  List<pointIndexHit> nearest;
89  surfPtr_().findNearest
90  (
91  cc,
92  scalarField(cc.size(), GREAT),
93  nearest
94  );
95 
96  if (signed_)
97  {
98  List<volumeType> volType;
99 
100  surfPtr_().getVolumeType(cc, volType);
101 
102  forAll(volType, i)
103  {
104  volumeType vT = volType[i];
105 
106  if (vT == volumeType::OUTSIDE)
107  {
108  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
109  }
110  else if (vT == volumeType::INSIDE)
111  {
112  fld[i] = -Foam::mag(cc[i] - nearest[i].hitPoint());
113  }
114  else
115  {
117  << "getVolumeType failure, neither INSIDE or OUTSIDE"
118  << exit(FatalError);
119  }
120  }
121  }
122  else
123  {
124  forAll(nearest, i)
125  {
126  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
127  }
128  }
129  }
130 
131  // Patch fields
132  {
133  forAll(fvm.C().boundaryField(), patchI)
134  {
135  const pointField& cc = fvm.C().boundaryField()[patchI];
136  fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
137 
138  List<pointIndexHit> nearest;
139  surfPtr_().findNearest
140  (
141  cc,
142  scalarField(cc.size(), GREAT),
143  nearest
144  );
145 
146  if (signed_)
147  {
148  List<volumeType> volType;
149 
150  surfPtr_().getVolumeType(cc, volType);
151 
152  forAll(volType, i)
153  {
154  volumeType vT = volType[i];
155 
156  if (vT == volumeType::OUTSIDE)
157  {
158  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
159  }
160  else if (vT == volumeType::INSIDE)
161  {
162  fld[i] = -Foam::mag(cc[i] - nearest[i].hitPoint());
163  }
164  else
165  {
167  << "getVolumeType failure, "
168  << "neither INSIDE or OUTSIDE"
169  << exit(FatalError);
170  }
171  }
172  }
173  else
174  {
175  forAll(nearest, i)
176  {
177  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
178  }
179  }
180  }
181  }
182 
183 
184  // On processor patches the mesh.C() will already be the cell centre
185  // on the opposite side so no need to swap cellDistance.
186 
187 
188  // Distance to points
189  pointDistance_.setSize(fvm.nPoints());
190  {
191  const pointField& pts = fvm.points();
192 
193  List<pointIndexHit> nearest;
194  surfPtr_().findNearest
195  (
196  pts,
197  scalarField(pts.size(), GREAT),
198  nearest
199  );
200 
201  if (signed_)
202  {
203  List<volumeType> volType;
204 
205  surfPtr_().getVolumeType(pts, volType);
206 
207  forAll(volType, i)
208  {
209  volumeType vT = volType[i];
210 
211  if (vT == volumeType::OUTSIDE)
212  {
213  pointDistance_[i] =
214  Foam::mag(pts[i] - nearest[i].hitPoint());
215  }
216  else if (vT == volumeType::INSIDE)
217  {
218  pointDistance_[i] =
219  -Foam::mag(pts[i] - nearest[i].hitPoint());
220  }
221  else
222  {
224  << "getVolumeType failure, neither INSIDE or OUTSIDE"
225  << exit(FatalError);
226  }
227  }
228  }
229  else
230  {
231  forAll(nearest, i)
232  {
233  pointDistance_[i] = Foam::mag(pts[i]-nearest[i].hitPoint());
234  }
235  }
236  }
237 
238 
239  if (debug)
240  {
241  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
242  cellDistance.write();
243  pointScalarField pDist
244  (
245  IOobject
246  (
247  "pointDistance",
248  fvm.time().timeName(),
249  fvm.time(),
252  false
253  ),
254  pointMesh::New(fvm),
255  dimensionedScalar("zero", dimLength, 0)
256  );
257  pDist.internalField() = pointDistance_;
258 
259  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
260  pDist.write();
261  }
262 
263 
264  //- Direct from cell field and point field.
265  if (cell_)
266  {
267  isoSurfCellPtr_.reset
268  (
269  new isoSurfaceCell
270  (
271  fvm,
272  cellDistance,
274  distance_,
275  regularise_,
276  bounds_
277  )
278  );
279  }
280  else
281  {
282  isoSurfPtr_.reset
283  (
284  new isoSurface
285  (
286  cellDistance,
288  distance_,
289  regularise_,
290  bounds_
291  )
292  );
293  }
294 
295  if (debug)
296  {
297  print(Pout);
298  Pout<< endl;
299  }
300 }
301 
302 
303 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
304 
306 (
307  const word& name,
308  const polyMesh& mesh,
309  const dictionary& dict
310 )
311 :
313  surfPtr_
314  (
316  (
317  dict.lookup("surfaceType"),
318  IOobject
319  (
320  dict.lookupOrDefault("surfaceName", name), // name
321  mesh.time().constant(), // directory
322  "triSurface", // instance
323  mesh.time(), // registry
326  ),
327  dict
328  )
329  ),
330  distance_(readScalar(dict.lookup("distance"))),
331  signed_(readBool(dict.lookup("signed"))),
332  cell_(dict.lookupOrDefault("cell", true)),
333  regularise_(dict.lookupOrDefault("regularise", true)),
334  average_(dict.lookupOrDefault("average", false)),
335  bounds_(dict.lookupOrDefault("bounds", boundBox::greatBox)),
336  zoneKey_(keyType::null),
337  needsUpdate_(true),
338  isoSurfCellPtr_(NULL),
339  isoSurfPtr_(NULL),
340  facesPtr_(NULL)
341 {
342 // dict.readIfPresent("zone", zoneKey_);
343 //
344 // if (debug && zoneKey_.size() && mesh.cellZones().findZoneID(zoneKey_) < 0)
345 // {
346 // Info<< "cellZone " << zoneKey_
347 // << " not found - using entire mesh" << endl;
348 // }
349 }
350 
351 
352 
354 (
355  const word& name,
356  const polyMesh& mesh,
357  const bool interpolate,
358  const word& surfaceType,
359  const word& surfaceName,
360  const scalar distance,
361  const bool signedDistance,
362  const bool cell,
363  const Switch regularise,
364  const Switch average,
365  const boundBox& bounds
366 )
367 :
369  surfPtr_
370  (
372  (
373  surfaceType,
374  IOobject
375  (
376  surfaceName, // name
377  mesh.time().constant(), // directory
378  "triSurface", // instance
379  mesh.time(), // registry
382  ),
383  dictionary()
384  )
385  ),
386  distance_(distance),
387  signed_(signedDistance),
388  cell_(cell),
389  regularise_(regularise),
390  average_(average),
391  bounds_(bounds),
392  zoneKey_(keyType::null),
393  needsUpdate_(true),
394  isoSurfCellPtr_(NULL),
395  isoSurfPtr_(NULL),
396  facesPtr_(NULL)
397 {}
398 
399 
400 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
401 
403 {}
404 
405 
406 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
407 
409 {
410  return needsUpdate_;
411 }
412 
413 
415 {
416  if (debug)
417  {
418  Pout<< "distanceSurface::expire :"
419  << " have-facesPtr_:" << facesPtr_.valid()
420  << " needsUpdate_:" << needsUpdate_ << endl;
421  }
422 
423  // Clear any stored topologies
424  facesPtr_.clear();
425 
426  // Clear derived data
427  clearGeom();
428 
429  // already marked as expired
430  if (needsUpdate_)
431  {
432  return false;
433  }
434 
435  needsUpdate_ = true;
436  return true;
437 }
438 
439 
441 {
442  if (debug)
443  {
444  Pout<< "distanceSurface::update :"
445  << " have-facesPtr_:" << facesPtr_.valid()
446  << " needsUpdate_:" << needsUpdate_ << endl;
447  }
448 
449  if (!needsUpdate_)
450  {
451  return false;
452  }
453 
454  createGeometry();
455 
456  needsUpdate_ = false;
457  return true;
458 }
459 
460 
462 (
463  const volScalarField& vField
464 ) const
465 {
466  return sampleField(vField);
467 }
468 
469 
471 (
472  const volVectorField& vField
473 ) const
474 {
475  return sampleField(vField);
476 }
477 
478 
480 (
481  const volSphericalTensorField& vField
482 ) const
483 {
484  return sampleField(vField);
485 }
486 
487 
489 (
490  const volSymmTensorField& vField
491 ) const
492 {
493  return sampleField(vField);
494 }
495 
496 
498 (
499  const volTensorField& vField
500 ) const
501 {
502  return sampleField(vField);
503 }
504 
505 
507 (
508  const interpolation<scalar>& interpolator
509 ) const
510 {
511  return interpolateField(interpolator);
512 }
513 
514 
516 (
517  const interpolation<vector>& interpolator
518 ) const
519 {
520  return interpolateField(interpolator);
521 }
522 
524 (
525  const interpolation<sphericalTensor>& interpolator
526 ) const
527 {
528  return interpolateField(interpolator);
529 }
530 
531 
533 (
534  const interpolation<symmTensor>& interpolator
535 ) const
536 {
537  return interpolateField(interpolator);
538 }
539 
540 
542 (
543  const interpolation<tensor>& interpolator
544 ) const
545 {
546  return interpolateField(interpolator);
547 }
548 
549 
551 {
552  os << "distanceSurface: " << name() << " :"
553  << " surface:" << surfPtr_().name()
554  << " distance:" << distance_
555  << " faces:" << faces().size()
556  << " points:" << points().size();
557 }
558 
559 
560 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::volumeType::OUTSIDE
@ OUTSIDE
Definition: volumeType.H:64
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
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::volumeType::INSIDE
@ INSIDE
Definition: volumeType.H:63
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::distanceSurface::distance_
const scalar distance_
Distance value.
Definition: distanceSurface.H:63
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::distanceSurface::isoSurfCellPtr_
autoPtr< isoSurfaceCell > isoSurfCellPtr_
Constructed iso surface.
Definition: distanceSurface.H:94
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::distanceSurface::update
virtual bool update()
Update the surface as required.
Definition: distanceSurface.C:440
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::distanceSurface::regularise_
const Switch regularise_
Whether to coarsen.
Definition: distanceSurface.H:72
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::distanceSurface::facesPtr_
autoPtr< faceList > facesPtr_
Triangles converted to faceList.
Definition: distanceSurface.H:100
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
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::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::distanceSurface::signed_
const bool signed_
Signed distance.
Definition: distanceSurface.H:66
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &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::distanceSurface::cellDistancePtr_
autoPtr< volScalarField > cellDistancePtr_
Distance to cell centres.
Definition: distanceSurface.H:88
Foam::distanceSurface::distanceSurface
distanceSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: distanceSurface.C:306
Foam::distanceSurface::cell_
const bool cell_
Whether to use isoSurfaceCell or isoSurface.
Definition: distanceSurface.H:69
Foam::distanceSurface::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: distanceSurface.C:414
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
volumeType.H
Foam::isoSurface
A surface formed by the iso value. After "Regularised Marching Tetrahedra: improved iso-surface extra...
Definition: isoSurface.H:85
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::volumeType
Definition: volumeType.H:54
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::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::distanceSurface::~distanceSurface
virtual ~distanceSurface()
Destructor.
Definition: distanceSurface.C:402
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 >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::distanceSurface::isoSurfPtr_
autoPtr< isoSurface > isoSurfPtr_
Constructed iso surface.
Definition: distanceSurface.H:97
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
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
distanceSurface.H
Foam::distanceSurface::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: distanceSurface.C:408
Foam::searchableSurface::New
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
Definition: searchableSurface.C:40
Foam::distanceSurface::print
virtual void print(Ostream &) const
Write.
Definition: distanceSurface.C:550
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
volPointInterpolation.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::distanceSurface::surfPtr_
const autoPtr< searchableSurface > surfPtr_
Surface.
Definition: distanceSurface.H:60
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::distanceSurface::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: distanceSurface.C:462
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::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
points
const pointField & points
Definition: gmvOutputHeader.H:1
dictionary.H
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Definition: sampledSurface.C:42
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::distanceSurface::createGeometry
void createGeometry()
Create iso surface.
Definition: distanceSurface.C:44
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:244
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::distanceSurface::pointDistance_
scalarField pointDistance_
Distance to points.
Definition: distanceSurface.H:91
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::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:335
Foam::distanceSurface::bounds_
const boundBox bounds_
Optional bounding box to trim triangles against.
Definition: distanceSurface.H:78