Ifpack2 Templated Preconditioning Package  Version 1.0
Ifpack2_Details_CrsArrays.hpp
Go to the documentation of this file.
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
47 
48 #ifndef __IFPACK2_CRSARRAYS_DECL_HPP__
49 #define __IFPACK2_CRSARRAYS_DECL_HPP__
50 
51 #include <Tpetra_RowMatrix.hpp>
52 #include <Tpetra_CrsMatrix.hpp>
53 #include <KokkosCompat_DefaultNode.hpp>
54 #include <KokkosSparse_CrsMatrix.hpp>
55 #include <Ifpack2_LocalFilter.hpp>
56 #include <Ifpack2_ReorderFilter.hpp>
57 
58 namespace Ifpack2
59 {
60 namespace Details
61 {
62 
63 //Utility for getting the local values, rowptrs and colinds (in Kokkos::Views) for any RowMatrix
64 //Used by Fic, Filu and Fildl but may also be useful in other classes
65 template<typename Scalar, typename ImplScalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
66 struct CrsArrayReader
67 {
68  typedef typename Node::device_type device_type;
69  typedef typename device_type::execution_space execution_space;
70  typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TRowMatrix;
71  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TCrsMatrix;
72  typedef Ifpack2::LocalFilter<TRowMatrix> Filter;
73  typedef Ifpack2::ReorderFilter<TRowMatrix> ReordFilter;
74  typedef KokkosSparse::CrsMatrix<ImplScalar, LocalOrdinal, execution_space> KCrsMatrix;
75  typedef Kokkos::View<LocalOrdinal*, execution_space> OrdinalArray;
76  typedef Kokkos::View<ImplScalar*, execution_space> ScalarArray;
77  typedef typename OrdinalArray::HostMirror OrdinalArrayHost;
79  typedef Kokkos::Serial functor_space;
80  typedef Kokkos::RangePolicy<functor_space, int> RangePol;
81 
83  // \param[in] A The matrix
84  // \param[out] vals The values array on device (allocated inside function)
85  // \param[in] rowptrs The rowptrs host view provided by getStructure()
86  static void getValues(const TRowMatrix* A, ScalarArray& vals, OrdinalArrayHost& rowptrs)
87  {
88  auto Acrs = dynamic_cast<const TCrsMatrix*>(A);
89  if(Acrs)
90  {
91  getValuesCrs(Acrs, vals);
92  return;
93  }
94  using range_type = Kokkos::pair<int, int>;
95  using local_inds_host_view_type = typename TRowMatrix::nonconst_local_inds_host_view_type;
96  using values_host_view_type = typename TRowMatrix::nonconst_values_host_view_type;
97  using scalar_type = typename values_host_view_type::value_type;
98 
99  LocalOrdinal nrows = A->getLocalNumRows();
100  size_t nnz = A->getLocalNumEntries();
101  size_t maxNnz = A->getLocalMaxNumRowEntries();
102 
103  vals = ScalarArray("Values", nnz);
104  auto valsHost = Kokkos::create_mirror(vals);
105  local_inds_host_view_type lclColInds ("lclColinds", maxNnz);
106 
107  nnz = 0;
108  for(LocalOrdinal i = 0; i < nrows; i++) {
109  size_t NumEntries = A->getNumEntriesInLocalRow(i);
110  auto constLclValues = Kokkos::subview (valsHost, range_type (nnz, nnz+NumEntries));
111  values_host_view_type lclValues (const_cast<scalar_type*>(constLclValues.data()), NumEntries);
112 
113  A->getLocalRowCopy (i, lclColInds, lclValues, NumEntries);
114  nnz += NumEntries;
115  }
116  Kokkos::deep_copy(vals, valsHost);
117  }
118 
120  // \param[in] A The matrix
121  // \param[out] rowptrsHost The rowptrs array, in host space (allocated inside function)
122  // \param[out] rowptrs The rowptrs host array, in device space (allocated inside function). Will have exactly the same values as rowptrsHost.
123  // \param[out] colinds The colinds array, in device space (allocated inside function)
124  static void getStructure(const TRowMatrix* A, OrdinalArrayHost& rowptrsHost, OrdinalArray& rowptrs, OrdinalArray& colinds)
125  {
126  auto Acrs = dynamic_cast<const TCrsMatrix*>(A);
127  if(Acrs)
128  {
129  getStructureCrs(Acrs, rowptrsHost, rowptrs, colinds);
130  return;
131  }
132  //Need to allocate new array, then copy in one row at a time
133  //note: actual rowptrs in the CrsMatrix implementation is size_t, but
134  //FastILU needs them as LocalOrdinal so this function provides an OrdinalArray
135  LocalOrdinal nrows = A->getLocalNumRows();
136  rowptrsHost = OrdinalArrayHost("RowPtrs (host)", nrows + 1);
137 
138  using range_type = Kokkos::pair<int, int>;
139  using values_host_view_type = typename TRowMatrix::nonconst_values_host_view_type;
140  using local_inds_host_view_type = typename TRowMatrix::nonconst_local_inds_host_view_type;
141  using local_ind_type = typename local_inds_host_view_type::value_type;
142  size_t nnz = A->getLocalNumEntries();
143  size_t maxNnz = A->getLocalMaxNumRowEntries();
144 
145  colinds = OrdinalArray("ColInds", nnz);
146  auto colindsHost = Kokkos::create_mirror(colinds);
147  values_host_view_type lclValues ("lclValues", maxNnz);
148 
149  nnz = 0;
150  rowptrsHost[0] = nnz;
151  for(LocalOrdinal i = 0; i < nrows; i++) {
152  size_t NumEntries = A->getNumEntriesInLocalRow(i);
153  auto constLclValues = Kokkos::subview (colindsHost, range_type (nnz, nnz+NumEntries));
154  local_inds_host_view_type lclColInds (const_cast<local_ind_type*>(constLclValues.data()), NumEntries);
155  A->getLocalRowCopy (i, lclColInds, lclValues, NumEntries);
156 
157  nnz += NumEntries;
158  rowptrsHost[i+1] = nnz;
159  }
160 
161  rowptrs = OrdinalArray("RowPtrs", nrows + 1);
162  Kokkos::deep_copy(rowptrs, rowptrsHost);
163  Kokkos::deep_copy(colinds, colindsHost);
164  }
165 
166  private:
167 
169  static void getValuesCrs(const TCrsMatrix* A, ScalarArray& values_)
170  {
171  auto localA = A->getLocalMatrixDevice();
172  auto values = localA.values;
173  auto nnz = values.extent(0);
174  values_ = ScalarArray("Values", nnz );
175  Kokkos::deep_copy(values_, values);
176  }
177 
179  static void getStructureCrs(const TCrsMatrix* A, OrdinalArrayHost& rowptrsHost_, OrdinalArray& rowptrs_, OrdinalArray& colinds_)
180  {
181  //rowptrs really have data type size_t, but need them as LocalOrdinal, so must convert manually
182  auto localA = A->getLocalMatrixDevice();
183  auto rowptrs = localA.graph.row_map;
184  auto colinds = localA.graph.entries;
185  auto numRows = A->getLocalNumRows();
186  auto nnz = colinds.extent(0);
187  //allocate rowptrs, it's a deep copy (colinds is a shallow copy so not necessary for it)
188  rowptrs_ = OrdinalArray("RowPtrs", numRows + 1);
189  colinds_ = OrdinalArray("ColInds", nnz );
190  Kokkos::deep_copy(rowptrs_, rowptrs);
191  Kokkos::deep_copy(colinds_, colinds);
192  // deep-copy to host
193  rowptrsHost_ = Kokkos::create_mirror(rowptrs_);
194  Kokkos::deep_copy(rowptrsHost_, rowptrs_);
195  }
196 };
197 
198 } //Details
199 } //Ifpack2
200 
201 #endif
202 
Wraps a Tpetra::RowMatrix in a filter that reorders local rows and columns.
Definition: Ifpack2_ReorderFilter_decl.hpp:69
Ifpack2 implementation details.
Access only local rows and columns of a sparse matrix.
Definition: Ifpack2_LocalFilter_decl.hpp:160
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:74