xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision c98fd787ce49cfe1bba58d59168f91d3fdf7ec65)
1421e10b8SBarry Smith 
2421e10b8SBarry Smith /*
3421e10b8SBarry Smith    This provides a matrix that consists of Mats
4421e10b8SBarry Smith */
5421e10b8SBarry Smith 
6b45d2f2cSJed Brown #include <petsc-private/matimpl.h>              /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>    /* use the common AIJ data-structure */
8421e10b8SBarry Smith 
9421e10b8SBarry Smith typedef struct {
10421e10b8SBarry Smith   SEQAIJHEADER(Mat);
11421e10b8SBarry Smith   SEQBAIJHEADER;
12ccb205f8SBarry Smith   Mat *diags;
13421e10b8SBarry Smith 
146e63c7a1SBarry Smith   Vec left,right,middle,workb;                 /* dummy vectors to perform local parts of product */
15421e10b8SBarry Smith } Mat_BlockMat;
16421e10b8SBarry Smith 
177087cfbeSBarry Smith extern PetscErrorCode  MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt*);
18a5973382SShri Abhyankar 
19421e10b8SBarry Smith #undef __FUNCT__
2041f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat_Symmetric"
2141f059aeSBarry Smith PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
22290bbb0aSBarry Smith {
23290bbb0aSBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
24290bbb0aSBarry Smith   PetscScalar       *x;
25a2ea699eSBarry Smith   const Mat         *v;
26290bbb0aSBarry Smith   const PetscScalar *b;
27290bbb0aSBarry Smith   PetscErrorCode    ierr;
28d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
29290bbb0aSBarry Smith   const PetscInt    *idx;
30290bbb0aSBarry Smith   IS                row,col;
31290bbb0aSBarry Smith   MatFactorInfo     info;
326e63c7a1SBarry Smith   Vec               left = a->left,right = a->right, middle = a->middle;
33290bbb0aSBarry Smith   Mat               *diag;
34290bbb0aSBarry Smith 
35290bbb0aSBarry Smith   PetscFunctionBegin;
36290bbb0aSBarry Smith   its = its*lits;
37e32f2f54SBarry Smith   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
38e32f2f54SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
39e32f2f54SBarry Smith   if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
40e32f2f54SBarry Smith   if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
4126fbe8dcSKarl Rupp   if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) {
42e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
4326fbe8dcSKarl Rupp   }
44290bbb0aSBarry Smith 
45290bbb0aSBarry Smith   if (!a->diags) {
46785e854fSJed Brown     ierr = PetscMalloc1(mbs,&a->diags);CHKERRQ(ierr);
47290bbb0aSBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
48290bbb0aSBarry Smith     for (i=0; i<mbs; i++) {
492692d6eeSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr);
50719d5645SBarry Smith       ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr);
51719d5645SBarry Smith       ierr = MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr);
526bf464f9SBarry Smith       ierr = ISDestroy(&row);CHKERRQ(ierr);
536bf464f9SBarry Smith       ierr = ISDestroy(&col);CHKERRQ(ierr);
54290bbb0aSBarry Smith     }
556e63c7a1SBarry Smith     ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr);
56290bbb0aSBarry Smith   }
57290bbb0aSBarry Smith   diag = a->diags;
58290bbb0aSBarry Smith 
59290bbb0aSBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
60290bbb0aSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
61290bbb0aSBarry Smith   /* copy right hand side because it must be modified during iteration */
626e63c7a1SBarry Smith   ierr = VecCopy(bb,a->workb);CHKERRQ(ierr);
633649974fSBarry Smith   ierr = VecGetArrayRead(a->workb,&b);CHKERRQ(ierr);
64290bbb0aSBarry Smith 
6541f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
66290bbb0aSBarry Smith   while (its--) {
67290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
68290bbb0aSBarry Smith 
69290bbb0aSBarry Smith       for (i=0; i<mbs; i++) {
706e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
716e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
726e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
73290bbb0aSBarry Smith 
74290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
75290bbb0aSBarry Smith         for (j=0; j<n; j++) {
76290bbb0aSBarry Smith           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
77290bbb0aSBarry Smith           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
78290bbb0aSBarry Smith           ierr = VecResetArray(right);CHKERRQ(ierr);
79290bbb0aSBarry Smith         }
80290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
81290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
82290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
83290bbb0aSBarry Smith 
84290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
85290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
866e63c7a1SBarry Smith 
8741f059aeSBarry Smith         /* now adjust right hand side, see MatSOR_SeqSBAIJ */
886e63c7a1SBarry Smith         for (j=0; j<n; j++) {
896e63c7a1SBarry Smith           ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr);
906e63c7a1SBarry Smith           ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr);
916e63c7a1SBarry Smith           ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr);
926e63c7a1SBarry Smith           ierr = VecResetArray(middle);CHKERRQ(ierr);
936e63c7a1SBarry Smith         }
94290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
956e63c7a1SBarry Smith 
96290bbb0aSBarry Smith       }
97290bbb0aSBarry Smith     }
98290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
99290bbb0aSBarry Smith 
100290bbb0aSBarry Smith       for (i=mbs-1; i>=0; i--) {
1016e63c7a1SBarry Smith         n   = a->i[i+1] - a->i[i] - 1;
1026e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
1036e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
104290bbb0aSBarry Smith 
105290bbb0aSBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
106290bbb0aSBarry Smith         for (j=0; j<n; j++) {
107290bbb0aSBarry Smith           ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
108290bbb0aSBarry Smith           ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
109290bbb0aSBarry Smith           ierr = VecResetArray(right);CHKERRQ(ierr);
110290bbb0aSBarry Smith         }
111290bbb0aSBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
112290bbb0aSBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
113290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
114290bbb0aSBarry Smith 
115290bbb0aSBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
116290bbb0aSBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
117290bbb0aSBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
118290bbb0aSBarry Smith 
119290bbb0aSBarry Smith       }
120290bbb0aSBarry Smith     }
121290bbb0aSBarry Smith   }
122290bbb0aSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1233649974fSBarry Smith   ierr = VecRestoreArrayRead(a->workb,&b);CHKERRQ(ierr);
124290bbb0aSBarry Smith   PetscFunctionReturn(0);
125290bbb0aSBarry Smith }
126290bbb0aSBarry Smith 
127290bbb0aSBarry Smith #undef __FUNCT__
12841f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat"
12941f059aeSBarry Smith PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
130ccb205f8SBarry Smith {
131ccb205f8SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
132ccb205f8SBarry Smith   PetscScalar       *x;
133a2ea699eSBarry Smith   const Mat         *v;
134ccb205f8SBarry Smith   const PetscScalar *b;
135ccb205f8SBarry Smith   PetscErrorCode    ierr;
136d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
137ccb205f8SBarry Smith   const PetscInt    *idx;
138ccb205f8SBarry Smith   IS                row,col;
139ccb205f8SBarry Smith   MatFactorInfo     info;
140ccb205f8SBarry Smith   Vec               left = a->left,right = a->right;
141ccb205f8SBarry Smith   Mat               *diag;
142ccb205f8SBarry Smith 
143ccb205f8SBarry Smith   PetscFunctionBegin;
144ccb205f8SBarry Smith   its = its*lits;
145e32f2f54SBarry Smith   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
146e32f2f54SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
147e32f2f54SBarry Smith   if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
148e32f2f54SBarry Smith   if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
149ccb205f8SBarry Smith 
150ccb205f8SBarry Smith   if (!a->diags) {
151785e854fSJed Brown     ierr = PetscMalloc1(mbs,&a->diags);CHKERRQ(ierr);
152ccb205f8SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
153ccb205f8SBarry Smith     for (i=0; i<mbs; i++) {
1542692d6eeSBarry Smith       ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr);
155719d5645SBarry Smith       ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr);
156719d5645SBarry Smith       ierr = MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr);
1576bf464f9SBarry Smith       ierr = ISDestroy(&row);CHKERRQ(ierr);
1586bf464f9SBarry Smith       ierr = ISDestroy(&col);CHKERRQ(ierr);
159ccb205f8SBarry Smith     }
160ccb205f8SBarry Smith   }
161ccb205f8SBarry Smith   diag = a->diags;
162ccb205f8SBarry Smith 
163ccb205f8SBarry Smith   ierr = VecSet(xx,0.0);CHKERRQ(ierr);
164ccb205f8SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1653649974fSBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
166ccb205f8SBarry Smith 
16741f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
168ccb205f8SBarry Smith   while (its--) {
169ccb205f8SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
170ccb205f8SBarry Smith 
171ccb205f8SBarry Smith       for (i=0; i<mbs; i++) {
172ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
173ccb205f8SBarry Smith         idx = a->j + a->i[i];
174ccb205f8SBarry Smith         v   = a->a + a->i[i];
175ccb205f8SBarry Smith 
176ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
177ccb205f8SBarry Smith         for (j=0; j<n; j++) {
178ccb205f8SBarry Smith           if (idx[j] != i) {
179ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
180ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
181ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
182ccb205f8SBarry Smith           }
183ccb205f8SBarry Smith         }
184ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
185ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
186ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
187ccb205f8SBarry Smith 
188ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
189ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
190ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
191ccb205f8SBarry Smith       }
192ccb205f8SBarry Smith     }
193ccb205f8SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
194ccb205f8SBarry Smith 
195ccb205f8SBarry Smith       for (i=mbs-1; i>=0; i--) {
196ccb205f8SBarry Smith         n   = a->i[i+1] - a->i[i];
197ccb205f8SBarry Smith         idx = a->j + a->i[i];
198ccb205f8SBarry Smith         v   = a->a + a->i[i];
199ccb205f8SBarry Smith 
200ccb205f8SBarry Smith         ierr = VecSet(left,0.0);CHKERRQ(ierr);
201ccb205f8SBarry Smith         for (j=0; j<n; j++) {
202ccb205f8SBarry Smith           if (idx[j] != i) {
203ccb205f8SBarry Smith             ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr);
204ccb205f8SBarry Smith             ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr);
205ccb205f8SBarry Smith             ierr = VecResetArray(right);CHKERRQ(ierr);
206ccb205f8SBarry Smith           }
207ccb205f8SBarry Smith         }
208ccb205f8SBarry Smith         ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr);
209ccb205f8SBarry Smith         ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr);
210ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
211ccb205f8SBarry Smith 
212ccb205f8SBarry Smith         ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr);
213ccb205f8SBarry Smith         ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr);
214ccb205f8SBarry Smith         ierr = VecResetArray(right);CHKERRQ(ierr);
215ccb205f8SBarry Smith 
216ccb205f8SBarry Smith       }
217ccb205f8SBarry Smith     }
218ccb205f8SBarry Smith   }
219ccb205f8SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2203649974fSBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
221ccb205f8SBarry Smith   PetscFunctionReturn(0);
222ccb205f8SBarry Smith }
223ccb205f8SBarry Smith 
224ccb205f8SBarry Smith #undef __FUNCT__
225421e10b8SBarry Smith #define __FUNCT__ "MatSetValues_BlockMat"
226421e10b8SBarry Smith PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
227421e10b8SBarry Smith {
228421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
229421e10b8SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
230421e10b8SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
231d0f46423SBarry Smith   PetscInt       *aj  =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
232421e10b8SBarry Smith   PetscErrorCode ierr;
233421e10b8SBarry Smith   PetscInt       ridx,cidx;
234ace3abfcSBarry Smith   PetscBool      roworiented=a->roworiented;
235421e10b8SBarry Smith   MatScalar      value;
236421e10b8SBarry Smith   Mat            *ap,*aa = a->a;
237421e10b8SBarry Smith 
238421e10b8SBarry Smith   PetscFunctionBegin;
239421e10b8SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
240421e10b8SBarry Smith     row  = im[k];
241421e10b8SBarry Smith     brow = row/bs;
242421e10b8SBarry Smith     if (row < 0) continue;
243421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
244e32f2f54SBarry Smith     if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
245421e10b8SBarry Smith #endif
246421e10b8SBarry Smith     rp   = aj + ai[brow];
247421e10b8SBarry Smith     ap   = aa + ai[brow];
248421e10b8SBarry Smith     rmax = imax[brow];
249421e10b8SBarry Smith     nrow = ailen[brow];
250421e10b8SBarry Smith     low  = 0;
251421e10b8SBarry Smith     high = nrow;
252421e10b8SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
253421e10b8SBarry Smith       if (in[l] < 0) continue;
254421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
255e32f2f54SBarry Smith       if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1);
256421e10b8SBarry Smith #endif
257421e10b8SBarry Smith       col = in[l]; bcol = col/bs;
2586e63c7a1SBarry Smith       if (A->symmetric && brow > bcol) continue;
259421e10b8SBarry Smith       ridx = row % bs; cidx = col % bs;
2602205254eSKarl Rupp       if (roworiented) value = v[l + k*n];
2612205254eSKarl Rupp       else value = v[k + l*m];
2622205254eSKarl Rupp 
2632205254eSKarl Rupp       if (col <= lastcol) low = 0;
2642205254eSKarl Rupp       else high = nrow;
265421e10b8SBarry Smith       lastcol = col;
266421e10b8SBarry Smith       while (high-low > 7) {
267421e10b8SBarry Smith         t = (low+high)/2;
268421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
269421e10b8SBarry Smith         else              low  = t;
270421e10b8SBarry Smith       }
271421e10b8SBarry Smith       for (i=low; i<high; i++) {
272421e10b8SBarry Smith         if (rp[i] > bcol) break;
2732205254eSKarl Rupp         if (rp[i] == bcol) goto noinsert1;
274421e10b8SBarry Smith       }
275421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
276e32f2f54SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
277fef13f97SBarry 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 */
280421e10b8SBarry Smith       for (ii=N; ii>=i; ii--) {
281421e10b8SBarry Smith         rp[ii+1] = rp[ii];
282421e10b8SBarry Smith         ap[ii+1] = ap[ii];
283421e10b8SBarry Smith       }
284421e10b8SBarry Smith       if (N>=i) ap[i] = 0;
285421e10b8SBarry Smith       rp[i] = bcol;
286421e10b8SBarry Smith       a->nz++;
287e56f5c9eSBarry Smith       A->nonzerostate++;
288421e10b8SBarry Smith noinsert1:;
289421e10b8SBarry Smith       if (!*(ap+i)) {
2906e63c7a1SBarry Smith         ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
291b0223207SBarry Smith       }
292421e10b8SBarry Smith       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
293421e10b8SBarry Smith       low  = i;
294421e10b8SBarry Smith     }
295421e10b8SBarry Smith     ailen[brow] = nrow;
296421e10b8SBarry Smith   }
297421e10b8SBarry Smith   PetscFunctionReturn(0);
298421e10b8SBarry Smith }
299421e10b8SBarry Smith 
300421e10b8SBarry Smith #undef __FUNCT__
3015bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_BlockMat"
302112444f4SShri Abhyankar PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
303a5973382SShri Abhyankar {
304a5973382SShri Abhyankar   PetscErrorCode    ierr;
305a5973382SShri Abhyankar   Mat               tmpA;
306a5973382SShri Abhyankar   PetscInt          i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0;
307a5973382SShri Abhyankar   const PetscInt    *cols;
308a5973382SShri Abhyankar   const PetscScalar *values;
309ace3abfcSBarry Smith   PetscBool         flg = PETSC_FALSE,notdone;
310a5973382SShri Abhyankar   Mat_SeqAIJ        *a;
311a5973382SShri Abhyankar   Mat_BlockMat      *amat;
312a5973382SShri Abhyankar 
313a5973382SShri Abhyankar   PetscFunctionBegin;
314*c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
315*c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
316a5973382SShri Abhyankar   ierr = MatCreate(PETSC_COMM_SELF,&tmpA);CHKERRQ(ierr);
317a5973382SShri Abhyankar   ierr = MatSetType(tmpA,MATSEQAIJ);CHKERRQ(ierr);
318112444f4SShri Abhyankar   ierr = MatLoad_SeqAIJ(tmpA,viewer);CHKERRQ(ierr);
319a5973382SShri Abhyankar 
320a5973382SShri Abhyankar   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
3210298fd71SBarry Smith   ierr = PetscOptionsBegin(PETSC_COMM_SELF,NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr);
3220298fd71SBarry Smith   ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);CHKERRQ(ierr);
3230298fd71SBarry Smith   ierr = PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,NULL);CHKERRQ(ierr);
324a5973382SShri Abhyankar   ierr = PetscOptionsEnd();CHKERRQ(ierr);
325a5973382SShri Abhyankar 
326a5973382SShri Abhyankar   /* Determine number of nonzero blocks for each block row */
327a5973382SShri Abhyankar   a    = (Mat_SeqAIJ*) tmpA->data;
328a5973382SShri Abhyankar   mbs  = m/bs;
329dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&lens,bs,&ii,bs,&ilens);CHKERRQ(ierr);
330a5973382SShri Abhyankar   ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr);
331a5973382SShri Abhyankar 
332a5973382SShri Abhyankar   for (i=0; i<mbs; i++) {
333a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
334a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
335a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
336a5973382SShri Abhyankar     }
337a5973382SShri Abhyankar 
338a5973382SShri Abhyankar     currentcol = -1;
339a5973382SShri Abhyankar     notdone    = PETSC_TRUE;
340a5973382SShri Abhyankar     while (PETSC_TRUE) {
341a5973382SShri Abhyankar       notdone = PETSC_FALSE;
342a5973382SShri Abhyankar       nextcol = 1000000000;
343a5973382SShri Abhyankar       for (j=0; j<bs; j++) {
344a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) {
345a5973382SShri Abhyankar           ii[j]++;
346a5973382SShri Abhyankar           ilens[j]--;
347a5973382SShri Abhyankar         }
348a5973382SShri Abhyankar         if (ilens[j] > 0) {
349a5973382SShri Abhyankar           notdone = PETSC_TRUE;
350a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
351a5973382SShri Abhyankar         }
352a5973382SShri Abhyankar       }
353a5973382SShri Abhyankar       if (!notdone) break;
354a5973382SShri Abhyankar       if (!flg || (nextcol >= i)) lens[i]++;
355a5973382SShri Abhyankar       currentcol = nextcol;
356a5973382SShri Abhyankar     }
357a5973382SShri Abhyankar   }
358a5973382SShri Abhyankar 
359a5973382SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
360a5973382SShri Abhyankar     ierr = MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
361a5973382SShri Abhyankar   }
362a5973382SShri Abhyankar   ierr = MatBlockMatSetPreallocation(newmat,bs,0,lens);CHKERRQ(ierr);
363a5973382SShri Abhyankar   if (flg) {
364a5973382SShri Abhyankar     ierr = MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
365a5973382SShri Abhyankar   }
366a5973382SShri Abhyankar   amat = (Mat_BlockMat*)(newmat)->data;
367a5973382SShri Abhyankar 
368a5973382SShri Abhyankar   /* preallocate the submatrices */
369785e854fSJed Brown   ierr = PetscMalloc1(bs,&llens);CHKERRQ(ierr);
370a5973382SShri Abhyankar   for (i=0; i<mbs; i++) { /* loops for block rows */
371a5973382SShri Abhyankar     for (j=0; j<bs; j++) {
372a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i*bs + j];
373a5973382SShri Abhyankar       ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
374a5973382SShri Abhyankar     }
375a5973382SShri Abhyankar 
376a5973382SShri Abhyankar     currentcol = 1000000000;
377a5973382SShri Abhyankar     for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */
378a5973382SShri Abhyankar       if (ilens[j] > 0) {
379a5973382SShri Abhyankar         currentcol = PetscMin(currentcol,ii[j][0]/bs);
380a5973382SShri Abhyankar       }
381a5973382SShri Abhyankar     }
382a5973382SShri Abhyankar 
383a5973382SShri Abhyankar     notdone = PETSC_TRUE;
384a5973382SShri Abhyankar     while (PETSC_TRUE) {  /* loops over blocks in block row */
385a5973382SShri Abhyankar 
386a5973382SShri Abhyankar       notdone = PETSC_FALSE;
387a5973382SShri Abhyankar       nextcol = 1000000000;
388a5973382SShri Abhyankar       ierr    = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr);
389a5973382SShri Abhyankar       for (j=0; j<bs; j++) { /* loop over rows in block */
390a5973382SShri Abhyankar         while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */
391a5973382SShri Abhyankar           ii[j]++;
392a5973382SShri Abhyankar           ilens[j]--;
393a5973382SShri Abhyankar           llens[j]++;
394a5973382SShri Abhyankar         }
395a5973382SShri Abhyankar         if (ilens[j] > 0) {
396a5973382SShri Abhyankar           notdone = PETSC_TRUE;
397a5973382SShri Abhyankar           nextcol = PetscMin(nextcol,ii[j][0]/bs);
398a5973382SShri Abhyankar         }
399a5973382SShri Abhyankar       }
400a5973382SShri Abhyankar       if (cnt >= amat->maxnz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt);
401a5973382SShri Abhyankar       if (!flg || currentcol >= i) {
402a5973382SShri Abhyankar         amat->j[cnt] = currentcol;
403a5973382SShri Abhyankar         ierr         = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr);
404a5973382SShri Abhyankar       }
405a5973382SShri Abhyankar 
406a5973382SShri Abhyankar       if (!notdone) break;
407a5973382SShri Abhyankar       currentcol = nextcol;
408a5973382SShri Abhyankar     }
409a5973382SShri Abhyankar     amat->ilen[i] = lens[i];
410a5973382SShri Abhyankar   }
411a5973382SShri Abhyankar 
412a5973382SShri Abhyankar   ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr);
413a5973382SShri Abhyankar   ierr = PetscFree(llens);CHKERRQ(ierr);
414a5973382SShri Abhyankar 
415a5973382SShri Abhyankar   /* copy over the matrix, one row at a time */
416a5973382SShri Abhyankar   for (i=0; i<m; i++) {
417a5973382SShri Abhyankar     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
418a5973382SShri Abhyankar     ierr = MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
419a5973382SShri Abhyankar     ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
420a5973382SShri Abhyankar   }
421a5973382SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
422a5973382SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
423a5973382SShri Abhyankar   PetscFunctionReturn(0);
424a5973382SShri Abhyankar }
425a5973382SShri Abhyankar 
426a5973382SShri Abhyankar #undef __FUNCT__
427ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat"
428ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
429ccb205f8SBarry Smith {
43064075487SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat*)A->data;
431ccb205f8SBarry Smith   PetscErrorCode    ierr;
432ccb205f8SBarry Smith   const char        *name;
433ccb205f8SBarry Smith   PetscViewerFormat format;
434ccb205f8SBarry Smith 
435ccb205f8SBarry Smith   PetscFunctionBegin;
436ccb205f8SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
437ccb205f8SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
438ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
439ccb205f8SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr);
44064075487SBarry Smith     if (A->symmetric) {
4418c1ad04dSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr);
44264075487SBarry Smith     }
443ccb205f8SBarry Smith   }
444ccb205f8SBarry Smith   PetscFunctionReturn(0);
445ccb205f8SBarry Smith }
446421e10b8SBarry Smith 
447421e10b8SBarry Smith #undef __FUNCT__
448421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat"
449421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat)
450421e10b8SBarry Smith {
451421e10b8SBarry Smith   PetscErrorCode ierr;
452421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
453ccb205f8SBarry Smith   PetscInt       i;
454421e10b8SBarry Smith 
455421e10b8SBarry Smith   PetscFunctionBegin;
4566bf464f9SBarry Smith   ierr = VecDestroy(&bmat->right);CHKERRQ(ierr);
4576bf464f9SBarry Smith   ierr = VecDestroy(&bmat->left);CHKERRQ(ierr);
4586bf464f9SBarry Smith   ierr = VecDestroy(&bmat->middle);CHKERRQ(ierr);
4596bf464f9SBarry Smith   ierr = VecDestroy(&bmat->workb);CHKERRQ(ierr);
460ccb205f8SBarry Smith   if (bmat->diags) {
461d0f46423SBarry Smith     for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) {
4626bf464f9SBarry Smith       ierr = MatDestroy(&bmat->diags[i]);CHKERRQ(ierr);
463ccb205f8SBarry Smith     }
464ccb205f8SBarry Smith   }
465ccb205f8SBarry Smith   if (bmat->a) {
466ccb205f8SBarry Smith     for (i=0; i<bmat->nz; i++) {
4676bf464f9SBarry Smith       ierr = MatDestroy(&bmat->a[i]);CHKERRQ(ierr);
468ccb205f8SBarry Smith     }
469ccb205f8SBarry Smith   }
470ccb205f8SBarry Smith   ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
471bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
472421e10b8SBarry Smith   PetscFunctionReturn(0);
473421e10b8SBarry Smith }
474421e10b8SBarry Smith 
475421e10b8SBarry Smith #undef __FUNCT__
476421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat"
477421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
478421e10b8SBarry Smith {
479421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
480421e10b8SBarry Smith   PetscErrorCode ierr;
481421e10b8SBarry Smith   PetscScalar    *xx,*yy;
482d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
483421e10b8SBarry Smith   Mat            *aa;
484421e10b8SBarry Smith 
485421e10b8SBarry Smith   PetscFunctionBegin;
486421e10b8SBarry Smith   /*
487421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
488421e10b8SBarry Smith   */
489421e10b8SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
490ccb205f8SBarry Smith 
491421e10b8SBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
492421e10b8SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
493421e10b8SBarry Smith   aj   = bmat->j;
494421e10b8SBarry Smith   aa   = bmat->a;
495421e10b8SBarry Smith   ii   = bmat->i;
496421e10b8SBarry Smith   for (i=0; i<m; i++) {
497421e10b8SBarry Smith     jrow = ii[i];
498ccb205f8SBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
499421e10b8SBarry Smith     n    = ii[i+1] - jrow;
500421e10b8SBarry Smith     for (j=0; j<n; j++) {
501421e10b8SBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
502ccb205f8SBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
503421e10b8SBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
504421e10b8SBarry Smith       jrow++;
505421e10b8SBarry Smith     }
506421e10b8SBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
507421e10b8SBarry Smith   }
508421e10b8SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
509421e10b8SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
510421e10b8SBarry Smith   PetscFunctionReturn(0);
511421e10b8SBarry Smith }
512421e10b8SBarry Smith 
513421e10b8SBarry Smith #undef __FUNCT__
514290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric"
515290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
516290bbb0aSBarry Smith {
517290bbb0aSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
518290bbb0aSBarry Smith   PetscErrorCode ierr;
519290bbb0aSBarry Smith   PetscScalar    *xx,*yy;
520d0f46423SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
521290bbb0aSBarry Smith   Mat            *aa;
522290bbb0aSBarry Smith 
523290bbb0aSBarry Smith   PetscFunctionBegin;
524290bbb0aSBarry Smith   /*
525290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
526290bbb0aSBarry Smith   */
527290bbb0aSBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
528290bbb0aSBarry Smith 
529290bbb0aSBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
530290bbb0aSBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
531290bbb0aSBarry Smith   aj   = bmat->j;
532290bbb0aSBarry Smith   aa   = bmat->a;
533290bbb0aSBarry Smith   ii   = bmat->i;
534290bbb0aSBarry Smith   for (i=0; i<m; i++) {
535290bbb0aSBarry Smith     jrow = ii[i];
536290bbb0aSBarry Smith     n    = ii[i+1] - jrow;
537290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr);
538290bbb0aSBarry Smith     ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr);
539290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
540290bbb0aSBarry Smith     if (aj[jrow] == i) {
541290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
542290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
543290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
544290bbb0aSBarry Smith       jrow++;
545290bbb0aSBarry Smith       n--;
546290bbb0aSBarry Smith     }
547290bbb0aSBarry Smith     for (j=0; j<n; j++) {
548290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);            /* upper triangular part */
549290bbb0aSBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr);
550290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
551290bbb0aSBarry Smith 
552290bbb0aSBarry Smith       ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr);            /* lower triangular part */
553290bbb0aSBarry Smith       ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr);
554290bbb0aSBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
555290bbb0aSBarry Smith       jrow++;
556290bbb0aSBarry Smith     }
557290bbb0aSBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
558290bbb0aSBarry Smith     ierr = VecResetArray(bmat->middle);CHKERRQ(ierr);
559290bbb0aSBarry Smith   }
560290bbb0aSBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
561290bbb0aSBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
562290bbb0aSBarry Smith   PetscFunctionReturn(0);
563290bbb0aSBarry Smith }
564290bbb0aSBarry Smith 
565290bbb0aSBarry Smith #undef __FUNCT__
566421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat"
567421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
568421e10b8SBarry Smith {
569421e10b8SBarry Smith   PetscFunctionBegin;
570421e10b8SBarry Smith   PetscFunctionReturn(0);
571421e10b8SBarry Smith }
572421e10b8SBarry Smith 
573421e10b8SBarry Smith #undef __FUNCT__
574421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat"
575421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
576421e10b8SBarry Smith {
577421e10b8SBarry Smith   PetscFunctionBegin;
578421e10b8SBarry Smith   PetscFunctionReturn(0);
579421e10b8SBarry Smith }
580421e10b8SBarry Smith 
581421e10b8SBarry Smith #undef __FUNCT__
582421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat"
583421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
584421e10b8SBarry Smith {
585421e10b8SBarry Smith   PetscFunctionBegin;
586421e10b8SBarry Smith   PetscFunctionReturn(0);
587421e10b8SBarry Smith }
588421e10b8SBarry Smith 
589ccb205f8SBarry Smith /*
590ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
591ccb205f8SBarry Smith */
592ccb205f8SBarry Smith #undef __FUNCT__
593ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat"
594ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
595ccb205f8SBarry Smith {
596ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
597ccb205f8SBarry Smith   PetscErrorCode ierr;
598d0f46423SBarry Smith   PetscInt       i,j,mbs = A->rmap->n/A->rmap->bs;
599ccb205f8SBarry Smith 
600ccb205f8SBarry Smith   PetscFunctionBegin;
601ccb205f8SBarry Smith   if (!a->diag) {
602785e854fSJed Brown     ierr = PetscMalloc1(mbs,&a->diag);CHKERRQ(ierr);
603ccb205f8SBarry Smith   }
604ccb205f8SBarry Smith   for (i=0; i<mbs; i++) {
605ccb205f8SBarry Smith     a->diag[i] = a->i[i+1];
606ccb205f8SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
607ccb205f8SBarry Smith       if (a->j[j] == i) {
608ccb205f8SBarry Smith         a->diag[i] = j;
609ccb205f8SBarry Smith         break;
610ccb205f8SBarry Smith       }
611ccb205f8SBarry Smith     }
612ccb205f8SBarry Smith   }
613ccb205f8SBarry Smith   PetscFunctionReturn(0);
614ccb205f8SBarry Smith }
615ccb205f8SBarry Smith 
616ccb205f8SBarry Smith #undef __FUNCT__
617ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat"
6184aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
619ccb205f8SBarry Smith {
620ccb205f8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
621ccb205f8SBarry Smith   Mat_SeqAIJ     *c;
622ccb205f8SBarry Smith   PetscErrorCode ierr;
623ccb205f8SBarry Smith   PetscInt       i,k,first,step,lensi,nrows,ncols;
624d2a212eaSBarry Smith   PetscInt       *j_new,*i_new,*aj = a->j,*ailen = a->ilen;
6251620fd73SBarry Smith   PetscScalar    *a_new;
626ccb205f8SBarry Smith   Mat            C,*aa = a->a;
627ace3abfcSBarry Smith   PetscBool      stride,equal;
628ccb205f8SBarry Smith 
629ccb205f8SBarry Smith   PetscFunctionBegin;
630ccb205f8SBarry Smith   ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr);
631e32f2f54SBarry Smith   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices");
632251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);CHKERRQ(ierr);
633e32f2f54SBarry Smith   if (!stride) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices");
634ccb205f8SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
635e32f2f54SBarry Smith   if (step != A->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block");
636ccb205f8SBarry Smith 
637ccb205f8SBarry Smith   ierr  = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
638ccb205f8SBarry Smith   ncols = nrows;
639ccb205f8SBarry Smith 
640ccb205f8SBarry Smith   /* create submatrix */
641ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
642ccb205f8SBarry Smith     PetscInt n_cols,n_rows;
643ccb205f8SBarry Smith     C    = *B;
644ccb205f8SBarry Smith     ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr);
645e32f2f54SBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
646ccb205f8SBarry Smith     ierr = MatZeroEntries(C);CHKERRQ(ierr);
647ccb205f8SBarry Smith   } else {
648ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
649ccb205f8SBarry Smith     ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
6506e63c7a1SBarry Smith     if (A->symmetric) {
6516e63c7a1SBarry Smith       ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
6526e63c7a1SBarry Smith     } else {
653ccb205f8SBarry Smith       ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr);
6546e63c7a1SBarry Smith     }
6556e63c7a1SBarry Smith     ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr);
6566e63c7a1SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr);
657ccb205f8SBarry Smith   }
658ccb205f8SBarry Smith   c = (Mat_SeqAIJ*)C->data;
659ccb205f8SBarry Smith 
660ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
661ccb205f8SBarry Smith   a_new = c->a;
662ccb205f8SBarry Smith   j_new = c->j;
663ccb205f8SBarry Smith   i_new = c->i;
664ccb205f8SBarry Smith 
665ccb205f8SBarry Smith   for (i=0; i<nrows; i++) {
666ccb205f8SBarry Smith     lensi = ailen[i];
667ccb205f8SBarry Smith     for (k=0; k<lensi; k++) {
668ccb205f8SBarry Smith       *j_new++ = *aj++;
6691620fd73SBarry Smith       ierr     = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr);
670ccb205f8SBarry Smith     }
671ccb205f8SBarry Smith     i_new[i+1] = i_new[i] + lensi;
672ccb205f8SBarry Smith     c->ilen[i] = lensi;
673ccb205f8SBarry Smith   }
674ccb205f8SBarry Smith 
675ccb205f8SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
676ccb205f8SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
677ccb205f8SBarry Smith   *B   = C;
678ccb205f8SBarry Smith   PetscFunctionReturn(0);
679ccb205f8SBarry Smith }
680ccb205f8SBarry Smith 
681421e10b8SBarry Smith #undef __FUNCT__
682421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat"
683421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
684421e10b8SBarry Smith {
685421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
686421e10b8SBarry Smith   PetscErrorCode ierr;
687421e10b8SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
688ccb205f8SBarry Smith   PetscInt       m      = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
689421e10b8SBarry Smith   Mat            *aa    = a->a,*ap;
690421e10b8SBarry Smith 
691421e10b8SBarry Smith   PetscFunctionBegin;
692421e10b8SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
693421e10b8SBarry Smith 
694421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
695421e10b8SBarry Smith   for (i=1; i<m; i++) {
696421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
697421e10b8SBarry Smith     fshift += imax[i-1] - ailen[i-1];
698421e10b8SBarry Smith     rmax    = PetscMax(rmax,ailen[i]);
699421e10b8SBarry Smith     if (fshift) {
700421e10b8SBarry Smith       ip = aj + ai[i];
701421e10b8SBarry Smith       ap = aa + ai[i];
702421e10b8SBarry Smith       N  = ailen[i];
703421e10b8SBarry Smith       for (j=0; j<N; j++) {
704421e10b8SBarry Smith         ip[j-fshift] = ip[j];
705421e10b8SBarry Smith         ap[j-fshift] = ap[j];
706421e10b8SBarry Smith       }
707421e10b8SBarry Smith     }
708421e10b8SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
709421e10b8SBarry Smith   }
710421e10b8SBarry Smith   if (m) {
711421e10b8SBarry Smith     fshift += imax[m-1] - ailen[m-1];
712421e10b8SBarry Smith     ai[m]   = ai[m-1] + ailen[m-1];
713421e10b8SBarry Smith   }
714421e10b8SBarry Smith   /* reset ilen and imax for each row */
715421e10b8SBarry Smith   for (i=0; i<m; i++) {
716421e10b8SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
717421e10b8SBarry Smith   }
718421e10b8SBarry Smith   a->nz = ai[m];
719ccb205f8SBarry Smith   for (i=0; i<a->nz; i++) {
720ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG)
721e32f2f54SBarry Smith     if (!aa[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz);
722ccb205f8SBarry Smith #endif
723ccb205f8SBarry Smith     ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
724ccb205f8SBarry Smith     ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
725ccb205f8SBarry Smith   }
726d0f46423SBarry Smith   ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap->n/A->cmap->bs,fshift,a->nz);CHKERRQ(ierr);
727421e10b8SBarry Smith   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr);
728421e10b8SBarry Smith   ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr);
7292205254eSKarl Rupp 
7308e58a170SBarry Smith   A->info.mallocs    += a->reallocs;
731421e10b8SBarry Smith   a->reallocs         = 0;
732421e10b8SBarry Smith   A->info.nz_unneeded = (double)fshift;
733421e10b8SBarry Smith   a->rmax             = rmax;
734ccb205f8SBarry Smith   ierr                = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr);
735421e10b8SBarry Smith   PetscFunctionReturn(0);
736421e10b8SBarry Smith }
737421e10b8SBarry Smith 
738290bbb0aSBarry Smith #undef __FUNCT__
739290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat"
740ace3abfcSBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg)
741290bbb0aSBarry Smith {
742290bbb0aSBarry Smith   PetscFunctionBegin;
7434e0d8c25SBarry Smith   if (opt == MAT_SYMMETRIC && flg) {
74441f059aeSBarry Smith     A->ops->sor  = MatSOR_BlockMat_Symmetric;
745290bbb0aSBarry Smith     A->ops->mult = MatMult_BlockMat_Symmetric;
746290bbb0aSBarry Smith   } else {
747290bbb0aSBarry Smith     PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]);
748290bbb0aSBarry Smith   }
749290bbb0aSBarry Smith   PetscFunctionReturn(0);
750290bbb0aSBarry Smith }
751290bbb0aSBarry Smith 
752290bbb0aSBarry Smith 
7533964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
754421e10b8SBarry Smith                                        0,
755421e10b8SBarry Smith                                        0,
756421e10b8SBarry Smith                                        MatMult_BlockMat,
757421e10b8SBarry Smith                                /*  4*/ MatMultAdd_BlockMat,
758421e10b8SBarry Smith                                        MatMultTranspose_BlockMat,
759421e10b8SBarry Smith                                        MatMultTransposeAdd_BlockMat,
760421e10b8SBarry Smith                                        0,
761421e10b8SBarry Smith                                        0,
762421e10b8SBarry Smith                                        0,
763421e10b8SBarry Smith                                /* 10*/ 0,
764421e10b8SBarry Smith                                        0,
765421e10b8SBarry Smith                                        0,
76641f059aeSBarry Smith                                        MatSOR_BlockMat,
767421e10b8SBarry Smith                                        0,
768421e10b8SBarry Smith                                /* 15*/ 0,
769421e10b8SBarry Smith                                        0,
770421e10b8SBarry Smith                                        0,
771421e10b8SBarry Smith                                        0,
772421e10b8SBarry Smith                                        0,
773421e10b8SBarry Smith                                /* 20*/ 0,
774421e10b8SBarry Smith                                        MatAssemblyEnd_BlockMat,
775290bbb0aSBarry Smith                                        MatSetOption_BlockMat,
776421e10b8SBarry Smith                                        0,
777d519adbfSMatthew Knepley                                /* 24*/ 0,
778421e10b8SBarry Smith                                        0,
779421e10b8SBarry Smith                                        0,
780421e10b8SBarry Smith                                        0,
781421e10b8SBarry Smith                                        0,
782d519adbfSMatthew Knepley                                /* 29*/ 0,
783421e10b8SBarry Smith                                        0,
784421e10b8SBarry Smith                                        0,
785421e10b8SBarry Smith                                        0,
786421e10b8SBarry Smith                                        0,
787d519adbfSMatthew Knepley                                /* 34*/ 0,
788421e10b8SBarry Smith                                        0,
789421e10b8SBarry Smith                                        0,
790421e10b8SBarry Smith                                        0,
791421e10b8SBarry Smith                                        0,
792d519adbfSMatthew Knepley                                /* 39*/ 0,
793421e10b8SBarry Smith                                        0,
794421e10b8SBarry Smith                                        0,
795421e10b8SBarry Smith                                        0,
796421e10b8SBarry Smith                                        0,
797d519adbfSMatthew Knepley                                /* 44*/ 0,
798421e10b8SBarry Smith                                        0,
799421e10b8SBarry Smith                                        0,
800421e10b8SBarry Smith                                        0,
801421e10b8SBarry Smith                                        0,
802d519adbfSMatthew Knepley                                /* 49*/ 0,
803421e10b8SBarry Smith                                        0,
804421e10b8SBarry Smith                                        0,
805421e10b8SBarry Smith                                        0,
806421e10b8SBarry Smith                                        0,
807d519adbfSMatthew Knepley                                /* 54*/ 0,
808421e10b8SBarry Smith                                        0,
809421e10b8SBarry Smith                                        0,
810421e10b8SBarry Smith                                        0,
811421e10b8SBarry Smith                                        0,
812d519adbfSMatthew Knepley                                /* 59*/ MatGetSubMatrix_BlockMat,
813421e10b8SBarry Smith                                        MatDestroy_BlockMat,
814ccb205f8SBarry Smith                                        MatView_BlockMat,
815421e10b8SBarry Smith                                        0,
816421e10b8SBarry Smith                                        0,
817d519adbfSMatthew Knepley                                /* 64*/ 0,
818421e10b8SBarry Smith                                        0,
819421e10b8SBarry Smith                                        0,
820421e10b8SBarry Smith                                        0,
821421e10b8SBarry Smith                                        0,
822d519adbfSMatthew Knepley                                /* 69*/ 0,
823421e10b8SBarry Smith                                        0,
824421e10b8SBarry Smith                                        0,
825421e10b8SBarry Smith                                        0,
826421e10b8SBarry Smith                                        0,
827d519adbfSMatthew Knepley                                /* 74*/ 0,
828421e10b8SBarry Smith                                        0,
829421e10b8SBarry Smith                                        0,
830421e10b8SBarry Smith                                        0,
831421e10b8SBarry Smith                                        0,
832d519adbfSMatthew Knepley                                /* 79*/ 0,
833421e10b8SBarry Smith                                        0,
834421e10b8SBarry Smith                                        0,
835421e10b8SBarry Smith                                        0,
8365bba2384SShri Abhyankar                                        MatLoad_BlockMat,
837d519adbfSMatthew Knepley                                /* 84*/ 0,
838421e10b8SBarry Smith                                        0,
839421e10b8SBarry Smith                                        0,
840421e10b8SBarry Smith                                        0,
841421e10b8SBarry Smith                                        0,
842d519adbfSMatthew Knepley                                /* 89*/ 0,
843421e10b8SBarry Smith                                        0,
844421e10b8SBarry Smith                                        0,
845421e10b8SBarry Smith                                        0,
846421e10b8SBarry Smith                                        0,
847d519adbfSMatthew Knepley                                /* 94*/ 0,
848421e10b8SBarry Smith                                        0,
849421e10b8SBarry Smith                                        0,
850a5973382SShri Abhyankar                                        0,
851a5973382SShri Abhyankar                                        0,
852a5973382SShri Abhyankar                                /* 99*/ 0,
853a5973382SShri Abhyankar                                        0,
854a5973382SShri Abhyankar                                        0,
855a5973382SShri Abhyankar                                        0,
856a5973382SShri Abhyankar                                        0,
857a5973382SShri Abhyankar                                /*104*/ 0,
858a5973382SShri Abhyankar                                        0,
859a5973382SShri Abhyankar                                        0,
860a5973382SShri Abhyankar                                        0,
861a5973382SShri Abhyankar                                        0,
862a5973382SShri Abhyankar                                /*109*/ 0,
863a5973382SShri Abhyankar                                        0,
864a5973382SShri Abhyankar                                        0,
865a5973382SShri Abhyankar                                        0,
866a5973382SShri Abhyankar                                        0,
867a5973382SShri Abhyankar                                /*114*/ 0,
868a5973382SShri Abhyankar                                        0,
869a5973382SShri Abhyankar                                        0,
870a5973382SShri Abhyankar                                        0,
871a5973382SShri Abhyankar                                        0,
872a5973382SShri Abhyankar                                /*119*/ 0,
873a5973382SShri Abhyankar                                        0,
874a5973382SShri Abhyankar                                        0,
8753964eb88SJed Brown                                        0,
8763964eb88SJed Brown                                        0,
8773964eb88SJed Brown                                /*124*/ 0,
8783964eb88SJed Brown                                        0,
8793964eb88SJed Brown                                        0,
8803964eb88SJed Brown                                        0,
8813964eb88SJed Brown                                        0,
8823964eb88SJed Brown                                /*129*/ 0,
8833964eb88SJed Brown                                        0,
8843964eb88SJed Brown                                        0,
8853964eb88SJed Brown                                        0,
8863964eb88SJed Brown                                        0,
8873964eb88SJed Brown                                /*134*/ 0,
8883964eb88SJed Brown                                        0,
8893964eb88SJed Brown                                        0,
8903964eb88SJed Brown                                        0,
8913964eb88SJed Brown                                        0,
8923964eb88SJed Brown                                /*139*/ 0,
893f9426fe0SMark Adams                                        0,
8945bba2384SShri Abhyankar                                        0
895a5973382SShri Abhyankar };
896421e10b8SBarry Smith 
8978cdf2d9bSBarry Smith #undef __FUNCT__
8988cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation"
8998cdf2d9bSBarry Smith /*@C
9008cdf2d9bSBarry Smith    MatBlockMatSetPreallocation - For good matrix assembly performance
9018cdf2d9bSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
9028cdf2d9bSBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
9038cdf2d9bSBarry Smith    during matrix assembly can be increased by more than a factor of 50.
9048cdf2d9bSBarry Smith 
9058cdf2d9bSBarry Smith    Collective on MPI_Comm
9068cdf2d9bSBarry Smith 
9078cdf2d9bSBarry Smith    Input Parameters:
9088cdf2d9bSBarry Smith +  B - The matrix
9098cdf2d9bSBarry Smith .  bs - size of each block in matrix
9108cdf2d9bSBarry Smith .  nz - number of nonzeros per block row (same for all rows)
9118cdf2d9bSBarry Smith -  nnz - array containing the number of nonzeros in the various block rows
9120298fd71SBarry Smith          (possibly different for each row) or NULL
9138cdf2d9bSBarry Smith 
9148cdf2d9bSBarry Smith    Notes:
9158cdf2d9bSBarry Smith      If nnz is given then nz is ignored
9168cdf2d9bSBarry Smith 
9178cdf2d9bSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
9180298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
9198cdf2d9bSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
9208cdf2d9bSBarry Smith    will get TERRIBLE performance, see the users' manual chapter on matrices.
9218cdf2d9bSBarry Smith 
9228cdf2d9bSBarry Smith    Level: intermediate
9238cdf2d9bSBarry Smith 
9248cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues()
9258cdf2d9bSBarry Smith 
9268cdf2d9bSBarry Smith @*/
9277087cfbeSBarry Smith PetscErrorCode  MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
9288cdf2d9bSBarry Smith {
9294ac538c5SBarry Smith   PetscErrorCode ierr;
9308cdf2d9bSBarry Smith 
9318cdf2d9bSBarry Smith   PetscFunctionBegin;
9324ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr);
9338cdf2d9bSBarry Smith   PetscFunctionReturn(0);
9348cdf2d9bSBarry Smith }
9358cdf2d9bSBarry Smith 
9368cdf2d9bSBarry Smith #undef __FUNCT__
9378cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat"
9387087cfbeSBarry Smith PetscErrorCode  MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz)
9398cdf2d9bSBarry Smith {
9408cdf2d9bSBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
9418cdf2d9bSBarry Smith   PetscErrorCode ierr;
9428cdf2d9bSBarry Smith   PetscInt       i;
9438cdf2d9bSBarry Smith 
9448cdf2d9bSBarry Smith   PetscFunctionBegin;
945e02043d6SBarry Smith   ierr = PetscLayoutSetBlockSize(A->rmap,bs);CHKERRQ(ierr);
946e02043d6SBarry Smith   ierr = PetscLayoutSetBlockSize(A->cmap,bs);CHKERRQ(ierr);
94734ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
94834ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
949e02043d6SBarry Smith   ierr = PetscLayoutGetBlockSize(A->rmap,&bs);CHKERRQ(ierr);
95034ef9618SShri Abhyankar 
9518cdf2d9bSBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
952e32f2f54SBarry Smith   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
9538cdf2d9bSBarry Smith   if (nnz) {
954d0f46423SBarry Smith     for (i=0; i<A->rmap->n/bs; i++) {
955e32f2f54SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
956e32f2f54SBarry Smith       if (nnz[i] > A->cmap->n/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],A->cmap->n/bs);
9578cdf2d9bSBarry Smith     }
9588cdf2d9bSBarry Smith   }
959d0f46423SBarry Smith   bmat->mbs = A->rmap->n/bs;
9608cdf2d9bSBarry Smith 
9610298fd71SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->right);CHKERRQ(ierr);
9620298fd71SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->middle);CHKERRQ(ierr);
9638cdf2d9bSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr);
9648cdf2d9bSBarry Smith 
9652ee49352SLisandro Dalcin   if (!bmat->imax) {
966dcca6d9dSJed Brown     ierr = PetscMalloc2(A->rmap->n,&bmat->imax,A->rmap->n,&bmat->ilen);CHKERRQ(ierr);
9673bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
9682ee49352SLisandro Dalcin   }
9698cdf2d9bSBarry Smith   if (nnz) {
9708cdf2d9bSBarry Smith     nz = 0;
971d0f46423SBarry Smith     for (i=0; i<A->rmap->n/A->rmap->bs; i++) {
9728cdf2d9bSBarry Smith       bmat->imax[i] = nnz[i];
9738cdf2d9bSBarry Smith       nz           += nnz[i];
9748cdf2d9bSBarry Smith     }
975f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation");
9768cdf2d9bSBarry Smith 
9778cdf2d9bSBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
9782205254eSKarl Rupp   for (i=0; i<bmat->mbs; i++) bmat->ilen[i] = 0;
9798cdf2d9bSBarry Smith 
9808cdf2d9bSBarry Smith   /* allocate the matrix space */
9812ee49352SLisandro Dalcin   ierr       = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr);
982dcca6d9dSJed Brown   ierr       = PetscMalloc3(nz,&bmat->a,nz,&bmat->j,A->rmap->n+1,&bmat->i);CHKERRQ(ierr);
9833bb1ff40SBarry Smith   ierr       = PetscLogObjectMemory((PetscObject)A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
9848cdf2d9bSBarry Smith   bmat->i[0] = 0;
9858cdf2d9bSBarry Smith   for (i=1; i<bmat->mbs+1; i++) {
9868cdf2d9bSBarry Smith     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
9878cdf2d9bSBarry Smith   }
9888cdf2d9bSBarry Smith   bmat->singlemalloc = PETSC_TRUE;
9898cdf2d9bSBarry Smith   bmat->free_a       = PETSC_TRUE;
9908cdf2d9bSBarry Smith   bmat->free_ij      = PETSC_TRUE;
9918cdf2d9bSBarry Smith 
9928cdf2d9bSBarry Smith   bmat->nz            = 0;
9938cdf2d9bSBarry Smith   bmat->maxnz         = nz;
9948cdf2d9bSBarry Smith   A->info.nz_unneeded = (double)bmat->maxnz;
9957827cd58SJed Brown   ierr                = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
9968cdf2d9bSBarry Smith   PetscFunctionReturn(0);
9978cdf2d9bSBarry Smith }
9988cdf2d9bSBarry Smith 
999421e10b8SBarry Smith /*MC
1000421e10b8SBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
1001421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
1002421e10b8SBarry Smith 
1003421e10b8SBarry Smith   Level: advanced
1004421e10b8SBarry Smith 
1005421e10b8SBarry Smith .seealso: MatCreateBlockMat()
1006421e10b8SBarry Smith 
1007421e10b8SBarry Smith M*/
1008421e10b8SBarry Smith 
1009421e10b8SBarry Smith #undef __FUNCT__
1010421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat"
10118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A)
1012421e10b8SBarry Smith {
1013421e10b8SBarry Smith   Mat_BlockMat   *b;
1014421e10b8SBarry Smith   PetscErrorCode ierr;
1015421e10b8SBarry Smith 
1016421e10b8SBarry Smith   PetscFunctionBegin;
1017b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1018421e10b8SBarry Smith   A->data = (void*)b;
101938f2d2fdSLisandro Dalcin   ierr    = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1020421e10b8SBarry Smith 
1021421e10b8SBarry Smith   A->assembled    = PETSC_TRUE;
1022421e10b8SBarry Smith   A->preallocated = PETSC_FALSE;
1023421e10b8SBarry Smith   ierr            = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
1024290bbb0aSBarry Smith 
1025bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatBlockMatSetPreallocation_C",MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr);
1026421e10b8SBarry Smith   PetscFunctionReturn(0);
1027421e10b8SBarry Smith }
1028421e10b8SBarry Smith 
1029421e10b8SBarry Smith #undef __FUNCT__
1030421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat"
1031421e10b8SBarry Smith /*@C
1032bbd3f336SJed Brown    MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential Mat object
1033421e10b8SBarry Smith 
1034421e10b8SBarry Smith   Collective on MPI_Comm
1035421e10b8SBarry Smith 
1036421e10b8SBarry Smith    Input Parameters:
1037421e10b8SBarry Smith +  comm - MPI communicator
1038421e10b8SBarry Smith .  m - number of rows
1039421e10b8SBarry Smith .  n  - number of columns
10408cdf2d9bSBarry Smith .  bs - size of each submatrix
10418cdf2d9bSBarry Smith .  nz  - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known)
10420298fd71SBarry Smith -  nnz - expected number of nonzers per block row if known (use NULL otherwise)
1043421e10b8SBarry Smith 
1044421e10b8SBarry Smith 
1045421e10b8SBarry Smith    Output Parameter:
1046421e10b8SBarry Smith .  A - the matrix
1047421e10b8SBarry Smith 
1048421e10b8SBarry Smith    Level: intermediate
1049421e10b8SBarry Smith 
1050bbd3f336SJed Brown    Notes: Matrices of this type are nominally-sparse matrices in which each "entry" is a Mat object.  Each Mat must
1051bbd3f336SJed Brown    have the same size and be sequential.  The local and global sizes must be compatible with this decomposition.
1052bbd3f336SJed Brown 
1053bbd3f336SJed Brown    For matrices containing parallel submatrices and variable block sizes, see MATNEST.
1054421e10b8SBarry Smith 
1055421e10b8SBarry Smith .keywords: matrix, bmat, create
1056421e10b8SBarry Smith 
1057bbd3f336SJed Brown .seealso: MATBLOCKMAT, MatCreateNest()
1058421e10b8SBarry Smith @*/
10597087cfbeSBarry Smith PetscErrorCode  MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A)
1060421e10b8SBarry Smith {
1061421e10b8SBarry Smith   PetscErrorCode ierr;
1062421e10b8SBarry Smith 
1063421e10b8SBarry Smith   PetscFunctionBegin;
1064421e10b8SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1065421e10b8SBarry Smith   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1066421e10b8SBarry Smith   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
10678cdf2d9bSBarry Smith   ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1068421e10b8SBarry Smith   PetscFunctionReturn(0);
1069421e10b8SBarry Smith }
1070