Tpetra parallel linear algebra  Version of the Day
Tpetra_Details_getDiagCopyWithoutOffsets_def.hpp
Go to the documentation of this file.
1 /*
2 // @HEADER
3 // ***********************************************************************
4 //
5 // Tpetra: Templated Linear Algebra Services Package
6 // Copyright (2008) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 */
43 
44 #ifndef TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DEF_HPP
45 #define TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DEF_HPP
46 
54 
56 #include "Tpetra_RowGraph.hpp"
57 #include "Tpetra_CrsGraph.hpp"
58 #include "Tpetra_RowMatrix.hpp"
59 #include "Tpetra_Vector.hpp"
60 
61 namespace Tpetra {
62 namespace Details {
63 
64 // Work-around for #499: Implementation of one-argument (no offsets)
65 // getLocalDiagCopy for the NOT fill-complete case.
66 //
67 // NOTE (mfh 18 Jul 2016) This calls functions that are NOT GPU device
68 // functions! Thus, we do NOT use KOKKOS_INLINE_FUNCTION or
69 // KOKKOS_FUNCTION here, because those attempt to mark the functions
70 // they modify as CUDA device functions. This functor is ONLY for
71 // non-CUDA execution spaces!
72 template<class SC, class LO, class GO, class NT>
73 class GetLocalDiagCopyWithoutOffsetsNotFillCompleteFunctor {
74 public:
75  typedef ::Tpetra::RowMatrix<SC, LO, GO, NT> row_matrix_type;
76  typedef ::Tpetra::Vector<SC, LO, GO, NT> vec_type;
77 
78  typedef typename vec_type::impl_scalar_type IST;
79  // The output Vector determines the execution space.
80  typedef typename vec_type::device_type device_type;
81 
82 private:
83  typedef typename vec_type::dual_view_type::t_host::execution_space host_execution_space;
84  typedef typename vec_type::map_type map_type;
85 
86  static bool
87  graphIsSorted (const row_matrix_type& A)
88  {
89  using Teuchos::RCP;
90  using Teuchos::rcp_dynamic_cast;
91  typedef Tpetra::CrsGraph<LO, GO, NT> crs_graph_type;
92  typedef Tpetra::RowGraph<LO, GO, NT> row_graph_type;
93 
94  // We conservatively assume not sorted. RowGraph lacks an
95  // "isSorted" predicate, so we can't know for sure unless the cast
96  // to CrsGraph succeeds.
97  bool sorted = false;
98 
99  RCP<const row_graph_type> G_row = A.getGraph ();
100  if (! G_row.is_null ()) {
101  RCP<const crs_graph_type> G_crs =
102  rcp_dynamic_cast<const crs_graph_type> (G_row);
103  if (! G_crs.is_null ()) {
104  sorted = G_crs->isSorted ();
105  }
106  }
107 
108  return sorted;
109  }
110 
111 public:
112  // lclNumErrs [out] Total count of errors on this process.
113  GetLocalDiagCopyWithoutOffsetsNotFillCompleteFunctor (LO& lclNumErrs,
114  vec_type& diag,
115  const row_matrix_type& A) :
116  A_ (A),
117  lclRowMap_ (*A.getRowMap ()),
118  lclColMap_ (*A.getColMap ()),
119  sorted_ (graphIsSorted (A))
120  {
121  const LO lclNumRows = static_cast<LO> (diag.getLocalLength ());
122  {
123  const LO matLclNumRows =
124  static_cast<LO> (lclRowMap_.getNodeNumElements ());
125  TEUCHOS_TEST_FOR_EXCEPTION
126  (lclNumRows != matLclNumRows, std::invalid_argument,
127  "diag.getLocalLength() = " << lclNumRows << " != "
128  "A.getRowMap()->getNodeNumElements() = " << matLclNumRows << ".");
129  }
130 
131  // Side effects start below this point.
132  D_lcl_ = diag.getLocalViewHost(Access::OverwriteAll);
133  D_lcl_1d_ = Kokkos::subview (D_lcl_, Kokkos::ALL (), 0);
134 
135  Kokkos::RangePolicy<host_execution_space, LO> range (0, lclNumRows);
136  lclNumErrs = 0;
137  Kokkos::parallel_reduce (range, *this, lclNumErrs);
138 
139  // sync changes back to device, since the user doesn't know that
140  // we had to run on host.
141  //diag.template sync<typename device_type::memory_space> ();
142  }
143 
144  void operator () (const LO& lclRowInd, LO& errCount) const {
145  using KokkosSparse::findRelOffset;
146 
147  D_lcl_1d_(lclRowInd) = Kokkos::Details::ArithTraits<IST>::zero ();
148  const GO gblInd = lclRowMap_.getGlobalElement (lclRowInd);
149  const LO lclColInd = lclColMap_.getLocalElement (gblInd);
150 
151  if (lclColInd == Tpetra::Details::OrdinalTraits<LO>::invalid ()) {
152  errCount++;
153  }
154  else { // row index is also in the column Map on this process
155  typename row_matrix_type::local_inds_host_view_type lclColInds;
156  typename row_matrix_type::values_host_view_type curVals;
157  A_.getLocalRowView(lclRowInd, lclColInds, curVals);
158  LO numEnt = lclColInds.extent(0);
159  // The search hint is always zero, since we only call this
160  // once per row of the matrix.
161  const LO hint = 0;
162  const LO offset =
163  findRelOffset (lclColInds, numEnt, lclColInd, hint, sorted_);
164  if (offset == numEnt) { // didn't find the diagonal column index
165  errCount++;
166  }
167  else {
168  D_lcl_1d_(lclRowInd) = curVals[offset];
169  }
170  }
171  }
172 
173 private:
174  const row_matrix_type& A_;
175  map_type lclRowMap_;
176  map_type lclColMap_;
177  typename vec_type::dual_view_type::t_host D_lcl_;
178  decltype (Kokkos::subview (D_lcl_, Kokkos::ALL (), 0)) D_lcl_1d_;
179  const bool sorted_;
180 };
181 
182 
183 template<class SC, class LO, class GO, class NT>
184 LO
186  const ::Tpetra::RowMatrix<SC, LO, GO, NT>& A,
187  const bool debug)
188 {
189  using ::Tpetra::Details::gathervPrint;
190  using Teuchos::outArg;
191  using Teuchos::REDUCE_MIN;
192  using Teuchos::reduceAll;
193  typedef GetLocalDiagCopyWithoutOffsetsNotFillCompleteFunctor<SC,
194  LO, GO, NT> functor_type;
195 
196  // The functor's constructor does error checking and executes the
197  // thread-parallel kernel.
198 
199  LO lclNumErrs = 0;
200 
201  if (debug) {
202  int lclSuccess = 1;
203  int gblSuccess = 0;
204  std::ostringstream errStrm;
205  Teuchos::RCP<const Teuchos::Comm<int> > commPtr = A.getComm ();
206  if (commPtr.is_null ()) {
207  return lclNumErrs; // this process does not participate
208  }
209  const Teuchos::Comm<int>& comm = *commPtr;
210 
211  try {
212  functor_type functor (lclNumErrs, diag, A);
213  }
214  catch (std::exception& e) {
215  lclSuccess = -1;
216  errStrm << "Process " << A.getComm ()->getRank () << ": "
217  << e.what () << std::endl;
218  }
219  if (lclNumErrs != 0) {
220  lclSuccess = 0;
221  }
222 
223  reduceAll<int, int> (comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
224  if (gblSuccess == -1) {
225  if (comm.getRank () == 0) {
226  // We gather into std::cerr, rather than using an
227  // std::ostringstream, because there might be a lot of MPI
228  // processes. It could take too much memory to gather all the
229  // messages to Process 0 before printing. gathervPrint gathers
230  // and prints one message at a time, thus saving memory. I
231  // don't want to run out of memory while trying to print an
232  // error message; that would hide the real problem.
233  std::cerr << "getLocalDiagCopyWithoutOffsetsNotFillComplete threw an "
234  "exception on one or more MPI processes in the matrix's comunicator."
235  << std::endl;
236  }
237  gathervPrint (std::cerr, errStrm.str (), comm);
238  // Don't need to print anything here, since we've already
239  // printed to std::cerr above.
240  TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "");
241  }
242  else if (gblSuccess == 0) {
243  TEUCHOS_TEST_FOR_EXCEPTION
244  (gblSuccess != 1, std::runtime_error,
245  "getLocalDiagCopyWithoutOffsetsNotFillComplete failed on "
246  "one or more MPI processes in the matrix's communicator.");
247  }
248  }
249  else { // ! debug
250  functor_type functor (lclNumErrs, diag, A);
251  }
252 
253  return lclNumErrs;
254 }
255 
256 } // namespace Details
257 } // namespace Tpetra
258 
259 // Explicit template instantiation macro for
260 // getLocalDiagCopyWithoutOffsetsNotFillComplete. NOT FOR USERS!!!
261 // Must be used inside the Tpetra namespace.
262 #define TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_INSTANT( SCALAR, LO, GO, NODE ) \
263  template LO \
264  Details::getLocalDiagCopyWithoutOffsetsNotFillComplete< SCALAR, LO, GO, NODE > \
265  ( ::Tpetra::Vector< SCALAR, LO, GO, NODE >& diag, \
266  const ::Tpetra::RowMatrix< SCALAR, LO, GO, NODE >& A, \
267  const bool debug);
268 
269 #endif // TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DEF_HPP
Declaration of a function that prints strings from each process.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
An abstract interface for graphs accessed by rows.
A read-only, row-oriented interface to a sparse matrix.
A distributed dense vector.
base_type::impl_scalar_type impl_scalar_type
The type used internally in place of Scalar.
Node::device_type device_type
The Kokkos device type.
base_type::map_type map_type
The type of the Map specialization used by this class.
Implementation details of Tpetra.
LO getLocalDiagCopyWithoutOffsetsNotFillComplete(::Tpetra::Vector< SC, LO, GO, NT > &diag, const ::Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool debug=false)
Given a locally indexed, global sparse matrix, extract the matrix's diagonal entries into a Tpetra::V...
void gathervPrint(std::ostream &out, const std::string &s, const Teuchos::Comm< int > &comm)
On Process 0 in the given communicator, print strings from each process in that communicator,...
Namespace Tpetra contains the class and methods constituting the Tpetra library.