MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_CoalesceDropFactory_decl.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_COALESCEDROPFACTORY_DECL_HPP
11#define MUELU_COALESCEDROPFACTORY_DECL_HPP
12
13#include <Xpetra_Matrix_fwd.hpp>
14#include <Xpetra_MultiVector_fwd.hpp>
15#include <Xpetra_VectorFactory_fwd.hpp>
16#include <Xpetra_Vector_fwd.hpp>
17#include <Xpetra_ImportFactory_fwd.hpp>
18#include <Xpetra_CrsGraph_fwd.hpp>
19#include <Xpetra_CrsGraphFactory_fwd.hpp> //TODO
20#include <Xpetra_StridedMap_fwd.hpp>
21#include <Xpetra_Map_fwd.hpp>
22
23#include "MueLu_ConfigDefs.hpp"
27
28#include "MueLu_Level_fwd.hpp"
29#include "MueLu_LWGraph.hpp"
30
31#include "MueLu_LWGraph_fwd.hpp"
36
37namespace MueLu {
38
94template <class Scalar = DefaultScalar,
97 class Node = DefaultNode>
99#undef MUELU_COALESCEDROPFACTORY_SHORT
101
102 public:
104
105
108
111
112 RCP<const ParameterList> GetValidParameterList() const;
113
115
117
118
119 void DeclareInput(Level& currentLevel) const;
120
123
125
126 void Build(Level& currentLevel) const; // Build
127
128 private:
129 // pre-drop function
130 mutable RCP<PreDropFunctionBaseClass> predrop_;
131
133 void MergeRows(const Matrix& A, const LO row, Array<LO>& cols, const Array<LO>& translation) const;
134 void MergeRowsWithDropping(const Matrix& A, const LO row, const ArrayRCP<const SC>& ghostedDiagVals, SC threshold, Array<LO>& cols, const Array<LO>& translation) const;
135
136 // When we want to decouple a block diagonal system (returns Teuchos::null if generate_matrix is false)
137 Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > BlockDiagonalize(Level& currentLevel, const RCP<Matrix>& A, bool generate_matrix) const;
138
139 // When we want to decouple a block diagonal system via a *graph*
140 void BlockDiagonalizeGraph(const RCP<LWGraph>& inputGraph, const RCP<LocalOrdinalVector>& ghostedBlockNumber, RCP<LWGraph>& outputGraph, RCP<const Import>& importer) const;
141
142}; // class CoalesceDropFactory
143
144} // namespace MueLu
145
146#define MUELU_COALESCEDROPFACTORY_SHORT
147#endif // MUELU_COALESCEDROPFACTORY_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Factory for creating a graph based on a given matrix.
void MergeRows(const Matrix &A, const LO row, Array< LO > &cols, const Array< LO > &translation) const
Method to merge rows of matrix for systems of PDEs.
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BlockDiagonalize(Level &currentLevel, const RCP< Matrix > &A, bool generate_matrix) const
void DeclareInput(Level &currentLevel) const
Input.
void BlockDiagonalizeGraph(const RCP< LWGraph > &inputGraph, const RCP< LocalOrdinalVector > &ghostedBlockNumber, RCP< LWGraph > &outputGraph, RCP< const Import > &importer) const
void MergeRowsWithDropping(const Matrix &A, const LO row, const ArrayRCP< const SC > &ghostedDiagVals, SC threshold, Array< LO > &cols, const Array< LO > &translation) const
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.
void SetPreDropFunction(const RCP< MueLu::PreDropFunctionBaseClass< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &predrop)
set predrop function
RCP< PreDropFunctionBaseClass > predrop_
Class that holds all level-specific information.
Base class for factories that use one level (currentLevel).
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar