NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
meshBase Class Referenceabstract

A brief description of meshBase. More...

Detailed Description

Note
Virtual methods are usually implemented in vtkMesh.C. We use that class for more VTK-specific functions that we want wrapped by meshBase
Try using smart pointer instances of meshBase in drivers and pass the underlying raw pointer with the get member function of the smart pointer when required (e.g., source->transfer(target.get(), ...)) for smart pointer instance of target. Ideally, driver classes should have empty destructors, which can be realized by using smart pointers. In general, smart pointer should be used and scoped when they can.

Definition at line 64 of file meshBase.H.

Public Member Functions

 meshBase ()
 
virtual ~meshBase ()
 
virtual void read (const std::string &fname)=0
 abstract read method reserved for derived classes More...
 
virtual std::vector< double > getPoint (nemId_t id) const =0
 get point with id More...
 
virtual std::vector< std::vector< double > > getVertCrds () const =0
 get 3 vecs with x,y and z coords More...
 
virtual std::map< nemId_t, std::vector< double > > getCell (nemId_t id) const =0
 get cell with id More...
 
virtual std::vector< std::vector< double > > getCellVec (nemId_t id) const =0
 get vector of coords of cell with id More...
 
virtual void inspectEdges (const std::string &ofname) const =0
 get edge lengths of dataSet More...
 
vtkSmartPointer< vtkDataSet > getDataSet () const
 get this meshes' dataSet More...
 
virtual vtkSmartPointer< vtkDataSet > extractSurface ()=0
 extract the surface mesh More...
 
virtual void setPointDataArray (const std::string &name, const std::vector< std::vector< double >> &data)
 register data to dataSet's point data More...
 
virtual void setPointDataArray (const std::string &name, const std::vector< double > &data)
 register data to dataSet's point data More...
 
virtual void setCellDataArray (const std::string &name, const std::vector< std::vector< double >> &data)
 register data to dataSet's cell data More...
 
virtual void setCellDataArray (const std::string &name, const std::vector< double > &data)
 register data to dataSet's cell data More...
 
virtual void getPointDataArray (const std::string &name, std::vector< double > &data)
 get scalar point or cell data array. More...
 
virtual void getPointDataArray (int arrayId, std::vector< double > &data)
 get scalar point or cell data array. More...
 
virtual int getCellDataIdx (const std::string &name)
 <> More...
 
virtual void getCellDataArray (const std::string &name, std::vector< double > &data)
 <> More...
 
virtual void getCellDataArray (int arrayId, std::vector< double > &data)
 <> More...
 
virtual void unsetPointDataArray (int arrayID)
 delete array with id from dataSet's point data More...
 
virtual void unsetPointDataArray (const std::string &name)
 <> More...
 
virtual void unsetCellDataArray (int arrayID)
 delete array with id from dataSet's cell data More...
 
virtual void unsetCellDataArray (const std::string &name)
 <> More...
 
virtual void unsetFieldDataArray (const std::string &name)
 delete array with id from dataSet's field data More...
 
virtual std::vector< double > getCellLengths () const =0
 get diameter of circumsphere of each cell More...
 
virtual std::vector< double > getCellCenter (nemId_t cellID) const =0
 get center of a cell More...
 
vtkSmartPointer< vtkStaticCellLocator > buildStaticCellLocator ()
 build locators for efficient search operations More...
 
vtkSmartPointer< vtkStaticPointLocator > buildStaticPointLocator ()
 build thread-safe point locator for efficient search operations More...
 
virtual int getCellType () const =0
 get cell type as an integer assumes all elements are the same type More...
 
virtual std::vector< nemId_tgetConnectivities () const =0
 get connectivities. More...
 
void setMetadata (vtkSmartPointer< vtkModelMetadata > _metadata)
 
vtkSmartPointer< vtkModelMetadata > getMetadata ()
 
std::vector< std::vector< double > > integrateOverMesh (const std::vector< int > &arrayIDs)
 integrate arrays in arrayIDs over the mesh. More...
 
void generateSizeField (const std::string &method, int arrayID, double dev_mlt, bool maxIsmin, double sizeFactor=1.0, int order=1)
 generate size field based on method and given a point data array. More...
 
int IsArrayName (const std::string &name, bool pointOrCell=false) const
 check for named array in vtk and return its integer id. More...
 
void refineMesh (const std::string &method, int arrayID, double dev_mult, bool maxIsmin, double edge_scale, const std::string &ofname, bool transferData, double sizeFactor=1., bool constrainBoundary=false)
 perform sizefield-based h-refinement. More...
 
void refineMesh (const std::string &method, const std::string &arrayName, double dev_mult, bool maxIsmin, double edge_scale, const std::string &ofname, bool transferData, double sizeFactor=1.)
 perform sizefield-based h-refinement. More...
 
void refineMesh (const std::string &method, double edge_scale, const std::string &ofname, bool transferData, bool constrainBoundary=false)
 added for uniform refinement by driver More...
 
void refineMesh (const std::string &method, int arrayID, int order, const std::string &ofname, bool transferData)
 <> More...
 
void refineMesh (const std::string &method, const std::string &arrayName, int order, const std::string &ofname, bool transferData)
 <> More...
 
virtual void report () const
 generate a report of the mesh More...
 
nemId_t getNumberOfPoints () const
 return the number of points More...
 
nemId_t getNumberOfCells () const
 return the number of cells More...
 
void checkMesh (const std::string &ofname) const
 <> More...
 
std::map< nemId_t, nemId_tgetGlobToPartNodeMap ()
 global to local mapping of nodes More...
 
std::map< nemId_t, nemId_tgetGlobToPartCellMap ()
 global to local mapping of cells More...
 
std::map< nemId_t, nemId_tgetPartToGlobNodeMap ()
 local to global mapping of nodes More...
 
std::map< nemId_t, nemId_tgetPartToGlobCellMap ()
 local to global mapping of cells More...
 
virtual void write () const
 write the mesh to file named after the private var 'filename'. More...
 
virtual void write (const std::string &fname) const =0
 write the mesh to file named fname More...
 
void writeMSH (std::ofstream &outputStream)
 convert to gmsh format without data More...
 
void writeMSH (const std::string &fname)
 convert to gmsh format without data More...
 
void writeMSH (std::ofstream &outputStream, const std::string &pointOrCell, int arrayID)
 convert to gmsh format with specified point or cell data More...
 
void writeMSH (const std::string &fname, const std::string &pointOrCell, int arrayID)
 convert to gmsh format without data More...
 
void writeMSH (std::ofstream &outputStream, const std::string &pointOrCell, int arrayID, bool onlyVol)
 convert to gmsh format with specified point or cell data for only volume elements (USE ONLY FOR MADLIB STUFF) More...
 
void writeMSH (const std::string &fname, const std::string &pointOrCell, int arrayID, bool onlyVol)
 convert to gmsh format with specified point or cell data for only volume elements (USE ONLY FOR MADLIB STUFF) More...
 
void writeCobalt (meshBase *surfWithPatch, const std::string &mapFile, std::ofstream &outputStream)
 surfWithPatch must have patchNo array More...
 
void writeCobalt (meshBase *surfWithPatch, const std::string &mapFile, const std::string &ofname)
 surfWithPatch must have patchNo array More...
 
void setFileName (const std::string &fname)
 set the file name. More...
 
const std::string & getFileName () const
 get the current file name More...
 
void setCheckQuality (bool x)
 set whether to check quality of transfer by back-transfer and rmse More...
 
void setContBool (bool x)
 set weighted averaging/smoothing for cell data transfer (default is off) More...
 
meshBaseconvertQuads ()
 
std::vector< std::string > getNewArrayNames ()
 get new array names for use in transfer More...
 
std::vector< int > getArrayIDs (std::vector< std::string > arrayNames, bool fromPointArrays=false)
 given array names, return corresponding ids More...
 
void convertHexToTetVTK (vtkSmartPointer< vtkDataSet > meshdataSet)
 Converts given hexahedral VTK dataset into tetrahedral mesh and stores it into dataSet variable. More...
 

Static Public Member Functions

static meshBaseCreate (const std::string &fname)
 Construct vtkMesh from filename. More...
 
static meshBaseCreate (vtkSmartPointer< vtkDataSet > other, const std::string &newname)
 Construct from existing vtkDataSet and assign newname as filename. More...
 
static meshBaseCreate (const std::vector< double > &xCrds, const std::vector< double > &yCrds, const std::vector< double > &zCrds, const std::vector< nemId_t > &elmConn, const int cellType, const std::string &newname)
 create from coordinates and connectivities. More...
 
static std::shared_ptr< meshBaseCreateShared (const std::string &fname)
 Create shared ptr from fname. More...
 
static std::shared_ptr< meshBaseCreateShared (meshBase *mesh)
 Create shared ptr from existing meshbase. More...
 
static std::shared_ptr< meshBaseCreateShared (vtkSmartPointer< vtkDataSet > other, const std::string &newname)
 Create shared ptr from existing vtkDataset and assign newname as filename. More...
 
static std::shared_ptr< meshBaseCreateShared (const std::vector< double > &xCrds, const std::vector< double > &yCrds, const std::vector< double > &zCrds, const std::vector< nemId_t > &elmConn, int cellType, const std::string &newname)
 Version of raw data mesh creation for memory managed shared_ptr instance. More...
 
static std::unique_ptr< meshBaseCreateUnique (const std::string &fname)
 create unique ptr from fname More...
 
static std::unique_ptr< meshBaseCreateUnique (const std::vector< double > &xCrds, const std::vector< double > &yCrds, const std::vector< double > &zCrds, const std::vector< nemId_t > &elmConn, int cellType, const std::string &newname)
 version of raw data mesh creation for memory managed unique ptr instance More...
 
static std::unique_ptr< meshBaseCreateUnique (vtkSmartPointer< vtkDataSet > other, const std::string &newname)
 construct from existing vtkDataSet and assign newname as filename More...
 
static std::unique_ptr< meshBaseCreateUnique (meshBase *mesh)
 construct from existing meshbase object More...
 
static meshBaseexportGmshToVtk (const std::string &fname)
 construct vtkMesh from gmsh msh file (called in Create methods) More...
 
static meshBaseexportVolToVtk (const std::string &fname)
 construct vtkMesh from netgen vol file (called in Create methods) More...
 
static meshBaseexportPntToVtk (const std::string &fname)
 construct vtkMesh from netgen vol file (called in Create methods) More...
 
static meshBaseexportExoToVtk (const std::string &fname)
 construct vtkMesh from exodusII files More...
 
static meshBasestitchMB (const std::vector< meshBase *> &mbObjs)
 stitch together several meshBases More...
 
static std::shared_ptr< meshBasestitchMB (const std::vector< std::shared_ptr< meshBase >> &_mbObjs)
 stitch together several meshBase More...
 
static std::vector< std::shared_ptr< meshBase > > partition (const meshBase *mbObj, int numPartitions)
 mesh partitioning (with METIS) More...
 
static meshBaseextractSelectedCells (meshBase *mesh, const std::vector< nemId_t > &cellIds)
 extract subset of mesh given list of cell ids and return meshBase obj More...
 
static meshBaseextractSelectedCells (vtkSmartPointer< vtkDataSet > mesh, vtkSmartPointer< vtkIdTypeArray > cellIds)
 helper wrapped by function above More...
 

Protected Attributes

nemId_t numPoints
 number of points in mesh More...
 
nemId_t numCells
 number of cells in mesh More...
 
vtkSmartPointer< vtkDataSet > dataSet
 mesh points, topology and data More...
 
std::string filename
 name of mesh file More...
 
bool checkQuality
 check transfer quality when on More...
 
bool continuous
 switch on / off weighted averaging for cell data transfer (default is off) More...
 
std::vector< std::string > newArrayNames
 new names to set for transferred data More...
 
std::map< nemId_t, nemId_tglobToPartNodeMap
 map between global and local node idx in partition for distributed data sets More...
 
std::map< nemId_t, nemId_tglobToPartCellMap
 map between global and local cell idx in partition More...
 
std::map< nemId_t, nemId_tpartToGlobNodeMap
 map between local and global node idx in partition More...
 
std::map< nemId_t, nemId_tpartToGlobCellMap
 map between local and global cell idx in partition More...
 
vtkSmartPointer< vtkModelMetadata > metadata
 

Inherited by FOAM::foamMesh, gmshMesh, meshSrch, and vtkMesh.

Constructor & Destructor Documentation

◆ meshBase()

meshBase::meshBase ( )
inline

Definition at line 67 of file meshBase.H.

68  : numPoints(0),
69  numCells(0),
70  dataSet(nullptr),
71  checkQuality(false),
72  continuous(false),
73  metadata(nullptr) {
74  std::cout << "meshBase constructed" << std::endl;
75  }
nemId_t numPoints
number of points in mesh
Definition: meshBase.H:714
bool checkQuality
check transfer quality when on
Definition: meshBase.H:730
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
nemId_t numCells
number of cells in mesh
Definition: meshBase.H:718
vtkSmartPointer< vtkModelMetadata > metadata
Definition: meshBase.H:760
bool continuous
switch on / off weighted averaging for cell data transfer (default is off)
Definition: meshBase.H:735

◆ ~meshBase()

virtual meshBase::~meshBase ( )
inlinevirtual

Definition at line 95 of file meshBase.H.

References cellType, and mesh.

95 { std::cout << "meshBase destroyed" << std::endl; }

Member Function Documentation

◆ buildStaticCellLocator()

vtkSmartPointer< vtkStaticCellLocator > meshBase::buildStaticCellLocator ( )
Returns
<>

Definition at line 1640 of file meshBase.C.

References NEM::MSH::New().

Referenced by writeCobalt().

1640  {
1641  vtkSmartPointer<vtkStaticCellLocator> cellLocator =
1643  cellLocator->SetDataSet(dataSet);
1644  cellLocator->BuildLocator();
1645  return cellLocator;
1646 }
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
geoMeshBase * New(MeshType meshType)
Create a new mesh object.

◆ buildStaticPointLocator()

vtkSmartPointer< vtkStaticPointLocator > meshBase::buildStaticPointLocator ( )
Returns
<>

Definition at line 1648 of file meshBase.C.

References NEM::MSH::New().

1648  {
1649  auto pointLocator = vtkSmartPointer<vtkStaticPointLocator>::New();
1650  pointLocator->SetDataSet(dataSet);
1651  pointLocator->BuildLocator();
1652  return pointLocator;
1653 }
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
geoMeshBase * New(MeshType meshType)
Create a new mesh object.

◆ checkMesh()

void meshBase::checkMesh ( const std::string &  ofname) const
Parameters
ofname<>

