Actual source code: mpiaijAssemble.cu

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: #define PETSC_SKIP_COMPLEX
  2: #define PETSC_SKIP_SPINLOCK

  4: #include <petscconf.h>
  5:  #include <../src/mat/impls/aij/seq/aij.h>
  6:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  7:  #include <petscbt.h>
  8:  #include <../src/vec/vec/impls/dvecimpl.h>
  9:  #include <petsc/private/vecimpl.h>
 10: #undef VecType
 11:  #include <../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h>

 13: #include <thrust/reduce.h>
 14: #include <thrust/inner_product.h>

 16: #include <cusp/array1d.h>
 17: #include <cusp/print.h>
 18: #include <cusp/coo_matrix.h>

 20: #include <cusp/io/matrix_market.h>

 22: #include <thrust/iterator/counting_iterator.h>
 23: #include <thrust/iterator/transform_iterator.h>
 24: #include <thrust/iterator/permutation_iterator.h>
 25: #include <thrust/functional.h>
 26: #include <thrust/partition.h>
 27: #include <thrust/remove.h>

 29: // this example illustrates how to make repeated access to a range of values
 30: // examples:
 31: //   repeated_range([0, 1, 2, 3], 1) -> [0, 1, 2, 3]
 32: //   repeated_range([0, 1, 2, 3], 2) -> [0, 0, 1, 1, 2, 2, 3, 3]
 33: //   repeated_range([0, 1, 2, 3], 3) -> [0, 0, 0, 1, 1, 1, 2, 2, 2, 3, 3, 3]
 34: //   ...

 36: template <typename Iterator>
 37: class repeated_range
 38: {
 39: public:

 41:   typedef typename thrust::iterator_difference<Iterator>::type difference_type;

 43:   struct repeat_functor : public thrust::unary_function<difference_type,difference_type>
 44:   {
 45:     difference_type repeats;

 47:     repeat_functor(difference_type repeats) : repeats(repeats) {}

 49:     __host__ __device__
 50:     difference_type operator()(const difference_type &i) const
 51:     {
 52:       return i / repeats;
 53:     }
 54:   };

 56:   typedef typename thrust::counting_iterator<difference_type>                   CountingIterator;
 57:   typedef typename thrust::transform_iterator<repeat_functor, CountingIterator> TransformIterator;
 58:   typedef typename thrust::permutation_iterator<Iterator,TransformIterator>     PermutationIterator;

 60:   // type of the repeated_range iterator
 61:   typedef PermutationIterator iterator;

 63:   // construct repeated_range for the range [first,last)
 64:   repeated_range(Iterator first, Iterator last, difference_type repeats)
 65:     : first(first), last(last), repeats(repeats) {}

 67:   iterator begin(void) const
 68:   {
 69:     return PermutationIterator(first, TransformIterator(CountingIterator(0), repeat_functor(repeats)));
 70:   }

 72:   iterator end(void) const
 73:   {
 74:     return begin() + repeats * (last - first);
 75:   }

 77: protected:
 78:   difference_type repeats;
 79:   Iterator        first;
 80:   Iterator        last;

 82: };

 84: // this example illustrates how to repeat blocks in a range multiple times
 85: // examples:
 86: //   tiled_range([0, 1, 2, 3], 2)    -> [0, 1, 2, 3, 0, 1, 2, 3]
 87: //   tiled_range([0, 1, 2, 3], 4, 2) -> [0, 1, 2, 3, 0, 1, 2, 3]
 88: //   tiled_range([0, 1, 2, 3], 2, 2) -> [0, 1, 0, 1, 2, 3, 2, 3]
 89: //   tiled_range([0, 1, 2, 3], 2, 3) -> [0, 1, 0, 1 0, 1, 2, 3, 2, 3, 2, 3]
 90: //   ...

 92: template <typename Iterator>
 93: class tiled_range
 94: {
 95: public:

 97:   typedef typename thrust::iterator_difference<Iterator>::type difference_type;

 99:   struct tile_functor : public thrust::unary_function<difference_type,difference_type>
100:   {
101:     difference_type repeats;
102:     difference_type tile_size;

104:     tile_functor(difference_type repeats, difference_type tile_size)
105:       : tile_size(tile_size), repeats(repeats) {}

107:     __host__ __device__
108:     difference_type operator() (const difference_type &i) const
109:     {
110:       return tile_size * (i / (tile_size * repeats)) + i % tile_size;
111:     }
112:   };

114:   typedef typename thrust::counting_iterator<difference_type>                   CountingIterator;
115:   typedef typename thrust::transform_iterator<tile_functor, CountingIterator>   TransformIterator;
116:   typedef typename thrust::permutation_iterator<Iterator,TransformIterator>     PermutationIterator;

118:   // type of the tiled_range iterator
119:   typedef PermutationIterator iterator;

121:   // construct repeated_range for the range [first,last)
122:   tiled_range(Iterator first, Iterator last, difference_type repeats)
123:     : first(first), last(last), repeats(repeats), tile_size(last - first) {}

125:   tiled_range(Iterator first, Iterator last, difference_type repeats, difference_type tile_size)
126:     : first(first), last(last), repeats(repeats), tile_size(tile_size)
127:   {
128:     // ASSERT((last - first) % tile_size == 0)
129:   }

131:   iterator begin(void) const
132:   {
133:     return PermutationIterator(first, TransformIterator(CountingIterator(0), tile_functor(repeats, tile_size)));
134:   }

136:   iterator end(void) const
137:   {
138:     return begin() + repeats * (last - first);
139:   }

141: protected:
142:   difference_type repeats;
143:   difference_type tile_size;
144:   Iterator        first;
145:   Iterator        last;
146: };

