Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
ExporterParaView_def.hpp
1#ifndef ExporterParaView_DEF_hpp
2#define ExporterParaView_DEF_hpp
3
12
13namespace FEDD {
14template<class SC,class LO,class GO,class NO>
15ExporterParaView<SC,LO,GO,NO>::ExporterParaView():
16hdf5exporter_(),
17comm_(),
18closingLinesPosition_(),
19closingLinesPositionTimes_(),
20closingLines_(),
21xmf_out_(),
22xmf_times_out_(),
23filename_(),
24outputFilename_(),
25postfix_(),
26FEType_(),
27variables_(0),
28uniqueMaps_(0),
29varNames_(0),
30varTypes_(0),
31varDofPerNode_(0),
32pointsHDF_(),
33elementsHDF_(),
34dim_(0),
35nmbElementsGlob_(0),
36nmbPointsGlob_(0),
37nmbExportValuesGlob_(0),
38nmbPointsPerElement_(0),
39timeIndex_(0),
40writeDt_(0),
41saveTimestep_(1),
42verbose_(false),
43parameterList_(),
44pointsUnique_()
45{
46
47}
48
49template<class SC,class LO,class GO,class NO>
50void ExporterParaView<SC,LO,GO,NO>::setup(std::string filename,
51 MeshPtr_Type mesh,
52 std::string FEType,
53 int saveTimestep,
54 ParameterListPtr_Type parameterList){
55
56
57 setup( filename, mesh, FEType, parameterList);
58 saveTimestep_ = saveTimestep;
59}
60
61template<class SC,class LO,class GO,class NO>
62void ExporterParaView<SC,LO,GO,NO>::setup(std::string filename,
63 MeshPtr_Type mesh,
64 std::string FEType,
65 ParameterListPtr_Type parameterList){
66
67
68 parameterList_ = parameterList;
69 comm_ = mesh->getComm();
70 verbose_ = (comm_->getRank() == 0);
71 Teuchos::RCP<const Teuchos::MpiComm<int> > mpiComm = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<int> >( mesh->getComm() );
72 filename_ = filename;
73 outputFilename_ = filename_ + ".h5";
74 FEType_ = FEType;
75 dim_ = mesh->getDimension();
76 nmbElementsGlob_ = mesh->getNumElementsGlobal();
77
78 pointsUnique_ = mesh->getPointsUnique();
79
80 if (FEType_ == "P0") {
81 if (dim_ == 2 ) {
82 nmbPointsPerElement_ = 3 ;
83
84 }
85 else if(dim_ == 3){
86 nmbPointsPerElement_ = 4 ;
87
88 }
89 }
90 else if (FEType_ == "P1") {
91 if (dim_ == 2 ) {
92 nmbPointsPerElement_ = 3 ;
93
94 }
95 else if(dim_ == 3){
96 nmbPointsPerElement_ = 4 ;
97
98 }
99 }
100 else if (FEType_ == "P1-disc") {
101 if (dim_ == 2 ) {
102 nmbPointsPerElement_ = 3 ;
103
104 }
105 else if(dim_ == 3){
106 nmbPointsPerElement_ = 4 ;
107
108 }
109 }
110 else if(FEType_ == "P2"){
111
112 if (dim_ == 2 ) {
113 nmbPointsPerElement_ = 6 ;
114
115 }
116 else if(dim_ == 3){
117 nmbPointsPerElement_ = 10 ;
118 }
119 }
120 else if(FEType_ == "P2-CR"){
121 if (dim_ == 2 ) {
122 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Wrong dimension for P2-CR.");
123 }
124 else if(dim_ == 3){
125 nmbPointsPerElement_ = 10 ; // only export P2 points. Paraview 5.6.0 and prior versions ignore face centered values.
126 }
127 }
128 else if(FEType_ == "Q1"){
129 if (dim_ == 2 ) {
130 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Wrong dimension for Q1.");
131 }
132 else if(dim_ == 3){
133 nmbPointsPerElement_ = 8 ;
134 }
135 }
136
137 else if(FEType_ == "Q2"){
138 if (dim_ == 2 ) {
139 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Wrong dimension for Q2.");
140 }
141 else if(dim_ == 3){
142 nmbPointsPerElement_ = 20 ; // only export Q20 points.
143 }
144 }
145 else if(FEType_ == "Q2-20"){
146 if (dim_ == 2 ) {
147 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Wrong dimension for Q2-20.");
148 }
149 else if(dim_ == 3){
150 nmbPointsPerElement_ = 20 ; // only export Q20 points.
151 }
152 }
153 else {
154 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Wrong FEType, choose either P0 or P1 or P1-disc or P2 or P2-CR or Q1 or Q2 or Q2-20. Subdomain export with P0 - Export stopped");
155 }
156
157 nmbPointsGlob_ = mesh->getMapUnique()->getGlobalNumElements();
158
159 closingLines_ = "\n </Grid>\n\n </Domain>\n</Xdmf>\n";
160
161 timeIndex_ = 0;
162
163 // Something different happens to the element List and the elements
164 // Probably they have the following form:
165 // ElementMap : 0 1 2 3 4 5
166 // -> 0 1 2 3 4 5 6 7 8 ...
167 // Element GIDs: 0 1 2 3 ...
168 // Where the map tells us which IDs belong to which element
169 MapConstPtr_Type elementMap = mesh->getElementMap();
170 Teuchos::ArrayView<const GO> nodeElementList = elementMap->getNodeElementList(); // Global Ids on this processor
171 vec_GO_Type nodeElementListInteger( nmbPointsPerElement_ * nodeElementList.size() );
172 int counter=0;
173 for (int i=0; i<nodeElementList.size(); i++) { // Number of elements
174 for (int j=0; j<nmbPointsPerElement_; j++){ // number of points per element
175 nodeElementListInteger[counter] = (int) nmbPointsPerElement_*nodeElementList[i] + j;
176 counter++; // from 0 ... to nmbPointsPerElement_*localElements
177 }
178 }
179 Teuchos::ArrayView<GO> globalMapIDs = Teuchos::arrayViewFromVector( nodeElementListInteger);
180 MapPtr_Type mapElements = Teuchos::rcp( new Map_Type((int) (nmbPointsPerElement_*nmbElementsGlob_), globalMapIDs,elementMap->getIndexBase()*nmbPointsPerElement_, comm_));
181
182 // They contain global IDs of nodes corresponding to 'elements'
183 elementsHDF_.reset(new MultiVector_Type(mapElements,1));
184
185 ElementsPtr_Type elements = mesh->getElementsC();
186 counter = 0;
187 for (int i=0; i<elements->numberElements(); i++) {
188 for (int j=0; j<nmbPointsPerElement_; j++) {
189 int globalIndex = (int) mesh->getMapRepeated()->getGlobalElement( elements->getElement(i).getNode(j) );
190 (elementsHDF_->getDataNonConst(0))[counter] = globalIndex;
191 counter++;
192 }
193 }
194
195
196 pointsHDF_.reset(new MultiVector_Type(mesh->getMapUnique(),dim_));
197
198 updatePoints();
199
200 this->initHDF5();
201
202 this->initXmf();
203
204 writeDt_ = false;
205
206 saveTimestep_ = 1;
207}
208
209template<class SC,class LO,class GO,class NO>
210void ExporterParaView<SC,LO,GO,NO>::addVariable(MultiVecConstPtr_Type &u,
211 std::string varName,
212 std::string varType,
213 int dofPerNode,
214 MapConstPtrConst_Type& mapUnique){
215
216 variables_.push_back(u);
217 varNames_.push_back(varName);
218 varTypes_.push_back(varType);
219 varDofPerNode_.push_back(dofPerNode);
220
221 nmbExportValuesGlob_ = mapUnique->getGlobalNumElements();
222
223 uniqueMaps_.push_back(mapUnique);
224
225}
226
227template<class SC,class LO,class GO,class NO>
228void ExporterParaView<SC,LO,GO,NO>::save(double time){
229
230 if (timeIndex_ % saveTimestep_ == 0) {
231 makePostfix();
232
233 if (!parameterList_.is_null()){
234 if ( parameterList_->sublist("Exporter").get("Write new mesh",false) )
235 writeMeshPointsHDF5();
236 }
237
238 writeVariablesHDF5();
239
240 writeXmf(time);
241 }
242 else{
243 if (this->verbose_)
244 std::cout << "\n \t ### Export criterion not satisfied for current time step - no export of time step! ###" << std::endl;
245 }
246
247 timeIndex_++;
248
249}
250
251template<class SC,class LO,class GO,class NO>
252void ExporterParaView<SC,LO,GO,NO>::save(double time, double dt){
253
254 if (timeIndex_ % saveTimestep_ == 0) {
255
256 makePostfix();
257 if (!parameterList_.is_null()){
258 if ( parameterList_->sublist("Exporter").get("Write new mesh",false) )
259 writeMeshPointsHDF5();
260 }
261 writeVariablesHDF5();
262
263 writeXmf(time);
264
265 if (timeIndex_==0) {
266 initXmfTimes();
267 }
268 writeXmfTime(time, dt);
269 }
270 else{
271 if (this->verbose_)
272 std::cout << "\n \t ### Export criterion not satisfied for current time step - no export of time step! ###" << std::endl;
273 }
274
275 timeIndex_++;
276
277}
278
279template<class SC,class LO,class GO,class NO>
280void ExporterParaView<SC,LO,GO,NO>::closeExporter(){
281
282 hdf5exporter_->close();
283 xmf_out_.close();
284 if (writeDt_) {
285 xmf_times_out_.close();
286 }
287
288}
289
290template<class SC,class LO,class GO,class NO>
291void ExporterParaView<SC,LO,GO,NO>::initHDF5(){
292
293 hdf5exporter_.reset( new HDF5_Type(comm_) );
294 hdf5exporter_->create(outputFilename_);
295 std::string nameConn = "Connections";
296
297 writeMeshElements(nameConn);
298 // Mesh is only written once. If we have a new mesh for a new export, we call writeMeshPointsHDF5() in function save(...), when all the data is updated.
299 if (parameterList_.is_null())
300 writeMeshPointsHDF5();
301 else if ( parameterList_->sublist("Exporter").get("Write new mesh",false) == false )
302 writeMeshPointsHDF5();
303
304}
305
306template<class SC,class LO,class GO,class NO>
307void ExporterParaView<SC,LO,GO,NO>::writeMeshPointsHDF5(){
308 if (!parameterList_.is_null()){
309 if ( parameterList_->sublist("Exporter").get("Write new mesh",false) ) {
310 updatePoints();
311 std::string nameP_X = "PointsX" + std::to_string(timeIndex_);
312 std::string nameP_Y = "PointsY" + std::to_string(timeIndex_);
313 std::string nameP_Z = "PointsZ" + std::to_string(timeIndex_);
314
315 writeMeshPoints(nameP_X, nameP_Y, nameP_Z );
316 }
317
318 else{
319 std::string nameP_X = "PointsX";
320 std::string nameP_Y = "PointsY";
321 std::string nameP_Z = "PointsZ";
322 writeMeshPoints( nameP_X, nameP_Y, nameP_Z );
323 }
324 }
325 else{
326 std::string nameP_X = "PointsX";
327 std::string nameP_Y = "PointsY";
328 std::string nameP_Z = "PointsZ";
329 writeMeshPoints( nameP_X, nameP_Y, nameP_Z );
330 }
331}
332
333template<class SC,class LO,class GO,class NO>
334void ExporterParaView<SC,LO,GO,NO>::writeMeshElements( std::string nameConn ){
335 //Triangle/Tetrahedron Connections
336 hdf5exporter_->createGroup(nameConn);
337
338 hdf5exporter_->write(nameConn,elementsHDF_);
339}
340
341template<class SC,class LO,class GO,class NO>
342void ExporterParaView<SC,LO,GO,NO>::writeMeshPoints(std::string nameP_X,
343 std::string nameP_Y,
344 std::string nameP_Z){
345
346
347
348 hdf5exporter_->createGroup(nameP_X);
349 hdf5exporter_->createGroup(nameP_Y);
350
351 hdf5exporter_->write(nameP_X,pointsHDF_->getVector(0));
352 hdf5exporter_->write(nameP_Y,pointsHDF_->getVector(1));
353
354 if (dim_ == 3) {
355 hdf5exporter_->createGroup(nameP_Z);
356 hdf5exporter_->write(nameP_Z,pointsHDF_->getVector(2));
357 }
358 hdf5exporter_->flush();
359
360}
361
362template<class SC,class LO,class GO,class NO>
363void ExporterParaView<SC,LO,GO,NO>::updatePoints(){
364 int dim = -1;
365 if (pointsUnique_->size()>0)
366 dim = pointsUnique_->at(0).size();
367
368 for (int i=0; i<pointsUnique_->size(); i++) {
369 for (int j = 0; j < dim; j++) {
370 pointsHDF_->getDataNonConst(j)[i] = (*pointsUnique_)[i][j];
371 }
372 }
373}
374
375template<class SC,class LO,class GO,class NO>
376void ExporterParaView<SC,LO,GO,NO>::writeVariablesHDF5(){
377
378 for (int i=0; i<variables_.size(); i++) {
379
380 if (varTypes_.at(i)=="Vector") {
381 MultiVectorPtr_Type u_export(new MultiVector_Type(uniqueMaps_.at(i),3)); // ParaView always uses 3D Data. for 2D Data the last entries (for the 3rd Dim) are all zero.
382
383 hdf5exporter_->createGroup(varNames_[i]+postfix_);
384
385 prepareVectorField(variables_.at(i),u_export, varDofPerNode_.at(i));
386
387 hdf5exporter_->write(varNames_[i]+postfix_,u_export,true);
388 }
389 else if(varTypes_.at(i)=="Scalar"){
390 MultiVectorPtr_Type u_export(new MultiVector_Type(uniqueMaps_.at(i),1)); // ParaView always uses 3D Data. for 2D Data the last entries (for the 3rd Dim) are all zero.
391
392 hdf5exporter_->createGroup(varNames_[i]+postfix_);
393
394 prepareScalar(variables_.at(i),u_export); //conversion to int
395
396 hdf5exporter_->write(varNames_[i]+postfix_,u_export);
397 }
398 hdf5exporter_->flush();
399 }
400
401}
402
403template<class SC,class LO,class GO,class NO>
404void ExporterParaView<SC,LO,GO,NO>::initXmf(){
405
406 if (comm_->getRank()==0) {
407
408 xmf_out_.open((filename_ + ".xmf").c_str(),std::ios_base::out);
409 xmf_out_ << "<?xml version=\"1.0\" ?>\n"
410 << "<!DOCTYPE Xdmf SYSTEM \""
411 << filename_
412 << ".xdmf\" [\n"
413 << "<!ENTITY DataFile \""
414 << filename_
415 << ".h5\">\n"
416 << "]>\n"
417 << "<!-- "
418 << filename_
419 << ".h5 -->\n"
420 << "<Xdmf>\n"
421 << " <Domain Name=\""
422 << filename_
423 << "\">\n"
424 << " <Grid Name=\""
425 << filename_
426 << "Grid\" GridType=\"Collection\" CollectionType=\"Temporal\">\n"
427 << "\n";
428
429 closingLinesPosition_ = xmf_out_.tellp();
430 // write closing lines
431 xmf_out_ << closingLines_;
432 xmf_out_.flush();
433 }
434
435}
436
437template<class SC,class LO,class GO,class NO>
438void ExporterParaView<SC,LO,GO,NO>::initXmfTimes(){
439 writeDt_ = true;
440 if (comm_->getRank()==0) {
441
442 xmf_times_out_.open((filename_ + "_times.xmf").c_str(),std::ios_base::out);
443 xmf_times_out_ << "<?xml version=\"1.0\" ?>\n"
444 << "<!DOCTYPE Xdmf SYSTEM \""
445 << filename_
446 << ".xdmf\" []>\n"
447 << "<Xdmf>\n"
448 << " <Domain Name=\""
449 << filename_
450 << "\">\n"
451 << " <Grid Name=\""
452 << filename_
453 << "Grid\" GridType=\"Collection\" CollectionType=\"Temporal\">\n"
454 << "\n";
455
456 closingLinesPositionTimes_ = xmf_times_out_.tellp();
457 // write closing lines
458 xmf_times_out_ << closingLines_;
459 xmf_times_out_.flush();
460 }
461
462}
463template<class SC,class LO,class GO,class NO>
464void ExporterParaView<SC,LO,GO,NO>::writeXmf(double time){
465
466 std::string nameP_X;
467 std::string nameP_Y;
468 std::string nameP_Z;
469 std::string nameConn = "Connections";
470 if(redo_ == true)
471 nameConn = "Connections" + std::to_string(timeIndex_);
472 if (!parameterList_.is_null()){
473 if ( parameterList_->sublist("Exporter").get("Write new mesh",false) ) {
474 nameP_X = "PointsX" + std::to_string(timeIndex_);
475 nameP_Y = "PointsY" + std::to_string(timeIndex_);
476 nameP_Z = "PointsZ" + std::to_string(timeIndex_);
477 }
478 else{
479 nameP_X = "PointsX";
480 nameP_Y = "PointsY";
481 nameP_Z = "PointsZ";
482 }
483 }
484 else {
485 nameP_X = "PointsX";
486 nameP_Y = "PointsY";
487 nameP_Z = "PointsZ";
488 }
489 writeXmfElements( nameConn, time);
490 writeXmfPoints( nameP_X, nameP_Y, nameP_Z );
491 writeXmfVariables();
492
493}
494
495template<class SC,class LO,class GO,class NO>
496void ExporterParaView<SC,LO,GO,NO>::writeXmfElements( std::string nameConn, double time ){
497
498 if (verbose_) {
499 xmf_out_.seekp (closingLinesPosition_);
500
501 xmf_out_ <<
502 "<!-- Time " << time << " Iteration " << postfix_.substr (1, 5) << " -->\n" <<
503 " <Grid Name=\"Mesh" << FEType_ << " " << time << "\">\n" <<
504 " <Time TimeType=\"Single\" Value=\"" << time << "\" />\n";
505
506 // writeTopology (M_xdmf);
507 std::string FEstring;
508 if (FEType_=="P0") {
509 if (dim_ == 2) {
510 FEstring = "Triangle";
511 }
512 else if(dim_ == 3){
513 FEstring = "Tetrahedron";
514 }
515 }
516 else if (FEType_=="P1-disc") {
517 if (dim_ == 2) {
518 FEstring = "Triangle";
519 }
520 else if(dim_ == 3){
521 FEstring = "Tetrahedron";
522 }
523 }
524 else if (FEType_=="P1") {
525 if (dim_ == 2) {
526 FEstring = "Triangle";
527 }
528 else if(dim_ == 3){
529 FEstring = "Tetrahedron";
530 }
531 }
532
533 else if (FEType_=="P2"){
534 if (dim_ == 2) {
535 FEstring = "Tri_6";
536 }
537 else if(dim_ == 3){
538 FEstring = "Tet_10";
539 }
540 }
541 else if (FEType_=="P2-CR"){
542 //there is no 2D version
543 if(dim_ == 3){
544 FEstring = "Tet_10";
545 }
546 }
547 else if (FEType_=="Q1"){
548 //there is no 2D version
549 if(dim_ == 3){
550 FEstring = "Hexahedron";
551 }
552 }
553 else if (FEType_=="Q2"){
554 //there is no 2D version
555 if(dim_ == 3){
556 FEstring = "Hex_20";
557 }
558 }
559 else if (FEType_=="Q2-20"){
560 //there is no 2D version
561 if(dim_ == 3){
562 FEstring = "Hex_20";
563 }
564 }
565
566 xmf_out_ << " <Topology\n"
567 << " Type=\""
568 << FEstring
569 << "\"\n"
570 << " NumberOfElements=\""
571 << nmbElementsGlob_
572 << "\"\n"
573 << " BaseOffset=\""
574 << 0
575 << "\">\n"
576 << " <DataStructure Format=\"HDF\"\n"
577 << " Dimensions=\""
578 << nmbElementsGlob_//this->M_mesh->numGlobalElements()
579 << " "
580 << nmbPointsPerElement_//this->M_mesh->numLocalVertices()
581 << "\"\n" << " DataType=\"Int\"\n"
582 << " Precision=\"8\">\n" << " "
583 << outputFilename_ << ":/"<< nameConn <<"/Values\n"
584 << " </DataStructure>\n" << " </Topology>\n";
585
586 closingLinesPosition_ = xmf_out_.tellp();
587
588 }
589}
590template<class SC,class LO,class GO,class NO>
591void ExporterParaView<SC,LO,GO,NO>::writeXmfPoints(std::string nameP_X,
592 std::string nameP_Y,
593 std::string nameP_Z ){
594
595 if (verbose_) {
596 xmf_out_.seekp (closingLinesPosition_);
597 // writeGeometry (M_xdmf);
598 if (dim_ == 2) {
599 xmf_out_ <<
600 " <Geometry Type=\"X_Y\">\n" <<
601 " <DataStructure Format=\"HDF\"\n" <<
602 " Dimensions=\"" << nmbPointsGlob_ << "\"\n" <<
603 " DataType=\"Float\"\n" <<
604 " Precision=\"8\">\n" <<
605 " " << outputFilename_ << ":/" << nameP_X << "/Values\n" <<
606 " </DataStructure>\n" <<
607 " <DataStructure Format=\"HDF\"\n" <<
608 " Dimensions=\"" << nmbPointsGlob_ << "\"\n" <<
609 " DataType=\"Float\"\n" <<
610 " Precision=\"8\">\n" <<
611 " " << outputFilename_ << ":/" << nameP_Y << "/Values\n" <<
612 " </DataStructure>\n" <<
613 " </Geometry>\n" <<
614 "\n";
615 }
616 else if(dim_ ==3){
617 xmf_out_ <<
618 " <Geometry Type=\"X_Y_Z\">\n" <<
619 " <DataStructure Format=\"HDF\"\n" <<
620 " Dimensions=\"" << nmbPointsGlob_ << "\"\n" <<
621 " DataType=\"Float\"\n" <<
622 " Precision=\"8\">\n" <<
623 " " << outputFilename_ << ":/" << nameP_X << "/Values\n" <<
624 " </DataStructure>\n" <<
625 " <DataStructure Format=\"HDF\"\n" <<
626 " Dimensions=\"" << nmbPointsGlob_ << "\"\n" <<
627 " DataType=\"Float\"\n" <<
628 " Precision=\"8\">\n" <<
629 " " << outputFilename_ << ":/" << nameP_Y << "/Values\n" <<
630 " </DataStructure>\n" <<
631 " <DataStructure Format=\"HDF\"\n" <<
632 " Dimensions=\"" << nmbPointsGlob_ << "\"\n" <<
633 " DataType=\"Float\"\n" <<
634 " Precision=\"8\">\n" <<
635 " " << outputFilename_ << ":/" << nameP_Z << "/Values\n" <<
636 " </DataStructure>\n" <<
637 " </Geometry>\n" <<
638 "\n";
639 }
640 closingLinesPosition_ = xmf_out_.tellp();
641
642 }
643
644}
645
646template<class SC,class LO,class GO,class NO>
647void ExporterParaView<SC,LO,GO,NO>::writeXmfVariables(){
648
649 std::string centerString;
650 if (FEType_=="P0") {
651 centerString = "Cell";
652 }
653 else if (FEType_=="P1-disc") {
654
655 }
656 else if (FEType_=="P1") {
657 centerString = "Node";
658 }
659
660 else if (FEType_=="P2"){
661 centerString = "Node";
662 }
663 else if (FEType_=="Q2"){
664 centerString = "Node";
665 }
666
667 if(verbose_){
668 xmf_out_.seekp (closingLinesPosition_);
669
670 for (int i=0; i<varNames_.size(); i++) {
671 int dof = varDofPerNode_.at(i);
672 if (dof == 2) {
673 dof=3; // ParaView always uses 3D Data. for 2D Data the last entries (for the 3rd Dim) are all zero.
674 }
675 xmf_out_ <<
676 "\n <Attribute\n" <<
677 " Type=\"" << varTypes_.at(i) << "\"\n" <<
678 " Center=\"" << centerString << "\"\n" <<
679 " Name=\"" << varNames_.at(i) << "\">\n";
680
681
682 xmf_out_ <<
683 " <DataStructure ItemType=\"HyperSlab\"\n" <<
684 " Dimensions=\"" << nmbExportValuesGlob_ << " " << dof << "\"\n" <<
685 " Type=\"HyperSlab\">\n" <<
686 " <DataStructure Dimensions=\"3 2\"\n" <<
687 " Format=\"XML\">\n" <<
688 " 0 0\n" <<
689 " 1 1\n" <<
690 " " << nmbExportValuesGlob_ << " " << dof << "\n" <<
691 " </DataStructure>\n" <<
692
693 " <DataStructure Format=\"HDF\"\n" <<
694 " Dimensions=\"" << nmbExportValuesGlob_ << " " << dof << "\"\n" <<
695 " DataType=\"Float\"\n" <<
696 " Precision=\"8\">\n" <<
697 " " << outputFilename_ << ":/" << varNames_.at(i)
698 << postfix_ << "/Values\n" << // see also in writeVector/scalar
699 " </DataStructure>\n" <<
700 " </DataStructure>\n";
701
702
703 xmf_out_ <<
704 " </Attribute>\n";
705
706 }
707
708 xmf_out_ << "\n"
709 " </Grid>\n\n";
710 closingLinesPosition_ = xmf_out_.tellp();
711 // write closing lines
712 xmf_out_ << closingLines_;
713 xmf_out_.flush();
714 }
715}
716template<class SC,class LO,class GO,class NO>
717void ExporterParaView<SC,LO,GO,NO>::writeXmfTime(double time, double dt){
718
719 if (comm_->getRank()==0) {
720
721 xmf_times_out_.seekp (closingLinesPositionTimes_);
722
723 xmf_times_out_ <<
724 "<!-- Time " << time << " Iteration " << postfix_.substr (1, 5) << " -->\n" <<
725 " <Grid Name=\"Mesh Times " << time << "\">\n" <<
726 " <Time TimeType=\"Single\" Value=\"" << time << "\" />\n";
727 xmf_times_out_ <<
728 " <Topology\n"
729 << " Type=\"Polyvertex\"\n"
730 << " NumberOfElements=\""
731 << 2
732 << "\"\n"
733 << " BaseOffset=\""
734 << 0
735 << "\">\n"
736 << " <DataStructure Format=\"XML\"\n"
737 << " Dimensions=\"2\"\n"
738 << " DataType=\"Int\"\n"
739 << " Precision=\"8\">\n"
740 << " 0 1 \n"
741 << " </DataStructure>\n" << " </Topology>\n";
742
743
744 xmf_times_out_ <<
745 " <Geometry Type=\"X_Y\">\n" <<
746 " <DataStructure Format=\"XML\"\n" <<
747 " Dimensions=\"" << 2 << "\"\n" <<
748 " DataType=\"Float\"\n" <<
749 " Precision=\"8\">\n" <<
750 " " << "0.0 1.0 \n" <<
751 " </DataStructure>\n" <<
752 " <DataStructure Format=\"XML\"\n" <<
753 " Dimensions=\"" << 2 << "\"\n" <<
754 " DataType=\"Float\"\n" <<
755 " Precision=\"8\">\n" <<
756 " " << "0.0 1.0 \n" <<
757 " </DataStructure>\n" << " </Geometry>\n" <<
758 "\n";
759
760 xmf_times_out_ <<
761 "\n <Attribute\n" <<
762 " Type=\"" << "Scalar" << "\"\n" <<
763 " Center=\"" << "Node" << "\"\n" <<
764 " Name=\"" << "t_dt" << "\">\n";
765
766
767 xmf_times_out_ <<
768 " <DataStructure Format=\"XML\"\n" <<
769 " Dimensions=\"" << 2 << "\"\n" <<
770 " DataType=\"Float\"\n" <<
771 " Precision=\"8\">\n" <<
772 " " << time << " " << dt << "\n" <<
773 " </DataStructure>\n";
774
775 xmf_times_out_ <<
776 " </Attribute>\n";
777
778 xmf_times_out_ << "\n"
779 " </Grid>\n\n";
780 closingLinesPositionTimes_ = xmf_times_out_.tellp();
781 // write closing lines
782 xmf_times_out_ << closingLines_;
783 xmf_times_out_.flush();
784 }
785
786}
787
788template<class SC,class LO,class GO,class NO>
789void ExporterParaView<SC,LO,GO,NO>::prepareVectorField(MultiVecConstPtr_Type &u,
790 MultiVectorPtr_Type &u_export,
791 int dof) const{
792
793 TEUCHOS_TEST_FOR_EXCEPTION(u->getNumVectors()!=1, std::logic_error, "Can only export single vector");
794
795 Teuchos::ArrayRCP<const SC> tmpData = u->getData(0);
796 for (int i=0; i<(u->getLocalLength()/dof); i++) {
797 for (int j=0; j < dof; j++) {
798 u_export->getDataNonConst(j)[i] = tmpData[ dof * i + j ];
799 }
800 }
801}
802
803template<class SC,class LO,class GO,class NO>
804void ExporterParaView<SC,LO,GO,NO>::prepareScalar(MultiVecConstPtr_Type &u,
805 MultiVectorPtr_Type &u_export) const{
806
807 TEUCHOS_TEST_FOR_EXCEPTION(u->getNumVectors()!=1, std::logic_error, "Can only export single vector");
808
809 Teuchos::ArrayRCP<const SC> tmpData = u->getData(0);
810 for (int i=0; i<tmpData.size(); i++) {
811 u_export->getDataNonConst(0)[i] = tmpData[ i];
812 }
813
814}
815
816template<class SC,class LO,class GO,class NO>
817void ExporterParaView<SC,LO,GO,NO>::makePostfix(){
818
819 int postfixLength = 5;
820 std::ostringstream index;
821 index.fill ( '0' );
822
823 if (timeIndex_ % saveTimestep_ == 0){
824 index << std::setw (postfixLength) << ( timeIndex_ / saveTimestep_ );
825 postfix_ = "." + index.str();
826 }
827}
828}
829#endif
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36