Nektar++
MeshGraphXmlCompressed.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: MeshGraphXmlCompressed.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#include <boost/iostreams/copy.hpp>
52#include <boost/iostreams/filter/zlib.hpp>
53#include <boost/iostreams/filtering_stream.hpp>
54#include <boost/make_shared.hpp>
55
56#include <tinyxml.h>
57using namespace std;
58
60{
61
64 "XmlCompressed", MeshGraphXmlCompressed::create,
65 "IO with Xml geometry");
66
68{
69 // Now read the vertices
70 TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
71 ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
72
73 NekDouble xscale, yscale, zscale;
74
75 // check to see if any scaling parameters are in
76 // attributes and determine these values
77 LibUtilities::Interpreter expEvaluator;
78 const char *xscal = element->Attribute("XSCALE");
79 if (!xscal)
80 {
81 xscale = 1.0;
82 }
83 else
84 {
85 std::string xscalstr = xscal;
86 int expr_id = expEvaluator.DefineFunction("", xscalstr);
87 xscale = expEvaluator.Evaluate(expr_id);
88 }
89
90 const char *yscal = element->Attribute("YSCALE");
91 if (!yscal)
92 {
93 yscale = 1.0;
94 }
95 else
96 {
97 std::string yscalstr = yscal;
98 int expr_id = expEvaluator.DefineFunction("", yscalstr);
99 yscale = expEvaluator.Evaluate(expr_id);
100 }
101
102 const char *zscal = element->Attribute("ZSCALE");
103 if (!zscal)
104 {
105 zscale = 1.0;
106 }
107 else
108 {
109 std::string zscalstr = zscal;
110 int expr_id = expEvaluator.DefineFunction("", zscalstr);
111 zscale = expEvaluator.Evaluate(expr_id);
112 }
113
114 NekDouble xmove, ymove, zmove;
115
116 // check to see if any moving parameters are in
117 // attributes and determine these values
118
119 const char *xmov = element->Attribute("XMOVE");
120 if (!xmov)
121 {
122 xmove = 0.0;
123 }
124 else
125 {
126 std::string xmovstr = xmov;
127 int expr_id = expEvaluator.DefineFunction("", xmovstr);
128 xmove = expEvaluator.Evaluate(expr_id);
129 }
130
131 const char *ymov = element->Attribute("YMOVE");
132 if (!ymov)
133 {
134 ymove = 0.0;
135 }
136 else
137 {
138 std::string ymovstr = ymov;
139 int expr_id = expEvaluator.DefineFunction("", ymovstr);
140 ymove = expEvaluator.Evaluate(expr_id);
141 }
142
143 const char *zmov = element->Attribute("ZMOVE");
144 if (!zmov)
145 {
146 zmove = 0.0;
147 }
148 else
149 {
150 std::string zmovstr = zmov;
151 int expr_id = expEvaluator.DefineFunction("", zmovstr);
152 zmove = expEvaluator.Evaluate(expr_id);
153 }
154
155 string IsCompressed;
156 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
157
158 if (boost::iequals(IsCompressed,
160 {
161 // Extract the vertex body
162 TiXmlNode *vertexChild = element->FirstChild();
163 ASSERTL0(vertexChild, "Unable to extract the data from the compressed "
164 "vertex tag.");
165
166 std::string vertexStr;
167 if (vertexChild->Type() == TiXmlNode::TINYXML_TEXT)
168 {
169 vertexStr += vertexChild->ToText()->ValueStr();
170 }
171
172 std::vector<SpatialDomains::MeshVertex> vertData;
174 vertData);
175
176 int indx;
177 NekDouble xval, yval, zval;
178 for (int i = 0; i < vertData.size(); ++i)
179 {
180 indx = vertData[i].id;
181 xval = vertData[i].x;
182 yval = vertData[i].y;
183 zval = vertData[i].z;
184
185 xval = xval * xscale + xmove;
186 yval = yval * yscale + ymove;
187 zval = zval * zscale + zmove;
188
190 m_spaceDimension, indx, xval, yval, zval));
191
192 vert->SetGlobalID(indx);
193 m_vertSet[indx] = vert;
194 }
195 }
196 else
197 {
198 ASSERTL0(false, "Compressed formats do not match. Expected :" +
200 " but got " + IsCompressed);
201 }
202}
203
205{
206 // check to see if any scaling parameters are in
207 // attributes and determine these values
208 TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
209 ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
210
211 NekDouble xscale, yscale, zscale;
212
213 LibUtilities::Interpreter expEvaluator;
214 const char *xscal = element->Attribute("XSCALE");
215 if (!xscal)
216 {
217 xscale = 1.0;
218 }
219 else
220 {
221 std::string xscalstr = xscal;
222 int expr_id = expEvaluator.DefineFunction("", xscalstr);
223 xscale = expEvaluator.Evaluate(expr_id);
224 }
225
226 const char *yscal = element->Attribute("YSCALE");
227 if (!yscal)
228 {
229 yscale = 1.0;
230 }
231 else
232 {
233 std::string yscalstr = yscal;
234 int expr_id = expEvaluator.DefineFunction("", yscalstr);
235 yscale = expEvaluator.Evaluate(expr_id);
236 }
237
238 const char *zscal = element->Attribute("ZSCALE");
239 if (!zscal)
240 {
241 zscale = 1.0;
242 }
243 else
244 {
245 std::string zscalstr = zscal;
246 int expr_id = expEvaluator.DefineFunction("", zscalstr);
247 zscale = expEvaluator.Evaluate(expr_id);
248 }
249
250 NekDouble xmove, ymove, zmove;
251
252 // check to see if any moving parameters are in
253 // attributes and determine these values
254
255 const char *xmov = element->Attribute("XMOVE");
256 if (!xmov)
257 {
258 xmove = 0.0;
259 }
260 else
261 {
262 std::string xmovstr = xmov;
263 int expr_id = expEvaluator.DefineFunction("", xmovstr);
264 xmove = expEvaluator.Evaluate(expr_id);
265 }
266
267 const char *ymov = element->Attribute("YMOVE");
268 if (!ymov)
269 {
270 ymove = 0.0;
271 }
272 else
273 {
274 std::string ymovstr = ymov;
275 int expr_id = expEvaluator.DefineFunction("", ymovstr);
276 ymove = expEvaluator.Evaluate(expr_id);
277 }
278
279 const char *zmov = element->Attribute("ZMOVE");
280 if (!zmov)
281 {
282 zmove = 0.0;
283 }
284 else
285 {
286 std::string zmovstr = zmov;
287 int expr_id = expEvaluator.DefineFunction("", zmovstr);
288 zmove = expEvaluator.Evaluate(expr_id);
289 }
290
291 /// Look for elements in CURVE block.
292 TiXmlElement *field = m_xmlGeom->FirstChildElement("CURVED");
293
294 if (!field) // return if no curved entities
295 {
296 return;
297 }
298
299 string IsCompressed;
300 field->QueryStringAttribute("COMPRESSED", &IsCompressed);
301
302 if (IsCompressed.size() == 0)
303 {
304 // this could be that the curved tag is empty
305 // in this case we dont want to read it
306 return;
307 }
308
309 ASSERTL0(boost::iequals(IsCompressed,
311 "Compressed formats do not match. Expected :" +
313 IsCompressed);
314
315 std::vector<SpatialDomains::MeshCurvedInfo> edginfo;
316 std::vector<SpatialDomains::MeshCurvedInfo> facinfo;
318
319 // read edge, face info and curved poitns.
320 TiXmlElement *x = field->FirstChildElement();
321 while (x)
322 {
323 const char *entitytype = x->Value();
324 // read in edge or face info
325 if (boost::iequals(entitytype, "E"))
326 {
327 // read in data
328 std::string elmtStr;
329 TiXmlNode *child = x->FirstChild();
330
331 if (child->Type() == TiXmlNode::TINYXML_TEXT)
332 {
333 elmtStr += child->ToText()->ValueStr();
334 }
335
337 edginfo);
338 }
339 else if (boost::iequals(entitytype, "F"))
340 {
341 // read in data
342 std::string elmtStr;
343 TiXmlNode *child = x->FirstChild();
344
345 if (child->Type() == TiXmlNode::TINYXML_TEXT)
346 {
347 elmtStr += child->ToText()->ValueStr();
348 }
349
351 facinfo);
352 }
353 else if (boost::iequals(entitytype, "DATAPOINTS"))
354 {
355 NekInt id;
356 ASSERTL0(x->Attribute("ID", &id),
357 "Failed to get ID from PTS section");
358 cpts.id = id;
359
360 // read in data
361 std::string elmtStr;
362
363 TiXmlElement *DataIdx = x->FirstChildElement("INDEX");
364 ASSERTL0(DataIdx, "Cannot read data index tag in compressed "
365 "curved section");
366
367 TiXmlNode *child = DataIdx->FirstChild();
368 if (child->Type() == TiXmlNode::TINYXML_TEXT)
369 {
370 elmtStr = child->ToText()->ValueStr();
371 }
372
374 cpts.index);
375
376 TiXmlElement *DataPts = x->FirstChildElement("POINTS");
377 ASSERTL0(DataPts, "Cannot read data pts tag in compressed "
378 "curved section");
379
380 child = DataPts->FirstChild();
381 if (child->Type() == TiXmlNode::TINYXML_TEXT)
382 {
383 elmtStr = child->ToText()->ValueStr();
384 }
385
387 cpts.pts);
388 }
389 else
390 {
391 ASSERTL0(false, "Unknown tag in curved section");
392 }
393 x = x->NextSiblingElement();
394 }
395
396 // rescale (x,y,z) points;
397 for (int i = 0; i < cpts.pts.size(); ++i)
398 {
399 cpts.pts[i].x = xscale * cpts.pts[i].x + xmove;
400 cpts.pts[i].y = yscale * cpts.pts[i].y + ymove;
401 cpts.pts[i].z = zscale * cpts.pts[i].z + zmove;
402 }
403
404 for (int i = 0; i < edginfo.size(); ++i)
405 {
406 int edgeid = edginfo[i].entityid;
408
410 edgeid, ptype = (LibUtilities::PointsType)edginfo[i].ptype));
411
412 // load points
413 int offset = edginfo[i].ptoffset;
414 for (int j = 0; j < edginfo[i].npoints; ++j)
415 {
416 int idx = cpts.index[offset + j];
417
419 m_spaceDimension, edginfo[i].id, cpts.pts[idx].x,
420 cpts.pts[idx].y, cpts.pts[idx].z));
421 curve->m_points.push_back(vert);
422 }
423
424 m_curvedEdges[edgeid] = curve;
425 }
426
427 for (int i = 0; i < facinfo.size(); ++i)
428 {
429 int faceid = facinfo[i].entityid;
431
433 faceid, ptype = (LibUtilities::PointsType)facinfo[i].ptype));
434
435 int offset = facinfo[i].ptoffset;
436 for (int j = 0; j < facinfo[i].npoints; ++j)
437 {
438 int idx = cpts.index[offset + j];
439
441 m_spaceDimension, facinfo[i].id, cpts.pts[idx].x,
442 cpts.pts[idx].y, cpts.pts[idx].z));
443 curve->m_points.push_back(vert);
444 }
445
446 m_curvedFaces[faceid] = curve;
447 }
448}
449
451{
452 CurveMap::iterator it;
453
454 /// Look for elements in ELEMENT block.
455 TiXmlElement *field = m_xmlGeom->FirstChildElement("EDGE");
456
457 ASSERTL0(field, "Unable to find EDGE tag in file.");
458
459 string IsCompressed;
460 field->QueryStringAttribute("COMPRESSED", &IsCompressed);
461
462 ASSERTL0(boost::iequals(IsCompressed,
464 "Compressed formats do not match. Expected :" +
466 IsCompressed);
467 // Extract the edge body
468 TiXmlNode *edgeChild = field->FirstChild();
469 ASSERTL0(edgeChild, "Unable to extract the data from "
470 "the compressed edge tag.");
471
472 std::string edgeStr;
473 if (edgeChild->Type() == TiXmlNode::TINYXML_TEXT)
474 {
475 edgeStr += edgeChild->ToText()->ValueStr();
476 }
477
478 std::vector<SpatialDomains::MeshEdge> edgeData;
480
481 int indx;
482 for (int i = 0; i < edgeData.size(); ++i)
483 {
484 indx = edgeData[i].id;
485 PointGeomSharedPtr vertices[2] = {GetVertex(edgeData[i].v0),
486 GetVertex(edgeData[i].v1)};
487 SegGeomSharedPtr edge;
488
489 it = m_curvedEdges.find(indx);
490 if (it == m_curvedEdges.end())
491 {
493 indx, m_spaceDimension, vertices);
494 }
495 else
496 {
498 indx, m_spaceDimension, vertices, it->second);
499 }
500 m_segGeoms[indx] = edge;
501 }
502}
503
505{
506 /// Look for elements in FACE block.
507 TiXmlElement *field = m_xmlGeom->FirstChildElement("FACE");
508
509 ASSERTL0(field, "Unable to find FACE tag in file.");
510
511 /// All faces are of the form: "<? ID="#"> ... </?>", with
512 /// ? being an element type (either Q or T).
513 /// They might be in compressed format and so then need upacking.
514
515 TiXmlElement *element = field->FirstChildElement();
516 CurveMap::iterator it;
517
518 while (element)
519 {
520 std::string elementType(element->ValueStr());
521
522 ASSERTL0(elementType == "Q" || elementType == "T",
523 (std::string("Unknown 3D face type: ") + elementType).c_str());
524
525 string IsCompressed;
526 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
527
528 ASSERTL0(
529 boost::iequals(IsCompressed,
531 "Compressed formats do not match. Expected :" +
533 IsCompressed);
534
535 // Extract the face body
536 TiXmlNode *faceChild = element->FirstChild();
537 ASSERTL0(faceChild, "Unable to extract the data from "
538 "the compressed face tag.");
539
540 std::string faceStr;
541 if (faceChild->Type() == TiXmlNode::TINYXML_TEXT)
542 {
543 faceStr += faceChild->ToText()->ValueStr();
544 }
545
546 int indx;
547 if (elementType == "T")
548 {
549 std::vector<SpatialDomains::MeshTri> faceData;
551 faceData);
552
553 for (int i = 0; i < faceData.size(); ++i)
554 {
555 indx = faceData[i].id;
556
557 /// See if this face has curves.
558 it = m_curvedFaces.find(indx);
559
560 /// Create a TriGeom to hold the new definition.
562 GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
563 GetSegGeom(faceData[i].e[2])};
564
565 TriGeomSharedPtr trigeom;
566 if (it == m_curvedFaces.end())
567 {
568 trigeom =
570 }
571 else
572 {
574 indx, edges, it->second);
575 }
576 trigeom->SetGlobalID(indx);
577 m_triGeoms[indx] = trigeom;
578 }
579 }
580 else if (elementType == "Q")
581 {
582 std::vector<SpatialDomains::MeshQuad> faceData;
584 faceData);
585
586 for (int i = 0; i < faceData.size(); ++i)
587 {
588 indx = faceData[i].id;
589
590 /// See if this face has curves.
591 it = m_curvedFaces.find(indx);
592
593 /// Create a QuadGeom to hold the new definition.
595 GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
596 GetSegGeom(faceData[i].e[2]), GetSegGeom(faceData[i].e[3])};
597
598 QuadGeomSharedPtr quadgeom;
599 if (it == m_curvedFaces.end())
600 {
601 quadgeom =
603 }
604 else
605 {
607 indx, edges, it->second);
608 }
609 quadgeom->SetGlobalID(indx);
610 m_quadGeoms[indx] = quadgeom;
611 }
612 }
613 /// Keep looking
614 element = element->NextSiblingElement();
615 }
616}
617
619{
620 TiXmlElement *field = nullptr;
621
622 /// Look for elements in ELEMENT block.
623 field = m_xmlGeom->FirstChildElement("ELEMENT");
624
625 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
626
627 /// All elements are of the form: "<S ID = n> ... </S>", with
628 /// ? being the element type.
629
630 TiXmlElement *segment = field->FirstChildElement("S");
631 CurveMap::iterator it;
632
633 while (segment)
634 {
635 string IsCompressed;
636 segment->QueryStringAttribute("COMPRESSED", &IsCompressed);
637 ASSERTL0(
638 boost::iequals(IsCompressed,
640 "Compressed formats do not match. Expected :" +
642 IsCompressed);
643
644 // Extract the face body
645 TiXmlNode *child = segment->FirstChild();
646 ASSERTL0(child, "Unable to extract the data from "
647 "the compressed face tag.");
648
649 std::string str;
650 if (child->Type() == TiXmlNode::TINYXML_TEXT)
651 {
652 str += child->ToText()->ValueStr();
653 }
654
655 int indx;
656
657 std::vector<SpatialDomains::MeshEdge> data;
659
660 for (int i = 0; i < data.size(); ++i)
661 {
662 indx = data[i].id;
663
664 /// See if this face has curves.
665 it = m_curvedEdges.find(indx);
666
667 PointGeomSharedPtr vertices[2] = {GetVertex(data[i].v0),
668 GetVertex(data[i].v1)};
670
671 if (it == m_curvedEdges.end())
672 {
674 indx, m_spaceDimension, vertices);
675 seg->SetGlobalID(indx); // Set global mesh id
676 }
677 else
678 {
680 indx, m_spaceDimension, vertices, it->second);
681 seg->SetGlobalID(indx); // Set global mesh id
682 }
683 seg->SetGlobalID(indx);
684 m_segGeoms[indx] = seg;
685 }
686 /// Keep looking for additional segments
687 segment = segment->NextSiblingElement("S");
688 }
689}
690
692{
693 /// Look for elements in ELEMENT block.
694 TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
695
696 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
697
698 // Set up curve map for curved elements on an embedded manifold.
699 CurveMap::iterator it;
700
701 /// All elements are of the form: "<? ID="#"> ... </?>", with
702 /// ? being the element type.
703
704 TiXmlElement *element = field->FirstChildElement();
705
706 while (element)
707 {
708 std::string elementType(element->ValueStr());
709
710 ASSERTL0(
711 elementType == "Q" || elementType == "T",
712 (std::string("Unknown 2D element type: ") + elementType).c_str());
713
714 string IsCompressed;
715 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
716
717 ASSERTL0(
718 boost::iequals(IsCompressed,
720 "Compressed formats do not match. Expected :" +
722 IsCompressed);
723
724 // Extract the face body
725 TiXmlNode *faceChild = element->FirstChild();
726 ASSERTL0(faceChild, "Unable to extract the data from "
727 "the compressed face tag.");
728
729 std::string faceStr;
730 if (faceChild->Type() == TiXmlNode::TINYXML_TEXT)
731 {
732 faceStr += faceChild->ToText()->ValueStr();
733 }
734
735 int indx;
736 if (elementType == "T")
737 {
738 std::vector<SpatialDomains::MeshTri> faceData;
740 faceData);
741
742 for (int i = 0; i < faceData.size(); ++i)
743 {
744 indx = faceData[i].id;
745
746 /// See if this face has curves.
747 it = m_curvedFaces.find(indx);
748
749 /// Create a TriGeom to hold the new definition.
751 GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
752 GetSegGeom(faceData[i].e[2])};
753
754 TriGeomSharedPtr trigeom;
755 if (it == m_curvedFaces.end())
756 {
757 trigeom =
759 }
760 else
761 {
763 indx, edges, it->second);
764 }
765 trigeom->SetGlobalID(indx);
766 m_triGeoms[indx] = trigeom;
767 }
768 }
769 else if (elementType == "Q")
770 {
771 std::vector<SpatialDomains::MeshQuad> faceData;
773 faceData);
774
775 for (int i = 0; i < faceData.size(); ++i)
776 {
777 indx = faceData[i].id;
778
779 /// See if this face has curves.
780 it = m_curvedFaces.find(indx);
781
782 /// Create a QuadGeom to hold the new definition.
784 GetSegGeom(faceData[i].e[0]), GetSegGeom(faceData[i].e[1]),
785 GetSegGeom(faceData[i].e[2]), GetSegGeom(faceData[i].e[3])};
786
787 QuadGeomSharedPtr quadgeom;
788 if (it == m_curvedFaces.end())
789 {
790 quadgeom =
792 }
793 else
794 {
796 indx, edges, it->second);
797 }
798 quadgeom->SetGlobalID(indx);
799 m_quadGeoms[indx] = quadgeom;
800 }
801 }
802 /// Keep looking
803 element = element->NextSiblingElement();
804 }
805}
806
808{
809 /// Look for elements in ELEMENT block.
810 TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
811
812 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
813
814 /// All elements are of the form: "<? ID="#"> ... </?>", with
815 /// ? being the element type.
816
817 TiXmlElement *element = field->FirstChildElement();
818
819 while (element)
820 {
821 std::string elementType(element->ValueStr());
822
823 // A - tet, P - pyramid, R - prism, H - hex
824 ASSERTL0(
825 elementType == "A" || elementType == "P" || elementType == "R" ||
826 elementType == "H",
827 (std::string("Unknown 3D element type: ") + elementType).c_str());
828
829 string IsCompressed;
830 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
831
832 ASSERTL0(
833 boost::iequals(IsCompressed,
835 "Compressed formats do not match. Expected :" +
837 IsCompressed);
838
839 // Extract the face body
840 TiXmlNode *child = element->FirstChild();
841 ASSERTL0(child, "Unable to extract the data from "
842 "the compressed face tag.");
843
844 std::string str;
845 if (child->Type() == TiXmlNode::TINYXML_TEXT)
846 {
847 str += child->ToText()->ValueStr();
848 }
849
850 int indx;
851 if (elementType == "A")
852 {
853 std::vector<SpatialDomains::MeshTet> data;
855 TriGeomSharedPtr tfaces[4];
856 for (int i = 0; i < data.size(); ++i)
857 {
858 indx = data[i].id;
859 for (int j = 0; j < 4; ++j)
860 {
861 Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
862 tfaces[j] = static_pointer_cast<TriGeom>(face);
863 }
864
865 TetGeomSharedPtr tetgeom(
867 m_tetGeoms[indx] = tetgeom;
868 PopulateFaceToElMap(tetgeom, 4);
869 }
870 }
871 else if (elementType == "P")
872 {
873 std::vector<SpatialDomains::MeshPyr> data;
875 Geometry2DSharedPtr faces[5];
876 for (int i = 0; i < data.size(); ++i)
877 {
878 indx = data[i].id;
879 int Ntfaces = 0;
880 int Nqfaces = 0;
881 for (int j = 0; j < 5; ++j)
882 {
883 Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
884
885 if (face == Geometry2DSharedPtr() ||
886 (face->GetShapeType() != LibUtilities::eTriangle &&
887 face->GetShapeType() != LibUtilities::eQuadrilateral))
888 {
889 std::stringstream errorstring;
890 errorstring << "Element " << indx
891 << " has invalid face: " << j;
892 ASSERTL0(false, errorstring.str().c_str());
893 }
894 else if (face->GetShapeType() == LibUtilities::eTriangle)
895 {
896 faces[j] = static_pointer_cast<TriGeom>(face);
897 Ntfaces++;
898 }
899 else if (face->GetShapeType() ==
901 {
902 faces[j] = static_pointer_cast<QuadGeom>(face);
903 Nqfaces++;
904 }
905 }
906 ASSERTL0((Ntfaces == 4) && (Nqfaces == 1),
907 "Did not identify the correct number of "
908 "triangular and quadrilateral faces for a "
909 "pyramid");
910
911 PyrGeomSharedPtr pyrgeom(
913
914 m_pyrGeoms[indx] = pyrgeom;
915 PopulateFaceToElMap(pyrgeom, 5);
916 }
917 }
918 else if (elementType == "R")
919 {
920 std::vector<SpatialDomains::MeshPrism> data;
922 Geometry2DSharedPtr faces[5];
923 for (int i = 0; i < data.size(); ++i)
924 {
925 indx = data[i].id;
926 int Ntfaces = 0;
927 int Nqfaces = 0;
928 for (int j = 0; j < 5; ++j)
929 {
930 Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
931 if (face == Geometry2DSharedPtr() ||
932 (face->GetShapeType() != LibUtilities::eTriangle &&
933 face->GetShapeType() != LibUtilities::eQuadrilateral))
934 {
935 std::stringstream errorstring;
936 errorstring << "Element " << indx
937 << " has invalid face: " << j;
938 ASSERTL0(false, errorstring.str().c_str());
939 }
940 else if (face->GetShapeType() == LibUtilities::eTriangle)
941 {
942 faces[j] = static_pointer_cast<TriGeom>(face);
943 Ntfaces++;
944 }
945 else if (face->GetShapeType() ==
947 {
948 faces[j] = static_pointer_cast<QuadGeom>(face);
949 Nqfaces++;
950 }
951 }
952 ASSERTL0((Ntfaces == 2) && (Nqfaces == 3),
953 "Did not identify the correct number of "
954 "triangular and quadrilateral faces for a "
955 "prism");
956
957 PrismGeomSharedPtr prismgeom(
959
960 m_prismGeoms[indx] = prismgeom;
961 PopulateFaceToElMap(prismgeom, 5);
962 }
963 }
964 else if (elementType == "H")
965 {
966 std::vector<SpatialDomains::MeshHex> data;
968
969 QuadGeomSharedPtr faces[6];
970 for (int i = 0; i < data.size(); ++i)
971 {
972 indx = data[i].id;
973 for (int j = 0; j < 6; ++j)
974 {
975 Geometry2DSharedPtr face = GetGeometry2D(data[i].f[j]);
976 faces[j] = static_pointer_cast<QuadGeom>(face);
977 }
978
979 HexGeomSharedPtr hexgeom(
981 m_hexGeoms[indx] = hexgeom;
982 PopulateFaceToElMap(hexgeom, 6);
983 }
984 }
985 /// Keep looking
986 element = element->NextSiblingElement();
987 }
988}
989
991 PointGeomMap &verts)
992{
993 if (verts.size() == 0)
994 {
995 return;
996 }
997
998 TiXmlElement *vertTag = new TiXmlElement("VERTEX");
999
1000 vector<MeshVertex> vertInfo;
1001
1002 for (auto &i : verts)
1003 {
1004 MeshVertex v;
1005 v.id = i.first;
1006 v.x = i.second->x();
1007 v.y = i.second->y();
1008 v.z = i.second->z();
1009 vertInfo.push_back(v);
1010 }
1011
1012 vertTag->SetAttribute("COMPRESSED",
1014 vertTag->SetAttribute("BITSIZE",
1016
1017 string vertStr;
1019
1020 vertTag->LinkEndChild(new TiXmlText(vertStr));
1021
1022 geomTag->LinkEndChild(vertTag);
1023}
1024
1025void MeshGraphXmlCompressed::v_WriteEdges(TiXmlElement *geomTag,
1026 SegGeomMap &edges)
1027{
1028 if (edges.size() == 0)
1029 {
1030 return;
1031 }
1032
1033 TiXmlElement *edgeTag =
1034 new TiXmlElement(m_meshDimension == 1 ? "S" : "EDGE");
1035
1036 vector<MeshEdge> edgeInfo;
1037
1038 for (auto &i : edges)
1039 {
1040 MeshEdge e;
1041 e.id = i.first;
1042 e.v0 = i.second->GetVid(0);
1043 e.v1 = i.second->GetVid(1);
1044 edgeInfo.push_back(e);
1045 }
1046
1047 string edgeStr;
1049
1050 edgeTag->SetAttribute("COMPRESSED",
1052 edgeTag->SetAttribute("BITSIZE",
1054
1055 edgeTag->LinkEndChild(new TiXmlText(edgeStr));
1056
1057 if (m_meshDimension == 1)
1058 {
1059 TiXmlElement *tmp = new TiXmlElement("ELEMENT");
1060 tmp->LinkEndChild(edgeTag);
1061 geomTag->LinkEndChild(tmp);
1062 }
1063 else
1064 {
1065 geomTag->LinkEndChild(edgeTag);
1066 }
1067}
1068
1069void MeshGraphXmlCompressed::v_WriteTris(TiXmlElement *faceTag,
1070 TriGeomMap &tris)
1071{
1072 if (tris.size() == 0)
1073 {
1074 return;
1075 }
1076
1077 string tag = "T";
1078
1079 vector<MeshTri> triInfo;
1080
1081 for (auto &i : tris)
1082 {
1083 MeshTri t;
1084 t.id = i.first;
1085 t.e[0] = i.second->GetEid(0);
1086 t.e[1] = i.second->GetEid(1);
1087 t.e[2] = i.second->GetEid(2);
1088 triInfo.push_back(t);
1089 }
1090
1091 TiXmlElement *x = new TiXmlElement(tag);
1092 string triStr;
1094
1095 x->SetAttribute("COMPRESSED",
1097 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1098
1099 x->LinkEndChild(new TiXmlText(triStr));
1100
1101 faceTag->LinkEndChild(x);
1102}
1103
1104void MeshGraphXmlCompressed::v_WriteQuads(TiXmlElement *faceTag,
1105 QuadGeomMap &quads)
1106{
1107 if (quads.size() == 0)
1108 {
1109 return;
1110 }
1111
1112 string tag = "Q";
1113
1114 vector<MeshQuad> quadInfo;
1115
1116 for (auto &i : quads)
1117 {
1118 MeshQuad q;
1119 q.id = i.first;
1120 q.e[0] = i.second->GetEid(0);
1121 q.e[1] = i.second->GetEid(1);
1122 q.e[2] = i.second->GetEid(2);
1123 q.e[3] = i.second->GetEid(3);
1124 quadInfo.push_back(q);
1125 }
1126
1127 TiXmlElement *x = new TiXmlElement(tag);
1128 string quadStr;
1130
1131 x->SetAttribute("COMPRESSED",
1133 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1134
1135 x->LinkEndChild(new TiXmlText(quadStr));
1136
1137 faceTag->LinkEndChild(x);
1138}
1139
1140void MeshGraphXmlCompressed::v_WriteHexs(TiXmlElement *elmtTag,
1141 HexGeomMap &hexs)
1142{
1143 if (hexs.size() == 0)
1144 {
1145 return;
1146 }
1147
1148 string tag = "H";
1149
1150 vector<MeshHex> elementInfo;
1151
1152 for (auto &i : hexs)
1153 {
1154 MeshHex e;
1155 e.id = i.first;
1156 e.f[0] = i.second->GetFid(0);
1157 e.f[1] = i.second->GetFid(1);
1158 e.f[2] = i.second->GetFid(2);
1159 e.f[3] = i.second->GetFid(3);
1160 e.f[4] = i.second->GetFid(4);
1161 e.f[5] = i.second->GetFid(5);
1162 elementInfo.push_back(e);
1163 }
1164
1165 TiXmlElement *x = new TiXmlElement(tag);
1166 string elStr;
1168
1169 x->SetAttribute("COMPRESSED",
1171 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1172
1173 x->LinkEndChild(new TiXmlText(elStr));
1174
1175 elmtTag->LinkEndChild(x);
1176}
1177
1178void MeshGraphXmlCompressed::v_WritePrisms(TiXmlElement *elmtTag,
1179 PrismGeomMap &pris)
1180{
1181 if (pris.size() == 0)
1182 {
1183 return;
1184 }
1185
1186 string tag = "R";
1187
1188 vector<MeshPrism> elementInfo;
1189
1190 for (auto &i : pris)
1191 {
1192 MeshPrism e;
1193 e.id = i.first;
1194 e.f[0] = i.second->GetFid(0);
1195 e.f[1] = i.second->GetFid(1);
1196 e.f[2] = i.second->GetFid(2);
1197 e.f[3] = i.second->GetFid(3);
1198 e.f[4] = i.second->GetFid(4);
1199 elementInfo.push_back(e);
1200 }
1201
1202 TiXmlElement *x = new TiXmlElement(tag);
1203 string elStr;
1205
1206 x->SetAttribute("COMPRESSED",
1208 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1209
1210 x->LinkEndChild(new TiXmlText(elStr));
1211
1212 elmtTag->LinkEndChild(x);
1213}
1214
1215void MeshGraphXmlCompressed::v_WritePyrs(TiXmlElement *elmtTag,
1216 PyrGeomMap &pyrs)
1217{
1218 if (pyrs.size() == 0)
1219 {
1220 return;
1221 }
1222
1223 string tag = "P";
1224
1225 vector<MeshPyr> elementInfo;
1226
1227 for (auto &i : pyrs)
1228 {
1229 MeshPyr e;
1230 e.id = i.first;
1231 e.f[0] = i.second->GetFid(0);
1232 e.f[1] = i.second->GetFid(1);
1233 e.f[2] = i.second->GetFid(2);
1234 e.f[3] = i.second->GetFid(3);
1235 e.f[4] = i.second->GetFid(4);
1236 elementInfo.push_back(e);
1237 }
1238
1239 TiXmlElement *x = new TiXmlElement(tag);
1240 string elStr;
1242
1243 x->SetAttribute("COMPRESSED",
1245 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1246
1247 x->LinkEndChild(new TiXmlText(elStr));
1248
1249 elmtTag->LinkEndChild(x);
1250}
1251
1252void MeshGraphXmlCompressed::v_WriteTets(TiXmlElement *elmtTag,
1253 TetGeomMap &tets)
1254{
1255 if (tets.size() == 0)
1256 {
1257 return;
1258 }
1259
1260 string tag = "A";
1261
1262 vector<MeshTet> elementInfo;
1263
1264 for (auto &i : tets)
1265 {
1266 MeshTet e;
1267 e.id = i.first;
1268 e.f[0] = i.second->GetFid(0);
1269 e.f[1] = i.second->GetFid(1);
1270 e.f[2] = i.second->GetFid(2);
1271 e.f[3] = i.second->GetFid(3);
1272 elementInfo.push_back(e);
1273 }
1274
1275 TiXmlElement *x = new TiXmlElement(tag);
1276 string elStr;
1278
1279 x->SetAttribute("COMPRESSED",
1281 x->SetAttribute("BITSIZE", LibUtilities::CompressData::GetBitSizeStr());
1282
1283 x->LinkEndChild(new TiXmlText(elStr));
1284
1285 elmtTag->LinkEndChild(x);
1286}
1287
1288void MeshGraphXmlCompressed::v_WriteCurves(TiXmlElement *geomTag,
1289 CurveMap &edges, CurveMap &faces)
1290{
1291 if (edges.size() == 0 && faces.size() == 0)
1292 {
1293 return;
1294 }
1295
1296 TiXmlElement *curveTag = new TiXmlElement("CURVED");
1297
1298 vector<MeshCurvedInfo> edgeInfo;
1299 vector<MeshCurvedInfo> faceInfo;
1300 MeshCurvedPts curvedPts;
1301 curvedPts.id = 0;
1302 int ptOffset = 0;
1303 int newIdx = 0;
1304 int edgeCnt = 0;
1305 int faceCnt = 0;
1306
1307 for (auto &i : edges)
1308 {
1309 MeshCurvedInfo cinfo;
1310 cinfo.id = edgeCnt++;
1311 cinfo.entityid = i.first;
1312 cinfo.npoints = i.second->m_points.size();
1313 cinfo.ptype = i.second->m_ptype;
1314 cinfo.ptid = 0;
1315 cinfo.ptoffset = ptOffset;
1316
1317 edgeInfo.push_back(cinfo);
1318
1319 for (int j = 0; j < i.second->m_points.size(); j++)
1320 {
1321 MeshVertex v;
1322 v.id = newIdx;
1323 v.x = i.second->m_points[j]->x();
1324 v.y = i.second->m_points[j]->y();
1325 v.z = i.second->m_points[j]->z();
1326 curvedPts.pts.push_back(v);
1327 curvedPts.index.push_back(newIdx);
1328 newIdx++;
1329 }
1330 ptOffset += cinfo.npoints;
1331 }
1332
1333 for (auto &i : faces)
1334 {
1335 MeshCurvedInfo cinfo;
1336 cinfo.id = faceCnt++;
1337 cinfo.entityid = i.first;
1338 cinfo.npoints = i.second->m_points.size();
1339 cinfo.ptype = i.second->m_ptype;
1340 cinfo.ptid = 0;
1341 cinfo.ptoffset = ptOffset;
1342
1343 faceInfo.push_back(cinfo);
1344
1345 for (int j = 0; j < i.second->m_points.size(); j++)
1346 {
1347 MeshVertex v;
1348 v.id = newIdx;
1349 v.x = i.second->m_points[j]->x();
1350 v.y = i.second->m_points[j]->y();
1351 v.z = i.second->m_points[j]->z();
1352 curvedPts.pts.push_back(v);
1353 curvedPts.index.push_back(newIdx);
1354 newIdx++;
1355 }
1356 ptOffset += cinfo.npoints;
1357 }
1358
1359 curveTag->SetAttribute("COMPRESSED",
1361 curveTag->SetAttribute("BITSIZE",
1363
1364 if (edgeInfo.size())
1365 {
1366 TiXmlElement *x = new TiXmlElement("E");
1367 string dataStr;
1369
1370 x->LinkEndChild(new TiXmlText(dataStr));
1371 curveTag->LinkEndChild(x);
1372 }
1373
1374 if (faceInfo.size())
1375 {
1376 TiXmlElement *x = new TiXmlElement("F");
1377 string dataStr;
1379
1380 x->LinkEndChild(new TiXmlText(dataStr));
1381 curveTag->LinkEndChild(x);
1382 }
1383
1384 if (edgeInfo.size() || faceInfo.size())
1385 {
1386 TiXmlElement *x = new TiXmlElement("DATAPOINTS");
1387 x->SetAttribute("ID", curvedPts.id);
1388 TiXmlElement *subx = new TiXmlElement("INDEX");
1389 string dataStr;
1391 dataStr);
1392 subx->LinkEndChild(new TiXmlText(dataStr));
1393 x->LinkEndChild(subx);
1394
1395 subx = new TiXmlElement("POINTS");
1397 dataStr);
1398 subx->LinkEndChild(new TiXmlText(dataStr));
1399 x->LinkEndChild(subx);
1400 curveTag->LinkEndChild(x);
1401 }
1402
1403 geomTag->LinkEndChild(curveTag);
1404}
1405} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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.
Definition: NekFactory.hpp:197
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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:3977
Geometry2DSharedPtr GetGeometry2D(int gID)
Definition: MeshGraph.h:416
SegGeomSharedPtr GetSegGeom(int id)
Definition: MeshGraph.h:362
PointGeomSharedPtr GetVertex(int id)
Definition: MeshGraph.h:357
void v_WriteHexs(TiXmlElement *elmtTag, HexGeomMap &hexs) override
void v_WriteQuads(TiXmlElement *faceTag, QuadGeomMap &quads) override
void v_WriteTris(TiXmlElement *faceTag, TriGeomMap &tris) override
void v_WriteTets(TiXmlElement *elmtTag, TetGeomMap &tets) override
void v_WritePyrs(TiXmlElement *elmtTag, PyrGeomMap &pyrs) override
void v_WriteEdges(TiXmlElement *geomTag, SegGeomMap &edges) override
void v_WriteCurves(TiXmlElement *geomTag, CurveMap &edges, CurveMap &faces) override
void v_WriteVertices(TiXmlElement *geomTag, PointGeomMap &verts) override
void v_WritePrisms(TiXmlElement *elmtTag, PrismGeomMap &pris) override
static const int kNedges
Definition: QuadGeom.h:74
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:70
int ZlibDecodeFromBase64Str(std::string &in64, std::vector< T > &out)
Definition: CompressData.h:218
int ZlibEncodeToBase64Str(std::vector< T > &in, std::string &out64)
Definition: CompressData.h:130
std::map< int, TriGeomSharedPtr > TriGeomMap
Definition: TriGeom.h:57
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: HexGeom.h:45
std::map< int, PyrGeomSharedPtr > PyrGeomMap
Definition: PyrGeom.h:76
std::map< int, QuadGeomSharedPtr > QuadGeomMap
Definition: QuadGeom.h:52
std::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:82
std::map< int, SegGeomSharedPtr > SegGeomMap
Definition: SegGeom.h:49
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:58
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:59
std::shared_ptr< HexGeom > HexGeomSharedPtr
Definition: HexGeom.h:84
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:59
std::map< int, TetGeomSharedPtr > TetGeomMap
Definition: TetGeom.h:84
std::shared_ptr< PyrGeom > PyrGeomSharedPtr
Definition: PyrGeom.h:75
std::shared_ptr< TetGeom > TetGeomSharedPtr
Definition: TetGeom.h:82
std::map< int, PrismGeomSharedPtr > PrismGeomMap
Definition: PrismGeom.h:83
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:62
std::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:56
MeshGraphFactory & GetMeshGraphFactory()
Definition: MeshGraph.cpp:77
std::map< int, HexGeomSharedPtr > HexGeomMap
Definition: HexGeom.h:85
std::map< int, PointGeomSharedPtr > PointGeomMap
Definition: PointGeom.h:51
std::vector< double > q(NPUPPER *NPUPPER)
std::int32_t NekInt
double NekDouble
NekInt64 ptid
The number of points in this curved entity.
NekInt64 npoints
The entity id corresponding to the global edge/curve.
NekInt64 ptype
point offset of data entry for this curve
NekInt64 ptoffset
the id of point data map (currently always 0 since we are using just one set).
NekInt64 entityid
Id of this curved information.
std::vector< MeshVertex > pts
mapping to access pts value.
std::vector< NekInt64 > index
id of this Point set