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