xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1421e10b8SBarry Smith 
2421e10b8SBarry Smith /*
3421e10b8SBarry Smith    This provides a matrix that consists of Mats
4421e10b8SBarry Smith */
5421e10b8SBarry Smith 
6b45d2f2cSJed Brown #include <petsc-private/matimpl.h>              /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>    /* use the common AIJ data-structure */
8421e10b8SBarry Smith 
9421e10b8SBarry Smith typedef struct {
10421e10b8SBarry Smith   SEQAIJHEADER(Mat);
11421e10b8SBarry Smith   SEQBAIJHEADER;
12ccb205f8SBarry Smith   Mat *diags;
13421e10b8SBarry Smith 
146e63c7a1SBarry Smith   Vec left,right,middle,workb;                 /* dummy vectors to perform local parts of product */
15421e10b8SBarry Smith } Mat_BlockMat;
16421e10b8SBarry Smith 
177087cfbeSBarry Smith extern PetscErrorCode  MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt*);
18a5973382SShri Abhyankar 
19421e10b8SBarry Smith #undef __FUNCT__
2041f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat_Symmetric"
2141f059aeSBarry Smith PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
22290bbb0aSBarry Smith {
23290bbb0aSBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
24290bbb0aSBarry Smith   PetscScalar       *x;
25a2ea699eSBarry Smith   const Mat         *v;
26290bbb0aSBarry Smith   const PetscScalar *b;
27290bbb0aSBarry Smith   PetscErrorCode    ierr;
28d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
29290bbb0aSBarry Smith   const PetscInt    *idx;
30290bbb0aSBarry Smith   IS                row,col;
31290bbb0aSBarry Smith   MatFactorInfo     info;
326e63c7a1SBarry Smith   Vec               left = a->left,right = a->right, middle = a->middle;
33290bbb0aSBarry Smith   Mat               *diag;
34290bbb0aSBarry Smith 
35290bbb0aSBarry Smith   PetscFunctionBegin;
36290bbb0aSBarry Smith   its = its*lits;
37e32f2f54SBarry Smith   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
38e32f2f54SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
39e32f2f54SBarry Smith   if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
40e32f2f54SBarry Smith   if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
41290bbb0aSBarry Smith   if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP))
42e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
43290bbb0aSBarry Smith 
44290bbb0aSBarry Smith   if (!a->diags) {
45290bbb0aSBarry Smith     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
46290bbb0aSBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
47290bbb0aSBarry Smith     for (i=0; i<mbs; i++) {
482692d6eeSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr);
49719d5645SBarry Smith       ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr);
50719d5645SBarry Smith       ierr = MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr);
516bf464f9SBarry Smith       ierr = ISDestroy(&row);CHKERRQ(ierr);
526bf464f9SBarry Smith       ierr = ISDestroy(&col);CHKERRQ(ierr);
53290bbb0aSBarry Smith     }
546e63c7a1SBarry Smith     ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr);
55290bbb0aSBarry Smith   }
56290bbb0aSBarry Smith   diag = a->diags;
57290bbb0aSBarry Smith 
58290bbb0aSBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
59290bbb0aSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
60290bbb0aSBarry Smith   /* copy right hand side because it must be modified during iteration */
616e63c7a1SBarry Smith   ierr = VecCopy(bb,a->workb);CHKERRQ(ierr);
623649974fSBarry Smith   ierr = VecGetArrayRead(a->workb,&b);CHKERRQ(ierr);
63290bbb0aSBarry Smith 
6441f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
65290bbb0aSBarry Smith   while (its--) {
66290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
67290bbb0aSBarry Smith 
68290bbb0aSBarry Smith       for (i=0; i<mbs; i++) {
696e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
706e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
716e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
72290bbb0aSBarry Smith 
73290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
74290bbb0aSBarry Smith         for (j=0; j<n; j++) {
75290bbb0aSBarry Smith           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
76290bbb0aSBarry Smith           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
77290bbb0aSBarry Smith           ierr = VecResetArray(right);CHKERRQ(ierr);
78290bbb0aSBarry Smith         }
79290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
80290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
81290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
82290bbb0aSBarry Smith 
83290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
84290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
856e63c7a1SBarry Smith 
8641f059aeSBarry Smith         /* now adjust right hand side, see MatSOR_SeqSBAIJ */
876e63c7a1SBarry Smith         for (j=0; j<n; j++) {
886e63c7a1SBarry Smith           ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr);
896e63c7a1SBarry Smith           ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr);
906e63c7a1SBarry Smith           ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr);
916e63c7a1SBarry Smith           ierr = VecResetArray(middle);CHKERRQ(ierr);
926e63c7a1SBarry Smith         }
93290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
946e63c7a1SBarry Smith 
95290bbb0aSBarry Smith       }
96290bbb0aSBarry Smith     }
97290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
98290bbb0aSBarry Smith 
99290bbb0aSBarry Smith       for (i=mbs-1; i>=0; i--) {
1006e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
1016e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
1026e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
103290bbb0aSBarry Smith 
104290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
105290bbb0aSBarry Smith         for (j=0; j<n; j++) {
106290bbb0aSBarry Smith           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
107290bbb0aSBarry Smith           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
108290bbb0aSBarry Smith           ierr = VecResetArray(right);CHKERRQ(ierr);
109290bbb0aSBarry Smith         }
110290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
111290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
112290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
113290bbb0aSBarry Smith 
114290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
115290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
116290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
117290bbb0aSBarry Smith 
118290bbb0aSBarry Smith       }
119290bbb0aSBarry Smith     }
120290bbb0aSBarry Smith   }
121290bbb0aSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1223649974fSBarry Smith   ierr = VecRestoreArrayRead(a->workb,&b);CHKERRQ(ierr);
123290bbb0aSBarry Smith   PetscFunctionReturn(0);
124290bbb0aSBarry Smith }
125290bbb0aSBarry Smith 
126290bbb0aSBarry Smith #undef __FUNCT__
12741f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat"
12841f059aeSBarry Smith PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
129ccb205f8SBarry Smith {
130ccb205f8SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
131ccb205f8SBarry Smith   PetscScalar       *x;
132a2ea699eSBarry Smith   const Mat         *v;
133ccb205f8SBarry Smith   const PetscScalar *b;
134ccb205f8SBarry Smith   PetscErrorCode    ierr;
135d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
136ccb205f8SBarry Smith   const PetscInt    *idx;
137ccb205f8SBarry Smith   IS                row,col;
138ccb205f8SBarry Smith   MatFactorInfo     info;
139ccb205f8SBarry Smith   Vec               left = a->left,right = a->right;
140ccb205f8SBarry Smith   Mat               *diag;
141ccb205f8SBarry Smith 
142ccb205f8SBarry Smith   PetscFunctionBegin;
143ccb205f8SBarry Smith   its = its*lits;
144e32f2f54SBarry Smith   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
145e32f2f54SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
146e32f2f54SBarry Smith   if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
147e32f2f54SBarry Smith   if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
148ccb205f8SBarry Smith 
149ccb205f8SBarry Smith   if (!a->diags) {
150ccb205f8SBarry Smith     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
151ccb205f8SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
152ccb205f8SBarry Smith     for (i=0; i<mbs; i++) {
1532692d6eeSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr);
154719d5645SBarry Smith       ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr);
155719d5645SBarry Smith       ierr = MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr);
1566bf464f9SBarry Smith       ierr = ISDestroy(&row);CHKERRQ(ierr);
1576bf464f9SBarry Smith       ierr = ISDestroy(&col);CHKERRQ(ierr);
158ccb205f8SBarry Smith     }
159ccb205f8SBarry Smith   }
160ccb205f8SBarry Smith   diag = a->diags;
161ccb205f8SBarry Smith 
162ccb205f8SBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
163ccb205f8SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1643649974fSBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
165ccb205f8SBarry Smith 
16641f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
167ccb205f8SBarry Smith   while (its--) {
168ccb205f8SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
169ccb205f8SBarry Smith 
170ccb205f8SBarry Smith       for (i=0; i<mbs; i++) {
171ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
172ccb205f8SBarry Smith         idx = a->j + a->i[i];
173ccb205f8SBarry Smith         v   = a->a + a->i[i];
174ccb205f8SBarry Smith 
175ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
176ccb205f8SBarry Smith         for (j=0; j<n; j++) {
177ccb205f8SBarry Smith           if (idx[j] != i) {
178ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
179ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
180ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
181ccb205f8SBarry Smith           }
182ccb205f8SBarry Smith         }
183ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
184ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
185ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
186ccb205f8SBarry Smith 
187ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
188ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
189ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
190ccb205f8SBarry Smith       }
191ccb205f8SBarry Smith     }
192ccb205f8SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
193ccb205f8SBarry Smith 
194ccb205f8SBarry Smith       for (i=mbs-1; i>=0; i--) {
195ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
196ccb205f8SBarry Smith         idx = a->j + a->i[i];
197ccb205f8SBarry Smith         v   = a->a + a->i[i];
198ccb205f8SBarry Smith 
199ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
200ccb205f8SBarry Smith         for (j=0; j<n; j++) {
201ccb205f8SBarry Smith           if (idx[j] != i) {
202ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
203ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
204ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
205ccb205f8SBarry Smith           }
206ccb205f8SBarry Smith         }
207ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
208ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
209ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
210ccb205f8SBarry Smith 
211ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
212ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
213ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
214ccb205f8SBarry Smith 
215ccb205f8SBarry Smith       }
216ccb205f8SBarry Smith     }
217ccb205f8SBarry Smith   }
218ccb205f8SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2193649974fSBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
220ccb205f8SBarry Smith   PetscFunctionReturn(0);
221ccb205f8SBarry Smith }
222ccb205f8SBarry Smith 
223ccb205f8SBarry Smith #undef __FUNCT__
224421e10b8SBarry Smith #define __FUNCT__ "MatSetValues_BlockMat"
225421e10b8SBarry Smith PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
226421e10b8SBarry Smith {
227421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
228421e10b8SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
229421e10b8SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
230d0f46423SBarry Smith   PetscInt       *aj  =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
231421e10b8SBarry Smith   PetscErrorCode ierr;
232421e10b8SBarry Smith   PetscInt       ridx,cidx;
233ace3abfcSBarry Smith   PetscBool      roworiented=a->roworiented;
234421e10b8SBarry Smith   MatScalar      value;
235421e10b8SBarry Smith   Mat            *ap,*aa = a->a;
236421e10b8SBarry Smith 
237421e10b8SBarry Smith   PetscFunctionBegin;
23871fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
239421e10b8SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
240421e10b8SBarry Smith     row  = im[k];
241421e10b8SBarry Smith     brow = row/bs;
242421e10b8SBarry Smith     if (row < 0) continue;
243421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
244e32f2f54SBarry Smith     if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
245421e10b8SBarry Smith #endif
246421e10b8SBarry Smith     rp   = aj + ai[brow];
247421e10b8SBarry Smith     ap   = aa + ai[brow];
248421e10b8SBarry Smith     rmax = imax[brow];
249421e10b8SBarry Smith     nrow = ailen[brow];
250421e10b8SBarry Smith     low  = 0;
251421e10b8SBarry Smith     high = nrow;
252421e10b8SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
253421e10b8SBarry Smith       if (in[l] < 0) continue;
254421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
255e32f2f54SBarry Smith       if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1);
256421e10b8SBarry Smith #endif
257421e10b8SBarry Smith       col = in[l]; bcol = col/bs;
2586e63c7a1SBarry Smith       if (A->symmetric && brow > bcol) continue;
259421e10b8SBarry Smith       ridx = row % bs; cidx = col % bs;
260*2205254eSKarl Rupp       if (roworiented) value = v[l + k*n];
261*2205254eSKarl Rupp       else value = v[k + l*m];
262*2205254eSKarl Rupp 
263*2205254eSKarl Rupp       if (col <= lastcol) low = 0;
264*2205254eSKarl Rupp       else high = nrow;
265421e10b8SBarry Smith       lastcol = col;
266421e10b8SBarry Smith       while (high-low > 7) {
267421e10b8SBarry Smith         t = (low+high)/2;
268421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
269421e10b8SBarry Smith         else              low  = t;
270421e10b8SBarry Smith       }
271421e10b8SBarry Smith       for (i=low; i<high; i++) {
272421e10b8SBarry Smith         if (rp[i] > bcol) break;
273*2205254eSKarl Rupp         if (rp[i] == bcol) goto noinsert1;
274421e10b8SBarry Smith       }
275421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
276e32f2f54SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
277fef13f97SBarry Smith       MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
278421e10b8SBarry Smith       N = nrow++ - 1; high++;
279421e10b8SBarry Smith       /* shift up all the later entries in this row */
280421e10b8SBarry Smith       for (ii=N; ii>=i; ii--) {
281421e10b8SBarry Smith         rp[ii+1] = rp[ii];
282421e10b8SBarry Smith         ap[ii+1] = ap[ii];
283421e10b8SBarry Smith       }
284421e10b8SBarry Smith       if (N>=i) ap[i] = 0;
285421e10b8SBarry Smith       rp[i] = bcol;
286421e10b8SBarry Smith       a->nz++;
287421e10b8SBarry Smith noinsert1:;
288421e10b8SBarry Smith       if (!*(ap+i)) {
2896e63c7a1SBarry Smith         ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
290b0223207SBarry Smith       }
291421e10b8SBarry Smith       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
292421e10b8SBarry Smith       low  = i;
293421e10b8SBarry Smith     }
294421e10b8SBarry Smith     ailen[brow] = nrow;
295421e10b8SBarry Smith   }
296421e10b8SBarry Smith   A->same_nonzero = PETSC_FALSE;
297421e10b8SBarry Smith   PetscFunctionReturn(0);
298421e10b8SBarry Smith }
299421e10b8SBarry Smith 
300421e10b8SBarry Smith #undef __FUNCT__
3015bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_BlockMat"
302112444f4SShri Abhyankar PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
303a5973382SShri Abhyankar {
304a5973382SShri Abhyankar   PetscErrorCode    ierr;
305a5973382SShri Abhyankar   Mat               tmpA;
306a5973382SShri Abhyankar   PetscInt          i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0;
307a5973382SShri Abhyankar   const PetscInt    *cols;
308a5973382SShri Abhyankar   const PetscScalar *values;
309ace3abfcSBarry Smith   PetscBool         flg = PETSC_FALSE,notdone;
310a5973382SShri Abhyankar   Mat_SeqAIJ        *a;
311a5973382SShri Abhyankar   Mat_BlockMat      *amat;
312a5973382SShri Abhyankar 
313a5973382SShri Abhyankar   PetscFunctionBegin;
314a5973382SShri Abhyankar   ierr = MatCreate(PETSC_COMM_SELF,&tmpA);CHKERRQ(ierr);
315a5973382SShri Abhyankar   ierr = MatSetType(tmpA,MATSEQAIJ);CHKERRQ(ierr);
316112444f4SShri Abhyankar   ierr = MatLoad_SeqAIJ(tmpA,viewer);CHKERRQ(ierr);
317a5973382SShri Abhyankar 
318a5973382SShri Abhyankar   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
319a5973382SShri Abhyankar   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr);
320a5973382SShri Abhyankar   ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr);
321acfcf0e5SJed Brown   ierr = PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
322a5973382SShri Abhyankar   ierr = PetscOptionsEnd();CHKERRQ(ierr);
323a5973382SShri Abhyankar 
324a5973382SShri Abhyankar   /* Determine number of nonzero blocks for each block row */
325a5973382SShri Abhyankar   a    = (Mat_SeqAIJ*) tmpA->data;
326a5973382SShri Abhyankar   mbs  = m/bs;
327a5973382SShri Abhyankar   ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt*,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr);
328a5973382SShri Abhyankar   ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr);
329a5973382SShri Abhyankar 
330a5973382SShri Abhyankar   for (i=0; i<mbs; i++) {
331a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
332a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
333a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
334a5973382SShri Abhyankar     }
335a5973382SShri Abhyankar 
336a5973382SShri Abhyankar     currentcol = -1;
337a5973382SShri Abhyankar     notdone    = PETSC_TRUE;
338a5973382SShri Abhyankar     while (PETSC_TRUE) {
339a5973382SShri Abhyankar       notdone = PETSC_FALSE;
340a5973382SShri Abhyankar       nextcol = 1000000000;
341a5973382SShri Abhyankar       for (j=0; j<bs; j++) {
342a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) {
343a5973382SShri Abhyankar           ii[j]++;
344a5973382SShri Abhyankar           ilens[j]--;
345a5973382SShri Abhyankar         }
346a5973382SShri Abhyankar         if (ilens[j] > 0) {
347a5973382SShri Abhyankar           notdone = PETSC_TRUE;
348a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
349a5973382SShri Abhyankar         }
350a5973382SShri Abhyankar       }
351a5973382SShri Abhyankar       if (!notdone) break;
352a5973382SShri Abhyankar       if (!flg || (nextcol >= i)) lens[i]++;
353a5973382SShri Abhyankar       currentcol = nextcol;
354a5973382SShri Abhyankar     }
355a5973382SShri Abhyankar   }
356a5973382SShri Abhyankar 
357a5973382SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
358a5973382SShri Abhyankar     ierr = MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
359a5973382SShri Abhyankar   }
360a5973382SShri Abhyankar   ierr = MatBlockMatSetPreallocation(newmat,bs,0,lens);CHKERRQ(ierr);
361a5973382SShri Abhyankar   if (flg) {
362a5973382SShri Abhyankar     ierr = MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
363a5973382SShri Abhyankar   }
364a5973382SShri Abhyankar   amat = (Mat_BlockMat*)(newmat)->data;
365a5973382SShri Abhyankar 
366a5973382SShri Abhyankar   /* preallocate the submatrices */
367a5973382SShri Abhyankar   ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr);
368a5973382SShri Abhyankar   for (i=0; i<mbs; i++) { /* loops for block rows */
369a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
370a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
371a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
372a5973382SShri Abhyankar     }
373a5973382SShri Abhyankar 
374a5973382SShri Abhyankar     currentcol = 1000000000;
375a5973382SShri Abhyankar     for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */
376a5973382SShri Abhyankar       if (ilens[j] > 0) {
377a5973382SShri Abhyankar         currentcol = PetscMin(currentcol,ii[j][0]/bs);
378a5973382SShri Abhyankar       }
379a5973382SShri Abhyankar     }
380a5973382SShri Abhyankar 
381a5973382SShri Abhyankar     notdone = PETSC_TRUE;
382a5973382SShri Abhyankar     while (PETSC_TRUE) {  /* loops over blocks in block row */
383a5973382SShri Abhyankar 
384a5973382SShri Abhyankar       notdone = PETSC_FALSE;
385a5973382SShri Abhyankar       nextcol = 1000000000;
386a5973382SShri Abhyankar       ierr    = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr);
387a5973382SShri Abhyankar       for (j=0; j<bs; j++) { /* loop over rows in block */
388a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */
389a5973382SShri Abhyankar           ii[j]++;
390a5973382SShri Abhyankar           ilens[j]--;
391a5973382SShri Abhyankar           llens[j]++;
392a5973382SShri Abhyankar         }
393a5973382SShri Abhyankar         if (ilens[j] > 0) {
394a5973382SShri Abhyankar           notdone = PETSC_TRUE;
395a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
396a5973382SShri Abhyankar         }
397a5973382SShri Abhyankar       }
398a5973382SShri Abhyankar       if (cnt >= amat->maxnz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt);
399a5973382SShri Abhyankar       if (!flg || currentcol >= i) {
400a5973382SShri Abhyankar         amat->j[cnt] = currentcol;
401a5973382SShri Abhyankar         ierr         = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr);
402a5973382SShri Abhyankar       }
403a5973382SShri Abhyankar 
404a5973382SShri Abhyankar       if (!notdone) break;
405a5973382SShri Abhyankar       currentcol = nextcol;
406a5973382SShri Abhyankar     }
407a5973382SShri Abhyankar     amat->ilen[i] = lens[i];
408a5973382SShri Abhyankar   }
409a5973382SShri Abhyankar   CHKMEMQ;
410a5973382SShri Abhyankar 
411a5973382SShri Abhyankar   ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr);
412a5973382SShri Abhyankar   ierr = PetscFree(llens);CHKERRQ(ierr);
413a5973382SShri Abhyankar 
414a5973382SShri Abhyankar   /* copy over the matrix, one row at a time */
415a5973382SShri Abhyankar   for (i=0; i<m; i++) {
416a5973382SShri Abhyankar     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
417a5973382SShri Abhyankar     ierr = MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
418a5973382SShri Abhyankar     ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
419a5973382SShri Abhyankar   }
420a5973382SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
421a5973382SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
422a5973382SShri Abhyankar   PetscFunctionReturn(0);
423a5973382SShri Abhyankar }
424a5973382SShri Abhyankar 
425a5973382SShri Abhyankar #undef __FUNCT__
426ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat"
427ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
428ccb205f8SBarry Smith {
42964075487SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
430ccb205f8SBarry Smith   PetscErrorCode    ierr;
431ccb205f8SBarry Smith   const char        *name;
432ccb205f8SBarry Smith   PetscViewerFormat format;
433ccb205f8SBarry Smith 
434ccb205f8SBarry Smith   PetscFunctionBegin;
435ccb205f8SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
436ccb205f8SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
437ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
438ccb205f8SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr);
43964075487SBarry Smith     if (A->symmetric) {
4408c1ad04dSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr);
44164075487SBarry Smith     }
442ccb205f8SBarry Smith   }
443ccb205f8SBarry Smith   PetscFunctionReturn(0);
444ccb205f8SBarry Smith }
445421e10b8SBarry Smith 
446421e10b8SBarry Smith #undef __FUNCT__
447421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat"
448421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat)
449421e10b8SBarry Smith {
450421e10b8SBarry Smith   PetscErrorCode ierr;
451421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
452ccb205f8SBarry Smith   PetscInt       i;
453421e10b8SBarry Smith 
454421e10b8SBarry Smith   PetscFunctionBegin;
4556bf464f9SBarry Smith   ierr = VecDestroy(&bmat->right);CHKERRQ(ierr);
4566bf464f9SBarry Smith   ierr = VecDestroy(&bmat->left);CHKERRQ(ierr);
4576bf464f9SBarry Smith   ierr = VecDestroy(&bmat->middle);CHKERRQ(ierr);
4586bf464f9SBarry Smith   ierr = VecDestroy(&bmat->workb);CHKERRQ(ierr);
459ccb205f8SBarry Smith   if (bmat->diags) {
460d0f46423SBarry Smith     for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) {
4616bf464f9SBarry Smith       ierr = MatDestroy(&bmat->diags[i]);CHKERRQ(ierr);
462ccb205f8SBarry Smith     }
463ccb205f8SBarry Smith   }
464ccb205f8SBarry Smith   if (bmat->a) {
465ccb205f8SBarry Smith     for (i=0; i<bmat->nz; i++) {
4666bf464f9SBarry Smith       ierr = MatDestroy(&bmat->a[i]);CHKERRQ(ierr);
467ccb205f8SBarry Smith     }
468ccb205f8SBarry Smith   }
469ccb205f8SBarry Smith   ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
470bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
471421e10b8SBarry Smith   PetscFunctionReturn(0);
472421e10b8SBarry Smith }
473421e10b8SBarry Smith 
474421e10b8SBarry Smith #undef __FUNCT__
475421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat"
476421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
477421e10b8SBarry Smith {
478421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
479421e10b8SBarry Smith   PetscErrorCode ierr;
480421e10b8SBarry Smith   PetscScalar    *xx,*yy;
481d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
482421e10b8SBarry Smith   Mat            *aa;
483421e10b8SBarry Smith 
484421e10b8SBarry Smith   PetscFunctionBegin;
485ccb205f8SBarry Smith   CHKMEMQ;
486421e10b8SBarry Smith   /*
487421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
488421e10b8SBarry Smith   */
489421e10b8SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
490ccb205f8SBarry Smith 
491421e10b8SBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
492421e10b8SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
493421e10b8SBarry Smith   aj   = bmat->j;
494421e10b8SBarry Smith   aa   = bmat->a;
495421e10b8SBarry Smith   ii   = bmat->i;
496421e10b8SBarry Smith   for (i=0; i<m; i++) {
497421e10b8SBarry Smith     jrow = ii[i];
498ccb205f8SBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
499421e10b8SBarry Smith     n    = ii[i+1] - jrow;
500421e10b8SBarry Smith     for (j=0; j<n; j++) {
501421e10b8SBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
502ccb205f8SBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
503421e10b8SBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
504421e10b8SBarry Smith       jrow++;
505421e10b8SBarry Smith     }
506421e10b8SBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
507421e10b8SBarry Smith   }
508421e10b8SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
509421e10b8SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
510ccb205f8SBarry Smith   CHKMEMQ;
511421e10b8SBarry Smith   PetscFunctionReturn(0);
512421e10b8SBarry Smith }
513421e10b8SBarry Smith 
514421e10b8SBarry Smith #undef __FUNCT__
515290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric"
516290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
517290bbb0aSBarry Smith {
518290bbb0aSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
519290bbb0aSBarry Smith   PetscErrorCode ierr;
520290bbb0aSBarry Smith   PetscScalar    *xx,*yy;
521d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
522290bbb0aSBarry Smith   Mat            *aa;
523290bbb0aSBarry Smith 
524290bbb0aSBarry Smith   PetscFunctionBegin;
525290bbb0aSBarry Smith   CHKMEMQ;
526290bbb0aSBarry Smith   /*
527290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
528290bbb0aSBarry Smith   */
529290bbb0aSBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
530290bbb0aSBarry Smith 
531290bbb0aSBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
532290bbb0aSBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
533290bbb0aSBarry Smith   aj   = bmat->j;
534290bbb0aSBarry Smith   aa   = bmat->a;
535290bbb0aSBarry Smith   ii   = bmat->i;
536290bbb0aSBarry Smith   for (i=0; i<m; i++) {
537290bbb0aSBarry Smith     jrow = ii[i];
538290bbb0aSBarry Smith     n    = ii[i+1] - jrow;
539290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
540290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr);
541290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
542290bbb0aSBarry Smith     if (aj[jrow] == i) {
543290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
544290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
545290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
546290bbb0aSBarry Smith       jrow++;
547290bbb0aSBarry Smith       n--;
548290bbb0aSBarry Smith     }
549290bbb0aSBarry Smith     for (j=0; j<n; j++) {
550290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);            /* upper triangular part */
551290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
552290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
553290bbb0aSBarry Smith 
554290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr);            /* lower triangular part */
555290bbb0aSBarry Smith       ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr);
556290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
557290bbb0aSBarry Smith       jrow++;
558290bbb0aSBarry Smith     }
559290bbb0aSBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
560290bbb0aSBarry Smith     ierr = VecResetArray(bmat->middle);CHKERRQ(ierr);
561290bbb0aSBarry Smith   }
562290bbb0aSBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
563290bbb0aSBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
564290bbb0aSBarry Smith   CHKMEMQ;
565290bbb0aSBarry Smith   PetscFunctionReturn(0);
566290bbb0aSBarry Smith }
567290bbb0aSBarry Smith 
568290bbb0aSBarry Smith #undef __FUNCT__
569421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat"
570421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
571421e10b8SBarry Smith {
572421e10b8SBarry Smith   PetscFunctionBegin;
573421e10b8SBarry Smith   PetscFunctionReturn(0);
574421e10b8SBarry Smith }
575421e10b8SBarry Smith 
576421e10b8SBarry Smith #undef __FUNCT__
577421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat"
578421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
579421e10b8SBarry Smith {
580421e10b8SBarry Smith   PetscFunctionBegin;
581421e10b8SBarry Smith   PetscFunctionReturn(0);
582421e10b8SBarry Smith }
583421e10b8SBarry Smith 
584421e10b8SBarry Smith #undef __FUNCT__
585421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat"
586421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
587421e10b8SBarry Smith {
588421e10b8SBarry Smith   PetscFunctionBegin;
589421e10b8SBarry Smith   PetscFunctionReturn(0);
590421e10b8SBarry Smith }
591421e10b8SBarry Smith 
592ccb205f8SBarry Smith /*
593ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
594ccb205f8SBarry Smith */
595ccb205f8SBarry Smith #undef __FUNCT__
596ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat"
597ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
598ccb205f8SBarry Smith {
599ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
600ccb205f8SBarry Smith   PetscErrorCode ierr;
601d0f46423SBarry Smith   PetscInt       i,j,mbs = A->rmap->n/A->rmap->bs;
602ccb205f8SBarry Smith 
603ccb205f8SBarry Smith   PetscFunctionBegin;
604ccb205f8SBarry Smith   if (!a->diag) {
605ccb205f8SBarry Smith     ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
606ccb205f8SBarry Smith   }
607ccb205f8SBarry Smith   for (i=0; i<mbs; i++) {
608ccb205f8SBarry Smith     a->diag[i] = a->i[i+1];
609ccb205f8SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
610ccb205f8SBarry Smith       if (a->j[j] == i) {
611ccb205f8SBarry Smith         a->diag[i] = j;
612ccb205f8SBarry Smith         break;
613ccb205f8SBarry Smith       }
614ccb205f8SBarry Smith     }
615ccb205f8SBarry Smith   }
616ccb205f8SBarry Smith   PetscFunctionReturn(0);
617ccb205f8SBarry Smith }
618ccb205f8SBarry Smith 
619ccb205f8SBarry Smith #undef __FUNCT__
620ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat"
6214aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
622ccb205f8SBarry Smith {
623ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
624ccb205f8SBarry Smith   Mat_SeqAIJ     *c;
625ccb205f8SBarry Smith   PetscErrorCode ierr;
626ccb205f8SBarry Smith   PetscInt       i,k,first,step,lensi,nrows,ncols;
627d2a212eaSBarry Smith   PetscInt       *j_new,*i_new,*aj = a->j,*ailen = a->ilen;
6281620fd73SBarry Smith   PetscScalar    *a_new;
629ccb205f8SBarry Smith   Mat            C,*aa = a->a;
630ace3abfcSBarry Smith   PetscBool      stride,equal;
631ccb205f8SBarry Smith 
632ccb205f8SBarry Smith   PetscFunctionBegin;
633ccb205f8SBarry Smith   ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr);
634e32f2f54SBarry Smith   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices");
635251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);CHKERRQ(ierr);
636e32f2f54SBarry Smith   if (!stride) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices");
637ccb205f8SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
638e32f2f54SBarry Smith   if (step != A->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block");
639ccb205f8SBarry Smith 
640ccb205f8SBarry Smith   ierr  = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
641ccb205f8SBarry Smith   ncols = nrows;
642ccb205f8SBarry Smith 
643ccb205f8SBarry Smith   /* create submatrix */
644ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
645ccb205f8SBarry Smith     PetscInt n_cols,n_rows;
646ccb205f8SBarry Smith     C    = *B;
647ccb205f8SBarry Smith     ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr);
648e32f2f54SBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
649ccb205f8SBarry Smith     ierr = MatZeroEntries(C);CHKERRQ(ierr);
650ccb205f8SBarry Smith   } else {
6517adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
652ccb205f8SBarry Smith     ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
6536e63c7a1SBarry Smith     if (A->symmetric) {
6546e63c7a1SBarry Smith       ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
6556e63c7a1SBarry Smith     } else {
656ccb205f8SBarry Smith       ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr);
6576e63c7a1SBarry Smith     }
6586e63c7a1SBarry Smith     ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr);
6596e63c7a1SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr);
660ccb205f8SBarry Smith   }
661ccb205f8SBarry Smith   c = (Mat_SeqAIJ*)C->data;
662ccb205f8SBarry Smith 
663ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
664ccb205f8SBarry Smith   a_new = c->a;
665ccb205f8SBarry Smith   j_new = c->j;
666ccb205f8SBarry Smith   i_new = c->i;
667ccb205f8SBarry Smith 
668ccb205f8SBarry Smith   for (i=0; i<nrows; i++) {
669ccb205f8SBarry Smith     lensi = ailen[i];
670ccb205f8SBarry Smith     for (k=0; k<lensi; k++) {
671ccb205f8SBarry Smith       *j_new++ = *aj++;
6721620fd73SBarry Smith       ierr     = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr);
673ccb205f8SBarry Smith     }
674ccb205f8SBarry Smith     i_new[i+1] = i_new[i] + lensi;
675ccb205f8SBarry Smith     c->ilen[i] = lensi;
676ccb205f8SBarry Smith   }
677ccb205f8SBarry Smith 
678ccb205f8SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
679ccb205f8SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
680ccb205f8SBarry Smith   *B   = C;
681ccb205f8SBarry Smith   PetscFunctionReturn(0);
682ccb205f8SBarry Smith }
683ccb205f8SBarry Smith 
684421e10b8SBarry Smith #undef __FUNCT__
685421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat"
686421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
687421e10b8SBarry Smith {
688421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
689421e10b8SBarry Smith   PetscErrorCode ierr;
690421e10b8SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
691ccb205f8SBarry Smith   PetscInt       m      = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
692421e10b8SBarry Smith   Mat            *aa    = a->a,*ap;
693421e10b8SBarry Smith 
694421e10b8SBarry Smith   PetscFunctionBegin;
695421e10b8SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
696421e10b8SBarry Smith 
697421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
698421e10b8SBarry Smith   for (i=1; i<m; i++) {
699421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
700421e10b8SBarry Smith     fshift += imax[i-1] - ailen[i-1];
701421e10b8SBarry Smith     rmax    = PetscMax(rmax,ailen[i]);
702421e10b8SBarry Smith     if (fshift) {
703421e10b8SBarry Smith       ip = aj + ai[i];
704421e10b8SBarry Smith       ap = aa + ai[i];
705421e10b8SBarry Smith       N  = ailen[i];
706421e10b8SBarry Smith       for (j=0; j<N; j++) {
707421e10b8SBarry Smith         ip[j-fshift] = ip[j];
708421e10b8SBarry Smith         ap[j-fshift] = ap[j];
709421e10b8SBarry Smith       }
710421e10b8SBarry Smith     }
711421e10b8SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
712421e10b8SBarry Smith   }
713421e10b8SBarry Smith   if (m) {
714421e10b8SBarry Smith     fshift += imax[m-1] - ailen[m-1];
715421e10b8SBarry Smith     ai[m]   = ai[m-1] + ailen[m-1];
716421e10b8SBarry Smith   }
717421e10b8SBarry Smith   /* reset ilen and imax for each row */
718421e10b8SBarry Smith   for (i=0; i<m; i++) {
719421e10b8SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
720421e10b8SBarry Smith   }
721421e10b8SBarry Smith   a->nz = ai[m];
722ccb205f8SBarry Smith   for (i=0; i<a->nz; i++) {
723ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG)
724e32f2f54SBarry Smith     if (!aa[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz);
725ccb205f8SBarry Smith #endif
726ccb205f8SBarry Smith     ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
727ccb205f8SBarry Smith     ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
728ccb205f8SBarry Smith   }
729ccb205f8SBarry Smith   CHKMEMQ;
730d0f46423SBarry Smith   ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap->n/A->cmap->bs,fshift,a->nz);CHKERRQ(ierr);
731421e10b8SBarry Smith   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr);
732421e10b8SBarry Smith   ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr);
733*2205254eSKarl Rupp 
7348e58a170SBarry Smith   A->info.mallocs    += a->reallocs;
735421e10b8SBarry Smith   a->reallocs         = 0;
736421e10b8SBarry Smith   A->info.nz_unneeded = (double)fshift;
737421e10b8SBarry Smith   a->rmax             = rmax;
738421e10b8SBarry Smith 
739421e10b8SBarry Smith   A->same_nonzero = PETSC_TRUE;
740ccb205f8SBarry Smith   ierr            = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr);
741421e10b8SBarry Smith   PetscFunctionReturn(0);
742421e10b8SBarry Smith }
743421e10b8SBarry Smith 
744290bbb0aSBarry Smith #undef __FUNCT__
745290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat"
746ace3abfcSBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg)
747290bbb0aSBarry Smith {
748290bbb0aSBarry Smith   PetscFunctionBegin;
7494e0d8c25SBarry Smith   if (opt == MAT_SYMMETRIC && flg) {
75041f059aeSBarry Smith     A->ops->sor  = MatSOR_BlockMat_Symmetric;
751290bbb0aSBarry Smith     A->ops->mult = MatMult_BlockMat_Symmetric;
752290bbb0aSBarry Smith   } else {
753290bbb0aSBarry Smith     PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]);
754290bbb0aSBarry Smith   }
755290bbb0aSBarry Smith   PetscFunctionReturn(0);
756290bbb0aSBarry Smith }
757290bbb0aSBarry Smith 
758290bbb0aSBarry Smith 
759be332245SKarl Rupp static struct _MatOps MatOps_Values = {
760be332245SKarl Rupp         MatSetValues_BlockMat,
761421e10b8SBarry Smith         0,
762421e10b8SBarry Smith         0,
763421e10b8SBarry Smith         MatMult_BlockMat,
764421e10b8SBarry Smith /*  4*/ MatMultAdd_BlockMat,
765421e10b8SBarry Smith         MatMultTranspose_BlockMat,
766421e10b8SBarry Smith         MatMultTransposeAdd_BlockMat,
767421e10b8SBarry Smith         0,
768421e10b8SBarry Smith         0,
769421e10b8SBarry Smith         0,
770421e10b8SBarry Smith /* 10*/ 0,
771421e10b8SBarry Smith         0,
772421e10b8SBarry Smith         0,
77341f059aeSBarry Smith         MatSOR_BlockMat,
774421e10b8SBarry Smith         0,
775421e10b8SBarry Smith /* 15*/ 0,
776421e10b8SBarry Smith         0,
777421e10b8SBarry Smith         0,
778421e10b8SBarry Smith         0,
779421e10b8SBarry Smith         0,
780421e10b8SBarry Smith /* 20*/ 0,
781421e10b8SBarry Smith         MatAssemblyEnd_BlockMat,
782290bbb0aSBarry Smith         MatSetOption_BlockMat,
783421e10b8SBarry Smith         0,
784d519adbfSMatthew Knepley /* 24*/ 0,
785421e10b8SBarry Smith         0,
786421e10b8SBarry Smith         0,
787421e10b8SBarry Smith         0,
788421e10b8SBarry Smith         0,
789d519adbfSMatthew Knepley /* 29*/ 0,
790421e10b8SBarry Smith         0,
791421e10b8SBarry Smith         0,
792421e10b8SBarry Smith         0,
793421e10b8SBarry Smith         0,
794d519adbfSMatthew Knepley /* 34*/ 0,
795421e10b8SBarry Smith         0,
796421e10b8SBarry Smith         0,
797421e10b8SBarry Smith         0,
798421e10b8SBarry Smith         0,
799d519adbfSMatthew Knepley /* 39*/ 0,
800421e10b8SBarry Smith         0,
801421e10b8SBarry Smith         0,
802421e10b8SBarry Smith         0,
803421e10b8SBarry Smith         0,
804d519adbfSMatthew Knepley /* 44*/ 0,
805421e10b8SBarry Smith         0,
806421e10b8SBarry Smith         0,
807421e10b8SBarry Smith         0,
808421e10b8SBarry Smith         0,
809d519adbfSMatthew Knepley /* 49*/ 0,
810421e10b8SBarry Smith         0,
811421e10b8SBarry Smith         0,
812421e10b8SBarry Smith         0,
813421e10b8SBarry Smith         0,
814d519adbfSMatthew Knepley /* 54*/ 0,
815421e10b8SBarry Smith         0,
816421e10b8SBarry Smith         0,
817421e10b8SBarry Smith         0,
818421e10b8SBarry Smith         0,
819d519adbfSMatthew Knepley /* 59*/ MatGetSubMatrix_BlockMat,
820421e10b8SBarry Smith         MatDestroy_BlockMat,
821ccb205f8SBarry Smith         MatView_BlockMat,
822421e10b8SBarry Smith         0,
823421e10b8SBarry Smith         0,
824d519adbfSMatthew Knepley /* 64*/ 0,
825421e10b8SBarry Smith         0,
826421e10b8SBarry Smith         0,
827421e10b8SBarry Smith         0,
828421e10b8SBarry Smith         0,
829d519adbfSMatthew Knepley /* 69*/ 0,
830421e10b8SBarry Smith         0,
831421e10b8SBarry Smith         0,
832421e10b8SBarry Smith         0,
833421e10b8SBarry Smith         0,
834d519adbfSMatthew Knepley /* 74*/ 0,
835421e10b8SBarry Smith         0,
836421e10b8SBarry Smith         0,
837421e10b8SBarry Smith         0,
838421e10b8SBarry Smith         0,
839d519adbfSMatthew Knepley /* 79*/ 0,
840421e10b8SBarry Smith         0,
841421e10b8SBarry Smith         0,
842421e10b8SBarry Smith         0,
8435bba2384SShri Abhyankar         MatLoad_BlockMat,
844d519adbfSMatthew Knepley /* 84*/ 0,
845421e10b8SBarry Smith         0,
846421e10b8SBarry Smith         0,
847421e10b8SBarry Smith         0,
848421e10b8SBarry Smith         0,
849d519adbfSMatthew Knepley /* 89*/ 0,
850421e10b8SBarry Smith         0,
851421e10b8SBarry Smith         0,
852421e10b8SBarry Smith         0,
853421e10b8SBarry Smith         0,
854d519adbfSMatthew Knepley /* 94*/ 0,
855421e10b8SBarry Smith         0,
856421e10b8SBarry Smith         0,
857a5973382SShri Abhyankar         0,
858a5973382SShri Abhyankar         0,
859a5973382SShri Abhyankar /* 99*/ 0,
860a5973382SShri Abhyankar         0,
861a5973382SShri Abhyankar         0,
862a5973382SShri Abhyankar         0,
863a5973382SShri Abhyankar         0,
864a5973382SShri Abhyankar /*104*/ 0,
865a5973382SShri Abhyankar         0,
866a5973382SShri Abhyankar         0,
867a5973382SShri Abhyankar         0,
868a5973382SShri Abhyankar         0,
869a5973382SShri Abhyankar /*109*/ 0,
870a5973382SShri Abhyankar         0,
871a5973382SShri Abhyankar         0,
872a5973382SShri Abhyankar         0,
873a5973382SShri Abhyankar         0,
874a5973382SShri Abhyankar /*114*/ 0,
875a5973382SShri Abhyankar         0,
876a5973382SShri Abhyankar         0,
877a5973382SShri Abhyankar         0,
878a5973382SShri Abhyankar         0,
879a5973382SShri Abhyankar /*119*/ 0,
880a5973382SShri Abhyankar         0,
881a5973382SShri Abhyankar         0,
8825bba2384SShri Abhyankar         0
883a5973382SShri Abhyankar };
884421e10b8SBarry Smith 
8858cdf2d9bSBarry Smith #undef __FUNCT__
8868cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation"
8878cdf2d9bSBarry Smith /*@C
8888cdf2d9bSBarry Smith    MatBlockMatSetPreallocation - For good matrix assembly performance
8898cdf2d9bSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
8908cdf2d9bSBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
8918cdf2d9bSBarry Smith    during matrix assembly can be increased by more than a factor of 50.
8928cdf2d9bSBarry Smith 
8938cdf2d9bSBarry Smith    Collective on MPI_Comm
8948cdf2d9bSBarry Smith 
8958cdf2d9bSBarry Smith    Input Parameters:
8968cdf2d9bSBarry Smith +  B - The matrix
8978cdf2d9bSBarry Smith .  bs - size of each block in matrix
8988cdf2d9bSBarry Smith .  nz - number of nonzeros per block row (same for all rows)
8998cdf2d9bSBarry Smith -  nnz - array containing the number of nonzeros in the various block rows
9008cdf2d9bSBarry Smith          (possibly different for each row) or PETSC_NULL
9018cdf2d9bSBarry Smith 
9028cdf2d9bSBarry Smith    Notes:
9038cdf2d9bSBarry Smith      If nnz is given then nz is ignored
9048cdf2d9bSBarry Smith 
9058cdf2d9bSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
9068cdf2d9bSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
9078cdf2d9bSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
9088cdf2d9bSBarry Smith    will get TERRIBLE performance, see the users' manual chapter on matrices.
9098cdf2d9bSBarry Smith 
9108cdf2d9bSBarry Smith    Level: intermediate
9118cdf2d9bSBarry Smith 
9128cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues()
9138cdf2d9bSBarry Smith 
9148cdf2d9bSBarry Smith @*/
9157087cfbeSBarry Smith PetscErrorCode  MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
9168cdf2d9bSBarry Smith {
9174ac538c5SBarry Smith   PetscErrorCode ierr;
9188cdf2d9bSBarry Smith 
9198cdf2d9bSBarry Smith   PetscFunctionBegin;
9204ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr);
9218cdf2d9bSBarry Smith   PetscFunctionReturn(0);
9228cdf2d9bSBarry Smith }
9238cdf2d9bSBarry Smith 
9248cdf2d9bSBarry Smith EXTERN_C_BEGIN
9258cdf2d9bSBarry Smith #undef __FUNCT__
9268cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat"
9277087cfbeSBarry Smith PetscErrorCode  MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz)
9288cdf2d9bSBarry Smith {
9298cdf2d9bSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
9308cdf2d9bSBarry Smith   PetscErrorCode ierr;
9318cdf2d9bSBarry Smith   PetscInt       i;
9328cdf2d9bSBarry Smith 
9338cdf2d9bSBarry Smith   PetscFunctionBegin;
934e02043d6SBarry Smith   ierr = PetscLayoutSetBlockSize(A->rmap,bs);CHKERRQ(ierr);
935e02043d6SBarry Smith   ierr = PetscLayoutSetBlockSize(A->cmap,bs);CHKERRQ(ierr);
93634ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
93734ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
938e02043d6SBarry Smith   ierr = PetscLayoutGetBlockSize(A->rmap,&bs);CHKERRQ(ierr);
93934ef9618SShri Abhyankar 
9408cdf2d9bSBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
941e32f2f54SBarry Smith   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
9428cdf2d9bSBarry Smith   if (nnz) {
943d0f46423SBarry Smith     for (i=0; i<A->rmap->n/bs; i++) {
944e32f2f54SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
945e32f2f54SBarry Smith       if (nnz[i] > A->cmap->n/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],A->cmap->n/bs);
9468cdf2d9bSBarry Smith     }
9478cdf2d9bSBarry Smith   }
948d0f46423SBarry Smith   bmat->mbs = A->rmap->n/bs;
9498cdf2d9bSBarry Smith 
950778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr);
951778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,PETSC_NULL,&bmat->middle);CHKERRQ(ierr);
9528cdf2d9bSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr);
9538cdf2d9bSBarry Smith 
9542ee49352SLisandro Dalcin   if (!bmat->imax) {
955d0f46423SBarry Smith     ierr = PetscMalloc2(A->rmap->n,PetscInt,&bmat->imax,A->rmap->n,PetscInt,&bmat->ilen);CHKERRQ(ierr);
956d0f46423SBarry Smith     ierr = PetscLogObjectMemory(A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
9572ee49352SLisandro Dalcin   }
9588cdf2d9bSBarry Smith   if (nnz) {
9598cdf2d9bSBarry Smith     nz = 0;
960d0f46423SBarry Smith     for (i=0; i<A->rmap->n/A->rmap->bs; i++) {
9618cdf2d9bSBarry Smith       bmat->imax[i] = nnz[i];
9628cdf2d9bSBarry Smith       nz           += nnz[i];
9638cdf2d9bSBarry Smith     }
964f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation");
9658cdf2d9bSBarry Smith 
9668cdf2d9bSBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
967*2205254eSKarl Rupp   for (i=0; i<bmat->mbs; i++) bmat->ilen[i] = 0;
9688cdf2d9bSBarry Smith 
9698cdf2d9bSBarry Smith   /* allocate the matrix space */
9702ee49352SLisandro Dalcin   ierr       = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
971d0f46423SBarry Smith   ierr       = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap->n+1,PetscInt,&bmat->i);CHKERRQ(ierr);
972d0f46423SBarry Smith   ierr       = PetscLogObjectMemory(A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
9738cdf2d9bSBarry Smith   bmat->i[0] = 0;
9748cdf2d9bSBarry Smith   for (i=1; i<bmat->mbs+1; i++) {
9758cdf2d9bSBarry Smith     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
9768cdf2d9bSBarry Smith   }
9778cdf2d9bSBarry Smith   bmat->singlemalloc = PETSC_TRUE;
9788cdf2d9bSBarry Smith   bmat->free_a       = PETSC_TRUE;
9798cdf2d9bSBarry Smith   bmat->free_ij      = PETSC_TRUE;
9808cdf2d9bSBarry Smith 
9818cdf2d9bSBarry Smith   bmat->nz            = 0;
9828cdf2d9bSBarry Smith   bmat->maxnz         = nz;
9838cdf2d9bSBarry Smith   A->info.nz_unneeded = (double)bmat->maxnz;
9847827cd58SJed Brown   ierr                = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
9858cdf2d9bSBarry Smith   PetscFunctionReturn(0);
9868cdf2d9bSBarry Smith }
9878cdf2d9bSBarry Smith EXTERN_C_END
9888cdf2d9bSBarry Smith 
989421e10b8SBarry Smith /*MC
990421e10b8SBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
991421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
992421e10b8SBarry Smith 
993421e10b8SBarry Smith   Level: advanced
994421e10b8SBarry Smith 
995421e10b8SBarry Smith .seealso: MatCreateBlockMat()
996421e10b8SBarry Smith 
997421e10b8SBarry Smith M*/
998421e10b8SBarry Smith 
999421e10b8SBarry Smith EXTERN_C_BEGIN
1000421e10b8SBarry Smith #undef __FUNCT__
1001421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat"
10027087cfbeSBarry Smith PetscErrorCode  MatCreate_BlockMat(Mat A)
1003421e10b8SBarry Smith {
1004421e10b8SBarry Smith   Mat_BlockMat   *b;
1005421e10b8SBarry Smith   PetscErrorCode ierr;
1006421e10b8SBarry Smith 
1007421e10b8SBarry Smith   PetscFunctionBegin;
100838f2d2fdSLisandro Dalcin   ierr    = PetscNewLog(A,Mat_BlockMat,&b);CHKERRQ(ierr);
1009421e10b8SBarry Smith   A->data = (void*)b;
101038f2d2fdSLisandro Dalcin   ierr    = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1011421e10b8SBarry Smith 
1012421e10b8SBarry Smith   A->assembled    = PETSC_TRUE;
1013421e10b8SBarry Smith   A->preallocated = PETSC_FALSE;
1014421e10b8SBarry Smith   ierr            = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
1015290bbb0aSBarry Smith 
10168cdf2d9bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatBlockMatSetPreallocation_C",
10178cdf2d9bSBarry Smith                                            "MatBlockMatSetPreallocation_BlockMat",
10188cdf2d9bSBarry Smith                                            MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr);
1019421e10b8SBarry Smith   PetscFunctionReturn(0);
1020421e10b8SBarry Smith }
1021421e10b8SBarry Smith EXTERN_C_END
1022421e10b8SBarry Smith 
1023421e10b8SBarry Smith #undef __FUNCT__
1024421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat"
1025421e10b8SBarry Smith /*@C
1026421e10b8SBarry Smith    MatCreateBlockMat - Creates a new matrix based sparse Mat storage
1027421e10b8SBarry Smith 
1028421e10b8SBarry Smith   Collective on MPI_Comm
1029421e10b8SBarry Smith 
1030421e10b8SBarry Smith    Input Parameters:
1031421e10b8SBarry Smith +  comm - MPI communicator
1032421e10b8SBarry Smith .  m - number of rows
1033421e10b8SBarry Smith .  n  - number of columns
10348cdf2d9bSBarry Smith .  bs - size of each submatrix
10358cdf2d9bSBarry Smith .  nz  - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known)
10368cdf2d9bSBarry Smith -  nnz - expected number of nonzers per block row if known (use PETSC_NULL otherwise)
1037421e10b8SBarry Smith 
1038421e10b8SBarry Smith 
1039421e10b8SBarry Smith    Output Parameter:
1040421e10b8SBarry Smith .  A - the matrix
1041421e10b8SBarry Smith 
1042421e10b8SBarry Smith    Level: intermediate
1043421e10b8SBarry Smith 
1044421e10b8SBarry Smith    PETSc requires that matrices and vectors being used for certain
1045421e10b8SBarry Smith    operations are partitioned accordingly.  For example, when
1046421e10b8SBarry Smith    creating a bmat matrix, A, that supports parallel matrix-vector
1047421e10b8SBarry Smith    products using MatMult(A,x,y) the user should set the number
1048421e10b8SBarry Smith    of local matrix rows to be the number of local elements of the
1049421e10b8SBarry Smith    corresponding result vector, y. Note that this is information is
1050421e10b8SBarry Smith    required for use of the matrix interface routines, even though
1051421e10b8SBarry Smith    the bmat matrix may not actually be physically partitioned.
1052421e10b8SBarry Smith    For example,
1053421e10b8SBarry Smith 
1054421e10b8SBarry Smith .keywords: matrix, bmat, create
1055421e10b8SBarry Smith 
1056421e10b8SBarry Smith .seealso: MATBLOCKMAT
1057421e10b8SBarry Smith @*/
10587087cfbeSBarry Smith PetscErrorCode  MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A)
1059421e10b8SBarry Smith {
1060421e10b8SBarry Smith   PetscErrorCode ierr;
1061421e10b8SBarry Smith 
1062421e10b8SBarry Smith   PetscFunctionBegin;
1063421e10b8SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1064421e10b8SBarry Smith   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1065421e10b8SBarry Smith   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
10668cdf2d9bSBarry Smith   ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1067421e10b8SBarry Smith   PetscFunctionReturn(0);
1068421e10b8SBarry Smith }
1069421e10b8SBarry Smith 
1070421e10b8SBarry Smith 
1071421e10b8SBarry Smith 
1072