Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::SpatialDomains::MeshGraphXmlCompressed Class Reference

#include <MeshGraphXmlCompressed.h>

Inheritance diagram for Nektar::SpatialDomains::MeshGraphXmlCompressed:
[legend]

Public Member Functions

 MeshGraphXmlCompressed ()
 
virtual ~MeshGraphXmlCompressed ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::MeshGraphXml
 MeshGraphXml ()
 
virtual ~MeshGraphXml ()
 
void WriteXMLGeometry (std::string outname, std::vector< std::set< unsigned int >> elements, std::vector< unsigned int > partitions)
 
- Public Member Functions inherited from Nektar::SpatialDomains::MeshGraph
 MeshGraph ()
 
virtual ~MeshGraph ()
 
void WriteGeometry (std::string &outfilename, bool defaultExp=false, const LibUtilities::FieldMetaDataMap &metadata=LibUtilities::NullFieldMetaDataMap)
 
void Empty (int dim, int space)
 
void FillGraph ()
 
void FillBoundingBoxTree ()
 
std::vector< int > GetElementsContainingPoint (PointGeomSharedPtr p)
 
void ReadExpansionInfo ()
 
int GetMeshDimension ()
 Dimension of the mesh (can be a 1D curve in 3D space). More...
 
int GetSpaceDimension ()
 Dimension of the space (can be a 1D curve in 3D space). More...
 
void SetDomainRange (NekDouble xmin, NekDouble xmax, NekDouble ymin=NekConstants::kNekUnsetDouble, NekDouble ymax=NekConstants::kNekUnsetDouble, NekDouble zmin=NekConstants::kNekUnsetDouble, NekDouble zmax=NekConstants::kNekUnsetDouble)
 
bool CheckRange (Geometry2D &geom)
 Check if goemetry is in range definition if activated. More...
 
bool CheckRange (Geometry3D &geom)
 Check if goemetry is in range definition if activated. More...
 
CompositeSharedPtr GetComposite (int whichComposite)
 
GeometrySharedPtr GetCompositeItem (int whichComposite, int whichItem)
 
void GetCompositeList (const std::string &compositeStr, CompositeMap &compositeVector) const
 
std::map< int, CompositeSharedPtr > & GetComposites ()
 
std::map< int, std::string > & GetCompositesLabels ()
 
std::map< int, std::map< int, CompositeSharedPtr > > & GetDomain ()
 
std::map< int, CompositeSharedPtr > & GetDomain (int domain)
 
const ExpansionInfoMapGetExpansionInfo (const std::string variable="DefaultVar")
 
ExpansionInfoShPtr GetExpansionInfo (GeometrySharedPtr geom, const std::string variable="DefaultVar")
 
