10 #ifndef TPETRA_DETAILS_PACKTRIPLES_HPP
11 #define TPETRA_DETAILS_PACKTRIPLES_HPP
13 #include "TpetraCore_config.h"
14 #include "Teuchos_Comm.hpp"
15 #ifdef HAVE_TPETRACORE_MPI
18 #endif // HAVE_TPETRACORE_MPI
30 #ifdef HAVE_TPETRACORE_MPI
37 countPackTriplesCountMpi (MPI_Comm comm,
39 std::ostream* errStrm = NULL);
42 packTriplesCountMpi (
const int numEnt,
47 std::ostream* errStrm = NULL);
49 template<
class ScalarType,
class OrdinalType>
51 countPackTriplesMpi (MPI_Datatype ordinalDt,
52 MPI_Datatype scalarDt,
56 std::ostream* errStrm = NULL)
59 int errCode = MPI_SUCCESS;
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;
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;
96 template<
class ScalarType,
class OrdinalType>
98 packTripleMpi (
const OrdinalType gblRowInd,
99 const OrdinalType gblColInd,
100 MPI_Datatype ordinalDt,
101 const ScalarType& val,
102 MPI_Datatype scalarDt,
104 const int outBufSize,
107 std::ostream* errStrm = NULL)
110 int errCode = MPI_SUCCESS;
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;
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
139 template<
class ScalarType,
class OrdinalType>
141 packTriplesMpi (
const OrdinalType gblRowInds[],
142 const OrdinalType gblColInds[],
143 MPI_Datatype ordinalDt,
144 const ScalarType vals[],
145 MPI_Datatype scalarDt,
148 const int outBufSize,
151 std::ostream* errStrm = NULL)
154 int errCode = MPI_SUCCESS;
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
171 template<
class ScalarType,
class OrdinalType>
173 unpackTripleMpi (
const char inBuf[],
176 OrdinalType& gblRowInd,
177 OrdinalType& gblColInd,
178 MPI_Datatype ordinalDt,
180 MPI_Datatype scalarDt,
182 std::ostream* errStrm = NULL)
185 int errCode = MPI_SUCCESS;
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;
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;
222 unpackTriplesCountMpi (
const char inBuf[],
227 std::ostream* errStrm = NULL);
229 template<
class ScalarType,
class OrdinalType>
231 unpackTriplesMpi (
const char inBuf[],
234 OrdinalType gblRowInds[],
235 OrdinalType gblColInds[],
236 MPI_Datatype ordinalDt,
238 MPI_Datatype scalarDt,
241 std::ostream* errStrm = NULL)
244 int errCode = MPI_SUCCESS;
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
263 #endif // HAVE_TPETRACORE_MPI
292 std::ostream* errStrm = NULL);
319 template<
class ScalarType,
class OrdinalType>
322 const ::Teuchos::Comm<int>& comm,
324 std::ostream* errStrm = NULL)
326 #ifdef HAVE_TPETRACORE_MPI
327 using ::Tpetra::Details::extractMpiCommFromTeuchos;
328 using ::Tpetra::Details::MpiTypeTraits;
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.");
335 MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
336 MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
337 MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();
340 countPackTriplesMpi<ScalarType, OrdinalType> (ordinalDt, scalarDt,
343 if (MpiTypeTraits<ScalarType>::needsFree) {
344 (void) MPI_Type_free (&scalarDt);
346 if (MpiTypeTraits<OrdinalType>::needsFree) {
347 (void) MPI_Type_free (&ordinalDt);
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;
357 #endif // HAVE_TPETRACORE_MPI
388 const int outBufSize,
390 const ::Teuchos::Comm<int>& comm,
391 std::ostream* errStrm = NULL);
418 template<
class ScalarType,
class OrdinalType>
420 #ifdef HAVE_TPETRACORE_MPI
422 const OrdinalType gblColInds[],
423 const ScalarType vals[],
426 const int outBufSize,
428 const ::Teuchos::Comm<int>& comm,
429 std::ostream* errStrm = NULL)
430 #else // NOT HAVE_TPETRACORE_MPI
432 const OrdinalType [],
438 const ::Teuchos::Comm<int>& ,
439 std::ostream* errStrm = NULL)
440 #endif // HAVE_TPETRACORE_MPI
442 #ifdef HAVE_TPETRACORE_MPI
443 using ::Tpetra::Details::extractMpiCommFromTeuchos;
444 using ::Tpetra::Details::MpiTypeTraits;
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.");
451 MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
452 MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
453 MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();
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);
461 if (MpiTypeTraits<OrdinalType>::needsFree) {
462 (void) MPI_Type_free (&ordinalDt);
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;
472 #endif // HAVE_TPETRACORE_MPI
503 const ::Teuchos::Comm<int>& comm,
504 std::ostream* errStrm = NULL);
534 template<
class ScalarType,
class OrdinalType>
536 #ifdef HAVE_TPETRACORE_MPI
540 OrdinalType gblRowInds[],
541 OrdinalType gblColInds[],
544 const ::Teuchos::Comm<int>& comm,
545 std::ostream* errStrm = NULL)
546 #else // NOT HAVE_TPETRACORE_MPI
554 const ::Teuchos::Comm<int>& ,
555 std::ostream* errStrm = NULL)
556 #endif // HAVE_TPETRACORE_MPI
558 #ifdef HAVE_TPETRACORE_MPI
559 using ::Tpetra::Details::extractMpiCommFromTeuchos;
560 using ::Tpetra::Details::MpiTypeTraits;
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.");
567 MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
568 MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
569 MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();
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);
578 if (MpiTypeTraits<OrdinalType>::needsFree) {
579 (void) MPI_Type_free (&ordinalDt);
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;
589 #endif // HAVE_TPETRACORE_MPI
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.
Add specializations of Teuchos::Details::MpiTypeTraits for Kokkos::complex<float> 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 ("triples")...
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 ("triples").
int unpackTriples(const char[], const int, int &, OrdinalType[], OrdinalType[], ScalarType[], const int, const ::Teuchos::Comm< int > &, std::ostream *errStrm=NULL)
Unpack matrix entries ("triples" (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 ("triples" (i, j, A(i,j))) into the given output buffer.