xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision ef413494761eee08f9a8ba7148086983ab493390)
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;   /* 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;
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 = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr);
52       ierr = MatLUFactorNumeric(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   }
57   middle  = a->middle;
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,middle);CHKERRQ(ierr);
64   ierr = VecGetArray(middle,(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];
72         idx  = a->j + a->i[i];
73         v    = a->a + a->i[i];
74 
75         ierr = VecSet(left,0.0);CHKERRQ(ierr);
76         for (j=0; j<n; j++) {
77           if (idx[j] != i) {
78             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
79             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
80             ierr = VecResetArray(right);CHKERRQ(ierr);
81           }
82         }
83         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
84         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
85         ierr = VecResetArray(right);CHKERRQ(ierr);
86 
87         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
88         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
89         ierr = VecResetArray(right);CHKERRQ(ierr);
90       }
91     }
92     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
93 
94       for (i=mbs-1; i>=0; i--) {
95         n    = a->i[i+1] - a->i[i];
96         idx  = a->j + a->i[i];
97         v    = a->a + a->i[i];
98 
99         ierr = VecSet(left,0.0);CHKERRQ(ierr);
100         for (j=0; j<n; j++) {
101           if (idx[j] != i) {
102             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
103             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
104             ierr = VecResetArray(right);CHKERRQ(ierr);
105           }
106         }
107         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
108         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
109         ierr = VecResetArray(right);CHKERRQ(ierr);
110 
111         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
112         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
113         ierr = VecResetArray(right);CHKERRQ(ierr);
114 
115       }
116     }
117   }
118   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
119   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 #undef __FUNCT__
124 #define __FUNCT__ "MatRelax_BlockMat"
125 PetscErrorCode MatRelax_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
126 {
127   Mat_BlockMat       *a = (Mat_BlockMat*)A->data;
128   PetscScalar        *x;
129   const Mat          *v = a->a;
130   const PetscScalar  *b;
131   PetscErrorCode     ierr;
132   PetscInt           n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs;
133   const PetscInt     *idx;
134   IS                 row,col;
135   MatFactorInfo      info;
136   Vec                left = a->left,right = a->right;
137   Mat                *diag;
138 
139   PetscFunctionBegin;
140   its = its*lits;
141   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
142   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
143   if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
144   if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift");
145 
146   if (!a->diags) {
147     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
148     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
149     for (i=0; i<mbs; i++) {
150       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr);
151       ierr = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr);
152       ierr = MatLUFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr);
153       ierr = ISDestroy(row);CHKERRQ(ierr);
154       ierr = ISDestroy(col);CHKERRQ(ierr);
155     }
156   }
157   diag = a->diags;
158 
159   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
160   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
161   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
162 
163   /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */
164   while (its--) {
165     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
166 
167       for (i=0; i<mbs; i++) {
168         n    = a->i[i+1] - a->i[i];
169         idx  = a->j + a->i[i];
170         v    = a->a + a->i[i];
171 
172         ierr = VecSet(left,0.0);CHKERRQ(ierr);
173         for (j=0; j<n; j++) {
174           if (idx[j] != i) {
175             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
176             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
177             ierr = VecResetArray(right);CHKERRQ(ierr);
178           }
179         }
180         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
181         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
182         ierr = VecResetArray(right);CHKERRQ(ierr);
183 
184         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
185         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
186         ierr = VecResetArray(right);CHKERRQ(ierr);
187       }
188     }
189     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
190 
191       for (i=mbs-1; i>=0; i--) {
192         n    = a->i[i+1] - a->i[i];
193         idx  = a->j + a->i[i];
194         v    = a->a + a->i[i];
195 
196         ierr = VecSet(left,0.0);CHKERRQ(ierr);
197         for (j=0; j<n; j++) {
198           if (idx[j] != i) {
199             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
200             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
201             ierr = VecResetArray(right);CHKERRQ(ierr);
202           }
203         }
204         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
205         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
206         ierr = VecResetArray(right);CHKERRQ(ierr);
207 
208         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
209         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
210         ierr = VecResetArray(right);CHKERRQ(ierr);
211 
212       }
213     }
214   }
215   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
216   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
217   PetscFunctionReturn(0);
218 }
219 
220 #undef __FUNCT__
221 #define __FUNCT__ "MatSetValues_BlockMat"
222 PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
223 {
224   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
225   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
226   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
227   PetscInt       *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol;
228   PetscErrorCode ierr;
229   PetscInt       ridx,cidx;
230   PetscTruth     roworiented=a->roworiented;
231   MatScalar      value;
232   Mat            *ap,*aa = a->a;
233 
234   PetscFunctionBegin;
235   for (k=0; k<m; k++) { /* loop over added rows */
236     row  = im[k];
237     brow = row/bs;
238     if (row < 0) continue;
239 #if defined(PETSC_USE_DEBUG)
240     if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1);
241 #endif
242     rp   = aj + ai[brow];
243     ap   = aa + ai[brow];
244     rmax = imax[brow];
245     nrow = ailen[brow];
246     low  = 0;
247     high = nrow;
248     for (l=0; l<n; l++) { /* loop over added columns */
249       if (in[l] < 0) continue;
250 #if defined(PETSC_USE_DEBUG)
251       if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
252 #endif
253       col = in[l]; bcol = col/bs;
254       if (A->symmetric && row > col) continue;
255       ridx = row % bs; cidx = col % bs;
256       if (roworiented) {
257         value = v[l + k*n];
258       } else {
259         value = v[k + l*m];
260       }
261       if (col <= lastcol) low = 0; else high = nrow;
262       lastcol = col;
263       while (high-low > 7) {
264         t = (low+high)/2;
265         if (rp[t] > bcol) high = t;
266         else              low  = t;
267       }
268       for (i=low; i<high; i++) {
269         if (rp[i] > bcol) break;
270         if (rp[i] == bcol) {
271 	  /*          printf("row %d col %d found i %d\n",brow,bcol,i);*/
272           goto noinsert1;
273         }
274       }
275       if (nonew == 1) goto noinsert1;
276       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
277       MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
278       N = nrow++ - 1; high++;
279       /* shift up all the later entries in this row */
280       /*      printf("N %d i %d\n",N,i);*/
281       for (ii=N; ii>=i; ii--) {
282         rp[ii+1] = rp[ii];
283         ap[ii+1] = ap[ii];
284       }
285       if (N>=i) ap[i] = 0;
286       rp[i]           = bcol;
287       a->nz++;
288       noinsert1:;
289       if (!*(ap+i)) {
290 	/*        printf("create matrix at i %d rw %d col %d\n",i,brow,bcol);*/
291         if (A->symmetric && brow == bcol) {
292           ierr = MatCreateSeqSBAIJ(PETSC_COMM_SELF,1,bs,bs,0,0,ap+i);CHKERRQ(ierr);
293         } else {
294           ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
295         }
296       }
297       /*      printf("numerical value at i %d row %d col %d cidx %d ridx %d value %g\n",i,brow,bcol,cidx,ridx,value);*/
298       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
299       low = i;
300     }
301     /*    printf("nrow for row %d %d\n",nrow,brow);*/
302     ailen[brow] = nrow;
303   }
304   A->same_nonzero = PETSC_FALSE;
305   PetscFunctionReturn(0);
306 }
307 
308 #undef __FUNCT__
309 #define __FUNCT__ "MatLoad_BlockMat"
310 PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A)
311 {
312   PetscErrorCode    ierr;
313   Mat               tmpA;
314   PetscInt          i,m,n,bs = 1,ncols;
315   const PetscInt    *cols;
316   const PetscScalar *values;
317   PetscTruth        flg;
318 
319   PetscFunctionBegin;
320   ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr);
321 
322   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
323   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr);
324     ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr);
325   ierr = PetscOptionsEnd();CHKERRQ(ierr);
326   ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,A);CHKERRQ(ierr);
327   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr);
328     ierr = PetscOptionsName("-matload_symmetric","Store the matrix as symmetric","MatLoad",&flg);CHKERRQ(ierr);
329     if (flg) {
330       ierr = MatSetOption(*A,MAT_SYMMETRIC);CHKERRQ(ierr);
331     }
332   ierr = PetscOptionsEnd();CHKERRQ(ierr);
333   for (i=0; i<m; i++) {
334     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
335     ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
336     ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
337   }
338   ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
339   ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
340   PetscFunctionReturn(0);
341 }
342 
343 #undef __FUNCT__
344 #define __FUNCT__ "MatView_BlockMat"
345 PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
346 {
347   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
348   PetscErrorCode    ierr;
349   const char        *name;
350   PetscViewerFormat format;
351 
352   PetscFunctionBegin;
353   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
354   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
355   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
356     ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr);
357   }
358   PetscFunctionReturn(0);
359 }
360 
361 #undef __FUNCT__
362 #define __FUNCT__ "MatDestroy_BlockMat"
363 PetscErrorCode MatDestroy_BlockMat(Mat mat)
364 {
365   PetscErrorCode ierr;
366   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
367   PetscInt       i;
368 
369   PetscFunctionBegin;
370   if (bmat->right) {
371     ierr = VecDestroy(bmat->right);CHKERRQ(ierr);
372   }
373   if (bmat->left) {
374     ierr = VecDestroy(bmat->left);CHKERRQ(ierr);
375   }
376   if (bmat->middle) {
377     ierr = VecDestroy(bmat->middle);CHKERRQ(ierr);
378   }
379   if (bmat->diags) {
380     for (i=0; i<mat->rmap.n/mat->rmap.bs; i++) {
381       if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);}
382     }
383   }
384   if (bmat->a) {
385     for (i=0; i<bmat->nz; i++) {
386       if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);}
387     }
388   }
389   ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
390   ierr = PetscFree(bmat);CHKERRQ(ierr);
391   PetscFunctionReturn(0);
392 }
393 
394 #undef __FUNCT__
395 #define __FUNCT__ "MatMult_BlockMat"
396 PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
397 {
398   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
399   PetscErrorCode ierr;
400   PetscScalar    *xx,*yy;
401   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
402   Mat            *aa;
403 
404   PetscFunctionBegin;
405   CHKMEMQ;
406   /*
407      Standard CSR multiply except each entry is a Mat
408   */
409   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
410 
411   ierr = VecSet(y,0.0);CHKERRQ(ierr);
412   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
413   aj  = bmat->j;
414   aa  = bmat->a;
415   ii  = bmat->i;
416   for (i=0; i<m; i++) {
417     jrow = ii[i];
418     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
419     n    = ii[i+1] - jrow;
420     for (j=0; j<n; j++) {
421       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
422       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
423       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
424       jrow++;
425     }
426     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
427   }
428   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
429   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
430   CHKMEMQ;
431   PetscFunctionReturn(0);
432 }
433 
434 #undef __FUNCT__
435 #define __FUNCT__ "MatMult_BlockMat_Symmetric"
436 PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
437 {
438   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
439   PetscErrorCode ierr;
440   PetscScalar    *xx,*yy;
441   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
442   Mat            *aa;
443 
444   PetscFunctionBegin;
445   CHKMEMQ;
446   /*
447      Standard CSR multiply except each entry is a Mat
448   */
449   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
450 
451   ierr = VecSet(y,0.0);CHKERRQ(ierr);
452   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
453   aj  = bmat->j;
454   aa  = bmat->a;
455   ii  = bmat->i;
456   for (i=0; i<m; i++) {
457     jrow = ii[i];
458     n    = ii[i+1] - jrow;
459     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
460     ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr);
461     /* if we ALWAYS required a diagonal entry then could remove this if test */
462     if (aj[jrow] == i) {
463       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
464       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
465       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
466       jrow++;
467       n--;
468     }
469     for (j=0; j<n; j++) {
470       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);            /* upper triangular part */
471       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
472       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
473 
474       ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr);            /* lower triangular part */
475       ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr);
476       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
477       jrow++;
478     }
479     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
480     ierr = VecResetArray(bmat->middle);CHKERRQ(ierr);
481   }
482   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
483   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
484   CHKMEMQ;
485   PetscFunctionReturn(0);
486 }
487 
488 #undef __FUNCT__
489 #define __FUNCT__ "MatMultAdd_BlockMat"
490 PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
491 {
492   PetscFunctionBegin;
493   PetscFunctionReturn(0);
494 }
495 
496 #undef __FUNCT__
497 #define __FUNCT__ "MatMultTranspose_BlockMat"
498 PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
499 {
500   PetscFunctionBegin;
501   PetscFunctionReturn(0);
502 }
503 
504 #undef __FUNCT__
505 #define __FUNCT__ "MatMultTransposeAdd_BlockMat"
506 PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
507 {
508   PetscFunctionBegin;
509   PetscFunctionReturn(0);
510 }
511 
512 #undef __FUNCT__
513 #define __FUNCT__ "MatSetBlockSize_BlockMat"
514 PetscErrorCode MatSetBlockSize_BlockMat(Mat A,PetscInt bs)
515 {
516   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
517   PetscErrorCode ierr;
518   PetscInt       nz = 10,i;
519 
520   PetscFunctionBegin;
521   if (A->rmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap.n);
522   if (A->cmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap.n);
523   A->rmap.bs = A->cmap.bs = bs;
524 
525   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr);
526   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->middle);CHKERRQ(ierr);
527   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr);
528 
529 
530   ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr);
531   for (i=0; i<A->rmap.n; i++) bmat->imax[i] = nz;
532   nz = nz*A->rmap.n;
533 
534   bmat->mbs = A->rmap.n/A->rmap.bs;
535 
536   /* bmat->ilen will count nonzeros in each row so far. */
537   for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;}
538 
539   /* allocate the matrix space */
540   ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr);
541   bmat->i[0] = 0;
542   for (i=1; i<bmat->mbs+1; i++) {
543     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
544   }
545   bmat->singlemalloc = PETSC_TRUE;
546   bmat->free_a       = PETSC_TRUE;
547   bmat->free_ij      = PETSC_TRUE;
548 
549   bmat->nz                = 0;
550   bmat->maxnz             = nz;
551   A->info.nz_unneeded  = (double)bmat->maxnz;
552 
553   PetscFunctionReturn(0);
554 }
555 
556 /*
557      Adds diagonal pointers to sparse matrix structure.
558 */
559 #undef __FUNCT__
560 #define __FUNCT__ "MatMarkDiagonal_BlockMat"
561 PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
562 {
563   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
564   PetscErrorCode ierr;
565   PetscInt       i,j,mbs = A->rmap.n/A->rmap.bs;
566 
567   PetscFunctionBegin;
568   if (!a->diag) {
569     ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
570   }
571   for (i=0; i<mbs; i++) {
572     a->diag[i] = a->i[i+1];
573     for (j=a->i[i]; j<a->i[i+1]; j++) {
574       if (a->j[j] == i) {
575         a->diag[i] = j;
576         break;
577       }
578     }
579   }
580   PetscFunctionReturn(0);
581 }
582 
583 #undef __FUNCT__
584 #define __FUNCT__ "MatGetSubMatrix_BlockMat"
585 PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
586 {
587   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
588   Mat_SeqAIJ     *c;
589   PetscErrorCode ierr;
590   PetscInt       i,k,first,step,lensi,nrows,ncols;
591   PetscInt       *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
592   PetscScalar    *a_new,value;
593   Mat            C,*aa = a->a;
594   PetscTruth     stride,equal;
595 
596   PetscFunctionBegin;
597   ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr);
598   if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices");
599   ierr = ISStride(iscol,&stride);CHKERRQ(ierr);
600   if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices");
601   ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
602   if (step != A->rmap.bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block");
603 
604   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
605   ncols = nrows;
606 
607   /* create submatrix */
608   if (scall == MAT_REUSE_MATRIX) {
609     PetscInt n_cols,n_rows;
610     C = *B;
611     ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr);
612     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
613     ierr = MatZeroEntries(C);CHKERRQ(ierr);
614   } else {
615     ierr = MatCreate(A->comm,&C);CHKERRQ(ierr);
616     ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
617     ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr);
618     ierr = MatSeqAIJSetPreallocation_SeqAIJ(C,0,ailen);CHKERRQ(ierr);
619   }
620   c = (Mat_SeqAIJ*)C->data;
621 
622   /* loop over rows inserting into submatrix */
623   a_new    = c->a;
624   j_new    = c->j;
625   i_new    = c->i;
626 
627   for (i=0; i<nrows; i++) {
628     ii    = ai[i];
629     lensi = ailen[i];
630     for (k=0; k<lensi; k++) {
631       *j_new++ = *aj++;
632       ierr     = MatGetValue(*aa++,first,first,value);CHKERRQ(ierr);
633       *a_new++ = value;
634     }
635     i_new[i+1]  = i_new[i] + lensi;
636     c->ilen[i]  = lensi;
637   }
638 
639   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
640   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
641   *B = C;
642   PetscFunctionReturn(0);
643 }
644 
645 #undef __FUNCT__
646 #define __FUNCT__ "MatAssemblyEnd_BlockMat"
647 PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
648 {
649   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
650   PetscErrorCode ierr;
651   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
652   PetscInt       m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
653   Mat            *aa = a->a,*ap;
654 
655   PetscFunctionBegin;
656   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
657 
658   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
659   for (i=1; i<m; i++) {
660     /* move each row back by the amount of empty slots (fshift) before it*/
661     fshift += imax[i-1] - ailen[i-1];
662     rmax   = PetscMax(rmax,ailen[i]);
663     if (fshift) {
664       ip = aj + ai[i] ;
665       ap = aa + ai[i] ;
666       N  = ailen[i];
667       for (j=0; j<N; j++) {
668         ip[j-fshift] = ip[j];
669         ap[j-fshift] = ap[j];
670       }
671     }
672     ai[i] = ai[i-1] + ailen[i-1];
673   }
674   if (m) {
675     fshift += imax[m-1] - ailen[m-1];
676     ai[m]  = ai[m-1] + ailen[m-1];
677   }
678   /* reset ilen and imax for each row */
679   for (i=0; i<m; i++) {
680     ailen[i] = imax[i] = ai[i+1] - ai[i];
681   }
682   a->nz = ai[m];
683   for (i=0; i<a->nz; i++) {
684 #if defined(PETSC_USE_DEBUG)
685     if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz);
686 #endif
687     ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
688     ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
689   }
690   CHKMEMQ;
691   ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);CHKERRQ(ierr);
692   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr);
693   ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr);
694   a->reallocs          = 0;
695   A->info.nz_unneeded  = (double)fshift;
696   a->rmax              = rmax;
697 
698   A->same_nonzero = PETSC_TRUE;
699   ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr);
700   PetscFunctionReturn(0);
701 }
702 
703 #undef __FUNCT__
704 #define __FUNCT__ "MatSetOption_BlockMat"
705 PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt)
706 {
707   PetscFunctionBegin;
708   if (opt == MAT_SYMMETRIC) {
709     A->ops->relax = MatRelax_BlockMat_Symmetric;
710     A->ops->mult  = MatMult_BlockMat_Symmetric;
711   } else {
712     PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]);
713   }
714   PetscFunctionReturn(0);
715 }
716 
717 
718 static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
719        0,
720        0,
721        MatMult_BlockMat,
722 /* 4*/ MatMultAdd_BlockMat,
723        MatMultTranspose_BlockMat,
724        MatMultTransposeAdd_BlockMat,
725        0,
726        0,
727        0,
728 /*10*/ 0,
729        0,
730        0,
731        MatRelax_BlockMat,
732        0,
733 /*15*/ 0,
734        0,
735        0,
736        0,
737        0,
738 /*20*/ 0,
739        MatAssemblyEnd_BlockMat,
740        0,
741        MatSetOption_BlockMat,
742        0,
743 /*25*/ 0,
744        0,
745        0,
746        0,
747        0,
748 /*30*/ 0,
749        0,
750        0,
751        0,
752        0,
753 /*35*/ 0,
754        0,
755        0,
756        0,
757        0,
758 /*40*/ 0,
759        0,
760        0,
761        0,
762        0,
763 /*45*/ 0,
764        0,
765        0,
766        0,
767        0,
768 /*50*/ MatSetBlockSize_BlockMat,
769        0,
770        0,
771        0,
772        0,
773 /*55*/ 0,
774        0,
775        0,
776        0,
777        0,
778 /*60*/ MatGetSubMatrix_BlockMat,
779        MatDestroy_BlockMat,
780        MatView_BlockMat,
781        0,
782        0,
783 /*65*/ 0,
784        0,
785        0,
786        0,
787        0,
788 /*70*/ 0,
789        0,
790        0,
791        0,
792        0,
793 /*75*/ 0,
794        0,
795        0,
796        0,
797        0,
798 /*80*/ 0,
799        0,
800        0,
801        0,
802        MatLoad_BlockMat,
803 /*85*/ 0,
804        0,
805        0,
806        0,
807        0,
808 /*90*/ 0,
809        0,
810        0,
811        0,
812        0,
813 /*95*/ 0,
814        0,
815        0,
816        0};
817 
818 /*MC
819    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
820                  consisting of (usually) sparse blocks.
821 
822   Level: advanced
823 
824 .seealso: MatCreateBlockMat()
825 
826 M*/
827 
828 EXTERN_C_BEGIN
829 #undef __FUNCT__
830 #define __FUNCT__ "MatCreate_BlockMat"
831 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A)
832 {
833   Mat_BlockMat   *b;
834   PetscErrorCode ierr;
835 
836   PetscFunctionBegin;
837   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
838   ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr);
839 
840   A->data = (void*)b;
841 
842   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
843   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
844 
845   A->assembled     = PETSC_TRUE;
846   A->preallocated  = PETSC_FALSE;
847   ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
848 
849   ierr = PetscOptionsBegin(A->comm,A->prefix,"Matrix Option","Mat");CHKERRQ(ierr);
850   ierr = PetscOptionsEnd();
851 
852   PetscFunctionReturn(0);
853 }
854 EXTERN_C_END
855 
856 #undef __FUNCT__
857 #define __FUNCT__ "MatCreateBlockMat"
858 /*@C
859    MatCreateBlockMat - Creates a new matrix based sparse Mat storage
860 
861   Collective on MPI_Comm
862 
863    Input Parameters:
864 +  comm - MPI communicator
865 .  m - number of rows
866 .  n  - number of columns
867 -  bs - size of each submatrix
868 
869 
870    Output Parameter:
871 .  A - the matrix
872 
873    Level: intermediate
874 
875    PETSc requires that matrices and vectors being used for certain
876    operations are partitioned accordingly.  For example, when
877    creating a bmat matrix, A, that supports parallel matrix-vector
878    products using MatMult(A,x,y) the user should set the number
879    of local matrix rows to be the number of local elements of the
880    corresponding result vector, y. Note that this is information is
881    required for use of the matrix interface routines, even though
882    the bmat matrix may not actually be physically partitioned.
883    For example,
884 
885 .keywords: matrix, bmat, create
886 
887 .seealso: MATBLOCKMAT
888 @*/
889 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,Mat *A)
890 {
891   PetscErrorCode ierr;
892 
893   PetscFunctionBegin;
894   ierr = MatCreate(comm,A);CHKERRQ(ierr);
895   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
896   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
897   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
898   PetscFunctionReturn(0);
899 }
900 
901 
902 
903