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