Definition at line 1657 of file meshBase.C.

References MeshQuality::checkMesh().

1657  {
1658  std::unique_ptr<MeshQuality> qualCheck =
1659  std::unique_ptr<MeshQuality>(new MeshQuality(this));
1660  qualCheck->checkMesh(ofname);
1661 }
void checkMesh(std::ostream &outputStream)
Definition: MeshQuality.C:64

◆ convertHexToTetVTK()

void meshBase::convertHexToTetVTK ( vtkSmartPointer< vtkDataSet >  meshdataSet)
Parameters
meshdataSetInput hexahedral mesh dataset

Definition at line 1826 of file meshBase.C.

References NEM::MSH::New().

1826  {
1827  vtkSmartPointer<vtkDataSetTriangleFilter> triFilter =
1829  triFilter->SetInputData(meshdataSet);
1830  triFilter->Update();
1831 
1832  dataSet = triFilter->GetOutput();
1833 }
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
geoMeshBase * New(MeshType meshType)
Create a new mesh object.

◆ convertQuads()

meshBase * meshBase::convertQuads ( )

Definition at line 1769 of file meshBase.C.

References NEM::MSH::New(), and points.

Referenced by proteusHdf5::proteusHdf5().

1769  {
1770  // Create new dataset
1771  vtkSmartPointer<vtkUnstructuredGrid> vtkDataSetNoQuads =
1773  // Accumulate vector of ID lists for quads removed
1774  std::vector<std::vector<int>> removedQuadIds;
1775  // Set points
1776  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
1777  points->SetNumberOfPoints(this->getDataSet()->GetNumberOfPoints());
1778  double pt[3];
1779  for (int iPoint = 0; iPoint < this->getDataSet()->GetNumberOfPoints();
1780  iPoint++) {
1781  this->getDataSet()->GetPoint(iPoint, pt);
1782  points->SetPoint(iPoint, pt[0], pt[1], pt[2]);
1783  }
1784  vtkDataSetNoQuads->SetPoints(points);
1785 
1786  // loop through elements and remove quads by splitting into tris
1787  for (int iCell = 0; iCell < this->getDataSet()->GetNumberOfCells(); iCell++) {
1788  // is quad?
1789  if (this->getDataSet()->GetCell(iCell)->GetCellType() == VTK_QUAD) {
1790  vtkSmartPointer<vtkIdList> elmIds =
1791  this->getDataSet()->GetCell(iCell)->GetPointIds();
1792  std::vector<int> elmIdsVec = {static_cast<int>(elmIds->GetId(0)),
1793  static_cast<int>(elmIds->GetId(1)),
1794  static_cast<int>(elmIds->GetId(2)),
1795  static_cast<int>(elmIds->GetId(3))};
1796  removedQuadIds.push_back(elmIdsVec);
1797  }
1798  // only add non-quads to vtkdataset
1799  else {
1800  vtkDataSetNoQuads->InsertNextCell(
1801  this->getDataSet()->GetCell(iCell)->GetCellType(),
1802  this->getDataSet()->GetCell(iCell)->GetPointIds());
1803  }
1804  }
1805 
1806  // Add in new tris based on accumulated quad Ids
1807  for (auto elemItr = removedQuadIds.begin(); elemItr != removedQuadIds.end();
1808  elemItr++) {
1809  vtkSmartPointer<vtkIdList> triList1 = vtkSmartPointer<vtkIdList>::New();
1810  triList1->InsertNextId((*elemItr)[0]);
1811  triList1->InsertNextId((*elemItr)[1]);
1812  triList1->InsertNextId((*elemItr)[2]);
1813  vtkDataSetNoQuads->InsertNextCell(VTK_TRIANGLE, triList1);
1814 
1815  vtkSmartPointer<vtkIdList> triList2 = vtkSmartPointer<vtkIdList>::New();
1816  triList2->InsertNextId((*elemItr)[0]);
1817  triList2->InsertNextId((*elemItr)[2]);
1818  triList2->InsertNextId((*elemItr)[3]);
1819 
1820  vtkDataSetNoQuads->InsertNextCell(VTK_TRIANGLE, triList2);
1821  }
1822 
1823  return Create(vtkDataSetNoQuads, "removedQuads.vtu");
1824 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133

◆ Create() [1/3]

meshBase * meshBase::Create ( const std::string &  fname)
static

This method calls the other factory methods based on extension.

Parameters
fnamename of mesh file
Returns
<>

Caller must delete object after use.

Definition at line 78 of file meshBase.C.

Referenced by OrderOfAccuracy::computeMeshWithResolution(), NEM::GEO::rocPack::createCohesiveElements(), CreateShared(), CreateUnique(), NEM::DRV::BlockMeshMeshGenDriver::execute(), NEM::DRV::SnappyMeshMeshGenDriver::execute(), NEM::DRV::VtkToPntConversionDriver::execute(), NEM::DRV::TransferDriver::execute(), NEM::DRV::CheckMeshQualDriver::execute(), NEM::DRV::GmshToExoConversionDriver::execute(), hdf5Reader::exportToMeshBase(), extractSelectedCells(), NEM::DRV::MeshGenDriver::MeshGenDriver(), NEM::DRV::ConversionDriver::procExo(), proteusHdf5::proteusHdf5(), and stitchMB().

78  {
79  if (fname.find(".vt") != std::string::npos ||
80  fname.find(".stl") != std::string::npos) {
81  auto *vtkmesh = new vtkMesh(fname);
82  vtkmesh->setFileName(fname);
83  return vtkmesh;
84  } else if (fname.find(".msh") != std::string::npos) {
85  std::cout << "Detected file in GMSH format" << std::endl;
86  std::cout << "Exporting to VTK ...." << std::endl;
87  return exportGmshToVtk(fname);
88  } else if (fname.find(".vol") != std::string::npos) {
89  std::cout << "Detected file in Netgen .vol format" << std::endl;
90  std::cout << "Exporting to VTK ...." << std::endl;
91  return exportVolToVtk(fname);
92  } else if (fname.find(".pntmesh") != std::string::npos) {
93  std::cout << "Detected file in PNTmesh format" << std::endl;
94  std::cout << "Processing the file ...." << std::endl;
95  return exportPntToVtk(fname);
96  } else if (fname.find(".g") != std::string::npos ||
97  fname.find(".exo") != std::string::npos ||
98  fname.find(".e") != std::string::npos) {
99  std::cout << "Detected file in Exodus II format" << std::endl;
100  std::cout << "Processing the file ...." << std::endl;
101  return exportExoToVtk(fname);
102  } else {
103  std::cout << "mesh files with extension "
104  << fname.substr(fname.find_last_of('.')) << " are not supported!"
105  << std::endl;
106  exit(1);
107  }
108 }
static meshBase * exportGmshToVtk(const std::string &fname)
construct vtkMesh from gmsh msh file (called in Create methods)
Definition: meshBase.C:409
static meshBase * exportVolToVtk(const std::string &fname)
construct vtkMesh from netgen vol file (called in Create methods)
Definition: meshBase.C:766
static meshBase * exportPntToVtk(const std::string &fname)
construct vtkMesh from netgen vol file (called in Create methods)
Definition: meshBase.C:848
static meshBase * exportExoToVtk(const std::string &fname)
construct vtkMesh from exodusII files
Definition: meshBase.C:918

◆ Create() [2/3]

meshBase * meshBase::Create ( vtkSmartPointer< vtkDataSet >  other,
const std::string &  newname 
)
static

Caller must delete object after use.

Parameters
otherThe vtkDataSet used to construct the mesh
fnamename of mesh file
Returns
<>

Definition at line 112 of file meshBase.C.

113  {
114  return new vtkMesh(other, newname);
115 }

◆ Create() [3/3]

meshBase * meshBase::Create ( const std::vector< double > &  xCrds,
const std::vector< double > &  yCrds,
const std::vector< double > &  zCrds,
const std::vector< nemId_t > &  elmConn,
const int  cellType,
const std::string &  newname 
)
static

Use of this is only valid when mesh has one cell type.

Parameters
xCrds<>
yCrds<>
zCrds<>
elmConn<>
cellTypeone of the vtkCellType enums. Currently, only VTK_TETRA and VTK_TRIANGLE are supported.
newnamename of mesh file
Returns
<>

Caller must delete object after use.

Definition at line 120 of file meshBase.C.

124  {
125  return new vtkMesh(xCrds, yCrds, zCrds, elmConn, cellType, newname);
126 }
VTKCellType cellType
Definition: inpGeoMesh.C:129

◆ CreateShared() [1/4]

◆ CreateShared() [2/4]

std::shared_ptr< meshBase > meshBase::CreateShared ( meshBase mesh)
static

Memory is managed by shared pointer, so do not call delete after use.

Parameters
meshthe existing meshbase
Returns
<>

(be careful with this one!)

Definition at line 154 of file meshBase.C.

154  {
155  std::shared_ptr<meshBase> sharedMesh;
156  sharedMesh.reset(mesh);
157  return sharedMesh;
158 }

◆ CreateShared() [3/4]

std::shared_ptr< meshBase > meshBase::CreateShared ( vtkSmartPointer< vtkDataSet >  other,
const std::string &  newname 
)
static

Memory is managed by shared pointer, so do not call delete after use.

Parameters
otherThe vtkDataSet used to construct the mesh
fnamename of mesh file
Returns
<>

Definition at line 162 of file meshBase.C.

References Create(), and mesh.

163  {
164  std::shared_ptr<meshBase> mesh;
165  mesh.reset(meshBase::Create(other, newname));
166  return mesh;
167 }
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
std::shared_ptr< meshBase > mesh

◆ CreateShared() [4/4]

std::shared_ptr< meshBase > meshBase::CreateShared ( const std::vector< double > &  xCrds,
const std::vector< double > &  yCrds,
const std::vector< double > &  zCrds,
const std::vector< nemId_t > &  elmConn,
int  cellType,
const std::string &  newname 
)
static

Memory is managed by shared pointer, so do not call delete after use.

Parameters
xCrds<>
yCrds<>
zCrds<>
elmConn<>
cellTypeone of the vtkCellType enums. Currently, only VTK_TETRA and VTK_TRIANGLE are supported.
fnamename of mesh file
Returns
<>

Definition at line 179 of file meshBase.C.

References Create(), and mesh.

182  {
183  std::shared_ptr<meshBase> mesh;
184  mesh.reset(meshBase::Create(xCrds, yCrds, zCrds, elmConn, cellType, newname));
185  return mesh;
186 }
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
VTKCellType cellType
Definition: inpGeoMesh.C:129
std::shared_ptr< meshBase > mesh

◆ CreateUnique() [1/4]

std::unique_ptr< meshBase > meshBase::CreateUnique ( const std::string &  fname)
static

Memory is managed by shared pointer, so do not call delete after use.

Parameters
fnamename of mesh file
Returns
<>

Definition at line 190 of file meshBase.C.

References Create().

190  {
191  return std::unique_ptr<meshBase>(meshBase::Create(fname));
192 }
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78

◆ CreateUnique() [2/4]

std::unique_ptr< meshBase > meshBase::CreateUnique ( const std::vector< double > &  xCrds,
const std::vector< double > &  yCrds,
const std::vector< double > &  zCrds,
const std::vector< nemId_t > &  elmConn,
int  cellType,
const std::string &  newname 
)
static

Memory is managed by shared pointer, so do not call delete after use.

Parameters
xCrds<>
yCrds<>
zCrds<>
elmConn<>
cellTypeone of the vtkCellType enums. Currently, only VTK_TETRA and VTK_TRIANGLE are supported.
fnamename of mesh file
Returns
<>

Definition at line 130 of file meshBase.C.

References Create().

133  {
134  return std::unique_ptr<meshBase>(
135  meshBase::Create(xCrds, yCrds, zCrds, elmConn, cellType, newname));
136 }
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
VTKCellType cellType
Definition: inpGeoMesh.C:129

◆ CreateUnique() [3/4]

std::unique_ptr< meshBase > meshBase::CreateUnique ( vtkSmartPointer< vtkDataSet >  other,
const std::string &  newname 
)
static

Memory is managed by shared pointer, so do not call delete after use.

Parameters
otherThe vtkDataSet used to construct the mesh
fnamename of mesh file
Returns
<>

Definition at line 140 of file meshBase.C.

References Create().

141  {
142  return std::unique_ptr<meshBase>(meshBase::Create(other, newname));
143 }
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78

◆ CreateUnique() [4/4]

std::unique_ptr< meshBase > meshBase::CreateUnique ( meshBase mesh)
static

Memory is managed by shared pointer, so do not call delete after use.

Parameters
meshthe existing meshbase
Returns
<>

Definition at line 147 of file meshBase.C.

References mesh.

147  {
148  return std::unique_ptr<meshBase>(mesh);
149 }
std::shared_ptr< meshBase > mesh

◆ exportExoToVtk()

meshBase * meshBase::exportExoToVtk ( const std::string &  fname)
static

exports exodusII to vtk format

Parameters
fnamename of mesh file
Returns
<>

Definition at line 918 of file meshBase.C.

References dataSet, NEM::MSH::EXOMesh::e2vEMap(), NEM::MSH::EXOMesh::elmTypeNum(), NEM::MSH::New(), numCells, numPoints, points, setFileName(), vtkMesh::write(), and NEM::MSH::EXOMesh::wrnErrMsg().

918  {
919  // opening the file
920  int CPU_word_size, IO_word_size;
921  int fid, _exErr;
922  float version;
923  CPU_word_size = sizeof(float);
924  IO_word_size = 0;
925  _exErr = 0;
926 
927  /* open EXODUS II files */
928  fid =
929  ex_open(fname.c_str(), EX_READ, &CPU_word_size, &IO_word_size, &version);
930  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem opening file " + fname + "\n");
931 
932  // declare points to be pushed into dataSet_tmp
933  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
934 
935  // declare dataSet_tmp which will be associated to output vtkMesh
936  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp =
938 
939  int numPoints;
940  int numVolCells;
941  int numElmBlk;
942  int numNdeSet;
943  int numSideSet;
944 
945  // parameter inquiry from Exodus file
946  int num_props;
947  float fdum;
948  char cdum;
949  _exErr = ex_inquire(fid, EX_INQ_API_VERS, &num_props, &fdum, &cdum);
950  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
951  std::cout << "Exodus II API version is " << fdum << std::endl;
952 
953  _exErr = ex_inquire(fid, EX_INQ_DB_VERS, &num_props, &fdum, &cdum);
954  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
955  std::cout << "Exodus II Database version is " << fdum << std::endl;
956 
957  _exErr = ex_inquire(fid, EX_INQ_DIM, &num_props, &fdum, &cdum);
958  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
959  std::cout << "Number of coordinate dimensions is " << num_props << std::endl;
960  if (num_props != 3)
961  NEM::MSH::EXOMesh::wrnErrMsg(-1, "Only 3D mesh data is supported!\n");
962 
963  _exErr = ex_inquire(fid, EX_INQ_NODES, &num_props, &fdum, &cdum);
964  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
965  numPoints = num_props;
966  std::cout << "Number of points " << numPoints << std::endl;
967 
968  _exErr = ex_inquire(fid, EX_INQ_ELEM, &num_props, &fdum, &cdum);
969  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
970  numVolCells = num_props;
971  std::cout << "Number of elements " << numVolCells << std::endl;
972 
973  _exErr = ex_inquire(fid, EX_INQ_ELEM_BLK, &num_props, &fdum, &cdum);
974  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
975  numElmBlk = num_props;
976  std::cout << "Number of element blocks " << numElmBlk << std::endl;
977 
978  _exErr = ex_inquire(fid, EX_INQ_NODE_SETS, &num_props, &fdum, &cdum);
979  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
980  numNdeSet = num_props;
981  std::cout << "Number of node sets " << numNdeSet << std::endl;
982 
983  _exErr = ex_inquire(fid, EX_INQ_SIDE_SETS, &num_props, &fdum, &cdum);
984  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading file contents.\n");
985  numSideSet = num_props;
986  std::cout << "Number of side sets " << numSideSet << std::endl;
987 
988  // nodal coordinates
989  std::vector<float> x, y, z;
990  x.resize(numPoints, 0);
991  y.resize(numPoints, 0);
992  z.resize(numPoints, 0);
993  _exErr = ex_get_coord(fid, &x[0], &y[0], &z[0]);
994  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem reading nodal coordinates.\n");
995 
996  // allocate memory for points
997  points->SetNumberOfPoints(numPoints);
998  for (int i = 0; i < numPoints; ++i) {
999  std::vector<double> pt = {x[i], y[i], z[i]};
1000  points->SetPoint(i, &pt[0]);
1001  }
1002  // inserting point array into dataSet_tmp
1003  dataSet_tmp->SetPoints(points);
1004 
1005  // Connectivities
1006  // allocating space for cell connectivities
1007  dataSet_tmp->Allocate(numVolCells);
1008 
1009  // VTK int array for block IDs
1010  vtkSmartPointer<vtkIntArray> blocks = vtkSmartPointer<vtkIntArray>::New();
1011  blocks->SetNumberOfValues(numVolCells);
1012  blocks->SetName("BlockId");
1013 
1014  int count = 0;
1015  // read element blocks
1016  for (int iEB = 1; iEB <= numElmBlk; iEB++) {
1017  int num_el_in_blk, num_nod_per_el, num_attr /*, *connect*/;
1018  // float *attrib;
1019  char elem_type[MAX_STR_LENGTH + 1];
1020  // read element block parameters
1021  _exErr = ex_get_elem_block(fid, iEB, elem_type, &num_el_in_blk,
1022  &num_nod_per_el, &num_attr);
1024  "Problem reading element block parameters.\n");
1025  // read element connectivity
1026  std::vector<int> conn;
1027  conn.resize(num_el_in_blk * num_nod_per_el, 0);
1028  _exErr = ex_get_elem_conn(fid, iEB, &conn[0]);
1030  _exErr, "Problem reading element block connectivities.\n");
1031  // read element block attributes
1032  // std::vector<float> attr;
1033  // attr.resize(0.,num_el_in_blk*num_nod_per_el);
1034  //_exErr = ex_get_elem_attr (fid, iEB, &attrib[0]);
1035  // EXOMesh::wrnErrMsg(_exErr, "Problem reading element block
1036  // attributes.\n");
1037 
1038  for (int iEl = 0; iEl < num_el_in_blk; ++iEl) {
1039  vtkSmartPointer<vtkIdList> vtkcellIds = vtkSmartPointer<vtkIdList>::New();
1040  VTKCellType vct =
1042  for (int jc = iEl * num_nod_per_el; jc < (iEl + 1) * num_nod_per_el;
1043  jc++) {
1044  // insert connectivities for cell into cellIds container
1045  vtkcellIds->InsertNextId(conn[jc] - 1);
1046  }
1047  // insert connectivities
1048  dataSet_tmp->InsertNextCell(vct, vtkcellIds);
1049 
1050  vtkIdType i = count;
1051  blocks->SetValue(i, iEB);
1052  count++;
1053  }
1054  }
1055  // std::cout << "Trimmed name = " << nemAux::trim_fname(fname, ".vtu")
1056  // << std::endl;
1057 
1058  vtkMesh *vtkmesh = new vtkMesh();
1059  vtkmesh->dataSet = dataSet_tmp;
1060  vtkmesh->numCells = vtkmesh->dataSet->GetNumberOfCells();
1061  vtkmesh->numPoints = vtkmesh->dataSet->GetNumberOfPoints();
1062 
1063  vtkmesh->dataSet->GetCellData()->AddArray(blocks);
1064 
1065  vtkmesh->setFileName("vtkWithIds.vtu");
1066  vtkmesh->write();
1067  std::cout << "vtkMesh constructed" << std::endl;
1068 
1069  // closing the file
1070  _exErr = ex_close(fid);
1071  NEM::MSH::EXOMesh::wrnErrMsg(_exErr, "Problem closing the exodusII file.");
1072 
1073  return vtkmesh;
1074 }
nemId_t numPoints
number of points in mesh
Definition: meshBase.H:714
void setFileName(const std::string &fname)
set the file name.
Definition: meshBase.H:675
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
elementType elmTypeNum(std::string tag)
Convert string to EXODUS element type.
Definition: exoMesh.C:97
nemId_t numCells
number of cells in mesh
Definition: meshBase.H:718
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133
VTKCellType e2vEMap(elementType et)
Convert EXODUS element type to VTK cell type.
Definition: exoMesh.C:50
void write() const override
write the mesh to file named after the private var &#39;filename&#39;.
Definition: vtkMesh.H:152
void wrnErrMsg(int errCode, const std::string &msg)
Logging method.
Definition: exoMesh.C:202

