Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 
24 namespace Ifpack2 {
25 
26 namespace Details {
27 
32 template <class graph_type>
34 computeDiagonalGraph(graph_type const& graph) {
35  typedef typename graph_type::local_ordinal_type LO;
36  typedef typename graph_type::global_ordinal_type GO;
37  typedef typename graph_type::node_type NO;
38  typedef Tpetra::Map<LO, GO, NO> map_type;
39  typedef Tpetra::CrsGraph<LO, GO, NO> crs_graph_type;
40 
41  const size_t maxDiagEntPerRow = 1;
42  // NOTE (mfh 12 Aug 2014) We could also pass in the column Map
43  // here. However, we still would have to do LID->GID lookups to
44  // make sure that we are using the correct diagonal column
45  // indices, so it probably wouldn't help much.
46  Teuchos::RCP<graph_type> diagonalGraph;
47  diagonalGraph = Teuchos::rcp(new crs_graph_type(graph.getRowMap(), maxDiagEntPerRow));
48  const map_type& meshRowMap = *(graph.getRowMap());
49 
50  Teuchos::Array<GO> diagGblColInds(maxDiagEntPerRow);
51 
52  for (LO lclRowInd = meshRowMap.getMinLocalIndex(); lclRowInd <= meshRowMap.getMaxLocalIndex(); ++lclRowInd) {
53  const GO gblRowInd = meshRowMap.getGlobalElement(lclRowInd);
54  diagGblColInds[0] = gblRowInd;
55  diagonalGraph->insertGlobalIndices(gblRowInd, diagGblColInds());
56  }
57 
58  diagonalGraph->fillComplete(graph.getDomainMap(), graph.getRangeMap());
59  return diagonalGraph;
60 }
61 
63 std::string canonicalize(const std::string& precType);
64 
65 } // namespace Details
66 
67 } // namespace Ifpack2
68 
69 #endif // IFPACK2_UTILITIES_HPP
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)