sampledCuttingPlane.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 "sampledCuttingPlane.H"
27 #include "dictionary.H"
28 #include "volFields.H"
29 #include "volPointInterpolation.H"
31 #include "fvMesh.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(sampledCuttingPlane, 0);
39  (
40  sampledSurface,
41  sampledCuttingPlane,
42  word,
43  cuttingPlane
44  );
45 }
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
50 {
51  if (debug)
52  {
53  Pout<< "sampledCuttingPlane::createGeometry :updating geometry."
54  << endl;
55  }
56 
57  // Clear any stored topologies
58  facesPtr_.clear();
59  isoSurfPtr_.ptr();
60  pointDistance_.clear();
61  cellDistancePtr_.clear();
62 
63  // Clear derived data
64  clearGeom();
65 
66  // Get any subMesh
67  if (zoneID_.index() != -1 && !subMeshPtr_.valid())
68  {
70 
71  // Patch to put exposed internal faces into
72  const label exposedPatchI = patches.findPatchID(exposedPatchName_);
73 
74  if (debug)
75  {
76  Info<< "Allocating subset of size "
77  << mesh().cellZones()[zoneID_.index()].size()
78  << " with exposed faces into patch "
79  << patches[exposedPatchI].name() << endl;
80  }
81 
82  subMeshPtr_.reset
83  (
84  new fvMeshSubset(static_cast<const fvMesh&>(mesh()))
85  );
86  subMeshPtr_().setLargeCellSubset
87  (
88  labelHashSet(mesh().cellZones()[zoneID_.index()]),
89  exposedPatchI
90  );
91  }
92 
93 
94  // Select either the submesh or the underlying mesh
95  const fvMesh& fvm =
96  (
97  subMeshPtr_.valid()
98  ? subMeshPtr_().subMesh()
99  : static_cast<const fvMesh&>(mesh())
100  );
101 
102 
103  // Distance to cell centres
104  // ~~~~~~~~~~~~~~~~~~~~~~~~
105 
106  cellDistancePtr_.reset
107  (
108  new volScalarField
109  (
110  IOobject
111  (
112  "cellDistance",
113  fvm.time().timeName(),
114  fvm.time(),
117  false
118  ),
119  fvm,
120  dimensionedScalar("zero", dimLength, 0)
121  )
122  );
123  volScalarField& cellDistance = cellDistancePtr_();
124 
125  // Internal field
126  {
127  const pointField& cc = fvm.cellCentres();
128  scalarField& fld = cellDistance.internalField();
129 
130  forAll(cc, i)
131  {
132  // Signed distance
133  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
134  }
135  }
136 
137  // Patch fields
138  {
139  forAll(cellDistance.boundaryField(), patchI)
140  {
141  if
142  (
143  isA<emptyFvPatchScalarField>
144  (
145  cellDistance.boundaryField()[patchI]
146  )
147  )
148  {
149  cellDistance.boundaryField().set
150  (
151  patchI,
152  new calculatedFvPatchScalarField
153  (
154  fvm.boundary()[patchI],
155  cellDistance
156  )
157  );
158 
159  const polyPatch& pp = fvm.boundary()[patchI].patch();
161 
162  fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
163  fld.setSize(pp.size());
164  forAll(fld, i)
165  {
166  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
167  }
168  }
169  else
170  {
171  const pointField& cc = fvm.C().boundaryField()[patchI];
172  fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
173 
174  forAll(fld, i)
175  {
176  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
177  }
178  }
179  }
180  }
181 
182 
183  // On processor patches the mesh.C() will already be the cell centre
184  // on the opposite side so no need to swap cellDistance.
185 
186 
187  // Distance to points
188  pointDistance_.setSize(fvm.nPoints());
189  {
190  const pointField& pts = fvm.points();
191 
193  {
194  pointDistance_[i] = (pts[i] - plane_.refPoint()) & plane_.normal();
195  }
196  }
197 
198 
199  if (debug)
200  {
201  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
202  cellDistance.write();
203  pointScalarField pDist
204  (
205  IOobject
206  (
207  "pointDistance",
208  fvm.time().timeName(),
209  fvm.time(),
212  false
213  ),
214  pointMesh::New(fvm),
215  dimensionedScalar("zero", dimLength, 0)
216  );
217  pDist.internalField() = pointDistance_;
218 
219  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
220  pDist.write();
221  }
222 
223 
224  //- Direct from cell field and point field.
225  isoSurfPtr_.reset
226  (
227  new isoSurface
228  (
229  cellDistance,
231  0.0,
232  regularise_,
233  bounds_,
234  mergeTol_
235  )
236  //new isoSurfaceCell
237  //(
238  // fvm,
239  // cellDistance,
240  // pointDistance_,
241  // 0.0,
242  // regularise_,
243  // mergeTol_
244  //)
245  );
246 
247  if (debug)
248  {
249  print(Pout);
250  Pout<< endl;
251  }
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
256 
258 (
259  const word& name,
260  const polyMesh& mesh,
261  const dictionary& dict
262 )
263 :
265  plane_(dict),
266  bounds_(dict.lookupOrDefault("bounds", boundBox::greatBox)),
267  mergeTol_(dict.lookupOrDefault("mergeTol", 1e-6)),
268  regularise_(dict.lookupOrDefault("regularise", true)),
269  average_(dict.lookupOrDefault("average", false)),
270  zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
271  exposedPatchName_(word::null),
272  needsUpdate_(true),
273  subMeshPtr_(NULL),
274  cellDistancePtr_(NULL),
275  isoSurfPtr_(NULL),
276  facesPtr_(NULL)
277 {
278  if (zoneID_.index() != -1)
279  {
280  dict.lookup("exposedPatchName") >> exposedPatchName_;
281 
282  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
283  {
285  << "Cannot find patch " << exposedPatchName_
286  << " in which to put exposed faces." << endl
287  << "Valid patches are " << mesh.boundaryMesh().names()
288  << exit(FatalError);
289  }
290 
291  if (debug && zoneID_.index() != -1)
292  {
293  Info<< "Restricting to cellZone " << zoneID_.name()
294  << " with exposed internal faces into patch "
295  << exposedPatchName_ << endl;
296  }
297  }
298 }
299 
300 
301 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
302 
304 {}
305 
306 
307 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
308 
310 {
311  return needsUpdate_;
312 }
313 
314 
316 {
317  if (debug)
318  {
319  Pout<< "sampledCuttingPlane::expire :"
320  << " have-facesPtr_:" << facesPtr_.valid()
321  << " needsUpdate_:" << needsUpdate_ << endl;
322  }
323 
324  // Clear any stored topologies
325  facesPtr_.clear();
326 
327  // Clear derived data
328  clearGeom();
329 
330  // already marked as expired
331  if (needsUpdate_)
332  {
333  return false;
334  }
335 
336  needsUpdate_ = true;
337  return true;
338 }
339 
340 
342 {
343  if (debug)
344  {
345  Pout<< "sampledCuttingPlane::update :"
346  << " have-facesPtr_:" << facesPtr_.valid()
347  << " needsUpdate_:" << needsUpdate_ << endl;
348  }
349 
350  if (!needsUpdate_)
351  {
352  return false;
353  }
354 
355  createGeometry();
356 
357  needsUpdate_ = false;
358  return true;
359 }
360 
361 
364 (
365  const volScalarField& vField
366 ) const
367 {
368  return sampleField(vField);
369 }
370 
371 
374 (
375  const volVectorField& vField
376 ) const
377 {
378  return sampleField(vField);
379 }
380 
381 
384 (
385  const volSphericalTensorField& vField
386 ) const
387 {
388  return sampleField(vField);
389 }
390 
391 
394 (
395  const volSymmTensorField& vField
396 ) const
397 {
398  return sampleField(vField);
399 }
400 
401 
404 (
405  const volTensorField& vField
406 ) const
407 {
408  return sampleField(vField);
409 }
410 
411 
414 (
415  const interpolation<scalar>& interpolator
416 ) const
417 {
418  return interpolateField(interpolator);
419 }
420 
421 
424 (
425  const interpolation<vector>& interpolator
426 ) const
427 {
428  return interpolateField(interpolator);
429 }
430 
433 (
434  const interpolation<sphericalTensor>& interpolator
435 ) const
436 {
437  return interpolateField(interpolator);
438 }
439 
440 
443 (
444  const interpolation<symmTensor>& interpolator
445 ) const
446 {
447  return interpolateField(interpolator);
448 }
449 
450 
453 (
454  const interpolation<tensor>& interpolator
455 ) const
456 {
457  return interpolateField(interpolator);
458 }
459 
460 
462 {
463  os << "sampledCuttingPlane: " << name() << " :"
464  << " plane:" << plane_
465  << " faces:" << faces().size()
466  << " points:" << points().size();
467 }
468 
469 
470 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::sampledCuttingPlane::cellDistancePtr_
autoPtr< volScalarField > cellDistancePtr_
Distance to cell centres.
Definition: sampledCuttingPlane.H:88
volFields.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
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::plane::normal
const vector & normal() const
Return plane normal.
Definition: plane.C:248
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch, cyclicPeriodicAMI)
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::sampledCuttingPlane::isoSurfPtr_
autoPtr< isoSurface > isoSurfPtr_
Constructed iso surface.
Definition: sampledCuttingPlane.H:95
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
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::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::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::sampledCuttingPlane::regularise_
const Switch regularise_
Whether to coarsen.
Definition: sampledCuttingPlane.H:69
Foam::sampledCuttingPlane::print
virtual void print(Ostream &) const
Write.
Definition: sampledCuttingPlane.C:461
Foam::sampledCuttingPlane::createGeometry
void createGeometry()
Create iso surface.
Definition: sampledCuttingPlane.C:49
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::sampledCuttingPlane::plane_
const plane plane_
Plane.
Definition: sampledCuttingPlane.H:60
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::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:527
Foam::SubField
Pre-declare related SubField type.
Definition: Field.H:61
sampledCuttingPlane.H
Foam::sampledCuttingPlane::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledCuttingPlane.C:309
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::sampledCuttingPlane::bounds_
const boundBox bounds_
Optional bounding box to trim triangles against.
Definition: sampledCuttingPlane.H:63
Foam::isoSurface
A surface formed by the iso value. After "Regularised Marching Tetrahedra: improved iso-surface extra...
Definition: isoSurface.H:85
Foam::label
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
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::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::sampledCuttingPlane::mergeTol_
const scalar mergeTol_
Merge tolerance.
Definition: sampledCuttingPlane.H:66
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::sampledCuttingPlane::update
virtual bool update()
Update the surface as required.
Definition: sampledCuttingPlane.C:341
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::sampledCuttingPlane::~sampledCuttingPlane
virtual ~sampledCuttingPlane()
Destructor.
Definition: sampledCuttingPlane.C:303
Foam::interpolation< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
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::sampledCuttingPlane::subMeshPtr_
autoPtr< fvMeshSubset > subMeshPtr_
Optional subsetted mesh.
Definition: sampledCuttingPlane.H:85
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
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::sampledCuttingPlane::sampledCuttingPlane
sampledCuttingPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledCuttingPlane.C:258
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::sampledCuttingPlane::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledCuttingPlane.C:364
volPointInterpolation.H
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sampledCuttingPlane::zoneID_
cellZoneID zoneID_
Zone name/index (if restricted to zones)
Definition: sampledCuttingPlane.H:75
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::polyPatch::patchSlice
const List< T >::subList patchSlice(const UList< T > &l) const
Slice list to patch.
Definition: polyPatch.H:345
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::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::sampledCuttingPlane::pointDistance_
scalarField pointDistance_
Distance to points.
Definition: sampledCuttingPlane.H:91
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:130
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::plane::refPoint
const point & refPoint() const
Return or return plane base point.
Definition: plane.C:255
dictionary.H
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Definition: sampledSurface.C:42
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:244
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::sampledCuttingPlane::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledCuttingPlane.C:315
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation requested for surface.
Definition: sampledSurface.H:256
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::sampledCuttingPlane::facesPtr_
autoPtr< faceList > facesPtr_
Triangles converted to faceList.
Definition: sampledCuttingPlane.H:98
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::DynamicID::index
label index() const
Return index of first matching zone.
Definition: DynamicID.H:107
Foam::sampledCuttingPlane::exposedPatchName_
word exposedPatchName_
For zones: patch to put exposed faces into.
Definition: sampledCuttingPlane.H:78