xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision 290bbb0a1dcfb34dbf94efcfcc44171581b0efea)
1421e10b8SBarry Smith #define PETSCMAT_DLL
2421e10b8SBarry Smith 
3421e10b8SBarry Smith /*
4421e10b8SBarry Smith    This provides a matrix that consists of Mats
5421e10b8SBarry Smith */
6421e10b8SBarry Smith 
7421e10b8SBarry Smith #include "src/mat/matimpl.h"              /*I "petscmat.h" I*/
8421e10b8SBarry Smith #include "src/mat/impls/baij/seq/baij.h"    /* use the common AIJ data-structure */
9ccb205f8SBarry Smith #include "petscksp.h"
10421e10b8SBarry Smith 
11421e10b8SBarry Smith #define CHUNKSIZE   15
12421e10b8SBarry Smith 
13421e10b8SBarry Smith typedef struct {
14421e10b8SBarry Smith   SEQAIJHEADER(Mat);
15421e10b8SBarry Smith   SEQBAIJHEADER;
16ccb205f8SBarry Smith   Mat               *diags;
17421e10b8SBarry Smith 
18*290bbb0aSBarry Smith   Vec               left,right,middle;   /* dummy vectors to perform local parts of product */
19421e10b8SBarry Smith } Mat_BlockMat;
20421e10b8SBarry Smith 
21421e10b8SBarry Smith #undef __FUNCT__
22*290bbb0aSBarry Smith #define __FUNCT__ "MatRelax_BlockMat_Symmetric"
23*290bbb0aSBarry Smith PetscErrorCode MatRelax_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
24*290bbb0aSBarry Smith {
25*290bbb0aSBarry Smith   Mat_BlockMat       *a = (Mat_BlockMat*)A->data;
26*290bbb0aSBarry Smith   PetscScalar        *x;
27*290bbb0aSBarry Smith   const Mat          *v = a->a;
28*290bbb0aSBarry Smith   const PetscScalar  *b;
29*290bbb0aSBarry Smith   PetscErrorCode     ierr;
30*290bbb0aSBarry Smith   PetscInt           n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs;
31*290bbb0aSBarry Smith   const PetscInt     *idx;
32*290bbb0aSBarry Smith   IS                 row,col;
33*290bbb0aSBarry Smith   MatFactorInfo      info;
34*290bbb0aSBarry Smith   Vec                left = a->left,right = a->right, middle;
35*290bbb0aSBarry Smith   Mat                *diag;
36*290bbb0aSBarry Smith 
37*290bbb0aSBarry Smith   PetscFunctionBegin;
38*290bbb0aSBarry Smith   its = its*lits;
39*290bbb0aSBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
40*290bbb0aSBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
41*290bbb0aSBarry Smith   if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
42*290bbb0aSBarry Smith   if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift");
43*290bbb0aSBarry Smith   if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP))
44*290bbb0aSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
45*290bbb0aSBarry Smith 
46*290bbb0aSBarry Smith   if (!a->diags) {
47*290bbb0aSBarry Smith     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
48*290bbb0aSBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
49*290bbb0aSBarry Smith     for (i=0; i<mbs; i++) {
50*290bbb0aSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr);
51*290bbb0aSBarry Smith       ierr = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr);
52*290bbb0aSBarry Smith       ierr = MatLUFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr);
53*290bbb0aSBarry Smith       ierr = ISDestroy(row);CHKERRQ(ierr);
54*290bbb0aSBarry Smith       ierr = ISDestroy(col);CHKERRQ(ierr);
55*290bbb0aSBarry Smith     }
56*290bbb0aSBarry Smith   }
57*290bbb0aSBarry Smith   middle  = a->middle;
58*290bbb0aSBarry Smith   diag    = a->diags;
59*290bbb0aSBarry Smith 
60*290bbb0aSBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
61*290bbb0aSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
62*290bbb0aSBarry Smith   /* copy right hand side because it must be modified during iteration */
63*290bbb0aSBarry Smith   ierr = VecCopy(bb,middle);CHKERRQ(ierr);
64*290bbb0aSBarry Smith   ierr = VecGetArray(middle,(PetscScalar**)&b);CHKERRQ(ierr);
65*290bbb0aSBarry Smith 
66*290bbb0aSBarry Smith   /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */
67*290bbb0aSBarry Smith   while (its--) {
68*290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
69*290bbb0aSBarry Smith 
70*290bbb0aSBarry Smith       for (i=0; i<mbs; i++) {
71*290bbb0aSBarry Smith         n    = a->i[i+1] - a->i[i];
72*290bbb0aSBarry Smith         idx  = a->j + a->i[i];
73*290bbb0aSBarry Smith         v    = a->a + a->i[i];
74*290bbb0aSBarry Smith 
75*290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
76*290bbb0aSBarry Smith         for (j=0; j<n; j++) {
77*290bbb0aSBarry Smith           if (idx[j] != i) {
78*290bbb0aSBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
79*290bbb0aSBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
80*290bbb0aSBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
81*290bbb0aSBarry Smith           }
82*290bbb0aSBarry Smith         }
83*290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
84*290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
85*290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
86*290bbb0aSBarry Smith 
87*290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
88*290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
89*290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
90*290bbb0aSBarry Smith       }
91*290bbb0aSBarry Smith     }
92*290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
93*290bbb0aSBarry Smith 
94*290bbb0aSBarry Smith       for (i=mbs-1; i>=0; i--) {
95*290bbb0aSBarry Smith         n    = a->i[i+1] - a->i[i];
96*290bbb0aSBarry Smith         idx  = a->j + a->i[i];
97*290bbb0aSBarry Smith         v    = a->a + a->i[i];
98*290bbb0aSBarry Smith 
99*290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
100*290bbb0aSBarry Smith         for (j=0; j<n; j++) {
101*290bbb0aSBarry Smith           if (idx[j] != i) {
102*290bbb0aSBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
103*290bbb0aSBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
104*290bbb0aSBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
105*290bbb0aSBarry Smith           }
106*290bbb0aSBarry Smith         }
107*290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
108*290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
109*290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
110*290bbb0aSBarry Smith 
111*290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
112*290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
113*290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
114*290bbb0aSBarry Smith 
115*290bbb0aSBarry Smith       }
116*290bbb0aSBarry Smith     }
117*290bbb0aSBarry Smith   }
118*290bbb0aSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
119*290bbb0aSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
120*290bbb0aSBarry Smith   PetscFunctionReturn(0);
121*290bbb0aSBarry Smith }
122*290bbb0aSBarry Smith 
123*290bbb0aSBarry Smith #undef __FUNCT__
124ccb205f8SBarry Smith #define __FUNCT__ "MatRelax_BlockMat"
125ccb205f8SBarry Smith PetscErrorCode MatRelax_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
126ccb205f8SBarry Smith {
127ccb205f8SBarry Smith   Mat_BlockMat       *a = (Mat_BlockMat*)A->data;
128ccb205f8SBarry Smith   PetscScalar        *x;
129ccb205f8SBarry Smith   const Mat          *v = a->a;
130ccb205f8SBarry Smith   const PetscScalar  *b;
131ccb205f8SBarry Smith   PetscErrorCode     ierr;
132ccb205f8SBarry Smith   PetscInt           n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs;
133ccb205f8SBarry Smith   const PetscInt     *idx;
134ccb205f8SBarry Smith   IS                 row,col;
135ccb205f8SBarry Smith   MatFactorInfo      info;
136ccb205f8SBarry Smith   Vec                left = a->left,right = a->right;
137ccb205f8SBarry Smith   Mat                *diag;
138ccb205f8SBarry Smith 
139ccb205f8SBarry Smith   PetscFunctionBegin;
140ccb205f8SBarry Smith   its = its*lits;
141ccb205f8SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
142ccb205f8SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
143ccb205f8SBarry Smith   if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
144ccb205f8SBarry Smith   if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift");
145ccb205f8SBarry Smith 
146ccb205f8SBarry Smith   if (!a->diags) {
147ccb205f8SBarry Smith     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
148ccb205f8SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
149ccb205f8SBarry Smith     for (i=0; i<mbs; i++) {
150ccb205f8SBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr);
151ccb205f8SBarry Smith       ierr = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr);
152ccb205f8SBarry Smith       ierr = MatLUFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr);
153ccb205f8SBarry Smith       ierr = ISDestroy(row);CHKERRQ(ierr);
154ccb205f8SBarry Smith       ierr = ISDestroy(col);CHKERRQ(ierr);
155ccb205f8SBarry Smith     }
156ccb205f8SBarry Smith   }
157ccb205f8SBarry Smith   diag = a->diags;
158ccb205f8SBarry Smith 
159ccb205f8SBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
160ccb205f8SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
161ccb205f8SBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
162ccb205f8SBarry Smith 
163ccb205f8SBarry Smith   /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */
164ccb205f8SBarry Smith   while (its--) {
165ccb205f8SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
166ccb205f8SBarry Smith 
167ccb205f8SBarry Smith       for (i=0; i<mbs; i++) {
168ccb205f8SBarry Smith         n    = a->i[i+1] - a->i[i];
169ccb205f8SBarry Smith         idx  = a->j + a->i[i];
170ccb205f8SBarry Smith         v    = a->a + a->i[i];
171ccb205f8SBarry Smith 
172ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
173ccb205f8SBarry Smith         for (j=0; j<n; j++) {
174ccb205f8SBarry Smith           if (idx[j] != i) {
175ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
176ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
177ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
178ccb205f8SBarry Smith           }
179ccb205f8SBarry Smith         }
180ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
181ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
182ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
183ccb205f8SBarry Smith 
184ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
185ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
186ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
187ccb205f8SBarry Smith       }
188ccb205f8SBarry Smith     }
189ccb205f8SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
190ccb205f8SBarry Smith 
191ccb205f8SBarry Smith       for (i=mbs-1; i>=0; i--) {
192ccb205f8SBarry Smith         n    = a->i[i+1] - a->i[i];
193ccb205f8SBarry Smith         idx  = a->j + a->i[i];
194ccb205f8SBarry Smith         v    = a->a + a->i[i];
195ccb205f8SBarry Smith 
196ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
197ccb205f8SBarry Smith         for (j=0; j<n; j++) {
198ccb205f8SBarry Smith           if (idx[j] != i) {
199ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
200ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
201ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
202ccb205f8SBarry Smith           }
203ccb205f8SBarry Smith         }
204ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
205ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
206ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
207ccb205f8SBarry Smith 
208ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
209ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
210ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
211ccb205f8SBarry Smith 
212ccb205f8SBarry Smith       }
213ccb205f8SBarry Smith     }
214ccb205f8SBarry Smith   }
215ccb205f8SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
216ccb205f8SBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
217ccb205f8SBarry Smith   PetscFunctionReturn(0);
218ccb205f8SBarry Smith }
219ccb205f8SBarry Smith 
220ccb205f8SBarry Smith #undef __FUNCT__
221421e10b8SBarry Smith #define __FUNCT__ "MatSetValues_BlockMat"
222421e10b8SBarry Smith PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
223421e10b8SBarry Smith {
224421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
225421e10b8SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
226421e10b8SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
227421e10b8SBarry Smith   PetscInt       *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol;
228421e10b8SBarry Smith   PetscErrorCode ierr;
229421e10b8SBarry Smith   PetscInt       ridx,cidx;
230421e10b8SBarry Smith   PetscTruth     roworiented=a->roworiented;
231421e10b8SBarry Smith   MatScalar      value;
232421e10b8SBarry Smith   Mat            *ap,*aa = a->a;
233421e10b8SBarry Smith 
234421e10b8SBarry Smith   PetscFunctionBegin;
235421e10b8SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
236421e10b8SBarry Smith     row  = im[k];
237421e10b8SBarry Smith     brow = row/bs;
238421e10b8SBarry Smith     if (row < 0) continue;
239421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
240421e10b8SBarry Smith     if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1);
241421e10b8SBarry Smith #endif
242421e10b8SBarry Smith     rp   = aj + ai[brow];
243421e10b8SBarry Smith     ap   = aa + ai[brow];
244421e10b8SBarry Smith     rmax = imax[brow];
245421e10b8SBarry Smith     nrow = ailen[brow];
246421e10b8SBarry Smith     low  = 0;
247421e10b8SBarry Smith     high = nrow;
248421e10b8SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
249421e10b8SBarry Smith       if (in[l] < 0) continue;
250421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
251421e10b8SBarry Smith       if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
252421e10b8SBarry Smith #endif
253421e10b8SBarry Smith       col = in[l]; bcol = col/bs;
254*290bbb0aSBarry Smith       if (A->symmetric && brow > bcol) continue;
255421e10b8SBarry Smith       ridx = row % bs; cidx = col % bs;
256421e10b8SBarry Smith       if (roworiented) {
257421e10b8SBarry Smith         value = v[l + k*n];
258421e10b8SBarry Smith       } else {
259421e10b8SBarry Smith         value = v[k + l*m];
260421e10b8SBarry Smith       }
261421e10b8SBarry Smith       if (col <= lastcol) low = 0; else high = nrow;
262421e10b8SBarry Smith       lastcol = col;
263421e10b8SBarry Smith       while (high-low > 7) {
264421e10b8SBarry Smith         t = (low+high)/2;
265421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
266421e10b8SBarry Smith         else              low  = t;
267421e10b8SBarry Smith       }
268421e10b8SBarry Smith       for (i=low; i<high; i++) {
269421e10b8SBarry Smith         if (rp[i] > bcol) break;
270421e10b8SBarry Smith         if (rp[i] == bcol) {
271ccb205f8SBarry Smith 	  /*          printf("row %d col %d found i %d\n",brow,bcol,i);*/
272421e10b8SBarry Smith           goto noinsert1;
273421e10b8SBarry Smith         }
274421e10b8SBarry Smith       }
275421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
276421e10b8SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
277421e10b8SBarry Smith       MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
278421e10b8SBarry Smith       N = nrow++ - 1; high++;
279421e10b8SBarry Smith       /* shift up all the later entries in this row */
280ccb205f8SBarry Smith       /*      printf("N %d i %d\n",N,i);*/
281421e10b8SBarry Smith       for (ii=N; ii>=i; ii--) {
282421e10b8SBarry Smith         rp[ii+1] = rp[ii];
283421e10b8SBarry Smith         ap[ii+1] = ap[ii];
284421e10b8SBarry Smith       }
285421e10b8SBarry Smith       if (N>=i) ap[i] = 0;
286421e10b8SBarry Smith       rp[i]           = bcol;
287421e10b8SBarry Smith       a->nz++;
288421e10b8SBarry Smith       noinsert1:;
289421e10b8SBarry Smith       if (!*(ap+i)) {
290ccb205f8SBarry Smith 	/*        printf("create matrix at i %d rw %d col %d\n",i,brow,bcol);*/
291421e10b8SBarry Smith         ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
292421e10b8SBarry Smith       }
293ccb205f8SBarry Smith       /*      printf("numerical value at i %d row %d col %d cidx %d ridx %d value %g\n",i,brow,bcol,cidx,ridx,value);*/
294421e10b8SBarry Smith       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
295421e10b8SBarry Smith       low = i;
296421e10b8SBarry Smith     }
297ccb205f8SBarry Smith     /*    printf("nrow for row %d %d\n",nrow,brow);*/
298421e10b8SBarry Smith     ailen[brow] = nrow;
299421e10b8SBarry Smith   }
300421e10b8SBarry Smith   A->same_nonzero = PETSC_FALSE;
301421e10b8SBarry Smith   PetscFunctionReturn(0);
302421e10b8SBarry Smith }
303421e10b8SBarry Smith 
304421e10b8SBarry Smith #undef __FUNCT__
305421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat"
306421e10b8SBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A)
307421e10b8SBarry Smith {
308421e10b8SBarry Smith   PetscErrorCode    ierr;
309421e10b8SBarry Smith   Mat               tmpA;
310421e10b8SBarry Smith   PetscInt          i,m,n,bs = 1,ncols;
311421e10b8SBarry Smith   const PetscInt    *cols;
312421e10b8SBarry Smith   const PetscScalar *values;
313*290bbb0aSBarry Smith   PetscTruth        flg;
314421e10b8SBarry Smith 
315421e10b8SBarry Smith   PetscFunctionBegin;
316421e10b8SBarry Smith   ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr);
317421e10b8SBarry Smith 
318421e10b8SBarry Smith   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
319421e10b8SBarry Smith   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr);
320421e10b8SBarry Smith     ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr);
321421e10b8SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
322421e10b8SBarry Smith   ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,A);CHKERRQ(ierr);
323*290bbb0aSBarry Smith   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr);
324*290bbb0aSBarry Smith     ierr = PetscOptionsName("-matload_symmetric","Store the matrix as symmetric","MatLoad",&flg);CHKERRQ(ierr);
325*290bbb0aSBarry Smith     if (flg) {
326*290bbb0aSBarry Smith       ierr = MatSetOption(*A,MAT_SYMMETRIC);CHKERRQ(ierr);
327*290bbb0aSBarry Smith     }
328*290bbb0aSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
329421e10b8SBarry Smith   for (i=0; i<m; i++) {
330421e10b8SBarry Smith     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
331421e10b8SBarry Smith     ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
332ccb205f8SBarry Smith     ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
333421e10b8SBarry Smith   }
334421e10b8SBarry Smith   ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
335421e10b8SBarry Smith   ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
336421e10b8SBarry Smith   PetscFunctionReturn(0);
337421e10b8SBarry Smith }
338421e10b8SBarry Smith 
339ccb205f8SBarry Smith #undef __FUNCT__
340ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat"
341ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
342ccb205f8SBarry Smith {
343ccb205f8SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
344ccb205f8SBarry Smith   PetscErrorCode    ierr;
345ccb205f8SBarry Smith   const char        *name;
346ccb205f8SBarry Smith   PetscViewerFormat format;
347ccb205f8SBarry Smith 
348ccb205f8SBarry Smith   PetscFunctionBegin;
349ccb205f8SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
350ccb205f8SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
351ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
352ccb205f8SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr);
353ccb205f8SBarry Smith   }
354ccb205f8SBarry Smith   PetscFunctionReturn(0);
355ccb205f8SBarry Smith }
356421e10b8SBarry Smith 
357421e10b8SBarry Smith #undef __FUNCT__
358421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat"
359421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat)
360421e10b8SBarry Smith {
361421e10b8SBarry Smith   PetscErrorCode ierr;
362421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
363ccb205f8SBarry Smith   PetscInt       i;
364421e10b8SBarry Smith 
365421e10b8SBarry Smith   PetscFunctionBegin;
366421e10b8SBarry Smith   if (bmat->right) {
367421e10b8SBarry Smith     ierr = VecDestroy(bmat->right);CHKERRQ(ierr);
368421e10b8SBarry Smith   }
369421e10b8SBarry Smith   if (bmat->left) {
370421e10b8SBarry Smith     ierr = VecDestroy(bmat->left);CHKERRQ(ierr);
371421e10b8SBarry Smith   }
372*290bbb0aSBarry Smith   if (bmat->middle) {
373*290bbb0aSBarry Smith     ierr = VecDestroy(bmat->middle);CHKERRQ(ierr);
374*290bbb0aSBarry Smith   }
375ccb205f8SBarry Smith   if (bmat->diags) {
376ccb205f8SBarry Smith     for (i=0; i<mat->rmap.n/mat->rmap.bs; i++) {
377ccb205f8SBarry Smith       if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);}
378ccb205f8SBarry Smith     }
379ccb205f8SBarry Smith   }
380ccb205f8SBarry Smith   if (bmat->a) {
381ccb205f8SBarry Smith     for (i=0; i<bmat->nz; i++) {
382ccb205f8SBarry Smith       if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);}
383ccb205f8SBarry Smith     }
384ccb205f8SBarry Smith   }
385ccb205f8SBarry Smith   ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
386421e10b8SBarry Smith   ierr = PetscFree(bmat);CHKERRQ(ierr);
387421e10b8SBarry Smith   PetscFunctionReturn(0);
388421e10b8SBarry Smith }
389421e10b8SBarry Smith 
390421e10b8SBarry Smith #undef __FUNCT__
391421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat"
392421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
393421e10b8SBarry Smith {
394421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
395421e10b8SBarry Smith   PetscErrorCode ierr;
396421e10b8SBarry Smith   PetscScalar    *xx,*yy;
397421e10b8SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
398421e10b8SBarry Smith   Mat            *aa;
399421e10b8SBarry Smith 
400421e10b8SBarry Smith   PetscFunctionBegin;
401ccb205f8SBarry Smith   CHKMEMQ;
402421e10b8SBarry Smith   /*
403421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
404421e10b8SBarry Smith   */
405421e10b8SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
406ccb205f8SBarry Smith 
407421e10b8SBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
408421e10b8SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
409421e10b8SBarry Smith   aj  = bmat->j;
410421e10b8SBarry Smith   aa  = bmat->a;
411421e10b8SBarry Smith   ii  = bmat->i;
412421e10b8SBarry Smith   for (i=0; i<m; i++) {
413421e10b8SBarry Smith     jrow = ii[i];
414ccb205f8SBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
415421e10b8SBarry Smith     n    = ii[i+1] - jrow;
416421e10b8SBarry Smith     for (j=0; j<n; j++) {
417421e10b8SBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
418ccb205f8SBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
419421e10b8SBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
420421e10b8SBarry Smith       jrow++;
421421e10b8SBarry Smith     }
422421e10b8SBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
423421e10b8SBarry Smith   }
424421e10b8SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
425421e10b8SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
426ccb205f8SBarry Smith   CHKMEMQ;
427421e10b8SBarry Smith   PetscFunctionReturn(0);
428421e10b8SBarry Smith }
429421e10b8SBarry Smith 
430421e10b8SBarry Smith #undef __FUNCT__
431*290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric"
432*290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
433*290bbb0aSBarry Smith {
434*290bbb0aSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
435*290bbb0aSBarry Smith   PetscErrorCode ierr;
436*290bbb0aSBarry Smith   PetscScalar    *xx,*yy;
437*290bbb0aSBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
438*290bbb0aSBarry Smith   Mat            *aa;
439*290bbb0aSBarry Smith 
440*290bbb0aSBarry Smith   PetscFunctionBegin;
441*290bbb0aSBarry Smith   CHKMEMQ;
442*290bbb0aSBarry Smith   /*
443*290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
444*290bbb0aSBarry Smith   */
445*290bbb0aSBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
446*290bbb0aSBarry Smith 
447*290bbb0aSBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
448*290bbb0aSBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
449*290bbb0aSBarry Smith   aj  = bmat->j;
450*290bbb0aSBarry Smith   aa  = bmat->a;
451*290bbb0aSBarry Smith   ii  = bmat->i;
452*290bbb0aSBarry Smith   for (i=0; i<m; i++) {
453*290bbb0aSBarry Smith     jrow = ii[i];
454*290bbb0aSBarry Smith     n    = ii[i+1] - jrow;
455*290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
456*290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr);
457*290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
458*290bbb0aSBarry Smith     if (aj[jrow] == i) {
459*290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
460*290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
461*290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
462*290bbb0aSBarry Smith       jrow++;
463*290bbb0aSBarry Smith       n--;
464*290bbb0aSBarry Smith     }
465*290bbb0aSBarry Smith     for (j=0; j<n; j++) {
466*290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);            /* upper triangular part */
467*290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
468*290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
469*290bbb0aSBarry Smith 
470*290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr);            /* lower triangular part */
471*290bbb0aSBarry Smith       ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr);
472*290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
473*290bbb0aSBarry Smith       jrow++;
474*290bbb0aSBarry Smith     }
475*290bbb0aSBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
476*290bbb0aSBarry Smith     ierr = VecResetArray(bmat->middle);CHKERRQ(ierr);
477*290bbb0aSBarry Smith   }
478*290bbb0aSBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
479*290bbb0aSBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
480*290bbb0aSBarry Smith   CHKMEMQ;
481*290bbb0aSBarry Smith   PetscFunctionReturn(0);
482*290bbb0aSBarry Smith }
483*290bbb0aSBarry Smith 
484*290bbb0aSBarry Smith #undef __FUNCT__
485421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat"
486421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
487421e10b8SBarry Smith {
488421e10b8SBarry Smith   PetscFunctionBegin;
489421e10b8SBarry Smith   PetscFunctionReturn(0);
490421e10b8SBarry Smith }
491421e10b8SBarry Smith 
492421e10b8SBarry Smith #undef __FUNCT__
493421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat"
494421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
495421e10b8SBarry Smith {
496421e10b8SBarry Smith   PetscFunctionBegin;
497421e10b8SBarry Smith   PetscFunctionReturn(0);
498421e10b8SBarry Smith }
499421e10b8SBarry Smith 
500421e10b8SBarry Smith #undef __FUNCT__
501421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat"
502421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
503421e10b8SBarry Smith {
504421e10b8SBarry Smith   PetscFunctionBegin;
505421e10b8SBarry Smith   PetscFunctionReturn(0);
506421e10b8SBarry Smith }
507421e10b8SBarry Smith 
508421e10b8SBarry Smith #undef __FUNCT__
509421e10b8SBarry Smith #define __FUNCT__ "MatSetBlockSize_BlockMat"
510421e10b8SBarry Smith PetscErrorCode MatSetBlockSize_BlockMat(Mat A,PetscInt bs)
511421e10b8SBarry Smith {
512421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
513421e10b8SBarry Smith   PetscErrorCode ierr;
514ccb205f8SBarry Smith   PetscInt       nz = 10,i;
515421e10b8SBarry Smith 
516421e10b8SBarry Smith   PetscFunctionBegin;
517ccb205f8SBarry Smith   if (A->rmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap.n);
518ccb205f8SBarry Smith   if (A->cmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap.n);
519421e10b8SBarry Smith   A->rmap.bs = A->cmap.bs = bs;
520421e10b8SBarry Smith 
521ccb205f8SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr);
522*290bbb0aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->middle);CHKERRQ(ierr);
523ccb205f8SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr);
524ccb205f8SBarry Smith 
525ccb205f8SBarry Smith 
526421e10b8SBarry Smith   ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr);
527421e10b8SBarry Smith   for (i=0; i<A->rmap.n; i++) bmat->imax[i] = nz;
528421e10b8SBarry Smith   nz = nz*A->rmap.n;
529421e10b8SBarry Smith 
530ccb205f8SBarry Smith   bmat->mbs = A->rmap.n/A->rmap.bs;
531421e10b8SBarry Smith 
532421e10b8SBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
533ccb205f8SBarry Smith   for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;}
534421e10b8SBarry Smith 
535421e10b8SBarry Smith   /* allocate the matrix space */
536421e10b8SBarry Smith   ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr);
537421e10b8SBarry Smith   bmat->i[0] = 0;
538ccb205f8SBarry Smith   for (i=1; i<bmat->mbs+1; i++) {
539421e10b8SBarry Smith     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
540421e10b8SBarry Smith   }
541421e10b8SBarry Smith   bmat->singlemalloc = PETSC_TRUE;
542421e10b8SBarry Smith   bmat->free_a       = PETSC_TRUE;
543421e10b8SBarry Smith   bmat->free_ij      = PETSC_TRUE;
544421e10b8SBarry Smith 
545421e10b8SBarry Smith   bmat->nz                = 0;
546421e10b8SBarry Smith   bmat->maxnz             = nz;
547421e10b8SBarry Smith   A->info.nz_unneeded  = (double)bmat->maxnz;
548421e10b8SBarry Smith 
549421e10b8SBarry Smith   PetscFunctionReturn(0);
550421e10b8SBarry Smith }
551421e10b8SBarry Smith 
552ccb205f8SBarry Smith /*
553ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
554ccb205f8SBarry Smith */
555ccb205f8SBarry Smith #undef __FUNCT__
556ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat"
557ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
558ccb205f8SBarry Smith {
559ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
560ccb205f8SBarry Smith   PetscErrorCode ierr;
561ccb205f8SBarry Smith   PetscInt       i,j,mbs = A->rmap.n/A->rmap.bs;
562ccb205f8SBarry Smith 
563ccb205f8SBarry Smith   PetscFunctionBegin;
564ccb205f8SBarry Smith   if (!a->diag) {
565ccb205f8SBarry Smith     ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
566ccb205f8SBarry Smith   }
567ccb205f8SBarry Smith   for (i=0; i<mbs; i++) {
568ccb205f8SBarry Smith     a->diag[i] = a->i[i+1];
569ccb205f8SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
570ccb205f8SBarry Smith       if (a->j[j] == i) {
571ccb205f8SBarry Smith         a->diag[i] = j;
572ccb205f8SBarry Smith         break;
573ccb205f8SBarry Smith       }
574ccb205f8SBarry Smith     }
575ccb205f8SBarry Smith   }
576ccb205f8SBarry Smith   PetscFunctionReturn(0);
577ccb205f8SBarry Smith }
578ccb205f8SBarry Smith 
579ccb205f8SBarry Smith #undef __FUNCT__
580ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat"
581ccb205f8SBarry Smith PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
582ccb205f8SBarry Smith {
583ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
584ccb205f8SBarry Smith   Mat_SeqAIJ     *c;
585ccb205f8SBarry Smith   PetscErrorCode ierr;
586ccb205f8SBarry Smith   PetscInt       i,k,first,step,lensi,nrows,ncols;
587ccb205f8SBarry Smith   PetscInt       *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
588ccb205f8SBarry Smith   PetscScalar    *a_new,value;
589ccb205f8SBarry Smith   Mat            C,*aa = a->a;
590ccb205f8SBarry Smith   PetscTruth     stride,equal;
591ccb205f8SBarry Smith 
592ccb205f8SBarry Smith   PetscFunctionBegin;
593ccb205f8SBarry Smith   ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr);
594ccb205f8SBarry Smith   if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices");
595ccb205f8SBarry Smith   ierr = ISStride(iscol,&stride);CHKERRQ(ierr);
596ccb205f8SBarry Smith   if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices");
597ccb205f8SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
598ccb205f8SBarry Smith   if (step != A->rmap.bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block");
599ccb205f8SBarry Smith 
600ccb205f8SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
601ccb205f8SBarry Smith   ncols = nrows;
602ccb205f8SBarry Smith 
603ccb205f8SBarry Smith   /* create submatrix */
604ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
605ccb205f8SBarry Smith     PetscInt n_cols,n_rows;
606ccb205f8SBarry Smith     C = *B;
607ccb205f8SBarry Smith     ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr);
608ccb205f8SBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
609ccb205f8SBarry Smith     ierr = MatZeroEntries(C);CHKERRQ(ierr);
610ccb205f8SBarry Smith   } else {
611ccb205f8SBarry Smith     ierr = MatCreate(A->comm,&C);CHKERRQ(ierr);
612ccb205f8SBarry Smith     ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
613ccb205f8SBarry Smith     ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr);
614ccb205f8SBarry Smith     ierr = MatSeqAIJSetPreallocation_SeqAIJ(C,0,ailen);CHKERRQ(ierr);
615ccb205f8SBarry Smith   }
616ccb205f8SBarry Smith   c = (Mat_SeqAIJ*)C->data;
617ccb205f8SBarry Smith 
618ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
619ccb205f8SBarry Smith   a_new    = c->a;
620ccb205f8SBarry Smith   j_new    = c->j;
621ccb205f8SBarry Smith   i_new    = c->i;
622ccb205f8SBarry Smith 
623ccb205f8SBarry Smith   for (i=0; i<nrows; i++) {
624ccb205f8SBarry Smith     ii    = ai[i];
625ccb205f8SBarry Smith     lensi = ailen[i];
626ccb205f8SBarry Smith     for (k=0; k<lensi; k++) {
627ccb205f8SBarry Smith       *j_new++ = *aj++;
628ccb205f8SBarry Smith       ierr     = MatGetValue(*aa++,first,first,value);CHKERRQ(ierr);
629ccb205f8SBarry Smith       *a_new++ = value;
630ccb205f8SBarry Smith     }
631ccb205f8SBarry Smith     i_new[i+1]  = i_new[i] + lensi;
632ccb205f8SBarry Smith     c->ilen[i]  = lensi;
633ccb205f8SBarry Smith   }
634ccb205f8SBarry Smith 
635ccb205f8SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
636ccb205f8SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
637ccb205f8SBarry Smith   *B = C;
638ccb205f8SBarry Smith   PetscFunctionReturn(0);
639ccb205f8SBarry Smith }
640ccb205f8SBarry Smith 
641421e10b8SBarry Smith #undef __FUNCT__
642421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat"
643421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
644421e10b8SBarry Smith {
645421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
646421e10b8SBarry Smith   PetscErrorCode ierr;
647421e10b8SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
648ccb205f8SBarry Smith   PetscInt       m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
649421e10b8SBarry Smith   Mat            *aa = a->a,*ap;
650421e10b8SBarry Smith 
651421e10b8SBarry Smith   PetscFunctionBegin;
652421e10b8SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
653421e10b8SBarry Smith 
654421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
655421e10b8SBarry Smith   for (i=1; i<m; i++) {
656421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
657421e10b8SBarry Smith     fshift += imax[i-1] - ailen[i-1];
658421e10b8SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
659421e10b8SBarry Smith     if (fshift) {
660421e10b8SBarry Smith       ip = aj + ai[i] ;
661421e10b8SBarry Smith       ap = aa + ai[i] ;
662421e10b8SBarry Smith       N  = ailen[i];
663421e10b8SBarry Smith       for (j=0; j<N; j++) {
664421e10b8SBarry Smith         ip[j-fshift] = ip[j];
665421e10b8SBarry Smith         ap[j-fshift] = ap[j];
666421e10b8SBarry Smith       }
667421e10b8SBarry Smith     }
668421e10b8SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
669421e10b8SBarry Smith   }
670421e10b8SBarry Smith   if (m) {
671421e10b8SBarry Smith     fshift += imax[m-1] - ailen[m-1];
672421e10b8SBarry Smith     ai[m]  = ai[m-1] + ailen[m-1];
673421e10b8SBarry Smith   }
674421e10b8SBarry Smith   /* reset ilen and imax for each row */
675421e10b8SBarry Smith   for (i=0; i<m; i++) {
676421e10b8SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
677421e10b8SBarry Smith   }
678421e10b8SBarry Smith   a->nz = ai[m];
679ccb205f8SBarry Smith   for (i=0; i<a->nz; i++) {
680ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG)
681ccb205f8SBarry Smith     if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz);
682ccb205f8SBarry Smith #endif
683ccb205f8SBarry Smith     ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
684ccb205f8SBarry Smith     ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
685ccb205f8SBarry Smith   }
686ccb205f8SBarry Smith   CHKMEMQ;
687421e10b8SBarry Smith   ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);CHKERRQ(ierr);
688421e10b8SBarry Smith   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr);
689421e10b8SBarry Smith   ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr);
690421e10b8SBarry Smith   a->reallocs          = 0;
691421e10b8SBarry Smith   A->info.nz_unneeded  = (double)fshift;
692421e10b8SBarry Smith   a->rmax              = rmax;
693421e10b8SBarry Smith 
694421e10b8SBarry Smith   A->same_nonzero = PETSC_TRUE;
695ccb205f8SBarry Smith   ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr);
696421e10b8SBarry Smith   PetscFunctionReturn(0);
697421e10b8SBarry Smith }
698421e10b8SBarry Smith 
699*290bbb0aSBarry Smith #undef __FUNCT__
700*290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat"
701*290bbb0aSBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt)
702*290bbb0aSBarry Smith {
703*290bbb0aSBarry Smith   PetscFunctionBegin;
704*290bbb0aSBarry Smith   if (opt == MAT_SYMMETRIC) {
705*290bbb0aSBarry Smith     A->ops->relax = MatRelax_BlockMat_Symmetric;
706*290bbb0aSBarry Smith     A->ops->mult  = MatMult_BlockMat_Symmetric;
707*290bbb0aSBarry Smith   } else {
708*290bbb0aSBarry Smith     PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]);
709*290bbb0aSBarry Smith   }
710*290bbb0aSBarry Smith   PetscFunctionReturn(0);
711*290bbb0aSBarry Smith }
712*290bbb0aSBarry Smith 
713*290bbb0aSBarry Smith 
714421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
715421e10b8SBarry Smith        0,
716421e10b8SBarry Smith        0,
717421e10b8SBarry Smith        MatMult_BlockMat,
718421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat,
719421e10b8SBarry Smith        MatMultTranspose_BlockMat,
720421e10b8SBarry Smith        MatMultTransposeAdd_BlockMat,
721421e10b8SBarry Smith        0,
722421e10b8SBarry Smith        0,
723421e10b8SBarry Smith        0,
724421e10b8SBarry Smith /*10*/ 0,
725421e10b8SBarry Smith        0,
726421e10b8SBarry Smith        0,
727ccb205f8SBarry Smith        MatRelax_BlockMat,
728421e10b8SBarry Smith        0,
729421e10b8SBarry Smith /*15*/ 0,
730421e10b8SBarry Smith        0,
731421e10b8SBarry Smith        0,
732421e10b8SBarry Smith        0,
733421e10b8SBarry Smith        0,
734421e10b8SBarry Smith /*20*/ 0,
735421e10b8SBarry Smith        MatAssemblyEnd_BlockMat,
736421e10b8SBarry Smith        0,
737*290bbb0aSBarry Smith        MatSetOption_BlockMat,
738421e10b8SBarry Smith        0,
739421e10b8SBarry Smith /*25*/ 0,
740421e10b8SBarry Smith        0,
741421e10b8SBarry Smith        0,
742421e10b8SBarry Smith        0,
743421e10b8SBarry Smith        0,
744421e10b8SBarry Smith /*30*/ 0,
745421e10b8SBarry Smith        0,
746421e10b8SBarry Smith        0,
747421e10b8SBarry Smith        0,
748421e10b8SBarry Smith        0,
749421e10b8SBarry Smith /*35*/ 0,
750421e10b8SBarry Smith        0,
751421e10b8SBarry Smith        0,
752421e10b8SBarry Smith        0,
753421e10b8SBarry Smith        0,
754421e10b8SBarry Smith /*40*/ 0,
755421e10b8SBarry Smith        0,
756421e10b8SBarry Smith        0,
757421e10b8SBarry Smith        0,
758421e10b8SBarry Smith        0,
759421e10b8SBarry Smith /*45*/ 0,
760421e10b8SBarry Smith        0,
761421e10b8SBarry Smith        0,
762421e10b8SBarry Smith        0,
763421e10b8SBarry Smith        0,
764421e10b8SBarry Smith /*50*/ MatSetBlockSize_BlockMat,
765421e10b8SBarry Smith        0,
766421e10b8SBarry Smith        0,
767421e10b8SBarry Smith        0,
768421e10b8SBarry Smith        0,
769421e10b8SBarry Smith /*55*/ 0,
770421e10b8SBarry Smith        0,
771421e10b8SBarry Smith        0,
772421e10b8SBarry Smith        0,
773421e10b8SBarry Smith        0,
774ccb205f8SBarry Smith /*60*/ MatGetSubMatrix_BlockMat,
775421e10b8SBarry Smith        MatDestroy_BlockMat,
776ccb205f8SBarry Smith        MatView_BlockMat,
777421e10b8SBarry Smith        0,
778421e10b8SBarry Smith        0,
779421e10b8SBarry Smith /*65*/ 0,
780421e10b8SBarry Smith        0,
781421e10b8SBarry Smith        0,
782421e10b8SBarry Smith        0,
783421e10b8SBarry Smith        0,
784421e10b8SBarry Smith /*70*/ 0,
785421e10b8SBarry Smith        0,
786421e10b8SBarry Smith        0,
787421e10b8SBarry Smith        0,
788421e10b8SBarry Smith        0,
789421e10b8SBarry Smith /*75*/ 0,
790421e10b8SBarry Smith        0,
791421e10b8SBarry Smith        0,
792421e10b8SBarry Smith        0,
793421e10b8SBarry Smith        0,
794421e10b8SBarry Smith /*80*/ 0,
795421e10b8SBarry Smith        0,
796421e10b8SBarry Smith        0,
797421e10b8SBarry Smith        0,
798421e10b8SBarry Smith        MatLoad_BlockMat,
799421e10b8SBarry Smith /*85*/ 0,
800421e10b8SBarry Smith        0,
801421e10b8SBarry Smith        0,
802421e10b8SBarry Smith        0,
803421e10b8SBarry Smith        0,
804421e10b8SBarry Smith /*90*/ 0,
805421e10b8SBarry Smith        0,
806421e10b8SBarry Smith        0,
807421e10b8SBarry Smith        0,
808421e10b8SBarry Smith        0,
809421e10b8SBarry Smith /*95*/ 0,
810421e10b8SBarry Smith        0,
811421e10b8SBarry Smith        0,
812421e10b8SBarry Smith        0};
813421e10b8SBarry Smith 
814421e10b8SBarry Smith /*MC
815421e10b8SBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
816421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
817421e10b8SBarry Smith 
818421e10b8SBarry Smith   Level: advanced
819421e10b8SBarry Smith 
820421e10b8SBarry Smith .seealso: MatCreateBlockMat()
821421e10b8SBarry Smith 
822421e10b8SBarry Smith M*/
823421e10b8SBarry Smith 
824421e10b8SBarry Smith EXTERN_C_BEGIN
825421e10b8SBarry Smith #undef __FUNCT__
826421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat"
827421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A)
828421e10b8SBarry Smith {
829421e10b8SBarry Smith   Mat_BlockMat   *b;
830421e10b8SBarry Smith   PetscErrorCode ierr;
831421e10b8SBarry Smith 
832421e10b8SBarry Smith   PetscFunctionBegin;
833421e10b8SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
834421e10b8SBarry Smith   ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr);
835421e10b8SBarry Smith 
836421e10b8SBarry Smith   A->data = (void*)b;
837421e10b8SBarry Smith 
838421e10b8SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
839421e10b8SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
840421e10b8SBarry Smith 
841421e10b8SBarry Smith   A->assembled     = PETSC_TRUE;
842421e10b8SBarry Smith   A->preallocated  = PETSC_FALSE;
843421e10b8SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
844*290bbb0aSBarry Smith 
845*290bbb0aSBarry Smith   ierr = PetscOptionsBegin(A->comm,A->prefix,"Matrix Option","Mat");CHKERRQ(ierr);
846*290bbb0aSBarry Smith   ierr = PetscOptionsEnd();
847*290bbb0aSBarry Smith 
848421e10b8SBarry Smith   PetscFunctionReturn(0);
849421e10b8SBarry Smith }
850421e10b8SBarry Smith EXTERN_C_END
851421e10b8SBarry Smith 
852421e10b8SBarry Smith #undef __FUNCT__
853421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat"
854421e10b8SBarry Smith /*@C
855421e10b8SBarry Smith    MatCreateBlockMat - Creates a new matrix based sparse Mat storage
856421e10b8SBarry Smith 
857421e10b8SBarry Smith   Collective on MPI_Comm
858421e10b8SBarry Smith 
859421e10b8SBarry Smith    Input Parameters:
860421e10b8SBarry Smith +  comm - MPI communicator
861421e10b8SBarry Smith .  m - number of rows
862421e10b8SBarry Smith .  n  - number of columns
863421e10b8SBarry Smith -  bs - size of each submatrix
864421e10b8SBarry Smith 
865421e10b8SBarry Smith 
866421e10b8SBarry Smith    Output Parameter:
867421e10b8SBarry Smith .  A - the matrix
868421e10b8SBarry Smith 
869421e10b8SBarry Smith    Level: intermediate
870421e10b8SBarry Smith 
871421e10b8SBarry Smith    PETSc requires that matrices and vectors being used for certain
872421e10b8SBarry Smith    operations are partitioned accordingly.  For example, when
873421e10b8SBarry Smith    creating a bmat matrix, A, that supports parallel matrix-vector
874421e10b8SBarry Smith    products using MatMult(A,x,y) the user should set the number
875421e10b8SBarry Smith    of local matrix rows to be the number of local elements of the
876421e10b8SBarry Smith    corresponding result vector, y. Note that this is information is
877421e10b8SBarry Smith    required for use of the matrix interface routines, even though
878421e10b8SBarry Smith    the bmat matrix may not actually be physically partitioned.
879421e10b8SBarry Smith    For example,
880421e10b8SBarry Smith 
881421e10b8SBarry Smith .keywords: matrix, bmat, create
882421e10b8SBarry Smith 
883421e10b8SBarry Smith .seealso: MATBLOCKMAT
884421e10b8SBarry Smith @*/
885421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,Mat *A)
886421e10b8SBarry Smith {
887421e10b8SBarry Smith   PetscErrorCode ierr;
888421e10b8SBarry Smith 
889421e10b8SBarry Smith   PetscFunctionBegin;
890421e10b8SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
891421e10b8SBarry Smith   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
892421e10b8SBarry Smith   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
893421e10b8SBarry Smith   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
894421e10b8SBarry Smith   PetscFunctionReturn(0);
895421e10b8SBarry Smith }
896421e10b8SBarry Smith 
897421e10b8SBarry Smith 
898421e10b8SBarry Smith 
899