MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_ProjectorSmoother_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_PROJECTORSMOOTHER_DEF_HPP
11#define MUELU_PROJECTORSMOOTHER_DEF_HPP
12
13#include <Xpetra_Matrix.hpp>
14#include <Xpetra_MultiVectorFactory.hpp>
15
16#include "MueLu_ConfigDefs.hpp"
17
19#include "MueLu_Level.hpp"
20#include "MueLu_Utilities.hpp"
21#include "MueLu_Monitor.hpp"
22
23namespace MueLu {
24
25template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
27 : coarseSolver_(coarseSolver) {}
28
29template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31
32template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
34 Factory::Input(currentLevel, "A");
35 Factory::Input(currentLevel, "Nullspace");
36
37 coarseSolver_->DeclareInput(currentLevel);
38}
39
40template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42 FactoryMonitor monitor(*this, "Projector Smoother", currentLevel);
43
44 coarseSolver_->Setup(currentLevel);
45
46 if (SmootherPrototype::IsSetup() == true)
47 this->GetOStream(Warnings0) << "MueLu::ProjectorSmoother::Setup(): Setup() has already been called" << std::endl;
48
49 RCP<Matrix> A = Factory::Get<RCP<Matrix> >(currentLevel, "A");
50 RCP<MultiVector> B = Factory::Get<RCP<MultiVector> >(currentLevel, "Nullspace");
51
52 int m = B->getNumVectors();
53
54 // Find which vectors we want to keep
55 Array<Scalar> br(m), bb(m);
56 RCP<MultiVector> R = MultiVectorFactory::Build(B->getMap(), m);
57 A->apply(*B, *R);
58 B->dot(*R, br);
59 B->dot(*B, bb);
60
61 Array<size_t> selectedIndices;
62 for (int i = 0; i < m; i++) {
63 Scalar rayleigh = br[i] / bb[i];
64
65 if (Teuchos::ScalarTraits<Scalar>::magnitude(rayleigh) < 1e-12)
66 selectedIndices.push_back(i);
67 }
68 this->GetOStream(Runtime0) << "Coarse level orth indices: " << selectedIndices << std::endl;
69
70#ifdef HAVE_MUELU_TPETRA_INST_INT_INT
71 // Orthonormalize
72 RCP<const Tpetra::MultiVector<SC, LO, GO, NO> > B_ = toTpetra(B);
73 // TAW: Oct 16 2015: subCopy is not part of Xpetra. One should either add it to Xpetra (with an emulator for Epetra)
74 // or replace this call by a local loop. I'm not motivated to do this now...
75 RCP<Tpetra::MultiVector<SC, LO, GO, NO> > Borth = B_->subCopy(selectedIndices); // copy
76 for (int i = 0; i < selectedIndices.size(); i++) {
77 RCP<Tpetra::Vector<SC, LO, GO, NO> > Bi = Borth->getVectorNonConst(i);
78
79 Scalar dot;
80 typename Teuchos::ScalarTraits<Scalar>::magnitudeType norm;
81 for (int k = 0; k < i; k++) { // orthogonalize
82 RCP<const Tpetra::Vector<SC, LO, GO, NO> > Bk = Borth->getVector(k);
83
84 dot = Bi->dot(*Bk);
85 Bi->update(-dot, *Bk, Teuchos::ScalarTraits<Scalar>::one());
86 }
87
88 norm = Bi->norm2();
89 Bi->scale(Teuchos::ScalarTraits<Scalar>::one() / norm); // normalize
90 }
91
92 Borth_ = rcp(static_cast<MultiVector *>(new TpetraMultiVector(Borth)));
93#else
94 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Tpetra with GO=int not available. The code in ProjectorSmoother should be rewritten!");
95#endif
96
98}
99
100template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
101void ProjectorSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero) const {
102 coarseSolver_->Apply(X, B, InitialGuessIsZero);
103
104 int m = Borth_->getNumVectors();
105 int n = X.getNumVectors();
106
107 RCP<Xpetra::MultiVector<SC, LO, GO, NO> > X_ = Teuchos::rcpFromRef(X);
108 for (int i = 0; i < n; i++) {
109 RCP<Xpetra::Vector<SC, LO, GO, NO> > Xi = X_->getVectorNonConst(i);
110
111 Array<Scalar> dot(1);
112 for (int k = 0; k < m; k++) { // orthogonalize
113 RCP<const Xpetra::Vector<SC, LO, GO, NO> > Bk = Borth_->getVector(k);
114
115 Xi->dot(*Bk, dot());
116 Xi->update(-dot[0], *Bk, Teuchos::ScalarTraits<Scalar>::one());
117 }
118 }
119}
120
121template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
122RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> > ProjectorSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Copy() const {
123 return rcp(new ProjectorSmoother(*this));
124}
125
126template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
128 std::ostringstream out;
130 return out.str();
131}
132
133template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134void ProjectorSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream &out, const VerbLevel verbLevel) const {
136 out0 << "";
137}
138
139} // namespace MueLu
140
141#endif // MUELU_PROJECTORSMOOTHER_DEF_HPP
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
MueLu::DefaultScalar Scalar
virtual std::string description() const
Return a simple one-line description of this object.
Exception throws to report errors in the internal logical of the program.
void Input(Level &level, const std::string &varName) const
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
This class enables the elimination of the nullspace component of the solution through the use of proj...
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.
RCP< SmootherPrototype > Copy() const
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
void Setup(Level &currentLevel)
Set up the direct solver.
std::string description() const
Return a simple one-line description of this object.
ProjectorSmoother(RCP< SmootherPrototype > coarseSolver)
Constructor.
void DeclareInput(Level &currentLevel) const
Input.
bool IsSetup() const
Get the state of a smoother prototype.
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line)
@ Runtime0
One-liner description of what is happening.