void SetExpansionInfo (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
 Sets expansions given field definitions. More...
 
void SetExpansionInfo (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, std::vector< std::vector< LibUtilities::PointsType >> &pointstype)
 Sets expansions given field definition, quadrature points. More...
 
void SetExpansionInfoToEvenlySpacedPoints (int npoints=0)
 Sets expansions to have equispaced points. More...
 
void SetExpansionInfoToNumModes (int nmodes)
 Reset expansion to have specified polynomial order nmodes. More...
 
void SetExpansionInfoToPointOrder (int npts)
 Reset expansion to have specified point order npts. More...
 
void SetExpansionInfo (const std::string variable, ExpansionInfoMapShPtr &exp)
 This function sets the expansion #exp in map with entry #variable. More...
 
void SetSession (LibUtilities::SessionReaderSharedPtr pSession)
 
void SetBasisKey (LibUtilities::ShapeType shape, LibUtilities::BasisKeyVector &keys, std::string var="DefaultVar")
 Sets the basis key for all expansions of the given shape. More...
 
void ResetExpansionInfoToBasisKey (ExpansionInfoMapShPtr &expansionMap, LibUtilities::ShapeType shape, LibUtilities::BasisKeyVector &keys)
 
bool SameExpansionInfo (const std::string var1, const std::string var2)
 
bool ExpansionInfoDefined (const std::string var)
 
bool CheckForGeomInfo (std::string parameter)
 
const std::string GetGeomInfo (std::string parameter)
 
LibUtilities::BasisKeyVector DefineBasisKeyFromExpansionTypeHomo (GeometrySharedPtr in, ExpansionType type_x, ExpansionType type_y, ExpansionType type_z, const int nummodes_x, const int nummodes_y, const int nummodes_z)
 
int GetNvertices ()
 
PointGeomSharedPtr GetVertex (int id)
 
SegGeomSharedPtr GetSegGeom (int id)
 
CurveMapGetCurvedEdges ()
 
CurveMapGetCurvedFaces ()
 
std::map< int, PointGeomSharedPtr > & GetAllPointGeoms ()
 
std::map< int, SegGeomSharedPtr > & GetAllSegGeoms ()
 
TriGeomMapGetAllTriGeoms ()
 
QuadGeomMapGetAllQuadGeoms ()
 
TetGeomMapGetAllTetGeoms ()
 
PyrGeomMapGetAllPyrGeoms ()
 
PrismGeomMapGetAllPrismGeoms ()
 
HexGeomMapGetAllHexGeoms ()
 
int GetNumElements ()
 
Geometry2DSharedPtr GetGeometry2D (int gID)
 
LibUtilities::BasisKey GetEdgeBasisKey (SegGeomSharedPtr edge, const std::string variable="DefaultVar")
 
GeometryLinkSharedPtr GetElementsFromEdge (Geometry1DSharedPtr edge)
 
GeometryLinkSharedPtr GetElementsFromFace (Geometry2DSharedPtr face)
 
LibUtilities::BasisKey GetFaceBasisKey (Geometry2DSharedPtr face, const int facedir, const std::string variable="DefaultVar")
 3D functions More...
 
CompositeOrderingGetCompositeOrdering ()
 
void SetCompositeOrdering (CompositeOrdering p_compOrder)
 
BndRegionOrderingGetBndRegionOrdering ()
 
void SetBndRegionOrdering (BndRegionOrdering p_bndRegOrder)
 
void ReadGeometry (LibUtilities::DomainRangeShPtr rng, bool fillGraph)
 
void PartitionMesh (LibUtilities::SessionReaderSharedPtr session)
 
std::map< int, MeshEntityCreateMeshEntities ()
 Create mesh entities for this graph. More...
 
CompositeDescriptor CreateCompositeDescriptor ()
 
MovementSharedPtrGetMovement ()
 

Static Public Member Functions

static MeshGraphSharedPtr create ()
 
- Static Public Member Functions inherited from Nektar::SpatialDomains::MeshGraphXml
static MeshGraphSharedPtr create ()
 
- Static Public Member Functions inherited from Nektar::SpatialDomains::MeshGraph
static MeshGraphSharedPtr Read (const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true)
 
static LibUtilities::BasisKeyVector DefineBasisKeyFromExpansionType (GeometrySharedPtr in, ExpansionType type, const int order)
 

Static Public Attributes

static std::string className
 
- Static Public Attributes inherited from Nektar::SpatialDomains::MeshGraphXml
static std::string className
 

Protected Member Functions

virtual void v_ReadVertices () override
 
virtual void v_ReadCurves () override
 
virtual void v_ReadEdges () override
 
virtual void v_ReadFaces () override
 
virtual void v_ReadElements1D () override
 
virtual void v_ReadElements2D () override
 
virtual void v_ReadElements3D () override
 
virtual void v_WriteVertices (TiXmlElement *geomTag, PointGeomMap &verts) override
 
virtual void v_WriteEdges (TiXmlElement *geomTag, SegGeomMap &edges) override
 
virtual void v_WriteTris (TiXmlElement *faceTag, TriGeomMap &tris) override
 
virtual void v_WriteQuads (TiXmlElement *faceTag, QuadGeomMap &quads) override
 
virtual void v_WriteHexs (TiXmlElement *elmtTag, HexGeomMap &hexs) override
 
virtual void v_WritePrisms (TiXmlElement *elmtTag, PrismGeomMap &pris) override
 
virtual void v_WritePyrs (TiXmlElement *elmtTag, PyrGeomMap &pyrs) override
 
virtual void v_WriteTets (TiXmlElement *elmtTag, TetGeomMap &tets) override
 
virtual void v_WriteCurves (TiXmlElement *geomTag, CurveMap &edges, CurveMap &faces) override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::MeshGraphXml
virtual void v_WriteGeometry (std::string &outfilename, bool defaultExp=false, const LibUtilities::FieldMetaDataMap &metadata=LibUtilities::NullFieldMetaDataMap) override
 Write out an XML file containing the GEOMETRY block representing this MeshGraph instance inside a NEKTAR tag. More...
 
virtual void v_ReadGeometry (LibUtilities::DomainRangeShPtr rng, bool fillGraph) override
 
virtual void v_PartitionMesh (LibUtilities::SessionReaderSharedPtr session) override
 
void ReadDomain ()
 
void ReadElements ()
 
void ReadComposites ()
 
void ResolveGeomRef (const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
 
void ResolveGeomRef1D (const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
 
void ResolveGeomRef2D (const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
 
void ResolveGeomRef3D (const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
 
void WriteComposites (TiXmlElement *geomTag, CompositeMap &comps, std::map< int, std::string > &compLabels)
 
void WriteDomain (TiXmlElement *geomTag, std::map< int, CompositeMap > &domain)
 
void WriteDefaultExpansion (TiXmlElement *root)
 
CompositeOrdering CreateCompositeOrdering ()
 
- Protected Member Functions inherited from Nektar::SpatialDomains::MeshGraph
void PopulateFaceToElMap (Geometry3DSharedPtr element, int kNfaces)
 Given a 3D geometry object #element, populate the face to element map m_faceToElMap which maps faces to their corresponding element(s). More...
 
ExpansionInfoMapShPtr SetUpExpansionInfoMap ()
 
std::string GetCompositeString (CompositeSharedPtr comp)
 Returns a string representation of a composite. More...
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::SpatialDomains::MeshGraph
LibUtilities::SessionReaderSharedPtr m_session
 
PointGeomMap m_vertSet
 
CurveMap m_curvedEdges
 
CurveMap m_curvedFaces
 
SegGeomMap m_segGeoms
 
TriGeomMap m_triGeoms
 
QuadGeomMap m_quadGeoms
 
TetGeomMap m_tetGeoms
 
PyrGeomMap m_pyrGeoms
 
PrismGeomMap m_prismGeoms
 
HexGeomMap m_hexGeoms
 
int m_meshDimension
 
int m_spaceDimension
 
int m_partition
 
bool m_meshPartitioned
 
CompositeMap m_meshComposites
 
std::map< int, std::string > m_compositesLabels
 
std::map< int, CompositeMapm_domain
 
LibUtilities::DomainRangeShPtr m_domainRange
 
ExpansionInfoMapShPtrMap m_expansionMapShPtrMap
 
GeomInfoMap m_geomInfo
 
std::unordered_map< int, GeometryLinkSharedPtrm_faceToElMap
 
TiXmlElement * m_xmlGeom
 
CompositeOrdering m_compOrder
 
BndRegionOrdering m_bndRegOrder
 
std::unique_ptr< GeomRTreem_boundingBoxTree
 
MovementSharedPtr m_movement = nullptr
 

Detailed Description

Definition at line 44 of file MeshGraphXmlCompressed.h.

Constructor & Destructor Documentation

◆ MeshGraphXmlCompressed()

Nektar::SpatialDomains::MeshGraphXmlCompressed::MeshGraphXmlCompressed ( )
inline

Definition at line 47 of file MeshGraphXmlCompressed.h.

48  {
49  }

◆ ~MeshGraphXmlCompressed()

virtual Nektar::SpatialDomains::MeshGraphXmlCompressed::~MeshGraphXmlCompressed ( )
inlinevirtual

Definition at line 51 of file MeshGraphXmlCompressed.h.

52  {
53  }

Member Function Documentation

◆ create()

static MeshGraphSharedPtr Nektar::SpatialDomains::MeshGraphXmlCompressed::create ( )
inlinestatic

Definition at line 55 of file MeshGraphXmlCompressed.h.

56  {
58  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ v_ReadCurves()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_ReadCurves ( )
overrideprotectedvirtual

Look for elements in CURVE block.

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 206 of file MeshGraphXmlCompressed.cpp.

207 {
208  // check to see if any scaling parameters are in
209  // attributes and determine these values
210  TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
211  ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
212 
213  NekDouble xscale, yscale, zscale;
214 
215  LibUtilities::Interpreter expEvaluator;
216  const char *xscal = element->Attribute("XSCALE");
217  if (!xscal)
218  {
219  xscale = 1.0;
220  }
221  else
222  {
223  std::string xscalstr = xscal;
224  int expr_id = expEvaluator.DefineFunction("", xscalstr);
225  xscale = expEvaluator.Evaluate(expr_id);
226  }
227 
228  const char *yscal = element->Attribute("YSCALE");
229  if (!yscal)
230  {
231  yscale = 1.0;
232  }
233  else
234  {
235  std::string yscalstr = yscal;
236  int expr_id = expEvaluator.DefineFunction("", yscalstr);
237  yscale = expEvaluator.Evaluate(expr_id);
238  }
239 
240  const char *zscal = element->Attribute("ZSCALE");
241  if (!zscal)
242  {
243  zscale = 1.0;
244  }
245  else
246  {
247  std::string zscalstr = zscal;
248  int expr_id = expEvaluator.DefineFunction("", zscalstr);
249  zscale = expEvaluator.Evaluate(expr_id);
250  }
251 
252  NekDouble xmove, ymove, zmove;
253 
254  // check to see if any moving parameters are in
255  // attributes and determine these values
256 
257  const char *xmov = element->Attribute("XMOVE");
258  if (!xmov)
259  {
260  xmove = 0.0;
261  }
262  else
263  {
264  std::string xmovstr = xmov;
265  int expr_id = expEvaluator.DefineFunction("", xmovstr);
266  xmove = expEvaluator.Evaluate(expr_id);
267  }
268 
269  const char *ymov = element->Attribute("YMOVE");
270  if (!ymov)
271  {
272  ymove = 0.0;
273  }
274  else
275  {
276  std::string ymovstr = ymov;
277  int expr_id = expEvaluator.DefineFunction("", ymovstr);
278  ymove = expEvaluator.Evaluate(expr_id);
279  }
280 
281  const char *zmov = element->Attribute("ZMOVE");
282  if (!zmov)
283  {
284  zmove = 0.0;
285  }
286  else
287  {
288  std::string zmovstr = zmov;
289  int expr_id = expEvaluator.DefineFunction("", zmovstr);
290  zmove = expEvaluator.Evaluate(expr_id);
291  }
292 
293  /// Look for elements in CURVE block.
294  TiXmlElement *field = m_xmlGeom->FirstChildElement("CURVED");
295 
296  if (!field) // return if no curved entities
297  {
298  return;
299  }
300 
301  string IsCompressed;
302  field->QueryStringAttribute("COMPRESSED", &IsCompressed);
303 
304  if (IsCompressed.size() == 0)
305  {
306  // this could be that the curved tag is empty
307  // in this case we dont want to read it
308  return;
309  }
310 
311  ASSERTL0(boost::iequals(IsCompressed,
313  "Compressed formats do not match. Expected :" +
315  boost::lexical_cast<std::string>(IsCompressed));
316 
317  std::vector<SpatialDomains::MeshCurvedInfo> edginfo;
318  std::vector<SpatialDomains::MeshCurvedInfo> facinfo;
319  SpatialDomains::MeshCurvedPts cpts;
320 
321  // read edge, face info and curved poitns.
322  TiXmlElement *x = field->FirstChildElement();
323  while (x)
324  {
325  const char *entitytype = x->Value();
326  // read in edge or face info
327  if (boost::iequals(entitytype, "E"))
328  {
329  // read in data
330  std::string elmtStr;
331  TiXmlNode *child = x->FirstChild();
332 
333  if (child->Type() == TiXmlNode::TINYXML_TEXT)
334  {
335  elmtStr += child->ToText()->ValueStr();
336  }
337 
339  edginfo);
340  }
341  else if (boost::iequals(entitytype, "F"))
342  {
343  // read in data
344  std::string elmtStr;
345  TiXmlNode *child = x->FirstChild();
346 
347  if (child->Type() == TiXmlNode::TINYXML_TEXT)
348  {
349  elmtStr += child->ToText()->ValueStr();
350  }
351 
353  facinfo);
354  }
355  else if (boost::iequals(entitytype, "DATAPOINTS"))
356  {
357  NekInt id;
358  ASSERTL0(x->Attribute("ID", &id),
359  "Failed to get ID from PTS section");
360  cpts.id = id;
361 
362  // read in data
363  std::string elmtStr;
364 
365  TiXmlElement *DataIdx = x->FirstChildElement("INDEX");
366  ASSERTL0(DataIdx, "Cannot read data index tag in compressed "
367  "curved section");
368 
369  TiXmlNode *child = DataIdx->FirstChild();
370  if (child->Type() == TiXmlNode::TINYXML_TEXT)
371  {
372  elmtStr = child->ToText()->ValueStr();
373  }
374 
376  cpts.index);
377 
378  TiXmlElement *DataPts = x->FirstChildElement("POINTS");
379  ASSERTL0(DataPts, "Cannot read data pts tag in compressed "
380  "curved section");
381 
382  child = DataPts->FirstChild();
383  if (child->Type() == TiXmlNode::TINYXML_TEXT)
384  {
385  elmtStr = child->ToText()->ValueStr();
386  }
387 
389  cpts.pts);
390  }
391  else
392  {
393  ASSERTL0(false, "Unknown tag in curved section");
394  }
395  x = x->NextSiblingElement();
396  }
397 
398  // rescale (x,y,z) points;
399  for (int i = 0; i < cpts.pts.size(); ++i)
400  {
401  cpts.pts[i].x = xscale * cpts.pts[i].x + xmove;
402  cpts.pts[i].y = yscale * cpts.pts[i].y + ymove;
403  cpts.pts[i].z = zscale * cpts.pts[i].z + zmove;
404  }
405 
406  for (int i = 0; i < edginfo.size(); ++i)
407  {
408  int edgeid = edginfo[i].entityid;
410 
412  edgeid, ptype = (LibUtilities::PointsType)edginfo[i].ptype));
413 
414  // load points
415  int offset = edginfo[i].ptoffset;
416  for (int j = 0; j < edginfo[i].npoints; ++j)
417  {
418  int idx = cpts.index[offset + j];
419 
421  m_spaceDimension, edginfo[i].id, cpts.pts[idx].x,
422  cpts.pts[idx].y, cpts.pts[idx].z));
423  curve->m_points.push_back(vert);
424  }
425 
426  m_curvedEdges[edgeid] = curve;
427  }
428 
429  for (int i = 0; i < facinfo.size(); ++i)
430  {
431  int faceid = facinfo[i].entityid;
433 
435  faceid, ptype = (LibUtilities::PointsType)facinfo[i].ptype));
436 
437  int offset = facinfo[i].ptoffset;
438  for (int j = 0; j < facinfo[i].npoints; ++j)
439  {
440  int idx = cpts.index[offset + j];
441 
443  m_spaceDimension, facinfo[i].id, cpts.pts[idx].x,
444  cpts.pts[idx].y, cpts.pts[idx].z));
445  curve->m_points.push_back(vert);
446  }
447 
448  m_curvedFaces[faceid] = curve;
449  }
450 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
int ZlibDecodeFromBase64Str(std::string &in64, std::vector< T > &out)
Definition: CompressData.h:220
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:60
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:59
std::int32_t NekInt
double NekDouble

