Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_checkView.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_DETAILS_CHECKVIEW_HPP
11#define TPETRA_DETAILS_CHECKVIEW_HPP
12
19
21#include "Tpetra_Details_WrappedDualView.hpp"
22#include "Kokkos_DualView.hpp"
23#include "Teuchos_TypeNameTraits.hpp"
24#include "Teuchos_Comm.hpp"
25#include "Teuchos_CommHelpers.hpp"
26#include <sstream>
27
28namespace Tpetra {
29namespace Details {
30
31std::string memorySpaceName(const void* ptr);
32
53template <class DataType, class... Properties>
55 const int myMpiProcessRank,
56 const Kokkos::View<DataType, Properties...>& view) {
57 using std::endl;
58 using Teuchos::TypeNameTraits;
59 using view_type = Kokkos::View<DataType, Properties...>;
60
61 if (view.size() == 0) {
62 // Kokkos::View can be zero size with a nonnull pointer.
63 // Even std::vector can have this behavior.
64 return true;
65 } else { // nonzero size
66 auto ptr = view.data();
67
68 if (ptr == nullptr) {
69 if (lclErrStrm != nullptr) {
70 const std::string viewName = TypeNameTraits<view_type>::name();
71 *lclErrStrm << "Proc " << myMpiProcessRank << ": Kokkos::View "
72 "of type "
73 << viewName << " has nonzero size " << view.size()
74 << " but a null pointer." << endl;
75 }
76 return false;
77 } else
78 return true;
79 }
80}
81
85template <class DataType,
86 class... Args>
87bool checkLocalDualViewValidity(std::ostream* const lclErrStrm,
88 const int myMpiProcessRank,
89 const Kokkos::DualView<DataType, Args...>& dv) {
90 const bool dev_good =
92 dv.view_device());
93 const bool host_good =
95 dv.view_host());
96 const bool good = dev_good && host_good;
97 if (!good && lclErrStrm != nullptr) {
98 using std::endl;
99 using Teuchos::TypeNameTraits;
100 using dv_type =
101 Kokkos::DualView<DataType, Args...>;
102
103 const std::string dvName = TypeNameTraits<dv_type>::name();
104 *lclErrStrm << "Proc " << myMpiProcessRank << ": Kokkos::DualView "
105 "of type "
106 << dvName << " has one or more invalid Views. See "
107 "above error messages from this MPI process for details."
108 << endl;
109 }
110 return good;
111}
112
113template <class DataType,
114 class... Args>
115bool checkGlobalDualViewValidity(std::ostream* const gblErrStrm,
116 const Kokkos::DualView<DataType, Args...>& dv,
117 const bool verbose,
118 const Teuchos::Comm<int>* const comm) {
119 using std::endl;
120 const int myRank = comm == nullptr ? 0 : comm->getRank();
121 std::ostringstream lclErrStrm;
122 int lclSuccess = 1;
123
124 try {
125 const bool lclValid =
127 lclSuccess = lclValid ? 1 : 0;
128 } catch (std::exception& e) {
129 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
130 "threw an exception: "
131 << e.what() << endl;
132 lclSuccess = 0;
133 } catch (...) {
134 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
135 "threw an exception not a subclass of std::exception."
136 << endl;
137 lclSuccess = 0;
138 }
139
140 int gblSuccess = 0; // output argument
141 if (comm == nullptr) {
142 gblSuccess = lclSuccess;
143 } else {
144 using Teuchos::outArg;
145 using Teuchos::REDUCE_MIN;
146 using Teuchos::reduceAll;
147 reduceAll(*comm, REDUCE_MIN, lclSuccess, outArg(gblSuccess));
148 }
149
150 if (gblSuccess != 1 && gblErrStrm != nullptr) {
151 *gblErrStrm << "On at least one (MPI) process, the "
152 "Kokkos::DualView has "
153 "either the device or host pointer in the "
154 "DualView equal to null, but the DualView has a nonzero number of "
155 "rows. For more detailed information, please rerun with the "
156 "TPETRA_VERBOSE environment variable set to 1. ";
157 if (verbose) {
158 *gblErrStrm << " Here are error messages from all "
159 "processes:"
160 << endl;
161 if (comm == nullptr) {
162 *gblErrStrm << lclErrStrm.str();
163 } else {
165 gathervPrint(*gblErrStrm, lclErrStrm.str(), *comm);
166 }
167 }
168 *gblErrStrm << endl;
169 }
170 return gblSuccess == 1;
171}
172
176template <class DataType,
177 class... Args>
179 const int myMpiProcessRank,
180 const Tpetra::Details::WrappedDualView<Kokkos::DualView<DataType, Args...> >& dv) {
181 const bool dev_good = dv.is_valid_device();
182 const bool host_good = dv.is_valid_host();
183 const bool good = dev_good && host_good;
184 if (!good && lclErrStrm != nullptr) {
185 using std::endl;
186 using Teuchos::TypeNameTraits;
187 using dv_type =
188 Tpetra::Details::WrappedDualView<Kokkos::DualView<DataType, Args...> >;
189
190 const std::string dvName = TypeNameTraits<dv_type>::name();
191 *lclErrStrm << "Proc " << myMpiProcessRank << ": Tpetra::WrappedDualView "
192 "of type "
193 << dvName << " has one or more invalid Views. See "
194 "above error messages from this MPI process for details."
195 << endl;
196 }
197 return good;
198}
199
200template <class DataType,
201 class... Args>
202bool checkGlobalWrappedDualViewValidity(std::ostream* const gblErrStrm,
203 const Tpetra::Details::WrappedDualView<Kokkos::DualView<DataType, Args...> >& dv,
204 const bool verbose,
205 const Teuchos::Comm<int>* const comm) {
206 using std::endl;
207 const int myRank = comm == nullptr ? 0 : comm->getRank();
208 std::ostringstream lclErrStrm;
209 int lclSuccess = 1;
210
211 try {
212 const bool lclValid =
214 lclSuccess = lclValid ? 1 : 0;
215 } catch (std::exception& e) {
216 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
217 "threw an exception: "
218 << e.what() << endl;
219 lclSuccess = 0;
220 } catch (...) {
221 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
222 "threw an exception not a subclass of std::exception."
223 << endl;
224 lclSuccess = 0;
225 }
226
227 int gblSuccess = 0; // output argument
228 if (comm == nullptr) {
229 gblSuccess = lclSuccess;
230 } else {
231 using Teuchos::outArg;
232 using Teuchos::REDUCE_MIN;
233 using Teuchos::reduceAll;
234 reduceAll(*comm, REDUCE_MIN, lclSuccess, outArg(gblSuccess));
235 }
236
237 if (gblSuccess != 1 && gblErrStrm != nullptr) {
238 *gblErrStrm << "On at least one (MPI) process, the "
239 "Kokkos::DualView has "
240 "either the device or host pointer in the "
241 "DualView equal to null, but the DualView has a nonzero number of "
242 "rows. For more detailed information, please rerun with the "
243 "TPETRA_VERBOSE environment variable set to 1. ";
244 if (verbose) {
245 *gblErrStrm << " Here are error messages from all "
246 "processes:"
247 << endl;
248 if (comm == nullptr) {
249 *gblErrStrm << lclErrStrm.str();
250 } else {
252 gathervPrint(*gblErrStrm, lclErrStrm.str(), *comm);
253 }
254 }
255 *gblErrStrm << endl;
256 }
257 return gblSuccess == 1;
258}
259
260} // namespace Details
261} // namespace Tpetra
262
263#endif // TPETRA_DETAILS_CHECKVIEW_HPP
Declaration of a function that prints strings from each process.
Struct that holds views of the contents of a CrsMatrix.
A wrapper around Kokkos::DualView to safely manage data that might be replicated between host and dev...
Implementation details of Tpetra.
bool checkLocalDualViewValidity(std::ostream *const lclErrStrm, const int myMpiProcessRank, const Kokkos::DualView< DataType, Args... > &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?
bool checkLocalWrappedDualViewValidity(std::ostream *const lclErrStrm, const int myMpiProcessRank, const Tpetra::Details::WrappedDualView< Kokkos::DualView< DataType, Args... > > &dv)
Is the given Tpetra::WrappedDualView 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.