vtkPVFoamPointFields.H
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-2025 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 #include "vtkPVFoam.H"
33 #include "vtkPVFoamReader.h"
34 #include "vtkOpenFOAMTupleRemap.H"
35 
36 // OpenFOAM includes
37 #include "domainDecomposition.H"
38 #include "pointFields.H"
39 #include "interpolatePointToCell.H"
41 
42 // VTK includes
43 #include "vtkFloatArray.h"
44 #include "vtkUnstructuredGrid.h"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 template<class Type>
49 void Foam::vtkPVFoam::convertPointFields
50 (
51  const IOobjectList& objects,
52  vtkMultiBlockDataSet* output
53 )
54 {
55  const polyBoundaryMesh& patches =
56  procMeshesPtr_->completeMesh().boundaryMesh();
57 
58  forAllConstIter(IOobjectList, objects, iter)
59  {
60  // restrict to this GeometricField<Type, ...>
61  if (iter()->headerClassName() != PointField<Type>::typeName)
62  {
63  continue;
64  }
65 
66  // Load the field
67  tmp<PointField<Type>> tptf;
68  if (reader_->GetDecomposedCase())
69  {
70  if (!pointReconstructorPtr_.valid())
71  {
72  pointReconstructorPtr_.set
73  (
74  new pointFieldReconstructor
75  (
76  pointMesh::New(procMeshesPtr_->completeMesh()),
77  procMeshesPtr_->procMeshes(),
78  procMeshesPtr_->procPointAddressing()
79  )
80  );
81  }
82 
83  tptf =
84  pointReconstructorPtr_
85  ->reconstructField<Type>(*iter());
86  }
87  else
88  {
89  tptf =
90  new PointField<Type>
91  (
92  *iter(),
93  pointMesh::New(procMeshesPtr_->completeMesh())
94  );
95  }
96  const PointField<Type>& ptf = tptf();
97 
98  // Convert activated internalMesh regions
99  convertPointFieldBlock
100  (
101  ptf,
102  output,
103  arrayRangeVolume_,
104  regionPolyDecomp_
105  );
106 
107  // Convert activated cellZones
108  convertPointFieldBlock
109  (
110  ptf,
111  output,
112  arrayRangeCellZones_,
113  zonePolyDecomp_
114  );
115 
116  // Convert activated cellSets
117  convertPointFieldBlock
118  (
119  ptf,
120  output,
121  arrayRangeCellSets_,
122  setPolyDecomp_
123  );
124 
125  // Convert patches - if activated
126  for
127  (
128  int partId = arrayRangePatches_.start();
129  partId < arrayRangePatches_.end();
130  ++partId
131  )
132  {
133  const word patchName = getPartName(partId);
134  const label datasetNo = partDataset_[partId];
135  const label patchId = patches.findIndex(patchName);
136 
137  if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
138  {
139  continue;
140  }
141 
142  convertPatchPointField
143  (
144  ptf.name(),
145  ptf.boundaryField()[patchId].patchInternalField()(),
146  output,
147  arrayRangePatches_,
148  datasetNo
149  );
150  }
151 
152  // Convert faceZones - if activated
153  for
154  (
155  int partId = arrayRangeFaceZones_.start();
156  partId < arrayRangeFaceZones_.end();
157  ++partId
158  )
159  {
160  const word zoneName = getPartName(partId);
161  const label datasetNo = partDataset_[partId];
162  const label zoneId =
163  ptf.mesh().mesh().faceZones().findIndex(zoneName);
164 
165  if (!partStatus_[partId] || datasetNo < 0 || zoneId < 0)
166  {
167  continue;
168  }
169 
170  // Extract the field on the zone
171  Field<Type> fld
172  (
173  ptf.primitiveField(),
174  ptf.mesh().mesh().faceZones()[zoneId].patch().meshPoints()
175  );
176 
177  convertPatchPointField
178  (
179  ptf.name(),
180  fld,
181  output,
182  arrayRangeFaceZones_,
183  datasetNo
184  );
185  }
186  }
187 }
188 
189 
190 template<class Type>
191 void Foam::vtkPVFoam::convertPointFieldBlock
192 (
193  const PointField<Type>& ptf,
194  vtkMultiBlockDataSet* output,
195  const arrayRange& range,
196  const List<polyDecomp>& decompLst
197 )
198 {
199  for (int partId = range.start(); partId < range.end(); ++partId)
200  {
201  const label datasetNo = partDataset_[partId];
202 
203  if (datasetNo >= 0 && partStatus_[partId])
204  {
205  convertPointField
206  (
207  ptf,
209  output,
210  range,
211  datasetNo,
212  decompLst[datasetNo]
213  );
214  }
215  }
216 }
217 
218 
219 template<class Type>
220 void Foam::vtkPVFoam::convertPointField
221 (
222  const PointField<Type>& ptf,
223  const VolField<Type>& tf,
224  vtkMultiBlockDataSet* output,
225  const arrayRange& range,
226  const label datasetNo,
227  const polyDecomp& decomp
228 )
229 {
230  const label nComp = pTraits<Type>::nComponents;
231  const labelList& addPointCellLabels = decomp.addPointCellLabels();
232  const labelList& pointMap = decomp.pointMap();
233 
234  // use a pointMap or address directly into mesh
235  label nPoints;
236  if (pointMap.size())
237  {
238  nPoints = pointMap.size();
239  }
240  else
241  {
242  nPoints = ptf.size();
243  }
244 
245  vtkFloatArray* pointData = vtkFloatArray::New();
246  pointData->SetNumberOfTuples(nPoints + addPointCellLabels.size());
247  pointData->SetNumberOfComponents(nComp);
248  pointData->Allocate(nComp*(nPoints + addPointCellLabels.size()));
249 
250  // Note: using the name of the original volField
251  // not the name generated by the interpolation "volPointInterpolate(<name>)"
252  if (&tf != &VolField<Type>::null())
253  {
254  pointData->SetName(tf.name().c_str());
255  }
256  else
257  {
258  pointData->SetName(ptf.name().c_str());
259  }
260 
262  << "Converting Point field: " << tf.name()
263  << " size=" << nPoints << " (" << nPoints + addPointCellLabels.size()
264  << "), nComp=" << nComp << endl;
265 
266  float vec[nComp];
267 
268  if (pointMap.size())
269  {
270  forAll(pointMap, i)
271  {
272  const Type& t = ptf[pointMap[i]];
273  for (direction d=0; d<nComp; ++d)
274  {
275  vec[d] = component(t, d);
276  }
277  vtkOpenFOAMTupleRemap<Type>(vec);
278 
279  pointData->InsertTuple(i, vec);
280  }
281  }
282  else
283  {
284  forAll(ptf, i)
285  {
286  const Type& t = ptf[i];
287  for (direction d=0; d<nComp; ++d)
288  {
289  vec[d] = component(t, d);
290  }
291  vtkOpenFOAMTupleRemap<Type>(vec);
292 
293  pointData->InsertTuple(i, vec);
294  }
295  }
296 
297  // continue insertion from here
298  label i = nPoints;
299 
300  if (&tf != &VolField<Type>::null())
301  {
302  forAll(addPointCellLabels, apI)
303  {
304  const Type& t = tf[addPointCellLabels[apI]];
305  for (direction d=0; d<nComp; ++d)
306  {
307  vec[d] = component(t, d);
308  }
309  vtkOpenFOAMTupleRemap<Type>(vec);
310 
311  pointData->InsertTuple(i++, vec);
312  }
313  }
314  else
315  {
316  forAll(addPointCellLabels, apI)
317  {
318  Type t = interpolatePointToCell(ptf, addPointCellLabels[apI]);
319  for (direction d=0; d<nComp; ++d)
320  {
321  vec[d] = component(t, d);
322  }
323  vtkOpenFOAMTupleRemap<Type>(vec);
324 
325  pointData->InsertTuple(i++, vec);
326  }
327  }
328 
329  vtkUnstructuredGrid::SafeDownCast
330  (
331  GetDataSetFromBlock(output, range, datasetNo)
332  ) ->GetPointData()
333  ->AddArray(pointData);
334 
335  pointData->Delete();
336 }
337 
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 #endif
342 
343 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static const GeometricField< Type, GeoMesh, PrimitiveField > & null()
Return a null geometric field.
const tensorField & tf
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Interpolates (averages) the vertex values to the cell center.
label patchId(-1)
const fvPatchList & patches
#define DebugInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
Type interpolatePointToCell(const PointField< Type > &ptf, const label celli)
uint8_t direction
Definition: direction.H:45
objects