xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision 351fe11bbcb6a8357b6c8fdaad7b5c40d80ea5ec)
1 #define PETSCMAT_DLL
2 
3 /*
4    This provides a matrix that consists of Mats
5 */
6 
7 #include "src/mat/matimpl.h"              /*I "petscmat.h" I*/
8 #include "src/mat/impls/baij/seq/baij.h"    /* use the common AIJ data-structure */
9 #include "petscksp.h"
10 
11 #define CHUNKSIZE   15
12 
13 typedef struct {
14   SEQAIJHEADER(Mat);
15   SEQBAIJHEADER;
16   Mat               *diags;
17 
18   Vec               left,right,middle,workb;   /* dummy vectors to perform local parts of product */
19 } Mat_BlockMat;
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "MatRelax_BlockMat_Symmetric"
23 PetscErrorCode MatRelax_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
24 {
25   Mat_BlockMat       *a = (Mat_BlockMat*)A->data;
26   PetscScalar        *x;
27   const Mat          *v = a->a;
28   const PetscScalar  *b;
29   PetscErrorCode     ierr;
30   PetscInt           n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs;
31   const PetscInt     *idx;
32   IS                 row,col;
33   MatFactorInfo      info;
34   Vec                left = a->left,right = a->right, middle = a->middle;
35   Mat                *diag;
36 
37   PetscFunctionBegin;
38   its = its*lits;
39   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
40   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
41   if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
42   if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift");
43   if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP))
44     SETERRQ(PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
45 
46   if (!a->diags) {
47     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
48     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
49     for (i=0; i<mbs; i++) {
50       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr);
51       ierr = MatCholeskyFactorSymbolic(a->a[a->diag[i]],row,&info,a->diags+i);CHKERRQ(ierr);
52       ierr = MatCholeskyFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr);
53       ierr = ISDestroy(row);CHKERRQ(ierr);
54       ierr = ISDestroy(col);CHKERRQ(ierr);
55     }
56     ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr);
57   }
58   diag    = a->diags;
59 
60   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
61   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
62   /* copy right hand side because it must be modified during iteration */
63   ierr = VecCopy(bb,a->workb);CHKERRQ(ierr);
64   ierr = VecGetArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr);
65 
66   /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */
67   while (its--) {
68     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
69 
70       for (i=0; i<mbs; i++) {
71         n    = a->i[i+1] - a->i[i] - 1;
72         idx  = a->j + a->i[i] + 1;
73         v    = a->a + a->i[i] + 1;
74 
75         ierr = VecSet(left,0.0);CHKERRQ(ierr);
76         for (j=0; j<n; j++) {
77           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
78           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
79           ierr = VecResetArray(right);CHKERRQ(ierr);
80         }
81         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
82         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
83         ierr = VecResetArray(right);CHKERRQ(ierr);
84 
85         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
86         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
87 
88         /* now adjust right hand side, see MatRelax_SeqSBAIJ */
89         for (j=0; j<n; j++) {
90           ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr);
91           ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr);
92           ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr);
93           ierr = VecResetArray(middle);CHKERRQ(ierr);
94         }
95         ierr = VecResetArray(right);CHKERRQ(ierr);
96 
97       }
98     }
99     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
100 
101       for (i=mbs-1; i>=0; i--) {
102         n    = a->i[i+1] - a->i[i] - 1;
103         idx  = a->j + a->i[i] + 1;
104         v    = a->a + a->i[i] + 1;
105 
106         ierr = VecSet(left,0.0);CHKERRQ(ierr);
107         for (j=0; j<n; j++) {
108           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
109           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
110           ierr = VecResetArray(right);CHKERRQ(ierr);
111         }
112         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
113         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
114         ierr = VecResetArray(right);CHKERRQ(ierr);
115 
116         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
117         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
118         ierr = VecResetArray(right);CHKERRQ(ierr);
119 
120       }
121     }
122   }
123   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
124   ierr = VecRestoreArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr);
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNCT__
129 #define __FUNCT__ "MatRelax_BlockMat"
130 PetscErrorCode MatRelax_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
131 {
132   Mat_BlockMat       *a = (Mat_BlockMat*)A->data;
133   PetscScalar        *x;
134   const Mat          *v = a->a;
135   const PetscScalar  *b;
136   PetscErrorCode     ierr;
137   PetscInt           n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs;
138   const PetscInt     *idx;
139   IS                 row,col;
140   MatFactorInfo      info;
141   Vec                left = a->left,right = a->right;
142   Mat                *diag;
143 
144   PetscFunctionBegin;
145   its = its*lits;
146   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
147   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
148   if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
149   if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift");
150 
151   if (!a->diags) {
152     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
153     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
154     for (i=0; i<mbs; i++) {
155       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr);
156       ierr = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr);
157       ierr = MatLUFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr);
158       ierr = ISDestroy(row);CHKERRQ(ierr);
159       ierr = ISDestroy(col);CHKERRQ(ierr);
160     }
161   }
162   diag = a->diags;
163 
164   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
165   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
166   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
167 
168   /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */
169   while (its--) {
170     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
171 
172       for (i=0; i<mbs; i++) {
173         n    = a->i[i+1] - a->i[i];
174         idx  = a->j + a->i[i];
175         v    = a->a + a->i[i];
176 
177         ierr = VecSet(left,0.0);CHKERRQ(ierr);
178         for (j=0; j<n; j++) {
179           if (idx[j] != i) {
180             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
181             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
182             ierr = VecResetArray(right);CHKERRQ(ierr);
183           }
184         }
185         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
186         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
187         ierr = VecResetArray(right);CHKERRQ(ierr);
188 
189         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
190         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
191         ierr = VecResetArray(right);CHKERRQ(ierr);
192       }
193     }
194     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
195 
196       for (i=mbs-1; i>=0; i--) {
197         n    = a->i[i+1] - a->i[i];
198         idx  = a->j + a->i[i];
199         v    = a->a + a->i[i];
200 
201         ierr = VecSet(left,0.0);CHKERRQ(ierr);
202         for (j=0; j<n; j++) {
203           if (idx[j] != i) {
204             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
205             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
206             ierr = VecResetArray(right);CHKERRQ(ierr);
207           }
208         }
209         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
210         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
211         ierr = VecResetArray(right);CHKERRQ(ierr);
212 
213         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
214         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
215         ierr = VecResetArray(right);CHKERRQ(ierr);
216 
217       }
218     }
219   }
220   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
221   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
222   PetscFunctionReturn(0);
223 }
224 
225 #undef __FUNCT__
226 #define __FUNCT__ "MatSetValues_BlockMat"
227 PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
228 {
229   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
230   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
231   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
232   PetscInt       *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol;
233   PetscErrorCode ierr;
234   PetscInt       ridx,cidx;
235   PetscTruth     roworiented=a->roworiented;
236   MatScalar      value;
237   Mat            *ap,*aa = a->a;
238 
239   PetscFunctionBegin;
240   for (k=0; k<m; k++) { /* loop over added rows */
241     row  = im[k];
242     brow = row/bs;
243     if (row < 0) continue;
244 #if defined(PETSC_USE_DEBUG)
245     if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1);
246 #endif
247     rp   = aj + ai[brow];
248     ap   = aa + ai[brow];
249     rmax = imax[brow];
250     nrow = ailen[brow];
251     low  = 0;
252     high = nrow;
253     for (l=0; l<n; l++) { /* loop over added columns */
254       if (in[l] < 0) continue;
255 #if defined(PETSC_USE_DEBUG)
256       if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
257 #endif
258       col = in[l]; bcol = col/bs;
259       if (A->symmetric && brow > bcol) continue;
260       ridx = row % bs; cidx = col % bs;
261       if (roworiented) {
262         value = v[l + k*n];
263       } else {
264         value = v[k + l*m];
265       }
266       /*      printf(" brow %d bcol %d\n",brow,bcol);*/
267       if (col <= lastcol) low = 0; else high = nrow;
268       lastcol = col;
269       while (high-low > 7) {
270         t = (low+high)/2;
271         if (rp[t] > bcol) high = t;
272         else              low  = t;
273       }
274       for (i=low; i<high; i++) {
275         if (rp[i] > bcol) break;
276         if (rp[i] == bcol) {
277 	  /*          printf("row %d col %d found i %d\n",brow,bcol,i);*/
278           goto noinsert1;
279         }
280       }
281       if (nonew == 1) goto noinsert1;
282       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
283       MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
284       N = nrow++ - 1; high++;
285       /* shift up all the later entries in this row */
286       for (ii=N; ii>=i; ii--) {
287 	/* printf("shiffting N %d i %d ii %d brow %d bcol %d\n",N,i,ii,brow,bcol);*/
288         rp[ii+1] = rp[ii];
289         ap[ii+1] = ap[ii];
290       }
291       if (N>=i) ap[i] = 0;
292       rp[i]           = bcol;
293       a->nz++;
294       noinsert1:;
295       if (!*(ap+i)) {
296 	/*printf("create matrix at i %d rw %d col %d\n",i,brow,bcol);*/
297         if (A->symmetric && brow == bcol) {
298           /* don't use SBAIJ since want to reorder in sparse factorization */
299           ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
300         } else {
301           ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
302         }
303       }
304       /*      printf("numerical value at i %d row %d col %d cidx %d ridx %d value %g\n",i,brow,bcol,cidx,ridx,value);*/
305       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
306       low = i;
307     }
308     /*    printf("nrow for row %d %d\n",nrow,brow);*/
309     ailen[brow] = nrow;
310   }
311   /*printf("nz %d\n",a->nz);*/
312   A->same_nonzero = PETSC_FALSE;
313   PetscFunctionReturn(0);
314 }
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "MatLoad_BlockMat"
318 PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A)
319 {
320   PetscErrorCode    ierr;
321   Mat               tmpA;
322   PetscInt          i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0;
323   const PetscInt    *cols;
324   const PetscScalar *values;
325   PetscTruth        flg,notdone;
326   Mat_SeqAIJ        *a;
327   Mat_BlockMat      *amat;
328 
329   PetscFunctionBegin;
330   ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr);
331 
332   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
333   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr);
334     ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr);
335     ierr = PetscOptionsName("-matload_symmetric","Store the matrix as symmetric","MatLoad",&flg);CHKERRQ(ierr);
336   ierr = PetscOptionsEnd();CHKERRQ(ierr);
337 
338   /* Determine number of nonzero blocks for each block row */
339   a    = (Mat_SeqAIJ*) tmpA->data;
340   mbs  = m/bs;
341   ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr);
342   ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr);
343 
344   for (i=0; i<mbs; i++) {
345     for (j=0; j<bs; j++) {
346       ii[j]         = a->j + a->i[i*bs + j];
347       ilens[j]      = a->i[i*bs + j + 1] - a->i[i*bs + j];
348       /*      printf("j %d length %d\n",j,ilens[j]);*/
349     }
350 
351     currentcol = -1;
352     notdone = PETSC_TRUE;
353     while (PETSC_TRUE) {
354       notdone = PETSC_FALSE;
355       nextcol = 1000000000;
356       for (j=0; j<bs; j++) {
357 	/*	printf("loop j %d length %d\n",j,ilens[j]); */
358         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) {
359           ii[j]++;
360           ilens[j]--;
361         }
362         if (ilens[j] > 0) {
363           notdone = PETSC_TRUE;
364           nextcol = PetscMin(nextcol,ii[j][0]/bs);
365         }
366       }
367       if (!notdone) break;
368       printf("i %d currentcol %d lens[i] %d\n",i,currentcol,lens[i]);
369       if (!flg || (nextcol >= i)) lens[i]++;
370       currentcol = nextcol;
371     }
372      printf("len of i %d %d\n",i,lens[i]);
373   }
374 
375   ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,0,lens,A);CHKERRQ(ierr);
376   if (flg) {
377     ierr = MatSetOption(*A,MAT_SYMMETRIC);CHKERRQ(ierr);
378   }
379   amat = (Mat_BlockMat*)(*A)->data;
380 
381   /* preallocate the submatrices */
382   ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr);
383   for (i=0; i<mbs; i++) { /* loops for block rows */
384     for (j=0; j<bs; j++) {
385       ii[j]         = a->j + a->i[i*bs + j];
386       ilens[j]      = a->i[i*bs + j + 1] - a->i[i*bs + j];
387       /* printf("j %d length %d\n",j,ilens[j]);*/
388     }
389 
390     currentcol = 1000000000;
391     for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */
392       if (ilens[j] > 0) {
393 	currentcol = PetscMin(currentcol,ii[j][0]/bs);
394       }
395     }
396 
397     notdone = PETSC_TRUE;
398     while (PETSC_TRUE) {  /* loops over blocks in block row */
399 
400       notdone = PETSC_FALSE;
401       nextcol = 1000000000;
402       ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr);
403       for (j=0; j<bs; j++) { /* loop over rows in block */
404 	/*printf("loop j %d length %d\n",j,ilens[j]); */
405 	/*printf("current col %d %d\n",currentcol,ii[j][0]/bs);*/
406         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */
407 	  /*printf("j %d in llens[j] %d\n",j,llens[j]);*/
408           ii[j]++;
409           ilens[j]--;
410           llens[j]++;
411         }
412         if (ilens[j] > 0) {
413           notdone = PETSC_TRUE;
414           nextcol = PetscMin(nextcol,ii[j][0]/bs);
415         }
416       }
417       printf("cnt %d llens %d %d %d %d\n",cnt,llens[0],llens[1],llens[2],llens[3]);
418       if (cnt >= amat->maxnz) SETERRQ1(PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt);
419       if (!flg || currentcol >= i) {
420         amat->j[cnt] = currentcol;
421         ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr);
422         /*printf("mat %p row %d col %d\n",amat->a[cnt-1],i,currentcol);*/
423       }
424 
425       if (!notdone) break;
426       currentcol = nextcol;
427     }
428     amat->ilen[i] = lens[i];
429   }
430   CHKMEMQ;
431   /*printf("total cnt %d\n",cnt);*/
432 
433   ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr);
434   ierr = PetscFree(llens);CHKERRQ(ierr);
435 
436   /* copy over the matrix, one row at a time */
437   for (i=0; i<m; i++) {
438     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
439     ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
440     ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
441   }
442   ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
443   ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
444   PetscFunctionReturn(0);
445 }
446 
447 #undef __FUNCT__
448 #define __FUNCT__ "MatView_BlockMat"
449 PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
450 {
451   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
452   PetscErrorCode    ierr;
453   const char        *name;
454   PetscViewerFormat format;
455 
456   PetscFunctionBegin;
457   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
458   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
459   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
460     ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr);
461   }
462   PetscFunctionReturn(0);
463 }
464 
465 #undef __FUNCT__
466 #define __FUNCT__ "MatDestroy_BlockMat"
467 PetscErrorCode MatDestroy_BlockMat(Mat mat)
468 {
469   PetscErrorCode ierr;
470   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
471   PetscInt       i;
472 
473   PetscFunctionBegin;
474   if (bmat->right) {
475     ierr = VecDestroy(bmat->right);CHKERRQ(ierr);
476   }
477   if (bmat->left) {
478     ierr = VecDestroy(bmat->left);CHKERRQ(ierr);
479   }
480   if (bmat->middle) {
481     ierr = VecDestroy(bmat->middle);CHKERRQ(ierr);
482   }
483   if (bmat->workb) {
484     ierr = VecDestroy(bmat->workb);CHKERRQ(ierr);
485   }
486   if (bmat->diags) {
487     for (i=0; i<mat->rmap.n/mat->rmap.bs; i++) {
488       if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);}
489     }
490   }
491   if (bmat->a) {
492     for (i=0; i<bmat->nz; i++) {
493       if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);}
494     }
495   }
496   ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
497   ierr = PetscFree(bmat);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 
501 #undef __FUNCT__
502 #define __FUNCT__ "MatMult_BlockMat"
503 PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
504 {
505   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
506   PetscErrorCode ierr;
507   PetscScalar    *xx,*yy;
508   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
509   Mat            *aa;
510 
511   PetscFunctionBegin;
512   CHKMEMQ;
513   /*
514      Standard CSR multiply except each entry is a Mat
515   */
516   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
517 
518   ierr = VecSet(y,0.0);CHKERRQ(ierr);
519   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
520   aj  = bmat->j;
521   aa  = bmat->a;
522   ii  = bmat->i;
523   for (i=0; i<m; i++) {
524     jrow = ii[i];
525     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
526     n    = ii[i+1] - jrow;
527     for (j=0; j<n; j++) {
528       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
529       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
530       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
531       jrow++;
532     }
533     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
534   }
535   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
536   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
537   CHKMEMQ;
538   PetscFunctionReturn(0);
539 }
540 
541 #undef __FUNCT__
542 #define __FUNCT__ "MatMult_BlockMat_Symmetric"
543 PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
544 {
545   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
546   PetscErrorCode ierr;
547   PetscScalar    *xx,*yy;
548   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
549   Mat            *aa;
550 
551   PetscFunctionBegin;
552   CHKMEMQ;
553   /*
554      Standard CSR multiply except each entry is a Mat
555   */
556   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
557 
558   ierr = VecSet(y,0.0);CHKERRQ(ierr);
559   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
560   aj  = bmat->j;
561   aa  = bmat->a;
562   ii  = bmat->i;
563   for (i=0; i<m; i++) {
564     jrow = ii[i];
565     n    = ii[i+1] - jrow;
566     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
567     ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr);
568     /* if we ALWAYS required a diagonal entry then could remove this if test */
569     if (aj[jrow] == i) {
570       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
571       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
572       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
573       jrow++;
574       n--;
575     }
576     for (j=0; j<n; j++) {
577       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);            /* upper triangular part */
578       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
579       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
580 
581       ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr);            /* lower triangular part */
582       ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr);
583       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
584       jrow++;
585     }
586     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
587     ierr = VecResetArray(bmat->middle);CHKERRQ(ierr);
588   }
589   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
590   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
591   CHKMEMQ;
592   PetscFunctionReturn(0);
593 }
594 
595 #undef __FUNCT__
596 #define __FUNCT__ "MatMultAdd_BlockMat"
597 PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
598 {
599   PetscFunctionBegin;
600   PetscFunctionReturn(0);
601 }
602 
603 #undef __FUNCT__
604 #define __FUNCT__ "MatMultTranspose_BlockMat"
605 PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
606 {
607   PetscFunctionBegin;
608   PetscFunctionReturn(0);
609 }
610 
611 #undef __FUNCT__
612 #define __FUNCT__ "MatMultTransposeAdd_BlockMat"
613 PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
614 {
615   PetscFunctionBegin;
616   PetscFunctionReturn(0);
617 }
618 
619 /*
620      Adds diagonal pointers to sparse matrix structure.
621 */
622 #undef __FUNCT__
623 #define __FUNCT__ "MatMarkDiagonal_BlockMat"
624 PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
625 {
626   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
627   PetscErrorCode ierr;
628   PetscInt       i,j,mbs = A->rmap.n/A->rmap.bs;
629 
630   PetscFunctionBegin;
631   if (!a->diag) {
632     ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
633   }
634   for (i=0; i<mbs; i++) {
635     a->diag[i] = a->i[i+1];
636     for (j=a->i[i]; j<a->i[i+1]; j++) {
637       if (a->j[j] == i) {
638         a->diag[i] = j;
639         break;
640       }
641     }
642   }
643   PetscFunctionReturn(0);
644 }
645 
646 #undef __FUNCT__
647 #define __FUNCT__ "MatGetSubMatrix_BlockMat"
648 PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
649 {
650   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
651   Mat_SeqAIJ     *c;
652   PetscErrorCode ierr;
653   PetscInt       i,k,first,step,lensi,nrows,ncols;
654   PetscInt       *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
655   PetscScalar    *a_new,value;
656   Mat            C,*aa = a->a;
657   PetscTruth     stride,equal;
658 
659   PetscFunctionBegin;
660   ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr);
661   if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices");
662   ierr = ISStride(iscol,&stride);CHKERRQ(ierr);
663   if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices");
664   ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
665   if (step != A->rmap.bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block");
666 
667   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
668   ncols = nrows;
669 
670   /* create submatrix */
671   if (scall == MAT_REUSE_MATRIX) {
672     PetscInt n_cols,n_rows;
673     C = *B;
674     ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr);
675     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
676     ierr = MatZeroEntries(C);CHKERRQ(ierr);
677   } else {
678     ierr = MatCreate(A->comm,&C);CHKERRQ(ierr);
679     ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
680     if (A->symmetric) {
681       ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
682     } else {
683       ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr);
684     }
685     ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr);
686     ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr);
687   }
688   c = (Mat_SeqAIJ*)C->data;
689 
690   /* loop over rows inserting into submatrix */
691   a_new    = c->a;
692   j_new    = c->j;
693   i_new    = c->i;
694 
695   for (i=0; i<nrows; i++) {
696     ii    = ai[i];
697     lensi = ailen[i];
698     for (k=0; k<lensi; k++) {
699       *j_new++ = *aj++;
700       ierr     = MatGetValue(*aa++,first,first,value);CHKERRQ(ierr);
701       *a_new++ = value;
702     }
703     i_new[i+1]  = i_new[i] + lensi;
704     c->ilen[i]  = lensi;
705   }
706 
707   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
708   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
709   *B = C;
710   PetscFunctionReturn(0);
711 }
712 
713 #undef __FUNCT__
714 #define __FUNCT__ "MatAssemblyEnd_BlockMat"
715 PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
716 {
717   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
718   PetscErrorCode ierr;
719   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
720   PetscInt       m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
721   Mat            *aa = a->a,*ap;
722 
723   PetscFunctionBegin;
724   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
725 
726   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
727   for (i=1; i<m; i++) {
728     /* move each row back by the amount of empty slots (fshift) before it*/
729     fshift += imax[i-1] - ailen[i-1];
730     rmax   = PetscMax(rmax,ailen[i]);
731     if (fshift) {
732       ip = aj + ai[i] ;
733       ap = aa + ai[i] ;
734       N  = ailen[i];
735       for (j=0; j<N; j++) {
736         ip[j-fshift] = ip[j];
737         ap[j-fshift] = ap[j];
738       }
739     }
740     ai[i] = ai[i-1] + ailen[i-1];
741   }
742   if (m) {
743     fshift += imax[m-1] - ailen[m-1];
744     ai[m]  = ai[m-1] + ailen[m-1];
745   }
746   /* reset ilen and imax for each row */
747   for (i=0; i<m; i++) {
748     ailen[i] = imax[i] = ai[i+1] - ai[i];
749   }
750   a->nz = ai[m];
751   for (i=0; i<a->nz; i++) {
752 #if defined(PETSC_USE_DEBUG)
753     if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz);
754 #endif
755     /*printf("mat assembly %p col %d\n",aa[i],aj[i]);*/
756     ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
757     ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
758   }
759   CHKMEMQ;
760   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);
761   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr);
762   ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr);
763   a->reallocs          = 0;
764   A->info.nz_unneeded  = (double)fshift;
765   a->rmax              = rmax;
766 
767   A->same_nonzero = PETSC_TRUE;
768   ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr);
769   PetscFunctionReturn(0);
770 }
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "MatSetOption_BlockMat"
774 PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt)
775 {
776   PetscFunctionBegin;
777   if (opt == MAT_SYMMETRIC) {
778     A->ops->relax = MatRelax_BlockMat_Symmetric;
779     A->ops->mult  = MatMult_BlockMat_Symmetric;
780   } else {
781     PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]);
782   }
783   PetscFunctionReturn(0);
784 }
785 
786 
787 static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
788        0,
789        0,
790        MatMult_BlockMat,
791 /* 4*/ MatMultAdd_BlockMat,
792        MatMultTranspose_BlockMat,
793        MatMultTransposeAdd_BlockMat,
794        0,
795        0,
796        0,
797 /*10*/ 0,
798        0,
799        0,
800        MatRelax_BlockMat,
801        0,
802 /*15*/ 0,
803        0,
804        0,
805        0,
806        0,
807 /*20*/ 0,
808        MatAssemblyEnd_BlockMat,
809        0,
810        MatSetOption_BlockMat,
811        0,
812 /*25*/ 0,
813        0,
814        0,
815        0,
816        0,
817 /*30*/ 0,
818        0,
819        0,
820        0,
821        0,
822 /*35*/ 0,
823        0,
824        0,
825        0,
826        0,
827 /*40*/ 0,
828        0,
829        0,
830        0,
831        0,
832 /*45*/ 0,
833        0,
834        0,
835        0,
836        0,
837 /*50*/ 0,
838        0,
839        0,
840        0,
841        0,
842 /*55*/ 0,
843        0,
844        0,
845        0,
846        0,
847 /*60*/ MatGetSubMatrix_BlockMat,
848        MatDestroy_BlockMat,
849        MatView_BlockMat,
850        0,
851        0,
852 /*65*/ 0,
853        0,
854        0,
855        0,
856        0,
857 /*70*/ 0,
858        0,
859        0,
860        0,
861        0,
862 /*75*/ 0,
863        0,
864        0,
865        0,
866        0,
867 /*80*/ 0,
868        0,
869        0,
870        0,
871        MatLoad_BlockMat,
872 /*85*/ 0,
873        0,
874        0,
875        0,
876        0,
877 /*90*/ 0,
878        0,
879        0,
880        0,
881        0,
882 /*95*/ 0,
883        0,
884        0,
885        0};
886 
887 #undef __FUNCT__
888 #define __FUNCT__ "MatBlockMatSetPreallocation"
889 /*@C
890    MatBlockMatSetPreallocation - For good matrix assembly performance
891    the user should preallocate the matrix storage by setting the parameter nz
892    (or the array nnz).  By setting these parameters accurately, performance
893    during matrix assembly can be increased by more than a factor of 50.
894 
895    Collective on MPI_Comm
896 
897    Input Parameters:
898 +  B - The matrix
899 .  bs - size of each block in matrix
900 .  nz - number of nonzeros per block row (same for all rows)
901 -  nnz - array containing the number of nonzeros in the various block rows
902          (possibly different for each row) or PETSC_NULL
903 
904    Notes:
905      If nnz is given then nz is ignored
906 
907    Specify the preallocated storage with either nz or nnz (not both).
908    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
909    allocation.  For large problems you MUST preallocate memory or you
910    will get TERRIBLE performance, see the users' manual chapter on matrices.
911 
912    Level: intermediate
913 
914 .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues()
915 
916 @*/
917 PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
918 {
919   PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]);
920 
921   PetscFunctionBegin;
922   ierr = PetscObjectQueryFunction((PetscObject)B,"MatBlockMatSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
923   if (f) {
924     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
925   }
926   PetscFunctionReturn(0);
927 }
928 
929 EXTERN_C_BEGIN
930 #undef __FUNCT__
931 #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat"
932 PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz)
933 {
934   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
935   PetscErrorCode ierr;
936   PetscInt       i;
937 
938   PetscFunctionBegin;
939   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Block size given %D must be great than zero",bs);
940   if (A->rmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap.n);
941   if (A->cmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap.n);
942   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
943   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
944   if (nnz) {
945     for (i=0; i<A->rmap.n/bs; i++) {
946       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
947       if (nnz[i] > A->cmap.n/bs) SETERRQ3(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);
948     }
949   }
950   A->rmap.bs = A->cmap.bs = bs;
951   bmat->mbs  = A->rmap.n/bs;
952 
953   /*printf("A->rmap.n%d %d\n",A->rmap.n, bs);*/
954   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr);
955   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->middle);CHKERRQ(ierr);
956   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr);
957 
958   ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr);
959   if (nnz) {
960     nz = 0;
961     for (i=0; i<A->rmap.n/A->rmap.bs; i++) {
962       bmat->imax[i] = nnz[i];
963       nz           += nnz[i];
964     }
965   } else {
966     SETERRQ(PETSC_ERR_SUP,"Currently requires block row by row preallocation");
967   }
968 
969   /*printf("nz in p;real %d\n",nz);*/
970 
971   /* bmat->ilen will count nonzeros in each row so far. */
972   for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;}
973 
974   /* allocate the matrix space */
975   ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr);
976   bmat->i[0] = 0;
977   for (i=1; i<bmat->mbs+1; i++) {
978     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
979   }
980   bmat->singlemalloc = PETSC_TRUE;
981   bmat->free_a       = PETSC_TRUE;
982   bmat->free_ij      = PETSC_TRUE;
983 
984   bmat->nz                = 0;
985   bmat->maxnz             = nz;
986   A->info.nz_unneeded  = (double)bmat->maxnz;
987 
988   PetscFunctionReturn(0);
989 }
990 EXTERN_C_END
991 
992 /*MC
993    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
994                  consisting of (usually) sparse blocks.
995 
996   Level: advanced
997 
998 .seealso: MatCreateBlockMat()
999 
1000 M*/
1001 
1002 EXTERN_C_BEGIN
1003 #undef __FUNCT__
1004 #define __FUNCT__ "MatCreate_BlockMat"
1005 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A)
1006 {
1007   Mat_BlockMat   *b;
1008   PetscErrorCode ierr;
1009 
1010   PetscFunctionBegin;
1011   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1012   ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr);
1013 
1014   A->data = (void*)b;
1015 
1016   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
1017   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
1018 
1019   A->assembled     = PETSC_TRUE;
1020   A->preallocated  = PETSC_FALSE;
1021   ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
1022 
1023   ierr = PetscOptionsBegin(A->comm,A->prefix,"Matrix Option","Mat");CHKERRQ(ierr);
1024   ierr = PetscOptionsEnd();
1025 
1026   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatBlockMatSetPreallocation_C",
1027                                      "MatBlockMatSetPreallocation_BlockMat",
1028                                       MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr);
1029 
1030   PetscFunctionReturn(0);
1031 }
1032 EXTERN_C_END
1033 
1034 #undef __FUNCT__
1035 #define __FUNCT__ "MatCreateBlockMat"
1036 /*@C
1037    MatCreateBlockMat - Creates a new matrix based sparse Mat storage
1038 
1039   Collective on MPI_Comm
1040 
1041    Input Parameters:
1042 +  comm - MPI communicator
1043 .  m - number of rows
1044 .  n  - number of columns
1045 .  bs - size of each submatrix
1046 .  nz  - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known)
1047 -  nnz - expected number of nonzers per block row if known (use PETSC_NULL otherwise)
1048 
1049 
1050    Output Parameter:
1051 .  A - the matrix
1052 
1053    Level: intermediate
1054 
1055    PETSc requires that matrices and vectors being used for certain
1056    operations are partitioned accordingly.  For example, when
1057    creating a bmat matrix, A, that supports parallel matrix-vector
1058    products using MatMult(A,x,y) the user should set the number
1059    of local matrix rows to be the number of local elements of the
1060    corresponding result vector, y. Note that this is information is
1061    required for use of the matrix interface routines, even though
1062    the bmat matrix may not actually be physically partitioned.
1063    For example,
1064 
1065 .keywords: matrix, bmat, create
1066 
1067 .seealso: MATBLOCKMAT
1068 @*/
1069 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A)
1070 {
1071   PetscErrorCode ierr;
1072 
1073   PetscFunctionBegin;
1074   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1075   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1076   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
1077   ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 
1082 
1083