Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
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
19namespace Ifpack2 {
20
37template <class GraphType>
38Teuchos::RCP<const GraphType>
39createOverlapGraph(const Teuchos::RCP<const GraphType>& inputGraph,
40 const int overlapLevel) {
41 using Teuchos::RCP;
42 using Teuchos::rcp;
43 typedef Tpetra::Map<typename GraphType::local_ordinal_type,
44 typename GraphType::global_ordinal_type,
45 typename GraphType::node_type>
46 map_type;
47 typedef Tpetra::Import<typename GraphType::local_ordinal_type,
48 typename GraphType::global_ordinal_type,
49 typename GraphType::node_type>
50 import_type;
51 TEUCHOS_TEST_FOR_EXCEPTION(
52 overlapLevel < 0, std::invalid_argument,
53 "Ifpack2::createOverlapGraph: overlapLevel must be >= 0, "
54 "but you specified overlapLevel = "
55 << overlapLevel << ".");
56
57 const int numProcs = inputGraph->getMap()->getComm()->getSize();
58 if (overlapLevel == 0 || numProcs < 2) {
59 return inputGraph;
60 }
61
62 RCP<const map_type> overlapRowMap = inputGraph->getRowMap();
63 RCP<const map_type> domainMap = inputGraph->getDomainMap();
64 RCP<const map_type> rangeMap = inputGraph->getRangeMap();
65
66 RCP<GraphType> overlapGraph;
67 RCP<const GraphType> oldGraph;
68 RCP<const map_type> oldRowMap;
69 for (int level = 0; level < overlapLevel; ++level) {
70 oldGraph = overlapGraph;
71 oldRowMap = overlapRowMap;
72
73 RCP<const import_type> overlapImporter;
74 if (level == 0) {
75 overlapImporter = inputGraph->getImporter();
76 } else {
77 overlapImporter = oldGraph->getImporter();
78 }
79
80 overlapRowMap = overlapImporter->getTargetMap();
81 if (level < overlapLevel - 1) {
82 overlapGraph = rcp(new GraphType(overlapRowMap, 0));
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
105template <class MatrixType>
106Teuchos::RCP<const MatrixType>
107createOverlapMatrix(const Teuchos::RCP<const MatrixType>& inputMatrix,
108 const int overlapLevel) {
109 using Teuchos::RCP;
110 using Teuchos::rcp;
111 typedef typename MatrixType::map_type map_type;
112 typedef Tpetra::Import<typename MatrixType::local_ordinal_type,
113 typename MatrixType::global_ordinal_type,
114 typename MatrixType::node_type>
115 import_type;
116
117 TEUCHOS_TEST_FOR_EXCEPTION(
118 overlapLevel < 0, std::invalid_argument,
119 "Ifpack2::createOverlapMatrix: overlapLevel must be >= 0, "
120 "but you specified overlapLevel = "
121 << overlapLevel << ".");
122
123 const int numProcs = inputMatrix->getMap()->getComm()->getSize();
124 if (overlapLevel == 0 || numProcs < 2) {
125 return inputMatrix;
126 }
127
128 RCP<const map_type> overlapRowMap = inputMatrix->getRowMap();
129 RCP<const map_type> domainMap = inputMatrix->getDomainMap();
130 RCP<const map_type> rangeMap = inputMatrix->getRangeMap();
131
132 RCP<MatrixType> overlapMatrix;
133 RCP<const MatrixType> oldMatrix;
134 RCP<const map_type> oldRowMap;
135 for (int level = 0; level < overlapLevel; ++level) {
136 oldMatrix = overlapMatrix;
137 oldRowMap = overlapRowMap;
138
139 RCP<const import_type> overlapImporter;
140 if (level == 0) {
141 overlapImporter = inputMatrix->getGraph()->getImporter();
142 } else {
143 overlapImporter = oldMatrix->getGraph()->getImporter();
144 }
145
146 overlapRowMap = overlapImporter->getTargetMap();
147 if (level < overlapLevel - 1) {
148 overlapMatrix = rcp(new MatrixType(overlapRowMap, 0));
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
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
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
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:39