MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_SaPFactory_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_SAPFACTORY_DECL_HPP
11#define MUELU_SAPFACTORY_DECL_HPP
12
13#include <string>
14
15#include "MueLu_ConfigDefs.hpp"
16
17#include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
18
20
21#include "MueLu_Level_fwd.hpp"
24#include "MueLu_PFactory.hpp"
26
27namespace MueLu {
28
58template <class Scalar = DefaultScalar,
61 class Node = DefaultNode>
62class SaPFactory : public PFactory {
63 public:
66 typedef typename Node::execution_space execution_space;
67 typedef Node node_type;
68
69 private:
70#undef MUELU_SAPFACTORY_SHORT
72
73 public:
75
76
81
83 virtual ~SaPFactory();
84
85 RCP<const ParameterList> GetValidParameterList() const;
86
88
90
91
92 void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
93
95
97
98
105 void Build(Level& fineLevel, Level& coarseLevel) const;
106
107 void BuildP(Level& fineLevel, Level& coarseLevel) const;
108
109 void SatisfyPConstraintsNonKokkos(RCP<Matrix> A, RCP<Matrix>& P) const;
110
111 void optimalSatisfyPConstraintsForScalarPDEsNonKokkos(RCP<Matrix>& P) const;
112
113 void SatisfyPConstraints(RCP<Matrix> A, RCP<Matrix>& P) const;
114
115 void optimalSatisfyPConstraintsForScalarPDEs(RCP<Matrix>& P) const;
116
117 bool constrainRow(Scalar* orig, LocalOrdinal nEntries, Scalar leftBound, Scalar rghtBound, Scalar rsumTarget, Scalar* fixedUnsorted, Scalar* scalarData) const;
118
120};
121
122} // namespace MueLu
123
124#define MUELU_SAPFACTORY_SHORT
125#endif // MUELU_SAPFACTORY_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Class that holds all level-specific information.
Factory that provides an interface for a concrete implementation of a prolongation operator.
Factory for building Smoothed Aggregation prolongators.
void SatisfyPConstraints(RCP< Matrix > A, RCP< Matrix > &P) const
void Build(Level &fineLevel, Level &coarseLevel) const
virtual ~SaPFactory()
Destructor.
void optimalSatisfyPConstraintsForScalarPDEsNonKokkos(RCP< Matrix > &P) const
GlobalOrdinal global_ordinal_type
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void optimalSatisfyPConstraintsForScalarPDEs(RCP< Matrix > &P) const
bool constrainRow(Scalar *orig, LocalOrdinal nEntries, Scalar leftBound, Scalar rghtBound, Scalar rsumTarget, Scalar *fixedUnsorted, Scalar *scalarData) const
Node::execution_space execution_space
SaPFactory()
Constructor. User can supply a factory for generating the tentative prolongator.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void SatisfyPConstraintsNonKokkos(RCP< Matrix > A, RCP< Matrix > &P) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar