Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_PackTriples.hpp
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_PACKTRIPLES_HPP
11 #define TPETRA_DETAILS_PACKTRIPLES_HPP
12 
13 #include "TpetraCore_config.h"
14 #include "Teuchos_Comm.hpp"
15 #ifdef HAVE_TPETRACORE_MPI
18 #endif // HAVE_TPETRACORE_MPI
19 #include <ostream>
20 
21 namespace Tpetra {
22 namespace Details {
23 
24 //
25 // Search for "SKIP DOWN TO HERE" (omit quotes) for the "public"
26 // interface. I put "public" in quotes because it's public only for
27 // Tpetra developers, NOT for Tpetra users.
28 //
29 
30 #ifdef HAVE_TPETRACORE_MPI
31 
32 // It's not useful to pack or unpack if not using MPI. Thus, we only
33 // need to define these with MPI. For a non-MPI build, we just need
34 // some outer stub interface that throws an exception.
35 
36 int
37 countPackTriplesCountMpi (MPI_Comm comm,
38  int& size,
39  std::ostream* errStrm = NULL);
40 
41 int
42 packTriplesCountMpi (const int numEnt,
43  char outBuf[],
44  const int outBufSize,
45  int& outBufCurPos,
46  MPI_Comm comm,
47  std::ostream* errStrm = NULL);
48 
49 template<class ScalarType, class OrdinalType>
50 int
51 countPackTriplesMpi (MPI_Datatype ordinalDt,
52  MPI_Datatype scalarDt,
53  const int numEnt,
54  MPI_Comm comm,
55  int& size,
56  std::ostream* errStrm = NULL)
57 {
58  using std::endl;
59  int errCode = MPI_SUCCESS;
60 
61  int totalSize = 0; // return value
62 
63  // Count the global row and column indices.
64  {
65  int curSize = 0;
66  // We're packing two ordinals, the row resp. column index.
67  errCode = MPI_Pack_size (2, ordinalDt, comm, &curSize);
68  if (errCode != MPI_SUCCESS) {
69  if (errStrm != NULL) {
70  *errStrm << "countPackTripleMpi: MPI_Pack_size failed on "
71  "ordinalDt call" << endl;
72  }
73  return errCode;
74  }
75  totalSize += curSize;
76  }
77  // Count the matrix value.
78  {
79  int curSize = 0;
80  errCode = MPI_Pack_size (1, scalarDt, comm, &curSize);
81  if (errCode != MPI_SUCCESS) {
82  if (errStrm != NULL) {
83  *errStrm << "countPackTripleMpi: MPI_Pack_size failed on "
84  "scalarDt call" << endl;
85  }
86  return errCode;
87  }
88  totalSize += curSize; // one Scalar
89  }
90  totalSize *= numEnt; // all the entries we want to pack
91 
92  size = totalSize; // "commit" the result
93  return errCode;
94 }
95 
96 template<class ScalarType, class OrdinalType>
97 int
98 packTripleMpi (const OrdinalType gblRowInd,
99  const OrdinalType gblColInd,
100  MPI_Datatype ordinalDt,
101  const ScalarType& val,
102  MPI_Datatype scalarDt,
103  char outBuf[],
104  const int outBufSize,
105  int& outBufCurPos,
106  MPI_Comm comm,
107  std::ostream* errStrm = NULL)
108 {
109  using std::endl;
110  int errCode = MPI_SUCCESS;
111 
112  // Combine the two indices into a single MPI_Pack call.
113  OrdinalType inds[2] = {gblRowInd, gblColInd};
114  errCode = MPI_Pack (inds, 2, ordinalDt,
115  outBuf, outBufSize, &outBufCurPos, comm);
116  if (errCode != MPI_SUCCESS) {
117  if (errStrm != NULL) {
118  *errStrm << "packTripleMpi: MPI_Pack failed for indices i=" << gblRowInd
119  << ", j=" << gblColInd << ", where outBufSize=" << outBufSize
120  << " and outBufCurPos=" << outBufCurPos << "." << endl;
121  }
122  return errCode;
123  }
124  // mfh 17,20 Jan 2017: Some (generally older) MPI implementations
125  // want the first argument to be a pointer to nonconst, even though
126  // MPI_Pack does not modify that argument.
127  errCode = MPI_Pack (const_cast<ScalarType*> (&val), 1, scalarDt,
128  outBuf, outBufSize, &outBufCurPos, comm);
129  if (errCode != MPI_SUCCESS) {
130  if (errStrm != NULL) {
131  *errStrm << "packTripleMpi: MPI_Pack failed on val = " << val
132  << endl;
133  }
134  return errCode;
135  }
136  return errCode;
137 }
138 
139 template<class ScalarType, class OrdinalType>
140 int
141 packTriplesMpi (const OrdinalType gblRowInds[],
142  const OrdinalType gblColInds[],
143  MPI_Datatype ordinalDt,
144  const ScalarType vals[],
145  MPI_Datatype scalarDt,
146  const int numEnt,
147  char outBuf[],
148  const int outBufSize,
149  int& outBufCurPos,
150  MPI_Comm comm,
151  std::ostream* errStrm = NULL)
152 {
153  using std::endl;
154  int errCode = MPI_SUCCESS;
155 
156  for (int k = 0; k < numEnt; ++k) {
157  errCode = packTripleMpi (gblRowInds[k], gblColInds[k], ordinalDt,
158  vals[k], scalarDt, outBuf, outBufSize,
159  outBufCurPos, comm, errStrm);
160  if (errCode != MPI_SUCCESS) {
161  if (errStrm != NULL) {
162  *errStrm << "packTriplesMpi: packTripleMpi failed at entry k=" << k
163  << endl;
164  }
165  return errCode;
166  }
167  }
168  return errCode;
169 }
170 
171 template<class ScalarType, class OrdinalType>
172 int
173 unpackTripleMpi (const char inBuf[],
174  const int inBufSize,
175  int& inBufCurPos,
176  OrdinalType& gblRowInd,
177  OrdinalType& gblColInd,
178  MPI_Datatype ordinalDt,
179  ScalarType& val,
180  MPI_Datatype scalarDt,
181  MPI_Comm comm,
182  std::ostream* errStrm = NULL)
183 {
184  using std::endl;
185  int errCode = MPI_SUCCESS;
186 
187  // mfh 17 Jan 2017: Some (generally older) MPI implementations want
188  // the input buffer argument to be a pointer to nonconst, even
189  // though MPI_Unpack does not modify that argument.
190 
191  // Combine the two indices into a single MPI_Pack call.
192  OrdinalType inds[2] = {static_cast<OrdinalType> (0),
193  static_cast<OrdinalType> (0)};
194  errCode = MPI_Unpack (const_cast<char*> (inBuf), inBufSize, &inBufCurPos,
195  inds, 2, ordinalDt, comm);
196  if (errCode != MPI_SUCCESS) {
197  if (errStrm != NULL) {
198  *errStrm << "unpackTripleMpi: MPI_Unpack failed when unpacking indices: "
199  "inBufSize=" << inBufSize << ", inBufCurPos=" << inBufCurPos << endl;
200  }
201  return errCode;
202  }
203  gblRowInd = inds[0];
204  gblColInd = inds[1];
205 
206  // mfh 17 Jan 2017: Some (generally older) MPI implementations want
207  // the input buffer argument to be a pointer to nonconst, even
208  // though MPI_Unpack does not modify that argument.
209  errCode = MPI_Unpack (const_cast<char*> (inBuf), inBufSize, &inBufCurPos,
210  &val, 1, scalarDt, comm);
211  if (errCode != MPI_SUCCESS) {
212  if (errStrm != NULL) {
213  *errStrm << "unpackTripleMpi: MPI_Unpack failed when unpacking value: "
214  "inBufSize=" << inBufSize << ", inBufCurPos=" << inBufCurPos << endl;
215  }
216  return errCode;
217  }
218  return errCode;
219 }
220 
221 int
222 unpackTriplesCountMpi (const char inBuf[],
223  const int inBufSize,
224  int& inBufCurPos,
225  int& numEnt,
226  MPI_Comm comm,
227  std::ostream* errStrm = NULL);
228 
229 template<class ScalarType, class OrdinalType>
230 int
231 unpackTriplesMpi (const char inBuf[],
232  const int inBufSize,
233  int& inBufCurPos,
234  OrdinalType gblRowInds[],
235  OrdinalType gblColInds[],
236  MPI_Datatype ordinalDt,
237  ScalarType vals[],
238  MPI_Datatype scalarDt,
239  const int numEnt, // input arg, from unpackTriplesCountMpi
240  MPI_Comm comm,
241  std::ostream* errStrm = NULL)
242 {
243  using std::endl;
244  int errCode = MPI_SUCCESS;
245 
246  for (int k = 0; k < numEnt; ++k) {
247  errCode = unpackTripleMpi (inBuf, inBufSize, inBufCurPos,
248  gblRowInds[k], gblColInds[k], ordinalDt,
249  vals[k], scalarDt, comm, errStrm);
250  if (errCode != MPI_SUCCESS) {
251  if (errStrm != NULL) {
252  *errStrm << "unpackTriplesMpi: packTripleMpi failed at entry k=" << k
253  << ": inBufSize=" << inBufSize
254  << ", inBufCurPos=" << inBufCurPos
255  << "." << endl;
256  }
257  return errCode;
258  }
259  }
260  return errCode;
261 }
262 
263 #endif // HAVE_TPETRACORE_MPI
264 
265 //
266 // SKIP DOWN TO HERE FOR "PUBLIC" INTERFACE
267 //
268 
289 int
290 countPackTriplesCount (const ::Teuchos::Comm<int>& comm,
291  int& size,
292  std::ostream* errStrm = NULL);
293 
319 template<class ScalarType, class OrdinalType>
320 int
321 countPackTriples (const int numEnt,
322  const ::Teuchos::Comm<int>& comm,
323  int& size, // output argument
324  std::ostream* errStrm = NULL)
325 {
326 #ifdef HAVE_TPETRACORE_MPI
327  using ::Tpetra::Details::extractMpiCommFromTeuchos;
328  using ::Tpetra::Details::MpiTypeTraits;
329 
330  static_assert (MpiTypeTraits<ScalarType>::isSpecialized, "countPackTriples: "
331  "ScalarType lacks an MpiTypeTraits specialization.");
332  static_assert (MpiTypeTraits<OrdinalType>::isSpecialized, "countPackTriples: "
333  "OrdinalType lacks an MpiTypeTraits specialization.");
334 
335  MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
336  MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
337  MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();
338 
339  const int errCode =
340  countPackTriplesMpi<ScalarType, OrdinalType> (ordinalDt, scalarDt,
341  numEnt, mpiComm,
342  size, errStrm);
343  if (MpiTypeTraits<ScalarType>::needsFree) {
344  (void) MPI_Type_free (&scalarDt);
345  }
346  if (MpiTypeTraits<OrdinalType>::needsFree) {
347  (void) MPI_Type_free (&ordinalDt);
348  }
349  return errCode;
350 
351 #else // NOT HAVE_TPETRACORE_MPI
352  if (errStrm != NULL) {
353  *errStrm << "countPackTriples: Not implemented (no need; there's no need "
354  "to pack or unpack anything if there's only one process)." << std::endl;
355  }
356  return -1;
357 #endif // HAVE_TPETRACORE_MPI
358 }
359 
385 int
386 packTriplesCount (const int numEnt,
387  char outBuf[],
388  const int outBufSize,
389  int& outBufCurPos,
390  const ::Teuchos::Comm<int>& comm,
391  std::ostream* errStrm = NULL);
392 
418 template<class ScalarType, class OrdinalType>
419 int
420 #ifdef HAVE_TPETRACORE_MPI
421 packTriples (const OrdinalType gblRowInds[],
422  const OrdinalType gblColInds[],
423  const ScalarType vals[],
424  const int numEnt,
425  char outBuf[],
426  const int outBufSize,
427  int& outBufCurPos,
428  const ::Teuchos::Comm<int>& comm,
429  std::ostream* errStrm = NULL)
430 #else // NOT HAVE_TPETRACORE_MPI
431 packTriples (const OrdinalType /* gblRowInds */ [],
432  const OrdinalType /* gblColInds */ [],
433  const ScalarType /* vals */ [],
434  const int /* numEnt */,
435  char /* outBuf */ [],
436  const int /* outBufSize */,
437  int& /* outBufCurPos */,
438  const ::Teuchos::Comm<int>& /* comm */,
439  std::ostream* errStrm = NULL)
440 #endif // HAVE_TPETRACORE_MPI
441 {
442 #ifdef HAVE_TPETRACORE_MPI
443  using ::Tpetra::Details::extractMpiCommFromTeuchos;
444  using ::Tpetra::Details::MpiTypeTraits;
445 
446  static_assert (MpiTypeTraits<ScalarType>::isSpecialized, "packTriples: "
447  "ScalarType lacks an MpiTypeTraits specialization.");
448  static_assert (MpiTypeTraits<OrdinalType>::isSpecialized, "packTriples: "
449  "OrdinalType lacks an MpiTypeTraits specialization.");
450 
451  MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
452  MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
453  MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();
454 
455  const int errCode =
456  packTriplesMpi (gblRowInds, gblColInds, ordinalDt, vals, scalarDt,
457  numEnt, outBuf, outBufSize, outBufCurPos, mpiComm, errStrm);
458  if (MpiTypeTraits<ScalarType>::needsFree) {
459  (void) MPI_Type_free (&scalarDt);
460  }
461  if (MpiTypeTraits<OrdinalType>::needsFree) {
462  (void) MPI_Type_free (&ordinalDt);
463  }
464  return errCode;
465 
466 #else // NOT HAVE_TPETRACORE_MPI
467  if (errStrm != NULL) {
468  *errStrm << "packTriples: Not implemented (no need; there's no need to "
469  "pack or unpack anything if there's only one process)." << std::endl;
470  }
471  return -1;
472 #endif // HAVE_TPETRACORE_MPI
473 }
474 
498 int
499 unpackTriplesCount (const char inBuf[],
500  const int inBufSize,
501  int& inBufCurPos,
502  int& numEnt, // output argument!
503  const ::Teuchos::Comm<int>& comm,
504  std::ostream* errStrm = NULL);
505 
534 template<class ScalarType, class OrdinalType>
535 int
536 #ifdef HAVE_TPETRACORE_MPI
537 unpackTriples (const char inBuf[],
538  const int inBufSize,
539  int& inBufCurPos,
540  OrdinalType gblRowInds[],
541  OrdinalType gblColInds[],
542  ScalarType vals[],
543  const int numEnt,
544  const ::Teuchos::Comm<int>& comm,
545  std::ostream* errStrm = NULL)
546 #else // NOT HAVE_TPETRACORE_MPI
547 unpackTriples (const char /* inBuf */ [],
548  const int /* inBufSize */,
549  int& /* inBufCurPos */,
550  OrdinalType /* gblRowInds */ [],
551  OrdinalType /* gblColInds */ [],
552  ScalarType /* vals */ [],
553  const int /* numEnt */,
554  const ::Teuchos::Comm<int>& /* comm */,
555  std::ostream* errStrm = NULL)
556 #endif // HAVE_TPETRACORE_MPI
557 {
558 #ifdef HAVE_TPETRACORE_MPI
559  using ::Tpetra::Details::extractMpiCommFromTeuchos;
560  using ::Tpetra::Details::MpiTypeTraits;
561 
562  static_assert (MpiTypeTraits<ScalarType>::isSpecialized, "unpackTriples: "
563  "ScalarType lacks an MpiTypeTraits specialization.");
564  static_assert (MpiTypeTraits<OrdinalType>::isSpecialized, "unpackTriples: "
565  "OrdinalType lacks an MpiTypeTraits specialization.");
566 
567  MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
568  MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
569  MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();
570 
571  const int errCode =
572  unpackTriplesMpi (inBuf, inBufSize, inBufCurPos,
573  gblRowInds, gblColInds, ordinalDt,
574  vals, scalarDt, numEnt, mpiComm, errStrm);
575  if (MpiTypeTraits<ScalarType>::needsFree) {
576  (void) MPI_Type_free (&scalarDt);
577  }
578  if (MpiTypeTraits<OrdinalType>::needsFree) {
579  (void) MPI_Type_free (&ordinalDt);
580  }
581  return errCode;
582 
583 #else // NOT HAVE_TPETRACORE_MPI
584  if (errStrm != NULL) {
585  *errStrm << "unpackTriples: Not implemented (no need; there's no need to "
586  "pack or unpack anything if there's only one process)." << std::endl;
587  }
588  return -1;
589 #endif // HAVE_TPETRACORE_MPI
590 }
591 
592 } // namespace Details
593 } // namespace Tpetra
594 
595 #endif // TPETRA_DETAILS_PACKTRIPLES_HPP
int packTriplesCount(const int, char[], const int, int &, const ::Teuchos::Comm< int > &, std::ostream *errStrm)
Pack the count (number) of matrix triples.
Declaration of Tpetra::Details::extractMpiCommFromTeuchos.
Add specializations of Teuchos::Details::MpiTypeTraits for Kokkos::complex&lt;float&gt; and Kokkos::complex...
int unpackTriplesCount(const char[], const int, int &, int &, const ::Teuchos::Comm< int > &, std::ostream *errStrm)
Unpack just the count of triples from the given input buffer.
int countPackTriplesCount(const ::Teuchos::Comm< int > &, int &size, std::ostream *errStrm)
Compute the buffer size required by packTriples for packing the number of matrix entries (&quot;triples&quot;)...
int countPackTriples(const int numEnt, const ::Teuchos::Comm< int > &comm, int &size, std::ostream *errStrm=NULL)
Compute the buffer size required by packTriples for packing numEnt number of (i,j,A(i,j)) matrix entries (&quot;triples&quot;).
int unpackTriples(const char[], const int, int &, OrdinalType[], OrdinalType[], ScalarType[], const int, const ::Teuchos::Comm< int > &, std::ostream *errStrm=NULL)
Unpack matrix entries (&quot;triples&quot; (i, j, A(i,j))) from the given input buffer.
int packTriples(const OrdinalType[], const OrdinalType[], const ScalarType[], const int, char[], const int, int &, const ::Teuchos::Comm< int > &, std::ostream *errStrm=NULL)
Pack matrix entries (&quot;triples&quot; (i, j, A(i,j))) into the given output buffer.