xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision 421e10b8f79bbed49dcc4e803c884835d979c6ea)
1*421e10b8SBarry Smith #define PETSCMAT_DLL
2*421e10b8SBarry Smith 
3*421e10b8SBarry Smith /*
4*421e10b8SBarry Smith    This provides a matrix that consists of Mats
5*421e10b8SBarry Smith */
6*421e10b8SBarry Smith 
7*421e10b8SBarry Smith #include "src/mat/matimpl.h"              /*I "petscmat.h" I*/
8*421e10b8SBarry Smith #include "src/mat/impls/baij/seq/baij.h"    /* use the common AIJ data-structure */
9*421e10b8SBarry Smith 
10*421e10b8SBarry Smith #define CHUNKSIZE   15
11*421e10b8SBarry Smith 
12*421e10b8SBarry Smith typedef struct {
13*421e10b8SBarry Smith   SEQAIJHEADER(Mat);
14*421e10b8SBarry Smith   SEQBAIJHEADER;
15*421e10b8SBarry Smith 
16*421e10b8SBarry Smith   Vec  left,right;   /* dummy vectors to perform local parts of product */
17*421e10b8SBarry Smith } Mat_BlockMat;
18*421e10b8SBarry Smith 
19*421e10b8SBarry Smith #undef __FUNCT__
20*421e10b8SBarry Smith #define __FUNCT__ "MatSetValues_BlockMat"
21*421e10b8SBarry Smith PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
22*421e10b8SBarry Smith {
23*421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
24*421e10b8SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
25*421e10b8SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen;
26*421e10b8SBarry Smith   PetscInt       *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol;
27*421e10b8SBarry Smith   PetscErrorCode ierr;
28*421e10b8SBarry Smith   PetscInt       ridx,cidx;
29*421e10b8SBarry Smith   PetscTruth     roworiented=a->roworiented;
30*421e10b8SBarry Smith   MatScalar      value;
31*421e10b8SBarry Smith   Mat            *ap,*aa = a->a;
32*421e10b8SBarry Smith 
33*421e10b8SBarry Smith   PetscFunctionBegin;
34*421e10b8SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
35*421e10b8SBarry Smith     row  = im[k];
36*421e10b8SBarry Smith     brow = row/bs;
37*421e10b8SBarry Smith     if (row < 0) continue;
38*421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
39*421e10b8SBarry Smith     if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1);
40*421e10b8SBarry Smith #endif
41*421e10b8SBarry Smith     rp   = aj + ai[brow];
42*421e10b8SBarry Smith     ap   = aa + ai[brow];
43*421e10b8SBarry Smith     rmax = imax[brow];
44*421e10b8SBarry Smith     nrow = ailen[brow];
45*421e10b8SBarry Smith     low  = 0;
46*421e10b8SBarry Smith     high = nrow;
47*421e10b8SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
48*421e10b8SBarry Smith       if (in[l] < 0) continue;
49*421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG)
50*421e10b8SBarry 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);
51*421e10b8SBarry Smith #endif
52*421e10b8SBarry Smith       col = in[l]; bcol = col/bs;
53*421e10b8SBarry Smith       ridx = row % bs; cidx = col % bs;
54*421e10b8SBarry Smith       if (roworiented) {
55*421e10b8SBarry Smith         value = v[l + k*n];
56*421e10b8SBarry Smith       } else {
57*421e10b8SBarry Smith         value = v[k + l*m];
58*421e10b8SBarry Smith       }
59*421e10b8SBarry Smith       if (col <= lastcol) low = 0; else high = nrow;
60*421e10b8SBarry Smith       lastcol = col;
61*421e10b8SBarry Smith       while (high-low > 7) {
62*421e10b8SBarry Smith         t = (low+high)/2;
63*421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
64*421e10b8SBarry Smith         else              low  = t;
65*421e10b8SBarry Smith       }
66*421e10b8SBarry Smith       for (i=low; i<high; i++) {
67*421e10b8SBarry Smith         if (rp[i] > bcol) break;
68*421e10b8SBarry Smith         if (rp[i] == bcol) {
69*421e10b8SBarry Smith           goto noinsert1;
70*421e10b8SBarry Smith         }
71*421e10b8SBarry Smith       }
72*421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
73*421e10b8SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
74*421e10b8SBarry Smith       MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
75*421e10b8SBarry Smith       N = nrow++ - 1; high++;
76*421e10b8SBarry Smith       /* shift up all the later entries in this row */
77*421e10b8SBarry Smith       for (ii=N; ii>=i; ii--) {
78*421e10b8SBarry Smith         rp[ii+1] = rp[ii];
79*421e10b8SBarry Smith         ap[ii+1] = ap[ii];
80*421e10b8SBarry Smith       }
81*421e10b8SBarry Smith       if (N>=i) ap[i] = 0;
82*421e10b8SBarry Smith       rp[i]           = bcol;
83*421e10b8SBarry Smith       a->nz++;
84*421e10b8SBarry Smith       noinsert1:;
85*421e10b8SBarry Smith       if (!*(ap+i)) {
86*421e10b8SBarry Smith         ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr);
87*421e10b8SBarry Smith       }
88*421e10b8SBarry Smith       ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr);
89*421e10b8SBarry Smith       low = i;
90*421e10b8SBarry Smith     }
91*421e10b8SBarry Smith     ailen[brow] = nrow;
92*421e10b8SBarry Smith   }
93*421e10b8SBarry Smith   A->same_nonzero = PETSC_FALSE;
94*421e10b8SBarry Smith   PetscFunctionReturn(0);
95*421e10b8SBarry Smith }
96*421e10b8SBarry Smith 
97*421e10b8SBarry Smith #undef __FUNCT__
98*421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat"
99*421e10b8SBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A)
100*421e10b8SBarry Smith {
101*421e10b8SBarry Smith   PetscErrorCode    ierr;
102*421e10b8SBarry Smith   Mat               tmpA;
103*421e10b8SBarry Smith   PetscInt          i,m,n,bs = 1,ncols;
104*421e10b8SBarry Smith   const PetscInt    *cols;
105*421e10b8SBarry Smith   const PetscScalar *values;
106*421e10b8SBarry Smith 
107*421e10b8SBarry Smith   PetscFunctionBegin;
108*421e10b8SBarry Smith   ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr);
109*421e10b8SBarry Smith 
110*421e10b8SBarry Smith   ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr);
111*421e10b8SBarry Smith   ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr);
112*421e10b8SBarry Smith     ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr);
113*421e10b8SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
114*421e10b8SBarry Smith   ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,A);CHKERRQ(ierr);
115*421e10b8SBarry Smith   for (i=0; i<m; i++) {
116*421e10b8SBarry Smith     ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr);
117*421e10b8SBarry Smith     ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr);
118*421e10b8SBarry Smith   }
119*421e10b8SBarry Smith   ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
120*421e10b8SBarry Smith   ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
121*421e10b8SBarry Smith   PetscFunctionReturn(0);
122*421e10b8SBarry Smith }
123*421e10b8SBarry Smith 
124*421e10b8SBarry Smith 
125*421e10b8SBarry Smith #undef __FUNCT__
126*421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat"
127*421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat)
128*421e10b8SBarry Smith {
129*421e10b8SBarry Smith   PetscErrorCode ierr;
130*421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)mat->data;
131*421e10b8SBarry Smith 
132*421e10b8SBarry Smith   PetscFunctionBegin;
133*421e10b8SBarry Smith   if (bmat->right) {
134*421e10b8SBarry Smith     ierr = VecDestroy(bmat->right);CHKERRQ(ierr);
135*421e10b8SBarry Smith   }
136*421e10b8SBarry Smith   if (bmat->left) {
137*421e10b8SBarry Smith     ierr = VecDestroy(bmat->left);CHKERRQ(ierr);
138*421e10b8SBarry Smith   }
139*421e10b8SBarry Smith   ierr = PetscFree(bmat);CHKERRQ(ierr);
140*421e10b8SBarry Smith   PetscFunctionReturn(0);
141*421e10b8SBarry Smith }
142*421e10b8SBarry Smith 
143*421e10b8SBarry Smith #undef __FUNCT__
144*421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat"
145*421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
146*421e10b8SBarry Smith {
147*421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
148*421e10b8SBarry Smith   PetscErrorCode ierr;
149*421e10b8SBarry Smith   PetscScalar    *xx,*yy;
150*421e10b8SBarry Smith   PetscInt       *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j;
151*421e10b8SBarry Smith   Mat            *aa;
152*421e10b8SBarry Smith 
153*421e10b8SBarry Smith   PetscFunctionBegin;
154*421e10b8SBarry Smith   /*
155*421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
156*421e10b8SBarry Smith   */
157*421e10b8SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
158*421e10b8SBarry Smith   ierr = VecSet(y,0.0);CHKERRQ(ierr);
159*421e10b8SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
160*421e10b8SBarry Smith   aj  = bmat->j;
161*421e10b8SBarry Smith   aa  = bmat->a;
162*421e10b8SBarry Smith   ii  = bmat->i;
163*421e10b8SBarry Smith   for (i=0; i<m; i++) {
164*421e10b8SBarry Smith     jrow = ii[i];
165*421e10b8SBarry Smith     ierr = VecPlaceArray(bmat->left,yy + bs*jrow);CHKERRQ(ierr);
166*421e10b8SBarry Smith     n    = ii[i+1] - jrow;
167*421e10b8SBarry Smith     ierr = VecSet(bmat->left,0.0);CHKERRQ(ierr);
168*421e10b8SBarry Smith     for (j=0; j<n; j++) {
169*421e10b8SBarry Smith       ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr);
170*421e10b8SBarry Smith       ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);
171*421e10b8SBarry Smith       ierr = VecResetArray(bmat->right);CHKERRQ(ierr);
172*421e10b8SBarry Smith       jrow++;
173*421e10b8SBarry Smith     }
174*421e10b8SBarry Smith     ierr = VecResetArray(bmat->left);CHKERRQ(ierr);
175*421e10b8SBarry Smith   }
176*421e10b8SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
177*421e10b8SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
178*421e10b8SBarry Smith   PetscFunctionReturn(0);
179*421e10b8SBarry Smith }
180*421e10b8SBarry Smith 
181*421e10b8SBarry Smith #undef __FUNCT__
182*421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat"
183*421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
184*421e10b8SBarry Smith {
185*421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
186*421e10b8SBarry Smith   PetscErrorCode ierr;
187*421e10b8SBarry Smith 
188*421e10b8SBarry Smith   PetscFunctionBegin;
189*421e10b8SBarry Smith   PetscFunctionReturn(0);
190*421e10b8SBarry Smith }
191*421e10b8SBarry Smith 
192*421e10b8SBarry Smith #undef __FUNCT__
193*421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat"
194*421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
195*421e10b8SBarry Smith {
196*421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
197*421e10b8SBarry Smith   PetscErrorCode ierr;
198*421e10b8SBarry Smith 
199*421e10b8SBarry Smith   PetscFunctionBegin;
200*421e10b8SBarry Smith   PetscFunctionReturn(0);
201*421e10b8SBarry Smith }
202*421e10b8SBarry Smith 
203*421e10b8SBarry Smith #undef __FUNCT__
204*421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat"
205*421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
206*421e10b8SBarry Smith {
207*421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
208*421e10b8SBarry Smith   PetscErrorCode ierr;
209*421e10b8SBarry Smith 
210*421e10b8SBarry Smith   PetscFunctionBegin;
211*421e10b8SBarry Smith   PetscFunctionReturn(0);
212*421e10b8SBarry Smith }
213*421e10b8SBarry Smith 
214*421e10b8SBarry Smith #undef __FUNCT__
215*421e10b8SBarry Smith #define __FUNCT__ "MatSetBlockSize_BlockMat"
216*421e10b8SBarry Smith PetscErrorCode MatSetBlockSize_BlockMat(Mat A,PetscInt bs)
217*421e10b8SBarry Smith {
218*421e10b8SBarry Smith   Mat_BlockMat   *bmat = (Mat_BlockMat*)A->data;
219*421e10b8SBarry Smith   PetscErrorCode ierr;
220*421e10b8SBarry Smith   PetscInt       nz = 10,i,m;
221*421e10b8SBarry Smith 
222*421e10b8SBarry Smith   PetscFunctionBegin;
223*421e10b8SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr);
224*421e10b8SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->left);CHKERRQ(ierr);
225*421e10b8SBarry Smith 
226*421e10b8SBarry Smith   A->rmap.bs = A->cmap.bs = bs;
227*421e10b8SBarry Smith 
228*421e10b8SBarry Smith   ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr);
229*421e10b8SBarry Smith   for (i=0; i<A->rmap.n; i++) bmat->imax[i] = nz;
230*421e10b8SBarry Smith   nz = nz*A->rmap.n;
231*421e10b8SBarry Smith 
232*421e10b8SBarry Smith 
233*421e10b8SBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
234*421e10b8SBarry Smith   for (i=0; i<A->rmap.n; i++) { bmat->ilen[i] = 0;}
235*421e10b8SBarry Smith 
236*421e10b8SBarry Smith   /* allocate the matrix space */
237*421e10b8SBarry Smith   ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr);
238*421e10b8SBarry Smith   bmat->i[0] = 0;
239*421e10b8SBarry Smith   for (i=1; i<A->rmap.n+1; i++) {
240*421e10b8SBarry Smith     bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
241*421e10b8SBarry Smith   }
242*421e10b8SBarry Smith   bmat->singlemalloc = PETSC_TRUE;
243*421e10b8SBarry Smith   bmat->free_a       = PETSC_TRUE;
244*421e10b8SBarry Smith   bmat->free_ij      = PETSC_TRUE;
245*421e10b8SBarry Smith 
246*421e10b8SBarry Smith   bmat->nz                = 0;
247*421e10b8SBarry Smith   bmat->maxnz             = nz;
248*421e10b8SBarry Smith   A->info.nz_unneeded  = (double)bmat->maxnz;
249*421e10b8SBarry Smith 
250*421e10b8SBarry Smith   PetscFunctionReturn(0);
251*421e10b8SBarry Smith }
252*421e10b8SBarry Smith 
253*421e10b8SBarry Smith #undef __FUNCT__
254*421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat"
255*421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
256*421e10b8SBarry Smith {
257*421e10b8SBarry Smith   Mat_BlockMat   *a = (Mat_BlockMat*)A->data;
258*421e10b8SBarry Smith   PetscErrorCode ierr;
259*421e10b8SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
260*421e10b8SBarry Smith   PetscInt       m = A->rmap.n/A->rmap.bs,*ip,N,*ailen = a->ilen,rmax = 0;
261*421e10b8SBarry Smith   Mat            *aa = a->a,*ap;
262*421e10b8SBarry Smith   PetscReal      ratio=0.6;
263*421e10b8SBarry Smith 
264*421e10b8SBarry Smith   PetscFunctionBegin;
265*421e10b8SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
266*421e10b8SBarry Smith 
267*421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
268*421e10b8SBarry Smith   for (i=1; i<m; i++) {
269*421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
270*421e10b8SBarry Smith     fshift += imax[i-1] - ailen[i-1];
271*421e10b8SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
272*421e10b8SBarry Smith     if (fshift) {
273*421e10b8SBarry Smith       ip = aj + ai[i] ;
274*421e10b8SBarry Smith       ap = aa + ai[i] ;
275*421e10b8SBarry Smith       N  = ailen[i];
276*421e10b8SBarry Smith       for (j=0; j<N; j++) {
277*421e10b8SBarry Smith         ip[j-fshift] = ip[j];
278*421e10b8SBarry Smith         ap[j-fshift] = ap[j];
279*421e10b8SBarry Smith       }
280*421e10b8SBarry Smith     }
281*421e10b8SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
282*421e10b8SBarry Smith   }
283*421e10b8SBarry Smith   if (m) {
284*421e10b8SBarry Smith     fshift += imax[m-1] - ailen[m-1];
285*421e10b8SBarry Smith     ai[m]  = ai[m-1] + ailen[m-1];
286*421e10b8SBarry Smith   }
287*421e10b8SBarry Smith   /* reset ilen and imax for each row */
288*421e10b8SBarry Smith   for (i=0; i<m; i++) {
289*421e10b8SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
290*421e10b8SBarry Smith   }
291*421e10b8SBarry Smith   a->nz = ai[m];
292*421e10b8SBarry Smith 
293*421e10b8SBarry 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);
294*421e10b8SBarry Smith   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr);
295*421e10b8SBarry Smith   ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr);
296*421e10b8SBarry Smith   a->reallocs          = 0;
297*421e10b8SBarry Smith   A->info.nz_unneeded  = (double)fshift;
298*421e10b8SBarry Smith   a->rmax              = rmax;
299*421e10b8SBarry Smith 
300*421e10b8SBarry Smith   A->same_nonzero = PETSC_TRUE;
301*421e10b8SBarry Smith   PetscFunctionReturn(0);
302*421e10b8SBarry Smith }
303*421e10b8SBarry Smith 
304*421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
305*421e10b8SBarry Smith        0,
306*421e10b8SBarry Smith        0,
307*421e10b8SBarry Smith        MatMult_BlockMat,
308*421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat,
309*421e10b8SBarry Smith        MatMultTranspose_BlockMat,
310*421e10b8SBarry Smith        MatMultTransposeAdd_BlockMat,
311*421e10b8SBarry Smith        0,
312*421e10b8SBarry Smith        0,
313*421e10b8SBarry Smith        0,
314*421e10b8SBarry Smith /*10*/ 0,
315*421e10b8SBarry Smith        0,
316*421e10b8SBarry Smith        0,
317*421e10b8SBarry Smith        0,
318*421e10b8SBarry Smith        0,
319*421e10b8SBarry Smith /*15*/ 0,
320*421e10b8SBarry Smith        0,
321*421e10b8SBarry Smith        0,
322*421e10b8SBarry Smith        0,
323*421e10b8SBarry Smith        0,
324*421e10b8SBarry Smith /*20*/ 0,
325*421e10b8SBarry Smith        MatAssemblyEnd_BlockMat,
326*421e10b8SBarry Smith        0,
327*421e10b8SBarry Smith        0,
328*421e10b8SBarry Smith        0,
329*421e10b8SBarry Smith /*25*/ 0,
330*421e10b8SBarry Smith        0,
331*421e10b8SBarry Smith        0,
332*421e10b8SBarry Smith        0,
333*421e10b8SBarry Smith        0,
334*421e10b8SBarry Smith /*30*/ 0,
335*421e10b8SBarry Smith        0,
336*421e10b8SBarry Smith        0,
337*421e10b8SBarry Smith        0,
338*421e10b8SBarry Smith        0,
339*421e10b8SBarry Smith /*35*/ 0,
340*421e10b8SBarry Smith        0,
341*421e10b8SBarry Smith        0,
342*421e10b8SBarry Smith        0,
343*421e10b8SBarry Smith        0,
344*421e10b8SBarry Smith /*40*/ 0,
345*421e10b8SBarry Smith        0,
346*421e10b8SBarry Smith        0,
347*421e10b8SBarry Smith        0,
348*421e10b8SBarry Smith        0,
349*421e10b8SBarry Smith /*45*/ 0,
350*421e10b8SBarry Smith        0,
351*421e10b8SBarry Smith        0,
352*421e10b8SBarry Smith        0,
353*421e10b8SBarry Smith        0,
354*421e10b8SBarry Smith /*50*/ MatSetBlockSize_BlockMat,
355*421e10b8SBarry Smith        0,
356*421e10b8SBarry Smith        0,
357*421e10b8SBarry Smith        0,
358*421e10b8SBarry Smith        0,
359*421e10b8SBarry Smith /*55*/ 0,
360*421e10b8SBarry Smith        0,
361*421e10b8SBarry Smith        0,
362*421e10b8SBarry Smith        0,
363*421e10b8SBarry Smith        0,
364*421e10b8SBarry Smith /*60*/ 0,
365*421e10b8SBarry Smith        MatDestroy_BlockMat,
366*421e10b8SBarry Smith        0,
367*421e10b8SBarry Smith        0,
368*421e10b8SBarry Smith        0,
369*421e10b8SBarry Smith /*65*/ 0,
370*421e10b8SBarry Smith        0,
371*421e10b8SBarry Smith        0,
372*421e10b8SBarry Smith        0,
373*421e10b8SBarry Smith        0,
374*421e10b8SBarry Smith /*70*/ 0,
375*421e10b8SBarry Smith        0,
376*421e10b8SBarry Smith        0,
377*421e10b8SBarry Smith        0,
378*421e10b8SBarry Smith        0,
379*421e10b8SBarry Smith /*75*/ 0,
380*421e10b8SBarry Smith        0,
381*421e10b8SBarry Smith        0,
382*421e10b8SBarry Smith        0,
383*421e10b8SBarry Smith        0,
384*421e10b8SBarry Smith /*80*/ 0,
385*421e10b8SBarry Smith        0,
386*421e10b8SBarry Smith        0,
387*421e10b8SBarry Smith        0,
388*421e10b8SBarry Smith        MatLoad_BlockMat,
389*421e10b8SBarry Smith /*85*/ 0,
390*421e10b8SBarry Smith        0,
391*421e10b8SBarry Smith        0,
392*421e10b8SBarry Smith        0,
393*421e10b8SBarry Smith        0,
394*421e10b8SBarry Smith /*90*/ 0,
395*421e10b8SBarry Smith        0,
396*421e10b8SBarry Smith        0,
397*421e10b8SBarry Smith        0,
398*421e10b8SBarry Smith        0,
399*421e10b8SBarry Smith /*95*/ 0,
400*421e10b8SBarry Smith        0,
401*421e10b8SBarry Smith        0,
402*421e10b8SBarry Smith        0};
403*421e10b8SBarry Smith 
404*421e10b8SBarry Smith /*MC
405*421e10b8SBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
406*421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
407*421e10b8SBarry Smith 
408*421e10b8SBarry Smith   Level: advanced
409*421e10b8SBarry Smith 
410*421e10b8SBarry Smith .seealso: MatCreateBlockMat()
411*421e10b8SBarry Smith 
412*421e10b8SBarry Smith M*/
413*421e10b8SBarry Smith 
414*421e10b8SBarry Smith EXTERN_C_BEGIN
415*421e10b8SBarry Smith #undef __FUNCT__
416*421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat"
417*421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A)
418*421e10b8SBarry Smith {
419*421e10b8SBarry Smith   Mat_BlockMat    *b;
420*421e10b8SBarry Smith   PetscErrorCode ierr;
421*421e10b8SBarry Smith 
422*421e10b8SBarry Smith   PetscFunctionBegin;
423*421e10b8SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
424*421e10b8SBarry Smith   ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr);
425*421e10b8SBarry Smith 
426*421e10b8SBarry Smith   A->data = (void*)b;
427*421e10b8SBarry Smith 
428*421e10b8SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
429*421e10b8SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
430*421e10b8SBarry Smith 
431*421e10b8SBarry Smith   A->assembled     = PETSC_TRUE;
432*421e10b8SBarry Smith   A->preallocated  = PETSC_FALSE;
433*421e10b8SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr);
434*421e10b8SBarry Smith   PetscFunctionReturn(0);
435*421e10b8SBarry Smith }
436*421e10b8SBarry Smith EXTERN_C_END
437*421e10b8SBarry Smith 
438*421e10b8SBarry Smith #undef __FUNCT__
439*421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat"
440*421e10b8SBarry Smith /*@C
441*421e10b8SBarry Smith    MatCreateBlockMat - Creates a new matrix based sparse Mat storage
442*421e10b8SBarry Smith 
443*421e10b8SBarry Smith   Collective on MPI_Comm
444*421e10b8SBarry Smith 
445*421e10b8SBarry Smith    Input Parameters:
446*421e10b8SBarry Smith +  comm - MPI communicator
447*421e10b8SBarry Smith .  m - number of rows
448*421e10b8SBarry Smith .  n  - number of columns
449*421e10b8SBarry Smith -  bs - size of each submatrix
450*421e10b8SBarry Smith 
451*421e10b8SBarry Smith 
452*421e10b8SBarry Smith    Output Parameter:
453*421e10b8SBarry Smith .  A - the matrix
454*421e10b8SBarry Smith 
455*421e10b8SBarry Smith    Level: intermediate
456*421e10b8SBarry Smith 
457*421e10b8SBarry Smith    PETSc requires that matrices and vectors being used for certain
458*421e10b8SBarry Smith    operations are partitioned accordingly.  For example, when
459*421e10b8SBarry Smith    creating a bmat matrix, A, that supports parallel matrix-vector
460*421e10b8SBarry Smith    products using MatMult(A,x,y) the user should set the number
461*421e10b8SBarry Smith    of local matrix rows to be the number of local elements of the
462*421e10b8SBarry Smith    corresponding result vector, y. Note that this is information is
463*421e10b8SBarry Smith    required for use of the matrix interface routines, even though
464*421e10b8SBarry Smith    the bmat matrix may not actually be physically partitioned.
465*421e10b8SBarry Smith    For example,
466*421e10b8SBarry Smith 
467*421e10b8SBarry Smith .keywords: matrix, bmat, create
468*421e10b8SBarry Smith 
469*421e10b8SBarry Smith .seealso: MATBLOCKMAT
470*421e10b8SBarry Smith @*/
471*421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,Mat *A)
472*421e10b8SBarry Smith {
473*421e10b8SBarry Smith   PetscErrorCode ierr;
474*421e10b8SBarry Smith 
475*421e10b8SBarry Smith   PetscFunctionBegin;
476*421e10b8SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
477*421e10b8SBarry Smith   ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
478*421e10b8SBarry Smith   ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr);
479*421e10b8SBarry Smith   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
480*421e10b8SBarry Smith   PetscFunctionReturn(0);
481*421e10b8SBarry Smith }
482*421e10b8SBarry Smith 
483*421e10b8SBarry Smith 
484*421e10b8SBarry Smith 
485