46 #include "Epetra_MpiComm.h" 
   49 #include "Epetra_SerialComm.h" 
   52 #include "Epetra_Map.h" 
   53 #include "Epetra_CrsGraph.h" 
   54 #include "Epetra_CrsMatrix.h" 
   55 #include "Epetra_MultiVector.h" 
   56 #include "Teuchos_ParameterList.hpp" 
   57 #include "Teuchos_XMLParameterListWriter.hpp" 
   58 #include "Teuchos_Assert.hpp" 
   60 using namespace Teuchos;
 
   74   if (Comm_.MyPID() == 0) 
 
   76     std::ofstream of(FileName_.c_str());
 
   77     of << 
"<ObjectCollection Label=\"" << Label << 
"\">" << std::endl;
 
   87   if (Comm_.MyPID() == 0) 
 
   89     std::ofstream of(FileName_.c_str(), std::ios::app);
 
   90     of << 
"</ObjectCollection>" << std::endl;
 
   99 Write(
const std::string& Label, 
const std::vector<std::string>& Content)
 
  101   TEUCHOS_TEST_FOR_EXCEPTION(IsOpen_ == 
false, std::logic_error,
 
  102                      "No file has been opened");
 
  104   if (Comm_.MyPID()) 
return;
 
  106   std::ofstream of(FileName_.c_str(), std::ios::app);
 
  108   of << 
"<Text Label=\"" << Label << 
"\">" << std::endl;
 
  109   int Csize = (int) Content.size();
 
  110   for (
int i = 0; i < Csize; ++i)
 
  111     of << Content[i] << std::endl;
 
  113   of << 
"</Text>" << std::endl;
 
  122   TEUCHOS_TEST_FOR_EXCEPTION(IsOpen_ == 
false, std::logic_error,
 
  123                      "No file has been opened");
 
  125   long long Rows = Matrix.NumGlobalRows64();
 
  126   long long Cols = Matrix.NumGlobalRows64();
 
  127   long long Nonzeros = Matrix.NumGlobalNonzeros64();
 
  129   if (Comm_.MyPID() == 0)
 
  131     std::ofstream of(FileName_.c_str(), std::ios::app);
 
  132     of << 
"<PointMatrix Label=\"" << Label << 
'"' 
  133       << 
" Rows=\"" << Rows << 
'"' 
  134       << 
" Columns=\"" << Cols<< 
'"' 
  135       << 
" Nonzeros=\"" << Nonzeros << 
'"' 
  136       << 
" Type=\"double\" StartingIndex=\"0\">" << std::endl;
 
  140   std::vector<int> Indices(Length);
 
  141   std::vector<double> Values(Length);
 
  143   for (
int iproc = 0; iproc < Comm_.NumProc(); iproc++)
 
  145     if (iproc == Comm_.MyPID())
 
  147       std::ofstream of(FileName_.c_str(), std::ios::app);
 
  150       for (
int i = 0; i < Matrix.
NumMyRows(); ++i)
 
  157         for (
int j = 0; j < NumMyEntries; ++j)
 
  159              << 
" " << std::setiosflags(std::ios::scientific) << Values[j] << std::endl;
 
  166   if (Comm_.MyPID() == 0)
 
  168     std::ofstream of(FileName_.c_str(), std::ios::app);
 
  169     of << 
"</PointMatrix>" << std::endl;
 
  178   TEUCHOS_TEST_FOR_EXCEPTION(IsOpen_ == 
false, std::logic_error,
 
  179                      "No file has been opened");
 
  181   long long Length = MultiVector.GlobalLength64();
 
  182   int NumVectors = MultiVector.NumVectors();
 
  184   if (Comm_.MyPID() == 0)
 
  186     std::ofstream of(FileName_.c_str(), std::ios::app);
 
  188     of << 
"<MultiVector Label=\"" << Label 
 
  189       << 
"\" Length=\"" << Length << 
'"' 
  190       << 
" NumVectors=\"" << NumVectors << 
'"' 
  191       << 
" Type=\"double\">" << std::endl;
 
  195   for (
int iproc = 0; iproc < Comm_.NumProc(); iproc++)
 
  197     if (iproc == Comm_.MyPID())
 
  199       std::ofstream of(FileName_.c_str(), std::ios::app);
 
  202       for (
int i = 0; i < MultiVector.MyLength(); ++i)
 
  204         for (
int j = 0; j < NumVectors; ++j)
 
  205           of << std::setiosflags(std::ios::scientific) << MultiVector[j][i] << 
" ";
 
  213   if (Comm_.MyPID() == 0)
 
  215     std::ofstream of(FileName_.c_str(), std::ios::app);
 
  216     of << 
"</MultiVector>" << std::endl;
 
  225   TEUCHOS_TEST_FOR_EXCEPTION(IsOpen_ == 
false, std::logic_error,
 
  226                      "No file has been opened");
 
  228   long long NumGlobalElements = Map.NumGlobalElements64();
 
  229   const int* MyGlobalElements_int = 0;
 
  230   const long long* MyGlobalElements_LL = 0;
 
  233   if(!MyGlobalElements_int || !MyGlobalElements_LL)
 
  234     throw "EpetraExt::XMLWriter::Write: ERROR, GlobalIndices type unknown.";
 
  236   if (Comm_.MyPID() == 0)
 
  238     std::ofstream of(FileName_.c_str(), std::ios::app);
 
  240     of << 
"<Map Label=\"" << Label 
 
  241       << 
"\" NumElements=\"" << NumGlobalElements << 
'"' 
  242       << 
" IndexBase=\"" << Map.IndexBase64() << 
'"' 
  243       << 
" NumProc=\"" << Comm_.NumProc() << 
'"';
 
  248   for (
int iproc = 0; iproc < Comm_.NumProc(); ++iproc)
 
  250     if (iproc == Comm_.MyPID())
 
  252       std::ofstream of(FileName_.c_str(), std::ios::app);
 
  254       of << 
" ElementsOnProc" << iproc << 
"=\"" << Map.
NumMyElements() << 
'"';
 
  260   if (Comm_.MyPID() == 0)
 
  262     std::ofstream of(FileName_.c_str(), std::ios::app);
 
  263     of << '>
' << std::endl; 
  267   for (int iproc = 0; iproc < Comm_.NumProc(); iproc++) 
  269     if (iproc == Comm_.MyPID()) 
  271       std::ofstream of(FileName_.c_str(), std::ios::app); 
  273       of << "<Proc ID=\"" << Comm_.MyPID() << "\">" << std::endl; 
  275       if(MyGlobalElements_int) 
  277         for (int i = 0; i < Map.NumMyElements(); ++i) 
  279           of << MyGlobalElements_int[i] << std::endl; 
  284         for (int i = 0; i < Map.NumMyElements(); ++i) 
  286           of << MyGlobalElements_LL[i] << std::endl; 
  290       of << "</Proc>" << std::endl; 
  296   if (Comm_.MyPID() == 0) 
  298     std::ofstream of(FileName_.c_str(), std::ios::app); 
  299     of << "</Map>" << std::endl; 
  304 // ============================================================================ 
  305 void EpetraExt::XMLWriter:: 
  306 Write(const std::string& Label, Teuchos::ParameterList& List) 
  308   TEUCHOS_TEST_FOR_EXCEPTION(IsOpen_ == false, std::logic_error, 
  309                      "No file has been opened"); 
  311   if (Comm_.MyPID()) return; 
  313   std::ofstream of(FileName_.c_str(), std::ios::app); 
  315   of << "<List Label=\"" << Label << "\">" << std::endl; 
  317   XMLParameterListWriter Writer; 
  318   XMLObject Obj = Writer.toXML(List); 
  320   of << Obj.toString(); 
  322   of << "</List>" << std::endl; 
virtual const Epetra_Map & RowMatrixRowMap() const =0
 
int MyGlobalElements(int *MyGlobalElementList) const 
 
XMLWriter(const Epetra_Comm &Comm, const std::string &FileName)
ctor 
 
int NumMyElements() const 
 
virtual int MaxNumEntries() const =0
 
void Create(const std::string &Label)
Creates the file, giving Label to the whole object. 
 
virtual int NumMyRows() const =0
 
virtual int ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const =0
 
virtual const Epetra_Map & RowMatrixColMap() const =0
 
void Close()
Closes the file. No Write operations can follow. 
 
void Write(const std::string &Label, const Epetra_Map &Map)
Writes an Epetra_Map using label Label.