References ASSERTL0, Nektar::LibUtilities::Interpreter::DefineFunction(), Nektar::LibUtilities::Interpreter::Evaluate(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshCurvedPts::id, Nektar::SpatialDomains::MeshCurvedPts::index, Nektar::SpatialDomains::MeshCurvedPts::pts, and Nektar::LibUtilities::CompressData::ZlibDecodeFromBase64Str().

◆ v_ReadEdges()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_ReadEdges ( )
overrideprotectedvirtual

Look for elements in ELEMENT block.

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 452 of file MeshGraphXmlCompressed.cpp.

453 {
454  CurveMap::iterator it;
455 
456  /// Look for elements in ELEMENT block.
457  TiXmlElement *field = m_xmlGeom->FirstChildElement("EDGE");
458 
459  ASSERTL0(field, "Unable to find EDGE tag in file.");
460 
461  string IsCompressed;
462  field->QueryStringAttribute("COMPRESSED", &IsCompressed);
463 
464  ASSERTL0(boost::iequals(IsCompressed,
466  "Compressed formats do not match. Expected :" +
468  std::string(IsCompressed));
469  // Extract the edge body
470  TiXmlNode *edgeChild = field->FirstChild();
471  ASSERTL0(edgeChild, "Unable to extract the data from "
472  "the compressed edge tag.");
473 
474  std::string edgeStr;
475  if (edgeChild->Type() == TiXmlNode::TINYXML_TEXT)
476  {
477  edgeStr += edgeChild->ToText()->ValueStr();
478  }
479 
480  std::vector<SpatialDomains::MeshEdge> edgeData;
482 
483  int indx;
484  for (int i = 0; i < edgeData.size(); ++i)
485  {
486  indx = edgeData[i].id;
487  PointGeomSharedPtr vertices[2] = {GetVertex(edgeData[i].v0),
488  GetVertex(edgeData[i].v1)};
489  SegGeomSharedPtr edge;
490 
491  it = m_curvedEdges.find(indx);
492  if (it == m_curvedEdges.end())
493  {
495  indx, m_spaceDimension, vertices);
496  }
497  else
498  {
500  indx, m_spaceDimension, vertices, it->second);
501  }
502  m_segGeoms[indx] = edge;
503  }
504 }
PointGeomSharedPtr GetVertex(int id)
Definition: MeshGraph.h:343
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:62