◆ exportGmshToVtk()

meshBase * meshBase::exportGmshToVtk ( const std::string &  fname)
static
Parameters
fnamename of mesh file
Returns
<>

Definition at line 409 of file meshBase.C.

References data, dataSet, id, NEM::MSH::New(), numCells, numPoints, points, vtkMesh::setCellDataArray(), and vtkMesh::setPointDataArray().

Referenced by NEM::DRV::GmshToVtkConversionDriver::execute(), NEM::DRV::GmshMeshGenDriver::execute(), NEM::GEO::rocPack::geomToVTK(), NEM::DRV::MeshGenDriver::MeshGenDriver(), and NEM::ADP::Refine::run().

409  {
410  std::ifstream meshStream(fname);
411  if (!meshStream.good()) {
412  std::cout << "Error opening file " << fname << std::endl;
413  exit(1);
414  }
415 
416  bool warning = true;
417 
418  std::string line;
419  int numPoints = 0, numCells = 0, numPhysGrps = 0;
420  bool fndPhyGrp = false;
421  std::vector<int> physGrpDims;
422  std::map<int, std::string> physGrpIdName;
423  std::vector<std::vector<std::vector<double>>> pointData;
424  std::vector<std::vector<std::vector<double>>> cellData;
425  std::vector<std::vector<double>> cellPhysGrpIds;
426  std::vector<std::string> pointDataNames;
427  std::vector<std::string> cellDataNames;
428  std::vector<std::string> fieldDataNames;
429  std::vector<std::string> fieldData;
430 
431  // declare points to be pushed into dataSet_tmp
432  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
433  points->SetDataTypeToDouble();
434  // declare dataSet_tmp which will be associated to output vtkMesh
435  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp =
437  // map to hold true index of points (gmsh allows non-contiguous ordering)
438  std::map<int, int> trueIndex;
439 
440  while (getline(meshStream, line)) {
441  if (line.find("$PhysicalNames") != -1) {
442  fndPhyGrp = true;
443  getline(meshStream, line);
444  std::stringstream ss(line);
445  ss >> numPhysGrps;
446  std::cout << "Found " << numPhysGrps << " physical groups!\n";
447  for (int i = 0; i < numPhysGrps; ++i) {
448  int grpDim, grpId;
449  std::string grpName;
450  getline(meshStream, line);
451  std::stringstream ss(line);
452  ss >> grpDim >> grpId >> grpName;
453  grpName.erase(std::remove(grpName.begin(), grpName.end(), '\"'),
454  grpName.end());
455  /*std::cout << "Group Dim = " << grpDim
456  << "\nGroup Id = " << grpId
457  << "\nGroup Name = " << grpName
458  << std::endl;*/
459  physGrpDims.push_back(grpDim);
460  physGrpIdName[grpId] = grpName;
461  }
462  }
463 
464  if (line.find("$Nodes") != -1) {
465  getline(meshStream, line);
466  std::stringstream ss(line);
467  ss >> numPoints;
468  int id;
469  double x, y, z;
470  // allocate memory for points
471  points->SetNumberOfPoints(numPoints);
472  for (int i = 0; i < numPoints; ++i) {
473  getline(meshStream, line);
474  std::stringstream ss(line);
475  ss >> id >> std::setprecision(16) >> x >> y >> z;
476  double point[3];
477  point[0] = x;
478  point[1] = y;
479  point[2] = z;
480  // insert point i
481  points->SetPoint(i, point);
482  trueIndex.insert(std::pair<int, int>(id, i));
483  }
484  // inserting point array into dataSet_tmp
485  dataSet_tmp->SetPoints(points);
486  }
487 
488  if (line.find("$Elements") != -1) {
489  getline(meshStream, line);
490  // std::cout << "line = " << line << std::endl;
491  std::stringstream ss(line);
492  ss >> numCells;
493  int id, type, numTags;
494  // double tmp2[1];
495  // allocate space for cell connectivities
496  dataSet_tmp->Allocate(numCells);
497  for (int i = 0; i < numCells; ++i) {
498  getline(meshStream, line);
499  std::stringstream ss(line);
500  ss >> id >> type >> numTags;
501  vtkSmartPointer<vtkIdList> vtkCellIds =
503  if (type == 2) {
504  int tmp;
505  if (!fndPhyGrp) {
506  for (int j = 0; j < numTags; ++j) ss >> tmp;
507  } else {
508  std::vector<double> physGrpId(1);
509  ss >> physGrpId[0];
510  cellPhysGrpIds.push_back(physGrpId);
511  for (int j = 0; j < numTags - 1; ++j) ss >> tmp;
512  }
513  for (int j = 0; j < 3; ++j) {
514  ss >> tmp;
515  // insert connectivities for cell into cellIds container
516  vtkCellIds->InsertNextId(trueIndex[tmp]); //-1);
517  }
518  // insert connectivities for triangle elements into dataSet
519  dataSet_tmp->InsertNextCell(VTK_TRIANGLE, vtkCellIds);
520  } else if (type == 4) {
521  int tmp;
522  if (!fndPhyGrp) {
523  for (int j = 0; j < numTags; ++j) ss >> tmp;
524  } else {
525  std::vector<double> physGrpId(1);
526  ss >> physGrpId[0];
527  cellPhysGrpIds.push_back(physGrpId);
528  for (int j = 0; j < numTags - 1; ++j) ss >> tmp;
529  }
530  for (int j = 0; j < 4; ++j) {
531  ss >> tmp;
532  // insert connectivities for cell into cellIds container
533  vtkCellIds->InsertNextId(trueIndex[tmp]); //-1);
534  }
535  // insert connectivities for tet elements into dataSet
536  dataSet_tmp->InsertNextCell(VTK_TETRA, vtkCellIds);
537  } else if (type == 3) {
538  int tmp;
539  if (!fndPhyGrp) {
540  for (int j = 0; j < numTags; ++j) ss >> tmp;
541  } else {
542  std::vector<double> physGrpId(1);
543  ss >> physGrpId[0];
544  cellPhysGrpIds.push_back(physGrpId);
545  for (int j = 0; j < numTags - 1; ++j) ss >> tmp;
546  }
547  for (int j = 0; j < 4; ++j) {
548  ss >> tmp;
549  // insert connectivities for cell into cellIds container
550  vtkCellIds->InsertNextId(trueIndex[tmp]); //-1);
551  }
552  // insert connectivities for tet elements into dataSet
553  dataSet_tmp->InsertNextCell(VTK_QUAD, vtkCellIds);
554  } else if (type == 5) {
555  int tmp;
556  if (!fndPhyGrp) {
557  for (int j = 0; j < numTags; ++j) ss >> tmp;
558  } else {
559  std::vector<double> physGrpId(1);
560  ss >> physGrpId[0];
561  cellPhysGrpIds.push_back(physGrpId);
562  for (int j = 0; j < numTags - 1; ++j) ss >> tmp;
563  }
564  for (int j = 0; j < 8; ++j) {
565  ss >> tmp;
566  // insert connectivities for cell into cellIds container
567  vtkCellIds->InsertNextId(trueIndex[tmp]); //-1);
568  }
569  // insert connectivities for hex elements into dataSet
570  dataSet_tmp->InsertNextCell(VTK_HEXAHEDRON, vtkCellIds);
571  } else if (type == 6) {
572  int tmp;
573  if (!fndPhyGrp) {
574  for (int j = 0; j < numTags; ++j) ss >> tmp;
575  } else {
576  std::vector<double> physGrpId(1);
577  ss >> physGrpId[0];
578  cellPhysGrpIds.push_back(physGrpId);
579  for (int j = 0; j < numTags - 1; ++j) ss >> tmp;
580  }
581  for (int j = 0; j < 6; ++j) {
582  ss >> tmp;
583  // insert connectivities for cell into cellIds container
584  vtkCellIds->InsertNextId(trueIndex[tmp]); //-1);
585  }
586  // insert connectivities for wedge/prism elements into dataSet
587  dataSet_tmp->InsertNextCell(VTK_WEDGE, vtkCellIds);
588  } else if (type == 11) {
589  int tmp;
590  if (!fndPhyGrp) {
591  for (int j = 0; j < numTags; ++j) ss >> tmp;
592  } else {
593  std::vector<double> physGrpId(1);
594  ss >> physGrpId[0];
595  cellPhysGrpIds.push_back(physGrpId);
596  for (int j = 0; j < numTags - 1; ++j) ss >> tmp;
597  }
598  for (int j = 0; j < 10; ++j) {
599  ss >> tmp;
600  // insert connectivities for cell into cellIds container
601  vtkCellIds->InsertNextId(trueIndex[tmp]); //-1);
602  }
603  // insert connectivities for tet elements into dataSet
604  dataSet_tmp->InsertNextCell(VTK_QUADRATIC_TETRA, vtkCellIds);
605  } else {
606  if (warning) {
607  std::cout
608  << "Warning: Only triangular, quadrilateral, "
609  "tetrahedral, hexahedral, wedge, and quadratic tetrahedral"
610  "elements are supported, "
611  << "everything else is ignored! " << std::endl;
612  warning = false;
613  // exit(1);
614  }
615  }
616  }
617  }
618 
619  if (line.find("$NodeData") != -1) {
620  std::vector<std::vector<double>> currPointData;
621  getline(meshStream, line); // moving to num_string_tags
622  std::stringstream ss(line);
623  int num_string_tags;
624  ss >> num_string_tags;
625  std::string dataname;
626  for (int i = 0; i < num_string_tags; ++i) {
627  getline(meshStream, line); // get string tag
628  if (i == 0) {
629  std::stringstream ss(line);
630  ss >> dataname;
631  }
632  }
633  dataname.erase(std::remove(dataname.begin(), dataname.end(), '\"'),
634  dataname.end());
635  pointDataNames.push_back(dataname);
636  getline(meshStream, line); // moving to num_real_tags
637  std::stringstream ss1(line);
638  int num_real_tags;
639  ss1 >> num_real_tags;
640  for (int i = 0; i < num_real_tags; ++i) getline(meshStream, line);
641 
642  getline(meshStream, line); // moving to num_int_tags
643  std::stringstream ss2(line);
644  int num_int_tags;
645  ss2 >> num_int_tags;
646  int dt, dim, numFields, tmp;
647  for (int i = 0; i < num_int_tags; ++i) {
648  getline(meshStream, line); // get int tag
649  std::stringstream ss(line);
650  if (i == 0)
651  ss >> dt;
652  else if (i == 1)
653  ss >> dim;
654  else if (i == 2)
655  ss >> numFields;
656  else
657  ss >> tmp;
658  }
659  for (int i = 0; i < numFields; ++i) {
660  std::vector<double> data(dim);
661  int id;
662  double val;
663  getline(meshStream, line);
664  std::stringstream ss(line);
665  ss >> id;
666  for (int j = 0; j < dim; ++j) {
667  ss >> val;
668  data[j] = val;
669  }
670  currPointData.push_back(data);
671  }
672  pointData.push_back(currPointData);
673  }
674 
675  if (line.find("$ElementData") != -1) {
676  std::vector<std::vector<double>> currCellData;
677  getline(meshStream, line); // moving to num_string_tags
678  std::stringstream ss(line);
679  int num_string_tags;
680  ss >> num_string_tags;
681  std::string dataname;
682  for (int i = 0; i < num_string_tags; ++i) {
683  getline(meshStream, line); // get string tag
684  if (i == 0) {
685  std::stringstream ss(line);
686  ss >> dataname;
687  }
688  }
689  dataname.erase(std::remove(dataname.begin(), dataname.end(), '\"'),
690  dataname.end());
691  cellDataNames.push_back(dataname);
692  getline(meshStream, line); // moving to num_real_tags
693  std::stringstream ss1(line);
694  int num_real_tags;
695  ss1 >> num_real_tags;
696  for (int i = 0; i < num_real_tags; ++i) getline(meshStream, line);
697 
698  getline(meshStream, line); // moving to num_int_tags
699  std::stringstream ss2(line);
700  int num_int_tags;
701  ss2 >> num_int_tags;
702  int dt, dim, numFields, tmp;
703  for (int i = 0; i < num_int_tags; ++i) {
704  getline(meshStream, line); // get int tag
705  std::stringstream ss(line);
706  if (i == 0)
707  ss >> dt;
708  else if (i == 1)
709  ss >> dim;
710  else if (i == 2)
711  ss >> numFields;
712  else
713  ss >> tmp;
714  }
715  for (int i = 0; i < numFields; ++i) {
716  std::vector<double> data(dim);
717  int id;
718  double val;
719  getline(meshStream, line);
720  std::stringstream ss(line);
721  ss >> id;
722  for (int j = 0; j < dim; ++j) {
723  ss >> val;
724  data[j] = val;
725  }
726  currCellData.push_back(data);
727  }
728  cellData.push_back(currCellData);
729  }
730  }
731 
732  if (fndPhyGrp) {
733  cellDataNames.push_back("PhysGrpId");
734  cellData.push_back(cellPhysGrpIds);
735  }
736 
737  vtkMesh *vtkmesh = new vtkMesh();
738  // vtkmesh->dataSet = dataSet_tmp->NewInstance();
739  // vtkmesh->dataSet->DeepCopy(dataSet_tmp);
740  vtkmesh->dataSet = dataSet_tmp;
741  vtkmesh->numCells = vtkmesh->dataSet->GetNumberOfCells();
742  if (numCells == 0) std::cerr << "Warning: No cells were found in the mesh!\n";
743  vtkmesh->numPoints = vtkmesh->dataSet->GetNumberOfPoints();
744  if (numPoints == 0)
745  std::cerr << "Warning: No points were found in the mesh!\n";
746 
747  if (numPoints > 0) {
748  for (int i = 0; i < pointData.size(); ++i)
749  vtkmesh->setPointDataArray(&(pointDataNames[i])[0u], pointData[i]);
750  }
751 
752  if (numCells > 0) {
753  for (int i = 0; i < cellData.size(); ++i)
754  vtkmesh->setCellDataArray(&(cellDataNames[i])[0u], cellData[i]);
755  }
756 
757  // vtkmesh->setFileName(nemAux::trim_fname(fname, ".vtu"));
758  // vtkmesh->write();
759  std::cout << "vtkMesh constructed" << std::endl;
760 
761  return vtkmesh;
762 }
void setPointDataArray(const std::string &name, const std::vector< std::vector< double >> &data) override
register data to dataSet&#39;s point data
Definition: vtkMesh.C:896
data_type data
Edge/face with sorted point ids (a, b, c, ...) is located at some index i in data[b], with data[b][i].first == [a, c] (for edges, third point id treated as -1).
nemId_t numPoints
number of points in mesh
Definition: meshBase.H:714
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
nemId_t numCells
number of cells in mesh
Definition: meshBase.H:718
void setCellDataArray(const std::string &name, const std::vector< std::vector< double >> &data) override
register data to dataSet&#39;s cell data
Definition: vtkMesh.C:1007
vtkIdType id
id in .inp file
Definition: inpGeoMesh.C:128
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133

