MueLu  Version of the Day
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 #include <unistd.h> //necessary for "sleep" function in debugging methods
50 #include <string>
51 
52 #include "MueLu_ConfigDefs.hpp"
53 
54 #include <Teuchos_DefaultComm.hpp>
55 #include <Teuchos_ScalarTraits.hpp>
56 #include <Teuchos_ParameterList.hpp>
57 
58 #include <Xpetra_BlockedCrsMatrix_fwd.hpp>
59 #include <Xpetra_CrsMatrix_fwd.hpp>
60 #include <Xpetra_CrsMatrixWrap_fwd.hpp>
61 #include <Xpetra_Map_fwd.hpp>
62 #include <Xpetra_BlockedMap_fwd.hpp>
63 #include <Xpetra_MapFactory_fwd.hpp>
64 #include <Xpetra_Matrix_fwd.hpp>
65 #include <Xpetra_MatrixFactory_fwd.hpp>
66 #include <Xpetra_MultiVector_fwd.hpp>
67 #include <Xpetra_MultiVectorFactory_fwd.hpp>
68 #include <Xpetra_Operator_fwd.hpp>
69 #include <Xpetra_Vector_fwd.hpp>
70 #include <Xpetra_BlockedMultiVector.hpp>
71 #include <Xpetra_BlockedVector.hpp>
72 #include <Xpetra_VectorFactory_fwd.hpp>
73 #include <Xpetra_ExportFactory.hpp>
74 
75 #include <Xpetra_Import.hpp>
76 #include <Xpetra_ImportFactory.hpp>
77 #include <Xpetra_MatrixMatrix.hpp>
78 #include <Xpetra_CrsMatrixWrap.hpp>
79 
80 #include "MueLu_Exceptions.hpp"
81 
82 
83 
84 
85 namespace MueLu {
86 
87 // MPI helpers
88 #define MueLu_sumAll(rcpComm, in, out) \
89  Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_SUM, in, Teuchos::outArg(out))
90 #define MueLu_minAll(rcpComm, in, out) \
91  Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MIN, in, Teuchos::outArg(out))
92 #define MueLu_maxAll(rcpComm, in, out) \
93  Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MAX, in, Teuchos::outArg(out))
94 
102  template <class Scalar,
103  class LocalOrdinal = int,
104  class GlobalOrdinal = LocalOrdinal,
105  class Node = KokkosClassic::DefaultNode::DefaultNodeType>
106  class UtilitiesBase {
107  public:
108 #undef MUELU_UTILITIESBASE_SHORT
109 //#include "MueLu_UseShortNames.hpp"
110  private:
111  typedef Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> CrsMatrixWrap;
112  typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> CrsMatrix;
113  typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> Matrix;
114  typedef Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> Vector;
115  typedef Xpetra::BlockedVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> BlockedVector;
116  typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MultiVector;
117  typedef Xpetra::BlockedMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> BlockedMultiVector;
118  typedef Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node> BlockedMap;
119  typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> Map;
120  public:
121  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
122 
123 
124  static RCP<Matrix> Crs2Op(RCP<CrsMatrix> Op) {
125  if (Op.is_null())
126  return Teuchos::null;
127  return rcp(new CrsMatrixWrap(Op));
128  }
129 
136  static Teuchos::ArrayRCP<Scalar> GetMatrixDiagonal(const Matrix& A) {
137  size_t numRows = A.getRowMap()->getNodeNumElements();
138  Teuchos::ArrayRCP<Scalar> diag(numRows);
139  Teuchos::ArrayView<const LocalOrdinal> cols;
140  Teuchos::ArrayView<const Scalar> vals;
141  for (size_t i = 0; i < numRows; ++i) {
142  A.getLocalRowView(i, cols, vals);
143  LocalOrdinal j = 0;
144  for (; j < cols.size(); ++j) {
145  if (Teuchos::as<size_t>(cols[j]) == i) {
146  diag[i] = vals[j];
147  break;
148  }
149  }
150  if (j == cols.size()) {
151  // Diagonal entry is absent
152  diag[i] = Teuchos::ScalarTraits<Scalar>::zero();
153  }
154  }
155  return diag;
156  }
157 
164  static RCP<Vector> GetMatrixDiagonalInverse(const Matrix& A, Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100) {
165 
166  RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
167 
168  RCP<const Map> rowMap = rcpA->getRowMap();
169  RCP<Vector> diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,true);
170 
171  rcpA->getLocalDiagCopy(*diag);
172 
173  RCP<Vector> inv = MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::GetInverse(diag, tol, Teuchos::ScalarTraits<Scalar>::zero());
174 
175  return inv;
176  }
177 
178 
179 
186  static Teuchos::ArrayRCP<Scalar> GetLumpedMatrixDiagonal(const Matrix& A) {
187 
188  size_t numRows = A.getRowMap()->getNodeNumElements();
189  Teuchos::ArrayRCP<Scalar> diag(numRows);
190  Teuchos::ArrayView<const LocalOrdinal> cols;
191  Teuchos::ArrayView<const Scalar> vals;
192  for (size_t i = 0; i < numRows; ++i) {
193  A.getLocalRowView(i, cols, vals);
194  diag[i] = Teuchos::ScalarTraits<Scalar>::zero();
195  for (LocalOrdinal j = 0; j < cols.size(); ++j) {
196  diag[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
197  }
198  }
199  return diag;
200  }
201 
208  static Teuchos::RCP<Vector> GetLumpedMatrixDiagonal(Teuchos::RCP<const Matrix> rcpA) {
209 
210  RCP<Vector> diag = Teuchos::null;
211 
212  RCP<const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > bA =
213  Teuchos::rcp_dynamic_cast<const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(rcpA);
214  if(bA == Teuchos::null) {
215  RCP<const Map> rowMap = rcpA->getRowMap();
216  diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,true);
217  ArrayRCP<Scalar> diagVals = diag->getDataNonConst(0);
218  Teuchos::ArrayView<const LocalOrdinal> cols;
219  Teuchos::ArrayView<const Scalar> vals;
220  for (size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
221  rcpA->getLocalRowView(i, cols, vals);
222  diagVals[i] = Teuchos::ScalarTraits<Scalar>::zero();
223  for (LocalOrdinal j = 0; j < cols.size(); ++j) {
224  diagVals[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
225  }
226  }
227 
228  } else {
229  //TEUCHOS_TEST_FOR_EXCEPTION(bA->Rows() != bA->Cols(), Xpetra::Exceptions::RuntimeError,
230  // "UtilitiesBase::GetLumpedMatrixDiagonal(): you cannot extract the diagonal of a "<< bA->Rows() << "x"<< bA->Cols() << " blocked matrix." );
231 
232  diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(bA->getRangeMapExtractor()->getFullMap(),true);
233 
234  for (size_t row = 0; row < bA->Rows(); ++row) {
235  for (size_t col = 0; col < bA->Cols(); ++col) {
236  if (!bA->getMatrix(row,col).is_null()) {
237  // 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!
238  bool bThyraMode = bA->getRangeMapExtractor()->getThyraMode() && (Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA->getMatrix(row,col)) == Teuchos::null);
239  RCP<Vector> ddtemp = bA->getRangeMapExtractor()->ExtractVector(diag,row,bThyraMode);
240  RCP<const Vector> dd = GetLumpedMatrixDiagonal(bA->getMatrix(row,col));
241  ddtemp->update(Teuchos::as<Scalar>(1.0),*dd,Teuchos::as<Scalar>(1.0));
242  bA->getRangeMapExtractor()->InsertVector(ddtemp,row,diag,bThyraMode);
243  }
244  }
245  }
246 
247  }
248 
249  // we should never get here...
250  return diag;
251  }
252 
260  static Teuchos::RCP<Vector> GetInverse(Teuchos::RCP<const Vector> v, Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100, Scalar tolReplacement = Teuchos::ScalarTraits<Scalar>::zero()) {
261 
262  RCP<Vector> ret = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(v->getMap(),true);
263 
264  // check whether input vector "v" is a BlockedVector
265  RCP<const BlockedVector> bv = Teuchos::rcp_dynamic_cast<const BlockedVector>(v);
266  if(bv.is_null() == false) {
267  RCP<BlockedVector> bret = Teuchos::rcp_dynamic_cast<BlockedVector>(ret);
268  TEUCHOS_TEST_FOR_EXCEPTION(bret.is_null() == true, MueLu::Exceptions::RuntimeError,"MueLu::UtilitiesBase::GetInverse: return vector should be of type BlockedVector");
269  RCP<const BlockedMap> bmap = bv->getBlockedMap();
270  for(size_t r = 0; r < bmap->getNumMaps(); ++r) {
271  RCP<const MultiVector> submvec = bv->getMultiVector(r,bmap->getThyraMode());
272  RCP<const Vector> subvec = submvec->getVector(0);
273  RCP<Vector> subvecinf = MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::GetInverse(subvec,tol,tolReplacement);
274  bret->setMultiVector(r, subvecinf, bmap->getThyraMode());
275  }
276  return ret;
277  }
278 
279  // v is an {Epetra,Tpetra}Vector: work with the underlying raw data
280  ArrayRCP<Scalar> retVals = ret->getDataNonConst(0);
281  ArrayRCP<const Scalar> inputVals = v->getData(0);
282  for (size_t i = 0; i < v->getMap()->getNodeNumElements(); ++i) {
283  if(Teuchos::ScalarTraits<Scalar>::magnitude(inputVals[i]) > tol)
284  retVals[i] = Teuchos::ScalarTraits<Scalar>::one() / inputVals[i];
285  else
286  retVals[i] = tolReplacement;
287  }
288  return ret;
289  }
290 
298  static RCP<Vector> GetMatrixOverlappedDiagonal(const Matrix& A) {
299  RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
300  RCP<Vector> localDiag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap);
301 
302  try {
303  const CrsMatrixWrap* crsOp = dynamic_cast<const CrsMatrixWrap*>(&A);
304  if (crsOp == NULL) {
305  throw Exceptions::RuntimeError("cast to CrsMatrixWrap failed");
306  }
307  Teuchos::ArrayRCP<size_t> offsets;
308  crsOp->getLocalDiagOffsets(offsets);
309  crsOp->getLocalDiagCopy(*localDiag,offsets());
310  }
311  catch (...) {
312  ArrayRCP<Scalar> localDiagVals = localDiag->getDataNonConst(0);
313  Teuchos::ArrayRCP<Scalar> diagVals = GetMatrixDiagonal(A);
314  for (LocalOrdinal i = 0; i < localDiagVals.size(); i++)
315  localDiagVals[i] = diagVals[i];
316  localDiagVals = diagVals = null;
317  }
318 
319  RCP<Vector> diagonal = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap);
320  RCP< const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer;
321  importer = A.getCrsGraph()->getImporter();
322  if (importer == Teuchos::null) {
323  importer = Xpetra::ImportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap, colMap);
324  }
325  diagonal->doImport(*localDiag, *(importer), Xpetra::INSERT);
326  return diagonal;
327  }
328 
329 
330  // TODO: should NOT return an Array. Definition must be changed to:
331  // - ArrayRCP<> ResidualNorm(Matrix const &Op, MultiVector const &X, MultiVector const &RHS)
332  // or
333  // - void ResidualNorm(Matrix const &Op, MultiVector const &X, MultiVector const &RHS, Array &)
334  static Teuchos::Array<Magnitude> ResidualNorm(const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const MultiVector& X, const MultiVector& RHS) {
335  TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
336  const size_t numVecs = X.getNumVectors();
337  RCP<MultiVector> RES = Residual(Op, X, RHS);
338  Teuchos::Array<Magnitude> norms(numVecs);
339  RES->norm2(norms);
340  return norms;
341  }
342 
343  static RCP<MultiVector> Residual(const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const MultiVector& X, const MultiVector& RHS) {
344  TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
345  const size_t numVecs = X.getNumVectors();
346  Scalar one = Teuchos::ScalarTraits<Scalar>::one(), negone = -one, zero = Teuchos::ScalarTraits<Scalar>::zero();
347  // TODO Op.getRangeMap should return a BlockedMap if it is a BlockedCrsOperator
348  RCP<MultiVector> RES = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(RHS.getMap(), numVecs, false); // no need to initialize to zero
349  Op.apply(X, *RES, Teuchos::NO_TRANS, one, zero);
350  RES->update(one, RHS, negone);
351  return RES;
352  }
353 
354 #ifndef _WIN32
355 #include <unistd.h>
356  static void PauseForDebugger() {
357  RCP<const Teuchos::Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();
358  int myPID = comm->getRank();
359  int pid = getpid();
360  char hostname[80];
361  for (int i = 0; i <comm->getSize(); i++) {
362  if (i == myPID) {
363  gethostname(hostname, sizeof(hostname));
364  std::cout << "Host: " << hostname << "\tMPI rank: " << myPID << ",\tPID: " << pid << "\n\tattach " << pid << std::endl;
365  sleep(1);
366  }
367  }
368  if (myPID == 0) {
369  std::cout << "** Enter a character to continue > " << std::endl;
370  char go = ' ';
371  int r = scanf("%c", &go);
372  (void)r;
373  assert(r > 0);
374  }
375  comm->barrier();
376  }
377 #else
378  static void PauseForDebugger() {
379  throw(Exceptions::RuntimeError("MueLu Utils: PauseForDebugger not implemented on Windows."));
380  }
381 #endif
382 
398  static Scalar PowerMethod(const Matrix& A, bool scaleByDiag = true,
399  LocalOrdinal niters = 10, Magnitude tolerance = 1e-2, bool verbose = false, unsigned int seed = 123) {
400  TEUCHOS_TEST_FOR_EXCEPTION(!(A.getRangeMap()->isSameAs(*(A.getDomainMap()))), Exceptions::Incompatible,
401  "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
402 
403  // Create three vectors, fill z with random numbers
404  RCP<Vector> q = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getDomainMap());
405  RCP<Vector> r = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
406  RCP<Vector> z = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
407 
408  z->setSeed(seed); // seed random number generator
409  z->randomize(true);// use Xpetra implementation: -> same results for Epetra and Tpetra
410 
411  Teuchos::Array<Magnitude> norms(1);
412 
413  typedef Teuchos::ScalarTraits<Scalar> STS;
414 
415  const Scalar zero = STS::zero(), one = STS::one();
416 
417  Scalar lambda = zero;
418  Magnitude residual = STS::magnitude(zero);
419 
420  // power iteration
421  RCP<Vector> diagInvVec;
422  if (scaleByDiag) {
423  RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
424  A.getLocalDiagCopy(*diagVec);
425  diagInvVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
426  diagInvVec->reciprocal(*diagVec);
427  }
428 
429  for (int iter = 0; iter < niters; ++iter) {
430  z->norm2(norms); // Compute 2-norm of z
431  q->update(one/norms[0], *z, zero); // Set q = z / normz
432  A.apply(*q, *z); // Compute z = A*q
433  if (scaleByDiag)
434  z->elementWiseMultiply(one, *diagInvVec, *z, zero);
435  lambda = q->dot(*z); // Approximate maximum eigenvalue: lamba = dot(q,z)
436 
437  if (iter % 100 == 0 || iter + 1 == niters) {
438  r->update(1.0, *z, -lambda, *q, zero); // Compute A*q - lambda*q
439  r->norm2(norms);
440  residual = STS::magnitude(norms[0] / lambda);
441  if (verbose) {
442  std::cout << "Iter = " << iter
443  << " Lambda = " << lambda
444  << " Residual of A*q - lambda*q = " << residual
445  << std::endl;
446  }
447  }
448  if (residual < tolerance)
449  break;
450  }
451  return lambda;
452  }
453 
454 
455 
456  static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) {
457  RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(os));
458  return fancy;
459  }
460 
465  static typename Teuchos::ScalarTraits<Scalar>::magnitudeType Distance2(const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& v, LocalOrdinal i0, LocalOrdinal i1) {
466  size_t numVectors = v.getNumVectors();
467 
468  Scalar d = Teuchos::ScalarTraits<Scalar>::zero();
469  for (size_t j = 0; j < numVectors; j++) {
470  Teuchos::ArrayRCP<const Scalar> vv = v.getData(j);
471  d += (vv[i0] - vv[i1])*(vv[i0] - vv[i1]);
472  }
473  return Teuchos::ScalarTraits<Scalar>::magnitude(d);
474  }
475 
488  static Teuchos::ArrayRCP<const bool> DetectDirichletRows(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero()) {
489  LocalOrdinal numRows = A.getNodeNumRows();
490  typedef Teuchos::ScalarTraits<Scalar> STS;
491  ArrayRCP<bool> boundaryNodes(numRows, true);
492  for (LocalOrdinal row = 0; row < numRows; row++) {
493  ArrayView<const LocalOrdinal> indices;
494  ArrayView<const Scalar> vals;
495  A.getLocalRowView(row, indices, vals);
496  size_t nnz = A.getNumEntriesInLocalRow(row);
497  if (nnz > 1)
498  for (size_t col = 0; col < nnz; col++)
499  if ( (indices[col] != row) && STS::magnitude(vals[col]) > tol) {
500  boundaryNodes[row] = false;
501  break;
502  }
503  }
504  return boundaryNodes;
505  }
506 
507 
520  static Teuchos::ArrayRCP<const bool> DetectDirichletRowsExt(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, bool & bHasZeroDiagonal, const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero()) {
521 
522  // assume that there is no zero diagonal in matrix
523  bHasZeroDiagonal = false;
524 
525  Teuchos::RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
526  A.getLocalDiagCopy(*diagVec);
527  Teuchos::ArrayRCP< const Scalar > diagVecData = diagVec->getData(0);
528 
529  LocalOrdinal numRows = A.getNodeNumRows();
530  typedef Teuchos::ScalarTraits<Scalar> STS;
531  ArrayRCP<bool> boundaryNodes(numRows, false);
532  for (LocalOrdinal row = 0; row < numRows; row++) {
533  ArrayView<const LocalOrdinal> indices;
534  ArrayView<const Scalar> vals;
535  A.getLocalRowView(row, indices, vals);
536  size_t nnz = 0; // collect nonzeros in row (excluding the diagonal)
537  bool bHasDiag = false;
538  for (decltype(indices.size()) col = 0; col < indices.size(); col++) {
539  if ( indices[col] != row) {
540  if (STS::magnitude(vals[col] / sqrt(STS::magnitude(diagVecData[row]) * STS::magnitude(diagVecData[col])) ) > tol) {
541  nnz++;
542  }
543  } else bHasDiag = true; // found a diagonal entry
544  }
545  if (bHasDiag == false) bHasZeroDiagonal = true; // we found at least one row without a diagonal
546  else if(nnz == 0) boundaryNodes[row] = true;
547  }
548  return boundaryNodes;
549  }
550 
555  static Scalar Frobenius(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B) {
556  // We check only row maps. Column may be different. One would hope that they are the same, as we typically
557  // calculate frobenius norm of the specified sparsity pattern with an updated matrix from the previous step,
558  // but matrix addition, even when one is submatrix of the other, changes column map (though change may be as
559  // simple as couple of elements swapped)
560  TEUCHOS_TEST_FOR_EXCEPTION(!A.getRowMap()->isSameAs(*B.getRowMap()), Exceptions::Incompatible, "MueLu::CGSolver::Frobenius: row maps are incompatible");
561  TEUCHOS_TEST_FOR_EXCEPTION(!A.isFillComplete() || !B.isFillComplete(), Exceptions::RuntimeError, "Matrices must be fill completed");
562 
563  const Map& AColMap = *A.getColMap();
564  const Map& BColMap = *B.getColMap();
565 
566  Teuchos::ArrayView<const LocalOrdinal> indA, indB;
567  Teuchos::ArrayView<const Scalar> valA, valB;
568  size_t nnzA = 0, nnzB = 0;
569 
570  // We use a simple algorithm
571  // for each row we fill valBAll array with the values in the corresponding row of B
572  // as such, it serves as both sorted array and as storage, so we don't need to do a
573  // tricky problem: "find a value in the row of B corresponding to the specific GID"
574  // Once we do that, we translate LID of entries of row of A to LID of B, and multiply
575  // corresponding entries.
576  // The algorithm should be reasonably cheap, as it does not sort anything, provided
577  // that getLocalElement and getGlobalElement functions are reasonably effective. It
578  // *is* possible that the costs are hidden in those functions, but if maps are close
579  // to linear maps, we should be fine
580  Teuchos::Array<Scalar> valBAll(BColMap.getNodeNumElements());
581 
582  LocalOrdinal invalid = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
583  Scalar zero = Teuchos::ScalarTraits<Scalar> ::zero(), f = zero, gf;
584  size_t numRows = A.getNodeNumRows();
585  for (size_t i = 0; i < numRows; i++) {
586  A.getLocalRowView(i, indA, valA);
587  B.getLocalRowView(i, indB, valB);
588  nnzA = indA.size();
589  nnzB = indB.size();
590 
591  // Set up array values
592  for (size_t j = 0; j < nnzB; j++)
593  valBAll[indB[j]] = valB[j];
594 
595  for (size_t j = 0; j < nnzA; j++) {
596  // The cost of the whole Frobenius dot product function depends on the
597  // cost of the getLocalElement and getGlobalElement functions here.
598  LocalOrdinal ind = BColMap.getLocalElement(AColMap.getGlobalElement(indA[j]));
599  if (ind != invalid)
600  f += valBAll[ind] * valA[j];
601  }
602 
603  // Clean up array values
604  for (size_t j = 0; j < nnzB; j++)
605  valBAll[indB[j]] = zero;
606  }
607 
608  MueLu_sumAll(AColMap.getComm(), f, gf);
609 
610  return gf;
611  }
612 
622  static void SetRandomSeed(const Teuchos::Comm<int> &comm) {
623  // Distribute the seeds evenly in [1,maxint-1]. This guarantees nothing
624  // about where in random number stream we are, but avoids overflow situations
625  // in parallel when multiplying by a PID. It would be better to use
626  // a good parallel random number generator.
627  double one = 1.0;
628  int maxint = INT_MAX; //= 2^31-1 = 2147483647 for 32-bit integers
629  int mySeed = Teuchos::as<int>((maxint-1) * (one -(comm.getRank()+1)/(comm.getSize()+one)) );
630  if (mySeed < 1 || mySeed == maxint) {
631  std::ostringstream errStr;
632  errStr << "Error detected with random seed = " << mySeed << ". It should be in the interval [1,2^31-2].";
633  throw Exceptions::RuntimeError(errStr.str());
634  }
635  std::srand(mySeed);
636  // For Tpetra, we could use Kokkos' random number generator here.
637  Teuchos::ScalarTraits<Scalar>::seedrandom(mySeed);
638  // Epetra
639  // MultiVector::Random() -> Epetra_Util::RandomDouble() -> Epetra_Utils::RandomInt()
640  // Its own random number generator, based on Seed_. Seed_ is initialized in Epetra_Util constructor with std::rand()
641  // So our setting std::srand() affects that too
642  }
643 
644 
645 
646  // Finds the OAZ Dirichlet rows for this matrix
647  // so far only used in IntrepidPCoarsenFactory
648  // TODO check whether we can use DetectDirichletRows instead
649  static void FindDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &A,
650  std::vector<LocalOrdinal>& dirichletRows, bool count_twos_as_dirichlet=false) {
651  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
652  dirichletRows.resize(0);
653  for(size_t i=0; i<A->getNodeNumRows(); i++) {
654  Teuchos::ArrayView<const LocalOrdinal> indices;
655  Teuchos::ArrayView<const Scalar> values;
656  A->getLocalRowView(i,indices,values);
657  int nnz=0;
658  for (size_t j=0; j<(size_t)indices.size(); j++) {
659  if (Teuchos::ScalarTraits<Scalar>::magnitude(values[j]) > Teuchos::ScalarTraits<MT>::eps()) {
660  nnz++;
661  }
662  }
663  if (nnz == 1 || (count_twos_as_dirichlet && nnz == 2)) {
664  dirichletRows.push_back(i);
665  }
666  }
667  }
668 
669  // Applies Ones-and-Zeros to matrix rows
670  static void ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
671  const std::vector<LocalOrdinal>& dirichletRows) {
672  RCP<const Map> Rmap = A->getColMap();
673  RCP<const Map> Cmap = A->getColMap();
674  Scalar one =Teuchos::ScalarTraits<Scalar>::one();
675  Scalar zero =Teuchos::ScalarTraits<Scalar>::zero();
676 
677  for(size_t i=0; i<dirichletRows.size(); i++) {
678  GlobalOrdinal row_gid = Rmap->getGlobalElement(dirichletRows[i]);
679 
680  Teuchos::ArrayView<const LocalOrdinal> indices;
681  Teuchos::ArrayView<const Scalar> values;
682  A->getLocalRowView(dirichletRows[i],indices,values);
683  // NOTE: This won't work with fancy node types.
684  Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
685  for(size_t j=0; j<(size_t)indices.size(); j++) {
686  if(Cmap->getGlobalElement(indices[j])==row_gid)
687  valuesNC[j]=one;
688  else
689  valuesNC[j]=zero;
690  }
691  }
692  }
693 
694  // Zeros out rows
695  static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
696  const std::vector<LocalOrdinal>& dirichletRows) {
697  Scalar zero =Teuchos::ScalarTraits<Scalar>::zero();
698 
699  for(size_t i=0; i<dirichletRows.size(); i++) {
700  Teuchos::ArrayView<const LocalOrdinal> indices;
701  Teuchos::ArrayView<const Scalar> values;
702  A->getLocalRowView(dirichletRows[i],indices,values);
703  // NOTE: This won't work with fancy node types.
704  Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
705  for(size_t j=0; j<(size_t)indices.size(); j++)
706  valuesNC[j]=zero;
707  }
708  }
709 
710  // Finds the OAZ Dirichlet rows for this matrix
711  static void FindDirichletRowsAndPropagateToCols(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &A,
712  Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletRow,
713  Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletCol) {
714 
715  // Make sure A's RowMap == DomainMap
716  if(!A->getRowMap()->isSameAs(*A->getDomainMap())) {
717  throw std::runtime_error("UtilitiesBase::FindDirichletRowsAndPropagateToCols row and domain maps must match.");
718  }
719  RCP<const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer = A->getCrsGraph()->getImporter();
720  bool has_import = !importer.is_null();
721 
722  // Find the Dirichlet rows
723  std::vector<LocalOrdinal> dirichletRows;
724  FindDirichletRows(A,dirichletRows);
725 
726 #if 0
727  printf("[%d] DirichletRow Ids = ",A->getRowMap()->getComm()->getRank());
728  for(size_t i=0; i<(size_t) dirichletRows.size(); i++)
729  printf("%d ",dirichletRows[i]);
730  printf("\n");
731  fflush(stdout);
732 #endif
733  // Allocate all as non-Dirichlet
734  isDirichletRow = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getRowMap(),true);
735  isDirichletCol = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getColMap(),true);
736 
737  Teuchos::ArrayRCP<int> dr_rcp = isDirichletRow->getDataNonConst(0);
738  Teuchos::ArrayView<int> dr = dr_rcp();
739  Teuchos::ArrayRCP<int> dc_rcp = isDirichletCol->getDataNonConst(0);
740  Teuchos::ArrayView<int> dc = dc_rcp();
741  for(size_t i=0; i<(size_t) dirichletRows.size(); i++) {
742  dr[dirichletRows[i]] = 1;
743  if(!has_import) dc[dirichletRows[i]] = 1;
744  }
745 
746  if(has_import)
747  isDirichletCol->doImport(*isDirichletRow,*importer,Xpetra::CombineMode::ADD);
748 
749  }
750 
751  }; // class Utils
752 
753 
754 
756 
757 } //namespace MueLu
758 
759 #define MUELU_UTILITIESBASE_SHORT
760 #endif // MUELU_UTILITIESBASE_DECL_HPP
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedMultiVector
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.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows)
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.
static void FindDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, std::vector< LocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
static Teuchos::ScalarTraits< Scalar >::magnitudeType Distance2(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &v, LocalOrdinal i0, LocalOrdinal i1)
Squared distance between two rows in a multivector.
Namespace for MueLu classes and methods.
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)
#define MueLu_sumAll(rcpComm, in, out)
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.
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Vector
Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > BlockedMap
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Matrix > Crs2Op(RCP< CrsMatrix > Op)
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixWrap
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100)
Extract Matrix Diagonal.
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero())
Detect Dirichlet rows.
Exception throws to report errors in the internal logical of the program.
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows)
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Teuchos::RCP< const Matrix > rcpA)
Extract Matrix Diagonal of lumped matrix.
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.