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 17e0877f53SBarry Smith static PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 18290bbb0aSBarry Smith { 19290bbb0aSBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 20290bbb0aSBarry Smith PetscScalar *x; 21a2ea699eSBarry Smith const Mat *v; 22290bbb0aSBarry Smith const PetscScalar *b; 23d0f46423SBarry Smith PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 24290bbb0aSBarry Smith const PetscInt *idx; 25290bbb0aSBarry Smith IS row,col; 26290bbb0aSBarry Smith MatFactorInfo info; 276e63c7a1SBarry Smith Vec left = a->left,right = a->right, middle = a->middle; 28290bbb0aSBarry Smith Mat *diag; 29290bbb0aSBarry Smith 30290bbb0aSBarry Smith PetscFunctionBegin; 31290bbb0aSBarry Smith its = its*lits; 3208401ef6SPierre Jolivet PetscCheck(its > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 33aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT),PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 3408401ef6SPierre Jolivet PetscCheck(omega == 1.0,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 3528b400f6SJacob Faibussowitsch PetscCheck(!fshift,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 36*7a46b595SBarry Smith PetscCheck(!((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)),PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep"); 37290bbb0aSBarry Smith 38290bbb0aSBarry Smith if (!a->diags) { 399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs,&a->diags)); 409566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 41290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 429566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col)); 439566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info)); 449566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info)); 459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 47290bbb0aSBarry Smith } 489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(bb,&a->workb)); 49290bbb0aSBarry Smith } 50290bbb0aSBarry Smith diag = a->diags; 51290bbb0aSBarry Smith 529566063dSJacob Faibussowitsch PetscCall(VecSet(xx,0.0)); 539566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx,&x)); 54290bbb0aSBarry Smith /* copy right hand side because it must be modified during iteration */ 559566063dSJacob Faibussowitsch PetscCall(VecCopy(bb,a->workb)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(a->workb,&b)); 57290bbb0aSBarry Smith 5841f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 59290bbb0aSBarry Smith while (its--) { 60290bbb0aSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 61290bbb0aSBarry Smith 62290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 636e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 646e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 656e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 66290bbb0aSBarry Smith 679566063dSJacob Faibussowitsch PetscCall(VecSet(left,0.0)); 68290bbb0aSBarry Smith for (j=0; j<n; j++) { 699566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + idx[j]*bs)); 709566063dSJacob Faibussowitsch PetscCall(MatMultAdd(v[j],right,left,left)); 719566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 72290bbb0aSBarry Smith } 739566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,b + i*bs)); 749566063dSJacob Faibussowitsch PetscCall(VecAYPX(left,-1.0,right)); 759566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 76290bbb0aSBarry Smith 779566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + i*bs)); 789566063dSJacob Faibussowitsch PetscCall(MatSolve(diag[i],left,right)); 796e63c7a1SBarry Smith 8041f059aeSBarry Smith /* now adjust right hand side, see MatSOR_SeqSBAIJ */ 816e63c7a1SBarry Smith for (j=0; j<n; j++) { 829566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(v[j],right,left)); 839566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(middle,b + idx[j]*bs)); 849566063dSJacob Faibussowitsch PetscCall(VecAXPY(middle,-1.0,left)); 859566063dSJacob Faibussowitsch PetscCall(VecResetArray(middle)); 866e63c7a1SBarry Smith } 879566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 886e63c7a1SBarry Smith 89290bbb0aSBarry Smith } 90290bbb0aSBarry Smith } 91290bbb0aSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 92290bbb0aSBarry Smith 93290bbb0aSBarry Smith for (i=mbs-1; i>=0; i--) { 946e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 956e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 966e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 97290bbb0aSBarry Smith 989566063dSJacob Faibussowitsch PetscCall(VecSet(left,0.0)); 99290bbb0aSBarry Smith for (j=0; j<n; j++) { 1009566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + idx[j]*bs)); 1019566063dSJacob Faibussowitsch PetscCall(MatMultAdd(v[j],right,left,left)); 1029566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 103290bbb0aSBarry Smith } 1049566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,b + i*bs)); 1059566063dSJacob Faibussowitsch PetscCall(VecAYPX(left,-1.0,right)); 1069566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 107290bbb0aSBarry Smith 1089566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + i*bs)); 1099566063dSJacob Faibussowitsch PetscCall(MatSolve(diag[i],left,right)); 1109566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 111290bbb0aSBarry Smith 112290bbb0aSBarry Smith } 113290bbb0aSBarry Smith } 114290bbb0aSBarry Smith } 1159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx,&x)); 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(a->workb,&b)); 117290bbb0aSBarry Smith PetscFunctionReturn(0); 118290bbb0aSBarry Smith } 119290bbb0aSBarry Smith 12081f0254dSBarry Smith static PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 121ccb205f8SBarry Smith { 122ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 123ccb205f8SBarry Smith PetscScalar *x; 124a2ea699eSBarry Smith const Mat *v; 125ccb205f8SBarry Smith const PetscScalar *b; 126d0f46423SBarry Smith PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 127ccb205f8SBarry Smith const PetscInt *idx; 128ccb205f8SBarry Smith IS row,col; 129ccb205f8SBarry Smith MatFactorInfo info; 130ccb205f8SBarry Smith Vec left = a->left,right = a->right; 131ccb205f8SBarry Smith Mat *diag; 132ccb205f8SBarry Smith 133ccb205f8SBarry Smith PetscFunctionBegin; 134ccb205f8SBarry Smith its = its*lits; 13508401ef6SPierre Jolivet PetscCheck(its > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 136aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT),PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 13708401ef6SPierre Jolivet PetscCheck(omega == 1.0,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 13828b400f6SJacob Faibussowitsch PetscCheck(!fshift,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 139ccb205f8SBarry Smith 140ccb205f8SBarry Smith if (!a->diags) { 1419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs,&a->diags)); 1429566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 143ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 1449566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col)); 1459566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info)); 1469566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info)); 1479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 149ccb205f8SBarry Smith } 150ccb205f8SBarry Smith } 151ccb205f8SBarry Smith diag = a->diags; 152ccb205f8SBarry Smith 1539566063dSJacob Faibussowitsch PetscCall(VecSet(xx,0.0)); 1549566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx,&x)); 1559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb,&b)); 156ccb205f8SBarry Smith 15741f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 158ccb205f8SBarry Smith while (its--) { 159ccb205f8SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 160ccb205f8SBarry Smith 161ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 162ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 163ccb205f8SBarry Smith idx = a->j + a->i[i]; 164ccb205f8SBarry Smith v = a->a + a->i[i]; 165ccb205f8SBarry Smith 1669566063dSJacob Faibussowitsch PetscCall(VecSet(left,0.0)); 167ccb205f8SBarry Smith for (j=0; j<n; j++) { 168ccb205f8SBarry Smith if (idx[j] != i) { 1699566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + idx[j]*bs)); 1709566063dSJacob Faibussowitsch PetscCall(MatMultAdd(v[j],right,left,left)); 1719566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 172ccb205f8SBarry Smith } 173ccb205f8SBarry Smith } 1749566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,b + i*bs)); 1759566063dSJacob Faibussowitsch PetscCall(VecAYPX(left,-1.0,right)); 1769566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 177ccb205f8SBarry Smith 1789566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + i*bs)); 1799566063dSJacob Faibussowitsch PetscCall(MatSolve(diag[i],left,right)); 1809566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 181ccb205f8SBarry Smith } 182ccb205f8SBarry Smith } 183ccb205f8SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 184ccb205f8SBarry Smith 185ccb205f8SBarry Smith for (i=mbs-1; i>=0; i--) { 186ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 187ccb205f8SBarry Smith idx = a->j + a->i[i]; 188ccb205f8SBarry Smith v = a->a + a->i[i]; 189ccb205f8SBarry Smith 1909566063dSJacob Faibussowitsch PetscCall(VecSet(left,0.0)); 191ccb205f8SBarry Smith for (j=0; j<n; j++) { 192ccb205f8SBarry Smith if (idx[j] != i) { 1939566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + idx[j]*bs)); 1949566063dSJacob Faibussowitsch PetscCall(MatMultAdd(v[j],right,left,left)); 1959566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 196ccb205f8SBarry Smith } 197ccb205f8SBarry Smith } 1989566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,b + i*bs)); 1999566063dSJacob Faibussowitsch PetscCall(VecAYPX(left,-1.0,right)); 2009566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 201ccb205f8SBarry Smith 2029566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right,x + i*bs)); 2039566063dSJacob Faibussowitsch PetscCall(MatSolve(diag[i],left,right)); 2049566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 205ccb205f8SBarry Smith 206ccb205f8SBarry Smith } 207ccb205f8SBarry Smith } 208ccb205f8SBarry Smith } 2099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx,&x)); 2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb,&b)); 211ccb205f8SBarry Smith PetscFunctionReturn(0); 212ccb205f8SBarry Smith } 213ccb205f8SBarry Smith 21481f0254dSBarry Smith static PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 215421e10b8SBarry Smith { 216421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 217421e10b8SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 218421e10b8SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 219d0f46423SBarry Smith PetscInt *aj =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 220421e10b8SBarry Smith PetscInt ridx,cidx; 221ace3abfcSBarry Smith PetscBool roworiented=a->roworiented; 222421e10b8SBarry Smith MatScalar value; 223421e10b8SBarry Smith Mat *ap,*aa = a->a; 224421e10b8SBarry Smith 225421e10b8SBarry Smith PetscFunctionBegin; 226421e10b8SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 227421e10b8SBarry Smith row = im[k]; 228421e10b8SBarry Smith brow = row/bs; 229421e10b8SBarry Smith if (row < 0) continue; 2306bdcaf15SBarry Smith PetscCheck(row < A->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,row,A->rmap->N-1); 231421e10b8SBarry Smith rp = aj + ai[brow]; 232421e10b8SBarry Smith ap = aa + ai[brow]; 233421e10b8SBarry Smith rmax = imax[brow]; 234421e10b8SBarry Smith nrow = ailen[brow]; 235421e10b8SBarry Smith low = 0; 236421e10b8SBarry Smith high = nrow; 237421e10b8SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 238421e10b8SBarry Smith if (in[l] < 0) continue; 2396bdcaf15SBarry Smith PetscCheck(in[l] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,in[l],A->cmap->n-1); 240421e10b8SBarry Smith col = in[l]; bcol = col/bs; 241b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE && brow > bcol) continue; 242421e10b8SBarry Smith ridx = row % bs; cidx = col % bs; 2432205254eSKarl Rupp if (roworiented) value = v[l + k*n]; 2442205254eSKarl Rupp else value = v[k + l*m]; 2452205254eSKarl Rupp 2462205254eSKarl Rupp if (col <= lastcol) low = 0; 2472205254eSKarl Rupp else high = nrow; 248421e10b8SBarry Smith lastcol = col; 249421e10b8SBarry Smith while (high-low > 7) { 250421e10b8SBarry Smith t = (low+high)/2; 251421e10b8SBarry Smith if (rp[t] > bcol) high = t; 252421e10b8SBarry Smith else low = t; 253421e10b8SBarry Smith } 254421e10b8SBarry Smith for (i=low; i<high; i++) { 255421e10b8SBarry Smith if (rp[i] > bcol) break; 2562205254eSKarl Rupp if (rp[i] == bcol) goto noinsert1; 257421e10b8SBarry Smith } 258421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 25908401ef6SPierre Jolivet PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col); 260fef13f97SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 261421e10b8SBarry Smith N = nrow++ - 1; high++; 262421e10b8SBarry Smith /* shift up all the later entries in this row */ 263421e10b8SBarry Smith for (ii=N; ii>=i; ii--) { 264421e10b8SBarry Smith rp[ii+1] = rp[ii]; 265421e10b8SBarry Smith ap[ii+1] = ap[ii]; 266421e10b8SBarry Smith } 267f4259b30SLisandro Dalcin if (N>=i) ap[i] = NULL; 268421e10b8SBarry Smith rp[i] = bcol; 269421e10b8SBarry Smith a->nz++; 270e56f5c9eSBarry Smith A->nonzerostate++; 271421e10b8SBarry Smith noinsert1:; 272421e10b8SBarry Smith if (!*(ap+i)) { 2739566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,NULL,ap+i)); 274b0223207SBarry Smith } 2759566063dSJacob Faibussowitsch PetscCall(MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is)); 276421e10b8SBarry Smith low = i; 277421e10b8SBarry Smith } 278421e10b8SBarry Smith ailen[brow] = nrow; 279421e10b8SBarry Smith } 280421e10b8SBarry Smith PetscFunctionReturn(0); 281421e10b8SBarry Smith } 282421e10b8SBarry Smith 28381f0254dSBarry Smith static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer) 284a5973382SShri Abhyankar { 285a5973382SShri Abhyankar Mat tmpA; 286a5973382SShri Abhyankar PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 287a5973382SShri Abhyankar const PetscInt *cols; 288a5973382SShri Abhyankar const PetscScalar *values; 289ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE,notdone; 290a5973382SShri Abhyankar Mat_SeqAIJ *a; 291a5973382SShri Abhyankar Mat_BlockMat *amat; 292a5973382SShri Abhyankar 293a5973382SShri Abhyankar PetscFunctionBegin; 294c98fd787SBarry Smith /* force binary viewer to load .info file if it has not yet done so */ 2959566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 2969566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&tmpA)); 2979566063dSJacob Faibussowitsch PetscCall(MatSetType(tmpA,MATSEQAIJ)); 2989566063dSJacob Faibussowitsch PetscCall(MatLoad_SeqAIJ(tmpA,viewer)); 299a5973382SShri Abhyankar 3009566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(tmpA,&m,&n)); 301d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_SELF,NULL,"Options for loading BlockMat matrix 1","Mat"); 3029566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL)); 3039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,NULL)); 304d0609cedSBarry Smith PetscOptionsEnd(); 305a5973382SShri Abhyankar 306a5973382SShri Abhyankar /* Determine number of nonzero blocks for each block row */ 307a5973382SShri Abhyankar a = (Mat_SeqAIJ*) tmpA->data; 308a5973382SShri Abhyankar mbs = m/bs; 3099566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(mbs,&lens,bs,&ii,bs,&ilens)); 3109566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(lens,mbs)); 311a5973382SShri Abhyankar 312a5973382SShri Abhyankar for (i=0; i<mbs; i++) { 313a5973382SShri Abhyankar for (j=0; j<bs; j++) { 314a5973382SShri Abhyankar ii[j] = a->j + a->i[i*bs + j]; 315a5973382SShri Abhyankar ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 316a5973382SShri Abhyankar } 317a5973382SShri Abhyankar 318a5973382SShri Abhyankar currentcol = -1; 319a5973382SShri Abhyankar while (PETSC_TRUE) { 320a5973382SShri Abhyankar notdone = PETSC_FALSE; 321a5973382SShri Abhyankar nextcol = 1000000000; 322a5973382SShri Abhyankar for (j=0; j<bs; j++) { 323a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 324a5973382SShri Abhyankar ii[j]++; 325a5973382SShri Abhyankar ilens[j]--; 326a5973382SShri Abhyankar } 327a5973382SShri Abhyankar if (ilens[j] > 0) { 328a5973382SShri Abhyankar notdone = PETSC_TRUE; 329a5973382SShri Abhyankar nextcol = PetscMin(nextcol,ii[j][0]/bs); 330a5973382SShri Abhyankar } 331a5973382SShri Abhyankar } 332a5973382SShri Abhyankar if (!notdone) break; 333a5973382SShri Abhyankar if (!flg || (nextcol >= i)) lens[i]++; 334a5973382SShri Abhyankar currentcol = nextcol; 335a5973382SShri Abhyankar } 336a5973382SShri Abhyankar } 337a5973382SShri Abhyankar 338a5973382SShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 3399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE)); 340a5973382SShri Abhyankar } 3419566063dSJacob Faibussowitsch PetscCall(MatBlockMatSetPreallocation(newmat,bs,0,lens)); 3421baa6e33SBarry Smith if (flg) PetscCall(MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE)); 343a5973382SShri Abhyankar amat = (Mat_BlockMat*)(newmat)->data; 344a5973382SShri Abhyankar 345a5973382SShri Abhyankar /* preallocate the submatrices */ 3469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs,&llens)); 347a5973382SShri Abhyankar for (i=0; i<mbs; i++) { /* loops for block rows */ 348a5973382SShri Abhyankar for (j=0; j<bs; j++) { 349a5973382SShri Abhyankar ii[j] = a->j + a->i[i*bs + j]; 350a5973382SShri Abhyankar ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 351a5973382SShri Abhyankar } 352a5973382SShri Abhyankar 353a5973382SShri Abhyankar currentcol = 1000000000; 354a5973382SShri Abhyankar for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 355a5973382SShri Abhyankar if (ilens[j] > 0) { 356a5973382SShri Abhyankar currentcol = PetscMin(currentcol,ii[j][0]/bs); 357a5973382SShri Abhyankar } 358a5973382SShri Abhyankar } 359a5973382SShri Abhyankar 360a5973382SShri Abhyankar while (PETSC_TRUE) { /* loops over blocks in block row */ 361a5973382SShri Abhyankar notdone = PETSC_FALSE; 362a5973382SShri Abhyankar nextcol = 1000000000; 3639566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(llens,bs)); 364a5973382SShri Abhyankar for (j=0; j<bs; j++) { /* loop over rows in block */ 365a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 366a5973382SShri Abhyankar ii[j]++; 367a5973382SShri Abhyankar ilens[j]--; 368a5973382SShri Abhyankar llens[j]++; 369a5973382SShri Abhyankar } 370a5973382SShri Abhyankar if (ilens[j] > 0) { 371a5973382SShri Abhyankar notdone = PETSC_TRUE; 372a5973382SShri Abhyankar nextcol = PetscMin(nextcol,ii[j][0]/bs); 373a5973382SShri Abhyankar } 374a5973382SShri Abhyankar } 37508401ef6SPierre Jolivet PetscCheck(cnt < amat->maxnz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %" PetscInt_FMT,cnt); 376a5973382SShri Abhyankar if (!flg || currentcol >= i) { 377a5973382SShri Abhyankar amat->j[cnt] = currentcol; 3789566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++)); 379a5973382SShri Abhyankar } 380a5973382SShri Abhyankar 381a5973382SShri Abhyankar if (!notdone) break; 382a5973382SShri Abhyankar currentcol = nextcol; 383a5973382SShri Abhyankar } 384a5973382SShri Abhyankar amat->ilen[i] = lens[i]; 385a5973382SShri Abhyankar } 386a5973382SShri Abhyankar 3879566063dSJacob Faibussowitsch PetscCall(PetscFree3(lens,ii,ilens)); 3889566063dSJacob Faibussowitsch PetscCall(PetscFree(llens)); 389a5973382SShri Abhyankar 390a5973382SShri Abhyankar /* copy over the matrix, one row at a time */ 391a5973382SShri Abhyankar for (i=0; i<m; i++) { 3929566063dSJacob Faibussowitsch PetscCall(MatGetRow(tmpA,i,&ncols,&cols,&values)); 3939566063dSJacob Faibussowitsch PetscCall(MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES)); 3949566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tmpA,i,&ncols,&cols,&values)); 395a5973382SShri Abhyankar } 3969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY)); 3979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY)); 398a5973382SShri Abhyankar PetscFunctionReturn(0); 399a5973382SShri Abhyankar } 400a5973382SShri Abhyankar 40181f0254dSBarry Smith static PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 402ccb205f8SBarry Smith { 40364075487SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 404ccb205f8SBarry Smith const char *name; 405ccb205f8SBarry Smith PetscViewerFormat format; 406ccb205f8SBarry Smith 407ccb205f8SBarry Smith PetscFunctionBegin; 4089566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)A,&name)); 4099566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 410ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 4119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %" PetscInt_FMT " \n",a->nz)); 412b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n")); 413ccb205f8SBarry Smith } 414ccb205f8SBarry Smith PetscFunctionReturn(0); 415ccb205f8SBarry Smith } 416421e10b8SBarry Smith 41781f0254dSBarry Smith static PetscErrorCode MatDestroy_BlockMat(Mat mat) 418421e10b8SBarry Smith { 419421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 420ccb205f8SBarry Smith PetscInt i; 421421e10b8SBarry Smith 422421e10b8SBarry Smith PetscFunctionBegin; 4239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->right)); 4249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->left)); 4259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->middle)); 4269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->workb)); 427ccb205f8SBarry Smith if (bmat->diags) { 428d0f46423SBarry Smith for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) { 4299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bmat->diags[i])); 430ccb205f8SBarry Smith } 431ccb205f8SBarry Smith } 432ccb205f8SBarry Smith if (bmat->a) { 433ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 4349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bmat->a[i])); 435ccb205f8SBarry Smith } 436ccb205f8SBarry Smith } 4379566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i)); 4389566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 439421e10b8SBarry Smith PetscFunctionReturn(0); 440421e10b8SBarry Smith } 441421e10b8SBarry Smith 44281f0254dSBarry Smith static PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 443421e10b8SBarry Smith { 444421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 445421e10b8SBarry Smith PetscScalar *xx,*yy; 446d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 447421e10b8SBarry Smith Mat *aa; 448421e10b8SBarry Smith 449421e10b8SBarry Smith PetscFunctionBegin; 450421e10b8SBarry Smith /* 451421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 452421e10b8SBarry Smith */ 4539566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&xx)); 454ccb205f8SBarry Smith 4559566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 4569566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 457421e10b8SBarry Smith aj = bmat->j; 458421e10b8SBarry Smith aa = bmat->a; 459421e10b8SBarry Smith ii = bmat->i; 460421e10b8SBarry Smith for (i=0; i<m; i++) { 461421e10b8SBarry Smith jrow = ii[i]; 4629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->left,yy + bs*i)); 463421e10b8SBarry Smith n = ii[i+1] - jrow; 464421e10b8SBarry Smith for (j=0; j<n; j++) { 4659566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right,xx + bs*aj[jrow])); 4669566063dSJacob Faibussowitsch PetscCall(MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left)); 4679566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 468421e10b8SBarry Smith jrow++; 469421e10b8SBarry Smith } 4709566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->left)); 471421e10b8SBarry Smith } 4729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xx)); 4739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 474421e10b8SBarry Smith PetscFunctionReturn(0); 475421e10b8SBarry Smith } 476421e10b8SBarry Smith 477290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 478290bbb0aSBarry Smith { 479290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 480290bbb0aSBarry Smith PetscScalar *xx,*yy; 481d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 482290bbb0aSBarry Smith Mat *aa; 483290bbb0aSBarry Smith 484290bbb0aSBarry Smith PetscFunctionBegin; 485290bbb0aSBarry Smith /* 486290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 487290bbb0aSBarry Smith */ 4889566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&xx)); 489290bbb0aSBarry Smith 4909566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 4919566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 492290bbb0aSBarry Smith aj = bmat->j; 493290bbb0aSBarry Smith aa = bmat->a; 494290bbb0aSBarry Smith ii = bmat->i; 495290bbb0aSBarry Smith for (i=0; i<m; i++) { 496290bbb0aSBarry Smith jrow = ii[i]; 497290bbb0aSBarry Smith n = ii[i+1] - jrow; 4989566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->left,yy + bs*i)); 4999566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->middle,xx + bs*i)); 500290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 501290bbb0aSBarry Smith if (aj[jrow] == i) { 5029566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right,xx + bs*aj[jrow])); 5039566063dSJacob Faibussowitsch PetscCall(MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left)); 5049566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 505290bbb0aSBarry Smith jrow++; 506290bbb0aSBarry Smith n--; 507290bbb0aSBarry Smith } 508290bbb0aSBarry Smith for (j=0; j<n; j++) { 5099566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right,xx + bs*aj[jrow])); /* upper triangular part */ 5109566063dSJacob Faibussowitsch PetscCall(MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left)); 5119566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 512290bbb0aSBarry Smith 5139566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right,yy + bs*aj[jrow])); /* lower triangular part */ 5149566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right)); 5159566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 516290bbb0aSBarry Smith jrow++; 517290bbb0aSBarry Smith } 5189566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->left)); 5199566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->middle)); 520290bbb0aSBarry Smith } 5219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xx)); 5229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 523290bbb0aSBarry Smith PetscFunctionReturn(0); 524290bbb0aSBarry Smith } 525290bbb0aSBarry Smith 52681f0254dSBarry Smith static PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 527421e10b8SBarry Smith { 528421e10b8SBarry Smith PetscFunctionBegin; 529421e10b8SBarry Smith PetscFunctionReturn(0); 530421e10b8SBarry Smith } 531421e10b8SBarry Smith 53281f0254dSBarry Smith static PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 533421e10b8SBarry Smith { 534421e10b8SBarry Smith PetscFunctionBegin; 535421e10b8SBarry Smith PetscFunctionReturn(0); 536421e10b8SBarry Smith } 537421e10b8SBarry Smith 53881f0254dSBarry Smith static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 539421e10b8SBarry Smith { 540421e10b8SBarry Smith PetscFunctionBegin; 541421e10b8SBarry Smith PetscFunctionReturn(0); 542421e10b8SBarry Smith } 543421e10b8SBarry Smith 544ccb205f8SBarry Smith /* 545ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 546ccb205f8SBarry Smith */ 54781f0254dSBarry Smith static PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 548ccb205f8SBarry Smith { 549ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 550d0f46423SBarry Smith PetscInt i,j,mbs = A->rmap->n/A->rmap->bs; 551ccb205f8SBarry Smith 552ccb205f8SBarry Smith PetscFunctionBegin; 553ccb205f8SBarry Smith if (!a->diag) { 5549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs,&a->diag)); 555ccb205f8SBarry Smith } 556ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 557ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 558ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 559ccb205f8SBarry Smith if (a->j[j] == i) { 560ccb205f8SBarry Smith a->diag[i] = j; 561ccb205f8SBarry Smith break; 562ccb205f8SBarry Smith } 563ccb205f8SBarry Smith } 564ccb205f8SBarry Smith } 565ccb205f8SBarry Smith PetscFunctionReturn(0); 566ccb205f8SBarry Smith } 567ccb205f8SBarry Smith 5687dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 569ccb205f8SBarry Smith { 570ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 571ccb205f8SBarry Smith Mat_SeqAIJ *c; 572ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 573d2a212eaSBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ailen = a->ilen; 5741620fd73SBarry Smith PetscScalar *a_new; 575ccb205f8SBarry Smith Mat C,*aa = a->a; 576ace3abfcSBarry Smith PetscBool stride,equal; 577ccb205f8SBarry Smith 578ccb205f8SBarry Smith PetscFunctionBegin; 5799566063dSJacob Faibussowitsch PetscCall(ISEqual(isrow,iscol,&equal)); 58028b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for identical column and row indices"); 5819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride)); 58228b400f6SJacob Faibussowitsch PetscCheck(stride,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices"); 5839566063dSJacob Faibussowitsch PetscCall(ISStrideGetInfo(iscol,&first,&step)); 58408401ef6SPierre Jolivet PetscCheck(step == A->rmap->bs,PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block"); 585ccb205f8SBarry Smith 5869566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow,&nrows)); 587ccb205f8SBarry Smith ncols = nrows; 588ccb205f8SBarry Smith 589ccb205f8SBarry Smith /* create submatrix */ 590ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 591ccb205f8SBarry Smith PetscInt n_cols,n_rows; 592ccb205f8SBarry Smith C = *B; 5939566063dSJacob Faibussowitsch PetscCall(MatGetSize(C,&n_rows,&n_cols)); 594aed4548fSBarry Smith PetscCheck(n_rows == nrows && n_cols == ncols,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 5959566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(C)); 596ccb205f8SBarry Smith } else { 5979566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 5989566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE)); 599b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetType(C,MATSEQSBAIJ)); 600b94d7dedSBarry Smith else PetscCall(MatSetType(C,MATSEQAIJ)); 6019566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(C,0,ailen)); 6029566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(C,1,0,ailen)); 603ccb205f8SBarry Smith } 604ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 605ccb205f8SBarry Smith 606ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 607ccb205f8SBarry Smith a_new = c->a; 608ccb205f8SBarry Smith j_new = c->j; 609ccb205f8SBarry Smith i_new = c->i; 610ccb205f8SBarry Smith 611ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 612ccb205f8SBarry Smith lensi = ailen[i]; 613ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 614ccb205f8SBarry Smith *j_new++ = *aj++; 6159566063dSJacob Faibussowitsch PetscCall(MatGetValue(*aa++,first,first,a_new++)); 616ccb205f8SBarry Smith } 617ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 618ccb205f8SBarry Smith c->ilen[i] = lensi; 619ccb205f8SBarry Smith } 620ccb205f8SBarry Smith 6219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 6229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 623ccb205f8SBarry Smith *B = C; 624ccb205f8SBarry Smith PetscFunctionReturn(0); 625ccb205f8SBarry Smith } 626ccb205f8SBarry Smith 62781f0254dSBarry Smith static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 628421e10b8SBarry Smith { 629421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 630421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 631ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 632421e10b8SBarry Smith Mat *aa = a->a,*ap; 633421e10b8SBarry Smith 634421e10b8SBarry Smith PetscFunctionBegin; 635421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 636421e10b8SBarry Smith 637421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 638421e10b8SBarry Smith for (i=1; i<m; i++) { 639421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 640421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 641421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 642421e10b8SBarry Smith if (fshift) { 643421e10b8SBarry Smith ip = aj + ai[i]; 644421e10b8SBarry Smith ap = aa + ai[i]; 645421e10b8SBarry Smith N = ailen[i]; 646421e10b8SBarry Smith for (j=0; j<N; j++) { 647421e10b8SBarry Smith ip[j-fshift] = ip[j]; 648421e10b8SBarry Smith ap[j-fshift] = ap[j]; 649421e10b8SBarry Smith } 650421e10b8SBarry Smith } 651421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 652421e10b8SBarry Smith } 653421e10b8SBarry Smith if (m) { 654421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 655421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 656421e10b8SBarry Smith } 657421e10b8SBarry Smith /* reset ilen and imax for each row */ 658421e10b8SBarry Smith for (i=0; i<m; i++) { 659421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 660421e10b8SBarry Smith } 661421e10b8SBarry Smith a->nz = ai[m]; 662ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 6636bdcaf15SBarry Smith PetscAssert(aa[i],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %" PetscInt_FMT " column %" PetscInt_FMT " nz %" PetscInt_FMT,i,aj[i],a->nz); 6649566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY)); 6659566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY)); 666ccb205f8SBarry Smith } 6679566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; storage space: %" PetscInt_FMT " unneeded,%" PetscInt_FMT " used\n",m,A->cmap->n/A->cmap->bs,fshift,a->nz)); 6689566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n",a->reallocs)); 6699566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Maximum nonzeros in any row is %" PetscInt_FMT "\n",rmax)); 6702205254eSKarl Rupp 6718e58a170SBarry Smith A->info.mallocs += a->reallocs; 672421e10b8SBarry Smith a->reallocs = 0; 673421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 674421e10b8SBarry Smith a->rmax = rmax; 6759566063dSJacob Faibussowitsch PetscCall(MatMarkDiagonal_BlockMat(A)); 676421e10b8SBarry Smith PetscFunctionReturn(0); 677421e10b8SBarry Smith } 678421e10b8SBarry Smith 67981f0254dSBarry Smith static PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg) 680290bbb0aSBarry Smith { 681290bbb0aSBarry Smith PetscFunctionBegin; 6824e0d8c25SBarry Smith if (opt == MAT_SYMMETRIC && flg) { 68341f059aeSBarry Smith A->ops->sor = MatSOR_BlockMat_Symmetric; 684290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 685290bbb0aSBarry Smith } else { 6869566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Unused matrix option %s\n",MatOptions[opt])); 687290bbb0aSBarry Smith } 688290bbb0aSBarry Smith PetscFunctionReturn(0); 689290bbb0aSBarry Smith } 690290bbb0aSBarry Smith 6913964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 692f4259b30SLisandro Dalcin NULL, 693f4259b30SLisandro Dalcin NULL, 694421e10b8SBarry Smith MatMult_BlockMat, 695421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 696421e10b8SBarry Smith MatMultTranspose_BlockMat, 697421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 698f4259b30SLisandro Dalcin NULL, 699f4259b30SLisandro Dalcin NULL, 700f4259b30SLisandro Dalcin NULL, 701f4259b30SLisandro Dalcin /* 10*/ NULL, 702f4259b30SLisandro Dalcin NULL, 703f4259b30SLisandro Dalcin NULL, 70441f059aeSBarry Smith MatSOR_BlockMat, 705f4259b30SLisandro Dalcin NULL, 706f4259b30SLisandro Dalcin /* 15*/ NULL, 707f4259b30SLisandro Dalcin NULL, 708f4259b30SLisandro Dalcin NULL, 709f4259b30SLisandro Dalcin NULL, 710f4259b30SLisandro Dalcin NULL, 711f4259b30SLisandro Dalcin /* 20*/ NULL, 712421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 713290bbb0aSBarry Smith MatSetOption_BlockMat, 714f4259b30SLisandro Dalcin NULL, 715f4259b30SLisandro Dalcin /* 24*/ NULL, 716f4259b30SLisandro Dalcin NULL, 717f4259b30SLisandro Dalcin NULL, 718f4259b30SLisandro Dalcin NULL, 719f4259b30SLisandro Dalcin NULL, 720f4259b30SLisandro Dalcin /* 29*/ NULL, 721f4259b30SLisandro Dalcin NULL, 722f4259b30SLisandro Dalcin NULL, 723f4259b30SLisandro Dalcin NULL, 724f4259b30SLisandro Dalcin NULL, 725f4259b30SLisandro Dalcin /* 34*/ NULL, 726f4259b30SLisandro Dalcin NULL, 727f4259b30SLisandro Dalcin NULL, 728f4259b30SLisandro Dalcin NULL, 729f4259b30SLisandro Dalcin NULL, 730f4259b30SLisandro Dalcin /* 39*/ NULL, 731f4259b30SLisandro Dalcin NULL, 732f4259b30SLisandro Dalcin NULL, 733f4259b30SLisandro Dalcin NULL, 734f4259b30SLisandro Dalcin NULL, 735f4259b30SLisandro Dalcin /* 44*/ NULL, 736f4259b30SLisandro Dalcin NULL, 7377d68702bSBarry Smith MatShift_Basic, 738f4259b30SLisandro Dalcin NULL, 739f4259b30SLisandro Dalcin NULL, 740f4259b30SLisandro Dalcin /* 49*/ NULL, 741f4259b30SLisandro Dalcin NULL, 742f4259b30SLisandro Dalcin NULL, 743f4259b30SLisandro Dalcin NULL, 744f4259b30SLisandro Dalcin NULL, 745f4259b30SLisandro Dalcin /* 54*/ NULL, 746f4259b30SLisandro Dalcin NULL, 747f4259b30SLisandro Dalcin NULL, 748f4259b30SLisandro Dalcin NULL, 749f4259b30SLisandro Dalcin NULL, 7507dae84e0SHong Zhang /* 59*/ MatCreateSubMatrix_BlockMat, 751421e10b8SBarry Smith MatDestroy_BlockMat, 752ccb205f8SBarry Smith MatView_BlockMat, 753f4259b30SLisandro Dalcin NULL, 754f4259b30SLisandro Dalcin NULL, 755f4259b30SLisandro Dalcin /* 64*/ NULL, 756f4259b30SLisandro Dalcin NULL, 757f4259b30SLisandro Dalcin NULL, 758f4259b30SLisandro Dalcin NULL, 759f4259b30SLisandro Dalcin NULL, 760f4259b30SLisandro Dalcin /* 69*/ NULL, 761f4259b30SLisandro Dalcin NULL, 762f4259b30SLisandro Dalcin NULL, 763f4259b30SLisandro Dalcin NULL, 764f4259b30SLisandro Dalcin NULL, 765f4259b30SLisandro Dalcin /* 74*/ NULL, 766f4259b30SLisandro Dalcin NULL, 767f4259b30SLisandro Dalcin NULL, 768f4259b30SLisandro Dalcin NULL, 769f4259b30SLisandro Dalcin NULL, 770f4259b30SLisandro Dalcin /* 79*/ NULL, 771f4259b30SLisandro Dalcin NULL, 772f4259b30SLisandro Dalcin NULL, 773f4259b30SLisandro Dalcin NULL, 7745bba2384SShri Abhyankar MatLoad_BlockMat, 775f4259b30SLisandro Dalcin /* 84*/ NULL, 776f4259b30SLisandro Dalcin NULL, 777f4259b30SLisandro Dalcin NULL, 778f4259b30SLisandro Dalcin NULL, 779f4259b30SLisandro Dalcin NULL, 780f4259b30SLisandro Dalcin /* 89*/ NULL, 781f4259b30SLisandro Dalcin NULL, 782f4259b30SLisandro Dalcin NULL, 783f4259b30SLisandro Dalcin NULL, 784f4259b30SLisandro Dalcin NULL, 785f4259b30SLisandro Dalcin /* 94*/ NULL, 786f4259b30SLisandro Dalcin NULL, 787f4259b30SLisandro Dalcin NULL, 788f4259b30SLisandro Dalcin NULL, 789f4259b30SLisandro Dalcin NULL, 790f4259b30SLisandro Dalcin /* 99*/ NULL, 791f4259b30SLisandro Dalcin NULL, 792f4259b30SLisandro Dalcin NULL, 793f4259b30SLisandro Dalcin NULL, 794f4259b30SLisandro Dalcin NULL, 795f4259b30SLisandro Dalcin /*104*/ NULL, 796f4259b30SLisandro Dalcin NULL, 797f4259b30SLisandro Dalcin NULL, 798f4259b30SLisandro Dalcin NULL, 799f4259b30SLisandro Dalcin NULL, 800f4259b30SLisandro Dalcin /*109*/ NULL, 801f4259b30SLisandro Dalcin NULL, 802f4259b30SLisandro Dalcin NULL, 803f4259b30SLisandro Dalcin NULL, 804f4259b30SLisandro Dalcin NULL, 805f4259b30SLisandro Dalcin /*114*/ NULL, 806f4259b30SLisandro Dalcin NULL, 807f4259b30SLisandro Dalcin NULL, 808f4259b30SLisandro Dalcin NULL, 809f4259b30SLisandro Dalcin NULL, 810f4259b30SLisandro Dalcin /*119*/ NULL, 811f4259b30SLisandro Dalcin NULL, 812f4259b30SLisandro Dalcin NULL, 813f4259b30SLisandro Dalcin NULL, 814f4259b30SLisandro Dalcin NULL, 815f4259b30SLisandro Dalcin /*124*/ NULL, 816f4259b30SLisandro Dalcin NULL, 817f4259b30SLisandro Dalcin NULL, 818f4259b30SLisandro Dalcin NULL, 819f4259b30SLisandro Dalcin NULL, 820f4259b30SLisandro Dalcin /*129*/ NULL, 821f4259b30SLisandro Dalcin NULL, 822f4259b30SLisandro Dalcin NULL, 823f4259b30SLisandro Dalcin NULL, 824f4259b30SLisandro Dalcin NULL, 825f4259b30SLisandro Dalcin /*134*/ NULL, 826f4259b30SLisandro Dalcin NULL, 827f4259b30SLisandro Dalcin NULL, 828f4259b30SLisandro Dalcin NULL, 829f4259b30SLisandro Dalcin NULL, 830f4259b30SLisandro Dalcin /*139*/ NULL, 831f4259b30SLisandro Dalcin NULL, 832d70f29a3SPierre Jolivet NULL, 833d70f29a3SPierre Jolivet NULL, 834d70f29a3SPierre Jolivet NULL, 835d70f29a3SPierre Jolivet /*144*/ NULL, 836d70f29a3SPierre Jolivet NULL, 837d70f29a3SPierre Jolivet NULL, 83899a7f59eSMark Adams NULL, 83999a7f59eSMark Adams NULL, 8407fb60732SBarry Smith NULL, 8417fb60732SBarry Smith /*150*/ NULL 842a5973382SShri Abhyankar }; 843421e10b8SBarry Smith 8448cdf2d9bSBarry Smith /*@C 8458cdf2d9bSBarry Smith MatBlockMatSetPreallocation - For good matrix assembly performance 8468cdf2d9bSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 8478cdf2d9bSBarry Smith (or the array nnz). By setting these parameters accurately, performance 8488cdf2d9bSBarry Smith during matrix assembly can be increased by more than a factor of 50. 8498cdf2d9bSBarry Smith 850d083f849SBarry Smith Collective 8518cdf2d9bSBarry Smith 8528cdf2d9bSBarry Smith Input Parameters: 8538cdf2d9bSBarry Smith + B - The matrix 8548cdf2d9bSBarry Smith . bs - size of each block in matrix 8558cdf2d9bSBarry Smith . nz - number of nonzeros per block row (same for all rows) 8568cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows 8570298fd71SBarry Smith (possibly different for each row) or NULL 8588cdf2d9bSBarry Smith 8598cdf2d9bSBarry Smith Notes: 8608cdf2d9bSBarry Smith If nnz is given then nz is ignored 8618cdf2d9bSBarry Smith 8628cdf2d9bSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 8630298fd71SBarry Smith Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 8648cdf2d9bSBarry Smith allocation. For large problems you MUST preallocate memory or you 8658cdf2d9bSBarry Smith will get TERRIBLE performance, see the users' manual chapter on matrices. 8668cdf2d9bSBarry Smith 8678cdf2d9bSBarry Smith Level: intermediate 8688cdf2d9bSBarry Smith 869db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateBlockMat()`, `MatSetValues()` 8708cdf2d9bSBarry Smith 8718cdf2d9bSBarry Smith @*/ 8727087cfbeSBarry Smith PetscErrorCode MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 8738cdf2d9bSBarry Smith { 8748cdf2d9bSBarry Smith PetscFunctionBegin; 875cac4c232SBarry Smith PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz)); 8768cdf2d9bSBarry Smith PetscFunctionReturn(0); 8778cdf2d9bSBarry Smith } 8788cdf2d9bSBarry Smith 87981f0254dSBarry Smith static PetscErrorCode MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz) 8808cdf2d9bSBarry Smith { 8818cdf2d9bSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 8828cdf2d9bSBarry Smith PetscInt i; 8838cdf2d9bSBarry Smith 8848cdf2d9bSBarry Smith PetscFunctionBegin; 8859566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->rmap,bs)); 8869566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->cmap,bs)); 8879566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 8889566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 8899566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetBlockSize(A->rmap,&bs)); 89034ef9618SShri Abhyankar 8918cdf2d9bSBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 89208401ef6SPierre Jolivet PetscCheck(nz >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %" PetscInt_FMT,nz); 8938cdf2d9bSBarry Smith if (nnz) { 894d0f46423SBarry Smith for (i=0; i<A->rmap->n/bs; i++) { 89508401ef6SPierre Jolivet PetscCheck(nnz[i] >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT,i,nnz[i]); 89608401ef6SPierre Jolivet PetscCheck(nnz[i] <= A->cmap->n/bs,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " rowlength %" PetscInt_FMT,i,nnz[i],A->cmap->n/bs); 8978cdf2d9bSBarry Smith } 8988cdf2d9bSBarry Smith } 899d0f46423SBarry Smith bmat->mbs = A->rmap->n/bs; 9008cdf2d9bSBarry Smith 9019566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->right)); 9029566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->middle)); 9039566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left)); 9048cdf2d9bSBarry Smith 9052ee49352SLisandro Dalcin if (!bmat->imax) { 9069566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->rmap->n,&bmat->imax,A->rmap->n,&bmat->ilen)); 9079566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)A,2*A->rmap->n*sizeof(PetscInt))); 9082ee49352SLisandro Dalcin } 909c0aa6a63SJacob Faibussowitsch if (PetscLikely(nnz)) { 9108cdf2d9bSBarry Smith nz = 0; 911d0f46423SBarry Smith for (i=0; i<A->rmap->n/A->rmap->bs; i++) { 9128cdf2d9bSBarry Smith bmat->imax[i] = nnz[i]; 9138cdf2d9bSBarry Smith nz += nnz[i]; 9148cdf2d9bSBarry Smith } 915f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation"); 9168cdf2d9bSBarry Smith 9178cdf2d9bSBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 9189566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(bmat->ilen,bmat->mbs)); 9198cdf2d9bSBarry Smith 9208cdf2d9bSBarry Smith /* allocate the matrix space */ 9219566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i)); 9229566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nz,&bmat->a,nz,&bmat->j,A->rmap->n+1,&bmat->i)); 9239566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)))); 9248cdf2d9bSBarry Smith bmat->i[0] = 0; 9258cdf2d9bSBarry Smith for (i=1; i<bmat->mbs+1; i++) { 9268cdf2d9bSBarry Smith bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 9278cdf2d9bSBarry Smith } 9288cdf2d9bSBarry Smith bmat->singlemalloc = PETSC_TRUE; 9298cdf2d9bSBarry Smith bmat->free_a = PETSC_TRUE; 9308cdf2d9bSBarry Smith bmat->free_ij = PETSC_TRUE; 9318cdf2d9bSBarry Smith 9328cdf2d9bSBarry Smith bmat->nz = 0; 9338cdf2d9bSBarry Smith bmat->maxnz = nz; 9348cdf2d9bSBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 9359566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 9368cdf2d9bSBarry Smith PetscFunctionReturn(0); 9378cdf2d9bSBarry Smith } 9388cdf2d9bSBarry Smith 939421e10b8SBarry Smith /*MC 940421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 941421e10b8SBarry Smith consisting of (usually) sparse blocks. 942421e10b8SBarry Smith 943421e10b8SBarry Smith Level: advanced 944421e10b8SBarry Smith 945db781477SPatrick Sanan .seealso: `MatCreateBlockMat()` 946421e10b8SBarry Smith 947421e10b8SBarry Smith M*/ 948421e10b8SBarry Smith 9498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A) 950421e10b8SBarry Smith { 951421e10b8SBarry Smith Mat_BlockMat *b; 952421e10b8SBarry Smith 953421e10b8SBarry Smith PetscFunctionBegin; 9549566063dSJacob Faibussowitsch PetscCall(PetscNewLog(A,&b)); 955421e10b8SBarry Smith A->data = (void*)b; 9569566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps))); 957421e10b8SBarry Smith 958421e10b8SBarry Smith A->assembled = PETSC_TRUE; 959421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 9609566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT)); 961290bbb0aSBarry Smith 9629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatBlockMatSetPreallocation_C",MatBlockMatSetPreallocation_BlockMat)); 963421e10b8SBarry Smith PetscFunctionReturn(0); 964421e10b8SBarry Smith } 965421e10b8SBarry Smith 966421e10b8SBarry Smith /*@C 967bbd3f336SJed Brown MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential Mat object 968421e10b8SBarry Smith 969d083f849SBarry Smith Collective 970421e10b8SBarry Smith 971421e10b8SBarry Smith Input Parameters: 972421e10b8SBarry Smith + comm - MPI communicator 973421e10b8SBarry Smith . m - number of rows 974421e10b8SBarry Smith . n - number of columns 9758cdf2d9bSBarry Smith . bs - size of each submatrix 9768cdf2d9bSBarry Smith . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 9770298fd71SBarry Smith - nnz - expected number of nonzers per block row if known (use NULL otherwise) 978421e10b8SBarry Smith 979421e10b8SBarry Smith Output Parameter: 980421e10b8SBarry Smith . A - the matrix 981421e10b8SBarry Smith 982421e10b8SBarry Smith Level: intermediate 983421e10b8SBarry Smith 98495452b02SPatrick Sanan Notes: 98595452b02SPatrick Sanan Matrices of this type are nominally-sparse matrices in which each "entry" is a Mat object. Each Mat must 986bbd3f336SJed Brown have the same size and be sequential. The local and global sizes must be compatible with this decomposition. 987bbd3f336SJed Brown 988bbd3f336SJed Brown For matrices containing parallel submatrices and variable block sizes, see MATNEST. 989421e10b8SBarry Smith 990db781477SPatrick Sanan .seealso: `MATBLOCKMAT`, `MatCreateNest()` 991421e10b8SBarry Smith @*/ 9927087cfbeSBarry Smith PetscErrorCode MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A) 993421e10b8SBarry Smith { 994421e10b8SBarry Smith PetscFunctionBegin; 9959566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,A)); 9969566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE)); 9979566063dSJacob Faibussowitsch PetscCall(MatSetType(*A,MATBLOCKMAT)); 9989566063dSJacob Faibussowitsch PetscCall(MatBlockMatSetPreallocation(*A,bs,nz,nnz)); 999421e10b8SBarry Smith PetscFunctionReturn(0); 1000421e10b8SBarry Smith } 1001