Actual source code: mpiaijAssemble.cu
petsc-3.5.4 2015-05-23
1: #include <petscconf.h>
2: PETSC_CUDA_EXTERN_C_BEGIN
3: #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
4: #include <../src/mat/impls/aij/mpi/mpiaij.h>
5: #include <petscbt.h>
6: #include <../src/vec/vec/impls/dvecimpl.h>
7: #include <petsc-private/vecimpl.h>
8: PETSC_CUDA_EXTERN_C_END
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>
19: #include <cusp/io/matrix_market.h>
21: #include <thrust/iterator/counting_iterator.h>
22: #include <thrust/iterator/transform_iterator.h>
23: #include <thrust/iterator/permutation_iterator.h>
24: #include <thrust/functional.h>
25: #include <thrust/partition.h>
26: #include <thrust/remove.h>
28: // this example illustrates how to make repeated access to a range of values
29: // examples:
30: // repeated_range([0, 1, 2, 3], 1) -> [0, 1, 2, 3]
31: // repeated_range([0, 1, 2, 3], 2) -> [0, 0, 1, 1, 2, 2, 3, 3]
32: // repeated_range([0, 1, 2, 3], 3) -> [0, 0, 0, 1, 1, 1, 2, 2, 2, 3, 3, 3]
33: // ...
35: template <typename Iterator>
36: class repeated_range
37: {
38: public:
40: typedef typename thrust::iterator_difference<Iterator>::type difference_type;
42: struct repeat_functor : public thrust::unary_function<difference_type,difference_type>
43: {
44: difference_type repeats;
46: repeat_functor(difference_type repeats) : repeats(repeats) {}
48: __host__ __device__
49: difference_type operator()(const difference_type &i) const
50: {
51: return i / repeats;
52: }
53: };
55: typedef typename thrust::counting_iterator<difference_type> CountingIterator;
56: typedef typename thrust::transform_iterator<repeat_functor, CountingIterator> TransformIterator;
57: typedef typename thrust::permutation_iterator<Iterator,TransformIterator> PermutationIterator;
59: // type of the repeated_range iterator
60: typedef PermutationIterator iterator;
62: // construct repeated_range for the range [first,last)
63: repeated_range(Iterator first, Iterator last, difference_type repeats)
64: : first(first), last(last), repeats(repeats) {}
66: iterator begin(void) const
67: {
68: return PermutationIterator(first, TransformIterator(CountingIterator(0), repeat_functor(repeats)));
69: }
71: iterator end(void) const
72: {
73: return begin() + repeats * (last - first);
74: }
76: protected:
77: difference_type repeats;
78: Iterator first;
79: Iterator last;
81: };
83: // this example illustrates how to repeat blocks in a range multiple times
84: // examples:
85: // tiled_range([0, 1, 2, 3], 2) -> [0, 1, 2, 3, 0, 1, 2, 3]
86: // tiled_range([0, 1, 2, 3], 4, 2) -> [0, 1, 2, 3, 0, 1, 2, 3]
87: // tiled_range([0, 1, 2, 3], 2, 2) -> [0, 1, 0, 1, 2, 3, 2, 3]
88: // tiled_range([0, 1, 2, 3], 2, 3) -> [0, 1, 0, 1 0, 1, 2, 3, 2, 3, 2, 3]
89: // ...
91: template <typename Iterator>
92: class tiled_range
93: {
94: public:
96: typedef typename thrust::iterator_difference<Iterator>::type difference_type;
98: struct tile_functor : public thrust::unary_function<difference_type,difference_type>
99: {
100: difference_type repeats;
101: difference_type tile_size;
103: tile_functor(difference_type repeats, difference_type tile_size)
104: : tile_size(tile_size), repeats(repeats) {}
106: __host__ __device__
107: difference_type operator() (const difference_type &i) const
108: {
109: return tile_size * (i / (tile_size * repeats)) + i % tile_size;
110: }
111: };
113: typedef typename thrust::counting_iterator<difference_type> CountingIterator;
114: typedef typename thrust::transform_iterator<tile_functor, CountingIterator> TransformIterator;
115: typedef typename thrust::permutation_iterator<Iterator,TransformIterator> PermutationIterator;
117: // type of the tiled_range iterator
118: typedef PermutationIterator iterator;
120: // construct repeated_range for the range [first,last)
121: tiled_range(Iterator first, Iterator last, difference_type repeats)
122: : first(first), last(last), repeats(repeats), tile_size(last - first) {}
124: tiled_range(Iterator first, Iterator last, difference_type repeats, difference_type tile_size)
125: : first(first), last(last), repeats(repeats), tile_size(tile_size)
126: {
127: // ASSERT((last - first) % tile_size == 0)
128: }
130: iterator begin(void) const
131: {
132: return PermutationIterator(first, TransformIterator(CountingIterator(0), tile_functor(repeats, tile_size)));
133: }
135: iterator end(void) const
136: {
137: return begin() + repeats * (last - first);
138: }
140: protected:
141: difference_type repeats;
142: difference_type tile_size;
143: Iterator first;
144: Iterator last;
145: };
147: typedef cusp::device_memory memSpace;
148: typedef int IndexType;
149: typedef PetscScalar ValueType;
150: typedef cusp::array1d<IndexType, memSpace> IndexArray;
151: typedef cusp::array1d<ValueType, memSpace> ValueArray;
152: typedef cusp::array1d<IndexType, cusp::host_memory> IndexHostArray;
153: typedef IndexArray::iterator IndexArrayIterator;
154: typedef ValueArray::iterator ValueArrayIterator;
156: struct is_diag
157: {
158: IndexType first, last;
160: is_diag(IndexType first, IndexType last) : first(first), last(last) {}
162: template <typename Tuple>
163: __host__ __device__
164: bool operator()(Tuple t)
165: {
166: // Check column
167: IndexType row = thrust::get<0>(t);
168: IndexType col = thrust::get<1>(t);
169: return (row >= first) && (row < last) && (col >= first) && (col < last);
170: }
171: };
173: struct is_nonlocal
174: {
175: IndexType first, last;
177: is_nonlocal(IndexType first, IndexType last) : first(first), last(last) {}
179: template <typename Tuple>
180: __host__ __device__
181: bool operator() (Tuple t)
182: {
183: // Check column
184: 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;
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;
240: // copy elemRows and elemMat to device
241: IndexArray d_elemRows(elemRows, elemRows + N);
242: ValueArray d_elemMats(elemMats, elemMats + No);
245: PetscObjectGetComm((PetscObject)J,&comm);
246: MPI_Comm_size(comm, &numProcs);
247: MPI_Comm_rank(comm, &rank);
248: // get matrix information
249: MatGetLocalSize(J, &Nr, NULL);
250: MatGetOwnershipRange(J, &firstRow, &lastRow);
251: MatGetOwnershipRanges(J, &rowRanges);
252: MatGetOwnershipRangeColumn(J, &firstCol, NULL);
253: PetscInfo3(J, "Assembling matrix of size %d (rows %d -- %d)\n", Nr, firstRow, lastRow);
255: // repeat elemRows entries Nl times
256: PetscInfo(J, "Making row indices\n");
257: repeated_range<IndexArrayIterator> rowInd(d_elemRows.begin(), d_elemRows.end(), Nl);
259: // tile rows of elemRows Nl times
260: PetscInfo(J, "Making column indices\n");
261: tiled_range<IndexArrayIterator> colInd(d_elemRows.begin(), d_elemRows.end(), Nl, Nl);
263: // Find number of nonlocal rows, convert nonlocal rows to procs, and send sizes of off-proc entries (could send diag and offdiag sizes)
264: // TODO: Ask Nathan how to do this on GPU
265: PetscLogEventBegin(MAT_SetValuesBatchI,0,0,0,0);
266: PetscMPIInt *procSendSizes, *procRecvSizes;
268: PetscCalloc2(numProcs, &procSendSizes, numProcs, &procRecvSizes);
270: numNonlocalRows = 0;
271: for (size_t i = 0; i < N; ++i) {
272: const PetscInt row = elemRows[i];
274: if ((row < firstRow) || (row >= lastRow)) {
275: numNonlocalRows++;
276: for (IndexType p = 0; p < numProcs; ++p) {
277: if ((row >= rowRanges[p]) && (row < rowRanges[p+1])) {
278: procSendSizes[p] += Nl;
279: break;
280: }
281: }
282: }
283: }
284: numSendEntries = numNonlocalRows*Nl;
286: PetscInfo2(J, "Nonlocal rows %d total entries %d\n", numNonlocalRows, No);
287: MPI_Alltoall(procSendSizes, 1, MPIU_INT, procRecvSizes, 1, MPIU_INT, comm);
289: numRecvEntries = 0;
290: for (PetscInt p = 0; p < numProcs; ++p) numRecvEntries += procRecvSizes[p];
291: PetscInfo2(j->A, "Send entries %d Recv Entries %d\n", numSendEntries, numRecvEntries);
292: PetscLogEventEnd(MAT_SetValuesBatchI,0,0,0,0);
293: // Allocate storage for "fat" COO representation of matrix
294: PetscLogEventBegin(MAT_SetValuesBatchII,0,0,0,0);
295: PetscInfo2(j->A, "Making COO matrices, diag entries %d, nondiag entries %d\n", No-numSendEntries+numRecvEntries, numSendEntries*2);
296: cusp::coo_matrix<IndexType,ValueType, memSpace> diagCOO(Nr, Nr, No-numSendEntries+numRecvEntries); // ALLOC: This is oversized because I also count offdiagonal entries
297: IndexArray nondiagonalRows(numSendEntries+numSendEntries); // ALLOC: This is oversized because numSendEntries > on-process offdiagonal entries
298: IndexArray nondiagonalCols(numSendEntries+numSendEntries); // ALLOC: This is oversized because numSendEntries > on-process offdiagonal entries
299: ValueArray nondiagonalVals(numSendEntries+numSendEntries); // ALLOC: This is oversized because numSendEntries > on-process offdiagonal entries
300: IndexArray nonlocalRows(numSendEntries);
301: IndexArray nonlocalCols(numSendEntries);
302: ValueArray nonlocalVals(numSendEntries);
303: // partition on-process entries into diagonal and off-diagonal+nonlocal
304: PetscInfo(J, "Splitting on-process entries into diagonal and off-diagonal+nonlocal\n");
305: thrust::fill(diagCOO.row_indices.begin(), diagCOO.row_indices.end(), -1);
306: thrust::fill(nondiagonalRows.begin(), nondiagonalRows.end(), -1);
307: thrust::partition_copy(thrust::make_zip_iterator(thrust::make_tuple(rowInd.begin(), colInd.begin(), d_elemMats.begin())),
308: thrust::make_zip_iterator(thrust::make_tuple(rowInd.end(), colInd.end(), d_elemMats.end())),
309: thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin(), diagCOO.values.begin())),
310: thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin())),
311: is_diag(firstRow, lastRow));
312: // Current size without off-proc entries
313: PetscInt diagonalSize = (diagCOO.row_indices.end() - diagCOO.row_indices.begin()) - thrust::count(diagCOO.row_indices.begin(), diagCOO.row_indices.end(), -1);
314: PetscInt nondiagonalSize = No - diagonalSize;
315: PetscInfo2(j->A, "Diagonal size %d Nondiagonal size %d\n", diagonalSize, nondiagonalSize);
316: ///cusp::print(diagCOO);
317: ///cusp::print(nondiagonalRows);
318: // partition on-process entries again into off-diagonal and nonlocal
319: PetscInfo(J, "Splitting on-process entries into off-diagonal and nonlocal\n");
320: thrust::stable_partition(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin())),
321: thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(), nondiagonalCols.end(), nondiagonalVals.end())),
322: is_nonlocal(firstRow, lastRow));
323: PetscInt nonlocalSize = numSendEntries;
324: PetscInt offdiagonalSize = nondiagonalSize - nonlocalSize;
325: PetscInfo2(j->A, "Nonlocal size %d Offdiagonal size %d\n", nonlocalSize, offdiagonalSize);
326: PetscLogEventEnd(MAT_SetValuesBatchII,0,0,0,0);
327: ///cusp::print(nondiagonalRows);
328: // send off-proc entries (pack this up later)
329: PetscLogEventBegin(MAT_SetValuesBatchIII,0,0,0,0);
330: PetscMPIInt *procSendDispls, *procRecvDispls;
331: PetscInt *sendRows, *recvRows;
332: PetscInt *sendCols, *recvCols;
333: PetscScalar *sendVals, *recvVals;
335: PetscMalloc2(numProcs, &procSendDispls, numProcs, &procRecvDispls);
336: PetscMalloc3(numSendEntries, &sendRows, numSendEntries, &sendCols, numSendEntries, &sendVals);
337: PetscMalloc3(numRecvEntries, &recvRows, numRecvEntries, &recvCols, numRecvEntries, &recvVals);
339: procSendDispls[0] = procRecvDispls[0] = 0;
340: for (PetscInt p = 1; p < numProcs; ++p) {
341: procSendDispls[p] = procSendDispls[p-1] + procSendSizes[p-1];
342: procRecvDispls[p] = procRecvDispls[p-1] + procRecvSizes[p-1];
343: }
344: #if 0
345: thrust::copy(nondiagonalRows.begin(), nondiagonalRows.begin()+nonlocalSize, sendRows);
346: thrust::copy(nondiagonalCols.begin(), nondiagonalCols.begin()+nonlocalSize, sendCols);
347: thrust::copy(nondiagonalVals.begin(), nondiagonalVals.begin()+nonlocalSize, sendVals);
348: #else
349: thrust::copy(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin())),
350: thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.begin(), nondiagonalCols.begin(), nondiagonalVals.begin()))+nonlocalSize,
351: thrust::make_zip_iterator(thrust::make_tuple(sendRows, sendCols, sendVals)));
352: #endif
353: // could pack into a struct and unpack
354: MPI_Alltoallv(sendRows, procSendSizes, procSendDispls, MPIU_INT, recvRows, procRecvSizes, procRecvDispls, MPIU_INT, comm);
355: MPI_Alltoallv(sendCols, procSendSizes, procSendDispls, MPIU_INT, recvCols, procRecvSizes, procRecvDispls, MPIU_INT, comm);
356: MPI_Alltoallv(sendVals, procSendSizes, procSendDispls, MPIU_SCALAR, recvVals, procRecvSizes, procRecvDispls, MPIU_SCALAR, comm);
357: PetscFree2(procSendSizes, procRecvSizes);
358: PetscFree2(procSendDispls, procRecvDispls);
359: PetscFree3(sendRows, sendCols, sendVals);
360: PetscLogEventEnd(MAT_SetValuesBatchIII,0,0,0,0);
362: PetscLogEventBegin(MAT_SetValuesBatchIV,0,0,0,0);
363: // Create off-diagonal matrix
364: cusp::coo_matrix<IndexType,ValueType, memSpace> offdiagCOO(Nr, Nr, offdiagonalSize+numRecvEntries); // ALLOC: This is oversizes because we count diagonal entries in numRecvEntries
365: // partition again into diagonal and off-diagonal
366: IndexArray d_recvRows(recvRows, recvRows+numRecvEntries);
367: IndexArray d_recvCols(recvCols, recvCols+numRecvEntries);
368: ValueArray d_recvVals(recvVals, recvVals+numRecvEntries);
369: #if 0
370: thrust::copy(nondiagonalRows.end()-offdiagonalSize, nondiagonalRows.end(), offdiagCOO.row_indices.begin());
371: thrust::copy(nondiagonalCols.end()-offdiagonalSize, nondiagonalCols.end(), offdiagCOO.column_indices.begin());
372: thrust::copy(nondiagonalVals.end()-offdiagonalSize, nondiagonalVals.end(), offdiagCOO.values.begin());
373: #else
374: thrust::copy(thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(), nondiagonalCols.end(), nondiagonalVals.end()))-offdiagonalSize,
375: thrust::make_zip_iterator(thrust::make_tuple(nondiagonalRows.end(), nondiagonalCols.end(), nondiagonalVals.end())),
376: thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin(), offdiagCOO.values.begin())));
377: #endif
378: thrust::fill(diagCOO.row_indices.begin()+diagonalSize, diagCOO.row_indices.end(), -1);
379: thrust::fill(offdiagCOO.row_indices.begin()+offdiagonalSize, offdiagCOO.row_indices.end(), -1);
380: thrust::partition_copy(thrust::make_zip_iterator(thrust::make_tuple(d_recvRows.begin(), d_recvCols.begin(), d_recvVals.begin())),
381: thrust::make_zip_iterator(thrust::make_tuple(d_recvRows.end(), d_recvCols.end(), d_recvVals.end())),
382: thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin()+diagonalSize, diagCOO.column_indices.begin()+diagonalSize, diagCOO.values.begin()+diagonalSize)),
383: thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin()+offdiagonalSize, offdiagCOO.column_indices.begin()+offdiagonalSize, offdiagCOO.values.begin()+offdiagonalSize)),
384: is_diag(firstRow, lastRow));
386: PetscFree3(recvRows, recvCols, recvVals);
388: diagonalSize = (diagCOO.row_indices.end() - diagCOO.row_indices.begin()) - thrust::count(diagCOO.row_indices.begin(), diagCOO.row_indices.end(), -1);
389: offdiagonalSize = (offdiagCOO.row_indices.end() - offdiagCOO.row_indices.begin()) - thrust::count(offdiagCOO.row_indices.begin(), offdiagCOO.row_indices.end(), -1);
391: // sort COO format by (i,j), this is the most costly step
392: PetscInfo(J, "Sorting rows and columns\n");
393: diagCOO.sort_by_row_and_column();
394: offdiagCOO.sort_by_row_and_column();
395: PetscInt diagonalOffset = (diagCOO.row_indices.end() - diagCOO.row_indices.begin()) - diagonalSize;
396: PetscInt offdiagonalOffset = (offdiagCOO.row_indices.end() - offdiagCOO.row_indices.begin()) - offdiagonalSize;
398: // print the "fat" COO representation
399: if (PetscLogPrintInfo) {
400: #if !defined(PETSC_USE_COMPLEX)
401: cusp::print(diagCOO);
402: cusp::print(offdiagCOO);
403: #endif
404: }
406: // Figure out the number of unique off-diagonal columns
407: Nc = thrust::inner_product(offdiagCOO.column_indices.begin()+offdiagonalOffset,
408: offdiagCOO.column_indices.end() - 1,
409: offdiagCOO.column_indices.begin()+offdiagonalOffset + 1,
410: size_t(1), thrust::plus<size_t>(), thrust::not_equal_to<IndexType>());
412: // compute number of unique (i,j) entries
413: // this counts the number of changes as we move along the (i,j) list
414: PetscInfo(J, "Computing number of unique entries\n");
415: size_t num_diag_entries = thrust::inner_product
416: (thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin())) + diagonalOffset,
417: thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.end(), diagCOO.column_indices.end())) - 1,
418: thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin())) + diagonalOffset + 1,
419: size_t(1),
420: thrust::plus<size_t>(),
421: thrust::not_equal_to< thrust::tuple<IndexType,IndexType> >());
422: size_t num_offdiag_entries = thrust::inner_product
423: (thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin())) + offdiagonalOffset,
424: thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.end(), offdiagCOO.column_indices.end())) - 1,
425: thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin())) + offdiagonalOffset + 1,
426: size_t(1),
427: thrust::plus<size_t>(),
428: thrust::not_equal_to< thrust::tuple<IndexType,IndexType> >());
430: // allocate COO storage for final matrix
431: PetscInfo(J, "Allocating compressed matrices\n");
432: cusp::coo_matrix<IndexType, ValueType, memSpace> A(Nr, Nr, num_diag_entries);
433: cusp::coo_matrix<IndexType, ValueType, memSpace> B(Nr, Nc, num_offdiag_entries);
435: // sum values with the same (i,j) index
436: // XXX thrust::reduce_by_key is unoptimized right now, so we provide a SpMV-based one in cusp::detail
437: // the Cusp one is 2x faster, but still not optimal
438: // This could possibly be done in-place
439: PetscInfo(J, "Compressing matrices\n");
440: thrust::reduce_by_key(thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.begin(), diagCOO.column_indices.begin())) + diagonalOffset,
441: thrust::make_zip_iterator(thrust::make_tuple(diagCOO.row_indices.end(), diagCOO.column_indices.end())),
442: diagCOO.values.begin() + diagonalOffset,
443: thrust::make_zip_iterator(thrust::make_tuple(A.row_indices.begin(), A.column_indices.begin())),
444: A.values.begin(),
445: thrust::equal_to< thrust::tuple<IndexType,IndexType> >(),
446: thrust::plus<ValueType>());
448: thrust::reduce_by_key(thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.begin(), offdiagCOO.column_indices.begin())) + offdiagonalOffset,
449: thrust::make_zip_iterator(thrust::make_tuple(offdiagCOO.row_indices.end(), offdiagCOO.column_indices.end())),
450: offdiagCOO.values.begin() + offdiagonalOffset,
451: thrust::make_zip_iterator(thrust::make_tuple(B.row_indices.begin(), B.column_indices.begin())),
452: B.values.begin(),
453: thrust::equal_to< thrust::tuple<IndexType,IndexType> >(),
454: thrust::plus<ValueType>());
456: // Convert row and column numbers
457: if (firstRow) {
458: thrust::transform(A.row_indices.begin(), A.row_indices.end(), thrust::make_constant_iterator(firstRow), A.row_indices.begin(), thrust::minus<IndexType>());
459: thrust::transform(B.row_indices.begin(), B.row_indices.end(), thrust::make_constant_iterator(firstRow), B.row_indices.begin(), thrust::minus<IndexType>());
460: }
461: if (firstCol) {
462: thrust::transform(A.column_indices.begin(), A.column_indices.end(), thrust::make_constant_iterator(firstCol), A.column_indices.begin(), thrust::minus<IndexType>());
463: }
464: #if 0 // This is done by MatSetUpMultiply_MPIAIJ()
465: // TODO: Get better code from Nathan
466: IndexArray d_colmap(Nc);
467: thrust::unique_copy(B.column_indices.begin(), B.column_indices.end(), d_colmap.begin());
468: IndexHostArray colmap(d_colmap.begin(), d_colmap.end());
469: IndexType newCol = 0;
470: for (IndexHostArray::iterator c_iter = colmap.begin(); c_iter != colmap.end(); ++c_iter, ++newCol) {
471: thrust::replace(B.column_indices.begin(), B.column_indices.end(), *c_iter, newCol);
472: }
473: #endif
475: // print the final matrix
476: if (PetscLogPrintInfo) {
477: #if !defined(PETSC_USE_COMPLEX)
478: cusp::print(A);
479: cusp::print(B);
480: #endif
481: }
483: PetscInfo(J, "Converting to PETSc matrix\n");
484: MatSetType(J, MATMPIAIJCUSP);
485: CUSPMATRIX *Agpu = new CUSPMATRIX;
486: CUSPMATRIX *Bgpu = new CUSPMATRIX;
487: cusp::convert(A, *Agpu);
488: cusp::convert(B, *Bgpu);
489: if (PetscLogPrintInfo) {
490: #if !defined(PETSC_USE_COMPLEX)
491: cusp::print(*Agpu);
492: cusp::print(*Bgpu);
493: #endif
494: }
495: {
496: PetscInfo(J, "Copying to CPU matrix");
497: MatCUSPCopyFromGPU(j->A, Agpu);
498: MatCUSPCopyFromGPU(j->B, Bgpu);
499: #if 0 // This is done by MatSetUpMultiply_MPIAIJ()
501: // Create the column map
502: PetscFree(j->garray);
503: PetscMalloc1(Nc, &j->garray);
504: PetscInt c = 0;
505: for (IndexHostArray::iterator c_iter = colmap.begin(); c_iter != colmap.end(); ++c_iter, ++c) {
506: j->garray[c] = *c_iter;
507: }
508: #endif
509: }
510: PetscLogEventEnd(MAT_SetValuesBatchIV,0,0,0,0);
511: return(0);
512: }