1 #define PETSCMAT_DLL 2 3 /* 4 This provides a matrix that consists of Mats 5 */ 6 7 #include "private/matimpl.h" /*I "petscmat.h" I*/ 8 #include "../src/mat/impls/baij/seq/baij.h" /* use the common AIJ data-structure */ 9 #include "petscksp.h" 10 11 //#define CHUNKSIZE 15 12 13 typedef struct { 14 SEQAIJHEADER(Mat); 15 SEQBAIJHEADER; 16 Mat *diags; 17 18 Vec left,right,middle,workb; /* dummy vectors to perform local parts of product */ 19 } Mat_BlockMat; 20 21 #undef __FUNCT__ 22 #define __FUNCT__ "MatSOR_BlockMat_Symmetric" 23 PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 24 { 25 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 26 PetscScalar *x; 27 const Mat *v = a->a; 28 const PetscScalar *b; 29 PetscErrorCode ierr; 30 PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 31 const PetscInt *idx; 32 IS row,col; 33 MatFactorInfo info; 34 Vec left = a->left,right = a->right, middle = a->middle; 35 Mat *diag; 36 37 PetscFunctionBegin; 38 its = its*lits; 39 if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 40 if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 41 if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 42 if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 43 if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) 44 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep"); 45 46 if (!a->diags) { 47 ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr); 48 ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 49 for (i=0; i<mbs; i++) { 50 ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 51 ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr); 52 ierr = MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr); 53 ierr = ISDestroy(row);CHKERRQ(ierr); 54 ierr = ISDestroy(col);CHKERRQ(ierr); 55 } 56 ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr); 57 } 58 diag = a->diags; 59 60 ierr = VecSet(xx,0.0);CHKERRQ(ierr); 61 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 62 /* copy right hand side because it must be modified during iteration */ 63 ierr = VecCopy(bb,a->workb);CHKERRQ(ierr); 64 ierr = VecGetArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr); 65 66 /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 67 while (its--) { 68 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 69 70 for (i=0; i<mbs; i++) { 71 n = a->i[i+1] - a->i[i] - 1; 72 idx = a->j + a->i[i] + 1; 73 v = a->a + a->i[i] + 1; 74 75 ierr = VecSet(left,0.0);CHKERRQ(ierr); 76 for (j=0; j<n; j++) { 77 ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 78 ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 79 ierr = VecResetArray(right);CHKERRQ(ierr); 80 } 81 ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 82 ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 83 ierr = VecResetArray(right);CHKERRQ(ierr); 84 85 ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 86 ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 87 88 /* now adjust right hand side, see MatSOR_SeqSBAIJ */ 89 for (j=0; j<n; j++) { 90 ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr); 91 ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr); 92 ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr); 93 ierr = VecResetArray(middle);CHKERRQ(ierr); 94 } 95 ierr = VecResetArray(right);CHKERRQ(ierr); 96 97 } 98 } 99 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 100 101 for (i=mbs-1; i>=0; i--) { 102 n = a->i[i+1] - a->i[i] - 1; 103 idx = a->j + a->i[i] + 1; 104 v = a->a + a->i[i] + 1; 105 106 ierr = VecSet(left,0.0);CHKERRQ(ierr); 107 for (j=0; j<n; j++) { 108 ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 109 ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 110 ierr = VecResetArray(right);CHKERRQ(ierr); 111 } 112 ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 113 ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 114 ierr = VecResetArray(right);CHKERRQ(ierr); 115 116 ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 117 ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 118 ierr = VecResetArray(right);CHKERRQ(ierr); 119 120 } 121 } 122 } 123 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 124 ierr = VecRestoreArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr); 125 PetscFunctionReturn(0); 126 } 127 128 #undef __FUNCT__ 129 #define __FUNCT__ "MatSOR_BlockMat" 130 PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 131 { 132 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 133 PetscScalar *x; 134 const Mat *v = a->a; 135 const PetscScalar *b; 136 PetscErrorCode ierr; 137 PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 138 const PetscInt *idx; 139 IS row,col; 140 MatFactorInfo info; 141 Vec left = a->left,right = a->right; 142 Mat *diag; 143 144 PetscFunctionBegin; 145 its = its*lits; 146 if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 147 if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 148 if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 149 if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 150 151 if (!a->diags) { 152 ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr); 153 ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 154 for (i=0; i<mbs; i++) { 155 ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 156 ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr); 157 ierr = MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr); 158 ierr = ISDestroy(row);CHKERRQ(ierr); 159 ierr = ISDestroy(col);CHKERRQ(ierr); 160 } 161 } 162 diag = a->diags; 163 164 ierr = VecSet(xx,0.0);CHKERRQ(ierr); 165 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 166 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 167 168 /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 169 while (its--) { 170 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 171 172 for (i=0; i<mbs; i++) { 173 n = a->i[i+1] - a->i[i]; 174 idx = a->j + a->i[i]; 175 v = a->a + a->i[i]; 176 177 ierr = VecSet(left,0.0);CHKERRQ(ierr); 178 for (j=0; j<n; j++) { 179 if (idx[j] != i) { 180 ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 181 ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 182 ierr = VecResetArray(right);CHKERRQ(ierr); 183 } 184 } 185 ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 186 ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 187 ierr = VecResetArray(right);CHKERRQ(ierr); 188 189 ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 190 ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 191 ierr = VecResetArray(right);CHKERRQ(ierr); 192 } 193 } 194 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 195 196 for (i=mbs-1; i>=0; i--) { 197 n = a->i[i+1] - a->i[i]; 198 idx = a->j + a->i[i]; 199 v = a->a + a->i[i]; 200 201 ierr = VecSet(left,0.0);CHKERRQ(ierr); 202 for (j=0; j<n; j++) { 203 if (idx[j] != i) { 204 ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 205 ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 206 ierr = VecResetArray(right);CHKERRQ(ierr); 207 } 208 } 209 ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 210 ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 211 ierr = VecResetArray(right);CHKERRQ(ierr); 212 213 ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 214 ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 215 ierr = VecResetArray(right);CHKERRQ(ierr); 216 217 } 218 } 219 } 220 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 221 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 222 PetscFunctionReturn(0); 223 } 224 225 /* 226 Same as MatSeqXAIJReallocateAIJ but uses Mat datatype instead of MatScalar 227 */ 228 #undef __FUNCT__ 229 #define __FUNCT__ "MatSeqBlockMatReallocateAIJ" 230 PETSC_STATIC_INLINE PetscErrorCode MatSeqBlockMatReallocateAIJ(Mat Amat,PetscInt AM,PetscInt BS2,PetscInt NROW,PetscInt ROW,PetscInt COL,PetscInt RMAX,Mat *AA,PetscInt *AI,PetscInt *AJ,PetscInt *RP,Mat *AP,PetscInt *AIMAX,PetscInt NONEW) 231 { 232 PetscErrorCode ierr; 233 if (NROW >= RMAX) { 234 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; 235 PetscInt CHUNKSIZE=15; 236 /* there is no extra room in row, therefore enlarge */ 237 PetscInt new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0,ii; 238 Mat *new_a; 239 240 if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); 241 /* malloc new storage space */ 242 ierr = PetscMalloc3(BS2*new_nz,Mat, &new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr); 243 /* copy over old data into new slots */ 244 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} 245 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} 246 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); 247 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); 248 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); 249 ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(Mat));CHKERRQ(ierr); 250 ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(Mat));CHKERRQ(ierr); 251 ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(Mat));CHKERRQ(ierr); 252 /* free up old matrix storage */ 253 ierr = MatSeqXAIJFreeAIJ(Amat,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); 254 AA = new_a; 255 Ain->a = (MatScalar*) new_a; 256 AI = Ain->i = new_i; AJ = Ain->j = new_j; 257 Ain->singlemalloc = PETSC_TRUE; 258 RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; 259 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; 260 Ain->maxnz += BS2*CHUNKSIZE; 261 Ain->reallocs++; 262 } 263 return(0); 264 } 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "MatSetValues_BlockMat" 268 PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 269 { 270 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 271 PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 272 PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 273 PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 274 PetscErrorCode ierr; 275 PetscInt ridx,cidx; 276 PetscTruth roworiented=a->roworiented; 277 MatScalar value; 278 Mat *ap,*aa = a->a; 279 280 PetscFunctionBegin; 281 if (v) PetscValidScalarPointer(v,6); 282 for (k=0; k<m; k++) { /* loop over added rows */ 283 row = im[k]; 284 brow = row/bs; 285 if (row < 0) continue; 286 #if defined(PETSC_USE_DEBUG) 287 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); 288 #endif 289 rp = aj + ai[brow]; 290 ap = aa + ai[brow]; 291 rmax = imax[brow]; 292 nrow = ailen[brow]; 293 low = 0; 294 high = nrow; 295 for (l=0; l<n; l++) { /* loop over added columns */ 296 if (in[l] < 0) continue; 297 #if defined(PETSC_USE_DEBUG) 298 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); 299 #endif 300 col = in[l]; bcol = col/bs; 301 if (A->symmetric && brow > bcol) continue; 302 ridx = row % bs; cidx = col % bs; 303 if (roworiented) { 304 value = v[l + k*n]; 305 } else { 306 value = v[k + l*m]; 307 } 308 if (col <= lastcol) low = 0; else high = nrow; 309 lastcol = col; 310 while (high-low > 7) { 311 t = (low+high)/2; 312 if (rp[t] > bcol) high = t; 313 else low = t; 314 } 315 for (i=low; i<high; i++) { 316 if (rp[i] > bcol) break; 317 if (rp[i] == bcol) { 318 goto noinsert1; 319 } 320 } 321 if (nonew == 1) goto noinsert1; 322 if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 323 ierr = MatSeqBlockMatReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew);CHKERRQ(ierr); 324 N = nrow++ - 1; high++; 325 /* shift up all the later entries in this row */ 326 for (ii=N; ii>=i; ii--) { 327 rp[ii+1] = rp[ii]; 328 ap[ii+1] = ap[ii]; 329 } 330 if (N>=i) ap[i] = 0; 331 rp[i] = bcol; 332 a->nz++; 333 noinsert1:; 334 if (!*(ap+i)) { 335 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 336 } 337 ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 338 low = i; 339 } 340 ailen[brow] = nrow; 341 } 342 A->same_nonzero = PETSC_FALSE; 343 PetscFunctionReturn(0); 344 } 345 346 #undef __FUNCT__ 347 #define __FUNCT__ "MatLoad_BlockMat" 348 PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, const MatType type,Mat *A) 349 { 350 PetscErrorCode ierr; 351 Mat tmpA; 352 PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 353 const PetscInt *cols; 354 const PetscScalar *values; 355 PetscTruth flg = PETSC_FALSE,notdone; 356 Mat_SeqAIJ *a; 357 Mat_BlockMat *amat; 358 359 PetscFunctionBegin; 360 ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr); 361 362 ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 363 ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr); 364 ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 365 ierr = PetscOptionsTruth("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 366 ierr = PetscOptionsEnd();CHKERRQ(ierr); 367 368 /* Determine number of nonzero blocks for each block row */ 369 a = (Mat_SeqAIJ*) tmpA->data; 370 mbs = m/bs; 371 ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt*,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr); 372 ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr); 373 374 for (i=0; i<mbs; i++) { 375 for (j=0; j<bs; j++) { 376 ii[j] = a->j + a->i[i*bs + j]; 377 ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 378 } 379 380 currentcol = -1; 381 notdone = PETSC_TRUE; 382 while (PETSC_TRUE) { 383 notdone = PETSC_FALSE; 384 nextcol = 1000000000; 385 for (j=0; j<bs; j++) { 386 while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 387 ii[j]++; 388 ilens[j]--; 389 } 390 if (ilens[j] > 0) { 391 notdone = PETSC_TRUE; 392 nextcol = PetscMin(nextcol,ii[j][0]/bs); 393 } 394 } 395 if (!notdone) break; 396 if (!flg || (nextcol >= i)) lens[i]++; 397 currentcol = nextcol; 398 } 399 } 400 401 ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,0,lens,A);CHKERRQ(ierr); 402 if (flg) { 403 ierr = MatSetOption(*A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 404 } 405 amat = (Mat_BlockMat*)(*A)->data; 406 407 /* preallocate the submatrices */ 408 ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr); 409 for (i=0; i<mbs; i++) { /* loops for block rows */ 410 for (j=0; j<bs; j++) { 411 ii[j] = a->j + a->i[i*bs + j]; 412 ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 413 } 414 415 currentcol = 1000000000; 416 for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 417 if (ilens[j] > 0) { 418 currentcol = PetscMin(currentcol,ii[j][0]/bs); 419 } 420 } 421 422 notdone = PETSC_TRUE; 423 while (PETSC_TRUE) { /* loops over blocks in block row */ 424 425 notdone = PETSC_FALSE; 426 nextcol = 1000000000; 427 ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr); 428 for (j=0; j<bs; j++) { /* loop over rows in block */ 429 while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 430 ii[j]++; 431 ilens[j]--; 432 llens[j]++; 433 } 434 if (ilens[j] > 0) { 435 notdone = PETSC_TRUE; 436 nextcol = PetscMin(nextcol,ii[j][0]/bs); 437 } 438 } 439 if (cnt >= amat->maxnz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt); 440 if (!flg || currentcol >= i) { 441 amat->j[cnt] = currentcol; 442 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr); 443 } 444 445 if (!notdone) break; 446 currentcol = nextcol; 447 } 448 amat->ilen[i] = lens[i]; 449 } 450 CHKMEMQ; 451 452 ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr); 453 ierr = PetscFree(llens);CHKERRQ(ierr); 454 455 /* copy over the matrix, one row at a time */ 456 for (i=0; i<m; i++) { 457 ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 458 ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 459 ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 460 } 461 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 462 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNCT__ 467 #define __FUNCT__ "MatView_BlockMat" 468 PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 469 { 470 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 471 PetscErrorCode ierr; 472 const char *name; 473 PetscViewerFormat format; 474 475 PetscFunctionBegin; 476 ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 477 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 478 if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 479 ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 480 if (A->symmetric) { 481 ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr); 482 } 483 } 484 PetscFunctionReturn(0); 485 } 486 487 #undef __FUNCT__ 488 #define __FUNCT__ "MatDestroy_BlockMat" 489 PetscErrorCode MatDestroy_BlockMat(Mat mat) 490 { 491 PetscErrorCode ierr; 492 Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 493 PetscInt i; 494 495 PetscFunctionBegin; 496 if (bmat->right) { 497 ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 498 } 499 if (bmat->left) { 500 ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 501 } 502 if (bmat->middle) { 503 ierr = VecDestroy(bmat->middle);CHKERRQ(ierr); 504 } 505 if (bmat->workb) { 506 ierr = VecDestroy(bmat->workb);CHKERRQ(ierr); 507 } 508 if (bmat->diags) { 509 for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) { 510 if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);} 511 } 512 } 513 if (bmat->a) { 514 for (i=0; i<bmat->nz; i++) { 515 if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);} 516 } 517 } 518 ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 519 ierr = PetscFree(bmat);CHKERRQ(ierr); 520 PetscFunctionReturn(0); 521 } 522 523 #undef __FUNCT__ 524 #define __FUNCT__ "MatMult_BlockMat" 525 PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 526 { 527 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 528 PetscErrorCode ierr; 529 PetscScalar *xx,*yy; 530 PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 531 Mat *aa; 532 533 PetscFunctionBegin; 534 CHKMEMQ; 535 /* 536 Standard CSR multiply except each entry is a Mat 537 */ 538 ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 539 540 ierr = VecSet(y,0.0);CHKERRQ(ierr); 541 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 542 aj = bmat->j; 543 aa = bmat->a; 544 ii = bmat->i; 545 for (i=0; i<m; i++) { 546 jrow = ii[i]; 547 ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 548 n = ii[i+1] - jrow; 549 for (j=0; j<n; j++) { 550 ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 551 ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 552 ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 553 jrow++; 554 } 555 ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 556 } 557 ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 558 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 559 CHKMEMQ; 560 PetscFunctionReturn(0); 561 } 562 563 #undef __FUNCT__ 564 #define __FUNCT__ "MatMult_BlockMat_Symmetric" 565 PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 566 { 567 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 568 PetscErrorCode ierr; 569 PetscScalar *xx,*yy; 570 PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 571 Mat *aa; 572 573 PetscFunctionBegin; 574 CHKMEMQ; 575 /* 576 Standard CSR multiply except each entry is a Mat 577 */ 578 ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 579 580 ierr = VecSet(y,0.0);CHKERRQ(ierr); 581 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 582 aj = bmat->j; 583 aa = bmat->a; 584 ii = bmat->i; 585 for (i=0; i<m; i++) { 586 jrow = ii[i]; 587 n = ii[i+1] - jrow; 588 ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 589 ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr); 590 /* if we ALWAYS required a diagonal entry then could remove this if test */ 591 if (aj[jrow] == i) { 592 ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 593 ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 594 ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 595 jrow++; 596 n--; 597 } 598 for (j=0; j<n; j++) { 599 ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); /* upper triangular part */ 600 ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 601 ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 602 603 ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr); /* lower triangular part */ 604 ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr); 605 ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 606 jrow++; 607 } 608 ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 609 ierr = VecResetArray(bmat->middle);CHKERRQ(ierr); 610 } 611 ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 612 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 613 CHKMEMQ; 614 PetscFunctionReturn(0); 615 } 616 617 #undef __FUNCT__ 618 #define __FUNCT__ "MatMultAdd_BlockMat" 619 PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 620 { 621 PetscFunctionBegin; 622 PetscFunctionReturn(0); 623 } 624 625 #undef __FUNCT__ 626 #define __FUNCT__ "MatMultTranspose_BlockMat" 627 PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 628 { 629 PetscFunctionBegin; 630 PetscFunctionReturn(0); 631 } 632 633 #undef __FUNCT__ 634 #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 635 PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 636 { 637 PetscFunctionBegin; 638 PetscFunctionReturn(0); 639 } 640 641 /* 642 Adds diagonal pointers to sparse matrix structure. 643 */ 644 #undef __FUNCT__ 645 #define __FUNCT__ "MatMarkDiagonal_BlockMat" 646 PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 647 { 648 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 649 PetscErrorCode ierr; 650 PetscInt i,j,mbs = A->rmap->n/A->rmap->bs; 651 652 PetscFunctionBegin; 653 if (!a->diag) { 654 ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 655 } 656 for (i=0; i<mbs; i++) { 657 a->diag[i] = a->i[i+1]; 658 for (j=a->i[i]; j<a->i[i+1]; j++) { 659 if (a->j[j] == i) { 660 a->diag[i] = j; 661 break; 662 } 663 } 664 } 665 PetscFunctionReturn(0); 666 } 667 668 #undef __FUNCT__ 669 #define __FUNCT__ "MatGetSubMatrix_BlockMat" 670 PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 671 { 672 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 673 Mat_SeqAIJ *c; 674 PetscErrorCode ierr; 675 PetscInt i,k,first,step,lensi,nrows,ncols; 676 PetscInt *j_new,*i_new,*aj = a->j,*ailen = a->ilen; 677 PetscScalar *a_new; 678 Mat C,*aa = a->a; 679 PetscTruth stride,equal; 680 681 PetscFunctionBegin; 682 ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 683 if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices"); 684 ierr = ISStride(iscol,&stride);CHKERRQ(ierr); 685 if (!stride) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices"); 686 ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 687 if (step != A->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block"); 688 689 ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 690 ncols = nrows; 691 692 /* create submatrix */ 693 if (scall == MAT_REUSE_MATRIX) { 694 PetscInt n_cols,n_rows; 695 C = *B; 696 ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 697 if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 698 ierr = MatZeroEntries(C);CHKERRQ(ierr); 699 } else { 700 ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 701 ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 702 if (A->symmetric) { 703 ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 704 } else { 705 ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 706 } 707 ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr); 708 ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr); 709 } 710 c = (Mat_SeqAIJ*)C->data; 711 712 /* loop over rows inserting into submatrix */ 713 a_new = c->a; 714 j_new = c->j; 715 i_new = c->i; 716 717 for (i=0; i<nrows; i++) { 718 lensi = ailen[i]; 719 for (k=0; k<lensi; k++) { 720 *j_new++ = *aj++; 721 ierr = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr); 722 } 723 i_new[i+1] = i_new[i] + lensi; 724 c->ilen[i] = lensi; 725 } 726 727 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 728 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 729 *B = C; 730 PetscFunctionReturn(0); 731 } 732 733 #undef __FUNCT__ 734 #define __FUNCT__ "MatAssemblyEnd_BlockMat" 735 PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 736 { 737 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 738 PetscErrorCode ierr; 739 PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 740 PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 741 Mat *aa = a->a,*ap; 742 743 PetscFunctionBegin; 744 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 745 746 if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 747 for (i=1; i<m; i++) { 748 /* move each row back by the amount of empty slots (fshift) before it*/ 749 fshift += imax[i-1] - ailen[i-1]; 750 rmax = PetscMax(rmax,ailen[i]); 751 if (fshift) { 752 ip = aj + ai[i] ; 753 ap = aa + ai[i] ; 754 N = ailen[i]; 755 for (j=0; j<N; j++) { 756 ip[j-fshift] = ip[j]; 757 ap[j-fshift] = ap[j]; 758 } 759 } 760 ai[i] = ai[i-1] + ailen[i-1]; 761 } 762 if (m) { 763 fshift += imax[m-1] - ailen[m-1]; 764 ai[m] = ai[m-1] + ailen[m-1]; 765 } 766 /* reset ilen and imax for each row */ 767 for (i=0; i<m; i++) { 768 ailen[i] = imax[i] = ai[i+1] - ai[i]; 769 } 770 a->nz = ai[m]; 771 for (i=0; i<a->nz; i++) { 772 #if defined(PETSC_USE_DEBUG) 773 if (!aa[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 774 #endif 775 ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 776 ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 777 } 778 CHKMEMQ; 779 ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap->n/A->cmap->bs,fshift,a->nz);CHKERRQ(ierr); 780 ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 781 ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 782 A->info.mallocs += a->reallocs; 783 a->reallocs = 0; 784 A->info.nz_unneeded = (double)fshift; 785 a->rmax = rmax; 786 787 A->same_nonzero = PETSC_TRUE; 788 ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 789 PetscFunctionReturn(0); 790 } 791 792 #undef __FUNCT__ 793 #define __FUNCT__ "MatSetOption_BlockMat" 794 PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscTruth flg) 795 { 796 PetscFunctionBegin; 797 if (opt == MAT_SYMMETRIC && flg) { 798 A->ops->sor = MatSOR_BlockMat_Symmetric; 799 A->ops->mult = MatMult_BlockMat_Symmetric; 800 } else { 801 PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]); 802 } 803 PetscFunctionReturn(0); 804 } 805 806 807 static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 808 0, 809 0, 810 MatMult_BlockMat, 811 /* 4*/ MatMultAdd_BlockMat, 812 MatMultTranspose_BlockMat, 813 MatMultTransposeAdd_BlockMat, 814 0, 815 0, 816 0, 817 /*10*/ 0, 818 0, 819 0, 820 MatSOR_BlockMat, 821 0, 822 /*15*/ 0, 823 0, 824 0, 825 0, 826 0, 827 /*20*/ 0, 828 MatAssemblyEnd_BlockMat, 829 MatSetOption_BlockMat, 830 0, 831 /*24*/ 0, 832 0, 833 0, 834 0, 835 0, 836 /*29*/ 0, 837 0, 838 0, 839 0, 840 0, 841 /*34*/ 0, 842 0, 843 0, 844 0, 845 0, 846 /*39*/ 0, 847 0, 848 0, 849 0, 850 0, 851 /*44*/ 0, 852 0, 853 0, 854 0, 855 0, 856 /*49*/ 0, 857 0, 858 0, 859 0, 860 0, 861 /*54*/ 0, 862 0, 863 0, 864 0, 865 0, 866 /*59*/ MatGetSubMatrix_BlockMat, 867 MatDestroy_BlockMat, 868 MatView_BlockMat, 869 0, 870 0, 871 /*64*/ 0, 872 0, 873 0, 874 0, 875 0, 876 /*69*/ 0, 877 0, 878 0, 879 0, 880 0, 881 /*74*/ 0, 882 0, 883 0, 884 0, 885 0, 886 /*79*/ 0, 887 0, 888 0, 889 0, 890 MatLoad_BlockMat, 891 /*84*/ 0, 892 0, 893 0, 894 0, 895 0, 896 /*89*/ 0, 897 0, 898 0, 899 0, 900 0, 901 /*94*/ 0, 902 0, 903 0, 904 0}; 905 906 #undef __FUNCT__ 907 #define __FUNCT__ "MatBlockMatSetPreallocation" 908 /*@C 909 MatBlockMatSetPreallocation - For good matrix assembly performance 910 the user should preallocate the matrix storage by setting the parameter nz 911 (or the array nnz). By setting these parameters accurately, performance 912 during matrix assembly can be increased by more than a factor of 50. 913 914 Collective on MPI_Comm 915 916 Input Parameters: 917 + B - The matrix 918 . bs - size of each block in matrix 919 . nz - number of nonzeros per block row (same for all rows) 920 - nnz - array containing the number of nonzeros in the various block rows 921 (possibly different for each row) or PETSC_NULL 922 923 Notes: 924 If nnz is given then nz is ignored 925 926 Specify the preallocated storage with either nz or nnz (not both). 927 Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 928 allocation. For large problems you MUST preallocate memory or you 929 will get TERRIBLE performance, see the users' manual chapter on matrices. 930 931 Level: intermediate 932 933 .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues() 934 935 @*/ 936 PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 937 { 938 PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 939 940 PetscFunctionBegin; 941 ierr = PetscObjectQueryFunction((PetscObject)B,"MatBlockMatSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 942 if (f) { 943 ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 944 } 945 PetscFunctionReturn(0); 946 } 947 948 EXTERN_C_BEGIN 949 #undef __FUNCT__ 950 #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat" 951 PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz) 952 { 953 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 954 PetscErrorCode ierr; 955 PetscInt i; 956 957 PetscFunctionBegin; 958 if (bs < 1) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Block size given %D must be great than zero",bs); 959 if (A->rmap->n % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap->n); 960 if (A->cmap->n % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap->n); 961 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 962 if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 963 if (nnz) { 964 for (i=0; i<A->rmap->n/bs; i++) { 965 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]); 966 if (nnz[i] > A->cmap->n/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],A->cmap->n/bs); 967 } 968 } 969 A->rmap->bs = A->cmap->bs = bs; 970 bmat->mbs = A->rmap->n/bs; 971 972 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr); 973 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->middle);CHKERRQ(ierr); 974 ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr); 975 976 if (!bmat->imax) { 977 ierr = PetscMalloc2(A->rmap->n,PetscInt,&bmat->imax,A->rmap->n,PetscInt,&bmat->ilen);CHKERRQ(ierr); 978 ierr = PetscLogObjectMemory(A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 979 } 980 if (nnz) { 981 nz = 0; 982 for (i=0; i<A->rmap->n/A->rmap->bs; i++) { 983 bmat->imax[i] = nnz[i]; 984 nz += nnz[i]; 985 } 986 } else { 987 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation"); 988 } 989 990 /* bmat->ilen will count nonzeros in each row so far. */ 991 for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;} 992 993 /* allocate the matrix space */ 994 ierr = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 995 ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap->n+1,PetscInt,&bmat->i);CHKERRQ(ierr); 996 ierr = PetscLogObjectMemory(A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 997 bmat->i[0] = 0; 998 for (i=1; i<bmat->mbs+1; i++) { 999 bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 1000 } 1001 bmat->singlemalloc = PETSC_TRUE; 1002 bmat->free_a = PETSC_TRUE; 1003 bmat->free_ij = PETSC_TRUE; 1004 1005 bmat->nz = 0; 1006 bmat->maxnz = nz; 1007 A->info.nz_unneeded = (double)bmat->maxnz; 1008 1009 PetscFunctionReturn(0); 1010 } 1011 EXTERN_C_END 1012 1013 /*MC 1014 MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 1015 consisting of (usually) sparse blocks. 1016 1017 Level: advanced 1018 1019 .seealso: MatCreateBlockMat() 1020 1021 M*/ 1022 1023 EXTERN_C_BEGIN 1024 #undef __FUNCT__ 1025 #define __FUNCT__ "MatCreate_BlockMat" 1026 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A) 1027 { 1028 Mat_BlockMat *b; 1029 PetscErrorCode ierr; 1030 1031 PetscFunctionBegin; 1032 ierr = PetscNewLog(A,Mat_BlockMat,&b);CHKERRQ(ierr); 1033 A->data = (void*)b; 1034 ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1035 1036 ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 1037 ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 1038 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1039 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1040 1041 A->assembled = PETSC_TRUE; 1042 A->preallocated = PETSC_FALSE; 1043 ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 1044 1045 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatBlockMatSetPreallocation_C", 1046 "MatBlockMatSetPreallocation_BlockMat", 1047 MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr); 1048 1049 PetscFunctionReturn(0); 1050 } 1051 EXTERN_C_END 1052 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "MatCreateBlockMat" 1055 /*@C 1056 MatCreateBlockMat - Creates a new matrix based sparse Mat storage 1057 1058 Collective on MPI_Comm 1059 1060 Input Parameters: 1061 + comm - MPI communicator 1062 . m - number of rows 1063 . n - number of columns 1064 . bs - size of each submatrix 1065 . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 1066 - nnz - expected number of nonzers per block row if known (use PETSC_NULL otherwise) 1067 1068 1069 Output Parameter: 1070 . A - the matrix 1071 1072 Level: intermediate 1073 1074 PETSc requires that matrices and vectors being used for certain 1075 operations are partitioned accordingly. For example, when 1076 creating a bmat matrix, A, that supports parallel matrix-vector 1077 products using MatMult(A,x,y) the user should set the number 1078 of local matrix rows to be the number of local elements of the 1079 corresponding result vector, y. Note that this is information is 1080 required for use of the matrix interface routines, even though 1081 the bmat matrix may not actually be physically partitioned. 1082 For example, 1083 1084 .keywords: matrix, bmat, create 1085 1086 .seealso: MATBLOCKMAT 1087 @*/ 1088 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A) 1089 { 1090 PetscErrorCode ierr; 1091 1092 PetscFunctionBegin; 1093 ierr = MatCreate(comm,A);CHKERRQ(ierr); 1094 ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1095 ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 1096 ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1097 PetscFunctionReturn(0); 1098 } 1099 1100 1101 1102