1#ifndef ExporterParaView_DEF_hpp
2#define ExporterParaView_DEF_hpp
14template<
class SC,
class LO,
class GO,
class NO>
15ExporterParaView<SC,LO,GO,NO>::ExporterParaView():
18closingLinesPosition_(),
19closingLinesPositionTimes_(),
37nmbExportValuesGlob_(0),
38nmbPointsPerElement_(0),
49template<
class SC,
class LO,
class GO,
class NO>
50void ExporterParaView<SC,LO,GO,NO>::setup(std::string filename,
54 ParameterListPtr_Type parameterList){
57 setup( filename, mesh, FEType, parameterList);
58 saveTimestep_ = saveTimestep;
61template<
class SC,
class LO,
class GO,
class NO>
62void ExporterParaView<SC,LO,GO,NO>::setup(std::string filename,
65 ParameterListPtr_Type parameterList){
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() );
73 outputFilename_ = filename_ +
".h5";
75 dim_ = mesh->getDimension();
76 nmbElementsGlob_ = mesh->getNumElementsGlobal();
78 pointsUnique_ = mesh->getPointsUnique();
80 if (FEType_ ==
"P0") {
82 nmbPointsPerElement_ = 3 ;
86 nmbPointsPerElement_ = 4 ;
90 else if (FEType_ ==
"P1") {
92 nmbPointsPerElement_ = 3 ;
96 nmbPointsPerElement_ = 4 ;
100 else if (FEType_ ==
"P1-disc") {
102 nmbPointsPerElement_ = 3 ;
106 nmbPointsPerElement_ = 4 ;
110 else if(FEType_ ==
"P2"){
113 nmbPointsPerElement_ = 6 ;
117 nmbPointsPerElement_ = 10 ;
120 else if(FEType_ ==
"P2-CR"){
122 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Wrong dimension for P2-CR.");
125 nmbPointsPerElement_ = 10 ;
128 else if(FEType_ ==
"Q1"){
130 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Wrong dimension for Q1.");
133 nmbPointsPerElement_ = 8 ;
137 else if(FEType_ ==
"Q2"){
139 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Wrong dimension for Q2.");
142 nmbPointsPerElement_ = 20 ;
145 else if(FEType_ ==
"Q2-20"){
147 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Wrong dimension for Q2-20.");
150 nmbPointsPerElement_ = 20 ;
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");
157 nmbPointsGlob_ = mesh->getMapUnique()->getGlobalNumElements();
159 closingLines_ =
"\n </Grid>\n\n </Domain>\n</Xdmf>\n";
169 MapConstPtr_Type elementMap = mesh->getElementMap();
170 Teuchos::ArrayView<const GO> nodeElementList = elementMap->getNodeElementList();
171 vec_GO_Type nodeElementListInteger( nmbPointsPerElement_ * nodeElementList.size() );
173 for (
int i=0; i<nodeElementList.size(); i++) {
174 for (
int j=0; j<nmbPointsPerElement_; j++){
175 nodeElementListInteger[counter] = (int) nmbPointsPerElement_*nodeElementList[i] + j;
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_));
183 elementsHDF_.reset(
new MultiVector_Type(mapElements,1));
185 ElementsPtr_Type elements = mesh->getElementsC();
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;
196 pointsHDF_.reset(
new MultiVector_Type(mesh->getMapUnique(),dim_));
209template<
class SC,
class LO,
class GO,
class NO>
210void ExporterParaView<SC,LO,GO,NO>::addVariable(MultiVecConstPtr_Type &u,
214 MapConstPtrConst_Type& mapUnique){
216 variables_.push_back(u);
217 varNames_.push_back(varName);
218 varTypes_.push_back(varType);
219 varDofPerNode_.push_back(dofPerNode);
221 nmbExportValuesGlob_ = mapUnique->getGlobalNumElements();
223 uniqueMaps_.push_back(mapUnique);
227template<
class SC,
class LO,
class GO,
class NO>
228void ExporterParaView<SC,LO,GO,NO>::save(
double time){
230 if (timeIndex_ % saveTimestep_ == 0) {
233 if (!parameterList_.is_null()){
234 if ( parameterList_->sublist(
"Exporter").get(
"Write new mesh",
false) )
235 writeMeshPointsHDF5();
238 writeVariablesHDF5();
244 std::cout <<
"\n \t ### Export criterion not satisfied for current time step - no export of time step! ###" << std::endl;
251template<
class SC,
class LO,
class GO,
class NO>
252void ExporterParaView<SC,LO,GO,NO>::save(
double time,
double dt){
254 if (timeIndex_ % saveTimestep_ == 0) {
257 if (!parameterList_.is_null()){
258 if ( parameterList_->sublist(
"Exporter").get(
"Write new mesh",
false) )
259 writeMeshPointsHDF5();
261 writeVariablesHDF5();
268 writeXmfTime(time, dt);
272 std::cout <<
"\n \t ### Export criterion not satisfied for current time step - no export of time step! ###" << std::endl;
279template<
class SC,
class LO,
class GO,
class NO>
280void ExporterParaView<SC,LO,GO,NO>::closeExporter(){
282 hdf5exporter_->close();
285 xmf_times_out_.close();
290template<
class SC,
class LO,
class GO,
class NO>
291void ExporterParaView<SC,LO,GO,NO>::initHDF5(){
293 hdf5exporter_.reset(
new HDF5_Type(comm_) );
294 hdf5exporter_->create(outputFilename_);
295 std::string nameConn =
"Connections";
297 writeMeshElements(nameConn);
299 if (parameterList_.is_null())
300 writeMeshPointsHDF5();
301 else if ( parameterList_->sublist(
"Exporter").get(
"Write new mesh",
false) ==
false )
302 writeMeshPointsHDF5();
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) ) {
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_);
315 writeMeshPoints(nameP_X, nameP_Y, nameP_Z );
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 );
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 );
333template<
class SC,
class LO,
class GO,
class NO>
334void ExporterParaView<SC,LO,GO,NO>::writeMeshElements( std::string nameConn ){
336 hdf5exporter_->createGroup(nameConn);
338 hdf5exporter_->write(nameConn,elementsHDF_);
341template<
class SC,
class LO,
class GO,
class NO>
342void ExporterParaView<SC,LO,GO,NO>::writeMeshPoints(std::string nameP_X,
344 std::string nameP_Z){
348 hdf5exporter_->createGroup(nameP_X);
349 hdf5exporter_->createGroup(nameP_Y);
351 hdf5exporter_->write(nameP_X,pointsHDF_->getVector(0));
352 hdf5exporter_->write(nameP_Y,pointsHDF_->getVector(1));
355 hdf5exporter_->createGroup(nameP_Z);
356 hdf5exporter_->write(nameP_Z,pointsHDF_->getVector(2));
358 hdf5exporter_->flush();
362template<
class SC,
class LO,
class GO,
class NO>
363void ExporterParaView<SC,LO,GO,NO>::updatePoints(){
365 if (pointsUnique_->size()>0)
366 dim = pointsUnique_->at(0).size();
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];
375template<
class SC,
class LO,
class GO,
class NO>
376void ExporterParaView<SC,LO,GO,NO>::writeVariablesHDF5(){
378 for (
int i=0; i<variables_.size(); i++) {
380 if (varTypes_.at(i)==
"Vector") {
381 MultiVectorPtr_Type u_export(
new MultiVector_Type(uniqueMaps_.at(i),3));
383 hdf5exporter_->createGroup(varNames_[i]+postfix_);
385 prepareVectorField(variables_.at(i),u_export, varDofPerNode_.at(i));
387 hdf5exporter_->write(varNames_[i]+postfix_,u_export,
true);
389 else if(varTypes_.at(i)==
"Scalar"){
390 MultiVectorPtr_Type u_export(
new MultiVector_Type(uniqueMaps_.at(i),1));
392 hdf5exporter_->createGroup(varNames_[i]+postfix_);
394 prepareScalar(variables_.at(i),u_export);
396 hdf5exporter_->write(varNames_[i]+postfix_,u_export);
398 hdf5exporter_->flush();
403template<
class SC,
class LO,
class GO,
class NO>
404void ExporterParaView<SC,LO,GO,NO>::initXmf(){
406 if (comm_->getRank()==0) {
408 xmf_out_.open((filename_ +
".xmf").c_str(),std::ios_base::out);
409 xmf_out_ <<
"<?xml version=\"1.0\" ?>\n"
410 <<
"<!DOCTYPE Xdmf SYSTEM \""
413 <<
"<!ENTITY DataFile \""
421 <<
" <Domain Name=\""
426 <<
"Grid\" GridType=\"Collection\" CollectionType=\"Temporal\">\n"
429 closingLinesPosition_ = xmf_out_.tellp();
431 xmf_out_ << closingLines_;
437template<
class SC,
class LO,
class GO,
class NO>
438void ExporterParaView<SC,LO,GO,NO>::initXmfTimes(){
440 if (comm_->getRank()==0) {
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 \""
448 <<
" <Domain Name=\""
453 <<
"Grid\" GridType=\"Collection\" CollectionType=\"Temporal\">\n"
456 closingLinesPositionTimes_ = xmf_times_out_.tellp();
458 xmf_times_out_ << closingLines_;
459 xmf_times_out_.flush();
463template<
class SC,
class LO,
class GO,
class NO>
464void ExporterParaView<SC,LO,GO,NO>::writeXmf(
double time){
469 std::string nameConn =
"Connections";
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_);
489 writeXmfElements( nameConn, time);
490 writeXmfPoints( nameP_X, nameP_Y, nameP_Z );
495template<
class SC,
class LO,
class GO,
class NO>
496void ExporterParaView<SC,LO,GO,NO>::writeXmfElements( std::string nameConn,
double time ){
499 xmf_out_.seekp (closingLinesPosition_);
502 "<!-- Time " << time <<
" Iteration " << postfix_.substr (1, 5) <<
" -->\n" <<
503 " <Grid Name=\"Mesh" << FEType_ <<
" " << time <<
"\">\n" <<
504 " <Time TimeType=\"Single\" Value=\"" << time <<
"\" />\n";
507 std::string FEstring;
510 FEstring =
"Triangle";
513 FEstring =
"Tetrahedron";
516 else if (FEType_==
"P1-disc") {
518 FEstring =
"Triangle";
521 FEstring =
"Tetrahedron";
524 else if (FEType_==
"P1") {
526 FEstring =
"Triangle";
529 FEstring =
"Tetrahedron";
533 else if (FEType_==
"P2"){
541 else if (FEType_==
"P2-CR"){
547 else if (FEType_==
"Q1"){
550 FEstring =
"Hexahedron";
553 else if (FEType_==
"Q2"){
559 else if (FEType_==
"Q2-20"){
566 xmf_out_ <<
" <Topology\n"
570 <<
" NumberOfElements=\""
576 <<
" <DataStructure Format=\"HDF\"\n"
580 << nmbPointsPerElement_
581 <<
"\"\n" <<
" DataType=\"Int\"\n"
582 <<
" Precision=\"8\">\n" <<
" "
583 << outputFilename_ <<
":/"<< nameConn <<
"/Values\n"
584 <<
" </DataStructure>\n" <<
" </Topology>\n";
586 closingLinesPosition_ = xmf_out_.tellp();
590template<
class SC,
class LO,
class GO,
class NO>
591void ExporterParaView<SC,LO,GO,NO>::writeXmfPoints(std::string nameP_X,
593 std::string nameP_Z ){
596 xmf_out_.seekp (closingLinesPosition_);
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" <<
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" <<
640 closingLinesPosition_ = xmf_out_.tellp();
646template<
class SC,
class LO,
class GO,
class NO>
647void ExporterParaView<SC,LO,GO,NO>::writeXmfVariables(){
649 std::string centerString;
651 centerString =
"Cell";
653 else if (FEType_==
"P1-disc") {
656 else if (FEType_==
"P1") {
657 centerString =
"Node";
660 else if (FEType_==
"P2"){
661 centerString =
"Node";
663 else if (FEType_==
"Q2"){
664 centerString =
"Node";
668 xmf_out_.seekp (closingLinesPosition_);
670 for (
int i=0; i<varNames_.size(); i++) {
671 int dof = varDofPerNode_.at(i);
677 " Type=\"" << varTypes_.at(i) <<
"\"\n" <<
678 " Center=\"" << centerString <<
"\"\n" <<
679 " Name=\"" << varNames_.at(i) <<
"\">\n";
683 " <DataStructure ItemType=\"HyperSlab\"\n" <<
684 " Dimensions=\"" << nmbExportValuesGlob_ <<
" " << dof <<
"\"\n" <<
685 " Type=\"HyperSlab\">\n" <<
686 " <DataStructure Dimensions=\"3 2\"\n" <<
687 " Format=\"XML\">\n" <<
690 " " << nmbExportValuesGlob_ <<
" " << dof <<
"\n" <<
691 " </DataStructure>\n" <<
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" <<
699 " </DataStructure>\n" <<
700 " </DataStructure>\n";
710 closingLinesPosition_ = xmf_out_.tellp();
712 xmf_out_ << closingLines_;
716template<
class SC,
class LO,
class GO,
class NO>
717void ExporterParaView<SC,LO,GO,NO>::writeXmfTime(
double time,
double dt){
719 if (comm_->getRank()==0) {
721 xmf_times_out_.seekp (closingLinesPositionTimes_);
724 "<!-- Time " << time <<
" Iteration " << postfix_.substr (1, 5) <<
" -->\n" <<
725 " <Grid Name=\"Mesh Times " << time <<
"\">\n" <<
726 " <Time TimeType=\"Single\" Value=\"" << time <<
"\" />\n";
729 <<
" Type=\"Polyvertex\"\n"
730 <<
" NumberOfElements=\""
736 <<
" <DataStructure Format=\"XML\"\n"
737 <<
" Dimensions=\"2\"\n"
738 <<
" DataType=\"Int\"\n"
739 <<
" Precision=\"8\">\n"
741 <<
" </DataStructure>\n" <<
" </Topology>\n";
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" <<
762 " Type=\"" <<
"Scalar" <<
"\"\n" <<
763 " Center=\"" <<
"Node" <<
"\"\n" <<
764 " Name=\"" <<
"t_dt" <<
"\">\n";
768 " <DataStructure Format=\"XML\"\n" <<
769 " Dimensions=\"" << 2 <<
"\"\n" <<
770 " DataType=\"Float\"\n" <<
771 " Precision=\"8\">\n" <<
772 " " << time <<
" " << dt <<
"\n" <<
773 " </DataStructure>\n";
778 xmf_times_out_ <<
"\n"
780 closingLinesPositionTimes_ = xmf_times_out_.tellp();
782 xmf_times_out_ << closingLines_;
783 xmf_times_out_.flush();
788template<
class SC,
class LO,
class GO,
class NO>
789void ExporterParaView<SC,LO,GO,NO>::prepareVectorField(MultiVecConstPtr_Type &u,
790 MultiVectorPtr_Type &u_export,
793 TEUCHOS_TEST_FOR_EXCEPTION(u->getNumVectors()!=1, std::logic_error,
"Can only export single vector");
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 ];
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{
807 TEUCHOS_TEST_FOR_EXCEPTION(u->getNumVectors()!=1, std::logic_error,
"Can only export single vector");
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];
816template<
class SC,
class LO,
class GO,
class NO>
817void ExporterParaView<SC,LO,GO,NO>::makePostfix(){
819 int postfixLength = 5;
820 std::ostringstream index;
823 if (timeIndex_ % saveTimestep_ == 0){
824 index << std::setw (postfixLength) << ( timeIndex_ / saveTimestep_ );
825 postfix_ =
"." + index.str();
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36