Line data Source code
1 : /*
2 : * Copyright (c) 2014-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Sebastian Reiter
4 : *
5 : * This file is part of UG4.
6 : *
7 : * UG4 is free software: you can redistribute it and/or modify it under the
8 : * terms of the GNU Lesser General Public License version 3 (as published by the
9 : * Free Software Foundation) with the following additional attribution
10 : * requirements (according to LGPL/GPL v3 §7):
11 : *
12 : * (1) The following notice must be displayed in the Appropriate Legal Notices
13 : * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 : *
15 : * (2) The following notice must be displayed at a prominent place in the
16 : * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 : *
18 : * (3) The following bibliography is recommended for citation and must be
19 : * preserved in all covered files:
20 : * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 : * parallel geometric multigrid solver on hierarchically distributed grids.
22 : * Computing and visualization in science 16, 4 (2013), 151-164"
23 : * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 : * flexible software system for simulating pde based models on high performance
25 : * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 : *
27 : * This program is distributed in the hope that it will be useful,
28 : * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : * GNU Lesser General Public License for more details.
31 : */
32 :
33 : #include "file_io_vtu.h"
34 : #include <string>
35 : #include <iostream>
36 :
37 : using namespace std;
38 : using namespace rapidxml;
39 :
40 : namespace ug{
41 :
42 : enum VTKCellTypes{
43 : VTK_VERTEX = 1,
44 : VTK_POLY_VERTEX = 2,
45 : VTK_LINE = 3,
46 : VTK_POLY_LINE = 4,
47 : VTK_TRIANGLE = 5,
48 : VTK_TRIANGLE_STRIP = 6,
49 : VTK_POLYGON = 7,
50 : VTK_PIXEL = 8,
51 : VTK_QUAD = 9,
52 : VTK_TETRA = 10,
53 : VTK_VOXEL = 11,
54 : VTK_HEXAHEDRON = 12,
55 : VTK_WEDGE = 13,
56 : VTK_PYRAMID = 14,
57 :
58 : VTK_QUADRATIC_EDGE = 21,
59 : VTK_QUADRATIC_TRIANGLE = 22,
60 : VTK_QUADRATIC_QUAD = 23,
61 : VTK_QUADRATIC_TETRA = 24,
62 : VTK_QUADRATIC_HEXAHEDRON = 25,
63 :
64 : VTK_NUM_TYPES
65 : };
66 :
67 : const char* VTKCellNames[] = { "UNDEFINED",
68 : "VERTEX",
69 : "POLY_VERTEX",
70 : "LINE",
71 : "POLY_LINE",
72 : "TRIANGLE",
73 : "TRIANGLE_STRIP",
74 : "POLYGON",
75 : "PIXEL",
76 : "QUAD",
77 : "TETRA",
78 : "VOXEL",
79 : "HEXAHEDRON",
80 : "WEDGE",
81 : "PYRAMID",
82 : "UNDEFINED",
83 : "UNDEFINED",
84 : "UNDEFINED",
85 : "UNDEFINED",
86 : "UNDEFINED",
87 : "UNDEFINED",
88 : "QUADRATIC_EDGE",
89 : "QUADRATIC_TRIANGLE",
90 : "QUADRATIC_QUAD",
91 : "QUADRATIC_TETRA",
92 : "QUADRATIC_HEXAHEDRON"};
93 :
94 : const int ugRefObjIdToVTKCellType[] = { VTK_VERTEX,
95 : VTK_LINE,
96 : VTK_TRIANGLE,
97 : VTK_QUAD,
98 : VTK_TETRA,
99 : VTK_HEXAHEDRON,
100 : VTK_WEDGE,
101 : VTK_PYRAMID};
102 :
103 0 : bool LoadGridFromVTU(Grid& grid, ISubsetHandler& sh,
104 : const char* filename)
105 : {
106 0 : if(grid.has_vertex_attachment(aPosition))
107 0 : return LoadGridFromVTU(grid, sh, filename, aPosition);
108 0 : else if(grid.has_vertex_attachment(aPosition2))
109 0 : return LoadGridFromVTU(grid, sh, filename, aPosition2);
110 0 : else if(grid.has_vertex_attachment(aPosition1))
111 0 : return LoadGridFromVTU(grid, sh, filename, aPosition1);
112 :
113 : // no standard position attachments are available.
114 : // Attach aPosition and use it.
115 : grid.attach_to_vertices(aPosition);
116 0 : return LoadGridFromVTU(grid, sh, filename, aPosition);
117 : }
118 :
119 :
120 : ////////////////////////////////////////////////////////////////////////////////
121 : ////////////////////////////////////////////////////////////////////////////////
122 : // GridWriterVTU
123 : ////////////////////////////////////////////////////////////////////////////////
124 0 : GridWriterVTU::
125 0 : GridWriterVTU() :
126 0 : m_pout(NULL),
127 0 : m_pieceMode(NONE),
128 0 : m_pointDataMode(NONE),
129 0 : m_cellDataMode(NONE),
130 0 : m_curGrid(NULL)
131 : {
132 :
133 0 : }
134 :
135 0 : GridWriterVTU::
136 0 : GridWriterVTU(std::ostream* pout) :
137 0 : m_pout(NULL),
138 0 : m_pieceMode(NONE),
139 0 : m_pointDataMode(NONE),
140 0 : m_cellDataMode(NONE),
141 0 : m_curGrid(NULL)
142 : {
143 0 : set_stream(pout);
144 0 : }
145 :
146 0 : GridWriterVTU::
147 0 : ~GridWriterVTU()
148 : {
149 :
150 0 : }
151 :
152 0 : void GridWriterVTU::
153 : set_stream(std::ostream* pout)
154 : {
155 0 : if(m_pout)
156 0 : finish();
157 0 : m_pout = pout;
158 :
159 0 : if(!m_pout)
160 : return;
161 :
162 : std::ostream& out = *m_pout;
163 0 : UG_COND_THROW(!out, "Invalid ostream specified fo GridWriterVTU!");
164 :
165 : // todo: check endianiess for binary writes!
166 : out << "<?xml version=\"1.0\"?>" << endl;
167 : out << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">" << endl;
168 : out << " <UnstructuredGrid>" << endl;
169 : }
170 :
171 0 : void GridWriterVTU::
172 : add_subset_handler(ISubsetHandler& sh, const std::string& name)
173 : {
174 0 : UG_COND_THROW(m_cellDataMode == CLOSED,
175 : "A subset handler has to be added before a call to 'end_cell_data'");
176 :
177 0 : m_pieceSubsetHandlers.push_back(make_pair(&sh, name));
178 0 : }
179 :
180 0 : void GridWriterVTU::
181 : write_data_array_header(const char* type, const char* name,
182 : int numberOfComponents)
183 : {
184 0 : ostream& out = out_stream();
185 :
186 0 : out << " <DataArray type=\"" << type << "\"";
187 0 : if(strlen(name) > 0)
188 0 : out << " Name=\"" << name << "\"";
189 0 : if(numberOfComponents > 0)
190 0 : out << " NumberOfComponents=\"" << numberOfComponents << "\"";
191 : out << " format=\"ascii\">" << endl;
192 0 : }
193 :
194 0 : void GridWriterVTU::
195 : write_data_array_footer()
196 : {
197 0 : out_stream() << " </DataArray>" << endl;
198 0 : }
199 :
200 0 : void GridWriterVTU::
201 : begin_point_data()
202 : {
203 0 : UG_COND_THROW(m_pointDataMode != NONE,
204 : "begin_point_data can only be called once per piece!");
205 :
206 0 : UG_COND_THROW(m_cellDataMode == OPEN,
207 : "begin_point_data may not be called between "
208 : "begin_cell_data and end_cell_data!");
209 :
210 0 : out_stream() << " <PointData>" << endl;
211 :
212 0 : m_pointDataMode = OPEN;
213 0 : }
214 :
215 0 : void GridWriterVTU::
216 : end_point_data()
217 : {
218 0 : UG_COND_THROW(m_pointDataMode != OPEN,
219 : "end_point_data has to be called after begin_point_data!");
220 :
221 0 : out_stream() << " </PointData>" << endl;
222 :
223 0 : m_pointDataMode = CLOSED;
224 0 : }
225 :
226 :
227 0 : void GridWriterVTU::
228 : begin_cell_data()
229 : {
230 0 : UG_COND_THROW(m_cellDataMode != NONE,
231 : "begin_cell_data can only be called once per piece!");
232 :
233 0 : UG_COND_THROW(m_pointDataMode == OPEN,
234 : "begin_cell_data may not be called between "
235 : "begin_point_data and end_point_data!");
236 :
237 0 : out_stream() << " <CellData>" << endl;
238 :
239 0 : m_cellDataMode = OPEN;
240 0 : }
241 :
242 :
243 0 : void GridWriterVTU::
244 : end_cell_data()
245 : {
246 0 : UG_COND_THROW(m_cellDataMode != OPEN,
247 : "end_cell_data has to be called after begin_cell_data!");
248 :
249 0 : ostream& out = out_stream();
250 :
251 0 : for(size_t i = 0; i < m_pieceSubsetHandlers.size(); ++i){
252 0 : ISubsetHandler* sh = m_pieceSubsetHandlers[i].first;
253 : const string& name = m_pieceSubsetHandlers[i].second;
254 0 : write_data_array_header("Int32", name.c_str(), 1);
255 :
256 0 : out << " ";
257 0 : for(size_t icell = 0; icell < m_cells.size(); ++icell){
258 0 : out << " " << sh->get_subset_index(m_cells[icell]);
259 : }
260 : out << endl;
261 :
262 0 : write_data_array_footer();
263 : }
264 :
265 : out << " </CellData>" << endl;
266 :
267 0 : m_cellDataMode = CLOSED;
268 :
269 : // finally write regionos-info
270 0 : for(size_t i = 0; i < m_pieceSubsetHandlers.size(); ++i){
271 0 : ISubsetHandler& sh = *m_pieceSubsetHandlers[i].first;
272 : const string& name = m_pieceSubsetHandlers[i].second;
273 : out << " <RegionInfo Name=\"" << name << "\">" << endl;
274 :
275 0 : for(int si = 0; si < sh.num_subsets(); ++si){
276 0 : out << " <Region Name=\"" << sh.subset_info(si).name << "\">";// << endl;
277 : out << "</Region>" << endl;
278 : //out << " </Region>" << endl;
279 : }
280 :
281 : out << " </RegionInfo>" << endl;
282 : }
283 0 : }
284 :
285 0 : void GridWriterVTU::
286 : write_cells(std::vector<GridObject*>& cells, Grid & grid,
287 : AAVrtIndex aaInd)
288 : {
289 0 : ostream& out = out_stream();
290 :
291 : out << " <Cells>" << endl;
292 :
293 0 : write_data_array_header("Int32", "connectivity", 0);
294 : const char* dataPrefix = " ";
295 0 : out << dataPrefix;
296 :
297 : Grid::vertex_traits::secure_container vrts;
298 0 : for(size_t icell = 0; icell < cells.size(); ++icell){
299 0 : GridObject* cell = cells[icell];
300 : grid.associated_elements(vrts, cell);
301 : // Prisms have a different vertex order in ug than in vtk.
302 0 : if(cell->reference_object_id() == ROID_PRISM){
303 0 : out << " " << aaInd[vrts[1]];
304 0 : out << " " << aaInd[vrts[0]];
305 0 : out << " " << aaInd[vrts[2]];
306 0 : out << " " << aaInd[vrts[4]];
307 0 : out << " " << aaInd[vrts[3]];
308 0 : out << " " << aaInd[vrts[5]];
309 : }
310 : else{
311 0 : for(size_t ivrt = 0; ivrt < vrts.size(); ++ivrt){
312 0 : out << " " << aaInd[vrts[ivrt]];
313 : }
314 : }
315 : }
316 :
317 : out << endl;
318 0 : write_data_array_footer();
319 :
320 :
321 0 : write_data_array_header("Int32", "offsets", 0);
322 0 : out << dataPrefix << " ";
323 :
324 : size_t offset = 0;
325 0 : for(size_t icell = 0; icell < cells.size(); ++icell){
326 0 : GridObject* cell = cells[icell];
327 : grid.associated_elements(vrts, cell);
328 0 : offset += vrts.size();
329 : out << " " << offset;
330 : }
331 :
332 : out << endl;
333 0 : write_data_array_footer();
334 :
335 :
336 0 : write_data_array_header("Int8", "types", 0);
337 0 : out << dataPrefix << " ";
338 :
339 0 : for(size_t icell = 0; icell < cells.size(); ++icell){
340 0 : GridObject* cell = cells[icell];
341 0 : int roid = cell->reference_object_id();
342 0 : if(roid >= 0 && roid <= ROID_PYRAMID){
343 0 : out << " " << ugRefObjIdToVTKCellType[roid];
344 : }
345 : else{
346 0 : UG_THROW("Can't map grid-object with ROID " << roid << " to vtk cell type!");
347 : }
348 : }
349 :
350 : out << endl;
351 0 : write_data_array_footer();
352 :
353 : out << " </Cells>" << endl;
354 0 : }
355 :
356 :
357 0 : void GridWriterVTU::
358 : end_piece()
359 : {
360 0 : if(m_pieceMode == OPEN){
361 0 : if(m_pointDataMode == NONE)
362 0 : begin_point_data();
363 0 : if(m_pointDataMode == OPEN)
364 0 : end_point_data();
365 :
366 0 : if(m_cellDataMode == NONE)
367 0 : begin_cell_data();
368 0 : if(m_cellDataMode == OPEN)
369 0 : end_cell_data();
370 :
371 0 : m_pieceMode = CLOSED;
372 0 : out_stream() << " </Piece>" << endl;
373 : }
374 :
375 : m_pieceSubsetHandlers.clear();
376 : m_cells.clear();
377 0 : }
378 :
379 0 : void GridWriterVTU::
380 : finish()
381 : {
382 : // make sure that point and cell data are present in the file
383 0 : if(m_pieceMode == NONE){
384 0 : m_pieceMode = OPEN;
385 : m_pieceSubsetHandlers.clear();
386 : m_cells.clear();
387 0 : out_stream() << " <Piece NumberOfPoints=\"0\" NumberOfCells=\"0\">" << endl;
388 0 : out_stream() << " <Points></Points>" << endl;
389 0 : out_stream() << " <Cells></Cells>" << endl;
390 : }
391 :
392 0 : if(m_pieceMode == OPEN)
393 0 : end_piece();
394 :
395 0 : out_stream() << " </UnstructuredGrid>" << endl;
396 0 : out_stream() << "</VTKFile>" << endl;
397 :
398 : // invalidate out-stream
399 0 : m_pout = NULL;
400 0 : }
401 :
402 :
403 :
404 : ////////////////////////////////////////////////////////////////////////////////
405 : ////////////////////////////////////////////////////////////////////////////////
406 : // GridReaderVTU
407 : ////////////////////////////////////////////////////////////////////////////////
408 0 : GridReaderVTU::GridReaderVTU()
409 : {
410 0 : }
411 :
412 0 : GridReaderVTU::~GridReaderVTU()
413 : {
414 0 : }
415 :
416 : #if 0
417 : // original function
418 : bool GridReaderVTU::
419 : parse_file(const char* filename)
420 : {
421 : ifstream in(filename, ios::binary);
422 : if(!in)
423 : return false;
424 :
425 : m_filename = filename;
426 :
427 : // get the length of the file
428 : streampos posStart = in.tellg();
429 : in.seekg(0, ios_base::end);
430 : streampos posEnd = in.tellg();
431 : streamsize size = posEnd - posStart;
432 :
433 : // go back to the start of the file
434 : in.seekg(posStart);
435 :
436 : // read the whole file en-block and terminate it with 0
437 : char* fileContent = m_doc.allocate_string(0, size + 1);
438 : in.read(fileContent, size);
439 : fileContent[size] = 0;
440 : in.close();
441 :
442 : // parse the xml-data
443 : m_doc.parse<0>(fileContent);
444 :
445 : // notify derived classes that a new document has been parsed.
446 : return new_document_parsed();
447 : }
448 : #else
449 :
450 0 : bool GridReaderVTU::
451 : parse_file(const char* filename)
452 : {
453 0 : ifstream in(filename, ios::binary);
454 0 : if(!in)
455 : return false;
456 :
457 0 : m_filename = filename;
458 :
459 : // get the length of the file
460 0 : streampos posStart = in.tellg();
461 0 : in.seekg(0, ios_base::end);
462 0 : streampos posEnd = in.tellg();
463 : streamsize size = posEnd - posStart;
464 :
465 : // go back to the start of the file
466 0 : in.seekg(posStart);
467 :
468 0 : char * fileContentOriginal = new char[size + 1];
469 :
470 0 : in.read(fileContentOriginal,size);
471 0 : fileContentOriginal[size] = 0;
472 0 : in.close();
473 :
474 0 : std::string fiCo2Str(fileContentOriginal);
475 :
476 0 : delete [] fileContentOriginal;
477 : fileContentOriginal = NULL;
478 :
479 0 : std::string regInf("RegionInfo");
480 :
481 0 : if ( fiCo2Str.find(regInf) == std::string::npos && fiCo2Str.find(m_regionOfInterest) != std::string::npos )
482 : {
483 : // we need to insert the additional string
484 :
485 : std::string regInfLines;
486 :
487 0 : regInfLines.append( "\n<RegionInfo Name=\"" );
488 :
489 : regInfLines.append( m_regionOfInterest );
490 :
491 0 : regInfLines.append("\">\n");
492 :
493 0 : regInfLines.append( "</RegionInfo>" );
494 :
495 0 : std::string insAft( "</CellData>" );
496 :
497 : size_t cedava = fiCo2Str.find( insAft );
498 :
499 0 : if( cedava == std::string::npos )
500 : return false;
501 :
502 0 : size_t insVal = cedava + insAft.size();
503 :
504 : fiCo2Str.insert( insVal, regInfLines );
505 :
506 : }
507 :
508 0 : char* fileContent = m_doc.allocate_string(0, fiCo2Str.size() );
509 :
510 : strcpy(fileContent, fiCo2Str.c_str());
511 :
512 :
513 : // parse the xml-data
514 0 : m_doc.parse<0>(fileContent);
515 :
516 : // notify derived classes that a new document has been parsed.
517 0 : return new_document_parsed();
518 0 : }
519 :
520 : #endif
521 :
522 :
523 :
524 0 : bool GridReaderVTU::
525 : new_document_parsed()
526 : {
527 : // update entries
528 : m_entries.clear();
529 :
530 0 : xml_node<>* vtkNode = m_doc.first_node("VTKFile");
531 0 : UG_COND_THROW(!vtkNode, "Specified file is not a valid VTKFile!");
532 :
533 0 : xml_node<>* ugridNode = vtkNode->first_node("UnstructuredGrid");
534 0 : UG_COND_THROW(!ugridNode, "Specified file does not contain an unstructured grid!");
535 :
536 : // iterate through all grids
537 0 : xml_node<>* curNode = ugridNode->first_node("Piece");
538 0 : while(curNode){
539 0 : m_entries.push_back(GridEntry(curNode));
540 : GridEntry& gridEntry = m_entries.back();
541 :
542 :
543 : // collect associated subset handlers
544 0 : xml_node<>* curSHNode = curNode->first_node("RegionInfo");
545 0 : while(curSHNode){
546 0 : gridEntry.subsetHandlerEntries.push_back(SubsetHandlerEntry(curSHNode));
547 0 : curSHNode = curSHNode->next_sibling("RegionInfo");
548 : }
549 :
550 0 : curNode = curNode->next_sibling("Piece");
551 : }
552 :
553 0 : return true;
554 : }
555 :
556 :
557 0 : const char* GridReaderVTU::
558 : get_grid_name(size_t index) const
559 : {
560 : assert(index < num_grids() && "Bad index!");
561 0 : xml_attribute<>* attrib = m_entries[index].node->first_attribute("name");
562 0 : if(attrib)
563 : return attrib->value();
564 : return "";
565 : }
566 :
567 0 : size_t GridReaderVTU::num_subset_handlers(size_t refGridIndex) const
568 : {
569 : // access the referred grid-entry
570 0 : if(refGridIndex >= m_entries.size()){
571 : UG_LOG("GridReaderVTU::num_subset_handlers: bad refGridIndex. Aborting.\n");
572 0 : return 0;
573 : }
574 :
575 0 : return m_entries[refGridIndex].subsetHandlerEntries.size();
576 : }
577 :
578 0 : const char* GridReaderVTU::
579 : get_subset_handler_name(size_t refGridIndex, size_t subsetHandlerIndex) const
580 : {
581 : assert(refGridIndex < num_grids() && "Bad refGridIndex!");
582 : const GridEntry& ge = m_entries[refGridIndex];
583 : assert(subsetHandlerIndex < ge.subsetHandlerEntries.size() && "Bad subsetHandlerIndex!");
584 :
585 0 : xml_attribute<>* attrib = ge.subsetHandlerEntries[subsetHandlerIndex].node->first_attribute("Name");
586 0 : if(attrib)
587 : return attrib->value();
588 : return "";
589 : }
590 :
591 0 : bool GridReaderVTU::
592 : subset_handler(ISubsetHandler& shOut,
593 : size_t refGridIndex,
594 : size_t subsetHandlerIndex)
595 : {
596 0 : if(refGridIndex >= m_entries.size()){
597 0 : UG_THROW("GridReaderVTU::subset_handler: bad refGridIndex " << refGridIndex
598 : << ". Only " << m_entries.size() << " grids available in the file "
599 : << m_filename);
600 : }
601 :
602 : GridEntry& gridEntry = m_entries[refGridIndex];
603 :
604 : // get the referenced subset-handler entry
605 0 : if(subsetHandlerIndex >= gridEntry.subsetHandlerEntries.size()){
606 : UG_LOG("GridReaderVTU::subset_handler: bad subsetHandlerIndex."
607 : << " Assigning all cells to subset 0.\n");
608 : vector<GridObject*>& cells = gridEntry.cells;
609 0 : for(size_t i = 0; i < cells.size(); ++i){
610 0 : shOut.assign_subset(cells[i], 0);
611 : }
612 : return false;
613 : }
614 :
615 : SubsetHandlerEntry& shEntry = gridEntry.subsetHandlerEntries[subsetHandlerIndex];
616 0 : shEntry.sh = &shOut;
617 0 : string shName = get_subset_handler_name(refGridIndex, subsetHandlerIndex);
618 :
619 : // read region-infos
620 0 : xml_node<>* regionNode = shEntry.node->first_node("Region");
621 : int subsetInd = 0;
622 0 : while(regionNode)
623 : {
624 : // set subset info
625 : // retrieve an initial subset-info from shOut, so that initialised values are kept.
626 0 : SubsetInfo& si = shOut.subset_info(subsetInd);
627 :
628 0 : xml_attribute<>* attrib = regionNode->first_attribute("Name");
629 0 : if(attrib)
630 0 : si.name = attrib->value();
631 :
632 0 : regionNode = regionNode->next_sibling("Region");
633 0 : ++subsetInd;
634 : }
635 :
636 : // read subset-index for each cell
637 0 : xml_node<>* cellDataNode = gridEntry.node->first_node("CellData");
638 0 : UG_COND_THROW(cellDataNode == NULL, "CellData has to be available if regions are defined!");
639 :
640 0 : xml_node<>* regionDataNode = find_child_node_by_argument_value(cellDataNode,
641 : "DataArray",
642 : "Name",
643 : shName.c_str());
644 0 : UG_COND_THROW(regionDataNode == NULL, "No Cell-Data-Array with name " << shName
645 : << " has been defined!");
646 :
647 : vector<GridObject*>& cells = gridEntry.cells;
648 : vector<int> subsetIndices;
649 0 : read_scalar_data(subsetIndices, regionDataNode, true);
650 :
651 0 : if( subsetIndices.size() != cells.size() )
652 : {
653 : vector<double> subsetIndicesDbl;
654 0 : read_scalar_data(subsetIndicesDbl, regionDataNode, true);
655 :
656 0 : trafoDblVec2Int( subsetIndicesDbl, subsetIndices );
657 0 : }
658 :
659 0 : UG_COND_THROW(subsetIndices.size() != cells.size(),
660 : "Mismatch regarding number of cells and number of region-indices!");
661 :
662 0 : for(size_t i = 0; i < cells.size(); ++i){
663 0 : shOut.assign_subset(cells[i], subsetIndices[i]);
664 : }
665 :
666 :
667 : // attrib = subsetNode->first_attribute("color");
668 : // if(attrib){
669 : // stringstream ss(attrib->value(), ios_base::in);
670 : // for(size_t i = 0; i < 4; ++i)
671 : // ss >> si.color[i];
672 : // }
673 :
674 : // attrib = subsetNode->first_attribute("state");
675 : // if(attrib){
676 : // stringstream ss(attrib->value(), ios_base::in);
677 : // size_t state;
678 : // ss >> state;
679 : // si.subsetState = (uint)state;
680 : // }
681 :
682 : // shOut.set_subset_info(subsetInd, si);
683 :
684 : // // read elements of this subset
685 : // if(shOut.elements_are_supported(SHE_VERTEX))
686 : // read_subset_handler_elements<Vertex>(shOut, "vertices",
687 : // subsetNode, subsetInd,
688 : // gridEntry.vertices);
689 : // if(shOut.elements_are_supported(SHE_EDGE))
690 : // read_subset_handler_elements<Edge>(shOut, "edges",
691 : // subsetNode, subsetInd,
692 : // gridEntry.edges);
693 : // if(shOut.elements_are_supported(SHE_FACE))
694 : // read_subset_handler_elements<Face>(shOut, "faces",
695 : // subsetNode, subsetInd,
696 : // gridEntry.faces);
697 : // if(shOut.elements_are_supported(SHE_VOLUME))
698 : // read_subset_handler_elements<Volume>(shOut, "volumes",
699 : // subsetNode, subsetInd,
700 : // gridEntry.volumes);
701 : // // next subset
702 : // subsetNode = subsetNode->next_sibling("subset");
703 : // ++subsetInd;
704 : // }
705 :
706 : return true;
707 0 : }
708 :
709 :
710 :
711 : // This macro should only be used during cell creation in 'create_cells'
712 : // locInd specified the vertex relative to the current offset
713 : #define VRT(locInd) vertices[connectivity[curOffset + (locInd)] + pieceVrtOffset]
714 :
715 0 : bool GridReaderVTU::
716 : create_cells(std::vector<GridObject*>& cellsOut,
717 : Grid& grid,
718 : rapidxml::xml_node<>* cellNode,
719 : std::vector<Vertex*> vertices,
720 : size_t pieceVrtOffset)
721 : {
722 0 : const size_t numPieceVrts = vertices.size() - pieceVrtOffset;
723 :
724 : std::vector<size_t> connectivity;
725 : std::vector<size_t> offsets;
726 : std::vector<size_t> types;
727 :
728 0 : xml_node<>* dataNode = cellNode->first_node("DataArray");
729 0 : for(; dataNode; dataNode = dataNode->next_sibling()){
730 0 : xml_attribute<>* attrib = dataNode->first_attribute("Name");
731 0 : if(!attrib)
732 0 : attrib = dataNode->first_attribute("name");
733 0 : if(!attrib)
734 0 : continue;
735 :
736 0 : if(strcmp(attrib->value(), "connectivity") == 0){
737 0 : read_scalar_data(connectivity, dataNode);
738 : }
739 0 : else if(strcmp(attrib->value(), "offsets") == 0){
740 0 : read_scalar_data(offsets, dataNode);
741 : }
742 0 : else if(strcmp(attrib->value(), "types") == 0){
743 0 : read_scalar_data(types, dataNode);
744 : }
745 : }
746 :
747 0 : UG_COND_THROW(offsets.size() != types.size(),
748 : "VTU parsing error in file " << m_filename
749 : << ": There have to be as many cell-offsets "
750 : "as there are cell-types in a 'cells' node.");
751 :
752 : // UG_LOG("connectivity:");
753 : // for(size_t i = 0; i < connectivity.size(); ++i){
754 : // UG_LOG(" " << connectivity[i]);
755 : // }
756 : // UG_LOG("\n");
757 :
758 : // UG_LOG("offsets:");
759 : // for(size_t i = 0; i < offsets.size(); ++i){
760 : // UG_LOG(" " << offsets[i]);
761 : // }
762 : // UG_LOG("\n");
763 :
764 : // UG_LOG("types:");
765 : // for(size_t i = 0; i < types.size(); ++i){
766 : // UG_LOG(" " << types[i]);
767 : // }
768 : // UG_LOG("\n");
769 :
770 : size_t curOffset = 0;
771 0 : for(size_t icell = 0; icell < types.size(); ++icell){
772 : try{
773 0 : size_t nextOffset = offsets[icell];
774 0 : switch(types[icell]){
775 0 : case VTK_VERTEX:{
776 0 : check_indices(connectivity, curOffset, 1, numPieceVrts);
777 0 : cellsOut.push_back(VRT(0));
778 0 : }break;
779 :
780 0 : case VTK_LINE:{
781 0 : check_indices(connectivity, curOffset, 2, numPieceVrts);
782 0 : RegularEdge* e = *grid.create<RegularEdge>(EdgeDescriptor(VRT(0), VRT(1)));
783 0 : cellsOut.push_back(e);
784 0 : }break;
785 :
786 0 : case VTK_TRIANGLE:{
787 0 : check_indices(connectivity, curOffset, 3, numPieceVrts);
788 : Triangle* f =
789 0 : *grid.create<Triangle>(
790 0 : TriangleDescriptor(VRT(0), VRT(1), VRT(2)));
791 0 : cellsOut.push_back(f);
792 0 : }break;
793 :
794 : // case VTK_TRIANGLE_STRIP:{
795 :
796 : // }break;
797 :
798 0 : case VTK_QUAD:{
799 0 : check_indices(connectivity, curOffset, 4, numPieceVrts);
800 : Quadrilateral* f =
801 0 : *grid.create<Quadrilateral>(
802 0 : QuadrilateralDescriptor(VRT(0), VRT(1), VRT(2), VRT(3)));
803 0 : cellsOut.push_back(f);
804 0 : }break;
805 :
806 0 : case VTK_TETRA:{
807 0 : check_indices(connectivity, curOffset, 4, numPieceVrts);
808 0 : Volume* v = *grid.create<Tetrahedron>(
809 0 : TetrahedronDescriptor(VRT(0), VRT(1),
810 0 : VRT(2), VRT(3)));
811 0 : cellsOut.push_back(v);
812 0 : }break;
813 :
814 0 : case VTK_HEXAHEDRON:{
815 0 : check_indices(connectivity, curOffset, 8, numPieceVrts);
816 0 : Volume* v = *grid.create<Hexahedron>(
817 0 : HexahedronDescriptor(VRT(0), VRT(1), VRT(2), VRT(3),
818 0 : VRT(4), VRT(5), VRT(6), VRT(7)));
819 0 : cellsOut.push_back(v);
820 0 : }break;
821 :
822 0 : case VTK_WEDGE:{
823 0 : check_indices(connectivity, curOffset, 6, numPieceVrts);
824 0 : Volume* v = *grid.create<Prism>(
825 0 : PrismDescriptor(VRT(1), VRT(0), VRT(2),
826 0 : VRT(4), VRT(3), VRT(5)));
827 0 : cellsOut.push_back(v);
828 0 : }break;
829 :
830 0 : case VTK_PYRAMID:{
831 0 : check_indices(connectivity, curOffset, 5, numPieceVrts);
832 0 : Volume* v = *grid.create<Pyramid>(
833 0 : PyramidDescriptor(VRT(0), VRT(1), VRT(2),
834 0 : VRT(3), VRT(4)));
835 0 : cellsOut.push_back(v);
836 0 : }break;
837 :
838 0 : default:{
839 0 : UG_THROW("VTU parsing error in file " << m_filename
840 : << ": Unsupported cell type encountered. ");
841 : }break;
842 : }
843 : curOffset = nextOffset;
844 : }
845 0 : catch(UGError& err){
846 0 : UG_LOG("icell: " << icell << ", types[icell]: " << types[icell] << endl);
847 :
848 0 : string cellName = "UNDEFINED";
849 0 : if(types[icell] > 0 && types[icell] < VTK_NUM_TYPES)
850 0 : cellName = VTKCellNames[types[icell]];
851 :
852 0 : std::stringstream ss;
853 : ss << "Parsed cell type: " << cellName;
854 0 : err.push_msg(ss.str(),__FILE__,__LINE__);
855 0 : throw(err);
856 0 : }
857 : }
858 :
859 0 : return true;
860 0 : }
861 :
862 0 : void GridReaderVTU::
863 : trafoDblVec2Int( std::vector<double> const & dblVec, std::vector<int> & intVec )
864 : {
865 0 : intVec = std::vector<int>();
866 :
867 0 : for( auto d : dblVec )
868 : {
869 0 : intVec.push_back( static_cast<int>( d ) );
870 : }
871 0 : }
872 :
873 :
874 0 : xml_node<>* GridReaderVTU::
875 : find_child_node_by_argument_value(rapidxml::xml_node<>* parent,
876 : const char* nodeName,
877 : const char* argName,
878 : const char* argValue)
879 : {
880 0 : xml_node<>* curNode = parent->first_node(nodeName);
881 0 : while(curNode){
882 0 : xml_attribute<>* attrib = curNode->first_attribute("Name");
883 0 : if(attrib){
884 0 : if(strcmp(attrib->value(), argValue) == 0)
885 0 : return curNode;
886 : }
887 0 : curNode = curNode->next_sibling(nodeName);
888 : }
889 : return NULL;
890 : }
891 :
892 : std::string GridReaderVTU::m_regionOfInterest = "regions";
893 :
894 : }// end of namespace
|