Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_CrsArrays.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
14
15#ifndef __IFPACK2_CRSARRAYS_DECL_HPP__
16#define __IFPACK2_CRSARRAYS_DECL_HPP__
17
18#include <Tpetra_RowMatrix.hpp>
19#include <Tpetra_CrsMatrix.hpp>
20#include <Tpetra_KokkosCompat_DefaultNode.hpp>
21#include <Tpetra_BlockCrsMatrix_Helpers_decl.hpp>
22#include <KokkosSparse_CrsMatrix.hpp>
23#include <Ifpack2_LocalFilter.hpp>
24#include <Ifpack2_ReorderFilter.hpp>
25
26namespace Ifpack2
27{
28namespace Details
29{
30
31//Utility for getting the local values, rowptrs and colinds (in Kokkos::Views) for any RowMatrix
32//Used by Fic, Filu and Fildl but may also be useful in other classes
33template<typename Scalar, typename ImplScalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
34struct CrsArrayReader
35{
36 typedef typename Node::device_type device_type;
37 typedef typename device_type::execution_space execution_space;
38 typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TRowMatrix;
39 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TCrsMatrix;
40 typedef Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TBcrsMatrix;
41 typedef Ifpack2::LocalFilter<TRowMatrix> Filter;
42 typedef Ifpack2::ReorderFilter<TRowMatrix> ReordFilter;
43 typedef KokkosSparse::CrsMatrix<ImplScalar, LocalOrdinal, execution_space> KCrsMatrix;
44 typedef Kokkos::View<LocalOrdinal*, execution_space> OrdinalArray;
45 typedef Kokkos::View<ImplScalar*, execution_space> ScalarArray;
46 typedef typename OrdinalArray::HostMirror OrdinalArrayHost;
48 typedef Kokkos::Serial functor_space;
49 typedef Kokkos::RangePolicy<functor_space, int> RangePol;
50
52 // \param[in] A The matrix
53 // \param[out] vals The values array on device (allocated inside function)
54 // \param[in] rowptrs The rowptrs host view provided by getStructure()
55 static void getValues(const TRowMatrix* A, ScalarArray& vals, OrdinalArrayHost& rowptrs)
56 {
57 auto Acrs = dynamic_cast<const TCrsMatrix*>(A);
58 auto Abcrs = dynamic_cast<const TBcrsMatrix*>(A);
59 if(Acrs)
60 {
61 getValuesCrs(Acrs, vals);
62 return;
63 }
64 if (Abcrs) {
65 getValuesBcrs(Abcrs, vals);
66 return;
67 }
68 using range_type = Kokkos::pair<int, int>;
69 using local_inds_host_view_type = typename TRowMatrix::nonconst_local_inds_host_view_type;
70 using values_host_view_type = typename TRowMatrix::nonconst_values_host_view_type;
71 using scalar_type = typename values_host_view_type::value_type;
72
73 LocalOrdinal nrows = A->getLocalNumRows();
74 size_t nnz = A->getLocalNumEntries();
75 size_t maxNnz = A->getLocalMaxNumRowEntries();
76
77 vals = ScalarArray("Values", nnz);
78 auto valsHost = Kokkos::create_mirror(vals);
79 local_inds_host_view_type lclColInds ("lclColinds", maxNnz);
80
81 nnz = 0;
82 for(LocalOrdinal i = 0; i < nrows; i++) {
83 size_t NumEntries = A->getNumEntriesInLocalRow(i);
84 auto constLclValues = Kokkos::subview (valsHost, range_type (nnz, nnz+NumEntries));
85 values_host_view_type lclValues (const_cast<scalar_type*>(constLclValues.data()), NumEntries);
86
87 A->getLocalRowCopy (i, lclColInds, lclValues, NumEntries);
88 nnz += NumEntries;
89 }
90 Kokkos::deep_copy(vals, valsHost);
91 }
92
94 // \param[in] A The matrix
95 // \param[out] rowptrsHost The rowptrs array, in host space (allocated inside function)
96 // \param[out] rowptrs The rowptrs host array, in device space (allocated inside function). Will have exactly the same values as rowptrsHost.
97 // \param[out] colinds The colinds array, in device space (allocated inside function)
98 static void getStructure(const TRowMatrix* A, OrdinalArrayHost& rowptrsHost, OrdinalArray& rowptrs, OrdinalArray& colinds)
99 {
100 auto Acrs = dynamic_cast<const TCrsMatrix*>(A);
101 auto Abcrs = dynamic_cast<const TBcrsMatrix*>(A);
102 if(Acrs)
103 {
104 getStructureCrs(Acrs, rowptrsHost, rowptrs, colinds);
105 return;
106 }
107 if (Abcrs) {
108 getStructureBcrs(Abcrs, rowptrsHost, rowptrs, colinds);
109 return;
110 }
111 //Need to allocate new array, then copy in one row at a time
112 //note: actual rowptrs in the CrsMatrix implementation is size_t, but
113 //FastILU needs them as LocalOrdinal so this function provides an OrdinalArray
114 LocalOrdinal nrows = A->getLocalNumRows();
115 rowptrsHost = OrdinalArrayHost("RowPtrs (host)", nrows + 1);
116
117 using range_type = Kokkos::pair<int, int>;
118 using values_host_view_type = typename TRowMatrix::nonconst_values_host_view_type;
119 using local_inds_host_view_type = typename TRowMatrix::nonconst_local_inds_host_view_type;
120 using local_ind_type = typename local_inds_host_view_type::value_type;
121 size_t nnz = A->getLocalNumEntries();
122 size_t maxNnz = A->getLocalMaxNumRowEntries();
123
124 colinds = OrdinalArray("ColInds", nnz);
125 auto colindsHost = Kokkos::create_mirror(colinds);
126 values_host_view_type lclValues ("lclValues", maxNnz);
127
128 nnz = 0;
129 rowptrsHost[0] = nnz;
130 for(LocalOrdinal i = 0; i < nrows; i++) {
131 size_t NumEntries = A->getNumEntriesInLocalRow(i);
132 auto constLclValues = Kokkos::subview (colindsHost, range_type (nnz, nnz+NumEntries));
133 local_inds_host_view_type lclColInds (const_cast<local_ind_type*>(constLclValues.data()), NumEntries);
134 A->getLocalRowCopy (i, lclColInds, lclValues, NumEntries);
135
136 nnz += NumEntries;
137 rowptrsHost[i+1] = nnz;
138 }
139
140 rowptrs = OrdinalArray("RowPtrs", nrows + 1);
141 Kokkos::deep_copy(rowptrs, rowptrsHost);
142 Kokkos::deep_copy(colinds, colindsHost);
143 }
144
145 private:
146
148 static void getValuesCrs(const TCrsMatrix* A, ScalarArray& values_)
149 {
150 auto localA = A->getLocalMatrixDevice();
151 auto values = localA.values;
152 auto nnz = values.extent(0);
153 values_ = ScalarArray("Values", nnz );
154 Kokkos::deep_copy(values_, values);
155 }
156
158 static void getStructureCrs(const TCrsMatrix* A, OrdinalArrayHost& rowptrsHost_, OrdinalArray& rowptrs_, OrdinalArray& colinds_)
159 {
160 //rowptrs really have data type size_t, but need them as LocalOrdinal, so must convert manually
161 auto localA = A->getLocalMatrixDevice();
162 auto rowptrs = localA.graph.row_map;
163 auto colinds = localA.graph.entries;
164 auto numRows = A->getLocalNumRows();
165 auto nnz = colinds.extent(0);
166 //allocate rowptrs, it's a deep copy (colinds is a shallow copy so not necessary for it)
167 rowptrs_ = OrdinalArray("RowPtrs", numRows + 1);
168 colinds_ = OrdinalArray("ColInds", nnz );
169 Kokkos::deep_copy(rowptrs_, rowptrs);
170 Kokkos::deep_copy(colinds_, colinds);
171 // deep-copy to host
172 rowptrsHost_ = Kokkos::create_mirror(rowptrs_);
173 Kokkos::deep_copy(rowptrsHost_, rowptrs_);
174 }
175
177 static void getValuesBcrs(const TBcrsMatrix* A, ScalarArray& values_)
178 {
179 auto localA = A->getLocalMatrixDevice();
180 auto values = localA.values;
181 auto nnz = values.extent(0);
182 values_ = ScalarArray("Values", nnz );
183 Kokkos::deep_copy(values_, values);
184 }
185
187 static void getStructureBcrs(const TBcrsMatrix* A, OrdinalArrayHost& rowptrsHost_, OrdinalArray& rowptrs_, OrdinalArray& colinds_)
188 {
189 //rowptrs really have data type size_t, but need them as LocalOrdinal, so must convert manually
190 auto localA = A->getLocalMatrixDevice();
191 auto rowptrs = localA.graph.row_map;
192 auto colinds = localA.graph.entries;
193 auto numRows = A->getLocalNumRows();
194 auto nnz = colinds.extent(0);
195 //allocate rowptrs, it's a deep copy (colinds is a shallow copy so not necessary for it)
196 rowptrs_ = OrdinalArray("RowPtrs", numRows + 1);
197 colinds_ = OrdinalArray("ColInds", nnz );
198 Kokkos::deep_copy(rowptrs_, rowptrs);
199 Kokkos::deep_copy(colinds_, colinds);
200 // deep-copy to host
201 rowptrsHost_ = Kokkos::create_mirror(rowptrs_);
202 Kokkos::deep_copy(rowptrsHost_, rowptrs_);
203 }
204
205};
206
207} //Details
208} //Ifpack2
209
210#endif
211
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41