Actual source code: mpiaijAssemble.cu

  1: #include "petscconf.h"
 3:  #include ../src/mat/impls/aij/seq/aij.h
 4:  #include ../src/mat/impls/aij/mpi/mpiaij.h
  5: #include "petscbt.h"
 6:  #include ../src/vec/vec/impls/dvecimpl.h
  7: #include "private/vecimpl.h"
  9: #undef VecType
 10:  #include ../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h

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

 15: #include <cusp/array1d.h>
 16: #include <cusp/print.h>
 17: #include <cusp/coo_matrix.h>
 18: #include <cusp/detail/device/reduce_by_key.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)
 48:             : repeats(repeats) {}

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

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

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

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

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

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

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

 83: };

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

162:   is_diag(IndexType first, IndexType last) : first(first), last(last) {};

164:   template <typename Tuple>
165:   __host__ __device__
166:   bool operator()(Tuple t) {
167:     // Check column
168:     const IndexType row = thrust::get<0>(t);
169:     const 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:     // Check column
184:     const IndexType row = thrust::get<0>(t);
185:     return (row < first) || (row >= last);
186:   }
187: };

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

192:   Not collective

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

201:   Level: advanced

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

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

245:   MPI_Comm_size(comm, &numProcs);
246:   MPI_Comm_rank(comm, &rank);
247:   // get matrix information
248:   MatGetLocalSize(J, &Nr, PETSC_NULL);
249:   MatGetOwnershipRange(J, &firstRow, &lastRow);
250:   MatGetOwnershipRanges(J, &rowRanges);
251:   MatGetOwnershipRangeColumn(J, &firstCol, PETSC_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;
266:   PetscMalloc2(numProcs, PetscInt, &procSendSizes, numProcs, PetscInt, &procRecvSizes);
267:   PetscMemzero(procSendSizes, numProcs * sizeof(PetscInt));
268:   PetscMemzero(procRecvSizes, numProcs * sizeof(PetscInt));
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 < numProcs; ++p) {
276:         if ((row >= rowRanges[p]) && (row < rowRanges[p+1])) {
277:           procSendSizes[p] += Nl;
278:           break;
279:         }
280:       }
281:     }
282:   }
283:   numSendEntries = numNonlocalRows*Nl;
284:   PetscInfo2(J, "Nonlocal rows %d total entries %d\n", numNonlocalRows, No);
285:   MPI_Alltoall(procSendSizes, 1, MPIU_INT, procRecvSizes, 1, MPIU_INT, comm);
286:   numRecvEntries = 0;
287:   for(PetscInt p = 0; p < numProcs; ++p) {
288:     numRecvEntries += procRecvSizes[p];
289:   }
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;
333:   PetscMalloc2(numProcs, PetscInt, &procSendDispls, numProcs, PetscInt, &procRecvDispls);
334:   PetscMalloc3(numSendEntries, PetscInt, &sendRows, numSendEntries, PetscInt, &sendCols, numSendEntries, PetscScalar, &sendVals);
335:   PetscMalloc3(numRecvEntries, PetscInt, &recvRows, numRecvEntries, PetscInt, &recvCols, numRecvEntries, PetscScalar, &recvVals);
336:   procSendDispls[0] = procRecvDispls[0] = 0;
337:   for(PetscInt p = 1; p < numProcs; ++p) {
338:     procSendDispls[p] = procSendDispls[p-1] + procSendSizes[p-1];
339:     procRecvDispls[p] = procRecvDispls[p-1] + procRecvSizes[p-1];
340:   }
341: #if 0
342:   thrust::copy(nondiagonalRows.begin(), nondiagonalRows.begin()+nonlocalSize, sendRows);
343:   thrust::copy(nondiagonalCols.begin(), nondiagonalCols.begin()+nonlocalSize, sendCols);
344:   thrust::copy(nondiagonalVals.begin(), nondiagonalVals.begin()+nonlocalSize, sendVals);
345: #else
346:   thrust::copy(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin())),
347:                thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin()))+nonlocalSize,
348:                thrust::make_zip_iterator(thrust::make_tuple(sendRows,                sendCols,                sendVals)));
349: #endif
350:   //   could pack into a struct and unpack
351:   MPI_Alltoallv(sendRows, procSendSizes, procSendDispls, MPIU_INT,    recvRows, procRecvSizes, procRecvDispls, MPIU_INT,    comm);
352:   MPI_Alltoallv(sendCols, procSendSizes, procSendDispls, MPIU_INT,    recvCols, procRecvSizes, procRecvDispls, MPIU_INT,    comm);
353:   MPI_Alltoallv(sendVals, procSendSizes, procSendDispls, MPIU_SCALAR, recvVals, procRecvSizes, procRecvDispls, MPIU_SCALAR, comm);
354:   PetscFree2(procSendSizes, procRecvSizes);
355:   PetscFree2(procSendDispls, procRecvDispls);
356:   PetscFree3(sendRows, sendCols, sendVals);
357:   PetscLogEventEnd(MAT_SetValuesBatchIII,0,0,0,0);

