Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_CreateOverlapGraph.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4 //
5 // Copyright 2009 NTESS and the Ifpack2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef IFPACK2_CREATEOVERLAPGRAPH_HPP
11 #define IFPACK2_CREATEOVERLAPGRAPH_HPP
12 
13 #include "Ifpack2_ConfigDefs.hpp"
14 #include "Tpetra_CrsGraph.hpp"
15 #include "Tpetra_CrsMatrix.hpp"
16 #include "Tpetra_Import.hpp"
17 #include "Teuchos_RefCountPtr.hpp"
18 
19 
20 namespace Ifpack2 {
21 
38 template<class GraphType>
41  const int overlapLevel)
42 {
43  using Teuchos::RCP;
44  using Teuchos::rcp;
45  typedef Tpetra::Map<typename GraphType::local_ordinal_type,
46  typename GraphType::global_ordinal_type,
47  typename GraphType::node_type> map_type;
48  typedef Tpetra::Import<typename GraphType::local_ordinal_type,
49  typename GraphType::global_ordinal_type,
50  typename GraphType::node_type> import_type;
52  overlapLevel < 0, std::invalid_argument,
53  "Ifpack2::createOverlapGraph: overlapLevel must be >= 0, "
54  "but you specified overlapLevel = " << overlapLevel << ".");
55 
56  const int numProcs = inputGraph->getMap ()->getComm ()->getSize ();
57  if (overlapLevel == 0 || numProcs < 2) {
58  return inputGraph;
59  }
60 
61  RCP<const map_type> overlapRowMap = inputGraph->getRowMap ();
62  RCP<const map_type> domainMap = inputGraph->getDomainMap ();
63  RCP<const map_type> rangeMap = inputGraph->getRangeMap ();
64 
65  RCP<GraphType> overlapGraph;
66  RCP<const GraphType> oldGraph;
67  RCP<const map_type> oldRowMap;
68  for (int level = 0; level < overlapLevel; ++level) {
69  oldGraph = overlapGraph;
70  oldRowMap = overlapRowMap;
71 
72  RCP<const import_type> overlapImporter;
73  if (level == 0) {
74  overlapImporter = inputGraph->getImporter ();
75  } else {
76  overlapImporter = oldGraph->getImporter ();
77  }
78 
79  overlapRowMap = overlapImporter->getTargetMap ();
80  if (level < overlapLevel - 1) {
81  overlapGraph = rcp (new GraphType (overlapRowMap, 0));
82  }
83  else {
84  // On last iteration, we want to filter out all columns except those that
85  // correspond to rows in the graph. This ensures that our graph is square
86  overlapGraph = rcp (new GraphType (overlapRowMap, overlapRowMap, 0));
87  }
88 
89  overlapGraph->doImport (*inputGraph, *overlapImporter, Tpetra::INSERT);
90  overlapGraph->fillComplete (domainMap, rangeMap);
91  }
92 
93  return overlapGraph;
94 }
95 
105 template<class MatrixType>
108  const int overlapLevel)
109 {
110  using Teuchos::RCP;
111  using Teuchos::rcp;
112  typedef typename MatrixType::map_type map_type;
113  typedef Tpetra::Import<typename MatrixType::local_ordinal_type,
114  typename MatrixType::global_ordinal_type,
115  typename MatrixType::node_type> import_type;
116 
118  overlapLevel < 0, std::invalid_argument,
119  "Ifpack2::createOverlapMatrix: overlapLevel must be >= 0, "
120  "but you specified overlapLevel = " << overlapLevel << ".");
121 
122  const int numProcs = inputMatrix->getMap ()->getComm ()->getSize ();
123  if (overlapLevel == 0 || numProcs < 2) {
124  return inputMatrix;
125  }
126 
127  RCP<const map_type> overlapRowMap = inputMatrix->getRowMap ();
128  RCP<const map_type> domainMap = inputMatrix->getDomainMap ();
129  RCP<const map_type> rangeMap = inputMatrix->getRangeMap ();
130 
131  RCP<MatrixType> overlapMatrix;
132  RCP<const MatrixType> oldMatrix;
133  RCP<const map_type> oldRowMap;
134  for (int level = 0; level < overlapLevel; ++level) {
135  oldMatrix = overlapMatrix;
136  oldRowMap = overlapRowMap;
137 
138  RCP<const import_type> overlapImporter;
139  if (level == 0) {
140  overlapImporter = inputMatrix->getGraph ()->getImporter ();
141  } else {
142  overlapImporter = oldMatrix->getGraph ()->getImporter ();
143  }
144 
145  overlapRowMap = overlapImporter->getTargetMap ();
146  if (level < overlapLevel - 1) {
147  overlapMatrix = rcp (new MatrixType (overlapRowMap, 0));
148  }
149  else {
150  // On last iteration, we want to filter out all columns except those that
151  // correspond to rows in the matrix. This assures that our matrix is square
152  overlapMatrix = rcp (new MatrixType (overlapRowMap, overlapRowMap, 0));
153  }
154 
155  overlapMatrix->doImport (*inputMatrix, *overlapImporter, Tpetra::INSERT);
156  overlapMatrix->fillComplete (domainMap, rangeMap);
157  }
158 
159  return overlapMatrix;
160 }
161 
162 } // namespace Ifpack2
163 
164 #endif // IFPACK2_CREATEOVERLAPGRAPH_HPP
Teuchos::RCP< const GraphType > createOverlapGraph(const Teuchos::RCP< const GraphType > &inputGraph, const int overlapLevel)
Construct an overlapped graph for use with Ifpack2 preconditioners.
Definition: Ifpack2_CreateOverlapGraph.hpp:40
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const MatrixType > createOverlapMatrix(const Teuchos::RCP< const MatrixType > &inputMatrix, const int overlapLevel)
Construct an overlapped matrix for use with Ifpack2 preconditioners.
Definition: Ifpack2_CreateOverlapGraph.hpp:107