40 #ifndef TPETRA_DETAILS_CHECKVIEW_HPP
41 #define TPETRA_DETAILS_CHECKVIEW_HPP
51 #include "Kokkos_DualView.hpp"
52 #include "Teuchos_TypeNameTraits.hpp"
53 #include "Teuchos_Comm.hpp"
54 #include "Teuchos_CommHelpers.hpp"
60 std::string memorySpaceName (
const void* ptr);
82 template<
class DataType,
class ... Properties>
85 (std::ostream* lclErrStrm,
86 const int myMpiProcessRank,
87 const Kokkos::View<DataType, Properties...>& view)
89 using Teuchos::TypeNameTraits;
91 using view_type = Kokkos::View<DataType, Properties...>;
93 if (view.size () == 0) {
99 auto ptr = view.data ();
101 if (ptr ==
nullptr) {
102 if (lclErrStrm !=
nullptr) {
103 const std::string viewName = TypeNameTraits<view_type>::name ();
104 *lclErrStrm <<
"Proc " << myMpiProcessRank <<
": Kokkos::View "
105 "of type " << viewName <<
" has nonzero size " << view.size ()
106 <<
" but a null pointer." << endl;
118 template<
class DataType ,
119 class Arg1Type = void ,
120 class Arg2Type = void ,
121 class Arg3Type =
void>
124 (std::ostream*
const lclErrStrm,
125 const int myMpiProcessRank,
126 const Kokkos::DualView<DataType, Arg1Type, Arg2Type, Arg3Type>& dv)
128 const bool dev_good =
131 const bool host_good =
134 const bool good = dev_good && host_good;
135 if (! good && lclErrStrm !=
nullptr) {
136 using Teuchos::TypeNameTraits;
139 Kokkos::DualView<DataType, Arg1Type, Arg2Type, Arg3Type>;
141 const std::string dvName = TypeNameTraits<dv_type>::name ();
142 *lclErrStrm <<
"Proc " << myMpiProcessRank <<
": Kokkos::DualView "
143 "of type " << dvName <<
" has one or more invalid Views. See "
144 "above error messages from this MPI process for details." << endl;
149 template<
class DataType ,
150 class Arg1Type = void ,
151 class Arg2Type = void ,
152 class Arg3Type =
void>
154 checkGlobalDualViewValidity
155 (std::ostream*
const gblErrStrm,
156 const Kokkos::DualView<DataType, Arg1Type, Arg2Type, Arg3Type>& dv,
158 const Teuchos::Comm<int>*
const comm)
161 const int myRank = comm ==
nullptr ? 0 : comm->getRank ();
162 std::ostringstream lclErrStrm;
166 const bool lclValid =
168 lclSuccess = lclValid ? 1 : 0;
170 catch (std::exception& e) {
171 lclErrStrm <<
"Proc " << myRank <<
": checkLocalDualViewValidity "
172 "threw an exception: " << e.what () << endl;
176 lclErrStrm <<
"Proc " << myRank <<
": checkLocalDualViewValidity "
177 "threw an exception not a subclass of std::exception." << endl;
182 if (comm ==
nullptr) {
183 gblSuccess = lclSuccess;
186 using Teuchos::outArg;
187 using Teuchos::REDUCE_MIN;
188 using Teuchos::reduceAll;
189 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
192 if (gblSuccess != 1 && gblErrStrm !=
nullptr) {
193 *gblErrStrm <<
"On at least one (MPI) process, the "
194 "Kokkos::DualView has "
195 "either the device or host pointer in the "
196 "DualView equal to null, but the DualView has a nonzero number of "
197 "rows. For more detailed information, please rerun with the "
198 "TPETRA_VERBOSE environment variable set to 1. ";
200 *gblErrStrm <<
" Here are error messages from all "
201 "processes:" << endl;
202 if (comm ==
nullptr) {
203 *gblErrStrm << lclErrStrm.str ();
212 return gblSuccess == 1;
Declaration of a function that prints strings from each process.
Implementation details of Tpetra.
bool checkLocalDualViewValidity(std::ostream *const lclErrStrm, const int myMpiProcessRank, const Kokkos::DualView< DataType, Arg1Type, Arg2Type, Arg3Type > &dv)
Is the given Kokkos::DualView valid?
bool checkLocalViewValidity(std::ostream *lclErrStrm, const int myMpiProcessRank, const Kokkos::View< DataType, Properties... > &view)
Is the given View valid?
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.