359:   PetscLogEventBegin(MAT_SetValuesBatchIV,0,0,0,0);
360:   // Create off-diagonal matrix
361:   cusp::coo_matrix<IndexType,ValueType, memSpace> offdiagCOO(Nr, Nr, offdiagonalSize+numRecvEntries); // ALLOC: This is oversizes because we count diagonal entries in numRecvEntries
362:   // partition again into diagonal and off-diagonal
363:   IndexArray d_recvRows(recvRows, recvRows+numRecvEntries);
364:   IndexArray d_recvCols(recvCols, recvCols+numRecvEntries);
365:   ValueArray d_recvVals(recvVals, recvVals+numRecvEntries);
366: #if 0
367:   thrust::copy(nondiagonalRows.end()-offdiagonalSize, nondiagonalRows.end(), offdiagCOO.row_indices.begin());
368:   thrust::copy(nondiagonalCols.end()-offdiagonalSize, nondiagonalCols.end(), offdiagCOO.column_indices.begin());
369:   thrust::copy(nondiagonalVals.end()-offdiagonalSize, nondiagonalVals.end(), offdiagCOO.values.begin());
370: #else
371:   thrust::copy(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(),          nondiagonalCols.end(),             nondiagonalVals.end()))-offdiagonalSize,
372:                thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(),          nondiagonalCols.end(),             nondiagonalVals.end())),
373:                thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin(), offdiagCOO.values.begin())));
374: #endif
375:   thrust::fill(diagCOO.row_indices.begin()+diagonalSize,       diagCOO.row_indices.end(),    -1);
376:   thrust::fill(offdiagCOO.row_indices.begin()+offdiagonalSize, offdiagCOO.row_indices.end(), -1);
377:   thrust::partition_copy(thrust::make_zip_iterator(thrust::make_tuple(d_recvRows.begin(), d_recvCols.begin(), d_recvVals.begin())),
378:                          thrust::make_zip_iterator(thrust::make_tuple(d_recvRows.end(),   d_recvCols.end(),   d_recvVals.end())),
379:                          thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin()+diagonalSize, diagCOO.column_indices.begin()+diagonalSize, diagCOO.values.begin()+diagonalSize)),
380:                          thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin()+offdiagonalSize, offdiagCOO.column_indices.begin()+offdiagonalSize, offdiagCOO.values.begin()+offdiagonalSize)),
381:                          is_diag(firstRow, lastRow));
382:   PetscFree3(recvRows, recvCols, recvVals);
383:   diagonalSize    = (diagCOO.row_indices.end()    - diagCOO.row_indices.begin())    - thrust::count(diagCOO.row_indices.begin(),    diagCOO.row_indices.end(),    -1);
384:   offdiagonalSize = (offdiagCOO.row_indices.end() - offdiagCOO.row_indices.begin()) - thrust::count(offdiagCOO.row_indices.begin(), offdiagCOO.row_indices.end(), -1);

386:   // sort COO format by (i,j), this is the most costly step
387:   PetscInfo(J, "Sorting rows and columns\n");
388:   diagCOO.sort_by_row_and_column();
389:   offdiagCOO.sort_by_row_and_column();
390:   PetscInt diagonalOffset    = (diagCOO.row_indices.end()    - diagCOO.row_indices.begin())    - diagonalSize;
391:   PetscInt offdiagonalOffset = (offdiagCOO.row_indices.end() - offdiagCOO.row_indices.begin()) - offdiagonalSize;

393:   // print the "fat" COO representation
394:   if (PetscLogPrintInfo) {
395:     cusp::print(diagCOO);
396:     cusp::print(offdiagCOO);
397:   }