References ASSERTL0, Nektar::LibUtilities::CompressData::GetCompressString(), and Nektar::LibUtilities::CompressData::ZlibDecodeFromBase64Str().

◆ v_ReadElements1D()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_ReadElements1D ( )
overrideprotectedvirtual

Look for elements in ELEMENT block.

All elements are of the form: "<S ID = n> ... </S>", with ? being the element type.

See if this face has curves.

Keep looking for additional segments

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 620 of file MeshGraphXmlCompressed.cpp.

621 {
622  TiXmlElement *field = NULL;
623 
624  /// Look for elements in ELEMENT block.
625  field = m_xmlGeom->FirstChildElement("ELEMENT");
626 
627  ASSERTL0(field, "Unable to find ELEMENT tag in file.");
628 
629  /// All elements are of the form: "<S ID = n> ... </S>", with
630  /// ? being the element type.
631 
632  TiXmlElement *segment = field->FirstChildElement("S");
633  CurveMap::iterator it;
634 
635  while (segment)
636  {
637  string IsCompressed;
638  segment->QueryStringAttribute("COMPRESSED", &IsCompressed);
639  ASSERTL0(
640  boost::iequals(IsCompressed,
642  "Compressed formats do not match. Expected :" +
644  std::string(IsCompressed));
645 
646  // Extract the face body
647  TiXmlNode *child = segment->FirstChild();
648  ASSERTL0(child, "Unable to extract the data from "
649  "the compressed face tag.");
650 
651  std::string str;
652  if (child->Type() == TiXmlNode::TINYXML_TEXT)
653  {
654  str += child->ToText()->ValueStr();
655  }
656 
657  int indx;
658 
659  std::vector<SpatialDomains::MeshEdge> data;
661 
662  for (int i = 0; i < data.size(); ++i)
663  {
664  indx = data[i].id;
665 
666  /// See if this face has curves.
667  it = m_curvedEdges.find(indx);
668 
669  PointGeomSharedPtr vertices[2] = {GetVertex(data[i].v0),
670  GetVertex(data[i].v1)};
671  SegGeomSharedPtr seg;
672 
673  if (it == m_curvedEdges.end())
674  {
676  indx, m_spaceDimension, vertices);
677  seg->SetGlobalID(indx); // Set global mesh id
678  }
679  else
680  {
682  indx, m_spaceDimension, vertices, it->second);
683  seg->SetGlobalID(indx); // Set global mesh id
684  }
685  seg->SetGlobalID(indx);
686  m_segGeoms[indx] = seg;
687  }
688  /// Keep looking for additional segments
689  segment = segment->NextSiblingElement("S");
690  }
691 }

References ASSERTL0, Nektar::LibUtilities::CompressData::GetCompressString(), and Nektar::LibUtilities::CompressData::ZlibDecodeFromBase64Str().

◆ v_ReadElements2D()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_ReadElements2D ( )
overrideprotectedvirtual

Look for elements in ELEMENT block.

All elements are of the form: "<? ID="#"> ... </?>", with ? being the element type.

See if this face has curves.

Create a TriGeom to hold the new definition.

See if this face has curves.

Create a QuadGeom to hold the new definition.

Keep looking

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 693 of file MeshGraphXmlCompressed.cpp.

694 {
695  /// Look for elements in ELEMENT block.
696  TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
697 
698  ASSERTL0(field, "Unable to find ELEMENT tag in file.");
699 
700  // Set up curve map for curved elements on an embedded manifold.
701  CurveMap::iterator it;
702 
703  /// All elements are of the form: "<? ID="#"> ... </?>", with
704  /// ? being the element type.
705 
706  TiXmlElement *element = field->FirstChildElement();
707 
708  while (element)
709  {
710  std::string elementType(element->ValueStr());
711 
712  ASSERTL0(
713  elementType == "Q" || elementType == "T",
714  (std::string("Unknown 2D element type: ") + elementType).c_str());
715 
716  string IsCompressed;
717  element->QueryStringAttribute("COMPRESSED", &IsCompressed);
718 
719  ASSERTL0(
720  boost::iequals(IsCompressed,
722  "Compressed formats do not match. Expected :" +
724  std::string(IsCompressed));
725 
726  // Extract the face body
727  TiXmlNode *faceChild = element->FirstChild();
728  ASSERTL0(faceChild, "Unable to extract the data from "
729  "the compressed face tag.");
730 
731  std::string faceStr;
732  if (faceChild->Type() == TiXmlNode::TINYXML_TEXT)
733  {
734  faceStr += faceChild->ToText()->ValueStr();
735  }
736 
737  int indx;
738  if (elementType == "T")
739  {
740  std::vector<SpatialDomains::MeshTri> faceData;
742  faceData);
743 
744  for (int i = 0; i < faceData.size(); ++i)
745  {
746  indx = faceData[i].id;
747 
748  /// See if this face has curves.
749  it = m_curvedFaces.find(indx);
750 
751  /// Create a TriGeom to hold the new definition.
753  GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
754  GetSegGeom(faceData[i].e[2])};
755 
756  TriGeomSharedPtr trigeom;
757  if (it == m_curvedFaces.end())
758  {
759  trigeom =
761  }
762  else
763  {
765  indx, edges, it->second);
766  }
767  trigeom->SetGlobalID(indx);
768  m_triGeoms[indx] = trigeom;
769  }
770  }
771  else if (elementType == "Q")
772  {
773  std::vector<SpatialDomains::MeshQuad> faceData;
775  faceData);
776 
777  for (int i = 0; i < faceData.size(); ++i)
778  {
779  indx = faceData[i].id;
780 
781  /// See if this face has curves.
782  it = m_curvedFaces.find(indx);
783 
784  /// Create a QuadGeom to hold the new definition.
786  GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
787  GetSegGeom(faceData[i].e[2]), GetSegGeom(faceData[i].e[3])};
788 
789  QuadGeomSharedPtr quadgeom;
790  if (it == m_curvedFaces.end())
791  {
792  quadgeom =
794  }
795  else
796  {
798  indx, edges, it->second);
799  }
800  quadgeom->SetGlobalID(indx);
801  m_quadGeoms[indx] = quadgeom;
802  }
803  }
804  /// Keep looking
805  element = element->NextSiblingElement();
806  }
807 }
SegGeomSharedPtr GetSegGeom(int id)
Definition: MeshGraph.h:348
static const int kNedges
Definition: QuadGeom.h:76
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:72
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: HexGeom.h:46
std::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58

References ASSERTL0, Nektar::LibUtilities::CompressData::GetCompressString(), and Nektar::LibUtilities::CompressData::ZlibDecodeFromBase64Str().

◆ v_ReadElements3D()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_ReadElements3D ( )
overrideprotectedvirtual

Look for elements in ELEMENT block.

All elements are of the form: "<? ID="#"> ... </?>", with ? being the element type.

Keep looking

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 809 of file MeshGraphXmlCompressed.cpp.