148: typedef cusp::device_memory memSpace;
149: typedef int IndexType;
150: typedef PetscScalar ValueType;
151: typedef cusp::array1d<IndexType, memSpace> IndexArray;
152: typedef cusp::array1d<ValueType, memSpace> ValueArray;
153: typedef cusp::array1d<IndexType, cusp::host_memory> IndexHostArray;
154: typedef IndexArray::iterator IndexArrayIterator;
155: typedef ValueArray::iterator ValueArrayIterator;

157: struct is_diag
158: {
159:   IndexType first, last;

161:   is_diag(IndexType first, IndexType last) : first(first), last(last) {}

163:   template <typename Tuple>
164:   __host__ __device__
165:   bool operator()(Tuple t)
166:   {
167:     // Check column
168:     IndexType row = thrust::get<0>(t);
169:     IndexType col = thrust::get<1>(t);
170:     return (row >= first) && (row < last) && (col >= first) && (col < last);
171:   }
172: };

174: struct is_nonlocal
175: {
176:   IndexType first, last;

178:   is_nonlocal(IndexType first, IndexType last) : first(first), last(last) {}

180:   template <typename Tuple>
181:   __host__ __device__
182:   bool operator() (Tuple t)
183:   {
184:     // Check column
185:     IndexType row = thrust::get<0>(t);
186:     return (row < first) || (row >= last);
187:   }
188: };

