xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision b02232071040b670e3b9229765e5d1cb65043fd5)
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 
18290bbb0aSBarry Smith   Vec               left,right,middle;   /* dummy vectors to perform local parts of product */
19421e10b8SBarry Smith } Mat_BlockMat;
20421e10b8SBarry Smith 
21421e10b8SBarry Smith #undef __FUNCT__
22290bbb0aSBarry Smith #define __FUNCT__ "MatRelax_BlockMat_Symmetric"
23290bbb0aSBarry Smith PetscErrorCode MatRelax_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
24290bbb0aSBarry Smith {
25290bbb0aSBarry Smith   Mat_BlockMat       *a = (Mat_BlockMat*)A->data;
26290bbb0aSBarry Smith   PetscScalar        *x;
27290bbb0aSBarry Smith   const Mat          *v = a->a;
28290bbb0aSBarry Smith   const PetscScalar  *b;
29290bbb0aSBarry Smith   PetscErrorCode     ierr;
30290bbb0aSBarry Smith   PetscInt           n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs;
31290bbb0aSBarry Smith   const PetscInt     *idx;
32290bbb0aSBarry Smith   IS                 row,col;
33290bbb0aSBarry Smith   MatFactorInfo      info;
34290bbb0aSBarry Smith   Vec                left = a->left,right = a->right, middle;
35290bbb0aSBarry Smith   Mat                *diag;
36290bbb0aSBarry Smith 
37290bbb0aSBarry Smith   PetscFunctionBegin;
38290bbb0aSBarry Smith   its = its*lits;
39290bbb0aSBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
40290bbb0aSBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
41290bbb0aSBarry Smith   if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
42290bbb0aSBarry Smith   if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift");
43290bbb0aSBarry Smith   if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP))
44290bbb0aSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
45290bbb0aSBarry Smith 
46290bbb0aSBarry Smith   if (!a->diags) {
47290bbb0aSBarry Smith     ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr);
48290bbb0aSBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
49290bbb0aSBarry Smith     for (i=0; i<mbs; i++) {
50290bbb0aSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr);
51290bbb0aSBarry Smith       ierr = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr);
52290bbb0aSBarry Smith       ierr = MatLUFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr);
53290bbb0aSBarry Smith       ierr = ISDestroy(row);CHKERRQ(ierr);
54290bbb0aSBarry Smith       ierr = ISDestroy(col);CHKERRQ(ierr);
55290bbb0aSBarry Smith     }
56290bbb0aSBarry Smith   }
57290bbb0aSBarry Smith   middle  = a->middle;
58290bbb0aSBarry Smith   diag    = a->diags;
59290bbb0aSBarry Smith 
60290bbb0aSBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
61290bbb0aSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
62290bbb0aSBarry Smith   /* copy right hand side because it must be modified during iteration */
63290bbb0aSBarry Smith   ierr = VecCopy(bb,middle);CHKERRQ(ierr);
64290bbb0aSBarry Smith   ierr = VecGetArray(middle,(PetscScalar**)&b);CHKERRQ(ierr);
65290bbb0aSBarry Smith 
66290bbb0aSBarry Smith   /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */
67290bbb0aSBarry Smith   while (its--) {
68290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
69290bbb0aSBarry Smith 
70290bbb0aSBarry Smith       for (i=0; i<mbs; i++) {
71290bbb0aSBarry Smith         n    = a->i[i+1] - a->i[i];
72290bbb0aSBarry Smith         idx  = a->j + a->i[i];
73290bbb0aSBarry Smith         v    = a->a + a->i[i];
74290bbb0aSBarry Smith 
75290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
76290bbb0aSBarry Smith         for (j=0; j<n; j++) {
77290bbb0aSBarry Smith           if (idx[j] != i) {
78290bbb0aSBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
79290bbb0aSBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
80290bbb0aSBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
81290bbb0aSBarry Smith           }
82290bbb0aSBarry Smith         }
83290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
84290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
85290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
86290bbb0aSBarry Smith 
87290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
88290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
89290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
90290bbb0aSBarry Smith       }
91290bbb0aSBarry Smith     }
92290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
93290bbb0aSBarry Smith 
94290bbb0aSBarry Smith       for (i=mbs-1; i>=0; i--) {
95290bbb0aSBarry Smith         n    = a->i[i+1] - a->i[i];
96290bbb0aSBarry Smith         idx  = a->j + a->i[i];
97290bbb0aSBarry Smith         v    = a->a + a->i[i];
98290bbb0aSBarry Smith 
99290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
100290bbb0aSBarry Smith         for (j=0; j<n; j++) {
101290bbb0aSBarry Smith           if (idx[j] != i) {
102290bbb0aSBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
103290bbb0aSBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
104290bbb0aSBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
105290bbb0aSBarry Smith           }
106290bbb0aSBarry Smith         }
107290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
108290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
109290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
110290bbb0aSBarry Smith 
111290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
112290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
113290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
114290bbb0aSBarry Smith 
115290bbb0aSBarry Smith       }
116290bbb0aSBarry Smith     }
117290bbb0aSBarry Smith   }
118290bbb0aSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
119290bbb0aSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
120290bbb0aSBarry Smith   PetscFunctionReturn(0);
121290bbb0aSBarry Smith }
122290bbb0aSBarry Smith 
123290bbb0aSBarry 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*b0223207SBarry Smith       if (A->symmetric && row > col) 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);*/
291*b0223207SBarry Smith         if (A->symmetric && brow == bcol) {
292*b0223207SBarry Smith           ierr = MatCreateSeqSBAIJ(PETSC_COMM_SELF,1,bs,bs,0,0,ap+i);CHKERRQ(ierr);
293*b0223207SBarry Smith         } else {
294421e10b8SBarry Smith           ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
295421e10b8SBarry Smith         }
296*b0223207SBarry Smith       }
297ccb205f8SBarry Smith       /*      printf("numerical value at i %d row %d col %d cidx %d ridx %d value %g\n",i,brow,bcol,cidx,ridx,value);*/
298421e10b8SBarry Smith       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
299421e10b8SBarry Smith       low = i;
300421e10b8SBarry Smith     }
301ccb205f8SBarry Smith     /*    printf("nrow for row %d %d\n",nrow,brow);*/
302421e10b8SBarry Smith     ailen[brow] = nrow;
303421e10b8SBarry Smith   }
304421e10b8SBarry Smith   A->same_nonzero = PETSC_FALSE;
305421e10b8SBarry Smith   PetscFunctionReturn(0);
306421e10b8SBarry Smith }
307421e10b8SBarry Smith 
308421e10b8SBarry Smith #undef __FUNCT__
309421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat"
310421e10b8SBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A)
311421e10b8SBarry Smith {
312421e10b8SBarry Smith   PetscErrorCode    ierr;
313421e10b8SBarry Smith   Mat               tmpA;
314421e10b8SBarry Smith   PetscInt          i,m,n,bs = 1,ncols;
315421e10b8SBarry Smith   const PetscInt    *cols;
316421e10b8SBarry Smith   const PetscScalar *values;
317290bbb0aSBarry Smith   PetscTruth        flg;
318421e10b8SBarry Smith 
319421e10b8SBarry Smith   PetscFunctionBegin;
320421e10b8SBarry Smith   ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr);
321421e10b8SBarry Smith 
322421e10b8SBarry Smith   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
323421e10b8SBarry Smith   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr);
324421e10b8SBarry Smith     ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr);
325421e10b8SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
326421e10b8SBarry Smith   ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,A);CHKERRQ(ierr);
327290bbb0aSBarry Smith   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr);
328290bbb0aSBarry Smith     ierr = PetscOptionsName("-matload_symmetric","Store the matrix as symmetric","MatLoad",&flg);CHKERRQ(ierr);
329290bbb0aSBarry Smith     if (flg) {
330290bbb0aSBarry Smith       ierr = MatSetOption(*A,MAT_SYMMETRIC);CHKERRQ(ierr);
331290bbb0aSBarry Smith     }
332290bbb0aSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
333421e10b8SBarry Smith   for (i=0; i<m; i++) {
334421e10b8SBarry Smith     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
335421e10b8SBarry Smith     ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
336ccb205f8SBarry Smith     ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
337421e10b8SBarry Smith   }
338421e10b8SBarry Smith   ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
339421e10b8SBarry Smith   ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
340421e10b8SBarry Smith   PetscFunctionReturn(0);
341421e10b8SBarry Smith }
342421e10b8SBarry Smith 
343ccb205f8SBarry Smith #undef __FUNCT__
344ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat"
345ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
346ccb205f8SBarry Smith {
347ccb205f8SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
348ccb205f8SBarry Smith   PetscErrorCode    ierr;
349ccb205f8SBarry Smith   const char        *name;
350ccb205f8SBarry Smith   PetscViewerFormat format;
351ccb205f8SBarry Smith 
352ccb205f8SBarry Smith   PetscFunctionBegin;
353ccb205f8SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
354ccb205f8SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
355ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
356ccb205f8SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr);
357ccb205f8SBarry Smith   }
358ccb205f8SBarry Smith   PetscFunctionReturn(0);
359ccb205f8SBarry Smith }
360421e10b8SBarry Smith 
361421e10b8SBarry Smith #undef __FUNCT__
362421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat"
363421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat)
364421e10b8SBarry Smith {
365421e10b8SBarry Smith   PetscErrorCode ierr;
366421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
367ccb205f8SBarry Smith   PetscInt       i;
368421e10b8SBarry Smith 
369421e10b8SBarry Smith   PetscFunctionBegin;
370421e10b8SBarry Smith   if (bmat->right) {
371421e10b8SBarry Smith     ierr = VecDestroy(bmat->right);CHKERRQ(ierr);
372421e10b8SBarry Smith   }
373421e10b8SBarry Smith   if (bmat->left) {
374421e10b8SBarry Smith     ierr = VecDestroy(bmat->left);CHKERRQ(ierr);
375421e10b8SBarry Smith   }
376290bbb0aSBarry Smith   if (bmat->middle) {
377290bbb0aSBarry Smith     ierr = VecDestroy(bmat->middle);CHKERRQ(ierr);
378290bbb0aSBarry Smith   }
379ccb205f8SBarry Smith   if (bmat->diags) {
380ccb205f8SBarry Smith     for (i=0; i<mat->rmap.n/mat->rmap.bs; i++) {
381ccb205f8SBarry Smith       if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);}
382ccb205f8SBarry Smith     }
383ccb205f8SBarry Smith   }
384ccb205f8SBarry Smith   if (bmat->a) {
385ccb205f8SBarry Smith     for (i=0; i<bmat->nz; i++) {
386ccb205f8SBarry Smith       if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);}
387ccb205f8SBarry Smith     }
388ccb205f8SBarry Smith   }
389ccb205f8SBarry Smith   ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
390421e10b8SBarry Smith   ierr = PetscFree(bmat);CHKERRQ(ierr);
391421e10b8SBarry Smith   PetscFunctionReturn(0);
392421e10b8SBarry Smith }
393421e10b8SBarry Smith 
394421e10b8SBarry Smith #undef __FUNCT__
395421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat"
396421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
397421e10b8SBarry Smith {
398421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
399421e10b8SBarry Smith   PetscErrorCode ierr;
400421e10b8SBarry Smith   PetscScalar    *xx,*yy;
401421e10b8SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
402421e10b8SBarry Smith   Mat            *aa;
403421e10b8SBarry Smith 
404421e10b8SBarry Smith   PetscFunctionBegin;
405ccb205f8SBarry Smith   CHKMEMQ;
406421e10b8SBarry Smith   /*
407421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
408421e10b8SBarry Smith   */
409421e10b8SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
410ccb205f8SBarry Smith 
411421e10b8SBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
412421e10b8SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
413421e10b8SBarry Smith   aj  = bmat->j;
414421e10b8SBarry Smith   aa  = bmat->a;
415421e10b8SBarry Smith   ii  = bmat->i;
416421e10b8SBarry Smith   for (i=0; i<m; i++) {
417421e10b8SBarry Smith     jrow = ii[i];
418ccb205f8SBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
419421e10b8SBarry Smith     n    = ii[i+1] - jrow;
420421e10b8SBarry Smith     for (j=0; j<n; j++) {
421421e10b8SBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
422ccb205f8SBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
423421e10b8SBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
424421e10b8SBarry Smith       jrow++;
425421e10b8SBarry Smith     }
426421e10b8SBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
427421e10b8SBarry Smith   }
428421e10b8SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
429421e10b8SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
430ccb205f8SBarry Smith   CHKMEMQ;
431421e10b8SBarry Smith   PetscFunctionReturn(0);
432421e10b8SBarry Smith }
433421e10b8SBarry Smith 
434421e10b8SBarry Smith #undef __FUNCT__
435290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric"
436290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
437290bbb0aSBarry Smith {
438290bbb0aSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
439290bbb0aSBarry Smith   PetscErrorCode ierr;
440290bbb0aSBarry Smith   PetscScalar    *xx,*yy;
441290bbb0aSBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
442290bbb0aSBarry Smith   Mat            *aa;
443290bbb0aSBarry Smith 
444290bbb0aSBarry Smith   PetscFunctionBegin;
445290bbb0aSBarry Smith   CHKMEMQ;
446290bbb0aSBarry Smith   /*
447290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
448290bbb0aSBarry Smith   */
449290bbb0aSBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
450290bbb0aSBarry Smith 
451290bbb0aSBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
452290bbb0aSBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
453290bbb0aSBarry Smith   aj  = bmat->j;
454290bbb0aSBarry Smith   aa  = bmat->a;
455290bbb0aSBarry Smith   ii  = bmat->i;
456290bbb0aSBarry Smith   for (i=0; i<m; i++) {
457290bbb0aSBarry Smith     jrow = ii[i];
458290bbb0aSBarry Smith     n    = ii[i+1] - jrow;
459290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
460290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr);
461290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
462290bbb0aSBarry Smith     if (aj[jrow] == i) {
463290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
464290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
465290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
466290bbb0aSBarry Smith       jrow++;
467290bbb0aSBarry Smith       n--;
468290bbb0aSBarry Smith     }
469290bbb0aSBarry Smith     for (j=0; j<n; j++) {
470290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);            /* upper triangular part */
471290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
472290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
473290bbb0aSBarry Smith 
474290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr);            /* lower triangular part */
475290bbb0aSBarry Smith       ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr);
476290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
477290bbb0aSBarry Smith       jrow++;
478290bbb0aSBarry Smith     }
479290bbb0aSBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
480290bbb0aSBarry Smith     ierr = VecResetArray(bmat->middle);CHKERRQ(ierr);
481290bbb0aSBarry Smith   }
482290bbb0aSBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
483290bbb0aSBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
484290bbb0aSBarry Smith   CHKMEMQ;
485290bbb0aSBarry Smith   PetscFunctionReturn(0);
486290bbb0aSBarry Smith }
487290bbb0aSBarry Smith 
488290bbb0aSBarry Smith #undef __FUNCT__
489421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat"
490421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
491421e10b8SBarry Smith {
492421e10b8SBarry Smith   PetscFunctionBegin;
493421e10b8SBarry Smith   PetscFunctionReturn(0);
494421e10b8SBarry Smith }
495421e10b8SBarry Smith 
496421e10b8SBarry Smith #undef __FUNCT__
497421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat"
498421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
499421e10b8SBarry Smith {
500421e10b8SBarry Smith   PetscFunctionBegin;
501421e10b8SBarry Smith   PetscFunctionReturn(0);
502421e10b8SBarry Smith }
503421e10b8SBarry Smith 
504421e10b8SBarry Smith #undef __FUNCT__
505421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat"
506421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
507421e10b8SBarry Smith {
508421e10b8SBarry Smith   PetscFunctionBegin;
509421e10b8SBarry Smith   PetscFunctionReturn(0);
510421e10b8SBarry Smith }
511421e10b8SBarry Smith 
512421e10b8SBarry Smith #undef __FUNCT__
513421e10b8SBarry Smith #define __FUNCT__ "MatSetBlockSize_BlockMat"
514421e10b8SBarry Smith PetscErrorCode MatSetBlockSize_BlockMat(Mat A,PetscInt bs)
515421e10b8SBarry Smith {
516421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
517421e10b8SBarry Smith   PetscErrorCode ierr;
518ccb205f8SBarry Smith   PetscInt       nz = 10,i;
519421e10b8SBarry Smith 
520421e10b8SBarry Smith   PetscFunctionBegin;
521ccb205f8SBarry Smith   if (A->rmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap.n);
522ccb205f8SBarry Smith   if (A->cmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap.n);
523421e10b8SBarry Smith   A->rmap.bs = A->cmap.bs = bs;
524421e10b8SBarry Smith 
525ccb205f8SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr);
526290bbb0aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->middle);CHKERRQ(ierr);
527ccb205f8SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr);
528ccb205f8SBarry Smith 
529ccb205f8SBarry Smith 
530421e10b8SBarry Smith   ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr);
531421e10b8SBarry Smith   for (i=0; i<A->rmap.n; i++) bmat->imax[i] = nz;
532421e10b8SBarry Smith   nz = nz*A->rmap.n;
533421e10b8SBarry Smith 
534ccb205f8SBarry Smith   bmat->mbs = A->rmap.n/A->rmap.bs;
535421e10b8SBarry Smith 
536421e10b8SBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
537ccb205f8SBarry Smith   for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;}
538421e10b8SBarry Smith 
539421e10b8SBarry Smith   /* allocate the matrix space */
540421e10b8SBarry Smith   ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr);
541421e10b8SBarry Smith   bmat->i[0] = 0;
542ccb205f8SBarry Smith   for (i=1; i<bmat->mbs+1; i++) {
543421e10b8SBarry Smith     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
544421e10b8SBarry Smith   }
545421e10b8SBarry Smith   bmat->singlemalloc = PETSC_TRUE;
546421e10b8SBarry Smith   bmat->free_a       = PETSC_TRUE;
547421e10b8SBarry Smith   bmat->free_ij      = PETSC_TRUE;
548421e10b8SBarry Smith 
549421e10b8SBarry Smith   bmat->nz                = 0;
550421e10b8SBarry Smith   bmat->maxnz             = nz;
551421e10b8SBarry Smith   A->info.nz_unneeded  = (double)bmat->maxnz;
552421e10b8SBarry Smith 
553421e10b8SBarry Smith   PetscFunctionReturn(0);
554421e10b8SBarry Smith }
555421e10b8SBarry Smith 
556ccb205f8SBarry Smith /*
557ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
558ccb205f8SBarry Smith */
559ccb205f8SBarry Smith #undef __FUNCT__
560ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat"
561ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
562ccb205f8SBarry Smith {
563ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
564ccb205f8SBarry Smith   PetscErrorCode ierr;
565ccb205f8SBarry Smith   PetscInt       i,j,mbs = A->rmap.n/A->rmap.bs;
566ccb205f8SBarry Smith 
567ccb205f8SBarry Smith   PetscFunctionBegin;
568ccb205f8SBarry Smith   if (!a->diag) {
569ccb205f8SBarry Smith     ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
570ccb205f8SBarry Smith   }
571ccb205f8SBarry Smith   for (i=0; i<mbs; i++) {
572ccb205f8SBarry Smith     a->diag[i] = a->i[i+1];
573ccb205f8SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
574ccb205f8SBarry Smith       if (a->j[j] == i) {
575ccb205f8SBarry Smith         a->diag[i] = j;
576ccb205f8SBarry Smith         break;
577ccb205f8SBarry Smith       }
578ccb205f8SBarry Smith     }
579ccb205f8SBarry Smith   }
580ccb205f8SBarry Smith   PetscFunctionReturn(0);
581ccb205f8SBarry Smith }
582ccb205f8SBarry Smith 
583ccb205f8SBarry Smith #undef __FUNCT__
584ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat"
585ccb205f8SBarry Smith PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
586ccb205f8SBarry Smith {
587ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
588ccb205f8SBarry Smith   Mat_SeqAIJ     *c;
589ccb205f8SBarry Smith   PetscErrorCode ierr;
590ccb205f8SBarry Smith   PetscInt       i,k,first,step,lensi,nrows,ncols;
591ccb205f8SBarry Smith   PetscInt       *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
592ccb205f8SBarry Smith   PetscScalar    *a_new,value;
593ccb205f8SBarry Smith   Mat            C,*aa = a->a;
594ccb205f8SBarry Smith   PetscTruth     stride,equal;
595ccb205f8SBarry Smith 
596ccb205f8SBarry Smith   PetscFunctionBegin;
597ccb205f8SBarry Smith   ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr);
598ccb205f8SBarry Smith   if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices");
599ccb205f8SBarry Smith   ierr = ISStride(iscol,&stride);CHKERRQ(ierr);
600ccb205f8SBarry Smith   if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices");
601ccb205f8SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
602ccb205f8SBarry Smith   if (step != A->rmap.bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block");
603ccb205f8SBarry Smith 
604ccb205f8SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
605ccb205f8SBarry Smith   ncols = nrows;
606ccb205f8SBarry Smith 
607ccb205f8SBarry Smith   /* create submatrix */
608ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
609ccb205f8SBarry Smith     PetscInt n_cols,n_rows;
610ccb205f8SBarry Smith     C = *B;
611ccb205f8SBarry Smith     ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr);
612ccb205f8SBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
613ccb205f8SBarry Smith     ierr = MatZeroEntries(C);CHKERRQ(ierr);
614ccb205f8SBarry Smith   } else {
615ccb205f8SBarry Smith     ierr = MatCreate(A->comm,&C);CHKERRQ(ierr);
616ccb205f8SBarry Smith     ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
617ccb205f8SBarry Smith     ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr);
618ccb205f8SBarry Smith     ierr = MatSeqAIJSetPreallocation_SeqAIJ(C,0,ailen);CHKERRQ(ierr);
619ccb205f8SBarry Smith   }
620ccb205f8SBarry Smith   c = (Mat_SeqAIJ*)C->data;
621ccb205f8SBarry Smith 
622ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
623ccb205f8SBarry Smith   a_new    = c->a;
624ccb205f8SBarry Smith   j_new    = c->j;
625ccb205f8SBarry Smith   i_new    = c->i;
626ccb205f8SBarry Smith 
627ccb205f8SBarry Smith   for (i=0; i<nrows; i++) {
628ccb205f8SBarry Smith     ii    = ai[i];
629ccb205f8SBarry Smith     lensi = ailen[i];
630ccb205f8SBarry Smith     for (k=0; k<lensi; k++) {
631ccb205f8SBarry Smith       *j_new++ = *aj++;
632ccb205f8SBarry Smith       ierr     = MatGetValue(*aa++,first,first,value);CHKERRQ(ierr);
633ccb205f8SBarry Smith       *a_new++ = value;
634ccb205f8SBarry Smith     }
635ccb205f8SBarry Smith     i_new[i+1]  = i_new[i] + lensi;
636ccb205f8SBarry Smith     c->ilen[i]  = lensi;
637ccb205f8SBarry Smith   }
638ccb205f8SBarry Smith 
639ccb205f8SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
640ccb205f8SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
641ccb205f8SBarry Smith   *B = C;
642ccb205f8SBarry Smith   PetscFunctionReturn(0);
643ccb205f8SBarry Smith }
644ccb205f8SBarry Smith 
645421e10b8SBarry Smith #undef __FUNCT__
646421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat"
647421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
648421e10b8SBarry Smith {
649421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
650421e10b8SBarry Smith   PetscErrorCode ierr;
651421e10b8SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
652ccb205f8SBarry Smith   PetscInt       m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
653421e10b8SBarry Smith   Mat            *aa = a->a,*ap;
654421e10b8SBarry Smith 
655421e10b8SBarry Smith   PetscFunctionBegin;
656421e10b8SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
657421e10b8SBarry Smith 
658421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
659421e10b8SBarry Smith   for (i=1; i<m; i++) {
660421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
661421e10b8SBarry Smith     fshift += imax[i-1] - ailen[i-1];
662421e10b8SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
663421e10b8SBarry Smith     if (fshift) {
664421e10b8SBarry Smith       ip = aj + ai[i] ;
665421e10b8SBarry Smith       ap = aa + ai[i] ;
666421e10b8SBarry Smith       N  = ailen[i];
667421e10b8SBarry Smith       for (j=0; j<N; j++) {
668421e10b8SBarry Smith         ip[j-fshift] = ip[j];
669421e10b8SBarry Smith         ap[j-fshift] = ap[j];
670421e10b8SBarry Smith       }
671421e10b8SBarry Smith     }
672421e10b8SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
673421e10b8SBarry Smith   }
674421e10b8SBarry Smith   if (m) {
675421e10b8SBarry Smith     fshift += imax[m-1] - ailen[m-1];
676421e10b8SBarry Smith     ai[m]  = ai[m-1] + ailen[m-1];
677421e10b8SBarry Smith   }
678421e10b8SBarry Smith   /* reset ilen and imax for each row */
679421e10b8SBarry Smith   for (i=0; i<m; i++) {
680421e10b8SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
681421e10b8SBarry Smith   }
682421e10b8SBarry Smith   a->nz = ai[m];
683ccb205f8SBarry Smith   for (i=0; i<a->nz; i++) {
684ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG)
685ccb205f8SBarry Smith     if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz);
686ccb205f8SBarry Smith #endif
687ccb205f8SBarry Smith     ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
688ccb205f8SBarry Smith     ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
689ccb205f8SBarry Smith   }
690ccb205f8SBarry Smith   CHKMEMQ;
691421e10b8SBarry 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);
692421e10b8SBarry Smith   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr);
693421e10b8SBarry Smith   ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr);
694421e10b8SBarry Smith   a->reallocs          = 0;
695421e10b8SBarry Smith   A->info.nz_unneeded  = (double)fshift;
696421e10b8SBarry Smith   a->rmax              = rmax;
697421e10b8SBarry Smith 
698421e10b8SBarry Smith   A->same_nonzero = PETSC_TRUE;
699ccb205f8SBarry Smith   ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr);
700421e10b8SBarry Smith   PetscFunctionReturn(0);
701421e10b8SBarry Smith }
702421e10b8SBarry Smith 
703290bbb0aSBarry Smith #undef __FUNCT__
704290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat"
705290bbb0aSBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt)
706290bbb0aSBarry Smith {
707290bbb0aSBarry Smith   PetscFunctionBegin;
708290bbb0aSBarry Smith   if (opt == MAT_SYMMETRIC) {
709290bbb0aSBarry Smith     A->ops->relax = MatRelax_BlockMat_Symmetric;
710290bbb0aSBarry Smith     A->ops->mult  = MatMult_BlockMat_Symmetric;
711290bbb0aSBarry Smith   } else {
712290bbb0aSBarry Smith     PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]);
713290bbb0aSBarry Smith   }
714290bbb0aSBarry Smith   PetscFunctionReturn(0);
715290bbb0aSBarry Smith }
716290bbb0aSBarry Smith 
717290bbb0aSBarry Smith 
718421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
719421e10b8SBarry Smith        0,
720421e10b8SBarry Smith        0,
721421e10b8SBarry Smith        MatMult_BlockMat,
722421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat,
723421e10b8SBarry Smith        MatMultTranspose_BlockMat,
724421e10b8SBarry Smith        MatMultTransposeAdd_BlockMat,
725421e10b8SBarry Smith        0,
726421e10b8SBarry Smith        0,
727421e10b8SBarry Smith        0,
728421e10b8SBarry Smith /*10*/ 0,
729421e10b8SBarry Smith        0,
730421e10b8SBarry Smith        0,
731ccb205f8SBarry Smith        MatRelax_BlockMat,
732421e10b8SBarry Smith        0,
733421e10b8SBarry Smith /*15*/ 0,
734421e10b8SBarry Smith        0,
735421e10b8SBarry Smith        0,
736421e10b8SBarry Smith        0,
737421e10b8SBarry Smith        0,
738421e10b8SBarry Smith /*20*/ 0,
739421e10b8SBarry Smith        MatAssemblyEnd_BlockMat,
740421e10b8SBarry Smith        0,
741290bbb0aSBarry Smith        MatSetOption_BlockMat,
742421e10b8SBarry Smith        0,
743421e10b8SBarry Smith /*25*/ 0,
744421e10b8SBarry Smith        0,
745421e10b8SBarry Smith        0,
746421e10b8SBarry Smith        0,
747421e10b8SBarry Smith        0,
748421e10b8SBarry Smith /*30*/ 0,
749421e10b8SBarry Smith        0,
750421e10b8SBarry Smith        0,
751421e10b8SBarry Smith        0,
752421e10b8SBarry Smith        0,
753421e10b8SBarry Smith /*35*/ 0,
754421e10b8SBarry Smith        0,
755421e10b8SBarry Smith        0,
756421e10b8SBarry Smith        0,
757421e10b8SBarry Smith        0,
758421e10b8SBarry Smith /*40*/ 0,
759421e10b8SBarry Smith        0,
760421e10b8SBarry Smith        0,
761421e10b8SBarry Smith        0,
762421e10b8SBarry Smith        0,
763421e10b8SBarry Smith /*45*/ 0,
764421e10b8SBarry Smith        0,
765421e10b8SBarry Smith        0,
766421e10b8SBarry Smith        0,
767421e10b8SBarry Smith        0,
768421e10b8SBarry Smith /*50*/ MatSetBlockSize_BlockMat,
769421e10b8SBarry Smith        0,
770421e10b8SBarry Smith        0,
771421e10b8SBarry Smith        0,
772421e10b8SBarry Smith        0,
773421e10b8SBarry Smith /*55*/ 0,
774421e10b8SBarry Smith        0,
775421e10b8SBarry Smith        0,
776421e10b8SBarry Smith        0,
777421e10b8SBarry Smith        0,
778ccb205f8SBarry Smith /*60*/ MatGetSubMatrix_BlockMat,
779421e10b8SBarry Smith        MatDestroy_BlockMat,
780ccb205f8SBarry Smith        MatView_BlockMat,
781421e10b8SBarry Smith        0,
782421e10b8SBarry Smith        0,
783421e10b8SBarry Smith /*65*/ 0,
784421e10b8SBarry Smith        0,
785421e10b8SBarry Smith        0,
786421e10b8SBarry Smith        0,
787421e10b8SBarry Smith        0,
788421e10b8SBarry Smith /*70*/ 0,
789421e10b8SBarry Smith        0,
790421e10b8SBarry Smith        0,
791421e10b8SBarry Smith        0,
792421e10b8SBarry Smith        0,
793421e10b8SBarry Smith /*75*/ 0,
794421e10b8SBarry Smith        0,
795421e10b8SBarry Smith        0,
796421e10b8SBarry Smith        0,
797421e10b8SBarry Smith        0,
798421e10b8SBarry Smith /*80*/ 0,
799421e10b8SBarry Smith        0,
800421e10b8SBarry Smith        0,
801421e10b8SBarry Smith        0,
802421e10b8SBarry Smith        MatLoad_BlockMat,
803421e10b8SBarry Smith /*85*/ 0,
804421e10b8SBarry Smith        0,
805421e10b8SBarry Smith        0,
806421e10b8SBarry Smith        0,
807421e10b8SBarry Smith        0,
808421e10b8SBarry Smith /*90*/ 0,
809421e10b8SBarry Smith        0,
810421e10b8SBarry Smith        0,
811421e10b8SBarry Smith        0,
812421e10b8SBarry Smith        0,
813421e10b8SBarry Smith /*95*/ 0,
814421e10b8SBarry Smith        0,
815421e10b8SBarry Smith        0,
816421e10b8SBarry Smith        0};
817421e10b8SBarry Smith 
818421e10b8SBarry Smith /*MC
819421e10b8SBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
820421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
821421e10b8SBarry Smith 
822421e10b8SBarry Smith   Level: advanced
823421e10b8SBarry Smith 
824421e10b8SBarry Smith .seealso: MatCreateBlockMat()
825421e10b8SBarry Smith 
826421e10b8SBarry Smith M*/
827421e10b8SBarry Smith 
828421e10b8SBarry Smith EXTERN_C_BEGIN
829421e10b8SBarry Smith #undef __FUNCT__
830421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat"
831421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A)
832421e10b8SBarry Smith {
833421e10b8SBarry Smith   Mat_BlockMat   *b;
834421e10b8SBarry Smith   PetscErrorCode ierr;
835421e10b8SBarry Smith 
836421e10b8SBarry Smith   PetscFunctionBegin;
837421e10b8SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
838421e10b8SBarry Smith   ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr);
839421e10b8SBarry Smith 
840421e10b8SBarry Smith   A->data = (void*)b;
841421e10b8SBarry Smith 
842421e10b8SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
843421e10b8SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
844421e10b8SBarry Smith 
845421e10b8SBarry Smith   A->assembled     = PETSC_TRUE;
846421e10b8SBarry Smith   A->preallocated  = PETSC_FALSE;
847421e10b8SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
848290bbb0aSBarry Smith 
849290bbb0aSBarry Smith   ierr = PetscOptionsBegin(A->comm,A->prefix,"Matrix Option","Mat");CHKERRQ(ierr);
850290bbb0aSBarry Smith   ierr = PetscOptionsEnd();
851290bbb0aSBarry Smith 
852421e10b8SBarry Smith   PetscFunctionReturn(0);
853421e10b8SBarry Smith }
854421e10b8SBarry Smith EXTERN_C_END
855421e10b8SBarry Smith 
856421e10b8SBarry Smith #undef __FUNCT__
857421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat"
858421e10b8SBarry Smith /*@C
859421e10b8SBarry Smith    MatCreateBlockMat - Creates a new matrix based sparse Mat storage
860421e10b8SBarry Smith 
861421e10b8SBarry Smith   Collective on MPI_Comm
862421e10b8SBarry Smith 
863421e10b8SBarry Smith    Input Parameters:
864421e10b8SBarry Smith +  comm - MPI communicator
865421e10b8SBarry Smith .  m - number of rows
866421e10b8SBarry Smith .  n  - number of columns
867421e10b8SBarry Smith -  bs - size of each submatrix
868421e10b8SBarry Smith 
869421e10b8SBarry Smith 
870421e10b8SBarry Smith    Output Parameter:
871421e10b8SBarry Smith .  A - the matrix
872421e10b8SBarry Smith 
873421e10b8SBarry Smith    Level: intermediate
874421e10b8SBarry Smith 
875421e10b8SBarry Smith    PETSc requires that matrices and vectors being used for certain
876421e10b8SBarry Smith    operations are partitioned accordingly.  For example, when
877421e10b8SBarry Smith    creating a bmat matrix, A, that supports parallel matrix-vector
878421e10b8SBarry Smith    products using MatMult(A,x,y) the user should set the number
879421e10b8SBarry Smith    of local matrix rows to be the number of local elements of the
880421e10b8SBarry Smith    corresponding result vector, y. Note that this is information is
881421e10b8SBarry Smith    required for use of the matrix interface routines, even though
882421e10b8SBarry Smith    the bmat matrix may not actually be physically partitioned.
883421e10b8SBarry Smith    For example,
884421e10b8SBarry Smith 
885421e10b8SBarry Smith .keywords: matrix, bmat, create
886421e10b8SBarry Smith 
887421e10b8SBarry Smith .seealso: MATBLOCKMAT
888421e10b8SBarry Smith @*/
889421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,Mat *A)
890421e10b8SBarry Smith {
891421e10b8SBarry Smith   PetscErrorCode ierr;
892421e10b8SBarry Smith 
893421e10b8SBarry Smith   PetscFunctionBegin;
894421e10b8SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
895421e10b8SBarry Smith   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
896421e10b8SBarry Smith   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
897421e10b8SBarry Smith   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
898421e10b8SBarry Smith   PetscFunctionReturn(0);
899421e10b8SBarry Smith }
900421e10b8SBarry Smith 
901421e10b8SBarry Smith 
902421e10b8SBarry Smith 
903