810 {
811  /// Look for elements in ELEMENT block.
812  TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
813 
814  ASSERTL0(field, "Unable to find ELEMENT tag in file.");
815 
816  /// All elements are of the form: "<? ID="#"> ... </?>", with
817  /// ? being the element type.
818 
819  TiXmlElement *element = field->FirstChildElement();
820 
821  while (element)
822  {
823  std::string elementType(element->ValueStr());
824 
825  // A - tet, P - pyramid, R - prism, H - hex
826  ASSERTL0(
827  elementType == "A" || elementType == "P" || elementType == "R" ||
828  elementType == "H",
829  (std::string("Unknown 3D element type: ") + elementType).c_str());
830 
831  string IsCompressed;
832  element->QueryStringAttribute("COMPRESSED", &IsCompressed);
833 
834  ASSERTL0(
835  boost::iequals(IsCompressed,
837  "Compressed formats do not match. Expected :" +
839  std::string(IsCompressed));
840 
841  // Extract the face body
842  TiXmlNode *child = element->FirstChild();
843  ASSERTL0(child, "Unable to extract the data from "
844  "the compressed face tag.");
845 
846  std::string str;
847  if (child->Type() == TiXmlNode::TINYXML_TEXT)
848  {
849  str += child->ToText()->ValueStr();
850  }
851 
852  int indx;
853  if (elementType == "A")
854  {
855  std::vector<SpatialDomains::MeshTet> data;
857  TriGeomSharedPtr tfaces[4];
858  for (int i = 0; i < data.size(); ++i)
859  {
860  indx = data[i].id;
861  for (int j = 0; j < 4; ++j)
862  {
863  Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
864  tfaces[j] = static_pointer_cast<TriGeom>(face);
865  }
866 
867  TetGeomSharedPtr tetgeom(
869  m_tetGeoms[indx] = tetgeom;
870  PopulateFaceToElMap(tetgeom, 4);
871  }
872  }
873  else if (elementType == "P")
874  {
875  std::vector<SpatialDomains::MeshPyr> data;
877  Geometry2DSharedPtr faces[5];
878  for (int i = 0; i < data.size(); ++i)
879  {
880  indx = data[i].id;
881  int Ntfaces = 0;
882  int Nqfaces = 0;
883  for (int j = 0; j < 5; ++j)
884  {
885  Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
886 
887  if (face == Geometry2DSharedPtr() ||
888  (face->GetShapeType() != LibUtilities::eTriangle &&
889  face->GetShapeType() != LibUtilities::eQuadrilateral))
890  {
891  std::stringstream errorstring;
892  errorstring << "Element " << indx
893  << " has invalid face: " << j;
894  ASSERTL0(false, errorstring.str().c_str());
895  }
896  else if (face->GetShapeType() == LibUtilities::eTriangle)
897  {
898  faces[j] = static_pointer_cast<TriGeom>(face);
899  Ntfaces++;
900  }
901  else if (face->GetShapeType() ==
903  {
904  faces[j] = static_pointer_cast<QuadGeom>(face);
905  Nqfaces++;
906  }
907  }
908  ASSERTL0((Ntfaces == 4) && (Nqfaces == 1),
909  "Did not identify the correct number of "
910  "triangular and quadrilateral faces for a "
911  "pyramid");
912 
913  PyrGeomSharedPtr pyrgeom(
915 
916  m_pyrGeoms[indx] = pyrgeom;
917  PopulateFaceToElMap(pyrgeom, 5);
918  }
919  }
920  else if (elementType == "R")
921  {
922  std::vector<SpatialDomains::MeshPrism> data;
924  Geometry2DSharedPtr faces[5];
925  for (int i = 0; i < data.size(); ++i)
926  {
927  indx = data[i].id;
928  int Ntfaces = 0;
929  int Nqfaces = 0;
930  for (int j = 0; j < 5; ++j)
931  {
932  Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
933  if (face == Geometry2DSharedPtr() ||
934  (face->GetShapeType() != LibUtilities::eTriangle &&
935  face->GetShapeType() != LibUtilities::eQuadrilateral))
936  {
937  std::stringstream errorstring;
938  errorstring << "Element " << indx
939  << " has invalid face: " << j;
940  ASSERTL0(false, errorstring.str().c_str());
941  }
942  else if (face->GetShapeType() == LibUtilities::eTriangle)
943  {
944  faces[j] = static_pointer_cast<TriGeom>(face);
945  Ntfaces++;
946  }
947  else if (face->GetShapeType() ==
949  {
950  faces[j] = static_pointer_cast<QuadGeom>(face);
951  Nqfaces++;
952  }
953  }
954  ASSERTL0((Ntfaces == 2) && (Nqfaces == 3),
955  "Did not identify the correct number of "
956  "triangular and quadrilateral faces for a "
957  "prism");
958 
959  PrismGeomSharedPtr prismgeom(
961 
962  m_prismGeoms[indx] = prismgeom;
963  PopulateFaceToElMap(prismgeom, 5);
964  }
965  }
966  else if (elementType == "H")
967  {
968  std::vector<SpatialDomains::MeshHex> data;
970 
971  QuadGeomSharedPtr faces[6];
972  for (int i = 0; i < data.size(); ++i)
973  {
974  indx = data[i].id;
975  for (int j = 0; j < 6; ++j)
976  {
977  Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
978  faces[j] = static_pointer_cast<QuadGeom>(face);
979  }
980 
981  HexGeomSharedPtr hexgeom(
983  m_hexGeoms[indx] = hexgeom;
984  PopulateFaceToElMap(hexgeom, 6);
985  }
986  }
987  /// Keep looking
988  element = element->NextSiblingElement();
989  }
990 }
void PopulateFaceToElMap(Geometry3DSharedPtr element, int kNfaces)
Given a 3D geometry object #element, populate the face to element map m_faceToElMap which maps faces ...
Definition: MeshGraph.cpp:3833
Geometry2DSharedPtr GetGeometry2D(int gID)
Definition: MeshGraph.h:397
std::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:85
std::shared_ptr< HexGeom > HexGeomSharedPtr
Definition: HexGeom.h:87
std::shared_ptr< PyrGeom > PyrGeomSharedPtr
Definition: PyrGeom.h:77
std::shared_ptr< TetGeom > TetGeomSharedPtr
Definition: TetGeom.h:85
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65

References ASSERTL0, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTriangle, Nektar::LibUtilities::CompressData::GetCompressString(), and Nektar::LibUtilities::CompressData::ZlibDecodeFromBase64Str().

◆ v_ReadFaces()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_ReadFaces ( )
overrideprotectedvirtual

Look for elements in FACE block.

All faces are of the form: "<? ID="#"> ... </?>", with ? being an element type (either Q or T). They might be in compressed format and so then need upacking.

See if this face has curves.

Create a TriGeom to hold the new definition.

See if this face has curves.

Create a QuadGeom to hold the new definition.

Keep looking

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 506 of file MeshGraphXmlCompressed.cpp.

507 {
508  /// Look for elements in FACE block.
509  TiXmlElement *field = m_xmlGeom->FirstChildElement("FACE");
510 
511  ASSERTL0(field, "Unable to find FACE tag in file.");
512 
513  /// All faces are of the form: "<? ID="#"> ... </?>", with
514  /// ? being an element type (either Q or T).
515  /// They might be in compressed format and so then need upacking.
516 
517  TiXmlElement *element = field->FirstChildElement();
518  CurveMap::iterator it;
519 
520  while (element)
521  {
522  std::string elementType(element->ValueStr());
523 
524  ASSERTL0(elementType == "Q" || elementType == "T",
525  (std::string("Unknown 3D face type: ") + elementType).c_str());
526 
527  string IsCompressed;
528  element->QueryStringAttribute("COMPRESSED", &IsCompressed);
529 
530  ASSERTL0(
531  boost::iequals(IsCompressed,
533  "Compressed formats do not match. Expected :" +
535  std::string(IsCompressed));
536 
537  // Extract the face body
538  TiXmlNode *faceChild = element->FirstChild();
539  ASSERTL0(faceChild, "Unable to extract the data from "
540  "the compressed face tag.");
541 
542  std::string faceStr;
543  if (faceChild->Type() == TiXmlNode::TINYXML_TEXT)
544  {
545  faceStr += faceChild->ToText()->ValueStr();
546  }
547 
548  int indx;
549  if (elementType == "T")
550  {
551  std::vector<SpatialDomains::MeshTri> faceData;
553  faceData);
554 
555  for (int i = 0; i < faceData.size(); ++i)
556  {
557  indx = faceData[i].id;
558 
559  /// See if this face has curves.
560  it = m_curvedFaces.find(indx);
561 
562  /// Create a TriGeom to hold the new definition.
564  GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
565  GetSegGeom(faceData[i].e[2])};
566 
567  TriGeomSharedPtr trigeom;
568  if (it == m_curvedFaces.end())
569  {
570  trigeom =
572  }
573  else
574  {
576  indx, edges, it->second);
577  }
578  trigeom->SetGlobalID(indx);
579  m_triGeoms[indx] = trigeom;
580  }
581  }
582  else if (elementType == "Q")
583  {
584  std::vector<SpatialDomains::MeshQuad> faceData;
586  faceData);
587 
588  for (int i = 0; i < faceData.size(); ++i)
589  {
590  indx = faceData[i].id;
591 
592  /// See if this face has curves.
593  it = m_curvedFaces.find(indx);
594 
595  /// Create a QuadGeom to hold the new definition.
597  GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
598  GetSegGeom(faceData[i].e[2]), GetSegGeom(faceData[i].e[3])};
599 
600  QuadGeomSharedPtr quadgeom;
601  if (it == m_curvedFaces.end())
602  {
603  quadgeom =
605  }
606  else
607  {
609  indx, edges, it->second);
610  }
611  quadgeom->SetGlobalID(indx);
612  m_quadGeoms[indx] = quadgeom;
613  }
614  }
615  /// Keep looking
616  element = element->NextSiblingElement();
617  }
618 }

References ASSERTL0, Nektar::LibUtilities::CompressData::GetCompressString(), and Nektar::LibUtilities::CompressData::ZlibDecodeFromBase64Str().

◆ v_ReadVertices()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_ReadVertices ( )
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 69 of file MeshGraphXmlCompressed.cpp.

70 {
71  // Now read the vertices
72  TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
73  ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
74 
75  NekDouble xscale, yscale, zscale;
76 
77  // check to see if any scaling parameters are in
78  // attributes and determine these values
79  LibUtilities::Interpreter expEvaluator;
80  const char *xscal = element->Attribute("XSCALE");
81  if (!xscal)
82  {
83  xscale = 1.0;
84  }
85  else
86  {
87  std::string xscalstr = xscal;
88  int expr_id = expEvaluator.DefineFunction("", xscalstr);
89  xscale = expEvaluator.Evaluate(expr_id);
90  }
91 
92  const char *yscal = element->Attribute("YSCALE");
93  if (!yscal)
94  {
95  yscale = 1.0;
96  }
97  else
98  {
99  std::string yscalstr = yscal;
100  int expr_id = expEvaluator.DefineFunction("", yscalstr);
101  yscale = expEvaluator.Evaluate(expr_id);
102  }
103 
104  const char *zscal = element->Attribute("ZSCALE");
105  if (!zscal)
106  {
107  zscale = 1.0;
108  }
109  else
110  {
111  std::string zscalstr = zscal;
112  int expr_id = expEvaluator.DefineFunction("", zscalstr);
113  zscale = expEvaluator.Evaluate(expr_id);
114  }
115 
116  NekDouble xmove, ymove, zmove;
117 
118  // check to see if any moving parameters are in
119  // attributes and determine these values
120 
121  const char *xmov = element->Attribute("XMOVE");
122  if (!xmov)
123  {
124  xmove = 0.0;
125  }
126  else
127  {
128  std::string xmovstr = xmov;
129  int expr_id = expEvaluator.DefineFunction("", xmovstr);
130  xmove = expEvaluator.Evaluate(expr_id);
131  }
132 
133  const char *ymov = element->Attribute("YMOVE");
134  if (!ymov)
135  {
136  ymove = 0.0;
137  }
138  else
139  {
140  std::string ymovstr = ymov;
141  int expr_id = expEvaluator.DefineFunction("", ymovstr);
142  ymove = expEvaluator.Evaluate(expr_id);
143  }
144 
145  const char *zmov = element->Attribute("ZMOVE");
146  if (!zmov)
147  {
148  zmove = 0.0;
149  }
150  else
151  {
152  std::string zmovstr = zmov;
153  int expr_id = expEvaluator.DefineFunction("", zmovstr);
154  zmove = expEvaluator.Evaluate(expr_id);
155  }
156 
157  string IsCompressed;
158  element->QueryStringAttribute("COMPRESSED", &IsCompressed);
159 
160  if (boost::iequals(IsCompressed,
162  {
163  // Extract the vertex body
164  TiXmlNode *vertexChild = element->FirstChild();
165  ASSERTL0(vertexChild, "Unable to extract the data from the compressed "
166  "vertex tag.");
167 
168  std::string vertexStr;
169  if (vertexChild->Type() == TiXmlNode::TINYXML_TEXT)
170  {
171  vertexStr += vertexChild->ToText()->ValueStr();
172  }
173 
174  std::vector<SpatialDomains::MeshVertex> vertData;
176  vertData);
177 
178  int indx;
179  NekDouble xval, yval, zval;
180  for (int i = 0; i < vertData.size(); ++i)
181  {
182  indx = vertData[i].id;
183  xval = vertData[i].x;
184  yval = vertData[i].y;
185  zval = vertData[i].z;
186 
187  xval = xval * xscale + xmove;
188  yval = yval * yscale + ymove;
189  zval = zval * zscale + zmove;
190 
192  m_spaceDimension, indx, xval, yval, zval));
193 
194  vert->SetGlobalID(indx);
195  m_vertSet[indx] = vert;
196  }
197  }
198  else
199  {
200  ASSERTL0(false, "Compressed formats do not match. Expected :" +
202  " but got " + std::string(IsCompressed));
203  }
204 }

References ASSERTL0, Nektar::LibUtilities::Interpreter::DefineFunction(), Nektar::LibUtilities::Interpreter::Evaluate(), Nektar::LibUtilities::CompressData::GetCompressString(), and Nektar::LibUtilities::CompressData::ZlibDecodeFromBase64Str().

◆ v_WriteCurves()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WriteCurves ( TiXmlElement *  geomTag,
CurveMap edges,
CurveMap faces 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1290 of file MeshGraphXmlCompressed.cpp.

1292 {
1293  if (edges.size() == 0 && faces.size() == 0)
1294  {
1295  return;
1296  }
1297 
1298  TiXmlElement *curveTag = new TiXmlElement("CURVED");
1299 
1300  vector<MeshCurvedInfo> edgeInfo;
1301  vector<MeshCurvedInfo> faceInfo;
1302  MeshCurvedPts curvedPts;
1303  curvedPts.id = 0;
1304  int ptOffset = 0;
1305  int newIdx = 0;
1306  int edgeCnt = 0;
1307  int faceCnt = 0;
1308 
1309  for (auto &i : edges)
1310  {
1311  MeshCurvedInfo cinfo;
1312  cinfo.id = edgeCnt++;
1313  cinfo.entityid = i.first;
1314  cinfo.npoints = i.second->m_points.size();
1315  cinfo.ptype = i.second->m_ptype;
1316  cinfo.ptid = 0;
1317  cinfo.ptoffset = ptOffset;
1318 
1319  edgeInfo.push_back(cinfo);
1320 
1321  for (int j = 0; j < i.second->m_points.size(); j++)
1322  {
1323  MeshVertex v;
1324  v.id = newIdx;
1325  v.x = i.second->m_points[j]->x();
1326  v.y = i.second->m_points[j]->y();
1327  v.z = i.second->m_points[j]->z();
1328  curvedPts.pts.push_back(v);
1329  curvedPts.index.push_back(newIdx);
1330  newIdx++;
1331  }
1332  ptOffset += cinfo.npoints;
1333  }
1334 
1335  for (auto &i : faces)
1336  {
1337  MeshCurvedInfo cinfo;
1338  cinfo.id = faceCnt++;
1339  cinfo.entityid = i.first;
1340  cinfo.npoints = i.second->m_points.size();
1341  cinfo.ptype = i.second->m_ptype;
1342  cinfo.ptid = 0;
1343  cinfo.ptoffset = ptOffset;
1344 
1345  faceInfo.push_back(cinfo);
1346 
1347  for (int j = 0; j < i.second->m_points.size(); j++)
1348  {
1349  MeshVertex v;
1350  v.id = newIdx;
1351  v.x = i.second->m_points[j]->x();
1352  v.y = i.second->m_points[j]->y();
1353  v.z = i.second->m_points[j]->z();
1354  curvedPts.pts.push_back(v);
1355  curvedPts.index.push_back(newIdx);
1356  newIdx++;
1357  }
1358  ptOffset += cinfo.npoints;
1359  }
1360 
1361  curveTag->SetAttribute("COMPRESSED",
1363  curveTag->SetAttribute("BITSIZE",
1365 
1366  if (edgeInfo.size())
1367  {
1368  TiXmlElement *x = new TiXmlElement("E");
1369  string dataStr;
1371 
1372  x->LinkEndChild(new TiXmlText(dataStr));
1373  curveTag->LinkEndChild(x);
1374  }
1375 
1376  if (faceInfo.size())
1377  {
1378  TiXmlElement *x = new TiXmlElement("F");
1379  string dataStr;
1381 
1382  x->LinkEndChild(new TiXmlText(dataStr));
1383  curveTag->LinkEndChild(x);
1384  }
1385 
1386  if (edgeInfo.size() || faceInfo.size())
1387  {
1388  TiXmlElement *x = new TiXmlElement("DATAPOINTS");
1389  x->SetAttribute("ID", curvedPts.id);
1390  TiXmlElement *subx = new TiXmlElement("INDEX");
1391  string dataStr;
1393  dataStr);
1394  subx->LinkEndChild(new TiXmlText(dataStr));
1395  x->LinkEndChild(subx);
1396 
1397  subx = new TiXmlElement("POINTS");
1399  dataStr);
1400  subx->LinkEndChild(new TiXmlText(dataStr));
1401  x->LinkEndChild(subx);
1402  curveTag->LinkEndChild(x);
1403  }
1404 
1405  geomTag->LinkEndChild(curveTag);
1406 }
int ZlibEncodeToBase64Str(std::vector< T > &in, std::string &out64)
Definition: CompressData.h:132

References Nektar::SpatialDomains::MeshCurvedInfo::entityid, Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshVertex::id, Nektar::SpatialDomains::MeshCurvedInfo::id, Nektar::SpatialDomains::MeshCurvedPts::id, Nektar::SpatialDomains::MeshCurvedPts::index, Nektar::SpatialDomains::MeshCurvedInfo::npoints, Nektar::SpatialDomains::MeshCurvedInfo::ptid, Nektar::SpatialDomains::MeshCurvedInfo::ptoffset, Nektar::SpatialDomains::MeshCurvedPts::pts, Nektar::SpatialDomains::MeshCurvedInfo::ptype, Nektar::SpatialDomains::MeshVertex::x, Nektar::SpatialDomains::MeshVertex::y, Nektar::SpatialDomains::MeshVertex::z, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WriteEdges()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WriteEdges ( TiXmlElement *  geomTag,
SegGeomMap edges 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1027 of file MeshGraphXmlCompressed.cpp.

1029 {
1030  if (edges.size() == 0)
1031  {
1032  return;
1033  }
1034 
1035  TiXmlElement *edgeTag =
1036  new TiXmlElement(m_meshDimension == 1 ? "S" : "EDGE");
1037 
1038  vector<MeshEdge> edgeInfo;
1039 
1040  for (auto &i : edges)
1041  {
1042  MeshEdge e;
1043  e.id = i.first;
1044  e.v0 = i.second->GetVid(0);
1045  e.v1 = i.second->GetVid(1);
1046  edgeInfo.push_back(e);
1047  }
1048 
1049  string edgeStr;
1051 
1052  edgeTag->SetAttribute("COMPRESSED",
1054  edgeTag->SetAttribute("BITSIZE",
1056 
1057  edgeTag->LinkEndChild(new TiXmlText(edgeStr));
1058 
1059  if (m_meshDimension == 1)
1060  {
1061  TiXmlElement *tmp = new TiXmlElement("ELEMENT");
1062  tmp->LinkEndChild(edgeTag);
1063  geomTag->LinkEndChild(tmp);
1064  }
1065  else
1066  {
1067  geomTag->LinkEndChild(edgeTag);
1068  }
1069 }

References Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshEdge::id, Nektar::SpatialDomains::MeshEdge::v0, Nektar::SpatialDomains::MeshEdge::v1, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WriteHexs()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WriteHexs ( TiXmlElement *  elmtTag,
HexGeomMap hexs 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1142 of file MeshGraphXmlCompressed.cpp.

1144 {
1145  if (hexs.size() == 0)
1146  {
1147  return;
1148  }
1149 
1150  string tag = "H";
1151 
1152  vector<MeshHex> elementInfo;
1153 
1154  for (auto &i : hexs)
1155  {
1156  MeshHex e;
1157  e.id = i.first;
1158  e.f[0] = i.second->GetFid(0);
1159  e.f[1] = i.second->GetFid(1);
1160  e.f[2] = i.second->GetFid(2);
1161  e.f[3] = i.second->GetFid(3);
1162  e.f[4] = i.second->GetFid(4);
1163  e.f[5] = i.second->GetFid(5);
1164  elementInfo.push_back(e);
1165  }
1166 
1167  TiXmlElement *x = new TiXmlElement(tag);
1168  string elStr;
1170 
1171  x->SetAttribute("COMPRESSED",
1173  x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1174 
1175  x->LinkEndChild(new TiXmlText(elStr));
1176 
1177  elmtTag->LinkEndChild(x);
1178 }

References Nektar::SpatialDomains::MeshHex::f, Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshHex::id, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WritePrisms()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WritePrisms ( TiXmlElement *  elmtTag,
PrismGeomMap pris 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1180 of file MeshGraphXmlCompressed.cpp.

1182 {
1183  if (pris.size() == 0)
1184  {
1185  return;
1186  }
1187 
1188  string tag = "R";
1189 
1190  vector<MeshPrism> elementInfo;
1191 
1192  for (auto &i : pris)
1193  {
1194  MeshPrism e;
1195  e.id = i.first;
1196  e.f[0] = i.second->GetFid(0);
1197  e.f[1] = i.second->GetFid(1);
1198  e.f[2] = i.second->GetFid(2);
1199  e.f[3] = i.second->GetFid(3);
1200  e.f[4] = i.second->GetFid(4);
1201  elementInfo.push_back(e);
1202  }
1203 
1204  TiXmlElement *x = new TiXmlElement(tag);
1205  string elStr;
1207 
1208  x->SetAttribute("COMPRESSED",
1210  x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1211 
1212  x->LinkEndChild(new TiXmlText(elStr));
1213 
1214  elmtTag->LinkEndChild(x);
1215 }

References Nektar::SpatialDomains::MeshPrism::f, Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshPrism::id, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WritePyrs()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WritePyrs ( TiXmlElement *  elmtTag,
PyrGeomMap pyrs 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1217 of file MeshGraphXmlCompressed.cpp.

1219 {
1220  if (pyrs.size() == 0)
1221  {
1222  return;
1223  }
1224 
1225  string tag = "P";
1226 
1227  vector<MeshPyr> elementInfo;
1228 
1229  for (auto &i : pyrs)
1230  {
1231  MeshPyr e;
1232  e.id = i.first;
1233  e.f[0] = i.second->GetFid(0);
1234  e.f[1] = i.second->GetFid(1);
1235  e.f[2] = i.second->GetFid(2);
1236  e.f[3] = i.second->GetFid(3);
1237  e.f[4] = i.second->GetFid(4);
1238  elementInfo.push_back(e);
1239  }
1240 
1241  TiXmlElement *x = new TiXmlElement(tag);
1242  string elStr;
1244 
1245  x->SetAttribute("COMPRESSED",
1247  x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1248 
1249  x->LinkEndChild(new TiXmlText(elStr));
1250 
1251  elmtTag->LinkEndChild(x);
1252 }

References Nektar::SpatialDomains::MeshPyr::f, Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshPyr::id, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WriteQuads()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WriteQuads ( TiXmlElement *  faceTag,
QuadGeomMap quads 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1106 of file MeshGraphXmlCompressed.cpp.

1108 {
1109  if (quads.size() == 0)
1110  {
1111  return;
1112  }
1113 
1114  string tag = "Q";
1115 
1116  vector<MeshQuad> quadInfo;
1117 
1118  for (auto &i : quads)
1119  {
1120  MeshQuad q;
1121  q.id = i.first;
1122  q.e[0] = i.second->GetEid(0);
1123  q.e[1] = i.second->GetEid(1);
1124  q.e[2] = i.second->GetEid(2);
1125  q.e[3] = i.second->GetEid(3);
1126  quadInfo.push_back(q);
1127  }
1128 
1129  TiXmlElement *x = new TiXmlElement(tag);
1130  string quadStr;
1132 
1133  x->SetAttribute("COMPRESSED",
1135  x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1136 
1137  x->LinkEndChild(new TiXmlText(quadStr));
1138 
1139  faceTag->LinkEndChild(x);
1140 }

References Nektar::SpatialDomains::MeshQuad::e, Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshQuad::id, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WriteTets()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WriteTets ( TiXmlElement *  elmtTag,
TetGeomMap tets 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1254 of file MeshGraphXmlCompressed.cpp.

1256 {
1257  if (tets.size() == 0)
1258  {
1259  return;
1260  }
1261 
1262  string tag = "A";
1263 
1264  vector<MeshTet> elementInfo;
1265 
1266  for (auto &i : tets)
1267  {
1268  MeshTet e;
1269  e.id = i.first;
1270  e.f[0] = i.second->GetFid(0);
1271  e.f[1] = i.second->GetFid(1);
1272  e.f[2] = i.second->GetFid(2);
1273  e.f[3] = i.second->GetFid(3);
1274  elementInfo.push_back(e);
1275  }
1276 
1277  TiXmlElement *x = new TiXmlElement(tag);
1278  string elStr;
1280 
1281  x->SetAttribute("COMPRESSED",
1283  x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1284 
1285  x->LinkEndChild(new TiXmlText(elStr));
1286 
1287  elmtTag->LinkEndChild(x);
1288 }

References Nektar::SpatialDomains::MeshTet::f, Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshTet::id, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WriteTris()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WriteTris ( TiXmlElement *  faceTag,
TriGeomMap tris 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 1071 of file MeshGraphXmlCompressed.cpp.

1073 {
1074  if (tris.size() == 0)
1075  {
1076  return;
1077  }
1078 
1079  string tag = "T";
1080 
1081  vector<MeshTri> triInfo;
1082 
1083  for (auto &i : tris)
1084  {
1085  MeshTri t;
1086  t.id = i.first;
1087  t.e[0] = i.second->GetEid(0);
1088  t.e[1] = i.second->GetEid(1);
1089  t.e[2] = i.second->GetEid(2);
1090  triInfo.push_back(t);
1091  }
1092 
1093  TiXmlElement *x = new TiXmlElement(tag);
1094  string triStr;
1096 
1097  x->SetAttribute("COMPRESSED",
1099  x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1100 
1101  x->LinkEndChild(new TiXmlText(triStr));
1102 
1103  faceTag->LinkEndChild(x);
1104 }

References Nektar::SpatialDomains::MeshTri::e, Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshTri::id, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

◆ v_WriteVertices()

void Nektar::SpatialDomains::MeshGraphXmlCompressed::v_WriteVertices ( TiXmlElement *  geomTag,
PointGeomMap verts 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::MeshGraphXml.

Definition at line 992 of file MeshGraphXmlCompressed.cpp.

994 {
995  if (verts.size() == 0)
996  {
997  return;
998  }
999 
1000  TiXmlElement *vertTag = new TiXmlElement("VERTEX");
1001 
1002  vector<MeshVertex> vertInfo;
1003 
1004  for (auto &i : verts)
1005  {
1006  MeshVertex v;
1007  v.id = i.first;
1008  v.x = i.second->x();
1009  v.y = i.second->y();
1010  v.z = i.second->z();
1011  vertInfo.push_back(v);
1012  }
1013 
1014  vertTag->SetAttribute("COMPRESSED",
1016  vertTag->SetAttribute("BITSIZE",
1018 
1019  string vertStr;
1021 
1022  vertTag->LinkEndChild(new TiXmlText(vertStr));
1023 
1024  geomTag->LinkEndChild(vertTag);
1025 }

References Nektar::LibUtilities::CompressData::GetBitSizeStr(), Nektar::LibUtilities::CompressData::GetCompressString(), Nektar::SpatialDomains::MeshVertex::id, Nektar::SpatialDomains::MeshVertex::x, Nektar::SpatialDomains::MeshVertex::y, Nektar::SpatialDomains::MeshVertex::z, and Nektar::LibUtilities::CompressData::ZlibEncodeToBase64Str().

Member Data Documentation

◆ className

std::string Nektar::SpatialDomains::MeshGraphXmlCompressed::className
static
Initial value:
=
"IO with Xml geometry")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
MeshGraphFactory & GetMeshGraphFactory()
Definition: MeshGraph.cpp:74

Definition at line 60 of file MeshGraphXmlCompressed.h.