For example, to mesh a sphere with a single block the geometry is defined in the blockMeshDict as a searchableSurface: geometry { sphere { type searchableSphere; centre (0 0 0); radius 1; } } The vertices, block topology and curved edges are defined in the usual way, for example v 0.5773502; mv -0.5773502; a 0.7071067; ma -0.7071067; vertices ( ($mv $mv $mv) ( $v $mv $mv) ( $v $v $mv) ($mv $v $mv) ($mv $mv $v) ( $v $mv $v) ( $v $v $v) ($mv $v $v) ); blocks ( hex (0 1 2 3 4 5 6 7) (10 10 10) simpleGrading (1 1 1) ); edges ( arc 0 1 (0 $ma $ma) arc 2 3 (0 $a $ma) arc 6 7 (0 $a $a) arc 4 5 (0 $ma $a) arc 0 3 ($ma 0 $ma) arc 1 2 ($a 0 $ma) arc 5 6 ($a 0 $a) arc 4 7 ($ma 0 $a) arc 0 4 ($ma $ma 0) arc 1 5 ($a $ma 0) arc 2 6 ($a $a 0) arc 3 7 ($ma $a 0) ); which will produce a mesh in which the block edges conform to the sphere but the faces of the block lie somewhere between the original cube and the spherical surface which is a consequence of the edge-based transfinite interpolation. Now the projection of the block faces to the geometry specified above can also be specified: faces ( project (0 4 7 3) sphere project (2 6 5 1) sphere project (1 5 4 0) sphere project (3 7 6 2) sphere project (0 3 2 1) sphere project (4 5 6 7) sphere ); which produces a mesh that actually conforms to the sphere. See OpenFOAM-dev/tutorials/mesh/blockMesh/sphere This functionality is experimental and will undergo further development and generalization in the future to support more complex surfaces, feature edge specification and extraction etc. Please get involved if you would like to see blockMesh become a more flexible block-structured mesher. Henry G. Weller, CFD Direct.
320 lines
7.8 KiB
C
320 lines
7.8 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
|
under the terms of the GNU General Public License as published by
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
(at your option) any later version.
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
for more details.
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "vtkPV3blockMesh.H"
|
|
#include "vtkPV3blockMeshReader.h"
|
|
|
|
// OpenFOAM includes
|
|
#include "blockMesh.H"
|
|
#include "Time.H"
|
|
|
|
#include "vtkOpenFOAMPoints.H"
|
|
|
|
// VTK includes
|
|
#include "vtkCellArray.h"
|
|
#include "vtkDataArraySelection.h"
|
|
#include "vtkMultiBlockDataSet.h"
|
|
#include "vtkPoints.h"
|
|
#include "vtkPolyData.h"
|
|
#include "vtkUnstructuredGrid.h"
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
void Foam::vtkPV3blockMesh::convertMeshBlocks
|
|
(
|
|
vtkMultiBlockDataSet* output,
|
|
int& blockNo
|
|
)
|
|
{
|
|
vtkDataArraySelection* selection = reader_->GetBlockSelection();
|
|
arrayRange& range = arrayRangeBlocks_;
|
|
range.block(blockNo); // set output block
|
|
label datasetNo = 0; // restart at dataset 0
|
|
|
|
const blockMesh& blkMesh = *meshPtr_;
|
|
const Foam::pointField& blockPoints = blkMesh.vertices();
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<beg> Foam::vtkPV3blockMesh::convertMeshBlocks" << endl;
|
|
}
|
|
|
|
int blockI = 0;
|
|
const scalar scaleFactor = blkMesh.scaleFactor();
|
|
|
|
for
|
|
(
|
|
int partId = range.start();
|
|
partId < range.end();
|
|
++partId, ++blockI
|
|
)
|
|
{
|
|
if (!blockStatus_[partId])
|
|
{
|
|
continue;
|
|
}
|
|
|
|
const blockDescriptor& blockDef = blkMesh[blockI];
|
|
|
|
vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();
|
|
|
|
// Convert OpenFOAM mesh vertices to VTK
|
|
vtkPoints *vtkpoints = vtkPoints::New();
|
|
vtkpoints->Allocate( blockDef.nPoints() );
|
|
const labelList& blockLabels = blockDef.blockShape();
|
|
|
|
vtkmesh->Allocate(1);
|
|
vtkIdType nodeIds[8];
|
|
|
|
forAll(blockLabels, ptI)
|
|
{
|
|
vtkInsertNextOpenFOAMPoint
|
|
(
|
|
vtkpoints,
|
|
blockPoints[blockLabels[ptI]],
|
|
scaleFactor
|
|
);
|
|
|
|
nodeIds[ptI] = ptI;
|
|
}
|
|
|
|
vtkmesh->InsertNextCell
|
|
(
|
|
VTK_HEXAHEDRON,
|
|
8,
|
|
nodeIds
|
|
);
|
|
|
|
vtkmesh->SetPoints(vtkpoints);
|
|
vtkpoints->Delete();
|
|
|
|
AddToBlock
|
|
(
|
|
output, vtkmesh, range, datasetNo,
|
|
selection->GetArrayName(partId)
|
|
);
|
|
|
|
vtkmesh->Delete();
|
|
datasetNo++;
|
|
}
|
|
|
|
|
|
// anything added?
|
|
if (datasetNo)
|
|
{
|
|
++blockNo;
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<end> Foam::vtkPV3blockMesh::convertMeshBlocks" << endl;
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::vtkPV3blockMesh::convertMeshEdges
|
|
(
|
|
vtkMultiBlockDataSet* output,
|
|
int& blockNo
|
|
)
|
|
{
|
|
vtkDataArraySelection* selection = reader_->GetCurvedEdgesSelection();
|
|
arrayRange& range = arrayRangeEdges_;
|
|
|
|
range.block(blockNo); // set output block
|
|
label datasetNo = 0; // restart at dataset 0
|
|
|
|
const blockMesh& blkMesh = *meshPtr_;
|
|
const blockEdgeList& edges = blkMesh.edges();
|
|
|
|
int edgeI = 0;
|
|
const scalar scaleFactor = blkMesh.scaleFactor();
|
|
|
|
for
|
|
(
|
|
int partId = range.start();
|
|
partId < range.end();
|
|
++partId, ++edgeI
|
|
)
|
|
{
|
|
if (!edgeStatus_[partId])
|
|
{
|
|
continue;
|
|
}
|
|
|
|
// search each block
|
|
forAll(blkMesh, blockI)
|
|
{
|
|
const blockDescriptor& blockDef = blkMesh[blockI];
|
|
|
|
edgeList blkEdges = blockDef.blockShape().edges();
|
|
|
|
// find the corresponding edge within the block
|
|
label foundEdgeI = -1;
|
|
forAll(blkEdges, blkEdgeI)
|
|
{
|
|
if (edges[edgeI].compare(blkEdges[blkEdgeI]))
|
|
{
|
|
foundEdgeI = blkEdgeI;
|
|
break;
|
|
}
|
|
}
|
|
|
|
if (foundEdgeI != -1)
|
|
{
|
|
const List<point>& edgePoints =
|
|
blockDef.blockEdgePoints()[foundEdgeI];
|
|
|
|
|
|
vtkPolyData* vtkmesh = vtkPolyData::New();
|
|
vtkPoints* vtkpoints = vtkPoints::New();
|
|
|
|
vtkpoints->Allocate( edgePoints.size() );
|
|
vtkmesh->Allocate(1);
|
|
|
|
vtkIdType pointIds[edgePoints.size()];
|
|
forAll(edgePoints, ptI)
|
|
{
|
|
vtkInsertNextOpenFOAMPoint
|
|
(
|
|
vtkpoints,
|
|
edgePoints[ptI],
|
|
scaleFactor
|
|
);
|
|
pointIds[ptI] = ptI;
|
|
}
|
|
|
|
vtkmesh->InsertNextCell
|
|
(
|
|
VTK_POLY_LINE,
|
|
edgePoints.size(),
|
|
pointIds
|
|
);
|
|
|
|
vtkmesh->SetPoints(vtkpoints);
|
|
vtkpoints->Delete();
|
|
|
|
AddToBlock
|
|
(
|
|
output, vtkmesh, range, datasetNo,
|
|
selection->GetArrayName(partId)
|
|
);
|
|
|
|
vtkmesh->Delete();
|
|
datasetNo++;
|
|
|
|
break;
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
// anything added?
|
|
if (datasetNo)
|
|
{
|
|
++blockNo;
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<end> Foam::vtkPV3blockMesh::convertMeshEdges" << endl;
|
|
}
|
|
|
|
}
|
|
|
|
|
|
void Foam::vtkPV3blockMesh::convertMeshCorners
|
|
(
|
|
vtkMultiBlockDataSet* output,
|
|
int& blockNo
|
|
)
|
|
{
|
|
arrayRange& range = arrayRangeCorners_;
|
|
range.block(blockNo); // set output block
|
|
label datasetNo = 0; // restart at dataset 0
|
|
|
|
const pointField& blockPoints = meshPtr_->vertices();
|
|
const scalar& scaleFactor = meshPtr_->scaleFactor();
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<beg> Foam::vtkPV3blockMesh::convertMeshCorners" << endl;
|
|
}
|
|
|
|
if (true) // or some flag or other condition
|
|
{
|
|
vtkPolyData* vtkmesh = vtkPolyData::New();
|
|
vtkPoints* vtkpoints = vtkPoints::New();
|
|
vtkCellArray* vtkcells = vtkCellArray::New();
|
|
|
|
vtkpoints->Allocate( blockPoints.size() );
|
|
vtkcells->Allocate( blockPoints.size() );
|
|
|
|
vtkIdType pointId = 0;
|
|
forAll(blockPoints, ptI)
|
|
{
|
|
vtkInsertNextOpenFOAMPoint
|
|
(
|
|
vtkpoints,
|
|
blockPoints[ptI],
|
|
scaleFactor
|
|
);
|
|
|
|
vtkcells->InsertNextCell(1, &pointId);
|
|
pointId++;
|
|
}
|
|
|
|
vtkmesh->SetPoints(vtkpoints);
|
|
vtkpoints->Delete();
|
|
|
|
vtkmesh->SetVerts(vtkcells);
|
|
vtkcells->Delete();
|
|
|
|
AddToBlock
|
|
(
|
|
output, vtkmesh, range, datasetNo,
|
|
arrayRangeCorners_.name()
|
|
);
|
|
vtkmesh->Delete();
|
|
|
|
datasetNo++;
|
|
}
|
|
|
|
// anything added?
|
|
if (datasetNo)
|
|
{
|
|
++blockNo;
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<end> Foam::vtkPV3blockMesh::convertMeshCorners" << endl;
|
|
}
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|