MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_RegionRFactory_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_REGIONRFACTORY_KOKKOS_DECL_HPP
11#define MUELU_REGIONRFACTORY_KOKKOS_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14#include "MueLu_Types.hpp"
15
16#include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
17#include <Kokkos_Core.hpp>
18
21
22namespace MueLu {
23
29template <class Scalar = DefaultScalar,
32 class Node = DefaultNode>
34#undef MUELU_REGIONRFACTORY_KOKKOS_SHORT
36
37 public:
38 using real_type = typename Teuchos::ScalarTraits<SC>::coordinateType;
39 using realvaluedmultivector_type = typename Xpetra::MultiVector<real_type, LO, GO, NO>;
40 using execution_space = typename Node::execution_space;
41 using memory_space = typename Node::memory_space;
42 using device_type = Kokkos::Device<execution_space, memory_space>;
43 using intTupleView = typename Kokkos::View<int[3], device_type>;
44 using LOTupleView = typename Kokkos::View<LO[3], device_type>;
45
47
48
51
53 virtual ~RegionRFactory_kokkos() = default;
55
57
58 RCP<const ParameterList> GetValidParameterList() const;
59
60 void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
61
63
65
66
67 void Build(Level& fineLevel, Level& coarseLevel) const;
68
69 void Build3D(const int numDimensions,
70 Array<LO>& lFineNodesPerDim,
71 const RCP<Matrix>& A,
72 const RCP<realvaluedmultivector_type>& fineCoordinates,
73 RCP<Matrix>& R,
74 RCP<realvaluedmultivector_type>& coarseCoordinates,
75 Array<LO>& lCoarseNodesPerDim) const;
76
78
79}; // class RegionRFactory_kokkos
80
81} // namespace MueLu
82
83#define MUELU_REGIONRFACTORY_KOKKOS_SHORT
84#endif // MUELU_REGIONRFACTORY_KOKKOS_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Class that holds all level-specific information.
Factory that builds a restriction operator for region multigrid.
typename Kokkos::View< LO[3], device_type > LOTupleView
RegionRFactory_kokkos()=default
Default Constructor.
typename Teuchos::ScalarTraits< SC >::coordinateType real_type
typename Xpetra::MultiVector< real_type, LO, GO, NO > realvaluedmultivector_type
typename Kokkos::View< int[3], device_type > intTupleView
Kokkos::Device< execution_space, memory_space > device_type
void Build3D(const int numDimensions, Array< LO > &lFineNodesPerDim, const RCP< Matrix > &A, const RCP< realvaluedmultivector_type > &fineCoordinates, RCP< Matrix > &R, RCP< realvaluedmultivector_type > &coarseCoordinates, Array< LO > &lCoarseNodesPerDim) const
typename Node::execution_space execution_space
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
virtual ~RegionRFactory_kokkos()=default
Destructor.
Base class for factories that use two levels (fineLevel and coarseLevel).
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar