MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_TentativePFactory_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_TENTATIVEPFACTORY_DECL_HPP
11#define MUELU_TENTATIVEPFACTORY_DECL_HPP
12
13#include <Teuchos_ScalarTraits.hpp>
14#include <Teuchos_SerialDenseMatrix.hpp>
15#include <Teuchos_SerialQRDenseSolver.hpp>
16
17#include <Xpetra_CrsGraphFactory_fwd.hpp>
18#include <Xpetra_CrsMatrix_fwd.hpp>
19#include <Xpetra_Matrix_fwd.hpp>
20#include <Xpetra_MultiVector_fwd.hpp>
21#include <Xpetra_MapFactory_fwd.hpp>
22#include <Xpetra_Map_fwd.hpp>
23#include <Xpetra_MultiVectorFactory_fwd.hpp>
24#include <Xpetra_Import_fwd.hpp>
25#include <Xpetra_ImportFactory_fwd.hpp>
26#include <Xpetra_CrsMatrixWrap_fwd.hpp>
27
28#include "MueLu_ConfigDefs.hpp"
30
35#include "MueLu_Level_fwd.hpp"
37#include "MueLu_PFactory.hpp"
39
40namespace MueLu {
41
81template <class Scalar = DefaultScalar,
84 class Node = DefaultNode>
86#undef MUELU_TENTATIVEPFACTORY_SHORT
88
89 public:
91
92
95
97 virtual ~TentativePFactory() {}
99
100 RCP<const ParameterList> GetValidParameterList() const;
101
103
104
105 void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
106
108
110
111
112 void Build(Level& fineLevel, Level& coarseLevel) const;
113 void BuildP(Level& fineLevel, Level& coarseLevel) const;
114
116
117 private:
118 void BuildPuncoupled(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
119 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const;
120 void BuildPcoupled(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
121 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace) const;
122 void BuildPuncoupledBlockCrs(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
123 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const;
124
125 mutable bool bTransferCoordinates_ = false;
126
127}; // class TentativePFactory
128
129} // namespace MueLu
130
131// TODO: noQR_
132
133#define MUELU_TENTATIVEPFACTORY_SHORT
134#endif // MUELU_TENTATIVEPFACTORY_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 tentative prolongator.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void BuildPuncoupled(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace, const int levelID) const
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void BuildPuncoupledBlockCrs(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace, const int levelID) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void BuildPcoupled(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace) const
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar