This file is indexed.

/usr/include/trilinos/Tpetra_Details_PackTriples.hpp is in libtrilinos-tpetra-dev 12.12.1-5.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
// @HEADER
// ***********************************************************************
//
//          Tpetra: Templated Linear Algebra Services Package
//                 Copyright (2008) Sandia Corporation
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ************************************************************************
// @HEADER

#ifndef TPETRA_DETAILS_PACKTRIPLES_HPP
#define TPETRA_DETAILS_PACKTRIPLES_HPP

#include "TpetraCore_config.h"
#include "Teuchos_Comm.hpp"
#ifdef HAVE_TPETRACORE_MPI
#  include "Tpetra_Details_extractMpiCommFromTeuchos.hpp"
#  include "Tpetra_Details_MpiTypeTraits.hpp"
#endif // HAVE_TPETRACORE_MPI
#include <ostream>

namespace Tpetra {
namespace Details {

//
// Search for "SKIP DOWN TO HERE" (omit quotes) for the "public"
// interface.  I put "public" in quotes because it's public only for
// Tpetra developers, NOT for Tpetra users.
//

#ifdef HAVE_TPETRACORE_MPI

// It's not useful to pack or unpack if not using MPI.  Thus, we only
// need to define these with MPI.  For a non-MPI build, we just need
// some outer stub interface that throws an exception.

int
countPackTriplesCountMpi (MPI_Comm comm,
                          int& size,
                          std::ostream* errStrm = NULL);

int
packTriplesCountMpi (const int numEnt,
                     char outBuf[],
                     const int outBufSize,
                     int& outBufCurPos,
                     MPI_Comm comm,
                     std::ostream* errStrm = NULL);

template<class ScalarType, class OrdinalType>
int
countPackTriplesMpi (MPI_Datatype ordinalDt,
                     MPI_Datatype scalarDt,
                     const int numEnt,
                     MPI_Comm comm,
                     int& size,
                     std::ostream* errStrm = NULL)
{
  using std::endl;
  int errCode = MPI_SUCCESS;

  int totalSize = 0; // return value

  // Count the global row and column indices.
  {
    int curSize = 0;
    // We're packing two ordinals, the row resp. column index.
    errCode = MPI_Pack_size (2, ordinalDt, comm, &curSize);
    if (errCode != MPI_SUCCESS) {
      if (errStrm != NULL) {
        *errStrm << "countPackTripleMpi: MPI_Pack_size failed on "
          "ordinalDt call" << endl;
      }
      return errCode;
    }
    totalSize += curSize;
  }
  // Count the matrix value.
  {
    int curSize = 0;
    errCode = MPI_Pack_size (1, scalarDt, comm, &curSize);
    if (errCode != MPI_SUCCESS) {
      if (errStrm != NULL) {
        *errStrm << "countPackTripleMpi: MPI_Pack_size failed on "
          "scalarDt call" << endl;
      }
      return errCode;
    }
    totalSize += curSize; // one Scalar
  }
  totalSize *= numEnt; // all the entries we want to pack

  size = totalSize; // "commit" the result
  return errCode;
}

template<class ScalarType, class OrdinalType>
int
packTripleMpi (const OrdinalType gblRowInd,
               const OrdinalType gblColInd,
               MPI_Datatype ordinalDt,
               const ScalarType& val,
               MPI_Datatype scalarDt,
               char outBuf[],
               const int outBufSize,
               int& outBufCurPos,
               MPI_Comm comm,
               std::ostream* errStrm = NULL)
{
  using std::endl;
  int errCode = MPI_SUCCESS;

  // Combine the two indices into a single MPI_Pack call.
  OrdinalType inds[2] = {gblRowInd, gblColInd};
  errCode = MPI_Pack (inds, 2, ordinalDt,
                      outBuf, outBufSize, &outBufCurPos, comm);
  if (errCode != MPI_SUCCESS) {
    if (errStrm != NULL) {
      *errStrm << "packTripleMpi: MPI_Pack failed for indices i=" << gblRowInd
               << ", j=" << gblColInd << ", where outBufSize=" << outBufSize
               << " and outBufCurPos=" << outBufCurPos << "." << endl;
    }
    return errCode;
  }
  // mfh 17,20 Jan 2017: Some (generally older) MPI implementations
  // want the first argument to be a pointer to nonconst, even though
  // MPI_Pack does not modify that argument.
  errCode = MPI_Pack (const_cast<ScalarType*> (&val), 1, scalarDt,
                      outBuf, outBufSize, &outBufCurPos, comm);
  if (errCode != MPI_SUCCESS) {
    if (errStrm != NULL) {
      *errStrm << "packTripleMpi: MPI_Pack failed on val = " << val
               << endl;
    }
    return errCode;
  }
  return errCode;
}

template<class ScalarType, class OrdinalType>
int
packTriplesMpi (const OrdinalType gblRowInds[],
                const OrdinalType gblColInds[],
                MPI_Datatype ordinalDt,
                const ScalarType vals[],
                MPI_Datatype scalarDt,
                const int numEnt,
                char outBuf[],
                const int outBufSize,
                int& outBufCurPos,
                MPI_Comm comm,
                std::ostream* errStrm = NULL)
{
  using std::endl;
  int errCode = MPI_SUCCESS;

  for (int k = 0; k < numEnt; ++k) {
    errCode = packTripleMpi (gblRowInds[k], gblColInds[k], ordinalDt,
                             vals[k], scalarDt, outBuf, outBufSize,
                             outBufCurPos, comm, errStrm);
    if (errCode != MPI_SUCCESS) {
      if (errStrm != NULL) {
        *errStrm << "packTriplesMpi: packTripleMpi failed at entry k=" << k
                 << endl;
      }
      return errCode;
    }
  }
  return errCode;
}

template<class ScalarType, class OrdinalType>
int
unpackTripleMpi (const char inBuf[],
                 const int inBufSize,
                 int& inBufCurPos,
                 OrdinalType& gblRowInd,
                 OrdinalType& gblColInd,
                 MPI_Datatype ordinalDt,
                 ScalarType& val,
                 MPI_Datatype scalarDt,
                 MPI_Comm comm,
                 std::ostream* errStrm = NULL)
{
  using std::endl;
  int errCode = MPI_SUCCESS;

  // mfh 17 Jan 2017: Some (generally older) MPI implementations want
  // the input buffer argument to be a pointer to nonconst, even
  // though MPI_Unpack does not modify that argument.

  // Combine the two indices into a single MPI_Pack call.
  OrdinalType inds[2] = {static_cast<OrdinalType> (0),
                         static_cast<OrdinalType> (0)};
  errCode = MPI_Unpack (const_cast<char*> (inBuf), inBufSize, &inBufCurPos,
                        inds, 2, ordinalDt, comm);
  if (errCode != MPI_SUCCESS) {
    if (errStrm != NULL) {
      *errStrm << "unpackTripleMpi: MPI_Unpack failed when unpacking indices: "
        "inBufSize=" << inBufSize << ", inBufCurPos=" << inBufCurPos << endl;
    }
    return errCode;
  }
  gblRowInd = inds[0];
  gblColInd = inds[1];

  // mfh 17 Jan 2017: Some (generally older) MPI implementations want
  // the input buffer argument to be a pointer to nonconst, even
  // though MPI_Unpack does not modify that argument.
  errCode = MPI_Unpack (const_cast<char*> (inBuf), inBufSize, &inBufCurPos,
                        &val, 1, scalarDt, comm);
  if (errCode != MPI_SUCCESS) {
    if (errStrm != NULL) {
      *errStrm << "unpackTripleMpi: MPI_Unpack failed when unpacking value: "
        "inBufSize=" << inBufSize << ", inBufCurPos=" << inBufCurPos << endl;
    }
    return errCode;
  }
  return errCode;
}

int
unpackTriplesCountMpi (const char inBuf[],
                       const int inBufSize,
                       int& inBufCurPos,
                       int& numEnt,
                       MPI_Comm comm,
                       std::ostream* errStrm = NULL);

template<class ScalarType, class OrdinalType>
int
unpackTriplesMpi (const char inBuf[],
                  const int inBufSize,
                  int& inBufCurPos,
                  OrdinalType gblRowInds[],
                  OrdinalType gblColInds[],
                  MPI_Datatype ordinalDt,
                  ScalarType vals[],
                  MPI_Datatype scalarDt,
                  const int numEnt, // input arg, from unpackTriplesCountMpi
                  MPI_Comm comm,
                  std::ostream* errStrm = NULL)
{
  using std::endl;
  int errCode = MPI_SUCCESS;

  for (int k = 0; k < numEnt; ++k) {
    errCode = unpackTripleMpi (inBuf, inBufSize, inBufCurPos,
                               gblRowInds[k], gblColInds[k], ordinalDt,
                               vals[k], scalarDt, comm, errStrm);
    if (errCode != MPI_SUCCESS) {
      if (errStrm != NULL) {
        *errStrm << "unpackTriplesMpi: packTripleMpi failed at entry k=" << k
                 << ": inBufSize=" << inBufSize
                 << ", inBufCurPos=" << inBufCurPos
                 << "." << endl;
      }
      return errCode;
    }
  }
  return errCode;
}

#endif // HAVE_TPETRACORE_MPI

//
// SKIP DOWN TO HERE FOR "PUBLIC" INTERFACE
//

/// \brief Compute the buffer size required by packTriples for packing
///   the number of matrix entries ("triples").
///
/// countPackTriples tells me an upper bound on how much buffer space
/// I need to hold numEnt triples.  packTriplesCount actually packs
/// numEnt, the number of triples.  countPackTriplesCount tells me an
/// upper bound on how much buffer space I need to hold the number of
/// triples, not the triples themselves.
///
/// \param comm [in] Communicator used in sending and receiving the
///   packed entries.  (MPI wants this, so we have to include it.).
/// \param size [out] Pack buffer size in bytes (sizeof(char)).
/// \param errStrm [out] If nonnull, print any error messages to this
///   stream, else don't print error messages.
///
/// \return Error code.  MPI_SUCCESS (0) if successful, else nonzero.
///
/// \warning It only makes sense to call this function if using MPI.
///   If <i>not</i> building with MPI, this function is a stub that
///   returns nonzero.
int
countPackTriplesCount (const ::Teuchos::Comm<int>& comm,
                       int& size,
                       std::ostream* errStrm = NULL);

/// \brief Compute the buffer size required by packTriples for packing
///   \c numEnt number of (i,j,A(i,j)) matrix entries ("triples").
///
/// This function is NOT the same thing as packTriplesCount.
/// countPackTriples tells me an upper bound on how much buffer space
/// I need to hold numEnt triples.  packTriplesCount actually packs
/// numEnt, the number of triples.  countPackTriplesCount tells me an
/// upper bound on how much buffer space I need to hold the number of
/// triples, not the triples themselves.
///
/// \tparam ScalarType Type of each matrix entry A(i,j).
/// \tparam OrdinalType Type of each matrix index i or j.
///
/// \param numEnt [in] Number of matrix entries ("triples") to pack.
/// \param comm [in] Communicator used in sending and receiving the
///   packed entries.  (MPI wants this, so we have to include it.).
/// \param size [out] Pack buffer size in bytes (sizeof(char)).
/// \param errStrm [out] If nonnull, print any error messages to this
///   stream, else don't print error messages.
///
/// \return Error code.  MPI_SUCCESS (0) if successful, else nonzero.
///
/// \warning It only makes sense to call this function if using MPI.
///   If <i>not</i> building with MPI, this function is a stub that
///   returns nonzero.
template<class ScalarType, class OrdinalType>
int
countPackTriples (const int numEnt,
                  const ::Teuchos::Comm<int>& comm,
                  int& size, // output argument
                  std::ostream* errStrm = NULL)
{
#ifdef HAVE_TPETRACORE_MPI
  using ::Tpetra::Details::extractMpiCommFromTeuchos;
  using ::Tpetra::Details::MpiTypeTraits;

  static_assert (MpiTypeTraits<ScalarType>::isSpecialized, "countPackTriples: "
                 "ScalarType lacks an MpiTypeTraits specialization.");
  static_assert (MpiTypeTraits<OrdinalType>::isSpecialized, "countPackTriples: "
                 "OrdinalType lacks an MpiTypeTraits specialization.");

  MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
  MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
  MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();

  const int errCode =
    countPackTriplesMpi<ScalarType, OrdinalType> (ordinalDt, scalarDt,
                                                  numEnt, mpiComm,
                                                  size, errStrm);
  if (MpiTypeTraits<ScalarType>::needsFree) {
    (void) MPI_Type_free (&scalarDt);
  }
  if (MpiTypeTraits<OrdinalType>::needsFree) {
    (void) MPI_Type_free (&ordinalDt);
  }
  return errCode;

#else // NOT HAVE_TPETRACORE_MPI
  if (errStrm != NULL) {
    *errStrm << "countPackTriples: Not implemented (no need; there's no need "
      "to pack or unpack anything if there's only one process)." << std::endl;
  }
  return -1;
#endif // HAVE_TPETRACORE_MPI
}

/// \brief Pack the count (number) of matrix triples.
///
/// This function is NOT the same thing as countPackTriples.
/// countPackTriples tells me an upper bound on how much buffer space
/// I need to hold numEnt triples.  packTriplesCount actually packs
/// numEnt, the number of triples.  countPackTriplesCount tells me an
/// upper bound on how much buffer space I need to hold the number of
/// triples, not the triples themselves.
///
/// \param numEnt [in] Number of matrix entries ("triples") to pack.
/// \param outBuf [out] Output buffer.
/// \param outBufSize [out] Total output buffer size in bytes.
/// \param outBufCurPos [in/out] Current position from which to start
///   writing to the output buffer.  This corresponds to the
///   'position' in/out argument of MPI_Pack.
/// \param comm [in] Communicator used in sending and receiving the
///   packed entries.  (MPI wants this, so we have to include it.).
/// \param errStrm [out] If nonnull, print any error messages to this
///   stream, else don't print error messages.
///
/// \return Error code.  MPI_SUCCESS (0) if successful, else nonzero.
///
/// \warning It only makes sense to call this function if using MPI.
///   If <i>not</i> building with MPI, this function is a stub that
///   returns nonzero.
int
packTriplesCount (const int numEnt,
                  char outBuf[],
                  const int outBufSize,
                  int& outBufCurPos,
                  const ::Teuchos::Comm<int>& comm,
                  std::ostream* errStrm = NULL);

/// \brief Pack matrix entries ("triples" (i, j, A(i,j))) into the
///   given output buffer.
///
/// \tparam ScalarType Type of each matrix entry A(i,j).
/// \tparam OrdinalType Type of each matrix index i or j.
///
/// \param gblRowInds [in] Row indices to pack.
/// \param gblColInds [in] Column indices to pack.
/// \param val [in] Matrix values A(i,j) to pack.
/// \param numEnt [in] Number of matrix entries ("triples") to pack.
/// \param outBuf [out] Output buffer.
/// \param outBufSize [out] Total output buffer size in bytes.
/// \param outBufCurPos [in/out] Current position from which to start
///   writing to the output buffer.  This corresponds to the
///   'position' in/out argument of MPI_Pack.
/// \param comm [in] Communicator used in sending and receiving the
///   packed entries.  (MPI wants this, so we have to include it.).
/// \param errStrm [out] If nonnull, print any error messages to this
///   stream, else don't print error messages.
///
/// \return Error code.  MPI_SUCCESS (0) if successful, else nonzero.
///
/// \warning It only makes sense to call this function if using MPI.
///   If <i>not</i> building with MPI, this function is a stub that
///   returns nonzero.
template<class ScalarType, class OrdinalType>
int
#ifdef HAVE_TPETRACORE_MPI
packTriples (const OrdinalType gblRowInds[],
             const OrdinalType gblColInds[],
             const ScalarType vals[],
             const int numEnt,
             char outBuf[],
             const int outBufSize,
             int& outBufCurPos,
             const ::Teuchos::Comm<int>& comm,
             std::ostream* errStrm = NULL)
#else // NOT HAVE_TPETRACORE_MPI
packTriples (const OrdinalType /* gblRowInds */ [],
             const OrdinalType /* gblColInds */ [],
             const ScalarType /* vals */ [],
             const int /* numEnt */,
             char /* outBuf */ [],
             const int /* outBufSize */,
             int& /* outBufCurPos */,
             const ::Teuchos::Comm<int>& /* comm */,
             std::ostream* errStrm = NULL)
#endif // HAVE_TPETRACORE_MPI
{
#ifdef HAVE_TPETRACORE_MPI
  using ::Tpetra::Details::extractMpiCommFromTeuchos;
  using ::Tpetra::Details::MpiTypeTraits;

  static_assert (MpiTypeTraits<ScalarType>::isSpecialized, "packTriples: "
                 "ScalarType lacks an MpiTypeTraits specialization.");
  static_assert (MpiTypeTraits<OrdinalType>::isSpecialized, "packTriples: "
                 "OrdinalType lacks an MpiTypeTraits specialization.");

  MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
  MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
  MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();

  const int errCode =
    packTriplesMpi (gblRowInds, gblColInds, ordinalDt, vals, scalarDt,
                    numEnt, outBuf, outBufSize, outBufCurPos, mpiComm, errStrm);
  if (MpiTypeTraits<ScalarType>::needsFree) {
    (void) MPI_Type_free (&scalarDt);
  }
  if (MpiTypeTraits<OrdinalType>::needsFree) {
    (void) MPI_Type_free (&ordinalDt);
  }
  return errCode;

#else // NOT HAVE_TPETRACORE_MPI
  if (errStrm != NULL) {
    *errStrm << "packTriples: Not implemented (no need; there's no need to "
      "pack or unpack anything if there's only one process)." << std::endl;
  }
  return -1;
#endif // HAVE_TPETRACORE_MPI
}

/// \brief Unpack just the count of triples from the given input
///   buffer.
///
/// We store the count of triples as an \c int, because MPI buffer
/// sizes are \c int.
///
/// \param inBuf [in] Input buffer.
/// \param inBufSize [out] Total input buffer size in bytes.
/// \param inBufCurPos [in/out] Current position from which to start
///   reading from the input buffer.  This corresponds to the
///   'position' in/out argument of MPI_Unpack.
/// \param numEnt [out] Number of matrix entries ("triples") that were
///   packed.
/// \param comm [in] Communicator used in sending and receiving the
///   packed entries.  (MPI wants this, so we have to include it.).
/// \param errStrm [out] If nonnull, print any error messages to this
///   stream, else don't print error messages.
///
/// \return Error code.  MPI_SUCCESS (0) if successful, else nonzero.
///
/// \warning It only makes sense to call this function if using MPI.
///   If <i>not</i> building with MPI, this function is a stub that
///   returns nonzero.
int
unpackTriplesCount (const char inBuf[],
                    const int inBufSize,
                    int& inBufCurPos,
                    int& numEnt, // output argument!
                    const ::Teuchos::Comm<int>& comm,
                    std::ostream* errStrm = NULL);

/// \brief Unpack matrix entries ("triples" (i, j, A(i,j))) from the
///   given input buffer.
///
/// \tparam ScalarType Type of each matrix entry A(i,j).
/// \tparam OrdinalType Type of each matrix index i or j.
///
/// \param inBuf [in] Input buffer.
/// \param inBufSize [out] Total pack buffer size in bytes (sizeof(char)).
/// \param inBufCurPos [in/out] Current position from which to start
///   reading from the input buffer.  This corresponds to the
///   'position' in/out argument of MPI_Unpack.
/// \param gblRowInds [out] Row indices unpacked.
/// \param gblColInds [out] Column indices unpacked.
/// \param val [out] Matrix values A(i,j) unpacked.
/// \param numEnt [in] Number of matrix entries ("triples") to unpack.
///   If you don't know it, then you should have senders pack the
///   triples count as the first thing in the buffer, and unpack it
///   first via unpackTriplesCount().
/// \param comm [in] Communicator used in sending and receiving the
///   packed entries.  (MPI wants this, so we have to include it.).
/// \param errStrm [out] If nonnull, print any error messages to this
///   stream, else don't print error messages.
///
/// \return Error code.  MPI_SUCCESS (0) if successful, else nonzero.
///
/// \warning It only makes sense to call this function if using MPI.
///   If <i>not</i> building with MPI, this function is a stub that
///   returns nonzero.
template<class ScalarType, class OrdinalType>
int
#ifdef HAVE_TPETRACORE_MPI
unpackTriples (const char inBuf[],
               const int inBufSize,
               int& inBufCurPos,
               OrdinalType gblRowInds[],
               OrdinalType gblColInds[],
               ScalarType vals[],
               const int numEnt,
               const ::Teuchos::Comm<int>& comm,
               std::ostream* errStrm = NULL)
#else // NOT HAVE_TPETRACORE_MPI
unpackTriples (const char /* inBuf */ [],
               const int /* inBufSize */,
               int& /* inBufCurPos */,
               OrdinalType /* gblRowInds */ [],
               OrdinalType /* gblColInds */ [],
               ScalarType /* vals */ [],
               const int /* numEnt */,
               const ::Teuchos::Comm<int>& /* comm */,
               std::ostream* errStrm = NULL)
#endif // HAVE_TPETRACORE_MPI
{
#ifdef HAVE_TPETRACORE_MPI
  using ::Tpetra::Details::extractMpiCommFromTeuchos;
  using ::Tpetra::Details::MpiTypeTraits;

  static_assert (MpiTypeTraits<ScalarType>::isSpecialized, "unpackTriples: "
                 "ScalarType lacks an MpiTypeTraits specialization.");
  static_assert (MpiTypeTraits<OrdinalType>::isSpecialized, "unpackTriples: "
                 "OrdinalType lacks an MpiTypeTraits specialization.");

  MPI_Comm mpiComm = extractMpiCommFromTeuchos (comm);
  MPI_Datatype ordinalDt = MpiTypeTraits<OrdinalType>::getType ();
  MPI_Datatype scalarDt = MpiTypeTraits<ScalarType>::getType ();

  const int errCode =
    unpackTriplesMpi (inBuf, inBufSize, inBufCurPos,
                      gblRowInds, gblColInds, ordinalDt,
                      vals, scalarDt, numEnt, mpiComm, errStrm);
  if (MpiTypeTraits<ScalarType>::needsFree) {
    (void) MPI_Type_free (&scalarDt);
  }
  if (MpiTypeTraits<OrdinalType>::needsFree) {
    (void) MPI_Type_free (&ordinalDt);
  }
  return errCode;

#else // NOT HAVE_TPETRACORE_MPI
  if (errStrm != NULL) {
    *errStrm << "unpackTriples: Not implemented (no need; there's no need to "
      "pack or unpack anything if there's only one process)." << std::endl;
  }
  return -1;
#endif // HAVE_TPETRACORE_MPI
}

} // namespace Details
} // namespace Tpetra

#endif // TPETRA_DETAILS_PACKTRIPLES_HPP