MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_MergedBlockedMatrixFactory_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_MERGEDBLOCKEDMATRIXFACTORY_DEF_HPP_
11#define MUELU_MERGEDBLOCKEDMATRIXFACTORY_DEF_HPP_
12
13#include <Xpetra_BlockedCrsMatrix.hpp>
14#include "MueLu_Level.hpp"
15#include "MueLu_Monitor.hpp"
16#include "MueLu_PerfUtils.hpp"
17
19
20namespace MueLu {
21
22template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
24
25template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
27 RCP<ParameterList> validParamList = rcp(new ParameterList());
28
29 validParamList->set<RCP<const FactoryBase> >("A", MueLu::NoFactory::getRCP() /*Teuchos::null*/, "Generating factory of the matrix A used for building SchurComplement (must be a 2x2 blocked operator, default = MueLu::NoFactory::getRCP())");
30
31 return validParamList;
32}
33
34template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36 Input(currentLevel, "A");
37}
38
39template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
41 FactoryMonitor m(*this, "MergedBlockedMatrix", currentLevel);
42 Teuchos::RCP<Matrix> A = Get<RCP<Matrix> >(currentLevel, "A");
43
44 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A);
45 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::MergedBlockedMatrixFactory::Build: input matrix A is not of type BlockedCrsMatrix! A generated by AFact_ must be a 2x2 block operator. error.");
46
47 Teuchos::RCP<Matrix> mergedA = bA->Merge();
48 {
49 GetOStream(Statistics1) << PerfUtils::PrintMatrixInfo(*mergedA, "A (merged)");
50
51 // note: variable "A" generated by this MergedBlockedMatrix factory is in fact the a merged version
52 // of the blocked matrix A (from the input)
53 Set(currentLevel, "A", mergedA);
54 }
55}
56
57} // namespace MueLu
58
59#endif /* MUELU_MERGEDBLOCKEDMATRIXFACTORY_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 DeclareInput(Level &currentLevel) const
Specifies the data that this class needs, and the factories that generate that data.
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 const RCP< const NoFactory > getRCP()
Static Get() functions.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
Namespace for MueLu classes and methods.
@ Statistics1
Print more statistics.