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