190: /*@C
191:   MatMPIAIJSetValuesBatch - Set multiple blocks of values into a matrix

193:   Not collective

195:   Input Parameters:
196: + J  - the assembled Mat object
197: . Ne -  the number of blocks (elements)
198: . Nl -  the block size (number of dof per element)
199: . elemRows - List of block row indices, in bunches of length Nl
200: - elemMats - List of block values, in bunches of Nl*Nl

202:   Level: advanced

204: .seealso MatSetValues()
205: @*/
206: PetscErrorCode MatSetValuesBatch_MPIAIJCUSP(Mat J, PetscInt Ne, PetscInt Nl, PetscInt *elemRows, const PetscScalar *elemMats)
207: {
208:   // Assumptions:
209:   //   1) Each elemMat is square, of size Nl x Nl
210:   //
211:   //      This means that any nonlocal entry (i,j) where i is owned by another process is matched to
212:   //        a) an offdiagonal entry (j,i) if j is diagonal, or
213:   //        b) another nonlocal entry (j,i) if j is offdiagonal
214:   //
215:   //      No - numSendEntries: The number of on-process  diagonal+offdiagonal entries
216:   //      numRecvEntries:      The number of off-process diagonal+offdiagonal entries
217:   //
218:   //  Glossary:
219:   //     diagonal: (i,j) such that i,j in [firstRow, lastRow)
220:   //  offdiagonal: (i,j) such that i in [firstRow, lastRow), and j not in [firstRow, lastRow)
221:   //     nonlocal: (i,j) such that i not in [firstRow, lastRow)
222:   //  nondiagonal: (i,j) such that i not in [firstRow, lastRow), or j not in [firstRow, lastRow)
223:   //   on-process: entries provided by elemMats
224:   //  off-process: entries received from other processes
225:   MPI_Comm       comm;
226:   Mat_MPIAIJ     *j   = (Mat_MPIAIJ*) J->data;
227:   size_t         N    = Ne * Nl;     // Length of elemRows (dimension of unassembled space)
228:   size_t         No   = Ne * Nl*Nl;  // Length of elemMats (total number of values)
229:   PetscInt       Nr;                 // Size of J          (dimension of assembled space)
230:   PetscInt       firstRow, lastRow, firstCol;
231:   const PetscInt *rowRanges;
232:   PetscInt       numNonlocalRows;    // Number of rows in elemRows not owned by this process
233:   PetscInt       numSendEntries;     // Number of (i,j,v) entries sent to other processes
234:   PetscInt       numRecvEntries;     // Number of (i,j,v) entries received from other processes
235:   PetscInt       Nc;
236:   PetscMPIInt    size, rank;

239:   // copy elemRows and elemMat to device
240:   IndexArray d_elemRows(elemRows, elemRows + N);
241:   ValueArray d_elemMats(elemMats, elemMats + No);

244:   PetscObjectGetComm((PetscObject)J,&comm);
245:   MPI_Comm_size(comm, &size);
246:   MPI_Comm_rank(comm, &rank);
247:   // get matrix information
248:   MatGetLocalSize(J, &Nr, NULL);
249:   MatGetOwnershipRange(J, &firstRow, &lastRow);
250:   MatGetOwnershipRanges(J, &rowRanges);
251:   MatGetOwnershipRangeColumn(J, &firstCol, NULL);
252:   PetscInfo3(J, "Assembling matrix of size %d (rows %d -- %d)\n", Nr, firstRow, lastRow);

254:   // repeat elemRows entries Nl times
255:   PetscInfo(J, "Making row indices\n");
256:   repeated_range<IndexArrayIterator> rowInd(d_elemRows.begin(), d_elemRows.end(), Nl);

258:   // tile rows of elemRows Nl times
259:   PetscInfo(J, "Making column indices\n");
260:   tiled_range<IndexArrayIterator> colInd(d_elemRows.begin(), d_elemRows.end(), Nl, Nl);

262:   // Find number of nonlocal rows, convert nonlocal rows to procs, and send sizes of off-proc entries (could send diag and offdiag sizes)
263:   // TODO: Ask Nathan how to do this on GPU
264:   PetscLogEventBegin(MAT_SetValuesBatchI,0,0,0,0);
265:   PetscMPIInt *procSendSizes, *procRecvSizes;

267:   PetscCalloc2(size, &procSendSizes, size, &procRecvSizes);

269:   numNonlocalRows = 0;
270:   for (size_t i = 0; i < N; ++i) {
271:     const PetscInt row = elemRows[i];

273:     if ((row < firstRow) || (row >= lastRow)) {
274:       numNonlocalRows++;
275:       for (IndexType p = 0; p < size; ++p) {
276:         if ((row >= rowRanges[p]) && (row < rowRanges[p+1])) {
277:           procSendSizes[p] += Nl;
278:           break;
279:         }
280:       }
281:     }
282:   }
283:   numSendEntries = numNonlocalRows*Nl;

285:   PetscInfo2(J, "Nonlocal rows %d total entries %d\n", numNonlocalRows, No);
286:   MPI_Alltoall(procSendSizes, 1, MPIU_INT, procRecvSizes, 1, MPIU_INT, comm);

288:   numRecvEntries = 0;
289:   for (PetscInt p = 0; p < size; ++p) numRecvEntries += procRecvSizes[p];
290:   PetscInfo2(j->A, "Send entries %d Recv Entries %d\n", numSendEntries, numRecvEntries);
291:   PetscLogEventEnd(MAT_SetValuesBatchI,0,0,0,0);
292:   // Allocate storage for "fat" COO representation of matrix
293:   PetscLogEventBegin(MAT_SetValuesBatchII,0,0,0,0);
294:   PetscInfo2(j->A, "Making COO matrices, diag entries %d, nondiag entries %d\n", No-numSendEntries+numRecvEntries, numSendEntries*2);
295:   cusp::coo_matrix<IndexType,ValueType, memSpace> diagCOO(Nr, Nr, No-numSendEntries+numRecvEntries); // ALLOC: This is oversized because I also count offdiagonal entries
296:   IndexArray nondiagonalRows(numSendEntries+numSendEntries); // ALLOC: This is oversized because numSendEntries > on-process offdiagonal entries
297:   IndexArray nondiagonalCols(numSendEntries+numSendEntries); // ALLOC: This is oversized because numSendEntries > on-process offdiagonal entries
298:   ValueArray nondiagonalVals(numSendEntries+numSendEntries); // ALLOC: This is oversized because numSendEntries > on-process offdiagonal entries
299:   IndexArray nonlocalRows(numSendEntries);
300:   IndexArray nonlocalCols(numSendEntries);
301:   ValueArray nonlocalVals(numSendEntries);
302:   // partition on-process entries into diagonal and off-diagonal+nonlocal
303:   PetscInfo(J, "Splitting on-process entries into diagonal and off-diagonal+nonlocal\n");
304:   thrust::fill(diagCOO.row_indices.begin(), diagCOO.row_indices.end(), -1);
305:   thrust::fill(nondiagonalRows.begin(),     nondiagonalRows.end(),     -1);
306:   thrust::partition_copy(thrust::make_zip_iterator(thrust::make_tuple(rowInd.begin(), colInd.begin(), d_elemMats.begin())),
307:                          thrust::make_zip_iterator(thrust::make_tuple(rowInd.end(),   colInd.end(),   d_elemMats.end())),
308:                          thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(),    diagCOO.column_indices.begin(), diagCOO.values.begin())),
309:                          thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(),        nondiagonalCols.begin(),        nondiagonalVals.begin())),
310:                          is_diag(firstRow, lastRow));
311:   // Current size without off-proc entries
312:   PetscInt diagonalSize    = (diagCOO.row_indices.end() - diagCOO.row_indices.begin()) - thrust::count(diagCOO.row_indices.begin(), diagCOO.row_indices.end(), -1);
313:   PetscInt nondiagonalSize = No - diagonalSize;
314:   PetscInfo2(j->A, "Diagonal size %d Nondiagonal size %d\n", diagonalSize, nondiagonalSize);
315:   ///cusp::print(diagCOO);
316:   ///cusp::print(nondiagonalRows);
317:   // partition on-process entries again into off-diagonal and nonlocal
318:   PetscInfo(J, "Splitting on-process entries into off-diagonal and nonlocal\n");
319:   thrust::stable_partition(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin())),
320:                            thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(),   nondiagonalCols.end(),   nondiagonalVals.end())),
321:                            is_nonlocal(firstRow, lastRow));
322:   PetscInt nonlocalSize    = numSendEntries;
323:   PetscInt offdiagonalSize = nondiagonalSize - nonlocalSize;
324:   PetscInfo2(j->A, "Nonlocal size %d Offdiagonal size %d\n", nonlocalSize, offdiagonalSize);
325:   PetscLogEventEnd(MAT_SetValuesBatchII,0,0,0,0);
326:   ///cusp::print(nondiagonalRows);
327:   // send off-proc entries (pack this up later)
328:   PetscLogEventBegin(MAT_SetValuesBatchIII,0,0,0,0);
329:   PetscMPIInt *procSendDispls, *procRecvDispls;
330:   PetscInt    *sendRows, *recvRows;
331:   PetscInt    *sendCols, *recvCols;
332:   PetscScalar *sendVals, *recvVals;

