vtkPVFoamPointFields.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-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 InClass
25  vtkPVFoam
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #ifndef vtkPVFoamPointFields_H
30 #define vtkPVFoamPointFields_H
31 
32 // OpenFOAM includes
33 #include "interpolatePointToCell.H"
34 
35 #include "vtkOpenFOAMTupleRemap.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 template<class Type>
40 void Foam::vtkPVFoam::convertPointFields
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::vtkPVFoam::convertPointFields : "
66  << fieldName << endl;
67  }
68 
69  GeometricField<Type, pointPatchField, pointMesh> ptf
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
153  Field<Type> fld
154  (
155  ptf.primitiveField(),
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>
173 void Foam::vtkPVFoam::convertPointFieldBlock
174 (
175  const GeometricField<Type, pointPatchField, pointMesh>& ptf,
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>
202 void Foam::vtkPVFoam::convertPointField
203 (
204  const GeometricField<Type, pointPatchField, pointMesh>& ptf,
205  const GeometricField<Type, fvPatchField, volMesh>& tf,
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 // ************************************************************************* //
label patchId(-1)
#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
uint8_t direction
Definition: direction.H:46
static const char *const typeName
Definition: Field.H:94
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
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))
label nPoints
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Type interpolatePointToCell(const GeometricField< Type, pointPatchField, pointMesh > &ptf, const label celli)
messageStream Info
Interpolates (averages) the vertex values to the cell center.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)