1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
27 \*---------------------------------------------------------------------------*/
29 #ifndef vtkPV3FoamPointFields_H
30 #define vtkPV3FoamPointFields_H
33 #include "interpolatePointToCell.H"
35 #include "vtkFOAMTupleRemap.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 void Foam::vtkPV3Foam::convertPointFields
43 const pointMesh& pMesh,
44 const IOobjectList& objects,
45 vtkMultiBlockDataSet* output
48 const polyBoundaryMesh& patches = mesh.boundaryMesh();
50 forAllConstIter(IOobjectList, objects, iter)
52 const word& fieldName = iter()->name();
53 // restrict to this GeometricField<Type, ...>
56 iter()->headerClassName()
57 != GeometricField<Type, pointPatchField, pointMesh>::typeName
65 Info<< "Foam::vtkPV3Foam::convertPointFields : "
69 GeometricField<Type, pointPatchField, pointMesh> ptf
76 // Convert activated internalMesh regions
77 convertPointFieldBlock
85 // Convert activated cellZones
86 convertPointFieldBlock
94 // Convert activated cellSets
95 convertPointFieldBlock
105 // Convert patches - if activated
109 int partId = partInfoPatches_.start();
110 partId < partInfoPatches_.end();
114 const word patchName = getPartName(partId);
115 const label datasetNo = partDataset_[partId];
116 const label patchId = patches.findPatchID(patchName);
118 if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
123 convertPatchPointField
126 ptf.boundaryField()[patchId].patchInternalField()(),
137 void Foam::vtkPV3Foam::convertPointFieldBlock
139 const GeometricField<Type, pointPatchField, pointMesh>& ptf,
140 vtkMultiBlockDataSet* output,
141 const partInfo& selector,
142 const List<polyDecomp>& decompLst
145 for (int partId = selector.start(); partId < selector.end(); ++partId)
147 const label datasetNo = partDataset_[partId];
149 if (datasetNo >= 0 && partStatus_[partId])
154 GeometricField<Type, fvPatchField, volMesh>::null(),
166 void Foam::vtkPV3Foam::convertPointField
168 const GeometricField<Type, pointPatchField, pointMesh>& ptf,
169 const GeometricField<Type, fvPatchField, volMesh>& tf,
170 vtkMultiBlockDataSet* output,
171 const partInfo& selector,
172 const label datasetNo,
173 const polyDecomp& decomp
176 const label nComp = pTraits<Type>::nComponents;
177 const labelList& addPointCellLabels = decomp.addPointCellLabels();
178 const labelList& pointMap = decomp.pointMap();
180 // use a pointMap or address directly into mesh
184 nPoints = pointMap.size();
188 nPoints = ptf.size();
191 vtkFloatArray *pointData = vtkFloatArray::New();
192 pointData->SetNumberOfTuples(nPoints + addPointCellLabels.size());
193 pointData->SetNumberOfComponents(nComp);
194 pointData->Allocate(nComp*(nPoints + addPointCellLabels.size()));
195 pointData->SetName(ptf.name().c_str());
199 Info<< "convert convertPointField: "
201 << " size = " << nPoints
202 << " nComp=" << nComp
203 << " nTuples = " << (nPoints + addPointCellLabels.size())
213 const Type& t = ptf[pointMap[i]];
214 for (direction d=0; d<nComp; d++)
216 vec[d] = component(t, d);
218 vtkFOAMTupleRemap<Type>(vec);
220 pointData->InsertTuple(i, vec);
227 const Type& t = ptf[i];
228 for (direction d=0; d<nComp; d++)
230 vec[d] = component(t, d);
232 vtkFOAMTupleRemap<Type>(vec);
234 pointData->InsertTuple(i, vec);
238 // continue insertion from here
241 if (&tf != &GeometricField<Type, fvPatchField, volMesh>::null())
243 forAll(addPointCellLabels, apI)
245 const Type& t = tf[addPointCellLabels[apI]];
246 for (direction d=0; d<nComp; d++)
248 vec[d] = component(t, d);
250 vtkFOAMTupleRemap<Type>(vec);
252 pointData->InsertTuple(i++, vec);
257 forAll(addPointCellLabels, apI)
259 Type t = interpolatePointToCell(ptf, addPointCellLabels[apI]);
260 for (direction d=0; d<nComp; d++)
262 vec[d] = component(t, d);
264 vtkFOAMTupleRemap<Type>(vec);
266 pointData->InsertTuple(i++, vec);
270 vtkUnstructuredGrid::SafeDownCast
272 GetDataSetFromBlock(output, selector, datasetNo)
274 ->AddArray(pointData);
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 // ************************************************************************* //