MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_RepartitionBlockDiagonalFactory_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_REPARTITIONBLOCKDIAGONALFACTORY_DEF_HPP_
11#define MUELU_REPARTITIONBLOCKDIAGONALFACTORY_DEF_HPP_
12
14
15#include <Teuchos_Utils.hpp>
16
17#include <Xpetra_BlockedCrsMatrix.hpp>
18
19#include "MueLu_Exceptions.hpp"
20#include "MueLu_Level.hpp"
21#include "MueLu_Monitor.hpp"
22
23namespace MueLu {
24
25template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
27 RCP<ParameterList> validParamList = rcp(new ParameterList());
28 validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Factory of the matrix A");
29
30 return validParamList;
31}
32
33template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35 Input(currentLevel, "A");
36} // DeclareInput()
37
38template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40 FactoryMonitor m(*this, "Build", currentLevel);
41 typedef Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> BlockCrs;
42
43 RCP<Matrix> originalA = Get<RCP<Matrix> >(currentLevel, "A");
44 RCP<BlockCrs> A = Teuchos::rcp_dynamic_cast<BlockCrs>(originalA);
45
46 // RCP<BlockCrs> A = Get< RCP<BlockCrs> >(currentLevel, "A");
47 TEUCHOS_TEST_FOR_EXCEPTION(A == Teuchos::null, Exceptions::BadCast, "MueLu::RepartitionBlockDiagonalFactory::Build: input matrix A is not of type BlockedCrsMatrix! error.");
48
49 // Build the block diagonal
50 RCP<BlockCrs> DiagonalMatrix = Teuchos::rcp(new BlockCrs(A->getBlockedRangeMap(), A->getBlockedDomainMap(), 0));
51 for (size_t i = 0; i < A->Rows(); i++)
52 DiagonalMatrix->setMatrix(i, i, A->getMatrix(i, i));
53
54 Set(currentLevel, "A", Teuchos::rcp_dynamic_cast<Matrix>(DiagonalMatrix));
55
56} // Build()
57
58} // end namespace MueLu
59
60#endif /* MUELU_REPARTITIONBLOCKDIAGONALFACTORY_DEF_HPP_ */
Exception indicating invalid cast attempted.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
void Build(Level &level) const
Build an object with this factory.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &level) const
Specifies the data that this class needs, and the factories that generate that data.
Namespace for MueLu classes and methods.