MueLu  Version of the Day
Thyra_MueLuTpetraQ2Q1PreconditionerFactory_def.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 THYRA_MUELU_TPETRA_Q2Q1PRECONDITIONER_FACTORY_DEF_HPP
47 #define THYRA_MUELU_TPETRA_Q2Q1PRECONDITIONER_FACTORY_DEF_HPP
48 
49 #ifdef HAVE_MUELU_EXPERIMENTAL
50 
52 
53 #include <Thyra_DefaultPreconditioner.hpp>
54 #include <Thyra_TpetraLinearOp.hpp>
55 #include <Thyra_TpetraThyraWrappers.hpp>
56 
57 #include <Teuchos_Ptr.hpp>
58 #include <Teuchos_TestForException.hpp>
59 #include <Teuchos_Assert.hpp>
60 #include <Teuchos_Time.hpp>
61 #include <Teuchos_FancyOStream.hpp>
62 #include <Teuchos_VerbosityLevel.hpp>
63 
64 #include <Teko_Utilities.hpp>
65 
66 #include <Xpetra_BlockedCrsMatrix.hpp>
67 #include <Xpetra_CrsMatrixWrap.hpp>
68 #include <Xpetra_IO.hpp>
69 #include <Xpetra_MapExtractorFactory.hpp>
70 #include <Xpetra_Matrix.hpp>
71 #include <Xpetra_MatrixMatrix.hpp>
72 
73 #include "MueLu.hpp"
74 
75 #include "../../research/q2q1/MueLu_Q2Q1PFactory.hpp"
76 #include "../../research/q2q1/MueLu_Q2Q1uPFactory.hpp"
77 
78 #include "MueLu_AmalgamationFactory.hpp"
79 #include "MueLu_BaseClass.hpp"
80 #include "MueLu_BlockedDirectSolver.hpp"
81 #include "MueLu_BlockedPFactory.hpp"
82 #include "MueLu_BlockedRAPFactory.hpp"
83 #include "MueLu_BraessSarazinSmoother.hpp"
84 #include "MueLu_CoalesceDropFactory.hpp"
85 #include "MueLu_ConstraintFactory.hpp"
87 #include "MueLu_DirectSolver.hpp"
88 #include "MueLu_EminPFactory.hpp"
89 #include "MueLu_FactoryManager.hpp"
90 #include "MueLu_FilteredAFactory.hpp"
91 #include "MueLu_GenericRFactory.hpp"
92 #include "MueLu_Level.hpp"
93 #include "MueLu_PatternFactory.hpp"
94 #include "MueLu_SchurComplementFactory.hpp"
95 #include "MueLu_SmootherFactory.hpp"
96 #include "MueLu_SmootherPrototype.hpp"
97 #include "MueLu_SubBlockAFactory.hpp"
98 #include "MueLu_TpetraOperator.hpp"
99 #include "MueLu_TrilinosSmoother.hpp"
100 
101 #include <string>
102 
103 namespace Thyra {
104 
105 #define MUELU_GPD(name, type, defaultValue) \
106  (paramList.isParameter(name) ? paramList.get<type>(name) : defaultValue)
107 
108  using Teuchos::RCP;
109  using Teuchos::rcp;
110  using Teuchos::rcp_const_cast;
111  using Teuchos::rcp_dynamic_cast;
112  using Teuchos::ParameterList;
113  using Teuchos::ArrayView;
114  using Teuchos::ArrayRCP;
115  using Teuchos::as;
116  using Teuchos::Array;
117 
118  // Constructors/initializers/accessors
119  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121 
122 
123  // Overridden from PreconditionerFactoryBase
124  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
126  typedef Thyra ::TpetraLinearOp<SC,LO,GO,NO> ThyraTpetraLinOp;
127  typedef Tpetra::Operator <SC,LO,GO,NO> TpetraLinOp;
128  typedef Tpetra::CrsMatrix <SC,LO,GO,NO> TpetraCrsMat;
129 
130  const RCP<const LinearOpBase<SC> > fwdOp = fwdOpSrc.getOp();
131  const RCP<const ThyraTpetraLinOp> thyraTpetraFwdOp = rcp_dynamic_cast<const ThyraTpetraLinOp>(fwdOp);
132  const RCP<const TpetraLinOp> tpetraFwdOp = Teuchos::nonnull(thyraTpetraFwdOp) ? thyraTpetraFwdOp->getConstTpetraOperator() : Teuchos::null;
133  const RCP<const TpetraCrsMat> tpetraFwdCrsMat = rcp_dynamic_cast<const TpetraCrsMat>(tpetraFwdOp);
134 
135  return Teuchos::nonnull(tpetraFwdCrsMat);
136  }
137 
138  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139  RCP<PreconditionerBase<Scalar> >
141  return rcp(new DefaultPreconditioner<SC>);
142  }
143 
144  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
146  initializePrec(const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc, PreconditionerBase<Scalar> *prec, const ESupportSolveUse supportSolveUse) const {
147  // Check precondition
148  TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
149  TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
150  TEUCHOS_ASSERT(prec);
151 
152  // Retrieve wrapped concrete Tpetra matrix from FwdOp
153  const RCP<const LinearOpBase<SC> > fwdOp = fwdOpSrc->getOp();
154  TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));
155 
156  typedef Thyra::TpetraLinearOp<SC,LO,GO,NO> ThyraTpetraLinOp;
157  const RCP<const ThyraTpetraLinOp> thyraTpetraFwdOp = rcp_dynamic_cast<const ThyraTpetraLinOp>(fwdOp);
158  TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraTpetraFwdOp));
159 
160  typedef Tpetra::Operator<SC,LO,GO,NO> TpetraLinOp;
161  const RCP<const TpetraLinOp> tpetraFwdOp = thyraTpetraFwdOp->getConstTpetraOperator();
162  TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(tpetraFwdOp));
163 
164  typedef Tpetra::CrsMatrix<SC,LO,GO,NO> TpetraCrsMat;
165  const RCP<const TpetraCrsMat> tpetraFwdCrsMat = rcp_dynamic_cast<const TpetraCrsMat>(tpetraFwdOp);
166  TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(tpetraFwdCrsMat));
167 
168  // Retrieve concrete preconditioner object
169  const Teuchos::Ptr<DefaultPreconditioner<SC> > defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<SC> *>(prec));
170  TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
171 
172  // Workaround since MueLu interface does not accept const matrix as input
173  const RCP<TpetraCrsMat> tpetraFwdCrsMatNonConst = rcp_const_cast<TpetraCrsMat>(tpetraFwdCrsMat);
174 
175  // Create and compute the initial preconditioner
176 
177  // Create a copy, as we may remove some things from the list
178  ParameterList paramList = *paramList_;
179 
180  typedef Tpetra::MultiVector<SC,LO,GO,NO> MultiVector;
181  RCP<MultiVector> coords, nullspace, velCoords, presCoords;
182  ArrayRCP<LO> p2vMap;
183  Teko::LinearOp thA11, thA12, thA21, thA11_9Pt;
184  if (paramList.isType<RCP<MultiVector> >("Coordinates")) { coords = paramList.get<RCP<MultiVector> >("Coordinates"); paramList.remove("Coordinates"); }
185  if (paramList.isType<RCP<MultiVector> >("Nullspace")) { nullspace = paramList.get<RCP<MultiVector> >("Nullspace"); paramList.remove("Nullspace"); }
186  if (paramList.isType<RCP<MultiVector> >("Velcoords")) { velCoords = paramList.get<RCP<MultiVector> >("Velcoords"); paramList.remove("Velcoords"); }
187  if (paramList.isType<RCP<MultiVector> >("Prescoords")) { presCoords = paramList.get<RCP<MultiVector> >("Prescoords"); paramList.remove("Prescoords"); }
188  if (paramList.isType<ArrayRCP<LO> > ("p2vMap")) { p2vMap = paramList.get<ArrayRCP<LO> > ("p2vMap"); paramList.remove("p2vMap"); }
189  if (paramList.isType<Teko::LinearOp> ("A11")) { thA11 = paramList.get<Teko::LinearOp> ("A11"); paramList.remove("A11"); }
190  if (paramList.isType<Teko::LinearOp> ("A12")) { thA12 = paramList.get<Teko::LinearOp> ("A12"); paramList.remove("A12"); }
191  if (paramList.isType<Teko::LinearOp> ("A21")) { thA21 = paramList.get<Teko::LinearOp> ("A21"); paramList.remove("A21"); }
192  if (paramList.isType<Teko::LinearOp> ("A11_9Pt")) { thA11_9Pt = paramList.get<Teko::LinearOp> ("A11_9Pt"); paramList.remove("A11_9Pt"); }
193 
194  typedef MueLu::TpetraOperator<SC,LO,GO,NO> MueLuOperator;
195  const RCP<MueLuOperator> mueluPrecOp = Q2Q1MkPrecond(paramList, velCoords, presCoords, p2vMap, thA11, thA12, thA21, thA11_9Pt);
196 
197  const RCP<LinearOpBase<SC> > thyraPrecOp = Thyra::createLinearOp(RCP<TpetraLinOp>(mueluPrecOp));
198  defaultPrec->initializeUnspecified(thyraPrecOp);
199  }
200 
201  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
203  uninitializePrec(PreconditionerBase<Scalar> *prec, RCP<const LinearOpSourceBase<Scalar> >* fwdOp, ESupportSolveUse* supportSolveUse) const {
204  // Check precondition
205  TEUCHOS_ASSERT(prec);
206 
207  // Retrieve concrete preconditioner object
208  const Teuchos::Ptr<DefaultPreconditioner<SC> > defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<SC> *>(prec));
209  TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
210 
211  if (fwdOp) {
212  // TODO: Implement properly instead of returning default value
213  *fwdOp = Teuchos::null;
214  }
215 
216  if (supportSolveUse) {
217  // TODO: Implement properly instead of returning default value
218  *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
219  }
220 
221  defaultPrec->uninitialize();
222  }
223 
224 
225  // Overridden from ParameterListAcceptor
226  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
228  TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(paramList));
229  paramList_ = paramList;
230  }
231 
232  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
233  RCP<ParameterList>
235  return paramList_;
236  }
237 
238 
239  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
240  RCP<ParameterList>
242  RCP<ParameterList> savedParamList = paramList_;
243  paramList_ = Teuchos::null;
244  return savedParamList;
245  }
246 
247  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
248  RCP<const ParameterList>
250  return paramList_;
251  }
252 
253  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
254  RCP<const ParameterList>
256  static RCP<const ParameterList> validPL;
257 
258  if (validPL.is_null())
259  validPL = rcp(new ParameterList());
260 
261  return validPL;
262  }
263 
264  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
265  RCP<MueLu::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
267  Q2Q1MkPrecond(const ParameterList& paramList,
268  const RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& velCoords,
269  const RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& presCoords,
270  const ArrayRCP<LocalOrdinal>& p2vMap,
271  const Teko::LinearOp& thA11, const Teko::LinearOp& thA12, const Teko::LinearOp& thA21, const Teko::LinearOp& thA11_9Pt) const
272  {
273  using Teuchos::null;
274 
275  typedef Tpetra::CrsMatrix <SC,LO,GO,NO> TP_Crs;
276  typedef Tpetra::Operator <SC,LO,GO,NO> TP_Op;
277 
278  typedef Xpetra::BlockedCrsMatrix <SC,LO,GO,NO> BlockedCrsMatrix;
279  typedef Xpetra::CrsMatrix <SC,LO,GO,NO> CrsMatrix;
280  typedef Xpetra::CrsMatrixWrap <SC,LO,GO,NO> CrsMatrixWrap;
281  typedef Xpetra::MapExtractorFactory <SC,LO,GO,NO> MapExtractorFactory;
282  typedef Xpetra::MapExtractor <SC,LO,GO,NO> MapExtractor;
283  typedef Xpetra::Map <LO,GO,NO> Map;
284  typedef Xpetra::MapFactory <LO,GO,NO> MapFactory;
285  typedef Xpetra::Matrix <SC,LO,GO,NO> Matrix;
286  typedef Xpetra::MatrixFactory <SC,LO,GO,NO> MatrixFactory;
287  typedef Xpetra::StridedMapFactory <LO,GO,NO> StridedMapFactory;
288 
289  typedef MueLu::Hierarchy <SC,LO,GO,NO> Hierarchy;
290 
291  const RCP<const Teuchos::Comm<int> > comm = velCoords->getMap()->getComm();
292 
293  // Pull out Tpetra matrices
294  RCP<Thyra::LinearOpBase<SC> > ThNonConstA11 = rcp_const_cast<Thyra::LinearOpBase<double> >(thA11);
295  RCP<Thyra::LinearOpBase<SC> > ThNonConstA21 = rcp_const_cast<Thyra::LinearOpBase<double> >(thA21);
296  RCP<Thyra::LinearOpBase<SC> > ThNonConstA12 = rcp_const_cast<Thyra::LinearOpBase<double> >(thA12);
297  RCP<Thyra::LinearOpBase<SC> > ThNonConstA11_9Pt = rcp_const_cast<Thyra::LinearOpBase<double> >(thA11_9Pt);
298 
299  RCP<TP_Op> TpetA11 = Thyra::TpetraOperatorVectorExtraction<SC,LO,GO,NO>::getTpetraOperator(ThNonConstA11);
300  RCP<TP_Op> TpetA21 = Thyra::TpetraOperatorVectorExtraction<SC,LO,GO,NO>::getTpetraOperator(ThNonConstA21);
301  RCP<TP_Op> TpetA12 = Thyra::TpetraOperatorVectorExtraction<SC,LO,GO,NO>::getTpetraOperator(ThNonConstA12);
302  RCP<TP_Op> TpetA11_9Pt = Thyra::TpetraOperatorVectorExtraction<SC,LO,GO,NO>::getTpetraOperator(ThNonConstA11_9Pt);
303 
304  RCP<TP_Crs> TpetCrsA11 = rcp_dynamic_cast<TP_Crs>(TpetA11);
305  RCP<TP_Crs> TpetCrsA21 = rcp_dynamic_cast<TP_Crs>(TpetA21);
306  RCP<TP_Crs> TpetCrsA12 = rcp_dynamic_cast<TP_Crs>(TpetA12);
307  RCP<TP_Crs> TpetCrsA11_9Pt = rcp_dynamic_cast<TP_Crs>(TpetA11_9Pt);
308 
309  RCP<Matrix> A_11 = MueLu::TpetraCrs_To_XpetraMatrix(TpetCrsA11);
310  RCP<Matrix> tmp_A_21 = MueLu::TpetraCrs_To_XpetraMatrix(TpetCrsA21); // needs map modification
311  RCP<Matrix> tmp_A_12 = MueLu::TpetraCrs_To_XpetraMatrix(TpetCrsA12); // needs map modification
312  RCP<Matrix> A_11_9Pt = MueLu::TpetraCrs_To_XpetraMatrix(TpetCrsA11_9Pt);
313 
314  Xpetra::global_size_t numVel = A_11->getRowMap()->getNodeNumElements();
315  Xpetra::global_size_t numPres = tmp_A_21->getRowMap()->getNodeNumElements();
316 
317  // Create new A21 with map so that the global indices of the row map starts
318  // from numVel+1 (where numVel is the number of rows in the A11 block)
319  RCP<const Map> domainMap2 = tmp_A_12->getDomainMap();
320  RCP<const Map> rangeMap2 = tmp_A_21->getRangeMap();
321  Xpetra::global_size_t numRows2 = rangeMap2->getNodeNumElements();
322  Xpetra::global_size_t numCols2 = domainMap2->getNodeNumElements();
323  ArrayView<const GO> rangeElem2 = rangeMap2->getNodeElementList();
324  ArrayView<const GO> domainElem2 = domainMap2->getNodeElementList();
325  ArrayView<const GO> rowElem1 = tmp_A_12->getRowMap()->getNodeElementList();
326  ArrayView<const GO> colElem1 = tmp_A_21->getColMap()->getNodeElementList();
327 
328  Xpetra::UnderlyingLib lib = domainMap2->lib();
329  GO indexBase = domainMap2->getIndexBase();
330 
331  Array<GO> newRowElem2(numRows2, 0);
332  for (Xpetra::global_size_t i = 0; i < numRows2; i++)
333  newRowElem2[i] = numVel + rangeElem2[i];
334 
335  RCP<const Map> newRangeMap2 = MapFactory::Build(lib, numRows2, newRowElem2, indexBase, comm);
336 
337  // maybe should be column map???
338  Array<GO> newColElem2(numCols2, 0);
339  for (Xpetra::global_size_t i = 0; i < numCols2; i++)
340  newColElem2[i] = numVel + domainElem2[i];
341 
342  RCP<const Map> newDomainMap2 = MapFactory::Build(lib, numCols2, newColElem2, indexBase, comm);
343 
344  RCP<Matrix> A_12 = MatrixFactory::Build(tmp_A_12->getRangeMap(), newDomainMap2, tmp_A_12->getNodeMaxNumRowEntries());
345  RCP<Matrix> A_21 = MatrixFactory::Build(newRangeMap2, tmp_A_21->getDomainMap(), tmp_A_21->getNodeMaxNumRowEntries());
346 
347  RCP<CrsMatrix> A_11_crs = rcp_dynamic_cast<CrsMatrixWrap>(A_11) ->getCrsMatrix();
348  RCP<CrsMatrix> A_12_crs = rcp_dynamic_cast<CrsMatrixWrap>(A_12) ->getCrsMatrix();
349  RCP<CrsMatrix> A_21_crs = rcp_dynamic_cast<CrsMatrixWrap>(A_21) ->getCrsMatrix();
350  RCP<CrsMatrix> A_11_crs_9Pt = rcp_dynamic_cast<CrsMatrixWrap>(A_11_9Pt)->getCrsMatrix();
351 
352 #if 0
353  RCP<Matrix> A_22 = MatrixFactory::Build(newRangeMap2, newDomainMap2, 1);
354  RCP<CrsMatrix> A_22_crs = rcp_dynamic_cast<CrsMatrixWrap>(A_22) ->getCrsMatrix();
355 
356  // FIXME: why do we need to perturb A_22?
357  Array<SC> smallVal(1, 1.0e-10);
358 
359  // FIXME: could this be sped up using expertStaticFillComplete?
360  // There was an attempt on doing it, but it did not do the proper thing
361  // with empty columns. See git history
362  ArrayView<const LO> inds;
363  ArrayView<const SC> vals;
364  for (LO row = 0; row < as<LO>(numRows2); ++row) {
365  tmp_A_21->getLocalRowView(row, inds, vals);
366 
367  size_t nnz = inds.size();
368  Array<GO> newInds(nnz, 0);
369  for (LO colID = 0; colID < as<LO>(nnz); colID++)
370  newInds[colID] = colElem1[inds[colID]];
371 
372  A_21_crs->insertGlobalValues(newRowElem2[row], newInds, vals);
373  A_22_crs->insertGlobalValues(newRowElem2[row], Array<LO>(1, newRowElem2[row]), smallVal);
374  }
375  A_21_crs->fillComplete(tmp_A_21->getDomainMap(), newRangeMap2);
376  A_22_crs->fillComplete(newDomainMap2, newRangeMap2);
377 #else
378  RCP<Matrix> A_22 = Teuchos::null;
379  RCP<CrsMatrix> A_22_crs = Teuchos::null;
380 
381  ArrayView<const LO> inds;
382  ArrayView<const SC> vals;
383  for (LO row = 0; row < as<LO>(numRows2); ++row) {
384  tmp_A_21->getLocalRowView(row, inds, vals);
385 
386  size_t nnz = inds.size();
387  Array<GO> newInds(nnz, 0);
388  for (LO colID = 0; colID < as<LO>(nnz); colID++)
389  newInds[colID] = colElem1[inds[colID]];
390 
391  A_21_crs->insertGlobalValues(newRowElem2[row], newInds, vals);
392  }
393  A_21_crs->fillComplete(tmp_A_21->getDomainMap(), newRangeMap2);
394 #endif
395 
396  // Create new A12 with map so that the global indices of the ColMap starts
397  // from numVel+1 (where numVel is the number of rows in the A11 block)
398  for (LO row = 0; row < as<LO>(tmp_A_12->getRowMap()->getNodeNumElements()); ++row) {
399  tmp_A_12->getLocalRowView(row, inds, vals);
400 
401  size_t nnz = inds.size();
402  Array<GO> newInds(nnz, 0);
403  for (LO colID = 0; colID < as<LO>(nnz); colID++)
404  newInds[colID] = newColElem2[inds[colID]];
405 
406  A_12_crs->insertGlobalValues(rowElem1[row], newInds, vals);
407  }
408  A_12_crs->fillComplete(newDomainMap2, tmp_A_12->getRangeMap());
409 
410  RCP<Matrix> A_12_abs = Absolute(*A_12);
411  RCP<Matrix> A_21_abs = Absolute(*A_21);
412 
413  // =========================================================================
414  // Preconditioner construction - I (block)
415  // =========================================================================
416  RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
417  Teuchos::FancyOStream& out = *fancy;
418  out.setOutputToRootOnly(0);
419  RCP<Matrix> BBt = Xpetra::MatrixMatrix<SC,LO,GO,NO>::Multiply(*A_21, false, *A_12, false, out);
420  RCP<Matrix> BBt_abs = Xpetra::MatrixMatrix<SC,LO,GO,NO>::Multiply(*A_21_abs, false, *A_12_abs, false, out);
421 
422  SC dropTol = (paramList.get<int>("useFilters") ? paramList.get<double>("tau_1") : 0.00);
423  RCP<Matrix> filteredA = FilterMatrix(*A_11, *A_11, dropTol);
424  RCP<Matrix> filteredB = FilterMatrix(*BBt, *BBt_abs, dropTol);
425 
426  RCP<Matrix> fA_11_crs = rcp_dynamic_cast<CrsMatrixWrap>(filteredA);
427  RCP<Matrix> fA_12_crs = Teuchos::null;
428  RCP<Matrix> fA_21_crs = Teuchos::null;
429  RCP<Matrix> fA_22_crs = rcp_dynamic_cast<CrsMatrixWrap>(filteredB);
430 
431  // Build the large filtered matrix which requires strided maps
432  std::vector<size_t> stridingInfo(1, 1);
433  int stridedBlockId = -1;
434 
435  Array<GO> elementList(numVel+numPres); // Not RCP ... does this get cleared ?
436  Array<GO> velElem = A_12_crs->getRangeMap()->getNodeElementList();
437  Array<GO> presElem = A_21_crs->getRangeMap()->getNodeElementList();
438 
439  for (Xpetra::global_size_t i = 0 ; i < numVel; i++) elementList[i] = velElem[i];
440  for (Xpetra::global_size_t i = numVel; i < numVel+numPres; i++) elementList[i] = presElem[i-numVel];
441  RCP<const Map> fullMap = StridedMapFactory::Build(Xpetra::UseTpetra, numVel+numPres, elementList(), indexBase, stridingInfo, comm);
442 
443  std::vector<RCP<const Map> > partMaps(2);
444  partMaps[0] = StridedMapFactory::Build(Xpetra::UseTpetra, numVel, velElem, indexBase, stridingInfo, comm);
445  partMaps[1] = StridedMapFactory::Build(Xpetra::UseTpetra, numPres, presElem, indexBase, stridingInfo, comm, stridedBlockId, numVel);
446 
447  // Map extractors are necessary for Xpetra's block operators
448  RCP<const MapExtractor> mapExtractor = MapExtractorFactory::Build(fullMap, partMaps);
449  RCP<BlockedCrsMatrix> fA = rcp(new BlockedCrsMatrix(mapExtractor, mapExtractor, 10));
450  fA->setMatrix(0, 0, fA_11_crs);
451  fA->setMatrix(0, 1, fA_12_crs);
452  fA->setMatrix(1, 0, fA_21_crs);
453  fA->setMatrix(1, 1, fA_22_crs);
454  fA->fillComplete();
455 
456  // -------------------------------------------------------------------------
457  // Preconditioner construction - I.a (filtered hierarchy)
458  // -------------------------------------------------------------------------
460  SetDependencyTree(M, paramList);
461 
462  RCP<Hierarchy> H = rcp(new Hierarchy);
463  RCP<MueLu::Level> finestLevel = H->GetLevel(0);
464  finestLevel->Set("A", rcp_dynamic_cast<Matrix>(fA));
465  finestLevel->Set("p2vMap", p2vMap);
466  finestLevel->Set("CoordinatesVelocity", Xpetra::toXpetra(velCoords));
467  finestLevel->Set("CoordinatesPressure", Xpetra::toXpetra(presCoords));
468  finestLevel->Set("AForPat", A_11_9Pt);
469  H->SetMaxCoarseSize(MUELU_GPD("coarse: max size", int, 1));
470 
471  // The first invocation of Setup() builds the hierarchy using the filtered
472  // matrix. This build includes the grid transfers but not the creation of the
473  // smoothers.
474  // NOTE: we need to indicate what should be kept from the first invocation
475  // for the second invocation, which then focuses on building the smoothers
476  // for the unfiltered matrix.
477  H->Keep("P", M.GetFactory("P") .get());
478  H->Keep("R", M.GetFactory("R") .get());
479  H->Keep("Ptent", M.GetFactory("Ptent").get());
480  H->Setup(M, 0, MUELU_GPD("max levels", int, 3));
481 
482 #if 0
483  for (int i = 1; i < H->GetNumLevels(); i++) {
484  RCP<Matrix> P = H->GetLevel(i)->template Get<RCP<Matrix> >("P");
485  RCP<BlockedCrsMatrix> Pcrs = rcp_dynamic_cast<BlockedCrsMatrix>(P);
486  RCP<Matrix> Pp = Pcrs->getMatrix(1,1);
487  RCP<Matrix> Pv = Pcrs->getMatrix(0,0);
488 
489  Xpetra::IO<SC,LO,GO,NO>::Write("Pp_l" + MueLu::toString(i) + ".mm", *Pp);
490  Xpetra::IO<SC,LO,GO,NO>::Write("Pv_l" + MueLu::toString(i) + ".mm", *Pv);
491  }
492 #endif
493 
494  // -------------------------------------------------------------------------
495  // Preconditioner construction - I.b (smoothers for unfiltered matrix)
496  // -------------------------------------------------------------------------
497  std::string smootherType = MUELU_GPD("smoother: type", std::string, "vanka");
498  ParameterList smootherParams;
499  if (paramList.isSublist("smoother: params"))
500  smootherParams = paramList.sublist("smoother: params");
501  M.SetFactory("Smoother", GetSmoother(smootherType, smootherParams, false/*coarseSolver?*/));
502 
503  std::string coarseType = MUELU_GPD("coarse: type", std::string, "direct");
504  ParameterList coarseParams;
505  if (paramList.isSublist("coarse: params"))
506  coarseParams = paramList.sublist("coarse: params");
507  M.SetFactory("CoarseSolver", GetSmoother(coarseType, coarseParams, true/*coarseSolver?*/));
508 
509 #ifdef HAVE_MUELU_DEBUG
510  M.ResetDebugData();
511 #endif
512 
513  RCP<BlockedCrsMatrix> A = rcp(new BlockedCrsMatrix(mapExtractor, mapExtractor, 10));
514  A->setMatrix(0, 0, A_11);
515  A->setMatrix(0, 1, A_12);
516  A->setMatrix(1, 0, A_21);
517  A->setMatrix(1, 1, A_22);
518  A->fillComplete();
519 
520  H->GetLevel(0)->Set("A", rcp_dynamic_cast<Matrix>(A));
521 
522  H->Setup(M, 0, H->GetNumLevels());
523 
524  return rcp(new MueLu::TpetraOperator<SC,LO,GO,NO>(H));
525  }
526 
527  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
528  RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
530  FilterMatrix(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Pattern, Scalar dropTol) const {
531  typedef Xpetra::Matrix<SC,LO,GO,NO> Matrix;
532  typedef MueLu::AmalgamationFactory<SC,LO,GO,NO> AmalgamationFactory;
533  typedef MueLu::CoalesceDropFactory<SC,LO,GO,NO> CoalesceDropFactory;
534  typedef MueLu::FactoryManager<SC,LO,GO,NO> FactoryManager;
535  typedef MueLu::FilteredAFactory<SC,LO,GO,NO> FilteredAFactory;
536  typedef MueLu::GraphBase<LO,GO,NO> GraphBase;
537 
538  RCP<GraphBase> filteredGraph;
539  {
540  // Get graph pattern for the pattern matrix
541  MueLu::Level level;
542  level.SetLevelID(1);
543 
544  level.Set<RCP<Matrix> >("A", rcpFromRef(Pattern));
545 
546  RCP<AmalgamationFactory> amalgFactory = rcp(new AmalgamationFactory());
547 
548  RCP<CoalesceDropFactory> dropFactory = rcp(new CoalesceDropFactory());
549  ParameterList dropParams = *(dropFactory->GetValidParameterList());
550  dropParams.set("lightweight wrap", true);
551  dropParams.set("aggregation: drop scheme", "classical");
552  dropParams.set("aggregation: drop tol", dropTol);
553  // dropParams.set("Dirichlet detection threshold", <>);
554  dropFactory->SetParameterList(dropParams);
555  dropFactory->SetFactory("UnAmalgamationInfo", amalgFactory);
556 
557  // Build
558  level.Request("Graph", dropFactory.get());
559  dropFactory->Build(level);
560 
561  level.Get("Graph", filteredGraph, dropFactory.get());
562  }
563 
564  RCP<Matrix> filteredA;
565  {
566  // Filter the original matrix, not the pattern one
567  MueLu::Level level;
568  level.SetLevelID(1);
569 
570  level.Set("A", rcpFromRef(A));
571  level.Set("Graph", filteredGraph);
572  level.Set("Filtering", true);
573 
574  RCP<FilteredAFactory> filterFactory = rcp(new FilteredAFactory());
575  ParameterList filterParams = *(filterFactory->GetValidParameterList());
576  // We need a graph that has proper structure in it. Therefore, we need to
577  // drop older pattern, i.e. not to reuse it
578  filterParams.set("filtered matrix: reuse graph", false);
579  filterParams.set("filtered matrix: use lumping", false);
580  filterFactory->SetParameterList(filterParams);
581 
582  // Build
583  level.Request("A", filterFactory.get());
584  filterFactory->Build(level);
585 
586  level.Get("A", filteredA, filterFactory.get());
587  }
588 
589  // Zero out row sums by fixing the diagonal
590  filteredA->resumeFill();
591  size_t numRows = filteredA->getRowMap()->getNodeNumElements();
592  for (size_t i = 0; i < numRows; i++) {
593  ArrayView<const LO> inds;
594  ArrayView<const SC> vals;
595  filteredA->getLocalRowView(i, inds, vals);
596 
597  size_t nnz = inds.size();
598 
599  Array<SC> valsNew = vals;
600 
601  LO diagIndex = -1;
602  SC diag = Teuchos::ScalarTraits<SC>::zero();
603  for (size_t j = 0; j < inds.size(); j++) {
604  diag += vals[j];
605  if (inds[j] == i)
606  diagIndex = j;
607  }
608  TEUCHOS_TEST_FOR_EXCEPTION(diagIndex == -1, MueLu::Exceptions::RuntimeError,
609  "No diagonal found");
610  if (nnz <= 1)
611  continue;
612 
613  valsNew[diagIndex] -= diag;
614 
615  filteredA->replaceLocalValues(i, inds, valsNew);
616  }
617  filteredA->fillComplete();
618 
619  return filteredA;
620  }
621 
622  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
623  void
626  typedef MueLu::BlockedPFactory <SC,LO,GO,NO> BlockedPFactory;
627  typedef MueLu::GenericRFactory <SC,LO,GO,NO> GenericRFactory;
628  typedef MueLu::BlockedRAPFactory <SC,LO,GO,NO> BlockedRAPFactory;
629  typedef MueLu::SmootherFactory <SC,LO,GO,NO> SmootherFactory;
630  typedef MueLu::BlockedDirectSolver <SC,LO,GO,NO> BlockedDirectSolver;
631  typedef MueLu::FactoryManager <SC,LO,GO,NO> FactoryManager;
632 
633  // Pressure and velocity dependency trees are identical. The only
634  // difference is that pressure has to go first, so that velocity can use
635  // some of pressure data
636  RCP<FactoryManager> M11 = rcp(new FactoryManager()), M22 = rcp(new FactoryManager());
637  SetBlockDependencyTree(*M11, 0, 0, "velocity", paramList);
638  SetBlockDependencyTree(*M22, 1, 1, "pressure", paramList);
639 
640  RCP<BlockedPFactory> PFact = rcp(new BlockedPFactory());
641  ParameterList pParamList = *(PFact->GetValidParameterList());
642  pParamList.set("backwards", true); // do pressure first
643  PFact->SetParameterList(pParamList);
644  PFact->AddFactoryManager(M11);
645  PFact->AddFactoryManager(M22);
646  M.SetFactory("P", PFact);
647 
648  RCP<GenericRFactory> RFact = rcp(new GenericRFactory());
649  RFact->SetFactory("P", PFact);
650  M.SetFactory("R", RFact);
651 
652  RCP<MueLu::Factory > AcFact = rcp(new BlockedRAPFactory());
653  AcFact->SetFactory("R", RFact);
654  AcFact->SetFactory("P", PFact);
655  M.SetFactory("A", AcFact);
656 
657  // Smoothers will be set later
658  M.SetFactory("Smoother", Teuchos::null);
659 
660  RCP<MueLu::Factory> coarseFact = rcp(new SmootherFactory(rcp(new BlockedDirectSolver()), Teuchos::null));
661  // M.SetFactory("CoarseSolver", coarseFact);
662  M.SetFactory("CoarseSolver", Teuchos::null);
663  }
664 
665  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
666  void
668  SetBlockDependencyTree(MueLu::FactoryManager<Scalar,LocalOrdinal,GlobalOrdinal,Node>& M, LocalOrdinal row, LocalOrdinal col, const std::string& mode, const ParameterList& paramList) const {
669  typedef MueLu::ConstraintFactory <SC,LO,GO,NO> ConstraintFactory;
670  typedef MueLu::EminPFactory <SC,LO,GO,NO> EminPFactory;
671  typedef MueLu::GenericRFactory <SC,LO,GO,NO> GenericRFactory;
672  typedef MueLu::PatternFactory <SC,LO,GO,NO> PatternFactory;
673  typedef MueLu::Q2Q1PFactory <SC,LO,GO,NO> Q2Q1PFactory;
674  typedef MueLu::Q2Q1uPFactory <SC,LO,GO,NO> Q2Q1uPFactory;
675  typedef MueLu::SubBlockAFactory <SC,LO,GO,NO> SubBlockAFactory;
676 
677  RCP<SubBlockAFactory> AFact = rcp(new SubBlockAFactory());
678  AFact->SetFactory ("A", MueLu::NoFactory::getRCP());
679  AFact->SetParameter("block row", Teuchos::ParameterEntry(row));
680  AFact->SetParameter("block col", Teuchos::ParameterEntry(col));
681  M.SetFactory("A", AFact);
682 
683  RCP<MueLu::Factory> Q2Q1Fact;
684 
685  const bool isStructured = false;
686 
687  if (isStructured) {
688  Q2Q1Fact = rcp(new Q2Q1PFactory);
689 
690  } else {
691  Q2Q1Fact = rcp(new Q2Q1uPFactory);
692  ParameterList q2q1ParamList = *(Q2Q1Fact->GetValidParameterList());
693  q2q1ParamList.set("mode", mode);
694  if (paramList.isParameter("dump status"))
695  q2q1ParamList.set("dump status", paramList.get<bool>("dump status"));
696  if (paramList.isParameter("phase2"))
697  q2q1ParamList.set("phase2", paramList.get<bool>("phase2"));
698  if (paramList.isParameter("tau_2"))
699  q2q1ParamList.set("tau_2", paramList.get<double>("tau_2"));
700  Q2Q1Fact->SetParameterList(q2q1ParamList);
701  }
702  Q2Q1Fact->SetFactory("A", AFact);
703  M.SetFactory("Ptent", Q2Q1Fact);
704 
705  RCP<PatternFactory> patternFact = rcp(new PatternFactory);
706  ParameterList patternParams = *(patternFact->GetValidParameterList());
707  // Our prolongator constructs the exact pattern we are going to use,
708  // therefore we do not expand it
709  patternParams.set("emin: pattern order", 0);
710  patternFact->SetParameterList(patternParams);
711  patternFact->SetFactory("A", AFact);
712  patternFact->SetFactory("P", Q2Q1Fact);
713  M.SetFactory("Ppattern", patternFact);
714 
715  RCP<ConstraintFactory> CFact = rcp(new ConstraintFactory);
716  CFact->SetFactory("Ppattern", patternFact);
717  M.SetFactory("Constraint", CFact);
718 
719  RCP<EminPFactory> EminPFact = rcp(new EminPFactory());
720  ParameterList eminParams = *(EminPFact->GetValidParameterList());
721  if (paramList.isParameter("emin: num iterations"))
722  eminParams.set("emin: num iterations", paramList.get<int>("emin: num iterations"));
723  if (mode == "pressure") {
724  eminParams.set("emin: iterative method", "cg");
725  } else {
726  eminParams.set("emin: iterative method", "gmres");
727  if (paramList.isParameter("emin: iterative method"))
728  eminParams.set("emin: iterative method", paramList.get<std::string>("emin: iterative method"));
729  }
730  EminPFact->SetParameterList(eminParams);
731  EminPFact->SetFactory("A", AFact);
732  EminPFact->SetFactory("Constraint", CFact);
733  EminPFact->SetFactory("P", Q2Q1Fact);
734  M.SetFactory("P", EminPFact);
735 
736  if (mode == "velocity" && (!paramList.isParameter("velocity: use transpose") || paramList.get<bool>("velocity: use transpose") == false)) {
737  // Pressure system is symmetric, so it does not matter
738  // Velocity system may benefit from running emin in restriction mode (with A^T)
739  RCP<GenericRFactory> RFact = rcp(new GenericRFactory());
740  RFact->SetFactory("P", EminPFact);
741  M.SetFactory("R", RFact);
742  }
743  }
744 
745  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
746  RCP<MueLu::FactoryBase>
748  GetSmoother(const std::string& type, const ParameterList& paramList, bool coarseSolver) const {
749  typedef Teuchos::ParameterEntry ParameterEntry;
750 
751  typedef MueLu::BlockedDirectSolver <SC,LO,GO,NO> BlockedDirectSolver;
752  typedef MueLu::BraessSarazinSmoother <SC,LO,GO,NO> BraessSarazinSmoother;
753  typedef MueLu::DirectSolver <SC,LO,GO,NO> DirectSolver;
754  typedef MueLu::FactoryManager <SC,LO,GO,NO> FactoryManager;
755  typedef MueLu::SchurComplementFactory<SC,LO,GO,NO> SchurComplementFactory;
756  typedef MueLu::SmootherFactory <SC,LO,GO,NO> SmootherFactory;
757  typedef MueLu::SmootherPrototype <SC,LO,GO,NO> SmootherPrototype;
758  typedef MueLu::TrilinosSmoother <SC,LO,GO,NO> TrilinosSmoother;
759 
760  RCP<SmootherPrototype> smootherPrototype;
761  if (type == "none") {
762  return Teuchos::null;
763 
764  } else if (type == "vanka") {
765  // Set up Vanka smoothing via a combination of Schwarz and block relaxation.
766  ParameterList schwarzList;
767  schwarzList.set("schwarz: overlap level", as<int>(0));
768  schwarzList.set("schwarz: zero starting solution", false);
769  schwarzList.set("subdomain solver name", "Block_Relaxation");
770 
771  ParameterList& innerSolverList = schwarzList.sublist("subdomain solver parameters");
772  innerSolverList.set("partitioner: type", "user");
773  innerSolverList.set("partitioner: overlap", MUELU_GPD("partitioner: overlap", int, 1));
774  innerSolverList.set("relaxation: type", MUELU_GPD("relaxation: type", std::string, "Gauss-Seidel"));
775  innerSolverList.set("relaxation: sweeps", MUELU_GPD("relaxation: sweeps", int, 1));
776  innerSolverList.set("relaxation: damping factor", MUELU_GPD("relaxation: damping factor", double, 0.5));
777  innerSolverList.set("relaxation: zero starting solution", false);
778  // innerSolverList.set("relaxation: backward mode", MUELU_GPD("relaxation: backward mode", bool, true); NOT SUPPORTED YET
779 
780  std::string ifpackType = "SCHWARZ";
781 
782  smootherPrototype = rcp(new TrilinosSmoother(ifpackType, schwarzList));
783 
784  } else if (type == "schwarz") {
785 
786  std::string ifpackType = "SCHWARZ";
787 
788  smootherPrototype = rcp(new TrilinosSmoother(ifpackType, paramList));
789 
790  } else if (type == "braess-sarazin") {
791  // Define smoother/solver for BraessSarazin
792  SC omega = MUELU_GPD("bs: omega", double, 1.0);
793  bool lumping = MUELU_GPD("bs: lumping", bool, false);
794 
795  RCP<SchurComplementFactory> schurFact = rcp(new SchurComplementFactory());
796  schurFact->SetParameter("omega", ParameterEntry(omega));
797  schurFact->SetParameter("lumping", ParameterEntry(lumping));
798  schurFact->SetFactory ("A", MueLu::NoFactory::getRCP());
799 
800  // Schur complement solver
801  RCP<SmootherPrototype> schurSmootherPrototype;
802  std::string schurSmootherType = (paramList.isParameter("schur smoother: type") ? paramList.get<std::string>("schur smoother: type") : "RELAXATION");
803  if (schurSmootherType == "RELAXATION") {
804  ParameterList schurSmootherParams = paramList.sublist("schur smoother: params");
805  // schurSmootherParams.set("relaxation: damping factor", omega);
806  schurSmootherPrototype = rcp(new TrilinosSmoother(schurSmootherType, schurSmootherParams));
807  } else {
808  schurSmootherPrototype = rcp(new DirectSolver());
809  }
810  schurSmootherPrototype->SetFactory("A", schurFact);
811 
812  RCP<SmootherFactory> schurSmootherFact = rcp(new SmootherFactory(schurSmootherPrototype));
813 
814  // Define temporary FactoryManager that is used as input for BraessSarazin smoother
815  RCP<FactoryManager> braessManager = rcp(new FactoryManager());
816  braessManager->SetFactory("A", schurFact); // SchurComplement operator for correction step (defined as "A")
817  braessManager->SetFactory("Smoother", schurSmootherFact); // solver/smoother for correction step
818  braessManager->SetFactory("PreSmoother", schurSmootherFact);
819  braessManager->SetFactory("PostSmoother", schurSmootherFact);
820  braessManager->SetIgnoreUserData(true); // always use data from factories defined in factory manager
821 
822  smootherPrototype = rcp(new BraessSarazinSmoother());
823  smootherPrototype->SetParameter("Sweeps", ParameterEntry(MUELU_GPD("bs: sweeps", int, 1)));
824  smootherPrototype->SetParameter("lumping", ParameterEntry(lumping));
825  smootherPrototype->SetParameter("Damping factor", ParameterEntry(omega));
826  smootherPrototype->SetParameter("q2q1 mode", ParameterEntry(true));
827  rcp_dynamic_cast<BraessSarazinSmoother>(smootherPrototype)->AddFactoryManager(braessManager, 0); // set temporary factory manager in BraessSarazin smoother
828 
829  } else if (type == "ilu") {
830  std::string ifpackType = "RILUK";
831 
832  smootherPrototype = rcp(new TrilinosSmoother(ifpackType, paramList));
833 
834  } else if (type == "direct") {
835  smootherPrototype = rcp(new BlockedDirectSolver());
836 
837  } else {
838  throw MueLu::Exceptions::RuntimeError("Unknown smoother type: \"" + type + "\"");
839  }
840 
841  return coarseSolver ? rcp(new SmootherFactory(smootherPrototype, Teuchos::null)) : rcp(new SmootherFactory(smootherPrototype));
842  }
843 
844  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
845  RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
846  MueLuTpetraQ2Q1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Absolute(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A) const {
847  typedef Xpetra::CrsMatrix <SC,LO,GO,NO> CrsMatrix;
848  typedef Xpetra::CrsMatrixWrap<SC,LO,GO,NO> CrsMatrixWrap;
849  typedef Xpetra::Matrix <SC,LO,GO,NO> Matrix;
850 
851  const CrsMatrixWrap& Awrap = dynamic_cast<const CrsMatrixWrap&>(A);
852 
853  ArrayRCP<const size_t> iaA;
854  ArrayRCP<const LO> jaA;
855  ArrayRCP<const SC> valA;
856  Awrap.getCrsMatrix()->getAllValues(iaA, jaA, valA);
857 
858  ArrayRCP<size_t> iaB (iaA .size());
859  ArrayRCP<LO> jaB (jaA .size());
860  ArrayRCP<SC> valB(valA.size());
861  for (int i = 0; i < iaA .size(); i++) iaB [i] = iaA[i];
862  for (int i = 0; i < jaA .size(); i++) jaB [i] = jaA[i];
863  for (int i = 0; i < valA.size(); i++) valB[i] = Teuchos::ScalarTraits<SC>::magnitude(valA[i]);
864 
865  RCP<Matrix> B = rcp(new CrsMatrixWrap(A.getRowMap(), A.getColMap(), 0, Xpetra::StaticProfile));
866  RCP<CrsMatrix> Bcrs = rcp_dynamic_cast<CrsMatrixWrap>(B)->getCrsMatrix();
867  Bcrs->setAllValues(iaB, jaB, valB);
868  Bcrs->expertStaticFillComplete(A.getDomainMap(), A.getRangeMap());
869 
870  return B;
871  }
872 
873  // Public functions overridden from Teuchos::Describable
874  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
876  return "Thyra::MueLuTpetraQ2Q1PreconditionerFactory";
877  }
878 
879 } // namespace Thyra
880 
881 #endif
882 #endif // ifdef THYRA_MUELU_TPETRA_Q2Q1PRECONDITIONER_FACTORY_DEF_HPP
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
This class specifies the default factory that should generate some data on a Level if the data does n...
Teuchos::RCP< MueLu::TpetraOperator< SC, LO, GO, NO > > Q2Q1MkPrecond(const ParameterList &paramList, const Teuchos::RCP< Tpetra::MultiVector< SC, LO, GO, NO > > &velCoords, const Teuchos::RCP< Tpetra::MultiVector< SC, LO, GO, NO > > &presCoords, const Teuchos::ArrayRCP< LO > &p2vMap, const Teko::LinearOp &thA11, const Teko::LinearOp &thA12, const Teko::LinearOp &thA21, const Teko::LinearOp &thA11_9Pt) const
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level->Get< RCP<Matrix> >("A", factory) if factory == NULL => use default factory.
BraessSarazin smoother for 2x2 block matrices.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Teuchos::RCP< Xpetra::Matrix< SC, LO, GO, NO > > FilterMatrix(Xpetra::Matrix< SC, LO, GO, NO > &A, Xpetra::Matrix< SC, LO, GO, NO > &Pattern, SC dropTol) const
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
Factory for building blocked, segregated prolongation operators.
bool isCompatible(const LinearOpSourceBase< SC > &fwdOp) const
Class that encapsulates external library smoothers.
Factory for building coarse matrices.
Base class for smoother prototypes.
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
void SetBlockDependencyTree(MueLu::FactoryManager< SC, LO, GO, NO > &M, LO row, LO col, const std::string &mode, const ParameterList &paramList) const
Factory for building restriction operators using a prolongator factory.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
Factory for building a thresholded operator.
AmalgamationFactory for subblocks of strided map based amalgamation data.
Factory for building the constraint operator.
direct solver for nxn blocked matrices
Teuchos::RCP< Xpetra::Matrix< SC, LO, GO, NO > > Absolute(const Xpetra::Matrix< SC, LO, GO, NO > &A) const
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
Various adapters that will create a MueLu preconditioner that is a Tpetra::Operator.
void SetFactory(const std::string &varName, const RCP< const FactoryBase > &factory)
Set Factory.
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Get factory associated with a particular data name.
MueLu representation of a graph.
Factory for building the Schur Complement for a 2x2 block matrix.
Factory for creating a graph base on a given matrix.
RCP< MueLu::FactoryBase > GetSmoother(const std::string &type, const ParameterList &paramList, bool coarseSolver) const
void SetLevelID(int levelID)
Set level number.
Definition: MueLu_Level.cpp:78
Factory for building nonzero patterns for energy minimization.
RCP< Xpetra::Matrix< SC, LO, GO, NO > > TpetraCrs_To_XpetraMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > &Atpetra)
Factory for building Energy Minimization prolongators.
Wraps an existing MueLu::Hierarchy as a Tpetra::Operator.
void SetDependencyTree(MueLu::FactoryManager< SC, LO, GO, NO > &M, const ParameterList &paramList) const
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
Exception throws to report errors in the internal logical of the program.
Factory for building filtered matrices using filtered graphs.
void initializePrec(const Teuchos::RCP< const LinearOpSourceBase< SC > > &fwdOp, PreconditionerBase< SC > *prec, const ESupportSolveUse supportSolveUse) const
#define MUELU_GPD(name, type, defaultValue)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
void Request(const FactoryBase &factory)
Increment the storage counter for all the inputs of a factory.
void uninitializePrec(PreconditionerBase< SC > *prec, Teuchos::RCP< const LinearOpSourceBase< SC > > *fwdOp, ESupportSolveUse *supportSolveUse) const
static const RCP< const NoFactory > getRCP()
Static Get() functions.