MueLu  Version of the Day
MueLu_AlgebraicPermutationStrategy_decl.hpp
Go to the documentation of this file.
1 /*
2  * MueLu_AlgebraicPermutationStrategy_decl.hpp
3  *
4  * Created on: Feb 20, 2013
5  * Author: tobias
6  */
7 
8 #ifndef MUELU_ALGEBRAICPERMUTATIONSTRATEGY_DECL_HPP_
9 #define MUELU_ALGEBRAICPERMUTATIONSTRATEGY_DECL_HPP_
10 
11 #include <Xpetra_MultiVector_fwd.hpp>
12 #include <Xpetra_Matrix_fwd.hpp>
13 #include <Xpetra_CrsGraph_fwd.hpp>
14 #include <Xpetra_Vector_fwd.hpp>
15 #include <Xpetra_MultiVectorFactory_fwd.hpp>
16 #include <Xpetra_VectorFactory_fwd.hpp>
17 #include <Xpetra_CrsMatrixWrap_fwd.hpp>
18 #include <Xpetra_Export_fwd.hpp>
19 #include <Xpetra_ExportFactory_fwd.hpp>
20 #include <Xpetra_Import_fwd.hpp>
21 #include <Xpetra_ImportFactory_fwd.hpp>
22 
23 #include "MueLu_ConfigDefs.hpp"
24 #include "MueLu_Level.hpp"
25 #include "MueLu_BaseClass.hpp"
26 
27 namespace MueLu {
28 
29 // template struct for comparing pairs
30 template<class Scalar = double, class LocalOrdinal = int>
31 struct CompPairs {
32  CompPairs(const std::vector<Scalar> & v) : vinternal_(v) {}
33  std::vector<Scalar> vinternal_;
34  bool operator()(LocalOrdinal a, LocalOrdinal b) {
35  //return vinternal_[a] < vinternal_[b];
36  return Teuchos::ScalarTraits<Scalar>::magnitude(vinternal_[a]) > Teuchos::ScalarTraits<Scalar>::magnitude(vinternal_[b]);
37  }
38 };
39 
40 // template function for comparison
41 template<class Scalar, class LocalOrdinal>
42 CompPairs<Scalar,LocalOrdinal> CreateCmpPairs(const std::vector<Scalar> & v) {
44 }
45 
46 // template function for sorting permutations
47 template<class Scalar, class LocalOrdinal>
48 void sortingPermutation(const std::vector<Scalar> & values, std::vector<LocalOrdinal> & v) {
49  size_t size = values.size();
50  v.clear(); v.reserve(size);
51  for(size_t i=0; i<size; ++i)
52  v.push_back(i);
53 
54  std::sort(v.begin(),v.end(), MueLu::CreateCmpPairs<Scalar,LocalOrdinal>(values));
55 }
56 
58 
64  template<class Scalar = double,
65  class LocalOrdinal = int,
66  class GlobalOrdinal = LocalOrdinal,
67  class Node = KokkosClassic::DefaultNode::DefaultNodeType>
68  class AlgebraicPermutationStrategy : public BaseClass {
69 #undef MUELU_ALGEBRAICPERMUTATIONSTRATEGY_SHORT
70 #include "MueLu_UseShortNames.hpp"
71  public:
72 
78 
80 
82 
95  void BuildPermutation(const Teuchos::RCP<Matrix> & A, const Teuchos::RCP<const Map> permRowMap, Level & currentLevel, const FactoryBase* genFactory) const;
96 
98 
99 
100  private:
101  };
102 
103 } // namespace MueLu
104 
105 #define MUELU_ALGEBRAICPERMUTATIONSTRATEGY_SHORT
106 
107 
108 #endif /* MUELU_ALGEBRAICPERMUTATIONSTRATEGY_DECL_HPP_ */
void sortingPermutation(const std::vector< Scalar > &values, std::vector< LocalOrdinal > &v)
void BuildPermutation(const Teuchos::RCP< Matrix > &A, const Teuchos::RCP< const Map > permRowMap, Level &currentLevel, const FactoryBase *genFactory) const
build permutation operators
Namespace for MueLu classes and methods.
bool operator()(LocalOrdinal a, LocalOrdinal b)
CompPairs(const std::vector< Scalar > &v)
Class which defines local permutations of matrix columns.
CompPairs< Scalar, LocalOrdinal > CreateCmpPairs(const std::vector< Scalar > &v)