distanceSurface.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 "distanceSurface.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace sampledSurfaces
34 {
37 }
38 }
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::sampledSurfaces::distanceSurface::createGeometry()
43 {
44  if (debug)
45  {
46  Pout<< "distanceSurface::createGeometry :updating geometry." << endl;
47  }
48 
49  // Clear any stored topologies
50  isoSurfPtr_.clear();
51 
52  // Clear derived data
53  clearGeom();
54 
55  // Get any subMesh
56  if (zoneKey_.size() && !subMeshPtr_.valid())
57  {
58  const polyBoundaryMesh& patches = mesh().boundaryMesh();
59 
60  // Patch to put exposed internal faces into
61  const label exposedPatchi = patches.findPatchID(exposedPatchName_);
62 
63  subMeshPtr_.reset
64  (
65  new fvMeshSubset(static_cast<const fvMesh&>(mesh()))
66  );
67  subMeshPtr_().setLargeCellSubset
68  (
69  labelHashSet(mesh().cellZones().findMatching(zoneKey_).used()),
70  exposedPatchi
71  );
72 
73  DebugInfo
74  << "Allocating subset of size " << subMeshPtr_().subMesh().nCells()
75  << " with exposed faces into patch "
76  << patches[exposedPatchi].name() << endl;
77  }
78 
79  // Select either the submesh or the underlying mesh
80  const fvMesh& mesh =
81  (
82  subMeshPtr_.valid()
83  ? subMeshPtr_().subMesh()
84  : static_cast<const fvMesh&>(this->mesh())
85  );
86 
87 
88  // Distance to cell centres
89  // ~~~~~~~~~~~~~~~~~~~~~~~~
90 
91  cellDistancePtr_.reset
92  (
93  new volScalarField
94  (
95  IOobject
96  (
97  "cellDistance",
98  mesh.time().timeName(),
99  mesh.time(),
102  false
103  ),
104  mesh,
106  )
107  );
108  volScalarField& cellDistance = cellDistancePtr_();
109 
110  // Internal field
111  {
112  const pointField& cc = mesh.C();
113  scalarField& fld = cellDistance.primitiveFieldRef();
114 
115  List<pointIndexHit> nearest;
116  surfPtr_().findNearest
117  (
118  cc,
119  scalarField(cc.size(), great),
120  nearest
121  );
122 
123  if (signed_)
124  {
125  List<volumeType> volType;
126 
127  surfPtr_().getVolumeType(cc, volType);
128 
129  forAll(volType, i)
130  {
131  volumeType vT = volType[i];
132 
133  if (vT == volumeType::outside)
134  {
135  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
136  }
137  else if (vT == volumeType::inside)
138  {
139  fld[i] = -Foam::mag(cc[i] - nearest[i].hitPoint());
140  }
141  else
142  {
144  << "getVolumeType failure, neither INSIDE or OUTSIDE"
145  << exit(FatalError);
146  }
147  }
148  }
149  else
150  {
151  forAll(nearest, i)
152  {
153  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
154  }
155  }
156  }
157 
158  volScalarField::Boundary& cellDistanceBf =
159  cellDistance.boundaryFieldRef();
160 
161  // Patch fields
162  {
163  forAll(mesh.C().boundaryField(), patchi)
164  {
165  const pointField& cc = mesh.C().boundaryField()[patchi];
166  fvPatchScalarField& fld = cellDistanceBf[patchi];
167 
168  List<pointIndexHit> nearest;
169  surfPtr_().findNearest
170  (
171  cc,
172  scalarField(cc.size(), great),
173  nearest
174  );
175 
176  if (signed_)
177  {
178  List<volumeType> volType;
179 
180  surfPtr_().getVolumeType(cc, volType);
181 
182  forAll(volType, i)
183  {
184  volumeType vT = volType[i];
185 
186  if (vT == volumeType::outside)
187  {
188  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
189  }
190  else if (vT == volumeType::inside)
191  {
192  fld[i] = -Foam::mag(cc[i] - nearest[i].hitPoint());
193  }
194  else
195  {
197  << "getVolumeType failure, "
198  << "neither INSIDE or OUTSIDE"
199  << exit(FatalError);
200  }
201  }
202  }
203  else
204  {
205  forAll(nearest, i)
206  {
207  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
208  }
209  }
210  }
211  }
212 
213 
214  // On processor patches the mesh.C() will already be the cell centre
215  // on the opposite side so no need to swap cellDistance.
216 
217 
218  // Distance to points
219  pointDistance_.setSize(mesh.nPoints());
220  {
221  const pointField& pts = mesh.points();
222 
223  List<pointIndexHit> nearest;
224  surfPtr_().findNearest
225  (
226  pts,
227  scalarField(pts.size(), great),
228  nearest
229  );
230 
231  if (signed_)
232  {
233  List<volumeType> volType;
234 
235  surfPtr_().getVolumeType(pts, volType);
236 
237  forAll(volType, i)
238  {
239  volumeType vT = volType[i];
240 
241  if (vT == volumeType::outside)
242  {
243  pointDistance_[i] =
244  Foam::mag(pts[i] - nearest[i].hitPoint());
245  }
246  else if (vT == volumeType::inside)
247  {
248  pointDistance_[i] =
249  -Foam::mag(pts[i] - nearest[i].hitPoint());
250  }
251  else
252  {
254  << "getVolumeType failure, neither INSIDE or OUTSIDE"
255  << exit(FatalError);
256  }
257  }
258  }
259  else
260  {
261  forAll(nearest, i)
262  {
263  pointDistance_[i] = Foam::mag(pts[i]-nearest[i].hitPoint());
264  }
265  }
266  }
267 
268 
269  if (debug)
270  {
271  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
272  cellDistance.write();
273  pointScalarField pDist
274  (
275  IOobject
276  (
277  "pointDistance",
278  mesh.time().timeName(),
279  mesh.time(),
282  false
283  ),
284  pointMesh::New(mesh),
286  );
287  pDist.primitiveFieldRef() = pointDistance_;
288 
289  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
290  pDist.write();
291  }
292 
293 
294  //- Direct from cell field and point field.
295  isoSurfPtr_.reset
296  (
297  new isoSurface
298  (
299  mesh,
300  cellDistance,
301  pointDistance_,
302  distance_,
303  filter_
304  )
305  );
306 
307  if (debug)
308  {
309  print(Pout);
310  Pout<< endl;
311  }
312 }
313 
314 
315 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
316 
318 (
319  const word& name,
320  const polyMesh& mesh,
321  const dictionary& dict
322 )
323 :
324  sampledSurface(name, mesh, dict),
325  surfPtr_
326  (
328  (
329  dict.lookup("surfaceType"),
330  IOobject
331  (
332  dict.lookupOrDefault("surfaceName", name),
333  mesh.time().constant(),
335  mesh.time(),
338  ),
339  dict
340  )
341  ),
342  distance_(dict.lookup<scalar>("distance")),
343  signed_(readBool(dict.lookup("signed"))),
344  filter_
345  (
346  dict.found("filtering")
347  ? isoSurface::filterTypeNames_.read(dict.lookup("filtering"))
349  ),
350  average_(dict.lookupOrDefault("average", false)),
351  zoneKey_(dict.lookupOrDefault("zone", wordRe::null)),
352  needsUpdate_(true),
353  subMeshPtr_(nullptr),
354  isoSurfPtr_(nullptr)
355 {
356  if (zoneKey_.size())
357  {
358  dict.lookup("exposedPatchName") >> exposedPatchName_;
359 
360  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
361  {
363  << "Cannot find patch " << exposedPatchName_
364  << " in which to put exposed faces." << endl
365  << "Valid patches are " << mesh.boundaryMesh().names()
366  << exit(FatalError);
367  }
368 
369  if (debug)
370  {
371  Info<< "Restricting to cellZone " << zoneKey_
372  << " with exposed internal faces into patch "
373  << exposedPatchName_ << endl;
374  }
375 
376  if (mesh.cellZones().findIndex(zoneKey_) < 0)
377  {
379  << "cellZone " << zoneKey_
380  << " not found - using entire mesh" << endl;
381  }
382  }
383 }
384 
385 
386 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
387 
389 {}
390 
391 
392 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
393 
395 {
396  return needsUpdate_;
397 }
398 
399 
401 {
402  if (debug)
403  {
404  Pout<< "distanceSurface::expire :"
405  << " needsUpdate_:" << needsUpdate_ << endl;
406  }
407 
408  // Clear derived data
409  clearGeom();
410 
411  // already marked as expired
412  if (needsUpdate_)
413  {
414  return false;
415  }
416 
417  needsUpdate_ = true;
418  return true;
419 }
420 
421 
423 {
424  if (debug)
425  {
426  Pout<< "distanceSurface::update :"
427  << " needsUpdate_:" << needsUpdate_ << endl;
428  }
429 
430  if (!needsUpdate_)
431  {
432  return false;
433  }
434 
435  createGeometry();
436 
437  needsUpdate_ = false;
438  return true;
439 }
440 
441 
444 (
445  const volScalarField& vField
446 ) const
447 {
448  return sampleField(vField);
449 }
450 
451 
454 (
455  const volVectorField& vField
456 ) const
457 {
458  return sampleField(vField);
459 }
460 
461 
464 (
465  const volSphericalTensorField& vField
466 ) const
467 {
468  return sampleField(vField);
469 }
470 
471 
474 (
475  const volSymmTensorField& vField
476 ) const
477 {
478  return sampleField(vField);
479 }
480 
481 
484 (
485  const volTensorField& vField
486 ) const
487 {
488  return sampleField(vField);
489 }
490 
491 
494 (
495  const interpolation<scalar>& interpolator
496 ) const
497 {
498  return interpolateField(interpolator);
499 }
500 
501 
504 (
505  const interpolation<vector>& interpolator
506 ) const
507 {
508  return interpolateField(interpolator);
509 }
510 
513 (
514  const interpolation<sphericalTensor>& interpolator
515 ) const
516 {
517  return interpolateField(interpolator);
518 }
519 
520 
523 (
524  const interpolation<symmTensor>& interpolator
525 ) const
526 {
527  return interpolateField(interpolator);
528 }
529 
530 
533 (
534  const interpolation<tensor>& interpolator
535 ) const
536 {
537  return interpolateField(interpolator);
538 }
539 
540 
542 {
543  os << "distanceSurface: " << name() << " :"
544  << " surface:" << surfPtr_().name()
545  << " distance:" << distance_
546  << " faces:" << faces().size()
547  << " points:" << points().size();
548 }
549 
550 
551 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
const word & name() const
Return name.
Definition: IOobject.H:303
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:482
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:323
An abstract class for surfaces with sampling.
defineTypeNameAndDebug(distanceSurface, 0)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const NamedEnum< filterType, 4 > filterTypeNames_
Definition: isoSurface.H:92
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
virtual bool update()
Update the surface as required.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
bool interpolate() const
Interpolation requested for surface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
label findPatchID(const word &patchName) const
Find patch index given a name.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
bool readBool(Istream &)
Definition: boolIO.C:60
patches[0]
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
A sampledSurface defined by a distance to a surface.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
virtual void print(Ostream &) const
Write.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
dynamicFvMesh & mesh
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))
const pointField & points
static const word & geometryDir()
Return the geometry directory name.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList names() const
Return a list of patch names.
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool expire()
Mark the surface as needing an update.
A sampledSurface defined by a surface of iso value.
#define DebugInfo
Report an information message using Foam::Info.
Foam::polyBoundaryMesh.
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const Time & time() const
Return time.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
distanceSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
label findIndex(const wordRe &) const
Return zone index for the first match, return -1 if not found.
Definition: MeshZones.C:306
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label nPoints() const
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const volVectorField & C() const
Return cell centres as volVectorField.
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:419
static const wordRe null
An empty wordRe.
Definition: wordRe.H:88
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844