1d4002b98SHong Zhang 26524c165SJacob Faibussowitsch #ifndef __SELL_H 3d4002b98SHong Zhang #define __SELL_H 4d4002b98SHong Zhang 5d4002b98SHong Zhang #include <petsc/private/matimpl.h> 6eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h> 7d4002b98SHong Zhang 8d4002b98SHong Zhang /* 9d4002b98SHong Zhang Struct header for SeqSELL matrix format 10d4002b98SHong Zhang */ 11d4002b98SHong Zhang #define SEQSELLHEADER(datatype) \ 12d4002b98SHong Zhang PetscBool roworiented; /* if true, row-oriented input, default */ \ 13d4002b98SHong Zhang PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \ 14d4002b98SHong Zhang PetscInt nounused; /* -1 generate error on unused space */ \ 15d4002b98SHong Zhang PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \ 16d4002b98SHong Zhang PetscInt maxallocmat; /* max allocated space for the matrix */ \ 17d4002b98SHong Zhang PetscInt maxallocrow; /* max allocated space for each row */ \ 18d4002b98SHong Zhang PetscInt nz; /* actual nonzeros */ \ 19d4002b98SHong Zhang PetscInt rlenmax; /* max actual row length, rmax cannot exceed maxallocrow */ \ 20d4002b98SHong Zhang PetscInt *rlen; /* actual length of each row (padding zeros excluded) */ \ 21d4002b98SHong Zhang PetscBool free_rlen; /* free rlen array ? */ \ 22d4002b98SHong Zhang PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 23d4002b98SHong Zhang as more values are set than were prealloced */ \ 24d4002b98SHong Zhang PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \ 25d4002b98SHong Zhang PetscBool ignorezeroentries; \ 26d4002b98SHong Zhang PetscBool free_colidx; /* free the column indices colidx when the matrix is destroyed */ \ 27d4002b98SHong Zhang PetscBool free_val; /* free the numerical values when matrix is destroy */ \ 28d4002b98SHong Zhang PetscInt *colidx; /* column index */ \ 29d4002b98SHong Zhang PetscInt *diag; /* pointers to diagonal elements */ \ 30d4002b98SHong Zhang PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 31d4002b98SHong Zhang PetscBool free_diag; /* free diag ? */ \ 32d4002b98SHong Zhang datatype *val; /* elements including nonzeros and padding zeros */ \ 33d4002b98SHong Zhang PetscScalar *solve_work; /* work space used in MatSolve */ \ 34d4002b98SHong Zhang IS row, col, icol; /* index sets, used for reorderings */ \ 35d4002b98SHong Zhang PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 36d4002b98SHong Zhang Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); 37d4002b98SHong Zhang means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \ 38d4002b98SHong Zhang PetscInt *sliidx; /* slice index */ \ 396108893eSStefano Zampini PetscInt totalslices; /* total number of slices */ \ 40*07e43b41SHong Zhang PetscInt sliceheight; /* slice height */ \ 41*07e43b41SHong Zhang PetscReal fillratio; /* ratio of number of padded zeros over total number of elements */ \ 42*07e43b41SHong Zhang PetscReal avgslicewidth; /* average slice width */ \ 43*07e43b41SHong Zhang PetscInt maxslicewidth; /* maximum slice width */ \ 44*07e43b41SHong Zhang PetscInt *sliperm; /* slice permutation array, CUDA only */ \ 45*07e43b41SHong Zhang PetscInt totalblocks; /* total number of blocks, CUDA only */ \ 46*07e43b41SHong Zhang PetscInt *blockidx; /* block index, CUDA only */ \ 47*07e43b41SHong Zhang PetscInt *block_row_map; /* starting row of the current block, CUDA only */ \ 486108893eSStefano Zampini PetscInt *getrowcols; /* workarray for MatGetRow_SeqSELL */ \ 499371c9d4SSatish Balay PetscScalar *getrowvals /* workarray for MatGetRow_SeqSELL */ 50d4002b98SHong Zhang 51d4002b98SHong Zhang typedef struct { 52d4002b98SHong Zhang SEQSELLHEADER(MatScalar); 53d4002b98SHong Zhang MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 54d4002b98SHong Zhang PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 55d4002b98SHong Zhang PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 56d4002b98SHong Zhang PetscScalar fshift, omega; /* last used omega and fshift */ 57d4002b98SHong Zhang ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 58d4002b98SHong Zhang } Mat_SeqSELL; 59d4002b98SHong Zhang 60d4002b98SHong Zhang /* 61d4002b98SHong Zhang Frees the arrays from the XSELLPACK matrix type 62d4002b98SHong Zhang */ 63d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSeqXSELLFreeSELL(Mat AA, MatScalar **val, PetscInt **colidx) 64d71ae5a4SJacob Faibussowitsch { 65d4002b98SHong Zhang Mat_SeqSELL *A = (Mat_SeqSELL *)AA->data; 66d4002b98SHong Zhang if (A->singlemalloc) { 679566063dSJacob Faibussowitsch PetscCall(PetscFree2(*val, *colidx)); 68d4002b98SHong Zhang } else { 699566063dSJacob Faibussowitsch if (A->free_val) PetscCall(PetscFree(*val)); 709566063dSJacob Faibussowitsch if (A->free_colidx) PetscCall(PetscFree(*colidx)); 71d4002b98SHong Zhang } 723ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 73d4002b98SHong Zhang } 74d4002b98SHong Zhang 75*07e43b41SHong Zhang #define MatSeqXSELLReallocateSELL(Amat, AM, BS2, WIDTH, SIDX, SH, SID, ROW, COL, COLIDX, VAL, CP, VP, NONEW, datatype) \ 76*07e43b41SHong Zhang if (WIDTH >= (SIDX[SID + 1] - SIDX[SID]) / SH) { \ 77d4002b98SHong Zhang Mat_SeqSELL *Ain = (Mat_SeqSELL *)Amat->data; \ 782d1451d4SHong Zhang /* there is no extra room in row, therefore enlarge 1 slice column */ \ 79*07e43b41SHong Zhang PetscInt new_size = Ain->maxallocmat + SH, *new_colidx; \ 80d4002b98SHong Zhang datatype *new_val; \ 81d4002b98SHong Zhang \ 8208401ef6SPierre Jolivet PetscCheck(NONEW != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", ROW, COL); \ 83d4002b98SHong Zhang /* malloc new storage space */ \ 849566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(BS2 *new_size, &new_val, BS2 *new_size, &new_colidx)); \ 85d4002b98SHong Zhang \ 86d4002b98SHong Zhang /* copy over old data into new slots by two steps: one step for data before the current slice and the other for the rest */ \ 879566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_val, VAL, SIDX[SID + 1])); \ 889566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_colidx, COLIDX, SIDX[SID + 1])); \ 89*07e43b41SHong Zhang PetscCall(PetscArraycpy(new_val + SIDX[SID + 1] + SH, VAL + SIDX[SID + 1], SIDX[Ain->totalslices] - SIDX[SID + 1])); \ 90*07e43b41SHong Zhang PetscCall(PetscArraycpy(new_colidx + SIDX[SID + 1] + SH, COLIDX + SIDX[SID + 1], SIDX[Ain->totalslices] - SIDX[SID + 1])); \ 91d4002b98SHong Zhang /* update slice_idx */ \ 92*07e43b41SHong Zhang for (ii = SID + 1; ii <= Ain->totalslices; ii++) { SIDX[ii] += SH; } \ 932d1451d4SHong Zhang /* update pointers. Notice that they point to the FIRST postion of the row */ \ 94*07e43b41SHong Zhang CP = new_colidx + SIDX[SID] + (ROW % SH); \ 95*07e43b41SHong Zhang VP = new_val + SIDX[SID] + (ROW % SH); \ 96d4002b98SHong Zhang /* free up old matrix storage */ \ 979566063dSJacob Faibussowitsch PetscCall(MatSeqXSELLFreeSELL(A, &Ain->val, &Ain->colidx)); \ 98d4002b98SHong Zhang Ain->val = (MatScalar *)new_val; \ 99d4002b98SHong Zhang Ain->colidx = new_colidx; \ 100d4002b98SHong Zhang Ain->singlemalloc = PETSC_TRUE; \ 101d4002b98SHong Zhang Ain->maxallocmat = new_size; \ 102d4002b98SHong Zhang Ain->reallocs++; \ 103d4002b98SHong Zhang if (WIDTH >= Ain->maxallocrow) Ain->maxallocrow++; \ 104d4002b98SHong Zhang if (WIDTH >= Ain->rlenmax) Ain->rlenmax++; \ 1059371c9d4SSatish Balay } 106d4002b98SHong Zhang 107d4002b98SHong Zhang #define MatSetValue_SeqSELL_Private(A, row, col, value, addv, orow, ocol, cp, vp, lastcol, low, high) \ 108d4002b98SHong Zhang { \ 109d4002b98SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; \ 110d4002b98SHong Zhang found = PETSC_FALSE; \ 111d4002b98SHong Zhang if (col <= lastcol) low = 0; \ 112d4002b98SHong Zhang else high = a->rlen[row]; \ 113d4002b98SHong Zhang lastcol = col; \ 114d4002b98SHong Zhang while (high - low > 5) { \ 115d4002b98SHong Zhang t = (low + high) / 2; \ 116*07e43b41SHong Zhang if (*(cp + a->sliceheight * t) > col) high = t; \ 117d4002b98SHong Zhang else low = t; \ 118d4002b98SHong Zhang } \ 119d4002b98SHong Zhang for (_i = low; _i < high; _i++) { \ 120*07e43b41SHong Zhang if (*(cp + a->sliceheight * _i) > col) break; \ 121*07e43b41SHong Zhang if (*(cp + a->sliceheight * _i) == col) { \ 122*07e43b41SHong Zhang if (addv == ADD_VALUES) *(vp + a->sliceheight * _i) += value; \ 123*07e43b41SHong Zhang else *(vp + a->sliceheight * _i) = value; \ 124d4002b98SHong Zhang found = PETSC_TRUE; \ 125d4002b98SHong Zhang break; \ 126d4002b98SHong Zhang } \ 127d4002b98SHong Zhang } \ 128d4002b98SHong Zhang if (!found) { \ 12908401ef6SPierre Jolivet PetscCheck(a->nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \ 130*07e43b41SHong Zhang if (a->nonew != 1 && !(value == 0.0 && a->ignorezeroentries) && a->rlen[row] >= (a->sliidx[row / a->sliceheight + 1] - a->sliidx[row / a->sliceheight]) / a->sliceheight) { \ 1312d1451d4SHong Zhang /* there is no extra room in row, therefore enlarge 1 slice column */ \ 132*07e43b41SHong Zhang if (a->maxallocmat < a->sliidx[a->totalslices] + a->sliceheight) { \ 133d4002b98SHong Zhang /* allocates a larger array for the XSELL matrix types; only extend the current slice by one more column. */ \ 134*07e43b41SHong Zhang PetscInt new_size = a->maxallocmat + a->sliceheight, *new_colidx; \ 135d4002b98SHong Zhang MatScalar *new_val; \ 13608401ef6SPierre Jolivet PetscCheck(a->nonew != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", orow, ocol); \ 137d4002b98SHong Zhang /* malloc new storage space */ \ 1389566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(new_size, &new_val, new_size, &new_colidx)); \ 139d4002b98SHong Zhang /* copy over old data into new slots by two steps: one step for data before the current slice and the other for the rest */ \ 140*07e43b41SHong Zhang PetscCall(PetscArraycpy(new_val, a->val, a->sliidx[row / a->sliceheight + 1])); \ 141*07e43b41SHong Zhang PetscCall(PetscArraycpy(new_colidx, a->colidx, a->sliidx[row / a->sliceheight + 1])); \ 142*07e43b41SHong Zhang PetscCall(PetscArraycpy(new_val + a->sliidx[row / a->sliceheight + 1] + a->sliceheight, a->val + a->sliidx[row / a->sliceheight + 1], a->sliidx[a->totalslices] - a->sliidx[row / a->sliceheight + 1])); \ 143*07e43b41SHong Zhang PetscCall(PetscArraycpy(new_colidx + a->sliidx[row / a->sliceheight + 1] + a->sliceheight, a->colidx + a->sliidx[row / a->sliceheight + 1], a->sliidx[a->totalslices] - a->sliidx[row / a->sliceheight + 1])); \ 1442d1451d4SHong Zhang /* update pointers. Notice that they point to the FIRST postion of the row */ \ 145*07e43b41SHong Zhang cp = new_colidx + a->sliidx[row / a->sliceheight] + (row % a->sliceheight); \ 146*07e43b41SHong Zhang vp = new_val + a->sliidx[row / a->sliceheight] + (row % a->sliceheight); \ 147d4002b98SHong Zhang /* free up old matrix storage */ \ 1489566063dSJacob Faibussowitsch PetscCall(MatSeqXSELLFreeSELL(A, &a->val, &a->colidx)); \ 149d4002b98SHong Zhang a->val = (MatScalar *)new_val; \ 150d4002b98SHong Zhang a->colidx = new_colidx; \ 151d4002b98SHong Zhang a->singlemalloc = PETSC_TRUE; \ 152d4002b98SHong Zhang a->maxallocmat = new_size; \ 153d4002b98SHong Zhang a->reallocs++; \ 154d4002b98SHong Zhang } else { \ 155d4002b98SHong Zhang /* no need to reallocate, just shift the following slices to create space for the added slice column */ \ 156*07e43b41SHong Zhang PetscCall(PetscArraymove(a->val + a->sliidx[row / a->sliceheight + 1] + a->sliceheight, a->val + a->sliidx[row / a->sliceheight + 1], a->sliidx[a->totalslices] - a->sliidx[row / a->sliceheight + 1])); \ 157*07e43b41SHong Zhang PetscCall(PetscArraymove(a->colidx + a->sliidx[row / a->sliceheight + 1] + a->sliceheight, a->colidx + a->sliidx[row / a->sliceheight + 1], a->sliidx[a->totalslices] - a->sliidx[row / a->sliceheight + 1])); \ 158d4002b98SHong Zhang } \ 159d4002b98SHong Zhang /* update slice_idx */ \ 160*07e43b41SHong Zhang for (ii = row / a->sliceheight + 1; ii <= a->totalslices; ii++) a->sliidx[ii] += a->sliceheight; \ 161d4002b98SHong Zhang if (a->rlen[row] >= a->maxallocrow) a->maxallocrow++; \ 162d4002b98SHong Zhang if (a->rlen[row] >= a->rlenmax) a->rlenmax++; \ 163d4002b98SHong Zhang } \ 164d4002b98SHong Zhang /* shift up all the later entries in this row */ \ 165d4002b98SHong Zhang for (ii = a->rlen[row] - 1; ii >= _i; ii--) { \ 166*07e43b41SHong Zhang *(cp + a->sliceheight * (ii + 1)) = *(cp + a->sliceheight * ii); \ 167*07e43b41SHong Zhang *(vp + a->sliceheight * (ii + 1)) = *(vp + a->sliceheight * ii); \ 168d4002b98SHong Zhang } \ 169*07e43b41SHong Zhang *(cp + a->sliceheight * _i) = col; \ 170*07e43b41SHong Zhang *(vp + a->sliceheight * _i) = value; \ 1719371c9d4SSatish Balay a->nz++; \ 1729371c9d4SSatish Balay a->rlen[row]++; \ 1739371c9d4SSatish Balay A->nonzerostate++; \ 1749371c9d4SSatish Balay low = _i + 1; \ 1759371c9d4SSatish Balay high++; \ 176d4002b98SHong Zhang } \ 1779371c9d4SSatish Balay } 178d4002b98SHong Zhang 179d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSeqSELLSetPreallocation_SeqSELL(Mat, PetscInt, const PetscInt[]); 180d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatMult_SeqSELL(Mat, Vec, Vec); 181d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatMultAdd_SeqSELL(Mat, Vec, Vec, Vec); 182d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatMultTranspose_SeqSELL(Mat, Vec, Vec); 183d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqSELL(Mat, Vec, Vec, Vec); 184d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqSELL(Mat, PetscBool *, PetscInt *); 185d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqSELL(Mat); 186d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatInvertDiagonal_SeqSELL(Mat, PetscScalar, PetscScalar); 187d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatZeroEntries_SeqSELL(Mat); 188d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqSELL(Mat); 189d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSetOption_SeqSELL(Mat, MatOption, PetscBool); 190d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatGetDiagonal_SeqSELL(Mat, Vec v); 191d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatGetValues_SeqSELL(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]); 192d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatView_SeqSELL(Mat, PetscViewer); 193d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqSELL(Mat, MatAssemblyType); 194d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatGetInfo_SeqSELL(Mat, MatInfoType, MatInfo *); 195d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSetValues_SeqSELL(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 196d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatCopy_SeqSELL(Mat, Mat, MatStructure); 197d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSetUp_SeqSELL(Mat); 198d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSeqSELLGetArray_SeqSELL(Mat, PetscScalar *[]); 199d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSeqSELLRestoreArray_SeqSELL(Mat, PetscScalar *[]); 200d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatShift_SeqSELL(Mat, PetscScalar); 201d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSOR_SeqSELL(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 202d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreate_SeqSELL(Mat); 203d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatDuplicate_SeqSELL(Mat, MatDuplicateOption, Mat *); 204d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqSELL(Mat, Mat, PetscBool *); 205d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatSeqSELLInvalidateDiagonal(Mat); 206d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqSELL_SeqAIJ(Mat, MatType, MatReuse, Mat *); 207d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat, MatType, MatReuse, Mat *); 208d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqSELL(Mat, ISColoring, MatFDColoring); 209d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqSELL(Mat, ISColoring, MatFDColoring); 210d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqSELL_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 211d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqSELL_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 212d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatConjugate_SeqSELL(Mat A); 213d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatScale_SeqSELL(Mat, PetscScalar); 214d4002b98SHong Zhang PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqSELL(Mat, Vec, Vec); 215d4002b98SHong Zhang #endif 216