MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_ProjectorSmoother_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_PROJECTORSMOOTHER_DEF_HPP
47 #define MUELU_PROJECTORSMOOTHER_DEF_HPP
48 
49 #include <Xpetra_Matrix.hpp>
50 #include <Xpetra_MultiVectorFactory.hpp>
51 
52 #include "MueLu_ConfigDefs.hpp"
53 
55 #include "MueLu_Level.hpp"
56 #include "MueLu_Utilities.hpp"
57 #include "MueLu_Monitor.hpp"
58 
59 namespace MueLu {
60 
61  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
63  : coarseSolver_(coarseSolver)
64  { }
65 
66  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
68 
69  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
71  Factory::Input(currentLevel, "A");
72  Factory::Input(currentLevel, "Nullspace");
73 
74  coarseSolver_->DeclareInput(currentLevel);
75  }
76 
77  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
79  FactoryMonitor monitor(*this, "Projector Smoother", currentLevel);
80 
81  coarseSolver_->Setup(currentLevel);
82 
83  if (SmootherPrototype::IsSetup() == true)
84  this->GetOStream(Warnings0) << "MueLu::ProjectorSmoother::Setup(): Setup() has already been called" << std::endl;
85 
86  RCP<Matrix> A = Factory::Get< RCP<Matrix> > (currentLevel, "A");
87  RCP<MultiVector> B = Factory::Get< RCP<MultiVector> >(currentLevel, "Nullspace");
88 
89  int m = B->getNumVectors();
90 
91  // Find which vectors we want to keep
92  Array<Scalar> br(m), bb(m);
93  RCP<MultiVector> R = MultiVectorFactory::Build(B->getMap(), m);
94  A->apply(*B, *R);
95  B-> dot(*R, br);
96  B-> dot(*B, bb);
97 
98  Array<size_t> selectedIndices;
99  for (int i = 0; i < m; i++) {
100  Scalar rayleigh = br[i] / bb[i];
101 
102  if (Teuchos::ScalarTraits<Scalar>::magnitude(rayleigh) < 1e-12)
103  selectedIndices.push_back(i);
104  }
105  this->GetOStream(Runtime0) << "Coarse level orth indices: " << selectedIndices << std::endl;
106 
107 #if defined(HAVE_MUELU_TPETRA) && defined(HAVE_XPETRA_TPETRA)
108 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
109  // Orthonormalize
111  // TAW: Oct 16 2015: subCopy is not part of Xpetra. One should either add it to Xpetra (with an emulator for Epetra)
112  // or replace this call by a local loop. I'm not motivated to do this now...
113  RCP<Tpetra::MultiVector<SC,LO,GO,NO> > Borth = B_->subCopy(selectedIndices); // copy
114  for (int i = 0; i < selectedIndices.size(); i++) {
115  RCP<Tpetra::Vector<SC,LO,GO,NO> > Bi = Borth->getVectorNonConst(i);
116 
117  Scalar dot;
119  for (int k = 0; k < i; k++) { // orthogonalize
120  RCP<const Tpetra::Vector<SC,LO,GO,NO> > Bk = Borth->getVector(k);
121 
122  dot = Bi->dot(*Bk);
123  Bi->update(-dot, *Bk, Teuchos::ScalarTraits<Scalar>::one());
124  }
125 
126  norm = Bi->norm2();
127  Bi->scale(Teuchos::ScalarTraits<Scalar>::one()/norm); // normalize
128  }
129 
130  Borth_ = rcp(static_cast<MultiVector*>(new TpetraMultiVector(Borth)));
131 #else
132  TEUCHOS_TEST_FOR_EXCEPTION(true,Exceptions::RuntimeError,"Tpetra with GO=int not available. The code in ProjectorSmoother should be rewritten!");
133 #endif
134 #endif
135 
137  }
138 
139  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
140  void ProjectorSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
141  coarseSolver_->Apply(X, B, InitialGuessIsZero);
142 
143  int m = Borth_->getNumVectors();
144  int n = X.getNumVectors();
145 
146  RCP<Xpetra::MultiVector<SC,LO,GO,NO> > X_ = Teuchos::rcpFromRef(X);
147  for (int i = 0; i < n; i++) {
148  RCP<Xpetra::Vector<SC,LO,GO,NO> > Xi = X_->getVectorNonConst(i);
149 
150  Array<Scalar> dot(1);
151  for (int k = 0; k < m; k++) { // orthogonalize
152  RCP<const Xpetra::Vector<SC,LO,GO,NO> > Bk = Borth_->getVector(k);
153 
154  Xi->dot(*Bk, dot());
155  Xi->update(-dot[0], *Bk, Teuchos::ScalarTraits<Scalar>::one());
156  }
157  }
158  }
159 
160  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
162  return rcp( new ProjectorSmoother(*this) );
163  }
164 
165  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
167  std::ostringstream out;
169  return out.str();
170  }
171 
172  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
175  out0 << "";
176  }
177 
178 } // namespace MueLu
179 
180 #endif // MUELU_PROJECTORSMOOTHER_DEF_HPP
Important warning messages (one line)
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
std::string description() const
Return a simple one-line description of this object.
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
One-liner description of what is happening.
virtual ~ProjectorSmoother()
Destructor.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
RCP< SmootherPrototype > Copy() const
bool IsSetup() const
Get the state of a smoother prototype.
void Setup(Level &currentLevel)
Set up the direct solver.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
This class enables the elimination of the nullspace component of the solution through the use of proj...
void push_back(const value_type &x)
size_type size() const
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the direct solver. Solves the linear system AX=B using the constructed solver.
Exception throws to report errors in the internal logical of the program.
ProjectorSmoother(RCP< SmootherPrototype > coarseSolver)
Constructor.
void DeclareInput(Level &currentLevel) const
Input.
virtual std::string description() const
Return a simple one-line description of this object.