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