cuttingPlaneFilter.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2019 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // OpenFOAM includes
29 #include "cuttingPlaneFilter.H"
30 #include "runTimePostProcessing.H"
32 
33 // VTK includes
34 #include "vtkActor.h"
35 #include "vtkCellDataToPointData.h"
36 #include "vtkCompositeDataGeometryFilter.h"
37 #include "vtkCompositeDataSet.h"
38 #include "vtkCompositePolyDataMapper.h"
39 #include "vtkCutter.h"
40 #include "vtkMultiPieceDataSet.h"
41 #include "vtkPlane.h"
42 #include "vtkPolyData.h"
43 #include "vtkPolyDataMapper.h"
44 #include "vtkRenderer.h"
45 #include "vtkSmartPointer.h"
46 
47 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 namespace functionObjects
52 {
53 namespace runTimePostPro
54 {
55  defineTypeName(cuttingPlaneFilter);
56  addToRunTimeSelectionTable(surface, cuttingPlaneFilter, dictionary);
57 }
58 }
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 (
66  const runTimePostProcessing& parent,
67  const dictionary& dict,
68  const HashPtrTable<Function1<vector>>& colours
69 )
70 :
71  volumeFilter(parent, dict, colours),
72  fieldVisualisationBase(dict, colours),
73  plane_(dict),
74  values_()
75 {
76  dict.readIfPresent("offsets", values_);
77 
78  if (values_.empty())
79  {
80  values_.resize(1);
81  values_.first() = Zero;
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
87 
90 (
91  const scalar position,
92  vtkRenderer* renderer
93 )
94 {
95  if (!visible_)
96  {
97  return false;
98  }
99 
100  if (needsCollective())
101  {
102  Info<< type() << " : Not available for collective operation" << endl;
103  return false;
104  }
105 
106  DebugInfo << " Adding cutting plane" << endl;
107 
108 
109  // Bookkeeping for vtkUnstructuredGrid
110  vtk::vtuAdaptor adaptor;
111  vtkSmartPointer<vtkMultiPieceDataSet> multiPiece = mesh(adaptor);
112 
113 
114  // Add (scalar/vector) field.
115  // - Need field(s) for glyphs or colourByField:
116 
117  int nCmpt = 0;
118  if (representation_ == rtGlyph || colourBy_ == cbField)
119  {
120  const auto* ioptr =
121  parent().mesh().cfindObject<regIOobject>(fieldName_);
122 
123  if (!nCmpt)
124  {
125  nCmpt = addDimField<scalar>
126  (
127  multiPiece, adaptor, ioptr, fieldName_
128  );
129  }
130  if (!nCmpt)
131  {
132  nCmpt = addDimField<vector>
133  (
134  multiPiece, adaptor, ioptr, fieldName_
135  );
136  }
137  }
138 
139 
140  // Now have a multi-piece dataset that is one of the following:
141  //
142  // - one-piece per processor (OpenFOAM = parallel, VTK=parallel)
143 
144 
145  // Re-query field information - we may have stored it differently
146  // than the original source.
147 
148  fieldSummary fieldInfo = queryFieldSummary(fieldName_, multiPiece);
149  fieldInfo.reduce();
150 
151 
152  // Not rendered on this processor?
153  // This is where we stop, but could also have an MPI barrier
154  if (!renderer)
155  {
156  return true;
157  }
158 
159 
160  // Rendering
161  {
162  // OpenFOAM plane -> vtkPlane definition
163 
164  auto pln = vtkSmartPointer<vtkPlane>::New();
165 
166  pln->SetNormal
167  (
168  plane_.normal().x(),
169  plane_.normal().y(),
170  plane_.normal().z()
171  );
172  pln->SetOrigin
173  (
174  plane_.origin().x(),
175  plane_.origin().y(),
176  plane_.origin().z()
177  );
178 
179 
180  // Plane cutting algorithm
181 
182  auto cutter = vtkSmartPointer<vtkCutter>::New();
183 
184  cutter->SetInputData(multiPiece);
185  cutter->SetCutFunction(pln);
186 
187  cutter->SetNumberOfContours(values_.size());
188 
189  forAll(values_, pointi)
190  {
191  cutter->SetValue(pointi, values_[pointi]);
192  }
193 
194  cutter->SetInputArrayToProcess
195  (
196  (nCmpt == 3 ? 1 : 0), // index: scalars(0), vectors(1)
197  0, // port
198  0, // connection
199  vtkDataObject::FIELD_ASSOCIATION_CELLS,
200  fieldName_.c_str()
201  );
202 
203  cutter->Modified();
204  cutter->Update();
205 
207 
208  polyData->SetInputConnection(cutter->GetOutputPort());
209  polyData->Update();
210 
212  mapper->SetInputConnection(polyData->GetOutputPort());
213 
214  if (representation_ == rtGlyph)
215  {
216  addGlyphs
217  (
218  position,
219  fieldName_, fieldInfo, // scaling
220  fieldName_, fieldInfo, // colouring
221  maxGlyphLength_,
222  polyData->GetOutput(),
223  surfaceActor_,
224  renderer
225  );
226  }
227  else
228  {
230 
231  // CellData - Need a cell->point filter
232  if (smooth_ && !fieldInfo.hasPointData())
233  {
235  cellToPoint->SetInputConnection(cutter->GetOutputPort());
236 
237  polyData->SetInputConnection(cellToPoint->GetOutputPort());
238  polyData->Update();
239  }
240 
241  setField
242  (
243  position,
244  fieldName_,
245  (
246  smooth_
248  : FieldAssociation(fieldInfo.association_)
249  ),
250  mapper,
251  renderer
252  );
253 
254  surfaceActor_->SetMapper(mapper);
255 
256  setRepresentation(surfaceActor_);
257 
258  renderer->AddActor(surfaceActor_);
259  }
260  }
261 
262  return true;
263 }
264 
265 
268 (
269  const scalar position,
270  vtkRenderer* renderer
271 )
272 {
273  if (visible_)
274  {
275  // Live source
276  if (addGeometry(position, renderer))
277  {
278  return;
279  }
280 
282  << "Unsupported for OpenFOAM parallel and VTK serial"
283  << endl;
284  }
285 }
286 
287 
289 {
290  return true;
291 }
292 
293 
294 // ************************************************************************* //
Foam::functionObjects::runTimePostPro::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surface, contourFilter, dictionary)
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::fieldSummary
General field characteristics.
Definition: fieldVisualisationBase.H:172
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::expressions::patchExpr::POINT_DATA
Point data.
Definition: patchExprFwd.H:60
Foam::cellToPoint
A topoSetPointSource to select points based on usage in cells.
Definition: cellToPoint.H:83
Foam::functionObjects::fieldInfo
Helper class to store a wordRe and label used by Foam::functionObjects::fieldSelection.
Definition: fieldInfo.H:56
cuttingPlaneFilter.H
Foam::functionObjects::runTimePostPro::cuttingPlaneFilter::clear
virtual bool clear()
Cleanup files etc.
Definition: cuttingPlaneFilter.C:288
vtkSmartPointer
Definition: runTimePostProcessing.H:148
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:56
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::expressions::patchExpr::FieldAssociation
FieldAssociation
The field association for patch expressions (mutually exclusive)
Definition: patchExprFwd.H:57
Foam::functionObjects::runTimePostPro::cuttingPlaneFilter::cuttingPlaneFilter
cuttingPlaneFilter(const cuttingPlaneFilter &)=delete
No copy construct.
Foam::functionObjects::runTimePostPro::cuttingPlaneFilter::addGeometry
bool addGeometry(const scalar position, vtkRenderer *renderer)
Add cutting planes to scene (using simulation source)
Definition: cuttingPlaneFilter.C:90
Foam::functionObjects::runTimePostProcessing
Generate images during run-time.
Definition: runTimePostProcessing.H:170
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
setField
surfacesMesh setField(triSurfaceToAgglom)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:67
Foam::functionObjects::runTimePostPro::fieldVisualisationBase
Base class for scene objects.
Definition: fieldVisualisationBase.H:125
runTimePostProcessing.H
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::objectRegistry::cfindObject
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:390
Foam::HashPtrTable
A HashTable of pointers to objects of type <T>.
Definition: HashPtrTable.H:54
Foam::functionObjects::runTimePostPro::defineTypeName
defineTypeName(contourFilter)
Foam::vtk::vtuAdaptor
Bookkeeping for vtkUnstructuredGrid.
Definition: foamVtkVtuAdaptor.H:78
Foam::functionObjects::runTimePostPro::cuttingPlaneFilter::addGeometryToScene
virtual void addGeometryToScene(const scalar position, vtkRenderer *renderer)
Add cutting planes to scene (using simulation source)
Definition: cuttingPlaneFilter.C:268
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::functionObjects::runTimePostProcessing::mesh
const fvMesh & mesh() const
Reference to the underlying OpenFOAM mesh.
Definition: runTimePostProcessing.H:259
Foam::functionObjects::runTimePostPro::volumeFilter
Visualisation of OpenFOAM volume fields as surface data using a VTK filter cascade.
Definition: volumeFilter.H:74