MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_UtilitiesBase_decl.hpp
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 #ifndef MUELU_UTILITIESBASE_DECL_HPP
47 #define MUELU_UTILITIESBASE_DECL_HPP
48 
49 #ifndef _WIN32
50 #include <unistd.h> //necessary for "sleep" function in debugging methods (PauseForDebugging)
51 #endif
52 
53 #include <string>
54 
55 #include "MueLu_ConfigDefs.hpp"
56 
57 #include <Teuchos_DefaultComm.hpp>
58 #include <Teuchos_ScalarTraits.hpp>
60 
61 #include <Xpetra_BlockedCrsMatrix_fwd.hpp>
62 #include <Xpetra_CrsMatrix_fwd.hpp>
63 #include <Xpetra_CrsMatrixWrap_fwd.hpp>
64 #include <Xpetra_Map_fwd.hpp>
65 #include <Xpetra_BlockedMap_fwd.hpp>
66 #include <Xpetra_MapFactory_fwd.hpp>
67 #include <Xpetra_Matrix_fwd.hpp>
68 #include <Xpetra_MatrixFactory_fwd.hpp>
69 #include <Xpetra_MultiVector_fwd.hpp>
70 #include <Xpetra_MultiVectorFactory_fwd.hpp>
71 #include <Xpetra_Operator_fwd.hpp>
72 #include <Xpetra_Vector_fwd.hpp>
73 #include <Xpetra_BlockedMultiVector.hpp>
74 #include <Xpetra_BlockedVector.hpp>
75 #include <Xpetra_VectorFactory_fwd.hpp>
76 #include <Xpetra_ExportFactory.hpp>
77 
78 #include <Xpetra_Import.hpp>
79 #include <Xpetra_ImportFactory.hpp>
80 #include <Xpetra_MatrixMatrix.hpp>
81 #include <Xpetra_CrsMatrixWrap.hpp>
82 #include <Xpetra_StridedMap.hpp>
83 
84 #include "MueLu_Exceptions.hpp"
85 
86 
87 namespace MueLu {
88 
89 // MPI helpers
90 #define MueLu_sumAll(rcpComm, in, out) \
91  Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_SUM, in, Teuchos::outArg(out))
92 #define MueLu_minAll(rcpComm, in, out) \
93  Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MIN, in, Teuchos::outArg(out))
94 #define MueLu_maxAll(rcpComm, in, out) \
95  Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MAX, in, Teuchos::outArg(out))
96 
104  template <class Scalar,
107  class Node = DefaultNode>
108  class UtilitiesBase {
109  public:
110 #undef MUELU_UTILITIESBASE_SHORT
111 //#include "MueLu_UseShortNames.hpp"
112  private:
113  typedef Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> CrsMatrixWrap;
114  typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> CrsMatrix;
115  typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> Matrix;
116  typedef Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> Vector;
117  typedef Xpetra::BlockedVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> BlockedVector;
118  typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MultiVector;
119  typedef Xpetra::BlockedMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> BlockedMultiVector;
120  typedef Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node> BlockedMap;
121  typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> Map;
122  public:
124 
125 
127  if (Op.is_null())
128  return Teuchos::null;
129  return rcp(new CrsMatrixWrap(Op));
130  }
131 
139  size_t numRows = A.getRowMap()->getNodeNumElements();
140  Teuchos::ArrayRCP<Scalar> diag(numRows);
143  for (size_t i = 0; i < numRows; ++i) {
144  A.getLocalRowView(i, cols, vals);
145  LocalOrdinal j = 0;
146  for (; j < cols.size(); ++j) {
147  if (Teuchos::as<size_t>(cols[j]) == i) {
148  diag[i] = vals[j];
149  break;
150  }
151  }
152  if (j == cols.size()) {
153  // Diagonal entry is absent
155  }
156  }
157  return diag;
158  }
159 
167  Teuchos::TimeMonitor MM = *Teuchos::TimeMonitor::getNewTimer("UtilitiesBase::GetMatrixDiagonalInverse");
168 
169  RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
170 
171  RCP<const Map> rowMap = rcpA->getRowMap();
172  RCP<Vector> diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,true);
173 
174  rcpA->getLocalDiagCopy(*diag);
175 
177 
178  return inv;
179  }
180 
181 
182 
190 
191  size_t numRows = A.getRowMap()->getNodeNumElements();
192  Teuchos::ArrayRCP<Scalar> diag(numRows);
195  for (size_t i = 0; i < numRows; ++i) {
196  A.getLocalRowView(i, cols, vals);
198  for (LocalOrdinal j = 0; j < cols.size(); ++j) {
199  diag[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
200  }
201  }
202  return diag;
203  }
204 
212 
213  RCP<Vector> diag = Teuchos::null;
214 
216  Teuchos::rcp_dynamic_cast<const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(rcpA);
217  if(bA == Teuchos::null) {
218  RCP<const Map> rowMap = rcpA->getRowMap();
219  diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,true);
220  ArrayRCP<Scalar> diagVals = diag->getDataNonConst(0);
223  for (size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
224  rcpA->getLocalRowView(i, cols, vals);
225  diagVals[i] = Teuchos::ScalarTraits<Scalar>::zero();
226  for (LocalOrdinal j = 0; j < cols.size(); ++j) {
227  diagVals[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
228  }
229  }
230 
231  } else {
232  //TEUCHOS_TEST_FOR_EXCEPTION(bA->Rows() != bA->Cols(), Xpetra::Exceptions::RuntimeError,
233  // "UtilitiesBase::GetLumpedMatrixDiagonal(): you cannot extract the diagonal of a "<< bA->Rows() << "x"<< bA->Cols() << " blocked matrix." );
234 
235  diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(bA->getRangeMapExtractor()->getFullMap(),true);
236 
237  for (size_t row = 0; row < bA->Rows(); ++row) {
238  for (size_t col = 0; col < bA->Cols(); ++col) {
239  if (!bA->getMatrix(row,col).is_null()) {
240  // if we are in Thyra mode, but the block (row,row) is again a blocked operator, we have to use (pseudo) Xpetra-style GIDs with offset!
241  bool bThyraMode = bA->getRangeMapExtractor()->getThyraMode() && (Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA->getMatrix(row,col)) == Teuchos::null);
242  RCP<Vector> ddtemp = bA->getRangeMapExtractor()->ExtractVector(diag,row,bThyraMode);
243  RCP<const Vector> dd = GetLumpedMatrixDiagonal(bA->getMatrix(row,col));
244  ddtemp->update(Teuchos::as<Scalar>(1.0),*dd,Teuchos::as<Scalar>(1.0));
245  bA->getRangeMapExtractor()->InsertVector(ddtemp,row,diag,bThyraMode);
246  }
247  }
248  }
249 
250  }
251 
252  // we should never get here...
253  return diag;
254  }
255 
264 
265  RCP<Vector> ret = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(v->getMap(),true);
266 
267  // check whether input vector "v" is a BlockedVector
268  RCP<const BlockedVector> bv = Teuchos::rcp_dynamic_cast<const BlockedVector>(v);
269  if(bv.is_null() == false) {
270  RCP<BlockedVector> bret = Teuchos::rcp_dynamic_cast<BlockedVector>(ret);
271  TEUCHOS_TEST_FOR_EXCEPTION(bret.is_null() == true, MueLu::Exceptions::RuntimeError,"MueLu::UtilitiesBase::GetInverse: return vector should be of type BlockedVector");
272  RCP<const BlockedMap> bmap = bv->getBlockedMap();
273  for(size_t r = 0; r < bmap->getNumMaps(); ++r) {
274  RCP<const MultiVector> submvec = bv->getMultiVector(r,bmap->getThyraMode());
275  RCP<const Vector> subvec = submvec->getVector(0);
277  bret->setMultiVector(r, subvecinf, bmap->getThyraMode());
278  }
279  return ret;
280  }
281 
282  // v is an {Epetra,Tpetra}Vector: work with the underlying raw data
283  ArrayRCP<Scalar> retVals = ret->getDataNonConst(0);
284  ArrayRCP<const Scalar> inputVals = v->getData(0);
285  for (size_t i = 0; i < v->getMap()->getNodeNumElements(); ++i) {
286  if(Teuchos::ScalarTraits<Scalar>::magnitude(inputVals[i]) > tol)
287  retVals[i] = Teuchos::ScalarTraits<Scalar>::one() / inputVals[i];
288  else
289  retVals[i] = tolReplacement;
290  }
291  return ret;
292  }
293 
302  RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
303 
304  // Undo block map (if we have one)
305  RCP<const BlockedMap> browMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
306  if(!browMap.is_null()) rowMap = browMap->getMap();
307 
308  RCP<Vector> localDiag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap);
309  try {
310  const CrsMatrixWrap* crsOp = dynamic_cast<const CrsMatrixWrap*>(&A);
311  if (crsOp == NULL) {
312  throw Exceptions::RuntimeError("cast to CrsMatrixWrap failed");
313  }
315  crsOp->getLocalDiagOffsets(offsets);
316  crsOp->getLocalDiagCopy(*localDiag,offsets());
317  }
318  catch (...) {
319  ArrayRCP<Scalar> localDiagVals = localDiag->getDataNonConst(0);
321  for (LocalOrdinal i = 0; i < localDiagVals.size(); i++)
322  localDiagVals[i] = diagVals[i];
323  localDiagVals = diagVals = null;
324  }
325 
326  RCP<Vector> diagonal = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap);
328  importer = A.getCrsGraph()->getImporter();
329  if (importer == Teuchos::null) {
330  importer = Xpetra::ImportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap, colMap);
331  }
332  diagonal->doImport(*localDiag, *(importer), Xpetra::INSERT);
333  return diagonal;
334  }
335 
336 
345  using LO = LocalOrdinal;
346  using GO = GlobalOrdinal;
347  using SC = Scalar;
348  using STS = typename Teuchos::ScalarTraits<SC>;
349 
350  // Undo block map (if we have one)
351  RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
352  RCP<const BlockedMap> browMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
353  if(!browMap.is_null()) rowMap = browMap->getMap();
354 
355  RCP<Vector> local = Xpetra::VectorFactory<SC,LO,GO,Node>::Build(rowMap);
356  RCP<Vector> ghosted = Xpetra::VectorFactory<SC,LO,GO,Node>::Build(colMap,true);
357  ArrayRCP<SC> localVals = local->getDataNonConst(0);
358 
359  for (LO row = 0; row < static_cast<LO>(A.getRowMap()->getNodeNumElements()); ++row) {
360  size_t nnz = A.getNumEntriesInLocalRow(row);
361  ArrayView<const LO> indices;
362  ArrayView<const SC> vals;
363  A.getLocalRowView(row, indices, vals);
364 
365  SC si = STS::zero();
366 
367  for (LO colID = 0; colID < static_cast<LO>(nnz); colID++) {
368  if(indices[colID] != row) {
369  si += vals[colID];
370  }
371  }
372  localVals[row] = si;
373  }
374 
376  importer = A.getCrsGraph()->getImporter();
377  if (importer == Teuchos::null) {
378  importer = Xpetra::ImportFactory<LO,GO,Node>::Build(rowMap, colMap);
379  }
380  ghosted->doImport(*local, *(importer), Xpetra::INSERT);
381  return ghosted;
382  }
383 
384 
385 
388  RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
389  using STS = typename Teuchos::ScalarTraits<Scalar>;
390  using MTS = typename Teuchos::ScalarTraits<Magnitude>;
391  using MT = Magnitude;
392  using LO = LocalOrdinal;
393  using GO = GlobalOrdinal;
394  using SC = Scalar;
395  using RealValuedVector = Xpetra::Vector<MT,LO,GO,Node>;
396 
397  // Undo block map (if we have one)
398  RCP<const BlockedMap> browMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
399  if(!browMap.is_null()) rowMap = browMap->getMap();
400 
401  RCP<RealValuedVector> local = Xpetra::VectorFactory<MT,LO,GO,Node>::Build(rowMap);
402  RCP<RealValuedVector> ghosted = Xpetra::VectorFactory<MT,LO,GO,Node>::Build(colMap,true);
403  ArrayRCP<MT> localVals = local->getDataNonConst(0);
404 
405  for (LO rowIdx = 0; rowIdx < static_cast<LO>(A.getRowMap()->getNodeNumElements()); ++rowIdx) {
406  size_t nnz = A.getNumEntriesInLocalRow(rowIdx);
407  ArrayView<const LO> indices;
408  ArrayView<const SC> vals;
409  A.getLocalRowView(rowIdx, indices, vals);
410 
411  MT si = MTS::zero();
412 
413  for (LO colID = 0; colID < static_cast<LO>(nnz); ++colID) {
414  if(indices[colID] != rowIdx) {
415  si += STS::magnitude(vals[colID]);
416  }
417  }
418  localVals[rowIdx] = si;
419  }
420 
422  importer = A.getCrsGraph()->getImporter();
423  if (importer == Teuchos::null) {
424  importer = Xpetra::ImportFactory<LO,GO,Node>::Build(rowMap, colMap);
425  }
426  ghosted->doImport(*local, *(importer), Xpetra::INSERT);
427  return ghosted;
428  }
429 
430 
431 
432  // TODO: should NOT return an Array. Definition must be changed to:
433  // - ArrayRCP<> ResidualNorm(Matrix const &Op, MultiVector const &X, MultiVector const &RHS)
434  // or
435  // - void ResidualNorm(Matrix const &Op, MultiVector const &X, MultiVector const &RHS, Array &)
436  static Teuchos::Array<Magnitude> ResidualNorm(const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const MultiVector& X, const MultiVector& RHS) {
437  TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
438  const size_t numVecs = X.getNumVectors();
439  RCP<MultiVector> RES = Residual(Op, X, RHS);
440  Teuchos::Array<Magnitude> norms(numVecs);
441  RES->norm2(norms);
442  return norms;
443  }
444 
445  static Teuchos::Array<Magnitude> ResidualNorm(const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const MultiVector& X, const MultiVector& RHS, MultiVector & Resid) {
446  TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
447  const size_t numVecs = X.getNumVectors();
448  Residual(Op,X,RHS,Resid);
449  Teuchos::Array<Magnitude> norms(numVecs);
450  Resid.norm2(norms);
451  return norms;
452  }
453 
454  static RCP<MultiVector> Residual(const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const MultiVector& X, const MultiVector& RHS) {
455  TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
456  const size_t numVecs = X.getNumVectors();
457  // TODO Op.getRangeMap should return a BlockedMap if it is a BlockedCrsOperator
458  RCP<MultiVector> RES = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(RHS.getMap(), numVecs, false); // no need to initialize to zero
459  Op.residual(X,RHS,*RES);
460  return RES;
461  }
462 
463 
464  static void Residual(const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const MultiVector& X, const MultiVector& RHS, MultiVector & Resid) {
465  TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides");
466  TEUCHOS_TEST_FOR_EXCEPTION(Resid.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of residual vectors != number of right-hand sides");
467  Op.residual(X,RHS,Resid);
468  }
469 
470 
471 #ifndef _WIN32
472  static void PauseForDebugger() {
474  int myPID = comm->getRank();
475  int pid = getpid();
476  char hostname[80];
477  for (int i = 0; i <comm->getSize(); i++) {
478  if (i == myPID) {
479  gethostname(hostname, sizeof(hostname));
480  std::cout << "Host: " << hostname << "\tMPI rank: " << myPID << ",\tPID: " << pid << "\n\tattach " << pid << std::endl;
481  sleep(1);
482  }
483  }
484  if (myPID == 0) {
485  std::cout << "** Enter a character to continue > " << std::endl;
486  char go = ' ';
487  int r = scanf("%c", &go);
488  (void)r;
489  assert(r > 0);
490  }
491  comm->barrier();
492  }
493 #else
494  static void PauseForDebugger() {
495  throw(Exceptions::RuntimeError("MueLu Utils: PauseForDebugger not implemented on Windows."));
496  }
497 #endif
498 
514  static Scalar PowerMethod(const Matrix& A, bool scaleByDiag = true,
515  LocalOrdinal niters = 10, Magnitude tolerance = 1e-2, bool verbose = false, unsigned int seed = 123) {
516  TEUCHOS_TEST_FOR_EXCEPTION(!(A.getRangeMap()->isSameAs(*(A.getDomainMap()))), Exceptions::Incompatible,
517  "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
518 
519  // Create three vectors, fill z with random numbers
520  RCP<Vector> q = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getDomainMap());
521  RCP<Vector> r = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
522  RCP<Vector> z = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
523 
524  z->setSeed(seed); // seed random number generator
525  z->randomize(true);// use Xpetra implementation: -> same results for Epetra and Tpetra
526 
527  Teuchos::Array<Magnitude> norms(1);
528 
529  typedef Teuchos::ScalarTraits<Scalar> STS;
530 
531  const Scalar zero = STS::zero(), one = STS::one();
532 
533  Scalar lambda = zero;
534  Magnitude residual = STS::magnitude(zero);
535 
536  // power iteration
537  RCP<Vector> diagInvVec;
538  if (scaleByDiag) {
539  RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
540  A.getLocalDiagCopy(*diagVec);
541  diagInvVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
542  diagInvVec->reciprocal(*diagVec);
543  }
544 
545  for (int iter = 0; iter < niters; ++iter) {
546  z->norm2(norms); // Compute 2-norm of z
547  q->update(one/norms[0], *z, zero); // Set q = z / normz
548  A.apply(*q, *z); // Compute z = A*q
549  if (scaleByDiag)
550  z->elementWiseMultiply(one, *diagInvVec, *z, zero);
551  lambda = q->dot(*z); // Approximate maximum eigenvalue: lamba = dot(q,z)
552 
553  if (iter % 100 == 0 || iter + 1 == niters) {
554  r->update(1.0, *z, -lambda, *q, zero); // Compute A*q - lambda*q
555  r->norm2(norms);
556  residual = STS::magnitude(norms[0] / lambda);
557  if (verbose) {
558  std::cout << "Iter = " << iter
559  << " Lambda = " << lambda
560  << " Residual of A*q - lambda*q = " << residual
561  << std::endl;
562  }
563  }
564  if (residual < tolerance)
565  break;
566  }
567  return lambda;
568  }
569 
570 
571 
572  static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) {
573  RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(os));
574  return fancy;
575  }
576 
582  const size_t numVectors = v.size();
583 
585  for (size_t j = 0; j < numVectors; j++) {
586  d += (v[j][i0] - v[j][i1])*(v[j][i0] - v[j][i1]);
587  }
589  }
590 
603  static Teuchos::ArrayRCP<const bool> DetectDirichletRows(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero(), bool count_twos_as_dirichlet=false) {
604  LocalOrdinal numRows = A.getNodeNumRows();
605  typedef Teuchos::ScalarTraits<Scalar> STS;
606  ArrayRCP<bool> boundaryNodes(numRows, true);
607  if (count_twos_as_dirichlet) {
608  for (LocalOrdinal row = 0; row < numRows; row++) {
611  A.getLocalRowView(row, indices, vals);
612  size_t nnz = A.getNumEntriesInLocalRow(row);
613  if (nnz > 2) {
614  size_t col;
615  for (col = 0; col < nnz; col++)
616  if ( (indices[col] != row) && STS::magnitude(vals[col]) > tol) {
617  if (!boundaryNodes[row])
618  break;
619  boundaryNodes[row] = false;
620  }
621  if (col == nnz)
622  boundaryNodes[row] = true;
623  }
624  }
625  } else {
626  for (LocalOrdinal row = 0; row < numRows; row++) {
629  A.getLocalRowView(row, indices, vals);
630  size_t nnz = A.getNumEntriesInLocalRow(row);
631  if (nnz > 1)
632  for (size_t col = 0; col < nnz; col++)
633  if ( (indices[col] != row) && STS::magnitude(vals[col]) > tol) {
634  boundaryNodes[row] = false;
635  break;
636  }
637  }
638  }
639  return boundaryNodes;
640  }
641 
642 
655  static Teuchos::ArrayRCP<const bool> DetectDirichletRowsExt(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, bool & bHasZeroDiagonal, const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero()) {
656 
657  // assume that there is no zero diagonal in matrix
658  bHasZeroDiagonal = false;
659 
660  Teuchos::RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
661  A.getLocalDiagCopy(*diagVec);
662  Teuchos::ArrayRCP< const Scalar > diagVecData = diagVec->getData(0);
663 
664  LocalOrdinal numRows = A.getNodeNumRows();
665  typedef Teuchos::ScalarTraits<Scalar> STS;
666  ArrayRCP<bool> boundaryNodes(numRows, false);
667  for (LocalOrdinal row = 0; row < numRows; row++) {
670  A.getLocalRowView(row, indices, vals);
671  size_t nnz = 0; // collect nonzeros in row (excluding the diagonal)
672  bool bHasDiag = false;
673  for (decltype(indices.size()) col = 0; col < indices.size(); col++) {
674  if ( indices[col] != row) {
675  if (STS::magnitude(vals[col] / STS::magnitude(sqrt(STS::magnitude(diagVecData[row]) * STS::magnitude(diagVecData[col]))) ) > tol) {
676  nnz++;
677  }
678  } else bHasDiag = true; // found a diagonal entry
679  }
680  if (bHasDiag == false) bHasZeroDiagonal = true; // we found at least one row without a diagonal
681  else if(nnz == 0) boundaryNodes[row] = true;
682  }
683  return boundaryNodes;
684  }
685 
686 
697  static Teuchos::ArrayRCP<const bool> DetectDirichletCols(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
698  const Teuchos::ArrayRCP<const bool>& dirichletRows) {
703  Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > myColsToZero = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap,1);
704  myColsToZero->putScalar(zero);
705  // Find all local column indices that are in Dirichlet rows, record in myColsToZero as 1.0
706  for(size_t i=0; i<(size_t) dirichletRows.size(); i++) {
707  if (dirichletRows[i]) {
710  A.getLocalRowView(i,indices,values);
711  for(size_t j=0; j<static_cast<size_t>(indices.size()); j++)
712  myColsToZero->replaceLocalValue(indices[j],0,one);
713  }
714  }
715 
716  Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > globalColsToZero = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(domMap,1);
717  globalColsToZero->putScalar(zero);
718  Teuchos::RCP<Xpetra::Export<LocalOrdinal,GlobalOrdinal,Node> > exporter = Xpetra::ExportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(colMap,domMap);
719  // export to domain map
720  globalColsToZero->doExport(*myColsToZero,*exporter,Xpetra::ADD);
721  // import to column map
722  myColsToZero->doImport(*globalColsToZero,*exporter,Xpetra::INSERT);
723  Teuchos::ArrayRCP<const Scalar> myCols = myColsToZero->getData(0);
724  Teuchos::ArrayRCP<bool> dirichletCols(colMap->getNodeNumElements(), true);
726  for(size_t i=0; i<colMap->getNodeNumElements(); i++) {
727  dirichletCols[i] = Teuchos::ScalarTraits<Scalar>::magnitude(myCols[i])>2.0*eps;
728  }
729  return dirichletCols;
730  }
731 
732 
737  static Scalar Frobenius(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B) {
738  // We check only row maps. Column may be different. One would hope that they are the same, as we typically
739  // calculate frobenius norm of the specified sparsity pattern with an updated matrix from the previous step,
740  // but matrix addition, even when one is submatrix of the other, changes column map (though change may be as
741  // simple as couple of elements swapped)
742  TEUCHOS_TEST_FOR_EXCEPTION(!A.getRowMap()->isSameAs(*B.getRowMap()), Exceptions::Incompatible, "MueLu::CGSolver::Frobenius: row maps are incompatible");
743  TEUCHOS_TEST_FOR_EXCEPTION(!A.isFillComplete() || !B.isFillComplete(), Exceptions::RuntimeError, "Matrices must be fill completed");
744 
745  const Map& AColMap = *A.getColMap();
746  const Map& BColMap = *B.getColMap();
747 
750  size_t nnzA = 0, nnzB = 0;
751 
752  // We use a simple algorithm
753  // for each row we fill valBAll array with the values in the corresponding row of B
754  // as such, it serves as both sorted array and as storage, so we don't need to do a
755  // tricky problem: "find a value in the row of B corresponding to the specific GID"
756  // Once we do that, we translate LID of entries of row of A to LID of B, and multiply
757  // corresponding entries.
758  // The algorithm should be reasonably cheap, as it does not sort anything, provided
759  // that getLocalElement and getGlobalElement functions are reasonably effective. It
760  // *is* possible that the costs are hidden in those functions, but if maps are close
761  // to linear maps, we should be fine
762  Teuchos::Array<Scalar> valBAll(BColMap.getNodeNumElements());
763 
765  Scalar zero = Teuchos::ScalarTraits<Scalar> ::zero(), f = zero, gf;
766  size_t numRows = A.getNodeNumRows();
767  for (size_t i = 0; i < numRows; i++) {
768  A.getLocalRowView(i, indA, valA);
769  B.getLocalRowView(i, indB, valB);
770  nnzA = indA.size();
771  nnzB = indB.size();
772 
773  // Set up array values
774  for (size_t j = 0; j < nnzB; j++)
775  valBAll[indB[j]] = valB[j];
776 
777  for (size_t j = 0; j < nnzA; j++) {
778  // The cost of the whole Frobenius dot product function depends on the
779  // cost of the getLocalElement and getGlobalElement functions here.
780  LocalOrdinal ind = BColMap.getLocalElement(AColMap.getGlobalElement(indA[j]));
781  if (ind != invalid)
782  f += valBAll[ind] * valA[j];
783  }
784 
785  // Clean up array values
786  for (size_t j = 0; j < nnzB; j++)
787  valBAll[indB[j]] = zero;
788  }
789 
790  MueLu_sumAll(AColMap.getComm(), f, gf);
791 
792  return gf;
793  }
794 
804  static void SetRandomSeed(const Teuchos::Comm<int> &comm) {
805  // Distribute the seeds evenly in [1,maxint-1]. This guarantees nothing
806  // about where in random number stream we are, but avoids overflow situations
807  // in parallel when multiplying by a PID. It would be better to use
808  // a good parallel random number generator.
809  double one = 1.0;
810  int maxint = INT_MAX; //= 2^31-1 = 2147483647 for 32-bit integers
811  int mySeed = Teuchos::as<int>((maxint-1) * (one -(comm.getRank()+1)/(comm.getSize()+one)) );
812  if (mySeed < 1 || mySeed == maxint) {
813  std::ostringstream errStr;
814  errStr << "Error detected with random seed = " << mySeed << ". It should be in the interval [1,2^31-2].";
815  throw Exceptions::RuntimeError(errStr.str());
816  }
817  std::srand(mySeed);
818  // For Tpetra, we could use Kokkos' random number generator here.
820  // Epetra
821  // MultiVector::Random() -> Epetra_Util::RandomDouble() -> Epetra_Utils::RandomInt()
822  // Its own random number generator, based on Seed_. Seed_ is initialized in Epetra_Util constructor with std::rand()
823  // So our setting std::srand() affects that too
824  }
825 
826 
827 
828  // Finds the OAZ Dirichlet rows for this matrix
829  // so far only used in IntrepidPCoarsenFactory
830  // TODO check whether we can use DetectDirichletRows instead
831  static void FindDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &A,
832  std::vector<LocalOrdinal>& dirichletRows, bool count_twos_as_dirichlet=false) {
834  dirichletRows.resize(0);
835  for(size_t i=0; i<A->getNodeNumRows(); i++) {
838  A->getLocalRowView(i,indices,values);
839  int nnz=0;
840  for (size_t j=0; j<(size_t)indices.size(); j++) {
842  nnz++;
843  }
844  }
845  if (nnz == 1 || (count_twos_as_dirichlet && nnz == 2)) {
846  dirichletRows.push_back(i);
847  }
848  }
849  }
850 
851  // Applies Ones-and-Zeros to matrix rows
852  // Takes a vector of row indices
853  static void ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
854  const std::vector<LocalOrdinal>& dirichletRows) {
855  RCP<const Map> Rmap = A->getRowMap();
856  RCP<const Map> Cmap = A->getColMap();
859 
860  for(size_t i=0; i<dirichletRows.size(); i++) {
861  GlobalOrdinal row_gid = Rmap->getGlobalElement(dirichletRows[i]);
862 
865  A->getLocalRowView(dirichletRows[i],indices,values);
866  // NOTE: This won't work with fancy node types.
867  Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
868  for(size_t j=0; j<(size_t)indices.size(); j++) {
869  if(Cmap->getGlobalElement(indices[j])==row_gid)
870  valuesNC[j]=one;
871  else
872  valuesNC[j]=zero;
873  }
874  }
875  }
876 
877  // Applies Ones-and-Zeros to matrix rows
878  // Takes a Boolean array.
879  static void ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
880  const Teuchos::ArrayRCP<const bool>& dirichletRows) {
881  TEUCHOS_ASSERT(A->isFillComplete());
882  RCP<const Map> domMap = A->getDomainMap();
883  RCP<const Map> ranMap = A->getRangeMap();
884  RCP<const Map> Rmap = A->getRowMap();
885  RCP<const Map> Cmap = A->getColMap();
886  TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == Rmap->getNodeNumElements());
889  A->resumeFill();
890  for(size_t i=0; i<(size_t) dirichletRows.size(); i++) {
891  if (dirichletRows[i]){
892  GlobalOrdinal row_gid = Rmap->getGlobalElement(i);
893 
896  A->getLocalRowView(i,indices,values);
897 
898  Teuchos::ArrayRCP<Scalar> valuesNC(values.size());
899  for(size_t j=0; j<(size_t)indices.size(); j++) {
900  if(Cmap->getGlobalElement(indices[j])==row_gid)
901  valuesNC[j]=one;
902  else
903  valuesNC[j]=zero;
904  }
905  A->replaceLocalValues(i,indices,valuesNC());
906  }
907  }
908  A->fillComplete(domMap, ranMap);
909  }
910 
911  // Zeros out rows
912  // Takes a vector containg Dirichlet row indices
913  static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
914  const std::vector<LocalOrdinal>& dirichletRows,
916  for(size_t i=0; i<dirichletRows.size(); i++) {
919  A->getLocalRowView(dirichletRows[i],indices,values);
920  // NOTE: This won't work with fancy node types.
921  Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
922  for(size_t j=0; j<(size_t)indices.size(); j++)
923  valuesNC[j]=replaceWith;
924  }
925  }
926 
927  // Zeros out rows
928  // Takes a Boolean ArrayRCP
929  static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
930  const Teuchos::ArrayRCP<const bool>& dirichletRows,
932  TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == A->getRowMap()->getNodeNumElements());
933  for(size_t i=0; i<(size_t) dirichletRows.size(); i++) {
934  if (dirichletRows[i]) {
937  A->getLocalRowView(i,indices,values);
938  // NOTE: This won't work with fancy node types.
939  Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
940  for(size_t j=0; j<(size_t)indices.size(); j++)
941  valuesNC[j]=replaceWith;
942  }
943  }
944  }
945 
946  // Zeros out rows
947  // Takes a Boolean ArrayRCP
948  static void ZeroDirichletRows(Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& X,
949  const Teuchos::ArrayRCP<const bool>& dirichletRows,
951  TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == X->getMap()->getNodeNumElements());
952  for(size_t i=0; i<(size_t) dirichletRows.size(); i++) {
953  if (dirichletRows[i]) {
954  for(size_t j=0; j<X->getNumVectors(); j++)
955  X->replaceLocalValue(i,j,replaceWith);
956  }
957  }
958  }
959 
960  // Zeros out columns
961  // Takes a Boolean vector
963  const Teuchos::ArrayRCP<const bool>& dirichletCols,
965  TEUCHOS_ASSERT(static_cast<size_t>(dirichletCols.size()) == A->getColMap()->getNodeNumElements());
966  for(size_t i=0; i<A->getNodeNumRows(); i++) {
969  A->getLocalRowView(i,indices,values);
970  // NOTE: This won't work with fancy node types.
971  Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
972  for(size_t j=0; j<static_cast<size_t>(indices.size()); j++)
973  if (dirichletCols[indices[j]])
974  valuesNC[j] = replaceWith;
975  }
976  }
977 
978  // Finds the OAZ Dirichlet rows for this matrix
979  static void FindDirichletRowsAndPropagateToCols(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &A,
980  Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletRow,
981  Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletCol) {
982 
983  // Make sure A's RowMap == DomainMap
984  if(!A->getRowMap()->isSameAs(*A->getDomainMap())) {
985  throw std::runtime_error("UtilitiesBase::FindDirichletRowsAndPropagateToCols row and domain maps must match.");
986  }
987  RCP<const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer = A->getCrsGraph()->getImporter();
988  bool has_import = !importer.is_null();
989 
990  // Find the Dirichlet rows
991  std::vector<LocalOrdinal> dirichletRows;
992  FindDirichletRows(A,dirichletRows);
993 
994 #if 0
995  printf("[%d] DirichletRow Ids = ",A->getRowMap()->getComm()->getRank());
996  for(size_t i=0; i<(size_t) dirichletRows.size(); i++)
997  printf("%d ",dirichletRows[i]);
998  printf("\n");
999  fflush(stdout);
1000 #endif
1001  // Allocate all as non-Dirichlet
1002  isDirichletRow = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getRowMap(),true);
1003  isDirichletCol = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getColMap(),true);
1004 
1005  Teuchos::ArrayRCP<int> dr_rcp = isDirichletRow->getDataNonConst(0);
1006  Teuchos::ArrayView<int> dr = dr_rcp();
1007  Teuchos::ArrayRCP<int> dc_rcp = isDirichletCol->getDataNonConst(0);
1008  Teuchos::ArrayView<int> dc = dc_rcp();
1009  for(size_t i=0; i<(size_t) dirichletRows.size(); i++) {
1010  dr[dirichletRows[i]] = 1;
1011  if(!has_import) dc[dirichletRows[i]] = 1;
1012  }
1013 
1014  if(has_import)
1015  isDirichletCol->doImport(*isDirichletRow,*importer,Xpetra::CombineMode::ADD);
1016 
1017  }
1018 
1019  // This routine takes a BlockedMap and an Importer (assuming that the BlockedMap matches the source of the importer) and generates a BlockedMap corresponding
1020  // to the Importer's target map. We assume that the targetMap is unique (which, is not a strict requirement of an Importer, but is here and no, we don't check)
1021  // This is largely intended to be used in repartitioning of blocked matrices
1022  static RCP<const Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node> > GeneratedBlockedTargetMap(const Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node> & sourceBlockedMap,
1023  const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> & Importer) {
1024  typedef Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> IntVector;
1025  Xpetra::UnderlyingLib lib = sourceBlockedMap.lib();
1026 
1027  // De-stride the map if we have to (might regret this later)
1028  RCP<const Map> fullMap = sourceBlockedMap.getMap();
1029  RCP<const Map> stridedMap = Teuchos::rcp_dynamic_cast<const Xpetra::StridedMap<LocalOrdinal,GlobalOrdinal,Node> >(fullMap);
1030  if(!stridedMap.is_null()) fullMap = stridedMap->getMap();
1031 
1032  // Initial sanity checking for map compatibil
1033  const size_t numSubMaps = sourceBlockedMap.getNumMaps();
1034  if(!Importer.getSourceMap()->isCompatible(*fullMap))
1035  throw std::runtime_error("GenerateBlockedTargetMap(): Map compatibility error");
1036 
1037  // Build an indicator vector
1038  RCP<IntVector> block_ids = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(fullMap);
1039 
1040  for(size_t i=0; i<numSubMaps; i++) {
1041  RCP<const Map> map = sourceBlockedMap.getMap(i);
1042 
1043  for(size_t j=0; j<map->getNodeNumElements(); j++) {
1044  LocalOrdinal jj = fullMap->getLocalElement(map->getGlobalElement(j));
1045  block_ids->replaceLocalValue(jj,(int)i);
1046  }
1047  }
1048 
1049  // Get the block ids for the new map
1050  RCP<const Map> targetMap = Importer.getTargetMap();
1051  RCP<IntVector> new_block_ids = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(targetMap);
1052  new_block_ids->doImport(*block_ids,Importer,Xpetra::CombineMode::ADD);
1053  Teuchos::ArrayRCP<const int> dataRCP = new_block_ids->getData(0);
1054  Teuchos::ArrayView<const int> data = dataRCP();
1055 
1056 
1057  // Get the GIDs for each subblock
1058  Teuchos::Array<Teuchos::Array<GlobalOrdinal> > elementsInSubMap(numSubMaps);
1059  for(size_t i=0; i<targetMap->getNodeNumElements(); i++) {
1060  elementsInSubMap[data[i]].push_back(targetMap->getGlobalElement(i));
1061  }
1062 
1063  // Generate the new submaps
1064  std::vector<RCP<const Map> > subMaps(numSubMaps);
1065  for(size_t i=0; i<numSubMaps; i++) {
1066  subMaps[i] = Xpetra::MapFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(lib,Teuchos::OrdinalTraits<GlobalOrdinal>::invalid(),elementsInSubMap[i](),targetMap->getIndexBase(),targetMap->getComm());
1067  }
1068 
1069  // Build the BlockedMap
1070  return rcp(new BlockedMap(targetMap,subMaps));
1071 
1072  }
1073 
1074  }; // class Utils
1075 
1076 
1078 
1079 } //namespace MueLu
1080 
1081 #define MUELU_UTILITIESBASE_SHORT
1082 #endif // MUELU_UTILITIESBASE_DECL_HPP
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedMultiVector
virtual int getSize() const =0
MueLu::DefaultLocalOrdinal LocalOrdinal
static Teuchos::ArrayRCP< const bool > DetectDirichletCols(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows)
Detect Dirichlet columns based on Dirichlet rows.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
virtual int getRank() const =0
static magnitudeType eps()
static Scalar PowerMethod(const Matrix &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
Simple transpose for Tpetra::CrsMatrix types.
static void Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS, MultiVector &Resid)
static Teuchos::RCP< Vector > GetInverse(Teuchos::RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero())
Return vector containing inverse of input vector.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > GeneratedBlockedTargetMap(const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > &sourceBlockedMap, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &Importer)
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows)
size_type size() const
static void FindDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, std::vector< LocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
static RCP< Vector > GetMatrixOverlappedDeletedRowsum(const Matrix &A)
Extract Overlapped Matrix Deleted Rowsum.
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
size_type size() const
Exception throws to report incompatible objects (like maps).
static Teuchos::ArrayRCP< const bool > DetectDirichletRowsExt(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool &bHasZeroDiagonal, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero())
Detect Dirichlet rows (extended version)
MueLu::DefaultNode Node
#define MueLu_sumAll(rcpComm, in, out)
static Teuchos::ScalarTraits< Scalar >::magnitudeType Distance2(const Teuchos::Array< Teuchos::ArrayRCP< const Scalar >> &v, LocalOrdinal i0, LocalOrdinal i1)
Squared distance between two rows in a multivector.
static RCP< Time > getNewTimer(const std::string &name)
static RCP< Xpetra::Vector< Magnitude, LocalOrdinal, GlobalOrdinal, Node > > GetMatrixOverlappedAbsDeletedRowsum(const Matrix &A)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
static void FindDirichletRowsAndPropagateToCols(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node > > &isDirichletRow, Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node > > &isDirichletCol)
Xpetra::BlockedVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedVector
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
static Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal of lumped matrix.
static Scalar Frobenius(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B)
Frobenius inner product of two matrices.
T * getRawPtr() const
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static void seedrandom(unsigned int s)
Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Vector
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows)
static magnitudeType magnitude(T a)
TransListIter iter
Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > BlockedMap
static void ZeroDirichletCols(Teuchos::RCP< Matrix > &A, const Teuchos::ArrayRCP< const bool > &dirichletCols, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
void push_back(const value_type &x)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero(), bool count_twos_as_dirichlet=false)
Detect Dirichlet rows.
static RCP< Matrix > Crs2Op(RCP< CrsMatrix > Op)
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixWrap
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
Exception throws to report errors in the internal logical of the program.
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
#define TEUCHOS_ASSERT(assertion_test)
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS, MultiVector &Resid)
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Teuchos::RCP< const Matrix > rcpA)
Extract Matrix Diagonal of lumped matrix.
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100)
Extract Matrix Diagonal.
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &X, const Teuchos::ArrayRCP< const bool > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.
bool is_null() const