◆ exportPntToVtk()

meshBase * meshBase::exportPntToVtk ( const std::string &  fname)
static

exports pntMesh to vtk format

Parameters
fnamename of mesh file
Returns
<>

Definition at line 848 of file meshBase.C.

References dataSet, PNTMesh::pntMesh::getElmConn(), PNTMesh::pntMesh::getElmOrder(), PNTMesh::pntMesh::getElmType(), PNTMesh::pntMesh::getNumberOfCells(), PNTMesh::pntMesh::getNumberOfPoints(), PNTMesh::pntMesh::getPointCrd(), PNTMesh::pntMesh::getVtkCellTag(), NEM::MSH::New(), numCells, numPoints, points, setFileName(), and nemAux::trim_fname().

848  {
849  PNTMesh::pntMesh *pMesh;
850  pMesh = new PNTMesh::pntMesh(fname);
851 
852  vtkMesh *vtkmesh = new vtkMesh();
853 
854  /*
855  if (!pMesh->isCompatible())
856  {
857  std::cerr << "Mesh contains unsupported element types.\n";
858  std::cerr << "Only 3-Node TRIs and 4-Node TETs are supported.\n";
859  vtkmesh->numCells = 0;
860  vtkmesh->numPoints = 0;
861  return vtkmesh;
862  }
863  */
864 
865  // declare points to be pushed into dataSet_tmp
866  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
867  // declare dataSet_tmp which will be associated to output vtkMesh
868  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp =
870 
871  int numPoints = pMesh->getNumberOfPoints();
872  int numVolCells = pMesh->getNumberOfCells();
873 
874  // allocate memory for points
875  points->SetNumberOfPoints(numPoints);
876  for (int i = 0; i < numPoints; ++i) {
877  std::vector<double> point;
878  point = pMesh->getPointCrd(i);
879  // insert point i
880  points->SetPoint(i, &point[0]);
881  }
882  // inserting point array into dataSet_tmp
883  dataSet_tmp->SetPoints(points);
884 
885  // allocating space for cell connectivities
886  dataSet_tmp->Allocate(numVolCells);
887  for (int i = 0; i < numVolCells; ++i) {
888  std::cout << "Element " << i << std::endl;
889  vtkSmartPointer<vtkIdList> vtkcellIds = vtkSmartPointer<vtkIdList>::New();
890  std::vector<int> conn;
891  VTKCellType vct =
892  pMesh->getVtkCellTag(pMesh->getElmType(i), pMesh->getElmOrder(i));
893  conn = pMesh->getElmConn(i, vct);
894  for (int j = 0; j < conn.size(); ++j) {
895  // insert connectivies for cell into cellIds container
896  vtkcellIds->InsertNextId(conn[j]);
897  }
898  // insert connectivies
899  dataSet_tmp->InsertNextCell(vct, vtkcellIds);
900  }
901 
902  vtkmesh->dataSet = dataSet_tmp;
903  vtkmesh->numCells = vtkmesh->dataSet->GetNumberOfCells();
904  vtkmesh->numPoints = vtkmesh->dataSet->GetNumberOfPoints();
905 
906  //
907  std::cout << "Trimmed name = " << nemAux::trim_fname(fname, ".vtu")
908  << std::endl;
909  vtkmesh->setFileName(nemAux::trim_fname(fname, ".vtu"));
910  // vtkmesh->write();
911  std::cout << "vtkMesh constructed" << std::endl;
912 
913  return vtkmesh;
914 }
nemId_t numPoints
number of points in mesh
Definition: meshBase.H:714
void setFileName(const std::string &fname)
set the file name.
Definition: meshBase.H:675
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
VTKCellType getVtkCellTag(elementType et, int order) const
Definition: pntMesh.C:558
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
int getElmOrder(int id) const
Definition: pntMesh.C:547
elementType getElmType(int id) const
Definition: pntMesh.C:537
nemId_t numCells
number of cells in mesh
Definition: meshBase.H:718
std::vector< double > getPointCrd(int id) const
Definition: pntMesh.H:113
std::string trim_fname(const std::string &name, const std::string &ext)
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133
std::vector< int > getElmConn(int id) const
Definition: pntMesh.H:114
int getNumberOfPoints() const
Definition: pntMesh.H:109
int getNumberOfCells() const
Definition: pntMesh.H:110

◆ exportVolToVtk()

meshBase * meshBase::exportVolToVtk ( const std::string &  fname)
static
Parameters
fnamename of mesh file
Returns
<>

Definition at line 766 of file meshBase.C.

References dataSet, NEM::MSH::New(), numCells, numPoints, points, setFileName(), and nemAux::trim_fname().

Referenced by NEM::DRV::NetgenMeshGenDriver::execute(), and NEM::DRV::MeshGenDriver::MeshGenDriver().

766  {
767 #ifdef HAVE_NGEN
768  nglib::Ng_Mesh *Ngmesh;
769  nglib::Ng_Init();
770  Ngmesh = nglib::Ng_NewMesh();
771 
772  int status = nglib::Ng_MergeMesh(Ngmesh, &fname[0u]);
773  if (status) {
774  std::cout << "Error: NetGen Returned: " << status << std::endl;
775  std::cout << "Could not load " << fname << " into netgen" << std::endl;
776  exit(1);
777  }
778 
779  // declare points to be pushed into dataSet_tmp
780  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
781  // declare dataSet_tmp which will be associated to output vtkMesh
782  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp =
784  int numNgPoints = nglib::Ng_GetNP(Ngmesh);
785  int numSurfCells = nglib::Ng_GetNSE(Ngmesh);
786  int numVolCells = nglib::Ng_GetNE(Ngmesh);
787 
788  // allocate memory for points
789  points->SetNumberOfPoints(numNgPoints);
790  for (int i = 1; i <= numNgPoints; ++i) {
791  double point[3];
792  nglib::Ng_GetPoint(Ngmesh, i, point);
793  // insert point i
794  points->SetPoint(i - 1, point);
795  }
796  // inserting point array into dataSet_tmp
797  dataSet_tmp->SetPoints(points);
798 
799  // allocating space for cell connectivities
800  dataSet_tmp->Allocate(numSurfCells + numVolCells);
801 
802  for (int i = 1; i <= numSurfCells; ++i) {
803  vtkSmartPointer<vtkIdList> vtkcellIds = vtkSmartPointer<vtkIdList>::New();
804  int tri[3];
805  nglib::Ng_GetSurfaceElement(Ngmesh, i, tri);
806  for (int j = 0; j < 3; ++j) {
807  // insert connectivies for cell into cellIds container
808  vtkcellIds->InsertNextId(tri[j] - 1);
809  }
810  // insert connectivies for triangle elements into dataSet
811  dataSet_tmp->InsertNextCell(VTK_TRIANGLE, vtkcellIds);
812  }
813 
814  for (int i = 1; i <= numVolCells; ++i) {
815  vtkSmartPointer<vtkIdList> vtkcellIds = vtkSmartPointer<vtkIdList>::New();
816  int tet[4];
817  nglib::Ng_GetVolumeElement(Ngmesh, i, tet);
818  for (int j = 0; j < 4; ++j) {
819  // insert connectivies for cell into cellIds container
820  vtkcellIds->InsertNextId(tet[j] - 1);
821  }
822  // insert connectivies for triangle elements into dataSet
823  dataSet_tmp->InsertNextCell(VTK_TETRA, vtkcellIds);
824  }
825 
826  vtkMesh *vtkmesh = new vtkMesh();
827  // vtkmesh->dataSet = dataSet_tmp->NewInstance();//
828  // vtkmesh->dataSet->DeepCopy(dataSet_tmp);//vtkDataSet::SafeDownCast(dataSet_tmp));
829  vtkmesh->dataSet = dataSet_tmp;
830  vtkmesh->numCells = vtkmesh->dataSet->GetNumberOfCells();
831  vtkmesh->numPoints = vtkmesh->dataSet->GetNumberOfPoints();
832 
833  vtkmesh->setFileName(nemAux::trim_fname(fname, ".vtu"));
834  std::cout << "vtkMesh constructed" << std::endl;
835 
836  if (Ngmesh) nglib::Ng_DeleteMesh(Ngmesh);
837  nglib::Ng_Exit();
838  return vtkmesh;
839 #else
840  std::cerr << "ENABLE_NETGEN is not used during build."
841  << " Build with ENABLE_NETGEN to use this function." << std::endl;
842  exit(1);
843 #endif
844 }
nemId_t numPoints
number of points in mesh
Definition: meshBase.H:714
void setFileName(const std::string &fname)
set the file name.
Definition: meshBase.H:675
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
nemId_t numCells
number of cells in mesh
Definition: meshBase.H:718
std::string trim_fname(const std::string &name, const std::string &ext)
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133

