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-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 "sampledCuttingPlane.H"
27 #include "emptyFvPatchFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace sampledSurfaces
35 {
36  defineTypeNameAndDebug(cuttingPlane, 0);
37  addToRunTimeSelectionTable(sampledSurface, cuttingPlane, word);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::sampledSurfaces::cuttingPlane::createGeometry()
45 {
46  if (debug)
47  {
48  Pout<< "cuttingPlane::createGeometry :updating geometry."
49  << endl;
50  }
51 
52  // Clear any stored topologies
53  facesPtr_.clear();
54  isoSurfPtr_.clear();
55  pointDistance_.clear();
56  cellDistancePtr_.clear();
57 
58  // Clear derived data
59  clearGeom();
60 
61  // Get any subMesh
62  if (zoneKey_.size() && !subMeshPtr_.valid())
63  {
64  const polyBoundaryMesh& patches = mesh().boundaryMesh();
65 
66  // Patch to put exposed internal faces into
67  const label exposedPatchi = patches.findPatchID(exposedPatchName_);
68 
69  subMeshPtr_.reset
70  (
71  new fvMeshSubset(static_cast<const fvMesh&>(mesh()))
72  );
73  subMeshPtr_().setLargeCellSubset
74  (
75  labelHashSet(mesh().cellZones().findMatching(zoneKey_).used()),
76  exposedPatchi
77  );
78 
79  DebugInfo
80  << "Allocating subset of size " << subMeshPtr_().subMesh().nCells()
81  << " with exposed faces into patch "
82  << patches[exposedPatchi].name() << endl;
83  }
84 
85  // Select either the submesh or the underlying mesh
86  const fvMesh& mesh =
87  (
88  subMeshPtr_.valid()
89  ? subMeshPtr_().subMesh()
90  : static_cast<const fvMesh&>(this->mesh())
91  );
92 
93 
94  // Distance to cell centres
95  // ~~~~~~~~~~~~~~~~~~~~~~~~
96 
97  cellDistancePtr_.reset
98  (
99  new volScalarField
100  (
101  IOobject
102  (
103  "cellDistance",
104  mesh.time().timeName(),
105  mesh.time(),
108  false
109  ),
110  mesh,
112  )
113  );
114  volScalarField& cellDistance = cellDistancePtr_();
115 
116  // Internal field
117  {
118  const pointField& cc = mesh.cellCentres();
119  scalarField& fld = cellDistance.primitiveFieldRef();
120 
121  forAll(cc, i)
122  {
123  // Signed distance
124  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
125  }
126  }
127 
128  volScalarField::Boundary& cellDistanceBf =
129  cellDistance.boundaryFieldRef();
130 
131  // Patch fields
132  {
133  forAll(cellDistanceBf, patchi)
134  {
135  if
136  (
137  isA<emptyFvPatchScalarField>
138  (
139  cellDistanceBf[patchi]
140  )
141  )
142  {
143  cellDistanceBf.set
144  (
145  patchi,
146  new calculatedFvPatchScalarField
147  (
148  mesh.boundary()[patchi],
149  cellDistance
150  )
151  );
152 
153  const polyPatch& pp = mesh.boundary()[patchi].patch();
154  pointField::subField cc = pp.patchSlice(mesh.faceCentres());
155 
156  fvPatchScalarField& fld = cellDistanceBf[patchi];
157  fld.setSize(pp.size());
158  forAll(fld, i)
159  {
160  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
161  }
162  }
163  else
164  {
165  const pointField& cc = mesh.C().boundaryField()[patchi];
166  fvPatchScalarField& fld = cellDistanceBf[patchi];
167 
168  forAll(fld, i)
169  {
170  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
171  }
172  }
173  }
174  }
175 
176 
177  // On processor patches the mesh.C() will already be the cell centre
178  // on the opposite side so no need to swap cellDistance.
179 
180 
181  // Distance to points
182  pointDistance_.setSize(mesh.nPoints());
183  {
184  const pointField& pts = mesh.points();
185 
186  forAll(pointDistance_, i)
187  {
188  pointDistance_[i] = (pts[i] - plane_.refPoint()) & plane_.normal();
189  }
190  }
191 
192 
193  if (debug)
194  {
195  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
196  cellDistance.write();
197  pointScalarField pDist
198  (
199  IOobject
200  (
201  "pointDistance",
202  mesh.time().timeName(),
203  mesh.time(),
206  false
207  ),
208  pointMesh::New(mesh),
210  );
211  pDist.primitiveFieldRef() = pointDistance_;
212 
213  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
214  pDist.write();
215  }
216 
217 
218  //- Direct from cell field and point field.
219  isoSurfPtr_.reset
220  (
221  new isoSurface
222  (
223  mesh,
224  cellDistance,
225  pointDistance_,
226  0,
227  filter_
228  )
229  );
230 
231  if (debug)
232  {
233  print(Pout);
234  Pout<< endl;
235  }
236 }
237 
238 
239 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
240 
242 (
243  const word& name,
244  const polyMesh& mesh,
245  const dictionary& dict
246 )
247 :
248  sampledSurface(name, mesh, dict),
249  plane_(dict),
250  filter_
251  (
252  dict.found("filtering")
253  ? isoSurface::filterTypeNames_.read(dict.lookup("filtering"))
255  ),
256  average_(dict.lookupOrDefault("average", false)),
257  zoneKey_(dict.lookupOrDefault("zone", wordRe::null)),
258  exposedPatchName_(word::null),
259  needsUpdate_(true),
260  subMeshPtr_(nullptr),
261  cellDistancePtr_(nullptr),
262  isoSurfPtr_(nullptr)
263 {
264  if (zoneKey_.size())
265  {
266  dict.lookup("exposedPatchName") >> exposedPatchName_;
267 
268  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
269  {
271  << "Cannot find patch " << exposedPatchName_
272  << " in which to put exposed faces." << endl
273  << "Valid patches are " << mesh.boundaryMesh().names()
274  << exit(FatalError);
275  }
276 
277  if (debug)
278  {
279  Info<< "Restricting to cellZone " << zoneKey_
280  << " with exposed internal faces into patch "
281  << exposedPatchName_ << endl;
282  }
283 
284  if (mesh.cellZones().findIndex(zoneKey_) < 0)
285  {
287  << "cellZone " << zoneKey_
288  << " not found - using entire mesh" << endl;
289  }
290  }
291 }
292 
293 
294 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
295 
297 {}
298 
299 
300 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
301 
303 {
304  return needsUpdate_;
305 }
306 
307 
309 {
310  if (debug)
311  {
312  Pout<< "cuttingPlane::expire :"
313  << " needsUpdate_:" << needsUpdate_ << endl;
314  }
315 
316  // Clear derived data
317  clearGeom();
318 
319  // Already marked as expired
320  if (needsUpdate_)
321  {
322  return false;
323  }
324 
325  needsUpdate_ = true;
326  return true;
327 }
328 
329 
331 {
332  if (debug)
333  {
334  Pout<< "cuttingPlane::update :"
335  << " needsUpdate_:" << needsUpdate_ << endl;
336  }
337 
338  if (!needsUpdate_)
339  {
340  return false;
341  }
342 
343  createGeometry();
344 
345  needsUpdate_ = false;
346  return true;
347 }
348 
349 
352 (
353  const volScalarField& vField
354 ) const
355 {
356  return sampleField(vField);
357 }
358 
359 
362 (
363  const volVectorField& vField
364 ) const
365 {
366  return sampleField(vField);
367 }
368 
369 
372 (
373  const volSphericalTensorField& vField
374 ) const
375 {
376  return sampleField(vField);
377 }
378 
379 
382 (
383  const volSymmTensorField& vField
384 ) const
385 {
386  return sampleField(vField);
387 }
388 
389 
392 (
393  const volTensorField& vField
394 ) const
395 {
396  return sampleField(vField);
397 }
398 
399 
402 (
403  const interpolation<scalar>& interpolator
404 ) const
405 {
406  return interpolateField(interpolator);
407 }
408 
409 
412 (
413  const interpolation<vector>& interpolator
414 ) const
415 {
416  return interpolateField(interpolator);
417 }
418 
421 (
422  const interpolation<sphericalTensor>& interpolator
423 ) const
424 {
425  return interpolateField(interpolator);
426 }
427 
428 
431 (
432  const interpolation<symmTensor>& interpolator
433 ) const
434 {
435  return interpolateField(interpolator);
436 }
437 
438 
441 (
442  const interpolation<tensor>& interpolator
443 ) const
444 {
445  return interpolateField(interpolator);
446 }
447 
448 
450 {
451  os << "cuttingPlane: " << name() << " :"
452  << " plane:" << plane_
453  << " faces:" << faces().size()
454  << " points:" << points().size();
455 }
456 
457 
458 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
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
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 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.
virtual void print(Ostream &) const
Write.
defineTypeNameAndDebug(distanceSurface, 0)
static const NamedEnum< filterType, 4 > filterTypeNames_
Definition: isoSurface.H:92
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
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:251
label findPatchID(const word &patchName) const
Find patch index given a name.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
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.
wordList names() const
Return a list of patch names.
static const word null
An empty word.
Definition: word.H:77
#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:54
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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,.
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const wordRe &) const
Return zone index for the first match, return -1 if not found.
Definition: MeshZones.C:306
messageStream Info
SubField< Type > subField
Declare type of subField.
Definition: Field.H:100
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
A class for managing temporary objects.
Definition: PtrList.H:53
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
virtual bool expire()
Mark the surface as needing an update.