Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tpetra_RowMatrixTransposer_def.hpp
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
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 #ifndef TPETRA_ROWMATRIXTRANSPOSER_DEF_HPP
43 #define TPETRA_ROWMATRIXTRANSPOSER_DEF_HPP
44 
45 #include "Tpetra_CrsMatrix.hpp"
46 #include "Tpetra_Export.hpp"
49 #include "Teuchos_ParameterList.hpp"
50 #include "Teuchos_TimeMonitor.hpp"
51 
52 namespace Tpetra {
53 
54 template<class Scalar,
55  class LocalOrdinal,
56  class GlobalOrdinal,
57  class Node>
59 RowMatrixTransposer (const Teuchos::RCP<const crs_matrix_type>& origMatrix,
60  const std::string& label)
61  : origMatrix_ (origMatrix), label_ (label)
62 {}
63 
64 template<class Scalar,
65  class LocalOrdinal,
66  class GlobalOrdinal,
67  class Node>
68 Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
70 createTranspose (const Teuchos::RCP<Teuchos::ParameterList> &params)
71 {
72  using Teuchos::RCP;
73  // Do the local transpose
74  RCP<crs_matrix_type> transMatrixWithSharedRows = createTransposeLocal (params);
75 
76 #ifdef HAVE_TPETRA_MMM_TIMINGS
77  const std::string prefix = std::string ("Tpetra ") + label_ + ": ";
78  using Teuchos::TimeMonitor;
79  TimeMonitor MM (*TimeMonitor::getNewTimer (prefix + "Transpose TAFC"));
80 #endif
81 
82  // If transMatrixWithSharedRows has an exporter, that's what we
83  // want. If it doesn't, the rows aren't actually shared, and we're
84  // done!
86  RCP<const export_type> exporter =
87  transMatrixWithSharedRows->getGraph ()->getExporter ();
88  if (exporter.is_null ()) {
89  return transMatrixWithSharedRows;
90  }
91  else {
92  Teuchos::ParameterList labelList;
93 #ifdef HAVE_TPETRA_MMM_TIMINGS
94  labelList.set("Timer Label", label_);
95 #endif
96  if(! params.is_null ()) {
97  const char paramName[] = "compute global constants";
98  labelList.set (paramName, params->get (paramName, true));
99  }
100  // Use the Export object to do a fused Export and fillComplete.
101  return exportAndFillCompleteCrsMatrix<crs_matrix_type>
102  (transMatrixWithSharedRows, *exporter, Teuchos::null,
103  Teuchos::null, Teuchos::rcpFromRef (labelList));
104  }
105 }
106 
107 
108 // mfh 03 Feb 2013: In a definition outside the class like this, the
109 // return value is considered outside the class scope (for things like
110 // resolving typedefs), but the arguments are considered inside the
111 // class scope.
112 template<class Scalar,
113  class LocalOrdinal,
114  class GlobalOrdinal,
115  class Node>
116 Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
118 createTransposeLocal (const Teuchos::RCP<Teuchos::ParameterList>& params)
119 {
120  using Teuchos::Array;
121  using Teuchos::ArrayRCP;
122  using Teuchos::ArrayView;
123  using Teuchos::RCP;
124  using Teuchos::rcp;
125  using Teuchos::rcp_dynamic_cast;
126  using LO = LocalOrdinal;
127  using GO = GlobalOrdinal;
129  using import_type = Tpetra::Import<LO, GO, Node>;
130  using export_type = Tpetra::Export<LO, GO, Node>;
131 
132 #ifdef HAVE_TPETRA_MMM_TIMINGS
133  std::string prefix = std::string("Tpetra ") + label_ + ": ";
134  using Teuchos::TimeMonitor;
135  TimeMonitor MM (*TimeMonitor::getNewTimer (prefix + "Transpose Local"));
136 #endif
137 
138  const bool sort = [&] () {
139  constexpr bool sortDefault = true; // see #4607 discussion
140  const char sortParamName[] = "sort";
141  return params.get () == nullptr ? sortDefault :
142  params->get (sortParamName, sortDefault);
143  } ();
144 
145  const LO lclNumCols (origMatrix_->getNodeNumCols ());
146  const LO lclNumRows (origMatrix_->getNodeNumRows ());
147  const size_t nnz (origMatrix_->getNodeNumEntries ());
148 
149  RCP<const crs_matrix_type> crsMatrix =
150  rcp_dynamic_cast<const crs_matrix_type> (origMatrix_);
151  if (crsMatrix.is_null ()) {
152  auto rowMap = origMatrix_->getRowMap ();
153  if (rowMap->isOneToOne ()) {
154  Teuchos::Array<size_t> numEntPerRow (lclNumRows);
155  for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
156  numEntPerRow[lclRow] = origMatrix_->getNumEntriesInLocalRow (lclRow);
157  }
158  auto colMap = origMatrix_->getColMap ();
159 
160  RCP<crs_matrix_type> crsMatrix_nc =
161 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
162  rcp (new crs_matrix_type (rowMap, colMap, numEntPerRow (),
163  Tpetra::StaticProfile));
164 #else
165  rcp (new crs_matrix_type (rowMap, colMap, numEntPerRow ()));
166 #endif // TPETRA_ENABLE_DEPRECATED_CODE
167 
168  // When source & target Maps are same, Import just copies.
169  import_type imp (rowMap, rowMap);
170  crsMatrix_nc->doImport (*origMatrix_, imp, Tpetra::REPLACE);
171  crsMatrix_nc->fillComplete (origMatrix_->getDomainMap (),
172  origMatrix_->getRangeMap ());
173  crsMatrix = crsMatrix_nc;
174  }
175  else {
176  TEUCHOS_ASSERT( false ); // not implemented (it wasn't before)
177  }
178  }
179 
180  using local_matrix_type = typename crs_matrix_type::local_matrix_type;
181  using local_graph_type = typename crs_matrix_type::local_graph_type;
182  using offset_type = typename local_graph_type::size_type;
183  using row_map_type = typename local_matrix_type::row_map_type::non_const_type;
184  using index_type = typename local_matrix_type::index_type::non_const_type;
185  using values_type = typename local_matrix_type::values_type::non_const_type;
186  using execution_space = typename local_matrix_type::execution_space;
187 
188  local_matrix_type lclMatrix = crsMatrix->getLocalMatrix ();
189  local_graph_type lclGraph = lclMatrix.graph;
190 
191  // Determine how many nonzeros there are per row in the transpose.
192  using DT = typename crs_matrix_type::device_type;
193  Kokkos::View<LO*, DT> t_counts ("transpose row counts", lclNumCols);
194  using range_type = Kokkos::RangePolicy<LO, execution_space>;
195  Kokkos::parallel_for
196  ("Compute row counts of local transpose",
197  range_type (0, lclNumRows),
198  KOKKOS_LAMBDA (const LO row) {
199  auto rowView = lclGraph.rowConst(row);
200  const LO length = rowView.length;
201 
202  for (LO colID = 0; colID < length; ++colID) {
203  const LO col = rowView(colID);
204  Kokkos::atomic_fetch_add (&t_counts[col], LO (1));
205  }
206  });
207 
208  using Kokkos::view_alloc;
209  using Kokkos::WithoutInitializing;
210  row_map_type t_offsets
211  (view_alloc ("transpose ptr", WithoutInitializing), lclNumCols + 1);
212 
213  // TODO (mfh 10 Jul 2019): This returns the sum of all counts,
214  // which could be useful for checking nnz.
215  (void) Details::computeOffsetsFromCounts (t_offsets, t_counts);
216 
217  index_type t_cols
218  (view_alloc ("transpose lcl ind", WithoutInitializing), nnz);
219  values_type t_vals
220  (view_alloc ("transpose val", WithoutInitializing), nnz);
221  Kokkos::parallel_for
222  ("Compute local transpose",
223  range_type (0, lclNumRows),
224  KOKKOS_LAMBDA (const LO row) {
225  auto rowView = lclMatrix.rowConst(row);
226  const LO length = rowView.length;
227 
228  for (LO colID = 0; colID < length; colID++) {
229  const LO col = rowView.colidx(colID);
230  const offset_type beg = t_offsets[col];
231  const LO old_count =
232  Kokkos::atomic_fetch_sub (&t_counts[col], LO (1));
233  const LO len (t_offsets[col+1] - beg);
234  const offset_type insert_pos = beg + (len - old_count);
235  t_cols[insert_pos] = row;
236  t_vals[insert_pos] = rowView.value(colID);
237  }
238  });
239 
240  // Invariant: At this point, all entries of t_counts are zero.
241  // This means we can use it to store new post-merge counts.
242 
243  if (sort) {
244  // NOTE (mfh 11 Jul 2019) Merging is unnecessary: above
245  // parallel_for visits each row of the original matrix once, so
246  // there can be no duplicate column indices in the transpose.
248  Kokkos::parallel_for
249  ("Sort rows of local transpose",
250  range_type (0, lclNumCols),
251  KOKKOS_LAMBDA (const LO lclCol) {
252  const offset_type beg = t_offsets[lclCol];
253  const LO len (t_offsets[lclCol+1] - t_offsets[lclCol]);
254 
255  LO* cols_beg = t_cols.data () + beg;
256  IST* vals_beg = t_vals.data () + beg;
257  shellSortKeysAndValues (cols_beg, vals_beg, len);
258  });
259  }
260 
261  local_matrix_type lclTransposeMatrix ("transpose", lclNumCols,
262  lclNumRows, nnz,
263  t_vals, t_offsets, t_cols);
264 
265  // Prebuild the importers and exporters the no-communication way,
266  // flipping the importers and exporters around.
267  const auto origExport = origMatrix_->getGraph ()->getExporter ();
268  RCP<const import_type> myImport = origExport.is_null () ?
269  Teuchos::null : rcp (new import_type (*origExport));
270  const auto origImport = origMatrix_->getGraph ()->getImporter ();
271  RCP<const export_type> myExport = origImport.is_null () ?
272  Teuchos::null : rcp (new export_type (*origImport));
273 
274  return rcp (new crs_matrix_type (lclTransposeMatrix,
275  origMatrix_->getColMap (),
276  origMatrix_->getRowMap (),
277  origMatrix_->getRangeMap (),
278  origMatrix_->getDomainMap (),
279  myImport, myExport));
280 }
281 //
282 // Explicit instantiation macro
283 //
284 // Must be expanded from within the Tpetra namespace!
285 //
286 
287 #define TPETRA_ROWMATRIXTRANSPOSER_INSTANT(SCALAR,LO,GO,NODE) \
288  template class RowMatrixTransposer< SCALAR, LO , GO , NODE >;
289 
290 } // namespace Tpetra
291 
292 #endif
RowMatrixTransposer(const Teuchos::RCP< const crs_matrix_type > &origMatrix, const std::string &label=std::string())
Constructor that takes the matrix to transpose.
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
Declaration and definition of functions for sorting &quot;short&quot; arrays of keys and corresponding values...
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
Declare and define the functions Tpetra::Details::computeOffsetsFromCounts and Tpetra::computeOffsets...
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
typename crs_graph_type::local_graph_type local_graph_type
The part of the sparse matrix&#39;s graph on each MPI process.
KOKKOS_FUNCTION void shellSortKeysAndValues(KeyType keys[], ValueType values[], const IndexType n)
Shellsort (yes, it&#39;s one word) the input array keys, and apply the resulting permutation to the input...
typename Node::device_type device_type
The Kokkos device type.
Teuchos::RCP< crs_matrix_type > createTransposeLocal(const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Compute and return the transpose of the matrix given to the constructor.
Replace existing values with new values.
Teuchos::RCP< crs_matrix_type > createTranspose(const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Compute and return the transpose of the matrix given to the constructor.
KokkosSparse::CrsMatrix< impl_scalar_type, local_ordinal_type, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
OffsetsViewType::non_const_value_type computeOffsetsFromCounts(const ExecutionSpace &execSpace, const OffsetsViewType &ptr, const CountsViewType &counts)
Compute offsets from counts.
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.