FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

Post-processing mesh subset tool. Given the original mesh and the list of selected cells, it creates the mesh consisting only of the desired cells, with the mapping list for points, faces, and cells. More...

#include <finiteVolume/fvMeshSubset.H>


Detailed Description

Post-processing mesh subset tool. Given the original mesh and the list of selected cells, it creates the mesh consisting only of the desired cells, with the mapping list for points, faces, and cells.

Puts all exposed internal faces into either

  • a user supplied patch
  • a newly created patch "oldInternalFaces"
  • setCellSubset is for small subsets. Uses Maps to minimize memory.
  • setLargeCellSubset is for largish subsets (>10% of mesh). Uses labelLists instead.
  • setLargeCellSubset does coupled patch subsetting as well. If it detects a face on a coupled patch 'losing' its neighbour it will move the face into the oldInternalFaces patch.
  • if a user supplied patch is used the mapping becomes a problem. Do the new faces get the value of the internal face they came from? What if e.g. the user supplied patch is a fixedValue 0? So for now they get the face of existing patch face 0.
Source files

Definition at line 74 of file fvMeshSubset.H.

Collaboration diagram for fvMeshSubset:

List of all members.

Classes

class  patchFieldSubset
 Patch-field subset interpolation class. More...
class  pointPatchFieldSubset
 Patch-field subset interpolation class. More...

Public Member Functions

 fvMeshSubset (const fvMesh &)
 Construct given a mesh to subset.
void  setCellSubset (const labelHashSet &globalCellMap, const label patchID=-1)
 Set the subset. Create "oldInternalFaces" patch for exposed.
void  setLargeCellSubset (const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
 Set the subset from all cells with region == currentRegion.
void  setLargeCellSubset (const labelHashSet &globalCellMap, const label patchID=-1, const bool syncPar=true)
 setLargeCellSubset but with labelHashSet.
const fvMesh &  baseMesh () const
 Original mesh.
const fvMesh &  subMesh () const
 Return reference to subset mesh.
fvMesh &  subMesh ()
const labelList &  pointMap () const
 Return point map.
const labelList &  faceMap () const
 Return face map.
const labelList &  cellMap () const
 Return cell map.
const labelList &  patchMap () const
 Return patch map.
template<class Type >
tmp< GeometricField< Type,
fvPatchField, volMesh > >  
interpolate (const GeometricField< Type, fvPatchField, volMesh > &) const
template<class Type >
tmp< GeometricField< Type,
fvsPatchField, surfaceMesh > >  
interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &) const
template<class Type >
tmp< GeometricField< Type,
pointPatchField, pointMesh > >  
interpolate (const GeometricField< Type, pointPatchField, pointMesh > &) const

Static Public Member Functions

template<class Type >
static tmp< GeometricField
< Type, fvPatchField, volMesh > >  
interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
 Map volume field.
template<class Type >
static tmp< GeometricField
< Type, fvsPatchField,
surfaceMesh > >  
interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &faceMap)
 Map surface field.
template<class Type >
static tmp< GeometricField
< Type, pointPatchField,
pointMesh > >  
interpolate (const GeometricField< Type, pointPatchField, pointMesh > &, const pointMesh &sMesh, const labelList &patchMap, const labelList &pointMap)
 Map point field.

Constructor & Destructor Documentation

fvMeshSubset ( const fvMesh &   baseMesh  ) [explicit]

Construct given a mesh to subset.

Definition at line 370 of file fvMeshSubset.C.


Member Function Documentation

void setCellSubset ( const labelHashSet &   globalCellMap,
const label   patchID = -1  
)
void setLargeCellSubset ( const labelList &   region,
const label   currentRegion,
const label   patchID = -1,
const bool   syncCouples = true  
)

Set the subset from all cells with region == currentRegion.

Create "oldInternalFaces" patch for exposed internal faces (patchID==-1) or use supplied patch. Handles coupled patches by if nessecary making coupled patch face part of patchID (so uncoupled)

Definition at line 767 of file fvMeshSubset.C.

References Foam::abort(), IOobject::clone(), Foam::endl(), Foam::FatalError, FatalErrorIn, polyBoundaryMesh::findPatchID(), forAll, Pstream::gatherList(), Foam::identity(), Pstream::listCombineGather(), Pstream::listCombineScatter(), Pstream::myProcNo(), Foam::name(), polyBoundaryMesh::names(), IOobject::NO_READ, IOobject::NO_WRITE, Pstream::nProcs(), Pstream::parRun(), patchNames(), Foam::reduce(), Pstream::scatterList(), List< T >::setSize(), PtrList< T >::size(), List< T >::size(), polyPatch::start(), timeName, and Foam::xferMove().

Referenced by fvMeshDistribute::distribute().

void setLargeCellSubset ( const labelHashSet &   globalCellMap,
const label   patchID = -1,
const bool   syncPar = true  
)

setLargeCellSubset but with labelHashSet.

Definition at line 1362 of file fvMeshSubset.C.

References forAllConstIter.

const fvMesh& baseMesh (  ) const [inline]

Original mesh.

Definition at line 268 of file fvMeshSubset.H.

const fvMesh & subMesh (  ) const

Return reference to subset mesh.

Definition at line 1378 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute(), and vtkMesh::mesh().

fvMesh & subMesh (  )

Definition at line 1386 of file fvMeshSubset.C.

const labelList & pointMap (  ) const

Return point map.

Definition at line 1394 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute().

const labelList & faceMap (  ) const

Return face map.

Definition at line 1402 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute().

const labelList & cellMap (  ) const

Return cell map.

Definition at line 1410 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute().

const labelList & patchMap (  ) const

Return patch map.

Definition at line 1418 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute().

tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate ( const GeometricField< Type, fvPatchField, volMesh > &   vf  ) const

Definition at line 140 of file fvMeshSubsetInterpolate.C.

References Foam::interpolate().

tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate ( const GeometricField< Type, fvsPatchField, surfaceMesh > &   sf  ) const

Definition at line 281 of file fvMeshSubsetInterpolate.C.

References Foam::interpolate().

tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate ( const GeometricField< Type, pointPatchField, pointMesh > &   sf  ) const

The documentation for this class was generated from the following files: