Actual source code: aijcusparseband.cu

  1: /*
  2:   AIJCUSPARSE methods implemented with Cuda kernels. Uses cuSparse/Thrust maps from AIJCUSPARSE
  3: */
  4: #define PETSC_SKIP_SPINLOCK
  5: #define PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND 1

  7: #include <petscconf.h>
  8: #include <../src/mat/impls/aij/seq/aij.h>
  9: #include <../src/mat/impls/sbaij/seq/sbaij.h>
 10: #undef VecType
 11: #include <../src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h>
 12: #if defined(__CUDA_ARCH__) && __CUDA_ARCH__ > 600 && PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)
 13:   #define AIJBANDUSEGROUPS 1
 14: #endif
 15: #if defined(AIJBANDUSEGROUPS)
 16:   #include <cooperative_groups.h>
 17: #endif

 19: /*
 20:   LU BAND factorization with optimization for block diagonal (Nf blocks) in natural order (-mat_no_inode -pc_factor_mat_ordering_type rcm with Nf>1 fields)

 22:   requires:
 23:      structurally symmetric: fix with transpose/column meta data
 24: */

 26: static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSEBAND(Mat, Mat, IS, IS, const MatFactorInfo *);
 27: static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSEBAND(Mat, Mat, const MatFactorInfo *);

 29: /*
 30:   The GPU LU factor kernel
 31: */
 32: __global__ void __launch_bounds__(1024, 1) mat_lu_factor_band_init_set_i(const PetscInt n, const int bw, int bi_csr[])
 33: {
 34:   const PetscInt Nf = gridDim.x, Nblk = gridDim.y, nloc = n / Nf;
 35:   const PetscInt field = blockIdx.x, blkIdx = blockIdx.y;
 36:   const PetscInt nloc_i = (nloc / Nblk + !!(nloc % Nblk)), start_i = field * nloc + blkIdx * nloc_i, end_i = (start_i + nloc_i) > (field + 1) * nloc ? (field + 1) * nloc : (start_i + nloc_i);

 38:   // set i (row+1)
 39:   if (threadIdx.x + threadIdx.y + blockIdx.x + blockIdx.y == 0) bi_csr[0] = 0; // dummy at zero
 40:   for (int rowb = start_i + threadIdx.y; rowb < end_i; rowb += blockDim.y) {   // rows in block by thread y
 41:     if (rowb < end_i && threadIdx.x == 0) {
 42:       PetscInt i = rowb + 1, ni = (rowb > bw) ? bw + 1 : i, n1L = ni * (ni - 1) / 2, nug = i * bw, n2L = bw * ((rowb > bw) ? (rowb - bw) : 0), mi = bw + rowb + 1 - n, clip = (mi > 0) ? mi * (mi - 1) / 2 + mi : 0;
 43:       bi_csr[rowb + 1] = n1L + nug - clip + n2L + i;
 44:     }
 45:   }
 46: }
 47: // copy AIJ to AIJ_BAND
 48: __global__ void __launch_bounds__(1024, 1) mat_lu_factor_band_copy_aij_aij(const PetscInt n, const int bw, const PetscInt r[], const PetscInt ic[], const int ai_d[], const int aj_d[], const PetscScalar aa_d[], const int bi_csr[], PetscScalar ba_csr[])
 49: {
 50:   const PetscInt Nf = gridDim.x, Nblk = gridDim.y, nloc = n / Nf;
 51:   const PetscInt field = blockIdx.x, blkIdx = blockIdx.y;
 52:   const PetscInt nloc_i = (nloc / Nblk + !!(nloc % Nblk)), start_i = field * nloc + blkIdx * nloc_i, end_i = (start_i + nloc_i) > (field + 1) * nloc ? (field + 1) * nloc : (start_i + nloc_i);

 54:   // zero B
 55:   if (threadIdx.x + threadIdx.y + blockIdx.x + blockIdx.y == 0) ba_csr[bi_csr[n]] = 0; // flop count at end
 56:   for (int rowb = start_i + threadIdx.y; rowb < end_i; rowb += blockDim.y) {           // rows in block by thread y
 57:     if (rowb < end_i) {
 58:       PetscScalar   *batmp = ba_csr + bi_csr[rowb];
 59:       const PetscInt nzb   = bi_csr[rowb + 1] - bi_csr[rowb];
 60:       for (int j = threadIdx.x; j < nzb; j += blockDim.x) {
 61:         if (j < nzb) batmp[j] = 0;
 62:       }
 63:     }
 64:   }

 66:   // copy A into B with CSR format -- these two loops can be fused
 67:   for (int rowb = start_i + threadIdx.y; rowb < end_i; rowb += blockDim.y) { // rows in block by thread y
 68:     if (rowb < end_i) {
 69:       const PetscInt     rowa = r[rowb], nza = ai_d[rowa + 1] - ai_d[rowa];
 70:       const int         *ajtmp = aj_d + ai_d[rowa], bjStart = (rowb > bw) ? rowb - bw : 0;
 71:       const PetscScalar *av    = aa_d + ai_d[rowa];
 72:       PetscScalar       *batmp = ba_csr + bi_csr[rowb];
 73:       /* load in initial (unfactored row) */
 74:       for (int j = threadIdx.x; j < nza; j += blockDim.x) {
 75:         if (j < nza) {
 76:           PetscInt    colb = ic[ajtmp[j]], idx = colb - bjStart;
 77:           PetscScalar vala = av[j];
 78:           batmp[idx]       = vala;
 79:         }
 80:       }
 81:     }
 82:   }
 83: }
 84: // print AIJ_BAND
 85: __global__ void print_mat_aij_band(const PetscInt n, const int bi_csr[], const PetscScalar ba_csr[])
 86: {
 87:   // debug
 88:   if (threadIdx.x + threadIdx.y + blockIdx.x + blockIdx.y == 0) {
 89:     printf("B (AIJ) n=%d:\n", (int)n);
 90:     for (int rowb = 0; rowb < n; rowb++) {
 91:       const PetscInt     nz    = bi_csr[rowb + 1] - bi_csr[rowb];
 92:       const PetscScalar *batmp = ba_csr + bi_csr[rowb];
 93:       for (int j = 0; j < nz; j++) printf("(%13.6e) ", PetscRealPart(batmp[j]));
 94:       printf(" bi=%d\n", bi_csr[rowb + 1]);
 95:     }
 96:   }
 97: }
 98: // Band LU kernel ---  ba_csr bi_csr
 99: __global__ void __launch_bounds__(1024, 1) mat_lu_factor_band(const PetscInt n, const PetscInt bw, const int bi_csr[], PetscScalar ba_csr[], int *use_group_sync)
100: {
101:   const PetscInt Nf = gridDim.x, Nblk = gridDim.y, nloc = n / Nf;
102:   const PetscInt field = blockIdx.x, blkIdx = blockIdx.y;
103:   const PetscInt start = field * nloc, end = start + nloc;
104: #if defined(AIJBANDUSEGROUPS)
105:   auto g = cooperative_groups::this_grid();
106: #endif
107:   // A22 panel update for each row A(1,:) and col A(:,1)
108:   for (int glbDD = start, locDD = 0; glbDD < end; glbDD++, locDD++) {
109:     PetscInt           tnzUd = bw, maxU = end - 1 - glbDD;                                        // we are chopping off the inter ears
110:     const PetscInt     nzUd = (tnzUd > maxU) ? maxU : tnzUd, dOffset = (glbDD > bw) ? bw : glbDD; // global to go past ears after first
111:     PetscScalar       *pBdd   = ba_csr + bi_csr[glbDD] + dOffset;
112:     const PetscScalar *baUd   = pBdd + 1; // vector of data  U(i,i+1:end)
113:     const PetscScalar  Bdd    = *pBdd;
114:     const PetscInt     offset = blkIdx * blockDim.y + threadIdx.y, inc = Nblk * blockDim.y;
115:     if (threadIdx.x == 0) {
116:       for (int idx = offset, myi = glbDD + offset + 1; idx < nzUd; idx += inc, myi += inc) { /* assuming symmetric structure */
117:         const PetscInt bwi = myi > bw ? bw : myi, kIdx = bwi - (myi - glbDD);                // cuts off just the first (global) block
118:         PetscScalar   *Aid = ba_csr + bi_csr[myi] + kIdx;
119:         *Aid               = *Aid / Bdd;
120:       }
121:     }
122:     __syncthreads(); // synch on threadIdx.x only
123:     for (int idx = offset, myi = glbDD + offset + 1; idx < nzUd; idx += inc, myi += inc) {
124:       const PetscInt    bwi = myi > bw ? bw : myi, kIdx = bwi - (myi - glbDD); // cuts off just the first (global) block
125:       PetscScalar      *Aid = ba_csr + bi_csr[myi] + kIdx;
126:       PetscScalar      *Aij = Aid + 1;
127:       const PetscScalar Lid = *Aid;
128:       for (int jIdx = threadIdx.x; jIdx < nzUd; jIdx += blockDim.x) Aij[jIdx] -= Lid * baUd[jIdx];
129:     }
130: #if defined(AIJBANDUSEGROUPS)
131:     if (use_group_sync) {
132:       g.sync();
133:     } else {
134:       __syncthreads();
135:     }
136: #else
137:     __syncthreads();
138: #endif
139:   } /* endof for (i=0; i<n; i++) { */
140: }

142: static PetscErrorCode MatSolve_SeqAIJCUSPARSEBAND(Mat, Vec, Vec);
143: static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSEBAND(Mat B, Mat A, const MatFactorInfo *info)
144: {
145:   Mat_SeqAIJ                   *b                  = (Mat_SeqAIJ *)B->data;
146:   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors *)B->spptr;
148:   Mat_SeqAIJCUSPARSE           *cusparsestructA = (Mat_SeqAIJCUSPARSE *)A->spptr;
149:   Mat_SeqAIJCUSPARSEMultStruct *matstructA;
150:   CsrMatrix                    *matrixA;
151:   const PetscInt                n = A->rmap->n, *ic, *r;
152:   const int                    *ai_d, *aj_d;
153:   const PetscScalar            *aa_d;
154:   PetscScalar                  *ba_t = cusparseTriFactors->a_band_d;
155:   int                          *bi_t = cusparseTriFactors->i_band_d;
156:   int                           Ni = 10, team_size = 9, Nf = 1, nVec = 56, nconcurrent = 1, nsm = -1; // Nf is batch size - not used

158:   if (A->rmap->n == 0) return 0;
159:   // cusparse setup
161:   matstructA = (Mat_SeqAIJCUSPARSEMultStruct *)cusparsestructA->mat; //  matstruct->cprowIndices
163:   matrixA = (CsrMatrix *)matstructA->mat;

166:   // get data
167:   ic   = thrust::raw_pointer_cast(cusparseTriFactors->cpermIndices->data());
168:   ai_d = thrust::raw_pointer_cast(matrixA->row_offsets->data());
169:   aj_d = thrust::raw_pointer_cast(matrixA->column_indices->data());
170:   aa_d = thrust::raw_pointer_cast(matrixA->values->data().get());
171:   r    = thrust::raw_pointer_cast(cusparseTriFactors->rpermIndices->data());

173:   WaitForCUDA();
174:   PetscLogGpuTimeBegin();
175:   {
176:     int bw = (int)(2. * (double)n - 1. - (double)(PetscSqrtReal(1. + 4. * ((double)n * (double)n - (double)b->nz)) + PETSC_MACHINE_EPSILON)) / 2, bm1 = bw - 1, nl = n / Nf;
177: #if !defined(AIJBANDUSEGROUPS)
178:     Ni = 1 / nconcurrent;
179:     Ni = 1;
180: #else
181:     if (!cusparseTriFactors->init_dev_prop) {
182:       int gpuid;
183:       cusparseTriFactors->init_dev_prop = PETSC_TRUE;
184:       cudaGetDevice(&gpuid);
185:       cudaGetDeviceProperties(&cusparseTriFactors->dev_prop, gpuid);
186:     }
187:     nsm = cusparseTriFactors->dev_prop.multiProcessorCount;
188:     Ni  = nsm / Nf / nconcurrent;
189: #endif
190:     team_size = bw / Ni + !!(bw % Ni);
191:     nVec      = PetscMin(bw, 1024 / team_size);
192:     PetscInfo(A, "Matrix Bandwidth = %d, number SMs/block = %d, num concurrency = %d, num fields = %d, numSMs/GPU = %d, thread group size = %d,%d\n", bw, Ni, nconcurrent, Nf, nsm, team_size, nVec);
193:     {
194:       dim3 dimBlockTeam(nVec, team_size);
195:       dim3 dimBlockLeague(Nf, Ni);
196:       mat_lu_factor_band_copy_aij_aij<<<dimBlockLeague, dimBlockTeam>>>(n, bw, r, ic, ai_d, aj_d, aa_d, bi_t, ba_t);
197:       PetscCUDACheckLaunch; // does a sync
198: #if defined(AIJBANDUSEGROUPS)
199:       if (Ni > 1) {
200:         void *kernelArgs[] = {(void *)&n, (void *)&bw, (void *)&bi_t, (void *)&ba_t, (void *)&nsm};
201:         cudaLaunchCooperativeKernel((void *)mat_lu_factor_band, dimBlockLeague, dimBlockTeam, kernelArgs, 0, NULL);
202:       } else {
203:         mat_lu_factor_band<<<dimBlockLeague, dimBlockTeam>>>(n, bw, bi_t, ba_t, NULL);
204:       }
205: #else
206:       mat_lu_factor_band<<<dimBlockLeague, dimBlockTeam>>>(n, bw, bi_t, ba_t, NULL);
207: #endif
208:       PetscCUDACheckLaunch; // does a sync
209: #if defined(PETSC_USE_LOG)
210:       PetscLogGpuFlops((PetscLogDouble)Nf * (bm1 * (bm1 + 1) * (PetscLogDouble)(2 * bm1 + 1) / 3 + (PetscLogDouble)2 * (nl - bw) * bw * bw + (PetscLogDouble)nl * (nl + 1) / 2));
211: #endif
212:     }
213:   }
214:   PetscLogGpuTimeEnd();
215:   /* determine which version of MatSolve needs to be used. from MatLUFactorNumeric_AIJ_SeqAIJCUSPARSE */
216:   B->ops->solve             = MatSolve_SeqAIJCUSPARSEBAND;
217:   B->ops->solvetranspose    = NULL; // need transpose
218:   B->ops->matsolve          = NULL;
219:   B->ops->matsolvetranspose = NULL;
220:   return 0;
221: }

223: PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSEBAND(Mat B, Mat A, IS isrow, IS iscol, const MatFactorInfo *info)
224: {
225:   Mat_SeqAIJ                   *a = (Mat_SeqAIJ *)A->data, *b;
226:   IS                            isicol;
227:   const PetscInt               *ic, *ai = a->i, *aj = a->j;
228:   PetscScalar                  *ba_t;
229:   int                          *bi_t;
230:   PetscInt                      i, n = A->rmap->n, Nf = 1; // Nf batch size - not used
231:   PetscInt                      nzBcsr, bwL, bwU;
232:   PetscBool                     missing;
233:   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors *)B->spptr;

236:   MatMissingDiagonal(A, &missing, &i);
239:   MatIsStructurallySymmetric(A, &missing);

242:   ISInvertPermutation(iscol, PETSC_DECIDE, &isicol);
243:   ISGetIndices(isicol, &ic);

245:   MatSeqAIJSetPreallocation_SeqAIJ(B, MAT_SKIP_ALLOCATION, NULL);
246:   b = (Mat_SeqAIJ *)(B)->data;

248:   /* get band widths, MatComputeBandwidth should take a reordering ic and do this */
249:   bwL = bwU = 0;
250:   for (int rwb = 0; rwb < n; rwb++) {
251:     const PetscInt rwa = ic[rwb], anz = ai[rwb + 1] - ai[rwb], *ajtmp = aj + ai[rwb];
252:     for (int j = 0; j < anz; j++) {
253:       PetscInt colb = ic[ajtmp[j]];
254:       if (colb < rwa) { // L
255:         if (rwa - colb > bwL) bwL = rwa - colb;
256:       } else {
257:         if (colb - rwa > bwU) bwU = colb - rwa;
258:       }
259:     }
260:   }
261:   ISRestoreIndices(isicol, &ic);
262:   /* only support structurally symmetric, but it might work */
264:   MatSeqAIJCUSPARSETriFactors_Reset(&cusparseTriFactors);
265:   nzBcsr   = n + (2 * n - 1) * bwU - bwU * bwU;
266:   b->maxnz = b->nz        = nzBcsr;
267:   cusparseTriFactors->nnz = b->nz; // only meta data needed: n & nz
268:   PetscInfo(A, "Matrix Bandwidth = %" PetscInt_FMT ", nnz = %" PetscInt_FMT "\n", bwL, b->nz);
269:   if (!cusparseTriFactors->workVector) cusparseTriFactors->workVector = new THRUSTARRAY(n);
270:   cudaMalloc(&ba_t, (b->nz + 1) * sizeof(PetscScalar)); // include a place for flops
271:   cudaMalloc(&bi_t, (n + 1) * sizeof(int));
272:   cusparseTriFactors->a_band_d = ba_t;
273:   cusparseTriFactors->i_band_d = bi_t;
274:   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
275:   {
276:     dim3 dimBlockTeam(1, 128);
277:     dim3 dimBlockLeague(Nf, 1);
278:     mat_lu_factor_band_init_set_i<<<dimBlockLeague, dimBlockTeam>>>(n, bwU, bi_t);
279:   }
280:   PetscCUDACheckLaunch; // does a sync

282:   // setup data
283:   if (!cusparseTriFactors->rpermIndices) {
284:     const PetscInt *r;

286:     ISGetIndices(isrow, &r);
287:     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
288:     cusparseTriFactors->rpermIndices->assign(r, r + n);
289:     ISRestoreIndices(isrow, &r);
290:     PetscLogCpuToGpu(n * sizeof(PetscInt));
291:   }
292:   /* upper triangular indices */
293:   if (!cusparseTriFactors->cpermIndices) {
294:     const PetscInt *c;

296:     ISGetIndices(isicol, &c);
297:     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
298:     cusparseTriFactors->cpermIndices->assign(c, c + n);
299:     ISRestoreIndices(isicol, &c);
300:     PetscLogCpuToGpu(n * sizeof(PetscInt));
301:   }

303:   /* put together the new matrix */
304:   b->free_a       = PETSC_FALSE;
305:   b->free_ij      = PETSC_FALSE;
306:   b->singlemalloc = PETSC_FALSE;
307:   b->ilen         = NULL;
308:   b->imax         = NULL;
309:   b->row          = isrow;
310:   b->col          = iscol;
311:   PetscObjectReference((PetscObject)isrow);
312:   PetscObjectReference((PetscObject)iscol);
313:   b->icol = isicol;
314:   PetscMalloc1(n + 1, &b->solve_work);

316:   B->factortype            = MAT_FACTOR_LU;
317:   B->info.factor_mallocs   = 0;
318:   B->info.fill_ratio_given = 0;

320:   if (ai[n]) {
321:     B->info.fill_ratio_needed = ((PetscReal)(nzBcsr)) / ((PetscReal)ai[n]);
322:   } else {
323:     B->info.fill_ratio_needed = 0.0;
324:   }
325: #if defined(PETSC_USE_INFO)
326:   if (ai[n] != 0) {
327:     PetscReal af = B->info.fill_ratio_needed;
328:     PetscInfo(A, "Band fill ratio %g\n", (double)af);
329:   } else {
330:     PetscInfo(A, "Empty matrix\n");
331:   }
332: #endif
333:   if (a->inode.size) PetscInfo(A, "Warning: using inodes in band solver.\n");
334:   MatSeqAIJCheckInode_FactorLU(B);
335:   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSEBAND;
336:   B->offloadmask          = PETSC_OFFLOAD_GPU;

338:   return 0;
339: }

341: /* Use -pc_factor_mat_solver_type cusparseband */
342: PetscErrorCode MatFactorGetSolverType_seqaij_cusparse_band(Mat A, MatSolverType *type)
343: {
344:   *type = MATSOLVERCUSPARSEBAND;
345:   return 0;
346: }

