1421e10b8SBarry Smith 2421e10b8SBarry Smith /* 3421e10b8SBarry Smith This provides a matrix that consists of Mats 4421e10b8SBarry Smith */ 5421e10b8SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> /* use the common AIJ data-structure */ 8421e10b8SBarry Smith 9421e10b8SBarry Smith typedef struct { 10421e10b8SBarry Smith SEQAIJHEADER(Mat); 11421e10b8SBarry Smith SEQBAIJHEADER; 12ccb205f8SBarry Smith Mat *diags; 13421e10b8SBarry Smith 146e63c7a1SBarry Smith Vec left,right,middle,workb; /* dummy vectors to perform local parts of product */ 15421e10b8SBarry Smith } Mat_BlockMat; 16421e10b8SBarry Smith 177087cfbeSBarry Smith extern PetscErrorCode MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt*); 18a5973382SShri Abhyankar 19e0877f53SBarry Smith static PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 20290bbb0aSBarry Smith { 21290bbb0aSBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 22290bbb0aSBarry Smith PetscScalar *x; 23a2ea699eSBarry Smith const Mat *v; 24290bbb0aSBarry Smith const PetscScalar *b; 25290bbb0aSBarry Smith PetscErrorCode ierr; 26d0f46423SBarry Smith PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 27290bbb0aSBarry Smith const PetscInt *idx; 28290bbb0aSBarry Smith IS row,col; 29290bbb0aSBarry Smith MatFactorInfo info; 306e63c7a1SBarry Smith Vec left = a->left,right = a->right, middle = a->middle; 31290bbb0aSBarry Smith Mat *diag; 32290bbb0aSBarry Smith 33290bbb0aSBarry Smith PetscFunctionBegin; 34290bbb0aSBarry Smith its = its*lits; 35e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 36e32f2f54SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 37e32f2f54SBarry Smith if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 38e32f2f54SBarry Smith if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 3926fbe8dcSKarl Rupp if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) { 40e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep"); 4126fbe8dcSKarl Rupp } 42290bbb0aSBarry Smith 43290bbb0aSBarry Smith if (!a->diags) { 44785e854fSJed Brown ierr = PetscMalloc1(mbs,&a->diags);CHKERRQ(ierr); 45290bbb0aSBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 46290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 472692d6eeSBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 48719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr); 49719d5645SBarry Smith ierr = MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr); 506bf464f9SBarry Smith ierr = ISDestroy(&row);CHKERRQ(ierr); 516bf464f9SBarry Smith ierr = ISDestroy(&col);CHKERRQ(ierr); 52290bbb0aSBarry Smith } 536e63c7a1SBarry Smith ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr); 54290bbb0aSBarry Smith } 55290bbb0aSBarry Smith diag = a->diags; 56290bbb0aSBarry Smith 57290bbb0aSBarry Smith ierr = VecSet(xx,0.0);CHKERRQ(ierr); 58290bbb0aSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 59290bbb0aSBarry Smith /* copy right hand side because it must be modified during iteration */ 606e63c7a1SBarry Smith ierr = VecCopy(bb,a->workb);CHKERRQ(ierr); 613649974fSBarry Smith ierr = VecGetArrayRead(a->workb,&b);CHKERRQ(ierr); 62290bbb0aSBarry Smith 6341f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 64290bbb0aSBarry Smith while (its--) { 65290bbb0aSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 66290bbb0aSBarry Smith 67290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 686e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 696e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 706e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 71290bbb0aSBarry Smith 72290bbb0aSBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 73290bbb0aSBarry Smith for (j=0; j<n; j++) { 74290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 75290bbb0aSBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 76290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 77290bbb0aSBarry Smith } 78290bbb0aSBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 79290bbb0aSBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 80290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 81290bbb0aSBarry Smith 82290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 83290bbb0aSBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 846e63c7a1SBarry Smith 8541f059aeSBarry Smith /* now adjust right hand side, see MatSOR_SeqSBAIJ */ 866e63c7a1SBarry Smith for (j=0; j<n; j++) { 876e63c7a1SBarry Smith ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr); 886e63c7a1SBarry Smith ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr); 896e63c7a1SBarry Smith ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr); 906e63c7a1SBarry Smith ierr = VecResetArray(middle);CHKERRQ(ierr); 916e63c7a1SBarry Smith } 92290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 936e63c7a1SBarry Smith 94290bbb0aSBarry Smith } 95290bbb0aSBarry Smith } 96290bbb0aSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 97290bbb0aSBarry Smith 98290bbb0aSBarry Smith for (i=mbs-1; i>=0; i--) { 996e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 1006e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 1016e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 102290bbb0aSBarry Smith 103290bbb0aSBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 104290bbb0aSBarry Smith for (j=0; j<n; j++) { 105290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 106290bbb0aSBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 107290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 108290bbb0aSBarry Smith } 109290bbb0aSBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 110290bbb0aSBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 111290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 112290bbb0aSBarry Smith 113290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 114290bbb0aSBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 115290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 116290bbb0aSBarry Smith 117290bbb0aSBarry Smith } 118290bbb0aSBarry Smith } 119290bbb0aSBarry Smith } 120290bbb0aSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1213649974fSBarry Smith ierr = VecRestoreArrayRead(a->workb,&b);CHKERRQ(ierr); 122290bbb0aSBarry Smith PetscFunctionReturn(0); 123290bbb0aSBarry Smith } 124290bbb0aSBarry Smith 12581f0254dSBarry Smith static PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 126ccb205f8SBarry Smith { 127ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 128ccb205f8SBarry Smith PetscScalar *x; 129a2ea699eSBarry Smith const Mat *v; 130ccb205f8SBarry Smith const PetscScalar *b; 131ccb205f8SBarry Smith PetscErrorCode ierr; 132d0f46423SBarry Smith PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 133ccb205f8SBarry Smith const PetscInt *idx; 134ccb205f8SBarry Smith IS row,col; 135ccb205f8SBarry Smith MatFactorInfo info; 136ccb205f8SBarry Smith Vec left = a->left,right = a->right; 137ccb205f8SBarry Smith Mat *diag; 138ccb205f8SBarry Smith 139ccb205f8SBarry Smith PetscFunctionBegin; 140ccb205f8SBarry Smith its = its*lits; 141e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 142e32f2f54SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 143e32f2f54SBarry Smith if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 144e32f2f54SBarry Smith if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 145ccb205f8SBarry Smith 146ccb205f8SBarry Smith if (!a->diags) { 147785e854fSJed Brown ierr = PetscMalloc1(mbs,&a->diags);CHKERRQ(ierr); 148ccb205f8SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 149ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 1502692d6eeSBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 151719d5645SBarry Smith ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr); 152719d5645SBarry Smith ierr = MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr); 1536bf464f9SBarry Smith ierr = ISDestroy(&row);CHKERRQ(ierr); 1546bf464f9SBarry Smith ierr = ISDestroy(&col);CHKERRQ(ierr); 155ccb205f8SBarry Smith } 156ccb205f8SBarry Smith } 157ccb205f8SBarry Smith diag = a->diags; 158ccb205f8SBarry Smith 159ccb205f8SBarry Smith ierr = VecSet(xx,0.0);CHKERRQ(ierr); 160ccb205f8SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1613649974fSBarry Smith ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 162ccb205f8SBarry Smith 16341f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 164ccb205f8SBarry Smith while (its--) { 165ccb205f8SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 166ccb205f8SBarry Smith 167ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 168ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 169ccb205f8SBarry Smith idx = a->j + a->i[i]; 170ccb205f8SBarry Smith v = a->a + a->i[i]; 171ccb205f8SBarry Smith 172ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 173ccb205f8SBarry Smith for (j=0; j<n; j++) { 174ccb205f8SBarry Smith if (idx[j] != i) { 175ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 176ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 177ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 178ccb205f8SBarry Smith } 179ccb205f8SBarry Smith } 180ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 181ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 182ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 183ccb205f8SBarry Smith 184ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 185ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 186ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 187ccb205f8SBarry Smith } 188ccb205f8SBarry Smith } 189ccb205f8SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 190ccb205f8SBarry Smith 191ccb205f8SBarry Smith for (i=mbs-1; i>=0; i--) { 192ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 193ccb205f8SBarry Smith idx = a->j + a->i[i]; 194ccb205f8SBarry Smith v = a->a + a->i[i]; 195ccb205f8SBarry Smith 196ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 197ccb205f8SBarry Smith for (j=0; j<n; j++) { 198ccb205f8SBarry Smith if (idx[j] != i) { 199ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 200ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 201ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 202ccb205f8SBarry Smith } 203ccb205f8SBarry Smith } 204ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 205ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 206ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 207ccb205f8SBarry Smith 208ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 209ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 210ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 211ccb205f8SBarry Smith 212ccb205f8SBarry Smith } 213ccb205f8SBarry Smith } 214ccb205f8SBarry Smith } 215ccb205f8SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2163649974fSBarry Smith ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 217ccb205f8SBarry Smith PetscFunctionReturn(0); 218ccb205f8SBarry Smith } 219ccb205f8SBarry Smith 22081f0254dSBarry Smith static PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 221421e10b8SBarry Smith { 222421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 223421e10b8SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 224421e10b8SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 225d0f46423SBarry Smith PetscInt *aj =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 226421e10b8SBarry Smith PetscErrorCode ierr; 227421e10b8SBarry Smith PetscInt ridx,cidx; 228ace3abfcSBarry Smith PetscBool roworiented=a->roworiented; 229421e10b8SBarry Smith MatScalar value; 230421e10b8SBarry Smith Mat *ap,*aa = a->a; 231421e10b8SBarry Smith 232421e10b8SBarry Smith PetscFunctionBegin; 233421e10b8SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 234421e10b8SBarry Smith row = im[k]; 235421e10b8SBarry Smith brow = row/bs; 236421e10b8SBarry Smith if (row < 0) continue; 237421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 238e32f2f54SBarry Smith 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); 239421e10b8SBarry Smith #endif 240421e10b8SBarry Smith rp = aj + ai[brow]; 241421e10b8SBarry Smith ap = aa + ai[brow]; 242421e10b8SBarry Smith rmax = imax[brow]; 243421e10b8SBarry Smith nrow = ailen[brow]; 244421e10b8SBarry Smith low = 0; 245421e10b8SBarry Smith high = nrow; 246421e10b8SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 247421e10b8SBarry Smith if (in[l] < 0) continue; 248421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 249e32f2f54SBarry Smith 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); 250421e10b8SBarry Smith #endif 251421e10b8SBarry Smith col = in[l]; bcol = col/bs; 2526e63c7a1SBarry Smith if (A->symmetric && brow > bcol) continue; 253421e10b8SBarry Smith ridx = row % bs; cidx = col % bs; 2542205254eSKarl Rupp if (roworiented) value = v[l + k*n]; 2552205254eSKarl Rupp else value = v[k + l*m]; 2562205254eSKarl Rupp 2572205254eSKarl Rupp if (col <= lastcol) low = 0; 2582205254eSKarl Rupp else high = nrow; 259421e10b8SBarry Smith lastcol = col; 260421e10b8SBarry Smith while (high-low > 7) { 261421e10b8SBarry Smith t = (low+high)/2; 262421e10b8SBarry Smith if (rp[t] > bcol) high = t; 263421e10b8SBarry Smith else low = t; 264421e10b8SBarry Smith } 265421e10b8SBarry Smith for (i=low; i<high; i++) { 266421e10b8SBarry Smith if (rp[i] > bcol) break; 2672205254eSKarl Rupp if (rp[i] == bcol) goto noinsert1; 268421e10b8SBarry Smith } 269421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 270e32f2f54SBarry Smith if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 271fef13f97SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 272421e10b8SBarry Smith N = nrow++ - 1; high++; 273421e10b8SBarry Smith /* shift up all the later entries in this row */ 274421e10b8SBarry Smith for (ii=N; ii>=i; ii--) { 275421e10b8SBarry Smith rp[ii+1] = rp[ii]; 276421e10b8SBarry Smith ap[ii+1] = ap[ii]; 277421e10b8SBarry Smith } 278421e10b8SBarry Smith if (N>=i) ap[i] = 0; 279421e10b8SBarry Smith rp[i] = bcol; 280421e10b8SBarry Smith a->nz++; 281e56f5c9eSBarry Smith A->nonzerostate++; 282421e10b8SBarry Smith noinsert1:; 283421e10b8SBarry Smith if (!*(ap+i)) { 2846e63c7a1SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 285b0223207SBarry Smith } 286421e10b8SBarry Smith ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 287421e10b8SBarry Smith low = i; 288421e10b8SBarry Smith } 289421e10b8SBarry Smith ailen[brow] = nrow; 290421e10b8SBarry Smith } 291421e10b8SBarry Smith PetscFunctionReturn(0); 292421e10b8SBarry Smith } 293421e10b8SBarry Smith 29481f0254dSBarry Smith static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer) 295a5973382SShri Abhyankar { 296a5973382SShri Abhyankar PetscErrorCode ierr; 297a5973382SShri Abhyankar Mat tmpA; 298a5973382SShri Abhyankar PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 299a5973382SShri Abhyankar const PetscInt *cols; 300a5973382SShri Abhyankar const PetscScalar *values; 301ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE,notdone; 302a5973382SShri Abhyankar Mat_SeqAIJ *a; 303a5973382SShri Abhyankar Mat_BlockMat *amat; 304a5973382SShri Abhyankar 305a5973382SShri Abhyankar PetscFunctionBegin; 306c98fd787SBarry Smith /* force binary viewer to load .info file if it has not yet done so */ 307c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 308a5973382SShri Abhyankar ierr = MatCreate(PETSC_COMM_SELF,&tmpA);CHKERRQ(ierr); 309a5973382SShri Abhyankar ierr = MatSetType(tmpA,MATSEQAIJ);CHKERRQ(ierr); 310112444f4SShri Abhyankar ierr = MatLoad_SeqAIJ(tmpA,viewer);CHKERRQ(ierr); 311a5973382SShri Abhyankar 312a5973382SShri Abhyankar ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 3130298fd71SBarry Smith ierr = PetscOptionsBegin(PETSC_COMM_SELF,NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr); 3140298fd71SBarry Smith ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);CHKERRQ(ierr); 3150298fd71SBarry Smith ierr = PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,NULL);CHKERRQ(ierr); 316a5973382SShri Abhyankar ierr = PetscOptionsEnd();CHKERRQ(ierr); 317a5973382SShri Abhyankar 318a5973382SShri Abhyankar /* Determine number of nonzero blocks for each block row */ 319a5973382SShri Abhyankar a = (Mat_SeqAIJ*) tmpA->data; 320a5973382SShri Abhyankar mbs = m/bs; 321dcca6d9dSJed Brown ierr = PetscMalloc3(mbs,&lens,bs,&ii,bs,&ilens);CHKERRQ(ierr); 322a5973382SShri Abhyankar ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr); 323a5973382SShri Abhyankar 324a5973382SShri Abhyankar for (i=0; i<mbs; i++) { 325a5973382SShri Abhyankar for (j=0; j<bs; j++) { 326a5973382SShri Abhyankar ii[j] = a->j + a->i[i*bs + j]; 327a5973382SShri Abhyankar ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 328a5973382SShri Abhyankar } 329a5973382SShri Abhyankar 330a5973382SShri Abhyankar currentcol = -1; 331a5973382SShri Abhyankar while (PETSC_TRUE) { 332a5973382SShri Abhyankar notdone = PETSC_FALSE; 333a5973382SShri Abhyankar nextcol = 1000000000; 334a5973382SShri Abhyankar for (j=0; j<bs; j++) { 335a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 336a5973382SShri Abhyankar ii[j]++; 337a5973382SShri Abhyankar ilens[j]--; 338a5973382SShri Abhyankar } 339a5973382SShri Abhyankar if (ilens[j] > 0) { 340a5973382SShri Abhyankar notdone = PETSC_TRUE; 341a5973382SShri Abhyankar nextcol = PetscMin(nextcol,ii[j][0]/bs); 342a5973382SShri Abhyankar } 343a5973382SShri Abhyankar } 344a5973382SShri Abhyankar if (!notdone) break; 345a5973382SShri Abhyankar if (!flg || (nextcol >= i)) lens[i]++; 346a5973382SShri Abhyankar currentcol = nextcol; 347a5973382SShri Abhyankar } 348a5973382SShri Abhyankar } 349a5973382SShri Abhyankar 350a5973382SShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 351a5973382SShri Abhyankar ierr = MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 352a5973382SShri Abhyankar } 353a5973382SShri Abhyankar ierr = MatBlockMatSetPreallocation(newmat,bs,0,lens);CHKERRQ(ierr); 354a5973382SShri Abhyankar if (flg) { 355a5973382SShri Abhyankar ierr = MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 356a5973382SShri Abhyankar } 357a5973382SShri Abhyankar amat = (Mat_BlockMat*)(newmat)->data; 358a5973382SShri Abhyankar 359a5973382SShri Abhyankar /* preallocate the submatrices */ 360785e854fSJed Brown ierr = PetscMalloc1(bs,&llens);CHKERRQ(ierr); 361a5973382SShri Abhyankar for (i=0; i<mbs; i++) { /* loops for block rows */ 362a5973382SShri Abhyankar for (j=0; j<bs; j++) { 363a5973382SShri Abhyankar ii[j] = a->j + a->i[i*bs + j]; 364a5973382SShri Abhyankar ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 365a5973382SShri Abhyankar } 366a5973382SShri Abhyankar 367a5973382SShri Abhyankar currentcol = 1000000000; 368a5973382SShri Abhyankar for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 369a5973382SShri Abhyankar if (ilens[j] > 0) { 370a5973382SShri Abhyankar currentcol = PetscMin(currentcol,ii[j][0]/bs); 371a5973382SShri Abhyankar } 372a5973382SShri Abhyankar } 373a5973382SShri Abhyankar 374a5973382SShri Abhyankar while (PETSC_TRUE) { /* loops over blocks in block row */ 375a5973382SShri Abhyankar notdone = PETSC_FALSE; 376a5973382SShri Abhyankar nextcol = 1000000000; 377a5973382SShri Abhyankar ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr); 378a5973382SShri Abhyankar for (j=0; j<bs; j++) { /* loop over rows in block */ 379a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 380a5973382SShri Abhyankar ii[j]++; 381a5973382SShri Abhyankar ilens[j]--; 382a5973382SShri Abhyankar llens[j]++; 383a5973382SShri Abhyankar } 384a5973382SShri Abhyankar if (ilens[j] > 0) { 385a5973382SShri Abhyankar notdone = PETSC_TRUE; 386a5973382SShri Abhyankar nextcol = PetscMin(nextcol,ii[j][0]/bs); 387a5973382SShri Abhyankar } 388a5973382SShri Abhyankar } 389a5973382SShri Abhyankar if (cnt >= amat->maxnz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt); 390a5973382SShri Abhyankar if (!flg || currentcol >= i) { 391a5973382SShri Abhyankar amat->j[cnt] = currentcol; 392a5973382SShri Abhyankar ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr); 393a5973382SShri Abhyankar } 394a5973382SShri Abhyankar 395a5973382SShri Abhyankar if (!notdone) break; 396a5973382SShri Abhyankar currentcol = nextcol; 397a5973382SShri Abhyankar } 398a5973382SShri Abhyankar amat->ilen[i] = lens[i]; 399a5973382SShri Abhyankar } 400a5973382SShri Abhyankar 401a5973382SShri Abhyankar ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr); 402a5973382SShri Abhyankar ierr = PetscFree(llens);CHKERRQ(ierr); 403a5973382SShri Abhyankar 404a5973382SShri Abhyankar /* copy over the matrix, one row at a time */ 405a5973382SShri Abhyankar for (i=0; i<m; i++) { 406a5973382SShri Abhyankar ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 407a5973382SShri Abhyankar ierr = MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 408a5973382SShri Abhyankar ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 409a5973382SShri Abhyankar } 410a5973382SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 411a5973382SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 412a5973382SShri Abhyankar PetscFunctionReturn(0); 413a5973382SShri Abhyankar } 414a5973382SShri Abhyankar 41581f0254dSBarry Smith static PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 416ccb205f8SBarry Smith { 41764075487SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 418ccb205f8SBarry Smith PetscErrorCode ierr; 419ccb205f8SBarry Smith const char *name; 420ccb205f8SBarry Smith PetscViewerFormat format; 421ccb205f8SBarry Smith 422ccb205f8SBarry Smith PetscFunctionBegin; 423ccb205f8SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 424ccb205f8SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 425ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 426ccb205f8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 42764075487SBarry Smith if (A->symmetric) { 4288c1ad04dSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr); 42964075487SBarry Smith } 430ccb205f8SBarry Smith } 431ccb205f8SBarry Smith PetscFunctionReturn(0); 432ccb205f8SBarry Smith } 433421e10b8SBarry Smith 43481f0254dSBarry Smith static PetscErrorCode MatDestroy_BlockMat(Mat mat) 435421e10b8SBarry Smith { 436421e10b8SBarry Smith PetscErrorCode ierr; 437421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 438ccb205f8SBarry Smith PetscInt i; 439421e10b8SBarry Smith 440421e10b8SBarry Smith PetscFunctionBegin; 4416bf464f9SBarry Smith ierr = VecDestroy(&bmat->right);CHKERRQ(ierr); 4426bf464f9SBarry Smith ierr = VecDestroy(&bmat->left);CHKERRQ(ierr); 4436bf464f9SBarry Smith ierr = VecDestroy(&bmat->middle);CHKERRQ(ierr); 4446bf464f9SBarry Smith ierr = VecDestroy(&bmat->workb);CHKERRQ(ierr); 445ccb205f8SBarry Smith if (bmat->diags) { 446d0f46423SBarry Smith for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) { 4476bf464f9SBarry Smith ierr = MatDestroy(&bmat->diags[i]);CHKERRQ(ierr); 448ccb205f8SBarry Smith } 449ccb205f8SBarry Smith } 450ccb205f8SBarry Smith if (bmat->a) { 451ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 4526bf464f9SBarry Smith ierr = MatDestroy(&bmat->a[i]);CHKERRQ(ierr); 453ccb205f8SBarry Smith } 454ccb205f8SBarry Smith } 455ccb205f8SBarry Smith ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 456bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 457421e10b8SBarry Smith PetscFunctionReturn(0); 458421e10b8SBarry Smith } 459421e10b8SBarry Smith 46081f0254dSBarry Smith static PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 461421e10b8SBarry Smith { 462421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 463421e10b8SBarry Smith PetscErrorCode ierr; 464421e10b8SBarry Smith PetscScalar *xx,*yy; 465d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 466421e10b8SBarry Smith Mat *aa; 467421e10b8SBarry Smith 468421e10b8SBarry Smith PetscFunctionBegin; 469421e10b8SBarry Smith /* 470421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 471421e10b8SBarry Smith */ 472421e10b8SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 473ccb205f8SBarry Smith 474421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 475421e10b8SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 476421e10b8SBarry Smith aj = bmat->j; 477421e10b8SBarry Smith aa = bmat->a; 478421e10b8SBarry Smith ii = bmat->i; 479421e10b8SBarry Smith for (i=0; i<m; i++) { 480421e10b8SBarry Smith jrow = ii[i]; 481ccb205f8SBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 482421e10b8SBarry Smith n = ii[i+1] - jrow; 483421e10b8SBarry Smith for (j=0; j<n; j++) { 484421e10b8SBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 485ccb205f8SBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 486421e10b8SBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 487421e10b8SBarry Smith jrow++; 488421e10b8SBarry Smith } 489421e10b8SBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 490421e10b8SBarry Smith } 491421e10b8SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 492421e10b8SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 493421e10b8SBarry Smith PetscFunctionReturn(0); 494421e10b8SBarry Smith } 495421e10b8SBarry Smith 496290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 497290bbb0aSBarry Smith { 498290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 499290bbb0aSBarry Smith PetscErrorCode ierr; 500290bbb0aSBarry Smith PetscScalar *xx,*yy; 501d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 502290bbb0aSBarry Smith Mat *aa; 503290bbb0aSBarry Smith 504290bbb0aSBarry Smith PetscFunctionBegin; 505290bbb0aSBarry Smith /* 506290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 507290bbb0aSBarry Smith */ 508290bbb0aSBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 509290bbb0aSBarry Smith 510290bbb0aSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 511290bbb0aSBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 512290bbb0aSBarry Smith aj = bmat->j; 513290bbb0aSBarry Smith aa = bmat->a; 514290bbb0aSBarry Smith ii = bmat->i; 515290bbb0aSBarry Smith for (i=0; i<m; i++) { 516290bbb0aSBarry Smith jrow = ii[i]; 517290bbb0aSBarry Smith n = ii[i+1] - jrow; 518290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 519290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr); 520290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 521290bbb0aSBarry Smith if (aj[jrow] == i) { 522290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 523290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 524290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 525290bbb0aSBarry Smith jrow++; 526290bbb0aSBarry Smith n--; 527290bbb0aSBarry Smith } 528290bbb0aSBarry Smith for (j=0; j<n; j++) { 529290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); /* upper triangular part */ 530290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 531290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 532290bbb0aSBarry Smith 533290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr); /* lower triangular part */ 534290bbb0aSBarry Smith ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr); 535290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 536290bbb0aSBarry Smith jrow++; 537290bbb0aSBarry Smith } 538290bbb0aSBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 539290bbb0aSBarry Smith ierr = VecResetArray(bmat->middle);CHKERRQ(ierr); 540290bbb0aSBarry Smith } 541290bbb0aSBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 542290bbb0aSBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 543290bbb0aSBarry Smith PetscFunctionReturn(0); 544290bbb0aSBarry Smith } 545290bbb0aSBarry Smith 54681f0254dSBarry Smith static PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 547421e10b8SBarry Smith { 548421e10b8SBarry Smith PetscFunctionBegin; 549421e10b8SBarry Smith PetscFunctionReturn(0); 550421e10b8SBarry Smith } 551421e10b8SBarry Smith 55281f0254dSBarry Smith static PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 553421e10b8SBarry Smith { 554421e10b8SBarry Smith PetscFunctionBegin; 555421e10b8SBarry Smith PetscFunctionReturn(0); 556421e10b8SBarry Smith } 557421e10b8SBarry Smith 55881f0254dSBarry Smith static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 559421e10b8SBarry Smith { 560421e10b8SBarry Smith PetscFunctionBegin; 561421e10b8SBarry Smith PetscFunctionReturn(0); 562421e10b8SBarry Smith } 563421e10b8SBarry Smith 564ccb205f8SBarry Smith /* 565ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 566ccb205f8SBarry Smith */ 56781f0254dSBarry Smith static PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 568ccb205f8SBarry Smith { 569ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 570ccb205f8SBarry Smith PetscErrorCode ierr; 571d0f46423SBarry Smith PetscInt i,j,mbs = A->rmap->n/A->rmap->bs; 572ccb205f8SBarry Smith 573ccb205f8SBarry Smith PetscFunctionBegin; 574ccb205f8SBarry Smith if (!a->diag) { 575785e854fSJed Brown ierr = PetscMalloc1(mbs,&a->diag);CHKERRQ(ierr); 576ccb205f8SBarry Smith } 577ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 578ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 579ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 580ccb205f8SBarry Smith if (a->j[j] == i) { 581ccb205f8SBarry Smith a->diag[i] = j; 582ccb205f8SBarry Smith break; 583ccb205f8SBarry Smith } 584ccb205f8SBarry Smith } 585ccb205f8SBarry Smith } 586ccb205f8SBarry Smith PetscFunctionReturn(0); 587ccb205f8SBarry Smith } 588ccb205f8SBarry Smith 5897dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 590ccb205f8SBarry Smith { 591ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 592ccb205f8SBarry Smith Mat_SeqAIJ *c; 593ccb205f8SBarry Smith PetscErrorCode ierr; 594ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 595d2a212eaSBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ailen = a->ilen; 5961620fd73SBarry Smith PetscScalar *a_new; 597ccb205f8SBarry Smith Mat C,*aa = a->a; 598ace3abfcSBarry Smith PetscBool stride,equal; 599ccb205f8SBarry Smith 600ccb205f8SBarry Smith PetscFunctionBegin; 601ccb205f8SBarry Smith ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 602e32f2f54SBarry Smith if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices"); 603251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);CHKERRQ(ierr); 604e32f2f54SBarry Smith if (!stride) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices"); 605ccb205f8SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 606e32f2f54SBarry Smith if (step != A->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block"); 607ccb205f8SBarry Smith 608ccb205f8SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 609ccb205f8SBarry Smith ncols = nrows; 610ccb205f8SBarry Smith 611ccb205f8SBarry Smith /* create submatrix */ 612ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 613ccb205f8SBarry Smith PetscInt n_cols,n_rows; 614ccb205f8SBarry Smith C = *B; 615ccb205f8SBarry Smith ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 616e32f2f54SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 617ccb205f8SBarry Smith ierr = MatZeroEntries(C);CHKERRQ(ierr); 618ccb205f8SBarry Smith } else { 619ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 620ccb205f8SBarry Smith ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 6216e63c7a1SBarry Smith if (A->symmetric) { 6226e63c7a1SBarry Smith ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 6236e63c7a1SBarry Smith } else { 624ccb205f8SBarry Smith ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 6256e63c7a1SBarry Smith } 6266e63c7a1SBarry Smith ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr); 6276e63c7a1SBarry Smith ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr); 628ccb205f8SBarry Smith } 629ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 630ccb205f8SBarry Smith 631ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 632ccb205f8SBarry Smith a_new = c->a; 633ccb205f8SBarry Smith j_new = c->j; 634ccb205f8SBarry Smith i_new = c->i; 635ccb205f8SBarry Smith 636ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 637ccb205f8SBarry Smith lensi = ailen[i]; 638ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 639ccb205f8SBarry Smith *j_new++ = *aj++; 6401620fd73SBarry Smith ierr = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr); 641ccb205f8SBarry Smith } 642ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 643ccb205f8SBarry Smith c->ilen[i] = lensi; 644ccb205f8SBarry Smith } 645ccb205f8SBarry Smith 646ccb205f8SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 647ccb205f8SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 648ccb205f8SBarry Smith *B = C; 649ccb205f8SBarry Smith PetscFunctionReturn(0); 650ccb205f8SBarry Smith } 651ccb205f8SBarry Smith 65281f0254dSBarry Smith static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 653421e10b8SBarry Smith { 654421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 655421e10b8SBarry Smith PetscErrorCode ierr; 656421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 657ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 658421e10b8SBarry Smith Mat *aa = a->a,*ap; 659421e10b8SBarry Smith 660421e10b8SBarry Smith PetscFunctionBegin; 661421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 662421e10b8SBarry Smith 663421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 664421e10b8SBarry Smith for (i=1; i<m; i++) { 665421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 666421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 667421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 668421e10b8SBarry Smith if (fshift) { 669421e10b8SBarry Smith ip = aj + ai[i]; 670421e10b8SBarry Smith ap = aa + ai[i]; 671421e10b8SBarry Smith N = ailen[i]; 672421e10b8SBarry Smith for (j=0; j<N; j++) { 673421e10b8SBarry Smith ip[j-fshift] = ip[j]; 674421e10b8SBarry Smith ap[j-fshift] = ap[j]; 675421e10b8SBarry Smith } 676421e10b8SBarry Smith } 677421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 678421e10b8SBarry Smith } 679421e10b8SBarry Smith if (m) { 680421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 681421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 682421e10b8SBarry Smith } 683421e10b8SBarry Smith /* reset ilen and imax for each row */ 684421e10b8SBarry Smith for (i=0; i<m; i++) { 685421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 686421e10b8SBarry Smith } 687421e10b8SBarry Smith a->nz = ai[m]; 688ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 689ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG) 690e32f2f54SBarry Smith if (!aa[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 691ccb205f8SBarry Smith #endif 692ccb205f8SBarry Smith ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 693ccb205f8SBarry Smith ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 694ccb205f8SBarry Smith } 695d0f46423SBarry Smith 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); 696421e10b8SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 697421e10b8SBarry Smith ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 6982205254eSKarl Rupp 6998e58a170SBarry Smith A->info.mallocs += a->reallocs; 700421e10b8SBarry Smith a->reallocs = 0; 701421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 702421e10b8SBarry Smith a->rmax = rmax; 703ccb205f8SBarry Smith ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 704421e10b8SBarry Smith PetscFunctionReturn(0); 705421e10b8SBarry Smith } 706421e10b8SBarry Smith 70781f0254dSBarry Smith static PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg) 708290bbb0aSBarry Smith { 709290bbb0aSBarry Smith PetscFunctionBegin; 7104e0d8c25SBarry Smith if (opt == MAT_SYMMETRIC && flg) { 71141f059aeSBarry Smith A->ops->sor = MatSOR_BlockMat_Symmetric; 712290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 713290bbb0aSBarry Smith } else { 714290bbb0aSBarry Smith PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]); 715290bbb0aSBarry Smith } 716290bbb0aSBarry Smith PetscFunctionReturn(0); 717290bbb0aSBarry Smith } 718290bbb0aSBarry Smith 719290bbb0aSBarry Smith 7203964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 721421e10b8SBarry Smith 0, 722421e10b8SBarry Smith 0, 723421e10b8SBarry Smith MatMult_BlockMat, 724421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 725421e10b8SBarry Smith MatMultTranspose_BlockMat, 726421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 727421e10b8SBarry Smith 0, 728421e10b8SBarry Smith 0, 729421e10b8SBarry Smith 0, 730421e10b8SBarry Smith /* 10*/ 0, 731421e10b8SBarry Smith 0, 732421e10b8SBarry Smith 0, 73341f059aeSBarry Smith MatSOR_BlockMat, 734421e10b8SBarry Smith 0, 735421e10b8SBarry Smith /* 15*/ 0, 736421e10b8SBarry Smith 0, 737421e10b8SBarry Smith 0, 738421e10b8SBarry Smith 0, 739421e10b8SBarry Smith 0, 740421e10b8SBarry Smith /* 20*/ 0, 741421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 742290bbb0aSBarry Smith MatSetOption_BlockMat, 743421e10b8SBarry Smith 0, 744d519adbfSMatthew Knepley /* 24*/ 0, 745421e10b8SBarry Smith 0, 746421e10b8SBarry Smith 0, 747421e10b8SBarry Smith 0, 748421e10b8SBarry Smith 0, 749d519adbfSMatthew Knepley /* 29*/ 0, 750421e10b8SBarry Smith 0, 751421e10b8SBarry Smith 0, 752421e10b8SBarry Smith 0, 753421e10b8SBarry Smith 0, 754d519adbfSMatthew Knepley /* 34*/ 0, 755421e10b8SBarry Smith 0, 756421e10b8SBarry Smith 0, 757421e10b8SBarry Smith 0, 758421e10b8SBarry Smith 0, 759d519adbfSMatthew Knepley /* 39*/ 0, 760421e10b8SBarry Smith 0, 761421e10b8SBarry Smith 0, 762421e10b8SBarry Smith 0, 763421e10b8SBarry Smith 0, 764d519adbfSMatthew Knepley /* 44*/ 0, 765421e10b8SBarry Smith 0, 7667d68702bSBarry Smith MatShift_Basic, 767421e10b8SBarry Smith 0, 768421e10b8SBarry Smith 0, 769d519adbfSMatthew Knepley /* 49*/ 0, 770421e10b8SBarry Smith 0, 771421e10b8SBarry Smith 0, 772421e10b8SBarry Smith 0, 773421e10b8SBarry Smith 0, 774d519adbfSMatthew Knepley /* 54*/ 0, 775421e10b8SBarry Smith 0, 776421e10b8SBarry Smith 0, 777421e10b8SBarry Smith 0, 778421e10b8SBarry Smith 0, 7797dae84e0SHong Zhang /* 59*/ MatCreateSubMatrix_BlockMat, 780421e10b8SBarry Smith MatDestroy_BlockMat, 781ccb205f8SBarry Smith MatView_BlockMat, 782421e10b8SBarry Smith 0, 783421e10b8SBarry Smith 0, 784d519adbfSMatthew Knepley /* 64*/ 0, 785421e10b8SBarry Smith 0, 786421e10b8SBarry Smith 0, 787421e10b8SBarry Smith 0, 788421e10b8SBarry Smith 0, 789d519adbfSMatthew Knepley /* 69*/ 0, 790421e10b8SBarry Smith 0, 791421e10b8SBarry Smith 0, 792421e10b8SBarry Smith 0, 793421e10b8SBarry Smith 0, 794d519adbfSMatthew Knepley /* 74*/ 0, 795421e10b8SBarry Smith 0, 796421e10b8SBarry Smith 0, 797421e10b8SBarry Smith 0, 798421e10b8SBarry Smith 0, 799d519adbfSMatthew Knepley /* 79*/ 0, 800421e10b8SBarry Smith 0, 801421e10b8SBarry Smith 0, 802421e10b8SBarry Smith 0, 8035bba2384SShri Abhyankar MatLoad_BlockMat, 804d519adbfSMatthew Knepley /* 84*/ 0, 805421e10b8SBarry Smith 0, 806421e10b8SBarry Smith 0, 807421e10b8SBarry Smith 0, 808421e10b8SBarry Smith 0, 809d519adbfSMatthew Knepley /* 89*/ 0, 810421e10b8SBarry Smith 0, 811421e10b8SBarry Smith 0, 812421e10b8SBarry Smith 0, 813421e10b8SBarry Smith 0, 814d519adbfSMatthew Knepley /* 94*/ 0, 815421e10b8SBarry Smith 0, 816421e10b8SBarry Smith 0, 817a5973382SShri Abhyankar 0, 818a5973382SShri Abhyankar 0, 819a5973382SShri Abhyankar /* 99*/ 0, 820a5973382SShri Abhyankar 0, 821a5973382SShri Abhyankar 0, 822a5973382SShri Abhyankar 0, 823a5973382SShri Abhyankar 0, 824a5973382SShri Abhyankar /*104*/ 0, 825a5973382SShri Abhyankar 0, 826a5973382SShri Abhyankar 0, 827a5973382SShri Abhyankar 0, 828a5973382SShri Abhyankar 0, 829a5973382SShri Abhyankar /*109*/ 0, 830a5973382SShri Abhyankar 0, 831a5973382SShri Abhyankar 0, 832a5973382SShri Abhyankar 0, 833a5973382SShri Abhyankar 0, 834a5973382SShri Abhyankar /*114*/ 0, 835a5973382SShri Abhyankar 0, 836a5973382SShri Abhyankar 0, 837a5973382SShri Abhyankar 0, 838a5973382SShri Abhyankar 0, 839a5973382SShri Abhyankar /*119*/ 0, 840a5973382SShri Abhyankar 0, 841a5973382SShri Abhyankar 0, 8423964eb88SJed Brown 0, 8433964eb88SJed Brown 0, 8443964eb88SJed Brown /*124*/ 0, 8453964eb88SJed Brown 0, 8463964eb88SJed Brown 0, 8473964eb88SJed Brown 0, 8483964eb88SJed Brown 0, 8493964eb88SJed Brown /*129*/ 0, 8503964eb88SJed Brown 0, 8513964eb88SJed Brown 0, 8523964eb88SJed Brown 0, 8533964eb88SJed Brown 0, 8543964eb88SJed Brown /*134*/ 0, 8553964eb88SJed Brown 0, 8563964eb88SJed Brown 0, 8573964eb88SJed Brown 0, 8583964eb88SJed Brown 0, 8593964eb88SJed Brown /*139*/ 0, 860f9426fe0SMark Adams 0, 8615bba2384SShri Abhyankar 0 862a5973382SShri Abhyankar }; 863421e10b8SBarry Smith 8648cdf2d9bSBarry Smith /*@C 8658cdf2d9bSBarry Smith MatBlockMatSetPreallocation - For good matrix assembly performance 8668cdf2d9bSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 8678cdf2d9bSBarry Smith (or the array nnz). By setting these parameters accurately, performance 8688cdf2d9bSBarry Smith during matrix assembly can be increased by more than a factor of 50. 8698cdf2d9bSBarry Smith 8708cdf2d9bSBarry Smith Collective on MPI_Comm 8718cdf2d9bSBarry Smith 8728cdf2d9bSBarry Smith Input Parameters: 8738cdf2d9bSBarry Smith + B - The matrix 8748cdf2d9bSBarry Smith . bs - size of each block in matrix 8758cdf2d9bSBarry Smith . nz - number of nonzeros per block row (same for all rows) 8768cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows 8770298fd71SBarry Smith (possibly different for each row) or NULL 8788cdf2d9bSBarry Smith 8798cdf2d9bSBarry Smith Notes: 8808cdf2d9bSBarry Smith If nnz is given then nz is ignored 8818cdf2d9bSBarry Smith 8828cdf2d9bSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 8830298fd71SBarry Smith Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 8848cdf2d9bSBarry Smith allocation. For large problems you MUST preallocate memory or you 8858cdf2d9bSBarry Smith will get TERRIBLE performance, see the users' manual chapter on matrices. 8868cdf2d9bSBarry Smith 8878cdf2d9bSBarry Smith Level: intermediate 8888cdf2d9bSBarry Smith 8898cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues() 8908cdf2d9bSBarry Smith 8918cdf2d9bSBarry Smith @*/ 8927087cfbeSBarry Smith PetscErrorCode MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 8938cdf2d9bSBarry Smith { 8944ac538c5SBarry Smith PetscErrorCode ierr; 8958cdf2d9bSBarry Smith 8968cdf2d9bSBarry Smith PetscFunctionBegin; 8974ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr); 8988cdf2d9bSBarry Smith PetscFunctionReturn(0); 8998cdf2d9bSBarry Smith } 9008cdf2d9bSBarry Smith 90181f0254dSBarry Smith static PetscErrorCode MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz) 9028cdf2d9bSBarry Smith { 9038cdf2d9bSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 9048cdf2d9bSBarry Smith PetscErrorCode ierr; 9058cdf2d9bSBarry Smith PetscInt i; 9068cdf2d9bSBarry Smith 9078cdf2d9bSBarry Smith PetscFunctionBegin; 908e02043d6SBarry Smith ierr = PetscLayoutSetBlockSize(A->rmap,bs);CHKERRQ(ierr); 909e02043d6SBarry Smith ierr = PetscLayoutSetBlockSize(A->cmap,bs);CHKERRQ(ierr); 91034ef9618SShri Abhyankar ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 91134ef9618SShri Abhyankar ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 912e02043d6SBarry Smith ierr = PetscLayoutGetBlockSize(A->rmap,&bs);CHKERRQ(ierr); 91334ef9618SShri Abhyankar 9148cdf2d9bSBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 915e32f2f54SBarry Smith if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 9168cdf2d9bSBarry Smith if (nnz) { 917d0f46423SBarry Smith for (i=0; i<A->rmap->n/bs; i++) { 918e32f2f54SBarry Smith 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]); 919e32f2f54SBarry Smith 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); 9208cdf2d9bSBarry Smith } 9218cdf2d9bSBarry Smith } 922d0f46423SBarry Smith bmat->mbs = A->rmap->n/bs; 9238cdf2d9bSBarry Smith 9240298fd71SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->right);CHKERRQ(ierr); 9250298fd71SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->middle);CHKERRQ(ierr); 9268cdf2d9bSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr); 9278cdf2d9bSBarry Smith 9282ee49352SLisandro Dalcin if (!bmat->imax) { 929dcca6d9dSJed Brown ierr = PetscMalloc2(A->rmap->n,&bmat->imax,A->rmap->n,&bmat->ilen);CHKERRQ(ierr); 9303bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 9312ee49352SLisandro Dalcin } 9328cdf2d9bSBarry Smith if (nnz) { 9338cdf2d9bSBarry Smith nz = 0; 934d0f46423SBarry Smith for (i=0; i<A->rmap->n/A->rmap->bs; i++) { 9358cdf2d9bSBarry Smith bmat->imax[i] = nnz[i]; 9368cdf2d9bSBarry Smith nz += nnz[i]; 9378cdf2d9bSBarry Smith } 938f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation"); 9398cdf2d9bSBarry Smith 9408cdf2d9bSBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 9412205254eSKarl Rupp for (i=0; i<bmat->mbs; i++) bmat->ilen[i] = 0; 9428cdf2d9bSBarry Smith 9438cdf2d9bSBarry Smith /* allocate the matrix space */ 9442ee49352SLisandro Dalcin ierr = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 945dcca6d9dSJed Brown ierr = PetscMalloc3(nz,&bmat->a,nz,&bmat->j,A->rmap->n+1,&bmat->i);CHKERRQ(ierr); 9463bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 9478cdf2d9bSBarry Smith bmat->i[0] = 0; 9488cdf2d9bSBarry Smith for (i=1; i<bmat->mbs+1; i++) { 9498cdf2d9bSBarry Smith bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 9508cdf2d9bSBarry Smith } 9518cdf2d9bSBarry Smith bmat->singlemalloc = PETSC_TRUE; 9528cdf2d9bSBarry Smith bmat->free_a = PETSC_TRUE; 9538cdf2d9bSBarry Smith bmat->free_ij = PETSC_TRUE; 9548cdf2d9bSBarry Smith 9558cdf2d9bSBarry Smith bmat->nz = 0; 9568cdf2d9bSBarry Smith bmat->maxnz = nz; 9578cdf2d9bSBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 9587827cd58SJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 9598cdf2d9bSBarry Smith PetscFunctionReturn(0); 9608cdf2d9bSBarry Smith } 9618cdf2d9bSBarry Smith 962421e10b8SBarry Smith /*MC 963421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 964421e10b8SBarry Smith consisting of (usually) sparse blocks. 965421e10b8SBarry Smith 966421e10b8SBarry Smith Level: advanced 967421e10b8SBarry Smith 968421e10b8SBarry Smith .seealso: MatCreateBlockMat() 969421e10b8SBarry Smith 970421e10b8SBarry Smith M*/ 971421e10b8SBarry Smith 9728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A) 973421e10b8SBarry Smith { 974421e10b8SBarry Smith Mat_BlockMat *b; 975421e10b8SBarry Smith PetscErrorCode ierr; 976421e10b8SBarry Smith 977421e10b8SBarry Smith PetscFunctionBegin; 978b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 979421e10b8SBarry Smith A->data = (void*)b; 98038f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 981421e10b8SBarry Smith 982421e10b8SBarry Smith A->assembled = PETSC_TRUE; 983421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 984421e10b8SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 985290bbb0aSBarry Smith 986bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatBlockMatSetPreallocation_C",MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr); 987421e10b8SBarry Smith PetscFunctionReturn(0); 988421e10b8SBarry Smith } 989421e10b8SBarry Smith 990421e10b8SBarry Smith /*@C 991bbd3f336SJed Brown MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential Mat object 992421e10b8SBarry Smith 993421e10b8SBarry Smith Collective on MPI_Comm 994421e10b8SBarry Smith 995421e10b8SBarry Smith Input Parameters: 996421e10b8SBarry Smith + comm - MPI communicator 997421e10b8SBarry Smith . m - number of rows 998421e10b8SBarry Smith . n - number of columns 9998cdf2d9bSBarry Smith . bs - size of each submatrix 10008cdf2d9bSBarry Smith . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 10010298fd71SBarry Smith - nnz - expected number of nonzers per block row if known (use NULL otherwise) 1002421e10b8SBarry Smith 1003421e10b8SBarry Smith 1004421e10b8SBarry Smith Output Parameter: 1005421e10b8SBarry Smith . A - the matrix 1006421e10b8SBarry Smith 1007421e10b8SBarry Smith Level: intermediate 1008421e10b8SBarry Smith 1009*95452b02SPatrick Sanan Notes: 1010*95452b02SPatrick Sanan Matrices of this type are nominally-sparse matrices in which each "entry" is a Mat object. Each Mat must 1011bbd3f336SJed Brown have the same size and be sequential. The local and global sizes must be compatible with this decomposition. 1012bbd3f336SJed Brown 1013bbd3f336SJed Brown For matrices containing parallel submatrices and variable block sizes, see MATNEST. 1014421e10b8SBarry Smith 1015421e10b8SBarry Smith .keywords: matrix, bmat, create 1016421e10b8SBarry Smith 1017bbd3f336SJed Brown .seealso: MATBLOCKMAT, MatCreateNest() 1018421e10b8SBarry Smith @*/ 10197087cfbeSBarry Smith PetscErrorCode MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A) 1020421e10b8SBarry Smith { 1021421e10b8SBarry Smith PetscErrorCode ierr; 1022421e10b8SBarry Smith 1023421e10b8SBarry Smith PetscFunctionBegin; 1024421e10b8SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1025421e10b8SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1026421e10b8SBarry Smith ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 10278cdf2d9bSBarry Smith ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1028421e10b8SBarry Smith PetscFunctionReturn(0); 1029421e10b8SBarry Smith } 1030