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-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 "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 
49 void Foam::sampledCuttingPlane::createGeometry()
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  {
69  const polyBoundaryMesh& patches = mesh().boundaryMesh();
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.primitiveFieldRef();
129 
130  forAll(cc, i)
131  {
132  // Signed distance
133  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
134  }
135  }
136 
137  volScalarField::Boundary& cellDistanceBf =
138  cellDistance.boundaryFieldRef();
139 
140  // Patch fields
141  {
142  forAll(cellDistanceBf, patchi)
143  {
144  if
145  (
146  isA<emptyFvPatchScalarField>
147  (
148  cellDistanceBf[patchi]
149  )
150  )
151  {
152  cellDistanceBf.set
153  (
154  patchi,
155  new calculatedFvPatchScalarField
156  (
157  fvm.boundary()[patchi],
158  cellDistance
159  )
160  );
161 
162  const polyPatch& pp = fvm.boundary()[patchi].patch();
163  pointField::subField cc = pp.patchSlice(fvm.faceCentres());
164 
165  fvPatchScalarField& fld = cellDistanceBf[patchi];
166  fld.setSize(pp.size());
167  forAll(fld, i)
168  {
169  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
170  }
171  }
172  else
173  {
174  const pointField& cc = fvm.C().boundaryField()[patchi];
175  fvPatchScalarField& fld = cellDistanceBf[patchi];
176 
177  forAll(fld, i)
178  {
179  fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
180  }
181  }
182  }
183  }
184 
185 
186  // On processor patches the mesh.C() will already be the cell centre
187  // on the opposite side so no need to swap cellDistance.
188 
189 
190  // Distance to points
191  pointDistance_.setSize(fvm.nPoints());
192  {
193  const pointField& pts = fvm.points();
194 
195  forAll(pointDistance_, i)
196  {
197  pointDistance_[i] = (pts[i] - plane_.refPoint()) & plane_.normal();
198  }
199  }
200 
201 
202  if (debug)
203  {
204  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
205  cellDistance.write();
206  pointScalarField pDist
207  (
208  IOobject
209  (
210  "pointDistance",
211  fvm.time().timeName(),
212  fvm.time(),
215  false
216  ),
217  pointMesh::New(fvm),
218  dimensionedScalar("zero", dimLength, 0)
219  );
220  pDist.primitiveFieldRef() = pointDistance_;
221 
222  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
223  pDist.write();
224  }
225 
226 
227  //- Direct from cell field and point field.
228  isoSurfPtr_.reset
229  (
230  new isoSurface
231  (
232  cellDistance,
233  pointDistance_,
234  0.0,
235  regularise_,
236  mergeTol_
237  )
238  //new isoSurfaceCell
239  //(
240  // fvm,
241  // cellDistance,
242  // pointDistance_,
243  // 0.0,
244  // regularise_,
245  // mergeTol_
246  //)
247  );
248 
249  if (debug)
250  {
251  print(Pout);
252  Pout<< endl;
253  }
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
258 
260 (
261  const word& name,
262  const polyMesh& mesh,
263  const dictionary& dict
264 )
265 :
266  sampledSurface(name, mesh, dict),
267  plane_(dict),
268  mergeTol_(dict.lookupOrDefault("mergeTol", 1e-6)),
269  regularise_(dict.lookupOrDefault("regularise", true)),
270  average_(dict.lookupOrDefault("average", false)),
271  zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
272  exposedPatchName_(word::null),
273  needsUpdate_(true),
274  subMeshPtr_(NULL),
275  cellDistancePtr_(NULL),
276  isoSurfPtr_(NULL),
277  facesPtr_(NULL)
278 {
279  if (zoneID_.index() != -1)
280  {
281  dict.lookup("exposedPatchName") >> exposedPatchName_;
282 
283  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
284  {
286  << "Cannot find patch " << exposedPatchName_
287  << " in which to put exposed faces." << endl
288  << "Valid patches are " << mesh.boundaryMesh().names()
289  << exit(FatalError);
290  }
291 
292  if (debug && zoneID_.index() != -1)
293  {
294  Info<< "Restricting to cellZone " << zoneID_.name()
295  << " with exposed internal faces into patch "
296  << exposedPatchName_ << endl;
297  }
298  }
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
303 
305 {}
306 
307 
308 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
309 
311 {
312  return needsUpdate_;
313 }
314 
315 
317 {
318  if (debug)
319  {
320  Pout<< "sampledCuttingPlane::expire :"
321  << " have-facesPtr_:" << facesPtr_.valid()
322  << " needsUpdate_:" << needsUpdate_ << endl;
323  }
324 
325  // Clear any stored topologies
326  facesPtr_.clear();
327 
328  // Clear derived data
329  clearGeom();
330 
331  // already marked as expired
332  if (needsUpdate_)
333  {
334  return false;
335  }
336 
337  needsUpdate_ = true;
338  return true;
339 }
340 
341 
343 {
344  if (debug)
345  {
346  Pout<< "sampledCuttingPlane::update :"
347  << " have-facesPtr_:" << facesPtr_.valid()
348  << " needsUpdate_:" << needsUpdate_ << endl;
349  }
350 
351  if (!needsUpdate_)
352  {
353  return false;
354  }
355 
356  createGeometry();
357 
358  needsUpdate_ = false;
359  return true;
360 }
361 
362 
365 (
366  const volScalarField& vField
367 ) const
368 {
369  return sampleField(vField);
370 }
371 
372 
375 (
376  const volVectorField& vField
377 ) const
378 {
379  return sampleField(vField);
380 }
381 
382 
385 (
386  const volSphericalTensorField& vField
387 ) const
388 {
389  return sampleField(vField);
390 }
391 
392 
395 (
396  const volSymmTensorField& vField
397 ) const
398 {
399  return sampleField(vField);
400 }
401 
402 
405 (
406  const volTensorField& vField
407 ) const
408 {
409  return sampleField(vField);
410 }
411 
412 
415 (
416  const interpolation<scalar>& interpolator
417 ) const
418 {
419  return interpolateField(interpolator);
420 }
421 
422 
425 (
426  const interpolation<vector>& interpolator
427 ) const
428 {
429  return interpolateField(interpolator);
430 }
431 
434 (
435  const interpolation<sphericalTensor>& interpolator
436 ) const
437 {
438  return interpolateField(interpolator);
439 }
440 
441 
444 (
445  const interpolation<symmTensor>& interpolator
446 ) const
447 {
448  return interpolateField(interpolator);
449 }
450 
451 
454 (
455  const interpolation<tensor>& interpolator
456 ) const
457 {
458  return interpolateField(interpolator);
459 }
460 
461 
463 {
464  os << "sampledCuttingPlane: " << name() << " :"
465  << " plane:" << plane_
466  << " faces:" << faces().size()
467  << " points:" << points().size();
468 }
469 
470 
471 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
virtual ~sampledCuttingPlane()
Destructor.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
wordList names() const
Return a list of patch names.
virtual bool expire()
Mark the surface as needing an update.
const double e
Elementary charge.
Definition: doubleFloat.H:78
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:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Macros for easy insertion into run-time selection tables.
static const pointMesh & New(const polyMesh &mesh)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
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.
virtual bool needsUpdate() const
Does the surface need an update?
static const word null
An empty word.
Definition: word.H:77
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addNamedToRunTimeSelectionTable(GAMGProcAgglomeration, noneGAMGProcAgglomeration, GAMGAgglomeration, none)
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
bool interpolate() const
Interpolation requested for surface.
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
label patchi
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
sampledCuttingPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
messageStream Info
SubField< vector > subField
Declare type of subField.
Definition: Field.H:89
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:54
label findPatchID(const word &patchName) const
Find patch index given a name.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual bool update()
Update the surface as required.
virtual void print(Ostream &) const
Write.
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451