sampledCuttingPlane.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-2019 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace sampledSurfaces
34 {
35  defineTypeNameAndDebug(cuttingPlane, 0);
36  addToRunTimeSelectionTable(sampledSurface, cuttingPlane, word);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::sampledSurfaces::cuttingPlane::createGeometry()
44 {
45  if (debug)
46  {
47  Pout<< "cuttingPlane::createGeometry :updating geometry."
48  << endl;
49  }
50 
51  // Clear any stored topologies
52  facesPtr_.clear();
53  isoSurfPtr_.ptr();
54  pointDistance_.clear();
55  cellDistancePtr_.clear();
56 
57  // Clear derived data
58  clearGeom();
59 
60  // Get any subMesh
61  if (zoneID_.index() != -1 && !subMeshPtr_.valid())
62  {
63  const polyBoundaryMesh& patches = mesh().boundaryMesh();
64 
65  // Patch to put exposed internal faces into
66  const label exposedPatchi = patches.findPatchID(exposedPatchName_);
67 
68  DebugInfo
69  << "Allocating subset of size "
70  << mesh().cellZones()[zoneID_.index()].size()
71  << " with exposed faces into patch "
72  << patches[exposedPatchi].name() << endl;
73 
74  subMeshPtr_.reset
75  (
76  new fvMeshSubset(static_cast<const fvMesh&>(mesh()))
77  );
78  subMeshPtr_().setLargeCellSubset
79  (
80  labelHashSet(mesh().cellZones()[zoneID_.index()]),
81  exposedPatchi
82  );
83  }
84 
85 
86  // Select either the submesh or the underlying mesh
87  const fvMesh& mesh =
88  (
89  subMeshPtr_.valid()
90  ? subMeshPtr_().subMesh()
91  : static_cast<const fvMesh&>(this->mesh())
92  );
93 
94 
95  // Distance to cell centres
96  // ~~~~~~~~~~~~~~~~~~~~~~~~
97 
98  cellDistancePtr_.reset
99  (
100  new volScalarField
101  (
102  IOobject
103  (
104  "cellDistance",
105  mesh.time().timeName(),
106  mesh.time(),
109  false
110  ),
111  mesh,
113  )
114  );
115  volScalarField& cellDistance = cellDistancePtr_();
116 
117  // Internal field
118  {
119  const pointField& cc = mesh.cellCentres();
120  scalarField& fld = cellDistance.primitiveFieldRef();
121 
122  forAll(cc, i)
123  {
124  // Signed distance
125  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
126  }
127  }
128 
129  volScalarField::Boundary& cellDistanceBf =
130  cellDistance.boundaryFieldRef();
131 
132  // Patch fields
133  {
134  forAll(cellDistanceBf, patchi)
135  {
136  if
137  (
138  isA<emptyFvPatchScalarField>
139  (
140  cellDistanceBf[patchi]
141  )
142  )
143  {
144  cellDistanceBf.set
145  (
146  patchi,
147  new calculatedFvPatchScalarField
148  (
149  mesh.boundary()[patchi],
150  cellDistance
151  )
152  );
153 
154  const polyPatch& pp = mesh.boundary()[patchi].patch();
155  pointField::subField cc = pp.patchSlice(mesh.faceCentres());
156 
157  fvPatchScalarField& fld = cellDistanceBf[patchi];
158  fld.setSize(pp.size());
159  forAll(fld, i)
160  {
161  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
162  }
163  }
164  else
165  {
166  const pointField& cc = mesh.C().boundaryField()[patchi];
167  fvPatchScalarField& fld = cellDistanceBf[patchi];
168 
169  forAll(fld, i)
170  {
171  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
172  }
173  }
174  }
175  }
176 
177 
178  // On processor patches the mesh.C() will already be the cell centre
179  // on the opposite side so no need to swap cellDistance.
180 
181 
182  // Distance to points
183  pointDistance_.setSize(mesh.nPoints());
184  {
185  const pointField& pts = mesh.points();
186 
187  forAll(pointDistance_, i)
188  {
189  pointDistance_[i] = (pts[i] - plane_.refPoint()) & plane_.normal();
190  }
191  }
192 
193 
194  if (debug)
195  {
196  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
197  cellDistance.write();
198  pointScalarField pDist
199  (
200  IOobject
201  (
202  "pointDistance",
203  mesh.time().timeName(),
204  mesh.time(),
207  false
208  ),
209  pointMesh::New(mesh),
211  );
212  pDist.primitiveFieldRef() = pointDistance_;
213 
214  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
215  pDist.write();
216  }
217 
218 
219  //- Direct from cell field and point field.
220  isoSurfPtr_.reset
221  (
222  new isoSurface
223  (
224  mesh,
225  cellDistance,
226  pointDistance_,
227  0,
228  regularise_ ? isoSurface::DIAGCELL : isoSurface::NONE
229  )
230  );
231 
232  if (debug)
233  {
234  print(Pout);
235  Pout<< endl;
236  }
237 }
238 
239 
240 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
241 
243 (
244  const word& name,
245  const polyMesh& mesh,
246  const dictionary& dict
247 )
248 :
249  sampledSurface(name, mesh, dict),
250  plane_(dict),
251  regularise_(dict.lookupOrDefault("regularise", true)),
252  average_(dict.lookupOrDefault("average", false)),
253  zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
254  exposedPatchName_(word::null),
255  needsUpdate_(true),
256  subMeshPtr_(nullptr),
257  cellDistancePtr_(nullptr),
258  isoSurfPtr_(nullptr)
259 {
260  if (zoneID_.index() != -1)
261  {
262  dict.lookup("exposedPatchName") >> exposedPatchName_;
263 
264  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
265  {
267  << "Cannot find patch " << exposedPatchName_
268  << " in which to put exposed faces." << endl
269  << "Valid patches are " << mesh.boundaryMesh().names()
270  << exit(FatalError);
271  }
272 
273  if (debug && zoneID_.index() != -1)
274  {
275  Info<< "Restricting to cellZone " << zoneID_.name()
276  << " with exposed internal faces into patch "
277  << exposedPatchName_ << endl;
278  }
279  }
280 }
281 
282 
283 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
284 
286 {}
287 
288 
289 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
290 
292 {
293  return needsUpdate_;
294 }
295 
296 
298 {
299  if (debug)
300  {
301  Pout<< "cuttingPlane::expire :"
302  << " needsUpdate_:" << needsUpdate_ << endl;
303  }
304 
305  // Clear derived data
306  clearGeom();
307 
308  // Already marked as expired
309  if (needsUpdate_)
310  {
311  return false;
312  }
313 
314  needsUpdate_ = true;
315  return true;
316 }
317 
318 
320 {
321  if (debug)
322  {
323  Pout<< "cuttingPlane::update :"
324  << " needsUpdate_:" << needsUpdate_ << endl;
325  }
326 
327  if (!needsUpdate_)
328  {
329  return false;
330  }
331 
332  createGeometry();
333 
334  needsUpdate_ = false;
335  return true;
336 }
337 
338 
341 (
342  const volScalarField& vField
343 ) const
344 {
345  return sampleField(vField);
346 }
347 
348 
351 (
352  const volVectorField& vField
353 ) const
354 {
355  return sampleField(vField);
356 }
357 
358 
361 (
362  const volSphericalTensorField& vField
363 ) const
364 {
365  return sampleField(vField);
366 }
367 
368 
371 (
372  const volSymmTensorField& vField
373 ) const
374 {
375  return sampleField(vField);
376 }
377 
378 
381 (
382  const volTensorField& vField
383 ) const
384 {
385  return sampleField(vField);
386 }
387 
388 
391 (
392  const interpolation<scalar>& interpolator
393 ) const
394 {
395  return interpolateField(interpolator);
396 }
397 
398 
401 (
402  const interpolation<vector>& interpolator
403 ) const
404 {
405  return interpolateField(interpolator);
406 }
407 
410 (
411  const interpolation<sphericalTensor>& interpolator
412 ) const
413 {
414  return interpolateField(interpolator);
415 }
416 
417 
420 (
421  const interpolation<symmTensor>& interpolator
422 ) const
423 {
424  return interpolateField(interpolator);
425 }
426 
427 
430 (
431  const interpolation<tensor>& interpolator
432 ) const
433 {
434  return interpolateField(interpolator);
435 }
436 
437 
439 {
440  os << "cuttingPlane: " << name() << " :"
441  << " plane:" << plane_
442  << " faces:" << faces().size()
443  << " points:" << points().size();
444 }
445 
446 
447 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
cuttingPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const word & name() const
Return the name of this functionObject.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
An abstract class for surfaces with sampling.
virtual void print(Ostream &) const
Write.
defineTypeNameAndDebug(distanceSurface, 0)
virtual bool update()
Update the surface as required.
bool interpolate() const
Interpolation requested for surface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
label findPatchID(const word &patchName) const
Find patch index given a name.
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
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:208
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 cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:482
wordList names() const
Return a list of patch names.
static const word null
An empty word.
Definition: word.H:77
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
#define DebugInfo
Report an information message using Foam::Info.
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:53
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
virtual bool needsUpdate() const
Does the surface need an update?
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
messageStream Info
SubField< vector > subField
Declare type of subField.
Definition: Field.H:100
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
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
virtual bool expire()
Mark the surface as needing an update.