348: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse_band(Mat A, MatFactorType ftype, Mat *B)
349: {
350:   PetscInt n = A->rmap->n;

352:   MatCreate(PetscObjectComm((PetscObject)A), B);
353:   MatSetSizes(*B, n, n, n, n);
354:   (*B)->factortype     = ftype;
355:   (*B)->canuseordering = PETSC_TRUE;
356:   MatSetType(*B, MATSEQAIJCUSPARSE);

358:   if (ftype == MAT_FACTOR_LU) {
359:     MatSetBlockSizesFromMats(*B, A, A);
360:     (*B)->ops->ilufactorsymbolic = NULL; // MatILUFactorSymbolic_SeqAIJCUSPARSE;
361:     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJCUSPARSEBAND;
362:     PetscStrallocpy(MATORDERINGRCM, (char **)&(*B)->preferredordering[MAT_FACTOR_LU]);
363:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported for CUSPARSEBAND Matrix Types");

365:   MatSeqAIJSetPreallocation(*B, MAT_SKIP_ALLOCATION, NULL);
366:   PetscObjectComposeFunction((PetscObject)(*B), "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_cusparse_band);
367:   return 0;
368: }

370: #define WARP_SIZE 32
371: template <typename T>
372: __forceinline__ __device__ T wreduce(T a)
373: {
374:   T b;
375: #pragma unroll
376:   for (int i = WARP_SIZE / 2; i >= 1; i = i >> 1) {
377:     b = __shfl_down_sync(0xffffffff, a, i);
378:     a += b;
379:   }
380:   return a;
381: }
382: // reduce in a block, returns result in thread 0
383: template <typename T, int BLOCK_SIZE>
384: __device__ T breduce(T a)
385: {
386:   constexpr int     NWARP = BLOCK_SIZE / WARP_SIZE;
387:   __shared__ double buf[NWARP];
388:   int               wid    = threadIdx.x / WARP_SIZE;
389:   int               laneid = threadIdx.x % WARP_SIZE;
390:   T                 b      = wreduce<T>(a);
391:   if (laneid == 0) buf[wid] = b;
392:   __syncthreads();
393:   if (wid == 0) {
394:     if (threadIdx.x < NWARP) a = buf[threadIdx.x];
395:     else a = 0;
396:     for (int i = (NWARP + 1) / 2; i >= 1; i = i >> 1) a += __shfl_down_sync(0xffffffff, a, i);
397:   }
398:   return a;
399: }

401: // Band LU kernel ---  ba_csr bi_csr
402: template <int BLOCK_SIZE>
403: __global__ void __launch_bounds__(256, 1) mat_solve_band(const PetscInt n, const PetscInt bw, const PetscScalar ba_csr[], PetscScalar x[])
404: {
405:   const PetscInt     Nf = gridDim.x, nloc = n / Nf, field = blockIdx.x, start = field * nloc, end = start + nloc, chopnz = bw * (bw + 1) / 2, blocknz = (2 * bw + 1) * nloc, blocknz_0 = blocknz - chopnz;
406:   const PetscScalar *pLi;
407:   const int          tid = threadIdx.x;

409:   /* Next, solve L */
410:   pLi = ba_csr + (field == 0 ? 0 : blocknz_0 + (field - 1) * blocknz + bw); // diagonal (0,0) in field
411:   for (int glbDD = start, locDD = 0; glbDD < end; glbDD++, locDD++) {
412:     const PetscInt col = locDD < bw ? start : (glbDD - bw);
413:     PetscScalar    t   = 0;
414:     for (int j = col + tid, idx = tid; j < glbDD; j += blockDim.x, idx += blockDim.x) t += pLi[idx] * x[j];
415: #if defined(PETSC_USE_COMPLEX)
416:     PetscReal   tr = PetscRealPartComplex(t), ti = PetscImaginaryPartComplex(t);
417:     PetscScalar tt(breduce<PetscReal, BLOCK_SIZE>(tr), breduce<PetscReal, BLOCK_SIZE>(ti));
418:     t = tt;
419: #else
420:     t = breduce<PetscReal, BLOCK_SIZE>(t);
421: #endif
422:     if (threadIdx.x == 0) x[glbDD] -= t; // /1.0
423:     __syncthreads();
424:     // inc
425:     pLi += glbDD - col;                           // get to diagonal
426:     if (glbDD > n - 1 - bw) pLi += n - 1 - glbDD; // skip over U, only last block has funny offset
427:     else pLi += bw;
428:     pLi += 1;                                                   // skip to next row
429:     if (field > 0 && (locDD + 1) < bw) pLi += bw - (locDD + 1); // skip padding at beginning (ear)
430:   }
431:   /* Then, solve U */
432:   pLi = ba_csr + Nf * blocknz - 2 * chopnz - 1;                            // end of real data on block (diagonal)
433:   if (field != Nf - 1) pLi -= blocknz_0 + (Nf - 2 - field) * blocknz + bw; // diagonal of last local row

435:   for (int glbDD = end - 1, locDD = 0; glbDD >= start; glbDD--, locDD++) {
436:     const PetscInt col = (locDD < bw) ? end - 1 : glbDD + bw; // end of row in U
437:     PetscScalar    t   = 0;
438:     for (int j = col - tid, idx = tid; j > glbDD; j -= blockDim.x, idx += blockDim.x) t += pLi[-idx] * x[j];
439: #if defined(PETSC_USE_COMPLEX)
440:     PetscReal   tr = PetscRealPartComplex(t), ti = PetscImaginaryPartComplex(t);
441:     PetscScalar tt(breduce<PetscReal, BLOCK_SIZE>(tr), breduce<PetscReal, BLOCK_SIZE>(ti));
442:     t = tt;
443: #else
444:     t = breduce<PetscReal, BLOCK_SIZE>(PetscRealPart(t));
445: #endif
446:     pLi -= col - glbDD; // diagonal
447:     if (threadIdx.x == 0) {
448:       x[glbDD] -= t;
449:       x[glbDD] /= pLi[0];
450:     }
451:     __syncthreads();
452:     // inc past L to start of previous U
453:     pLi -= bw + 1;
454:     if (glbDD < bw) pLi += bw - glbDD;                                    // overshot in top left corner
455:     if (((locDD + 1) < bw) && field != Nf - 1) pLi -= (bw - (locDD + 1)); // skip past right corner
456:   }
457: }

459: static PetscErrorCode MatSolve_SeqAIJCUSPARSEBAND(Mat A, Vec bb, Vec xx)
460: {
461:   const PetscScalar                    *barray;
462:   PetscScalar                          *xarray;
463:   thrust::device_ptr<const PetscScalar> bGPU;
464:   thrust::device_ptr<PetscScalar>       xGPU;
465:   Mat_SeqAIJCUSPARSETriFactors         *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors *)A->spptr;
466:   THRUSTARRAY                          *tempGPU            = (THRUSTARRAY *)cusparseTriFactors->workVector;
467:   PetscInt                              n = A->rmap->n, nz = cusparseTriFactors->nnz, Nf = 1;                                                                                  // Nf is batch size - not used
468:   PetscInt                              bw = (int)(2. * (double)n - 1. - (double)(PetscSqrtReal(1. + 4. * ((double)n * (double)n - (double)nz)) + PETSC_MACHINE_EPSILON)) / 2; // quadric formula for bandwidth

470:   if (A->rmap->n == 0) return 0;

472:   /* Get the GPU pointers */
473:   VecCUDAGetArrayWrite(xx, &xarray);
474:   VecCUDAGetArrayRead(bb, &barray);
475:   xGPU = thrust::device_pointer_cast(xarray);
476:   bGPU = thrust::device_pointer_cast(barray);

478:   PetscLogGpuTimeBegin();
479:   /* First, reorder with the row permutation */
480:   thrust::copy(thrust::cuda::par.on(PetscDefaultCudaStream), thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()), thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->end()), tempGPU->begin());
481:   constexpr int block = 128;
482:   mat_solve_band<block><<<Nf, block>>>(n, bw, cusparseTriFactors->a_band_d, tempGPU->data().get());
483:   PetscCUDACheckLaunch; // does a sync

485:   /* Last, reorder with the column permutation */
486:   thrust::copy(thrust::cuda::par.on(PetscDefaultCudaStream), thrust::make_permutation_iterator(tempGPU->begin(), cusparseTriFactors->cpermIndices->begin()), thrust::make_permutation_iterator(tempGPU->begin(), cusparseTriFactors->cpermIndices->end()), xGPU);

488:   VecCUDARestoreArrayRead(bb, &barray);
489:   VecCUDARestoreArrayWrite(xx, &xarray);
490:   PetscLogGpuFlops(2.0 * cusparseTriFactors->nnz - A->cmap->n);
491:   PetscLogGpuTimeEnd();

493:   return 0;
494: }