MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_SemiCoarsenPFactory_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_SEMICOARSENPFACTORY_DECL_HPP
11#define MUELU_SEMICOARSENPFACTORY_DECL_HPP
12
13#include <Xpetra_Matrix_fwd.hpp>
14
15#include "MueLu_ConfigDefs.hpp"
16#include "MueLu_PFactory.hpp"
18
19#include "MueLu_Level_fwd.hpp"
20
21#define VERTICAL 1
22#define HORIZONTAL 2
23#define GRID_SUPPLIED -1
24#define NUM_ZPTS 0
25#define ORIENTATION 1
26
27namespace MueLu {
28
70template <class Scalar = DefaultScalar,
73 class Node = DefaultNode>
75#undef MUELU_SEMICOARSENPFACTORY_SHORT
77
78 public:
80
81
85
89
90 RCP<const ParameterList> GetValidParameterList() const;
91
93
94
95 void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
96
98
100
101
102 void Build(Level& fineLevel, Level& coarseLevel) const;
103 void BuildP(Level& fineLevel, Level& coarseLevel) const;
104
106
107 private:
108 LO FindCpts(LO const PtsPerLine, LO const CoarsenRate, LO const Thin, LO** LayerCpts) const;
109 LO MakeSemiCoarsenP(LO const Ntotal, LO const nz, LO const CoarsenRate, LO const LayerId[],
110 LO const VertLineId[], LO const DofsPerNode, RCP<Matrix>& Amat,
111 RCP<Matrix>& P, RCP<const Map>& coarseMap,
112 const RCP<MultiVector> fineNullspace, RCP<MultiVector>& coarseNullspace, RCP<Matrix>& R, bool buildRestriction, bool doLinear) const;
113 void RevertToPieceWiseConstant(RCP<Matrix>& P, LO BlkSize) const;
114
115 mutable bool bTransferCoordinates_; //< boolean which is true if coordinate information is available to be transferred to coarse coordinate information
116}; // class SemiCoarsenPFactory
117
118} // namespace MueLu
119
120#define MUELU_SEMICOARSENPFACTORY_SHORT
121#endif // MUELU_SEMICOARSENPFACTORY_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.
Prolongator factory performing semi-coarsening.
LO FindCpts(LO const PtsPerLine, LO const CoarsenRate, LO const Thin, LO **LayerCpts) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
LO MakeSemiCoarsenP(LO const Ntotal, LO const nz, LO const CoarsenRate, LO const LayerId[], LO const VertLineId[], LO const DofsPerNode, RCP< Matrix > &Amat, RCP< Matrix > &P, RCP< const Map > &coarseMap, const RCP< MultiVector > fineNullspace, RCP< MultiVector > &coarseNullspace, RCP< Matrix > &R, bool buildRestriction, bool doLinear) const
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void RevertToPieceWiseConstant(RCP< Matrix > &P, LO BlkSize) const
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar