MueLu  Version of the Day
MueLu_Utilities.cpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #include "MueLu_Utilities_def.hpp"
47 
48 #include <string>
49 
50 #ifdef HAVE_MUELU_EPETRAEXT
52 #endif
53 
54 namespace MueLu {
55 
56  RCP<Xpetra::Matrix<double, int, int> > Utils2<double, int, int>::Transpose(Matrix& Op, bool optimizeTranspose, const std::string & label) {
57  typedef double Scalar;
58  typedef int LocalOrdinal;
59  typedef int GlobalOrdinal;
60  typedef KokkosClassic::DefaultNode::DefaultNodeType Node;
61 
62 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
63  std::string TorE = "epetra";
64 #else
65  std::string TorE = "tpetra";
66 #endif
67 
68 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
69  try {
71  }
72  catch (...) {
73  TorE = "tpetra";
74  }
75 #endif
76 
77 #ifdef HAVE_MUELU_TPETRA
78  if (TorE == "tpetra") {
79  try {
80  const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = Utils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraCrs(Op);
81 
82  // Compute the transpose A of the Tpetra matrix tpetraOp.
83  RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A;
84  Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp),label);
85  A = transposer.createTranspose();
86  RCP<Xpetra::TpetraCrsMatrix<SC> > AA = rcp(new Xpetra::TpetraCrsMatrix<SC>(A));
87  RCP<Xpetra::CrsMatrix<SC> > AAA = rcp_implicit_cast<Xpetra::CrsMatrix<SC> >(AA);
88  RCP<Xpetra::CrsMatrixWrap<SC> > AAAA = rcp( new Xpetra::CrsMatrixWrap<SC> (AAA));
89 
90  return AAAA;
91  }
92  catch (std::exception& e) {
93  std::cout << "threw exception '" << e.what() << "'" << std::endl;
94  throw Exceptions::RuntimeError("Utils::Transpose failed, perhaps because matrix is not a Crs matrix");
95  }
96  } //if
97 #endif
98 
99  if (TorE == "tpetra") {
100 #ifdef HAVE_MUELU_TPETRA
101 #else
102  throw Exceptions::RuntimeError("Tpetra");
103 #endif // HAVE_MUELU_TPETRA
104 
105  } else {
106 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
107  Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer("ZZ Entire Transpose"));
108  // Epetra case
111  Epetra_CrsMatrix * A = dynamic_cast<Epetra_CrsMatrix*>(&transposer(epetraOp));
112  transposer.ReleaseTranspose(); // So we can keep A in Muelu...
113 
114  RCP<Epetra_CrsMatrix> rcpA(A);
115  RCP<EpetraCrsMatrix> AA = rcp(new EpetraCrsMatrix(rcpA));
116  RCP<Xpetra::CrsMatrix<SC> > AAA = rcp_implicit_cast<Xpetra::CrsMatrix<SC> >(AA);
117  RCP<Xpetra::CrsMatrixWrap<SC> > AAAA = rcp( new Xpetra::CrsMatrixWrap<SC>(AAA));
118  AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
119 
120  return AAAA;
121 #else
122  throw Exceptions::RuntimeError("Epetra (Err. 2)");
123 #endif
124  }
125  return Teuchos::null;
126  } //Transpose
127 
128  // -- ------------------------------------------------------- --
129 
130  void Utils2<double,int,int>::MyOldScaleMatrix_Epetra(Matrix& Op, const Teuchos::ArrayRCP<SC>& scalingVector,
131  bool doFillComplete,
132  bool doOptimizeStorage) {
133 #ifdef HAVE_MUELU_EPETRA
134  try {
135  const Epetra_CrsMatrix& epOp = Utils<double,int,int>::Op2NonConstEpetraCrs(Op);
136 
137  Epetra_Map const &rowMap = epOp.RowMap();
138  int nnz;
139  double *vals;
140  int *cols;
141 
142  for (int i = 0; i < rowMap.NumMyElements(); ++i) {
143  epOp.ExtractMyRowView(i, nnz, vals, cols);
144  for (int j = 0; j < nnz; ++j)
145  vals[j] *= scalingVector[i];
146  }
147 
148  } catch (...){
149  throw Exceptions::RuntimeError("Only Epetra_CrsMatrix types can be scaled");
150  }
151 #else
152  throw Exceptions::RuntimeError("Matrix scaling is not possible because Epetra has not been enabled.");
153 #endif // HAVE_MUELU_EPETRA
154  } //Utils2::MyOldScaleMatrix_Epetra()
155 
156  // -- ------------------------------------------------------- --
157 
158  RCP<Xpetra::MultiVector<double,int,int> > Utils2<double,int,int>::ReadMultiVector(const std::string& fileName, const RCP<const Map>& map) {
159  Xpetra::UnderlyingLib lib = map->lib();
160 
161  if (lib == Xpetra::UseEpetra) {
162 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
163  Epetra_MultiVector * MV;
164  EpetraExt::MatrixMarketFileToMultiVector(fileName.c_str(), toEpetra(map), MV);
165  return Xpetra::toXpetra<int>(rcp(MV));
166 #else
167  throw Exceptions::RuntimeError("MueLu has not been compiled with Epetra and EpetraExt support.");
168 #endif
169  } else if (lib == Xpetra::UseTpetra) {
170 #ifdef HAVE_MUELU_TPETRA
171  typedef Tpetra::CrsMatrix<SC,LO,GO,NO> sparse_matrix_type;
172  typedef Tpetra::MatrixMarket::Reader<sparse_matrix_type> reader_type;
173  typedef Tpetra::Map<LO,GO,NO> map_type;
174  typedef Tpetra::MultiVector<SC,LO,GO,NO> multivector_type;
175 
176  RCP<const map_type> temp = toTpetra(map);
177  RCP<multivector_type> TMV = reader_type::readDenseFile(fileName,map->getComm(),map->getNode(),temp);
178  RCP<MultiVector> rmv = Xpetra::toXpetra(TMV);
179  return rmv;
180 #else
181  throw Exceptions::RuntimeError("MueLu has not been compiled with Tpetra support.");
182 #endif
183  } else {
184  throw Exceptions::RuntimeError("Utils::Read : you must specify Xpetra::UseEpetra or Xpetra::UseTpetra.");
185  }
186 
187  return Teuchos::null;
188  }
189 
190  RCP<const Xpetra::Map<int,int> > Utils2<double,int,int>::ReadMap(const std::string& fileName, Xpetra::UnderlyingLib lib, const RCP<const Teuchos::Comm<int> >& comm) {
191  if (lib == Xpetra::UseEpetra) {
192 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
193  Epetra_Map *eMap;
194  int rv = EpetraExt::MatrixMarketFileToMap(fileName.c_str(), *(Xpetra::toEpetra(comm)), eMap);
195  if (rv != 0)
196  throw Exceptions::RuntimeError("Error reading matrix with EpetraExt::MatrixMarketToMap (returned " + toString(rv) + ")");
197 
198  RCP<Epetra_Map> eMap1 = rcp(new Epetra_Map(*eMap));
199  return Xpetra::toXpetra<int>(*eMap1);
200 #else
201  throw Exceptions::RuntimeError("MueLu has not been compiled with Epetra and EpetraExt support.");
202 #endif
203  } else if (lib == Xpetra::UseTpetra) {
204 #ifdef HAVE_MUELU_TPETRA
205  typedef Tpetra::CrsMatrix<double,int,int,NO> sparse_matrix_type;
206  typedef Tpetra::MatrixMarket::Reader<sparse_matrix_type> reader_type;
207 
208  RCP<NO> node = rcp(new NO());
209 
210  RCP<const Tpetra::Map<int,int,NO> > tMap = reader_type::readMapFile(fileName, comm, node);
211  if (tMap.is_null())
212  throw Exceptions::RuntimeError("The Tpetra::Map returned from readSparseFile() is null.");
213 
214  return Xpetra::toXpetra(tMap);
215 #else
216  throw Exceptions::RuntimeError("MueLu has not been compiled with Tpetra support.");
217 #endif
218  } else {
219  throw Exceptions::RuntimeError("Utils::Read : you must specify Xpetra::UseEpetra or Xpetra::UseTpetra.");
220  }
221  }
222 
223 
224  /* Removes the following non-serializable data (A,P,R,Nullspace,Coordinates)
225  from level-specific sublists from inList
226  and moves it to nonSerialList. Everything else is copied to serialList.
227  This function returns the level number of the highest level for which
228  non-serializable data was provided.
229  */
230  long ExtractNonSerializableData(const Teuchos::ParameterList& inList, Teuchos::ParameterList& serialList, Teuchos::ParameterList& nonSerialList) {
231  using Teuchos::ParameterList;
232 
233  ParameterList dummy;
234  long maxLevel = 0;
235 
236  for (ParameterList::ConstIterator it = inList.begin(); it != inList.end(); it++) {
237  const std::string& levelName = it->first;
238 
239  // Check for mach of the form "level X" where X is a positive integer
240  if (inList.isSublist(levelName) && levelName.find("level ") == 0 && levelName.size() > 6) {
241  int levelID = strtol(levelName.substr(6).c_str(), 0, 0);
242  if (maxLevel < levelID)
243  maxLevel = levelID;
244 
245  // Split the sublist
246  const ParameterList& levelList = inList.sublist(levelName);
247  for (ParameterList::ConstIterator it2 = levelList.begin(); it2 != levelList.end(); it2++) {
248  const std::string& name = it2->first;
249  if (name == "A" || name == "P" || name == "R" || name == "Nullspace" || name == "Coordinates")
250  nonSerialList.sublist(levelName).setEntry(name, it2->second);
251  #ifdef HAVE_MUELU_MATLAB
252  else if(IsParamMuemexVariable(name))
253  {
254  nonSerialList.sublist(levelName).setEntry(name, it2->second);
255  }
256  #endif
257  else
258  serialList.sublist(levelName).setEntry(name, it2->second);
259  }
260 
261  } else {
262  serialList.setEntry(it->first, it->second);
263  }
264  }
265 
266  return maxLevel;
267  }
268 
269  void TokenizeStringAndStripWhiteSpace(const std::string& stream, std::vector<std::string>& tokenList, const char* delimChars)
270  {
271  //note: default delimiter string is ","
272  // Take a comma-separated list and tokenize it, stripping out leading & trailing whitespace. Then add to tokenList
273  char* buf = (char*) malloc(stream.size() + 1);
274  strcpy(buf, stream.c_str());
275  char* token = strtok(buf, delimChars);
276  if(token == NULL)
277  {
278  free(buf);
279  return;
280  }
281  while(token)
282  {
283  //token points to start of string to add to tokenList
284  //remove front whitespace...
285  char* tokStart = token;
286  char* tokEnd = token + strlen(token) - 1;
287  while(*tokStart == ' ' && tokStart < tokEnd)
288  tokStart++;
289  while(*tokEnd == ' ' && tokStart < tokEnd)
290  tokEnd--;
291  tokEnd++;
292  if(tokStart < tokEnd)
293  {
294  std::string finishedToken(tokStart, tokEnd - tokStart); //use the constructor that takes a certain # of chars
295  tokenList.push_back(finishedToken);
296  }
297  token = strtok(NULL, delimChars);
298  }
299  free(buf);
300  }
301 
302  bool IsParamMuemexVariable(const std::string& name)
303  {
304  //see if paramName is exactly two "words" - like "OrdinalVector myNullspace" or something
305  char* str = (char*) malloc(name.length() + 1);
306  strcpy(str, name.c_str());
307  //Strip leading and trailing whitespace
308  char* firstWord = strtok(str, " ");
309  if(!firstWord)
310  return false;
311  char* secondWord = strtok(NULL, " ");
312  if(!secondWord)
313  return false;
314  char* thirdWord = strtok(NULL, " ");
315  if(thirdWord)
316  return false;
317  //convert first word to all lowercase for case insensitive compare
318  char* tolowerIt = firstWord;
319  while(*tolowerIt)
320  {
321  *tolowerIt = (char) tolower(*tolowerIt);
322  tolowerIt++;
323  }
324  //See if the first word is one of the custom variable names
325  if(strstr(firstWord, "matrix") ||
326  strstr(firstWord, "multivector") ||
327  strstr(firstWord, "map") ||
328  strstr(firstWord, "ordinalvector") ||
329  strstr(firstWord, "int") ||
330  strstr(firstWord, "scalar") ||
331  strstr(firstWord, "double") ||
332  strstr(firstWord, "complex") ||
333  strstr(firstWord, "string"))
334  //Add name to list of keys to remove
335  {
336  free(str);
337  return true;
338  }
339  else
340  {
341  free(str);
342  return false;
343  }
344  }
345 
346 } // namespace MueLu
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
int MatrixMarketFileToMultiVector(const char *filename, const Epetra_BlockMap &map, Epetra_MultiVector *&A)
static RCP< Matrix > Transpose(Matrix &Op, bool optimizeTranspose=false, const std::string &label=std::string())
Transpose a Xpetra::Matrix.
KokkosClassic::DefaultNode::DefaultNodeType NO
Namespace for MueLu classes and methods.
bool IsParamMuemexVariable(const std::string &name)
static void MyOldScaleMatrix_Epetra(Matrix &Op, const Teuchos::ArrayRCP< SC > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
Scale an Epetra matrix.
int MatrixMarketFileToMap(const char *filename, const Epetra_Comm &comm, Epetra_Map *&map)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Matrix > Op)
static RCP< const Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2TpetraCrs(RCP< const Matrix > Op)
Xpetra::Matrix< double, int, int, NO > Matrix
static RCP< const Map > ReadMap(const std::string &fileName, Xpetra::UnderlyingLib lib, const RCP< const Teuchos::Comm< int > > &comm)
static RCP< MultiVector > ReadMultiVector(const std::string &fileName, const RCP< const Map > &map)
Exception throws to report errors in the internal logical of the program.
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Matrix > Op)
void TokenizeStringAndStripWhiteSpace(const std::string &stream, std::vector< std::string > &tokenList, const char *delimChars)
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)