MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_IndefBlockedDiagonalSmoother_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/*
47 * MueLu_IndefBlockedDiagonalSmoother_def.hpp
48 *
49 * Created on: 13 May 2014
50 * Author: wiesner
51 */
52
53#ifndef MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_DEF_HPP_
54#define MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_DEF_HPP_
55
56#include "Teuchos_ArrayViewDecl.hpp"
57#include "Teuchos_ScalarTraits.hpp"
58
59#include "MueLu_ConfigDefs.hpp"
60
61#include <Xpetra_Matrix.hpp>
62#include <Xpetra_CrsMatrixWrap.hpp>
63#include <Xpetra_BlockedCrsMatrix.hpp>
64#include <Xpetra_MultiVectorFactory.hpp>
65#include <Xpetra_VectorFactory.hpp>
66#include <Xpetra_ReorderedBlockedCrsMatrix.hpp>
67
69#include "MueLu_Level.hpp"
70#include "MueLu_Monitor.hpp"
71#include "MueLu_HierarchyUtils.hpp"
73
74// include files for default FactoryManager
75#include "MueLu_SchurComplementFactory.hpp"
76#include "MueLu_FactoryManager.hpp"
77
78namespace MueLu {
79
80 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
85
86 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
88
89 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91 RCP<ParameterList> validParamList = rcp(new ParameterList());
92
93 validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A (must be a 2x2 block matrix)");
94 validParamList->set< Scalar > ("Damping factor", 1.0, "Damping/Scaling factor");
95 validParamList->set< LocalOrdinal > ("Sweeps", 1, "Number of SIMPLE sweeps (default = 1)");
96 //validParamList->set< bool > ("UseSIMPLEC", false, "Use SIMPLEC instead of SIMPLE (default = false)");
97
98 return validParamList;
99 }
100
102 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
104 TEUCHOS_TEST_FOR_EXCEPTION(pos < 0, Exceptions::RuntimeError, "MueLu::IndefBlockedDiagonalSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
105
106 size_t myPos = Teuchos::as<size_t>(pos);
107
108 if (myPos < FactManager_.size()) {
109 // replace existing entris in FactManager_ vector
110 FactManager_.at(myPos) = FactManager;
111 } else if( myPos == FactManager_.size()) {
112 // add new Factory manager in the end of the vector
113 FactManager_.push_back(FactManager);
114 } else { // if(myPos > FactManager_.size())
115 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
116 *out << "Warning: cannot add new FactoryManager at proper position " << pos << ". The FactoryManager is just appended to the end. Check this!" << std::endl;
117
118 // add new Factory manager in the end of the vector
119 FactManager_.push_back(FactManager);
120 }
121 }
122
123 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
125 currentLevel.DeclareInput("A",this->GetFactory("A").get());
126
127 TEUCHOS_TEST_FOR_EXCEPTION(FactManager_.size() != 2, Exceptions::RuntimeError,"MueLu::IndefBlockedDiagonalSmoother::DeclareInput: You have to declare two FactoryManagers with a \"Smoother\" object: One for predicting the primary variable and one for the SchurComplement system. The smoother for the SchurComplement system needs a SchurComplementFactory as input for variable \"A\"!");
128
129 // loop over all factory managers for the subblocks of blocked operator A
130 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
131 for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
132 SetFactoryManager currentSFM (rcpFromRef(currentLevel), *it);
133
134 // request "Smoother" for current subblock row.
135 currentLevel.DeclareInput("PreSmoother",(*it)->GetFactory("Smoother").get());
136 }
137
138 }
139
140 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
142 FactoryMonitor m(*this, "Setup for indefinite blocked diagonal smoother", currentLevel);
143
144 if (SmootherPrototype::IsSetup() == true)
145 this->GetOStream(Warnings0) << "MueLu::IndefBlockedDiagonalSmoother::Setup(): Setup() has already been called";
146
147 // extract blocked operator A from current level
148 A_ = Factory::Get<RCP<Matrix> > (currentLevel, "A");
149
150 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
151 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::IndefBlockedDiagonalSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
152
153 // store map extractors
154 rangeMapExtractor_ = bA->getRangeMapExtractor();
155 domainMapExtractor_ = bA->getDomainMapExtractor();
156
157 // Store the blocks in local member variables
158 Teuchos::RCP<Matrix> A00 = bA->getMatrix(0, 0);
159 Teuchos::RCP<Matrix> A01 = bA->getMatrix(0, 1);
160 Teuchos::RCP<Matrix> A10 = bA->getMatrix(1, 0);
161 Teuchos::RCP<Matrix> A11 = bA->getMatrix(1, 1);
162
163 F_ = A00;
164 Z_ = A11;
165
166 /*const ParameterList & pL = Factory::GetParameterList();
167 bool bSIMPLEC = pL.get<bool>("UseSIMPLEC");
168
169 // Create the inverse of the diagonal of F
170 RCP<Vector> diagFVector = VectorFactory::Build(F_->getRowMap());
171 if(!bSIMPLEC) {
172 F_->getLocalDiagCopy(*diagFVector); // extract diagonal of F
173 diagFVector->reciprocal(*diagFVector); // build reciprocal
174 } else {
175 const RCP<const Map> rowmap = F_->getRowMap();
176 size_t locSize = rowmap->getLocalNumElements();
177 Teuchos::ArrayRCP<SC> diag = diagFVector->getDataNonConst(0);
178 Teuchos::ArrayView<const LO> cols;
179 Teuchos::ArrayView<const SC> vals;
180 for (size_t i=0; i<locSize; ++i) { // loop over rows
181 F_->getLocalRowView(i,cols,vals);
182 Scalar absRowSum = Teuchos::ScalarTraits<Scalar>::zero();
183 for (LO j=0; j<cols.size(); ++j) { // loop over cols
184 absRowSum += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
185 }
186 diag[i] = absRowSum;
187 }
188 diagFVector->reciprocal(*diagFVector); // build reciprocal
189 }
190 diagFinv_ = diagFVector;*/
191
192 // Set the Smoother
193 // carefully switch to the SubFactoryManagers (defined by the users)
194 {
195 RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
196 SetFactoryManager currentSFM (rcpFromRef(currentLevel), velpredictFactManager);
197 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >("PreSmoother",velpredictFactManager->GetFactory("Smoother").get());
198 }
199 {
200 RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
201 SetFactoryManager currentSFM (rcpFromRef(currentLevel), schurFactManager);
202 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >("PreSmoother", schurFactManager->GetFactory("Smoother").get());
203 }
205 }
206
207 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
208 void IndefBlockedDiagonalSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool /* InitialGuessIsZero */) const
209 {
210 TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::IndefBlockedDiagonalSmoother::Apply(): Setup() has not been called");
211
212 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
213
214 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
215
216 // extract parameters from internal parameter list
217 const ParameterList & pL = Factory::GetParameterList();
218 LocalOrdinal nSweeps = pL.get<LocalOrdinal>("Sweeps");
219 Scalar omega = pL.get<Scalar>("Damping factor");
220
221 bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
222 bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
223
224 // wrap current solution vector in RCP
225 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
226 RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
227
228 // make sure that both rcpX and rcpB are BlockedMultiVector objects
229 bool bCopyResultX = false;
230 bool bReorderX = false;
231 bool bReorderB = false;
232 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
233 MUELU_TEST_FOR_EXCEPTION(bA.is_null() == true, Exceptions::RuntimeError, "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
234 RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
235 RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
236
237 // check the type of operator
238 RCP<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > rbA =
239 Teuchos::rcp_dynamic_cast<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA);
240
241 if(rbA.is_null() == false) {
242 // A is a ReorderedBlockedCrsMatrix and thus uses nested maps, retrieve BlockedCrsMatrix and use plain blocked
243 // map for the construction of the blocked multivectors
244
245 // check type of X vector
246 if (bX.is_null() == true) {
247 RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedDomainMap(), rcpX));
248 rcpX.swap(vectorWithBlockedMap);
249 bCopyResultX = true;
250 bReorderX = true;
251 }
252
253 // check type of B vector
254 if (bB.is_null() == true) {
255 RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedRangeMap(), rcpB));
256 rcpB.swap(vectorWithBlockedMap);
257 bReorderB = true;
258 }
259 }
260 else {
261 // A is a BlockedCrsMatrix and uses a plain blocked map
262 if (bX.is_null() == true) {
263 RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(), rcpX));
264 rcpX.swap(vectorWithBlockedMap);
265 bCopyResultX = true;
266 }
267
268 if(bB.is_null() == true) {
269 RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(),rcpB));
270 rcpB.swap(vectorWithBlockedMap);
271 }
272 }
273
274 // we now can guarantee that X and B are blocked multi vectors
275 bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
276 bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
277
278 // Finally we can do a reordering of the blocked multivectors if A is a ReorderedBlockedCrsMatrix
279 if(rbA.is_null() == false) {
280
281 Teuchos::RCP<const Xpetra::BlockReorderManager > brm = rbA->getBlockReorderManager();
282
283 // check type of X vector
284 if(bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
285 // X is a blocked multi vector but incompatible to the reordered blocked operator A
286 Teuchos::RCP<MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bX);
287 rcpX.swap(reorderedBlockedVector);
288 }
289
290 if(bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
291 // B is a blocked multi vector but incompatible to the reordered blocked operator A
292 Teuchos::RCP<const MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bB);
293 rcpB.swap(reorderedBlockedVector);
294 }
295 }
296
297 // Throughout the rest of the algorithm rcpX and rcpB are used for solution vector and RHS
298
299 // create residual vector
300 // contains current residual of current solution X with rhs B
301 RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
302 RCP<BlockedMultiVector> bresidual = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(residual);
303 Teuchos::RCP<MultiVector> r1 = bresidual->getMultiVector(0,bRangeThyraMode);
304 Teuchos::RCP<MultiVector> r2 = bresidual->getMultiVector(1,bRangeThyraMode);
305
306 // helper vector 1
307 RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
308 RCP<BlockedMultiVector> bxtilde = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xtilde);
309 RCP<MultiVector> xtilde1 = bxtilde->getMultiVector(0,bDomainThyraMode);
310 RCP<MultiVector> xtilde2 = bxtilde->getMultiVector(1,bDomainThyraMode);
311
312 // incrementally improve solution vector X
313 for (LocalOrdinal run = 0; run < nSweeps; ++run) {
314 // 1) calculate current residual
315 residual->update(one,*rcpB,zero); // residual = B
316 A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
317
318 // split residual vector
319 //Teuchos::RCP<MultiVector> r1 = rangeMapExtractor_->ExtractVector(residual, 0, bRangeThyraMode);
320 //Teuchos::RCP<MultiVector> r2 = rangeMapExtractor_->ExtractVector(residual, 1, bRangeThyraMode);
321
322 // 2) solve F * \Delta \tilde{x}_1 = r_1
323 // start with zero guess \Delta \tilde{x}_1
324 //RCP<MultiVector> xtilde1 = MultiVectorFactory::Build(r1->getMap(),rcpX->getNumVectors(),true);
325 //RCP<MultiVector> xtilde2 = MultiVectorFactory::Build(r2->getMap(),rcpX->getNumVectors(),true);
326 bxtilde->putScalar(zero);
327 velPredictSmoo_->Apply(*xtilde1,*r1);
328
329 // 3) solve SchurComp equation
330 // start with zero guess \Delta \tilde{x}_2
331 schurCompSmoo_->Apply(*xtilde2,*r2);
332#if 1
333 // 4) update solution vector
334 rcpX->update(omega,*bxtilde,one);
335#else
336 Teuchos::RCP<MultiVector> x1 = domainMapExtractor_->ExtractVector(rcpX, 0, bDomainThyraMode);
337 Teuchos::RCP<MultiVector> x2 = domainMapExtractor_->ExtractVector(rcpX, 1, bDomainThyraMode);
338
339 // 5) update solution vector with increments xhat1 and xhat2
340 // rescale increment for x2 with omega_
341 x1->update(omega,*xtilde1,one); // x1 = x1_old + omega xtilde1
342 x2->update(omega,*xtilde2,one); // x2 = x2_old + omega xtilde2
343
344 // write back solution in global vector X
345 domainMapExtractor_->InsertVector(x1, 0, rcpX, bDomainThyraMode);
346 domainMapExtractor_->InsertVector(x2, 1, rcpX, bDomainThyraMode);
347#endif
348 }
349
350 if (bCopyResultX == true) {
351 RCP<MultiVector> Xmerged = bX->Merge();
352 X.update(one, *Xmerged, zero);
353 }
354 }
355
356 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
357 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
361
362 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
364 std::ostringstream out;
366 out << "{type = " << type_ << "}";
367 return out.str();
368 }
369
370 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
373
374 if (verbLevel & Parameters0) {
375 out0 << "Prec. type: " << type_ << /*" Sweeps: " << nSweeps_ << " damping: " << omega_ <<*/ std::endl;
376 }
377
378 if (verbLevel & Debug) {
379 out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
380 }
381 }
382 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
384 // FIXME: This is a placeholder
385 return Teuchos::OrdinalTraits<size_t>::invalid();
386 }
387
388
389
390} // namespace MueLu
391
392#endif /* MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_DEF_HPP_ */
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
#define MUELU_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
virtual std::string description() const
Return a simple one-line description of this object.
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
T Get(Level &level, const std::string &varName) const
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory().
RCP< Matrix > Z_
pressure stabilization term or null block
std::string description() const
Return a simple one-line description of this object.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos=0)
Add a factory manager for BraessSarazin internal SchurComplement handling.
RCP< const MapExtractorClass > domainMapExtractor_
domain map extractor (from A_ generated by AFact)
RCP< const MapExtractorClass > rangeMapExtractor_
range map extractor (from A_ generated by AFact)
Teuchos::RCP< SmootherBase > schurCompSmoo_
smoother for SchurComplement equation
Teuchos::RCP< SmootherBase > velPredictSmoo_
Block smoothers.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
RCP< const ParameterList > GetValidParameterList() const
Input.
std::vector< Teuchos::RCP< const FactoryManagerBase > > FactManager_
vector of factory managers
Class that holds all level-specific information.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput().
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
virtual const Teuchos::ParameterList & GetParameterList() const =0
An exception safe way to call the method 'Level::SetFactoryManager()'.
bool IsSetup() const
Get the state of a smoother prototype.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line).
@ Debug
Print additional debugging information.
@ Parameters0
Print class parameters.