◆ extractSelectedCells() [1/2]

meshBase * meshBase::extractSelectedCells ( meshBase mesh,
const std::vector< nemId_t > &  cellIds 
)
static
Parameters
meshThe meshBase object to extract the subset from.
cellIds<>
Returns
meshBase object representing the subset

Definition at line 226 of file meshBase.C.

References getDataSet(), and NEM::MSH::New().

227  {
228  vtkSmartPointer<vtkIdTypeArray> selectionIds =
230  selectionIds->SetNumberOfComponents(1);
231  for (const auto &cellId : cellIds) {
232  selectionIds->InsertNextValue(cellId);
233  }
234  return extractSelectedCells(mesh->getDataSet(), selectionIds);
235 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308
static meshBase * extractSelectedCells(meshBase *mesh, const std::vector< nemId_t > &cellIds)
extract subset of mesh given list of cell ids and return meshBase obj
Definition: meshBase.C:226

◆ extractSelectedCells() [2/2]

meshBase * meshBase::extractSelectedCells ( vtkSmartPointer< vtkDataSet >  mesh,
vtkSmartPointer< vtkIdTypeArray >  cellIds 
)
static
Parameters
meshThe meshBase object to extract the subset from.
cellIds<>
Returns
meshBase object representing the subset

Definition at line 239 of file meshBase.C.

References Create(), and NEM::MSH::New().

240  {
241  vtkSmartPointer<vtkSelectionNode> selectionNode =
243  selectionNode->SetFieldType(vtkSelectionNode::CELL);
244  selectionNode->SetContentType(vtkSelectionNode::INDICES);
245  selectionNode->SetSelectionList(cellIds);
246  // create the selection
247  vtkSmartPointer<vtkSelection> selection =
249  selection->AddNode(selectionNode);
250  // creating extractor
251  vtkSmartPointer<vtkExtractSelection> extractSelection =
253  // set stitch surf as input on first port
254  extractSelection->SetInputData(0, mesh);
255  // set selectionNode as input on second port
256  extractSelection->SetInputData(1, selection);
257  extractSelection->Update();
258  vtkSmartPointer<vtkDataSet> extractedCellMesh =
259  vtkDataSet::SafeDownCast(extractSelection->GetOutput());
260  meshBase *selectedCellMesh =
261  meshBase::Create(extractedCellMesh, "extracted.vtu");
262  return selectedCellMesh;
263 }
A brief description of meshBase.
Definition: meshBase.H:64
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
std::shared_ptr< meshBase > mesh

◆ extractSurface()

virtual vtkSmartPointer<vtkDataSet> meshBase::extractSurface ( )
pure virtual
Returns
the surface mesh for this mesh.

Implemented in FOAM::foamMesh, vtkMesh, meshSrch, and gmshMesh.

Referenced by vtkMesh::read().

◆ generateSizeField()

void meshBase::generateSizeField ( const std::string &  method,
int  arrayID,
double  dev_mlt,
bool  maxIsmin,
double  sizeFactor = 1.0,
int  order = 1 
)
Parameters
method(e.g., "gradient", "value", "error estimator")
arrayID<>
dev_mltused to determine which cells to consider for refinement
maxIsminused to determine which cells to consider for refinement
sizeFactor<>
order<>

Definition at line 396 of file meshBase.C.

References NEM::ADP::SizeFieldBase::computeSizeField(), NEM::ADP::SizeFieldBase::CreateUnique(), and NEM::ADP::SizeFieldBase::setSizeFactor().

Referenced by NEM::ADP::Refine::Refine().

398  {
399  std::cout << "Size Factor = " << sizeFactor << std::endl;
400  std::unique_ptr<NEM::ADP::SizeFieldBase> sfobj =
401  NEM::ADP::SizeFieldBase::CreateUnique(dataSet, method, arrayID, dev_mult,
402  maxIsmin, sizeFactor, order);
403  sfobj->setSizeFactor(sizeFactor);
404  sfobj->computeSizeField(dataSet->GetPointData()->GetArray(arrayID));
405 }
virtual void computeSizeField(vtkDataArray *da)=0
static std::unique_ptr< SizeFieldBase > CreateUnique(vtkDataSet *_dataSet, const std::string &method, int arrayID, double _dev_mult, bool _maxIsmin, double _sizeFactor=1.0, int _order=1)
Definition: SizeFieldGen.C:103
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
void setSizeFactor(double sf)
Definition: SizeFieldGen.H:73

◆ getArrayIDs()

std::vector< int > meshBase::getArrayIDs ( std::vector< std::string >  arrayNames,
bool  fromPointArrays = false 
)

Definition at line 1835 of file meshBase.C.

References id.

1836  {
1837  std::vector<int> arrayIDs(arrayNames.size());
1838  for (int i = 0; i < arrayNames.size(); ++i) {
1839  int id = IsArrayName(arrayNames[i], fromPointArrays);
1840  if (id == -1) {
1841  std::cerr << "Array " << arrayNames[i]
1842  << " not found in set of data arrays." << std::endl;
1843  exit(1);
1844  }
1845  arrayIDs[i] = id;
1846  }
1847  return arrayIDs;
1848 }
int IsArrayName(const std::string &name, bool pointOrCell=false) const
check for named array in vtk and return its integer id.
Definition: meshBase.C:267
vtkIdType id
id in .inp file
Definition: inpGeoMesh.C:128

◆ getCell()

virtual std::map<nemId_t, std::vector<double> > meshBase::getCell ( nemId_t  id) const
pure virtual
Parameters
idThe id of the cell.
Returns
point indices and respective coordinates

Implemented in vtkMesh, FOAM::foamMesh, meshSrch, and gmshMesh.

Referenced by vtkMesh::read().

◆ getCellCenter()

virtual std::vector<double> meshBase::getCellCenter ( nemId_t  cellID) const
pure virtual
Parameters
cellID<>
Returns
<>

Implemented in vtkMesh, FOAM::foamMesh, gmshMesh, and meshSrch.

Referenced by meshSrch::CreateUnique(), vtkMesh::read(), and FETransfer::transferCellData().

◆ getCellDataArray() [1/2]

virtual void meshBase::getCellDataArray ( const std::string &  name,
std::vector< double > &  data 
)
inlinevirtual
Parameters
name<>
data<>

Reimplemented in vtkMesh.

Definition at line 369 of file meshBase.H.

Referenced by NEM::GEO::rocPack::createCohesiveElements(), NEM::DRV::ConversionDriver::genExo(), cgnsAnalyzer::overwriteSolData(), and vtkMesh::write().

370  {}

◆ getCellDataArray() [2/2]

virtual void meshBase::getCellDataArray ( int  arrayId,
std::vector< double > &  data 
)
inlinevirtual
Parameters
arrayId<>
data<>

Reimplemented in vtkMesh.

Definition at line 376 of file meshBase.H.

376 {}

◆ getCellDataIdx()

virtual int meshBase::getCellDataIdx ( const std::string &  name)
inlinevirtual
Parameters
name<>
Returns
<>

Reimplemented in vtkMesh, and meshSrch.

Definition at line 363 of file meshBase.H.

Referenced by vtkMesh::write().

363 { return 0; }

◆ getCellLengths()

virtual std::vector<double> meshBase::getCellLengths ( ) const
pure virtual
Returns
<>

Implemented in vtkMesh, FOAM::foamMesh, meshSrch, and gmshMesh.

Referenced by vtkMesh::read().

◆ getCellType()

virtual int meshBase::getCellType ( ) const
pure virtual
Returns
<>

Implemented in vtkMesh, FOAM::foamMesh, gmshMesh, and meshSrch.

Referenced by vtkMesh::read().

◆ getCellVec()

virtual std::vector<std::vector<double> > meshBase::getCellVec ( nemId_t  id) const
pure virtual
Parameters
idThe id of the cell.
Returns
vector of coords of cell

Implemented in vtkMesh, FOAM::foamMesh, gmshMesh, and meshSrch.

Referenced by meshSrch::CreateUnique(), diffMesh(), and vtkMesh::read().

◆ getConnectivities()

virtual std::vector<nemId_t> meshBase::getConnectivities ( ) const
pure virtual

This is only safe to use if mesh has cells of the same type or you have information on the number of cells of each type and the order in which they appear (for look up in resulting vector)

Returns
<>

Implemented in FOAM::foamMesh, vtkMesh, gmshMesh, and meshSrch.

Referenced by meshPartitioner::meshPartitioner(), and vtkMesh::read().

◆ getDataSet()

◆ getFileName()

const std::string& meshBase::getFileName ( ) const
inline
Returns
The current value of the private variable "filename"

Definition at line 680 of file meshBase.H.

Referenced by MeshQuality::checkMesh(), and partition().

680 { return filename; }
std::string filename
name of mesh file
Definition: meshBase.H:726

◆ getGlobToPartCellMap()

std::map<nemId_t, nemId_t> meshBase::getGlobToPartCellMap ( )
inline
Note
These are only generated if mesh is one of the partitions returned from a call to meshBase::partition

Definition at line 571 of file meshBase.H.

References globToPartCellMap.

571  {
572  return globToPartCellMap;
573  }
std::map< nemId_t, nemId_t > globToPartCellMap
map between global and local cell idx in partition
Definition: meshBase.H:750

◆ getGlobToPartNodeMap()

std::map<nemId_t, nemId_t> meshBase::getGlobToPartNodeMap ( )
inline
Note
These are only generated if mesh is one of the partitions returned from a call to meshBase::partition

Definition at line 563 of file meshBase.H.

References globToPartNodeMap.

563  {
564  return globToPartNodeMap;
565  }
std::map< nemId_t, nemId_t > globToPartNodeMap
map between global and local node idx in partition for distributed data sets
Definition: meshBase.H:746

◆ getMetadata()

vtkSmartPointer<vtkModelMetadata> meshBase::getMetadata ( )
inline

Definition at line 443 of file meshBase.H.

Referenced by NEM::DRV::ConversionDriver::genExo().

443 { return metadata; }
vtkSmartPointer< vtkModelMetadata > metadata
Definition: meshBase.H:760

◆ getNewArrayNames()

std::vector<std::string> meshBase::getNewArrayNames ( )
inline

Definition at line 698 of file meshBase.H.

Referenced by OrderOfAccuracy::computeRichardsonExtrapolation(), NEM::ADP::Refine::run(), and vtkMesh::vtkMesh().

698 { return newArrayNames; }
std::vector< std::string > newArrayNames
new names to set for transferred data
Definition: meshBase.H:739

◆ getNumberOfCells()

◆ getNumberOfPoints()

◆ getPartToGlobCellMap()

std::map<nemId_t, nemId_t> meshBase::getPartToGlobCellMap ( )
inline
Note
These are only generated if mesh is one of the partitions returned from a call to meshBase::partition

Definition at line 587 of file meshBase.H.

References partToGlobCellMap.

587  {
588  return partToGlobCellMap;
589  }
std::map< nemId_t, nemId_t > partToGlobCellMap
map between local and global cell idx in partition
Definition: meshBase.H:758

◆ getPartToGlobNodeMap()

std::map<nemId_t, nemId_t> meshBase::getPartToGlobNodeMap ( )
inline
Note
These are only generated if mesh is one of the partitions returned from a call to meshBase::partition

Definition at line 579 of file meshBase.H.

References partToGlobNodeMap.

579  {
580  return partToGlobNodeMap;
581  }
std::map< nemId_t, nemId_t > partToGlobNodeMap
map between local and global node idx in partition
Definition: meshBase.H:754

◆ getPoint()

virtual std::vector<double> meshBase::getPoint ( nemId_t  id) const
pure virtual

◆ getPointDataArray() [1/2]

virtual void meshBase::getPointDataArray ( const std::string &  name,
std::vector< double > &  data 
)
inlinevirtual

assumes data is not allocated prior to calling

Parameters
name<>
data<>

Reimplemented in vtkMesh.

Definition at line 349 of file meshBase.H.

Referenced by cgnsAnalyzer::overwriteSolData(), and vtkMesh::write().

350  {}

◆ getPointDataArray() [2/2]

virtual void meshBase::getPointDataArray ( int  arrayId,
std::vector< double > &  data 
)
inlinevirtual

assumes data is not allocated prior to calling

Parameters
arrayId<>
data<>

Reimplemented in vtkMesh.

Definition at line 357 of file meshBase.H.

357 {}

◆ getVertCrds()

virtual std::vector<std::vector<double> > meshBase::getVertCrds ( ) const
pure virtual
Returns
<>

Implemented in vtkMesh, FOAM::foamMesh, meshSrch, and gmshMesh.

Referenced by partition(), and vtkMesh::read().

◆ inspectEdges()

virtual void meshBase::inspectEdges ( const std::string &  ofname) const
pure virtual
Parameters
ofname<>

Implemented in vtkMesh, FOAM::foamMesh, meshSrch, and gmshMesh.

Referenced by vtkMesh::read().

◆ integrateOverMesh()

std::vector< std::vector< double > > meshBase::integrateOverMesh ( const std::vector< int > &  arrayIDs)

integrated data is available per cell as vtk cell data after operation

Parameters
arrayIDs<>
Returns
total integral for each datum is returned

Definition at line 387 of file meshBase.C.

References GaussCubature::CreateUnique(), and GaussCubature::integrateOverAllCells().

Referenced by OrderOfAccuracy::computeDiff(), OrderOfAccuracy::computeGCI_21(), and OrderOfAccuracy::computeGCI_32().

388  {
389  std::unique_ptr<GaussCubature> cubature =
391  return cubature->integrateOverAllCells();
392 }
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
std::vector< std::vector< double > > integrateOverAllCells()
Integrate arrays specified by arrayIDs over all cells.
Definition: Cubature.C:697
static std::unique_ptr< GaussCubature > CreateUnique(vtkDataSet *_dataSet)
Definition: Cubature.C:135

◆ IsArrayName()

int meshBase::IsArrayName ( const std::string &  name,
bool  pointOrCell = false 
) const

check for named array in vtk

Parameters
pointOrCellboolean that tells the method whether to transfer point (False) or cell (True) data.
Returns
<>

Definition at line 267 of file meshBase.C.

Referenced by writeCobalt().

268  {
269  if (!pointOrCell) {
270  vtkPointData *pd = dataSet->GetPointData();
271  if (pd->GetNumberOfArrays()) {
272  for (int i = 0; i < pd->GetNumberOfArrays(); ++i) {
273  if (name == (pd->GetArrayName(i) ? pd->GetArrayName(i) : "NULL")) {
274  return i;
275  }
276  }
277  }
278  } else {
279  vtkCellData *cd = dataSet->GetCellData();
280  if (cd->GetNumberOfArrays()) {
281  for (int i = 0; i < cd->GetNumberOfArrays(); ++i) {
282  if (name == (cd->GetArrayName(i) ? cd->GetArrayName(i) : "Null")) {
283  return i;
284  }
285  }
286  }
287  }
288  return -1;
289 }
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722

◆ partition()

std::vector< std::shared_ptr< meshBase > > meshBase::partition ( const meshBase mbObj,
int  numPartitions 
)
static

partition mesh into numPartition pieces (static fcn)

Parameters
mbObjThe meshBase object to partition.
numPartitionsThe number of partitions to partition the mesh into
Returns
<>

Memory is managed by shared pointer, so do not call delete after use.

Definition at line 295 of file meshBase.C.

References CreateShared(), getFileName(), getVertCrds(), globalNodeIds, NEM::MSH::New(), partToGlobCellMap, partToGlobNodeMap, and nemAux::trim_fname().

296  {
297  // construct partitioner with meshBase object
298  auto *mPart = new meshPartitioner(mbObj);
299  if (mPart->partition(numPartitions)) {
300  exit(1);
301  }
302  // initialize vector of meshBase partitions
303  std::vector<std::shared_ptr<meshBase>> mbParts(numPartitions);
304  for (int i = 0; i < numPartitions; ++i) {
305  // define coordinates
306  std::vector<std::vector<double>> comp_crds(mbObj->getVertCrds());
307  // get partition connectivity and zero index it
308  std::vector<int> vtkConn_int(mPart->getConns(i));
309  std::vector<nemId_t> vtkConn(vtkConn_int.begin(), vtkConn_int.end());
310  for (auto &&it : vtkConn) it -= 1;
311 
312  std::string basename(nemAux::trim_fname(mbObj->getFileName(), ""));
313  basename += std::to_string(i);
314  basename += ".vtu";
315  // construct meshBase partition from coordinates and connectivities
316  // from partitioner
317  mbParts[i] = meshBase::CreateShared(
318  mPart->getCrds(i, comp_crds[0]), mPart->getCrds(i, comp_crds[1]),
319  mPart->getCrds(i, comp_crds[2]), vtkConn, VTK_TETRA, basename);
320  // add partition id to node and cell data of mbPart
321  vtkSmartPointer<vtkIntArray> nodePartitionIds =
323  nodePartitionIds->SetName("NodePartitionIds");
324  nodePartitionIds->SetNumberOfComponents(1);
325  nodePartitionIds->SetNumberOfTuples(mbParts[i]->getNumberOfPoints());
326  nodePartitionIds->FillComponent(0, i);
327  mbParts[i]->getDataSet()->GetPointData()->AddArray(nodePartitionIds);
328 
329  vtkSmartPointer<vtkIntArray> cellPartitionIds =
331  cellPartitionIds->SetName("CellPartitionIds");
332  cellPartitionIds->SetNumberOfComponents(1);
333  cellPartitionIds->SetNumberOfTuples(mbParts[i]->getNumberOfCells());
334  cellPartitionIds->FillComponent(0, i);
335  mbParts[i]->getDataSet()->GetCellData()->AddArray(cellPartitionIds);
336 
337  // add global node index array to partition
338  std::map<int, int> partToGlobNodeMap(mPart->getPartToGlobNodeMap(i));
339  vtkSmartPointer<vtkIdTypeArray> globalNodeIds =
341  globalNodeIds->SetName("GlobalNodeIds");
342  globalNodeIds->SetNumberOfComponents(1);
343  globalNodeIds->SetNumberOfTuples(mbParts[i]->getNumberOfPoints());
344  globalNodeIds->SetNumberOfValues(mbParts[i]->getNumberOfPoints());
345  auto it = partToGlobNodeMap.begin();
346  int idx = 0;
347  while (it != partToGlobNodeMap.end()) {
348  int globidx = it->second - 1;
349  int locidx = it->first - 1;
350  globalNodeIds->SetTuple1(idx, globidx);
351  mbParts[i]->globToPartNodeMap[globidx] = locidx;
352  mbParts[i]->partToGlobNodeMap[locidx] = globidx;
353  ++idx;
354  ++it;
355  }
356  mbParts[i]->getDataSet()->GetPointData()->AddArray(globalNodeIds);
357  // add global cell index array to partition
358  std::map<int, int> partToGlobCellMap(mPart->getPartToGlobElmMap(i));
359  // vtkSmartPointer<vtkIdTypeArray> globalCellIds =
360  // vtkSmartPointer<vtkIdTypeArray>::New();
361  // globalCellIds->SetName("GlobalCellIds");
362  // globalCellIds->SetNumberOfComponents(1);
363  // globalCellIds->SetNumberOfTuples(mbPart->getNumberOfCells());
364  // globalCellIds->SetNumberOfValues(mbPart->getNumberOfCells());
365  it = partToGlobCellMap.begin();
366  idx = 0;
367  while (it != partToGlobCellMap.end()) {
368  int globidx = it->second - 1;
369  int locidx = it->first - 1;
370  // globalCellIds->SetTuple1(idx,globidx);
371  mbParts[i]->globToPartCellMap[globidx] = locidx;
372  mbParts[i]->partToGlobCellMap[locidx] = globidx;
373  ++idx;
374  ++it;
375  }
376  // mbPart->getDataSet()->GetCellData()->AddArray(globalCellIds);
377  mbParts[i]->write();
378  // mbParts[i] = mbPart;
379  }
380  delete mPart;
381  mPart = nullptr;
382  return mbParts;
383 }
std::vector< std::vector< nemId_t > > globalNodeIds
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
std::string trim_fname(const std::string &name, const std::string &ext)
const std::string & getFileName() const
get the current file name
Definition: meshBase.H:680
nemId_t getNumberOfPoints() const
return the number of points
Definition: meshBase.H:545
std::map< nemId_t, nemId_t > partToGlobNodeMap
map between local and global node idx in partition
Definition: meshBase.H:754
static std::shared_ptr< meshBase > CreateShared(const std::string &fname)
Create shared ptr from fname.
Definition: meshBase.C:171
std::map< nemId_t, nemId_t > partToGlobCellMap
map between local and global cell idx in partition
Definition: meshBase.H:758
virtual std::vector< std::vector< double > > getVertCrds() const =0
get 3 vecs with x,y and z coords
nemId_t getNumberOfCells() const
return the number of cells
Definition: meshBase.H:550

◆ read()

virtual void meshBase::read ( const std::string &  fname)
pure virtual

◆ refineMesh() [1/5]

void meshBase::refineMesh ( const std::string &  method,
int  arrayID,
double  dev_mult,
bool  maxIsmin,
double  edge_scale,
const std::string &  ofname,
bool  transferData,
double  sizeFactor = 1.,
bool  constrainBoundary = false 
)

edge_scale is for uniform refinement and is ignored in calls where method is "gradient", "value", etc.

Parameters
method<>
arrayID<>
dev_mult<>
maxIsmin<>
edge_scale<>
ofname<>
transferData<>
sizeFactor<>

instead of "uniform"

Definition at line 1565 of file meshBase.C.

References NEM::ADP::Refine::run().

Referenced by OrderOfAccuracy::computeMeshWithResolution(), and proteusHdf5::proteusHdf5().

1568  {
1569 #ifdef HAVE_GMSH
1570  std::unique_ptr<NEM::ADP::Refine> refineobj =
1571  std::unique_ptr<NEM::ADP::Refine>(
1572  new NEM::ADP::Refine(this, method, arrayID, dev_mult, maxIsmin,
1573  edge_scale, ofname, sizeFactor));
1574  refineobj->run(transferData, constrainBoundary);
1575 #else
1576  std::cerr << "Cannot use meshBase::refineMesh without Gmsh. Please configure "
1577  "with ENABLE_GMSH=ON.\n";
1578 #endif
1579 }
void run(bool transferData, bool bndryConstraint=false)
Definition: Refine.C:161

◆ refineMesh() [2/5]

void meshBase::refineMesh ( const std::string &  method,
const std::string &  arrayName,
double  dev_mult,
bool  maxIsmin,
double  edge_scale,
const std::string &  ofname,
bool  transferData,
double  sizeFactor = 1. 
)

edge_scale is for uniform refinement and is ignored in calls where method is "gradient", "value", etc.

Parameters
method<>
arrayName<>
dev_mult<>
maxIsmin<>
edge_scale<>
ofname<>
transferData<>
sizeFactor<>

instead of "uniform"

Definition at line 1599 of file meshBase.C.

1603  {
1604  int arrayID = IsArrayName(arrayName);
1605  if (arrayID == -1) {
1606  std::cout << "Array " << arrayName
1607  << " not found in set of point data arrays" << std::endl;
1608  exit(1);
1609  }
1610  refineMesh(method, arrayID, dev_mult, maxIsmin, edge_scale, ofname,
1611  transferData, sizeFactor);
1612 }
void refineMesh(const std::string &method, int arrayID, double dev_mult, bool maxIsmin, double edge_scale, const std::string &ofname, bool transferData, double sizeFactor=1., bool constrainBoundary=false)
perform sizefield-based h-refinement.
Definition: meshBase.C:1565
int IsArrayName(const std::string &name, bool pointOrCell=false) const
check for named array in vtk and return its integer id.
Definition: meshBase.C:267

◆ refineMesh() [3/5]

void meshBase::refineMesh ( const std::string &  method,
double  edge_scale,
const std::string &  ofname,
bool  transferData,
bool  constrainBoundary = false 
)
Parameters
method<>
edge_scale<>
ofname<>
transferData<>

Definition at line 1631 of file meshBase.C.

1633  {
1634  refineMesh(method, 0, 0, false, edge_scale, ofname, transferData, 1.0,
1635  constrainBoundary);
1636 }
void refineMesh(const std::string &method, int arrayID, double dev_mult, bool maxIsmin, double edge_scale, const std::string &ofname, bool transferData, double sizeFactor=1., bool constrainBoundary=false)
perform sizefield-based h-refinement.
Definition: meshBase.C:1565

◆ refineMesh() [4/5]

void meshBase::refineMesh ( const std::string &  method,
int  arrayID,
int  order,
const std::string &  ofname,
bool  transferData 
)
Parameters
method<>
arrayID<>
order<>
ofname<>
transferData<>

Definition at line 1583 of file meshBase.C.

References NEM::ADP::Refine::run().

1584  {
1585 #ifdef HAVE_GMSH
1586  std::unique_ptr<NEM::ADP::Refine> refineobj =
1587  std::unique_ptr<NEM::ADP::Refine>(new NEM::ADP::Refine(
1588  this, method, arrayID, 0.0, false, 0.0, ofname, 1.0, _order));
1589  refineobj->run(transferData);
1590 #else
1591  std::cerr << "Cannot use meshBase::refineMesh without Gmsh. Please configure "
1592  "with ENABLE_GMSH=ON.\n";
1593 #endif
1594 }
void run(bool transferData, bool bndryConstraint=false)
Definition: Refine.C:161

◆ refineMesh() [5/5]

void meshBase::refineMesh ( const std::string &  method,
const std::string &  arrayName,
int  order,
const std::string &  ofname,
bool  transferData 
)
Parameters
method<>
arrayName<>
order<>
ofname<>
transferData<>

Definition at line 1616 of file meshBase.C.

1618  {
1619  int arrayID = IsArrayName(arrayName);
1620  if (arrayID == -1) {
1621  std::cout << "Array " << arrayName
1622  << " not found in set of point data arrays" << std::endl;
1623  exit(1);
1624  }
1625 
1626  refineMesh(method, arrayID, order, ofname, transferData);
1627 }
void refineMesh(const std::string &method, int arrayID, double dev_mult, bool maxIsmin, double edge_scale, const std::string &ofname, bool transferData, double sizeFactor=1., bool constrainBoundary=false)
perform sizefield-based h-refinement.
Definition: meshBase.C:1565
int IsArrayName(const std::string &name, bool pointOrCell=false) const
check for named array in vtk and return its integer id.
Definition: meshBase.C:267

◆ report()

virtual void meshBase::report ( ) const
inlinevirtual

◆ setCellDataArray() [1/2]

virtual void meshBase::setCellDataArray ( const std::string &  name,
const std::vector< std::vector< double >> &  data 
)
inlinevirtual
Parameters
name<>
data<>

Reimplemented in vtkMesh.

Definition at line 334 of file meshBase.H.

Referenced by hdf5Reader::setFields(), and vtkMesh::write().

335  {}

◆ setCellDataArray() [2/2]

virtual void meshBase::setCellDataArray ( const std::string &  name,
const std::vector< double > &  data 
)
inlinevirtual
Parameters
name<>
data<>

Reimplemented in vtkMesh.

Definition at line 341 of file meshBase.H.

342  {}

◆ setCheckQuality()

void meshBase::setCheckQuality ( bool  x)
inline
Parameters
x<>

Definition at line 685 of file meshBase.H.

685 { checkQuality = x; }
bool checkQuality
check transfer quality when on
Definition: meshBase.H:730

◆ setContBool()

void meshBase::setContBool ( bool  x)
inline
Parameters
x<>

Definition at line 691 of file meshBase.H.

691 { continuous = x; }
bool continuous
switch on / off weighted averaging for cell data transfer (default is off)
Definition: meshBase.H:735

◆ setFileName()

void meshBase::setFileName ( const std::string &  fname)
inline

This will allow vtk to dispatch appropriate writers based on the extension and whether it is supported by vtk.

Parameters
fnameThe name to set the private variable "filename" to

Definition at line 675 of file meshBase.H.

Referenced by exportExoToVtk(), exportPntToVtk(), exportVolToVtk(), NEM::ADP::Refine::run(), and vtkMesh::vtkMesh().

675 { filename = fname; }
std::string filename
name of mesh file
Definition: meshBase.H:726

◆ setMetadata()

void meshBase::setMetadata ( vtkSmartPointer< vtkModelMetadata >  _metadata)
inline

Definition at line 440 of file meshBase.H.

Referenced by proteusHdf5::proteusHdf5().

440  {
441  metadata = _metadata;
442  }
vtkSmartPointer< vtkModelMetadata > metadata
Definition: meshBase.H:760

◆ setPointDataArray() [1/2]

virtual void meshBase::setPointDataArray ( const std::string &  name,
const std::vector< std::vector< double >> &  data 
)
inlinevirtual
Parameters
name<>
data<>

Reimplemented in vtkMesh.

Definition at line 319 of file meshBase.H.

Referenced by hdf5Reader::setFields(), and vtkMesh::write().

320  {
321  }

◆ setPointDataArray() [2/2]

virtual void meshBase::setPointDataArray ( const std::string &  name,
const std::vector< double > &  data 
)
inlinevirtual
Parameters
name<>
data<>

Reimplemented in vtkMesh.

Definition at line 327 of file meshBase.H.

328  {}

◆ stitchMB() [1/2]

meshBase * meshBase::stitchMB ( const std::vector< meshBase *> &  mbObjs)
static

Caller must delete object after use.

Parameters
mbObjsa vector of meshBase objects to stich together
Returns
<>

Definition at line 196 of file meshBase.C.

References Create(), and NEM::MSH::New().

Referenced by stitchMB().

196  {
197  if (!mbObjs.empty()) {
198  vtkSmartPointer<vtkAppendFilter> appender =
200  appender->MergePointsOn();
201  for (int i = 0; i < mbObjs.size(); ++i) {
202  appender->AddInputData(mbObjs[i]->getDataSet());
203  }
204  appender->Update();
205  return meshBase::Create(appender->GetOutput(), "stitched.vtu");
206  } else {
207  std::cerr << "Nothing to stitch!" << std::endl;
208  exit(1);
209  }
210 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
static meshBase * Create(const std::string &fname)
Construct vtkMesh from filename.
Definition: meshBase.C:78
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308

◆ stitchMB() [2/2]

std::shared_ptr< meshBase > meshBase::stitchMB ( const std::vector< std::shared_ptr< meshBase >> &  _mbObjs)
static

This is the shared pointer version of stitchMB.

Parameters
mbObjsa vector of meshBase objects to stich together
Returns
<>

Memory is managed by shared pointer, so do not call delete after use.

Definition at line 215 of file meshBase.C.

References CreateShared(), and stitchMB().

216  {
217  std::vector<meshBase *> mbObjs(_mbObjs.size());
218  for (int i = 0; i < mbObjs.size(); ++i) {
219  mbObjs[i] = _mbObjs[i].get();
220  }
222 }
static meshBase * stitchMB(const std::vector< meshBase *> &mbObjs)
stitch together several meshBases
Definition: meshBase.C:196
static std::shared_ptr< meshBase > CreateShared(const std::string &fname)
Create shared ptr from fname.
Definition: meshBase.C:171

◆ unsetCellDataArray() [1/2]

virtual void meshBase::unsetCellDataArray ( int  arrayID)
inlinevirtual
Parameters
arrayID<>

Reimplemented in vtkMesh.

Definition at line 391 of file meshBase.H.

Referenced by NEM::ADP::Refine::initAdaptive(), FETransfer::transferCellData(), and vtkMesh::write().

391 {}

◆ unsetCellDataArray() [2/2]

virtual void meshBase::unsetCellDataArray ( const std::string &  name)
inlinevirtual
Parameters
name<>

Reimplemented in vtkMesh.

Definition at line 396 of file meshBase.H.

396 {}

◆ unsetFieldDataArray()

virtual void meshBase::unsetFieldDataArray ( const std::string &  name)
inlinevirtual
Parameters
name<>

Reimplemented in vtkMesh.

Definition at line 401 of file meshBase.H.

Referenced by vtkMesh::write().

401 {}

◆ unsetPointDataArray() [1/2]

virtual void meshBase::unsetPointDataArray ( int  arrayID)
inlinevirtual
Parameters
arrayID<>

Reimplemented in vtkMesh.

Definition at line 381 of file meshBase.H.

Referenced by FETransfer::transferPointData(), and vtkMesh::write().

381 {}

◆ unsetPointDataArray() [2/2]

virtual void meshBase::unsetPointDataArray ( const std::string &  name)
inlinevirtual
Parameters
name<>

Reimplemented in vtkMesh.

Definition at line 386 of file meshBase.H.

386 {}

◆ write() [1/2]

virtual void meshBase::write ( ) const
inlinevirtual

The file extension of the private var "filename" determines the format of the output file

Reimplemented in vtkMesh.

Definition at line 598 of file meshBase.H.

References write().

Referenced by MeshQuality::checkMesh(), OrderOfAccuracy::computeMeshWithResolution(), NEM::DRV::GmshToVtkConversionDriver::execute(), NEM::DRV::GmshToExoConversionDriver::execute(), NEM::GEO::rocPack::geomToVTK(), proteusHdf5::proteusHdf5(), meshSrch::read(), NEM::ADP::Refine::run(), vtkMesh::write(), write(), and hdf5Reader::writeVTK().

598 { write(filename); }
virtual void write() const
write the mesh to file named after the private var &#39;filename&#39;.
Definition: meshBase.H:598
std::string filename
name of mesh file
Definition: meshBase.H:726

◆ write() [2/2]

virtual void meshBase::write ( const std::string &  fname) const
pure virtual
Parameters
fnameThe name of the file to write to

Implemented in FOAM::foamMesh, vtkMesh, meshSrch, and gmshMesh.

◆ writeCobalt() [1/2]

void meshBase::writeCobalt ( meshBase surfWithPatch,
const std::string &  mapFile,
std::ofstream &  outputStream 
)

for rocstar restart hack through rflupart/prep

Parameters
surfWithPatch<>
mapFile<>
outputStream<>

Definition at line 1407 of file meshBase.C.

References buildStaticCellLocator(), getDataSet(), getNumberOfCells(), IsArrayName(), NEM::MSH::New(), and writePatchMap().

1409  {
1410  if (!surfWithPatches) {
1411  std::cout << "surface mesh is empty!" << std::endl;
1412  exit(1);
1413  }
1414  if (surfWithPatches->IsArrayName("patchNo", 1) == -1) {
1415  std::cout << "surface mesh must have patchNo cell array" << std::endl;
1416  exit(1);
1417  }
1418  vtkSmartPointer<vtkIdList> facePtIds;
1419  vtkSmartPointer<vtkIdList> sharedCellPtIds =
1421  vtkSmartPointer<vtkGenericCell> genCell1 =
1423  vtkSmartPointer<vtkGenericCell> genCell2 =
1425  std::map<std::vector<nemId_t>, std::pair<nemId_t, nemId_t>,
1427  faceMap;
1428  // building cell locator for looking up patch number in remeshed surface mesh
1429  vtkSmartPointer<vtkStaticCellLocator> surfCellLocator =
1430  surfWithPatches->buildStaticCellLocator();
1431  // maximum number of vertices per face (to be found in proceeding loop)
1432  int nVerticesPerFaceMax = 0;
1433  // maximum number of faces per cell (to be found in proceeding loop)
1434  int nFacesPerCellMax = 0;
1435 
1436  for (nemId_t i = 0; i < this->getNumberOfCells(); ++i) {
1437  // get cell i
1438  dataSet->GetCell(i, genCell1);
1439  // get faces, find cells sharing it. if no cell shares it,
1440  // use the locator of the surfWithPatches to find the patch number
1441  int numFaces = genCell1->GetNumberOfFaces();
1442  nFacesPerCellMax =
1443  (nFacesPerCellMax < numFaces ? numFaces : nFacesPerCellMax);
1444  for (int j = 0; j < numFaces; ++j) {
1445  vtkCell *face = genCell1->GetFace(j);
1446  // bool shared = false;
1447  vtkIdType numVerts = face->GetNumberOfPoints();
1448  nVerticesPerFaceMax =
1449  (nVerticesPerFaceMax < numVerts ? numVerts : nVerticesPerFaceMax);
1450  facePtIds = face->GetPointIds();
1451  dataSet->GetCellNeighbors(i, facePtIds, sharedCellPtIds);
1452  std::vector<nemId_t> facePntIds(numVerts);
1453  for (vtkIdType k = 0; k < numVerts; ++k) {
1454  facePntIds[k] = face->GetPointId(k) + 1;
1455  }
1456  // std::cout << "sharedCellPtIds->GetNumberOfIds() = " <<
1457  // sharedCellPtIds->GetNumberOfIds() << std::endl;
1458  if (sharedCellPtIds->GetNumberOfIds()) {
1459  faceMap.insert(
1460  std::pair<std::vector<nemId_t>, std::pair<nemId_t, nemId_t>>(
1461  facePntIds,
1462  std::make_pair(i + 1, sharedCellPtIds->GetId(0) + 1)));
1463  } else {
1464  double p1[3], p2[3], p3[3];
1465  face->GetPoints()->GetPoint(0, p1);
1466  face->GetPoints()->GetPoint(1, p2);
1467  face->GetPoints()->GetPoint(2, p3);
1468  double faceCenter[3];
1469  for (vtkIdType k = 0; k < numVerts; ++k) {
1470  faceCenter[k] = (p1[k] + p2[k] + p3[k]) / 3.0;
1471  }
1472  vtkIdType closestCellId;
1473  int subId;
1474  double minDist2;
1475  double closestPoint[3];
1476  // find closest point and closest cell to faceCenter
1477  surfCellLocator->FindClosestPoint(faceCenter, closestPoint, genCell2,
1478  closestCellId, subId, minDist2);
1479  double patchNo[1];
1480  surfWithPatches->getDataSet()
1481  ->GetCellData()
1482  ->GetArray("patchNo")
1483  ->GetTuple(closestCellId, patchNo);
1484  faceMap.insert(
1485  std::pair<std::vector<nemId_t>, std::pair<nemId_t, nemId_t>>(
1486  facePntIds, std::make_pair(i + 1, -1 * patchNo[0])));
1487  }
1488  }
1489  }
1490 
1491  std::map<int, int> patchMap;
1492  for (int i = 0; i < surfWithPatches->getNumberOfCells(); ++i) {
1493  double patchNo[1];
1494  surfWithPatches->getDataSet()->GetCellData()->GetArray("patchNo")->GetTuple(
1495  i, patchNo);
1496  patchMap.insert(std::pair<int, int>(patchNo[0], i));
1497  }
1498 
1499  // write patch mapping file
1500  writePatchMap(mapFile, patchMap);
1501  // write cobalt file
1502  outputStream << 3 << " " << 1 << " " << patchMap.size() << std::endl;
1503  outputStream << this->getNumberOfPoints() << " " << faceMap.size() << " "
1504  << this->getNumberOfCells() << " " << nVerticesPerFaceMax << " "
1505  << nFacesPerCellMax << std::endl;
1506  for (int i = 0; i < this->getNumberOfPoints(); ++i) {
1507  std::vector<double> pnt(this->getPoint(i));
1508  outputStream << std::setw(21) << std::fixed << std::setprecision(15)
1509  << pnt[0] << " " << pnt[1] << " " << pnt[2] << std::endl;
1510  }
1511 
1512  auto it = faceMap.begin();
1513  while (it != faceMap.end()) {
1514  outputStream << it->first.size() << " ";
1515  auto faceIdIter = it->first.begin();
1516  while (faceIdIter != it->first.end()) {
1517  outputStream << *faceIdIter << " ";
1518  ++faceIdIter;
1519  }
1520  outputStream << it->second.first << " " << it->second.second << std::endl;
1521  ++it;
1522  }
1523 }
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
std::size_t nemId_t
Definition: meshBase.H:51
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
void writePatchMap(const std::string &mapFile, const std::map< int, int > &patchMap)
write patch map file for roc prep (trivial identity mapping)
Definition: meshBase.C:1376
nemId_t getNumberOfPoints() const
return the number of points
Definition: meshBase.H:545
sum comparison for vectors representing faces inserted into map
Definition: meshBase.H:765
nemId_t getNumberOfCells() const
return the number of cells
Definition: meshBase.H:550
virtual std::vector< double > getPoint(nemId_t id) const =0
get point with id

◆ writeCobalt() [2/2]

void meshBase::writeCobalt ( meshBase surfWithPatch,
const std::string &  mapFile,
const std::string &  ofname 
)

for rocstar restart hack through rflupart/prep

Parameters
surfWithPatch<>
mapFile<>
ofname<>

Definition at line 1527 of file meshBase.C.

1529  {
1530  std::ofstream outputStream(ofname);
1531  if (!outputStream.good()) {
1532  std::cout << "Cannot open file " << ofname << std::endl;
1533  exit(1);
1534  }
1535  writeCobalt(surfWithPatches, mapFile, outputStream);
1536 }
void writeCobalt(meshBase *surfWithPatch, const std::string &mapFile, std::ofstream &outputStream)
surfWithPatch must have patchNo array
Definition: meshBase.C:1407

◆ writeMSH() [1/6]

void meshBase::writeMSH ( std::ofstream &  outputStream)
Parameters
outputStream<>

Definition at line 1078 of file meshBase.C.

Referenced by NEM::ADP::Refine::initAdaptive().

1078  {
1079  if (!outputStream.good()) {
1080  std::cout << "Output file stream is bad" << std::endl;
1081  exit(1);
1082  }
1083 
1084  if (!dataSet) {
1085  std::cout << "No data to write" << std::endl;
1086  exit(1);
1087  }
1088  // --------- writing gmsh header ----------- //
1089  outputStream << "$MeshFormat" << std::endl
1090  << "2.2 0 8" << std::endl
1091  << "$EndMeshFormat" << std::endl;
1092 
1093  // -------- ensure all cell types are tri/tet or below -------------- //
1094  int type_id;
1095  for (int i = 0; i < numCells; i++) {
1096  type_id = dataSet->GetCellType(i);
1097  if (!(type_id == 3 || type_id == 5 || type_id == 10 || type_id == 9)) {
1098  std::cout << "Error: Only tetrahedral, triangular, and quadrilateral"
1099  << " meshes can be written to gmsh format" << std::endl;
1100  exit(3);
1101  }
1102  }
1103 
1104  // ------------------------ write point coords -------------------------- //
1105  outputStream << "$Nodes" << std::endl << numPoints << std::endl;
1106  for (int i = 0; i < numPoints; ++i) {
1107  std::vector<double> pntcrds = getPoint(i);
1108  outputStream << i + 1 << " " << pntcrds[0] << " " << pntcrds[1] << " "
1109  << pntcrds[2] << " " << std::endl;
1110  }
1111  outputStream << "$EndNodes" << std::endl;
1112 
1113  // ------------- write element type and connectivity --------------------- //
1114  outputStream << "$Elements" << std::endl << numCells << std::endl;
1115  for (int i = 0; i < numCells; ++i) {
1116  vtkIdList *point_ids = dataSet->GetCell(i)->GetPointIds();
1117  int numComponent = point_ids->GetNumberOfIds();
1118  type_id = dataSet->GetCell(i)->GetCellType();
1119  outputStream << i + 1 << " ";
1120  switch (numComponent) {
1121  case 2: {
1122  outputStream << 1 << " " << 2 << " " << 1 << " " << 1 << " ";
1123  break;
1124  }
1125  case 3: {
1126  outputStream << 2 << " " << 2 << " " << 1 << " " << 1 << " ";
1127  break;
1128  }
1129  case 4: {
1130  if (type_id == 10) // tetra
1131  {
1132  outputStream << 4 << " " << 2 << " " << 1 << " " << 1 << " ";
1133  break;
1134  } else if (type_id == 9) // quad
1135  {
1136  outputStream << 3 << " " << 2 << " " << 1 << " " << 1 << " ";
1137  break;
1138  } else
1139  {
1140  std::cerr << "Cells with 4 components must be tet or quad."
1141  << std::endl;
1142  exit(1);
1143  }
1144  }
1145 
1146  default: {
1147  std::cerr << "Components in cell should be <= 4" << std::endl;
1148  exit(1);
1149  }
1150  }
1151  for (int j = 0; j < numComponent; ++j)
1152  outputStream << point_ids->GetId(j) + 1 << " ";
1153  outputStream << std::endl;
1154  }
1155  outputStream << "$EndElements" << std::endl;
1156 }
nemId_t numPoints
number of points in mesh
Definition: meshBase.H:714
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
nemId_t numCells
number of cells in mesh
Definition: meshBase.H:718
virtual std::vector< double > getPoint(nemId_t id) const =0
get point with id

◆ writeMSH() [2/6]

void meshBase::writeMSH ( const std::string &  fname)
Parameters
fnameThe name of the file to write to

Definition at line 1549 of file meshBase.C.

1549  {
1550  std::ofstream outputStream(fname.c_str());
1551  writeMSH(outputStream);
1552 }
void writeMSH(std::ofstream &outputStream)
convert to gmsh format without data
Definition: meshBase.C:1078

◆ writeMSH() [3/6]

void meshBase::writeMSH ( std::ofstream &  outputStream,
const std::string &  pointOrCell,
int  arrayID 
)
Parameters
outputStream<>
pointOrCell<>
arrayID<>

Definition at line 1160 of file meshBase.C.

References data.

1161  {
1162  // write points and cells
1163  writeMSH(outputStream);
1164 
1165  if (!outputStream.good()) {
1166  std::cout << "Output file stream is bad" << std::endl;
1167  exit(1);
1168  }
1169 
1170  if (!pointOrCell.compare("point")) {
1171  // ---------------------------- write point data -------------------------
1172  // //
1173  vtkPointData *pointData = dataSet->GetPointData();
1174  if (pointData) {
1175  int numArr = pointData->GetNumberOfArrays();
1176  if (arrayID >= numArr) {
1177  std::cout << "ERROR: arrayID is out of bounds" << std::endl;
1178  std::cout << "There are " << numArr << " point data arrays"
1179  << std::endl;
1180  exit(1);
1181  } else if (numArr < 1) {
1182  std::cout << "no point data found" << std::endl;
1183  exit(1);
1184  }
1185  vtkDataArray *da = pointData->GetArray(arrayID);
1186  if (da) {
1187  int numComponent = da->GetNumberOfComponents();
1188  int numTuple = da->GetNumberOfTuples();
1189  std::string tmpname = "PointArray";
1190  tmpname += std::to_string(arrayID);
1191  outputStream << "$NodeData" << std::endl
1192  << 1 << std::endl // 1 string tag
1193  << "\""
1194  << (pointData->GetArrayName(arrayID)
1195  ? pointData->GetArrayName(arrayID)
1196  : tmpname) // name of view
1197  << "\"" << std::endl
1198  << 0 << std::endl // 0 real tag
1199  << 3 << std::endl // 3 int tags (dt index, dim of field,
1200  // number of fields)
1201  << 0 << std::endl // dt index
1202  << numComponent << std::endl // dim of field
1203  << numTuple << std::endl; // number of fields
1204  for (int j = 0; j < numTuple; ++j) {
1205  double *data = da->GetTuple(j);
1206  outputStream << j + 1 << " ";
1207  for (int k = 0; k < numComponent; ++k) {
1208  outputStream << data[k] << " ";
1209  }
1210  outputStream << std::endl;
1211  }
1212  outputStream << "$EndNodeData" << std::endl;
1213  }
1214  }
1215  } else if (!pointOrCell.compare("cell")) {
1216  // -------------------------- write cell data ----------------------------
1217  // //
1218 
1219  vtkCellData *cellData = dataSet->GetCellData();
1220  if (cellData) {
1221  int numArr = cellData->GetNumberOfArrays();
1222  if (arrayID >= numArr) {
1223  std::cout << "ERROR: arrayID is out of bounds" << std::endl;
1224  std::cout << "There are " << numArr << " cell data arrays" << std::endl;
1225  exit(1);
1226  } else if (numArr < 1) {
1227  std::cout << "no cell data found" << std::endl;
1228  exit(1);
1229  }
1230  vtkDataArray *da = cellData->GetArray(arrayID);
1231  if (da) {
1232  int numComponent = da->GetNumberOfComponents();
1233  int numTuple = da->GetNumberOfTuples();
1234  std::string tmpname = "CellArray";
1235  tmpname += std::to_string(arrayID);
1236  outputStream << "$ElementData" << std::endl
1237  << 1 << std::endl // 1 string tag
1238  << "\""
1239  << (cellData->GetArrayName(arrayID)
1240  ? cellData->GetArrayName(arrayID)
1241  : tmpname) // name of view
1242  << "\"" << std::endl
1243  << 0 << std::endl // 0 real tag
1244  << 3 << std::endl // 3 int tags (dt index, dim of field,
1245  // number of fields)
1246  << 0 << std::endl // dt index
1247  << numComponent << std::endl // dim of field
1248  << numTuple << std::endl; // number of fields
1249  for (int j = 0; j < numTuple; ++j) {
1250  double *data = da->GetTuple(j);
1251  outputStream << j + 1 << " ";
1252  for (int k = 0; k < numComponent; ++k) {
1253  outputStream << data[k] << " ";
1254  }
1255  outputStream << std::endl;
1256  }
1257  outputStream << "$EndElementData" << std::endl;
1258  }
1259  }
1260  }
1261 }
data_type data
Edge/face with sorted point ids (a, b, c, ...) is located at some index i in data[b], with data[b][i].first == [a, c] (for edges, third point id treated as -1).
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
void writeMSH(std::ofstream &outputStream)
convert to gmsh format without data
Definition: meshBase.C:1078

◆ writeMSH() [4/6]

void meshBase::writeMSH ( const std::string &  fname,
const std::string &  pointOrCell,
int  arrayID 
)
Parameters
fnameThe name of the file to write to
pointOrCell<>
arrayID<>

Definition at line 1556 of file meshBase.C.

1557  {
1558  std::ofstream outputStream(fname.c_str());
1559  writeMSH(outputStream, pointOrCell, arrayID);
1560 }
void writeMSH(std::ofstream &outputStream)
convert to gmsh format without data
Definition: meshBase.C:1078

◆ writeMSH() [5/6]

void meshBase::writeMSH ( std::ofstream &  outputStream,
const std::string &  pointOrCell,
int  arrayID,
bool  onlyVol 
)

convert to gmsh format with specified point or cell data for

Parameters
outputStream<>
pointOrCell<>
arrayID<>
onlyVol<>

Definition at line 1265 of file meshBase.C.

References data.

1267  {
1268  if (!outputStream.good()) {
1269  std::cout << "Output file stream is bad" << std::endl;
1270  exit(1);
1271  }
1272 
1273  if (!dataSet) {
1274  std::cout << "No data to write" << std::endl;
1275  exit(2);
1276  }
1277  // --------- writing gmsh header ----------- //
1278  outputStream << "$MeshFormat" << std::endl
1279  << "2.2 0 8" << std::endl
1280  << "$EndMeshFormat" << std::endl;
1281 
1282  // ---- get number of points and number of elements ---- //
1283  getNumberOfCells();
1285 
1286  // -------- ensure all cell types are tri/tet or below -------------- //
1287  int num_bad = 0;
1288  for (int i = 0; i < numCells; i++) {
1289  int type_id = dataSet->GetCellType(i);
1290  if (!(type_id == 3 || type_id == 5 || type_id == 10 || type_id == 9)) {
1291  std::cout << "Error: Only tetrahedral, triangular, and quadrilateral"
1292  << " meshes can be written to gmsh format" << std::endl;
1293  exit(3);
1294  }
1295  if (!(type_id == 10)) num_bad += 1;
1296  }
1297 
1298  // ------------------------ write point coords -------------------------- //
1299  outputStream << "$Nodes" << std::endl << numPoints << std::endl;
1300  for (int i = 0; i < numPoints; ++i) {
1301  std::vector<double> pntcrds = getPoint(i);
1302  outputStream << i + 1 << " ";
1303  outputStream << std::setprecision(16) << pntcrds[0] << " " << pntcrds[1]
1304  << " " << pntcrds[2] << " " << std::endl;
1305  }
1306  outputStream << "$EndNodes" << std::endl;
1307 
1308  // ------------- write element type and connectivity --------------------- //
1309  outputStream << "$Elements" << std::endl << numCells - num_bad << std::endl;
1310  // int k = 0;
1311  for (int i = 0; i < numCells; ++i) {
1312  vtkIdList *point_ids = dataSet->GetCell(i)->GetPointIds();
1313  int numComponent = point_ids->GetNumberOfIds();
1314  int type_id = dataSet->GetCellType(i);
1315  if (type_id == 10) {
1316  outputStream << i + 1 << " ";
1317  switch (numComponent) {
1318  case 2: {
1319  break;
1320  }
1321  case 3: {
1322  outputStream << 2 << " " << 2 << " " << 1 << " " << 1 << " ";
1323  break;
1324  }
1325  case 4: {
1326  outputStream << 4 << " " << 2 << " " << 1 << " " << 1 << " ";
1327  break;
1328  }
1329 
1330  default: {
1331  std::cerr << "Components in cell should be less than 4" << std::endl;
1332  exit(1);
1333  }
1334  }
1335  for (int j = 0; j < numComponent; ++j)
1336  outputStream << point_ids->GetId(j) + 1 << " ";
1337  outputStream << std::endl;
1338  // k+=1;
1339  }
1340  }
1341  outputStream << "$EndElements" << std::endl;
1342  // -------------------------- write cell data ---------------------------- //
1343  vtkCellData *cellData = dataSet->GetCellData();
1344  vtkDataArray *da = cellData->GetArray(arrayID);
1345  if (da) {
1346  std::string tmpname = "CellArray";
1347  tmpname += std::to_string(arrayID);
1348  outputStream
1349  << "$ElementData" << std::endl
1350  << 1 << std::endl // 1 string tag
1351  << "\""
1352  << (cellData->GetArrayName(arrayID) ? cellData->GetArrayName(arrayID)
1353  : tmpname) // name of view
1354  << "\"" << std::endl
1355  << 0 << std::endl // 0 real tag
1356  << 3
1357  << std::endl // 3 int tags (dt index, dim of field, number of fields)
1358  << 0 << std::endl // dt index
1359  << 1 << std::endl // dim of field
1360  << numCells - num_bad << std::endl; // number of fields
1361  for (int j = 0; j < numCells; ++j) {
1362  int type_id = dataSet->GetCellType(j);
1363  if (type_id == 10) {
1364  double *data = da->GetTuple(j);
1365  outputStream << j + 1 << " ";
1366  outputStream << data[0] << " ";
1367  outputStream << std::endl;
1368  }
1369  }
1370  outputStream << "$EndElementData" << std::endl;
1371  }
1372 }
data_type data
Edge/face with sorted point ids (a, b, c, ...) is located at some index i in data[b], with data[b][i].first == [a, c] (for edges, third point id treated as -1).
nemId_t numPoints
number of points in mesh
Definition: meshBase.H:714
vtkSmartPointer< vtkDataSet > dataSet
mesh points, topology and data
Definition: meshBase.H:722
nemId_t numCells
number of cells in mesh
Definition: meshBase.H:718
nemId_t getNumberOfPoints() const
return the number of points
Definition: meshBase.H:545
nemId_t getNumberOfCells() const
return the number of cells
Definition: meshBase.H:550
virtual std::vector< double > getPoint(nemId_t id) const =0
get point with id

◆ writeMSH() [6/6]

void meshBase::writeMSH ( const std::string &  fname,
const std::string &  pointOrCell,
int  arrayID,
bool  onlyVol 
)
Parameters
fnameThe name of the file to write to
pointOrCell<>
arrayID<>
onlyVol<>

Definition at line 1540 of file meshBase.C.

1542  {
1543  std::ofstream outputStream(fname.c_str());
1544  writeMSH(outputStream, pointOrCell, arrayID, onlyVol);
1545 }
void writeMSH(std::ofstream &outputStream)
convert to gmsh format without data
Definition: meshBase.C:1078

Member Data Documentation

◆ checkQuality

bool meshBase::checkQuality
protected

Definition at line 730 of file meshBase.H.

◆ continuous

bool meshBase::continuous
protected

Definition at line 735 of file meshBase.H.

◆ dataSet

◆ filename

std::string meshBase::filename
protected

Definition at line 726 of file meshBase.H.

Referenced by FOAM::foamMesh::report(), vtkMesh::report(), and vtkMesh::vtkMesh().

◆ globToPartCellMap

std::map<nemId_t, nemId_t> meshBase::globToPartCellMap
protected

Definition at line 750 of file meshBase.H.

◆ globToPartNodeMap

std::map<nemId_t, nemId_t> meshBase::globToPartNodeMap
protected

Only populated for mesh resulting from call to meshBase::partition

Definition at line 746 of file meshBase.H.

◆ metadata

vtkSmartPointer<vtkModelMetadata> meshBase::metadata
protected

Definition at line 760 of file meshBase.H.

◆ newArrayNames

std::vector<std::string> meshBase::newArrayNames
protected

Definition at line 739 of file meshBase.H.

◆ numCells

◆ numPoints

◆ partToGlobCellMap

std::map<nemId_t, nemId_t> meshBase::partToGlobCellMap
protected

Definition at line 758 of file meshBase.H.

◆ partToGlobNodeMap

std::map<nemId_t, nemId_t> meshBase::partToGlobNodeMap
protected

Definition at line 754 of file meshBase.H.


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