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//
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// ************************************************************************
38// @HEADER
39
40#ifndef TPETRA_DETAILS_CHECKVIEW_HPP
41#define TPETRA_DETAILS_CHECKVIEW_HPP
42
49
51#include "Tpetra_Details_WrappedDualView.hpp"
52#include "Kokkos_DualView.hpp"
53#include "Teuchos_TypeNameTraits.hpp"
54#include "Teuchos_Comm.hpp"
55#include "Teuchos_CommHelpers.hpp"
56#include <sstream>
57
58namespace Tpetra {
59namespace Details {
60
61std::string memorySpaceName (const void* ptr);
62
83template<class DataType, class ... Properties>
84bool
86 (std::ostream* lclErrStrm,
87 const int myMpiProcessRank,
88 const Kokkos::View<DataType, Properties...>& view)
89{
90 using Teuchos::TypeNameTraits;
91 using std::endl;
92 using view_type = Kokkos::View<DataType, Properties...>;
93
94 if (view.size () == 0) {
95 // Kokkos::View can be zero size with a nonnull pointer.
96 // Even std::vector can have this behavior.
97 return true;
98 }
99 else { // nonzero size
100 auto ptr = view.data ();
101
102 if (ptr == nullptr) {
103 if (lclErrStrm != nullptr) {
104 const std::string viewName = TypeNameTraits<view_type>::name ();
105 *lclErrStrm << "Proc " << myMpiProcessRank << ": Kokkos::View "
106 "of type " << viewName << " has nonzero size " << view.size ()
107 << " but a null pointer." << endl;
108 }
109 return false;
110 }
111 else
112 return true;
113 }
114}
115
119template<class DataType ,
120 class... Args>
121bool
123 (std::ostream* const lclErrStrm,
124 const int myMpiProcessRank,
125 const Kokkos::DualView<DataType, Args...>& dv)
126{
127 const bool dev_good =
128 checkLocalViewValidity (lclErrStrm, myMpiProcessRank,
129 dv.view_device ());
130 const bool host_good =
131 checkLocalViewValidity (lclErrStrm, myMpiProcessRank,
132 dv.view_host ());
133 const bool good = dev_good && host_good;
134 if (! good && lclErrStrm != nullptr) {
135 using Teuchos::TypeNameTraits;
136 using std::endl;
137 using dv_type =
138 Kokkos::DualView<DataType, Args...>;
139
140 const std::string dvName = TypeNameTraits<dv_type>::name ();
141 *lclErrStrm << "Proc " << myMpiProcessRank << ": Kokkos::DualView "
142 "of type " << dvName << " has one or more invalid Views. See "
143 "above error messages from this MPI process for details." << endl;
144 }
145 return good;
146}
147
148template<class DataType ,
149 class... Args>
150bool
151checkGlobalDualViewValidity
152(std::ostream* const gblErrStrm,
153 const Kokkos::DualView<DataType, Args...>& dv,
154 const bool verbose,
155 const Teuchos::Comm<int>* const comm)
156{
157 using std::endl;
158 const int myRank = comm == nullptr ? 0 : comm->getRank ();
159 std::ostringstream lclErrStrm;
160 int lclSuccess = 1;
161
162 try {
163 const bool lclValid =
164 checkLocalDualViewValidity (&lclErrStrm, myRank, dv);
165 lclSuccess = lclValid ? 1 : 0;
166 }
167 catch (std::exception& e) {
168 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
169 "threw an exception: " << e.what () << endl;
170 lclSuccess = 0;
171 }
172 catch (...) {
173 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
174 "threw an exception not a subclass of std::exception." << endl;
175 lclSuccess = 0;
176 }
177
178 int gblSuccess = 0; // output argument
179 if (comm == nullptr) {
180 gblSuccess = lclSuccess;
181 }
182 else {
183 using Teuchos::outArg;
184 using Teuchos::REDUCE_MIN;
185 using Teuchos::reduceAll;
186 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
187 }
188
189 if (gblSuccess != 1 && gblErrStrm != nullptr) {
190 *gblErrStrm << "On at least one (MPI) process, the "
191 "Kokkos::DualView has "
192 "either the device or host pointer in the "
193 "DualView equal to null, but the DualView has a nonzero number of "
194 "rows. For more detailed information, please rerun with the "
195 "TPETRA_VERBOSE environment variable set to 1. ";
196 if (verbose) {
197 *gblErrStrm << " Here are error messages from all "
198 "processes:" << endl;
199 if (comm == nullptr) {
200 *gblErrStrm << lclErrStrm.str ();
201 }
202 else {
204 gathervPrint (*gblErrStrm, lclErrStrm.str (), *comm);
205 }
206 }
207 *gblErrStrm << endl;
208 }
209 return gblSuccess == 1;
210}
211
212
216template<class DataType ,
217 class... Args>
218bool
220 (std::ostream* const lclErrStrm,
221 const int myMpiProcessRank,
222 const Tpetra::Details::WrappedDualView<Kokkos::DualView<DataType, Args...> >& dv)
223{
224 const bool dev_good = dv.is_valid_device();
225 const bool host_good = dv. is_valid_host();
226 const bool good = dev_good && host_good;
227 if (! good && lclErrStrm != nullptr) {
228 using Teuchos::TypeNameTraits;
229 using std::endl;
230 using dv_type =
231 Tpetra::Details::WrappedDualView<Kokkos::DualView<DataType, Args...> >;
232
233 const std::string dvName = TypeNameTraits<dv_type>::name ();
234 *lclErrStrm << "Proc " << myMpiProcessRank << ": Tpetra::WrappedDualView "
235 "of type " << dvName << " has one or more invalid Views. See "
236 "above error messages from this MPI process for details." << endl;
237 }
238 return good;
239}
240
241template<class DataType ,
242 class... Args>
243bool
244checkGlobalWrappedDualViewValidity
245(std::ostream* const gblErrStrm,
246 const Tpetra::Details::WrappedDualView<Kokkos::DualView<DataType, Args...> >& dv,
247 const bool verbose,
248 const Teuchos::Comm<int>* const comm)
249{
250 using std::endl;
251 const int myRank = comm == nullptr ? 0 : comm->getRank ();
252 std::ostringstream lclErrStrm;
253 int lclSuccess = 1;
254
255 try {
256 const bool lclValid =
257 checkLocalWrappedDualViewValidity (&lclErrStrm, myRank, dv);
258 lclSuccess = lclValid ? 1 : 0;
259 }
260 catch (std::exception& e) {
261 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
262 "threw an exception: " << e.what () << endl;
263 lclSuccess = 0;
264 }
265 catch (...) {
266 lclErrStrm << "Proc " << myRank << ": checkLocalDualViewValidity "
267 "threw an exception not a subclass of std::exception." << endl;
268 lclSuccess = 0;
269 }
270
271 int gblSuccess = 0; // output argument
272 if (comm == nullptr) {
273 gblSuccess = lclSuccess;
274 }
275 else {
276 using Teuchos::outArg;
277 using Teuchos::REDUCE_MIN;
278 using Teuchos::reduceAll;
279 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
280 }
281
282 if (gblSuccess != 1 && gblErrStrm != nullptr) {
283 *gblErrStrm << "On at least one (MPI) process, the "
284 "Kokkos::DualView has "
285 "either the device or host pointer in the "
286 "DualView equal to null, but the DualView has a nonzero number of "
287 "rows. For more detailed information, please rerun with the "
288 "TPETRA_VERBOSE environment variable set to 1. ";
289 if (verbose) {
290 *gblErrStrm << " Here are error messages from all "
291 "processes:" << endl;
292 if (comm == nullptr) {
293 *gblErrStrm << lclErrStrm.str ();
294 }
295 else {
297 gathervPrint (*gblErrStrm, lclErrStrm.str (), *comm);
298 }
299 }
300 *gblErrStrm << endl;
301 }
302 return gblSuccess == 1;
303}
304
305
306} // namespace Details
307} // namespace Tpetra
308
309#endif // TPETRA_DETAILS_CHECKVIEW_HPP
Declaration of a function that prints strings from each process.
A wrapper around Kokkos::DualView to safely manage data that might be replicated between host and dev...
Nonmember function that computes a residual Computes R = B - A * X.
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.