1 2 /* 3 Defines the basic matrix operations for the SBAIJ (compressed row) 4 matrix storage format. 5 */ 6 #include <../src/mat/impls/baij/seq/baij.h> /*I "petscmat.h" I*/ 7 #include <../src/mat/impls/sbaij/seq/sbaij.h> 8 #include <petscblaslapack.h> 9 10 #include <../src/mat/impls/sbaij/seq/relax.h> 11 #define USESHORT 12 #include <../src/mat/impls/sbaij/seq/relax.h> 13 14 #if defined(PETSC_HAVE_ELEMENTAL) 15 PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat,MatType,MatReuse,Mat*); 16 #endif 17 18 /* 19 Checks for missing diagonals 20 */ 21 PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscBool *missing,PetscInt *dd) 22 { 23 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 24 PetscErrorCode ierr; 25 PetscInt *diag,*ii = a->i,i; 26 27 PetscFunctionBegin; 28 ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 29 *missing = PETSC_FALSE; 30 if (A->rmap->n > 0 && !ii) { 31 *missing = PETSC_TRUE; 32 if (dd) *dd = 0; 33 ierr = PetscInfo(A,"Matrix has no entries therefore is missing diagonal\n");CHKERRQ(ierr); 34 } else { 35 diag = a->diag; 36 for (i=0; i<a->mbs; i++) { 37 if (diag[i] >= ii[i+1]) { 38 *missing = PETSC_TRUE; 39 if (dd) *dd = i; 40 break; 41 } 42 } 43 } 44 PetscFunctionReturn(0); 45 } 46 47 PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 48 { 49 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 50 PetscErrorCode ierr; 51 PetscInt i,j; 52 53 PetscFunctionBegin; 54 if (!a->diag) { 55 ierr = PetscMalloc1(a->mbs,&a->diag);CHKERRQ(ierr); 56 ierr = PetscLogObjectMemory((PetscObject)A,a->mbs*sizeof(PetscInt));CHKERRQ(ierr); 57 a->free_diag = PETSC_TRUE; 58 } 59 for (i=0; i<a->mbs; i++) { 60 a->diag[i] = a->i[i+1]; 61 for (j=a->i[i]; j<a->i[i+1]; j++) { 62 if (a->j[j] == i) { 63 a->diag[i] = j; 64 break; 65 } 66 } 67 } 68 PetscFunctionReturn(0); 69 } 70 71 static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 72 { 73 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 74 PetscErrorCode ierr; 75 PetscInt i,j,n = a->mbs,nz = a->i[n],*tia,*tja,bs = A->rmap->bs,k,l,cnt; 76 PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 77 78 PetscFunctionBegin; 79 *nn = n; 80 if (!ia) PetscFunctionReturn(0); 81 if (symmetric) { 82 ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,PETSC_FALSE,0,0,&tia,&tja);CHKERRQ(ierr); 83 nz = tia[n]; 84 } else { 85 tia = a->i; tja = a->j; 86 } 87 88 if (!blockcompressed && bs > 1) { 89 (*nn) *= bs; 90 /* malloc & create the natural set of indices */ 91 ierr = PetscMalloc1((n+1)*bs,ia);CHKERRQ(ierr); 92 if (n) { 93 (*ia)[0] = oshift; 94 for (j=1; j<bs; j++) { 95 (*ia)[j] = (tia[1]-tia[0])*bs+(*ia)[j-1]; 96 } 97 } 98 99 for (i=1; i<n; i++) { 100 (*ia)[i*bs] = (tia[i]-tia[i-1])*bs + (*ia)[i*bs-1]; 101 for (j=1; j<bs; j++) { 102 (*ia)[i*bs+j] = (tia[i+1]-tia[i])*bs + (*ia)[i*bs+j-1]; 103 } 104 } 105 if (n) { 106 (*ia)[n*bs] = (tia[n]-tia[n-1])*bs + (*ia)[n*bs-1]; 107 } 108 109 if (inja) { 110 ierr = PetscMalloc1(nz*bs*bs,ja);CHKERRQ(ierr); 111 cnt = 0; 112 for (i=0; i<n; i++) { 113 for (j=0; j<bs; j++) { 114 for (k=tia[i]; k<tia[i+1]; k++) { 115 for (l=0; l<bs; l++) { 116 (*ja)[cnt++] = bs*tja[k] + l; 117 } 118 } 119 } 120 } 121 } 122 123 if (symmetric) { /* deallocate memory allocated in MatToSymmetricIJ_SeqAIJ() */ 124 ierr = PetscFree(tia);CHKERRQ(ierr); 125 ierr = PetscFree(tja);CHKERRQ(ierr); 126 } 127 } else if (oshift == 1) { 128 if (symmetric) { 129 nz = tia[A->rmap->n/bs]; 130 /* add 1 to i and j indices */ 131 for (i=0; i<A->rmap->n/bs+1; i++) tia[i] = tia[i] + 1; 132 *ia = tia; 133 if (ja) { 134 for (i=0; i<nz; i++) tja[i] = tja[i] + 1; 135 *ja = tja; 136 } 137 } else { 138 nz = a->i[A->rmap->n/bs]; 139 /* malloc space and add 1 to i and j indices */ 140 ierr = PetscMalloc1(A->rmap->n/bs+1,ia);CHKERRQ(ierr); 141 for (i=0; i<A->rmap->n/bs+1; i++) (*ia)[i] = a->i[i] + 1; 142 if (ja) { 143 ierr = PetscMalloc1(nz,ja);CHKERRQ(ierr); 144 for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1; 145 } 146 } 147 } else { 148 *ia = tia; 149 if (ja) *ja = tja; 150 } 151 PetscFunctionReturn(0); 152 } 153 154 static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done) 155 { 156 PetscErrorCode ierr; 157 158 PetscFunctionBegin; 159 if (!ia) PetscFunctionReturn(0); 160 if ((!blockcompressed && A->rmap->bs > 1) || (symmetric || oshift == 1)) { 161 ierr = PetscFree(*ia);CHKERRQ(ierr); 162 if (ja) {ierr = PetscFree(*ja);CHKERRQ(ierr);} 163 } 164 PetscFunctionReturn(0); 165 } 166 167 PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 168 { 169 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 170 PetscErrorCode ierr; 171 172 PetscFunctionBegin; 173 #if defined(PETSC_USE_LOG) 174 PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz); 175 #endif 176 ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr); 177 if (a->free_diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 178 ierr = ISDestroy(&a->row);CHKERRQ(ierr); 179 ierr = ISDestroy(&a->col);CHKERRQ(ierr); 180 ierr = ISDestroy(&a->icol);CHKERRQ(ierr); 181 ierr = PetscFree(a->idiag);CHKERRQ(ierr); 182 ierr = PetscFree(a->inode.size);CHKERRQ(ierr); 183 if (a->free_imax_ilen) {ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr);} 184 ierr = PetscFree(a->solve_work);CHKERRQ(ierr); 185 ierr = PetscFree(a->sor_work);CHKERRQ(ierr); 186 ierr = PetscFree(a->solves_work);CHKERRQ(ierr); 187 ierr = PetscFree(a->mult_work);CHKERRQ(ierr); 188 ierr = PetscFree(a->saved_values);CHKERRQ(ierr); 189 if (a->free_jshort) {ierr = PetscFree(a->jshort);CHKERRQ(ierr);} 190 ierr = PetscFree(a->inew);CHKERRQ(ierr); 191 ierr = MatDestroy(&a->parent);CHKERRQ(ierr); 192 ierr = PetscFree(A->data);CHKERRQ(ierr); 193 194 ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 195 ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);CHKERRQ(ierr); 196 ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);CHKERRQ(ierr); 197 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C",NULL);CHKERRQ(ierr); 198 ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C",NULL);CHKERRQ(ierr); 199 ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C",NULL);CHKERRQ(ierr); 200 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",NULL);CHKERRQ(ierr); 201 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocationCSR_C",NULL);CHKERRQ(ierr); 202 #if defined(PETSC_HAVE_ELEMENTAL) 203 ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_elemental_C",NULL);CHKERRQ(ierr); 204 #endif 205 PetscFunctionReturn(0); 206 } 207 208 PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg) 209 { 210 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 211 PetscErrorCode ierr; 212 213 PetscFunctionBegin; 214 switch (op) { 215 case MAT_ROW_ORIENTED: 216 a->roworiented = flg; 217 break; 218 case MAT_KEEP_NONZERO_PATTERN: 219 a->keepnonzeropattern = flg; 220 break; 221 case MAT_NEW_NONZERO_LOCATIONS: 222 a->nonew = (flg ? 0 : 1); 223 break; 224 case MAT_NEW_NONZERO_LOCATION_ERR: 225 a->nonew = (flg ? -1 : 0); 226 break; 227 case MAT_NEW_NONZERO_ALLOCATION_ERR: 228 a->nonew = (flg ? -2 : 0); 229 break; 230 case MAT_UNUSED_NONZERO_LOCATION_ERR: 231 a->nounused = (flg ? -1 : 0); 232 break; 233 case MAT_NEW_DIAGONALS: 234 case MAT_IGNORE_OFF_PROC_ENTRIES: 235 case MAT_USE_HASH_TABLE: 236 ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 237 break; 238 case MAT_HERMITIAN: 239 #if defined(PETSC_USE_COMPLEX) /* MAT_HERMITIAN is a synonym for MAT_SYMMETRIC with reals */ 240 if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first"); 241 if (A->cmap->n < 65536 && A->cmap->bs == 1) { 242 A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort; 243 } else if (A->cmap->bs == 1) { 244 A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian; 245 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1"); 246 #endif 247 break; 248 case MAT_SPD: 249 /* These options are handled directly by MatSetOption() */ 250 break; 251 case MAT_SYMMETRIC: 252 case MAT_STRUCTURALLY_SYMMETRIC: 253 case MAT_SYMMETRY_ETERNAL: 254 case MAT_STRUCTURE_ONLY: 255 /* These options are handled directly by MatSetOption() */ 256 break; 257 case MAT_IGNORE_LOWER_TRIANGULAR: 258 a->ignore_ltriangular = flg; 259 break; 260 case MAT_ERROR_LOWER_TRIANGULAR: 261 a->ignore_ltriangular = flg; 262 break; 263 case MAT_GETROW_UPPERTRIANGULAR: 264 a->getrow_utriangular = flg; 265 break; 266 case MAT_SUBMAT_SINGLEIS: 267 break; 268 default: 269 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op); 270 } 271 PetscFunctionReturn(0); 272 } 273 274 PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 275 { 276 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 277 PetscErrorCode ierr; 278 279 PetscFunctionBegin; 280 if (A && !a->getrow_utriangular) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()"); 281 282 /* Get the upper triangular part of the row */ 283 ierr = MatGetRow_SeqBAIJ_private(A,row,nz,idx,v,a->i,a->j,a->a);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 288 { 289 PetscErrorCode ierr; 290 291 PetscFunctionBegin; 292 if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 293 if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 294 PetscFunctionReturn(0); 295 } 296 297 PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A) 298 { 299 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 300 301 PetscFunctionBegin; 302 a->getrow_utriangular = PETSC_TRUE; 303 PetscFunctionReturn(0); 304 } 305 306 PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A) 307 { 308 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 309 310 PetscFunctionBegin; 311 a->getrow_utriangular = PETSC_FALSE; 312 PetscFunctionReturn(0); 313 } 314 315 PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B) 316 { 317 PetscErrorCode ierr; 318 319 PetscFunctionBegin; 320 if (reuse == MAT_INITIAL_MATRIX) { 321 ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 322 } else if (reuse == MAT_REUSE_MATRIX) { 323 ierr = MatCopy(A,*B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 324 } 325 PetscFunctionReturn(0); 326 } 327 328 PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 329 { 330 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 331 PetscErrorCode ierr; 332 PetscInt i,j,bs = A->rmap->bs,k,l,bs2=a->bs2; 333 PetscViewerFormat format; 334 PetscInt *diag; 335 336 PetscFunctionBegin; 337 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 338 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 339 ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 340 } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 341 Mat aij; 342 const char *matname; 343 344 if (A->factortype && bs>1) { 345 ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");CHKERRQ(ierr); 346 PetscFunctionReturn(0); 347 } 348 ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr); 349 ierr = PetscObjectGetName((PetscObject)A,&matname);CHKERRQ(ierr); 350 ierr = PetscObjectSetName((PetscObject)aij,matname);CHKERRQ(ierr); 351 ierr = MatView(aij,viewer);CHKERRQ(ierr); 352 ierr = MatDestroy(&aij);CHKERRQ(ierr); 353 } else if (format == PETSC_VIEWER_ASCII_COMMON) { 354 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 355 for (i=0; i<a->mbs; i++) { 356 for (j=0; j<bs; j++) { 357 ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 358 for (k=a->i[i]; k<a->i[i+1]; k++) { 359 for (l=0; l<bs; l++) { 360 #if defined(PETSC_USE_COMPLEX) 361 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 362 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 363 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 364 } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 365 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 366 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 367 } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 368 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 369 } 370 #else 371 if (a->a[bs2*k + l*bs + j] != 0.0) { 372 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 373 } 374 #endif 375 } 376 } 377 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 378 } 379 } 380 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 381 } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 382 PetscFunctionReturn(0); 383 } else { 384 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 385 if (A->factortype) { /* for factored matrix */ 386 if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet"); 387 388 diag=a->diag; 389 for (i=0; i<a->mbs; i++) { /* for row block i */ 390 ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 391 /* diagonal entry */ 392 #if defined(PETSC_USE_COMPLEX) 393 if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) { 394 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),(double)PetscImaginaryPart(1.0/a->a[diag[i]]));CHKERRQ(ierr); 395 } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) { 396 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),-(double)PetscImaginaryPart(1.0/a->a[diag[i]]));CHKERRQ(ierr); 397 } else { 398 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]));CHKERRQ(ierr); 399 } 400 #else 401 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)(1.0/a->a[diag[i]]));CHKERRQ(ierr); 402 #endif 403 /* off-diagonal entries */ 404 for (k=a->i[i]; k<a->i[i+1]-1; k++) { 405 #if defined(PETSC_USE_COMPLEX) 406 if (PetscImaginaryPart(a->a[k]) > 0.0) { 407 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),(double)PetscImaginaryPart(a->a[k]));CHKERRQ(ierr); 408 } else if (PetscImaginaryPart(a->a[k]) < 0.0) { 409 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),-(double)PetscImaginaryPart(a->a[k]));CHKERRQ(ierr); 410 } else { 411 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k],(double)PetscRealPart(a->a[k]));CHKERRQ(ierr); 412 } 413 #else 414 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[k],(double)a->a[k]);CHKERRQ(ierr); 415 #endif 416 } 417 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 418 } 419 420 } else { /* for non-factored matrix */ 421 for (i=0; i<a->mbs; i++) { /* for row block i */ 422 for (j=0; j<bs; j++) { /* for row bs*i + j */ 423 ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 424 for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */ 425 for (l=0; l<bs; l++) { /* for column */ 426 #if defined(PETSC_USE_COMPLEX) 427 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 428 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 429 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 430 } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 431 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 432 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 433 } else { 434 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 435 } 436 #else 437 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 438 #endif 439 } 440 } 441 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 442 } 443 } 444 } 445 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 446 } 447 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 448 PetscFunctionReturn(0); 449 } 450 451 #include <petscdraw.h> 452 static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 453 { 454 Mat A = (Mat) Aa; 455 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 456 PetscErrorCode ierr; 457 PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2; 458 PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 459 MatScalar *aa; 460 PetscViewer viewer; 461 462 PetscFunctionBegin; 463 ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 464 ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 465 466 /* loop over matrix elements drawing boxes */ 467 468 ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 469 ierr = PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");CHKERRQ(ierr); 470 /* Blue for negative, Cyan for zero and Red for positive */ 471 color = PETSC_DRAW_BLUE; 472 for (i=0,row=0; i<mbs; i++,row+=bs) { 473 for (j=a->i[i]; j<a->i[i+1]; j++) { 474 y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 475 x_l = a->j[j]*bs; x_r = x_l + 1.0; 476 aa = a->a + j*bs2; 477 for (k=0; k<bs; k++) { 478 for (l=0; l<bs; l++) { 479 if (PetscRealPart(*aa++) >= 0.) continue; 480 ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 481 } 482 } 483 } 484 } 485 color = PETSC_DRAW_CYAN; 486 for (i=0,row=0; i<mbs; i++,row+=bs) { 487 for (j=a->i[i]; j<a->i[i+1]; j++) { 488 y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 489 x_l = a->j[j]*bs; x_r = x_l + 1.0; 490 aa = a->a + j*bs2; 491 for (k=0; k<bs; k++) { 492 for (l=0; l<bs; l++) { 493 if (PetscRealPart(*aa++) != 0.) continue; 494 ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 495 } 496 } 497 } 498 } 499 color = PETSC_DRAW_RED; 500 for (i=0,row=0; i<mbs; i++,row+=bs) { 501 for (j=a->i[i]; j<a->i[i+1]; j++) { 502 y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 503 x_l = a->j[j]*bs; x_r = x_l + 1.0; 504 aa = a->a + j*bs2; 505 for (k=0; k<bs; k++) { 506 for (l=0; l<bs; l++) { 507 if (PetscRealPart(*aa++) <= 0.) continue; 508 ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 509 } 510 } 511 } 512 } 513 ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 514 PetscFunctionReturn(0); 515 } 516 517 static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 518 { 519 PetscErrorCode ierr; 520 PetscReal xl,yl,xr,yr,w,h; 521 PetscDraw draw; 522 PetscBool isnull; 523 524 PetscFunctionBegin; 525 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 526 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 527 if (isnull) PetscFunctionReturn(0); 528 529 xr = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0; 530 xr += w; yr += h; xl = -w; yl = -h; 531 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 532 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 533 ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 534 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);CHKERRQ(ierr); 535 ierr = PetscDrawSave(draw);CHKERRQ(ierr); 536 PetscFunctionReturn(0); 537 } 538 539 PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 540 { 541 PetscErrorCode ierr; 542 PetscBool iascii,isdraw; 543 FILE *file = 0; 544 545 PetscFunctionBegin; 546 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 547 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 548 if (iascii) { 549 ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 550 } else if (isdraw) { 551 ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 552 } else { 553 Mat B; 554 const char *matname; 555 ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 556 ierr = PetscObjectGetName((PetscObject)A,&matname);CHKERRQ(ierr); 557 ierr = PetscObjectSetName((PetscObject)B,matname);CHKERRQ(ierr); 558 ierr = MatView(B,viewer);CHKERRQ(ierr); 559 ierr = MatDestroy(&B);CHKERRQ(ierr); 560 ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 561 if (file) { 562 fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs); 563 } 564 } 565 PetscFunctionReturn(0); 566 } 567 568 569 PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 570 { 571 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 572 PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 573 PetscInt *ai = a->i,*ailen = a->ilen; 574 PetscInt brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2; 575 MatScalar *ap,*aa = a->a; 576 577 PetscFunctionBegin; 578 for (k=0; k<m; k++) { /* loop over rows */ 579 row = im[k]; brow = row/bs; 580 if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */ 581 if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1); 582 rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 583 nrow = ailen[brow]; 584 for (l=0; l<n; l++) { /* loop over columns */ 585 if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */ 586 if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1); 587 col = in[l]; 588 bcol = col/bs; 589 cidx = col%bs; 590 ridx = row%bs; 591 high = nrow; 592 low = 0; /* assume unsorted */ 593 while (high-low > 5) { 594 t = (low+high)/2; 595 if (rp[t] > bcol) high = t; 596 else low = t; 597 } 598 for (i=low; i<high; i++) { 599 if (rp[i] > bcol) break; 600 if (rp[i] == bcol) { 601 *v++ = ap[bs2*i+bs*cidx+ridx]; 602 goto finished; 603 } 604 } 605 *v++ = 0.0; 606 finished:; 607 } 608 } 609 PetscFunctionReturn(0); 610 } 611 612 613 PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 614 { 615 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 616 PetscErrorCode ierr; 617 PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 618 PetscInt *imax =a->imax,*ai=a->i,*ailen=a->ilen; 619 PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval; 620 PetscBool roworiented=a->roworiented; 621 const PetscScalar *value = v; 622 MatScalar *ap,*aa = a->a,*bap; 623 624 PetscFunctionBegin; 625 if (roworiented) stepval = (n-1)*bs; 626 else stepval = (m-1)*bs; 627 628 for (k=0; k<m; k++) { /* loop over added rows */ 629 row = im[k]; 630 if (row < 0) continue; 631 #if defined(PETSC_USE_DEBUG) 632 if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index row too large %D max %D",row,a->mbs-1); 633 #endif 634 rp = aj + ai[row]; 635 ap = aa + bs2*ai[row]; 636 rmax = imax[row]; 637 nrow = ailen[row]; 638 low = 0; 639 high = nrow; 640 for (l=0; l<n; l++) { /* loop over added columns */ 641 if (in[l] < 0) continue; 642 col = in[l]; 643 #if defined(PETSC_USE_DEBUG) 644 if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index column too large %D max %D",col,a->nbs-1); 645 #endif 646 if (col < row) { 647 if (a->ignore_ltriangular) continue; /* ignore lower triangular block */ 648 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)"); 649 } 650 if (roworiented) value = v + k*(stepval+bs)*bs + l*bs; 651 else value = v + l*(stepval+bs)*bs + k*bs; 652 653 if (col <= lastcol) low = 0; 654 else high = nrow; 655 656 lastcol = col; 657 while (high-low > 7) { 658 t = (low+high)/2; 659 if (rp[t] > col) high = t; 660 else low = t; 661 } 662 for (i=low; i<high; i++) { 663 if (rp[i] > col) break; 664 if (rp[i] == col) { 665 bap = ap + bs2*i; 666 if (roworiented) { 667 if (is == ADD_VALUES) { 668 for (ii=0; ii<bs; ii++,value+=stepval) { 669 for (jj=ii; jj<bs2; jj+=bs) { 670 bap[jj] += *value++; 671 } 672 } 673 } else { 674 for (ii=0; ii<bs; ii++,value+=stepval) { 675 for (jj=ii; jj<bs2; jj+=bs) { 676 bap[jj] = *value++; 677 } 678 } 679 } 680 } else { 681 if (is == ADD_VALUES) { 682 for (ii=0; ii<bs; ii++,value+=stepval) { 683 for (jj=0; jj<bs; jj++) { 684 *bap++ += *value++; 685 } 686 } 687 } else { 688 for (ii=0; ii<bs; ii++,value+=stepval) { 689 for (jj=0; jj<bs; jj++) { 690 *bap++ = *value++; 691 } 692 } 693 } 694 } 695 goto noinsert2; 696 } 697 } 698 if (nonew == 1) goto noinsert2; 699 if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new block index nonzero block (%D, %D) in the matrix", row, col); 700 MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 701 N = nrow++ - 1; high++; 702 /* shift up all the later entries in this row */ 703 for (ii=N; ii>=i; ii--) { 704 rp[ii+1] = rp[ii]; 705 ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 706 } 707 if (N >= i) { 708 ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 709 } 710 rp[i] = col; 711 bap = ap + bs2*i; 712 if (roworiented) { 713 for (ii=0; ii<bs; ii++,value+=stepval) { 714 for (jj=ii; jj<bs2; jj+=bs) { 715 bap[jj] = *value++; 716 } 717 } 718 } else { 719 for (ii=0; ii<bs; ii++,value+=stepval) { 720 for (jj=0; jj<bs; jj++) { 721 *bap++ = *value++; 722 } 723 } 724 } 725 noinsert2:; 726 low = i; 727 } 728 ailen[row] = nrow; 729 } 730 PetscFunctionReturn(0); 731 } 732 733 /* 734 This is not yet used 735 */ 736 PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A) 737 { 738 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 739 PetscErrorCode ierr; 740 const PetscInt *ai = a->i, *aj = a->j,*cols; 741 PetscInt i = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts; 742 PetscBool flag; 743 744 PetscFunctionBegin; 745 ierr = PetscMalloc1(m,&ns);CHKERRQ(ierr); 746 while (i < m) { 747 nzx = ai[i+1] - ai[i]; /* Number of nonzeros */ 748 /* Limits the number of elements in a node to 'a->inode.limit' */ 749 for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) { 750 nzy = ai[j+1] - ai[j]; 751 if (nzy != (nzx - j + i)) break; 752 ierr = PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);CHKERRQ(ierr); 753 if (!flag) break; 754 } 755 ns[node_count++] = blk_size; 756 757 i = j; 758 } 759 if (!a->inode.size && m && node_count > .9*m) { 760 ierr = PetscFree(ns);CHKERRQ(ierr); 761 ierr = PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);CHKERRQ(ierr); 762 } else { 763 a->inode.node_count = node_count; 764 765 ierr = PetscMalloc1(node_count,&a->inode.size);CHKERRQ(ierr); 766 ierr = PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));CHKERRQ(ierr); 767 ierr = PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));CHKERRQ(ierr); 768 ierr = PetscFree(ns);CHKERRQ(ierr); 769 ierr = PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);CHKERRQ(ierr); 770 771 /* count collections of adjacent columns in each inode */ 772 row = 0; 773 cnt = 0; 774 for (i=0; i<node_count; i++) { 775 cols = aj + ai[row] + a->inode.size[i]; 776 nz = ai[row+1] - ai[row] - a->inode.size[i]; 777 for (j=1; j<nz; j++) { 778 if (cols[j] != cols[j-1]+1) cnt++; 779 } 780 cnt++; 781 row += a->inode.size[i]; 782 } 783 ierr = PetscMalloc1(2*cnt,&counts);CHKERRQ(ierr); 784 cnt = 0; 785 row = 0; 786 for (i=0; i<node_count; i++) { 787 cols = aj + ai[row] + a->inode.size[i]; 788 counts[2*cnt] = cols[0]; 789 nz = ai[row+1] - ai[row] - a->inode.size[i]; 790 cnt2 = 1; 791 for (j=1; j<nz; j++) { 792 if (cols[j] != cols[j-1]+1) { 793 counts[2*(cnt++)+1] = cnt2; 794 counts[2*cnt] = cols[j]; 795 cnt2 = 1; 796 } else cnt2++; 797 } 798 counts[2*(cnt++)+1] = cnt2; 799 row += a->inode.size[i]; 800 } 801 ierr = PetscIntView(2*cnt,counts,0);CHKERRQ(ierr); 802 } 803 PetscFunctionReturn(0); 804 } 805 806 PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 807 { 808 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 809 PetscErrorCode ierr; 810 PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 811 PetscInt m = A->rmap->N,*ip,N,*ailen = a->ilen; 812 PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 813 MatScalar *aa = a->a,*ap; 814 815 PetscFunctionBegin; 816 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 817 818 if (m) rmax = ailen[0]; 819 for (i=1; i<mbs; i++) { 820 /* move each row back by the amount of empty slots (fshift) before it*/ 821 fshift += imax[i-1] - ailen[i-1]; 822 rmax = PetscMax(rmax,ailen[i]); 823 if (fshift) { 824 ip = aj + ai[i]; ap = aa + bs2*ai[i]; 825 N = ailen[i]; 826 for (j=0; j<N; j++) { 827 ip[j-fshift] = ip[j]; 828 ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 829 } 830 } 831 ai[i] = ai[i-1] + ailen[i-1]; 832 } 833 if (mbs) { 834 fshift += imax[mbs-1] - ailen[mbs-1]; 835 ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 836 } 837 /* reset ilen and imax for each row */ 838 for (i=0; i<mbs; i++) { 839 ailen[i] = imax[i] = ai[i+1] - ai[i]; 840 } 841 a->nz = ai[mbs]; 842 843 /* diagonals may have moved, reset it */ 844 if (a->diag) { 845 ierr = PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));CHKERRQ(ierr); 846 } 847 if (fshift && a->nounused == -1) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Unused space detected in matrix: %D X %D block size %D, %D unneeded", m, A->cmap->n, A->rmap->bs, fshift*bs2); 848 849 ierr = PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap->N,A->rmap->bs,fshift*bs2,a->nz*bs2);CHKERRQ(ierr); 850 ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 851 ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 852 853 A->info.mallocs += a->reallocs; 854 a->reallocs = 0; 855 A->info.nz_unneeded = (PetscReal)fshift*bs2; 856 a->idiagvalid = PETSC_FALSE; 857 a->rmax = rmax; 858 859 if (A->cmap->n < 65536 && A->cmap->bs == 1) { 860 if (a->jshort && a->free_jshort) { 861 /* when matrix data structure is changed, previous jshort must be replaced */ 862 ierr = PetscFree(a->jshort);CHKERRQ(ierr); 863 } 864 ierr = PetscMalloc1(a->i[A->rmap->n],&a->jshort);CHKERRQ(ierr); 865 ierr = PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));CHKERRQ(ierr); 866 for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i]; 867 A->ops->mult = MatMult_SeqSBAIJ_1_ushort; 868 A->ops->sor = MatSOR_SeqSBAIJ_ushort; 869 a->free_jshort = PETSC_TRUE; 870 } 871 PetscFunctionReturn(0); 872 } 873 874 /* 875 This function returns an array of flags which indicate the locations of contiguous 876 blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 877 then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 878 Assume: sizes should be long enough to hold all the values. 879 */ 880 PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 881 { 882 PetscInt i,j,k,row; 883 PetscBool flg; 884 885 PetscFunctionBegin; 886 for (i=0,j=0; i<n; j++) { 887 row = idx[i]; 888 if (row%bs!=0) { /* Not the begining of a block */ 889 sizes[j] = 1; 890 i++; 891 } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 892 sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 893 i++; 894 } else { /* Begining of the block, so check if the complete block exists */ 895 flg = PETSC_TRUE; 896 for (k=1; k<bs; k++) { 897 if (row+k != idx[i+k]) { /* break in the block */ 898 flg = PETSC_FALSE; 899 break; 900 } 901 } 902 if (flg) { /* No break in the bs */ 903 sizes[j] = bs; 904 i += bs; 905 } else { 906 sizes[j] = 1; 907 i++; 908 } 909 } 910 } 911 *bs_max = j; 912 PetscFunctionReturn(0); 913 } 914 915 916 /* Only add/insert a(i,j) with i<=j (blocks). 917 Any a(i,j) with i>j input by user is ingored. 918 */ 919 920 PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 921 { 922 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 923 PetscErrorCode ierr; 924 PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 925 PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 926 PetscInt *aj =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 927 PetscInt ridx,cidx,bs2=a->bs2; 928 MatScalar *ap,value,*aa=a->a,*bap; 929 930 PetscFunctionBegin; 931 for (k=0; k<m; k++) { /* loop over added rows */ 932 row = im[k]; /* row number */ 933 brow = row/bs; /* block row number */ 934 if (row < 0) continue; 935 #if defined(PETSC_USE_DEBUG) 936 if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1); 937 #endif 938 rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 939 ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 940 rmax = imax[brow]; /* maximum space allocated for this row */ 941 nrow = ailen[brow]; /* actual length of this row */ 942 low = 0; 943 high = nrow; 944 for (l=0; l<n; l++) { /* loop over added columns */ 945 if (in[l] < 0) continue; 946 #if defined(PETSC_USE_DEBUG) 947 if (in[l] >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->rmap->N-1); 948 #endif 949 col = in[l]; 950 bcol = col/bs; /* block col number */ 951 952 if (brow > bcol) { 953 if (a->ignore_ltriangular) continue; /* ignore lower triangular values */ 954 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)"); 955 } 956 957 ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 958 if ((brow==bcol && ridx<=cidx) || (brow<bcol)) { 959 /* element value a(k,l) */ 960 if (roworiented) value = v[l + k*n]; 961 else value = v[k + l*m]; 962 963 /* move pointer bap to a(k,l) quickly and add/insert value */ 964 if (col <= lastcol) low = 0; 965 else high = nrow; 966 967 lastcol = col; 968 while (high-low > 7) { 969 t = (low+high)/2; 970 if (rp[t] > bcol) high = t; 971 else low = t; 972 } 973 for (i=low; i<high; i++) { 974 if (rp[i] > bcol) break; 975 if (rp[i] == bcol) { 976 bap = ap + bs2*i + bs*cidx + ridx; 977 if (is == ADD_VALUES) *bap += value; 978 else *bap = value; 979 /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 980 if (brow == bcol && ridx < cidx) { 981 bap = ap + bs2*i + bs*ridx + cidx; 982 if (is == ADD_VALUES) *bap += value; 983 else *bap = value; 984 } 985 goto noinsert1; 986 } 987 } 988 989 if (nonew == 1) goto noinsert1; 990 if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 991 MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 992 993 N = nrow++ - 1; high++; 994 /* shift up all the later entries in this row */ 995 for (ii=N; ii>=i; ii--) { 996 rp[ii+1] = rp[ii]; 997 ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 998 } 999 if (N>=i) { 1000 ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1001 } 1002 rp[i] = bcol; 1003 ap[bs2*i + bs*cidx + ridx] = value; 1004 /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 1005 if (brow == bcol && ridx < cidx) { 1006 ap[bs2*i + bs*ridx + cidx] = value; 1007 } 1008 A->nonzerostate++; 1009 noinsert1:; 1010 low = i; 1011 } 1012 } /* end of loop over added columns */ 1013 ailen[brow] = nrow; 1014 } /* end of loop over added rows */ 1015 PetscFunctionReturn(0); 1016 } 1017 1018 PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info) 1019 { 1020 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 1021 Mat outA; 1022 PetscErrorCode ierr; 1023 PetscBool row_identity; 1024 1025 PetscFunctionBegin; 1026 if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc"); 1027 ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1028 if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported"); 1029 if (inA->rmap->bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */ 1030 1031 outA = inA; 1032 inA->factortype = MAT_FACTOR_ICC; 1033 ierr = PetscFree(inA->solvertype);CHKERRQ(ierr); 1034 ierr = PetscStrallocpy(MATSOLVERPETSC,&inA->solvertype);CHKERRQ(ierr); 1035 1036 ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 1037 ierr = MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);CHKERRQ(ierr); 1038 1039 ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1040 ierr = ISDestroy(&a->row);CHKERRQ(ierr); 1041 a->row = row; 1042 ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1043 ierr = ISDestroy(&a->col);CHKERRQ(ierr); 1044 a->col = row; 1045 1046 /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */ 1047 if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);} 1048 ierr = PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);CHKERRQ(ierr); 1049 1050 if (!a->solve_work) { 1051 ierr = PetscMalloc1(inA->rmap->N+inA->rmap->bs,&a->solve_work);CHKERRQ(ierr); 1052 ierr = PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));CHKERRQ(ierr); 1053 } 1054 1055 ierr = MatCholeskyFactorNumeric(outA,inA,info);CHKERRQ(ierr); 1056 PetscFunctionReturn(0); 1057 } 1058 1059 PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 1060 { 1061 Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ*)mat->data; 1062 PetscInt i,nz,n; 1063 PetscErrorCode ierr; 1064 1065 PetscFunctionBegin; 1066 nz = baij->maxnz; 1067 n = mat->cmap->n; 1068 for (i=0; i<nz; i++) baij->j[i] = indices[i]; 1069 1070 baij->nz = nz; 1071 for (i=0; i<n; i++) baij->ilen[i] = baij->imax[i]; 1072 1073 ierr = MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1074 PetscFunctionReturn(0); 1075 } 1076 1077 /*@ 1078 MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 1079 in the matrix. 1080 1081 Input Parameters: 1082 + mat - the SeqSBAIJ matrix 1083 - indices - the column indices 1084 1085 Level: advanced 1086 1087 Notes: 1088 This can be called if you have precomputed the nonzero structure of the 1089 matrix and want to provide it to the matrix object to improve the performance 1090 of the MatSetValues() operation. 1091 1092 You MUST have set the correct numbers of nonzeros per row in the call to 1093 MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 1094 1095 MUST be called before any calls to MatSetValues() 1096 1097 .seealso: MatCreateSeqSBAIJ 1098 @*/ 1099 PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 1100 { 1101 PetscErrorCode ierr; 1102 1103 PetscFunctionBegin; 1104 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1105 PetscValidPointer(indices,2); 1106 ierr = PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));CHKERRQ(ierr); 1107 PetscFunctionReturn(0); 1108 } 1109 1110 PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 1111 { 1112 PetscErrorCode ierr; 1113 1114 PetscFunctionBegin; 1115 /* If the two matrices have the same copy implementation, use fast copy. */ 1116 if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 1117 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1118 Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 1119 1120 if (a->i[A->rmap->N] != b->i[B->rmap->N]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 1121 ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));CHKERRQ(ierr); 1122 ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr); 1123 } else { 1124 ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 1125 ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1126 ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 1127 } 1128 PetscFunctionReturn(0); 1129 } 1130 1131 PetscErrorCode MatSetUp_SeqSBAIJ(Mat A) 1132 { 1133 PetscErrorCode ierr; 1134 1135 PetscFunctionBegin; 1136 ierr = MatSeqSBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0);CHKERRQ(ierr); 1137 PetscFunctionReturn(0); 1138 } 1139 1140 PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1141 { 1142 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1143 1144 PetscFunctionBegin; 1145 *array = a->a; 1146 PetscFunctionReturn(0); 1147 } 1148 1149 PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1150 { 1151 PetscFunctionBegin; 1152 PetscFunctionReturn(0); 1153 } 1154 1155 PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y,Mat X,PetscInt *nnz) 1156 { 1157 PetscInt bs = Y->rmap->bs,mbs = Y->rmap->N/bs; 1158 Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ*)X->data; 1159 Mat_SeqSBAIJ *y = (Mat_SeqSBAIJ*)Y->data; 1160 PetscErrorCode ierr; 1161 1162 PetscFunctionBegin; 1163 /* Set the number of nonzeros in the new matrix */ 1164 ierr = MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz);CHKERRQ(ierr); 1165 PetscFunctionReturn(0); 1166 } 1167 1168 PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 1169 { 1170 Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data; 1171 PetscErrorCode ierr; 1172 PetscInt bs=Y->rmap->bs,bs2=bs*bs; 1173 PetscBLASInt one = 1; 1174 1175 PetscFunctionBegin; 1176 if (str == SAME_NONZERO_PATTERN) { 1177 PetscScalar alpha = a; 1178 PetscBLASInt bnz; 1179 ierr = PetscBLASIntCast(x->nz*bs2,&bnz);CHKERRQ(ierr); 1180 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one)); 1181 ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 1182 } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1183 ierr = MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 1184 ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1185 ierr = MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);CHKERRQ(ierr); 1186 } else { 1187 Mat B; 1188 PetscInt *nnz; 1189 if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size"); 1190 ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1191 ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 1192 ierr = PetscMalloc1(Y->rmap->N,&nnz);CHKERRQ(ierr); 1193 ierr = MatCreate(PetscObjectComm((PetscObject)Y),&B);CHKERRQ(ierr); 1194 ierr = PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);CHKERRQ(ierr); 1195 ierr = MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);CHKERRQ(ierr); 1196 ierr = MatSetBlockSizesFromMats(B,Y,Y);CHKERRQ(ierr); 1197 ierr = MatSetType(B,(MatType) ((PetscObject)Y)->type_name);CHKERRQ(ierr); 1198 ierr = MatAXPYGetPreallocation_SeqSBAIJ(Y,X,nnz);CHKERRQ(ierr); 1199 ierr = MatSeqSBAIJSetPreallocation(B,bs,0,nnz);CHKERRQ(ierr); 1200 1201 ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr); 1202 1203 ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr); 1204 ierr = PetscFree(nnz);CHKERRQ(ierr); 1205 ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 1206 ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 1207 } 1208 PetscFunctionReturn(0); 1209 } 1210 1211 PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg) 1212 { 1213 PetscFunctionBegin; 1214 *flg = PETSC_TRUE; 1215 PetscFunctionReturn(0); 1216 } 1217 1218 PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool *flg) 1219 { 1220 PetscFunctionBegin; 1221 *flg = PETSC_TRUE; 1222 PetscFunctionReturn(0); 1223 } 1224 1225 PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg) 1226 { 1227 PetscFunctionBegin; 1228 *flg = PETSC_FALSE; 1229 PetscFunctionReturn(0); 1230 } 1231 1232 PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 1233 { 1234 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1235 PetscInt i,nz = a->bs2*a->i[a->mbs]; 1236 MatScalar *aa = a->a; 1237 1238 PetscFunctionBegin; 1239 for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1240 PetscFunctionReturn(0); 1241 } 1242 1243 PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 1244 { 1245 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1246 PetscInt i,nz = a->bs2*a->i[a->mbs]; 1247 MatScalar *aa = a->a; 1248 1249 PetscFunctionBegin; 1250 for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1251 PetscFunctionReturn(0); 1252 } 1253 1254 PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b) 1255 { 1256 Mat_SeqSBAIJ *baij=(Mat_SeqSBAIJ*)A->data; 1257 PetscErrorCode ierr; 1258 PetscInt i,j,k,count; 1259 PetscInt bs =A->rmap->bs,bs2=baij->bs2,row,col; 1260 PetscScalar zero = 0.0; 1261 MatScalar *aa; 1262 const PetscScalar *xx; 1263 PetscScalar *bb; 1264 PetscBool *zeroed,vecs = PETSC_FALSE; 1265 1266 PetscFunctionBegin; 1267 /* fix right hand side if needed */ 1268 if (x && b) { 1269 ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 1270 ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 1271 vecs = PETSC_TRUE; 1272 } 1273 1274 /* zero the columns */ 1275 ierr = PetscCalloc1(A->rmap->n,&zeroed);CHKERRQ(ierr); 1276 for (i=0; i<is_n; i++) { 1277 if (is_idx[i] < 0 || is_idx[i] >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range",is_idx[i]); 1278 zeroed[is_idx[i]] = PETSC_TRUE; 1279 } 1280 if (vecs) { 1281 for (i=0; i<A->rmap->N; i++) { 1282 row = i/bs; 1283 for (j=baij->i[row]; j<baij->i[row+1]; j++) { 1284 for (k=0; k<bs; k++) { 1285 col = bs*baij->j[j] + k; 1286 if (col <= i) continue; 1287 aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k; 1288 if (!zeroed[i] && zeroed[col]) bb[i] -= aa[0]*xx[col]; 1289 if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i]; 1290 } 1291 } 1292 } 1293 for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]]; 1294 } 1295 1296 for (i=0; i<A->rmap->N; i++) { 1297 if (!zeroed[i]) { 1298 row = i/bs; 1299 for (j=baij->i[row]; j<baij->i[row+1]; j++) { 1300 for (k=0; k<bs; k++) { 1301 col = bs*baij->j[j] + k; 1302 if (zeroed[col]) { 1303 aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k; 1304 aa[0] = 0.0; 1305 } 1306 } 1307 } 1308 } 1309 } 1310 ierr = PetscFree(zeroed);CHKERRQ(ierr); 1311 if (vecs) { 1312 ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 1313 ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 1314 } 1315 1316 /* zero the rows */ 1317 for (i=0; i<is_n; i++) { 1318 row = is_idx[i]; 1319 count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 1320 aa = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs); 1321 for (k=0; k<count; k++) { 1322 aa[0] = zero; 1323 aa += bs; 1324 } 1325 if (diag != 0.0) { 1326 ierr = (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);CHKERRQ(ierr); 1327 } 1328 } 1329 ierr = MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1330 PetscFunctionReturn(0); 1331 } 1332 1333 PetscErrorCode MatShift_SeqSBAIJ(Mat Y,PetscScalar a) 1334 { 1335 PetscErrorCode ierr; 1336 Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)Y->data; 1337 1338 PetscFunctionBegin; 1339 if (!Y->preallocated || !aij->nz) { 1340 ierr = MatSeqSBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL);CHKERRQ(ierr); 1341 } 1342 ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 1343 PetscFunctionReturn(0); 1344 } 1345 1346 /* -------------------------------------------------------------------*/ 1347 static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 1348 MatGetRow_SeqSBAIJ, 1349 MatRestoreRow_SeqSBAIJ, 1350 MatMult_SeqSBAIJ_N, 1351 /* 4*/ MatMultAdd_SeqSBAIJ_N, 1352 MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1353 MatMultAdd_SeqSBAIJ_N, 1354 0, 1355 0, 1356 0, 1357 /* 10*/ 0, 1358 0, 1359 MatCholeskyFactor_SeqSBAIJ, 1360 MatSOR_SeqSBAIJ, 1361 MatTranspose_SeqSBAIJ, 1362 /* 15*/ MatGetInfo_SeqSBAIJ, 1363 MatEqual_SeqSBAIJ, 1364 MatGetDiagonal_SeqSBAIJ, 1365 MatDiagonalScale_SeqSBAIJ, 1366 MatNorm_SeqSBAIJ, 1367 /* 20*/ 0, 1368 MatAssemblyEnd_SeqSBAIJ, 1369 MatSetOption_SeqSBAIJ, 1370 MatZeroEntries_SeqSBAIJ, 1371 /* 24*/ 0, 1372 0, 1373 0, 1374 0, 1375 0, 1376 /* 29*/ MatSetUp_SeqSBAIJ, 1377 0, 1378 0, 1379 0, 1380 0, 1381 /* 34*/ MatDuplicate_SeqSBAIJ, 1382 0, 1383 0, 1384 0, 1385 MatICCFactor_SeqSBAIJ, 1386 /* 39*/ MatAXPY_SeqSBAIJ, 1387 MatCreateSubMatrices_SeqSBAIJ, 1388 MatIncreaseOverlap_SeqSBAIJ, 1389 MatGetValues_SeqSBAIJ, 1390 MatCopy_SeqSBAIJ, 1391 /* 44*/ 0, 1392 MatScale_SeqSBAIJ, 1393 MatShift_SeqSBAIJ, 1394 0, 1395 MatZeroRowsColumns_SeqSBAIJ, 1396 /* 49*/ 0, 1397 MatGetRowIJ_SeqSBAIJ, 1398 MatRestoreRowIJ_SeqSBAIJ, 1399 0, 1400 0, 1401 /* 54*/ 0, 1402 0, 1403 0, 1404 0, 1405 MatSetValuesBlocked_SeqSBAIJ, 1406 /* 59*/ MatCreateSubMatrix_SeqSBAIJ, 1407 0, 1408 0, 1409 0, 1410 0, 1411 /* 64*/ 0, 1412 0, 1413 0, 1414 0, 1415 0, 1416 /* 69*/ MatGetRowMaxAbs_SeqSBAIJ, 1417 0, 1418 0, 1419 0, 1420 0, 1421 /* 74*/ 0, 1422 0, 1423 0, 1424 0, 1425 0, 1426 /* 79*/ 0, 1427 0, 1428 0, 1429 MatGetInertia_SeqSBAIJ, 1430 MatLoad_SeqSBAIJ, 1431 /* 84*/ MatIsSymmetric_SeqSBAIJ, 1432 MatIsHermitian_SeqSBAIJ, 1433 MatIsStructurallySymmetric_SeqSBAIJ, 1434 0, 1435 0, 1436 /* 89*/ 0, 1437 0, 1438 0, 1439 0, 1440 0, 1441 /* 94*/ 0, 1442 0, 1443 0, 1444 0, 1445 0, 1446 /* 99*/ 0, 1447 0, 1448 0, 1449 0, 1450 0, 1451 /*104*/ 0, 1452 MatRealPart_SeqSBAIJ, 1453 MatImaginaryPart_SeqSBAIJ, 1454 MatGetRowUpperTriangular_SeqSBAIJ, 1455 MatRestoreRowUpperTriangular_SeqSBAIJ, 1456 /*109*/ 0, 1457 0, 1458 0, 1459 0, 1460 MatMissingDiagonal_SeqSBAIJ, 1461 /*114*/ 0, 1462 0, 1463 0, 1464 0, 1465 0, 1466 /*119*/ 0, 1467 0, 1468 0, 1469 0, 1470 0, 1471 /*124*/ 0, 1472 0, 1473 0, 1474 0, 1475 0, 1476 /*129*/ 0, 1477 0, 1478 0, 1479 0, 1480 0, 1481 /*134*/ 0, 1482 0, 1483 0, 1484 0, 1485 0, 1486 /*139*/ MatSetBlockSizes_Default, 1487 0, 1488 0, 1489 0, 1490 0, 1491 /*144*/MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ 1492 }; 1493 1494 PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat) 1495 { 1496 Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)mat->data; 1497 PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2; 1498 PetscErrorCode ierr; 1499 1500 PetscFunctionBegin; 1501 if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 1502 1503 /* allocate space for values if not already there */ 1504 if (!aij->saved_values) { 1505 ierr = PetscMalloc1(nz+1,&aij->saved_values);CHKERRQ(ierr); 1506 } 1507 1508 /* copy values over */ 1509 ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 1510 PetscFunctionReturn(0); 1511 } 1512 1513 PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat) 1514 { 1515 Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)mat->data; 1516 PetscErrorCode ierr; 1517 PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2; 1518 1519 PetscFunctionBegin; 1520 if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 1521 if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 1522 1523 /* copy values over */ 1524 ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 1525 PetscFunctionReturn(0); 1526 } 1527 1528 static PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 1529 { 1530 Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 1531 PetscErrorCode ierr; 1532 PetscInt i,mbs,nbs,bs2; 1533 PetscBool skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE; 1534 1535 PetscFunctionBegin; 1536 if (nz >= 0 || nnz) realalloc = PETSC_TRUE; 1537 1538 ierr = MatSetBlockSize(B,PetscAbs(bs));CHKERRQ(ierr); 1539 ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 1540 ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 1541 ierr = PetscLayoutGetBlockSize(B->rmap,&bs);CHKERRQ(ierr); 1542 1543 B->preallocated = PETSC_TRUE; 1544 1545 mbs = B->rmap->N/bs; 1546 nbs = B->cmap->n/bs; 1547 bs2 = bs*bs; 1548 1549 if (mbs*bs != B->rmap->N || nbs*bs!=B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 1550 1551 if (nz == MAT_SKIP_ALLOCATION) { 1552 skipallocation = PETSC_TRUE; 1553 nz = 0; 1554 } 1555 1556 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1557 if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 1558 if (nnz) { 1559 for (i=0; i<mbs; i++) { 1560 if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 1561 if (nnz[i] > nbs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D block rowlength %D",i,nnz[i],nbs); 1562 } 1563 } 1564 1565 B->ops->mult = MatMult_SeqSBAIJ_N; 1566 B->ops->multadd = MatMultAdd_SeqSBAIJ_N; 1567 B->ops->multtranspose = MatMult_SeqSBAIJ_N; 1568 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N; 1569 1570 ierr = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);CHKERRQ(ierr); 1571 if (!flg) { 1572 switch (bs) { 1573 case 1: 1574 B->ops->mult = MatMult_SeqSBAIJ_1; 1575 B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1576 B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1577 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 1578 break; 1579 case 2: 1580 B->ops->mult = MatMult_SeqSBAIJ_2; 1581 B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1582 B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1583 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 1584 break; 1585 case 3: 1586 B->ops->mult = MatMult_SeqSBAIJ_3; 1587 B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1588 B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1589 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 1590 break; 1591 case 4: 1592 B->ops->mult = MatMult_SeqSBAIJ_4; 1593 B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1594 B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1595 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 1596 break; 1597 case 5: 1598 B->ops->mult = MatMult_SeqSBAIJ_5; 1599 B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1600 B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1601 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 1602 break; 1603 case 6: 1604 B->ops->mult = MatMult_SeqSBAIJ_6; 1605 B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1606 B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1607 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 1608 break; 1609 case 7: 1610 B->ops->mult = MatMult_SeqSBAIJ_7; 1611 B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1612 B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1613 B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 1614 break; 1615 } 1616 } 1617 1618 b->mbs = mbs; 1619 b->nbs = nbs; 1620 if (!skipallocation) { 1621 if (!b->imax) { 1622 ierr = PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);CHKERRQ(ierr); 1623 1624 b->free_imax_ilen = PETSC_TRUE; 1625 1626 ierr = PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt));CHKERRQ(ierr); 1627 } 1628 if (!nnz) { 1629 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 1630 else if (nz <= 0) nz = 1; 1631 for (i=0; i<mbs; i++) b->imax[i] = nz; 1632 nz = nz*mbs; /* total nz */ 1633 } else { 1634 nz = 0; 1635 for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 1636 } 1637 /* b->ilen will count nonzeros in each block row so far. */ 1638 for (i=0; i<mbs; i++) b->ilen[i] = 0; 1639 /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 1640 1641 /* allocate the matrix space */ 1642 ierr = MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);CHKERRQ(ierr); 1643 ierr = PetscMalloc3(bs2*nz,&b->a,nz,&b->j,B->rmap->N+1,&b->i);CHKERRQ(ierr); 1644 ierr = PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 1645 ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1646 ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 1647 1648 b->singlemalloc = PETSC_TRUE; 1649 1650 /* pointer to beginning of each row */ 1651 b->i[0] = 0; 1652 for (i=1; i<mbs+1; i++) b->i[i] = b->i[i-1] + b->imax[i-1]; 1653 1654 b->free_a = PETSC_TRUE; 1655 b->free_ij = PETSC_TRUE; 1656 } else { 1657 b->free_a = PETSC_FALSE; 1658 b->free_ij = PETSC_FALSE; 1659 } 1660 1661 b->bs2 = bs2; 1662 b->nz = 0; 1663 b->maxnz = nz; 1664 b->inew = 0; 1665 b->jnew = 0; 1666 b->anew = 0; 1667 b->a2anew = 0; 1668 b->permute = PETSC_FALSE; 1669 1670 B->was_assembled = PETSC_FALSE; 1671 B->assembled = PETSC_FALSE; 1672 if (realalloc) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);} 1673 PetscFunctionReturn(0); 1674 } 1675 1676 PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[], const PetscScalar V[]) 1677 { 1678 PetscInt i,j,m,nz,nz_max=0,*nnz; 1679 PetscScalar *values=0; 1680 PetscBool roworiented = ((Mat_SeqSBAIJ*)B->data)->roworiented; 1681 PetscErrorCode ierr; 1682 PetscFunctionBegin; 1683 if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs); 1684 ierr = PetscLayoutSetBlockSize(B->rmap,bs);CHKERRQ(ierr); 1685 ierr = PetscLayoutSetBlockSize(B->cmap,bs);CHKERRQ(ierr); 1686 ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 1687 ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 1688 ierr = PetscLayoutGetBlockSize(B->rmap,&bs);CHKERRQ(ierr); 1689 m = B->rmap->n/bs; 1690 1691 if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]); 1692 ierr = PetscMalloc1(m+1,&nnz);CHKERRQ(ierr); 1693 for (i=0; i<m; i++) { 1694 nz = ii[i+1] - ii[i]; 1695 if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D has a negative number of columns %D",i,nz); 1696 nz_max = PetscMax(nz_max,nz); 1697 nnz[i] = nz; 1698 } 1699 ierr = MatSeqSBAIJSetPreallocation(B,bs,0,nnz);CHKERRQ(ierr); 1700 ierr = PetscFree(nnz);CHKERRQ(ierr); 1701 1702 values = (PetscScalar*)V; 1703 if (!values) { 1704 ierr = PetscCalloc1(bs*bs*nz_max,&values);CHKERRQ(ierr); 1705 } 1706 for (i=0; i<m; i++) { 1707 PetscInt ncols = ii[i+1] - ii[i]; 1708 const PetscInt *icols = jj + ii[i]; 1709 if (!roworiented || bs == 1) { 1710 const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0); 1711 ierr = MatSetValuesBlocked_SeqSBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES);CHKERRQ(ierr); 1712 } else { 1713 for (j=0; j<ncols; j++) { 1714 const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0); 1715 ierr = MatSetValuesBlocked_SeqSBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES);CHKERRQ(ierr); 1716 } 1717 } 1718 } 1719 if (!V) { ierr = PetscFree(values);CHKERRQ(ierr); } 1720 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1721 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1722 ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1723 PetscFunctionReturn(0); 1724 } 1725 1726 /* 1727 This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization 1728 */ 1729 PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural) 1730 { 1731 PetscErrorCode ierr; 1732 PetscBool flg = PETSC_FALSE; 1733 PetscInt bs = B->rmap->bs; 1734 1735 PetscFunctionBegin; 1736 ierr = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);CHKERRQ(ierr); 1737 if (flg) bs = 8; 1738 1739 if (!natural) { 1740 switch (bs) { 1741 case 1: 1742 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace; 1743 break; 1744 case 2: 1745 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 1746 break; 1747 case 3: 1748 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 1749 break; 1750 case 4: 1751 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 1752 break; 1753 case 5: 1754 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 1755 break; 1756 case 6: 1757 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 1758 break; 1759 case 7: 1760 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 1761 break; 1762 default: 1763 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N; 1764 break; 1765 } 1766 } else { 1767 switch (bs) { 1768 case 1: 1769 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace; 1770 break; 1771 case 2: 1772 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 1773 break; 1774 case 3: 1775 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 1776 break; 1777 case 4: 1778 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 1779 break; 1780 case 5: 1781 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 1782 break; 1783 case 6: 1784 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 1785 break; 1786 case 7: 1787 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 1788 break; 1789 default: 1790 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 1791 break; 1792 } 1793 } 1794 PetscFunctionReturn(0); 1795 } 1796 1797 PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1798 PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1799 1800 PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 1801 { 1802 PetscInt n = A->rmap->n; 1803 PetscErrorCode ierr; 1804 1805 PetscFunctionBegin; 1806 #if defined(PETSC_USE_COMPLEX) 1807 if (A->hermitian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported"); 1808 #endif 1809 ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr); 1810 ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 1811 if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 1812 ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 1813 ierr = MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 1814 1815 (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; 1816 (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; 1817 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 1818 1819 (*B)->factortype = ftype; 1820 ierr = PetscFree((*B)->solvertype);CHKERRQ(ierr); 1821 ierr = PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);CHKERRQ(ierr); 1822 PetscFunctionReturn(0); 1823 } 1824 1825 /*@C 1826 MatSeqSBAIJGetArray - gives access to the array where the data for a MATSEQSBAIJ matrix is stored 1827 1828 Not Collective 1829 1830 Input Parameter: 1831 . mat - a MATSEQSBAIJ matrix 1832 1833 Output Parameter: 1834 . array - pointer to the data 1835 1836 Level: intermediate 1837 1838 .seealso: MatSeqSBAIJRestoreArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray() 1839 @*/ 1840 PetscErrorCode MatSeqSBAIJGetArray(Mat A,PetscScalar **array) 1841 { 1842 PetscErrorCode ierr; 1843 1844 PetscFunctionBegin; 1845 ierr = PetscUseMethod(A,"MatSeqSBAIJGetArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 1846 PetscFunctionReturn(0); 1847 } 1848 1849 /*@C 1850 MatSeqSBAIJRestoreArray - returns access to the array where the data for a MATSEQSBAIJ matrix is stored obtained by MatSeqSBAIJGetArray() 1851 1852 Not Collective 1853 1854 Input Parameters: 1855 . mat - a MATSEQSBAIJ matrix 1856 . array - pointer to the data 1857 1858 Level: intermediate 1859 1860 .seealso: MatSeqSBAIJGetArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray() 1861 @*/ 1862 PetscErrorCode MatSeqSBAIJRestoreArray(Mat A,PetscScalar **array) 1863 { 1864 PetscErrorCode ierr; 1865 1866 PetscFunctionBegin; 1867 ierr = PetscUseMethod(A,"MatSeqSBAIJRestoreArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 1868 PetscFunctionReturn(0); 1869 } 1870 1871 /*MC 1872 MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 1873 based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 1874 1875 For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you 1876 can call MatSetOption(Mat, MAT_HERMITIAN); after MatAssemblyEnd() 1877 1878 Options Database Keys: 1879 . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 1880 1881 Notes: 1882 By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not 1883 stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use 1884 the options database -mat_ignore_lower_triangular false it will generate an error if you try to set a value in the lower triangular portion. 1885 1886 1887 Level: beginner 1888 1889 .seealso: MatCreateSeqSBAIJ 1890 M*/ 1891 PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B) 1892 { 1893 Mat_SeqSBAIJ *b; 1894 PetscErrorCode ierr; 1895 PetscMPIInt size; 1896 PetscBool no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE; 1897 1898 PetscFunctionBegin; 1899 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);CHKERRQ(ierr); 1900 if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1901 1902 ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 1903 B->data = (void*)b; 1904 ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1905 1906 B->ops->destroy = MatDestroy_SeqSBAIJ; 1907 B->ops->view = MatView_SeqSBAIJ; 1908 b->row = 0; 1909 b->icol = 0; 1910 b->reallocs = 0; 1911 b->saved_values = 0; 1912 b->inode.limit = 5; 1913 b->inode.max_limit = 5; 1914 1915 b->roworiented = PETSC_TRUE; 1916 b->nonew = 0; 1917 b->diag = 0; 1918 b->solve_work = 0; 1919 b->mult_work = 0; 1920 B->spptr = 0; 1921 B->info.nz_unneeded = (PetscReal)b->maxnz*b->bs2; 1922 b->keepnonzeropattern = PETSC_FALSE; 1923 1924 b->inew = 0; 1925 b->jnew = 0; 1926 b->anew = 0; 1927 b->a2anew = 0; 1928 b->permute = PETSC_FALSE; 1929 1930 b->ignore_ltriangular = PETSC_TRUE; 1931 1932 ierr = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,NULL);CHKERRQ(ierr); 1933 1934 b->getrow_utriangular = PETSC_FALSE; 1935 1936 ierr = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,NULL);CHKERRQ(ierr); 1937 1938 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJGetArray_C",MatSeqSBAIJGetArray_SeqSBAIJ);CHKERRQ(ierr); 1939 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJRestoreArray_C",MatSeqSBAIJRestoreArray_SeqSBAIJ);CHKERRQ(ierr); 1940 ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1941 ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1942 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1943 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1944 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1945 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 1946 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);CHKERRQ(ierr); 1947 #if defined(PETSC_HAVE_ELEMENTAL) 1948 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_elemental_C",MatConvert_SeqSBAIJ_Elemental);CHKERRQ(ierr); 1949 #endif 1950 1951 B->symmetric = PETSC_TRUE; 1952 B->structurally_symmetric = PETSC_TRUE; 1953 B->symmetric_set = PETSC_TRUE; 1954 B->structurally_symmetric_set = PETSC_TRUE; 1955 B->symmetric_eternal = PETSC_TRUE; 1956 1957 B->hermitian = PETSC_FALSE; 1958 B->hermitian_set = PETSC_FALSE; 1959 1960 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr); 1961 1962 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");CHKERRQ(ierr); 1963 ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);CHKERRQ(ierr); 1964 if (no_unroll) { 1965 ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr); 1966 } 1967 ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);CHKERRQ(ierr); 1968 if (no_inode) { 1969 ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr); 1970 } 1971 ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);CHKERRQ(ierr); 1972 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1973 b->inode.use = (PetscBool)(!(no_unroll || no_inode)); 1974 if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit; 1975 PetscFunctionReturn(0); 1976 } 1977 1978 /*@C 1979 MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1980 compressed row) format. For good matrix assembly performance the 1981 user should preallocate the matrix storage by setting the parameter nz 1982 (or the array nnz). By setting these parameters accurately, performance 1983 during matrix assembly can be increased by more than a factor of 50. 1984 1985 Collective on Mat 1986 1987 Input Parameters: 1988 + B - the symmetric matrix 1989 . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row 1990 blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs() 1991 . nz - number of block nonzeros per block row (same for all rows) 1992 - nnz - array containing the number of block nonzeros in the upper triangular plus 1993 diagonal portion of each block (possibly different for each block row) or NULL 1994 1995 Options Database Keys: 1996 . -mat_no_unroll - uses code that does not unroll the loops in the 1997 block calculations (much slower) 1998 . -mat_block_size - size of the blocks to use (only works if a negative bs is passed in 1999 2000 Level: intermediate 2001 2002 Notes: 2003 Specify the preallocated storage with either nz or nnz (not both). 2004 Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 2005 allocation. See Users-Manual: ch_mat for details. 2006 2007 You can call MatGetInfo() to get information on how effective the preallocation was; 2008 for example the fields mallocs,nz_allocated,nz_used,nz_unneeded; 2009 You can also run with the option -info and look for messages with the string 2010 malloc in them to see if additional memory allocation was needed. 2011 2012 If the nnz parameter is given then the nz parameter is ignored 2013 2014 2015 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ() 2016 @*/ 2017 PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 2018 { 2019 PetscErrorCode ierr; 2020 2021 PetscFunctionBegin; 2022 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 2023 PetscValidType(B,1); 2024 PetscValidLogicalCollectiveInt(B,bs,2); 2025 ierr = PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr); 2026 PetscFunctionReturn(0); 2027 } 2028 2029 /*@C 2030 MatSeqSBAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in symmetric block AIJ format. 2031 2032 Input Parameters: 2033 + B - the matrix 2034 . bs - size of block, the blocks are ALWAYS square. 2035 . i - the indices into j for the start of each local row (starts with zero) 2036 . j - the column indices for each local row (starts with zero) these must be sorted for each row 2037 - v - optional values in the matrix 2038 2039 Level: developer 2040 2041 Notes: 2042 The order of the entries in values is specified by the MatOption MAT_ROW_ORIENTED. For example, C programs 2043 may want to use the default MAT_ROW_ORIENTED=PETSC_TRUE and use an array v[nnz][bs][bs] where the second index is 2044 over rows within a block and the last index is over columns within a block row. Fortran programs will likely set 2045 MAT_ROW_ORIENTED=PETSC_FALSE and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a 2046 block column and the second index is over columns within a block. 2047 2048 .keywords: matrix, block, aij, compressed row, sparse 2049 2050 .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ 2051 @*/ 2052 PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[]) 2053 { 2054 PetscErrorCode ierr; 2055 2056 PetscFunctionBegin; 2057 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 2058 PetscValidType(B,1); 2059 PetscValidLogicalCollectiveInt(B,bs,2); 2060 ierr = PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));CHKERRQ(ierr); 2061 PetscFunctionReturn(0); 2062 } 2063 2064 /*@C 2065 MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 2066 compressed row) format. For good matrix assembly performance the 2067 user should preallocate the matrix storage by setting the parameter nz 2068 (or the array nnz). By setting these parameters accurately, performance 2069 during matrix assembly can be increased by more than a factor of 50. 2070 2071 Collective on MPI_Comm 2072 2073 Input Parameters: 2074 + comm - MPI communicator, set to PETSC_COMM_SELF 2075 . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row 2076 blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs() 2077 . m - number of rows, or number of columns 2078 . nz - number of block nonzeros per block row (same for all rows) 2079 - nnz - array containing the number of block nonzeros in the upper triangular plus 2080 diagonal portion of each block (possibly different for each block row) or NULL 2081 2082 Output Parameter: 2083 . A - the symmetric matrix 2084 2085 Options Database Keys: 2086 . -mat_no_unroll - uses code that does not unroll the loops in the 2087 block calculations (much slower) 2088 . -mat_block_size - size of the blocks to use 2089 2090 Level: intermediate 2091 2092 It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 2093 MatXXXXSetPreallocation() paradgm instead of this routine directly. 2094 [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 2095 2096 Notes: 2097 The number of rows and columns must be divisible by blocksize. 2098 This matrix type does not support complex Hermitian operation. 2099 2100 Specify the preallocated storage with either nz or nnz (not both). 2101 Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 2102 allocation. See Users-Manual: ch_mat for details. 2103 2104 If the nnz parameter is given then the nz parameter is ignored 2105 2106 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ() 2107 @*/ 2108 PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 2109 { 2110 PetscErrorCode ierr; 2111 2112 PetscFunctionBegin; 2113 ierr = MatCreate(comm,A);CHKERRQ(ierr); 2114 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 2115 ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 2116 ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 2117 PetscFunctionReturn(0); 2118 } 2119 2120 PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 2121 { 2122 Mat C; 2123 Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 2124 PetscErrorCode ierr; 2125 PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 2126 2127 PetscFunctionBegin; 2128 if (a->i[mbs] != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix"); 2129 2130 *B = 0; 2131 ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 2132 ierr = MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);CHKERRQ(ierr); 2133 ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 2134 c = (Mat_SeqSBAIJ*)C->data; 2135 2136 C->preallocated = PETSC_TRUE; 2137 C->factortype = A->factortype; 2138 c->row = 0; 2139 c->icol = 0; 2140 c->saved_values = 0; 2141 c->keepnonzeropattern = a->keepnonzeropattern; 2142 C->assembled = PETSC_TRUE; 2143 2144 ierr = PetscLayoutReference(A->rmap,&C->rmap);CHKERRQ(ierr); 2145 ierr = PetscLayoutReference(A->cmap,&C->cmap);CHKERRQ(ierr); 2146 c->bs2 = a->bs2; 2147 c->mbs = a->mbs; 2148 c->nbs = a->nbs; 2149 2150 if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 2151 c->imax = a->imax; 2152 c->ilen = a->ilen; 2153 c->free_imax_ilen = PETSC_FALSE; 2154 } else { 2155 ierr = PetscMalloc2((mbs+1),&c->imax,(mbs+1),&c->ilen);CHKERRQ(ierr); 2156 ierr = PetscLogObjectMemory((PetscObject)C,2*(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 2157 for (i=0; i<mbs; i++) { 2158 c->imax[i] = a->imax[i]; 2159 c->ilen[i] = a->ilen[i]; 2160 } 2161 c->free_imax_ilen = PETSC_TRUE; 2162 } 2163 2164 /* allocate the matrix space */ 2165 if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 2166 ierr = PetscMalloc1(bs2*nz,&c->a);CHKERRQ(ierr); 2167 ierr = PetscLogObjectMemory((PetscObject)C,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 2168 c->i = a->i; 2169 c->j = a->j; 2170 c->singlemalloc = PETSC_FALSE; 2171 c->free_a = PETSC_TRUE; 2172 c->free_ij = PETSC_FALSE; 2173 c->parent = A; 2174 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 2175 ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 2176 ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 2177 } else { 2178 ierr = PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i);CHKERRQ(ierr); 2179 ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 2180 ierr = PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 2181 c->singlemalloc = PETSC_TRUE; 2182 c->free_a = PETSC_TRUE; 2183 c->free_ij = PETSC_TRUE; 2184 } 2185 if (mbs > 0) { 2186 if (cpvalues != MAT_SHARE_NONZERO_PATTERN) { 2187 ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 2188 } 2189 if (cpvalues == MAT_COPY_VALUES) { 2190 ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 2191 } else { 2192 ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 2193 } 2194 if (a->jshort) { 2195 /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */ 2196 /* if the parent matrix is reassembled, this child matrix will never notice */ 2197 ierr = PetscMalloc1(nz,&c->jshort);CHKERRQ(ierr); 2198 ierr = PetscLogObjectMemory((PetscObject)C,nz*sizeof(unsigned short));CHKERRQ(ierr); 2199 ierr = PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));CHKERRQ(ierr); 2200 2201 c->free_jshort = PETSC_TRUE; 2202 } 2203 } 2204 2205 c->roworiented = a->roworiented; 2206 c->nonew = a->nonew; 2207 2208 if (a->diag) { 2209 if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 2210 c->diag = a->diag; 2211 c->free_diag = PETSC_FALSE; 2212 } else { 2213 ierr = PetscMalloc1(mbs,&c->diag);CHKERRQ(ierr); 2214 ierr = PetscLogObjectMemory((PetscObject)C,mbs*sizeof(PetscInt));CHKERRQ(ierr); 2215 for (i=0; i<mbs; i++) c->diag[i] = a->diag[i]; 2216 c->free_diag = PETSC_TRUE; 2217 } 2218 } 2219 c->nz = a->nz; 2220 c->maxnz = a->nz; /* Since we allocate exactly the right amount */ 2221 c->solve_work = 0; 2222 c->mult_work = 0; 2223 2224 *B = C; 2225 ierr = PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);CHKERRQ(ierr); 2226 PetscFunctionReturn(0); 2227 } 2228 2229 PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer) 2230 { 2231 Mat_SeqSBAIJ *a; 2232 PetscErrorCode ierr; 2233 int fd; 2234 PetscMPIInt size; 2235 PetscInt i,nz,header[4],*rowlengths=0,M,N,bs = newmat->rmap->bs; 2236 PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 2237 PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols; 2238 PetscInt *masked,nmask,tmp,bs2,ishift; 2239 PetscScalar *aa; 2240 MPI_Comm comm; 2241 PetscBool isbinary; 2242 2243 PetscFunctionBegin; 2244 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 2245 if (!isbinary) SETERRQ2(PetscObjectComm((PetscObject)newmat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newmat)->type_name); 2246 2247 /* force binary viewer to load .info file if it has not yet done so */ 2248 ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 2249 ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 2250 ierr = PetscOptionsGetInt(((PetscObject)newmat)->options,((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);CHKERRQ(ierr); 2251 if (bs < 0) bs = 1; 2252 bs2 = bs*bs; 2253 2254 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2255 if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 2256 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 2257 ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 2258 if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 2259 M = header[1]; N = header[2]; nz = header[3]; 2260 2261 if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 2262 2263 if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices"); 2264 2265 /* 2266 This code adds extra rows to make sure the number of rows is 2267 divisible by the blocksize 2268 */ 2269 mbs = M/bs; 2270 extra_rows = bs - M + bs*(mbs); 2271 if (extra_rows == bs) extra_rows = 0; 2272 else mbs++; 2273 if (extra_rows) { 2274 ierr = PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 2275 } 2276 2277 /* Set global sizes if not already set */ 2278 if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 2279 ierr = MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 2280 } else { /* Check if the matrix global sizes are correct */ 2281 ierr = MatGetSize(newmat,&rows,&cols);CHKERRQ(ierr); 2282 if (M != rows || N != cols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,rows,cols); 2283 } 2284 2285 /* read in row lengths */ 2286 ierr = PetscMalloc1(M+extra_rows,&rowlengths);CHKERRQ(ierr); 2287 ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 2288 for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 2289 2290 /* read in column indices */ 2291 ierr = PetscMalloc1(nz+extra_rows,&jj);CHKERRQ(ierr); 2292 ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 2293 for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 2294 2295 /* loop over row lengths determining block row lengths */ 2296 ierr = PetscCalloc1(mbs,&s_browlengths);CHKERRQ(ierr); 2297 ierr = PetscMalloc2(mbs,&mask,mbs,&masked);CHKERRQ(ierr); 2298 ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 2299 rowcount = 0; 2300 nzcount = 0; 2301 for (i=0; i<mbs; i++) { 2302 nmask = 0; 2303 for (j=0; j<bs; j++) { 2304 kmax = rowlengths[rowcount]; 2305 for (k=0; k<kmax; k++) { 2306 tmp = jj[nzcount++]/bs; /* block col. index */ 2307 if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 2308 } 2309 rowcount++; 2310 } 2311 s_browlengths[i] += nmask; 2312 2313 /* zero out the mask elements we set */ 2314 for (j=0; j<nmask; j++) mask[masked[j]] = 0; 2315 } 2316 2317 /* Do preallocation */ 2318 ierr = MatSeqSBAIJSetPreallocation(newmat,bs,0,s_browlengths);CHKERRQ(ierr); 2319 a = (Mat_SeqSBAIJ*)newmat->data; 2320 2321 /* set matrix "i" values */ 2322 a->i[0] = 0; 2323 for (i=1; i<= mbs; i++) { 2324 a->i[i] = a->i[i-1] + s_browlengths[i-1]; 2325 a->ilen[i-1] = s_browlengths[i-1]; 2326 } 2327 a->nz = a->i[mbs]; 2328 2329 /* read in nonzero values */ 2330 ierr = PetscMalloc1(nz+extra_rows,&aa);CHKERRQ(ierr); 2331 ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 2332 for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 2333 2334 /* set "a" and "j" values into matrix */ 2335 nzcount = 0; jcount = 0; 2336 for (i=0; i<mbs; i++) { 2337 nzcountb = nzcount; 2338 nmask = 0; 2339 for (j=0; j<bs; j++) { 2340 kmax = rowlengths[i*bs+j]; 2341 for (k=0; k<kmax; k++) { 2342 tmp = jj[nzcount++]/bs; /* block col. index */ 2343 if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 2344 } 2345 } 2346 /* sort the masked values */ 2347 ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 2348 2349 /* set "j" values into matrix */ 2350 maskcount = 1; 2351 for (j=0; j<nmask; j++) { 2352 a->j[jcount++] = masked[j]; 2353 mask[masked[j]] = maskcount++; 2354 } 2355 2356 /* set "a" values into matrix */ 2357 ishift = bs2*a->i[i]; 2358 for (j=0; j<bs; j++) { 2359 kmax = rowlengths[i*bs+j]; 2360 for (k=0; k<kmax; k++) { 2361 tmp = jj[nzcountb]/bs; /* block col. index */ 2362 if (tmp >= i) { 2363 block = mask[tmp] - 1; 2364 point = jj[nzcountb] - bs*tmp; 2365 idx = ishift + bs2*block + j + bs*point; 2366 a->a[idx] = aa[nzcountb]; 2367 } 2368 nzcountb++; 2369 } 2370 } 2371 /* zero out the mask elements we set */ 2372 for (j=0; j<nmask; j++) mask[masked[j]] = 0; 2373 } 2374 if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 2375 2376 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2377 ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 2378 ierr = PetscFree(aa);CHKERRQ(ierr); 2379 ierr = PetscFree(jj);CHKERRQ(ierr); 2380 ierr = PetscFree2(mask,masked);CHKERRQ(ierr); 2381 2382 ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2383 ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2384 PetscFunctionReturn(0); 2385 } 2386 2387 /*@ 2388 MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2389 (upper triangular entries in CSR format) provided by the user. 2390 2391 Collective on MPI_Comm 2392 2393 Input Parameters: 2394 + comm - must be an MPI communicator of size 1 2395 . bs - size of block 2396 . m - number of rows 2397 . n - number of columns 2398 . i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that row block row of the matrix 2399 . j - column indices 2400 - a - matrix values 2401 2402 Output Parameter: 2403 . mat - the matrix 2404 2405 Level: advanced 2406 2407 Notes: 2408 The i, j, and a arrays are not copied by this routine, the user must free these arrays 2409 once the matrix is destroyed 2410 2411 You cannot set new nonzero locations into this matrix, that will generate an error. 2412 2413 The i and j indices are 0 based 2414 2415 When block size is greater than 1 the matrix values must be stored using the SBAIJ storage format (see the SBAIJ code to determine this). For block size of 1 2416 it is the regular CSR format excluding the lower triangular elements. 2417 2418 .seealso: MatCreate(), MatCreateSBAIJ(), MatCreateSeqSBAIJ() 2419 2420 @*/ 2421 PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],Mat *mat) 2422 { 2423 PetscErrorCode ierr; 2424 PetscInt ii; 2425 Mat_SeqSBAIJ *sbaij; 2426 2427 PetscFunctionBegin; 2428 if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2429 if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2430 2431 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2432 ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2433 ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2434 ierr = MatSeqSBAIJSetPreallocation(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2435 sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2436 ierr = PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);CHKERRQ(ierr); 2437 ierr = PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));CHKERRQ(ierr); 2438 2439 sbaij->i = i; 2440 sbaij->j = j; 2441 sbaij->a = a; 2442 2443 sbaij->singlemalloc = PETSC_FALSE; 2444 sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2445 sbaij->free_a = PETSC_FALSE; 2446 sbaij->free_ij = PETSC_FALSE; 2447 sbaij->free_imax_ilen = PETSC_TRUE; 2448 2449 for (ii=0; ii<m; ii++) { 2450 sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2451 #if defined(PETSC_USE_DEBUG) 2452 if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]); 2453 #endif 2454 } 2455 #if defined(PETSC_USE_DEBUG) 2456 for (ii=0; ii<sbaij->i[m]; ii++) { 2457 if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2458 if (j[ii] > n - 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2459 } 2460 #endif 2461 2462 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2463 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2464 PetscFunctionReturn(0); 2465 } 2466 2467 PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat) 2468 { 2469 PetscErrorCode ierr; 2470 PetscMPIInt size; 2471 2472 PetscFunctionBegin; 2473 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2474 if (size == 1 && scall == MAT_REUSE_MATRIX) { 2475 ierr = MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2476 } else { 2477 ierr = MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(comm,inmat,n,scall,outmat);CHKERRQ(ierr); 2478 } 2479 PetscFunctionReturn(0); 2480 } 2481