334:   PetscMalloc2(size, &procSendDispls, size, &procRecvDispls);
335:   PetscMalloc3(numSendEntries, &sendRows, numSendEntries, &sendCols, numSendEntries, &sendVals);
336:   PetscMalloc3(numRecvEntries, &recvRows, numRecvEntries, &recvCols, numRecvEntries, &recvVals);

338:   procSendDispls[0] = procRecvDispls[0] = 0;
339:   for (PetscInt p = 1; p < size; ++p) {
340:     procSendDispls[p] = procSendDispls[p-1] + procSendSizes[p-1];
341:     procRecvDispls[p] = procRecvDispls[p-1] + procRecvSizes[p-1];
342:   }
343:   thrust::copy(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin())),
344:                thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin()))+nonlocalSize,
345:                thrust::make_zip_iterator(thrust::make_tuple(sendRows,                sendCols,                sendVals)));
346:   //   could pack into a struct and unpack
347:   MPI_Alltoallv(sendRows, procSendSizes, procSendDispls, MPIU_INT,    recvRows, procRecvSizes, procRecvDispls, MPIU_INT,    comm);
348:   MPI_Alltoallv(sendCols, procSendSizes, procSendDispls, MPIU_INT,    recvCols, procRecvSizes, procRecvDispls, MPIU_INT,    comm);
349:   MPI_Alltoallv(sendVals, procSendSizes, procSendDispls, MPIU_SCALAR, recvVals, procRecvSizes, procRecvDispls, MPIU_SCALAR, comm);
350:   PetscFree2(procSendSizes, procRecvSizes);
351:   PetscFree2(procSendDispls, procRecvDispls);
352:   PetscFree3(sendRows, sendCols, sendVals);
353:   PetscLogEventEnd(MAT_SetValuesBatchIII,0,0,0,0);

355:   PetscLogEventBegin(MAT_SetValuesBatchIV,0,0,0,0);
356:   // Create off-diagonal matrix
357:   cusp::coo_matrix<IndexType,ValueType, memSpace> offdiagCOO(Nr, Nr, offdiagonalSize+numRecvEntries); // ALLOC: This is oversizes because we count diagonal entries in numRecvEntries
358:   // partition again into diagonal and off-diagonal
359:   IndexArray d_recvRows(recvRows, recvRows+numRecvEntries);
360:   IndexArray d_recvCols(recvCols, recvCols+numRecvEntries);
361:   ValueArray d_recvVals(recvVals, recvVals+numRecvEntries);
362:   thrust::copy(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(),          nondiagonalCols.end(),             nondiagonalVals.end()))-offdiagonalSize,
363:                thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(),          nondiagonalCols.end(),             nondiagonalVals.end())),
364:                thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin(), offdiagCOO.values.begin())));
365:   thrust::fill(diagCOO.row_indices.begin()+diagonalSize,       diagCOO.row_indices.end(),    -1);
366:   thrust::fill(offdiagCOO.row_indices.begin()+offdiagonalSize, offdiagCOO.row_indices.end(), -1);
367:   thrust::partition_copy(thrust::make_zip_iterator(thrust::make_tuple(d_recvRows.begin(), d_recvCols.begin(), d_recvVals.begin())),
368:                          thrust::make_zip_iterator(thrust::make_tuple(d_recvRows.end(),   d_recvCols.end(),   d_recvVals.end())),
369:                          thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin()+diagonalSize, diagCOO.column_indices.begin()+diagonalSize, diagCOO.values.begin()+diagonalSize)),
370:                          thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin()+offdiagonalSize, offdiagCOO.column_indices.begin()+offdiagonalSize, offdiagCOO.values.begin()+offdiagonalSize)),
371:                          is_diag(firstRow, lastRow));

373:   PetscFree3(recvRows, recvCols, recvVals);

375:   diagonalSize    = (diagCOO.row_indices.end()    - diagCOO.row_indices.begin())    - thrust::count(diagCOO.row_indices.begin(),    diagCOO.row_indices.end(),    -1);
376:   offdiagonalSize = (offdiagCOO.row_indices.end() - offdiagCOO.row_indices.begin()) - thrust::count(offdiagCOO.row_indices.begin(), offdiagCOO.row_indices.end(), -1);

378:   // sort COO format by (i,j), this is the most costly step
379:   PetscInfo(J, "Sorting rows and columns\n");
380:   diagCOO.sort_by_row_and_column();
381:   offdiagCOO.sort_by_row_and_column();
382:   PetscInt diagonalOffset    = (diagCOO.row_indices.end()    - diagCOO.row_indices.begin())    - diagonalSize;
383:   PetscInt offdiagonalOffset = (offdiagCOO.row_indices.end() - offdiagCOO.row_indices.begin()) - offdiagonalSize;

385:   // print the "fat" COO representation
386:   if (PetscLogPrintInfo) {
387: #if !defined(PETSC_USE_COMPLEX)
388:     cusp::print(diagCOO);
389:     cusp::print(offdiagCOO);
390: #endif
391:   }

393:   // Figure out the number of unique off-diagonal columns
394:   Nc = thrust::inner_product(offdiagCOO.column_indices.begin()+offdiagonalOffset,
395:                              offdiagCOO.column_indices.end()   - 1,
396:                              offdiagCOO.column_indices.begin()+offdiagonalOffset + 1,
397:                              size_t(1), thrust::plus<size_t>(), thrust::not_equal_to<IndexType>());

