44#ifndef TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DEF_HPP
45#define TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DEF_HPP
56#include "Tpetra_RowGraph.hpp"
57#include "Tpetra_CrsGraph.hpp"
58#include "Tpetra_RowMatrix.hpp"
59#include "Tpetra_Vector.hpp"
72template<
class SC,
class LO,
class GO,
class NT>
73class GetLocalDiagCopyWithoutOffsetsNotFillCompleteFunctor {
75 typedef ::Tpetra::RowMatrix<SC, LO, GO, NT> row_matrix_type;
76 typedef ::Tpetra::Vector<SC, LO, GO, NT> vec_type;
83 typedef typename vec_type::dual_view_type::t_host::execution_space host_execution_space;
87 graphIsSorted (
const row_matrix_type& A)
90 using Teuchos::rcp_dynamic_cast;
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 ();
113 GetLocalDiagCopyWithoutOffsetsNotFillCompleteFunctor (LO& lclNumErrs,
115 const row_matrix_type& A) :
117 lclRowMap_ (*A.getRowMap ()),
118 lclColMap_ (*A.getColMap ()),
119 sorted_ (graphIsSorted (A))
121 const LO lclNumRows =
static_cast<LO
> (diag.getLocalLength ());
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 <<
".");
132 D_lcl_ = diag.getLocalViewHost(Access::OverwriteAll);
133 D_lcl_1d_ = Kokkos::subview (D_lcl_, Kokkos::ALL (), 0);
135 Kokkos::RangePolicy<host_execution_space, LO> range (0, lclNumRows);
137 Kokkos::parallel_reduce (range, *
this, lclNumErrs);
144 void operator () (
const LO& lclRowInd, LO& errCount)
const {
145 using KokkosSparse::findRelOffset;
147 D_lcl_1d_(lclRowInd) = Kokkos::Details::ArithTraits<IST>::zero ();
148 const GO gblInd = lclRowMap_.getGlobalElement (lclRowInd);
149 const LO lclColInd = lclColMap_.getLocalElement (gblInd);
151 if (lclColInd == Tpetra::Details::OrdinalTraits<LO>::invalid ()) {
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);
163 findRelOffset (lclColInds, numEnt, lclColInd, hint, sorted_);
164 if (offset == numEnt) {
168 D_lcl_1d_(lclRowInd) = curVals[offset];
174 const row_matrix_type& A_;
177 typename vec_type::dual_view_type::t_host D_lcl_;
178 decltype (Kokkos::subview (D_lcl_, Kokkos::ALL (), 0)) D_lcl_1d_;
183template<
class SC,
class LO,
class GO,
class NT>
186 const ::Tpetra::RowMatrix<SC, LO, GO, NT>& A,
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;
204 std::ostringstream errStrm;
205 Teuchos::RCP<const Teuchos::Comm<int> > commPtr = A.getComm ();
206 if (commPtr.is_null ()) {
209 const Teuchos::Comm<int>& comm = *commPtr;
212 functor_type functor (lclNumErrs, diag, A);
214 catch (std::exception& e) {
216 errStrm <<
"Process " << A.getComm ()->getRank () <<
": "
217 << e.what () << std::endl;
219 if (lclNumErrs != 0) {
223 reduceAll<int, int> (comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
224 if (gblSuccess == -1) {
225 if (comm.getRank () == 0) {
233 std::cerr <<
"getLocalDiagCopyWithoutOffsetsNotFillComplete threw an "
234 "exception on one or more MPI processes in the matrix's comunicator."
240 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::runtime_error,
"");
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.");
250 functor_type functor (lclNumErrs, diag, A);
262#define TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_INSTANT( SCALAR, LO, GO, NODE ) \
264 Details::getLocalDiagCopyWithoutOffsetsNotFillComplete< SCALAR, LO, GO, NODE > \
265 ( ::Tpetra::Vector< SCALAR, LO, GO, NODE >& diag, \
266 const ::Tpetra::RowMatrix< SCALAR, LO, GO, NODE >& A, \
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 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.