MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_TentativePFactory_kokkos_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_KOKKOS_DECL_HPP
11 #define MUELU_TENTATIVEPFACTORY_KOKKOS_DECL_HPP
12 
13 #include "MueLu_ConfigDefs.hpp"
14 
16 
17 #include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
18 
19 #include "Teuchos_ScalarTraits.hpp"
20 
22 
23 #include "MueLu_Aggregates_fwd.hpp"
25 #include "MueLu_Level_fwd.hpp"
26 #include "MueLu_PerfUtils_fwd.hpp"
27 #include "MueLu_PFactory.hpp"
28 
29 namespace MueLu {
30 
69 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71  public:
74  typedef typename Node::execution_space execution_space;
75  typedef Kokkos::RangePolicy<local_ordinal_type, execution_space> range_type;
76  typedef typename Node::device_type DeviceType;
77  typedef Node node_type;
80 
81  private:
82 #undef MUELU_TENTATIVEPFACTORY_KOKKOS_SHORT
83 #include "MueLu_UseShortNames.hpp"
84 
85  public:
87 
88 
91 
95 
97 
99 
100 
101  void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
102 
104 
106 
107 
108  void Build(Level& fineLevel, Level& coarseLevel) const;
109  void BuildP(Level& fineLevel, Level& coarseLevel) const;
110 
112 
113  // NOTE: All of thess should really be private, but CUDA doesn't like that
114 
115  void BuildPuncoupledBlockCrs(Level& coarseLevel, RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo,
116  RCP<MultiVector> fineNullspace, RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const;
117 
118  bool isGoodMap(const Map& rowMap, const Map& colMap) const;
119 
120  void BuildPcoupled(RCP<Matrix> A, RCP<Aggregates> aggregates,
121  RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
122  RCP<const Map> coarseMap, RCP<Matrix>& Ptentative,
123  RCP<MultiVector>& coarseNullspace) const;
124 
125  void BuildPuncoupled(Level& coarseLevel, RCP<Matrix> A, RCP<Aggregates> aggregates,
126  RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
127  RCP<const Map> coarseMap, RCP<Matrix>& Ptentative,
128  RCP<MultiVector>& coarseNullspace, const int levelID) const;
129 
130  mutable bool bTransferCoordinates_ = false;
131 };
132 
133 } // namespace MueLu
134 
135 #define MUELU_TENTATIVEPFACTORY_KOKKOS_SHORT
136 #endif // MUELU_TENTATIVEPFACTORY_KOKKOS_DECL_HPP
Teuchos::ScalarTraits< Scalar >::coordinateType real_type
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
MueLu::DefaultLocalOrdinal LocalOrdinal
void BuildPuncoupledBlockCrs(Level &coarseLevel, 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.
Xpetra::MultiVector< real_type, LocalOrdinal, GlobalOrdinal, node_type > RealValuedMultiVector
MueLu::DefaultNode Node
void BuildPuncoupled(Level &coarseLevel, 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
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Kokkos::RangePolicy< local_ordinal_type, execution_space > range_type
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
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.
Factory that provides an interface for a concrete implementation of a prolongation operator...
bool isGoodMap(const Map &rowMap, const Map &colMap) const