399:   // Figure out the number of unique off-diagonal columns
400:   Nc = thrust::inner_product(offdiagCOO.column_indices.begin()+offdiagonalOffset,
401:                              offdiagCOO.column_indices.end()   - 1,
402:                              offdiagCOO.column_indices.begin()+offdiagonalOffset + 1,
403:                              size_t(1), thrust::plus<size_t>(), thrust::not_equal_to<IndexType>());

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

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

428:   // sum values with the same (i,j) index
429:   // XXX thrust::reduce_by_key is unoptimized right now, so we provide a SpMV-based one in cusp::detail
430:   //     the Cusp one is 2x faster, but still not optimal
431:   // This could possibly be done in-place
432:   PetscInfo(J, "Compressing matrices\n");
433:   cusp::detail::device::reduce_by_key
434:     (thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin())) + diagonalOffset,
435:      thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.end(),   diagCOO.column_indices.end())),
436:      diagCOO.values.begin() + diagonalOffset,
437:      thrust::make_zip_iterator(thrust::make_tuple(A.row_indices.begin(), A.column_indices.begin())),
438:      A.values.begin(),
439:      thrust::equal_to< thrust::tuple<IndexType,IndexType> >(),
440:      thrust::plus<ValueType>());
441:   cusp::detail::device::reduce_by_key
442:     (thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin())) + offdiagonalOffset,
443:      thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.end(),   offdiagCOO.column_indices.end())),
444:      offdiagCOO.values.begin() + offdiagonalOffset,
445:      thrust::make_zip_iterator(thrust::make_tuple(B.row_indices.begin(), B.column_indices.begin())),
446:      B.values.begin(),
447:      thrust::equal_to< thrust::tuple<IndexType,IndexType> >(),
448:      thrust::plus<ValueType>());

450:   // Convert row and column numbers
451:   if (firstRow) {
452:     thrust::transform(A.row_indices.begin(), A.row_indices.end(), thrust::make_constant_iterator(firstRow), A.row_indices.begin(), thrust::minus<IndexType>());
453:     thrust::transform(B.row_indices.begin(), B.row_indices.end(), thrust::make_constant_iterator(firstRow), B.row_indices.begin(), thrust::minus<IndexType>());
454:   }
455:   if (firstCol) {
456:     thrust::transform(A.column_indices.begin(), A.column_indices.end(), thrust::make_constant_iterator(firstCol), A.column_indices.begin(), thrust::minus<IndexType>());
457:   }
458: #if 0 // This is done by MatSetUpMultiply_MPIAIJ()
459:   //   TODO: Get better code from Nathan
460:   IndexArray d_colmap(Nc);
461:   thrust::unique_copy(B.column_indices.begin(), B.column_indices.end(), d_colmap.begin());
462:   IndexHostArray colmap(d_colmap.begin(), d_colmap.end());
463:   IndexType      newCol = 0;
464:   for(IndexHostArray::iterator c_iter = colmap.begin(); c_iter != colmap.end(); ++c_iter, ++newCol) {
465:     thrust::replace(B.column_indices.begin(), B.column_indices.end(), *c_iter, newCol);
466:   }
467: #endif

469:   // print the final matrix
470:   if (PetscLogPrintInfo) {
471:     cusp::print(A);
472:     cusp::print(B);
473:   }

475:   PetscInfo(J, "Converting to PETSc matrix\n");
476:   MatSetType(J, MATMPIAIJCUSP);
477:   CUSPMATRIX *Agpu = new CUSPMATRIX;
478:   CUSPMATRIX *Bgpu = new CUSPMATRIX;
479:   cusp::convert(A, *Agpu);
480:   cusp::convert(B, *Bgpu);
481:   if (PetscLogPrintInfo) {
482:     cusp::print(*Agpu);
483:     cusp::print(*Bgpu);
484:   }
485:   {
486:     PetscInfo(J, "Copying to CPU matrix");
487:     MatCUSPCopyFromGPU(j->A, Agpu);
488:     MatCUSPCopyFromGPU(j->B, Bgpu);
489: #if 0 // This is done by MatSetUpMultiply_MPIAIJ()
490:     // Create the column map
491:     PetscFree(j->garray);
492:     PetscMalloc(Nc * sizeof(PetscInt), &j->garray);
493:     PetscInt c = 0;
494:     for(IndexHostArray::iterator c_iter = colmap.begin(); c_iter != colmap.end(); ++c_iter, ++c) {
495:       j->garray[c] = *c_iter;
496:     }
497: #endif
498:   }
499:   PetscLogEventEnd(MAT_SetValuesBatchIV,0,0,0,0);
500:   return(0);
501: }