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-2016 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 
44 void Foam::distanceSurface::createGeometry()
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.primitiveFieldRef();
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  volScalarField::Boundary& cellDistanceBf =
132  cellDistance.boundaryFieldRef();
133 
134  // Patch fields
135  {
136  forAll(fvm.C().boundaryField(), patchi)
137  {
138  const pointField& cc = fvm.C().boundaryField()[patchi];
139  fvPatchScalarField& fld = cellDistanceBf[patchi];
140 
141  List<pointIndexHit> nearest;
142  surfPtr_().findNearest
143  (
144  cc,
145  scalarField(cc.size(), GREAT),
146  nearest
147  );
148 
149  if (signed_)
150  {
151  List<volumeType> volType;
152 
153  surfPtr_().getVolumeType(cc, volType);
154 
155  forAll(volType, i)
156  {
157  volumeType vT = volType[i];
158 
159  if (vT == volumeType::OUTSIDE)
160  {
161  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
162  }
163  else if (vT == volumeType::INSIDE)
164  {
165  fld[i] = -Foam::mag(cc[i] - nearest[i].hitPoint());
166  }
167  else
168  {
170  << "getVolumeType failure, "
171  << "neither INSIDE or OUTSIDE"
172  << exit(FatalError);
173  }
174  }
175  }
176  else
177  {
178  forAll(nearest, i)
179  {
180  fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
181  }
182  }
183  }
184  }
185 
186 
187  // On processor patches the mesh.C() will already be the cell centre
188  // on the opposite side so no need to swap cellDistance.
189 
190 
191  // Distance to points
192  pointDistance_.setSize(fvm.nPoints());
193  {
194  const pointField& pts = fvm.points();
195 
196  List<pointIndexHit> nearest;
197  surfPtr_().findNearest
198  (
199  pts,
200  scalarField(pts.size(), GREAT),
201  nearest
202  );
203 
204  if (signed_)
205  {
206  List<volumeType> volType;
207 
208  surfPtr_().getVolumeType(pts, volType);
209 
210  forAll(volType, i)
211  {
212  volumeType vT = volType[i];
213 
214  if (vT == volumeType::OUTSIDE)
215  {
216  pointDistance_[i] =
217  Foam::mag(pts[i] - nearest[i].hitPoint());
218  }
219  else if (vT == volumeType::INSIDE)
220  {
221  pointDistance_[i] =
222  -Foam::mag(pts[i] - nearest[i].hitPoint());
223  }
224  else
225  {
227  << "getVolumeType failure, neither INSIDE or OUTSIDE"
228  << exit(FatalError);
229  }
230  }
231  }
232  else
233  {
234  forAll(nearest, i)
235  {
236  pointDistance_[i] = Foam::mag(pts[i]-nearest[i].hitPoint());
237  }
238  }
239  }
240 
241 
242  if (debug)
243  {
244  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
245  cellDistance.write();
246  pointScalarField pDist
247  (
248  IOobject
249  (
250  "pointDistance",
251  fvm.time().timeName(),
252  fvm.time(),
255  false
256  ),
257  pointMesh::New(fvm),
258  dimensionedScalar("zero", dimLength, 0)
259  );
260  pDist.primitiveFieldRef() = pointDistance_;
261 
262  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
263  pDist.write();
264  }
265 
266 
267  //- Direct from cell field and point field.
268  if (cell_)
269  {
270  isoSurfCellPtr_.reset
271  (
272  new isoSurfaceCell
273  (
274  fvm,
275  cellDistance,
276  pointDistance_,
277  distance_,
278  regularise_
279  )
280  );
281  }
282  else
283  {
284  isoSurfPtr_.reset
285  (
286  new isoSurface
287  (
288  cellDistance,
289  pointDistance_,
290  distance_,
291  regularise_
292  )
293  );
294  }
295 
296  if (debug)
297  {
298  print(Pout);
299  Pout<< endl;
300  }
301 }
302 
303 
304 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
305 
307 (
308  const word& name,
309  const polyMesh& mesh,
310  const dictionary& dict
311 )
312 :
313  sampledSurface(name, mesh, dict),
314  surfPtr_
315  (
317  (
318  dict.lookup("surfaceType"),
319  IOobject
320  (
321  dict.lookupOrDefault("surfaceName", name), // name
322  mesh.time().constant(), // directory
323  "triSurface", // instance
324  mesh.time(), // registry
327  ),
328  dict
329  )
330  ),
331  distance_(readScalar(dict.lookup("distance"))),
332  signed_(readBool(dict.lookup("signed"))),
333  cell_(dict.lookupOrDefault("cell", true)),
334  regularise_(dict.lookupOrDefault("regularise", true)),
335  average_(dict.lookupOrDefault("average", false)),
336  zoneKey_(keyType::null),
337  needsUpdate_(true),
338  isoSurfCellPtr_(nullptr),
339  isoSurfPtr_(nullptr),
340  facesPtr_(nullptr)
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 )
366 :
367  sampledSurface(name, mesh, interpolate),
368  surfPtr_
369  (
371  (
372  surfaceType,
373  IOobject
374  (
375  surfaceName, // name
376  mesh.time().constant(), // directory
377  "triSurface", // instance
378  mesh.time(), // registry
381  ),
382  dictionary()
383  )
384  ),
385  distance_(distance),
386  signed_(signedDistance),
387  cell_(cell),
388  regularise_(regularise),
389  average_(average),
390  zoneKey_(keyType::null),
391  needsUpdate_(true),
392  isoSurfCellPtr_(nullptr),
393  isoSurfPtr_(nullptr),
394  facesPtr_(nullptr)
395 {}
396 
397 
398 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
399 
401 {}
402 
403 
404 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
405 
407 {
408  return needsUpdate_;
409 }
410 
411 
413 {
414  if (debug)
415  {
416  Pout<< "distanceSurface::expire :"
417  << " have-facesPtr_:" << facesPtr_.valid()
418  << " needsUpdate_:" << needsUpdate_ << endl;
419  }
420 
421  // Clear any stored topologies
422  facesPtr_.clear();
423 
424  // Clear derived data
425  clearGeom();
426 
427  // already marked as expired
428  if (needsUpdate_)
429  {
430  return false;
431  }
432 
433  needsUpdate_ = true;
434  return true;
435 }
436 
437 
439 {
440  if (debug)
441  {
442  Pout<< "distanceSurface::update :"
443  << " have-facesPtr_:" << facesPtr_.valid()
444  << " needsUpdate_:" << needsUpdate_ << endl;
445  }
446 
447  if (!needsUpdate_)
448  {
449  return false;
450  }
451 
452  createGeometry();
453 
454  needsUpdate_ = false;
455  return true;
456 }
457 
458 
460 (
461  const volScalarField& vField
462 ) const
463 {
464  return sampleField(vField);
465 }
466 
467 
469 (
470  const volVectorField& vField
471 ) const
472 {
473  return sampleField(vField);
474 }
475 
476 
478 (
479  const volSphericalTensorField& vField
480 ) const
481 {
482  return sampleField(vField);
483 }
484 
485 
487 (
488  const volSymmTensorField& vField
489 ) const
490 {
491  return sampleField(vField);
492 }
493 
494 
496 (
497  const volTensorField& vField
498 ) const
499 {
500  return sampleField(vField);
501 }
502 
503 
505 (
506  const interpolation<scalar>& interpolator
507 ) const
508 {
509  return interpolateField(interpolator);
510 }
511 
512 
514 (
515  const interpolation<vector>& interpolator
516 ) const
517 {
518  return interpolateField(interpolator);
519 }
520 
522 (
523  const interpolation<sphericalTensor>& interpolator
524 ) const
525 {
526  return interpolateField(interpolator);
527 }
528 
529 
531 (
532  const interpolation<symmTensor>& interpolator
533 ) const
534 {
535  return interpolateField(interpolator);
536 }
537 
538 
540 (
541  const interpolation<tensor>& interpolator
542 ) const
543 {
544  return interpolateField(interpolator);
545 }
546 
547 
549 {
550  os << "distanceSurface: " << name() << " :"
551  << " surface:" << surfPtr_().name()
552  << " distance:" << distance_
553  << " faces:" << faces().size()
554  << " points:" << points().size();
555 }
556 
557 
558 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
An abstract class for surfaces with sampling.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
bool interpolate() const
Interpolation requested for surface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
bool readBool(Istream &)
Definition: boolIO.C:60
virtual ~distanceSurface()
Destructor.
virtual bool update()
Update the surface as required.
Macros for easy insertion into run-time selection tables.
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void print(Ostream &) const
Write.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
const pointField & points
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
distanceSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
virtual bool expire()
Mark the surface as needing an update.
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual bool needsUpdate() const
Does the surface need an update?
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
static const keyType null
An empty keyType.
Definition: keyType.H:83
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.