vtkPV4FoamPointFields.H
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-2013 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 InClass
25  vtkPV4Foam
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #ifndef vtkPV4FoamPointFields_H
30 #define vtkPV4FoamPointFields_H
31 
32 // OpenFOAM includes
33 #include "interpolatePointToCell.H"
34 
35 #include "vtkOpenFOAMTupleRemap.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 template<class Type>
41 (
42  const fvMesh& mesh,
43  const pointMesh& pMesh,
44  const IOobjectList& objects,
45  vtkMultiBlockDataSet* output
46 )
47 {
48  const polyBoundaryMesh& patches = mesh.boundaryMesh();
49 
50  forAllConstIter(IOobjectList, objects, iter)
51  {
52  const word& fieldName = iter()->name();
53  // restrict to this GeometricField<Type, ...>
54  if
55  (
56  iter()->headerClassName()
58  )
59  {
60  continue;
61  }
62 
63  if (debug)
64  {
65  Info<< "Foam::vtkPV4Foam::convertPointFields : "
66  << fieldName << endl;
67  }
68 
70  (
71  *iter(),
72  pMesh
73  );
74 
75 
76  // Convert activated internalMesh regions
77  convertPointFieldBlock
78  (
79  ptf,
80  output,
81  arrayRangeVolume_,
82  regionPolyDecomp_
83  );
84 
85  // Convert activated cellZones
86  convertPointFieldBlock
87  (
88  ptf,
89  output,
90  arrayRangeCellZones_,
91  zonePolyDecomp_
92  );
93 
94  // Convert activated cellSets
95  convertPointFieldBlock
96  (
97  ptf,
98  output,
99  arrayRangeCellSets_,
100  csetPolyDecomp_
101  );
102 
103 
104  //
105  // Convert patches - if activated
106  //
107  for
108  (
109  int partId = arrayRangePatches_.start();
110  partId < arrayRangePatches_.end();
111  ++partId
112  )
113  {
114  const word patchName = getPartName(partId);
115  const label datasetNo = partDataset_[partId];
116  const label patchId = patches.findPatchID(patchName);
117 
118  if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
119  {
120  continue;
121  }
122 
123  convertPatchPointField
124  (
125  fieldName,
126  ptf.boundaryField()[patchId].patchInternalField()(),
127  output,
128  arrayRangePatches_,
129  datasetNo
130  );
131  }
132 
133  //
134  // Convert faceZones - if activated
135  //
136  for
137  (
138  int partId = arrayRangeFaceZones_.start();
139  partId < arrayRangeFaceZones_.end();
140  ++partId
141  )
142  {
143  const word zoneName = getPartName(partId);
144  const label datasetNo = partDataset_[partId];
145  const label zoneId = mesh.faceZones().findZoneID(zoneName);
146 
147  if (!partStatus_[partId] || datasetNo < 0 || zoneId < 0)
148  {
149  continue;
150  }
151 
152  // Extract the field on the zone
154  (
155  ptf.internalField(),
156  mesh.faceZones()[zoneId]().meshPoints()
157  );
158 
159  convertPatchPointField
160  (
161  fieldName,
162  fld,
163  output,
164  arrayRangeFaceZones_,
165  datasetNo
166  );
167  }
168  }
169 }
170 
171 
172 template<class Type>
174 (
176  vtkMultiBlockDataSet* output,
177  const arrayRange& range,
178  const List<polyDecomp>& decompLst
179 )
180 {
181  for (int partId = range.start(); partId < range.end(); ++partId)
182  {
183  const label datasetNo = partDataset_[partId];
184 
185  if (datasetNo >= 0 && partStatus_[partId])
186  {
187  convertPointField
188  (
189  ptf,
191  output,
192  range,
193  datasetNo,
194  decompLst[datasetNo]
195  );
196  }
197  }
198 }
199 
200 
201 template<class Type>
203 (
206  vtkMultiBlockDataSet* output,
207  const arrayRange& range,
208  const label datasetNo,
209  const polyDecomp& decomp
210 )
211 {
212  const label nComp = pTraits<Type>::nComponents;
213  const labelList& addPointCellLabels = decomp.addPointCellLabels();
214  const labelList& pointMap = decomp.pointMap();
215 
216  // use a pointMap or address directly into mesh
217  label nPoints;
218  if (pointMap.size())
219  {
220  nPoints = pointMap.size();
221  }
222  else
223  {
224  nPoints = ptf.size();
225  }
226 
227  vtkFloatArray* pointData = vtkFloatArray::New();
228  pointData->SetNumberOfTuples(nPoints + addPointCellLabels.size());
229  pointData->SetNumberOfComponents(nComp);
230  pointData->Allocate(nComp*(nPoints + addPointCellLabels.size()));
231 
232  // Note: using the name of the original volField
233  // not the name generated by the interpolation "volPointInterpolate(<name>)"
234 
236  {
237  pointData->SetName(tf.name().c_str());
238  }
239  else
240  {
241  pointData->SetName(ptf.name().c_str());
242  }
243 
244  if (debug)
245  {
246  Info<< "convert convertPointField: "
247  << ptf.name()
248  << " size = " << nPoints
249  << " nComp=" << nComp
250  << " nTuples = " << (nPoints + addPointCellLabels.size())
251  << endl;
252  }
253 
254  float vec[nComp];
255 
256  if (pointMap.size())
257  {
258  forAll(pointMap, i)
259  {
260  const Type& t = ptf[pointMap[i]];
261  for (direction d=0; d<nComp; ++d)
262  {
263  vec[d] = component(t, d);
264  }
265  vtkOpenFOAMTupleRemap<Type>(vec);
266 
267  pointData->InsertTuple(i, vec);
268  }
269  }
270  else
271  {
272  forAll(ptf, i)
273  {
274  const Type& t = ptf[i];
275  for (direction d=0; d<nComp; ++d)
276  {
277  vec[d] = component(t, d);
278  }
279  vtkOpenFOAMTupleRemap<Type>(vec);
280 
281  pointData->InsertTuple(i, vec);
282  }
283  }
284 
285  // continue insertion from here
286  label i = nPoints;
287 
289  {
290  forAll(addPointCellLabels, apI)
291  {
292  const Type& t = tf[addPointCellLabels[apI]];
293  for (direction d=0; d<nComp; ++d)
294  {
295  vec[d] = component(t, d);
296  }
297  vtkOpenFOAMTupleRemap<Type>(vec);
298 
299  pointData->InsertTuple(i++, vec);
300  }
301  }
302  else
303  {
304  forAll(addPointCellLabels, apI)
305  {
306  Type t = interpolatePointToCell(ptf, addPointCellLabels[apI]);
307  for (direction d=0; d<nComp; ++d)
308  {
309  vec[d] = component(t, d);
310  }
311  vtkOpenFOAMTupleRemap<Type>(vec);
312 
313  pointData->InsertTuple(i++, vec);
314  }
315  }
316 
317  vtkUnstructuredGrid::SafeDownCast
318  (
319  GetDataSetFromBlock(output, range, datasetNo)
320  ) ->GetPointData()
321  ->AddArray(pointData);
322 
323  pointData->Delete();
324 }
325 
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 #endif
330 
331 // ************************************************************************* //
Foam::vtkPV4Foam::arrayRange
Bookkeeping for GUI checklists and the multi-block organization.
Definition: vtkPV4Foam.H:112
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:41
Foam::vtkPV4Foam::polyDecomp
Bookkeeping for polyhedral cell decomposition.
Definition: vtkPV4Foam.H:188
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::interpolatePointToCell
Type interpolatePointToCell(const GeometricField< Type, pointPatchField, pointMesh > &ptf, const label celli)
Definition: interpolatePointToCell.C:32
interpolatePointToCell.H
Interpolates (averages) the vertex values to the cell center.
Foam::vtkPV4Foam::convertPointField
void convertPointField(const GeometricField< Type, pointPatchField, pointMesh > &, const GeometricField< Type, fvPatchField, volMesh > &, vtkMultiBlockDataSet *output, const arrayRange &, const label datasetNo, const polyDecomp &)
Point fields.
Definition: vtkPV4FoamPointFields.H:203
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::pointData
Variant of pointEdgePoint with some transported additional data. WIP - should be templated on data li...
Definition: pointData.H:53
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
tf
const tensorField & tf
Definition: getPatchFieldTensor.H:36
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< Type >
Foam::Info
messageStream Info
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::vtkPV4Foam::polyDecomp::pointMap
labelList & pointMap()
Point labels for subsetted meshes.
Definition: vtkPV4Foam.H:224
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
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::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::vtkPV4Foam::polyDecomp::addPointCellLabels
labelList & addPointCellLabels()
Cell-centre labels for additional points of decomposed cells.
Definition: vtkPV4Foam.H:212
Foam::vtkPV4Foam::convertPointFieldBlock
void convertPointFieldBlock(const GeometricField< Type, pointPatchField, pointMesh > &, vtkMultiBlockDataSet *output, const arrayRange &, const List< polyDecomp > &)
Point field - all selected parts.
Definition: vtkPV4FoamPointFields.H:174
range
scalar range
Definition: LISASMDCalcMethod1.H:12
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::vtkPV3Readers::GetDataSetFromBlock
vtkDataSet * GetDataSetFromBlock(vtkMultiBlockDataSet *output, const partInfo &selector, const label datasetNo)
Convenience method use to convert the readers from VTK 5.
Definition: vtkPV3Readers.C:140
patches
patches[0]
Definition: createSingleCellMesh.H:36
patchId
label patchId(-1)
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::vtkPV4Foam::convertPointFields
void convertPointFields(vtkMultiBlockDataSet *)
Convert point fields.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
vtkOpenFOAMTupleRemap.H