xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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;
23290bbb0aSBarry Smith   PetscErrorCode    ierr;
24d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
25290bbb0aSBarry Smith   const PetscInt    *idx;
26290bbb0aSBarry Smith   IS                row,col;
27290bbb0aSBarry Smith   MatFactorInfo     info;
286e63c7a1SBarry Smith   Vec               left = a->left,right = a->right, middle = a->middle;
29290bbb0aSBarry Smith   Mat               *diag;
30290bbb0aSBarry Smith 
31290bbb0aSBarry Smith   PetscFunctionBegin;
32290bbb0aSBarry Smith   its = its*lits;
33*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(its <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits);
34*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(flag & SOR_EISENSTAT,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
35*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(omega != 1.0,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
36*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(fshift,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
3726fbe8dcSKarl Rupp   if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) {
38e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
3926fbe8dcSKarl Rupp   }
40290bbb0aSBarry Smith 
41290bbb0aSBarry Smith   if (!a->diags) {
42785e854fSJed Brown     ierr = PetscMalloc1(mbs,&a->diags);CHKERRQ(ierr);
43290bbb0aSBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
44290bbb0aSBarry Smith     for (i=0; i<mbs; i++) {
452692d6eeSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr);
46719d5645SBarry Smith       ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr);
47719d5645SBarry Smith       ierr = MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr);
486bf464f9SBarry Smith       ierr = ISDestroy(&row);CHKERRQ(ierr);
496bf464f9SBarry Smith       ierr = ISDestroy(&col);CHKERRQ(ierr);
50290bbb0aSBarry Smith     }
516e63c7a1SBarry Smith     ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr);
52290bbb0aSBarry Smith   }
53290bbb0aSBarry Smith   diag = a->diags;
54290bbb0aSBarry Smith 
55290bbb0aSBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
56290bbb0aSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
57290bbb0aSBarry Smith   /* copy right hand side because it must be modified during iteration */
586e63c7a1SBarry Smith   ierr = VecCopy(bb,a->workb);CHKERRQ(ierr);
593649974fSBarry Smith   ierr = VecGetArrayRead(a->workb,&b);CHKERRQ(ierr);
60290bbb0aSBarry Smith 
6141f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
62290bbb0aSBarry Smith   while (its--) {
63290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
64290bbb0aSBarry Smith 
65290bbb0aSBarry Smith       for (i=0; i<mbs; i++) {
666e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
676e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
686e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
69290bbb0aSBarry Smith 
70290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
71290bbb0aSBarry Smith         for (j=0; j<n; j++) {
72290bbb0aSBarry Smith           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
73290bbb0aSBarry Smith           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
74290bbb0aSBarry Smith           ierr = VecResetArray(right);CHKERRQ(ierr);
75290bbb0aSBarry Smith         }
76290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
77290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
78290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
79290bbb0aSBarry Smith 
80290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
81290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
826e63c7a1SBarry Smith 
8341f059aeSBarry Smith         /* now adjust right hand side, see MatSOR_SeqSBAIJ */
846e63c7a1SBarry Smith         for (j=0; j<n; j++) {
856e63c7a1SBarry Smith           ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr);
866e63c7a1SBarry Smith           ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr);
876e63c7a1SBarry Smith           ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr);
886e63c7a1SBarry Smith           ierr = VecResetArray(middle);CHKERRQ(ierr);
896e63c7a1SBarry Smith         }
90290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
916e63c7a1SBarry Smith 
92290bbb0aSBarry Smith       }
93290bbb0aSBarry Smith     }
94290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
95290bbb0aSBarry Smith 
96290bbb0aSBarry Smith       for (i=mbs-1; i>=0; i--) {
976e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
986e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
996e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
100290bbb0aSBarry Smith 
101290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
102290bbb0aSBarry Smith         for (j=0; j<n; j++) {
103290bbb0aSBarry Smith           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
104290bbb0aSBarry Smith           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
105290bbb0aSBarry Smith           ierr = VecResetArray(right);CHKERRQ(ierr);
106290bbb0aSBarry Smith         }
107290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
108290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
109290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
110290bbb0aSBarry Smith 
111290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
112290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
113290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
114290bbb0aSBarry Smith 
115290bbb0aSBarry Smith       }
116290bbb0aSBarry Smith     }
117290bbb0aSBarry Smith   }
118290bbb0aSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1193649974fSBarry Smith   ierr = VecRestoreArrayRead(a->workb,&b);CHKERRQ(ierr);
120290bbb0aSBarry Smith   PetscFunctionReturn(0);
121290bbb0aSBarry Smith }
122290bbb0aSBarry Smith 
12381f0254dSBarry Smith static PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
124ccb205f8SBarry Smith {
125ccb205f8SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
126ccb205f8SBarry Smith   PetscScalar       *x;
127a2ea699eSBarry Smith   const Mat         *v;
128ccb205f8SBarry Smith   const PetscScalar *b;
129ccb205f8SBarry Smith   PetscErrorCode    ierr;
130d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
131ccb205f8SBarry Smith   const PetscInt    *idx;
132ccb205f8SBarry Smith   IS                row,col;
133ccb205f8SBarry Smith   MatFactorInfo     info;
134ccb205f8SBarry Smith   Vec               left = a->left,right = a->right;
135ccb205f8SBarry Smith   Mat               *diag;
136ccb205f8SBarry Smith 
137ccb205f8SBarry Smith   PetscFunctionBegin;
138ccb205f8SBarry Smith   its = its*lits;
139*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(its <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits);
140*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(flag & SOR_EISENSTAT,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
141*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(omega != 1.0,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
142*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(fshift,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
143ccb205f8SBarry Smith 
144ccb205f8SBarry Smith   if (!a->diags) {
145785e854fSJed Brown     ierr = PetscMalloc1(mbs,&a->diags);CHKERRQ(ierr);
146ccb205f8SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
147ccb205f8SBarry Smith     for (i=0; i<mbs; i++) {
1482692d6eeSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr);
149719d5645SBarry Smith       ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr);
150719d5645SBarry Smith       ierr = MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr);
1516bf464f9SBarry Smith       ierr = ISDestroy(&row);CHKERRQ(ierr);
1526bf464f9SBarry Smith       ierr = ISDestroy(&col);CHKERRQ(ierr);
153ccb205f8SBarry Smith     }
154ccb205f8SBarry Smith   }
155ccb205f8SBarry Smith   diag = a->diags;
156ccb205f8SBarry Smith 
157ccb205f8SBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
158ccb205f8SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1593649974fSBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
160ccb205f8SBarry Smith 
16141f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
162ccb205f8SBarry Smith   while (its--) {
163ccb205f8SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
164ccb205f8SBarry Smith 
165ccb205f8SBarry Smith       for (i=0; i<mbs; i++) {
166ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
167ccb205f8SBarry Smith         idx = a->j + a->i[i];
168ccb205f8SBarry Smith         v   = a->a + a->i[i];
169ccb205f8SBarry Smith 
170ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
171ccb205f8SBarry Smith         for (j=0; j<n; j++) {
172ccb205f8SBarry Smith           if (idx[j] != i) {
173ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
174ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
175ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
176ccb205f8SBarry Smith           }
177ccb205f8SBarry Smith         }
178ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
179ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
180ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
181ccb205f8SBarry Smith 
182ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
183ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
184ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
185ccb205f8SBarry Smith       }
186ccb205f8SBarry Smith     }
187ccb205f8SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
188ccb205f8SBarry Smith 
189ccb205f8SBarry Smith       for (i=mbs-1; i>=0; i--) {
190ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
191ccb205f8SBarry Smith         idx = a->j + a->i[i];
192ccb205f8SBarry Smith         v   = a->a + a->i[i];
193ccb205f8SBarry Smith 
194ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
195ccb205f8SBarry Smith         for (j=0; j<n; j++) {
196ccb205f8SBarry Smith           if (idx[j] != i) {
197ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
198ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
199ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
200ccb205f8SBarry Smith           }
201ccb205f8SBarry Smith         }
202ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
203ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
204ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
205ccb205f8SBarry Smith 
206ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
207ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
208ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
209ccb205f8SBarry Smith 
210ccb205f8SBarry Smith       }
211ccb205f8SBarry Smith     }
212ccb205f8SBarry Smith   }
213ccb205f8SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2143649974fSBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
215ccb205f8SBarry Smith   PetscFunctionReturn(0);
216ccb205f8SBarry Smith }
217ccb205f8SBarry Smith 
21881f0254dSBarry Smith static PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
219421e10b8SBarry Smith {
220421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
221421e10b8SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
222421e10b8SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
223d0f46423SBarry Smith   PetscInt       *aj  =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
224421e10b8SBarry Smith   PetscErrorCode ierr;
225421e10b8SBarry Smith   PetscInt       ridx,cidx;
226ace3abfcSBarry Smith   PetscBool      roworiented=a->roworiented;
227421e10b8SBarry Smith   MatScalar      value;
228421e10b8SBarry Smith   Mat            *ap,*aa = a->a;
229421e10b8SBarry Smith 
230421e10b8SBarry Smith   PetscFunctionBegin;
231421e10b8SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
232421e10b8SBarry Smith     row  = im[k];
233421e10b8SBarry Smith     brow = row/bs;
234421e10b8SBarry Smith     if (row < 0) continue;
235*2c71b3e2SJacob Faibussowitsch     PetscAssertFalse(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);
236421e10b8SBarry Smith     rp   = aj + ai[brow];
237421e10b8SBarry Smith     ap   = aa + ai[brow];
238421e10b8SBarry Smith     rmax = imax[brow];
239421e10b8SBarry Smith     nrow = ailen[brow];
240421e10b8SBarry Smith     low  = 0;
241421e10b8SBarry Smith     high = nrow;
242421e10b8SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
243421e10b8SBarry Smith       if (in[l] < 0) continue;
244*2c71b3e2SJacob Faibussowitsch       PetscAssertFalse(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);
245421e10b8SBarry Smith       col = in[l]; bcol = col/bs;
2466e63c7a1SBarry Smith       if (A->symmetric && brow > bcol) continue;
247421e10b8SBarry Smith       ridx = row % bs; cidx = col % bs;
2482205254eSKarl Rupp       if (roworiented) value = v[l + k*n];
2492205254eSKarl Rupp       else value = v[k + l*m];
2502205254eSKarl Rupp 
2512205254eSKarl Rupp       if (col <= lastcol) low = 0;
2522205254eSKarl Rupp       else high = nrow;
253421e10b8SBarry Smith       lastcol = col;
254421e10b8SBarry Smith       while (high-low > 7) {
255421e10b8SBarry Smith         t = (low+high)/2;
256421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
257421e10b8SBarry Smith         else              low  = t;
258421e10b8SBarry Smith       }
259421e10b8SBarry Smith       for (i=low; i<high; i++) {
260421e10b8SBarry Smith         if (rp[i] > bcol) break;
2612205254eSKarl Rupp         if (rp[i] == bcol) goto noinsert1;
262421e10b8SBarry Smith       }
263421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
264*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(nonew == -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
265fef13f97SBarry Smith       MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
266421e10b8SBarry Smith       N = nrow++ - 1; high++;
267421e10b8SBarry Smith       /* shift up all the later entries in this row */
268421e10b8SBarry Smith       for (ii=N; ii>=i; ii--) {
269421e10b8SBarry Smith         rp[ii+1] = rp[ii];
270421e10b8SBarry Smith         ap[ii+1] = ap[ii];
271421e10b8SBarry Smith       }
272f4259b30SLisandro Dalcin       if (N>=i) ap[i] = NULL;
273421e10b8SBarry Smith       rp[i] = bcol;
274421e10b8SBarry Smith       a->nz++;
275e56f5c9eSBarry Smith       A->nonzerostate++;
276421e10b8SBarry Smith noinsert1:;
277421e10b8SBarry Smith       if (!*(ap+i)) {
278f4259b30SLisandro Dalcin         ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,NULL,ap+i);CHKERRQ(ierr);
279b0223207SBarry Smith       }
280421e10b8SBarry Smith       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
281421e10b8SBarry Smith       low  = i;
282421e10b8SBarry Smith     }
283421e10b8SBarry Smith     ailen[brow] = nrow;
284421e10b8SBarry Smith   }
285421e10b8SBarry Smith   PetscFunctionReturn(0);
286421e10b8SBarry Smith }
287421e10b8SBarry Smith 
28881f0254dSBarry Smith static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
289a5973382SShri Abhyankar {
290a5973382SShri Abhyankar   PetscErrorCode    ierr;
291a5973382SShri Abhyankar   Mat               tmpA;
292a5973382SShri Abhyankar   PetscInt          i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0;
293a5973382SShri Abhyankar   const PetscInt    *cols;
294a5973382SShri Abhyankar   const PetscScalar *values;
295ace3abfcSBarry Smith   PetscBool         flg = PETSC_FALSE,notdone;
296a5973382SShri Abhyankar   Mat_SeqAIJ        *a;
297a5973382SShri Abhyankar   Mat_BlockMat      *amat;
298a5973382SShri Abhyankar 
299a5973382SShri Abhyankar   PetscFunctionBegin;
300c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
301c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
302a5973382SShri Abhyankar   ierr = MatCreate(PETSC_COMM_SELF,&tmpA);CHKERRQ(ierr);
303a5973382SShri Abhyankar   ierr = MatSetType(tmpA,MATSEQAIJ);CHKERRQ(ierr);
304112444f4SShri Abhyankar   ierr = MatLoad_SeqAIJ(tmpA,viewer);CHKERRQ(ierr);
305a5973382SShri Abhyankar 
306a5973382SShri Abhyankar   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
3070298fd71SBarry Smith   ierr = PetscOptionsBegin(PETSC_COMM_SELF,NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr);
3080298fd71SBarry Smith   ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);CHKERRQ(ierr);
3090298fd71SBarry Smith   ierr = PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,NULL);CHKERRQ(ierr);
310a5973382SShri Abhyankar   ierr = PetscOptionsEnd();CHKERRQ(ierr);
311a5973382SShri Abhyankar 
312a5973382SShri Abhyankar   /* Determine number of nonzero blocks for each block row */
313a5973382SShri Abhyankar   a    = (Mat_SeqAIJ*) tmpA->data;
314a5973382SShri Abhyankar   mbs  = m/bs;
315dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&lens,bs,&ii,bs,&ilens);CHKERRQ(ierr);
316580bdb30SBarry Smith   ierr = PetscArrayzero(lens,mbs);CHKERRQ(ierr);
317a5973382SShri Abhyankar 
318a5973382SShri Abhyankar   for (i=0; i<mbs; i++) {
319a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
320a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
321a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
322a5973382SShri Abhyankar     }
323a5973382SShri Abhyankar 
324a5973382SShri Abhyankar     currentcol = -1;
325a5973382SShri Abhyankar     while (PETSC_TRUE) {
326a5973382SShri Abhyankar       notdone = PETSC_FALSE;
327a5973382SShri Abhyankar       nextcol = 1000000000;
328a5973382SShri Abhyankar       for (j=0; j<bs; j++) {
329a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) {
330a5973382SShri Abhyankar           ii[j]++;
331a5973382SShri Abhyankar           ilens[j]--;
332a5973382SShri Abhyankar         }
333a5973382SShri Abhyankar         if (ilens[j] > 0) {
334a5973382SShri Abhyankar           notdone = PETSC_TRUE;
335a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
336a5973382SShri Abhyankar         }
337a5973382SShri Abhyankar       }
338a5973382SShri Abhyankar       if (!notdone) break;
339a5973382SShri Abhyankar       if (!flg || (nextcol >= i)) lens[i]++;
340a5973382SShri Abhyankar       currentcol = nextcol;
341a5973382SShri Abhyankar     }
342a5973382SShri Abhyankar   }
343a5973382SShri Abhyankar 
344a5973382SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
345a5973382SShri Abhyankar     ierr = MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
346a5973382SShri Abhyankar   }
347a5973382SShri Abhyankar   ierr = MatBlockMatSetPreallocation(newmat,bs,0,lens);CHKERRQ(ierr);
348a5973382SShri Abhyankar   if (flg) {
349a5973382SShri Abhyankar     ierr = MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
350a5973382SShri Abhyankar   }
351a5973382SShri Abhyankar   amat = (Mat_BlockMat*)(newmat)->data;
352a5973382SShri Abhyankar 
353a5973382SShri Abhyankar   /* preallocate the submatrices */
354785e854fSJed Brown   ierr = PetscMalloc1(bs,&llens);CHKERRQ(ierr);
355a5973382SShri Abhyankar   for (i=0; i<mbs; i++) { /* loops for block rows */
356a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
357a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
358a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
359a5973382SShri Abhyankar     }
360a5973382SShri Abhyankar 
361a5973382SShri Abhyankar     currentcol = 1000000000;
362a5973382SShri Abhyankar     for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */
363a5973382SShri Abhyankar       if (ilens[j] > 0) {
364a5973382SShri Abhyankar         currentcol = PetscMin(currentcol,ii[j][0]/bs);
365a5973382SShri Abhyankar       }
366a5973382SShri Abhyankar     }
367a5973382SShri Abhyankar 
368a5973382SShri Abhyankar     while (PETSC_TRUE) {  /* loops over blocks in block row */
369a5973382SShri Abhyankar       notdone = PETSC_FALSE;
370a5973382SShri Abhyankar       nextcol = 1000000000;
371580bdb30SBarry Smith       ierr    = PetscArrayzero(llens,bs);CHKERRQ(ierr);
372a5973382SShri Abhyankar       for (j=0; j<bs; j++) { /* loop over rows in block */
373a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */
374a5973382SShri Abhyankar           ii[j]++;
375a5973382SShri Abhyankar           ilens[j]--;
376a5973382SShri Abhyankar           llens[j]++;
377a5973382SShri Abhyankar         }
378a5973382SShri Abhyankar         if (ilens[j] > 0) {
379a5973382SShri Abhyankar           notdone = PETSC_TRUE;
380a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
381a5973382SShri Abhyankar         }
382a5973382SShri Abhyankar       }
383*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(cnt >= amat->maxnz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %" PetscInt_FMT,cnt);
384a5973382SShri Abhyankar       if (!flg || currentcol >= i) {
385a5973382SShri Abhyankar         amat->j[cnt] = currentcol;
386a5973382SShri Abhyankar         ierr         = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr);
387a5973382SShri Abhyankar       }
388a5973382SShri Abhyankar 
389a5973382SShri Abhyankar       if (!notdone) break;
390a5973382SShri Abhyankar       currentcol = nextcol;
391a5973382SShri Abhyankar     }
392a5973382SShri Abhyankar     amat->ilen[i] = lens[i];
393a5973382SShri Abhyankar   }
394a5973382SShri Abhyankar 
395a5973382SShri Abhyankar   ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr);
396a5973382SShri Abhyankar   ierr = PetscFree(llens);CHKERRQ(ierr);
397a5973382SShri Abhyankar 
398a5973382SShri Abhyankar   /* copy over the matrix, one row at a time */
399a5973382SShri Abhyankar   for (i=0; i<m; i++) {
400a5973382SShri Abhyankar     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
401a5973382SShri Abhyankar     ierr = MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
402a5973382SShri Abhyankar     ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
403a5973382SShri Abhyankar   }
404a5973382SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
405a5973382SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
406a5973382SShri Abhyankar   PetscFunctionReturn(0);
407a5973382SShri Abhyankar }
408a5973382SShri Abhyankar 
40981f0254dSBarry Smith static PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
410ccb205f8SBarry Smith {
41164075487SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
412ccb205f8SBarry Smith   PetscErrorCode    ierr;
413ccb205f8SBarry Smith   const char        *name;
414ccb205f8SBarry Smith   PetscViewerFormat format;
415ccb205f8SBarry Smith 
416ccb205f8SBarry Smith   PetscFunctionBegin;
417ccb205f8SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
418ccb205f8SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
419ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
420c0aa6a63SJacob Faibussowitsch     ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %" PetscInt_FMT " \n",a->nz);CHKERRQ(ierr);
42164075487SBarry Smith     if (A->symmetric) {
4228c1ad04dSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr);
42364075487SBarry Smith     }
424ccb205f8SBarry Smith   }
425ccb205f8SBarry Smith   PetscFunctionReturn(0);
426ccb205f8SBarry Smith }
427421e10b8SBarry Smith 
42881f0254dSBarry Smith static PetscErrorCode MatDestroy_BlockMat(Mat mat)
429421e10b8SBarry Smith {
430421e10b8SBarry Smith   PetscErrorCode ierr;
431421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
432ccb205f8SBarry Smith   PetscInt       i;
433421e10b8SBarry Smith 
434421e10b8SBarry Smith   PetscFunctionBegin;
4356bf464f9SBarry Smith   ierr = VecDestroy(&bmat->right);CHKERRQ(ierr);
4366bf464f9SBarry Smith   ierr = VecDestroy(&bmat->left);CHKERRQ(ierr);
4376bf464f9SBarry Smith   ierr = VecDestroy(&bmat->middle);CHKERRQ(ierr);
4386bf464f9SBarry Smith   ierr = VecDestroy(&bmat->workb);CHKERRQ(ierr);
439ccb205f8SBarry Smith   if (bmat->diags) {
440d0f46423SBarry Smith     for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) {
4416bf464f9SBarry Smith       ierr = MatDestroy(&bmat->diags[i]);CHKERRQ(ierr);
442ccb205f8SBarry Smith     }
443ccb205f8SBarry Smith   }
444ccb205f8SBarry Smith   if (bmat->a) {
445ccb205f8SBarry Smith     for (i=0; i<bmat->nz; i++) {
4466bf464f9SBarry Smith       ierr = MatDestroy(&bmat->a[i]);CHKERRQ(ierr);
447ccb205f8SBarry Smith     }
448ccb205f8SBarry Smith   }
449ccb205f8SBarry Smith   ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
450bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
451421e10b8SBarry Smith   PetscFunctionReturn(0);
452421e10b8SBarry Smith }
453421e10b8SBarry Smith 
45481f0254dSBarry Smith static PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
455421e10b8SBarry Smith {
456421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
457421e10b8SBarry Smith   PetscErrorCode ierr;
458421e10b8SBarry Smith   PetscScalar    *xx,*yy;
459d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
460421e10b8SBarry Smith   Mat            *aa;
461421e10b8SBarry Smith 
462421e10b8SBarry Smith   PetscFunctionBegin;
463421e10b8SBarry Smith   /*
464421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
465421e10b8SBarry Smith   */
466421e10b8SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
467ccb205f8SBarry Smith 
468421e10b8SBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
469421e10b8SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
470421e10b8SBarry Smith   aj   = bmat->j;
471421e10b8SBarry Smith   aa   = bmat->a;
472421e10b8SBarry Smith   ii   = bmat->i;
473421e10b8SBarry Smith   for (i=0; i<m; i++) {
474421e10b8SBarry Smith     jrow = ii[i];
475ccb205f8SBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
476421e10b8SBarry Smith     n    = ii[i+1] - jrow;
477421e10b8SBarry Smith     for (j=0; j<n; j++) {
478421e10b8SBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
479ccb205f8SBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
480421e10b8SBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
481421e10b8SBarry Smith       jrow++;
482421e10b8SBarry Smith     }
483421e10b8SBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
484421e10b8SBarry Smith   }
485421e10b8SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
486421e10b8SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
487421e10b8SBarry Smith   PetscFunctionReturn(0);
488421e10b8SBarry Smith }
489421e10b8SBarry Smith 
490290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
491290bbb0aSBarry Smith {
492290bbb0aSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
493290bbb0aSBarry Smith   PetscErrorCode ierr;
494290bbb0aSBarry Smith   PetscScalar    *xx,*yy;
495d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
496290bbb0aSBarry Smith   Mat            *aa;
497290bbb0aSBarry Smith 
498290bbb0aSBarry Smith   PetscFunctionBegin;
499290bbb0aSBarry Smith   /*
500290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
501290bbb0aSBarry Smith   */
502290bbb0aSBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
503290bbb0aSBarry Smith 
504290bbb0aSBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
505290bbb0aSBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
506290bbb0aSBarry Smith   aj   = bmat->j;
507290bbb0aSBarry Smith   aa   = bmat->a;
508290bbb0aSBarry Smith   ii   = bmat->i;
509290bbb0aSBarry Smith   for (i=0; i<m; i++) {
510290bbb0aSBarry Smith     jrow = ii[i];
511290bbb0aSBarry Smith     n    = ii[i+1] - jrow;
512290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
513290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr);
514290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
515290bbb0aSBarry Smith     if (aj[jrow] == i) {
516290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
517290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
518290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
519290bbb0aSBarry Smith       jrow++;
520290bbb0aSBarry Smith       n--;
521290bbb0aSBarry Smith     }
522290bbb0aSBarry Smith     for (j=0; j<n; j++) {
523290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);            /* upper triangular part */
524290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
525290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
526290bbb0aSBarry Smith 
527290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr);            /* lower triangular part */
528290bbb0aSBarry Smith       ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr);
529290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
530290bbb0aSBarry Smith       jrow++;
531290bbb0aSBarry Smith     }
532290bbb0aSBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
533290bbb0aSBarry Smith     ierr = VecResetArray(bmat->middle);CHKERRQ(ierr);
534290bbb0aSBarry Smith   }
535290bbb0aSBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
536290bbb0aSBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
537290bbb0aSBarry Smith   PetscFunctionReturn(0);
538290bbb0aSBarry Smith }
539290bbb0aSBarry Smith 
54081f0254dSBarry Smith static PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
541421e10b8SBarry Smith {
542421e10b8SBarry Smith   PetscFunctionBegin;
543421e10b8SBarry Smith   PetscFunctionReturn(0);
544421e10b8SBarry Smith }
545421e10b8SBarry Smith 
54681f0254dSBarry Smith static PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
547421e10b8SBarry Smith {
548421e10b8SBarry Smith   PetscFunctionBegin;
549421e10b8SBarry Smith   PetscFunctionReturn(0);
550421e10b8SBarry Smith }
551421e10b8SBarry Smith 
55281f0254dSBarry Smith static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
553421e10b8SBarry Smith {
554421e10b8SBarry Smith   PetscFunctionBegin;
555421e10b8SBarry Smith   PetscFunctionReturn(0);
556421e10b8SBarry Smith }
557421e10b8SBarry Smith 
558ccb205f8SBarry Smith /*
559ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
560ccb205f8SBarry Smith */
56181f0254dSBarry Smith static PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
562ccb205f8SBarry Smith {
563ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
564ccb205f8SBarry Smith   PetscErrorCode ierr;
565d0f46423SBarry Smith   PetscInt       i,j,mbs = A->rmap->n/A->rmap->bs;
566ccb205f8SBarry Smith 
567ccb205f8SBarry Smith   PetscFunctionBegin;
568ccb205f8SBarry Smith   if (!a->diag) {
569785e854fSJed Brown     ierr = PetscMalloc1(mbs,&a->diag);CHKERRQ(ierr);
570ccb205f8SBarry Smith   }
571ccb205f8SBarry Smith   for (i=0; i<mbs; i++) {
572ccb205f8SBarry Smith     a->diag[i] = a->i[i+1];
573ccb205f8SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
574ccb205f8SBarry Smith       if (a->j[j] == i) {
575ccb205f8SBarry Smith         a->diag[i] = j;
576ccb205f8SBarry Smith         break;
577ccb205f8SBarry Smith       }
578ccb205f8SBarry Smith     }
579ccb205f8SBarry Smith   }
580ccb205f8SBarry Smith   PetscFunctionReturn(0);
581ccb205f8SBarry Smith }
582ccb205f8SBarry Smith 
5837dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
584ccb205f8SBarry Smith {
585ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
586ccb205f8SBarry Smith   Mat_SeqAIJ     *c;
587ccb205f8SBarry Smith   PetscErrorCode ierr;
588ccb205f8SBarry Smith   PetscInt       i,k,first,step,lensi,nrows,ncols;
589d2a212eaSBarry Smith   PetscInt       *j_new,*i_new,*aj = a->j,*ailen = a->ilen;
5901620fd73SBarry Smith   PetscScalar    *a_new;
591ccb205f8SBarry Smith   Mat            C,*aa = a->a;
592ace3abfcSBarry Smith   PetscBool      stride,equal;
593ccb205f8SBarry Smith 
594ccb205f8SBarry Smith   PetscFunctionBegin;
595ccb205f8SBarry Smith   ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr);
596*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for identical column and row indices");
597251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);CHKERRQ(ierr);
598*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!stride,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices");
599ccb205f8SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
600*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(step != A->rmap->bs,PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block");
601ccb205f8SBarry Smith 
602ccb205f8SBarry Smith   ierr  = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
603ccb205f8SBarry Smith   ncols = nrows;
604ccb205f8SBarry Smith 
605ccb205f8SBarry Smith   /* create submatrix */
606ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
607ccb205f8SBarry Smith     PetscInt n_cols,n_rows;
608ccb205f8SBarry Smith     C    = *B;
609ccb205f8SBarry Smith     ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr);
610*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(n_rows != nrows || n_cols != ncols,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
611ccb205f8SBarry Smith     ierr = MatZeroEntries(C);CHKERRQ(ierr);
612ccb205f8SBarry Smith   } else {
613ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
614ccb205f8SBarry Smith     ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
6156e63c7a1SBarry Smith     if (A->symmetric) {
6166e63c7a1SBarry Smith       ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
6176e63c7a1SBarry Smith     } else {
618ccb205f8SBarry Smith       ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr);
6196e63c7a1SBarry Smith     }
6206e63c7a1SBarry Smith     ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr);
6216e63c7a1SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr);
622ccb205f8SBarry Smith   }
623ccb205f8SBarry Smith   c = (Mat_SeqAIJ*)C->data;
624ccb205f8SBarry Smith 
625ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
626ccb205f8SBarry Smith   a_new = c->a;
627ccb205f8SBarry Smith   j_new = c->j;
628ccb205f8SBarry Smith   i_new = c->i;
629ccb205f8SBarry Smith 
630ccb205f8SBarry Smith   for (i=0; i<nrows; i++) {
631ccb205f8SBarry Smith     lensi = ailen[i];
632ccb205f8SBarry Smith     for (k=0; k<lensi; k++) {
633ccb205f8SBarry Smith       *j_new++ = *aj++;
6341620fd73SBarry Smith       ierr     = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr);
635ccb205f8SBarry Smith     }
636ccb205f8SBarry Smith     i_new[i+1] = i_new[i] + lensi;
637ccb205f8SBarry Smith     c->ilen[i] = lensi;
638ccb205f8SBarry Smith   }
639ccb205f8SBarry Smith 
640ccb205f8SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
641ccb205f8SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
642ccb205f8SBarry Smith   *B   = C;
643ccb205f8SBarry Smith   PetscFunctionReturn(0);
644ccb205f8SBarry Smith }
645ccb205f8SBarry Smith 
64681f0254dSBarry Smith static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
647421e10b8SBarry Smith {
648421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
649421e10b8SBarry Smith   PetscErrorCode ierr;
650421e10b8SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
651ccb205f8SBarry Smith   PetscInt       m      = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
652421e10b8SBarry Smith   Mat            *aa    = a->a,*ap;
653421e10b8SBarry Smith 
654421e10b8SBarry Smith   PetscFunctionBegin;
655421e10b8SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
656421e10b8SBarry Smith 
657421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
658421e10b8SBarry Smith   for (i=1; i<m; i++) {
659421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
660421e10b8SBarry Smith     fshift += imax[i-1] - ailen[i-1];
661421e10b8SBarry Smith     rmax    = PetscMax(rmax,ailen[i]);
662421e10b8SBarry Smith     if (fshift) {
663421e10b8SBarry Smith       ip = aj + ai[i];
664421e10b8SBarry Smith       ap = aa + ai[i];
665421e10b8SBarry Smith       N  = ailen[i];
666421e10b8SBarry Smith       for (j=0; j<N; j++) {
667421e10b8SBarry Smith         ip[j-fshift] = ip[j];
668421e10b8SBarry Smith         ap[j-fshift] = ap[j];
669421e10b8SBarry Smith       }
670421e10b8SBarry Smith     }
671421e10b8SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
672421e10b8SBarry Smith   }
673421e10b8SBarry Smith   if (m) {
674421e10b8SBarry Smith     fshift += imax[m-1] - ailen[m-1];
675421e10b8SBarry Smith     ai[m]   = ai[m-1] + ailen[m-1];
676421e10b8SBarry Smith   }
677421e10b8SBarry Smith   /* reset ilen and imax for each row */
678421e10b8SBarry Smith   for (i=0; i<m; i++) {
679421e10b8SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
680421e10b8SBarry Smith   }
681421e10b8SBarry Smith   a->nz = ai[m];
682ccb205f8SBarry Smith   for (i=0; i<a->nz; i++) {
683*2c71b3e2SJacob Faibussowitsch     PetscAssertFalse(!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);
684ccb205f8SBarry Smith     ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
685ccb205f8SBarry Smith     ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
686ccb205f8SBarry Smith   }
6877d3de750SJacob Faibussowitsch   ierr = 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);CHKERRQ(ierr);
6887d3de750SJacob Faibussowitsch   ierr = PetscInfo(A,"Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n",a->reallocs);CHKERRQ(ierr);
6897d3de750SJacob Faibussowitsch   ierr = PetscInfo(A,"Maximum nonzeros in any row is %" PetscInt_FMT "\n",rmax);CHKERRQ(ierr);
6902205254eSKarl Rupp 
6918e58a170SBarry Smith   A->info.mallocs    += a->reallocs;
692421e10b8SBarry Smith   a->reallocs         = 0;
693421e10b8SBarry Smith   A->info.nz_unneeded = (double)fshift;
694421e10b8SBarry Smith   a->rmax             = rmax;
695ccb205f8SBarry Smith   ierr                = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr);
696421e10b8SBarry Smith   PetscFunctionReturn(0);
697421e10b8SBarry Smith }
698421e10b8SBarry Smith 
69981f0254dSBarry Smith static PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg)
700290bbb0aSBarry Smith {
701994fe344SLisandro Dalcin   PetscErrorCode ierr;
702290bbb0aSBarry Smith   PetscFunctionBegin;
7034e0d8c25SBarry Smith   if (opt == MAT_SYMMETRIC && flg) {
70441f059aeSBarry Smith     A->ops->sor  = MatSOR_BlockMat_Symmetric;
705290bbb0aSBarry Smith     A->ops->mult = MatMult_BlockMat_Symmetric;
706290bbb0aSBarry Smith   } else {
7077d3de750SJacob Faibussowitsch     ierr = PetscInfo(A,"Unused matrix option %s\n",MatOptions[opt]);CHKERRQ(ierr);
708290bbb0aSBarry Smith   }
709290bbb0aSBarry Smith   PetscFunctionReturn(0);
710290bbb0aSBarry Smith }
711290bbb0aSBarry Smith 
7123964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
713f4259b30SLisandro Dalcin                                        NULL,
714f4259b30SLisandro Dalcin                                        NULL,
715421e10b8SBarry Smith                                        MatMult_BlockMat,
716421e10b8SBarry Smith                                /*  4*/ MatMultAdd_BlockMat,
717421e10b8SBarry Smith                                        MatMultTranspose_BlockMat,
718421e10b8SBarry Smith                                        MatMultTransposeAdd_BlockMat,
719f4259b30SLisandro Dalcin                                        NULL,
720f4259b30SLisandro Dalcin                                        NULL,
721f4259b30SLisandro Dalcin                                        NULL,
722f4259b30SLisandro Dalcin                                /* 10*/ NULL,
723f4259b30SLisandro Dalcin                                        NULL,
724f4259b30SLisandro Dalcin                                        NULL,
72541f059aeSBarry Smith                                        MatSOR_BlockMat,
726f4259b30SLisandro Dalcin                                        NULL,
727f4259b30SLisandro Dalcin                                /* 15*/ NULL,
728f4259b30SLisandro Dalcin                                        NULL,
729f4259b30SLisandro Dalcin                                        NULL,
730f4259b30SLisandro Dalcin                                        NULL,
731f4259b30SLisandro Dalcin                                        NULL,
732f4259b30SLisandro Dalcin                                /* 20*/ NULL,
733421e10b8SBarry Smith                                        MatAssemblyEnd_BlockMat,
734290bbb0aSBarry Smith                                        MatSetOption_BlockMat,
735f4259b30SLisandro Dalcin                                        NULL,
736f4259b30SLisandro Dalcin                                /* 24*/ NULL,
737f4259b30SLisandro Dalcin                                        NULL,
738f4259b30SLisandro Dalcin                                        NULL,
739f4259b30SLisandro Dalcin                                        NULL,
740f4259b30SLisandro Dalcin                                        NULL,
741f4259b30SLisandro Dalcin                                /* 29*/ NULL,
742f4259b30SLisandro Dalcin                                        NULL,
743f4259b30SLisandro Dalcin                                        NULL,
744f4259b30SLisandro Dalcin                                        NULL,
745f4259b30SLisandro Dalcin                                        NULL,
746f4259b30SLisandro Dalcin                                /* 34*/ NULL,
747f4259b30SLisandro Dalcin                                        NULL,
748f4259b30SLisandro Dalcin                                        NULL,
749f4259b30SLisandro Dalcin                                        NULL,
750f4259b30SLisandro Dalcin                                        NULL,
751f4259b30SLisandro Dalcin                                /* 39*/ NULL,
752f4259b30SLisandro Dalcin                                        NULL,
753f4259b30SLisandro Dalcin                                        NULL,
754f4259b30SLisandro Dalcin                                        NULL,
755f4259b30SLisandro Dalcin                                        NULL,
756f4259b30SLisandro Dalcin                                /* 44*/ NULL,
757f4259b30SLisandro Dalcin                                        NULL,
7587d68702bSBarry Smith                                        MatShift_Basic,
759f4259b30SLisandro Dalcin                                        NULL,
760f4259b30SLisandro Dalcin                                        NULL,
761f4259b30SLisandro Dalcin                                /* 49*/ NULL,
762f4259b30SLisandro Dalcin                                        NULL,
763f4259b30SLisandro Dalcin                                        NULL,
764f4259b30SLisandro Dalcin                                        NULL,
765f4259b30SLisandro Dalcin                                        NULL,
766f4259b30SLisandro Dalcin                                /* 54*/ NULL,
767f4259b30SLisandro Dalcin                                        NULL,
768f4259b30SLisandro Dalcin                                        NULL,
769f4259b30SLisandro Dalcin                                        NULL,
770f4259b30SLisandro Dalcin                                        NULL,
7717dae84e0SHong Zhang                                /* 59*/ MatCreateSubMatrix_BlockMat,
772421e10b8SBarry Smith                                        MatDestroy_BlockMat,
773ccb205f8SBarry Smith                                        MatView_BlockMat,
774f4259b30SLisandro Dalcin                                        NULL,
775f4259b30SLisandro Dalcin                                        NULL,
776f4259b30SLisandro Dalcin                                /* 64*/ NULL,
777f4259b30SLisandro Dalcin                                        NULL,
778f4259b30SLisandro Dalcin                                        NULL,
779f4259b30SLisandro Dalcin                                        NULL,
780f4259b30SLisandro Dalcin                                        NULL,
781f4259b30SLisandro Dalcin                                /* 69*/ NULL,
782f4259b30SLisandro Dalcin                                        NULL,
783f4259b30SLisandro Dalcin                                        NULL,
784f4259b30SLisandro Dalcin                                        NULL,
785f4259b30SLisandro Dalcin                                        NULL,
786f4259b30SLisandro Dalcin                                /* 74*/ NULL,
787f4259b30SLisandro Dalcin                                        NULL,
788f4259b30SLisandro Dalcin                                        NULL,
789f4259b30SLisandro Dalcin                                        NULL,
790f4259b30SLisandro Dalcin                                        NULL,
791f4259b30SLisandro Dalcin                                /* 79*/ NULL,
792f4259b30SLisandro Dalcin                                        NULL,
793f4259b30SLisandro Dalcin                                        NULL,
794f4259b30SLisandro Dalcin                                        NULL,
7955bba2384SShri Abhyankar                                        MatLoad_BlockMat,
796f4259b30SLisandro Dalcin                                /* 84*/ NULL,
797f4259b30SLisandro Dalcin                                        NULL,
798f4259b30SLisandro Dalcin                                        NULL,
799f4259b30SLisandro Dalcin                                        NULL,
800f4259b30SLisandro Dalcin                                        NULL,
801f4259b30SLisandro Dalcin                                /* 89*/ NULL,
802f4259b30SLisandro Dalcin                                        NULL,
803f4259b30SLisandro Dalcin                                        NULL,
804f4259b30SLisandro Dalcin                                        NULL,
805f4259b30SLisandro Dalcin                                        NULL,
806f4259b30SLisandro Dalcin                                /* 94*/ NULL,
807f4259b30SLisandro Dalcin                                        NULL,
808f4259b30SLisandro Dalcin                                        NULL,
809f4259b30SLisandro Dalcin                                        NULL,
810f4259b30SLisandro Dalcin                                        NULL,
811f4259b30SLisandro Dalcin                                /* 99*/ NULL,
812f4259b30SLisandro Dalcin                                        NULL,
813f4259b30SLisandro Dalcin                                        NULL,
814f4259b30SLisandro Dalcin                                        NULL,
815f4259b30SLisandro Dalcin                                        NULL,
816f4259b30SLisandro Dalcin                                /*104*/ NULL,
817f4259b30SLisandro Dalcin                                        NULL,
818f4259b30SLisandro Dalcin                                        NULL,
819f4259b30SLisandro Dalcin                                        NULL,
820f4259b30SLisandro Dalcin                                        NULL,
821f4259b30SLisandro Dalcin                                /*109*/ NULL,
822f4259b30SLisandro Dalcin                                        NULL,
823f4259b30SLisandro Dalcin                                        NULL,
824f4259b30SLisandro Dalcin                                        NULL,
825f4259b30SLisandro Dalcin                                        NULL,
826f4259b30SLisandro Dalcin                                /*114*/ NULL,
827f4259b30SLisandro Dalcin                                        NULL,
828f4259b30SLisandro Dalcin                                        NULL,
829f4259b30SLisandro Dalcin                                        NULL,
830f4259b30SLisandro Dalcin                                        NULL,
831f4259b30SLisandro Dalcin                                /*119*/ NULL,
832f4259b30SLisandro Dalcin                                        NULL,
833f4259b30SLisandro Dalcin                                        NULL,
834f4259b30SLisandro Dalcin                                        NULL,
835f4259b30SLisandro Dalcin                                        NULL,
836f4259b30SLisandro Dalcin                                /*124*/ NULL,
837f4259b30SLisandro Dalcin                                        NULL,
838f4259b30SLisandro Dalcin                                        NULL,
839f4259b30SLisandro Dalcin                                        NULL,
840f4259b30SLisandro Dalcin                                        NULL,
841f4259b30SLisandro Dalcin                                /*129*/ NULL,
842f4259b30SLisandro Dalcin                                        NULL,
843f4259b30SLisandro Dalcin                                        NULL,
844f4259b30SLisandro Dalcin                                        NULL,
845f4259b30SLisandro Dalcin                                        NULL,
846f4259b30SLisandro Dalcin                                /*134*/ NULL,
847f4259b30SLisandro Dalcin                                        NULL,
848f4259b30SLisandro Dalcin                                        NULL,
849f4259b30SLisandro Dalcin                                        NULL,
850f4259b30SLisandro Dalcin                                        NULL,
851f4259b30SLisandro Dalcin                                /*139*/ NULL,
852f4259b30SLisandro Dalcin                                        NULL,
853f4259b30SLisandro Dalcin                                        NULL
854a5973382SShri Abhyankar };
855421e10b8SBarry Smith 
8568cdf2d9bSBarry Smith /*@C
8578cdf2d9bSBarry Smith    MatBlockMatSetPreallocation - For good matrix assembly performance
8588cdf2d9bSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
8598cdf2d9bSBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
8608cdf2d9bSBarry Smith    during matrix assembly can be increased by more than a factor of 50.
8618cdf2d9bSBarry Smith 
862d083f849SBarry Smith    Collective
8638cdf2d9bSBarry Smith 
8648cdf2d9bSBarry Smith    Input Parameters:
8658cdf2d9bSBarry Smith +  B - The matrix
8668cdf2d9bSBarry Smith .  bs - size of each block in matrix
8678cdf2d9bSBarry Smith .  nz - number of nonzeros per block row (same for all rows)
8688cdf2d9bSBarry Smith -  nnz - array containing the number of nonzeros in the various block rows
8690298fd71SBarry Smith          (possibly different for each row) or NULL
8708cdf2d9bSBarry Smith 
8718cdf2d9bSBarry Smith    Notes:
8728cdf2d9bSBarry Smith      If nnz is given then nz is ignored
8738cdf2d9bSBarry Smith 
8748cdf2d9bSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
8750298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
8768cdf2d9bSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
8778cdf2d9bSBarry Smith    will get TERRIBLE performance, see the users' manual chapter on matrices.
8788cdf2d9bSBarry Smith 
8798cdf2d9bSBarry Smith    Level: intermediate
8808cdf2d9bSBarry Smith 
8818cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues()
8828cdf2d9bSBarry Smith 
8838cdf2d9bSBarry Smith @*/
8847087cfbeSBarry Smith PetscErrorCode  MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
8858cdf2d9bSBarry Smith {
8864ac538c5SBarry Smith   PetscErrorCode ierr;
8878cdf2d9bSBarry Smith 
8888cdf2d9bSBarry Smith   PetscFunctionBegin;
8894ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr);
8908cdf2d9bSBarry Smith   PetscFunctionReturn(0);
8918cdf2d9bSBarry Smith }
8928cdf2d9bSBarry Smith 
89381f0254dSBarry Smith static PetscErrorCode  MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz)
8948cdf2d9bSBarry Smith {
8958cdf2d9bSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
8968cdf2d9bSBarry Smith   PetscErrorCode ierr;
8978cdf2d9bSBarry Smith   PetscInt       i;
8988cdf2d9bSBarry Smith 
8998cdf2d9bSBarry Smith   PetscFunctionBegin;
900e02043d6SBarry Smith   ierr = PetscLayoutSetBlockSize(A->rmap,bs);CHKERRQ(ierr);
901e02043d6SBarry Smith   ierr = PetscLayoutSetBlockSize(A->cmap,bs);CHKERRQ(ierr);
90234ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
90334ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
904e02043d6SBarry Smith   ierr = PetscLayoutGetBlockSize(A->rmap,&bs);CHKERRQ(ierr);
90534ef9618SShri Abhyankar 
9068cdf2d9bSBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
907*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nz < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %" PetscInt_FMT,nz);
9088cdf2d9bSBarry Smith   if (nnz) {
909d0f46423SBarry Smith     for (i=0; i<A->rmap->n/bs; i++) {
910*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(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]);
911*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(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);
9128cdf2d9bSBarry Smith     }
9138cdf2d9bSBarry Smith   }
914d0f46423SBarry Smith   bmat->mbs = A->rmap->n/bs;
9158cdf2d9bSBarry Smith 
9160298fd71SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->right);CHKERRQ(ierr);
9170298fd71SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->middle);CHKERRQ(ierr);
9188cdf2d9bSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr);
9198cdf2d9bSBarry Smith 
9202ee49352SLisandro Dalcin   if (!bmat->imax) {
921dcca6d9dSJed Brown     ierr = PetscMalloc2(A->rmap->n,&bmat->imax,A->rmap->n,&bmat->ilen);CHKERRQ(ierr);
9223bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
9232ee49352SLisandro Dalcin   }
924c0aa6a63SJacob Faibussowitsch   if (PetscLikely(nnz)) {
9258cdf2d9bSBarry Smith     nz = 0;
926d0f46423SBarry Smith     for (i=0; i<A->rmap->n/A->rmap->bs; i++) {
9278cdf2d9bSBarry Smith       bmat->imax[i] = nnz[i];
9288cdf2d9bSBarry Smith       nz           += nnz[i];
9298cdf2d9bSBarry Smith     }
930f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation");
9318cdf2d9bSBarry Smith 
9328cdf2d9bSBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
933c0aa6a63SJacob Faibussowitsch   ierr = PetscArrayzero(bmat->ilen,bmat->mbs);CHKERRQ(ierr);
9348cdf2d9bSBarry Smith 
9358cdf2d9bSBarry Smith   /* allocate the matrix space */
9362ee49352SLisandro Dalcin   ierr       = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
937dcca6d9dSJed Brown   ierr       = PetscMalloc3(nz,&bmat->a,nz,&bmat->j,A->rmap->n+1,&bmat->i);CHKERRQ(ierr);
9383bb1ff40SBarry Smith   ierr       = PetscLogObjectMemory((PetscObject)A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
9398cdf2d9bSBarry Smith   bmat->i[0] = 0;
9408cdf2d9bSBarry Smith   for (i=1; i<bmat->mbs+1; i++) {
9418cdf2d9bSBarry Smith     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
9428cdf2d9bSBarry Smith   }
9438cdf2d9bSBarry Smith   bmat->singlemalloc = PETSC_TRUE;
9448cdf2d9bSBarry Smith   bmat->free_a       = PETSC_TRUE;
9458cdf2d9bSBarry Smith   bmat->free_ij      = PETSC_TRUE;
9468cdf2d9bSBarry Smith 
9478cdf2d9bSBarry Smith   bmat->nz            = 0;
9488cdf2d9bSBarry Smith   bmat->maxnz         = nz;
9498cdf2d9bSBarry Smith   A->info.nz_unneeded = (double)bmat->maxnz;
9507827cd58SJed Brown   ierr                = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
9518cdf2d9bSBarry Smith   PetscFunctionReturn(0);
9528cdf2d9bSBarry Smith }
9538cdf2d9bSBarry Smith 
954421e10b8SBarry Smith /*MC
955421e10b8SBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
956421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
957421e10b8SBarry Smith 
958421e10b8SBarry Smith   Level: advanced
959421e10b8SBarry Smith 
960421e10b8SBarry Smith .seealso: MatCreateBlockMat()
961421e10b8SBarry Smith 
962421e10b8SBarry Smith M*/
963421e10b8SBarry Smith 
9648cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A)
965421e10b8SBarry Smith {
966421e10b8SBarry Smith   Mat_BlockMat   *b;
967421e10b8SBarry Smith   PetscErrorCode ierr;
968421e10b8SBarry Smith 
969421e10b8SBarry Smith   PetscFunctionBegin;
970b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
971421e10b8SBarry Smith   A->data = (void*)b;
97238f2d2fdSLisandro Dalcin   ierr    = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
973421e10b8SBarry Smith 
974421e10b8SBarry Smith   A->assembled    = PETSC_TRUE;
975421e10b8SBarry Smith   A->preallocated = PETSC_FALSE;
976421e10b8SBarry Smith   ierr            = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
977290bbb0aSBarry Smith 
978bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatBlockMatSetPreallocation_C",MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr);
979421e10b8SBarry Smith   PetscFunctionReturn(0);
980421e10b8SBarry Smith }
981421e10b8SBarry Smith 
982421e10b8SBarry Smith /*@C
983bbd3f336SJed Brown    MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential Mat object
984421e10b8SBarry Smith 
985d083f849SBarry Smith   Collective
986421e10b8SBarry Smith 
987421e10b8SBarry Smith    Input Parameters:
988421e10b8SBarry Smith +  comm - MPI communicator
989421e10b8SBarry Smith .  m - number of rows
990421e10b8SBarry Smith .  n  - number of columns
9918cdf2d9bSBarry Smith .  bs - size of each submatrix
9928cdf2d9bSBarry Smith .  nz  - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known)
9930298fd71SBarry Smith -  nnz - expected number of nonzers per block row if known (use NULL otherwise)
994421e10b8SBarry Smith 
995421e10b8SBarry Smith    Output Parameter:
996421e10b8SBarry Smith .  A - the matrix
997421e10b8SBarry Smith 
998421e10b8SBarry Smith    Level: intermediate
999421e10b8SBarry Smith 
100095452b02SPatrick Sanan    Notes:
100195452b02SPatrick Sanan     Matrices of this type are nominally-sparse matrices in which each "entry" is a Mat object.  Each Mat must
1002bbd3f336SJed Brown    have the same size and be sequential.  The local and global sizes must be compatible with this decomposition.
1003bbd3f336SJed Brown 
1004bbd3f336SJed Brown    For matrices containing parallel submatrices and variable block sizes, see MATNEST.
1005421e10b8SBarry Smith 
1006bbd3f336SJed Brown .seealso: MATBLOCKMAT, MatCreateNest()
1007421e10b8SBarry Smith @*/
10087087cfbeSBarry Smith PetscErrorCode  MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A)
1009421e10b8SBarry Smith {
1010421e10b8SBarry Smith   PetscErrorCode ierr;
1011421e10b8SBarry Smith 
1012421e10b8SBarry Smith   PetscFunctionBegin;
1013421e10b8SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1014421e10b8SBarry Smith   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1015421e10b8SBarry Smith   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
10168cdf2d9bSBarry Smith   ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1017421e10b8SBarry Smith   PetscFunctionReturn(0);
1018421e10b8SBarry Smith }
1019