Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Utilities.hpp
Go to the documentation of this file.
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_UTILITIES_HPP
11#define IFPACK2_UTILITIES_HPP
12
13#include "Ifpack2_ConfigDefs.hpp"
14
15#include "Teuchos_RefCountPtr.hpp"
16#include "Teuchos_ScalarTraits.hpp"
17
18#include "Tpetra_ConfigDefs.hpp"
19#include "Tpetra_CrsGraph.hpp"
20
23
24namespace Ifpack2 {
25
26namespace Details {
27
28
33 template <class graph_type>
34 Teuchos::RCP<Tpetra::CrsGraph<typename graph_type::local_ordinal_type, typename graph_type::global_ordinal_type, typename graph_type::node_type> >
35 computeDiagonalGraph (graph_type const &graph)
36 {
37 typedef typename graph_type::local_ordinal_type LO;
38 typedef typename graph_type::global_ordinal_type GO;
39 typedef typename graph_type::node_type NO;
40 typedef Tpetra::Map<LO, GO, NO> map_type;
41 typedef Tpetra::CrsGraph<LO, GO, NO> crs_graph_type;
42
43 const size_t maxDiagEntPerRow = 1;
44 // NOTE (mfh 12 Aug 2014) We could also pass in the column Map
45 // here. However, we still would have to do LID->GID lookups to
46 // make sure that we are using the correct diagonal column
47 // indices, so it probably wouldn't help much.
48 Teuchos::RCP<graph_type> diagonalGraph;
49 diagonalGraph = Teuchos::rcp(new crs_graph_type(graph.getRowMap(), maxDiagEntPerRow));
50 const map_type& meshRowMap = *(graph.getRowMap());
51
52 Teuchos::Array<GO> diagGblColInds(maxDiagEntPerRow);
53
54 for (LO lclRowInd = meshRowMap.getMinLocalIndex(); lclRowInd <= meshRowMap.getMaxLocalIndex(); ++lclRowInd) {
55 const GO gblRowInd = meshRowMap.getGlobalElement(lclRowInd);
56 diagGblColInds[0] = gblRowInd;
57 diagonalGraph->insertGlobalIndices(gblRowInd, diagGblColInds());
58 }
59
60 diagonalGraph->fillComplete(graph.getDomainMap(), graph.getRangeMap());
61 return diagonalGraph;
62 }
63
64
66 std::string canonicalize(const std::string& precType);
67
68
69}// namespace Details
70
71}// namespace Ifpack2
72
73#endif //IFPACK2_UTILITIES_HPP
Teuchos::RCP< Tpetra::CrsGraph< typename graph_type::local_ordinal_type, typename graph_type::global_ordinal_type, typename graph_type::node_type > > computeDiagonalGraph(graph_type const &graph)
Compute and return the graph of the diagonal of the input graph.
Definition Ifpack2_Utilities.hpp:35
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41