1 /* 2 Defines the basic matrix operations for the AIJ (compressed row) 3 matrix storage format. 4 */ 5 6 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 7 #include <petscblaslapack.h> 8 #include <petscbt.h> 9 #include <petsc/private/kernels/blocktranspose.h> 10 11 /* defines MatSetValues_Seq_Hash(), MatAssemblyEnd_Seq_Hash(), MatSetUp_Seq_Hash() */ 12 #define TYPE AIJ 13 #define TYPE_BS 14 #include "../src/mat/impls/aij/seq/seqhashmatsetvalues.h" 15 #include "../src/mat/impls/aij/seq/seqhashmat.h" 16 #undef TYPE 17 #undef TYPE_BS 18 19 PetscErrorCode MatSeqAIJSetTypeFromOptions(Mat A) 20 { 21 PetscBool flg; 22 char type[256]; 23 24 PetscFunctionBegin; 25 PetscObjectOptionsBegin((PetscObject)A); 26 PetscCall(PetscOptionsFList("-mat_seqaij_type", "Matrix SeqAIJ type", "MatSeqAIJSetType", MatSeqAIJList, "seqaij", type, 256, &flg)); 27 if (flg) PetscCall(MatSeqAIJSetType(A, type)); 28 PetscOptionsEnd(); 29 PetscFunctionReturn(PETSC_SUCCESS); 30 } 31 32 PetscErrorCode MatGetColumnReductions_SeqAIJ(Mat A, PetscInt type, PetscReal *reductions) 33 { 34 PetscInt i, m, n; 35 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 36 37 PetscFunctionBegin; 38 PetscCall(MatGetSize(A, &m, &n)); 39 PetscCall(PetscArrayzero(reductions, n)); 40 if (type == NORM_2) { 41 for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscAbsScalar(aij->a[i] * aij->a[i]); 42 } else if (type == NORM_1) { 43 for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscAbsScalar(aij->a[i]); 44 } else if (type == NORM_INFINITY) { 45 for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] = PetscMax(PetscAbsScalar(aij->a[i]), reductions[aij->j[i]]); 46 } else if (type == REDUCTION_SUM_REALPART || type == REDUCTION_MEAN_REALPART) { 47 for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscRealPart(aij->a[i]); 48 } else if (type == REDUCTION_SUM_IMAGINARYPART || type == REDUCTION_MEAN_IMAGINARYPART) { 49 for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscImaginaryPart(aij->a[i]); 50 } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Unknown reduction type"); 51 52 if (type == NORM_2) { 53 for (i = 0; i < n; i++) reductions[i] = PetscSqrtReal(reductions[i]); 54 } else if (type == REDUCTION_MEAN_REALPART || type == REDUCTION_MEAN_IMAGINARYPART) { 55 for (i = 0; i < n; i++) reductions[i] /= m; 56 } 57 PetscFunctionReturn(PETSC_SUCCESS); 58 } 59 60 PetscErrorCode MatFindOffBlockDiagonalEntries_SeqAIJ(Mat A, IS *is) 61 { 62 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 63 PetscInt i, m = A->rmap->n, cnt = 0, bs = A->rmap->bs; 64 const PetscInt *jj = a->j, *ii = a->i; 65 PetscInt *rows; 66 67 PetscFunctionBegin; 68 for (i = 0; i < m; i++) { 69 if ((ii[i] != ii[i + 1]) && ((jj[ii[i]] < bs * (i / bs)) || (jj[ii[i + 1] - 1] > bs * ((i + bs) / bs) - 1))) cnt++; 70 } 71 PetscCall(PetscMalloc1(cnt, &rows)); 72 cnt = 0; 73 for (i = 0; i < m; i++) { 74 if ((ii[i] != ii[i + 1]) && ((jj[ii[i]] < bs * (i / bs)) || (jj[ii[i + 1] - 1] > bs * ((i + bs) / bs) - 1))) { 75 rows[cnt] = i; 76 cnt++; 77 } 78 } 79 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cnt, rows, PETSC_OWN_POINTER, is)); 80 PetscFunctionReturn(PETSC_SUCCESS); 81 } 82 83 PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat A, PetscInt *nrows, PetscInt **zrows) 84 { 85 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 86 const MatScalar *aa; 87 PetscInt i, m = A->rmap->n, cnt = 0; 88 const PetscInt *ii = a->i, *jj = a->j, *diag; 89 PetscInt *rows; 90 91 PetscFunctionBegin; 92 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 93 PetscCall(MatMarkDiagonal_SeqAIJ(A)); 94 diag = a->diag; 95 for (i = 0; i < m; i++) { 96 if ((diag[i] >= ii[i + 1]) || (jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) cnt++; 97 } 98 PetscCall(PetscMalloc1(cnt, &rows)); 99 cnt = 0; 100 for (i = 0; i < m; i++) { 101 if ((diag[i] >= ii[i + 1]) || (jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) rows[cnt++] = i; 102 } 103 *nrows = cnt; 104 *zrows = rows; 105 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 106 PetscFunctionReturn(PETSC_SUCCESS); 107 } 108 109 PetscErrorCode MatFindZeroDiagonals_SeqAIJ(Mat A, IS *zrows) 110 { 111 PetscInt nrows, *rows; 112 113 PetscFunctionBegin; 114 *zrows = NULL; 115 PetscCall(MatFindZeroDiagonals_SeqAIJ_Private(A, &nrows, &rows)); 116 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), nrows, rows, PETSC_OWN_POINTER, zrows)); 117 PetscFunctionReturn(PETSC_SUCCESS); 118 } 119 120 PetscErrorCode MatFindNonzeroRows_SeqAIJ(Mat A, IS *keptrows) 121 { 122 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 123 const MatScalar *aa; 124 PetscInt m = A->rmap->n, cnt = 0; 125 const PetscInt *ii; 126 PetscInt n, i, j, *rows; 127 128 PetscFunctionBegin; 129 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 130 *keptrows = NULL; 131 ii = a->i; 132 for (i = 0; i < m; i++) { 133 n = ii[i + 1] - ii[i]; 134 if (!n) { 135 cnt++; 136 goto ok1; 137 } 138 for (j = ii[i]; j < ii[i + 1]; j++) { 139 if (aa[j] != 0.0) goto ok1; 140 } 141 cnt++; 142 ok1:; 143 } 144 if (!cnt) { 145 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 146 PetscFunctionReturn(PETSC_SUCCESS); 147 } 148 PetscCall(PetscMalloc1(A->rmap->n - cnt, &rows)); 149 cnt = 0; 150 for (i = 0; i < m; i++) { 151 n = ii[i + 1] - ii[i]; 152 if (!n) continue; 153 for (j = ii[i]; j < ii[i + 1]; j++) { 154 if (aa[j] != 0.0) { 155 rows[cnt++] = i; 156 break; 157 } 158 } 159 } 160 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 161 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cnt, rows, PETSC_OWN_POINTER, keptrows)); 162 PetscFunctionReturn(PETSC_SUCCESS); 163 } 164 165 PetscErrorCode MatDiagonalSet_SeqAIJ(Mat Y, Vec D, InsertMode is) 166 { 167 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)Y->data; 168 PetscInt i, m = Y->rmap->n; 169 const PetscInt *diag; 170 MatScalar *aa; 171 const PetscScalar *v; 172 PetscBool missing; 173 174 PetscFunctionBegin; 175 if (Y->assembled) { 176 PetscCall(MatMissingDiagonal_SeqAIJ(Y, &missing, NULL)); 177 if (!missing) { 178 diag = aij->diag; 179 PetscCall(VecGetArrayRead(D, &v)); 180 PetscCall(MatSeqAIJGetArray(Y, &aa)); 181 if (is == INSERT_VALUES) { 182 for (i = 0; i < m; i++) aa[diag[i]] = v[i]; 183 } else { 184 for (i = 0; i < m; i++) aa[diag[i]] += v[i]; 185 } 186 PetscCall(MatSeqAIJRestoreArray(Y, &aa)); 187 PetscCall(VecRestoreArrayRead(D, &v)); 188 PetscFunctionReturn(PETSC_SUCCESS); 189 } 190 PetscCall(MatSeqAIJInvalidateDiagonal(Y)); 191 } 192 PetscCall(MatDiagonalSet_Default(Y, D, is)); 193 PetscFunctionReturn(PETSC_SUCCESS); 194 } 195 196 PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *m, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 197 { 198 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 199 PetscInt i, ishift; 200 201 PetscFunctionBegin; 202 if (m) *m = A->rmap->n; 203 if (!ia) PetscFunctionReturn(PETSC_SUCCESS); 204 ishift = 0; 205 if (symmetric && A->structurally_symmetric != PETSC_BOOL3_TRUE) { 206 PetscCall(MatToSymmetricIJ_SeqAIJ(A->rmap->n, a->i, a->j, PETSC_TRUE, ishift, oshift, (PetscInt **)ia, (PetscInt **)ja)); 207 } else if (oshift == 1) { 208 PetscInt *tia; 209 PetscInt nz = a->i[A->rmap->n]; 210 /* malloc space and add 1 to i and j indices */ 211 PetscCall(PetscMalloc1(A->rmap->n + 1, &tia)); 212 for (i = 0; i < A->rmap->n + 1; i++) tia[i] = a->i[i] + 1; 213 *ia = tia; 214 if (ja) { 215 PetscInt *tja; 216 PetscCall(PetscMalloc1(nz + 1, &tja)); 217 for (i = 0; i < nz; i++) tja[i] = a->j[i] + 1; 218 *ja = tja; 219 } 220 } else { 221 *ia = a->i; 222 if (ja) *ja = a->j; 223 } 224 PetscFunctionReturn(PETSC_SUCCESS); 225 } 226 227 PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 228 { 229 PetscFunctionBegin; 230 if (!ia) PetscFunctionReturn(PETSC_SUCCESS); 231 if ((symmetric && A->structurally_symmetric != PETSC_BOOL3_TRUE) || oshift == 1) { 232 PetscCall(PetscFree(*ia)); 233 if (ja) PetscCall(PetscFree(*ja)); 234 } 235 PetscFunctionReturn(PETSC_SUCCESS); 236 } 237 238 PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *nn, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 239 { 240 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 241 PetscInt i, *collengths, *cia, *cja, n = A->cmap->n, m = A->rmap->n; 242 PetscInt nz = a->i[m], row, *jj, mr, col; 243 244 PetscFunctionBegin; 245 *nn = n; 246 if (!ia) PetscFunctionReturn(PETSC_SUCCESS); 247 if (symmetric) { 248 PetscCall(MatToSymmetricIJ_SeqAIJ(A->rmap->n, a->i, a->j, PETSC_TRUE, 0, oshift, (PetscInt **)ia, (PetscInt **)ja)); 249 } else { 250 PetscCall(PetscCalloc1(n, &collengths)); 251 PetscCall(PetscMalloc1(n + 1, &cia)); 252 PetscCall(PetscMalloc1(nz, &cja)); 253 jj = a->j; 254 for (i = 0; i < nz; i++) collengths[jj[i]]++; 255 cia[0] = oshift; 256 for (i = 0; i < n; i++) cia[i + 1] = cia[i] + collengths[i]; 257 PetscCall(PetscArrayzero(collengths, n)); 258 jj = a->j; 259 for (row = 0; row < m; row++) { 260 mr = a->i[row + 1] - a->i[row]; 261 for (i = 0; i < mr; i++) { 262 col = *jj++; 263 264 cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 265 } 266 } 267 PetscCall(PetscFree(collengths)); 268 *ia = cia; 269 *ja = cja; 270 } 271 PetscFunctionReturn(PETSC_SUCCESS); 272 } 273 274 PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 275 { 276 PetscFunctionBegin; 277 if (!ia) PetscFunctionReturn(PETSC_SUCCESS); 278 279 PetscCall(PetscFree(*ia)); 280 PetscCall(PetscFree(*ja)); 281 PetscFunctionReturn(PETSC_SUCCESS); 282 } 283 284 /* 285 MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 286 MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 287 spidx[], index of a->a, to be used in MatTransposeColoringCreate_SeqAIJ() and MatFDColoringCreate_SeqXAIJ() 288 */ 289 PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *nn, const PetscInt *ia[], const PetscInt *ja[], PetscInt *spidx[], PetscBool *done) 290 { 291 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 292 PetscInt i, *collengths, *cia, *cja, n = A->cmap->n, m = A->rmap->n; 293 PetscInt nz = a->i[m], row, mr, col, tmp; 294 PetscInt *cspidx; 295 const PetscInt *jj; 296 297 PetscFunctionBegin; 298 *nn = n; 299 if (!ia) PetscFunctionReturn(PETSC_SUCCESS); 300 301 PetscCall(PetscCalloc1(n, &collengths)); 302 PetscCall(PetscMalloc1(n + 1, &cia)); 303 PetscCall(PetscMalloc1(nz, &cja)); 304 PetscCall(PetscMalloc1(nz, &cspidx)); 305 jj = a->j; 306 for (i = 0; i < nz; i++) collengths[jj[i]]++; 307 cia[0] = oshift; 308 for (i = 0; i < n; i++) cia[i + 1] = cia[i] + collengths[i]; 309 PetscCall(PetscArrayzero(collengths, n)); 310 jj = a->j; 311 for (row = 0; row < m; row++) { 312 mr = a->i[row + 1] - a->i[row]; 313 for (i = 0; i < mr; i++) { 314 col = *jj++; 315 tmp = cia[col] + collengths[col]++ - oshift; 316 cspidx[tmp] = a->i[row] + i; /* index of a->j */ 317 cja[tmp] = row + oshift; 318 } 319 } 320 PetscCall(PetscFree(collengths)); 321 *ia = cia; 322 *ja = cja; 323 *spidx = cspidx; 324 PetscFunctionReturn(PETSC_SUCCESS); 325 } 326 327 PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscInt *spidx[], PetscBool *done) 328 { 329 PetscFunctionBegin; 330 PetscCall(MatRestoreColumnIJ_SeqAIJ(A, oshift, symmetric, inodecompressed, n, ia, ja, done)); 331 PetscCall(PetscFree(*spidx)); 332 PetscFunctionReturn(PETSC_SUCCESS); 333 } 334 335 PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A, PetscInt row, const PetscScalar v[]) 336 { 337 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 338 PetscInt *ai = a->i; 339 PetscScalar *aa; 340 341 PetscFunctionBegin; 342 PetscCall(MatSeqAIJGetArray(A, &aa)); 343 PetscCall(PetscArraycpy(aa + ai[row], v, ai[row + 1] - ai[row])); 344 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 345 PetscFunctionReturn(PETSC_SUCCESS); 346 } 347 348 /* 349 MatSeqAIJSetValuesLocalFast - An optimized version of MatSetValuesLocal() for SeqAIJ matrices with several assumptions 350 351 - a single row of values is set with each call 352 - no row or column indices are negative or (in error) larger than the number of rows or columns 353 - the values are always added to the matrix, not set 354 - no new locations are introduced in the nonzero structure of the matrix 355 356 This does NOT assume the global column indices are sorted 357 358 */ 359 360 #include <petsc/private/isimpl.h> 361 PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is) 362 { 363 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 364 PetscInt low, high, t, row, nrow, i, col, l; 365 const PetscInt *rp, *ai = a->i, *ailen = a->ilen, *aj = a->j; 366 PetscInt lastcol = -1; 367 MatScalar *ap, value, *aa; 368 const PetscInt *ridx = A->rmap->mapping->indices, *cidx = A->cmap->mapping->indices; 369 370 PetscFunctionBegin; 371 PetscCall(MatSeqAIJGetArray(A, &aa)); 372 row = ridx[im[0]]; 373 rp = aj + ai[row]; 374 ap = aa + ai[row]; 375 nrow = ailen[row]; 376 low = 0; 377 high = nrow; 378 for (l = 0; l < n; l++) { /* loop over added columns */ 379 col = cidx[in[l]]; 380 value = v[l]; 381 382 if (col <= lastcol) low = 0; 383 else high = nrow; 384 lastcol = col; 385 while (high - low > 5) { 386 t = (low + high) / 2; 387 if (rp[t] > col) high = t; 388 else low = t; 389 } 390 for (i = low; i < high; i++) { 391 if (rp[i] == col) { 392 ap[i] += value; 393 low = i + 1; 394 break; 395 } 396 } 397 } 398 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 399 return PETSC_SUCCESS; 400 } 401 402 PetscErrorCode MatSetValues_SeqAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is) 403 { 404 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 405 PetscInt *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N; 406 PetscInt *imax = a->imax, *ai = a->i, *ailen = a->ilen; 407 PetscInt *aj = a->j, nonew = a->nonew, lastcol = -1; 408 MatScalar *ap = NULL, value = 0.0, *aa; 409 PetscBool ignorezeroentries = a->ignorezeroentries; 410 PetscBool roworiented = a->roworiented; 411 412 PetscFunctionBegin; 413 PetscCall(MatSeqAIJGetArray(A, &aa)); 414 for (k = 0; k < m; k++) { /* loop over added rows */ 415 row = im[k]; 416 if (row < 0) continue; 417 PetscCheck(row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->n - 1); 418 rp = aj + ai[row]; 419 if (!A->structure_only) ap = aa + ai[row]; 420 rmax = imax[row]; 421 nrow = ailen[row]; 422 low = 0; 423 high = nrow; 424 for (l = 0; l < n; l++) { /* loop over added columns */ 425 if (in[l] < 0) continue; 426 PetscCheck(in[l] < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[l], A->cmap->n - 1); 427 col = in[l]; 428 if (v && !A->structure_only) value = roworiented ? v[l + k * n] : v[k + l * m]; 429 if (!A->structure_only && value == 0.0 && ignorezeroentries && is == ADD_VALUES && row != col) continue; 430 431 if (col <= lastcol) low = 0; 432 else high = nrow; 433 lastcol = col; 434 while (high - low > 5) { 435 t = (low + high) / 2; 436 if (rp[t] > col) high = t; 437 else low = t; 438 } 439 for (i = low; i < high; i++) { 440 if (rp[i] > col) break; 441 if (rp[i] == col) { 442 if (!A->structure_only) { 443 if (is == ADD_VALUES) { 444 ap[i] += value; 445 (void)PetscLogFlops(1.0); 446 } else ap[i] = value; 447 } 448 low = i + 1; 449 goto noinsert; 450 } 451 } 452 if (value == 0.0 && ignorezeroentries && row != col) goto noinsert; 453 if (nonew == 1) goto noinsert; 454 PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") in the matrix", row, col); 455 if (A->structure_only) { 456 MatSeqXAIJReallocateAIJ_structure_only(A, A->rmap->n, 1, nrow, row, col, rmax, ai, aj, rp, imax, nonew, MatScalar); 457 } else { 458 MatSeqXAIJReallocateAIJ(A, A->rmap->n, 1, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar); 459 } 460 N = nrow++ - 1; 461 a->nz++; 462 high++; 463 /* shift up all the later entries in this row */ 464 PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1)); 465 rp[i] = col; 466 if (!A->structure_only) { 467 PetscCall(PetscArraymove(ap + i + 1, ap + i, N - i + 1)); 468 ap[i] = value; 469 } 470 low = i + 1; 471 A->nonzerostate++; 472 noinsert:; 473 } 474 ailen[row] = nrow; 475 } 476 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 477 PetscFunctionReturn(PETSC_SUCCESS); 478 } 479 480 PetscErrorCode MatSetValues_SeqAIJ_SortedFullNoPreallocation(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is) 481 { 482 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 483 PetscInt *rp, k, row; 484 PetscInt *ai = a->i; 485 PetscInt *aj = a->j; 486 MatScalar *aa, *ap; 487 488 PetscFunctionBegin; 489 PetscCheck(!A->was_assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot call on assembled matrix."); 490 PetscCheck(m * n + a->nz <= a->maxnz, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of entries in matrix will be larger than maximum nonzeros allocated for %" PetscInt_FMT " in MatSeqAIJSetTotalPreallocation()", a->maxnz); 491 492 PetscCall(MatSeqAIJGetArray(A, &aa)); 493 for (k = 0; k < m; k++) { /* loop over added rows */ 494 row = im[k]; 495 rp = aj + ai[row]; 496 ap = aa + ai[row]; 497 498 PetscCall(PetscMemcpy(rp, in, n * sizeof(PetscInt))); 499 if (!A->structure_only) { 500 if (v) { 501 PetscCall(PetscMemcpy(ap, v, n * sizeof(PetscScalar))); 502 v += n; 503 } else { 504 PetscCall(PetscMemzero(ap, n * sizeof(PetscScalar))); 505 } 506 } 507 a->ilen[row] = n; 508 a->imax[row] = n; 509 a->i[row + 1] = a->i[row] + n; 510 a->nz += n; 511 } 512 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 513 PetscFunctionReturn(PETSC_SUCCESS); 514 } 515 516 /*@ 517 MatSeqAIJSetTotalPreallocation - Sets an upper bound on the total number of expected nonzeros in the matrix. 518 519 Input Parameters: 520 + A - the `MATSEQAIJ` matrix 521 - nztotal - bound on the number of nonzeros 522 523 Level: advanced 524 525 Notes: 526 This can be called if you will be provided the matrix row by row (from row zero) with sorted column indices for each row. 527 Simply call `MatSetValues()` after this call to provide the matrix entries in the usual manner. This matrix may be used 528 as always with multiple matrix assemblies. 529 530 .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MAT_SORTED_FULL`, `MatSetValues()`, `MatSeqAIJSetPreallocation()` 531 @*/ 532 PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat A, PetscInt nztotal) 533 { 534 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 535 536 PetscFunctionBegin; 537 PetscCall(PetscLayoutSetUp(A->rmap)); 538 PetscCall(PetscLayoutSetUp(A->cmap)); 539 a->maxnz = nztotal; 540 if (!a->imax) { PetscCall(PetscMalloc1(A->rmap->n, &a->imax)); } 541 if (!a->ilen) { 542 PetscCall(PetscMalloc1(A->rmap->n, &a->ilen)); 543 } else { 544 PetscCall(PetscMemzero(a->ilen, A->rmap->n * sizeof(PetscInt))); 545 } 546 547 /* allocate the matrix space */ 548 if (A->structure_only) { 549 PetscCall(PetscMalloc1(nztotal, &a->j)); 550 PetscCall(PetscMalloc1(A->rmap->n + 1, &a->i)); 551 } else { 552 PetscCall(PetscMalloc3(nztotal, &a->a, nztotal, &a->j, A->rmap->n + 1, &a->i)); 553 } 554 a->i[0] = 0; 555 if (A->structure_only) { 556 a->singlemalloc = PETSC_FALSE; 557 a->free_a = PETSC_FALSE; 558 } else { 559 a->singlemalloc = PETSC_TRUE; 560 a->free_a = PETSC_TRUE; 561 } 562 a->free_ij = PETSC_TRUE; 563 A->ops->setvalues = MatSetValues_SeqAIJ_SortedFullNoPreallocation; 564 A->preallocated = PETSC_TRUE; 565 PetscFunctionReturn(PETSC_SUCCESS); 566 } 567 568 PetscErrorCode MatSetValues_SeqAIJ_SortedFull(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is) 569 { 570 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 571 PetscInt *rp, k, row; 572 PetscInt *ai = a->i, *ailen = a->ilen; 573 PetscInt *aj = a->j; 574 MatScalar *aa, *ap; 575 576 PetscFunctionBegin; 577 PetscCall(MatSeqAIJGetArray(A, &aa)); 578 for (k = 0; k < m; k++) { /* loop over added rows */ 579 row = im[k]; 580 PetscCheck(n <= a->imax[row], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Preallocation for row %" PetscInt_FMT " does not match number of columns provided", n); 581 rp = aj + ai[row]; 582 ap = aa + ai[row]; 583 if (!A->was_assembled) PetscCall(PetscMemcpy(rp, in, n * sizeof(PetscInt))); 584 if (!A->structure_only) { 585 if (v) { 586 PetscCall(PetscMemcpy(ap, v, n * sizeof(PetscScalar))); 587 v += n; 588 } else { 589 PetscCall(PetscMemzero(ap, n * sizeof(PetscScalar))); 590 } 591 } 592 ailen[row] = n; 593 a->nz += n; 594 } 595 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 596 PetscFunctionReturn(PETSC_SUCCESS); 597 } 598 599 PetscErrorCode MatGetValues_SeqAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[]) 600 { 601 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 602 PetscInt *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 603 PetscInt *ai = a->i, *ailen = a->ilen; 604 const MatScalar *ap, *aa; 605 606 PetscFunctionBegin; 607 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 608 for (k = 0; k < m; k++) { /* loop over rows */ 609 row = im[k]; 610 if (row < 0) { 611 v += n; 612 continue; 613 } /* negative row */ 614 PetscCheck(row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->n - 1); 615 rp = aj + ai[row]; 616 ap = aa + ai[row]; 617 nrow = ailen[row]; 618 for (l = 0; l < n; l++) { /* loop over columns */ 619 if (in[l] < 0) { 620 v++; 621 continue; 622 } /* negative column */ 623 PetscCheck(in[l] < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[l], A->cmap->n - 1); 624 col = in[l]; 625 high = nrow; 626 low = 0; /* assume unsorted */ 627 while (high - low > 5) { 628 t = (low + high) / 2; 629 if (rp[t] > col) high = t; 630 else low = t; 631 } 632 for (i = low; i < high; i++) { 633 if (rp[i] > col) break; 634 if (rp[i] == col) { 635 *v++ = ap[i]; 636 goto finished; 637 } 638 } 639 *v++ = 0.0; 640 finished:; 641 } 642 } 643 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 644 PetscFunctionReturn(PETSC_SUCCESS); 645 } 646 647 PetscErrorCode MatView_SeqAIJ_Binary(Mat mat, PetscViewer viewer) 648 { 649 Mat_SeqAIJ *A = (Mat_SeqAIJ *)mat->data; 650 const PetscScalar *av; 651 PetscInt header[4], M, N, m, nz, i; 652 PetscInt *rowlens; 653 654 PetscFunctionBegin; 655 PetscCall(PetscViewerSetUp(viewer)); 656 657 M = mat->rmap->N; 658 N = mat->cmap->N; 659 m = mat->rmap->n; 660 nz = A->nz; 661 662 /* write matrix header */ 663 header[0] = MAT_FILE_CLASSID; 664 header[1] = M; 665 header[2] = N; 666 header[3] = nz; 667 PetscCall(PetscViewerBinaryWrite(viewer, header, 4, PETSC_INT)); 668 669 /* fill in and store row lengths */ 670 PetscCall(PetscMalloc1(m, &rowlens)); 671 for (i = 0; i < m; i++) rowlens[i] = A->i[i + 1] - A->i[i]; 672 PetscCall(PetscViewerBinaryWrite(viewer, rowlens, m, PETSC_INT)); 673 PetscCall(PetscFree(rowlens)); 674 /* store column indices */ 675 PetscCall(PetscViewerBinaryWrite(viewer, A->j, nz, PETSC_INT)); 676 /* store nonzero values */ 677 PetscCall(MatSeqAIJGetArrayRead(mat, &av)); 678 PetscCall(PetscViewerBinaryWrite(viewer, av, nz, PETSC_SCALAR)); 679 PetscCall(MatSeqAIJRestoreArrayRead(mat, &av)); 680 681 /* write block size option to the viewer's .info file */ 682 PetscCall(MatView_Binary_BlockSizes(mat, viewer)); 683 PetscFunctionReturn(PETSC_SUCCESS); 684 } 685 686 static PetscErrorCode MatView_SeqAIJ_ASCII_structonly(Mat A, PetscViewer viewer) 687 { 688 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 689 PetscInt i, k, m = A->rmap->N; 690 691 PetscFunctionBegin; 692 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 693 for (i = 0; i < m; i++) { 694 PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i)); 695 for (k = a->i[i]; k < a->i[i + 1]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ") ", a->j[k])); 696 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 697 } 698 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 699 PetscFunctionReturn(PETSC_SUCCESS); 700 } 701 702 extern PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat, PetscViewer); 703 704 PetscErrorCode MatView_SeqAIJ_ASCII(Mat A, PetscViewer viewer) 705 { 706 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 707 const PetscScalar *av; 708 PetscInt i, j, m = A->rmap->n; 709 const char *name; 710 PetscViewerFormat format; 711 712 PetscFunctionBegin; 713 if (A->structure_only) { 714 PetscCall(MatView_SeqAIJ_ASCII_structonly(A, viewer)); 715 PetscFunctionReturn(PETSC_SUCCESS); 716 } 717 718 PetscCall(PetscViewerGetFormat(viewer, &format)); 719 if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS); 720 721 /* trigger copy to CPU if needed */ 722 PetscCall(MatSeqAIJGetArrayRead(A, &av)); 723 PetscCall(MatSeqAIJRestoreArrayRead(A, &av)); 724 if (format == PETSC_VIEWER_ASCII_MATLAB) { 725 PetscInt nofinalvalue = 0; 726 if (m && ((a->i[m] == a->i[m - 1]) || (a->j[a->nz - 1] != A->cmap->n - 1))) { 727 /* Need a dummy value to ensure the dimension of the matrix. */ 728 nofinalvalue = 1; 729 } 730 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 731 PetscCall(PetscViewerASCIIPrintf(viewer, "%% Size = %" PetscInt_FMT " %" PetscInt_FMT " \n", m, A->cmap->n)); 732 PetscCall(PetscViewerASCIIPrintf(viewer, "%% Nonzeros = %" PetscInt_FMT " \n", a->nz)); 733 #if defined(PETSC_USE_COMPLEX) 734 PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",4);\n", a->nz + nofinalvalue)); 735 #else 736 PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",3);\n", a->nz + nofinalvalue)); 737 #endif 738 PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = [\n")); 739 740 for (i = 0; i < m; i++) { 741 for (j = a->i[i]; j < a->i[i + 1]; j++) { 742 #if defined(PETSC_USE_COMPLEX) 743 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %18.16e %18.16e\n", i + 1, a->j[j] + 1, (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j]))); 744 #else 745 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n", i + 1, a->j[j] + 1, (double)a->a[j])); 746 #endif 747 } 748 } 749 if (nofinalvalue) { 750 #if defined(PETSC_USE_COMPLEX) 751 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %18.16e %18.16e\n", m, A->cmap->n, 0., 0.)); 752 #else 753 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n", m, A->cmap->n, 0.0)); 754 #endif 755 } 756 PetscCall(PetscObjectGetName((PetscObject)A, &name)); 757 PetscCall(PetscViewerASCIIPrintf(viewer, "];\n %s = spconvert(zzz);\n", name)); 758 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 759 } else if (format == PETSC_VIEWER_ASCII_COMMON) { 760 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 761 for (i = 0; i < m; i++) { 762 PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i)); 763 for (j = a->i[i]; j < a->i[i + 1]; j++) { 764 #if defined(PETSC_USE_COMPLEX) 765 if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) { 766 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j]))); 767 } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) { 768 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)-PetscImaginaryPart(a->a[j]))); 769 } else if (PetscRealPart(a->a[j]) != 0.0) { 770 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j]))); 771 } 772 #else 773 if (a->a[j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j])); 774 #endif 775 } 776 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 777 } 778 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 779 } else if (format == PETSC_VIEWER_ASCII_SYMMODU) { 780 PetscInt nzd = 0, fshift = 1, *sptr; 781 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 782 PetscCall(PetscMalloc1(m + 1, &sptr)); 783 for (i = 0; i < m; i++) { 784 sptr[i] = nzd + 1; 785 for (j = a->i[i]; j < a->i[i + 1]; j++) { 786 if (a->j[j] >= i) { 787 #if defined(PETSC_USE_COMPLEX) 788 if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++; 789 #else 790 if (a->a[j] != 0.0) nzd++; 791 #endif 792 } 793 } 794 } 795 sptr[m] = nzd + 1; 796 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT "\n\n", m, nzd)); 797 for (i = 0; i < m + 1; i += 6) { 798 if (i + 4 < m) { 799 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2], sptr[i + 3], sptr[i + 4], sptr[i + 5])); 800 } else if (i + 3 < m) { 801 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2], sptr[i + 3], sptr[i + 4])); 802 } else if (i + 2 < m) { 803 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2], sptr[i + 3])); 804 } else if (i + 1 < m) { 805 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2])); 806 } else if (i < m) { 807 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1])); 808 } else { 809 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "\n", sptr[i])); 810 } 811 } 812 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 813 PetscCall(PetscFree(sptr)); 814 for (i = 0; i < m; i++) { 815 for (j = a->i[i]; j < a->i[i + 1]; j++) { 816 if (a->j[j] >= i) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " ", a->j[j] + fshift)); 817 } 818 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 819 } 820 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 821 for (i = 0; i < m; i++) { 822 for (j = a->i[i]; j < a->i[i + 1]; j++) { 823 if (a->j[j] >= i) { 824 #if defined(PETSC_USE_COMPLEX) 825 if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " %18.16e %18.16e ", (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j]))); 826 #else 827 if (a->a[j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " %18.16e ", (double)a->a[j])); 828 #endif 829 } 830 } 831 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 832 } 833 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 834 } else if (format == PETSC_VIEWER_ASCII_DENSE) { 835 PetscInt cnt = 0, jcnt; 836 PetscScalar value; 837 #if defined(PETSC_USE_COMPLEX) 838 PetscBool realonly = PETSC_TRUE; 839 840 for (i = 0; i < a->i[m]; i++) { 841 if (PetscImaginaryPart(a->a[i]) != 0.0) { 842 realonly = PETSC_FALSE; 843 break; 844 } 845 } 846 #endif 847 848 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 849 for (i = 0; i < m; i++) { 850 jcnt = 0; 851 for (j = 0; j < A->cmap->n; j++) { 852 if (jcnt < a->i[i + 1] - a->i[i] && j == a->j[cnt]) { 853 value = a->a[cnt++]; 854 jcnt++; 855 } else { 856 value = 0.0; 857 } 858 #if defined(PETSC_USE_COMPLEX) 859 if (realonly) { 860 PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)PetscRealPart(value))); 861 } else { 862 PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e+%7.5e i ", (double)PetscRealPart(value), (double)PetscImaginaryPart(value))); 863 } 864 #else 865 PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)value)); 866 #endif 867 } 868 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 869 } 870 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 871 } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) { 872 PetscInt fshift = 1; 873 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 874 #if defined(PETSC_USE_COMPLEX) 875 PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate complex general\n")); 876 #else 877 PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate real general\n")); 878 #endif 879 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, A->cmap->n, a->nz)); 880 for (i = 0; i < m; i++) { 881 for (j = a->i[i]; j < a->i[i + 1]; j++) { 882 #if defined(PETSC_USE_COMPLEX) 883 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g %g\n", i + fshift, a->j[j] + fshift, (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j]))); 884 #else 885 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g\n", i + fshift, a->j[j] + fshift, (double)a->a[j])); 886 #endif 887 } 888 } 889 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 890 } else { 891 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 892 if (A->factortype) { 893 for (i = 0; i < m; i++) { 894 PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i)); 895 /* L part */ 896 for (j = a->i[i]; j < a->i[i + 1]; j++) { 897 #if defined(PETSC_USE_COMPLEX) 898 if (PetscImaginaryPart(a->a[j]) > 0.0) { 899 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j]))); 900 } else if (PetscImaginaryPart(a->a[j]) < 0.0) { 901 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)(-PetscImaginaryPart(a->a[j])))); 902 } else { 903 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j]))); 904 } 905 #else 906 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j])); 907 #endif 908 } 909 /* diagonal */ 910 j = a->diag[i]; 911 #if defined(PETSC_USE_COMPLEX) 912 if (PetscImaginaryPart(a->a[j]) > 0.0) { 913 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(1.0 / a->a[j]), (double)PetscImaginaryPart(1.0 / a->a[j]))); 914 } else if (PetscImaginaryPart(a->a[j]) < 0.0) { 915 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(1.0 / a->a[j]), (double)(-PetscImaginaryPart(1.0 / a->a[j])))); 916 } else { 917 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(1.0 / a->a[j]))); 918 } 919 #else 920 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)(1.0 / a->a[j]))); 921 #endif 922 923 /* U part */ 924 for (j = a->diag[i + 1] + 1; j < a->diag[i]; j++) { 925 #if defined(PETSC_USE_COMPLEX) 926 if (PetscImaginaryPart(a->a[j]) > 0.0) { 927 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j]))); 928 } else if (PetscImaginaryPart(a->a[j]) < 0.0) { 929 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)(-PetscImaginaryPart(a->a[j])))); 930 } else { 931 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j]))); 932 } 933 #else 934 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j])); 935 #endif 936 } 937 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 938 } 939 } else { 940 for (i = 0; i < m; i++) { 941 PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i)); 942 for (j = a->i[i]; j < a->i[i + 1]; j++) { 943 #if defined(PETSC_USE_COMPLEX) 944 if (PetscImaginaryPart(a->a[j]) > 0.0) { 945 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j]))); 946 } else if (PetscImaginaryPart(a->a[j]) < 0.0) { 947 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)-PetscImaginaryPart(a->a[j]))); 948 } else { 949 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j]))); 950 } 951 #else 952 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j])); 953 #endif 954 } 955 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 956 } 957 } 958 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 959 } 960 PetscCall(PetscViewerFlush(viewer)); 961 PetscFunctionReturn(PETSC_SUCCESS); 962 } 963 964 #include <petscdraw.h> 965 PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw, void *Aa) 966 { 967 Mat A = (Mat)Aa; 968 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 969 PetscInt i, j, m = A->rmap->n; 970 int color; 971 PetscReal xl, yl, xr, yr, x_l, x_r, y_l, y_r; 972 PetscViewer viewer; 973 PetscViewerFormat format; 974 const PetscScalar *aa; 975 976 PetscFunctionBegin; 977 PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer)); 978 PetscCall(PetscViewerGetFormat(viewer, &format)); 979 PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 980 981 /* loop over matrix elements drawing boxes */ 982 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 983 if (format != PETSC_VIEWER_DRAW_CONTOUR) { 984 PetscDrawCollectiveBegin(draw); 985 /* Blue for negative, Cyan for zero and Red for positive */ 986 color = PETSC_DRAW_BLUE; 987 for (i = 0; i < m; i++) { 988 y_l = m - i - 1.0; 989 y_r = y_l + 1.0; 990 for (j = a->i[i]; j < a->i[i + 1]; j++) { 991 x_l = a->j[j]; 992 x_r = x_l + 1.0; 993 if (PetscRealPart(aa[j]) >= 0.) continue; 994 PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color)); 995 } 996 } 997 color = PETSC_DRAW_CYAN; 998 for (i = 0; i < m; i++) { 999 y_l = m - i - 1.0; 1000 y_r = y_l + 1.0; 1001 for (j = a->i[i]; j < a->i[i + 1]; j++) { 1002 x_l = a->j[j]; 1003 x_r = x_l + 1.0; 1004 if (aa[j] != 0.) continue; 1005 PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color)); 1006 } 1007 } 1008 color = PETSC_DRAW_RED; 1009 for (i = 0; i < m; i++) { 1010 y_l = m - i - 1.0; 1011 y_r = y_l + 1.0; 1012 for (j = a->i[i]; j < a->i[i + 1]; j++) { 1013 x_l = a->j[j]; 1014 x_r = x_l + 1.0; 1015 if (PetscRealPart(aa[j]) <= 0.) continue; 1016 PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color)); 1017 } 1018 } 1019 PetscDrawCollectiveEnd(draw); 1020 } else { 1021 /* use contour shading to indicate magnitude of values */ 1022 /* first determine max of all nonzero values */ 1023 PetscReal minv = 0.0, maxv = 0.0; 1024 PetscInt nz = a->nz, count = 0; 1025 PetscDraw popup; 1026 1027 for (i = 0; i < nz; i++) { 1028 if (PetscAbsScalar(aa[i]) > maxv) maxv = PetscAbsScalar(aa[i]); 1029 } 1030 if (minv >= maxv) maxv = minv + PETSC_SMALL; 1031 PetscCall(PetscDrawGetPopup(draw, &popup)); 1032 PetscCall(PetscDrawScalePopup(popup, minv, maxv)); 1033 1034 PetscDrawCollectiveBegin(draw); 1035 for (i = 0; i < m; i++) { 1036 y_l = m - i - 1.0; 1037 y_r = y_l + 1.0; 1038 for (j = a->i[i]; j < a->i[i + 1]; j++) { 1039 x_l = a->j[j]; 1040 x_r = x_l + 1.0; 1041 color = PetscDrawRealToColor(PetscAbsScalar(aa[count]), minv, maxv); 1042 PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color)); 1043 count++; 1044 } 1045 } 1046 PetscDrawCollectiveEnd(draw); 1047 } 1048 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 1049 PetscFunctionReturn(PETSC_SUCCESS); 1050 } 1051 1052 #include <petscdraw.h> 1053 PetscErrorCode MatView_SeqAIJ_Draw(Mat A, PetscViewer viewer) 1054 { 1055 PetscDraw draw; 1056 PetscReal xr, yr, xl, yl, h, w; 1057 PetscBool isnull; 1058 1059 PetscFunctionBegin; 1060 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 1061 PetscCall(PetscDrawIsNull(draw, &isnull)); 1062 if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 1063 1064 xr = A->cmap->n; 1065 yr = A->rmap->n; 1066 h = yr / 10.0; 1067 w = xr / 10.0; 1068 xr += w; 1069 yr += h; 1070 xl = -w; 1071 yl = -h; 1072 PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr)); 1073 PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer)); 1074 PetscCall(PetscDrawZoom(draw, MatView_SeqAIJ_Draw_Zoom, A)); 1075 PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL)); 1076 PetscCall(PetscDrawSave(draw)); 1077 PetscFunctionReturn(PETSC_SUCCESS); 1078 } 1079 1080 PetscErrorCode MatView_SeqAIJ(Mat A, PetscViewer viewer) 1081 { 1082 PetscBool iascii, isbinary, isdraw; 1083 1084 PetscFunctionBegin; 1085 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1086 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 1087 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1088 if (iascii) PetscCall(MatView_SeqAIJ_ASCII(A, viewer)); 1089 else if (isbinary) PetscCall(MatView_SeqAIJ_Binary(A, viewer)); 1090 else if (isdraw) PetscCall(MatView_SeqAIJ_Draw(A, viewer)); 1091 PetscCall(MatView_SeqAIJ_Inode(A, viewer)); 1092 PetscFunctionReturn(PETSC_SUCCESS); 1093 } 1094 1095 PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A, MatAssemblyType mode) 1096 { 1097 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1098 PetscInt fshift = 0, i, *ai = a->i, *aj = a->j, *imax = a->imax; 1099 PetscInt m = A->rmap->n, *ip, N, *ailen = a->ilen, rmax = 0; 1100 MatScalar *aa = a->a, *ap; 1101 PetscReal ratio = 0.6; 1102 1103 PetscFunctionBegin; 1104 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS); 1105 PetscCall(MatSeqAIJInvalidateDiagonal(A)); 1106 if (A->was_assembled && A->ass_nonzerostate == A->nonzerostate) { 1107 /* we need to respect users asking to use or not the inodes routine in between matrix assemblies */ 1108 PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode)); 1109 PetscFunctionReturn(PETSC_SUCCESS); 1110 } 1111 1112 if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 1113 for (i = 1; i < m; i++) { 1114 /* move each row back by the amount of empty slots (fshift) before it*/ 1115 fshift += imax[i - 1] - ailen[i - 1]; 1116 rmax = PetscMax(rmax, ailen[i]); 1117 if (fshift) { 1118 ip = aj + ai[i]; 1119 ap = aa + ai[i]; 1120 N = ailen[i]; 1121 PetscCall(PetscArraymove(ip - fshift, ip, N)); 1122 if (!A->structure_only) PetscCall(PetscArraymove(ap - fshift, ap, N)); 1123 } 1124 ai[i] = ai[i - 1] + ailen[i - 1]; 1125 } 1126 if (m) { 1127 fshift += imax[m - 1] - ailen[m - 1]; 1128 ai[m] = ai[m - 1] + ailen[m - 1]; 1129 } 1130 /* reset ilen and imax for each row */ 1131 a->nonzerorowcnt = 0; 1132 if (A->structure_only) { 1133 PetscCall(PetscFree(a->imax)); 1134 PetscCall(PetscFree(a->ilen)); 1135 } else { /* !A->structure_only */ 1136 for (i = 0; i < m; i++) { 1137 ailen[i] = imax[i] = ai[i + 1] - ai[i]; 1138 a->nonzerorowcnt += ((ai[i + 1] - ai[i]) > 0); 1139 } 1140 } 1141 a->nz = ai[m]; 1142 PetscCheck(!fshift || a->nounused != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unused space detected in matrix: %" PetscInt_FMT " X %" PetscInt_FMT ", %" PetscInt_FMT " unneeded", m, A->cmap->n, fshift); 1143 1144 PetscCall(MatMarkDiagonal_SeqAIJ(A)); 1145 PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; storage space: %" PetscInt_FMT " unneeded,%" PetscInt_FMT " used\n", m, A->cmap->n, fshift, a->nz)); 1146 PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs)); 1147 PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", rmax)); 1148 1149 A->info.mallocs += a->reallocs; 1150 a->reallocs = 0; 1151 A->info.nz_unneeded = (PetscReal)fshift; 1152 a->rmax = rmax; 1153 1154 if (!A->structure_only) PetscCall(MatCheckCompressedRow(A, a->nonzerorowcnt, &a->compressedrow, a->i, m, ratio)); 1155 PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode)); 1156 PetscFunctionReturn(PETSC_SUCCESS); 1157 } 1158 1159 PetscErrorCode MatRealPart_SeqAIJ(Mat A) 1160 { 1161 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1162 PetscInt i, nz = a->nz; 1163 MatScalar *aa; 1164 1165 PetscFunctionBegin; 1166 PetscCall(MatSeqAIJGetArray(A, &aa)); 1167 for (i = 0; i < nz; i++) aa[i] = PetscRealPart(aa[i]); 1168 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 1169 PetscCall(MatSeqAIJInvalidateDiagonal(A)); 1170 PetscFunctionReturn(PETSC_SUCCESS); 1171 } 1172 1173 PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A) 1174 { 1175 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1176 PetscInt i, nz = a->nz; 1177 MatScalar *aa; 1178 1179 PetscFunctionBegin; 1180 PetscCall(MatSeqAIJGetArray(A, &aa)); 1181 for (i = 0; i < nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1182 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 1183 PetscCall(MatSeqAIJInvalidateDiagonal(A)); 1184 PetscFunctionReturn(PETSC_SUCCESS); 1185 } 1186 1187 PetscErrorCode MatZeroEntries_SeqAIJ(Mat A) 1188 { 1189 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1190 MatScalar *aa; 1191 1192 PetscFunctionBegin; 1193 PetscCall(MatSeqAIJGetArrayWrite(A, &aa)); 1194 PetscCall(PetscArrayzero(aa, a->i[A->rmap->n])); 1195 PetscCall(MatSeqAIJRestoreArrayWrite(A, &aa)); 1196 PetscCall(MatSeqAIJInvalidateDiagonal(A)); 1197 PetscFunctionReturn(PETSC_SUCCESS); 1198 } 1199 1200 PetscErrorCode MatDestroy_SeqAIJ(Mat A) 1201 { 1202 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1203 1204 PetscFunctionBegin; 1205 if (A->hash_active) { 1206 A->ops[0] = a->cops; 1207 PetscCall(PetscHMapIJVDestroy(&a->ht)); 1208 PetscCall(PetscFree(a->dnz)); 1209 A->hash_active = PETSC_FALSE; 1210 } 1211 1212 PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->n, A->cmap->n, a->nz)); 1213 PetscCall(MatSeqXAIJFreeAIJ(A, &a->a, &a->j, &a->i)); 1214 PetscCall(ISDestroy(&a->row)); 1215 PetscCall(ISDestroy(&a->col)); 1216 PetscCall(PetscFree(a->diag)); 1217 PetscCall(PetscFree(a->ibdiag)); 1218 PetscCall(PetscFree(a->imax)); 1219 PetscCall(PetscFree(a->ilen)); 1220 PetscCall(PetscFree(a->ipre)); 1221 PetscCall(PetscFree3(a->idiag, a->mdiag, a->ssor_work)); 1222 PetscCall(PetscFree(a->solve_work)); 1223 PetscCall(ISDestroy(&a->icol)); 1224 PetscCall(PetscFree(a->saved_values)); 1225 PetscCall(PetscFree2(a->compressedrow.i, a->compressedrow.rindex)); 1226 PetscCall(MatDestroy_SeqAIJ_Inode(A)); 1227 PetscCall(PetscFree(A->data)); 1228 1229 /* MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted may allocate this. 1230 That function is so heavily used (sometimes in an hidden way through multnumeric function pointers) 1231 that is hard to properly add this data to the MatProduct data. We free it here to avoid 1232 users reusing the matrix object with different data to incur in obscure segmentation faults 1233 due to different matrix sizes */ 1234 PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc__ab_dense", NULL)); 1235 1236 PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL)); 1237 PetscCall(PetscObjectComposeFunction((PetscObject)A, "PetscMatlabEnginePut_C", NULL)); 1238 PetscCall(PetscObjectComposeFunction((PetscObject)A, "PetscMatlabEngineGet_C", NULL)); 1239 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetColumnIndices_C", NULL)); 1240 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL)); 1241 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL)); 1242 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqsbaij_C", NULL)); 1243 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqbaij_C", NULL)); 1244 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijperm_C", NULL)); 1245 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijsell_C", NULL)); 1246 #if defined(PETSC_HAVE_MKL_SPARSE) 1247 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijmkl_C", NULL)); 1248 #endif 1249 #if defined(PETSC_HAVE_CUDA) 1250 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijcusparse_C", NULL)); 1251 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijcusparse_seqaij_C", NULL)); 1252 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaijcusparse_C", NULL)); 1253 #endif 1254 #if defined(PETSC_HAVE_HIP) 1255 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijhipsparse_C", NULL)); 1256 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijhipsparse_seqaij_C", NULL)); 1257 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaijhipsparse_C", NULL)); 1258 #endif 1259 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 1260 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijkokkos_C", NULL)); 1261 #endif 1262 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijcrl_C", NULL)); 1263 #if defined(PETSC_HAVE_ELEMENTAL) 1264 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_elemental_C", NULL)); 1265 #endif 1266 #if defined(PETSC_HAVE_SCALAPACK) 1267 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_scalapack_C", NULL)); 1268 #endif 1269 #if defined(PETSC_HAVE_HYPRE) 1270 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_hypre_C", NULL)); 1271 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_transpose_seqaij_seqaij_C", NULL)); 1272 #endif 1273 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqdense_C", NULL)); 1274 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqsell_C", NULL)); 1275 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_is_C", NULL)); 1276 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatIsTranspose_C", NULL)); 1277 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatIsHermitianTranspose_C", NULL)); 1278 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", NULL)); 1279 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatResetPreallocation_C", NULL)); 1280 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetPreallocationCSR_C", NULL)); 1281 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatReorderForNonzeroDiagonal_C", NULL)); 1282 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_is_seqaij_C", NULL)); 1283 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqdense_seqaij_C", NULL)); 1284 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaij_C", NULL)); 1285 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJKron_C", NULL)); 1286 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetPreallocationCOO_C", NULL)); 1287 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetValuesCOO_C", NULL)); 1288 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL)); 1289 /* these calls do not belong here: the subclasses Duplicate/Destroy are wrong */ 1290 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijsell_seqaij_C", NULL)); 1291 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijperm_seqaij_C", NULL)); 1292 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijviennacl_C", NULL)); 1293 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijviennacl_seqdense_C", NULL)); 1294 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijviennacl_seqaij_C", NULL)); 1295 PetscFunctionReturn(PETSC_SUCCESS); 1296 } 1297 1298 PetscErrorCode MatSetOption_SeqAIJ(Mat A, MatOption op, PetscBool flg) 1299 { 1300 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1301 1302 PetscFunctionBegin; 1303 switch (op) { 1304 case MAT_ROW_ORIENTED: 1305 a->roworiented = flg; 1306 break; 1307 case MAT_KEEP_NONZERO_PATTERN: 1308 a->keepnonzeropattern = flg; 1309 break; 1310 case MAT_NEW_NONZERO_LOCATIONS: 1311 a->nonew = (flg ? 0 : 1); 1312 break; 1313 case MAT_NEW_NONZERO_LOCATION_ERR: 1314 a->nonew = (flg ? -1 : 0); 1315 break; 1316 case MAT_NEW_NONZERO_ALLOCATION_ERR: 1317 a->nonew = (flg ? -2 : 0); 1318 break; 1319 case MAT_UNUSED_NONZERO_LOCATION_ERR: 1320 a->nounused = (flg ? -1 : 0); 1321 break; 1322 case MAT_IGNORE_ZERO_ENTRIES: 1323 a->ignorezeroentries = flg; 1324 break; 1325 case MAT_SPD: 1326 case MAT_SYMMETRIC: 1327 case MAT_STRUCTURALLY_SYMMETRIC: 1328 case MAT_HERMITIAN: 1329 case MAT_SYMMETRY_ETERNAL: 1330 case MAT_STRUCTURE_ONLY: 1331 case MAT_STRUCTURAL_SYMMETRY_ETERNAL: 1332 case MAT_SPD_ETERNAL: 1333 /* if the diagonal matrix is square it inherits some of the properties above */ 1334 break; 1335 case MAT_FORCE_DIAGONAL_ENTRIES: 1336 case MAT_IGNORE_OFF_PROC_ENTRIES: 1337 case MAT_USE_HASH_TABLE: 1338 PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op])); 1339 break; 1340 case MAT_USE_INODES: 1341 PetscCall(MatSetOption_SeqAIJ_Inode(A, MAT_USE_INODES, flg)); 1342 break; 1343 case MAT_SUBMAT_SINGLEIS: 1344 A->submat_singleis = flg; 1345 break; 1346 case MAT_SORTED_FULL: 1347 if (flg) A->ops->setvalues = MatSetValues_SeqAIJ_SortedFull; 1348 else A->ops->setvalues = MatSetValues_SeqAIJ; 1349 break; 1350 case MAT_FORM_EXPLICIT_TRANSPOSE: 1351 A->form_explicit_transpose = flg; 1352 break; 1353 default: 1354 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op); 1355 } 1356 PetscFunctionReturn(PETSC_SUCCESS); 1357 } 1358 1359 PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A, Vec v) 1360 { 1361 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1362 PetscInt i, j, n, *ai = a->i, *aj = a->j; 1363 PetscScalar *x; 1364 const PetscScalar *aa; 1365 1366 PetscFunctionBegin; 1367 PetscCall(VecGetLocalSize(v, &n)); 1368 PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector"); 1369 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 1370 if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) { 1371 PetscInt *diag = a->diag; 1372 PetscCall(VecGetArrayWrite(v, &x)); 1373 for (i = 0; i < n; i++) x[i] = 1.0 / aa[diag[i]]; 1374 PetscCall(VecRestoreArrayWrite(v, &x)); 1375 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 1376 PetscFunctionReturn(PETSC_SUCCESS); 1377 } 1378 1379 PetscCall(VecGetArrayWrite(v, &x)); 1380 for (i = 0; i < n; i++) { 1381 x[i] = 0.0; 1382 for (j = ai[i]; j < ai[i + 1]; j++) { 1383 if (aj[j] == i) { 1384 x[i] = aa[j]; 1385 break; 1386 } 1387 } 1388 } 1389 PetscCall(VecRestoreArrayWrite(v, &x)); 1390 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 1391 PetscFunctionReturn(PETSC_SUCCESS); 1392 } 1393 1394 #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h> 1395 PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A, Vec xx, Vec zz, Vec yy) 1396 { 1397 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1398 const MatScalar *aa; 1399 PetscScalar *y; 1400 const PetscScalar *x; 1401 PetscInt m = A->rmap->n; 1402 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) 1403 const MatScalar *v; 1404 PetscScalar alpha; 1405 PetscInt n, i, j; 1406 const PetscInt *idx, *ii, *ridx = NULL; 1407 Mat_CompressedRow cprow = a->compressedrow; 1408 PetscBool usecprow = cprow.use; 1409 #endif 1410 1411 PetscFunctionBegin; 1412 if (zz != yy) PetscCall(VecCopy(zz, yy)); 1413 PetscCall(VecGetArrayRead(xx, &x)); 1414 PetscCall(VecGetArray(yy, &y)); 1415 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 1416 1417 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) 1418 fortranmulttransposeaddaij_(&m, x, a->i, a->j, aa, y); 1419 #else 1420 if (usecprow) { 1421 m = cprow.nrows; 1422 ii = cprow.i; 1423 ridx = cprow.rindex; 1424 } else { 1425 ii = a->i; 1426 } 1427 for (i = 0; i < m; i++) { 1428 idx = a->j + ii[i]; 1429 v = aa + ii[i]; 1430 n = ii[i + 1] - ii[i]; 1431 if (usecprow) { 1432 alpha = x[ridx[i]]; 1433 } else { 1434 alpha = x[i]; 1435 } 1436 for (j = 0; j < n; j++) y[idx[j]] += alpha * v[j]; 1437 } 1438 #endif 1439 PetscCall(PetscLogFlops(2.0 * a->nz)); 1440 PetscCall(VecRestoreArrayRead(xx, &x)); 1441 PetscCall(VecRestoreArray(yy, &y)); 1442 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 1443 PetscFunctionReturn(PETSC_SUCCESS); 1444 } 1445 1446 PetscErrorCode MatMultTranspose_SeqAIJ(Mat A, Vec xx, Vec yy) 1447 { 1448 PetscFunctionBegin; 1449 PetscCall(VecSet(yy, 0.0)); 1450 PetscCall(MatMultTransposeAdd_SeqAIJ(A, xx, yy, yy)); 1451 PetscFunctionReturn(PETSC_SUCCESS); 1452 } 1453 1454 #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h> 1455 1456 PetscErrorCode MatMult_SeqAIJ(Mat A, Vec xx, Vec yy) 1457 { 1458 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1459 PetscScalar *y; 1460 const PetscScalar *x; 1461 const MatScalar *aa, *a_a; 1462 PetscInt m = A->rmap->n; 1463 const PetscInt *aj, *ii, *ridx = NULL; 1464 PetscInt n, i; 1465 PetscScalar sum; 1466 PetscBool usecprow = a->compressedrow.use; 1467 1468 #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 1469 #pragma disjoint(*x, *y, *aa) 1470 #endif 1471 1472 PetscFunctionBegin; 1473 if (a->inode.use && a->inode.checked) { 1474 PetscCall(MatMult_SeqAIJ_Inode(A, xx, yy)); 1475 PetscFunctionReturn(PETSC_SUCCESS); 1476 } 1477 PetscCall(MatSeqAIJGetArrayRead(A, &a_a)); 1478 PetscCall(VecGetArrayRead(xx, &x)); 1479 PetscCall(VecGetArray(yy, &y)); 1480 ii = a->i; 1481 if (usecprow) { /* use compressed row format */ 1482 PetscCall(PetscArrayzero(y, m)); 1483 m = a->compressedrow.nrows; 1484 ii = a->compressedrow.i; 1485 ridx = a->compressedrow.rindex; 1486 for (i = 0; i < m; i++) { 1487 n = ii[i + 1] - ii[i]; 1488 aj = a->j + ii[i]; 1489 aa = a_a + ii[i]; 1490 sum = 0.0; 1491 PetscSparseDensePlusDot(sum, x, aa, aj, n); 1492 /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */ 1493 y[*ridx++] = sum; 1494 } 1495 } else { /* do not use compressed row format */ 1496 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ) 1497 aj = a->j; 1498 aa = a_a; 1499 fortranmultaij_(&m, x, ii, aj, aa, y); 1500 #else 1501 for (i = 0; i < m; i++) { 1502 n = ii[i + 1] - ii[i]; 1503 aj = a->j + ii[i]; 1504 aa = a_a + ii[i]; 1505 sum = 0.0; 1506 PetscSparseDensePlusDot(sum, x, aa, aj, n); 1507 y[i] = sum; 1508 } 1509 #endif 1510 } 1511 PetscCall(PetscLogFlops(2.0 * a->nz - a->nonzerorowcnt)); 1512 PetscCall(VecRestoreArrayRead(xx, &x)); 1513 PetscCall(VecRestoreArray(yy, &y)); 1514 PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a)); 1515 PetscFunctionReturn(PETSC_SUCCESS); 1516 } 1517 1518 PetscErrorCode MatMultMax_SeqAIJ(Mat A, Vec xx, Vec yy) 1519 { 1520 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1521 PetscScalar *y; 1522 const PetscScalar *x; 1523 const MatScalar *aa, *a_a; 1524 PetscInt m = A->rmap->n; 1525 const PetscInt *aj, *ii, *ridx = NULL; 1526 PetscInt n, i, nonzerorow = 0; 1527 PetscScalar sum; 1528 PetscBool usecprow = a->compressedrow.use; 1529 1530 #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 1531 #pragma disjoint(*x, *y, *aa) 1532 #endif 1533 1534 PetscFunctionBegin; 1535 PetscCall(MatSeqAIJGetArrayRead(A, &a_a)); 1536 PetscCall(VecGetArrayRead(xx, &x)); 1537 PetscCall(VecGetArray(yy, &y)); 1538 if (usecprow) { /* use compressed row format */ 1539 m = a->compressedrow.nrows; 1540 ii = a->compressedrow.i; 1541 ridx = a->compressedrow.rindex; 1542 for (i = 0; i < m; i++) { 1543 n = ii[i + 1] - ii[i]; 1544 aj = a->j + ii[i]; 1545 aa = a_a + ii[i]; 1546 sum = 0.0; 1547 nonzerorow += (n > 0); 1548 PetscSparseDenseMaxDot(sum, x, aa, aj, n); 1549 /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */ 1550 y[*ridx++] = sum; 1551 } 1552 } else { /* do not use compressed row format */ 1553 ii = a->i; 1554 for (i = 0; i < m; i++) { 1555 n = ii[i + 1] - ii[i]; 1556 aj = a->j + ii[i]; 1557 aa = a_a + ii[i]; 1558 sum = 0.0; 1559 nonzerorow += (n > 0); 1560 PetscSparseDenseMaxDot(sum, x, aa, aj, n); 1561 y[i] = sum; 1562 } 1563 } 1564 PetscCall(PetscLogFlops(2.0 * a->nz - nonzerorow)); 1565 PetscCall(VecRestoreArrayRead(xx, &x)); 1566 PetscCall(VecRestoreArray(yy, &y)); 1567 PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a)); 1568 PetscFunctionReturn(PETSC_SUCCESS); 1569 } 1570 1571 PetscErrorCode MatMultAddMax_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz) 1572 { 1573 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1574 PetscScalar *y, *z; 1575 const PetscScalar *x; 1576 const MatScalar *aa, *a_a; 1577 PetscInt m = A->rmap->n, *aj, *ii; 1578 PetscInt n, i, *ridx = NULL; 1579 PetscScalar sum; 1580 PetscBool usecprow = a->compressedrow.use; 1581 1582 PetscFunctionBegin; 1583 PetscCall(MatSeqAIJGetArrayRead(A, &a_a)); 1584 PetscCall(VecGetArrayRead(xx, &x)); 1585 PetscCall(VecGetArrayPair(yy, zz, &y, &z)); 1586 if (usecprow) { /* use compressed row format */ 1587 if (zz != yy) PetscCall(PetscArraycpy(z, y, m)); 1588 m = a->compressedrow.nrows; 1589 ii = a->compressedrow.i; 1590 ridx = a->compressedrow.rindex; 1591 for (i = 0; i < m; i++) { 1592 n = ii[i + 1] - ii[i]; 1593 aj = a->j + ii[i]; 1594 aa = a_a + ii[i]; 1595 sum = y[*ridx]; 1596 PetscSparseDenseMaxDot(sum, x, aa, aj, n); 1597 z[*ridx++] = sum; 1598 } 1599 } else { /* do not use compressed row format */ 1600 ii = a->i; 1601 for (i = 0; i < m; i++) { 1602 n = ii[i + 1] - ii[i]; 1603 aj = a->j + ii[i]; 1604 aa = a_a + ii[i]; 1605 sum = y[i]; 1606 PetscSparseDenseMaxDot(sum, x, aa, aj, n); 1607 z[i] = sum; 1608 } 1609 } 1610 PetscCall(PetscLogFlops(2.0 * a->nz)); 1611 PetscCall(VecRestoreArrayRead(xx, &x)); 1612 PetscCall(VecRestoreArrayPair(yy, zz, &y, &z)); 1613 PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a)); 1614 PetscFunctionReturn(PETSC_SUCCESS); 1615 } 1616 1617 #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h> 1618 PetscErrorCode MatMultAdd_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz) 1619 { 1620 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1621 PetscScalar *y, *z; 1622 const PetscScalar *x; 1623 const MatScalar *aa, *a_a; 1624 const PetscInt *aj, *ii, *ridx = NULL; 1625 PetscInt m = A->rmap->n, n, i; 1626 PetscScalar sum; 1627 PetscBool usecprow = a->compressedrow.use; 1628 1629 PetscFunctionBegin; 1630 if (a->inode.use && a->inode.checked) { 1631 PetscCall(MatMultAdd_SeqAIJ_Inode(A, xx, yy, zz)); 1632 PetscFunctionReturn(PETSC_SUCCESS); 1633 } 1634 PetscCall(MatSeqAIJGetArrayRead(A, &a_a)); 1635 PetscCall(VecGetArrayRead(xx, &x)); 1636 PetscCall(VecGetArrayPair(yy, zz, &y, &z)); 1637 if (usecprow) { /* use compressed row format */ 1638 if (zz != yy) PetscCall(PetscArraycpy(z, y, m)); 1639 m = a->compressedrow.nrows; 1640 ii = a->compressedrow.i; 1641 ridx = a->compressedrow.rindex; 1642 for (i = 0; i < m; i++) { 1643 n = ii[i + 1] - ii[i]; 1644 aj = a->j + ii[i]; 1645 aa = a_a + ii[i]; 1646 sum = y[*ridx]; 1647 PetscSparseDensePlusDot(sum, x, aa, aj, n); 1648 z[*ridx++] = sum; 1649 } 1650 } else { /* do not use compressed row format */ 1651 ii = a->i; 1652 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) 1653 aj = a->j; 1654 aa = a_a; 1655 fortranmultaddaij_(&m, x, ii, aj, aa, y, z); 1656 #else 1657 for (i = 0; i < m; i++) { 1658 n = ii[i + 1] - ii[i]; 1659 aj = a->j + ii[i]; 1660 aa = a_a + ii[i]; 1661 sum = y[i]; 1662 PetscSparseDensePlusDot(sum, x, aa, aj, n); 1663 z[i] = sum; 1664 } 1665 #endif 1666 } 1667 PetscCall(PetscLogFlops(2.0 * a->nz)); 1668 PetscCall(VecRestoreArrayRead(xx, &x)); 1669 PetscCall(VecRestoreArrayPair(yy, zz, &y, &z)); 1670 PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a)); 1671 PetscFunctionReturn(PETSC_SUCCESS); 1672 } 1673 1674 /* 1675 Adds diagonal pointers to sparse matrix structure. 1676 */ 1677 PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A) 1678 { 1679 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1680 PetscInt i, j, m = A->rmap->n; 1681 PetscBool alreadySet = PETSC_TRUE; 1682 1683 PetscFunctionBegin; 1684 if (!a->diag) { 1685 PetscCall(PetscMalloc1(m, &a->diag)); 1686 alreadySet = PETSC_FALSE; 1687 } 1688 for (i = 0; i < A->rmap->n; i++) { 1689 /* If A's diagonal is already correctly set, this fast track enables cheap and repeated MatMarkDiagonal_SeqAIJ() calls */ 1690 if (alreadySet) { 1691 PetscInt pos = a->diag[i]; 1692 if (pos >= a->i[i] && pos < a->i[i + 1] && a->j[pos] == i) continue; 1693 } 1694 1695 a->diag[i] = a->i[i + 1]; 1696 for (j = a->i[i]; j < a->i[i + 1]; j++) { 1697 if (a->j[j] == i) { 1698 a->diag[i] = j; 1699 break; 1700 } 1701 } 1702 } 1703 PetscFunctionReturn(PETSC_SUCCESS); 1704 } 1705 1706 PetscErrorCode MatShift_SeqAIJ(Mat A, PetscScalar v) 1707 { 1708 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1709 const PetscInt *diag = (const PetscInt *)a->diag; 1710 const PetscInt *ii = (const PetscInt *)a->i; 1711 PetscInt i, *mdiag = NULL; 1712 PetscInt cnt = 0; /* how many diagonals are missing */ 1713 1714 PetscFunctionBegin; 1715 if (!A->preallocated || !a->nz) { 1716 PetscCall(MatSeqAIJSetPreallocation(A, 1, NULL)); 1717 PetscCall(MatShift_Basic(A, v)); 1718 PetscFunctionReturn(PETSC_SUCCESS); 1719 } 1720 1721 if (a->diagonaldense) { 1722 cnt = 0; 1723 } else { 1724 PetscCall(PetscCalloc1(A->rmap->n, &mdiag)); 1725 for (i = 0; i < A->rmap->n; i++) { 1726 if (i < A->cmap->n && diag[i] >= ii[i + 1]) { /* 'out of range' rows never have diagonals */ 1727 cnt++; 1728 mdiag[i] = 1; 1729 } 1730 } 1731 } 1732 if (!cnt) { 1733 PetscCall(MatShift_Basic(A, v)); 1734 } else { 1735 PetscScalar *olda = a->a; /* preserve pointers to current matrix nonzeros structure and values */ 1736 PetscInt *oldj = a->j, *oldi = a->i; 1737 PetscBool singlemalloc = a->singlemalloc, free_a = a->free_a, free_ij = a->free_ij; 1738 1739 a->a = NULL; 1740 a->j = NULL; 1741 a->i = NULL; 1742 /* increase the values in imax for each row where a diagonal is being inserted then reallocate the matrix data structures */ 1743 for (i = 0; i < PetscMin(A->rmap->n, A->cmap->n); i++) a->imax[i] += mdiag[i]; 1744 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(A, 0, a->imax)); 1745 1746 /* copy old values into new matrix data structure */ 1747 for (i = 0; i < A->rmap->n; i++) { 1748 PetscCall(MatSetValues(A, 1, &i, a->imax[i] - mdiag[i], &oldj[oldi[i]], &olda[oldi[i]], ADD_VALUES)); 1749 if (i < A->cmap->n) PetscCall(MatSetValue(A, i, i, v, ADD_VALUES)); 1750 } 1751 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1752 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1753 if (singlemalloc) { 1754 PetscCall(PetscFree3(olda, oldj, oldi)); 1755 } else { 1756 if (free_a) PetscCall(PetscFree(olda)); 1757 if (free_ij) PetscCall(PetscFree(oldj)); 1758 if (free_ij) PetscCall(PetscFree(oldi)); 1759 } 1760 } 1761 PetscCall(PetscFree(mdiag)); 1762 a->diagonaldense = PETSC_TRUE; 1763 PetscFunctionReturn(PETSC_SUCCESS); 1764 } 1765 1766 /* 1767 Checks for missing diagonals 1768 */ 1769 PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A, PetscBool *missing, PetscInt *d) 1770 { 1771 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1772 PetscInt *diag, *ii = a->i, i; 1773 1774 PetscFunctionBegin; 1775 *missing = PETSC_FALSE; 1776 if (A->rmap->n > 0 && !ii) { 1777 *missing = PETSC_TRUE; 1778 if (d) *d = 0; 1779 PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n")); 1780 } else { 1781 PetscInt n; 1782 n = PetscMin(A->rmap->n, A->cmap->n); 1783 diag = a->diag; 1784 for (i = 0; i < n; i++) { 1785 if (diag[i] >= ii[i + 1]) { 1786 *missing = PETSC_TRUE; 1787 if (d) *d = i; 1788 PetscCall(PetscInfo(A, "Matrix is missing diagonal number %" PetscInt_FMT "\n", i)); 1789 break; 1790 } 1791 } 1792 } 1793 PetscFunctionReturn(PETSC_SUCCESS); 1794 } 1795 1796 #include <petscblaslapack.h> 1797 #include <petsc/private/kernels/blockinvert.h> 1798 1799 /* 1800 Note that values is allocated externally by the PC and then passed into this routine 1801 */ 1802 PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJ(Mat A, PetscInt nblocks, const PetscInt *bsizes, PetscScalar *diag) 1803 { 1804 PetscInt n = A->rmap->n, i, ncnt = 0, *indx, j, bsizemax = 0, *v_pivots; 1805 PetscBool allowzeropivot, zeropivotdetected = PETSC_FALSE; 1806 const PetscReal shift = 0.0; 1807 PetscInt ipvt[5]; 1808 PetscCount flops = 0; 1809 PetscScalar work[25], *v_work; 1810 1811 PetscFunctionBegin; 1812 allowzeropivot = PetscNot(A->erroriffailure); 1813 for (i = 0; i < nblocks; i++) ncnt += bsizes[i]; 1814 PetscCheck(ncnt == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Total blocksizes %" PetscInt_FMT " doesn't match number matrix rows %" PetscInt_FMT, ncnt, n); 1815 for (i = 0; i < nblocks; i++) bsizemax = PetscMax(bsizemax, bsizes[i]); 1816 PetscCall(PetscMalloc1(bsizemax, &indx)); 1817 if (bsizemax > 7) PetscCall(PetscMalloc2(bsizemax, &v_work, bsizemax, &v_pivots)); 1818 ncnt = 0; 1819 for (i = 0; i < nblocks; i++) { 1820 for (j = 0; j < bsizes[i]; j++) indx[j] = ncnt + j; 1821 PetscCall(MatGetValues(A, bsizes[i], indx, bsizes[i], indx, diag)); 1822 switch (bsizes[i]) { 1823 case 1: 1824 *diag = 1.0 / (*diag); 1825 break; 1826 case 2: 1827 PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected)); 1828 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1829 PetscCall(PetscKernel_A_gets_transpose_A_2(diag)); 1830 break; 1831 case 3: 1832 PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected)); 1833 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1834 PetscCall(PetscKernel_A_gets_transpose_A_3(diag)); 1835 break; 1836 case 4: 1837 PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected)); 1838 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1839 PetscCall(PetscKernel_A_gets_transpose_A_4(diag)); 1840 break; 1841 case 5: 1842 PetscCall(PetscKernel_A_gets_inverse_A_5(diag, ipvt, work, shift, allowzeropivot, &zeropivotdetected)); 1843 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1844 PetscCall(PetscKernel_A_gets_transpose_A_5(diag)); 1845 break; 1846 case 6: 1847 PetscCall(PetscKernel_A_gets_inverse_A_6(diag, shift, allowzeropivot, &zeropivotdetected)); 1848 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1849 PetscCall(PetscKernel_A_gets_transpose_A_6(diag)); 1850 break; 1851 case 7: 1852 PetscCall(PetscKernel_A_gets_inverse_A_7(diag, shift, allowzeropivot, &zeropivotdetected)); 1853 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1854 PetscCall(PetscKernel_A_gets_transpose_A_7(diag)); 1855 break; 1856 default: 1857 PetscCall(PetscKernel_A_gets_inverse_A(bsizes[i], diag, v_pivots, v_work, allowzeropivot, &zeropivotdetected)); 1858 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1859 PetscCall(PetscKernel_A_gets_transpose_A_N(diag, bsizes[i])); 1860 } 1861 ncnt += bsizes[i]; 1862 diag += bsizes[i] * bsizes[i]; 1863 flops += 2 * PetscPowInt(bsizes[i], 3) / 3; 1864 } 1865 PetscCall(PetscLogFlops(flops)); 1866 if (bsizemax > 7) PetscCall(PetscFree2(v_work, v_pivots)); 1867 PetscCall(PetscFree(indx)); 1868 PetscFunctionReturn(PETSC_SUCCESS); 1869 } 1870 1871 /* 1872 Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways 1873 */ 1874 PetscErrorCode MatInvertDiagonal_SeqAIJ(Mat A, PetscScalar omega, PetscScalar fshift) 1875 { 1876 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1877 PetscInt i, *diag, m = A->rmap->n; 1878 const MatScalar *v; 1879 PetscScalar *idiag, *mdiag; 1880 1881 PetscFunctionBegin; 1882 if (a->idiagvalid) PetscFunctionReturn(PETSC_SUCCESS); 1883 PetscCall(MatMarkDiagonal_SeqAIJ(A)); 1884 diag = a->diag; 1885 if (!a->idiag) { PetscCall(PetscMalloc3(m, &a->idiag, m, &a->mdiag, m, &a->ssor_work)); } 1886 1887 mdiag = a->mdiag; 1888 idiag = a->idiag; 1889 PetscCall(MatSeqAIJGetArrayRead(A, &v)); 1890 if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) { 1891 for (i = 0; i < m; i++) { 1892 mdiag[i] = v[diag[i]]; 1893 if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */ 1894 if (PetscRealPart(fshift)) { 1895 PetscCall(PetscInfo(A, "Zero diagonal on row %" PetscInt_FMT "\n", i)); 1896 A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1897 A->factorerror_zeropivot_value = 0.0; 1898 A->factorerror_zeropivot_row = i; 1899 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Zero diagonal on row %" PetscInt_FMT, i); 1900 } 1901 idiag[i] = 1.0 / v[diag[i]]; 1902 } 1903 PetscCall(PetscLogFlops(m)); 1904 } else { 1905 for (i = 0; i < m; i++) { 1906 mdiag[i] = v[diag[i]]; 1907 idiag[i] = omega / (fshift + v[diag[i]]); 1908 } 1909 PetscCall(PetscLogFlops(2.0 * m)); 1910 } 1911 a->idiagvalid = PETSC_TRUE; 1912 PetscCall(MatSeqAIJRestoreArrayRead(A, &v)); 1913 PetscFunctionReturn(PETSC_SUCCESS); 1914 } 1915 1916 #include <../src/mat/impls/aij/seq/ftn-kernels/frelax.h> 1917 PetscErrorCode MatSOR_SeqAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) 1918 { 1919 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1920 PetscScalar *x, d, sum, *t, scale; 1921 const MatScalar *v, *idiag = NULL, *mdiag, *aa; 1922 const PetscScalar *b, *bs, *xb, *ts; 1923 PetscInt n, m = A->rmap->n, i; 1924 const PetscInt *idx, *diag; 1925 1926 PetscFunctionBegin; 1927 if (a->inode.use && a->inode.checked && omega == 1.0 && fshift == 0.0) { 1928 PetscCall(MatSOR_SeqAIJ_Inode(A, bb, omega, flag, fshift, its, lits, xx)); 1929 PetscFunctionReturn(PETSC_SUCCESS); 1930 } 1931 its = its * lits; 1932 1933 if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */ 1934 if (!a->idiagvalid) PetscCall(MatInvertDiagonal_SeqAIJ(A, omega, fshift)); 1935 a->fshift = fshift; 1936 a->omega = omega; 1937 1938 diag = a->diag; 1939 t = a->ssor_work; 1940 idiag = a->idiag; 1941 mdiag = a->mdiag; 1942 1943 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 1944 PetscCall(VecGetArray(xx, &x)); 1945 PetscCall(VecGetArrayRead(bb, &b)); 1946 /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */ 1947 if (flag == SOR_APPLY_UPPER) { 1948 /* apply (U + D/omega) to the vector */ 1949 bs = b; 1950 for (i = 0; i < m; i++) { 1951 d = fshift + mdiag[i]; 1952 n = a->i[i + 1] - diag[i] - 1; 1953 idx = a->j + diag[i] + 1; 1954 v = aa + diag[i] + 1; 1955 sum = b[i] * d / omega; 1956 PetscSparseDensePlusDot(sum, bs, v, idx, n); 1957 x[i] = sum; 1958 } 1959 PetscCall(VecRestoreArray(xx, &x)); 1960 PetscCall(VecRestoreArrayRead(bb, &b)); 1961 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 1962 PetscCall(PetscLogFlops(a->nz)); 1963 PetscFunctionReturn(PETSC_SUCCESS); 1964 } 1965 1966 PetscCheck(flag != SOR_APPLY_LOWER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_LOWER is not implemented"); 1967 if (flag & SOR_EISENSTAT) { 1968 /* Let A = L + U + D; where L is lower triangular, 1969 U is upper triangular, E = D/omega; This routine applies 1970 1971 (L + E)^{-1} A (U + E)^{-1} 1972 1973 to a vector efficiently using Eisenstat's trick. 1974 */ 1975 scale = (2.0 / omega) - 1.0; 1976 1977 /* x = (E + U)^{-1} b */ 1978 for (i = m - 1; i >= 0; i--) { 1979 n = a->i[i + 1] - diag[i] - 1; 1980 idx = a->j + diag[i] + 1; 1981 v = aa + diag[i] + 1; 1982 sum = b[i]; 1983 PetscSparseDenseMinusDot(sum, x, v, idx, n); 1984 x[i] = sum * idiag[i]; 1985 } 1986 1987 /* t = b - (2*E - D)x */ 1988 v = aa; 1989 for (i = 0; i < m; i++) t[i] = b[i] - scale * (v[*diag++]) * x[i]; 1990 1991 /* t = (E + L)^{-1}t */ 1992 ts = t; 1993 diag = a->diag; 1994 for (i = 0; i < m; i++) { 1995 n = diag[i] - a->i[i]; 1996 idx = a->j + a->i[i]; 1997 v = aa + a->i[i]; 1998 sum = t[i]; 1999 PetscSparseDenseMinusDot(sum, ts, v, idx, n); 2000 t[i] = sum * idiag[i]; 2001 /* x = x + t */ 2002 x[i] += t[i]; 2003 } 2004 2005 PetscCall(PetscLogFlops(6.0 * m - 1 + 2.0 * a->nz)); 2006 PetscCall(VecRestoreArray(xx, &x)); 2007 PetscCall(VecRestoreArrayRead(bb, &b)); 2008 PetscFunctionReturn(PETSC_SUCCESS); 2009 } 2010 if (flag & SOR_ZERO_INITIAL_GUESS) { 2011 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 2012 for (i = 0; i < m; i++) { 2013 n = diag[i] - a->i[i]; 2014 idx = a->j + a->i[i]; 2015 v = aa + a->i[i]; 2016 sum = b[i]; 2017 PetscSparseDenseMinusDot(sum, x, v, idx, n); 2018 t[i] = sum; 2019 x[i] = sum * idiag[i]; 2020 } 2021 xb = t; 2022 PetscCall(PetscLogFlops(a->nz)); 2023 } else xb = b; 2024 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 2025 for (i = m - 1; i >= 0; i--) { 2026 n = a->i[i + 1] - diag[i] - 1; 2027 idx = a->j + diag[i] + 1; 2028 v = aa + diag[i] + 1; 2029 sum = xb[i]; 2030 PetscSparseDenseMinusDot(sum, x, v, idx, n); 2031 if (xb == b) { 2032 x[i] = sum * idiag[i]; 2033 } else { 2034 x[i] = (1 - omega) * x[i] + sum * idiag[i]; /* omega in idiag */ 2035 } 2036 } 2037 PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */ 2038 } 2039 its--; 2040 } 2041 while (its--) { 2042 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 2043 for (i = 0; i < m; i++) { 2044 /* lower */ 2045 n = diag[i] - a->i[i]; 2046 idx = a->j + a->i[i]; 2047 v = aa + a->i[i]; 2048 sum = b[i]; 2049 PetscSparseDenseMinusDot(sum, x, v, idx, n); 2050 t[i] = sum; /* save application of the lower-triangular part */ 2051 /* upper */ 2052 n = a->i[i + 1] - diag[i] - 1; 2053 idx = a->j + diag[i] + 1; 2054 v = aa + diag[i] + 1; 2055 PetscSparseDenseMinusDot(sum, x, v, idx, n); 2056 x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */ 2057 } 2058 xb = t; 2059 PetscCall(PetscLogFlops(2.0 * a->nz)); 2060 } else xb = b; 2061 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 2062 for (i = m - 1; i >= 0; i--) { 2063 sum = xb[i]; 2064 if (xb == b) { 2065 /* whole matrix (no checkpointing available) */ 2066 n = a->i[i + 1] - a->i[i]; 2067 idx = a->j + a->i[i]; 2068 v = aa + a->i[i]; 2069 PetscSparseDenseMinusDot(sum, x, v, idx, n); 2070 x[i] = (1. - omega) * x[i] + (sum + mdiag[i] * x[i]) * idiag[i]; 2071 } else { /* lower-triangular part has been saved, so only apply upper-triangular */ 2072 n = a->i[i + 1] - diag[i] - 1; 2073 idx = a->j + diag[i] + 1; 2074 v = aa + diag[i] + 1; 2075 PetscSparseDenseMinusDot(sum, x, v, idx, n); 2076 x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */ 2077 } 2078 } 2079 if (xb == b) { 2080 PetscCall(PetscLogFlops(2.0 * a->nz)); 2081 } else { 2082 PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */ 2083 } 2084 } 2085 } 2086 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 2087 PetscCall(VecRestoreArray(xx, &x)); 2088 PetscCall(VecRestoreArrayRead(bb, &b)); 2089 PetscFunctionReturn(PETSC_SUCCESS); 2090 } 2091 2092 PetscErrorCode MatGetInfo_SeqAIJ(Mat A, MatInfoType flag, MatInfo *info) 2093 { 2094 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2095 2096 PetscFunctionBegin; 2097 info->block_size = 1.0; 2098 info->nz_allocated = a->maxnz; 2099 info->nz_used = a->nz; 2100 info->nz_unneeded = (a->maxnz - a->nz); 2101 info->assemblies = A->num_ass; 2102 info->mallocs = A->info.mallocs; 2103 info->memory = 0; /* REVIEW ME */ 2104 if (A->factortype) { 2105 info->fill_ratio_given = A->info.fill_ratio_given; 2106 info->fill_ratio_needed = A->info.fill_ratio_needed; 2107 info->factor_mallocs = A->info.factor_mallocs; 2108 } else { 2109 info->fill_ratio_given = 0; 2110 info->fill_ratio_needed = 0; 2111 info->factor_mallocs = 0; 2112 } 2113 PetscFunctionReturn(PETSC_SUCCESS); 2114 } 2115 2116 PetscErrorCode MatZeroRows_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 2117 { 2118 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2119 PetscInt i, m = A->rmap->n - 1; 2120 const PetscScalar *xx; 2121 PetscScalar *bb, *aa; 2122 PetscInt d = 0; 2123 2124 PetscFunctionBegin; 2125 if (x && b) { 2126 PetscCall(VecGetArrayRead(x, &xx)); 2127 PetscCall(VecGetArray(b, &bb)); 2128 for (i = 0; i < N; i++) { 2129 PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]); 2130 if (rows[i] >= A->cmap->n) continue; 2131 bb[rows[i]] = diag * xx[rows[i]]; 2132 } 2133 PetscCall(VecRestoreArrayRead(x, &xx)); 2134 PetscCall(VecRestoreArray(b, &bb)); 2135 } 2136 2137 PetscCall(MatSeqAIJGetArray(A, &aa)); 2138 if (a->keepnonzeropattern) { 2139 for (i = 0; i < N; i++) { 2140 PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]); 2141 PetscCall(PetscArrayzero(&aa[a->i[rows[i]]], a->ilen[rows[i]])); 2142 } 2143 if (diag != 0.0) { 2144 for (i = 0; i < N; i++) { 2145 d = rows[i]; 2146 if (rows[i] >= A->cmap->n) continue; 2147 PetscCheck(a->diag[d] < a->i[d + 1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry in the zeroed row %" PetscInt_FMT, d); 2148 } 2149 for (i = 0; i < N; i++) { 2150 if (rows[i] >= A->cmap->n) continue; 2151 aa[a->diag[rows[i]]] = diag; 2152 } 2153 } 2154 } else { 2155 if (diag != 0.0) { 2156 for (i = 0; i < N; i++) { 2157 PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]); 2158 if (a->ilen[rows[i]] > 0) { 2159 if (rows[i] >= A->cmap->n) { 2160 a->ilen[rows[i]] = 0; 2161 } else { 2162 a->ilen[rows[i]] = 1; 2163 aa[a->i[rows[i]]] = diag; 2164 a->j[a->i[rows[i]]] = rows[i]; 2165 } 2166 } else if (rows[i] < A->cmap->n) { /* in case row was completely empty */ 2167 PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES)); 2168 } 2169 } 2170 } else { 2171 for (i = 0; i < N; i++) { 2172 PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]); 2173 a->ilen[rows[i]] = 0; 2174 } 2175 } 2176 A->nonzerostate++; 2177 } 2178 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 2179 PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY); 2180 PetscFunctionReturn(PETSC_SUCCESS); 2181 } 2182 2183 PetscErrorCode MatZeroRowsColumns_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 2184 { 2185 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2186 PetscInt i, j, m = A->rmap->n - 1, d = 0; 2187 PetscBool missing, *zeroed, vecs = PETSC_FALSE; 2188 const PetscScalar *xx; 2189 PetscScalar *bb, *aa; 2190 2191 PetscFunctionBegin; 2192 if (!N) PetscFunctionReturn(PETSC_SUCCESS); 2193 PetscCall(MatSeqAIJGetArray(A, &aa)); 2194 if (x && b) { 2195 PetscCall(VecGetArrayRead(x, &xx)); 2196 PetscCall(VecGetArray(b, &bb)); 2197 vecs = PETSC_TRUE; 2198 } 2199 PetscCall(PetscCalloc1(A->rmap->n, &zeroed)); 2200 for (i = 0; i < N; i++) { 2201 PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]); 2202 PetscCall(PetscArrayzero(&aa[a->i[rows[i]]], a->ilen[rows[i]])); 2203 2204 zeroed[rows[i]] = PETSC_TRUE; 2205 } 2206 for (i = 0; i < A->rmap->n; i++) { 2207 if (!zeroed[i]) { 2208 for (j = a->i[i]; j < a->i[i + 1]; j++) { 2209 if (a->j[j] < A->rmap->n && zeroed[a->j[j]]) { 2210 if (vecs) bb[i] -= aa[j] * xx[a->j[j]]; 2211 aa[j] = 0.0; 2212 } 2213 } 2214 } else if (vecs && i < A->cmap->N) bb[i] = diag * xx[i]; 2215 } 2216 if (x && b) { 2217 PetscCall(VecRestoreArrayRead(x, &xx)); 2218 PetscCall(VecRestoreArray(b, &bb)); 2219 } 2220 PetscCall(PetscFree(zeroed)); 2221 if (diag != 0.0) { 2222 PetscCall(MatMissingDiagonal_SeqAIJ(A, &missing, &d)); 2223 if (missing) { 2224 for (i = 0; i < N; i++) { 2225 if (rows[i] >= A->cmap->N) continue; 2226 PetscCheck(!a->nonew || rows[i] < d, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry in row %" PetscInt_FMT " (%" PetscInt_FMT ")", d, rows[i]); 2227 PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES)); 2228 } 2229 } else { 2230 for (i = 0; i < N; i++) aa[a->diag[rows[i]]] = diag; 2231 } 2232 } 2233 PetscCall(MatSeqAIJRestoreArray(A, &aa)); 2234 PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY); 2235 PetscFunctionReturn(PETSC_SUCCESS); 2236 } 2237 2238 PetscErrorCode MatGetRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 2239 { 2240 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2241 const PetscScalar *aa; 2242 PetscInt *itmp; 2243 2244 PetscFunctionBegin; 2245 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 2246 *nz = a->i[row + 1] - a->i[row]; 2247 if (v) *v = (PetscScalar *)(aa + a->i[row]); 2248 if (idx) { 2249 itmp = a->j + a->i[row]; 2250 if (*nz) *idx = itmp; 2251 else *idx = NULL; 2252 } 2253 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 2254 PetscFunctionReturn(PETSC_SUCCESS); 2255 } 2256 2257 PetscErrorCode MatRestoreRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 2258 { 2259 PetscFunctionBegin; 2260 if (nz) *nz = 0; 2261 if (idx) *idx = NULL; 2262 if (v) *v = NULL; 2263 PetscFunctionReturn(PETSC_SUCCESS); 2264 } 2265 2266 PetscErrorCode MatNorm_SeqAIJ(Mat A, NormType type, PetscReal *nrm) 2267 { 2268 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2269 const MatScalar *v; 2270 PetscReal sum = 0.0; 2271 PetscInt i, j; 2272 2273 PetscFunctionBegin; 2274 PetscCall(MatSeqAIJGetArrayRead(A, &v)); 2275 if (type == NORM_FROBENIUS) { 2276 #if defined(PETSC_USE_REAL___FP16) 2277 PetscBLASInt one = 1, nz = a->nz; 2278 PetscCallBLAS("BLASnrm2", *nrm = BLASnrm2_(&nz, v, &one)); 2279 #else 2280 for (i = 0; i < a->nz; i++) { 2281 sum += PetscRealPart(PetscConj(*v) * (*v)); 2282 v++; 2283 } 2284 *nrm = PetscSqrtReal(sum); 2285 #endif 2286 PetscCall(PetscLogFlops(2.0 * a->nz)); 2287 } else if (type == NORM_1) { 2288 PetscReal *tmp; 2289 PetscInt *jj = a->j; 2290 PetscCall(PetscCalloc1(A->cmap->n + 1, &tmp)); 2291 *nrm = 0.0; 2292 for (j = 0; j < a->nz; j++) { 2293 tmp[*jj++] += PetscAbsScalar(*v); 2294 v++; 2295 } 2296 for (j = 0; j < A->cmap->n; j++) { 2297 if (tmp[j] > *nrm) *nrm = tmp[j]; 2298 } 2299 PetscCall(PetscFree(tmp)); 2300 PetscCall(PetscLogFlops(PetscMax(a->nz - 1, 0))); 2301 } else if (type == NORM_INFINITY) { 2302 *nrm = 0.0; 2303 for (j = 0; j < A->rmap->n; j++) { 2304 const PetscScalar *v2 = v + a->i[j]; 2305 sum = 0.0; 2306 for (i = 0; i < a->i[j + 1] - a->i[j]; i++) { 2307 sum += PetscAbsScalar(*v2); 2308 v2++; 2309 } 2310 if (sum > *nrm) *nrm = sum; 2311 } 2312 PetscCall(PetscLogFlops(PetscMax(a->nz - 1, 0))); 2313 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for two norm"); 2314 PetscCall(MatSeqAIJRestoreArrayRead(A, &v)); 2315 PetscFunctionReturn(PETSC_SUCCESS); 2316 } 2317 2318 PetscErrorCode MatIsTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f) 2319 { 2320 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data; 2321 PetscInt *adx, *bdx, *aii, *bii, *aptr, *bptr; 2322 const MatScalar *va, *vb; 2323 PetscInt ma, na, mb, nb, i; 2324 2325 PetscFunctionBegin; 2326 PetscCall(MatGetSize(A, &ma, &na)); 2327 PetscCall(MatGetSize(B, &mb, &nb)); 2328 if (ma != nb || na != mb) { 2329 *f = PETSC_FALSE; 2330 PetscFunctionReturn(PETSC_SUCCESS); 2331 } 2332 PetscCall(MatSeqAIJGetArrayRead(A, &va)); 2333 PetscCall(MatSeqAIJGetArrayRead(B, &vb)); 2334 aii = aij->i; 2335 bii = bij->i; 2336 adx = aij->j; 2337 bdx = bij->j; 2338 PetscCall(PetscMalloc1(ma, &aptr)); 2339 PetscCall(PetscMalloc1(mb, &bptr)); 2340 for (i = 0; i < ma; i++) aptr[i] = aii[i]; 2341 for (i = 0; i < mb; i++) bptr[i] = bii[i]; 2342 2343 *f = PETSC_TRUE; 2344 for (i = 0; i < ma; i++) { 2345 while (aptr[i] < aii[i + 1]) { 2346 PetscInt idc, idr; 2347 PetscScalar vc, vr; 2348 /* column/row index/value */ 2349 idc = adx[aptr[i]]; 2350 idr = bdx[bptr[idc]]; 2351 vc = va[aptr[i]]; 2352 vr = vb[bptr[idc]]; 2353 if (i != idr || PetscAbsScalar(vc - vr) > tol) { 2354 *f = PETSC_FALSE; 2355 goto done; 2356 } else { 2357 aptr[i]++; 2358 if (B || i != idc) bptr[idc]++; 2359 } 2360 } 2361 } 2362 done: 2363 PetscCall(PetscFree(aptr)); 2364 PetscCall(PetscFree(bptr)); 2365 PetscCall(MatSeqAIJRestoreArrayRead(A, &va)); 2366 PetscCall(MatSeqAIJRestoreArrayRead(B, &vb)); 2367 PetscFunctionReturn(PETSC_SUCCESS); 2368 } 2369 2370 PetscErrorCode MatIsHermitianTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f) 2371 { 2372 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data; 2373 PetscInt *adx, *bdx, *aii, *bii, *aptr, *bptr; 2374 MatScalar *va, *vb; 2375 PetscInt ma, na, mb, nb, i; 2376 2377 PetscFunctionBegin; 2378 PetscCall(MatGetSize(A, &ma, &na)); 2379 PetscCall(MatGetSize(B, &mb, &nb)); 2380 if (ma != nb || na != mb) { 2381 *f = PETSC_FALSE; 2382 PetscFunctionReturn(PETSC_SUCCESS); 2383 } 2384 aii = aij->i; 2385 bii = bij->i; 2386 adx = aij->j; 2387 bdx = bij->j; 2388 va = aij->a; 2389 vb = bij->a; 2390 PetscCall(PetscMalloc1(ma, &aptr)); 2391 PetscCall(PetscMalloc1(mb, &bptr)); 2392 for (i = 0; i < ma; i++) aptr[i] = aii[i]; 2393 for (i = 0; i < mb; i++) bptr[i] = bii[i]; 2394 2395 *f = PETSC_TRUE; 2396 for (i = 0; i < ma; i++) { 2397 while (aptr[i] < aii[i + 1]) { 2398 PetscInt idc, idr; 2399 PetscScalar vc, vr; 2400 /* column/row index/value */ 2401 idc = adx[aptr[i]]; 2402 idr = bdx[bptr[idc]]; 2403 vc = va[aptr[i]]; 2404 vr = vb[bptr[idc]]; 2405 if (i != idr || PetscAbsScalar(vc - PetscConj(vr)) > tol) { 2406 *f = PETSC_FALSE; 2407 goto done; 2408 } else { 2409 aptr[i]++; 2410 if (B || i != idc) bptr[idc]++; 2411 } 2412 } 2413 } 2414 done: 2415 PetscCall(PetscFree(aptr)); 2416 PetscCall(PetscFree(bptr)); 2417 PetscFunctionReturn(PETSC_SUCCESS); 2418 } 2419 2420 PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A, PetscReal tol, PetscBool *f) 2421 { 2422 PetscFunctionBegin; 2423 PetscCall(MatIsTranspose_SeqAIJ(A, A, tol, f)); 2424 PetscFunctionReturn(PETSC_SUCCESS); 2425 } 2426 2427 PetscErrorCode MatIsHermitian_SeqAIJ(Mat A, PetscReal tol, PetscBool *f) 2428 { 2429 PetscFunctionBegin; 2430 PetscCall(MatIsHermitianTranspose_SeqAIJ(A, A, tol, f)); 2431 PetscFunctionReturn(PETSC_SUCCESS); 2432 } 2433 2434 PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A, Vec ll, Vec rr) 2435 { 2436 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2437 const PetscScalar *l, *r; 2438 PetscScalar x; 2439 MatScalar *v; 2440 PetscInt i, j, m = A->rmap->n, n = A->cmap->n, M, nz = a->nz; 2441 const PetscInt *jj; 2442 2443 PetscFunctionBegin; 2444 if (ll) { 2445 /* The local size is used so that VecMPI can be passed to this routine 2446 by MatDiagonalScale_MPIAIJ */ 2447 PetscCall(VecGetLocalSize(ll, &m)); 2448 PetscCheck(m == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector wrong length"); 2449 PetscCall(VecGetArrayRead(ll, &l)); 2450 PetscCall(MatSeqAIJGetArray(A, &v)); 2451 for (i = 0; i < m; i++) { 2452 x = l[i]; 2453 M = a->i[i + 1] - a->i[i]; 2454 for (j = 0; j < M; j++) (*v++) *= x; 2455 } 2456 PetscCall(VecRestoreArrayRead(ll, &l)); 2457 PetscCall(PetscLogFlops(nz)); 2458 PetscCall(MatSeqAIJRestoreArray(A, &v)); 2459 } 2460 if (rr) { 2461 PetscCall(VecGetLocalSize(rr, &n)); 2462 PetscCheck(n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vector wrong length"); 2463 PetscCall(VecGetArrayRead(rr, &r)); 2464 PetscCall(MatSeqAIJGetArray(A, &v)); 2465 jj = a->j; 2466 for (i = 0; i < nz; i++) (*v++) *= r[*jj++]; 2467 PetscCall(MatSeqAIJRestoreArray(A, &v)); 2468 PetscCall(VecRestoreArrayRead(rr, &r)); 2469 PetscCall(PetscLogFlops(nz)); 2470 } 2471 PetscCall(MatSeqAIJInvalidateDiagonal(A)); 2472 PetscFunctionReturn(PETSC_SUCCESS); 2473 } 2474 2475 PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat A, IS isrow, IS iscol, PetscInt csize, MatReuse scall, Mat *B) 2476 { 2477 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *c; 2478 PetscInt *smap, i, k, kstart, kend, oldcols = A->cmap->n, *lens; 2479 PetscInt row, mat_i, *mat_j, tcol, first, step, *mat_ilen, sum, lensi; 2480 const PetscInt *irow, *icol; 2481 const PetscScalar *aa; 2482 PetscInt nrows, ncols; 2483 PetscInt *starts, *j_new, *i_new, *aj = a->j, *ai = a->i, ii, *ailen = a->ilen; 2484 MatScalar *a_new, *mat_a, *c_a; 2485 Mat C; 2486 PetscBool stride; 2487 2488 PetscFunctionBegin; 2489 PetscCall(ISGetIndices(isrow, &irow)); 2490 PetscCall(ISGetLocalSize(isrow, &nrows)); 2491 PetscCall(ISGetLocalSize(iscol, &ncols)); 2492 2493 PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISSTRIDE, &stride)); 2494 if (stride) { 2495 PetscCall(ISStrideGetInfo(iscol, &first, &step)); 2496 } else { 2497 first = 0; 2498 step = 0; 2499 } 2500 if (stride && step == 1) { 2501 /* special case of contiguous rows */ 2502 PetscCall(PetscMalloc2(nrows, &lens, nrows, &starts)); 2503 /* loop over new rows determining lens and starting points */ 2504 for (i = 0; i < nrows; i++) { 2505 kstart = ai[irow[i]]; 2506 kend = kstart + ailen[irow[i]]; 2507 starts[i] = kstart; 2508 for (k = kstart; k < kend; k++) { 2509 if (aj[k] >= first) { 2510 starts[i] = k; 2511 break; 2512 } 2513 } 2514 sum = 0; 2515 while (k < kend) { 2516 if (aj[k++] >= first + ncols) break; 2517 sum++; 2518 } 2519 lens[i] = sum; 2520 } 2521 /* create submatrix */ 2522 if (scall == MAT_REUSE_MATRIX) { 2523 PetscInt n_cols, n_rows; 2524 PetscCall(MatGetSize(*B, &n_rows, &n_cols)); 2525 PetscCheck(n_rows == nrows && n_cols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reused submatrix wrong size"); 2526 PetscCall(MatZeroEntries(*B)); 2527 C = *B; 2528 } else { 2529 PetscInt rbs, cbs; 2530 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 2531 PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE)); 2532 PetscCall(ISGetBlockSize(isrow, &rbs)); 2533 PetscCall(ISGetBlockSize(iscol, &cbs)); 2534 PetscCall(MatSetBlockSizes(C, rbs, cbs)); 2535 PetscCall(MatSetType(C, ((PetscObject)A)->type_name)); 2536 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(C, 0, lens)); 2537 } 2538 c = (Mat_SeqAIJ *)C->data; 2539 2540 /* loop over rows inserting into submatrix */ 2541 PetscCall(MatSeqAIJGetArrayWrite(C, &a_new)); // Not 'a_new = c->a-new', since that raw usage ignores offload state of C 2542 j_new = c->j; 2543 i_new = c->i; 2544 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 2545 for (i = 0; i < nrows; i++) { 2546 ii = starts[i]; 2547 lensi = lens[i]; 2548 for (k = 0; k < lensi; k++) *j_new++ = aj[ii + k] - first; 2549 PetscCall(PetscArraycpy(a_new, aa + starts[i], lensi)); 2550 a_new += lensi; 2551 i_new[i + 1] = i_new[i] + lensi; 2552 c->ilen[i] = lensi; 2553 } 2554 PetscCall(MatSeqAIJRestoreArrayWrite(C, &a_new)); // Set C's offload state properly 2555 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 2556 PetscCall(PetscFree2(lens, starts)); 2557 } else { 2558 PetscCall(ISGetIndices(iscol, &icol)); 2559 PetscCall(PetscCalloc1(oldcols, &smap)); 2560 PetscCall(PetscMalloc1(1 + nrows, &lens)); 2561 for (i = 0; i < ncols; i++) { 2562 PetscCheck(icol[i] < oldcols, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Requesting column beyond largest column icol[%" PetscInt_FMT "] %" PetscInt_FMT " >= A->cmap->n %" PetscInt_FMT, i, icol[i], oldcols); 2563 smap[icol[i]] = i + 1; 2564 } 2565 2566 /* determine lens of each row */ 2567 for (i = 0; i < nrows; i++) { 2568 kstart = ai[irow[i]]; 2569 kend = kstart + a->ilen[irow[i]]; 2570 lens[i] = 0; 2571 for (k = kstart; k < kend; k++) { 2572 if (smap[aj[k]]) lens[i]++; 2573 } 2574 } 2575 /* Create and fill new matrix */ 2576 if (scall == MAT_REUSE_MATRIX) { 2577 PetscBool equal; 2578 2579 c = (Mat_SeqAIJ *)((*B)->data); 2580 PetscCheck((*B)->rmap->n == nrows && (*B)->cmap->n == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size"); 2581 PetscCall(PetscArraycmp(c->ilen, lens, (*B)->rmap->n, &equal)); 2582 PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong no of nonzeros"); 2583 PetscCall(PetscArrayzero(c->ilen, (*B)->rmap->n)); 2584 C = *B; 2585 } else { 2586 PetscInt rbs, cbs; 2587 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 2588 PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE)); 2589 PetscCall(ISGetBlockSize(isrow, &rbs)); 2590 PetscCall(ISGetBlockSize(iscol, &cbs)); 2591 PetscCall(MatSetBlockSizes(C, rbs, cbs)); 2592 PetscCall(MatSetType(C, ((PetscObject)A)->type_name)); 2593 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(C, 0, lens)); 2594 } 2595 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 2596 2597 c = (Mat_SeqAIJ *)(C->data); 2598 PetscCall(MatSeqAIJGetArrayWrite(C, &c_a)); // Not 'c->a', since that raw usage ignores offload state of C 2599 for (i = 0; i < nrows; i++) { 2600 row = irow[i]; 2601 kstart = ai[row]; 2602 kend = kstart + a->ilen[row]; 2603 mat_i = c->i[i]; 2604 mat_j = c->j + mat_i; 2605 mat_a = c_a + mat_i; 2606 mat_ilen = c->ilen + i; 2607 for (k = kstart; k < kend; k++) { 2608 if ((tcol = smap[a->j[k]])) { 2609 *mat_j++ = tcol - 1; 2610 *mat_a++ = aa[k]; 2611 (*mat_ilen)++; 2612 } 2613 } 2614 } 2615 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 2616 /* Free work space */ 2617 PetscCall(ISRestoreIndices(iscol, &icol)); 2618 PetscCall(PetscFree(smap)); 2619 PetscCall(PetscFree(lens)); 2620 /* sort */ 2621 for (i = 0; i < nrows; i++) { 2622 PetscInt ilen; 2623 2624 mat_i = c->i[i]; 2625 mat_j = c->j + mat_i; 2626 mat_a = c_a + mat_i; 2627 ilen = c->ilen[i]; 2628 PetscCall(PetscSortIntWithScalarArray(ilen, mat_j, mat_a)); 2629 } 2630 PetscCall(MatSeqAIJRestoreArrayWrite(C, &c_a)); 2631 } 2632 #if defined(PETSC_HAVE_DEVICE) 2633 PetscCall(MatBindToCPU(C, A->boundtocpu)); 2634 #endif 2635 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 2636 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 2637 2638 PetscCall(ISRestoreIndices(isrow, &irow)); 2639 *B = C; 2640 PetscFunctionReturn(PETSC_SUCCESS); 2641 } 2642 2643 PetscErrorCode MatGetMultiProcBlock_SeqAIJ(Mat mat, MPI_Comm subComm, MatReuse scall, Mat *subMat) 2644 { 2645 Mat B; 2646 2647 PetscFunctionBegin; 2648 if (scall == MAT_INITIAL_MATRIX) { 2649 PetscCall(MatCreate(subComm, &B)); 2650 PetscCall(MatSetSizes(B, mat->rmap->n, mat->cmap->n, mat->rmap->n, mat->cmap->n)); 2651 PetscCall(MatSetBlockSizesFromMats(B, mat, mat)); 2652 PetscCall(MatSetType(B, MATSEQAIJ)); 2653 PetscCall(MatDuplicateNoCreate_SeqAIJ(B, mat, MAT_COPY_VALUES, PETSC_TRUE)); 2654 *subMat = B; 2655 } else { 2656 PetscCall(MatCopy_SeqAIJ(mat, *subMat, SAME_NONZERO_PATTERN)); 2657 } 2658 PetscFunctionReturn(PETSC_SUCCESS); 2659 } 2660 2661 PetscErrorCode MatILUFactor_SeqAIJ(Mat inA, IS row, IS col, const MatFactorInfo *info) 2662 { 2663 Mat_SeqAIJ *a = (Mat_SeqAIJ *)inA->data; 2664 Mat outA; 2665 PetscBool row_identity, col_identity; 2666 2667 PetscFunctionBegin; 2668 PetscCheck(info->levels == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only levels=0 supported for in-place ilu"); 2669 2670 PetscCall(ISIdentity(row, &row_identity)); 2671 PetscCall(ISIdentity(col, &col_identity)); 2672 2673 outA = inA; 2674 outA->factortype = MAT_FACTOR_LU; 2675 PetscCall(PetscFree(inA->solvertype)); 2676 PetscCall(PetscStrallocpy(MATSOLVERPETSC, &inA->solvertype)); 2677 2678 PetscCall(PetscObjectReference((PetscObject)row)); 2679 PetscCall(ISDestroy(&a->row)); 2680 2681 a->row = row; 2682 2683 PetscCall(PetscObjectReference((PetscObject)col)); 2684 PetscCall(ISDestroy(&a->col)); 2685 2686 a->col = col; 2687 2688 /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */ 2689 PetscCall(ISDestroy(&a->icol)); 2690 PetscCall(ISInvertPermutation(col, PETSC_DECIDE, &a->icol)); 2691 2692 if (!a->solve_work) { /* this matrix may have been factored before */ 2693 PetscCall(PetscMalloc1(inA->rmap->n + 1, &a->solve_work)); 2694 } 2695 2696 PetscCall(MatMarkDiagonal_SeqAIJ(inA)); 2697 if (row_identity && col_identity) { 2698 PetscCall(MatLUFactorNumeric_SeqAIJ_inplace(outA, inA, info)); 2699 } else { 2700 PetscCall(MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(outA, inA, info)); 2701 } 2702 PetscFunctionReturn(PETSC_SUCCESS); 2703 } 2704 2705 PetscErrorCode MatScale_SeqAIJ(Mat inA, PetscScalar alpha) 2706 { 2707 Mat_SeqAIJ *a = (Mat_SeqAIJ *)inA->data; 2708 PetscScalar *v; 2709 PetscBLASInt one = 1, bnz; 2710 2711 PetscFunctionBegin; 2712 PetscCall(MatSeqAIJGetArray(inA, &v)); 2713 PetscCall(PetscBLASIntCast(a->nz, &bnz)); 2714 PetscCallBLAS("BLASscal", BLASscal_(&bnz, &alpha, v, &one)); 2715 PetscCall(PetscLogFlops(a->nz)); 2716 PetscCall(MatSeqAIJRestoreArray(inA, &v)); 2717 PetscCall(MatSeqAIJInvalidateDiagonal(inA)); 2718 PetscFunctionReturn(PETSC_SUCCESS); 2719 } 2720 2721 PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *submatj) 2722 { 2723 PetscInt i; 2724 2725 PetscFunctionBegin; 2726 if (!submatj->id) { /* delete data that are linked only to submats[id=0] */ 2727 PetscCall(PetscFree4(submatj->sbuf1, submatj->ptr, submatj->tmp, submatj->ctr)); 2728 2729 for (i = 0; i < submatj->nrqr; ++i) PetscCall(PetscFree(submatj->sbuf2[i])); 2730 PetscCall(PetscFree3(submatj->sbuf2, submatj->req_size, submatj->req_source1)); 2731 2732 if (submatj->rbuf1) { 2733 PetscCall(PetscFree(submatj->rbuf1[0])); 2734 PetscCall(PetscFree(submatj->rbuf1)); 2735 } 2736 2737 for (i = 0; i < submatj->nrqs; ++i) PetscCall(PetscFree(submatj->rbuf3[i])); 2738 PetscCall(PetscFree3(submatj->req_source2, submatj->rbuf2, submatj->rbuf3)); 2739 PetscCall(PetscFree(submatj->pa)); 2740 } 2741 2742 #if defined(PETSC_USE_CTABLE) 2743 PetscCall(PetscHMapIDestroy(&submatj->rmap)); 2744 if (submatj->cmap_loc) PetscCall(PetscFree(submatj->cmap_loc)); 2745 PetscCall(PetscFree(submatj->rmap_loc)); 2746 #else 2747 PetscCall(PetscFree(submatj->rmap)); 2748 #endif 2749 2750 if (!submatj->allcolumns) { 2751 #if defined(PETSC_USE_CTABLE) 2752 PetscCall(PetscHMapIDestroy((PetscHMapI *)&submatj->cmap)); 2753 #else 2754 PetscCall(PetscFree(submatj->cmap)); 2755 #endif 2756 } 2757 PetscCall(PetscFree(submatj->row2proc)); 2758 2759 PetscCall(PetscFree(submatj)); 2760 PetscFunctionReturn(PETSC_SUCCESS); 2761 } 2762 2763 PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat C) 2764 { 2765 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 2766 Mat_SubSppt *submatj = c->submatis1; 2767 2768 PetscFunctionBegin; 2769 PetscCall((*submatj->destroy)(C)); 2770 PetscCall(MatDestroySubMatrix_Private(submatj)); 2771 PetscFunctionReturn(PETSC_SUCCESS); 2772 } 2773 2774 /* Note this has code duplication with MatDestroySubMatrices_SeqBAIJ() */ 2775 PetscErrorCode MatDestroySubMatrices_SeqAIJ(PetscInt n, Mat *mat[]) 2776 { 2777 PetscInt i; 2778 Mat C; 2779 Mat_SeqAIJ *c; 2780 Mat_SubSppt *submatj; 2781 2782 PetscFunctionBegin; 2783 for (i = 0; i < n; i++) { 2784 C = (*mat)[i]; 2785 c = (Mat_SeqAIJ *)C->data; 2786 submatj = c->submatis1; 2787 if (submatj) { 2788 if (--((PetscObject)C)->refct <= 0) { 2789 PetscCall(PetscFree(C->factorprefix)); 2790 PetscCall((*submatj->destroy)(C)); 2791 PetscCall(MatDestroySubMatrix_Private(submatj)); 2792 PetscCall(PetscFree(C->defaultvectype)); 2793 PetscCall(PetscFree(C->defaultrandtype)); 2794 PetscCall(PetscLayoutDestroy(&C->rmap)); 2795 PetscCall(PetscLayoutDestroy(&C->cmap)); 2796 PetscCall(PetscHeaderDestroy(&C)); 2797 } 2798 } else { 2799 PetscCall(MatDestroy(&C)); 2800 } 2801 } 2802 2803 /* Destroy Dummy submatrices created for reuse */ 2804 PetscCall(MatDestroySubMatrices_Dummy(n, mat)); 2805 2806 PetscCall(PetscFree(*mat)); 2807 PetscFunctionReturn(PETSC_SUCCESS); 2808 } 2809 2810 PetscErrorCode MatCreateSubMatrices_SeqAIJ(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *B[]) 2811 { 2812 PetscInt i; 2813 2814 PetscFunctionBegin; 2815 if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B)); 2816 2817 for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_SeqAIJ(A, irow[i], icol[i], PETSC_DECIDE, scall, &(*B)[i])); 2818 PetscFunctionReturn(PETSC_SUCCESS); 2819 } 2820 2821 PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A, PetscInt is_max, IS is[], PetscInt ov) 2822 { 2823 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2824 PetscInt row, i, j, k, l, ll, m, n, *nidx, isz, val; 2825 const PetscInt *idx; 2826 PetscInt start, end, *ai, *aj, bs = (A->rmap->bs > 0 && A->rmap->bs == A->cmap->bs) ? A->rmap->bs : 1; 2827 PetscBT table; 2828 2829 PetscFunctionBegin; 2830 m = A->rmap->n / bs; 2831 ai = a->i; 2832 aj = a->j; 2833 2834 PetscCheck(ov >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "illegal negative overlap value used"); 2835 2836 PetscCall(PetscMalloc1(m + 1, &nidx)); 2837 PetscCall(PetscBTCreate(m, &table)); 2838 2839 for (i = 0; i < is_max; i++) { 2840 /* Initialize the two local arrays */ 2841 isz = 0; 2842 PetscCall(PetscBTMemzero(m, table)); 2843 2844 /* Extract the indices, assume there can be duplicate entries */ 2845 PetscCall(ISGetIndices(is[i], &idx)); 2846 PetscCall(ISGetLocalSize(is[i], &n)); 2847 2848 if (bs > 1) { 2849 /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */ 2850 for (j = 0; j < n; ++j) { 2851 if (!PetscBTLookupSet(table, idx[j] / bs)) nidx[isz++] = idx[j] / bs; 2852 } 2853 PetscCall(ISRestoreIndices(is[i], &idx)); 2854 PetscCall(ISDestroy(&is[i])); 2855 2856 k = 0; 2857 for (j = 0; j < ov; j++) { /* for each overlap */ 2858 n = isz; 2859 for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */ 2860 for (ll = 0; ll < bs; ll++) { 2861 row = bs * nidx[k] + ll; 2862 start = ai[row]; 2863 end = ai[row + 1]; 2864 for (l = start; l < end; l++) { 2865 val = aj[l] / bs; 2866 if (!PetscBTLookupSet(table, val)) nidx[isz++] = val; 2867 } 2868 } 2869 } 2870 } 2871 PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, isz, nidx, PETSC_COPY_VALUES, (is + i))); 2872 } else { 2873 /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */ 2874 for (j = 0; j < n; ++j) { 2875 if (!PetscBTLookupSet(table, idx[j])) nidx[isz++] = idx[j]; 2876 } 2877 PetscCall(ISRestoreIndices(is[i], &idx)); 2878 PetscCall(ISDestroy(&is[i])); 2879 2880 k = 0; 2881 for (j = 0; j < ov; j++) { /* for each overlap */ 2882 n = isz; 2883 for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */ 2884 row = nidx[k]; 2885 start = ai[row]; 2886 end = ai[row + 1]; 2887 for (l = start; l < end; l++) { 2888 val = aj[l]; 2889 if (!PetscBTLookupSet(table, val)) nidx[isz++] = val; 2890 } 2891 } 2892 } 2893 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, PETSC_COPY_VALUES, (is + i))); 2894 } 2895 } 2896 PetscCall(PetscBTDestroy(&table)); 2897 PetscCall(PetscFree(nidx)); 2898 PetscFunctionReturn(PETSC_SUCCESS); 2899 } 2900 2901 PetscErrorCode MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B) 2902 { 2903 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2904 PetscInt i, nz = 0, m = A->rmap->n, n = A->cmap->n; 2905 const PetscInt *row, *col; 2906 PetscInt *cnew, j, *lens; 2907 IS icolp, irowp; 2908 PetscInt *cwork = NULL; 2909 PetscScalar *vwork = NULL; 2910 2911 PetscFunctionBegin; 2912 PetscCall(ISInvertPermutation(rowp, PETSC_DECIDE, &irowp)); 2913 PetscCall(ISGetIndices(irowp, &row)); 2914 PetscCall(ISInvertPermutation(colp, PETSC_DECIDE, &icolp)); 2915 PetscCall(ISGetIndices(icolp, &col)); 2916 2917 /* determine lengths of permuted rows */ 2918 PetscCall(PetscMalloc1(m + 1, &lens)); 2919 for (i = 0; i < m; i++) lens[row[i]] = a->i[i + 1] - a->i[i]; 2920 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 2921 PetscCall(MatSetSizes(*B, m, n, m, n)); 2922 PetscCall(MatSetBlockSizesFromMats(*B, A, A)); 2923 PetscCall(MatSetType(*B, ((PetscObject)A)->type_name)); 2924 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*B, 0, lens)); 2925 PetscCall(PetscFree(lens)); 2926 2927 PetscCall(PetscMalloc1(n, &cnew)); 2928 for (i = 0; i < m; i++) { 2929 PetscCall(MatGetRow_SeqAIJ(A, i, &nz, &cwork, &vwork)); 2930 for (j = 0; j < nz; j++) cnew[j] = col[cwork[j]]; 2931 PetscCall(MatSetValues_SeqAIJ(*B, 1, &row[i], nz, cnew, vwork, INSERT_VALUES)); 2932 PetscCall(MatRestoreRow_SeqAIJ(A, i, &nz, &cwork, &vwork)); 2933 } 2934 PetscCall(PetscFree(cnew)); 2935 2936 (*B)->assembled = PETSC_FALSE; 2937 2938 #if defined(PETSC_HAVE_DEVICE) 2939 PetscCall(MatBindToCPU(*B, A->boundtocpu)); 2940 #endif 2941 PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY)); 2942 PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY)); 2943 PetscCall(ISRestoreIndices(irowp, &row)); 2944 PetscCall(ISRestoreIndices(icolp, &col)); 2945 PetscCall(ISDestroy(&irowp)); 2946 PetscCall(ISDestroy(&icolp)); 2947 if (rowp == colp) PetscCall(MatPropagateSymmetryOptions(A, *B)); 2948 PetscFunctionReturn(PETSC_SUCCESS); 2949 } 2950 2951 PetscErrorCode MatCopy_SeqAIJ(Mat A, Mat B, MatStructure str) 2952 { 2953 PetscFunctionBegin; 2954 /* If the two matrices have the same copy implementation, use fast copy. */ 2955 if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 2956 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2957 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 2958 const PetscScalar *aa; 2959 2960 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 2961 PetscCheck(a->i[A->rmap->n] == b->i[B->rmap->n], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of nonzeros in two matrices are different %" PetscInt_FMT " != %" PetscInt_FMT, a->i[A->rmap->n], b->i[B->rmap->n]); 2962 PetscCall(PetscArraycpy(b->a, aa, a->i[A->rmap->n])); 2963 PetscCall(PetscObjectStateIncrease((PetscObject)B)); 2964 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 2965 } else { 2966 PetscCall(MatCopy_Basic(A, B, str)); 2967 } 2968 PetscFunctionReturn(PETSC_SUCCESS); 2969 } 2970 2971 PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat A, PetscScalar *array[]) 2972 { 2973 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2974 2975 PetscFunctionBegin; 2976 *array = a->a; 2977 PetscFunctionReturn(PETSC_SUCCESS); 2978 } 2979 2980 PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat A, PetscScalar *array[]) 2981 { 2982 PetscFunctionBegin; 2983 *array = NULL; 2984 PetscFunctionReturn(PETSC_SUCCESS); 2985 } 2986 2987 /* 2988 Computes the number of nonzeros per row needed for preallocation when X and Y 2989 have different nonzero structure. 2990 */ 2991 PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt m, const PetscInt *xi, const PetscInt *xj, const PetscInt *yi, const PetscInt *yj, PetscInt *nnz) 2992 { 2993 PetscInt i, j, k, nzx, nzy; 2994 2995 PetscFunctionBegin; 2996 /* Set the number of nonzeros in the new matrix */ 2997 for (i = 0; i < m; i++) { 2998 const PetscInt *xjj = xj + xi[i], *yjj = yj + yi[i]; 2999 nzx = xi[i + 1] - xi[i]; 3000 nzy = yi[i + 1] - yi[i]; 3001 nnz[i] = 0; 3002 for (j = 0, k = 0; j < nzx; j++) { /* Point in X */ 3003 for (; k < nzy && yjj[k] < xjj[j]; k++) nnz[i]++; /* Catch up to X */ 3004 if (k < nzy && yjj[k] == xjj[j]) k++; /* Skip duplicate */ 3005 nnz[i]++; 3006 } 3007 for (; k < nzy; k++) nnz[i]++; 3008 } 3009 PetscFunctionReturn(PETSC_SUCCESS); 3010 } 3011 3012 PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat Y, Mat X, PetscInt *nnz) 3013 { 3014 PetscInt m = Y->rmap->N; 3015 Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data; 3016 Mat_SeqAIJ *y = (Mat_SeqAIJ *)Y->data; 3017 3018 PetscFunctionBegin; 3019 /* Set the number of nonzeros in the new matrix */ 3020 PetscCall(MatAXPYGetPreallocation_SeqX_private(m, x->i, x->j, y->i, y->j, nnz)); 3021 PetscFunctionReturn(PETSC_SUCCESS); 3022 } 3023 3024 PetscErrorCode MatAXPY_SeqAIJ(Mat Y, PetscScalar a, Mat X, MatStructure str) 3025 { 3026 Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data, *y = (Mat_SeqAIJ *)Y->data; 3027 3028 PetscFunctionBegin; 3029 if (str == UNKNOWN_NONZERO_PATTERN || (PetscDefined(USE_DEBUG) && str == SAME_NONZERO_PATTERN)) { 3030 PetscBool e = x->nz == y->nz ? PETSC_TRUE : PETSC_FALSE; 3031 if (e) { 3032 PetscCall(PetscArraycmp(x->i, y->i, Y->rmap->n + 1, &e)); 3033 if (e) { 3034 PetscCall(PetscArraycmp(x->j, y->j, y->nz, &e)); 3035 if (e) str = SAME_NONZERO_PATTERN; 3036 } 3037 } 3038 if (!e) PetscCheck(str != SAME_NONZERO_PATTERN, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatStructure is not SAME_NONZERO_PATTERN"); 3039 } 3040 if (str == SAME_NONZERO_PATTERN) { 3041 const PetscScalar *xa; 3042 PetscScalar *ya, alpha = a; 3043 PetscBLASInt one = 1, bnz; 3044 3045 PetscCall(PetscBLASIntCast(x->nz, &bnz)); 3046 PetscCall(MatSeqAIJGetArray(Y, &ya)); 3047 PetscCall(MatSeqAIJGetArrayRead(X, &xa)); 3048 PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, xa, &one, ya, &one)); 3049 PetscCall(MatSeqAIJRestoreArrayRead(X, &xa)); 3050 PetscCall(MatSeqAIJRestoreArray(Y, &ya)); 3051 PetscCall(PetscLogFlops(2.0 * bnz)); 3052 PetscCall(MatSeqAIJInvalidateDiagonal(Y)); 3053 PetscCall(PetscObjectStateIncrease((PetscObject)Y)); 3054 } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 3055 PetscCall(MatAXPY_Basic(Y, a, X, str)); 3056 } else { 3057 Mat B; 3058 PetscInt *nnz; 3059 PetscCall(PetscMalloc1(Y->rmap->N, &nnz)); 3060 PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &B)); 3061 PetscCall(PetscObjectSetName((PetscObject)B, ((PetscObject)Y)->name)); 3062 PetscCall(MatSetLayouts(B, Y->rmap, Y->cmap)); 3063 PetscCall(MatSetType(B, ((PetscObject)Y)->type_name)); 3064 PetscCall(MatAXPYGetPreallocation_SeqAIJ(Y, X, nnz)); 3065 PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz)); 3066 PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str)); 3067 PetscCall(MatHeaderMerge(Y, &B)); 3068 PetscCall(MatSeqAIJCheckInode(Y)); 3069 PetscCall(PetscFree(nnz)); 3070 } 3071 PetscFunctionReturn(PETSC_SUCCESS); 3072 } 3073 3074 PETSC_INTERN PetscErrorCode MatConjugate_SeqAIJ(Mat mat) 3075 { 3076 #if defined(PETSC_USE_COMPLEX) 3077 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 3078 PetscInt i, nz; 3079 PetscScalar *a; 3080 3081 PetscFunctionBegin; 3082 nz = aij->nz; 3083 PetscCall(MatSeqAIJGetArray(mat, &a)); 3084 for (i = 0; i < nz; i++) a[i] = PetscConj(a[i]); 3085 PetscCall(MatSeqAIJRestoreArray(mat, &a)); 3086 #else 3087 PetscFunctionBegin; 3088 #endif 3089 PetscFunctionReturn(PETSC_SUCCESS); 3090 } 3091 3092 PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[]) 3093 { 3094 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 3095 PetscInt i, j, m = A->rmap->n, *ai, *aj, ncols, n; 3096 PetscReal atmp; 3097 PetscScalar *x; 3098 const MatScalar *aa, *av; 3099 3100 PetscFunctionBegin; 3101 PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3102 PetscCall(MatSeqAIJGetArrayRead(A, &av)); 3103 aa = av; 3104 ai = a->i; 3105 aj = a->j; 3106 3107 PetscCall(VecSet(v, 0.0)); 3108 PetscCall(VecGetArrayWrite(v, &x)); 3109 PetscCall(VecGetLocalSize(v, &n)); 3110 PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector"); 3111 for (i = 0; i < m; i++) { 3112 ncols = ai[1] - ai[0]; 3113 ai++; 3114 for (j = 0; j < ncols; j++) { 3115 atmp = PetscAbsScalar(*aa); 3116 if (PetscAbsScalar(x[i]) < atmp) { 3117 x[i] = atmp; 3118 if (idx) idx[i] = *aj; 3119 } 3120 aa++; 3121 aj++; 3122 } 3123 } 3124 PetscCall(VecRestoreArrayWrite(v, &x)); 3125 PetscCall(MatSeqAIJRestoreArrayRead(A, &av)); 3126 PetscFunctionReturn(PETSC_SUCCESS); 3127 } 3128 3129 PetscErrorCode MatGetRowMax_SeqAIJ(Mat A, Vec v, PetscInt idx[]) 3130 { 3131 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 3132 PetscInt i, j, m = A->rmap->n, *ai, *aj, ncols, n; 3133 PetscScalar *x; 3134 const MatScalar *aa, *av; 3135 3136 PetscFunctionBegin; 3137 PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3138 PetscCall(MatSeqAIJGetArrayRead(A, &av)); 3139 aa = av; 3140 ai = a->i; 3141 aj = a->j; 3142 3143 PetscCall(VecSet(v, 0.0)); 3144 PetscCall(VecGetArrayWrite(v, &x)); 3145 PetscCall(VecGetLocalSize(v, &n)); 3146 PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector"); 3147 for (i = 0; i < m; i++) { 3148 ncols = ai[1] - ai[0]; 3149 ai++; 3150 if (ncols == A->cmap->n) { /* row is dense */ 3151 x[i] = *aa; 3152 if (idx) idx[i] = 0; 3153 } else { /* row is sparse so already KNOW maximum is 0.0 or higher */ 3154 x[i] = 0.0; 3155 if (idx) { 3156 for (j = 0; j < ncols; j++) { /* find first implicit 0.0 in the row */ 3157 if (aj[j] > j) { 3158 idx[i] = j; 3159 break; 3160 } 3161 } 3162 /* in case first implicit 0.0 in the row occurs at ncols-th column */ 3163 if (j == ncols && j < A->cmap->n) idx[i] = j; 3164 } 3165 } 3166 for (j = 0; j < ncols; j++) { 3167 if (PetscRealPart(x[i]) < PetscRealPart(*aa)) { 3168 x[i] = *aa; 3169 if (idx) idx[i] = *aj; 3170 } 3171 aa++; 3172 aj++; 3173 } 3174 } 3175 PetscCall(VecRestoreArrayWrite(v, &x)); 3176 PetscCall(MatSeqAIJRestoreArrayRead(A, &av)); 3177 PetscFunctionReturn(PETSC_SUCCESS); 3178 } 3179 3180 PetscErrorCode MatGetRowMinAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[]) 3181 { 3182 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 3183 PetscInt i, j, m = A->rmap->n, *ai, *aj, ncols, n; 3184 PetscScalar *x; 3185 const MatScalar *aa, *av; 3186 3187 PetscFunctionBegin; 3188 PetscCall(MatSeqAIJGetArrayRead(A, &av)); 3189 aa = av; 3190 ai = a->i; 3191 aj = a->j; 3192 3193 PetscCall(VecSet(v, 0.0)); 3194 PetscCall(VecGetArrayWrite(v, &x)); 3195 PetscCall(VecGetLocalSize(v, &n)); 3196 PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector, %" PetscInt_FMT " vs. %" PetscInt_FMT " rows", m, n); 3197 for (i = 0; i < m; i++) { 3198 ncols = ai[1] - ai[0]; 3199 ai++; 3200 if (ncols == A->cmap->n) { /* row is dense */ 3201 x[i] = *aa; 3202 if (idx) idx[i] = 0; 3203 } else { /* row is sparse so already KNOW minimum is 0.0 or higher */ 3204 x[i] = 0.0; 3205 if (idx) { /* find first implicit 0.0 in the row */ 3206 for (j = 0; j < ncols; j++) { 3207 if (aj[j] > j) { 3208 idx[i] = j; 3209 break; 3210 } 3211 } 3212 /* in case first implicit 0.0 in the row occurs at ncols-th column */ 3213 if (j == ncols && j < A->cmap->n) idx[i] = j; 3214 } 3215 } 3216 for (j = 0; j < ncols; j++) { 3217 if (PetscAbsScalar(x[i]) > PetscAbsScalar(*aa)) { 3218 x[i] = *aa; 3219 if (idx) idx[i] = *aj; 3220 } 3221 aa++; 3222 aj++; 3223 } 3224 } 3225 PetscCall(VecRestoreArrayWrite(v, &x)); 3226 PetscCall(MatSeqAIJRestoreArrayRead(A, &av)); 3227 PetscFunctionReturn(PETSC_SUCCESS); 3228 } 3229 3230 PetscErrorCode MatGetRowMin_SeqAIJ(Mat A, Vec v, PetscInt idx[]) 3231 { 3232 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 3233 PetscInt i, j, m = A->rmap->n, ncols, n; 3234 const PetscInt *ai, *aj; 3235 PetscScalar *x; 3236 const MatScalar *aa, *av; 3237 3238 PetscFunctionBegin; 3239 PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3240 PetscCall(MatSeqAIJGetArrayRead(A, &av)); 3241 aa = av; 3242 ai = a->i; 3243 aj = a->j; 3244 3245 PetscCall(VecSet(v, 0.0)); 3246 PetscCall(VecGetArrayWrite(v, &x)); 3247 PetscCall(VecGetLocalSize(v, &n)); 3248 PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector"); 3249 for (i = 0; i < m; i++) { 3250 ncols = ai[1] - ai[0]; 3251 ai++; 3252 if (ncols == A->cmap->n) { /* row is dense */ 3253 x[i] = *aa; 3254 if (idx) idx[i] = 0; 3255 } else { /* row is sparse so already KNOW minimum is 0.0 or lower */ 3256 x[i] = 0.0; 3257 if (idx) { /* find first implicit 0.0 in the row */ 3258 for (j = 0; j < ncols; j++) { 3259 if (aj[j] > j) { 3260 idx[i] = j; 3261 break; 3262 } 3263 } 3264 /* in case first implicit 0.0 in the row occurs at ncols-th column */ 3265 if (j == ncols && j < A->cmap->n) idx[i] = j; 3266 } 3267 } 3268 for (j = 0; j < ncols; j++) { 3269 if (PetscRealPart(x[i]) > PetscRealPart(*aa)) { 3270 x[i] = *aa; 3271 if (idx) idx[i] = *aj; 3272 } 3273 aa++; 3274 aj++; 3275 } 3276 } 3277 PetscCall(VecRestoreArrayWrite(v, &x)); 3278 PetscCall(MatSeqAIJRestoreArrayRead(A, &av)); 3279 PetscFunctionReturn(PETSC_SUCCESS); 3280 } 3281 3282 PetscErrorCode MatInvertBlockDiagonal_SeqAIJ(Mat A, const PetscScalar **values) 3283 { 3284 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 3285 PetscInt i, bs = PetscAbs(A->rmap->bs), mbs = A->rmap->n / bs, ipvt[5], bs2 = bs * bs, *v_pivots, ij[7], *IJ, j; 3286 MatScalar *diag, work[25], *v_work; 3287 const PetscReal shift = 0.0; 3288 PetscBool allowzeropivot, zeropivotdetected = PETSC_FALSE; 3289 3290 PetscFunctionBegin; 3291 allowzeropivot = PetscNot(A->erroriffailure); 3292 if (a->ibdiagvalid) { 3293 if (values) *values = a->ibdiag; 3294 PetscFunctionReturn(PETSC_SUCCESS); 3295 } 3296 PetscCall(MatMarkDiagonal_SeqAIJ(A)); 3297 if (!a->ibdiag) { PetscCall(PetscMalloc1(bs2 * mbs, &a->ibdiag)); } 3298 diag = a->ibdiag; 3299 if (values) *values = a->ibdiag; 3300 /* factor and invert each block */ 3301 switch (bs) { 3302 case 1: 3303 for (i = 0; i < mbs; i++) { 3304 PetscCall(MatGetValues(A, 1, &i, 1, &i, diag + i)); 3305 if (PetscAbsScalar(diag[i] + shift) < PETSC_MACHINE_EPSILON) { 3306 if (allowzeropivot) { 3307 A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3308 A->factorerror_zeropivot_value = PetscAbsScalar(diag[i]); 3309 A->factorerror_zeropivot_row = i; 3310 PetscCall(PetscInfo(A, "Zero pivot, row %" PetscInt_FMT " pivot %g tolerance %g\n", i, (double)PetscAbsScalar(diag[i]), (double)PETSC_MACHINE_EPSILON)); 3311 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT " pivot %g tolerance %g", i, (double)PetscAbsScalar(diag[i]), (double)PETSC_MACHINE_EPSILON); 3312 } 3313 diag[i] = (PetscScalar)1.0 / (diag[i] + shift); 3314 } 3315 break; 3316 case 2: 3317 for (i = 0; i < mbs; i++) { 3318 ij[0] = 2 * i; 3319 ij[1] = 2 * i + 1; 3320 PetscCall(MatGetValues(A, 2, ij, 2, ij, diag)); 3321 PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected)); 3322 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3323 PetscCall(PetscKernel_A_gets_transpose_A_2(diag)); 3324 diag += 4; 3325 } 3326 break; 3327 case 3: 3328 for (i = 0; i < mbs; i++) { 3329 ij[0] = 3 * i; 3330 ij[1] = 3 * i + 1; 3331 ij[2] = 3 * i + 2; 3332 PetscCall(MatGetValues(A, 3, ij, 3, ij, diag)); 3333 PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected)); 3334 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3335 PetscCall(PetscKernel_A_gets_transpose_A_3(diag)); 3336 diag += 9; 3337 } 3338 break; 3339 case 4: 3340 for (i = 0; i < mbs; i++) { 3341 ij[0] = 4 * i; 3342 ij[1] = 4 * i + 1; 3343 ij[2] = 4 * i + 2; 3344 ij[3] = 4 * i + 3; 3345 PetscCall(MatGetValues(A, 4, ij, 4, ij, diag)); 3346 PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected)); 3347 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3348 PetscCall(PetscKernel_A_gets_transpose_A_4(diag)); 3349 diag += 16; 3350 } 3351 break; 3352 case 5: 3353 for (i = 0; i < mbs; i++) { 3354 ij[0] = 5 * i; 3355 ij[1] = 5 * i + 1; 3356 ij[2] = 5 * i + 2; 3357 ij[3] = 5 * i + 3; 3358 ij[4] = 5 * i + 4; 3359 PetscCall(MatGetValues(A, 5, ij, 5, ij, diag)); 3360 PetscCall(PetscKernel_A_gets_inverse_A_5(diag, ipvt, work, shift, allowzeropivot, &zeropivotdetected)); 3361 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3362 PetscCall(PetscKernel_A_gets_transpose_A_5(diag)); 3363 diag += 25; 3364 } 3365 break; 3366 case 6: 3367 for (i = 0; i < mbs; i++) { 3368 ij[0] = 6 * i; 3369 ij[1] = 6 * i + 1; 3370 ij[2] = 6 * i + 2; 3371 ij[3] = 6 * i + 3; 3372 ij[4] = 6 * i + 4; 3373 ij[5] = 6 * i + 5; 3374 PetscCall(MatGetValues(A, 6, ij, 6, ij, diag)); 3375 PetscCall(PetscKernel_A_gets_inverse_A_6(diag, shift, allowzeropivot, &zeropivotdetected)); 3376 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3377 PetscCall(PetscKernel_A_gets_transpose_A_6(diag)); 3378 diag += 36; 3379 } 3380 break; 3381 case 7: 3382 for (i = 0; i < mbs; i++) { 3383 ij[0] = 7 * i; 3384 ij[1] = 7 * i + 1; 3385 ij[2] = 7 * i + 2; 3386 ij[3] = 7 * i + 3; 3387 ij[4] = 7 * i + 4; 3388 ij[5] = 7 * i + 5; 3389 ij[6] = 7 * i + 6; 3390 PetscCall(MatGetValues(A, 7, ij, 7, ij, diag)); 3391 PetscCall(PetscKernel_A_gets_inverse_A_7(diag, shift, allowzeropivot, &zeropivotdetected)); 3392 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3393 PetscCall(PetscKernel_A_gets_transpose_A_7(diag)); 3394 diag += 49; 3395 } 3396 break; 3397 default: 3398 PetscCall(PetscMalloc3(bs, &v_work, bs, &v_pivots, bs, &IJ)); 3399 for (i = 0; i < mbs; i++) { 3400 for (j = 0; j < bs; j++) IJ[j] = bs * i + j; 3401 PetscCall(MatGetValues(A, bs, IJ, bs, IJ, diag)); 3402 PetscCall(PetscKernel_A_gets_inverse_A(bs, diag, v_pivots, v_work, allowzeropivot, &zeropivotdetected)); 3403 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3404 PetscCall(PetscKernel_A_gets_transpose_A_N(diag, bs)); 3405 diag += bs2; 3406 } 3407 PetscCall(PetscFree3(v_work, v_pivots, IJ)); 3408 } 3409 a->ibdiagvalid = PETSC_TRUE; 3410 PetscFunctionReturn(PETSC_SUCCESS); 3411 } 3412 3413 static PetscErrorCode MatSetRandom_SeqAIJ(Mat x, PetscRandom rctx) 3414 { 3415 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data; 3416 PetscScalar a, *aa; 3417 PetscInt m, n, i, j, col; 3418 3419 PetscFunctionBegin; 3420 if (!x->assembled) { 3421 PetscCall(MatGetSize(x, &m, &n)); 3422 for (i = 0; i < m; i++) { 3423 for (j = 0; j < aij->imax[i]; j++) { 3424 PetscCall(PetscRandomGetValue(rctx, &a)); 3425 col = (PetscInt)(n * PetscRealPart(a)); 3426 PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES)); 3427 } 3428 } 3429 } else { 3430 PetscCall(MatSeqAIJGetArrayWrite(x, &aa)); 3431 for (i = 0; i < aij->nz; i++) PetscCall(PetscRandomGetValue(rctx, aa + i)); 3432 PetscCall(MatSeqAIJRestoreArrayWrite(x, &aa)); 3433 } 3434 PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY)); 3435 PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY)); 3436 PetscFunctionReturn(PETSC_SUCCESS); 3437 } 3438 3439 /* Like MatSetRandom_SeqAIJ, but do not set values on columns in range of [low, high) */ 3440 PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat x, PetscInt low, PetscInt high, PetscRandom rctx) 3441 { 3442 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data; 3443 PetscScalar a; 3444 PetscInt m, n, i, j, col, nskip; 3445 3446 PetscFunctionBegin; 3447 nskip = high - low; 3448 PetscCall(MatGetSize(x, &m, &n)); 3449 n -= nskip; /* shrink number of columns where nonzeros can be set */ 3450 for (i = 0; i < m; i++) { 3451 for (j = 0; j < aij->imax[i]; j++) { 3452 PetscCall(PetscRandomGetValue(rctx, &a)); 3453 col = (PetscInt)(n * PetscRealPart(a)); 3454 if (col >= low) col += nskip; /* shift col rightward to skip the hole */ 3455 PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES)); 3456 } 3457 } 3458 PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY)); 3459 PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY)); 3460 PetscFunctionReturn(PETSC_SUCCESS); 3461 } 3462 3463 static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ, 3464 MatGetRow_SeqAIJ, 3465 MatRestoreRow_SeqAIJ, 3466 MatMult_SeqAIJ, 3467 /* 4*/ MatMultAdd_SeqAIJ, 3468 MatMultTranspose_SeqAIJ, 3469 MatMultTransposeAdd_SeqAIJ, 3470 NULL, 3471 NULL, 3472 NULL, 3473 /* 10*/ NULL, 3474 MatLUFactor_SeqAIJ, 3475 NULL, 3476 MatSOR_SeqAIJ, 3477 MatTranspose_SeqAIJ, 3478 /*1 5*/ MatGetInfo_SeqAIJ, 3479 MatEqual_SeqAIJ, 3480 MatGetDiagonal_SeqAIJ, 3481 MatDiagonalScale_SeqAIJ, 3482 MatNorm_SeqAIJ, 3483 /* 20*/ NULL, 3484 MatAssemblyEnd_SeqAIJ, 3485 MatSetOption_SeqAIJ, 3486 MatZeroEntries_SeqAIJ, 3487 /* 24*/ MatZeroRows_SeqAIJ, 3488 NULL, 3489 NULL, 3490 NULL, 3491 NULL, 3492 /* 29*/ MatSetUp_Seq_Hash, 3493 NULL, 3494 NULL, 3495 NULL, 3496 NULL, 3497 /* 34*/ MatDuplicate_SeqAIJ, 3498 NULL, 3499 NULL, 3500 MatILUFactor_SeqAIJ, 3501 NULL, 3502 /* 39*/ MatAXPY_SeqAIJ, 3503 MatCreateSubMatrices_SeqAIJ, 3504 MatIncreaseOverlap_SeqAIJ, 3505 MatGetValues_SeqAIJ, 3506 MatCopy_SeqAIJ, 3507 /* 44*/ MatGetRowMax_SeqAIJ, 3508 MatScale_SeqAIJ, 3509 MatShift_SeqAIJ, 3510 MatDiagonalSet_SeqAIJ, 3511 MatZeroRowsColumns_SeqAIJ, 3512 /* 49*/ MatSetRandom_SeqAIJ, 3513 MatGetRowIJ_SeqAIJ, 3514 MatRestoreRowIJ_SeqAIJ, 3515 MatGetColumnIJ_SeqAIJ, 3516 MatRestoreColumnIJ_SeqAIJ, 3517 /* 54*/ MatFDColoringCreate_SeqXAIJ, 3518 NULL, 3519 NULL, 3520 MatPermute_SeqAIJ, 3521 NULL, 3522 /* 59*/ NULL, 3523 MatDestroy_SeqAIJ, 3524 MatView_SeqAIJ, 3525 NULL, 3526 NULL, 3527 /* 64*/ NULL, 3528 MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ, 3529 NULL, 3530 NULL, 3531 NULL, 3532 /* 69*/ MatGetRowMaxAbs_SeqAIJ, 3533 MatGetRowMinAbs_SeqAIJ, 3534 NULL, 3535 NULL, 3536 NULL, 3537 /* 74*/ NULL, 3538 MatFDColoringApply_AIJ, 3539 NULL, 3540 NULL, 3541 NULL, 3542 /* 79*/ MatFindZeroDiagonals_SeqAIJ, 3543 NULL, 3544 NULL, 3545 NULL, 3546 MatLoad_SeqAIJ, 3547 /* 84*/ MatIsSymmetric_SeqAIJ, 3548 MatIsHermitian_SeqAIJ, 3549 NULL, 3550 NULL, 3551 NULL, 3552 /* 89*/ NULL, 3553 NULL, 3554 MatMatMultNumeric_SeqAIJ_SeqAIJ, 3555 NULL, 3556 NULL, 3557 /* 94*/ MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy, 3558 NULL, 3559 NULL, 3560 MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ, 3561 NULL, 3562 /* 99*/ MatProductSetFromOptions_SeqAIJ, 3563 NULL, 3564 NULL, 3565 MatConjugate_SeqAIJ, 3566 NULL, 3567 /*104*/ MatSetValuesRow_SeqAIJ, 3568 MatRealPart_SeqAIJ, 3569 MatImaginaryPart_SeqAIJ, 3570 NULL, 3571 NULL, 3572 /*109*/ MatMatSolve_SeqAIJ, 3573 NULL, 3574 MatGetRowMin_SeqAIJ, 3575 NULL, 3576 MatMissingDiagonal_SeqAIJ, 3577 /*114*/ NULL, 3578 NULL, 3579 NULL, 3580 NULL, 3581 NULL, 3582 /*119*/ NULL, 3583 NULL, 3584 NULL, 3585 NULL, 3586 MatGetMultiProcBlock_SeqAIJ, 3587 /*124*/ MatFindNonzeroRows_SeqAIJ, 3588 MatGetColumnReductions_SeqAIJ, 3589 MatInvertBlockDiagonal_SeqAIJ, 3590 MatInvertVariableBlockDiagonal_SeqAIJ, 3591 NULL, 3592 /*129*/ NULL, 3593 NULL, 3594 NULL, 3595 MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ, 3596 MatTransposeColoringCreate_SeqAIJ, 3597 /*134*/ MatTransColoringApplySpToDen_SeqAIJ, 3598 MatTransColoringApplyDenToSp_SeqAIJ, 3599 NULL, 3600 NULL, 3601 MatRARtNumeric_SeqAIJ_SeqAIJ, 3602 /*139*/ NULL, 3603 NULL, 3604 NULL, 3605 MatFDColoringSetUp_SeqXAIJ, 3606 MatFindOffBlockDiagonalEntries_SeqAIJ, 3607 MatCreateMPIMatConcatenateSeqMat_SeqAIJ, 3608 /*145*/ MatDestroySubMatrices_SeqAIJ, 3609 NULL, 3610 NULL, 3611 MatCreateGraph_Simple_AIJ, 3612 NULL, 3613 /*150*/ MatTransposeSymbolic_SeqAIJ, 3614 MatEliminateZeros_SeqAIJ}; 3615 3616 PetscErrorCode MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat, PetscInt *indices) 3617 { 3618 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 3619 PetscInt i, nz, n; 3620 3621 PetscFunctionBegin; 3622 nz = aij->maxnz; 3623 n = mat->rmap->n; 3624 for (i = 0; i < nz; i++) aij->j[i] = indices[i]; 3625 aij->nz = nz; 3626 for (i = 0; i < n; i++) aij->ilen[i] = aij->imax[i]; 3627 PetscFunctionReturn(PETSC_SUCCESS); 3628 } 3629 3630 /* 3631 * Given a sparse matrix with global column indices, compact it by using a local column space. 3632 * The result matrix helps saving memory in other algorithms, such as MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable() 3633 */ 3634 PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat mat, ISLocalToGlobalMapping *mapping) 3635 { 3636 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 3637 PetscHMapI gid1_lid1; 3638 PetscHashIter tpos; 3639 PetscInt gid, lid, i, ec, nz = aij->nz; 3640 PetscInt *garray, *jj = aij->j; 3641 3642 PetscFunctionBegin; 3643 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3644 PetscAssertPointer(mapping, 2); 3645 /* use a table */ 3646 PetscCall(PetscHMapICreateWithSize(mat->rmap->n, &gid1_lid1)); 3647 ec = 0; 3648 for (i = 0; i < nz; i++) { 3649 PetscInt data, gid1 = jj[i] + 1; 3650 PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &data)); 3651 if (!data) { 3652 /* one based table */ 3653 PetscCall(PetscHMapISet(gid1_lid1, gid1, ++ec)); 3654 } 3655 } 3656 /* form array of columns we need */ 3657 PetscCall(PetscMalloc1(ec, &garray)); 3658 PetscHashIterBegin(gid1_lid1, tpos); 3659 while (!PetscHashIterAtEnd(gid1_lid1, tpos)) { 3660 PetscHashIterGetKey(gid1_lid1, tpos, gid); 3661 PetscHashIterGetVal(gid1_lid1, tpos, lid); 3662 PetscHashIterNext(gid1_lid1, tpos); 3663 gid--; 3664 lid--; 3665 garray[lid] = gid; 3666 } 3667 PetscCall(PetscSortInt(ec, garray)); /* sort, and rebuild */ 3668 PetscCall(PetscHMapIClear(gid1_lid1)); 3669 for (i = 0; i < ec; i++) PetscCall(PetscHMapISet(gid1_lid1, garray[i] + 1, i + 1)); 3670 /* compact out the extra columns in B */ 3671 for (i = 0; i < nz; i++) { 3672 PetscInt gid1 = jj[i] + 1; 3673 PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &lid)); 3674 lid--; 3675 jj[i] = lid; 3676 } 3677 PetscCall(PetscLayoutDestroy(&mat->cmap)); 3678 PetscCall(PetscHMapIDestroy(&gid1_lid1)); 3679 PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)mat), ec, ec, 1, &mat->cmap)); 3680 PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, mat->cmap->bs, mat->cmap->n, garray, PETSC_OWN_POINTER, mapping)); 3681 PetscCall(ISLocalToGlobalMappingSetType(*mapping, ISLOCALTOGLOBALMAPPINGHASH)); 3682 PetscFunctionReturn(PETSC_SUCCESS); 3683 } 3684 3685 /*@ 3686 MatSeqAIJSetColumnIndices - Set the column indices for all the rows 3687 in the matrix. 3688 3689 Input Parameters: 3690 + mat - the `MATSEQAIJ` matrix 3691 - indices - the column indices 3692 3693 Level: advanced 3694 3695 Notes: 3696 This can be called if you have precomputed the nonzero structure of the 3697 matrix and want to provide it to the matrix object to improve the performance 3698 of the `MatSetValues()` operation. 3699 3700 You MUST have set the correct numbers of nonzeros per row in the call to 3701 `MatCreateSeqAIJ()`, and the columns indices MUST be sorted. 3702 3703 MUST be called before any calls to `MatSetValues()` 3704 3705 The indices should start with zero, not one. 3706 3707 .seealso: [](ch_matrices), `Mat`, `MATSEQAIJ` 3708 @*/ 3709 PetscErrorCode MatSeqAIJSetColumnIndices(Mat mat, PetscInt *indices) 3710 { 3711 PetscFunctionBegin; 3712 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3713 PetscAssertPointer(indices, 2); 3714 PetscUseMethod(mat, "MatSeqAIJSetColumnIndices_C", (Mat, PetscInt *), (mat, indices)); 3715 PetscFunctionReturn(PETSC_SUCCESS); 3716 } 3717 3718 PetscErrorCode MatStoreValues_SeqAIJ(Mat mat) 3719 { 3720 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 3721 size_t nz = aij->i[mat->rmap->n]; 3722 3723 PetscFunctionBegin; 3724 PetscCheck(aij->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 3725 3726 /* allocate space for values if not already there */ 3727 if (!aij->saved_values) { PetscCall(PetscMalloc1(nz + 1, &aij->saved_values)); } 3728 3729 /* copy values over */ 3730 PetscCall(PetscArraycpy(aij->saved_values, aij->a, nz)); 3731 PetscFunctionReturn(PETSC_SUCCESS); 3732 } 3733 3734 /*@ 3735 MatStoreValues - Stashes a copy of the matrix values; this allows reusing of the linear part of a Jacobian, while recomputing only the 3736 nonlinear portion. 3737 3738 Logically Collect 3739 3740 Input Parameter: 3741 . mat - the matrix (currently only `MATAIJ` matrices support this option) 3742 3743 Level: advanced 3744 3745 Example Usage: 3746 .vb 3747 Using SNES 3748 Create Jacobian matrix 3749 Set linear terms into matrix 3750 Apply boundary conditions to matrix, at this time matrix must have 3751 final nonzero structure (i.e. setting the nonlinear terms and applying 3752 boundary conditions again will not change the nonzero structure 3753 MatSetOption(mat, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE); 3754 MatStoreValues(mat); 3755 Call SNESSetJacobian() with matrix 3756 In your Jacobian routine 3757 MatRetrieveValues(mat); 3758 Set nonlinear terms in matrix 3759 3760 Without `SNESSolve()`, i.e. when you handle nonlinear solve yourself: 3761 // build linear portion of Jacobian 3762 MatSetOption(mat, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE); 3763 MatStoreValues(mat); 3764 loop over nonlinear iterations 3765 MatRetrieveValues(mat); 3766 // call MatSetValues(mat,...) to set nonliner portion of Jacobian 3767 // call MatAssemblyBegin/End() on matrix 3768 Solve linear system with Jacobian 3769 endloop 3770 .ve 3771 3772 Notes: 3773 Matrix must already be assembled before calling this routine 3774 Must set the matrix option `MatSetOption`(mat,`MAT_NEW_NONZERO_LOCATIONS`,`PETSC_FALSE`); before 3775 calling this routine. 3776 3777 When this is called multiple times it overwrites the previous set of stored values 3778 and does not allocated additional space. 3779 3780 .seealso: [](ch_matrices), `Mat`, `MatRetrieveValues()` 3781 @*/ 3782 PetscErrorCode MatStoreValues(Mat mat) 3783 { 3784 PetscFunctionBegin; 3785 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3786 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3787 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3788 PetscUseMethod(mat, "MatStoreValues_C", (Mat), (mat)); 3789 PetscFunctionReturn(PETSC_SUCCESS); 3790 } 3791 3792 PetscErrorCode MatRetrieveValues_SeqAIJ(Mat mat) 3793 { 3794 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 3795 PetscInt nz = aij->i[mat->rmap->n]; 3796 3797 PetscFunctionBegin; 3798 PetscCheck(aij->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 3799 PetscCheck(aij->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first"); 3800 /* copy values over */ 3801 PetscCall(PetscArraycpy(aij->a, aij->saved_values, nz)); 3802 PetscFunctionReturn(PETSC_SUCCESS); 3803 } 3804 3805 /*@ 3806 MatRetrieveValues - Retrieves the copy of the matrix values that was stored with `MatStoreValues()` 3807 3808 Logically Collect 3809 3810 Input Parameter: 3811 . mat - the matrix (currently only `MATAIJ` matrices support this option) 3812 3813 Level: advanced 3814 3815 .seealso: [](ch_matrices), `Mat`, `MatStoreValues()` 3816 @*/ 3817 PetscErrorCode MatRetrieveValues(Mat mat) 3818 { 3819 PetscFunctionBegin; 3820 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3821 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3822 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3823 PetscUseMethod(mat, "MatRetrieveValues_C", (Mat), (mat)); 3824 PetscFunctionReturn(PETSC_SUCCESS); 3825 } 3826 3827 /*@C 3828 MatCreateSeqAIJ - Creates a sparse matrix in `MATSEQAIJ` (compressed row) format 3829 (the default parallel PETSc format). For good matrix assembly performance 3830 the user should preallocate the matrix storage by setting the parameter `nz` 3831 (or the array `nnz`). 3832 3833 Collective 3834 3835 Input Parameters: 3836 + comm - MPI communicator, set to `PETSC_COMM_SELF` 3837 . m - number of rows 3838 . n - number of columns 3839 . nz - number of nonzeros per row (same for all rows) 3840 - nnz - array containing the number of nonzeros in the various rows 3841 (possibly different for each row) or NULL 3842 3843 Output Parameter: 3844 . A - the matrix 3845 3846 Options Database Keys: 3847 + -mat_no_inode - Do not use inodes 3848 - -mat_inode_limit <limit> - Sets inode limit (max limit=5) 3849 3850 Level: intermediate 3851 3852 Notes: 3853 If `nnz` is given then `nz` is ignored 3854 3855 The `MATSEQAIJ` format, also called 3856 compressed row storage, is fully compatible with standard Fortran 3857 storage. That is, the stored row and column indices can begin at 3858 either one (as in Fortran) or zero. 3859 3860 Specify the preallocated storage with either `nz` or `nnz` (not both). 3861 Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory 3862 allocation. 3863 3864 By default, this format uses inodes (identical nodes) when possible, to 3865 improve numerical efficiency of matrix-vector products and solves. We 3866 search for consecutive rows with the same nonzero structure, thereby 3867 reusing matrix information to achieve increased efficiency. 3868 3869 .seealso: [](ch_matrices), `Mat`, [Sparse Matrix Creation](sec_matsparse), `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()` 3870 @*/ 3871 PetscErrorCode MatCreateSeqAIJ(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A) 3872 { 3873 PetscFunctionBegin; 3874 PetscCall(MatCreate(comm, A)); 3875 PetscCall(MatSetSizes(*A, m, n, m, n)); 3876 PetscCall(MatSetType(*A, MATSEQAIJ)); 3877 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*A, nz, nnz)); 3878 PetscFunctionReturn(PETSC_SUCCESS); 3879 } 3880 3881 /*@C 3882 MatSeqAIJSetPreallocation - For good matrix assembly performance 3883 the user should preallocate the matrix storage by setting the parameter nz 3884 (or the array nnz). By setting these parameters accurately, performance 3885 during matrix assembly can be increased by more than a factor of 50. 3886 3887 Collective 3888 3889 Input Parameters: 3890 + B - The matrix 3891 . nz - number of nonzeros per row (same for all rows) 3892 - nnz - array containing the number of nonzeros in the various rows 3893 (possibly different for each row) or NULL 3894 3895 Options Database Keys: 3896 + -mat_no_inode - Do not use inodes 3897 - -mat_inode_limit <limit> - Sets inode limit (max limit=5) 3898 3899 Level: intermediate 3900 3901 Notes: 3902 If `nnz` is given then `nz` is ignored 3903 3904 The `MATSEQAIJ` format also called 3905 compressed row storage, is fully compatible with standard Fortran 3906 storage. That is, the stored row and column indices can begin at 3907 either one (as in Fortran) or zero. See the users' manual for details. 3908 3909 Specify the preallocated storage with either `nz` or `nnz` (not both). 3910 Set nz = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory 3911 allocation. 3912 3913 You can call `MatGetInfo()` to get information on how effective the preallocation was; 3914 for example the fields mallocs,nz_allocated,nz_used,nz_unneeded; 3915 You can also run with the option -info and look for messages with the string 3916 malloc in them to see if additional memory allocation was needed. 3917 3918 Developer Notes: 3919 Use nz of `MAT_SKIP_ALLOCATION` to not allocate any space for the matrix 3920 entries or columns indices 3921 3922 By default, this format uses inodes (identical nodes) when possible, to 3923 improve numerical efficiency of matrix-vector products and solves. We 3924 search for consecutive rows with the same nonzero structure, thereby 3925 reusing matrix information to achieve increased efficiency. 3926 3927 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`, `MatGetInfo()`, 3928 `MatSeqAIJSetTotalPreallocation()` 3929 @*/ 3930 PetscErrorCode MatSeqAIJSetPreallocation(Mat B, PetscInt nz, const PetscInt nnz[]) 3931 { 3932 PetscFunctionBegin; 3933 PetscValidHeaderSpecific(B, MAT_CLASSID, 1); 3934 PetscValidType(B, 1); 3935 PetscTryMethod(B, "MatSeqAIJSetPreallocation_C", (Mat, PetscInt, const PetscInt[]), (B, nz, nnz)); 3936 PetscFunctionReturn(PETSC_SUCCESS); 3937 } 3938 3939 PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat B, PetscInt nz, const PetscInt *nnz) 3940 { 3941 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 3942 PetscBool skipallocation = PETSC_FALSE, realalloc = PETSC_FALSE; 3943 PetscInt i; 3944 3945 PetscFunctionBegin; 3946 if (B->hash_active) { 3947 B->ops[0] = b->cops; 3948 PetscCall(PetscHMapIJVDestroy(&b->ht)); 3949 PetscCall(PetscFree(b->dnz)); 3950 B->hash_active = PETSC_FALSE; 3951 } 3952 if (nz >= 0 || nnz) realalloc = PETSC_TRUE; 3953 if (nz == MAT_SKIP_ALLOCATION) { 3954 skipallocation = PETSC_TRUE; 3955 nz = 0; 3956 } 3957 PetscCall(PetscLayoutSetUp(B->rmap)); 3958 PetscCall(PetscLayoutSetUp(B->cmap)); 3959 3960 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 3961 PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz); 3962 if (PetscUnlikelyDebug(nnz)) { 3963 for (i = 0; i < B->rmap->n; i++) { 3964 PetscCheck(nnz[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT, i, nnz[i]); 3965 PetscCheck(nnz[i] <= B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be greater than row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " rowlength %" PetscInt_FMT, i, nnz[i], B->cmap->n); 3966 } 3967 } 3968 3969 B->preallocated = PETSC_TRUE; 3970 if (!skipallocation) { 3971 if (!b->imax) { PetscCall(PetscMalloc1(B->rmap->n, &b->imax)); } 3972 if (!b->ilen) { 3973 /* b->ilen will count nonzeros in each row so far. */ 3974 PetscCall(PetscCalloc1(B->rmap->n, &b->ilen)); 3975 } else { 3976 PetscCall(PetscMemzero(b->ilen, B->rmap->n * sizeof(PetscInt))); 3977 } 3978 if (!b->ipre) PetscCall(PetscMalloc1(B->rmap->n, &b->ipre)); 3979 if (!nnz) { 3980 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10; 3981 else if (nz < 0) nz = 1; 3982 nz = PetscMin(nz, B->cmap->n); 3983 for (i = 0; i < B->rmap->n; i++) b->imax[i] = nz; 3984 nz = nz * B->rmap->n; 3985 } else { 3986 PetscInt64 nz64 = 0; 3987 for (i = 0; i < B->rmap->n; i++) { 3988 b->imax[i] = nnz[i]; 3989 nz64 += nnz[i]; 3990 } 3991 PetscCall(PetscIntCast(nz64, &nz)); 3992 } 3993 3994 /* allocate the matrix space */ 3995 /* FIXME: should B's old memory be unlogged? */ 3996 PetscCall(MatSeqXAIJFreeAIJ(B, &b->a, &b->j, &b->i)); 3997 if (B->structure_only) { 3998 PetscCall(PetscMalloc1(nz, &b->j)); 3999 PetscCall(PetscMalloc1(B->rmap->n + 1, &b->i)); 4000 } else { 4001 PetscCall(PetscMalloc3(nz, &b->a, nz, &b->j, B->rmap->n + 1, &b->i)); 4002 } 4003 b->i[0] = 0; 4004 for (i = 1; i < B->rmap->n + 1; i++) b->i[i] = b->i[i - 1] + b->imax[i - 1]; 4005 if (B->structure_only) { 4006 b->singlemalloc = PETSC_FALSE; 4007 b->free_a = PETSC_FALSE; 4008 } else { 4009 b->singlemalloc = PETSC_TRUE; 4010 b->free_a = PETSC_TRUE; 4011 } 4012 b->free_ij = PETSC_TRUE; 4013 } else { 4014 b->free_a = PETSC_FALSE; 4015 b->free_ij = PETSC_FALSE; 4016 } 4017 4018 if (b->ipre && nnz != b->ipre && b->imax) { 4019 /* reserve user-requested sparsity */ 4020 PetscCall(PetscArraycpy(b->ipre, b->imax, B->rmap->n)); 4021 } 4022 4023 b->nz = 0; 4024 b->maxnz = nz; 4025 B->info.nz_unneeded = (double)b->maxnz; 4026 if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 4027 B->was_assembled = PETSC_FALSE; 4028 B->assembled = PETSC_FALSE; 4029 /* We simply deem preallocation has changed nonzero state. Updating the state 4030 will give clients (like AIJKokkos) a chance to know something has happened. 4031 */ 4032 B->nonzerostate++; 4033 PetscFunctionReturn(PETSC_SUCCESS); 4034 } 4035 4036 PetscErrorCode MatResetPreallocation_SeqAIJ(Mat A) 4037 { 4038 Mat_SeqAIJ *a; 4039 PetscInt i; 4040 PetscBool skipreset; 4041 4042 PetscFunctionBegin; 4043 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 4044 4045 /* Check local size. If zero, then return */ 4046 if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS); 4047 4048 a = (Mat_SeqAIJ *)A->data; 4049 /* if no saved info, we error out */ 4050 PetscCheck(a->ipre, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "No saved preallocation info "); 4051 4052 PetscCheck(a->i && a->imax && a->ilen, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Memory info is incomplete, and can not reset preallocation "); 4053 4054 PetscCall(PetscArraycmp(a->ipre, a->ilen, A->rmap->n, &skipreset)); 4055 if (!skipreset) { 4056 PetscCall(PetscArraycpy(a->imax, a->ipre, A->rmap->n)); 4057 PetscCall(PetscArrayzero(a->ilen, A->rmap->n)); 4058 a->i[0] = 0; 4059 for (i = 1; i < A->rmap->n + 1; i++) a->i[i] = a->i[i - 1] + a->imax[i - 1]; 4060 A->preallocated = PETSC_TRUE; 4061 a->nz = 0; 4062 a->maxnz = a->i[A->rmap->n]; 4063 A->info.nz_unneeded = (double)a->maxnz; 4064 A->was_assembled = PETSC_FALSE; 4065 A->assembled = PETSC_FALSE; 4066 } 4067 PetscFunctionReturn(PETSC_SUCCESS); 4068 } 4069 4070 /*@ 4071 MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in `MATSEQAIJ` format. 4072 4073 Input Parameters: 4074 + B - the matrix 4075 . i - the indices into j for the start of each row (starts with zero) 4076 . j - the column indices for each row (starts with zero) these must be sorted for each row 4077 - v - optional values in the matrix 4078 4079 Level: developer 4080 4081 Notes: 4082 The `i`,`j`,`v` values are COPIED with this routine; to avoid the copy use `MatCreateSeqAIJWithArrays()` 4083 4084 This routine may be called multiple times with different nonzero patterns (or the same nonzero pattern). The nonzero 4085 structure will be the union of all the previous nonzero structures. 4086 4087 Developer Notes: 4088 An optimization could be added to the implementation where it checks if the `i`, and `j` are identical to the current `i` and `j` and 4089 then just copies the `v` values directly with `PetscMemcpy()`. 4090 4091 This routine could also take a `PetscCopyMode` argument to allow sharing the values instead of always copying them. 4092 4093 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatSeqAIJSetPreallocation()`, `MATSEQAIJ`, `MatResetPreallocation()` 4094 @*/ 4095 PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B, const PetscInt i[], const PetscInt j[], const PetscScalar v[]) 4096 { 4097 PetscFunctionBegin; 4098 PetscValidHeaderSpecific(B, MAT_CLASSID, 1); 4099 PetscValidType(B, 1); 4100 PetscTryMethod(B, "MatSeqAIJSetPreallocationCSR_C", (Mat, const PetscInt[], const PetscInt[], const PetscScalar[]), (B, i, j, v)); 4101 PetscFunctionReturn(PETSC_SUCCESS); 4102 } 4103 4104 PetscErrorCode MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B, const PetscInt Ii[], const PetscInt J[], const PetscScalar v[]) 4105 { 4106 PetscInt i; 4107 PetscInt m, n; 4108 PetscInt nz; 4109 PetscInt *nnz; 4110 4111 PetscFunctionBegin; 4112 PetscCheck(Ii[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %" PetscInt_FMT, Ii[0]); 4113 4114 PetscCall(PetscLayoutSetUp(B->rmap)); 4115 PetscCall(PetscLayoutSetUp(B->cmap)); 4116 4117 PetscCall(MatGetSize(B, &m, &n)); 4118 PetscCall(PetscMalloc1(m + 1, &nnz)); 4119 for (i = 0; i < m; i++) { 4120 nz = Ii[i + 1] - Ii[i]; 4121 PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT, i, nz); 4122 nnz[i] = nz; 4123 } 4124 PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz)); 4125 PetscCall(PetscFree(nnz)); 4126 4127 for (i = 0; i < m; i++) PetscCall(MatSetValues_SeqAIJ(B, 1, &i, Ii[i + 1] - Ii[i], J + Ii[i], v ? v + Ii[i] : NULL, INSERT_VALUES)); 4128 4129 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 4130 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 4131 4132 PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 4133 PetscFunctionReturn(PETSC_SUCCESS); 4134 } 4135 4136 /*@ 4137 MatSeqAIJKron - Computes `C`, the Kronecker product of `A` and `B`. 4138 4139 Input Parameters: 4140 + A - left-hand side matrix 4141 . B - right-hand side matrix 4142 - reuse - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 4143 4144 Output Parameter: 4145 . C - Kronecker product of `A` and `B` 4146 4147 Level: intermediate 4148 4149 Note: 4150 `MAT_REUSE_MATRIX` can only be used when the nonzero structure of the product matrix has not changed from that last call to `MatSeqAIJKron()`. 4151 4152 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MATKAIJ`, `MatReuse` 4153 @*/ 4154 PetscErrorCode MatSeqAIJKron(Mat A, Mat B, MatReuse reuse, Mat *C) 4155 { 4156 PetscFunctionBegin; 4157 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 4158 PetscValidType(A, 1); 4159 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 4160 PetscValidType(B, 2); 4161 PetscAssertPointer(C, 4); 4162 if (reuse == MAT_REUSE_MATRIX) { 4163 PetscValidHeaderSpecific(*C, MAT_CLASSID, 4); 4164 PetscValidType(*C, 4); 4165 } 4166 PetscTryMethod(A, "MatSeqAIJKron_C", (Mat, Mat, MatReuse, Mat *), (A, B, reuse, C)); 4167 PetscFunctionReturn(PETSC_SUCCESS); 4168 } 4169 4170 PetscErrorCode MatSeqAIJKron_SeqAIJ(Mat A, Mat B, MatReuse reuse, Mat *C) 4171 { 4172 Mat newmat; 4173 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 4174 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 4175 PetscScalar *v; 4176 const PetscScalar *aa, *ba; 4177 PetscInt *i, *j, m, n, p, q, nnz = 0, am = A->rmap->n, bm = B->rmap->n, an = A->cmap->n, bn = B->cmap->n; 4178 PetscBool flg; 4179 4180 PetscFunctionBegin; 4181 PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4182 PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4183 PetscCheck(!B->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4184 PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4185 PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &flg)); 4186 PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatType %s", ((PetscObject)B)->type_name); 4187 PetscCheck(reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatReuse %d", (int)reuse); 4188 if (reuse == MAT_INITIAL_MATRIX) { 4189 PetscCall(PetscMalloc2(am * bm + 1, &i, a->i[am] * b->i[bm], &j)); 4190 PetscCall(MatCreate(PETSC_COMM_SELF, &newmat)); 4191 PetscCall(MatSetSizes(newmat, am * bm, an * bn, am * bm, an * bn)); 4192 PetscCall(MatSetType(newmat, MATAIJ)); 4193 i[0] = 0; 4194 for (m = 0; m < am; ++m) { 4195 for (p = 0; p < bm; ++p) { 4196 i[m * bm + p + 1] = i[m * bm + p] + (a->i[m + 1] - a->i[m]) * (b->i[p + 1] - b->i[p]); 4197 for (n = a->i[m]; n < a->i[m + 1]; ++n) { 4198 for (q = b->i[p]; q < b->i[p + 1]; ++q) j[nnz++] = a->j[n] * bn + b->j[q]; 4199 } 4200 } 4201 } 4202 PetscCall(MatSeqAIJSetPreallocationCSR(newmat, i, j, NULL)); 4203 *C = newmat; 4204 PetscCall(PetscFree2(i, j)); 4205 nnz = 0; 4206 } 4207 PetscCall(MatSeqAIJGetArray(*C, &v)); 4208 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 4209 PetscCall(MatSeqAIJGetArrayRead(B, &ba)); 4210 for (m = 0; m < am; ++m) { 4211 for (p = 0; p < bm; ++p) { 4212 for (n = a->i[m]; n < a->i[m + 1]; ++n) { 4213 for (q = b->i[p]; q < b->i[p + 1]; ++q) v[nnz++] = aa[n] * ba[q]; 4214 } 4215 } 4216 } 4217 PetscCall(MatSeqAIJRestoreArray(*C, &v)); 4218 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 4219 PetscCall(MatSeqAIJRestoreArrayRead(B, &ba)); 4220 PetscFunctionReturn(PETSC_SUCCESS); 4221 } 4222 4223 #include <../src/mat/impls/dense/seq/dense.h> 4224 #include <petsc/private/kernels/petscaxpy.h> 4225 4226 /* 4227 Computes (B'*A')' since computing B*A directly is untenable 4228 4229 n p p 4230 [ ] [ ] [ ] 4231 m [ A ] * n [ B ] = m [ C ] 4232 [ ] [ ] [ ] 4233 4234 */ 4235 PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat A, Mat B, Mat C) 4236 { 4237 Mat_SeqDense *sub_a = (Mat_SeqDense *)A->data; 4238 Mat_SeqAIJ *sub_b = (Mat_SeqAIJ *)B->data; 4239 Mat_SeqDense *sub_c = (Mat_SeqDense *)C->data; 4240 PetscInt i, j, n, m, q, p; 4241 const PetscInt *ii, *idx; 4242 const PetscScalar *b, *a, *a_q; 4243 PetscScalar *c, *c_q; 4244 PetscInt clda = sub_c->lda; 4245 PetscInt alda = sub_a->lda; 4246 4247 PetscFunctionBegin; 4248 m = A->rmap->n; 4249 n = A->cmap->n; 4250 p = B->cmap->n; 4251 a = sub_a->v; 4252 b = sub_b->a; 4253 c = sub_c->v; 4254 if (clda == m) { 4255 PetscCall(PetscArrayzero(c, m * p)); 4256 } else { 4257 for (j = 0; j < p; j++) 4258 for (i = 0; i < m; i++) c[j * clda + i] = 0.0; 4259 } 4260 ii = sub_b->i; 4261 idx = sub_b->j; 4262 for (i = 0; i < n; i++) { 4263 q = ii[i + 1] - ii[i]; 4264 while (q-- > 0) { 4265 c_q = c + clda * (*idx); 4266 a_q = a + alda * i; 4267 PetscKernelAXPY(c_q, *b, a_q, m); 4268 idx++; 4269 b++; 4270 } 4271 } 4272 PetscFunctionReturn(PETSC_SUCCESS); 4273 } 4274 4275 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat A, Mat B, PetscReal fill, Mat C) 4276 { 4277 PetscInt m = A->rmap->n, n = B->cmap->n; 4278 PetscBool cisdense; 4279 4280 PetscFunctionBegin; 4281 PetscCheck(A->cmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "A->cmap->n %" PetscInt_FMT " != B->rmap->n %" PetscInt_FMT, A->cmap->n, B->rmap->n); 4282 PetscCall(MatSetSizes(C, m, n, m, n)); 4283 PetscCall(MatSetBlockSizesFromMats(C, A, B)); 4284 PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATSEQDENSECUDA, MATSEQDENSEHIP, "")); 4285 if (!cisdense) PetscCall(MatSetType(C, MATDENSE)); 4286 PetscCall(MatSetUp(C)); 4287 4288 C->ops->matmultnumeric = MatMatMultNumeric_SeqDense_SeqAIJ; 4289 PetscFunctionReturn(PETSC_SUCCESS); 4290 } 4291 4292 /*MC 4293 MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices, 4294 based on compressed sparse row format. 4295 4296 Options Database Key: 4297 . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions() 4298 4299 Level: beginner 4300 4301 Notes: 4302 `MatSetValues()` may be called for this matrix type with a `NULL` argument for the numerical values, 4303 in this case the values associated with the rows and columns one passes in are set to zero 4304 in the matrix 4305 4306 `MatSetOptions`(,`MAT_STRUCTURE_ONLY`,`PETSC_TRUE`) may be called for this matrix type. In this no 4307 space is allocated for the nonzero entries and any entries passed with `MatSetValues()` are ignored 4308 4309 Developer Note: 4310 It would be nice if all matrix formats supported passing `NULL` in for the numerical values 4311 4312 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatSetFromOptions()`, `MatSetType()`, `MatCreate()`, `MatType`, `MATSELL`, `MATSEQSELL`, `MATMPISELL` 4313 M*/ 4314 4315 /*MC 4316 MATAIJ - MATAIJ = "aij" - A matrix type to be used for sparse matrices. 4317 4318 This matrix type is identical to `MATSEQAIJ` when constructed with a single process communicator, 4319 and `MATMPIAIJ` otherwise. As a result, for single process communicators, 4320 `MatSeqAIJSetPreallocation()` is supported, and similarly `MatMPIAIJSetPreallocation()` is supported 4321 for communicators controlling multiple processes. It is recommended that you call both of 4322 the above preallocation routines for simplicity. 4323 4324 Options Database Key: 4325 . -mat_type aij - sets the matrix type to "aij" during a call to `MatSetFromOptions()` 4326 4327 Level: beginner 4328 4329 Note: 4330 Subclasses include `MATAIJCUSPARSE`, `MATAIJPERM`, `MATAIJSELL`, `MATAIJMKL`, `MATAIJCRL`, and also automatically switches over to use inodes when 4331 enough exist. 4332 4333 .seealso: [](ch_matrices), `Mat`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MATMPIAIJ`, `MATSELL`, `MATSEQSELL`, `MATMPISELL` 4334 M*/ 4335 4336 /*MC 4337 MATAIJCRL - MATAIJCRL = "aijcrl" - A matrix type to be used for sparse matrices. 4338 4339 Options Database Key: 4340 . -mat_type aijcrl - sets the matrix type to "aijcrl" during a call to `MatSetFromOptions()` 4341 4342 Level: beginner 4343 4344 Note: 4345 This matrix type is identical to `MATSEQAIJCRL` when constructed with a single process communicator, 4346 and `MATMPIAIJCRL` otherwise. As a result, for single process communicators, 4347 `MatSeqAIJSetPreallocation()` is supported, and similarly `MatMPIAIJSetPreallocation()` is supported 4348 for communicators controlling multiple processes. It is recommended that you call both of 4349 the above preallocation routines for simplicity. 4350 4351 .seealso: [](ch_matrices), `Mat`, `MatCreateMPIAIJCRL`, `MATSEQAIJCRL`, `MATMPIAIJCRL`, `MATSEQAIJCRL`, `MATMPIAIJCRL` 4352 M*/ 4353 4354 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *); 4355 #if defined(PETSC_HAVE_ELEMENTAL) 4356 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *); 4357 #endif 4358 #if defined(PETSC_HAVE_SCALAPACK) 4359 PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *); 4360 #endif 4361 #if defined(PETSC_HAVE_HYPRE) 4362 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType, MatReuse, Mat *); 4363 #endif 4364 4365 PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat, MatType, MatReuse, Mat *); 4366 PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat, MatType, MatReuse, Mat *); 4367 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_IS_XAIJ(Mat); 4368 4369 /*@C 4370 MatSeqAIJGetArray - gives read/write access to the array where the data for a `MATSEQAIJ` matrix is stored 4371 4372 Not Collective 4373 4374 Input Parameter: 4375 . A - a `MATSEQAIJ` matrix 4376 4377 Output Parameter: 4378 . array - pointer to the data 4379 4380 Level: intermediate 4381 4382 Fortran Notes: 4383 `MatSeqAIJGetArray()` Fortran binding is deprecated (since PETSc 3.19), use `MatSeqAIJGetArrayF90()` 4384 4385 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()` 4386 @*/ 4387 PetscErrorCode MatSeqAIJGetArray(Mat A, PetscScalar **array) 4388 { 4389 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 4390 4391 PetscFunctionBegin; 4392 if (aij->ops->getarray) { 4393 PetscCall((*aij->ops->getarray)(A, array)); 4394 } else { 4395 *array = aij->a; 4396 } 4397 PetscFunctionReturn(PETSC_SUCCESS); 4398 } 4399 4400 /*@C 4401 MatSeqAIJRestoreArray - returns access to the array where the data for a `MATSEQAIJ` matrix is stored obtained by `MatSeqAIJGetArray()` 4402 4403 Not Collective 4404 4405 Input Parameters: 4406 + A - a `MATSEQAIJ` matrix 4407 - array - pointer to the data 4408 4409 Level: intermediate 4410 4411 Fortran Notes: 4412 `MatSeqAIJRestoreArray()` Fortran binding is deprecated (since PETSc 3.19), use `MatSeqAIJRestoreArrayF90()` 4413 4414 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayF90()` 4415 @*/ 4416 PetscErrorCode MatSeqAIJRestoreArray(Mat A, PetscScalar **array) 4417 { 4418 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 4419 4420 PetscFunctionBegin; 4421 if (aij->ops->restorearray) { 4422 PetscCall((*aij->ops->restorearray)(A, array)); 4423 } else { 4424 *array = NULL; 4425 } 4426 PetscCall(MatSeqAIJInvalidateDiagonal(A)); 4427 PetscCall(PetscObjectStateIncrease((PetscObject)A)); 4428 PetscFunctionReturn(PETSC_SUCCESS); 4429 } 4430 4431 /*@C 4432 MatSeqAIJGetArrayRead - gives read-only access to the array where the data for a `MATSEQAIJ` matrix is stored 4433 4434 Not Collective; No Fortran Support 4435 4436 Input Parameter: 4437 . A - a `MATSEQAIJ` matrix 4438 4439 Output Parameter: 4440 . array - pointer to the data 4441 4442 Level: intermediate 4443 4444 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()` 4445 @*/ 4446 PetscErrorCode MatSeqAIJGetArrayRead(Mat A, const PetscScalar **array) 4447 { 4448 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 4449 4450 PetscFunctionBegin; 4451 if (aij->ops->getarrayread) { 4452 PetscCall((*aij->ops->getarrayread)(A, array)); 4453 } else { 4454 *array = aij->a; 4455 } 4456 PetscFunctionReturn(PETSC_SUCCESS); 4457 } 4458 4459 /*@C 4460 MatSeqAIJRestoreArrayRead - restore the read-only access array obtained from `MatSeqAIJGetArrayRead()` 4461 4462 Not Collective; No Fortran Support 4463 4464 Input Parameter: 4465 . A - a `MATSEQAIJ` matrix 4466 4467 Output Parameter: 4468 . array - pointer to the data 4469 4470 Level: intermediate 4471 4472 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()` 4473 @*/ 4474 PetscErrorCode MatSeqAIJRestoreArrayRead(Mat A, const PetscScalar **array) 4475 { 4476 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 4477 4478 PetscFunctionBegin; 4479 if (aij->ops->restorearrayread) { 4480 PetscCall((*aij->ops->restorearrayread)(A, array)); 4481 } else { 4482 *array = NULL; 4483 } 4484 PetscFunctionReturn(PETSC_SUCCESS); 4485 } 4486 4487 /*@C 4488 MatSeqAIJGetArrayWrite - gives write-only access to the array where the data for a `MATSEQAIJ` matrix is stored 4489 4490 Not Collective; No Fortran Support 4491 4492 Input Parameter: 4493 . A - a `MATSEQAIJ` matrix 4494 4495 Output Parameter: 4496 . array - pointer to the data 4497 4498 Level: intermediate 4499 4500 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()` 4501 @*/ 4502 PetscErrorCode MatSeqAIJGetArrayWrite(Mat A, PetscScalar **array) 4503 { 4504 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 4505 4506 PetscFunctionBegin; 4507 if (aij->ops->getarraywrite) { 4508 PetscCall((*aij->ops->getarraywrite)(A, array)); 4509 } else { 4510 *array = aij->a; 4511 } 4512 PetscCall(MatSeqAIJInvalidateDiagonal(A)); 4513 PetscCall(PetscObjectStateIncrease((PetscObject)A)); 4514 PetscFunctionReturn(PETSC_SUCCESS); 4515 } 4516 4517 /*@C 4518 MatSeqAIJRestoreArrayWrite - restore the read-only access array obtained from MatSeqAIJGetArrayRead 4519 4520 Not Collective; No Fortran Support 4521 4522 Input Parameter: 4523 . A - a MATSEQAIJ matrix 4524 4525 Output Parameter: 4526 . array - pointer to the data 4527 4528 Level: intermediate 4529 4530 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()` 4531 @*/ 4532 PetscErrorCode MatSeqAIJRestoreArrayWrite(Mat A, PetscScalar **array) 4533 { 4534 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 4535 4536 PetscFunctionBegin; 4537 if (aij->ops->restorearraywrite) { 4538 PetscCall((*aij->ops->restorearraywrite)(A, array)); 4539 } else { 4540 *array = NULL; 4541 } 4542 PetscFunctionReturn(PETSC_SUCCESS); 4543 } 4544 4545 /*@C 4546 MatSeqAIJGetCSRAndMemType - Get the CSR arrays and the memory type of the `MATSEQAIJ` matrix 4547 4548 Not Collective; No Fortran Support 4549 4550 Input Parameter: 4551 . mat - a matrix of type `MATSEQAIJ` or its subclasses 4552 4553 Output Parameters: 4554 + i - row map array of the matrix 4555 . j - column index array of the matrix 4556 . a - data array of the matrix 4557 - mtype - memory type of the arrays 4558 4559 Level: developer 4560 4561 Notes: 4562 Any of the output parameters can be `NULL`, in which case the corresponding value is not returned. 4563 If mat is a device matrix, the arrays are on the device. Otherwise, they are on the host. 4564 4565 One can call this routine on a preallocated but not assembled matrix to just get the memory of the CSR underneath the matrix. 4566 If the matrix is assembled, the data array `a` is guaranteed to have the latest values of the matrix. 4567 4568 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()` 4569 @*/ 4570 PetscErrorCode MatSeqAIJGetCSRAndMemType(Mat mat, const PetscInt **i, const PetscInt **j, PetscScalar **a, PetscMemType *mtype) 4571 { 4572 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 4573 4574 PetscFunctionBegin; 4575 PetscCheck(mat->preallocated, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "matrix is not preallocated"); 4576 if (aij->ops->getcsrandmemtype) { 4577 PetscCall((*aij->ops->getcsrandmemtype)(mat, i, j, a, mtype)); 4578 } else { 4579 if (i) *i = aij->i; 4580 if (j) *j = aij->j; 4581 if (a) *a = aij->a; 4582 if (mtype) *mtype = PETSC_MEMTYPE_HOST; 4583 } 4584 PetscFunctionReturn(PETSC_SUCCESS); 4585 } 4586 4587 /*@C 4588 MatSeqAIJGetMaxRowNonzeros - returns the maximum number of nonzeros in any row 4589 4590 Not Collective 4591 4592 Input Parameter: 4593 . A - a `MATSEQAIJ` matrix 4594 4595 Output Parameter: 4596 . nz - the maximum number of nonzeros in any row 4597 4598 Level: intermediate 4599 4600 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()` 4601 @*/ 4602 PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat A, PetscInt *nz) 4603 { 4604 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 4605 4606 PetscFunctionBegin; 4607 *nz = aij->rmax; 4608 PetscFunctionReturn(PETSC_SUCCESS); 4609 } 4610 4611 static PetscErrorCode MatCOOStructDestroy_SeqAIJ(void *data) 4612 { 4613 MatCOOStruct_SeqAIJ *coo = (MatCOOStruct_SeqAIJ *)data; 4614 PetscFunctionBegin; 4615 PetscCall(PetscFree(coo->perm)); 4616 PetscCall(PetscFree(coo->jmap)); 4617 PetscCall(PetscFree(coo)); 4618 PetscFunctionReturn(PETSC_SUCCESS); 4619 } 4620 4621 PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat mat, PetscCount coo_n, PetscInt coo_i[], PetscInt coo_j[]) 4622 { 4623 MPI_Comm comm; 4624 PetscInt *i, *j; 4625 PetscInt M, N, row; 4626 PetscCount k, p, q, nneg, nnz, start, end; /* Index the coo array, so use PetscCount as their type */ 4627 PetscInt *Ai; /* Change to PetscCount once we use it for row pointers */ 4628 PetscInt *Aj; 4629 PetscScalar *Aa; 4630 Mat_SeqAIJ *seqaij = (Mat_SeqAIJ *)(mat->data); 4631 MatType rtype; 4632 PetscCount *perm, *jmap; 4633 PetscContainer container; 4634 MatCOOStruct_SeqAIJ *coo; 4635 4636 PetscFunctionBegin; 4637 PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 4638 PetscCall(MatGetSize(mat, &M, &N)); 4639 i = coo_i; 4640 j = coo_j; 4641 PetscCall(PetscMalloc1(coo_n, &perm)); 4642 for (k = 0; k < coo_n; k++) { /* Ignore entries with negative row or col indices */ 4643 if (j[k] < 0) i[k] = -1; 4644 perm[k] = k; 4645 } 4646 4647 /* Sort by row */ 4648 PetscCall(PetscSortIntWithIntCountArrayPair(coo_n, i, j, perm)); 4649 4650 /* Advance k to the first row with a non-negative index */ 4651 for (k = 0; k < coo_n; k++) 4652 if (i[k] >= 0) break; 4653 nneg = k; 4654 PetscCall(PetscMalloc1(coo_n - nneg + 1, &jmap)); /* +1 to make a CSR-like data structure. jmap[i] originally is the number of repeats for i-th nonzero */ 4655 nnz = 0; /* Total number of unique nonzeros to be counted */ 4656 jmap++; /* Inc jmap by 1 for convenience */ 4657 4658 PetscCall(PetscCalloc1(M + 1, &Ai)); /* CSR of A */ 4659 PetscCall(PetscMalloc1(coo_n - nneg, &Aj)); /* We have at most coo_n-nneg unique nonzeros */ 4660 4661 /* Support for HYPRE */ 4662 PetscBool hypre; 4663 const char *name; 4664 PetscCall(PetscObjectGetName((PetscObject)mat, &name)); 4665 PetscCall(PetscStrcmp("_internal_COO_mat_for_hypre", name, &hypre)); 4666 4667 /* In each row, sort by column, then unique column indices to get row length */ 4668 Ai++; /* Inc by 1 for convenience */ 4669 q = 0; /* q-th unique nonzero, with q starting from 0 */ 4670 while (k < coo_n) { 4671 row = i[k]; 4672 start = k; /* [start,end) indices for this row */ 4673 while (k < coo_n && i[k] == row) k++; 4674 end = k; 4675 /* hack for HYPRE: swap min column to diag so that diagonal values will go first */ 4676 if (hypre) { 4677 PetscInt minj = PETSC_MAX_INT; 4678 PetscBool hasdiag = PETSC_FALSE; 4679 for (p = start; p < end; p++) { 4680 hasdiag = (PetscBool)(hasdiag || (j[p] == row)); 4681 minj = PetscMin(minj, j[p]); 4682 } 4683 if (hasdiag) { 4684 for (p = start; p < end; p++) { 4685 if (j[p] == minj) j[p] = row; 4686 else if (j[p] == row) j[p] = minj; 4687 } 4688 } 4689 } 4690 PetscCall(PetscSortIntWithCountArray(end - start, j + start, perm + start)); 4691 4692 /* Find number of unique col entries in this row */ 4693 Aj[q] = j[start]; /* Log the first nonzero in this row */ 4694 jmap[q] = 1; /* Number of repeats of this nonzero entry */ 4695 Ai[row] = 1; 4696 nnz++; 4697 4698 for (p = start + 1; p < end; p++) { /* Scan remaining nonzero in this row */ 4699 if (j[p] != j[p - 1]) { /* Meet a new nonzero */ 4700 q++; 4701 jmap[q] = 1; 4702 Aj[q] = j[p]; 4703 Ai[row]++; 4704 nnz++; 4705 } else { 4706 jmap[q]++; 4707 } 4708 } 4709 q++; /* Move to next row and thus next unique nonzero */ 4710 } 4711 Ai--; /* Back to the beginning of Ai[] */ 4712 for (k = 0; k < M; k++) Ai[k + 1] += Ai[k]; 4713 jmap--; /* Back to the beginning of jmap[] */ 4714 jmap[0] = 0; 4715 for (k = 0; k < nnz; k++) jmap[k + 1] += jmap[k]; 4716 if (nnz < coo_n - nneg) { /* Realloc with actual number of unique nonzeros */ 4717 PetscCount *jmap_new; 4718 PetscInt *Aj_new; 4719 4720 PetscCall(PetscMalloc1(nnz + 1, &jmap_new)); 4721 PetscCall(PetscArraycpy(jmap_new, jmap, nnz + 1)); 4722 PetscCall(PetscFree(jmap)); 4723 jmap = jmap_new; 4724 4725 PetscCall(PetscMalloc1(nnz, &Aj_new)); 4726 PetscCall(PetscArraycpy(Aj_new, Aj, nnz)); 4727 PetscCall(PetscFree(Aj)); 4728 Aj = Aj_new; 4729 } 4730 4731 if (nneg) { /* Discard heading entries with negative indices in perm[], as we'll access it from index 0 in MatSetValuesCOO */ 4732 PetscCount *perm_new; 4733 4734 PetscCall(PetscMalloc1(coo_n - nneg, &perm_new)); 4735 PetscCall(PetscArraycpy(perm_new, perm + nneg, coo_n - nneg)); 4736 PetscCall(PetscFree(perm)); 4737 perm = perm_new; 4738 } 4739 4740 PetscCall(MatGetRootType_Private(mat, &rtype)); 4741 PetscCall(PetscCalloc1(nnz, &Aa)); /* Zero the matrix */ 4742 PetscCall(MatSetSeqAIJWithArrays_private(PETSC_COMM_SELF, M, N, Ai, Aj, Aa, rtype, mat)); 4743 4744 seqaij->singlemalloc = PETSC_FALSE; /* Ai, Aj and Aa are not allocated in one big malloc */ 4745 seqaij->free_a = seqaij->free_ij = PETSC_TRUE; /* Let newmat own Ai, Aj and Aa */ 4746 4747 // Put the COO struct in a container and then attach that to the matrix 4748 PetscCall(PetscMalloc1(1, &coo)); 4749 coo->nz = nnz; 4750 coo->n = coo_n; 4751 coo->Atot = coo_n - nneg; // Annz is seqaij->nz, so no need to record that again 4752 coo->jmap = jmap; // of length nnz+1 4753 coo->perm = perm; 4754 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 4755 PetscCall(PetscContainerSetPointer(container, coo)); 4756 PetscCall(PetscContainerSetUserDestroy(container, MatCOOStructDestroy_SeqAIJ)); 4757 PetscCall(PetscObjectCompose((PetscObject)mat, "__PETSc_MatCOOStruct_Host", (PetscObject)container)); 4758 PetscCall(PetscContainerDestroy(&container)); 4759 PetscFunctionReturn(PETSC_SUCCESS); 4760 } 4761 4762 static PetscErrorCode MatSetValuesCOO_SeqAIJ(Mat A, const PetscScalar v[], InsertMode imode) 4763 { 4764 Mat_SeqAIJ *aseq = (Mat_SeqAIJ *)A->data; 4765 PetscCount i, j, Annz = aseq->nz; 4766 PetscCount *perm, *jmap; 4767 PetscScalar *Aa; 4768 PetscContainer container; 4769 MatCOOStruct_SeqAIJ *coo; 4770 4771 PetscFunctionBegin; 4772 PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_MatCOOStruct_Host", (PetscObject *)&container)); 4773 PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not found MatCOOStruct on this matrix"); 4774 PetscCall(PetscContainerGetPointer(container, (void **)&coo)); 4775 perm = coo->perm; 4776 jmap = coo->jmap; 4777 PetscCall(MatSeqAIJGetArray(A, &Aa)); 4778 for (i = 0; i < Annz; i++) { 4779 PetscScalar sum = 0.0; 4780 for (j = jmap[i]; j < jmap[i + 1]; j++) sum += v[perm[j]]; 4781 Aa[i] = (imode == INSERT_VALUES ? 0.0 : Aa[i]) + sum; 4782 } 4783 PetscCall(MatSeqAIJRestoreArray(A, &Aa)); 4784 PetscFunctionReturn(PETSC_SUCCESS); 4785 } 4786 4787 #if defined(PETSC_HAVE_CUDA) 4788 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat, MatType, MatReuse, Mat *); 4789 #endif 4790 #if defined(PETSC_HAVE_HIP) 4791 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJHIPSPARSE(Mat, MatType, MatReuse, Mat *); 4792 #endif 4793 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 4794 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJKokkos(Mat, MatType, MatReuse, Mat *); 4795 #endif 4796 4797 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat B) 4798 { 4799 Mat_SeqAIJ *b; 4800 PetscMPIInt size; 4801 4802 PetscFunctionBegin; 4803 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size)); 4804 PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Comm must be of size 1"); 4805 4806 PetscCall(PetscNew(&b)); 4807 4808 B->data = (void *)b; 4809 B->ops[0] = MatOps_Values; 4810 if (B->sortedfull) B->ops->setvalues = MatSetValues_SeqAIJ_SortedFull; 4811 4812 b->row = NULL; 4813 b->col = NULL; 4814 b->icol = NULL; 4815 b->reallocs = 0; 4816 b->ignorezeroentries = PETSC_FALSE; 4817 b->roworiented = PETSC_TRUE; 4818 b->nonew = 0; 4819 b->diag = NULL; 4820 b->solve_work = NULL; 4821 B->spptr = NULL; 4822 b->saved_values = NULL; 4823 b->idiag = NULL; 4824 b->mdiag = NULL; 4825 b->ssor_work = NULL; 4826 b->omega = 1.0; 4827 b->fshift = 0.0; 4828 b->idiagvalid = PETSC_FALSE; 4829 b->ibdiagvalid = PETSC_FALSE; 4830 b->keepnonzeropattern = PETSC_FALSE; 4831 4832 PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ)); 4833 #if defined(PETSC_HAVE_MATLAB) 4834 PetscCall(PetscObjectComposeFunction((PetscObject)B, "PetscMatlabEnginePut_C", MatlabEnginePut_SeqAIJ)); 4835 PetscCall(PetscObjectComposeFunction((PetscObject)B, "PetscMatlabEngineGet_C", MatlabEngineGet_SeqAIJ)); 4836 #endif 4837 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetColumnIndices_C", MatSeqAIJSetColumnIndices_SeqAIJ)); 4838 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqAIJ)); 4839 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqAIJ)); 4840 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqsbaij_C", MatConvert_SeqAIJ_SeqSBAIJ)); 4841 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqbaij_C", MatConvert_SeqAIJ_SeqBAIJ)); 4842 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijperm_C", MatConvert_SeqAIJ_SeqAIJPERM)); 4843 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijsell_C", MatConvert_SeqAIJ_SeqAIJSELL)); 4844 #if defined(PETSC_HAVE_MKL_SPARSE) 4845 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijmkl_C", MatConvert_SeqAIJ_SeqAIJMKL)); 4846 #endif 4847 #if defined(PETSC_HAVE_CUDA) 4848 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijcusparse_C", MatConvert_SeqAIJ_SeqAIJCUSPARSE)); 4849 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaijcusparse_seqaij_C", MatProductSetFromOptions_SeqAIJ)); 4850 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaijcusparse_C", MatProductSetFromOptions_SeqAIJ)); 4851 #endif 4852 #if defined(PETSC_HAVE_HIP) 4853 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijhipsparse_C", MatConvert_SeqAIJ_SeqAIJHIPSPARSE)); 4854 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaijhipsparse_seqaij_C", MatProductSetFromOptions_SeqAIJ)); 4855 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaijhipsparse_C", MatProductSetFromOptions_SeqAIJ)); 4856 #endif 4857 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 4858 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijkokkos_C", MatConvert_SeqAIJ_SeqAIJKokkos)); 4859 #endif 4860 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijcrl_C", MatConvert_SeqAIJ_SeqAIJCRL)); 4861 #if defined(PETSC_HAVE_ELEMENTAL) 4862 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_elemental_C", MatConvert_SeqAIJ_Elemental)); 4863 #endif 4864 #if defined(PETSC_HAVE_SCALAPACK) 4865 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_scalapack_C", MatConvert_AIJ_ScaLAPACK)); 4866 #endif 4867 #if defined(PETSC_HAVE_HYPRE) 4868 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_hypre_C", MatConvert_AIJ_HYPRE)); 4869 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_transpose_seqaij_seqaij_C", MatProductSetFromOptions_Transpose_AIJ_AIJ)); 4870 #endif 4871 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqdense_C", MatConvert_SeqAIJ_SeqDense)); 4872 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqsell_C", MatConvert_SeqAIJ_SeqSELL)); 4873 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_is_C", MatConvert_XAIJ_IS)); 4874 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsTranspose_C", MatIsTranspose_SeqAIJ)); 4875 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsHermitianTranspose_C", MatIsTranspose_SeqAIJ)); 4876 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetPreallocation_C", MatSeqAIJSetPreallocation_SeqAIJ)); 4877 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatResetPreallocation_C", MatResetPreallocation_SeqAIJ)); 4878 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetPreallocationCSR_C", MatSeqAIJSetPreallocationCSR_SeqAIJ)); 4879 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatReorderForNonzeroDiagonal_C", MatReorderForNonzeroDiagonal_SeqAIJ)); 4880 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_is_seqaij_C", MatProductSetFromOptions_IS_XAIJ)); 4881 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqdense_seqaij_C", MatProductSetFromOptions_SeqDense_SeqAIJ)); 4882 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaij_C", MatProductSetFromOptions_SeqAIJ)); 4883 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJKron_C", MatSeqAIJKron_SeqAIJ)); 4884 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSetPreallocationCOO_C", MatSetPreallocationCOO_SeqAIJ)); 4885 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSetValuesCOO_C", MatSetValuesCOO_SeqAIJ)); 4886 PetscCall(MatCreate_SeqAIJ_Inode(B)); 4887 PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ)); 4888 PetscCall(MatSeqAIJSetTypeFromOptions(B)); /* this allows changing the matrix subtype to say MATSEQAIJPERM */ 4889 PetscFunctionReturn(PETSC_SUCCESS); 4890 } 4891 4892 /* 4893 Given a matrix generated with MatGetFactor() duplicates all the information in A into C 4894 */ 4895 PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat C, Mat A, MatDuplicateOption cpvalues, PetscBool mallocmatspace) 4896 { 4897 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data, *a = (Mat_SeqAIJ *)A->data; 4898 PetscInt m = A->rmap->n, i; 4899 4900 PetscFunctionBegin; 4901 PetscCheck(A->assembled || cpvalues == MAT_DO_NOT_COPY_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot duplicate unassembled matrix"); 4902 4903 C->factortype = A->factortype; 4904 c->row = NULL; 4905 c->col = NULL; 4906 c->icol = NULL; 4907 c->reallocs = 0; 4908 4909 C->assembled = A->assembled; 4910 4911 if (A->preallocated) { 4912 PetscCall(PetscLayoutReference(A->rmap, &C->rmap)); 4913 PetscCall(PetscLayoutReference(A->cmap, &C->cmap)); 4914 4915 if (!A->hash_active) { 4916 PetscCall(PetscMalloc1(m, &c->imax)); 4917 PetscCall(PetscMemcpy(c->imax, a->imax, m * sizeof(PetscInt))); 4918 PetscCall(PetscMalloc1(m, &c->ilen)); 4919 PetscCall(PetscMemcpy(c->ilen, a->ilen, m * sizeof(PetscInt))); 4920 4921 /* allocate the matrix space */ 4922 if (mallocmatspace) { 4923 PetscCall(PetscMalloc3(a->i[m], &c->a, a->i[m], &c->j, m + 1, &c->i)); 4924 4925 c->singlemalloc = PETSC_TRUE; 4926 4927 PetscCall(PetscArraycpy(c->i, a->i, m + 1)); 4928 if (m > 0) { 4929 PetscCall(PetscArraycpy(c->j, a->j, a->i[m])); 4930 if (cpvalues == MAT_COPY_VALUES) { 4931 const PetscScalar *aa; 4932 4933 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 4934 PetscCall(PetscArraycpy(c->a, aa, a->i[m])); 4935 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 4936 } else { 4937 PetscCall(PetscArrayzero(c->a, a->i[m])); 4938 } 4939 } 4940 } 4941 C->preallocated = PETSC_TRUE; 4942 } else { 4943 PetscCheck(mallocmatspace, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Cannot malloc matrix memory from a non-preallocated matrix"); 4944 PetscCall(MatSetUp(C)); 4945 } 4946 4947 c->ignorezeroentries = a->ignorezeroentries; 4948 c->roworiented = a->roworiented; 4949 c->nonew = a->nonew; 4950 if (a->diag) { 4951 PetscCall(PetscMalloc1(m + 1, &c->diag)); 4952 PetscCall(PetscMemcpy(c->diag, a->diag, m * sizeof(PetscInt))); 4953 } else c->diag = NULL; 4954 4955 c->solve_work = NULL; 4956 c->saved_values = NULL; 4957 c->idiag = NULL; 4958 c->ssor_work = NULL; 4959 c->keepnonzeropattern = a->keepnonzeropattern; 4960 c->free_a = PETSC_TRUE; 4961 c->free_ij = PETSC_TRUE; 4962 4963 c->rmax = a->rmax; 4964 c->nz = a->nz; 4965 c->maxnz = a->nz; /* Since we allocate exactly the right amount */ 4966 4967 c->compressedrow.use = a->compressedrow.use; 4968 c->compressedrow.nrows = a->compressedrow.nrows; 4969 if (a->compressedrow.use) { 4970 i = a->compressedrow.nrows; 4971 PetscCall(PetscMalloc2(i + 1, &c->compressedrow.i, i, &c->compressedrow.rindex)); 4972 PetscCall(PetscArraycpy(c->compressedrow.i, a->compressedrow.i, i + 1)); 4973 PetscCall(PetscArraycpy(c->compressedrow.rindex, a->compressedrow.rindex, i)); 4974 } else { 4975 c->compressedrow.use = PETSC_FALSE; 4976 c->compressedrow.i = NULL; 4977 c->compressedrow.rindex = NULL; 4978 } 4979 c->nonzerorowcnt = a->nonzerorowcnt; 4980 C->nonzerostate = A->nonzerostate; 4981 4982 PetscCall(MatDuplicate_SeqAIJ_Inode(A, cpvalues, &C)); 4983 } 4984 PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist)); 4985 PetscFunctionReturn(PETSC_SUCCESS); 4986 } 4987 4988 PetscErrorCode MatDuplicate_SeqAIJ(Mat A, MatDuplicateOption cpvalues, Mat *B) 4989 { 4990 PetscFunctionBegin; 4991 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 4992 PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n)); 4993 if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) PetscCall(MatSetBlockSizesFromMats(*B, A, A)); 4994 PetscCall(MatSetType(*B, ((PetscObject)A)->type_name)); 4995 PetscCall(MatDuplicateNoCreate_SeqAIJ(*B, A, cpvalues, PETSC_TRUE)); 4996 PetscFunctionReturn(PETSC_SUCCESS); 4997 } 4998 4999 PetscErrorCode MatLoad_SeqAIJ(Mat newMat, PetscViewer viewer) 5000 { 5001 PetscBool isbinary, ishdf5; 5002 5003 PetscFunctionBegin; 5004 PetscValidHeaderSpecific(newMat, MAT_CLASSID, 1); 5005 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 5006 /* force binary viewer to load .info file if it has not yet done so */ 5007 PetscCall(PetscViewerSetUp(viewer)); 5008 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 5009 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 5010 if (isbinary) { 5011 PetscCall(MatLoad_SeqAIJ_Binary(newMat, viewer)); 5012 } else if (ishdf5) { 5013 #if defined(PETSC_HAVE_HDF5) 5014 PetscCall(MatLoad_AIJ_HDF5(newMat, viewer)); 5015 #else 5016 SETERRQ(PetscObjectComm((PetscObject)newMat), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5"); 5017 #endif 5018 } else { 5019 SETERRQ(PetscObjectComm((PetscObject)newMat), PETSC_ERR_SUP, "Viewer type %s not yet supported for reading %s matrices", ((PetscObject)viewer)->type_name, ((PetscObject)newMat)->type_name); 5020 } 5021 PetscFunctionReturn(PETSC_SUCCESS); 5022 } 5023 5024 PetscErrorCode MatLoad_SeqAIJ_Binary(Mat mat, PetscViewer viewer) 5025 { 5026 Mat_SeqAIJ *a = (Mat_SeqAIJ *)mat->data; 5027 PetscInt header[4], *rowlens, M, N, nz, sum, rows, cols, i; 5028 5029 PetscFunctionBegin; 5030 PetscCall(PetscViewerSetUp(viewer)); 5031 5032 /* read in matrix header */ 5033 PetscCall(PetscViewerBinaryRead(viewer, header, 4, NULL, PETSC_INT)); 5034 PetscCheck(header[0] == MAT_FILE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a matrix object in file"); 5035 M = header[1]; 5036 N = header[2]; 5037 nz = header[3]; 5038 PetscCheck(M >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix row size (%" PetscInt_FMT ") in file is negative", M); 5039 PetscCheck(N >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix column size (%" PetscInt_FMT ") in file is negative", N); 5040 PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix stored in special format on disk, cannot load as SeqAIJ"); 5041 5042 /* set block sizes from the viewer's .info file */ 5043 PetscCall(MatLoad_Binary_BlockSizes(mat, viewer)); 5044 /* set local and global sizes if not set already */ 5045 if (mat->rmap->n < 0) mat->rmap->n = M; 5046 if (mat->cmap->n < 0) mat->cmap->n = N; 5047 if (mat->rmap->N < 0) mat->rmap->N = M; 5048 if (mat->cmap->N < 0) mat->cmap->N = N; 5049 PetscCall(PetscLayoutSetUp(mat->rmap)); 5050 PetscCall(PetscLayoutSetUp(mat->cmap)); 5051 5052 /* check if the matrix sizes are correct */ 5053 PetscCall(MatGetSize(mat, &rows, &cols)); 5054 PetscCheck(M == rows && N == cols, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different sizes (%" PetscInt_FMT ", %" PetscInt_FMT ") than the input matrix (%" PetscInt_FMT ", %" PetscInt_FMT ")", M, N, rows, cols); 5055 5056 /* read in row lengths */ 5057 PetscCall(PetscMalloc1(M, &rowlens)); 5058 PetscCall(PetscViewerBinaryRead(viewer, rowlens, M, NULL, PETSC_INT)); 5059 /* check if sum(rowlens) is same as nz */ 5060 sum = 0; 5061 for (i = 0; i < M; i++) sum += rowlens[i]; 5062 PetscCheck(sum == nz, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Inconsistent matrix data in file: nonzeros = %" PetscInt_FMT ", sum-row-lengths = %" PetscInt_FMT, nz, sum); 5063 /* preallocate and check sizes */ 5064 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(mat, 0, rowlens)); 5065 PetscCall(MatGetSize(mat, &rows, &cols)); 5066 PetscCheck(M == rows && N == cols, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%" PetscInt_FMT ", %" PetscInt_FMT ") than the input matrix (%" PetscInt_FMT ", %" PetscInt_FMT ")", M, N, rows, cols); 5067 /* store row lengths */ 5068 PetscCall(PetscArraycpy(a->ilen, rowlens, M)); 5069 PetscCall(PetscFree(rowlens)); 5070 5071 /* fill in "i" row pointers */ 5072 a->i[0] = 0; 5073 for (i = 0; i < M; i++) a->i[i + 1] = a->i[i] + a->ilen[i]; 5074 /* read in "j" column indices */ 5075 PetscCall(PetscViewerBinaryRead(viewer, a->j, nz, NULL, PETSC_INT)); 5076 /* read in "a" nonzero values */ 5077 PetscCall(PetscViewerBinaryRead(viewer, a->a, nz, NULL, PETSC_SCALAR)); 5078 5079 PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 5080 PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 5081 PetscFunctionReturn(PETSC_SUCCESS); 5082 } 5083 5084 PetscErrorCode MatEqual_SeqAIJ(Mat A, Mat B, PetscBool *flg) 5085 { 5086 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 5087 const PetscScalar *aa, *ba; 5088 #if defined(PETSC_USE_COMPLEX) 5089 PetscInt k; 5090 #endif 5091 5092 PetscFunctionBegin; 5093 /* If the matrix dimensions are not equal,or no of nonzeros */ 5094 if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) || (a->nz != b->nz)) { 5095 *flg = PETSC_FALSE; 5096 PetscFunctionReturn(PETSC_SUCCESS); 5097 } 5098 5099 /* if the a->i are the same */ 5100 PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, flg)); 5101 if (!*flg) PetscFunctionReturn(PETSC_SUCCESS); 5102 5103 /* if a->j are the same */ 5104 PetscCall(PetscArraycmp(a->j, b->j, a->nz, flg)); 5105 if (!*flg) PetscFunctionReturn(PETSC_SUCCESS); 5106 5107 PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 5108 PetscCall(MatSeqAIJGetArrayRead(B, &ba)); 5109 /* if a->a are the same */ 5110 #if defined(PETSC_USE_COMPLEX) 5111 for (k = 0; k < a->nz; k++) { 5112 if (PetscRealPart(aa[k]) != PetscRealPart(ba[k]) || PetscImaginaryPart(aa[k]) != PetscImaginaryPart(ba[k])) { 5113 *flg = PETSC_FALSE; 5114 PetscFunctionReturn(PETSC_SUCCESS); 5115 } 5116 } 5117 #else 5118 PetscCall(PetscArraycmp(aa, ba, a->nz, flg)); 5119 #endif 5120 PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 5121 PetscCall(MatSeqAIJRestoreArrayRead(B, &ba)); 5122 PetscFunctionReturn(PETSC_SUCCESS); 5123 } 5124 5125 /*@ 5126 MatCreateSeqAIJWithArrays - Creates an sequential `MATSEQAIJ` matrix using matrix elements (in CSR format) 5127 provided by the user. 5128 5129 Collective 5130 5131 Input Parameters: 5132 + comm - must be an MPI communicator of size 1 5133 . m - number of rows 5134 . n - number of columns 5135 . i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix 5136 . j - column indices 5137 - a - matrix values 5138 5139 Output Parameter: 5140 . mat - the matrix 5141 5142 Level: intermediate 5143 5144 Notes: 5145 The `i`, `j`, and `a` arrays are not copied by this routine, the user must free these arrays 5146 once the matrix is destroyed and not before 5147 5148 You cannot set new nonzero locations into this matrix, that will generate an error. 5149 5150 The `i` and `j` indices are 0 based 5151 5152 The format which is used for the sparse matrix input, is equivalent to a 5153 row-major ordering.. i.e for the following matrix, the input data expected is 5154 as shown 5155 .vb 5156 1 0 0 5157 2 0 3 5158 4 5 6 5159 5160 i = {0,1,3,6} [size = nrow+1 = 3+1] 5161 j = {0,0,2,0,1,2} [size = 6]; values must be sorted for each row 5162 v = {1,2,3,4,5,6} [size = 6] 5163 .ve 5164 5165 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MatCreateMPIAIJWithArrays()`, `MatMPIAIJSetPreallocationCSR()` 5166 @*/ 5167 PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat) 5168 { 5169 PetscInt ii; 5170 Mat_SeqAIJ *aij; 5171 PetscInt jj; 5172 5173 PetscFunctionBegin; 5174 PetscCheck(m <= 0 || i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i (row indices) must start with 0"); 5175 PetscCall(MatCreate(comm, mat)); 5176 PetscCall(MatSetSizes(*mat, m, n, m, n)); 5177 /* PetscCall(MatSetBlockSizes(*mat,,)); */ 5178 PetscCall(MatSetType(*mat, MATSEQAIJ)); 5179 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, MAT_SKIP_ALLOCATION, NULL)); 5180 aij = (Mat_SeqAIJ *)(*mat)->data; 5181 PetscCall(PetscMalloc1(m, &aij->imax)); 5182 PetscCall(PetscMalloc1(m, &aij->ilen)); 5183 5184 aij->i = i; 5185 aij->j = j; 5186 aij->a = a; 5187 aij->singlemalloc = PETSC_FALSE; 5188 aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 5189 aij->free_a = PETSC_FALSE; 5190 aij->free_ij = PETSC_FALSE; 5191 5192 for (ii = 0, aij->nonzerorowcnt = 0, aij->rmax = 0; ii < m; ii++) { 5193 aij->ilen[ii] = aij->imax[ii] = i[ii + 1] - i[ii]; 5194 if (PetscDefined(USE_DEBUG)) { 5195 PetscCheck(i[ii + 1] - i[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative row length in i (row indices) row = %" PetscInt_FMT " length = %" PetscInt_FMT, ii, i[ii + 1] - i[ii]); 5196 for (jj = i[ii] + 1; jj < i[ii + 1]; jj++) { 5197 PetscCheck(j[jj] >= j[jj - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column entry number %" PetscInt_FMT " (actual column %" PetscInt_FMT ") in row %" PetscInt_FMT " is not sorted", jj - i[ii], j[jj], ii); 5198 PetscCheck(j[jj] != j[jj - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column entry number %" PetscInt_FMT " (actual column %" PetscInt_FMT ") in row %" PetscInt_FMT " is identical to previous entry", jj - i[ii], j[jj], ii); 5199 } 5200 } 5201 } 5202 if (PetscDefined(USE_DEBUG)) { 5203 for (ii = 0; ii < aij->i[m]; ii++) { 5204 PetscCheck(j[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative column index at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]); 5205 PetscCheck(j[ii] <= n - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column index to large at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]); 5206 } 5207 } 5208 5209 PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY)); 5210 PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY)); 5211 PetscFunctionReturn(PETSC_SUCCESS); 5212 } 5213 5214 /*@ 5215 MatCreateSeqAIJFromTriple - Creates an sequential `MATSEQAIJ` matrix using matrix elements (in COO format) 5216 provided by the user. 5217 5218 Collective 5219 5220 Input Parameters: 5221 + comm - must be an MPI communicator of size 1 5222 . m - number of rows 5223 . n - number of columns 5224 . i - row indices 5225 . j - column indices 5226 . a - matrix values 5227 . nz - number of nonzeros 5228 - idx - if the `i` and `j` indices start with 1 use `PETSC_TRUE` otherwise use `PETSC_FALSE` 5229 5230 Output Parameter: 5231 . mat - the matrix 5232 5233 Level: intermediate 5234 5235 Example: 5236 For the following matrix, the input data expected is as shown (using 0 based indexing) 5237 .vb 5238 1 0 0 5239 2 0 3 5240 4 5 6 5241 5242 i = {0,1,1,2,2,2} 5243 j = {0,0,2,0,1,2} 5244 v = {1,2,3,4,5,6} 5245 .ve 5246 5247 Note: 5248 Instead of using this function, users should also consider `MatSetPreallocationCOO()` and `MatSetValuesCOO()`, which allow repeated or remote entries, 5249 and are particularly useful in iterative applications. 5250 5251 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MatCreateSeqAIJWithArrays()`, `MatMPIAIJSetPreallocationCSR()`, `MatSetValuesCOO()`, `MatSetPreallocationCOO()` 5252 @*/ 5253 PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat, PetscInt nz, PetscBool idx) 5254 { 5255 PetscInt ii, *nnz, one = 1, row, col; 5256 5257 PetscFunctionBegin; 5258 PetscCall(PetscCalloc1(m, &nnz)); 5259 for (ii = 0; ii < nz; ii++) nnz[i[ii] - !!idx] += 1; 5260 PetscCall(MatCreate(comm, mat)); 5261 PetscCall(MatSetSizes(*mat, m, n, m, n)); 5262 PetscCall(MatSetType(*mat, MATSEQAIJ)); 5263 PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, 0, nnz)); 5264 for (ii = 0; ii < nz; ii++) { 5265 if (idx) { 5266 row = i[ii] - 1; 5267 col = j[ii] - 1; 5268 } else { 5269 row = i[ii]; 5270 col = j[ii]; 5271 } 5272 PetscCall(MatSetValues(*mat, one, &row, one, &col, &a[ii], ADD_VALUES)); 5273 } 5274 PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY)); 5275 PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY)); 5276 PetscCall(PetscFree(nnz)); 5277 PetscFunctionReturn(PETSC_SUCCESS); 5278 } 5279 5280 PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat A) 5281 { 5282 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 5283 5284 PetscFunctionBegin; 5285 a->idiagvalid = PETSC_FALSE; 5286 a->ibdiagvalid = PETSC_FALSE; 5287 5288 PetscCall(MatSeqAIJInvalidateDiagonal_Inode(A)); 5289 PetscFunctionReturn(PETSC_SUCCESS); 5290 } 5291 5292 PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat) 5293 { 5294 PetscFunctionBegin; 5295 PetscCall(MatCreateMPIMatConcatenateSeqMat_MPIAIJ(comm, inmat, n, scall, outmat)); 5296 PetscFunctionReturn(PETSC_SUCCESS); 5297 } 5298 5299 /* 5300 Permute A into C's *local* index space using rowemb,colemb. 5301 The embedding are supposed to be injections and the above implies that the range of rowemb is a subset 5302 of [0,m), colemb is in [0,n). 5303 If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A. 5304 */ 5305 PetscErrorCode MatSetSeqMat_SeqAIJ(Mat C, IS rowemb, IS colemb, MatStructure pattern, Mat B) 5306 { 5307 /* If making this function public, change the error returned in this function away from _PLIB. */ 5308 Mat_SeqAIJ *Baij; 5309 PetscBool seqaij; 5310 PetscInt m, n, *nz, i, j, count; 5311 PetscScalar v; 5312 const PetscInt *rowindices, *colindices; 5313 5314 PetscFunctionBegin; 5315 if (!B) PetscFunctionReturn(PETSC_SUCCESS); 5316 /* Check to make sure the target matrix (and embeddings) are compatible with C and each other. */ 5317 PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATSEQAIJ, &seqaij)); 5318 PetscCheck(seqaij, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is of wrong type"); 5319 if (rowemb) { 5320 PetscCall(ISGetLocalSize(rowemb, &m)); 5321 PetscCheck(m == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Row IS of size %" PetscInt_FMT " is incompatible with matrix row size %" PetscInt_FMT, m, B->rmap->n); 5322 } else { 5323 PetscCheck(C->rmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is row-incompatible with the target matrix"); 5324 } 5325 if (colemb) { 5326 PetscCall(ISGetLocalSize(colemb, &n)); 5327 PetscCheck(n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Diag col IS of size %" PetscInt_FMT " is incompatible with input matrix col size %" PetscInt_FMT, n, B->cmap->n); 5328 } else { 5329 PetscCheck(C->cmap->n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is col-incompatible with the target matrix"); 5330 } 5331 5332 Baij = (Mat_SeqAIJ *)(B->data); 5333 if (pattern == DIFFERENT_NONZERO_PATTERN) { 5334 PetscCall(PetscMalloc1(B->rmap->n, &nz)); 5335 for (i = 0; i < B->rmap->n; i++) nz[i] = Baij->i[i + 1] - Baij->i[i]; 5336 PetscCall(MatSeqAIJSetPreallocation(C, 0, nz)); 5337 PetscCall(PetscFree(nz)); 5338 } 5339 if (pattern == SUBSET_NONZERO_PATTERN) PetscCall(MatZeroEntries(C)); 5340 count = 0; 5341 rowindices = NULL; 5342 colindices = NULL; 5343 if (rowemb) PetscCall(ISGetIndices(rowemb, &rowindices)); 5344 if (colemb) PetscCall(ISGetIndices(colemb, &colindices)); 5345 for (i = 0; i < B->rmap->n; i++) { 5346 PetscInt row; 5347 row = i; 5348 if (rowindices) row = rowindices[i]; 5349 for (j = Baij->i[i]; j < Baij->i[i + 1]; j++) { 5350 PetscInt col; 5351 col = Baij->j[count]; 5352 if (colindices) col = colindices[col]; 5353 v = Baij->a[count]; 5354 PetscCall(MatSetValues(C, 1, &row, 1, &col, &v, INSERT_VALUES)); 5355 ++count; 5356 } 5357 } 5358 /* FIXME: set C's nonzerostate correctly. */ 5359 /* Assembly for C is necessary. */ 5360 C->preallocated = PETSC_TRUE; 5361 C->assembled = PETSC_TRUE; 5362 C->was_assembled = PETSC_FALSE; 5363 PetscFunctionReturn(PETSC_SUCCESS); 5364 } 5365 5366 PetscErrorCode MatEliminateZeros_SeqAIJ(Mat A, PetscBool keep) 5367 { 5368 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 5369 MatScalar *aa = a->a; 5370 PetscInt m = A->rmap->n, fshift = 0, fshift_prev = 0, i, k; 5371 PetscInt *ailen = a->ilen, *imax = a->imax, *ai = a->i, *aj = a->j, rmax = 0; 5372 5373 PetscFunctionBegin; 5374 PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot eliminate zeros for unassembled matrix"); 5375 if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 5376 for (i = 1; i <= m; i++) { 5377 /* move each nonzero entry back by the amount of zero slots (fshift) before it*/ 5378 for (k = ai[i - 1]; k < ai[i]; k++) { 5379 if (aa[k] == 0 && (aj[k] != i - 1 || !keep)) fshift++; 5380 else { 5381 if (aa[k] == 0 && aj[k] == i - 1) PetscCall(PetscInfo(A, "Keep the diagonal zero at row %" PetscInt_FMT "\n", i - 1)); 5382 aa[k - fshift] = aa[k]; 5383 aj[k - fshift] = aj[k]; 5384 } 5385 } 5386 ai[i - 1] -= fshift_prev; // safe to update ai[i-1] now since it will not be used in the next iteration 5387 fshift_prev = fshift; 5388 /* reset ilen and imax for each row */ 5389 ailen[i - 1] = imax[i - 1] = ai[i] - fshift - ai[i - 1]; 5390 a->nonzerorowcnt += ((ai[i] - fshift - ai[i - 1]) > 0); 5391 rmax = PetscMax(rmax, ailen[i - 1]); 5392 } 5393 if (fshift) { 5394 if (m) { 5395 ai[m] -= fshift; 5396 a->nz = ai[m]; 5397 } 5398 PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; zeros eliminated: %" PetscInt_FMT "; nonzeros left: %" PetscInt_FMT "\n", m, A->cmap->n, fshift, a->nz)); 5399 A->nonzerostate++; 5400 A->info.nz_unneeded += (PetscReal)fshift; 5401 a->rmax = rmax; 5402 if (a->inode.use && a->inode.checked) PetscCall(MatSeqAIJCheckInode(A)); 5403 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5404 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 5405 } 5406 PetscFunctionReturn(PETSC_SUCCESS); 5407 } 5408 5409 PetscFunctionList MatSeqAIJList = NULL; 5410 5411 /*@C 5412 MatSeqAIJSetType - Converts a `MATSEQAIJ` matrix to a subtype 5413 5414 Collective 5415 5416 Input Parameters: 5417 + mat - the matrix object 5418 - matype - matrix type 5419 5420 Options Database Key: 5421 . -mat_seqaij_type <method> - for example seqaijcrl 5422 5423 Level: intermediate 5424 5425 .seealso: [](ch_matrices), `Mat`, `PCSetType()`, `VecSetType()`, `MatCreate()`, `MatType` 5426 @*/ 5427 PetscErrorCode MatSeqAIJSetType(Mat mat, MatType matype) 5428 { 5429 PetscBool sametype; 5430 PetscErrorCode (*r)(Mat, MatType, MatReuse, Mat *); 5431 5432 PetscFunctionBegin; 5433 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5434 PetscCall(PetscObjectTypeCompare((PetscObject)mat, matype, &sametype)); 5435 if (sametype) PetscFunctionReturn(PETSC_SUCCESS); 5436 5437 PetscCall(PetscFunctionListFind(MatSeqAIJList, matype, &r)); 5438 PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Mat type given: %s", matype); 5439 PetscCall((*r)(mat, matype, MAT_INPLACE_MATRIX, &mat)); 5440 PetscFunctionReturn(PETSC_SUCCESS); 5441 } 5442 5443 /*@C 5444 MatSeqAIJRegister - - Adds a new sub-matrix type for sequential `MATSEQAIJ` matrices 5445 5446 Not Collective 5447 5448 Input Parameters: 5449 + sname - name of a new user-defined matrix type, for example `MATSEQAIJCRL` 5450 - function - routine to convert to subtype 5451 5452 Level: advanced 5453 5454 Notes: 5455 `MatSeqAIJRegister()` may be called multiple times to add several user-defined solvers. 5456 5457 Then, your matrix can be chosen with the procedural interface at runtime via the option 5458 $ -mat_seqaij_type my_mat 5459 5460 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRegisterAll()` 5461 @*/ 5462 PetscErrorCode MatSeqAIJRegister(const char sname[], PetscErrorCode (*function)(Mat, MatType, MatReuse, Mat *)) 5463 { 5464 PetscFunctionBegin; 5465 PetscCall(MatInitializePackage()); 5466 PetscCall(PetscFunctionListAdd(&MatSeqAIJList, sname, function)); 5467 PetscFunctionReturn(PETSC_SUCCESS); 5468 } 5469 5470 PetscBool MatSeqAIJRegisterAllCalled = PETSC_FALSE; 5471 5472 /*@C 5473 MatSeqAIJRegisterAll - Registers all of the matrix subtypes of `MATSSEQAIJ` 5474 5475 Not Collective 5476 5477 Level: advanced 5478 5479 Note: 5480 This registers the versions of `MATSEQAIJ` for GPUs 5481 5482 .seealso: [](ch_matrices), `Mat`, `MatRegisterAll()`, `MatSeqAIJRegister()` 5483 @*/ 5484 PetscErrorCode MatSeqAIJRegisterAll(void) 5485 { 5486 PetscFunctionBegin; 5487 if (MatSeqAIJRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 5488 MatSeqAIJRegisterAllCalled = PETSC_TRUE; 5489 5490 PetscCall(MatSeqAIJRegister(MATSEQAIJCRL, MatConvert_SeqAIJ_SeqAIJCRL)); 5491 PetscCall(MatSeqAIJRegister(MATSEQAIJPERM, MatConvert_SeqAIJ_SeqAIJPERM)); 5492 PetscCall(MatSeqAIJRegister(MATSEQAIJSELL, MatConvert_SeqAIJ_SeqAIJSELL)); 5493 #if defined(PETSC_HAVE_MKL_SPARSE) 5494 PetscCall(MatSeqAIJRegister(MATSEQAIJMKL, MatConvert_SeqAIJ_SeqAIJMKL)); 5495 #endif 5496 #if defined(PETSC_HAVE_CUDA) 5497 PetscCall(MatSeqAIJRegister(MATSEQAIJCUSPARSE, MatConvert_SeqAIJ_SeqAIJCUSPARSE)); 5498 #endif 5499 #if defined(PETSC_HAVE_HIP) 5500 PetscCall(MatSeqAIJRegister(MATSEQAIJHIPSPARSE, MatConvert_SeqAIJ_SeqAIJHIPSPARSE)); 5501 #endif 5502 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 5503 PetscCall(MatSeqAIJRegister(MATSEQAIJKOKKOS, MatConvert_SeqAIJ_SeqAIJKokkos)); 5504 #endif 5505 #if defined(PETSC_HAVE_VIENNACL) && defined(PETSC_HAVE_VIENNACL_NO_CUDA) 5506 PetscCall(MatSeqAIJRegister(MATMPIAIJVIENNACL, MatConvert_SeqAIJ_SeqAIJViennaCL)); 5507 #endif 5508 PetscFunctionReturn(PETSC_SUCCESS); 5509 } 5510 5511 /* 5512 Special version for direct calls from Fortran 5513 */ 5514 #include <petsc/private/fortranimpl.h> 5515 #if defined(PETSC_HAVE_FORTRAN_CAPS) 5516 #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ 5517 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 5518 #define matsetvaluesseqaij_ matsetvaluesseqaij 5519 #endif 5520 5521 /* Change these macros so can be used in void function */ 5522 5523 /* Change these macros so can be used in void function */ 5524 /* Identical to PetscCallVoid, except it assigns to *_ierr */ 5525 #undef PetscCall 5526 #define PetscCall(...) \ 5527 do { \ 5528 PetscErrorCode ierr_msv_mpiaij = __VA_ARGS__; \ 5529 if (PetscUnlikely(ierr_msv_mpiaij)) { \ 5530 *_ierr = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_msv_mpiaij, PETSC_ERROR_REPEAT, " "); \ 5531 return; \ 5532 } \ 5533 } while (0) 5534 5535 #undef SETERRQ 5536 #define SETERRQ(comm, ierr, ...) \ 5537 do { \ 5538 *_ierr = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 5539 return; \ 5540 } while (0) 5541 5542 PETSC_EXTERN void matsetvaluesseqaij_(Mat *AA, PetscInt *mm, const PetscInt im[], PetscInt *nn, const PetscInt in[], const PetscScalar v[], InsertMode *isis, PetscErrorCode *_ierr) 5543 { 5544 Mat A = *AA; 5545 PetscInt m = *mm, n = *nn; 5546 InsertMode is = *isis; 5547 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 5548 PetscInt *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N; 5549 PetscInt *imax, *ai, *ailen; 5550 PetscInt *aj, nonew = a->nonew, lastcol = -1; 5551 MatScalar *ap, value, *aa; 5552 PetscBool ignorezeroentries = a->ignorezeroentries; 5553 PetscBool roworiented = a->roworiented; 5554 5555 PetscFunctionBegin; 5556 MatCheckPreallocated(A, 1); 5557 imax = a->imax; 5558 ai = a->i; 5559 ailen = a->ilen; 5560 aj = a->j; 5561 aa = a->a; 5562 5563 for (k = 0; k < m; k++) { /* loop over added rows */ 5564 row = im[k]; 5565 if (row < 0) continue; 5566 PetscCheck(row < A->rmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large"); 5567 rp = aj + ai[row]; 5568 ap = aa + ai[row]; 5569 rmax = imax[row]; 5570 nrow = ailen[row]; 5571 low = 0; 5572 high = nrow; 5573 for (l = 0; l < n; l++) { /* loop over added columns */ 5574 if (in[l] < 0) continue; 5575 PetscCheck(in[l] < A->cmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Column too large"); 5576 col = in[l]; 5577 if (roworiented) value = v[l + k * n]; 5578 else value = v[k + l * m]; 5579 5580 if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue; 5581 5582 if (col <= lastcol) low = 0; 5583 else high = nrow; 5584 lastcol = col; 5585 while (high - low > 5) { 5586 t = (low + high) / 2; 5587 if (rp[t] > col) high = t; 5588 else low = t; 5589 } 5590 for (i = low; i < high; i++) { 5591 if (rp[i] > col) break; 5592 if (rp[i] == col) { 5593 if (is == ADD_VALUES) ap[i] += value; 5594 else ap[i] = value; 5595 goto noinsert; 5596 } 5597 } 5598 if (value == 0.0 && ignorezeroentries) goto noinsert; 5599 if (nonew == 1) goto noinsert; 5600 PetscCheck(nonew != -1, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero in the matrix"); 5601 MatSeqXAIJReallocateAIJ(A, A->rmap->n, 1, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar); 5602 N = nrow++ - 1; 5603 a->nz++; 5604 high++; 5605 /* shift up all the later entries in this row */ 5606 for (ii = N; ii >= i; ii--) { 5607 rp[ii + 1] = rp[ii]; 5608 ap[ii + 1] = ap[ii]; 5609 } 5610 rp[i] = col; 5611 ap[i] = value; 5612 A->nonzerostate++; 5613 noinsert:; 5614 low = i + 1; 5615 } 5616 ailen[row] = nrow; 5617 } 5618 PetscFunctionReturnVoid(); 5619 } 5620 /* Undefining these here since they were redefined from their original definition above! No 5621 * other PETSc functions should be defined past this point, as it is impossible to recover the 5622 * original definitions */ 5623 #undef PetscCall 5624 #undef SETERRQ 5625