Nektar++
Loading...
Searching...
No Matches
MeshGraphIOXmlCompressed.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: MeshGraphIOXmlCompressed.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description:
32//
33//
34////////////////////////////////////////////////////////////////////////////////
35
37
41
44
45// These are required for the Write(...) and Import(...) functions.
46#include <boost/algorithm/string/predicate.hpp>
47#include <boost/archive/iterators/base64_from_binary.hpp>
48#include <boost/archive/iterators/binary_from_base64.hpp>
49#include <boost/archive/iterators/transform_width.hpp>
50#include <boost/format.hpp>
51
52#include <tinyxml.h>
53
55{
56
59 "XmlCompressed", MeshGraphIOXmlCompressed::create,
60 "IO with compressed Xml geometry");
61
63{
64 int spaceDimension = m_meshGraph->GetSpaceDimension();
65
66 // Now read the vertices
67 TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
68 ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
69
70 NekDouble xscale, yscale, zscale;
71
72 // check to see if any scaling parameters are in
73 // attributes and determine these values
74 LibUtilities::Interpreter expEvaluator;
75 const char *xscal = element->Attribute("XSCALE");
76 if (!xscal)
77 {
78 xscale = 1.0;
79 }
80 else
81 {
82 std::string xscalstr = xscal;
83 int expr_id = expEvaluator.DefineFunction("", xscalstr);
84 xscale = expEvaluator.Evaluate(expr_id);
85 }
86
87 const char *yscal = element->Attribute("YSCALE");
88 if (!yscal)
89 {
90 yscale = 1.0;
91 }
92 else
93 {
94 std::string yscalstr = yscal;
95 int expr_id = expEvaluator.DefineFunction("", yscalstr);
96 yscale = expEvaluator.Evaluate(expr_id);
97 }
98
99 const char *zscal = element->Attribute("ZSCALE");
100 if (!zscal)
101 {
102 zscale = 1.0;
103 }
104 else
105 {
106 std::string zscalstr = zscal;
107 int expr_id = expEvaluator.DefineFunction("", zscalstr);
108 zscale = expEvaluator.Evaluate(expr_id);
109 }
110
111 NekDouble xmove, ymove, zmove;
112
113 // check to see if any moving parameters are in
114 // attributes and determine these values
115
116 const char *xmov = element->Attribute("XMOVE");
117 if (!xmov)
118 {
119 xmove = 0.0;
120 }
121 else
122 {
123 std::string xmovstr = xmov;
124 int expr_id = expEvaluator.DefineFunction("", xmovstr);
125 xmove = expEvaluator.Evaluate(expr_id);
126 }
127
128 const char *ymov = element->Attribute("YMOVE");
129 if (!ymov)
130 {
131 ymove = 0.0;
132 }
133 else
134 {
135 std::string ymovstr = ymov;
136 int expr_id = expEvaluator.DefineFunction("", ymovstr);
137 ymove = expEvaluator.Evaluate(expr_id);
138 }
139
140 const char *zmov = element->Attribute("ZMOVE");
141 if (!zmov)
142 {
143 zmove = 0.0;
144 }
145 else
146 {
147 std::string zmovstr = zmov;
148 int expr_id = expEvaluator.DefineFunction("", zmovstr);
149 zmove = expEvaluator.Evaluate(expr_id);
150 }
151
152 NekDouble zrotate;
153
154 const char *zrot = element->Attribute("ZROT");
155 if (!zrot)
156 {
157 zrotate = 0.0;
158 }
159 else
160 {
161 std::string zrotstr = zrot;
162 int expr_id = expEvaluator.DefineFunction("", zrotstr);
163 zrotate = expEvaluator.Evaluate(expr_id);
164 }
165
166 std::string IsCompressed;
167 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
168
169 if (boost::iequals(IsCompressed,
171 {
172 // Extract the vertex body
173 TiXmlNode *vertexChild = element->FirstChild();
174 ASSERTL0(vertexChild, "Unable to extract the data from the compressed "
175 "vertex tag.");
176
177 std::string vertexStr;
178 if (vertexChild->Type() == TiXmlNode::TINYXML_TEXT)
179 {
180 vertexStr += vertexChild->ToText()->ValueStr();
181 }
182
183 std::vector<SpatialDomains::MeshVertex> vertData;
185 vertData);
186
187 int indx;
188 NekDouble xval, yval, zval;
189 for (int i = 0; i < vertData.size(); ++i)
190 {
191 indx = vertData[i].id;
192 xval = vertData[i].x;
193 yval = vertData[i].y;
194 zval = vertData[i].z;
195
196 xval = xval * xscale + xmove;
197 yval = yval * yscale + ymove;
198 zval = zval * zscale + zmove;
199
200 if (zrotate != 0.0)
201 {
202 NekDouble xval_tmp = xval * cos(zrotate) - yval * sin(zrotate);
203 yval = xval * sin(zrotate) + yval * cos(zrotate);
204 xval = xval_tmp;
205 }
206
207 m_meshGraph->CreatePointGeom(spaceDimension, indx, xval, yval,
208 zval);
209 }
210 }
211 else
212 {
213 ASSERTL0(false, "Compressed formats do not match. Expected :" +
215 " but got " + IsCompressed);
216 }
217}
218
220{
221 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
222 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
223 auto &curveNodes = m_meshGraph->GetAllCurveNodes();
224 int spaceDimension = m_meshGraph->GetSpaceDimension();
225
226 // check to see if any scaling parameters are in
227 // attributes and determine these values
228 TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
229 ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
230
231 NekDouble xscale, yscale, zscale;
232
233 LibUtilities::Interpreter expEvaluator;
234 const char *xscal = element->Attribute("XSCALE");
235 if (!xscal)
236 {
237 xscale = 1.0;
238 }
239 else
240 {
241 std::string xscalstr = xscal;
242 int expr_id = expEvaluator.DefineFunction("", xscalstr);
243 xscale = expEvaluator.Evaluate(expr_id);
244 }
245
246 const char *yscal = element->Attribute("YSCALE");
247 if (!yscal)
248 {
249 yscale = 1.0;
250 }
251 else
252 {
253 std::string yscalstr = yscal;
254 int expr_id = expEvaluator.DefineFunction("", yscalstr);
255 yscale = expEvaluator.Evaluate(expr_id);
256 }
257
258 const char *zscal = element->Attribute("ZSCALE");
259 if (!zscal)
260 {
261 zscale = 1.0;
262 }
263 else
264 {
265 std::string zscalstr = zscal;
266 int expr_id = expEvaluator.DefineFunction("", zscalstr);
267 zscale = expEvaluator.Evaluate(expr_id);
268 }
269
270 NekDouble xmove, ymove, zmove;
271
272 // check to see if any moving parameters are in
273 // attributes and determine these values
274
275 const char *xmov = element->Attribute("XMOVE");
276 if (!xmov)
277 {
278 xmove = 0.0;
279 }
280 else
281 {
282 std::string xmovstr = xmov;
283 int expr_id = expEvaluator.DefineFunction("", xmovstr);
284 xmove = expEvaluator.Evaluate(expr_id);
285 }
286
287 const char *ymov = element->Attribute("YMOVE");
288 if (!ymov)
289 {
290 ymove = 0.0;
291 }
292 else
293 {
294 std::string ymovstr = ymov;
295 int expr_id = expEvaluator.DefineFunction("", ymovstr);
296 ymove = expEvaluator.Evaluate(expr_id);
297 }
298
299 const char *zmov = element->Attribute("ZMOVE");
300 if (!zmov)
301 {
302 zmove = 0.0;
303 }
304 else
305 {
306 std::string zmovstr = zmov;
307 int expr_id = expEvaluator.DefineFunction("", zmovstr);
308 zmove = expEvaluator.Evaluate(expr_id);
309 }
310
311 NekDouble zrotate;
312
313 const char *zrot = element->Attribute("ZROT");
314 if (!zrot)
315 {
316 zrotate = 0.0;
317 }
318 else
319 {
320 std::string zrotstr = zrot;
321 int expr_id = expEvaluator.DefineFunction("", zrotstr);
322 zrotate = expEvaluator.Evaluate(expr_id);
323 }
324
325 /// Look for elements in CURVE block.
326 TiXmlElement *field = m_xmlGeom->FirstChildElement("CURVED");
327
328 if (!field) // return if no curved entities
329 {
330 return;
331 }
332
333 std::string IsCompressed;
334 field->QueryStringAttribute("COMPRESSED", &IsCompressed);
335
336 if (IsCompressed.size() == 0)
337 {
338 // this could be that the curved tag is empty
339 // in this case we dont want to read it
340 return;
341 }
342
343 ASSERTL0(boost::iequals(IsCompressed,
345 "Compressed formats do not match. Expected :" +
347 IsCompressed);
348
349 std::vector<SpatialDomains::MeshCurvedInfo> edginfo;
350 std::vector<SpatialDomains::MeshCurvedInfo> facinfo;
352
353 // read edge, face info and curved poitns.
354 TiXmlElement *x = field->FirstChildElement();
355 while (x)
356 {
357 const char *entitytype = x->Value();
358 // read in edge or face info
359 if (boost::iequals(entitytype, "E"))
360 {
361 // read in data
362 std::string elmtStr;
363 TiXmlNode *child = x->FirstChild();
364
365 if (child->Type() == TiXmlNode::TINYXML_TEXT)
366 {
367 elmtStr += child->ToText()->ValueStr();
368 }
369
371 edginfo);
372 }
373 else if (boost::iequals(entitytype, "F"))
374 {
375 // read in data
376 std::string elmtStr;
377 TiXmlNode *child = x->FirstChild();
378
379 if (child->Type() == TiXmlNode::TINYXML_TEXT)
380 {
381 elmtStr += child->ToText()->ValueStr();
382 }
383
385 facinfo);
386 }
387 else if (boost::iequals(entitytype, "DATAPOINTS"))
388 {
389 int32_t id;
390 ASSERTL0(x->Attribute("ID", &id),
391 "Failed to get ID from PTS section");
392 cpts.id = id;
393
394 // read in data
395 std::string elmtStr;
396
397 TiXmlElement *DataIdx = x->FirstChildElement("INDEX");
398 ASSERTL0(DataIdx, "Cannot read data index tag in compressed "
399 "curved section");
400
401 TiXmlNode *child = DataIdx->FirstChild();
402 if (child->Type() == TiXmlNode::TINYXML_TEXT)
403 {
404 elmtStr = child->ToText()->ValueStr();
405 }
406
408 cpts.index);
409
410 TiXmlElement *DataPts = x->FirstChildElement("POINTS");
411 ASSERTL0(DataPts, "Cannot read data pts tag in compressed "
412 "curved section");
413
414 child = DataPts->FirstChild();
415 if (child->Type() == TiXmlNode::TINYXML_TEXT)
416 {
417 elmtStr = child->ToText()->ValueStr();
418 }
419
421 cpts.pts);
422 }
423 else
424 {
425 ASSERTL0(false, "Unknown tag in curved section");
426 }
427 x = x->NextSiblingElement();
428 }
429
430 // rescale (x,y,z) points;
431 for (int i = 0; i < cpts.pts.size(); ++i)
432 {
433 cpts.pts[i].x = xscale * cpts.pts[i].x + xmove;
434 cpts.pts[i].y = yscale * cpts.pts[i].y + ymove;
435 cpts.pts[i].z = zscale * cpts.pts[i].z + zmove;
436
437 if (zrotate != 0.0)
438 {
439 NekDouble xval_tmp =
440 cpts.pts[i].x * cos(zrotate) - cpts.pts[i].y * sin(zrotate);
441 cpts.pts[i].y =
442 cpts.pts[i].x * sin(zrotate) + cpts.pts[i].y * cos(zrotate);
443 cpts.pts[i].x = xval_tmp;
444 }
445 }
446
447 for (int i = 0; i < edginfo.size(); ++i)
448 {
449 int edgeid = edginfo[i].entityid;
451
452 curvedEdges[edgeid] = MemoryManager<Curve>::AllocateSharedPtr(
453 edgeid, ptype = (LibUtilities::PointsType)edginfo[i].ptype);
454
455 // load points
456 int offset = edginfo[i].ptoffset;
457 for (int j = 0; j < edginfo[i].npoints; ++j)
458 {
459 int idx = cpts.index[offset + j];
460 curveNodes.emplace_back(
462 spaceDimension, edginfo[i].id, cpts.pts[idx].x,
463 cpts.pts[idx].y, cpts.pts[idx].z));
464 curvedEdges[edgeid]->m_points.emplace_back(curveNodes.back().get());
465 }
466 }
467
468 for (int i = 0; i < facinfo.size(); ++i)
469 {
470 int faceid = facinfo[i].entityid;
472
473 curvedFaces[faceid] = MemoryManager<Curve>::AllocateSharedPtr(
474 faceid, ptype = (LibUtilities::PointsType)facinfo[i].ptype);
475
476 int offset = facinfo[i].ptoffset;
477 for (int j = 0; j < facinfo[i].npoints; ++j)
478 {
479 int idx = cpts.index[offset + j];
480 curveNodes.emplace_back(
482 spaceDimension, facinfo[i].id, cpts.pts[idx].x,
483 cpts.pts[idx].y, cpts.pts[idx].z));
484 curvedFaces[faceid]->m_points.emplace_back(curveNodes.back().get());
485 }
486 }
487}
488
490{
491 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
492 int spaceDimension = m_meshGraph->GetSpaceDimension();
493
494 CurveMap::iterator it;
495
496 /// Look for elements in ELEMENT block.
497 TiXmlElement *field = m_xmlGeom->FirstChildElement("EDGE");
498
499 ASSERTL0(field, "Unable to find EDGE tag in file.");
500
501 std::string IsCompressed;
502 field->QueryStringAttribute("COMPRESSED", &IsCompressed);
503
504 ASSERTL0(boost::iequals(IsCompressed,
506 "Compressed formats do not match. Expected :" +
508 IsCompressed);
509 // Extract the edge body
510 TiXmlNode *edgeChild = field->FirstChild();
511 ASSERTL0(edgeChild, "Unable to extract the data from "
512 "the compressed edge tag.");
513
514 std::string edgeStr;
515 if (edgeChild->Type() == TiXmlNode::TINYXML_TEXT)
516 {
517 edgeStr += edgeChild->ToText()->ValueStr();
518 }
519
520 std::vector<SpatialDomains::MeshEdge> edgeData;
522
523 int indx;
524 for (int i = 0; i < edgeData.size(); ++i)
525 {
526 indx = edgeData[i].id;
527 std::array<PointGeom *, 2> vertices = {
528 m_meshGraph->GetPointGeom(edgeData[i].v0),
529 m_meshGraph->GetPointGeom(edgeData[i].v1)};
530 SegGeomUniquePtr edge;
531 it = curvedEdges.find(indx);
532 if (it == curvedEdges.end())
533 {
534 m_meshGraph->AddGeom(indx,
536 indx, spaceDimension, vertices));
537 }
538 else
539 {
540 m_meshGraph->AddGeom(
542 indx, spaceDimension, vertices, it->second));
543 }
544 }
545}
546
548{
549 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
550
551 /// Look for elements in FACE block.
552 TiXmlElement *field = m_xmlGeom->FirstChildElement("FACE");
553
554 ASSERTL0(field, "Unable to find FACE tag in file.");
555
556 /// All faces are of the form: "<? ID="#"> ... </?>", with
557 /// ? being an element type (either Q or T).
558 /// They might be in compressed format and so then need upacking.
559
560 TiXmlElement *element = field->FirstChildElement();
561 CurveMap::iterator it;
562
563 while (element)
564 {
565 std::string elementType(element->ValueStr());
566
567 ASSERTL0(elementType == "Q" || elementType == "T",
568 (std::string("Unknown 3D face type: ") + elementType).c_str());
569
570 std::string IsCompressed;
571 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
572
573 ASSERTL0(
574 boost::iequals(IsCompressed,
576 "Compressed formats do not match. Expected :" +
578 IsCompressed);
579
580 // Extract the face body
581 TiXmlNode *faceChild = element->FirstChild();
582 ASSERTL0(faceChild, "Unable to extract the data from "
583 "the compressed face tag.");
584
585 std::string faceStr;
586 if (faceChild->Type() == TiXmlNode::TINYXML_TEXT)
587 {
588 faceStr += faceChild->ToText()->ValueStr();
589 }
590
591 int indx;
592 if (elementType == "T")
593 {
594 std::vector<SpatialDomains::MeshTri> faceData;
596 faceData);
597
598 for (int i = 0; i < faceData.size(); ++i)
599 {
600 indx = faceData[i].id;
601
602 /// See if this face has curves.
603 it = curvedFaces.find(indx);
604
605 /// Create a TriGeom to hold the new definition.
606 std::array<SegGeom *, TriGeom::kNedges> edges = {
607 m_meshGraph->GetSegGeom(faceData[i].e[0]),
608 m_meshGraph->GetSegGeom(faceData[i].e[1]),
609 m_meshGraph->GetSegGeom(faceData[i].e[2])};
610
611 if (it == curvedFaces.end())
612 {
613 m_meshGraph->AddGeom(
615 indx, edges));
616 }
617 else
618 {
619 m_meshGraph->AddGeom(
621 indx, edges, it->second));
622 }
623 }
624 }
625 else if (elementType == "Q")
626 {
627 std::vector<SpatialDomains::MeshQuad> faceData;
629 faceData);
630
631 for (int i = 0; i < faceData.size(); ++i)
632 {
633 indx = faceData[i].id;
634
635 /// See if this face has curves.
636 it = curvedFaces.find(indx);
637
638 /// Create a QuadGeom to hold the new definition.
639 std::array<SegGeom *, QuadGeom::kNedges> edges = {
640 m_meshGraph->GetSegGeom(faceData[i].e[0]),
641 m_meshGraph->GetSegGeom(faceData[i].e[1]),
642 m_meshGraph->GetSegGeom(faceData[i].e[2]),
643 m_meshGraph->GetSegGeom(faceData[i].e[3])};
644
645 if (it == curvedFaces.end())
646 {
647 m_meshGraph->AddGeom(
649 indx, edges));
650 }
651 else
652 {
653 m_meshGraph->AddGeom(
655 indx, edges, it->second));
656 }
657 }
658 }
659 /// Keep looking
660 element = element->NextSiblingElement();
661 }
662}
663
665{
666 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
667 int spaceDimension = m_meshGraph->GetSpaceDimension();
668
669 TiXmlElement *field = nullptr;
670
671 /// Look for elements in ELEMENT block.
672 field = m_xmlGeom->FirstChildElement("ELEMENT");
673
674 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
675
676 /// All elements are of the form: "<S ID = n> ... </S>", with
677 /// ? being the element type.
678
679 TiXmlElement *segment = field->FirstChildElement("S");
680 CurveMap::iterator it;
681
682 while (segment)
683 {
684 std::string IsCompressed;
685 segment->QueryStringAttribute("COMPRESSED", &IsCompressed);
686 ASSERTL0(
687 boost::iequals(IsCompressed,
689 "Compressed formats do not match. Expected :" +
691 IsCompressed);
692
693 // Extract the face body
694 TiXmlNode *child = segment->FirstChild();
695 ASSERTL0(child, "Unable to extract the data from "
696 "the compressed face tag.");
697
698 std::string str;
699 if (child->Type() == TiXmlNode::TINYXML_TEXT)
700 {
701 str += child->ToText()->ValueStr();
702 }
703
704 int indx;
705
706 std::vector<SpatialDomains::MeshEdge> data;
708
709 for (int i = 0; i < data.size(); ++i)
710 {
711 indx = data[i].id;
712
713 /// See if this face has curves.
714 it = curvedEdges.find(indx);
715
716 std::array<PointGeom *, 2> vertices = {
717 m_meshGraph->GetPointGeom(data[i].v0),
718 m_meshGraph->GetPointGeom(data[i].v1)};
720
721 if (it == curvedEdges.end())
722 {
723 m_meshGraph->AddGeom(indx,
725 indx, spaceDimension, vertices));
726 }
727 else
728 {
729 m_meshGraph->AddGeom(
731 indx, spaceDimension, vertices, it->second));
732 }
733 }
734 /// Keep looking for additional segments
735 segment = segment->NextSiblingElement("S");
736 }
737}
738
740{
741 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
742
743 /// Look for elements in ELEMENT block.
744 TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
745
746 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
747
748 // Set up curve map for curved elements on an embedded manifold.
749 CurveMap::iterator it;
750
751 /// All elements are of the form: "<? ID="#"> ... </?>", with
752 /// ? being the element type.
753
754 TiXmlElement *element = field->FirstChildElement();
755
756 while (element)
757 {
758 std::string elementType(element->ValueStr());
759
760 ASSERTL0(
761 elementType == "Q" || elementType == "T",
762 (std::string("Unknown 2D element type: ") + elementType).c_str());
763
764 std::string IsCompressed;
765 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
766
767 ASSERTL0(
768 boost::iequals(IsCompressed,
770 "Compressed formats do not match. Expected :" +
772 IsCompressed);
773
774 // Extract the face body
775 TiXmlNode *faceChild = element->FirstChild();
776 ASSERTL0(faceChild, "Unable to extract the data from "
777 "the compressed face tag.");
778
779 std::string faceStr;
780 if (faceChild->Type() == TiXmlNode::TINYXML_TEXT)
781 {
782 faceStr += faceChild->ToText()->ValueStr();
783 }
784
785 int indx;
786 if (elementType == "T")
787 {
788 std::vector<SpatialDomains::MeshTri> faceData;
790 faceData);
791
792 for (int i = 0; i < faceData.size(); ++i)
793 {
794 indx = faceData[i].id;
795
796 /// See if this face has curves.
797 it = curvedFaces.find(indx);
798
799 /// Create a TriGeom to hold the new definition.
800 std::array<SegGeom *, TriGeom::kNedges> edges = {
801 m_meshGraph->GetSegGeom(faceData[i].e[0]),
802 m_meshGraph->GetSegGeom(faceData[i].e[1]),
803 m_meshGraph->GetSegGeom(faceData[i].e[2])};
804
805 if (it == curvedFaces.end())
806 {
807 m_meshGraph->AddGeom(
809 indx, edges));
810 }
811 else
812 {
813 m_meshGraph->AddGeom(
815 indx, edges, it->second));
816 }
817 }
818 }
819 else if (elementType == "Q")
820 {
821 std::vector<SpatialDomains::MeshQuad> faceData;
823 faceData);
824
825 for (int i = 0; i < faceData.size(); ++i)
826 {
827 indx = faceData[i].id;
828
829 /// See if this face has curves.
830 it = curvedFaces.find(indx);
831
832 /// Create a QuadGeom to hold the new definition.
833 std::array<SegGeom *, QuadGeom::kNedges> edges = {
834 m_meshGraph->GetSegGeom(faceData[i].e[0]),
835 m_meshGraph->GetSegGeom(faceData[i].e[1]),
836 m_meshGraph->GetSegGeom(faceData[i].e[2]),
837 m_meshGraph->GetSegGeom(faceData[i].e[3])};
838
839 if (it == curvedFaces.end())
840 {
841 m_meshGraph->AddGeom(
843 indx, edges));
844 }
845 else
846 {
847 m_meshGraph->AddGeom(
849 indx, edges, it->second));
850 }
851 }
852 }
853 /// Keep looking
854 element = element->NextSiblingElement();
855 }
856}
857
859{
860 /// Look for elements in ELEMENT block.
861 TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
862
863 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
864
865 /// All elements are of the form: "<? ID="#"> ... </?>", with
866 /// ? being the element type.
867
868 TiXmlElement *element = field->FirstChildElement();
869
870 while (element)
871 {
872 std::string elementType(element->ValueStr());
873
874 // A - tet, P - pyramid, R - prism, H - hex
875 ASSERTL0(
876 elementType == "A" || elementType == "P" || elementType == "R" ||
877 elementType == "H",
878 (std::string("Unknown 3D element type: ") + elementType).c_str());
879
880 std::string IsCompressed;
881 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
882
883 ASSERTL0(
884 boost::iequals(IsCompressed,
886 "Compressed formats do not match. Expected :" +
888 IsCompressed);
889
890 // Extract the face body
891 TiXmlNode *child = element->FirstChild();
892 ASSERTL0(child, "Unable to extract the data from "
893 "the compressed face tag.");
894
895 std::string str;
896 if (child->Type() == TiXmlNode::TINYXML_TEXT)
897 {
898 str += child->ToText()->ValueStr();
899 }
900
901 int indx;
902 if (elementType == "A")
903 {
904 std::vector<SpatialDomains::MeshTet> data;
906 std::array<TriGeom *, 4> tfaces;
907 for (int i = 0; i < data.size(); ++i)
908 {
909 indx = data[i].id;
910 for (int j = 0; j < 4; ++j)
911 {
912 Geometry2D *face = m_meshGraph->GetGeometry2D(data[i].f[j]);
913 tfaces[j] = static_cast<TriGeom *>(face);
914 }
915
916 auto tetGeom =
918 m_meshGraph->PopulateFaceToElMap(tetGeom.get(), 4);
919 m_meshGraph->AddGeom(indx, std::move(tetGeom));
920 }
921 }
922 else if (elementType == "P")
923 {
924 std::vector<SpatialDomains::MeshPyr> data;
926 std::array<Geometry2D *, 5> faces;
927 for (int i = 0; i < data.size(); ++i)
928 {
929 indx = data[i].id;
930 int Ntfaces = 0;
931 int Nqfaces = 0;
932 for (int j = 0; j < 5; ++j)
933 {
934 Geometry2D *face = m_meshGraph->GetGeometry2D(data[i].f[j]);
935
936 if (face == nullptr ||
939 {
940 std::stringstream errorstring;
941 errorstring << "Element " << indx
942 << " has invalid face: " << j;
943 ASSERTL0(false, errorstring.str().c_str());
944 }
945 else if (face->GetShapeType() == LibUtilities::eTriangle)
946 {
947 faces[j] = static_cast<TriGeom *>(face);
948 Ntfaces++;
949 }
950 else if (face->GetShapeType() ==
952 {
953 faces[j] = static_cast<QuadGeom *>(face);
954 Nqfaces++;
955 }
956 }
957 ASSERTL0((Ntfaces == 4) && (Nqfaces == 1),
958 "Did not identify the correct number of "
959 "triangular and quadrilateral faces for a "
960 "pyramid");
961
962 auto pyrGeom =
964 m_meshGraph->PopulateFaceToElMap(pyrGeom.get(), 5);
965 m_meshGraph->AddGeom(indx, std::move(pyrGeom));
966 }
967 }
968 else if (elementType == "R")
969 {
970 std::vector<SpatialDomains::MeshPrism> data;
972 std::array<Geometry2D *, 5> faces;
973 for (int i = 0; i < data.size(); ++i)
974 {
975 indx = data[i].id;
976 int Ntfaces = 0;
977 int Nqfaces = 0;
978 for (int j = 0; j < 5; ++j)
979 {
980 Geometry2D *face = m_meshGraph->GetGeometry2D(data[i].f[j]);
981 if (face == nullptr ||
984 {
985 std::stringstream errorstring;
986 errorstring << "Element " << indx
987 << " has invalid face: " << j;
988 ASSERTL0(false, errorstring.str().c_str());
989 }
990 else if (face->GetShapeType() == LibUtilities::eTriangle)
991 {
992 faces[j] = static_cast<TriGeom *>(face);
993 Ntfaces++;
994 }
995 else if (face->GetShapeType() ==
997 {
998 faces[j] = static_cast<QuadGeom *>(face);
999 Nqfaces++;
1000 }
1001 }
1002 ASSERTL0((Ntfaces == 2) && (Nqfaces == 3),
1003 "Did not identify the correct number of "
1004 "triangular and quadrilateral faces for a "
1005 "prism");
1006
1007 auto prismGeom =
1009 m_meshGraph->PopulateFaceToElMap(prismGeom.get(), 5);
1010 m_meshGraph->AddGeom(indx, std::move(prismGeom));
1011 }
1012 }
1013 else if (elementType == "H")
1014 {
1015 std::vector<SpatialDomains::MeshHex> data;
1017
1018 std::array<QuadGeom *, 6> faces;
1019 for (int i = 0; i < data.size(); ++i)
1020 {
1021 indx = data[i].id;
1022 for (int j = 0; j < 6; ++j)
1023 {
1024 Geometry2D *face = m_meshGraph->GetGeometry2D(data[i].f[j]);
1025 faces[j] = static_cast<QuadGeom *>(face);
1026 }
1027
1028 auto hexGeom =
1030 m_meshGraph->PopulateFaceToElMap(hexGeom.get(), 6);
1031 m_meshGraph->AddGeom(indx, std::move(hexGeom));
1032 }
1033 }
1034 /// Keep looking
1035 element = element->NextSiblingElement();
1036 }
1037}
1038
1039void WriteVert(PointGeom *vert, std::vector<MeshVertex> &vertInfo, int vertID)
1040{
1041 MeshVertex v;
1042 v.id = vertID;
1043 v.x = vert->x();
1044 v.y = vert->y();
1045 v.z = vert->z();
1046 vertInfo.push_back(v);
1047}
1049 std::vector<int> keysToWrite)
1050{
1051 auto &verts = m_meshGraph->GetGeomMap<PointGeom>();
1052 if (m_meshGraph->GetGeomMap<PointGeom>().size() == 0)
1053 {
1054 return;
1055 }
1056
1057 TiXmlElement *vertTag = new TiXmlElement("VERTEX");
1058
1059 std::vector<MeshVertex> vertInfo;
1060
1061 if (keysToWrite.empty())
1062 {
1063 for (auto [id, vert] : verts)
1064 {
1065 WriteVert(vert, vertInfo, id);
1066 }
1067 }
1068 else
1069 {
1070 for (int id : keysToWrite)
1071 {
1072 WriteVert(verts.at(id), vertInfo, id);
1073 }
1074 }
1075
1076 vertTag->SetAttribute("COMPRESSED",
1078 vertTag->SetAttribute("BITSIZE",
1080
1081 std::string vertStr;
1083
1084 vertTag->LinkEndChild(new TiXmlText(vertStr));
1085
1086 geomTag->LinkEndChild(vertTag);
1087}
1088
1089void WriteEdge(SegGeom *seg, std::vector<MeshEdge> &edgeInfo, int edgeID)
1090{
1091 MeshEdge e;
1092 e.id = edgeID;
1093 e.v0 = seg->GetVid(0);
1094 e.v1 = seg->GetVid(1);
1095 edgeInfo.push_back(e);
1096}
1098 std::vector<int> keysToWrite)
1099{
1100 auto &edges = m_meshGraph->GetGeomMap<SegGeom>();
1101 if (m_meshGraph->GetGeomMap<SegGeom>().size() == 0)
1102 {
1103 return;
1104 }
1105
1106 int meshDimension = m_meshGraph->GetMeshDimension();
1107
1108 TiXmlElement *edgeTag = new TiXmlElement(meshDimension == 1 ? "S" : "EDGE");
1109
1110 std::vector<MeshEdge> edgeInfo;
1111
1112 if (keysToWrite.empty())
1113 {
1114 for (auto [id, edge] : edges)
1115 {
1116 WriteEdge(edge, edgeInfo, id);
1117 }
1118 }
1119 else
1120 {
1121 for (int id : keysToWrite)
1122 {
1123 WriteEdge(edges.at(id), edgeInfo, id);
1124 }
1125 }
1126
1127 std::string edgeStr;
1129
1130 edgeTag->SetAttribute("COMPRESSED",
1132 edgeTag->SetAttribute("BITSIZE",
1134
1135 edgeTag->LinkEndChild(new TiXmlText(edgeStr));
1136
1137 if (meshDimension == 1)
1138 {
1139 TiXmlElement *tmp = new TiXmlElement("ELEMENT");
1140 tmp->LinkEndChild(edgeTag);
1141 geomTag->LinkEndChild(tmp);
1142 }
1143 else
1144 {
1145 geomTag->LinkEndChild(edgeTag);
1146 }
1147}
1148
1149void WriteTri(TriGeom *tri, std::vector<MeshTri> &triInfo, int triID)
1150{
1151 MeshTri t;
1152 t.id = triID;
1153 t.e[0] = tri->GetEid(0);
1154 t.e[1] = tri->GetEid(1);
1155 t.e[2] = tri->GetEid(2);
1156 triInfo.push_back(t);
1157}
1158void MeshGraphIOXmlCompressed::v_WriteTris(TiXmlElement *faceTag,
1159 std::vector<int> keysToWrite)
1160{
1161 auto &tris = m_meshGraph->GetGeomMap<TriGeom>();
1162
1163 if (tris.size() == 0)
1164 {
1165 return;
1166 }
1167
1168 std::string tag = "T";
1169
1170 std::vector<MeshTri> triInfo;
1171
1172 if (keysToWrite.empty())
1173 {
1174 for (auto [id, tri] : tris)
1175 {
1176 WriteTri(tri, triInfo, id);
1177 }
1178 }
1179 else
1180 {
1181 for (int id : keysToWrite)
1182 {
1183 WriteTri(tris.at(id), triInfo, id);
1184 }
1185 }
1186
1187 TiXmlElement *x = new TiXmlElement(tag);
1188 std::string triStr;
1190
1191 x->SetAttribute("COMPRESSED",
1193 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1194
1195 x->LinkEndChild(new TiXmlText(triStr));
1196
1197 faceTag->LinkEndChild(x);
1198}
1199
1200void WriteQuad(QuadGeom *quad, std::vector<MeshQuad> &quadInfo, int quadID)
1201{
1202 MeshQuad q;
1203 q.id = quadID;
1204 q.e[0] = quad->GetEid(0);
1205 q.e[1] = quad->GetEid(1);
1206 q.e[2] = quad->GetEid(2);
1207 q.e[3] = quad->GetEid(3);
1208 quadInfo.push_back(q);
1209}
1211 std::vector<int> keysToWrite)
1212{
1213 auto &quads = m_meshGraph->GetGeomMap<QuadGeom>();
1214
1215 if (quads.size() == 0)
1216 {
1217 return;
1218 }
1219
1220 std::string tag = "Q";
1221
1222 std::vector<MeshQuad> quadInfo;
1223
1224 if (keysToWrite.empty())
1225 {
1226 for (auto [id, quad] : quads)
1227 {
1228 WriteQuad(quad, quadInfo, id);
1229 }
1230 }
1231 else
1232 {
1233 for (int id : keysToWrite)
1234 {
1235 WriteQuad(quads.at(id), quadInfo, id);
1236 }
1237 }
1238
1239 TiXmlElement *x = new TiXmlElement(tag);
1240 std::string quadStr;
1242
1243 x->SetAttribute("COMPRESSED",
1245 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1246
1247 x->LinkEndChild(new TiXmlText(quadStr));
1248
1249 faceTag->LinkEndChild(x);
1250}
1251
1252void WriteHex(HexGeom *hex, std::vector<MeshHex> &elementInfo, int hexID)
1253{
1254 MeshHex e;
1255 e.id = hexID;
1256 e.f[0] = hex->GetFid(0);
1257 e.f[1] = hex->GetFid(1);
1258 e.f[2] = hex->GetFid(2);
1259 e.f[3] = hex->GetFid(3);
1260 e.f[4] = hex->GetFid(4);
1261 e.f[5] = hex->GetFid(5);
1262 elementInfo.push_back(e);
1263}
1264void MeshGraphIOXmlCompressed::v_WriteHexs(TiXmlElement *elmtTag,
1265 std::vector<int> keysToWrite)
1266{
1267 auto &hexes = m_meshGraph->GetGeomMap<HexGeom>();
1268
1269 if (hexes.size() == 0)
1270 {
1271 return;
1272 }
1273
1274 std::string tag = "H";
1275
1276 std::vector<MeshHex> elementInfo;
1277
1278 if (keysToWrite.empty())
1279 {
1280 for (auto [id, hex] : hexes)
1281 {
1282 WriteHex(hex, elementInfo, id);
1283 }
1284 }
1285 else
1286 {
1287 for (int id : keysToWrite)
1288 {
1289 WriteHex(hexes.at(id), elementInfo, id);
1290 }
1291 }
1292
1293 TiXmlElement *x = new TiXmlElement(tag);
1294 std::string elStr;
1296
1297 x->SetAttribute("COMPRESSED",
1299 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1300
1301 x->LinkEndChild(new TiXmlText(elStr));
1302
1303 elmtTag->LinkEndChild(x);
1304}
1305
1306void WritePrism(PrismGeom *prism, std::vector<MeshPrism> &elementInfo,
1307 int prismID)
1308{
1309 MeshPrism e;
1310 e.id = prismID;
1311 e.f[0] = prism->GetFid(0);
1312 e.f[1] = prism->GetFid(1);
1313 e.f[2] = prism->GetFid(2);
1314 e.f[3] = prism->GetFid(3);
1315 e.f[4] = prism->GetFid(4);
1316 elementInfo.push_back(e);
1317}
1319 std::vector<int> keysToWrite)
1320{
1321 auto &prisms = m_meshGraph->GetGeomMap<PrismGeom>();
1322
1323 if (prisms.size() == 0)
1324 {
1325 return;
1326 }
1327
1328 std::string tag = "R";
1329
1330 std::vector<MeshPrism> elementInfo;
1331
1332 if (keysToWrite.empty())
1333 {
1334 for (auto [id, prism] : prisms)
1335 {
1336 WritePrism(prism, elementInfo, id);
1337 }
1338 }
1339 else
1340 {
1341 for (int id : keysToWrite)
1342 {
1343 WritePrism(prisms.at(id), elementInfo, id);
1344 }
1345 }
1346
1347 TiXmlElement *x = new TiXmlElement(tag);
1348 std::string elStr;
1350
1351 x->SetAttribute("COMPRESSED",
1353 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1354
1355 x->LinkEndChild(new TiXmlText(elStr));
1356
1357 elmtTag->LinkEndChild(x);
1358}
1359
1360void WritePyr(PyrGeom *pyr, std::vector<MeshPyr> &elementInfo, int pyrID)
1361{
1362 MeshPyr e;
1363 e.id = pyrID;
1364 e.f[0] = pyr->GetFid(0);
1365 e.f[1] = pyr->GetFid(1);
1366 e.f[2] = pyr->GetFid(2);
1367 e.f[3] = pyr->GetFid(3);
1368 e.f[4] = pyr->GetFid(4);
1369 elementInfo.push_back(e);
1370}
1371void MeshGraphIOXmlCompressed::v_WritePyrs(TiXmlElement *elmtTag,
1372 std::vector<int> keysToWrite)
1373{
1374 auto &pyrs = m_meshGraph->GetGeomMap<PyrGeom>();
1375
1376 if (pyrs.size() == 0)
1377 {
1378 return;
1379 }
1380
1381 std::string tag = "P";
1382
1383 std::vector<MeshPyr> elementInfo;
1384
1385 if (keysToWrite.empty())
1386 {
1387 for (auto [id, pyr] : pyrs)
1388 {
1389 WritePyr(pyr, elementInfo, id);
1390 }
1391 }
1392 else
1393 {
1394 for (int id : keysToWrite)
1395 {
1396 WritePyr(pyrs.at(id), elementInfo, id);
1397 }
1398 }
1399
1400 TiXmlElement *x = new TiXmlElement(tag);
1401 std::string elStr;
1403
1404 x->SetAttribute("COMPRESSED",
1406 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1407
1408 x->LinkEndChild(new TiXmlText(elStr));
1409
1410 elmtTag->LinkEndChild(x);
1411}
1412
1413void WriteTet(TetGeom *tet, std::vector<MeshTet> &elementInfo, int tetID)
1414{
1415 MeshTet e;
1416 e.id = tetID;
1417 e.f[0] = tet->GetFid(0);
1418 e.f[1] = tet->GetFid(1);
1419 e.f[2] = tet->GetFid(2);
1420 e.f[3] = tet->GetFid(3);
1421 elementInfo.push_back(e);
1422}
1423void MeshGraphIOXmlCompressed::v_WriteTets(TiXmlElement *elmtTag,
1424 std::vector<int> keysToWrite)
1425{
1426 auto &tets = m_meshGraph->GetGeomMap<TetGeom>();
1427
1428 if (tets.size() == 0)
1429 {
1430 return;
1431 }
1432
1433 std::string tag = "A";
1434
1435 std::vector<MeshTet> elementInfo;
1436
1437 if (keysToWrite.empty())
1438 {
1439 for (auto [id, tet] : tets)
1440 {
1441 WriteTet(tet, elementInfo, id);
1442 }
1443 }
1444 else
1445 {
1446 for (int id : keysToWrite)
1447 {
1448 WriteTet(tets.at(id), elementInfo, id);
1449 }
1450 }
1451
1452 TiXmlElement *x = new TiXmlElement(tag);
1453 std::string elStr;
1455
1456 x->SetAttribute("COMPRESSED",
1458 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1459
1460 x->LinkEndChild(new TiXmlText(elStr));
1461
1462 elmtTag->LinkEndChild(x);
1463}
1464
1466 CurveMap &edges, CurveMap &faces)
1467{
1468 if (edges.size() == 0 && faces.size() == 0)
1469 {
1470 return;
1471 }
1472
1473 TiXmlElement *curveTag = new TiXmlElement("CURVED");
1474
1475 std::vector<MeshCurvedInfo> edgeInfo;
1476 std::vector<MeshCurvedInfo> faceInfo;
1477 MeshCurvedPts curvedPts;
1478 curvedPts.id = 0;
1479 int ptOffset = 0;
1480 int newIdx = 0;
1481 int edgeCnt = 0;
1482 int faceCnt = 0;
1483
1484 for (auto &i : edges)
1485 {
1486 MeshCurvedInfo cinfo;
1487 cinfo.id = edgeCnt++;
1488 cinfo.entityid = i.first;
1489 cinfo.npoints = i.second->m_points.size();
1490 cinfo.ptype = i.second->m_ptype;
1491 cinfo.ptid = 0;
1492 cinfo.ptoffset = ptOffset;
1493
1494 edgeInfo.push_back(cinfo);
1495
1496 for (int j = 0; j < i.second->m_points.size(); j++)
1497 {
1498 MeshVertex v;
1499 v.id = newIdx;
1500 v.x = i.second->m_points[j]->x();
1501 v.y = i.second->m_points[j]->y();
1502 v.z = i.second->m_points[j]->z();
1503 curvedPts.pts.push_back(v);
1504 curvedPts.index.push_back(newIdx);
1505 newIdx++;
1506 }
1507 ptOffset += cinfo.npoints;
1508 }
1509
1510 for (auto &i : faces)
1511 {
1512 MeshCurvedInfo cinfo;
1513 cinfo.id = faceCnt++;
1514 cinfo.entityid = i.first;
1515 cinfo.npoints = i.second->m_points.size();
1516 cinfo.ptype = i.second->m_ptype;
1517 cinfo.ptid = 0;
1518 cinfo.ptoffset = ptOffset;
1519
1520 faceInfo.push_back(cinfo);
1521
1522 for (int j = 0; j < i.second->m_points.size(); j++)
1523 {
1524 MeshVertex v;
1525 v.id = newIdx;
1526 v.x = i.second->m_points[j]->x();
1527 v.y = i.second->m_points[j]->y();
1528 v.z = i.second->m_points[j]->z();
1529 curvedPts.pts.push_back(v);
1530 curvedPts.index.push_back(newIdx);
1531 newIdx++;
1532 }
1533 ptOffset += cinfo.npoints;
1534 }
1535
1536 curveTag->SetAttribute("COMPRESSED",
1538 curveTag->SetAttribute("BITSIZE",
1540
1541 if (edgeInfo.size())
1542 {
1543 TiXmlElement *x = new TiXmlElement("E");
1544 std::string dataStr;
1546
1547 x->LinkEndChild(new TiXmlText(dataStr));
1548 curveTag->LinkEndChild(x);
1549 }
1550
1551 if (faceInfo.size())
1552 {
1553 TiXmlElement *x = new TiXmlElement("F");
1554 std::string dataStr;
1556
1557 x->LinkEndChild(new TiXmlText(dataStr));
1558 curveTag->LinkEndChild(x);
1559 }
1560
1561 if (edgeInfo.size() || faceInfo.size())
1562 {
1563 TiXmlElement *x = new TiXmlElement("DATAPOINTS");
1564 x->SetAttribute("ID", curvedPts.id);
1565 TiXmlElement *subx = new TiXmlElement("INDEX");
1566 std::string dataStr;
1568 dataStr);
1569 subx->LinkEndChild(new TiXmlText(dataStr));
1570 x->LinkEndChild(subx);
1571
1572 subx = new TiXmlElement("POINTS");
1574 dataStr);
1575 subx->LinkEndChild(new TiXmlText(dataStr));
1576 x->LinkEndChild(subx);
1577 curveTag->LinkEndChild(x);
1578 }
1579
1580 geomTag->LinkEndChild(curveTag);
1581}
1582} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
Interpreter class for the evaluation of mathematical expressions.
Definition Interpreter.h:76
int DefineFunction(const std::string &vlist, const std::string &expr)
Defines a function for the purposes of evaluation.
NekDouble Evaluate(const int id)
Evaluate a function which depends only on constants and/or parameters.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
boost::call_traits< DataType >::const_reference x() const
Definition NekPoint.hpp:160
boost::call_traits< DataType >::const_reference z() const
Definition NekPoint.hpp:172
boost::call_traits< DataType >::const_reference y() const
Definition NekPoint.hpp:166
Generic object pool allocator/deallocator.
static std::unique_ptr< DataType, UniquePtrDeleter > AllocateUniquePtr(const Args &...args)
2D geometry information
Definition Geometry2D.h:50
LibUtilities::ShapeType GetShapeType(void)
Get the geometric shape type of this object.
Definition Geometry.h:314
int GetVid(int i) const
Returns global id of vertex i of this object.
Definition Geometry.h:353
int GetFid(int i) const
Get the ID of face i of this object.
Definition Geometry.cpp:118
int GetEid(int i) const
Get the ID of edge i of this object.
Definition Geometry.cpp:110
void v_WriteTris(TiXmlElement *faceTag, std::vector< int > keysToWrite=std::vector< int >()) override
void v_WritePyrs(TiXmlElement *elmtTag, std::vector< int > keysToWrite=std::vector< int >()) override
void v_WriteVertices(TiXmlElement *geomTag, std::vector< int > keysToWrite=std::vector< int >()) override
void v_WritePrisms(TiXmlElement *elmtTag, std::vector< int > keysToWrite=std::vector< int >()) override
void v_WriteHexs(TiXmlElement *elmtTag, std::vector< int > keysToWrite=std::vector< int >()) override
void v_WriteCurves(TiXmlElement *geomTag, CurveMap &edges, CurveMap &faces) override
void v_WriteTets(TiXmlElement *elmtTag, std::vector< int > keysToWrite=std::vector< int >()) override
void v_WriteEdges(TiXmlElement *geomTag, std::vector< int > keysToWrite=std::vector< int >()) override
void v_WriteQuads(TiXmlElement *faceTag, std::vector< int > keysToWrite=std::vector< int >()) override
int ZlibDecodeFromBase64Str(std::string &in64, std::vector< T > &out)
int ZlibEncodeToBase64Str(std::vector< T > &in, std::string &out64)
void WriteVert(PointGeom *vert, TiXmlElement *vertTag)
void WritePrism(PrismGeom *pri, TiXmlElement *elmtTag, std::string &tag, int priID)
void WriteTri(TriGeom *tri, TiXmlElement *faceTag, std::string &tag, int triID)
void WriteEdge(SegGeom *seg, TiXmlElement *edgeTag, std::string &tag, int edgeID)
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition Curve.hpp:61
void WriteQuad(QuadGeom *quad, TiXmlElement *faceTag, std::string &tag, int quadID)
void WriteHex(HexGeom *hex, TiXmlElement *elmtTag, std::string &tag, int hexID)
MeshGraphIOFactory & GetMeshGraphIOFactory()
void WritePyr(PyrGeom *pyr, TiXmlElement *elmtTag, std::string &tag, int pyrID)
unique_ptr_objpool< SegGeom > SegGeomUniquePtr
Definition MeshGraph.h:98
void WriteTet(TetGeom *tet, TiXmlElement *elmtTag, std::string &tag, int tetID)
std::int32_t int32_t
int64_t ptype
point offset of data entry for this curve
int64_t npoints
The entity id corresponding to the global edge/curve.
int64_t ptoffset
the id of point data map (currently always 0 since we are using just one set).
int64_t entityid
Id of this curved information.
int64_t ptid
The number of points in this curved entity.
std::vector< int64_t > index
Mapping to access the pts value. Given a 'ptoffset' value the npoints subsquent values provide the in...
std::vector< MeshVertex > pts
mapping to access pts value.