mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: paraview reader module with internal caching of the vtk geometries
- The reader module allows two levels of caching. The OpenFOAM fvMesh can be cached in memory, for faster loading of fields. Additionally, the translated VTK geometries are held in a local cache. The cached VTK geometries should incur no additional overhead since they use the VTK reference counting for their storage management.
This commit is contained in:
@ -40,6 +40,7 @@ License
|
||||
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
|
||||
|
||||
// file-scope
|
||||
// Widget properties
|
||||
static QWidget* setWidgetProperties
|
||||
(
|
||||
QWidget* widget,
|
||||
@ -67,6 +68,7 @@ static QWidget* setWidgetProperties
|
||||
|
||||
|
||||
// file-scope
|
||||
// Button properties
|
||||
static QAbstractButton* setButtonProperties
|
||||
(
|
||||
QAbstractButton* b,
|
||||
|
||||
@ -60,8 +60,7 @@ namespace Foam
|
||||
const Foam::point& pt
|
||||
)
|
||||
{
|
||||
vtkSmartPointer<vtkTextActor> txt =
|
||||
vtkSmartPointer<vtkTextActor>::New();
|
||||
auto txt = vtkSmartPointer<vtkTextActor>::New();
|
||||
|
||||
txt->SetInput(s.c_str());
|
||||
|
||||
@ -310,8 +309,8 @@ void Foam::vtkPVblockMesh::updateInfo()
|
||||
HashSet<string> enabledEdges;
|
||||
if (!firstTime)
|
||||
{
|
||||
enabledParts = getSelectedArrayEntries(blockSelection);
|
||||
enabledEdges = getSelectedArrayEntries(edgeSelection);
|
||||
enabledParts = getSelectedArraySet(blockSelection);
|
||||
enabledEdges = getSelectedArraySet(edgeSelection);
|
||||
}
|
||||
|
||||
// Clear current mesh parts list
|
||||
@ -407,14 +406,6 @@ void Foam::vtkPVblockMesh::Update
|
||||
{
|
||||
reader_->UpdateProgress(0.1);
|
||||
|
||||
// Set up mesh parts selection(s)
|
||||
getSelected(blockStatus_, reader_->GetBlockSelection());
|
||||
|
||||
// Set up curved edges selection(s)
|
||||
getSelected(edgeStatus_, reader_->GetCurvedEdgesSelection());
|
||||
|
||||
reader_->UpdateProgress(0.2);
|
||||
|
||||
// Update the OpenFOAM mesh
|
||||
updateFoamMesh();
|
||||
reader_->UpdateProgress(0.5);
|
||||
@ -427,7 +418,6 @@ void Foam::vtkPVblockMesh::Update
|
||||
convertMeshEdges(output, blockNo);
|
||||
|
||||
reader_->UpdateProgress(0.8);
|
||||
|
||||
}
|
||||
|
||||
|
||||
|
||||
@ -25,13 +25,18 @@ Class
|
||||
Foam::vtkPVblockMesh
|
||||
|
||||
Description
|
||||
Provides a reader interface for OpenFOAM blockMesh to VTK interaction
|
||||
The backend for the vtkPVblockMeshReader reader module -
|
||||
providing a paraview reader interface for OpenFOAM blockMesh.
|
||||
|
||||
The block reader module can assist when creating a blockMeshDict
|
||||
for use with the blockMesh utility. As well as blocks, it can be
|
||||
used to visualize edges,corners and patch names.
|
||||
|
||||
There is no native VTK equivalent for this functionality.
|
||||
|
||||
SourceFiles
|
||||
vtkPVblockMesh.C
|
||||
vtkPVblockMeshConvert.C
|
||||
vtkPVblockMeshUpdate.C
|
||||
vtkPVblockMeshUtils.C
|
||||
|
||||
// Needed by VTK:
|
||||
vtkDataArrayTemplateImplicit.txx
|
||||
@ -95,12 +100,6 @@ class vtkPVblockMesh
|
||||
//- The mesh directory for the region
|
||||
fileName meshDir_;
|
||||
|
||||
//- Selected geometrical parts
|
||||
boolList blockStatus_;
|
||||
|
||||
//- Selected curved edges
|
||||
boolList edgeStatus_;
|
||||
|
||||
//- First instance and size of bleckMesh blocks
|
||||
// used to index into blockStatus_
|
||||
arrayRange rangeBlocks_;
|
||||
|
||||
@ -53,7 +53,11 @@ void Foam::vtkPVblockMesh::convertMeshBlocks
|
||||
Info<< "<beg> convertMeshBlocks" << endl;
|
||||
}
|
||||
|
||||
vtkDataArraySelection* selection = reader_->GetBlockSelection();
|
||||
const Map<string> blockStatus = getSelectedArrayMap
|
||||
(
|
||||
reader_->GetBlockSelection()
|
||||
);
|
||||
|
||||
arrayRange& range = rangeBlocks_;
|
||||
range.block(blockNo); // set output block
|
||||
label datasetNo = 0; // restart at dataset 0
|
||||
@ -61,29 +65,23 @@ void Foam::vtkPVblockMesh::convertMeshBlocks
|
||||
const blockMesh& blkMesh = *meshPtr_;
|
||||
const pointField blkPoints(blkMesh.vertices() * blkMesh.scaleFactor());
|
||||
|
||||
int blockI = 0;
|
||||
for
|
||||
(
|
||||
auto iter = range.cbegin();
|
||||
iter != range.cend();
|
||||
++iter, ++blockI
|
||||
)
|
||||
vtkIdType nodeIds[8]; // Space for VTK_HEXAHEDRON vertices
|
||||
int blockId = -1;
|
||||
for (auto partId : range)
|
||||
{
|
||||
const auto partId = *iter;
|
||||
if (!blockStatus_[partId])
|
||||
++blockId; // Increment first
|
||||
if (!blockStatus.found(partId))
|
||||
{
|
||||
continue;
|
||||
}
|
||||
const auto& longName = blockStatus[partId];
|
||||
|
||||
const blockDescriptor& blockDef = blkMesh[blockI];
|
||||
const blockDescriptor& blockDef = blkMesh[blockId];
|
||||
const labelList& blockLabels = blockDef.blockShape();
|
||||
|
||||
vtkSmartPointer<vtkPoints> vtkpoints =
|
||||
vtkSmartPointer<vtkPoints>::New();
|
||||
|
||||
auto vtkpoints = vtkSmartPointer<vtkPoints>::New();
|
||||
vtkpoints->SetNumberOfPoints(blockLabels.size());
|
||||
|
||||
vtkIdType nodeIds[8];
|
||||
forAll(blockLabels, pointi)
|
||||
{
|
||||
vtkpoints->SetPoint
|
||||
@ -94,8 +92,7 @@ void Foam::vtkPVblockMesh::convertMeshBlocks
|
||||
nodeIds[pointi] = pointi;
|
||||
}
|
||||
|
||||
vtkSmartPointer<vtkUnstructuredGrid> vtkmesh =
|
||||
vtkSmartPointer<vtkUnstructuredGrid>::New();
|
||||
auto vtkmesh = vtkSmartPointer<vtkUnstructuredGrid>::New();
|
||||
|
||||
vtkmesh->Allocate(1);
|
||||
vtkmesh->InsertNextCell
|
||||
@ -107,11 +104,7 @@ void Foam::vtkPVblockMesh::convertMeshBlocks
|
||||
|
||||
vtkmesh->SetPoints(vtkpoints);
|
||||
|
||||
addToBlock
|
||||
(
|
||||
output, vtkmesh, range, datasetNo,
|
||||
selection->GetArrayName(partId)
|
||||
);
|
||||
addToBlock(output, vtkmesh, range, datasetNo, longName);
|
||||
++datasetNo;
|
||||
}
|
||||
|
||||
@ -135,7 +128,11 @@ void Foam::vtkPVblockMesh::convertMeshEdges
|
||||
int& blockNo
|
||||
)
|
||||
{
|
||||
vtkDataArraySelection* selection = reader_->GetCurvedEdgesSelection();
|
||||
const Map<string> edgeStatus = getSelectedArrayMap
|
||||
(
|
||||
reader_->GetCurvedEdgesSelection()
|
||||
);
|
||||
|
||||
arrayRange& range = rangeEdges_;
|
||||
|
||||
range.block(blockNo); // set output block
|
||||
@ -145,24 +142,20 @@ void Foam::vtkPVblockMesh::convertMeshEdges
|
||||
const blockEdgeList& edges = blkMesh.edges();
|
||||
const scalar scaleFactor = blkMesh.scaleFactor();
|
||||
|
||||
int edgeI = 0;
|
||||
for
|
||||
(
|
||||
auto iter = range.cbegin();
|
||||
iter != range.cend();
|
||||
++iter, ++edgeI
|
||||
)
|
||||
int edgeId = -1;
|
||||
for (auto partId : range)
|
||||
{
|
||||
const auto partId = *iter;
|
||||
if (!edgeStatus_[partId])
|
||||
++edgeId; // Increment first
|
||||
if (!edgeStatus.found(partId))
|
||||
{
|
||||
continue;
|
||||
}
|
||||
const auto& longName = edgeStatus[partId];
|
||||
|
||||
// search each block
|
||||
forAll(blkMesh, blockI)
|
||||
// Search each block
|
||||
forAll(blkMesh, blockId)
|
||||
{
|
||||
const blockDescriptor& blockDef = blkMesh[blockI];
|
||||
const blockDescriptor& blockDef = blkMesh[blockId];
|
||||
|
||||
edgeList blkEdges = blockDef.blockShape().edges();
|
||||
|
||||
@ -175,7 +168,7 @@ void Foam::vtkPVblockMesh::convertMeshEdges
|
||||
label foundEdgeI = -1;
|
||||
forAll(blkEdges, blkEdgeI)
|
||||
{
|
||||
if (edges[edgeI].compare(blkEdges[blkEdgeI]))
|
||||
if (edges[edgeId].compare(blkEdges[blkEdgeI]))
|
||||
{
|
||||
foundEdgeI = blkEdgeI;
|
||||
break;
|
||||
@ -186,9 +179,7 @@ void Foam::vtkPVblockMesh::convertMeshEdges
|
||||
{
|
||||
const List<point>& edgePoints = edgesPoints[foundEdgeI];
|
||||
|
||||
vtkSmartPointer<vtkPoints> vtkpoints =
|
||||
vtkSmartPointer<vtkPoints>::New();
|
||||
|
||||
auto vtkpoints = vtkSmartPointer<vtkPoints>::New();
|
||||
vtkpoints->SetNumberOfPoints(edgePoints.size());
|
||||
|
||||
vtkIdType pointIds[edgePoints.size()];
|
||||
@ -200,8 +191,7 @@ void Foam::vtkPVblockMesh::convertMeshEdges
|
||||
pointIds[pointi] = pointi;
|
||||
}
|
||||
|
||||
vtkSmartPointer<vtkPolyData> vtkmesh =
|
||||
vtkSmartPointer<vtkPolyData>::New();
|
||||
auto vtkmesh = vtkSmartPointer<vtkPolyData>::New();
|
||||
|
||||
vtkmesh->Allocate(1);
|
||||
vtkmesh->InsertNextCell
|
||||
@ -213,11 +203,7 @@ void Foam::vtkPVblockMesh::convertMeshEdges
|
||||
|
||||
vtkmesh->SetPoints(vtkpoints);
|
||||
|
||||
addToBlock
|
||||
(
|
||||
output, vtkmesh, range, datasetNo,
|
||||
selection->GetArrayName(partId)
|
||||
);
|
||||
addToBlock(output, vtkmesh, range, datasetNo, longName);
|
||||
++datasetNo;
|
||||
|
||||
break;
|
||||
@ -253,35 +239,23 @@ void Foam::vtkPVblockMesh::convertMeshCorners
|
||||
|
||||
if (debug)
|
||||
{
|
||||
Info<< "<beg> convertMeshCorners" << endl;
|
||||
Info<< "<beg> " << FUNCTION_NAME << endl;
|
||||
}
|
||||
|
||||
if (true) // or some flag or other condition
|
||||
if (true) // Or some flag or other condition
|
||||
{
|
||||
vtkSmartPointer<vtkPolyData> vtkmesh =
|
||||
vtkSmartPointer<vtkPolyData>::New();
|
||||
|
||||
vtkSmartPointer<vtkPoints> vtkpoints =
|
||||
vtkSmartPointer<vtkPoints>::New();
|
||||
|
||||
vtkSmartPointer<vtkCellArray> vtkcells =
|
||||
vtkSmartPointer<vtkCellArray>::New();
|
||||
|
||||
auto vtkpoints = vtkSmartPointer<vtkPoints>::New();
|
||||
vtkpoints->SetNumberOfPoints(blkPoints.size());
|
||||
vtkcells->Allocate(2*blkPoints.size());
|
||||
// If reusing memory, ensure insert always starts from 0
|
||||
vtkcells->Reset();
|
||||
|
||||
vtkIdType pointId = 0;
|
||||
forAll(blkPoints, pointi)
|
||||
{
|
||||
vtkpoints->SetPoint(pointi, blkPoints[pointi].v_);
|
||||
vtkcells->InsertNextCell(1, &pointId); // VTK_VERTEX
|
||||
pointId++;
|
||||
}
|
||||
|
||||
auto vtkmesh = vtkSmartPointer<vtkPolyData>::New();
|
||||
|
||||
vtkmesh->SetPoints(vtkpoints);
|
||||
vtkmesh->SetVerts(vtkcells);
|
||||
vtkmesh->SetVerts(foamPvCore::identityVertices(blkPoints.size()));
|
||||
|
||||
addToBlock(output, vtkmesh, range, datasetNo, range.name());
|
||||
++datasetNo;
|
||||
@ -295,7 +269,7 @@ void Foam::vtkPVblockMesh::convertMeshCorners
|
||||
|
||||
if (debug)
|
||||
{
|
||||
Info<< "<end> convertMeshCorners" << endl;
|
||||
Info<< "<end> " << FUNCTION_NAME << endl;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
Reference in New Issue
Block a user