MueLu  Version of the Day
MueLu_SchurComplementFactory_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 MUELU_SCHURCOMPLEMENTFACTORY_DEF_HPP_
47 #define MUELU_SCHURCOMPLEMENTFACTORY_DEF_HPP_
48 
49 #ifdef HAVE_MUELU_EXPERIMENTAL
50 
51 #include <Xpetra_BlockedCrsMatrix.hpp>
52 #include <Xpetra_MultiVectorFactory.hpp>
53 #include <Xpetra_VectorFactory.hpp>
54 #include <Xpetra_MatrixFactory.hpp>
55 #include <Xpetra_Matrix.hpp>
56 #include <Xpetra_MatrixMatrix.hpp>
57 #include <Xpetra_CrsMatrixWrap.hpp>
58 #include <Xpetra_BlockedCrsMatrix.hpp>
59 #include <Xpetra_CrsMatrix.hpp>
60 #include "MueLu_Level.hpp"
61 #include "MueLu_Monitor.hpp"
62 #include "MueLu_Utilities.hpp"
63 #include "MueLu_HierarchyHelpers.hpp"
64 
65 #include "MueLu_SchurComplementFactory.hpp"
66 
67 namespace MueLu {
68 
69  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71  RCP<ParameterList> validParamList = rcp(new ParameterList());
72 
73  SC one = Teuchos::ScalarTraits<SC>::one();
74 
75  validParamList->set<RCP<const FactoryBase> >("A", NoFactory::getRCP()/*null*/, "Generating factory of the matrix A used for building Schur complement\n"
76  "(must be a 2x2 blocked operator)");
77  validParamList->set<SC> ("omega", one, "Scaling parameter in S = A(1,1) - 1/omega A(1,0) diag{A(0,0)}^{-1} A(0,1)");
78  validParamList->set<bool> ("lumping", false, "Use lumping to construct diag(A(0,0)), i.e. use row sum of the abs values on the diagonal "
79  "as approximation of A00 (and A00^{-1})");
80  validParamList->set<bool> ("fixing", false, "Fix diagonal by replacing small entries with 1.0");
81 
82  return validParamList;
83  }
84 
85  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87  Input(currentLevel, "A");
88  }
89 
90  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
92  FactoryMonitor m(*this, "Build", currentLevel);
93 
94  typedef Teuchos::ScalarTraits<SC> STS;
95  SC zero = STS::zero(), one = STS::one();
96 
97  RCP<Matrix> A = Get<RCP<Matrix> >(currentLevel, "A");
98  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A);
99  TEUCHOS_TEST_FOR_EXCEPTION(bA.is_null(), Exceptions::BadCast,
100  "MueLu::SchurComplementFactory::Build: input matrix A is not of type BlockedCrsMatrix!");
101 
102  RCP<Matrix> A00 = Utilities::Crs2Op(bA->getMatrix(0,0));
103  RCP<Matrix> A01 = Utilities::Crs2Op(bA->getMatrix(0,1));
104  RCP<Matrix> A10 = Utilities::Crs2Op(bA->getMatrix(1,0));
105  RCP<Matrix> A11 = Utilities::Crs2Op(bA->getMatrix(1,1));
106 
107  // TODO move this to BlockedCrsMatrix->getMatrix routine...
108  A00->CreateView("stridedMaps", bA->getRangeMap(0), bA->getDomainMap(0));
109  A01->CreateView("stridedMaps", bA->getRangeMap(0), bA->getDomainMap(1));
110  A10->CreateView("stridedMaps", bA->getRangeMap(1), bA->getDomainMap(0));
111  if (!A11.is_null())
112  A11->CreateView("stridedMaps", bA->getRangeMap(1), bA->getDomainMap(1));
113 
114  const ParameterList& pL = GetParameterList();
115  SC omega = pL.get<Scalar>("omega");
116 
117  // Copy the value of A01 so we can do the left scale.
118  RCP<Matrix> T = MatrixFactory2::BuildCopy(A01);
119 
120  bool lumping = pL.get<bool>("lumping");
121  bool fixing = pL.get<bool>("fixing");
122  ArrayRCP<SC> D;
123  if (!lumping) {
125 
126  if (fixing) {
127  for (size_t k = 0; k < as<size_t>(D.size()); k++)
128  if (STS::magnitude(D[k]) < 1e-4)
129  D[k] = STS::one();
130  }
131 
132  } else {
134  }
135  // Update D to use omega
136  // As D is going to be used as D^{-1}, we must use -omega instead of -one/omega
137  for (size_t k = 0; k < as<size_t>(D.size()); k++)
138  D[k] *= -omega;
139  Utilities::MyOldScaleMatrix(*T, D, true/*doInverse*/, true/*doFillComplete*/, false/*doOptimizeStorage*/);
140 
141  RCP<Matrix> S = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*A10, false, *T, false, GetOStream(Statistics2));
142 
143  if (!A11.is_null()) {
144  T = Teuchos::null;
145  Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TwoMatrixAdd(*A11, false, one, *S, false, one, T, GetOStream(Statistics2));
146  T->fillComplete();
147  S.swap(T);
148  }
149 
150  // NOTE: "A" generated by this factory is actually the Schur complement
151  // matrix, but it is required as all smoothers expect "A"
152  Set(currentLevel, "A", S);
153  }
154 
155 } // namespace MueLu
156 
157 #endif /* HAVE_MUELU_EXPERIMENTAL */
158 #endif /* MUELU_SCHURCOMPLEMENTFACTORY_DEF_HPP_ */
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Exception indicating invalid cast attempted.
static void MyOldScaleMatrix(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
Timer to be used in factories. Similar to Monitor but with additional timers.
static Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Namespace for MueLu classes and methods.
void DeclareInput(Level &currentLevel) const
Input.
Print even more statistics.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void Build(Level &currentLevel) const
Build an object with this factory.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Crs2Op(RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static const RCP< const NoFactory > getRCP()
Static Get() functions.