xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
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");
3626fbe8dcSKarl Rupp   if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) {
37e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
3826fbe8dcSKarl Rupp   }
39290bbb0aSBarry Smith 
40290bbb0aSBarry Smith   if (!a->diags) {
419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mbs,&a->diags));
429566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
43290bbb0aSBarry Smith     for (i=0; i<mbs; i++) {
449566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col));
459566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info));
469566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info));
479566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&row));
489566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&col));
49290bbb0aSBarry Smith     }
509566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(bb,&a->workb));
51290bbb0aSBarry Smith   }
52290bbb0aSBarry Smith   diag = a->diags;
53290bbb0aSBarry Smith 
549566063dSJacob Faibussowitsch   PetscCall(VecSet(xx,0.0));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx,&x));
56290bbb0aSBarry Smith   /* copy right hand side because it must be modified during iteration */
579566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb,a->workb));
589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(a->workb,&b));
59290bbb0aSBarry Smith 
6041f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
61290bbb0aSBarry Smith   while (its--) {
62290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
63290bbb0aSBarry Smith 
64290bbb0aSBarry Smith       for (i=0; i<mbs; i++) {
656e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
666e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
676e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
68290bbb0aSBarry Smith 
699566063dSJacob Faibussowitsch         PetscCall(VecSet(left,0.0));
70290bbb0aSBarry Smith         for (j=0; j<n; j++) {
719566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(right,x + idx[j]*bs));
729566063dSJacob Faibussowitsch           PetscCall(MatMultAdd(v[j],right,left,left));
739566063dSJacob Faibussowitsch           PetscCall(VecResetArray(right));
74290bbb0aSBarry Smith         }
759566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,b + i*bs));
769566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left,-1.0,right));
779566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
78290bbb0aSBarry Smith 
799566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,x + i*bs));
809566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i],left,right));
816e63c7a1SBarry Smith 
8241f059aeSBarry Smith         /* now adjust right hand side, see MatSOR_SeqSBAIJ */
836e63c7a1SBarry Smith         for (j=0; j<n; j++) {
849566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(v[j],right,left));
859566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(middle,b + idx[j]*bs));
869566063dSJacob Faibussowitsch           PetscCall(VecAXPY(middle,-1.0,left));
879566063dSJacob Faibussowitsch           PetscCall(VecResetArray(middle));
886e63c7a1SBarry Smith         }
899566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
906e63c7a1SBarry Smith 
91290bbb0aSBarry Smith       }
92290bbb0aSBarry Smith     }
93290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
94290bbb0aSBarry Smith 
95290bbb0aSBarry Smith       for (i=mbs-1; i>=0; i--) {
966e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
976e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
986e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
99290bbb0aSBarry Smith 
1009566063dSJacob Faibussowitsch         PetscCall(VecSet(left,0.0));
101290bbb0aSBarry Smith         for (j=0; j<n; j++) {
1029566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(right,x + idx[j]*bs));
1039566063dSJacob Faibussowitsch           PetscCall(MatMultAdd(v[j],right,left,left));
1049566063dSJacob Faibussowitsch           PetscCall(VecResetArray(right));
105290bbb0aSBarry Smith         }
1069566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,b + i*bs));
1079566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left,-1.0,right));
1089566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
109290bbb0aSBarry Smith 
1109566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,x + i*bs));
1119566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i],left,right));
1129566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
113290bbb0aSBarry Smith 
114290bbb0aSBarry Smith       }
115290bbb0aSBarry Smith     }
116290bbb0aSBarry Smith   }
1179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx,&x));
1189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(a->workb,&b));
119290bbb0aSBarry Smith   PetscFunctionReturn(0);
120290bbb0aSBarry Smith }
121290bbb0aSBarry Smith 
12281f0254dSBarry Smith static PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
123ccb205f8SBarry Smith {
124ccb205f8SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
125ccb205f8SBarry Smith   PetscScalar       *x;
126a2ea699eSBarry Smith   const Mat         *v;
127ccb205f8SBarry Smith   const PetscScalar *b;
128d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
129ccb205f8SBarry Smith   const PetscInt    *idx;
130ccb205f8SBarry Smith   IS                row,col;
131ccb205f8SBarry Smith   MatFactorInfo     info;
132ccb205f8SBarry Smith   Vec               left = a->left,right = a->right;
133ccb205f8SBarry Smith   Mat               *diag;
134ccb205f8SBarry Smith 
135ccb205f8SBarry Smith   PetscFunctionBegin;
136ccb205f8SBarry Smith   its = its*lits;
13708401ef6SPierre 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);
138aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT),PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
13908401ef6SPierre Jolivet   PetscCheck(omega == 1.0,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
14028b400f6SJacob Faibussowitsch   PetscCheck(!fshift,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
141ccb205f8SBarry Smith 
142ccb205f8SBarry Smith   if (!a->diags) {
1439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mbs,&a->diags));
1449566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
145ccb205f8SBarry Smith     for (i=0; i<mbs; i++) {
1469566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col));
1479566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info));
1489566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info));
1499566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&row));
1509566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&col));
151ccb205f8SBarry Smith     }
152ccb205f8SBarry Smith   }
153ccb205f8SBarry Smith   diag = a->diags;
154ccb205f8SBarry Smith 
1559566063dSJacob Faibussowitsch   PetscCall(VecSet(xx,0.0));
1569566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx,&x));
1579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb,&b));
158ccb205f8SBarry Smith 
15941f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
160ccb205f8SBarry Smith   while (its--) {
161ccb205f8SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
162ccb205f8SBarry Smith 
163ccb205f8SBarry Smith       for (i=0; i<mbs; i++) {
164ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
165ccb205f8SBarry Smith         idx = a->j + a->i[i];
166ccb205f8SBarry Smith         v   = a->a + a->i[i];
167ccb205f8SBarry Smith 
1689566063dSJacob Faibussowitsch         PetscCall(VecSet(left,0.0));
169ccb205f8SBarry Smith         for (j=0; j<n; j++) {
170ccb205f8SBarry Smith           if (idx[j] != i) {
1719566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(right,x + idx[j]*bs));
1729566063dSJacob Faibussowitsch             PetscCall(MatMultAdd(v[j],right,left,left));
1739566063dSJacob Faibussowitsch             PetscCall(VecResetArray(right));
174ccb205f8SBarry Smith           }
175ccb205f8SBarry Smith         }
1769566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,b + i*bs));
1779566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left,-1.0,right));
1789566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
179ccb205f8SBarry Smith 
1809566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,x + i*bs));
1819566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i],left,right));
1829566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
183ccb205f8SBarry Smith       }
184ccb205f8SBarry Smith     }
185ccb205f8SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
186ccb205f8SBarry Smith 
187ccb205f8SBarry Smith       for (i=mbs-1; i>=0; i--) {
188ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
189ccb205f8SBarry Smith         idx = a->j + a->i[i];
190ccb205f8SBarry Smith         v   = a->a + a->i[i];
191ccb205f8SBarry Smith 
1929566063dSJacob Faibussowitsch         PetscCall(VecSet(left,0.0));
193ccb205f8SBarry Smith         for (j=0; j<n; j++) {
194ccb205f8SBarry Smith           if (idx[j] != i) {
1959566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(right,x + idx[j]*bs));
1969566063dSJacob Faibussowitsch             PetscCall(MatMultAdd(v[j],right,left,left));
1979566063dSJacob Faibussowitsch             PetscCall(VecResetArray(right));
198ccb205f8SBarry Smith           }
199ccb205f8SBarry Smith         }
2009566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,b + i*bs));
2019566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left,-1.0,right));
2029566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
203ccb205f8SBarry Smith 
2049566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right,x + i*bs));
2059566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i],left,right));
2069566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
207ccb205f8SBarry Smith 
208ccb205f8SBarry Smith       }
209ccb205f8SBarry Smith     }
210ccb205f8SBarry Smith   }
2119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx,&x));
2129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb,&b));
213ccb205f8SBarry Smith   PetscFunctionReturn(0);
214ccb205f8SBarry Smith }
215ccb205f8SBarry Smith 
21681f0254dSBarry Smith static PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
217421e10b8SBarry Smith {
218421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
219421e10b8SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
220421e10b8SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
221d0f46423SBarry Smith   PetscInt       *aj  =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
222421e10b8SBarry Smith   PetscInt       ridx,cidx;
223ace3abfcSBarry Smith   PetscBool      roworiented=a->roworiented;
224421e10b8SBarry Smith   MatScalar      value;
225421e10b8SBarry Smith   Mat            *ap,*aa = a->a;
226421e10b8SBarry Smith 
227421e10b8SBarry Smith   PetscFunctionBegin;
228421e10b8SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
229421e10b8SBarry Smith     row  = im[k];
230421e10b8SBarry Smith     brow = row/bs;
231421e10b8SBarry Smith     if (row < 0) continue;
2326bdcaf15SBarry 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);
233421e10b8SBarry Smith     rp   = aj + ai[brow];
234421e10b8SBarry Smith     ap   = aa + ai[brow];
235421e10b8SBarry Smith     rmax = imax[brow];
236421e10b8SBarry Smith     nrow = ailen[brow];
237421e10b8SBarry Smith     low  = 0;
238421e10b8SBarry Smith     high = nrow;
239421e10b8SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
240421e10b8SBarry Smith       if (in[l] < 0) continue;
2416bdcaf15SBarry 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);
242421e10b8SBarry Smith       col = in[l]; bcol = col/bs;
243*b94d7dedSBarry Smith       if (A->symmetric == PETSC_BOOL3_TRUE && brow > bcol) continue;
244421e10b8SBarry Smith       ridx = row % bs; cidx = col % bs;
2452205254eSKarl Rupp       if (roworiented) value = v[l + k*n];
2462205254eSKarl Rupp       else value = v[k + l*m];
2472205254eSKarl Rupp 
2482205254eSKarl Rupp       if (col <= lastcol) low = 0;
2492205254eSKarl Rupp       else high = nrow;
250421e10b8SBarry Smith       lastcol = col;
251421e10b8SBarry Smith       while (high-low > 7) {
252421e10b8SBarry Smith         t = (low+high)/2;
253421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
254421e10b8SBarry Smith         else              low  = t;
255421e10b8SBarry Smith       }
256421e10b8SBarry Smith       for (i=low; i<high; i++) {
257421e10b8SBarry Smith         if (rp[i] > bcol) break;
2582205254eSKarl Rupp         if (rp[i] == bcol) goto noinsert1;
259421e10b8SBarry Smith       }
260421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
26108401ef6SPierre Jolivet       PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
262fef13f97SBarry Smith       MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
263421e10b8SBarry Smith       N = nrow++ - 1; high++;
264421e10b8SBarry Smith       /* shift up all the later entries in this row */
265421e10b8SBarry Smith       for (ii=N; ii>=i; ii--) {
266421e10b8SBarry Smith         rp[ii+1] = rp[ii];
267421e10b8SBarry Smith         ap[ii+1] = ap[ii];
268421e10b8SBarry Smith       }
269f4259b30SLisandro Dalcin       if (N>=i) ap[i] = NULL;
270421e10b8SBarry Smith       rp[i] = bcol;
271421e10b8SBarry Smith       a->nz++;
272e56f5c9eSBarry Smith       A->nonzerostate++;
273421e10b8SBarry Smith noinsert1:;
274421e10b8SBarry Smith       if (!*(ap+i)) {
2759566063dSJacob Faibussowitsch         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,NULL,ap+i));
276b0223207SBarry Smith       }
2779566063dSJacob Faibussowitsch       PetscCall(MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is));
278421e10b8SBarry Smith       low  = i;
279421e10b8SBarry Smith     }
280421e10b8SBarry Smith     ailen[brow] = nrow;
281421e10b8SBarry Smith   }
282421e10b8SBarry Smith   PetscFunctionReturn(0);
283421e10b8SBarry Smith }
284421e10b8SBarry Smith 
28581f0254dSBarry Smith static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
286a5973382SShri Abhyankar {
287a5973382SShri Abhyankar   Mat               tmpA;
288a5973382SShri Abhyankar   PetscInt          i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0;
289a5973382SShri Abhyankar   const PetscInt    *cols;
290a5973382SShri Abhyankar   const PetscScalar *values;
291ace3abfcSBarry Smith   PetscBool         flg = PETSC_FALSE,notdone;
292a5973382SShri Abhyankar   Mat_SeqAIJ        *a;
293a5973382SShri Abhyankar   Mat_BlockMat      *amat;
294a5973382SShri Abhyankar 
295a5973382SShri Abhyankar   PetscFunctionBegin;
296c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
2979566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
2989566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&tmpA));
2999566063dSJacob Faibussowitsch   PetscCall(MatSetType(tmpA,MATSEQAIJ));
3009566063dSJacob Faibussowitsch   PetscCall(MatLoad_SeqAIJ(tmpA,viewer));
301a5973382SShri Abhyankar 
3029566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(tmpA,&m,&n));
303d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF,NULL,"Options for loading BlockMat matrix 1","Mat");
3049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL));
3059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,NULL));
306d0609cedSBarry Smith   PetscOptionsEnd();
307a5973382SShri Abhyankar 
308a5973382SShri Abhyankar   /* Determine number of nonzero blocks for each block row */
309a5973382SShri Abhyankar   a    = (Mat_SeqAIJ*) tmpA->data;
310a5973382SShri Abhyankar   mbs  = m/bs;
3119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs,&lens,bs,&ii,bs,&ilens));
3129566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lens,mbs));
313a5973382SShri Abhyankar 
314a5973382SShri Abhyankar   for (i=0; i<mbs; i++) {
315a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
316a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
317a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
318a5973382SShri Abhyankar     }
319a5973382SShri Abhyankar 
320a5973382SShri Abhyankar     currentcol = -1;
321a5973382SShri Abhyankar     while (PETSC_TRUE) {
322a5973382SShri Abhyankar       notdone = PETSC_FALSE;
323a5973382SShri Abhyankar       nextcol = 1000000000;
324a5973382SShri Abhyankar       for (j=0; j<bs; j++) {
325a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) {
326a5973382SShri Abhyankar           ii[j]++;
327a5973382SShri Abhyankar           ilens[j]--;
328a5973382SShri Abhyankar         }
329a5973382SShri Abhyankar         if (ilens[j] > 0) {
330a5973382SShri Abhyankar           notdone = PETSC_TRUE;
331a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
332a5973382SShri Abhyankar         }
333a5973382SShri Abhyankar       }
334a5973382SShri Abhyankar       if (!notdone) break;
335a5973382SShri Abhyankar       if (!flg || (nextcol >= i)) lens[i]++;
336a5973382SShri Abhyankar       currentcol = nextcol;
337a5973382SShri Abhyankar     }
338a5973382SShri Abhyankar   }
339a5973382SShri Abhyankar 
340a5973382SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
3419566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE));
342a5973382SShri Abhyankar   }
3439566063dSJacob Faibussowitsch   PetscCall(MatBlockMatSetPreallocation(newmat,bs,0,lens));
3441baa6e33SBarry Smith   if (flg) PetscCall(MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE));
345a5973382SShri Abhyankar   amat = (Mat_BlockMat*)(newmat)->data;
346a5973382SShri Abhyankar 
347a5973382SShri Abhyankar   /* preallocate the submatrices */
3489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs,&llens));
349a5973382SShri Abhyankar   for (i=0; i<mbs; i++) { /* loops for block rows */
350a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
351a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
352a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
353a5973382SShri Abhyankar     }
354a5973382SShri Abhyankar 
355a5973382SShri Abhyankar     currentcol = 1000000000;
356a5973382SShri Abhyankar     for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */
357a5973382SShri Abhyankar       if (ilens[j] > 0) {
358a5973382SShri Abhyankar         currentcol = PetscMin(currentcol,ii[j][0]/bs);
359a5973382SShri Abhyankar       }
360a5973382SShri Abhyankar     }
361a5973382SShri Abhyankar 
362a5973382SShri Abhyankar     while (PETSC_TRUE) {  /* loops over blocks in block row */
363a5973382SShri Abhyankar       notdone = PETSC_FALSE;
364a5973382SShri Abhyankar       nextcol = 1000000000;
3659566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(llens,bs));
366a5973382SShri Abhyankar       for (j=0; j<bs; j++) { /* loop over rows in block */
367a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */
368a5973382SShri Abhyankar           ii[j]++;
369a5973382SShri Abhyankar           ilens[j]--;
370a5973382SShri Abhyankar           llens[j]++;
371a5973382SShri Abhyankar         }
372a5973382SShri Abhyankar         if (ilens[j] > 0) {
373a5973382SShri Abhyankar           notdone = PETSC_TRUE;
374a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
375a5973382SShri Abhyankar         }
376a5973382SShri Abhyankar       }
37708401ef6SPierre Jolivet       PetscCheck(cnt < amat->maxnz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %" PetscInt_FMT,cnt);
378a5973382SShri Abhyankar       if (!flg || currentcol >= i) {
379a5973382SShri Abhyankar         amat->j[cnt] = currentcol;
3809566063dSJacob Faibussowitsch         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++));
381a5973382SShri Abhyankar       }
382a5973382SShri Abhyankar 
383a5973382SShri Abhyankar       if (!notdone) break;
384a5973382SShri Abhyankar       currentcol = nextcol;
385a5973382SShri Abhyankar     }
386a5973382SShri Abhyankar     amat->ilen[i] = lens[i];
387a5973382SShri Abhyankar   }
388a5973382SShri Abhyankar 
3899566063dSJacob Faibussowitsch   PetscCall(PetscFree3(lens,ii,ilens));
3909566063dSJacob Faibussowitsch   PetscCall(PetscFree(llens));
391a5973382SShri Abhyankar 
392a5973382SShri Abhyankar   /* copy over the matrix, one row at a time */
393a5973382SShri Abhyankar   for (i=0; i<m; i++) {
3949566063dSJacob Faibussowitsch     PetscCall(MatGetRow(tmpA,i,&ncols,&cols,&values));
3959566063dSJacob Faibussowitsch     PetscCall(MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES));
3969566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(tmpA,i,&ncols,&cols,&values));
397a5973382SShri Abhyankar   }
3989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY));
3999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY));
400a5973382SShri Abhyankar   PetscFunctionReturn(0);
401a5973382SShri Abhyankar }
402a5973382SShri Abhyankar 
40381f0254dSBarry Smith static PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
404ccb205f8SBarry Smith {
40564075487SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
406ccb205f8SBarry Smith   const char        *name;
407ccb205f8SBarry Smith   PetscViewerFormat format;
408ccb205f8SBarry Smith 
409ccb205f8SBarry Smith   PetscFunctionBegin;
4109566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)A,&name));
4119566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer,&format));
412ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
4139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %" PetscInt_FMT " \n",a->nz));
414*b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n"));
415ccb205f8SBarry Smith   }
416ccb205f8SBarry Smith   PetscFunctionReturn(0);
417ccb205f8SBarry Smith }
418421e10b8SBarry Smith 
41981f0254dSBarry Smith static PetscErrorCode MatDestroy_BlockMat(Mat mat)
420421e10b8SBarry Smith {
421421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
422ccb205f8SBarry Smith   PetscInt       i;
423421e10b8SBarry Smith 
424421e10b8SBarry Smith   PetscFunctionBegin;
4259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->right));
4269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->left));
4279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->middle));
4289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->workb));
429ccb205f8SBarry Smith   if (bmat->diags) {
430d0f46423SBarry Smith     for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) {
4319566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&bmat->diags[i]));
432ccb205f8SBarry Smith     }
433ccb205f8SBarry Smith   }
434ccb205f8SBarry Smith   if (bmat->a) {
435ccb205f8SBarry Smith     for (i=0; i<bmat->nz; i++) {
4369566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&bmat->a[i]));
437ccb205f8SBarry Smith     }
438ccb205f8SBarry Smith   }
4399566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i));
4409566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
441421e10b8SBarry Smith   PetscFunctionReturn(0);
442421e10b8SBarry Smith }
443421e10b8SBarry Smith 
44481f0254dSBarry Smith static PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
445421e10b8SBarry Smith {
446421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
447421e10b8SBarry Smith   PetscScalar    *xx,*yy;
448d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
449421e10b8SBarry Smith   Mat            *aa;
450421e10b8SBarry Smith 
451421e10b8SBarry Smith   PetscFunctionBegin;
452421e10b8SBarry Smith   /*
453421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
454421e10b8SBarry Smith   */
4559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&xx));
456ccb205f8SBarry Smith 
4579566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
4589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y,&yy));
459421e10b8SBarry Smith   aj   = bmat->j;
460421e10b8SBarry Smith   aa   = bmat->a;
461421e10b8SBarry Smith   ii   = bmat->i;
462421e10b8SBarry Smith   for (i=0; i<m; i++) {
463421e10b8SBarry Smith     jrow = ii[i];
4649566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left,yy + bs*i));
465421e10b8SBarry Smith     n    = ii[i+1] - jrow;
466421e10b8SBarry Smith     for (j=0; j<n; j++) {
4679566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right,xx + bs*aj[jrow]));
4689566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left));
4699566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
470421e10b8SBarry Smith       jrow++;
471421e10b8SBarry Smith     }
4729566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
473421e10b8SBarry Smith   }
4749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&xx));
4759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y,&yy));
476421e10b8SBarry Smith   PetscFunctionReturn(0);
477421e10b8SBarry Smith }
478421e10b8SBarry Smith 
479290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
480290bbb0aSBarry Smith {
481290bbb0aSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
482290bbb0aSBarry Smith   PetscScalar    *xx,*yy;
483d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
484290bbb0aSBarry Smith   Mat            *aa;
485290bbb0aSBarry Smith 
486290bbb0aSBarry Smith   PetscFunctionBegin;
487290bbb0aSBarry Smith   /*
488290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
489290bbb0aSBarry Smith   */
4909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&xx));
491290bbb0aSBarry Smith 
4929566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
4939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y,&yy));
494290bbb0aSBarry Smith   aj   = bmat->j;
495290bbb0aSBarry Smith   aa   = bmat->a;
496290bbb0aSBarry Smith   ii   = bmat->i;
497290bbb0aSBarry Smith   for (i=0; i<m; i++) {
498290bbb0aSBarry Smith     jrow = ii[i];
499290bbb0aSBarry Smith     n    = ii[i+1] - jrow;
5009566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left,yy + bs*i));
5019566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->middle,xx + bs*i));
502290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
503290bbb0aSBarry Smith     if (aj[jrow] == i) {
5049566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right,xx + bs*aj[jrow]));
5059566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left));
5069566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
507290bbb0aSBarry Smith       jrow++;
508290bbb0aSBarry Smith       n--;
509290bbb0aSBarry Smith     }
510290bbb0aSBarry Smith     for (j=0; j<n; j++) {
5119566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right,xx + bs*aj[jrow]));            /* upper triangular part */
5129566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left));
5139566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
514290bbb0aSBarry Smith 
5159566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right,yy + bs*aj[jrow]));            /* lower triangular part */
5169566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right));
5179566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
518290bbb0aSBarry Smith       jrow++;
519290bbb0aSBarry Smith     }
5209566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
5219566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->middle));
522290bbb0aSBarry Smith   }
5239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&xx));
5249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y,&yy));
525290bbb0aSBarry Smith   PetscFunctionReturn(0);
526290bbb0aSBarry Smith }
527290bbb0aSBarry Smith 
52881f0254dSBarry Smith static PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
529421e10b8SBarry Smith {
530421e10b8SBarry Smith   PetscFunctionBegin;
531421e10b8SBarry Smith   PetscFunctionReturn(0);
532421e10b8SBarry Smith }
533421e10b8SBarry Smith 
53481f0254dSBarry Smith static PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
535421e10b8SBarry Smith {
536421e10b8SBarry Smith   PetscFunctionBegin;
537421e10b8SBarry Smith   PetscFunctionReturn(0);
538421e10b8SBarry Smith }
539421e10b8SBarry Smith 
54081f0254dSBarry Smith static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
541421e10b8SBarry Smith {
542421e10b8SBarry Smith   PetscFunctionBegin;
543421e10b8SBarry Smith   PetscFunctionReturn(0);
544421e10b8SBarry Smith }
545421e10b8SBarry Smith 
546ccb205f8SBarry Smith /*
547ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
548ccb205f8SBarry Smith */
54981f0254dSBarry Smith static PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
550ccb205f8SBarry Smith {
551ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
552d0f46423SBarry Smith   PetscInt       i,j,mbs = A->rmap->n/A->rmap->bs;
553ccb205f8SBarry Smith 
554ccb205f8SBarry Smith   PetscFunctionBegin;
555ccb205f8SBarry Smith   if (!a->diag) {
5569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mbs,&a->diag));
557ccb205f8SBarry Smith   }
558ccb205f8SBarry Smith   for (i=0; i<mbs; i++) {
559ccb205f8SBarry Smith     a->diag[i] = a->i[i+1];
560ccb205f8SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
561ccb205f8SBarry Smith       if (a->j[j] == i) {
562ccb205f8SBarry Smith         a->diag[i] = j;
563ccb205f8SBarry Smith         break;
564ccb205f8SBarry Smith       }
565ccb205f8SBarry Smith     }
566ccb205f8SBarry Smith   }
567ccb205f8SBarry Smith   PetscFunctionReturn(0);
568ccb205f8SBarry Smith }
569ccb205f8SBarry Smith 
5707dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
571ccb205f8SBarry Smith {
572ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
573ccb205f8SBarry Smith   Mat_SeqAIJ     *c;
574ccb205f8SBarry Smith   PetscInt       i,k,first,step,lensi,nrows,ncols;
575d2a212eaSBarry Smith   PetscInt       *j_new,*i_new,*aj = a->j,*ailen = a->ilen;
5761620fd73SBarry Smith   PetscScalar    *a_new;
577ccb205f8SBarry Smith   Mat            C,*aa = a->a;
578ace3abfcSBarry Smith   PetscBool      stride,equal;
579ccb205f8SBarry Smith 
580ccb205f8SBarry Smith   PetscFunctionBegin;
5819566063dSJacob Faibussowitsch   PetscCall(ISEqual(isrow,iscol,&equal));
58228b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for identical column and row indices");
5839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride));
58428b400f6SJacob Faibussowitsch   PetscCheck(stride,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices");
5859566063dSJacob Faibussowitsch   PetscCall(ISStrideGetInfo(iscol,&first,&step));
58608401ef6SPierre Jolivet   PetscCheck(step == A->rmap->bs,PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block");
587ccb205f8SBarry Smith 
5889566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow,&nrows));
589ccb205f8SBarry Smith   ncols = nrows;
590ccb205f8SBarry Smith 
591ccb205f8SBarry Smith   /* create submatrix */
592ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
593ccb205f8SBarry Smith     PetscInt n_cols,n_rows;
594ccb205f8SBarry Smith     C    = *B;
5959566063dSJacob Faibussowitsch     PetscCall(MatGetSize(C,&n_rows,&n_cols));
596aed4548fSBarry Smith     PetscCheck(n_rows == nrows && n_cols == ncols,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
5979566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(C));
598ccb205f8SBarry Smith   } else {
5999566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C));
6009566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE));
601*b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetType(C,MATSEQSBAIJ));
602*b94d7dedSBarry Smith     else PetscCall(MatSetType(C,MATSEQAIJ));
6039566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(C,0,ailen));
6049566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(C,1,0,ailen));
605ccb205f8SBarry Smith   }
606ccb205f8SBarry Smith   c = (Mat_SeqAIJ*)C->data;
607ccb205f8SBarry Smith 
608ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
609ccb205f8SBarry Smith   a_new = c->a;
610ccb205f8SBarry Smith   j_new = c->j;
611ccb205f8SBarry Smith   i_new = c->i;
612ccb205f8SBarry Smith 
613ccb205f8SBarry Smith   for (i=0; i<nrows; i++) {
614ccb205f8SBarry Smith     lensi = ailen[i];
615ccb205f8SBarry Smith     for (k=0; k<lensi; k++) {
616ccb205f8SBarry Smith       *j_new++ = *aj++;
6179566063dSJacob Faibussowitsch       PetscCall(MatGetValue(*aa++,first,first,a_new++));
618ccb205f8SBarry Smith     }
619ccb205f8SBarry Smith     i_new[i+1] = i_new[i] + lensi;
620ccb205f8SBarry Smith     c->ilen[i] = lensi;
621ccb205f8SBarry Smith   }
622ccb205f8SBarry Smith 
6239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
6249566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
625ccb205f8SBarry Smith   *B   = C;
626ccb205f8SBarry Smith   PetscFunctionReturn(0);
627ccb205f8SBarry Smith }
628ccb205f8SBarry Smith 
62981f0254dSBarry Smith static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
630421e10b8SBarry Smith {
631421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
632421e10b8SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
633ccb205f8SBarry Smith   PetscInt       m      = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
634421e10b8SBarry Smith   Mat            *aa    = a->a,*ap;
635421e10b8SBarry Smith 
636421e10b8SBarry Smith   PetscFunctionBegin;
637421e10b8SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
638421e10b8SBarry Smith 
639421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
640421e10b8SBarry Smith   for (i=1; i<m; i++) {
641421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
642421e10b8SBarry Smith     fshift += imax[i-1] - ailen[i-1];
643421e10b8SBarry Smith     rmax    = PetscMax(rmax,ailen[i]);
644421e10b8SBarry Smith     if (fshift) {
645421e10b8SBarry Smith       ip = aj + ai[i];
646421e10b8SBarry Smith       ap = aa + ai[i];
647421e10b8SBarry Smith       N  = ailen[i];
648421e10b8SBarry Smith       for (j=0; j<N; j++) {
649421e10b8SBarry Smith         ip[j-fshift] = ip[j];
650421e10b8SBarry Smith         ap[j-fshift] = ap[j];
651421e10b8SBarry Smith       }
652421e10b8SBarry Smith     }
653421e10b8SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
654421e10b8SBarry Smith   }
655421e10b8SBarry Smith   if (m) {
656421e10b8SBarry Smith     fshift += imax[m-1] - ailen[m-1];
657421e10b8SBarry Smith     ai[m]   = ai[m-1] + ailen[m-1];
658421e10b8SBarry Smith   }
659421e10b8SBarry Smith   /* reset ilen and imax for each row */
660421e10b8SBarry Smith   for (i=0; i<m; i++) {
661421e10b8SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
662421e10b8SBarry Smith   }
663421e10b8SBarry Smith   a->nz = ai[m];
664ccb205f8SBarry Smith   for (i=0; i<a->nz; i++) {
6656bdcaf15SBarry 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);
6669566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY));
6679566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY));
668ccb205f8SBarry Smith   }
6699566063dSJacob 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));
6709566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A,"Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n",a->reallocs));
6719566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A,"Maximum nonzeros in any row is %" PetscInt_FMT "\n",rmax));
6722205254eSKarl Rupp 
6738e58a170SBarry Smith   A->info.mallocs    += a->reallocs;
674421e10b8SBarry Smith   a->reallocs         = 0;
675421e10b8SBarry Smith   A->info.nz_unneeded = (double)fshift;
676421e10b8SBarry Smith   a->rmax             = rmax;
6779566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_BlockMat(A));
678421e10b8SBarry Smith   PetscFunctionReturn(0);
679421e10b8SBarry Smith }
680421e10b8SBarry Smith 
68181f0254dSBarry Smith static PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg)
682290bbb0aSBarry Smith {
683290bbb0aSBarry Smith   PetscFunctionBegin;
6844e0d8c25SBarry Smith   if (opt == MAT_SYMMETRIC && flg) {
68541f059aeSBarry Smith     A->ops->sor  = MatSOR_BlockMat_Symmetric;
686290bbb0aSBarry Smith     A->ops->mult = MatMult_BlockMat_Symmetric;
687290bbb0aSBarry Smith   } else {
6889566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Unused matrix option %s\n",MatOptions[opt]));
689290bbb0aSBarry Smith   }
690290bbb0aSBarry Smith   PetscFunctionReturn(0);
691290bbb0aSBarry Smith }
692290bbb0aSBarry Smith 
6933964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
694f4259b30SLisandro Dalcin                                        NULL,
695f4259b30SLisandro Dalcin                                        NULL,
696421e10b8SBarry Smith                                        MatMult_BlockMat,
697421e10b8SBarry Smith                                /*  4*/ MatMultAdd_BlockMat,
698421e10b8SBarry Smith                                        MatMultTranspose_BlockMat,
699421e10b8SBarry Smith                                        MatMultTransposeAdd_BlockMat,
700f4259b30SLisandro Dalcin                                        NULL,
701f4259b30SLisandro Dalcin                                        NULL,
702f4259b30SLisandro Dalcin                                        NULL,
703f4259b30SLisandro Dalcin                                /* 10*/ NULL,
704f4259b30SLisandro Dalcin                                        NULL,
705f4259b30SLisandro Dalcin                                        NULL,
70641f059aeSBarry Smith                                        MatSOR_BlockMat,
707f4259b30SLisandro Dalcin                                        NULL,
708f4259b30SLisandro Dalcin                                /* 15*/ NULL,
709f4259b30SLisandro Dalcin                                        NULL,
710f4259b30SLisandro Dalcin                                        NULL,
711f4259b30SLisandro Dalcin                                        NULL,
712f4259b30SLisandro Dalcin                                        NULL,
713f4259b30SLisandro Dalcin                                /* 20*/ NULL,
714421e10b8SBarry Smith                                        MatAssemblyEnd_BlockMat,
715290bbb0aSBarry Smith                                        MatSetOption_BlockMat,
716f4259b30SLisandro Dalcin                                        NULL,
717f4259b30SLisandro Dalcin                                /* 24*/ NULL,
718f4259b30SLisandro Dalcin                                        NULL,
719f4259b30SLisandro Dalcin                                        NULL,
720f4259b30SLisandro Dalcin                                        NULL,
721f4259b30SLisandro Dalcin                                        NULL,
722f4259b30SLisandro Dalcin                                /* 29*/ NULL,
723f4259b30SLisandro Dalcin                                        NULL,
724f4259b30SLisandro Dalcin                                        NULL,
725f4259b30SLisandro Dalcin                                        NULL,
726f4259b30SLisandro Dalcin                                        NULL,
727f4259b30SLisandro Dalcin                                /* 34*/ NULL,
728f4259b30SLisandro Dalcin                                        NULL,
729f4259b30SLisandro Dalcin                                        NULL,
730f4259b30SLisandro Dalcin                                        NULL,
731f4259b30SLisandro Dalcin                                        NULL,
732f4259b30SLisandro Dalcin                                /* 39*/ NULL,
733f4259b30SLisandro Dalcin                                        NULL,
734f4259b30SLisandro Dalcin                                        NULL,
735f4259b30SLisandro Dalcin                                        NULL,
736f4259b30SLisandro Dalcin                                        NULL,
737f4259b30SLisandro Dalcin                                /* 44*/ NULL,
738f4259b30SLisandro Dalcin                                        NULL,
7397d68702bSBarry Smith                                        MatShift_Basic,
740f4259b30SLisandro Dalcin                                        NULL,
741f4259b30SLisandro Dalcin                                        NULL,
742f4259b30SLisandro Dalcin                                /* 49*/ NULL,
743f4259b30SLisandro Dalcin                                        NULL,
744f4259b30SLisandro Dalcin                                        NULL,
745f4259b30SLisandro Dalcin                                        NULL,
746f4259b30SLisandro Dalcin                                        NULL,
747f4259b30SLisandro Dalcin                                /* 54*/ NULL,
748f4259b30SLisandro Dalcin                                        NULL,
749f4259b30SLisandro Dalcin                                        NULL,
750f4259b30SLisandro Dalcin                                        NULL,
751f4259b30SLisandro Dalcin                                        NULL,
7527dae84e0SHong Zhang                                /* 59*/ MatCreateSubMatrix_BlockMat,
753421e10b8SBarry Smith                                        MatDestroy_BlockMat,
754ccb205f8SBarry Smith                                        MatView_BlockMat,
755f4259b30SLisandro Dalcin                                        NULL,
756f4259b30SLisandro Dalcin                                        NULL,
757f4259b30SLisandro Dalcin                                /* 64*/ NULL,
758f4259b30SLisandro Dalcin                                        NULL,
759f4259b30SLisandro Dalcin                                        NULL,
760f4259b30SLisandro Dalcin                                        NULL,
761f4259b30SLisandro Dalcin                                        NULL,
762f4259b30SLisandro Dalcin                                /* 69*/ NULL,
763f4259b30SLisandro Dalcin                                        NULL,
764f4259b30SLisandro Dalcin                                        NULL,
765f4259b30SLisandro Dalcin                                        NULL,
766f4259b30SLisandro Dalcin                                        NULL,
767f4259b30SLisandro Dalcin                                /* 74*/ NULL,
768f4259b30SLisandro Dalcin                                        NULL,
769f4259b30SLisandro Dalcin                                        NULL,
770f4259b30SLisandro Dalcin                                        NULL,
771f4259b30SLisandro Dalcin                                        NULL,
772f4259b30SLisandro Dalcin                                /* 79*/ NULL,
773f4259b30SLisandro Dalcin                                        NULL,
774f4259b30SLisandro Dalcin                                        NULL,
775f4259b30SLisandro Dalcin                                        NULL,
7765bba2384SShri Abhyankar                                        MatLoad_BlockMat,
777f4259b30SLisandro Dalcin                                /* 84*/ NULL,
778f4259b30SLisandro Dalcin                                        NULL,
779f4259b30SLisandro Dalcin                                        NULL,
780f4259b30SLisandro Dalcin                                        NULL,
781f4259b30SLisandro Dalcin                                        NULL,
782f4259b30SLisandro Dalcin                                /* 89*/ NULL,
783f4259b30SLisandro Dalcin                                        NULL,
784f4259b30SLisandro Dalcin                                        NULL,
785f4259b30SLisandro Dalcin                                        NULL,
786f4259b30SLisandro Dalcin                                        NULL,
787f4259b30SLisandro Dalcin                                /* 94*/ NULL,
788f4259b30SLisandro Dalcin                                        NULL,
789f4259b30SLisandro Dalcin                                        NULL,
790f4259b30SLisandro Dalcin                                        NULL,
791f4259b30SLisandro Dalcin                                        NULL,
792f4259b30SLisandro Dalcin                                /* 99*/ NULL,
793f4259b30SLisandro Dalcin                                        NULL,
794f4259b30SLisandro Dalcin                                        NULL,
795f4259b30SLisandro Dalcin                                        NULL,
796f4259b30SLisandro Dalcin                                        NULL,
797f4259b30SLisandro Dalcin                                /*104*/ NULL,
798f4259b30SLisandro Dalcin                                        NULL,
799f4259b30SLisandro Dalcin                                        NULL,
800f4259b30SLisandro Dalcin                                        NULL,
801f4259b30SLisandro Dalcin                                        NULL,
802f4259b30SLisandro Dalcin                                /*109*/ NULL,
803f4259b30SLisandro Dalcin                                        NULL,
804f4259b30SLisandro Dalcin                                        NULL,
805f4259b30SLisandro Dalcin                                        NULL,
806f4259b30SLisandro Dalcin                                        NULL,
807f4259b30SLisandro Dalcin                                /*114*/ NULL,
808f4259b30SLisandro Dalcin                                        NULL,
809f4259b30SLisandro Dalcin                                        NULL,
810f4259b30SLisandro Dalcin                                        NULL,
811f4259b30SLisandro Dalcin                                        NULL,
812f4259b30SLisandro Dalcin                                /*119*/ NULL,
813f4259b30SLisandro Dalcin                                        NULL,
814f4259b30SLisandro Dalcin                                        NULL,
815f4259b30SLisandro Dalcin                                        NULL,
816f4259b30SLisandro Dalcin                                        NULL,
817f4259b30SLisandro Dalcin                                /*124*/ NULL,
818f4259b30SLisandro Dalcin                                        NULL,
819f4259b30SLisandro Dalcin                                        NULL,
820f4259b30SLisandro Dalcin                                        NULL,
821f4259b30SLisandro Dalcin                                        NULL,
822f4259b30SLisandro Dalcin                                /*129*/ NULL,
823f4259b30SLisandro Dalcin                                        NULL,
824f4259b30SLisandro Dalcin                                        NULL,
825f4259b30SLisandro Dalcin                                        NULL,
826f4259b30SLisandro Dalcin                                        NULL,
827f4259b30SLisandro Dalcin                                /*134*/ NULL,
828f4259b30SLisandro Dalcin                                        NULL,
829f4259b30SLisandro Dalcin                                        NULL,
830f4259b30SLisandro Dalcin                                        NULL,
831f4259b30SLisandro Dalcin                                        NULL,
832f4259b30SLisandro Dalcin                                /*139*/ NULL,
833f4259b30SLisandro Dalcin                                        NULL,
834d70f29a3SPierre Jolivet                                        NULL,
835d70f29a3SPierre Jolivet                                        NULL,
836d70f29a3SPierre Jolivet                                        NULL,
837d70f29a3SPierre Jolivet                                /*144*/ NULL,
838d70f29a3SPierre Jolivet                                        NULL,
839d70f29a3SPierre Jolivet                                        NULL,
84099a7f59eSMark Adams                                        NULL,
84199a7f59eSMark Adams                                        NULL,
842f4259b30SLisandro Dalcin                                        NULL
843a5973382SShri Abhyankar };
844421e10b8SBarry Smith 
8458cdf2d9bSBarry Smith /*@C
8468cdf2d9bSBarry Smith    MatBlockMatSetPreallocation - For good matrix assembly performance
8478cdf2d9bSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
8488cdf2d9bSBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
8498cdf2d9bSBarry Smith    during matrix assembly can be increased by more than a factor of 50.
8508cdf2d9bSBarry Smith 
851d083f849SBarry Smith    Collective
8528cdf2d9bSBarry Smith 
8538cdf2d9bSBarry Smith    Input Parameters:
8548cdf2d9bSBarry Smith +  B - The matrix
8558cdf2d9bSBarry Smith .  bs - size of each block in matrix
8568cdf2d9bSBarry Smith .  nz - number of nonzeros per block row (same for all rows)
8578cdf2d9bSBarry Smith -  nnz - array containing the number of nonzeros in the various block rows
8580298fd71SBarry Smith          (possibly different for each row) or NULL
8598cdf2d9bSBarry Smith 
8608cdf2d9bSBarry Smith    Notes:
8618cdf2d9bSBarry Smith      If nnz is given then nz is ignored
8628cdf2d9bSBarry Smith 
8638cdf2d9bSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
8640298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
8658cdf2d9bSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
8668cdf2d9bSBarry Smith    will get TERRIBLE performance, see the users' manual chapter on matrices.
8678cdf2d9bSBarry Smith 
8688cdf2d9bSBarry Smith    Level: intermediate
8698cdf2d9bSBarry Smith 
870db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateBlockMat()`, `MatSetValues()`
8718cdf2d9bSBarry Smith 
8728cdf2d9bSBarry Smith @*/
8737087cfbeSBarry Smith PetscErrorCode  MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
8748cdf2d9bSBarry Smith {
8758cdf2d9bSBarry Smith   PetscFunctionBegin;
876cac4c232SBarry Smith   PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
8778cdf2d9bSBarry Smith   PetscFunctionReturn(0);
8788cdf2d9bSBarry Smith }
8798cdf2d9bSBarry Smith 
88081f0254dSBarry Smith static PetscErrorCode  MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz)
8818cdf2d9bSBarry Smith {
8828cdf2d9bSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
8838cdf2d9bSBarry Smith   PetscInt       i;
8848cdf2d9bSBarry Smith 
8858cdf2d9bSBarry Smith   PetscFunctionBegin;
8869566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(A->rmap,bs));
8879566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(A->cmap,bs));
8889566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
8899566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
8909566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetBlockSize(A->rmap,&bs));
89134ef9618SShri Abhyankar 
8928cdf2d9bSBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
89308401ef6SPierre Jolivet   PetscCheck(nz >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %" PetscInt_FMT,nz);
8948cdf2d9bSBarry Smith   if (nnz) {
895d0f46423SBarry Smith     for (i=0; i<A->rmap->n/bs; i++) {
89608401ef6SPierre 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]);
89708401ef6SPierre 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);
8988cdf2d9bSBarry Smith     }
8998cdf2d9bSBarry Smith   }
900d0f46423SBarry Smith   bmat->mbs = A->rmap->n/bs;
9018cdf2d9bSBarry Smith 
9029566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->right));
9039566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->middle));
9049566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left));
9058cdf2d9bSBarry Smith 
9062ee49352SLisandro Dalcin   if (!bmat->imax) {
9079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(A->rmap->n,&bmat->imax,A->rmap->n,&bmat->ilen));
9089566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)A,2*A->rmap->n*sizeof(PetscInt)));
9092ee49352SLisandro Dalcin   }
910c0aa6a63SJacob Faibussowitsch   if (PetscLikely(nnz)) {
9118cdf2d9bSBarry Smith     nz = 0;
912d0f46423SBarry Smith     for (i=0; i<A->rmap->n/A->rmap->bs; i++) {
9138cdf2d9bSBarry Smith       bmat->imax[i] = nnz[i];
9148cdf2d9bSBarry Smith       nz           += nnz[i];
9158cdf2d9bSBarry Smith     }
916f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation");
9178cdf2d9bSBarry Smith 
9188cdf2d9bSBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
9199566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(bmat->ilen,bmat->mbs));
9208cdf2d9bSBarry Smith 
9218cdf2d9bSBarry Smith   /* allocate the matrix space */
9229566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i));
9239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nz,&bmat->a,nz,&bmat->j,A->rmap->n+1,&bmat->i));
9249566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt))));
9258cdf2d9bSBarry Smith   bmat->i[0] = 0;
9268cdf2d9bSBarry Smith   for (i=1; i<bmat->mbs+1; i++) {
9278cdf2d9bSBarry Smith     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
9288cdf2d9bSBarry Smith   }
9298cdf2d9bSBarry Smith   bmat->singlemalloc = PETSC_TRUE;
9308cdf2d9bSBarry Smith   bmat->free_a       = PETSC_TRUE;
9318cdf2d9bSBarry Smith   bmat->free_ij      = PETSC_TRUE;
9328cdf2d9bSBarry Smith 
9338cdf2d9bSBarry Smith   bmat->nz            = 0;
9348cdf2d9bSBarry Smith   bmat->maxnz         = nz;
9358cdf2d9bSBarry Smith   A->info.nz_unneeded = (double)bmat->maxnz;
9369566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
9378cdf2d9bSBarry Smith   PetscFunctionReturn(0);
9388cdf2d9bSBarry Smith }
9398cdf2d9bSBarry Smith 
940421e10b8SBarry Smith /*MC
941421e10b8SBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
942421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
943421e10b8SBarry Smith 
944421e10b8SBarry Smith   Level: advanced
945421e10b8SBarry Smith 
946db781477SPatrick Sanan .seealso: `MatCreateBlockMat()`
947421e10b8SBarry Smith 
948421e10b8SBarry Smith M*/
949421e10b8SBarry Smith 
9508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A)
951421e10b8SBarry Smith {
952421e10b8SBarry Smith   Mat_BlockMat   *b;
953421e10b8SBarry Smith 
954421e10b8SBarry Smith   PetscFunctionBegin;
9559566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(A,&b));
956421e10b8SBarry Smith   A->data = (void*)b;
9579566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps)));
958421e10b8SBarry Smith 
959421e10b8SBarry Smith   A->assembled    = PETSC_TRUE;
960421e10b8SBarry Smith   A->preallocated = PETSC_FALSE;
9619566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT));
962290bbb0aSBarry Smith 
9639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatBlockMatSetPreallocation_C",MatBlockMatSetPreallocation_BlockMat));
964421e10b8SBarry Smith   PetscFunctionReturn(0);
965421e10b8SBarry Smith }
966421e10b8SBarry Smith 
967421e10b8SBarry Smith /*@C
968bbd3f336SJed Brown    MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential Mat object
969421e10b8SBarry Smith 
970d083f849SBarry Smith   Collective
971421e10b8SBarry Smith 
972421e10b8SBarry Smith    Input Parameters:
973421e10b8SBarry Smith +  comm - MPI communicator
974421e10b8SBarry Smith .  m - number of rows
975421e10b8SBarry Smith .  n  - number of columns
9768cdf2d9bSBarry Smith .  bs - size of each submatrix
9778cdf2d9bSBarry Smith .  nz  - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known)
9780298fd71SBarry Smith -  nnz - expected number of nonzers per block row if known (use NULL otherwise)
979421e10b8SBarry Smith 
980421e10b8SBarry Smith    Output Parameter:
981421e10b8SBarry Smith .  A - the matrix
982421e10b8SBarry Smith 
983421e10b8SBarry Smith    Level: intermediate
984421e10b8SBarry Smith 
98595452b02SPatrick Sanan    Notes:
98695452b02SPatrick Sanan     Matrices of this type are nominally-sparse matrices in which each "entry" is a Mat object.  Each Mat must
987bbd3f336SJed Brown    have the same size and be sequential.  The local and global sizes must be compatible with this decomposition.
988bbd3f336SJed Brown 
989bbd3f336SJed Brown    For matrices containing parallel submatrices and variable block sizes, see MATNEST.
990421e10b8SBarry Smith 
991db781477SPatrick Sanan .seealso: `MATBLOCKMAT`, `MatCreateNest()`
992421e10b8SBarry Smith @*/
9937087cfbeSBarry Smith PetscErrorCode  MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A)
994421e10b8SBarry Smith {
995421e10b8SBarry Smith   PetscFunctionBegin;
9969566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,A));
9979566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE));
9989566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A,MATBLOCKMAT));
9999566063dSJacob Faibussowitsch   PetscCall(MatBlockMatSetPreallocation(*A,bs,nz,nnz));
1000421e10b8SBarry Smith   PetscFunctionReturn(0);
1001421e10b8SBarry Smith }
1002