399:   // compute number of unique (i,j) entries
400:   //   this counts the number of changes as we move along the (i,j) list
401:   PetscInfo(J, "Computing number of unique entries\n");
402:   size_t num_diag_entries = thrust::inner_product
403:                               (thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin())) + diagonalOffset,
404:                               thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.end(),   diagCOO.column_indices.end())) - 1,
405:                               thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin())) + diagonalOffset + 1,
406:                               size_t(1),
407:                               thrust::plus<size_t>(),
408:                               thrust::not_equal_to< thrust::tuple<IndexType,IndexType> >());
409:   size_t num_offdiag_entries = thrust::inner_product
410:                                  (thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin())) + offdiagonalOffset,
411:                                  thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.end(),   offdiagCOO.column_indices.end())) - 1,
412:                                  thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin())) + offdiagonalOffset + 1,
413:                                  size_t(1),
414:                                  thrust::plus<size_t>(),
415:                                  thrust::not_equal_to< thrust::tuple<IndexType,IndexType> >());

417:   // allocate COO storage for final matrix
418:   PetscInfo(J, "Allocating compressed matrices\n");
419:   cusp::coo_matrix<IndexType, ValueType, memSpace> A(Nr, Nr, num_diag_entries);
420:   cusp::coo_matrix<IndexType, ValueType, memSpace> B(Nr, Nc, num_offdiag_entries);

422:   // sum values with the same (i,j) index
423:   // XXX thrust::reduce_by_key is unoptimized right now, so we provide a SpMV-based one in cusp::detail
424:   //     the Cusp one is 2x faster, but still not optimal
425:   // This could possibly be done in-place
426:   PetscInfo(J, "Compressing matrices\n");
427:   thrust::reduce_by_key(thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin())) + diagonalOffset,
428:                         thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.end(),   diagCOO.column_indices.end())),
429:                         diagCOO.values.begin() + diagonalOffset,
430:                         thrust::make_zip_iterator(thrust::make_tuple(A.row_indices.begin(), A.column_indices.begin())),
431:                         A.values.begin(),
432:                         thrust::equal_to< thrust::tuple<IndexType,IndexType> >(),
433:                         thrust::plus<ValueType>());

435:   thrust::reduce_by_key(thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin())) + offdiagonalOffset,
436:                         thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.end(),   offdiagCOO.column_indices.end())),
437:                         offdiagCOO.values.begin() + offdiagonalOffset,
438:                         thrust::make_zip_iterator(thrust::make_tuple(B.row_indices.begin(), B.column_indices.begin())),
439:                         B.values.begin(),
440:                         thrust::equal_to< thrust::tuple<IndexType,IndexType> >(),
441:                         thrust::plus<ValueType>());

443:   // Convert row and column numbers
444:   if (firstRow) {
445:     thrust::transform(A.row_indices.begin(), A.row_indices.end(), thrust::make_constant_iterator(firstRow), A.row_indices.begin(), thrust::minus<IndexType>());
446:     thrust::transform(B.row_indices.begin(), B.row_indices.end(), thrust::make_constant_iterator(firstRow), B.row_indices.begin(), thrust::minus<IndexType>());
447:   }
448:   if (firstCol) {
449:     thrust::transform(A.column_indices.begin(), A.column_indices.end(), thrust::make_constant_iterator(firstCol), A.column_indices.begin(), thrust::minus<IndexType>());
450:   }

452:   // print the final matrix
453:   if (PetscLogPrintInfo) {
454: #if !defined(PETSC_USE_COMPLEX)
455:     cusp::print(A);
456:     cusp::print(B);
457: #endif
458:   }

460:   PetscInfo(J, "Converting to PETSc matrix\n");
461:   MatSetType(J, MATMPIAIJCUSP);
462:   CUSPMATRIX *Agpu = new CUSPMATRIX;
463:   CUSPMATRIX *Bgpu = new CUSPMATRIX;
464:   cusp::convert(A, *Agpu);
465:   cusp::convert(B, *Bgpu);
466:   if (PetscLogPrintInfo) {
467: #if !defined(PETSC_USE_COMPLEX)
468:     cusp::print(*Agpu);
469:     cusp::print(*Bgpu);
470: #endif
471:   }
472:   {
473:     PetscInfo(J, "Copying to CPU matrix");
474:     MatCUSPCopyFromGPU(j->A, Agpu);
475:     MatCUSPCopyFromGPU(j->B, Bgpu);
476:   }
477:   PetscLogEventEnd(MAT_SetValuesBatchIV,0,0,0,0);
478:   return(0);
479: }