xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 61e22dc2a79bdea5a615071430fbf27c599d8ba0)
1*61e22dc2SBarry Smith /*$Id: baij.c,v 1.209 2000/05/10 16:40:51 bsmith Exp bsmith $*/
2*61e22dc2SBarry Smith 
3*61e22dc2SBarry Smith /*
4*61e22dc2SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
5*61e22dc2SBarry Smith   matrix storage format.
6*61e22dc2SBarry Smith */
7*61e22dc2SBarry Smith #include "petscsys.h"
8*61e22dc2SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
9*61e22dc2SBarry Smith #include "src/vec/vecimpl.h"
10*61e22dc2SBarry Smith #include "src/inline/spops.h"
11*61e22dc2SBarry Smith 
12*61e22dc2SBarry Smith /*  UGLY, ugly, ugly
13*61e22dc2SBarry Smith    When MatScalar == Scalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does
14*61e22dc2SBarry Smith    not exist. Otherwise ..._MatScalar() takes matrix elements in single precision and
15*61e22dc2SBarry Smith    inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ()
16*61e22dc2SBarry Smith    converts the entries into single precision and then calls ..._MatScalar() to put them
17*61e22dc2SBarry Smith    into the single precision data structures.
18*61e22dc2SBarry Smith */
19*61e22dc2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
20*61e22dc2SBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode);
21*61e22dc2SBarry Smith #else
22*61e22dc2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ
23*61e22dc2SBarry Smith #endif
24*61e22dc2SBarry Smith 
25*61e22dc2SBarry Smith #define CHUNKSIZE  10
26*61e22dc2SBarry Smith 
27*61e22dc2SBarry Smith /*
28*61e22dc2SBarry Smith      Checks for missing diagonals
29*61e22dc2SBarry Smith */
30*61e22dc2SBarry Smith #undef __FUNC__
31*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMissingDiagonal_SeqBAIJ"
32*61e22dc2SBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A)
33*61e22dc2SBarry Smith {
34*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
35*61e22dc2SBarry Smith   int         *diag,*jj = a->j,i,ierr;
36*61e22dc2SBarry Smith 
37*61e22dc2SBarry Smith   PetscFunctionBegin;
38*61e22dc2SBarry Smith   ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
39*61e22dc2SBarry Smith   diag = a->diag;
40*61e22dc2SBarry Smith   for (i=0; i<a->mbs; i++) {
41*61e22dc2SBarry Smith     if (jj[diag[i]] != i) {
42*61e22dc2SBarry Smith       SETERRQ1(1,1,"Matrix is missing diagonal number %d",i);
43*61e22dc2SBarry Smith     }
44*61e22dc2SBarry Smith   }
45*61e22dc2SBarry Smith   PetscFunctionReturn(0);
46*61e22dc2SBarry Smith }
47*61e22dc2SBarry Smith 
48*61e22dc2SBarry Smith #undef __FUNC__
49*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMarkDiagonal_SeqBAIJ"
50*61e22dc2SBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A)
51*61e22dc2SBarry Smith {
52*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
53*61e22dc2SBarry Smith   int         i,j,*diag,m = a->mbs;
54*61e22dc2SBarry Smith 
55*61e22dc2SBarry Smith   PetscFunctionBegin;
56*61e22dc2SBarry Smith   if (a->diag) PetscFunctionReturn(0);
57*61e22dc2SBarry Smith 
58*61e22dc2SBarry Smith   diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag);
59*61e22dc2SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
60*61e22dc2SBarry Smith   for (i=0; i<m; i++) {
61*61e22dc2SBarry Smith     diag[i] = a->i[i+1];
62*61e22dc2SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
63*61e22dc2SBarry Smith       if (a->j[j] == i) {
64*61e22dc2SBarry Smith         diag[i] = j;
65*61e22dc2SBarry Smith         break;
66*61e22dc2SBarry Smith       }
67*61e22dc2SBarry Smith     }
68*61e22dc2SBarry Smith   }
69*61e22dc2SBarry Smith   a->diag = diag;
70*61e22dc2SBarry Smith   PetscFunctionReturn(0);
71*61e22dc2SBarry Smith }
72*61e22dc2SBarry Smith 
73*61e22dc2SBarry Smith 
74*61e22dc2SBarry Smith EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
75*61e22dc2SBarry Smith 
76*61e22dc2SBarry Smith #undef __FUNC__
77*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRowIJ_SeqBAIJ"
78*61e22dc2SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
79*61e22dc2SBarry Smith {
80*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
81*61e22dc2SBarry Smith   int         ierr,n = a->mbs,i;
82*61e22dc2SBarry Smith 
83*61e22dc2SBarry Smith   PetscFunctionBegin;
84*61e22dc2SBarry Smith   *nn = n;
85*61e22dc2SBarry Smith   if (!ia) PetscFunctionReturn(0);
86*61e22dc2SBarry Smith   if (symmetric) {
87*61e22dc2SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
88*61e22dc2SBarry Smith   } else if (oshift == 1) {
89*61e22dc2SBarry Smith     /* temporarily add 1 to i and j indices */
90*61e22dc2SBarry Smith     int nz = a->i[n] + 1;
91*61e22dc2SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
92*61e22dc2SBarry Smith     for (i=0; i<n+1; i++) a->i[i]++;
93*61e22dc2SBarry Smith     *ia = a->i; *ja = a->j;
94*61e22dc2SBarry Smith   } else {
95*61e22dc2SBarry Smith     *ia = a->i; *ja = a->j;
96*61e22dc2SBarry Smith   }
97*61e22dc2SBarry Smith 
98*61e22dc2SBarry Smith   PetscFunctionReturn(0);
99*61e22dc2SBarry Smith }
100*61e22dc2SBarry Smith 
101*61e22dc2SBarry Smith #undef __FUNC__
102*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRowIJ_SeqBAIJ"
103*61e22dc2SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
104*61e22dc2SBarry Smith                                 PetscTruth *done)
105*61e22dc2SBarry Smith {
106*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
107*61e22dc2SBarry Smith   int         i,n = a->mbs,ierr;
108*61e22dc2SBarry Smith 
109*61e22dc2SBarry Smith   PetscFunctionBegin;
110*61e22dc2SBarry Smith   if (!ia) PetscFunctionReturn(0);
111*61e22dc2SBarry Smith   if (symmetric) {
112*61e22dc2SBarry Smith     ierr = PetscFree(*ia);CHKERRQ(ierr);
113*61e22dc2SBarry Smith     ierr = PetscFree(*ja);CHKERRQ(ierr);
114*61e22dc2SBarry Smith   } else if (oshift == 1) {
115*61e22dc2SBarry Smith     int nz = a->i[n];
116*61e22dc2SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
117*61e22dc2SBarry Smith     for (i=0; i<n+1; i++) a->i[i]--;
118*61e22dc2SBarry Smith   }
119*61e22dc2SBarry Smith   PetscFunctionReturn(0);
120*61e22dc2SBarry Smith }
121*61e22dc2SBarry Smith 
122*61e22dc2SBarry Smith #undef __FUNC__
123*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_SeqBAIJ"
124*61e22dc2SBarry Smith int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs)
125*61e22dc2SBarry Smith {
126*61e22dc2SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data;
127*61e22dc2SBarry Smith 
128*61e22dc2SBarry Smith   PetscFunctionBegin;
129*61e22dc2SBarry Smith   *bs = baij->bs;
130*61e22dc2SBarry Smith   PetscFunctionReturn(0);
131*61e22dc2SBarry Smith }
132*61e22dc2SBarry Smith 
133*61e22dc2SBarry Smith 
134*61e22dc2SBarry Smith #undef __FUNC__
135*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_SeqBAIJ"
136*61e22dc2SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
137*61e22dc2SBarry Smith {
138*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
139*61e22dc2SBarry Smith   int         ierr;
140*61e22dc2SBarry Smith 
141*61e22dc2SBarry Smith   PetscFunctionBegin;
142*61e22dc2SBarry Smith   if (A->mapping) {
143*61e22dc2SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr);
144*61e22dc2SBarry Smith   }
145*61e22dc2SBarry Smith   if (A->bmapping) {
146*61e22dc2SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr);
147*61e22dc2SBarry Smith   }
148*61e22dc2SBarry Smith   if (A->rmap) {
149*61e22dc2SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
150*61e22dc2SBarry Smith   }
151*61e22dc2SBarry Smith   if (A->cmap) {
152*61e22dc2SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
153*61e22dc2SBarry Smith   }
154*61e22dc2SBarry Smith #if defined(PETSC_USE_LOG)
155*61e22dc2SBarry Smith   PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
156*61e22dc2SBarry Smith #endif
157*61e22dc2SBarry Smith   ierr = PetscFree(a->a);CHKERRQ(ierr);
158*61e22dc2SBarry Smith   if (!a->singlemalloc) {
159*61e22dc2SBarry Smith     ierr = PetscFree(a->i);CHKERRQ(ierr);
160*61e22dc2SBarry Smith     ierr = PetscFree(a->j);CHKERRQ(ierr);
161*61e22dc2SBarry Smith   }
162*61e22dc2SBarry Smith   if (a->row) {
163*61e22dc2SBarry Smith     ierr = ISDestroy(a->row);CHKERRQ(ierr);
164*61e22dc2SBarry Smith   }
165*61e22dc2SBarry Smith   if (a->col) {
166*61e22dc2SBarry Smith     ierr = ISDestroy(a->col);CHKERRQ(ierr);
167*61e22dc2SBarry Smith   }
168*61e22dc2SBarry Smith   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
169*61e22dc2SBarry Smith   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
170*61e22dc2SBarry Smith   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
171*61e22dc2SBarry Smith   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
172*61e22dc2SBarry Smith   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
173*61e22dc2SBarry Smith   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
174*61e22dc2SBarry Smith   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
175*61e22dc2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
176*61e22dc2SBarry Smith   if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);}
177*61e22dc2SBarry Smith #endif
178*61e22dc2SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
179*61e22dc2SBarry Smith   PLogObjectDestroy(A);
180*61e22dc2SBarry Smith   PetscHeaderDestroy(A);
181*61e22dc2SBarry Smith   PetscFunctionReturn(0);
182*61e22dc2SBarry Smith }
183*61e22dc2SBarry Smith 
184*61e22dc2SBarry Smith #undef __FUNC__
185*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_SeqBAIJ"
186*61e22dc2SBarry Smith int MatSetOption_SeqBAIJ(Mat A,MatOption op)
187*61e22dc2SBarry Smith {
188*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
189*61e22dc2SBarry Smith 
190*61e22dc2SBarry Smith   PetscFunctionBegin;
191*61e22dc2SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented    = PETSC_TRUE;
192*61e22dc2SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented    = PETSC_FALSE;
193*61e22dc2SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted         = PETSC_TRUE;
194*61e22dc2SBarry Smith   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted         = PETSC_FALSE;
195*61e22dc2SBarry Smith   else if (op == MAT_KEEP_ZEROED_ROWS)             a->keepzeroedrows = PETSC_TRUE;
196*61e22dc2SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew          = 1;
197*61e22dc2SBarry Smith   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew          = -1;
198*61e22dc2SBarry Smith   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew          = -2;
199*61e22dc2SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew          = 0;
200*61e22dc2SBarry Smith   else if (op == MAT_ROWS_SORTED ||
201*61e22dc2SBarry Smith            op == MAT_ROWS_UNSORTED ||
202*61e22dc2SBarry Smith            op == MAT_SYMMETRIC ||
203*61e22dc2SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
204*61e22dc2SBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
205*61e22dc2SBarry Smith            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
206*61e22dc2SBarry Smith            op == MAT_USE_HASH_TABLE) {
207*61e22dc2SBarry Smith     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
208*61e22dc2SBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
209*61e22dc2SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
210*61e22dc2SBarry Smith   } else {
211*61e22dc2SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
212*61e22dc2SBarry Smith   }
213*61e22dc2SBarry Smith   PetscFunctionReturn(0);
214*61e22dc2SBarry Smith }
215*61e22dc2SBarry Smith 
216*61e22dc2SBarry Smith 
217*61e22dc2SBarry Smith #undef __FUNC__
218*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSize_SeqBAIJ"
219*61e22dc2SBarry Smith int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
220*61e22dc2SBarry Smith {
221*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
222*61e22dc2SBarry Smith 
223*61e22dc2SBarry Smith   PetscFunctionBegin;
224*61e22dc2SBarry Smith   if (m) *m = a->m;
225*61e22dc2SBarry Smith   if (n) *n = a->n;
226*61e22dc2SBarry Smith   PetscFunctionReturn(0);
227*61e22dc2SBarry Smith }
228*61e22dc2SBarry Smith 
229*61e22dc2SBarry Smith #undef __FUNC__
230*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_SeqBAIJ"
231*61e22dc2SBarry Smith int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
232*61e22dc2SBarry Smith {
233*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
234*61e22dc2SBarry Smith 
235*61e22dc2SBarry Smith   PetscFunctionBegin;
236*61e22dc2SBarry Smith   if (m) *m = 0;
237*61e22dc2SBarry Smith   if (n) *n = a->m;
238*61e22dc2SBarry Smith   PetscFunctionReturn(0);
239*61e22dc2SBarry Smith }
240*61e22dc2SBarry Smith 
241*61e22dc2SBarry Smith #undef __FUNC__
242*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_SeqBAIJ"
243*61e22dc2SBarry Smith int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
244*61e22dc2SBarry Smith {
245*61e22dc2SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
246*61e22dc2SBarry Smith   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
247*61e22dc2SBarry Smith   MatScalar    *aa,*aa_i;
248*61e22dc2SBarry Smith   Scalar       *v_i;
249*61e22dc2SBarry Smith 
250*61e22dc2SBarry Smith   PetscFunctionBegin;
251*61e22dc2SBarry Smith   bs  = a->bs;
252*61e22dc2SBarry Smith   ai  = a->i;
253*61e22dc2SBarry Smith   aj  = a->j;
254*61e22dc2SBarry Smith   aa  = a->a;
255*61e22dc2SBarry Smith   bs2 = a->bs2;
256*61e22dc2SBarry Smith 
257*61e22dc2SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
258*61e22dc2SBarry Smith 
259*61e22dc2SBarry Smith   bn  = row/bs;   /* Block number */
260*61e22dc2SBarry Smith   bp  = row % bs; /* Block Position */
261*61e22dc2SBarry Smith   M   = ai[bn+1] - ai[bn];
262*61e22dc2SBarry Smith   *nz = bs*M;
263*61e22dc2SBarry Smith 
264*61e22dc2SBarry Smith   if (v) {
265*61e22dc2SBarry Smith     *v = 0;
266*61e22dc2SBarry Smith     if (*nz) {
267*61e22dc2SBarry Smith       *v = (Scalar*)PetscMalloc((*nz)*sizeof(Scalar));CHKPTRQ(*v);
268*61e22dc2SBarry Smith       for (i=0; i<M; i++) { /* for each block in the block row */
269*61e22dc2SBarry Smith         v_i  = *v + i*bs;
270*61e22dc2SBarry Smith         aa_i = aa + bs2*(ai[bn] + i);
271*61e22dc2SBarry Smith         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
272*61e22dc2SBarry Smith       }
273*61e22dc2SBarry Smith     }
274*61e22dc2SBarry Smith   }
275*61e22dc2SBarry Smith 
276*61e22dc2SBarry Smith   if (idx) {
277*61e22dc2SBarry Smith     *idx = 0;
278*61e22dc2SBarry Smith     if (*nz) {
279*61e22dc2SBarry Smith       *idx = (int*)PetscMalloc((*nz)*sizeof(int));CHKPTRQ(*idx);
280*61e22dc2SBarry Smith       for (i=0; i<M; i++) { /* for each block in the block row */
281*61e22dc2SBarry Smith         idx_i = *idx + i*bs;
282*61e22dc2SBarry Smith         itmp  = bs*aj[ai[bn] + i];
283*61e22dc2SBarry Smith         for (j=0; j<bs; j++) {idx_i[j] = itmp++;}
284*61e22dc2SBarry Smith       }
285*61e22dc2SBarry Smith     }
286*61e22dc2SBarry Smith   }
287*61e22dc2SBarry Smith   PetscFunctionReturn(0);
288*61e22dc2SBarry Smith }
289*61e22dc2SBarry Smith 
290*61e22dc2SBarry Smith #undef __FUNC__
291*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_SeqBAIJ"
292*61e22dc2SBarry Smith int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
293*61e22dc2SBarry Smith {
294*61e22dc2SBarry Smith   int ierr;
295*61e22dc2SBarry Smith 
296*61e22dc2SBarry Smith   PetscFunctionBegin;
297*61e22dc2SBarry Smith   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
298*61e22dc2SBarry Smith   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
299*61e22dc2SBarry Smith   PetscFunctionReturn(0);
300*61e22dc2SBarry Smith }
301*61e22dc2SBarry Smith 
302*61e22dc2SBarry Smith #undef __FUNC__
303*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatTranspose_SeqBAIJ"
304*61e22dc2SBarry Smith int MatTranspose_SeqBAIJ(Mat A,Mat *B)
305*61e22dc2SBarry Smith {
306*61e22dc2SBarry Smith   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
307*61e22dc2SBarry Smith   Mat         C;
308*61e22dc2SBarry Smith   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
309*61e22dc2SBarry Smith   int         *rows,*cols,bs2=a->bs2,refct;
310*61e22dc2SBarry Smith   Scalar      *array;
311*61e22dc2SBarry Smith 
312*61e22dc2SBarry Smith   PetscFunctionBegin;
313*61e22dc2SBarry Smith   if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
314*61e22dc2SBarry Smith   col  = (int*)PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col);
315*61e22dc2SBarry Smith   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
316*61e22dc2SBarry Smith 
317*61e22dc2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
318*61e22dc2SBarry Smith   array = (Scalar *)PetscMalloc(a->bs2*a->nz*sizeof(Scalar));CHKPTRQ(array);
319*61e22dc2SBarry Smith   for (i=0; i<a->bs2*a->nz; i++) array[i] = (Scalar)a->a[i];
320*61e22dc2SBarry Smith #else
321*61e22dc2SBarry Smith   array = a->a;
322*61e22dc2SBarry Smith #endif
323*61e22dc2SBarry Smith 
324*61e22dc2SBarry Smith   for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1;
325*61e22dc2SBarry Smith   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
326*61e22dc2SBarry Smith   ierr = PetscFree(col);CHKERRQ(ierr);
327*61e22dc2SBarry Smith   rows = (int*)PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows);
328*61e22dc2SBarry Smith   cols = rows + bs;
329*61e22dc2SBarry Smith   for (i=0; i<mbs; i++) {
330*61e22dc2SBarry Smith     cols[0] = i*bs;
331*61e22dc2SBarry Smith     for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1;
332*61e22dc2SBarry Smith     len = ai[i+1] - ai[i];
333*61e22dc2SBarry Smith     for (j=0; j<len; j++) {
334*61e22dc2SBarry Smith       rows[0] = (*aj++)*bs;
335*61e22dc2SBarry Smith       for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1;
336*61e22dc2SBarry Smith       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
337*61e22dc2SBarry Smith       array += bs2;
338*61e22dc2SBarry Smith     }
339*61e22dc2SBarry Smith   }
340*61e22dc2SBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
341*61e22dc2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
342*61e22dc2SBarry Smith   ierr = PetscFree(array);
343*61e22dc2SBarry Smith #endif
344*61e22dc2SBarry Smith 
345*61e22dc2SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
346*61e22dc2SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
347*61e22dc2SBarry Smith 
348*61e22dc2SBarry Smith   if (B) {
349*61e22dc2SBarry Smith     *B = C;
350*61e22dc2SBarry Smith   } else {
351*61e22dc2SBarry Smith     PetscOps *Abops;
352*61e22dc2SBarry Smith     MatOps   Aops;
353*61e22dc2SBarry Smith 
354*61e22dc2SBarry Smith     /* This isn't really an in-place transpose */
355*61e22dc2SBarry Smith     ierr = PetscFree(a->a);CHKERRQ(ierr);
356*61e22dc2SBarry Smith     if (!a->singlemalloc) {
357*61e22dc2SBarry Smith       ierr = PetscFree(a->i);CHKERRQ(ierr);
358*61e22dc2SBarry Smith       ierr = PetscFree(a->j);CHKERRQ(ierr);
359*61e22dc2SBarry Smith     }
360*61e22dc2SBarry Smith     if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
361*61e22dc2SBarry Smith     if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
362*61e22dc2SBarry Smith     if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
363*61e22dc2SBarry Smith     if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
364*61e22dc2SBarry Smith     ierr = PetscFree(a);CHKERRQ(ierr);
365*61e22dc2SBarry Smith 
366*61e22dc2SBarry Smith 
367*61e22dc2SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
368*61e22dc2SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
369*61e22dc2SBarry Smith 
370*61e22dc2SBarry Smith     /*
371*61e22dc2SBarry Smith        This is horrible, horrible code. We need to keep the
372*61e22dc2SBarry Smith       A pointers for the bops and ops but copy everything
373*61e22dc2SBarry Smith       else from C.
374*61e22dc2SBarry Smith     */
375*61e22dc2SBarry Smith     Abops    = A->bops;
376*61e22dc2SBarry Smith     Aops     = A->ops;
377*61e22dc2SBarry Smith     refct    = A->refct;
378*61e22dc2SBarry Smith     ierr     = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
379*61e22dc2SBarry Smith     A->bops  = Abops;
380*61e22dc2SBarry Smith     A->ops   = Aops;
381*61e22dc2SBarry Smith     A->qlist = 0;
382*61e22dc2SBarry Smith     A->refct = refct;
383*61e22dc2SBarry Smith 
384*61e22dc2SBarry Smith     PetscHeaderDestroy(C);
385*61e22dc2SBarry Smith   }
386*61e22dc2SBarry Smith   PetscFunctionReturn(0);
387*61e22dc2SBarry Smith }
388*61e22dc2SBarry Smith 
389*61e22dc2SBarry Smith #undef __FUNC__
390*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Binary"
391*61e22dc2SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
392*61e22dc2SBarry Smith {
393*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
394*61e22dc2SBarry Smith   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
395*61e22dc2SBarry Smith   Scalar      *aa;
396*61e22dc2SBarry Smith   FILE        *file;
397*61e22dc2SBarry Smith 
398*61e22dc2SBarry Smith   PetscFunctionBegin;
399*61e22dc2SBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
400*61e22dc2SBarry Smith   col_lens = (int*)PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
401*61e22dc2SBarry Smith   col_lens[0] = MAT_COOKIE;
402*61e22dc2SBarry Smith 
403*61e22dc2SBarry Smith   col_lens[1] = a->m;
404*61e22dc2SBarry Smith   col_lens[2] = a->n;
405*61e22dc2SBarry Smith   col_lens[3] = a->nz*bs2;
406*61e22dc2SBarry Smith 
407*61e22dc2SBarry Smith   /* store lengths of each row and write (including header) to file */
408*61e22dc2SBarry Smith   count = 0;
409*61e22dc2SBarry Smith   for (i=0; i<a->mbs; i++) {
410*61e22dc2SBarry Smith     for (j=0; j<bs; j++) {
411*61e22dc2SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
412*61e22dc2SBarry Smith     }
413*61e22dc2SBarry Smith   }
414*61e22dc2SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr);
415*61e22dc2SBarry Smith   ierr = PetscFree(col_lens);CHKERRQ(ierr);
416*61e22dc2SBarry Smith 
417*61e22dc2SBarry Smith   /* store column indices (zero start index) */
418*61e22dc2SBarry Smith   jj = (int*)PetscMalloc((a->nz+1)*bs2*sizeof(int));CHKPTRQ(jj);
419*61e22dc2SBarry Smith   count = 0;
420*61e22dc2SBarry Smith   for (i=0; i<a->mbs; i++) {
421*61e22dc2SBarry Smith     for (j=0; j<bs; j++) {
422*61e22dc2SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
423*61e22dc2SBarry Smith         for (l=0; l<bs; l++) {
424*61e22dc2SBarry Smith           jj[count++] = bs*a->j[k] + l;
425*61e22dc2SBarry Smith         }
426*61e22dc2SBarry Smith       }
427*61e22dc2SBarry Smith     }
428*61e22dc2SBarry Smith   }
429*61e22dc2SBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
430*61e22dc2SBarry Smith   ierr = PetscFree(jj);CHKERRQ(ierr);
431*61e22dc2SBarry Smith 
432*61e22dc2SBarry Smith   /* store nonzero values */
433*61e22dc2SBarry Smith   aa = (Scalar*)PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa);
434*61e22dc2SBarry Smith   count = 0;
435*61e22dc2SBarry Smith   for (i=0; i<a->mbs; i++) {
436*61e22dc2SBarry Smith     for (j=0; j<bs; j++) {
437*61e22dc2SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
438*61e22dc2SBarry Smith         for (l=0; l<bs; l++) {
439*61e22dc2SBarry Smith           aa[count++] = a->a[bs2*k + l*bs + j];
440*61e22dc2SBarry Smith         }
441*61e22dc2SBarry Smith       }
442*61e22dc2SBarry Smith     }
443*61e22dc2SBarry Smith   }
444*61e22dc2SBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
445*61e22dc2SBarry Smith   ierr = PetscFree(aa);CHKERRQ(ierr);
446*61e22dc2SBarry Smith 
447*61e22dc2SBarry Smith   ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
448*61e22dc2SBarry Smith   if (file) {
449*61e22dc2SBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
450*61e22dc2SBarry Smith   }
451*61e22dc2SBarry Smith   PetscFunctionReturn(0);
452*61e22dc2SBarry Smith }
453*61e22dc2SBarry Smith 
454*61e22dc2SBarry Smith #undef __FUNC__
455*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_ASCII"
456*61e22dc2SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
457*61e22dc2SBarry Smith {
458*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
459*61e22dc2SBarry Smith   int         ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2;
460*61e22dc2SBarry Smith   char        *outputname;
461*61e22dc2SBarry Smith 
462*61e22dc2SBarry Smith   PetscFunctionBegin;
463*61e22dc2SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
464*61e22dc2SBarry Smith   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
465*61e22dc2SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
466*61e22dc2SBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
467*61e22dc2SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
468*61e22dc2SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
469*61e22dc2SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
470*61e22dc2SBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
471*61e22dc2SBarry Smith     for (i=0; i<a->mbs; i++) {
472*61e22dc2SBarry Smith       for (j=0; j<bs; j++) {
473*61e22dc2SBarry Smith         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
474*61e22dc2SBarry Smith         for (k=a->i[i]; k<a->i[i+1]; k++) {
475*61e22dc2SBarry Smith           for (l=0; l<bs; l++) {
476*61e22dc2SBarry Smith #if defined(PETSC_USE_COMPLEX)
477*61e22dc2SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
478*61e22dc2SBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g + %gi",bs*a->j[k]+l,
479*61e22dc2SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
480*61e22dc2SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
481*61e22dc2SBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g - %gi",bs*a->j[k]+l,
482*61e22dc2SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
483*61e22dc2SBarry Smith             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
484*61e22dc2SBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
485*61e22dc2SBarry Smith             }
486*61e22dc2SBarry Smith #else
487*61e22dc2SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
488*61e22dc2SBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
489*61e22dc2SBarry Smith             }
490*61e22dc2SBarry Smith #endif
491*61e22dc2SBarry Smith           }
492*61e22dc2SBarry Smith         }
493*61e22dc2SBarry Smith         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
494*61e22dc2SBarry Smith       }
495*61e22dc2SBarry Smith     }
496*61e22dc2SBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
497*61e22dc2SBarry Smith   } else {
498*61e22dc2SBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
499*61e22dc2SBarry Smith     for (i=0; i<a->mbs; i++) {
500*61e22dc2SBarry Smith       for (j=0; j<bs; j++) {
501*61e22dc2SBarry Smith         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
502*61e22dc2SBarry Smith         for (k=a->i[i]; k<a->i[i+1]; k++) {
503*61e22dc2SBarry Smith           for (l=0; l<bs; l++) {
504*61e22dc2SBarry Smith #if defined(PETSC_USE_COMPLEX)
505*61e22dc2SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
506*61e22dc2SBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
507*61e22dc2SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
508*61e22dc2SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
509*61e22dc2SBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
510*61e22dc2SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
511*61e22dc2SBarry Smith             } else {
512*61e22dc2SBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
513*61e22dc2SBarry Smith             }
514*61e22dc2SBarry Smith #else
515*61e22dc2SBarry Smith             ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
516*61e22dc2SBarry Smith #endif
517*61e22dc2SBarry Smith           }
518*61e22dc2SBarry Smith         }
519*61e22dc2SBarry Smith         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
520*61e22dc2SBarry Smith       }
521*61e22dc2SBarry Smith     }
522*61e22dc2SBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
523*61e22dc2SBarry Smith   }
524*61e22dc2SBarry Smith   ierr = ViewerFlush(viewer);CHKERRQ(ierr);
525*61e22dc2SBarry Smith   PetscFunctionReturn(0);
526*61e22dc2SBarry Smith }
527*61e22dc2SBarry Smith 
528*61e22dc2SBarry Smith #undef __FUNC__
529*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Draw_Zoom"
530*61e22dc2SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa)
531*61e22dc2SBarry Smith {
532*61e22dc2SBarry Smith   Mat          A = (Mat) Aa;
533*61e22dc2SBarry Smith   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ*)A->data;
534*61e22dc2SBarry Smith   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2;
535*61e22dc2SBarry Smith   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
536*61e22dc2SBarry Smith   MatScalar    *aa;
537*61e22dc2SBarry Smith   Viewer       viewer;
538*61e22dc2SBarry Smith 
539*61e22dc2SBarry Smith   PetscFunctionBegin;
540*61e22dc2SBarry Smith 
541*61e22dc2SBarry Smith   /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/
542*61e22dc2SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
543*61e22dc2SBarry Smith 
544*61e22dc2SBarry Smith   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
545*61e22dc2SBarry Smith 
546*61e22dc2SBarry Smith   /* loop over matrix elements drawing boxes */
547*61e22dc2SBarry Smith   color = DRAW_BLUE;
548*61e22dc2SBarry Smith   for (i=0,row=0; i<mbs; i++,row+=bs) {
549*61e22dc2SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
550*61e22dc2SBarry Smith       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
551*61e22dc2SBarry Smith       x_l = a->j[j]*bs; x_r = x_l + 1.0;
552*61e22dc2SBarry Smith       aa = a->a + j*bs2;
553*61e22dc2SBarry Smith       for (k=0; k<bs; k++) {
554*61e22dc2SBarry Smith         for (l=0; l<bs; l++) {
555*61e22dc2SBarry Smith           if (PetscRealPart(*aa++) >=  0.) continue;
556*61e22dc2SBarry Smith           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
557*61e22dc2SBarry Smith         }
558*61e22dc2SBarry Smith       }
559*61e22dc2SBarry Smith     }
560*61e22dc2SBarry Smith   }
561*61e22dc2SBarry Smith   color = DRAW_CYAN;
562*61e22dc2SBarry Smith   for (i=0,row=0; i<mbs; i++,row+=bs) {
563*61e22dc2SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
564*61e22dc2SBarry Smith       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
565*61e22dc2SBarry Smith       x_l = a->j[j]*bs; x_r = x_l + 1.0;
566*61e22dc2SBarry Smith       aa = a->a + j*bs2;
567*61e22dc2SBarry Smith       for (k=0; k<bs; k++) {
568*61e22dc2SBarry Smith         for (l=0; l<bs; l++) {
569*61e22dc2SBarry Smith           if (PetscRealPart(*aa++) != 0.) continue;
570*61e22dc2SBarry Smith           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
571*61e22dc2SBarry Smith         }
572*61e22dc2SBarry Smith       }
573*61e22dc2SBarry Smith     }
574*61e22dc2SBarry Smith   }
575*61e22dc2SBarry Smith 
576*61e22dc2SBarry Smith   color = DRAW_RED;
577*61e22dc2SBarry Smith   for (i=0,row=0; i<mbs; i++,row+=bs) {
578*61e22dc2SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
579*61e22dc2SBarry Smith       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
580*61e22dc2SBarry Smith       x_l = a->j[j]*bs; x_r = x_l + 1.0;
581*61e22dc2SBarry Smith       aa = a->a + j*bs2;
582*61e22dc2SBarry Smith       for (k=0; k<bs; k++) {
583*61e22dc2SBarry Smith         for (l=0; l<bs; l++) {
584*61e22dc2SBarry Smith           if (PetscRealPart(*aa++) <= 0.) continue;
585*61e22dc2SBarry Smith           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
586*61e22dc2SBarry Smith         }
587*61e22dc2SBarry Smith       }
588*61e22dc2SBarry Smith     }
589*61e22dc2SBarry Smith   }
590*61e22dc2SBarry Smith   PetscFunctionReturn(0);
591*61e22dc2SBarry Smith }
592*61e22dc2SBarry Smith 
593*61e22dc2SBarry Smith #undef __FUNC__
594*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Draw"
595*61e22dc2SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
596*61e22dc2SBarry Smith {
597*61e22dc2SBarry Smith   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ*)A->data;
598*61e22dc2SBarry Smith   int          ierr;
599*61e22dc2SBarry Smith   PetscReal    xl,yl,xr,yr,w,h;
600*61e22dc2SBarry Smith   Draw         draw;
601*61e22dc2SBarry Smith   PetscTruth   isnull;
602*61e22dc2SBarry Smith 
603*61e22dc2SBarry Smith   PetscFunctionBegin;
604*61e22dc2SBarry Smith 
605*61e22dc2SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
606*61e22dc2SBarry Smith   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
607*61e22dc2SBarry Smith 
608*61e22dc2SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
609*61e22dc2SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
610*61e22dc2SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
611*61e22dc2SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
612*61e22dc2SBarry Smith   ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
613*61e22dc2SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
614*61e22dc2SBarry Smith   PetscFunctionReturn(0);
615*61e22dc2SBarry Smith }
616*61e22dc2SBarry Smith 
617*61e22dc2SBarry Smith #undef __FUNC__
618*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ"
619*61e22dc2SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer)
620*61e22dc2SBarry Smith {
621*61e22dc2SBarry Smith   int        ierr;
622*61e22dc2SBarry Smith   PetscTruth issocket,isascii,isbinary,isdraw;
623*61e22dc2SBarry Smith 
624*61e22dc2SBarry Smith   PetscFunctionBegin;
625*61e22dc2SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
626*61e22dc2SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
627*61e22dc2SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
628*61e22dc2SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
629*61e22dc2SBarry Smith   if (issocket) {
630*61e22dc2SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported");
631*61e22dc2SBarry Smith   } else if (isascii){
632*61e22dc2SBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
633*61e22dc2SBarry Smith   } else if (isbinary) {
634*61e22dc2SBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
635*61e22dc2SBarry Smith   } else if (isdraw) {
636*61e22dc2SBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
637*61e22dc2SBarry Smith   } else {
638*61e22dc2SBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
639*61e22dc2SBarry Smith   }
640*61e22dc2SBarry Smith   PetscFunctionReturn(0);
641*61e22dc2SBarry Smith }
642*61e22dc2SBarry Smith 
643*61e22dc2SBarry Smith 
644*61e22dc2SBarry Smith #undef __FUNC__
645*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetValues_SeqBAIJ"
646*61e22dc2SBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
647*61e22dc2SBarry Smith {
648*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
649*61e22dc2SBarry Smith   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
650*61e22dc2SBarry Smith   int        *ai = a->i,*ailen = a->ilen;
651*61e22dc2SBarry Smith   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
652*61e22dc2SBarry Smith   MatScalar  *ap,*aa = a->a,zero = 0.0;
653*61e22dc2SBarry Smith 
654*61e22dc2SBarry Smith   PetscFunctionBegin;
655*61e22dc2SBarry Smith   for (k=0; k<m; k++) { /* loop over rows */
656*61e22dc2SBarry Smith     row  = im[k]; brow = row/bs;
657*61e22dc2SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
658*61e22dc2SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
659*61e22dc2SBarry Smith     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
660*61e22dc2SBarry Smith     nrow = ailen[brow];
661*61e22dc2SBarry Smith     for (l=0; l<n; l++) { /* loop over columns */
662*61e22dc2SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
663*61e22dc2SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
664*61e22dc2SBarry Smith       col  = in[l] ;
665*61e22dc2SBarry Smith       bcol = col/bs;
666*61e22dc2SBarry Smith       cidx = col%bs;
667*61e22dc2SBarry Smith       ridx = row%bs;
668*61e22dc2SBarry Smith       high = nrow;
669*61e22dc2SBarry Smith       low  = 0; /* assume unsorted */
670*61e22dc2SBarry Smith       while (high-low > 5) {
671*61e22dc2SBarry Smith         t = (low+high)/2;
672*61e22dc2SBarry Smith         if (rp[t] > bcol) high = t;
673*61e22dc2SBarry Smith         else             low  = t;
674*61e22dc2SBarry Smith       }
675*61e22dc2SBarry Smith       for (i=low; i<high; i++) {
676*61e22dc2SBarry Smith         if (rp[i] > bcol) break;
677*61e22dc2SBarry Smith         if (rp[i] == bcol) {
678*61e22dc2SBarry Smith           *v++ = ap[bs2*i+bs*cidx+ridx];
679*61e22dc2SBarry Smith           goto finished;
680*61e22dc2SBarry Smith         }
681*61e22dc2SBarry Smith       }
682*61e22dc2SBarry Smith       *v++ = zero;
683*61e22dc2SBarry Smith       finished:;
684*61e22dc2SBarry Smith     }
685*61e22dc2SBarry Smith   }
686*61e22dc2SBarry Smith   PetscFunctionReturn(0);
687*61e22dc2SBarry Smith }
688*61e22dc2SBarry Smith 
689*61e22dc2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
690*61e22dc2SBarry Smith #undef __FUNC__
691*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValuesBlocked_SeqBAIJ"
692*61e22dc2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv)
693*61e22dc2SBarry Smith {
694*61e22dc2SBarry Smith   Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data;
695*61e22dc2SBarry Smith   int         ierr,i,N = m*n*b->bs2;
696*61e22dc2SBarry Smith   MatScalar   *vsingle;
697*61e22dc2SBarry Smith 
698*61e22dc2SBarry Smith   PetscFunctionBegin;
699*61e22dc2SBarry Smith   if (N > b->setvalueslen) {
700*61e22dc2SBarry Smith     if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);}
701*61e22dc2SBarry Smith     b->setvaluescopy = (MatScalar*)PetscMalloc(N*sizeof(MatScalar));CHKPTRQ(b->setvaluescopy);
702*61e22dc2SBarry Smith     b->setvalueslen  = N;
703*61e22dc2SBarry Smith   }
704*61e22dc2SBarry Smith   vsingle = b->setvaluescopy;
705*61e22dc2SBarry Smith   for (i=0; i<N; i++) {
706*61e22dc2SBarry Smith     vsingle[i] = v[i];
707*61e22dc2SBarry Smith   }
708*61e22dc2SBarry Smith   ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr);
709*61e22dc2SBarry Smith   PetscFunctionReturn(0);
710*61e22dc2SBarry Smith }
711*61e22dc2SBarry Smith #endif
712*61e22dc2SBarry Smith 
713*61e22dc2SBarry Smith 
714*61e22dc2SBarry Smith #undef __FUNC__
715*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValuesBlocked_SeqBAIJ"
716*61e22dc2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is)
717*61e22dc2SBarry Smith {
718*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
719*61e22dc2SBarry Smith   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
720*61e22dc2SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
721*61e22dc2SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
722*61e22dc2SBarry Smith   MatScalar      *value = v,*ap,*aa = a->a,*bap;
723*61e22dc2SBarry Smith 
724*61e22dc2SBarry Smith   PetscFunctionBegin;
725*61e22dc2SBarry Smith   if (roworiented) {
726*61e22dc2SBarry Smith     stepval = (n-1)*bs;
727*61e22dc2SBarry Smith   } else {
728*61e22dc2SBarry Smith     stepval = (m-1)*bs;
729*61e22dc2SBarry Smith   }
730*61e22dc2SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
731*61e22dc2SBarry Smith     row  = im[k];
732*61e22dc2SBarry Smith     if (row < 0) continue;
733*61e22dc2SBarry Smith #if defined(PETSC_USE_BOPT_g)
734*61e22dc2SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
735*61e22dc2SBarry Smith #endif
736*61e22dc2SBarry Smith     rp   = aj + ai[row];
737*61e22dc2SBarry Smith     ap   = aa + bs2*ai[row];
738*61e22dc2SBarry Smith     rmax = imax[row];
739*61e22dc2SBarry Smith     nrow = ailen[row];
740*61e22dc2SBarry Smith     low  = 0;
741*61e22dc2SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
742*61e22dc2SBarry Smith       if (in[l] < 0) continue;
743*61e22dc2SBarry Smith #if defined(PETSC_USE_BOPT_g)
744*61e22dc2SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
745*61e22dc2SBarry Smith #endif
746*61e22dc2SBarry Smith       col = in[l];
747*61e22dc2SBarry Smith       if (roworiented) {
748*61e22dc2SBarry Smith         value = v + k*(stepval+bs)*bs + l*bs;
749*61e22dc2SBarry Smith       } else {
750*61e22dc2SBarry Smith         value = v + l*(stepval+bs)*bs + k*bs;
751*61e22dc2SBarry Smith       }
752*61e22dc2SBarry Smith       if (!sorted) low = 0; high = nrow;
753*61e22dc2SBarry Smith       while (high-low > 7) {
754*61e22dc2SBarry Smith         t = (low+high)/2;
755*61e22dc2SBarry Smith         if (rp[t] > col) high = t;
756*61e22dc2SBarry Smith         else             low  = t;
757*61e22dc2SBarry Smith       }
758*61e22dc2SBarry Smith       for (i=low; i<high; i++) {
759*61e22dc2SBarry Smith         if (rp[i] > col) break;
760*61e22dc2SBarry Smith         if (rp[i] == col) {
761*61e22dc2SBarry Smith           bap  = ap +  bs2*i;
762*61e22dc2SBarry Smith           if (roworiented) {
763*61e22dc2SBarry Smith             if (is == ADD_VALUES) {
764*61e22dc2SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
765*61e22dc2SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
766*61e22dc2SBarry Smith                   bap[jj] += *value++;
767*61e22dc2SBarry Smith                 }
768*61e22dc2SBarry Smith               }
769*61e22dc2SBarry Smith             } else {
770*61e22dc2SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
771*61e22dc2SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
772*61e22dc2SBarry Smith                   bap[jj] = *value++;
773*61e22dc2SBarry Smith                 }
774*61e22dc2SBarry Smith               }
775*61e22dc2SBarry Smith             }
776*61e22dc2SBarry Smith           } else {
777*61e22dc2SBarry Smith             if (is == ADD_VALUES) {
778*61e22dc2SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
779*61e22dc2SBarry Smith                 for (jj=0; jj<bs; jj++) {
780*61e22dc2SBarry Smith                   *bap++ += *value++;
781*61e22dc2SBarry Smith                 }
782*61e22dc2SBarry Smith               }
783*61e22dc2SBarry Smith             } else {
784*61e22dc2SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
785*61e22dc2SBarry Smith                 for (jj=0; jj<bs; jj++) {
786*61e22dc2SBarry Smith                   *bap++  = *value++;
787*61e22dc2SBarry Smith                 }
788*61e22dc2SBarry Smith               }
789*61e22dc2SBarry Smith             }
790*61e22dc2SBarry Smith           }
791*61e22dc2SBarry Smith           goto noinsert2;
792*61e22dc2SBarry Smith         }
793*61e22dc2SBarry Smith       }
794*61e22dc2SBarry Smith       if (nonew == 1) goto noinsert2;
795*61e22dc2SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
796*61e22dc2SBarry Smith       if (nrow >= rmax) {
797*61e22dc2SBarry Smith         /* there is no extra room in row, therefore enlarge */
798*61e22dc2SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
799*61e22dc2SBarry Smith         MatScalar *new_a;
800*61e22dc2SBarry Smith 
801*61e22dc2SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
802*61e22dc2SBarry Smith 
803*61e22dc2SBarry Smith         /* malloc new storage space */
804*61e22dc2SBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
805*61e22dc2SBarry Smith         new_a   = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a);
806*61e22dc2SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
807*61e22dc2SBarry Smith         new_i   = new_j + new_nz;
808*61e22dc2SBarry Smith 
809*61e22dc2SBarry Smith         /* copy over old data into new slots */
810*61e22dc2SBarry Smith         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
811*61e22dc2SBarry Smith         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
812*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
813*61e22dc2SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
814*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
815*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
816*61e22dc2SBarry Smith         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
817*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
818*61e22dc2SBarry Smith         /* free up old matrix storage */
819*61e22dc2SBarry Smith         ierr = PetscFree(a->a);CHKERRQ(ierr);
820*61e22dc2SBarry Smith         if (!a->singlemalloc) {
821*61e22dc2SBarry Smith           ierr = PetscFree(a->i);CHKERRQ(ierr);
822*61e22dc2SBarry Smith           ierr = PetscFree(a->j);CHKERRQ(ierr);
823*61e22dc2SBarry Smith         }
824*61e22dc2SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
825*61e22dc2SBarry Smith         a->singlemalloc = PETSC_TRUE;
826*61e22dc2SBarry Smith 
827*61e22dc2SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
828*61e22dc2SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
829*61e22dc2SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
830*61e22dc2SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
831*61e22dc2SBarry Smith         a->reallocs++;
832*61e22dc2SBarry Smith         a->nz++;
833*61e22dc2SBarry Smith       }
834*61e22dc2SBarry Smith       N = nrow++ - 1;
835*61e22dc2SBarry Smith       /* shift up all the later entries in this row */
836*61e22dc2SBarry Smith       for (ii=N; ii>=i; ii--) {
837*61e22dc2SBarry Smith         rp[ii+1] = rp[ii];
838*61e22dc2SBarry Smith         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
839*61e22dc2SBarry Smith       }
840*61e22dc2SBarry Smith       if (N >= i) {
841*61e22dc2SBarry Smith         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
842*61e22dc2SBarry Smith       }
843*61e22dc2SBarry Smith       rp[i] = col;
844*61e22dc2SBarry Smith       bap   = ap +  bs2*i;
845*61e22dc2SBarry Smith       if (roworiented) {
846*61e22dc2SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
847*61e22dc2SBarry Smith           for (jj=ii; jj<bs2; jj+=bs) {
848*61e22dc2SBarry Smith             bap[jj] = *value++;
849*61e22dc2SBarry Smith           }
850*61e22dc2SBarry Smith         }
851*61e22dc2SBarry Smith       } else {
852*61e22dc2SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
853*61e22dc2SBarry Smith           for (jj=0; jj<bs; jj++) {
854*61e22dc2SBarry Smith             *bap++  = *value++;
855*61e22dc2SBarry Smith           }
856*61e22dc2SBarry Smith         }
857*61e22dc2SBarry Smith       }
858*61e22dc2SBarry Smith       noinsert2:;
859*61e22dc2SBarry Smith       low = i;
860*61e22dc2SBarry Smith     }
861*61e22dc2SBarry Smith     ailen[row] = nrow;
862*61e22dc2SBarry Smith   }
863*61e22dc2SBarry Smith   PetscFunctionReturn(0);
864*61e22dc2SBarry Smith }
865*61e22dc2SBarry Smith 
866*61e22dc2SBarry Smith 
867*61e22dc2SBarry Smith #undef __FUNC__
868*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatAssemblyEnd_SeqBAIJ"
869*61e22dc2SBarry Smith int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
870*61e22dc2SBarry Smith {
871*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
872*61e22dc2SBarry Smith   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
873*61e22dc2SBarry Smith   int        m = a->m,*ip,N,*ailen = a->ilen;
874*61e22dc2SBarry Smith   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
875*61e22dc2SBarry Smith   MatScalar  *aa = a->a,*ap;
876*61e22dc2SBarry Smith 
877*61e22dc2SBarry Smith   PetscFunctionBegin;
878*61e22dc2SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
879*61e22dc2SBarry Smith 
880*61e22dc2SBarry Smith   if (m) rmax = ailen[0];
881*61e22dc2SBarry Smith   for (i=1; i<mbs; i++) {
882*61e22dc2SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
883*61e22dc2SBarry Smith     fshift += imax[i-1] - ailen[i-1];
884*61e22dc2SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
885*61e22dc2SBarry Smith     if (fshift) {
886*61e22dc2SBarry Smith       ip = aj + ai[i]; ap = aa + bs2*ai[i];
887*61e22dc2SBarry Smith       N = ailen[i];
888*61e22dc2SBarry Smith       for (j=0; j<N; j++) {
889*61e22dc2SBarry Smith         ip[j-fshift] = ip[j];
890*61e22dc2SBarry Smith         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
891*61e22dc2SBarry Smith       }
892*61e22dc2SBarry Smith     }
893*61e22dc2SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
894*61e22dc2SBarry Smith   }
895*61e22dc2SBarry Smith   if (mbs) {
896*61e22dc2SBarry Smith     fshift += imax[mbs-1] - ailen[mbs-1];
897*61e22dc2SBarry Smith     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
898*61e22dc2SBarry Smith   }
899*61e22dc2SBarry Smith   /* reset ilen and imax for each row */
900*61e22dc2SBarry Smith   for (i=0; i<mbs; i++) {
901*61e22dc2SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
902*61e22dc2SBarry Smith   }
903*61e22dc2SBarry Smith   a->nz = ai[mbs];
904*61e22dc2SBarry Smith 
905*61e22dc2SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
906*61e22dc2SBarry Smith   if (fshift && a->diag) {
907*61e22dc2SBarry Smith     ierr = PetscFree(a->diag);CHKERRQ(ierr);
908*61e22dc2SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
909*61e22dc2SBarry Smith     a->diag = 0;
910*61e22dc2SBarry Smith   }
911*61e22dc2SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
912*61e22dc2SBarry Smith            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
913*61e22dc2SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
914*61e22dc2SBarry Smith            a->reallocs);
915*61e22dc2SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
916*61e22dc2SBarry Smith   a->reallocs          = 0;
917*61e22dc2SBarry Smith   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
918*61e22dc2SBarry Smith 
919*61e22dc2SBarry Smith   PetscFunctionReturn(0);
920*61e22dc2SBarry Smith }
921*61e22dc2SBarry Smith 
922*61e22dc2SBarry Smith 
923*61e22dc2SBarry Smith 
924*61e22dc2SBarry Smith /*
925*61e22dc2SBarry Smith    This function returns an array of flags which indicate the locations of contiguous
926*61e22dc2SBarry Smith    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
927*61e22dc2SBarry Smith    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
928*61e22dc2SBarry Smith    Assume: sizes should be long enough to hold all the values.
929*61e22dc2SBarry Smith */
930*61e22dc2SBarry Smith #undef __FUNC__
931*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqBAIJ_Check_Blocks"
932*61e22dc2SBarry Smith static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
933*61e22dc2SBarry Smith {
934*61e22dc2SBarry Smith   int        i,j,k,row;
935*61e22dc2SBarry Smith   PetscTruth flg;
936*61e22dc2SBarry Smith 
937*61e22dc2SBarry Smith   PetscFunctionBegin;
938*61e22dc2SBarry Smith   for (i=0,j=0; i<n; j++) {
939*61e22dc2SBarry Smith     row = idx[i];
940*61e22dc2SBarry Smith     if (row%bs!=0) { /* Not the begining of a block */
941*61e22dc2SBarry Smith       sizes[j] = 1;
942*61e22dc2SBarry Smith       i++;
943*61e22dc2SBarry Smith     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
944*61e22dc2SBarry Smith       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
945*61e22dc2SBarry Smith       i++;
946*61e22dc2SBarry Smith     } else { /* Begining of the block, so check if the complete block exists */
947*61e22dc2SBarry Smith       flg = PETSC_TRUE;
948*61e22dc2SBarry Smith       for (k=1; k<bs; k++) {
949*61e22dc2SBarry Smith         if (row+k != idx[i+k]) { /* break in the block */
950*61e22dc2SBarry Smith           flg = PETSC_FALSE;
951*61e22dc2SBarry Smith           break;
952*61e22dc2SBarry Smith         }
953*61e22dc2SBarry Smith       }
954*61e22dc2SBarry Smith       if (flg == PETSC_TRUE) { /* No break in the bs */
955*61e22dc2SBarry Smith         sizes[j] = bs;
956*61e22dc2SBarry Smith         i+= bs;
957*61e22dc2SBarry Smith       } else {
958*61e22dc2SBarry Smith         sizes[j] = 1;
959*61e22dc2SBarry Smith         i++;
960*61e22dc2SBarry Smith       }
961*61e22dc2SBarry Smith     }
962*61e22dc2SBarry Smith   }
963*61e22dc2SBarry Smith   *bs_max = j;
964*61e22dc2SBarry Smith   PetscFunctionReturn(0);
965*61e22dc2SBarry Smith }
966*61e22dc2SBarry Smith 
967*61e22dc2SBarry Smith #undef __FUNC__
968*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqBAIJ"
969*61e22dc2SBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,Scalar *diag)
970*61e22dc2SBarry Smith {
971*61e22dc2SBarry Smith   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
972*61e22dc2SBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
973*61e22dc2SBarry Smith   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
974*61e22dc2SBarry Smith   Scalar      zero = 0.0;
975*61e22dc2SBarry Smith   MatScalar   *aa;
976*61e22dc2SBarry Smith 
977*61e22dc2SBarry Smith   PetscFunctionBegin;
978*61e22dc2SBarry Smith   /* Make a copy of the IS and  sort it */
979*61e22dc2SBarry Smith   ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr);
980*61e22dc2SBarry Smith   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
981*61e22dc2SBarry Smith 
982*61e22dc2SBarry Smith   /* allocate memory for rows,sizes */
983*61e22dc2SBarry Smith   rows  = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows);
984*61e22dc2SBarry Smith   sizes = rows + is_n;
985*61e22dc2SBarry Smith 
986*61e22dc2SBarry Smith   /* copy IS values to rows, and sort them */
987*61e22dc2SBarry Smith   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
988*61e22dc2SBarry Smith   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
989*61e22dc2SBarry Smith   if (baij->keepzeroedrows) {
990*61e22dc2SBarry Smith     for (i=0; i<is_n; i++) { sizes[i] = 1; }
991*61e22dc2SBarry Smith     bs_max = is_n;
992*61e22dc2SBarry Smith   } else {
993*61e22dc2SBarry Smith     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
994*61e22dc2SBarry Smith   }
995*61e22dc2SBarry Smith   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
996*61e22dc2SBarry Smith 
997*61e22dc2SBarry Smith   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
998*61e22dc2SBarry Smith     row   = rows[j];
999*61e22dc2SBarry Smith     if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row);
1000*61e22dc2SBarry Smith     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1001*61e22dc2SBarry Smith     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
1002*61e22dc2SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
1003*61e22dc2SBarry Smith       if (diag) {
1004*61e22dc2SBarry Smith         if (baij->ilen[row/bs] > 0) {
1005*61e22dc2SBarry Smith           baij->ilen[row/bs]       = 1;
1006*61e22dc2SBarry Smith           baij->j[baij->i[row/bs]] = row/bs;
1007*61e22dc2SBarry Smith           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
1008*61e22dc2SBarry Smith         }
1009*61e22dc2SBarry Smith         /* Now insert all the diagonal values for this bs */
1010*61e22dc2SBarry Smith         for (k=0; k<bs; k++) {
1011*61e22dc2SBarry Smith           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
1012*61e22dc2SBarry Smith         }
1013*61e22dc2SBarry Smith       } else { /* (!diag) */
1014*61e22dc2SBarry Smith         baij->ilen[row/bs] = 0;
1015*61e22dc2SBarry Smith       } /* end (!diag) */
1016*61e22dc2SBarry Smith     } else { /* (sizes[i] != bs) */
1017*61e22dc2SBarry Smith #if defined (PETSC_USE_DEBUG)
1018*61e22dc2SBarry Smith       if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1");
1019*61e22dc2SBarry Smith #endif
1020*61e22dc2SBarry Smith       for (k=0; k<count; k++) {
1021*61e22dc2SBarry Smith         aa[0] =  zero;
1022*61e22dc2SBarry Smith         aa    += bs;
1023*61e22dc2SBarry Smith       }
1024*61e22dc2SBarry Smith       if (diag) {
1025*61e22dc2SBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
1026*61e22dc2SBarry Smith       }
1027*61e22dc2SBarry Smith     }
1028*61e22dc2SBarry Smith   }
1029*61e22dc2SBarry Smith 
1030*61e22dc2SBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
1031*61e22dc2SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1032*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1033*61e22dc2SBarry Smith }
1034*61e22dc2SBarry Smith 
1035*61e22dc2SBarry Smith #undef __FUNC__
1036*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValues_SeqBAIJ"
1037*61e22dc2SBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
1038*61e22dc2SBarry Smith {
1039*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1040*61e22dc2SBarry Smith   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
1041*61e22dc2SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
1042*61e22dc2SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1043*61e22dc2SBarry Smith   int         ridx,cidx,bs2=a->bs2,ierr;
1044*61e22dc2SBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
1045*61e22dc2SBarry Smith 
1046*61e22dc2SBarry Smith   PetscFunctionBegin;
1047*61e22dc2SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
1048*61e22dc2SBarry Smith     row  = im[k]; brow = row/bs;
1049*61e22dc2SBarry Smith     if (row < 0) continue;
1050*61e22dc2SBarry Smith #if defined(PETSC_USE_BOPT_g)
1051*61e22dc2SBarry Smith     if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m);
1052*61e22dc2SBarry Smith #endif
1053*61e22dc2SBarry Smith     rp   = aj + ai[brow];
1054*61e22dc2SBarry Smith     ap   = aa + bs2*ai[brow];
1055*61e22dc2SBarry Smith     rmax = imax[brow];
1056*61e22dc2SBarry Smith     nrow = ailen[brow];
1057*61e22dc2SBarry Smith     low  = 0;
1058*61e22dc2SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
1059*61e22dc2SBarry Smith       if (in[l] < 0) continue;
1060*61e22dc2SBarry Smith #if defined(PETSC_USE_BOPT_g)
1061*61e22dc2SBarry Smith       if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n);
1062*61e22dc2SBarry Smith #endif
1063*61e22dc2SBarry Smith       col = in[l]; bcol = col/bs;
1064*61e22dc2SBarry Smith       ridx = row % bs; cidx = col % bs;
1065*61e22dc2SBarry Smith       if (roworiented) {
1066*61e22dc2SBarry Smith         value = v[l + k*n];
1067*61e22dc2SBarry Smith       } else {
1068*61e22dc2SBarry Smith         value = v[k + l*m];
1069*61e22dc2SBarry Smith       }
1070*61e22dc2SBarry Smith       if (!sorted) low = 0; high = nrow;
1071*61e22dc2SBarry Smith       while (high-low > 7) {
1072*61e22dc2SBarry Smith         t = (low+high)/2;
1073*61e22dc2SBarry Smith         if (rp[t] > bcol) high = t;
1074*61e22dc2SBarry Smith         else              low  = t;
1075*61e22dc2SBarry Smith       }
1076*61e22dc2SBarry Smith       for (i=low; i<high; i++) {
1077*61e22dc2SBarry Smith         if (rp[i] > bcol) break;
1078*61e22dc2SBarry Smith         if (rp[i] == bcol) {
1079*61e22dc2SBarry Smith           bap  = ap +  bs2*i + bs*cidx + ridx;
1080*61e22dc2SBarry Smith           if (is == ADD_VALUES) *bap += value;
1081*61e22dc2SBarry Smith           else                  *bap  = value;
1082*61e22dc2SBarry Smith           goto noinsert1;
1083*61e22dc2SBarry Smith         }
1084*61e22dc2SBarry Smith       }
1085*61e22dc2SBarry Smith       if (nonew == 1) goto noinsert1;
1086*61e22dc2SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
1087*61e22dc2SBarry Smith       if (nrow >= rmax) {
1088*61e22dc2SBarry Smith         /* there is no extra room in row, therefore enlarge */
1089*61e22dc2SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
1090*61e22dc2SBarry Smith         MatScalar *new_a;
1091*61e22dc2SBarry Smith 
1092*61e22dc2SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
1093*61e22dc2SBarry Smith 
1094*61e22dc2SBarry Smith         /* Malloc new storage space */
1095*61e22dc2SBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
1096*61e22dc2SBarry Smith         new_a   = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a);
1097*61e22dc2SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
1098*61e22dc2SBarry Smith         new_i   = new_j + new_nz;
1099*61e22dc2SBarry Smith 
1100*61e22dc2SBarry Smith         /* copy over old data into new slots */
1101*61e22dc2SBarry Smith         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
1102*61e22dc2SBarry Smith         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1103*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
1104*61e22dc2SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1105*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1106*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1107*61e22dc2SBarry Smith         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1108*61e22dc2SBarry Smith         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
1109*61e22dc2SBarry Smith         /* free up old matrix storage */
1110*61e22dc2SBarry Smith         ierr = PetscFree(a->a);CHKERRQ(ierr);
1111*61e22dc2SBarry Smith         if (!a->singlemalloc) {
1112*61e22dc2SBarry Smith           ierr = PetscFree(a->i);CHKERRQ(ierr);
1113*61e22dc2SBarry Smith           ierr = PetscFree(a->j);CHKERRQ(ierr);
1114*61e22dc2SBarry Smith         }
1115*61e22dc2SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1116*61e22dc2SBarry Smith         a->singlemalloc = PETSC_TRUE;
1117*61e22dc2SBarry Smith 
1118*61e22dc2SBarry Smith         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
1119*61e22dc2SBarry Smith         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1120*61e22dc2SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
1121*61e22dc2SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
1122*61e22dc2SBarry Smith         a->reallocs++;
1123*61e22dc2SBarry Smith         a->nz++;
1124*61e22dc2SBarry Smith       }
1125*61e22dc2SBarry Smith       N = nrow++ - 1;
1126*61e22dc2SBarry Smith       /* shift up all the later entries in this row */
1127*61e22dc2SBarry Smith       for (ii=N; ii>=i; ii--) {
1128*61e22dc2SBarry Smith         rp[ii+1] = rp[ii];
1129*61e22dc2SBarry Smith         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
1130*61e22dc2SBarry Smith       }
1131*61e22dc2SBarry Smith       if (N>=i) {
1132*61e22dc2SBarry Smith         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1133*61e22dc2SBarry Smith       }
1134*61e22dc2SBarry Smith       rp[i]                      = bcol;
1135*61e22dc2SBarry Smith       ap[bs2*i + bs*cidx + ridx] = value;
1136*61e22dc2SBarry Smith       noinsert1:;
1137*61e22dc2SBarry Smith       low = i;
1138*61e22dc2SBarry Smith     }
1139*61e22dc2SBarry Smith     ailen[brow] = nrow;
1140*61e22dc2SBarry Smith   }
1141*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1142*61e22dc2SBarry Smith }
1143*61e22dc2SBarry Smith 
1144*61e22dc2SBarry Smith EXTERN int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatLUInfo*,Mat*);
1145*61e22dc2SBarry Smith EXTERN int MatLUFactor_SeqBAIJ(Mat,IS,IS,MatLUInfo*);
1146*61e22dc2SBarry Smith EXTERN int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
1147*61e22dc2SBarry Smith EXTERN int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
1148*61e22dc2SBarry Smith EXTERN int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
1149*61e22dc2SBarry Smith EXTERN int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
1150*61e22dc2SBarry Smith EXTERN int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
1151*61e22dc2SBarry Smith EXTERN int MatScale_SeqBAIJ(Scalar*,Mat);
1152*61e22dc2SBarry Smith EXTERN int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
1153*61e22dc2SBarry Smith EXTERN int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
1154*61e22dc2SBarry Smith EXTERN int MatGetDiagonal_SeqBAIJ(Mat,Vec);
1155*61e22dc2SBarry Smith EXTERN int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
1156*61e22dc2SBarry Smith EXTERN int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
1157*61e22dc2SBarry Smith EXTERN int MatZeroEntries_SeqBAIJ(Mat);
1158*61e22dc2SBarry Smith 
1159*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
1160*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
1161*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
1162*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
1163*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
1164*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
1165*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
1166*61e22dc2SBarry Smith EXTERN int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
1167*61e22dc2SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
1168*61e22dc2SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
1169*61e22dc2SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
1170*61e22dc2SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
1171*61e22dc2SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
1172*61e22dc2SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
1173*61e22dc2SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
1174*61e22dc2SBarry Smith 
1175*61e22dc2SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
1176*61e22dc2SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
1177*61e22dc2SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
1178*61e22dc2SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
1179*61e22dc2SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
1180*61e22dc2SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
1181*61e22dc2SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
1182*61e22dc2SBarry Smith 
1183*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
1184*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
1185*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
1186*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
1187*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
1188*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_6(Mat,Vec,Vec);
1189*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
1190*61e22dc2SBarry Smith EXTERN int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
1191*61e22dc2SBarry Smith 
1192*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
1193*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
1194*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
1195*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
1196*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
1197*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
1198*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
1199*61e22dc2SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
1200*61e22dc2SBarry Smith 
1201*61e22dc2SBarry Smith #undef __FUNC__
1202*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatILUFactor_SeqBAIJ"
1203*61e22dc2SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
1204*61e22dc2SBarry Smith {
1205*61e22dc2SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
1206*61e22dc2SBarry Smith   Mat         outA;
1207*61e22dc2SBarry Smith   int         ierr;
1208*61e22dc2SBarry Smith   PetscTruth  row_identity,col_identity;
1209*61e22dc2SBarry Smith 
1210*61e22dc2SBarry Smith   PetscFunctionBegin;
1211*61e22dc2SBarry Smith   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU");
1212*61e22dc2SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1213*61e22dc2SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1214*61e22dc2SBarry Smith   if (!row_identity || !col_identity) {
1215*61e22dc2SBarry Smith     SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU");
1216*61e22dc2SBarry Smith   }
1217*61e22dc2SBarry Smith 
1218*61e22dc2SBarry Smith   outA          = inA;
1219*61e22dc2SBarry Smith   inA->factor   = FACTOR_LU;
1220*61e22dc2SBarry Smith 
1221*61e22dc2SBarry Smith   if (!a->diag) {
1222*61e22dc2SBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
1223*61e22dc2SBarry Smith   }
1224*61e22dc2SBarry Smith   /*
1225*61e22dc2SBarry Smith       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
1226*61e22dc2SBarry Smith       for ILU(0) factorization with natural ordering
1227*61e22dc2SBarry Smith   */
1228*61e22dc2SBarry Smith   switch (a->bs) {
1229*61e22dc2SBarry Smith   case 1:
1230*61e22dc2SBarry Smith     inA->ops->solvetranspose   = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
1231*61e22dc2SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
1232*61e22dc2SBarry Smith   case 2:
1233*61e22dc2SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
1234*61e22dc2SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_2_NaturalOrdering;
1235*61e22dc2SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
1236*61e22dc2SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
1237*61e22dc2SBarry Smith     break;
1238*61e22dc2SBarry Smith   case 3:
1239*61e22dc2SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
1240*61e22dc2SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_3_NaturalOrdering;
1241*61e22dc2SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
1242*61e22dc2SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
1243*61e22dc2SBarry Smith     break;
1244*61e22dc2SBarry Smith   case 4:
1245*61e22dc2SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
1246*61e22dc2SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
1247*61e22dc2SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering;
1248*61e22dc2SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
1249*61e22dc2SBarry Smith     break;
1250*61e22dc2SBarry Smith   case 5:
1251*61e22dc2SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
1252*61e22dc2SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_5_NaturalOrdering;
1253*61e22dc2SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
1254*61e22dc2SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
1255*61e22dc2SBarry Smith     break;
1256*61e22dc2SBarry Smith   case 6:
1257*61e22dc2SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
1258*61e22dc2SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_6_NaturalOrdering;
1259*61e22dc2SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
1260*61e22dc2SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
1261*61e22dc2SBarry Smith     break;
1262*61e22dc2SBarry Smith   case 7:
1263*61e22dc2SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
1264*61e22dc2SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering;
1265*61e22dc2SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_7_NaturalOrdering;
1266*61e22dc2SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
1267*61e22dc2SBarry Smith     break;
1268*61e22dc2SBarry Smith   default:
1269*61e22dc2SBarry Smith     a->row        = row;
1270*61e22dc2SBarry Smith     a->col        = col;
1271*61e22dc2SBarry Smith     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1272*61e22dc2SBarry Smith     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1273*61e22dc2SBarry Smith 
1274*61e22dc2SBarry Smith     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1275*61e22dc2SBarry Smith     ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr);
1276*61e22dc2SBarry Smith     PLogObjectParent(inA,a->icol);
1277*61e22dc2SBarry Smith 
1278*61e22dc2SBarry Smith     if (!a->solve_work) {
1279*61e22dc2SBarry Smith       a->solve_work = (Scalar*)PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
1280*61e22dc2SBarry Smith       PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
1281*61e22dc2SBarry Smith     }
1282*61e22dc2SBarry Smith   }
1283*61e22dc2SBarry Smith 
1284*61e22dc2SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1285*61e22dc2SBarry Smith 
1286*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1287*61e22dc2SBarry Smith }
1288*61e22dc2SBarry Smith #undef __FUNC__
1289*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatPrintHelp_SeqBAIJ"
1290*61e22dc2SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A)
1291*61e22dc2SBarry Smith {
1292*61e22dc2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1293*61e22dc2SBarry Smith   MPI_Comm          comm = A->comm;
1294*61e22dc2SBarry Smith   int               ierr;
1295*61e22dc2SBarry Smith 
1296*61e22dc2SBarry Smith   PetscFunctionBegin;
1297*61e22dc2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1298*61e22dc2SBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1299*61e22dc2SBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
1300*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1301*61e22dc2SBarry Smith }
1302*61e22dc2SBarry Smith 
1303*61e22dc2SBarry Smith EXTERN_C_BEGIN
1304*61e22dc2SBarry Smith #undef __FUNC__
1305*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqBAIJSetColumnIndices_SeqBAIJ"
1306*61e22dc2SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
1307*61e22dc2SBarry Smith {
1308*61e22dc2SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
1309*61e22dc2SBarry Smith   int         i,nz,n;
1310*61e22dc2SBarry Smith 
1311*61e22dc2SBarry Smith   PetscFunctionBegin;
1312*61e22dc2SBarry Smith   nz = baij->maxnz;
1313*61e22dc2SBarry Smith   n  = baij->n;
1314*61e22dc2SBarry Smith   for (i=0; i<nz; i++) {
1315*61e22dc2SBarry Smith     baij->j[i] = indices[i];
1316*61e22dc2SBarry Smith   }
1317*61e22dc2SBarry Smith   baij->nz = nz;
1318*61e22dc2SBarry Smith   for (i=0; i<n; i++) {
1319*61e22dc2SBarry Smith     baij->ilen[i] = baij->imax[i];
1320*61e22dc2SBarry Smith   }
1321*61e22dc2SBarry Smith 
1322*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1323*61e22dc2SBarry Smith }
1324*61e22dc2SBarry Smith EXTERN_C_END
1325*61e22dc2SBarry Smith 
1326*61e22dc2SBarry Smith #undef __FUNC__
1327*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqBAIJSetColumnIndices"
1328*61e22dc2SBarry Smith /*@
1329*61e22dc2SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
1330*61e22dc2SBarry Smith        in the matrix.
1331*61e22dc2SBarry Smith 
1332*61e22dc2SBarry Smith   Input Parameters:
1333*61e22dc2SBarry Smith +  mat - the SeqBAIJ matrix
1334*61e22dc2SBarry Smith -  indices - the column indices
1335*61e22dc2SBarry Smith 
1336*61e22dc2SBarry Smith   Level: advanced
1337*61e22dc2SBarry Smith 
1338*61e22dc2SBarry Smith   Notes:
1339*61e22dc2SBarry Smith     This can be called if you have precomputed the nonzero structure of the
1340*61e22dc2SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
1341*61e22dc2SBarry Smith   of the MatSetValues() operation.
1342*61e22dc2SBarry Smith 
1343*61e22dc2SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
1344*61e22dc2SBarry Smith   MatCreateSeqBAIJ().
1345*61e22dc2SBarry Smith 
1346*61e22dc2SBarry Smith     MUST be called before any calls to MatSetValues();
1347*61e22dc2SBarry Smith 
1348*61e22dc2SBarry Smith @*/
1349*61e22dc2SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
1350*61e22dc2SBarry Smith {
1351*61e22dc2SBarry Smith   int ierr,(*f)(Mat,int *);
1352*61e22dc2SBarry Smith 
1353*61e22dc2SBarry Smith   PetscFunctionBegin;
1354*61e22dc2SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1355*61e22dc2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr);
1356*61e22dc2SBarry Smith   if (f) {
1357*61e22dc2SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
1358*61e22dc2SBarry Smith   } else {
1359*61e22dc2SBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
1360*61e22dc2SBarry Smith   }
1361*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1362*61e22dc2SBarry Smith }
1363*61e22dc2SBarry Smith 
1364*61e22dc2SBarry Smith /* -------------------------------------------------------------------*/
1365*61e22dc2SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1366*61e22dc2SBarry Smith        MatGetRow_SeqBAIJ,
1367*61e22dc2SBarry Smith        MatRestoreRow_SeqBAIJ,
1368*61e22dc2SBarry Smith        MatMult_SeqBAIJ_N,
1369*61e22dc2SBarry Smith        MatMultAdd_SeqBAIJ_N,
1370*61e22dc2SBarry Smith        MatMultTranspose_SeqBAIJ,
1371*61e22dc2SBarry Smith        MatMultTransposeAdd_SeqBAIJ,
1372*61e22dc2SBarry Smith        MatSolve_SeqBAIJ_N,
1373*61e22dc2SBarry Smith        0,
1374*61e22dc2SBarry Smith        0,
1375*61e22dc2SBarry Smith        0,
1376*61e22dc2SBarry Smith        MatLUFactor_SeqBAIJ,
1377*61e22dc2SBarry Smith        0,
1378*61e22dc2SBarry Smith        0,
1379*61e22dc2SBarry Smith        MatTranspose_SeqBAIJ,
1380*61e22dc2SBarry Smith        MatGetInfo_SeqBAIJ,
1381*61e22dc2SBarry Smith        MatEqual_SeqBAIJ,
1382*61e22dc2SBarry Smith        MatGetDiagonal_SeqBAIJ,
1383*61e22dc2SBarry Smith        MatDiagonalScale_SeqBAIJ,
1384*61e22dc2SBarry Smith        MatNorm_SeqBAIJ,
1385*61e22dc2SBarry Smith        0,
1386*61e22dc2SBarry Smith        MatAssemblyEnd_SeqBAIJ,
1387*61e22dc2SBarry Smith        0,
1388*61e22dc2SBarry Smith        MatSetOption_SeqBAIJ,
1389*61e22dc2SBarry Smith        MatZeroEntries_SeqBAIJ,
1390*61e22dc2SBarry Smith        MatZeroRows_SeqBAIJ,
1391*61e22dc2SBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1392*61e22dc2SBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1393*61e22dc2SBarry Smith        0,
1394*61e22dc2SBarry Smith        0,
1395*61e22dc2SBarry Smith        MatGetSize_SeqBAIJ,
1396*61e22dc2SBarry Smith        MatGetSize_SeqBAIJ,
1397*61e22dc2SBarry Smith        MatGetOwnershipRange_SeqBAIJ,
1398*61e22dc2SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1399*61e22dc2SBarry Smith        0,
1400*61e22dc2SBarry Smith        0,
1401*61e22dc2SBarry Smith        0,
1402*61e22dc2SBarry Smith        MatDuplicate_SeqBAIJ,
1403*61e22dc2SBarry Smith        0,
1404*61e22dc2SBarry Smith        0,
1405*61e22dc2SBarry Smith        MatILUFactor_SeqBAIJ,
1406*61e22dc2SBarry Smith        0,
1407*61e22dc2SBarry Smith        0,
1408*61e22dc2SBarry Smith        MatGetSubMatrices_SeqBAIJ,
1409*61e22dc2SBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1410*61e22dc2SBarry Smith        MatGetValues_SeqBAIJ,
1411*61e22dc2SBarry Smith        0,
1412*61e22dc2SBarry Smith        MatPrintHelp_SeqBAIJ,
1413*61e22dc2SBarry Smith        MatScale_SeqBAIJ,
1414*61e22dc2SBarry Smith        0,
1415*61e22dc2SBarry Smith        0,
1416*61e22dc2SBarry Smith        0,
1417*61e22dc2SBarry Smith        MatGetBlockSize_SeqBAIJ,
1418*61e22dc2SBarry Smith        MatGetRowIJ_SeqBAIJ,
1419*61e22dc2SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1420*61e22dc2SBarry Smith        0,
1421*61e22dc2SBarry Smith        0,
1422*61e22dc2SBarry Smith        0,
1423*61e22dc2SBarry Smith        0,
1424*61e22dc2SBarry Smith        0,
1425*61e22dc2SBarry Smith        0,
1426*61e22dc2SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1427*61e22dc2SBarry Smith        MatGetSubMatrix_SeqBAIJ,
1428*61e22dc2SBarry Smith        MatDestroy_SeqBAIJ,
1429*61e22dc2SBarry Smith        MatView_SeqBAIJ,
1430*61e22dc2SBarry Smith        MatGetMaps_Petsc};
1431*61e22dc2SBarry Smith 
1432*61e22dc2SBarry Smith EXTERN_C_BEGIN
1433*61e22dc2SBarry Smith #undef __FUNC__
1434*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatStoreValues_SeqBAIJ"
1435*61e22dc2SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
1436*61e22dc2SBarry Smith {
1437*61e22dc2SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1438*61e22dc2SBarry Smith   int         nz = aij->i[aij->m]*aij->bs*aij->bs2;
1439*61e22dc2SBarry Smith   int         ierr;
1440*61e22dc2SBarry Smith 
1441*61e22dc2SBarry Smith   PetscFunctionBegin;
1442*61e22dc2SBarry Smith   if (aij->nonew != 1) {
1443*61e22dc2SBarry Smith     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
1444*61e22dc2SBarry Smith   }
1445*61e22dc2SBarry Smith 
1446*61e22dc2SBarry Smith   /* allocate space for values if not already there */
1447*61e22dc2SBarry Smith   if (!aij->saved_values) {
1448*61e22dc2SBarry Smith     aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values);
1449*61e22dc2SBarry Smith   }
1450*61e22dc2SBarry Smith 
1451*61e22dc2SBarry Smith   /* copy values over */
1452*61e22dc2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr);
1453*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1454*61e22dc2SBarry Smith }
1455*61e22dc2SBarry Smith EXTERN_C_END
1456*61e22dc2SBarry Smith 
1457*61e22dc2SBarry Smith EXTERN_C_BEGIN
1458*61e22dc2SBarry Smith #undef __FUNC__
1459*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRetrieveValues_SeqBAIJ"
1460*61e22dc2SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
1461*61e22dc2SBarry Smith {
1462*61e22dc2SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1463*61e22dc2SBarry Smith   int         nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr;
1464*61e22dc2SBarry Smith 
1465*61e22dc2SBarry Smith   PetscFunctionBegin;
1466*61e22dc2SBarry Smith   if (aij->nonew != 1) {
1467*61e22dc2SBarry Smith     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
1468*61e22dc2SBarry Smith   }
1469*61e22dc2SBarry Smith   if (!aij->saved_values) {
1470*61e22dc2SBarry Smith     SETERRQ(1,1,"Must call MatStoreValues(A);first");
1471*61e22dc2SBarry Smith   }
1472*61e22dc2SBarry Smith 
1473*61e22dc2SBarry Smith   /* copy values over */
1474*61e22dc2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr);
1475*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1476*61e22dc2SBarry Smith }
1477*61e22dc2SBarry Smith EXTERN_C_END
1478*61e22dc2SBarry Smith 
1479*61e22dc2SBarry Smith #undef __FUNC__
1480*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateSeqBAIJ"
1481*61e22dc2SBarry Smith /*@C
1482*61e22dc2SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
1483*61e22dc2SBarry Smith    compressed row) format.  For good matrix assembly performance the
1484*61e22dc2SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1485*61e22dc2SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1486*61e22dc2SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
1487*61e22dc2SBarry Smith 
1488*61e22dc2SBarry Smith    Collective on MPI_Comm
1489*61e22dc2SBarry Smith 
1490*61e22dc2SBarry Smith    Input Parameters:
1491*61e22dc2SBarry Smith +  comm - MPI communicator, set to PETSC_COMM_SELF
1492*61e22dc2SBarry Smith .  bs - size of block
1493*61e22dc2SBarry Smith .  m - number of rows
1494*61e22dc2SBarry Smith .  n - number of columns
1495*61e22dc2SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1496*61e22dc2SBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
1497*61e22dc2SBarry Smith          (possibly different for each block row) or PETSC_NULL
1498*61e22dc2SBarry Smith 
1499*61e22dc2SBarry Smith    Output Parameter:
1500*61e22dc2SBarry Smith .  A - the matrix
1501*61e22dc2SBarry Smith 
1502*61e22dc2SBarry Smith    Options Database Keys:
1503*61e22dc2SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1504*61e22dc2SBarry Smith                      block calculations (much slower)
1505*61e22dc2SBarry Smith .    -mat_block_size - size of the blocks to use
1506*61e22dc2SBarry Smith 
1507*61e22dc2SBarry Smith    Level: intermediate
1508*61e22dc2SBarry Smith 
1509*61e22dc2SBarry Smith    Notes:
1510*61e22dc2SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1511*61e22dc2SBarry Smith    storage.  That is, the stored row and column indices can begin at
1512*61e22dc2SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1513*61e22dc2SBarry Smith 
1514*61e22dc2SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1515*61e22dc2SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1516*61e22dc2SBarry Smith    allocation.  For additional details, see the users manual chapter on
1517*61e22dc2SBarry Smith    matrices.
1518*61e22dc2SBarry Smith 
1519*61e22dc2SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1520*61e22dc2SBarry Smith @*/
1521*61e22dc2SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1522*61e22dc2SBarry Smith {
1523*61e22dc2SBarry Smith   Mat         B;
1524*61e22dc2SBarry Smith   Mat_SeqBAIJ *b;
1525*61e22dc2SBarry Smith   int         i,len,ierr,mbs,nbs,bs2,size,newbs = bs;
1526*61e22dc2SBarry Smith   PetscTruth  flg;
1527*61e22dc2SBarry Smith 
1528*61e22dc2SBarry Smith   PetscFunctionBegin;
1529*61e22dc2SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1530*61e22dc2SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1531*61e22dc2SBarry Smith 
1532*61e22dc2SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr);
1533*61e22dc2SBarry Smith   if (nnz && newbs != bs) {
1534*61e22dc2SBarry Smith     SETERRQ(1,1,"Cannot change blocksize from command line if setting nnz");
1535*61e22dc2SBarry Smith   }
1536*61e22dc2SBarry Smith 
1537*61e22dc2SBarry Smith   mbs  = m/bs;
1538*61e22dc2SBarry Smith   nbs  = n/bs;
1539*61e22dc2SBarry Smith   bs2  = bs*bs;
1540*61e22dc2SBarry Smith 
1541*61e22dc2SBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1542*61e22dc2SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
1543*61e22dc2SBarry Smith   }
1544*61e22dc2SBarry Smith 
1545*61e22dc2SBarry Smith   if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz);
1546*61e22dc2SBarry Smith   if (nnz) {
1547*61e22dc2SBarry Smith     for (i=0; i<mbs; i++) {
1548*61e22dc2SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
1549*61e22dc2SBarry Smith     }
1550*61e22dc2SBarry Smith   }
1551*61e22dc2SBarry Smith 
1552*61e22dc2SBarry Smith   *A = 0;
1553*61e22dc2SBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView);
1554*61e22dc2SBarry Smith   PLogObjectCreate(B);
1555*61e22dc2SBarry Smith   B->data = (void*)(b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b);
1556*61e22dc2SBarry Smith   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1557*61e22dc2SBarry Smith   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1558*61e22dc2SBarry Smith   ierr    = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
1559*61e22dc2SBarry Smith   if (!flg) {
1560*61e22dc2SBarry Smith     switch (bs) {
1561*61e22dc2SBarry Smith     case 1:
1562*61e22dc2SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1563*61e22dc2SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1564*61e22dc2SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_1;
1565*61e22dc2SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1566*61e22dc2SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
1567*61e22dc2SBarry Smith       break;
1568*61e22dc2SBarry Smith     case 2:
1569*61e22dc2SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1570*61e22dc2SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1571*61e22dc2SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_2;
1572*61e22dc2SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1573*61e22dc2SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
1574*61e22dc2SBarry Smith       break;
1575*61e22dc2SBarry Smith     case 3:
1576*61e22dc2SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1577*61e22dc2SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1578*61e22dc2SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_3;
1579*61e22dc2SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1580*61e22dc2SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
1581*61e22dc2SBarry Smith       break;
1582*61e22dc2SBarry Smith     case 4:
1583*61e22dc2SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1584*61e22dc2SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1585*61e22dc2SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_4;
1586*61e22dc2SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1587*61e22dc2SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
1588*61e22dc2SBarry Smith       break;
1589*61e22dc2SBarry Smith     case 5:
1590*61e22dc2SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1591*61e22dc2SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1592*61e22dc2SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_5;
1593*61e22dc2SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1594*61e22dc2SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
1595*61e22dc2SBarry Smith       break;
1596*61e22dc2SBarry Smith     case 6:
1597*61e22dc2SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
1598*61e22dc2SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_6;
1599*61e22dc2SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_6;
1600*61e22dc2SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
1601*61e22dc2SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
1602*61e22dc2SBarry Smith       break;
1603*61e22dc2SBarry Smith     case 7:
1604*61e22dc2SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1605*61e22dc2SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1606*61e22dc2SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_7;
1607*61e22dc2SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
1608*61e22dc2SBarry Smith       break;
1609*61e22dc2SBarry Smith     }
1610*61e22dc2SBarry Smith   }
1611*61e22dc2SBarry Smith   B->factor           = 0;
1612*61e22dc2SBarry Smith   B->lupivotthreshold = 1.0;
1613*61e22dc2SBarry Smith   B->mapping          = 0;
1614*61e22dc2SBarry Smith   b->row              = 0;
1615*61e22dc2SBarry Smith   b->col              = 0;
1616*61e22dc2SBarry Smith   b->icol             = 0;
1617*61e22dc2SBarry Smith   b->reallocs         = 0;
1618*61e22dc2SBarry Smith   b->saved_values     = 0;
1619*61e22dc2SBarry Smith #if defined(PEYSC_USE_MAT_SINGLE)
1620*61e22dc2SBarry Smith   b->setvalueslen     = 0;
1621*61e22dc2SBarry Smith   b->setvaluescopy    = PETSC_NULL;
1622*61e22dc2SBarry Smith #endif
1623*61e22dc2SBarry Smith 
1624*61e22dc2SBarry Smith   b->m       = m; B->m = m; B->M = m;
1625*61e22dc2SBarry Smith   b->n       = n; B->n = n; B->N = n;
1626*61e22dc2SBarry Smith 
1627*61e22dc2SBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
1628*61e22dc2SBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1629*61e22dc2SBarry Smith 
1630*61e22dc2SBarry Smith   b->mbs     = mbs;
1631*61e22dc2SBarry Smith   b->nbs     = nbs;
1632*61e22dc2SBarry Smith   b->imax    = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax);
1633*61e22dc2SBarry Smith   if (!nnz) {
1634*61e22dc2SBarry Smith     if (nz == PETSC_DEFAULT) nz = 5;
1635*61e22dc2SBarry Smith     else if (nz <= 0)        nz = 1;
1636*61e22dc2SBarry Smith     for (i=0; i<mbs; i++) b->imax[i] = nz;
1637*61e22dc2SBarry Smith     nz = nz*mbs;
1638*61e22dc2SBarry Smith   } else {
1639*61e22dc2SBarry Smith     nz = 0;
1640*61e22dc2SBarry Smith     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1641*61e22dc2SBarry Smith   }
1642*61e22dc2SBarry Smith 
1643*61e22dc2SBarry Smith   /* allocate the matrix space */
1644*61e22dc2SBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int);
1645*61e22dc2SBarry Smith   b->a  = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a);
1646*61e22dc2SBarry Smith   ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1647*61e22dc2SBarry Smith   b->j  = (int*)(b->a + nz*bs2);
1648*61e22dc2SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
1649*61e22dc2SBarry Smith   b->i  = b->j + nz;
1650*61e22dc2SBarry Smith   b->singlemalloc = PETSC_TRUE;
1651*61e22dc2SBarry Smith 
1652*61e22dc2SBarry Smith   b->i[0] = 0;
1653*61e22dc2SBarry Smith   for (i=1; i<mbs+1; i++) {
1654*61e22dc2SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
1655*61e22dc2SBarry Smith   }
1656*61e22dc2SBarry Smith 
1657*61e22dc2SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1658*61e22dc2SBarry Smith   b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));
1659*61e22dc2SBarry Smith   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1660*61e22dc2SBarry Smith   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
1661*61e22dc2SBarry Smith 
1662*61e22dc2SBarry Smith   b->bs               = bs;
1663*61e22dc2SBarry Smith   b->bs2              = bs2;
1664*61e22dc2SBarry Smith   b->mbs              = mbs;
1665*61e22dc2SBarry Smith   b->nz               = 0;
1666*61e22dc2SBarry Smith   b->maxnz            = nz*bs2;
1667*61e22dc2SBarry Smith   b->sorted           = PETSC_FALSE;
1668*61e22dc2SBarry Smith   b->roworiented      = PETSC_TRUE;
1669*61e22dc2SBarry Smith   b->nonew            = 0;
1670*61e22dc2SBarry Smith   b->diag             = 0;
1671*61e22dc2SBarry Smith   b->solve_work       = 0;
1672*61e22dc2SBarry Smith   b->mult_work        = 0;
1673*61e22dc2SBarry Smith   b->spptr            = 0;
1674*61e22dc2SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1675*61e22dc2SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
1676*61e22dc2SBarry Smith 
1677*61e22dc2SBarry Smith   *A = B;
1678*61e22dc2SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
1679*61e22dc2SBarry Smith   if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); }
1680*61e22dc2SBarry Smith 
1681*61e22dc2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1682*61e22dc2SBarry Smith                                      "MatStoreValues_SeqBAIJ",
1683*61e22dc2SBarry Smith                                       MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1684*61e22dc2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1685*61e22dc2SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
1686*61e22dc2SBarry Smith                                       MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1687*61e22dc2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
1688*61e22dc2SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
1689*61e22dc2SBarry Smith                                       MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
1690*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1691*61e22dc2SBarry Smith }
1692*61e22dc2SBarry Smith 
1693*61e22dc2SBarry Smith #undef __FUNC__
1694*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDuplicate_SeqBAIJ"
1695*61e22dc2SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
1696*61e22dc2SBarry Smith {
1697*61e22dc2SBarry Smith   Mat         C;
1698*61e22dc2SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data;
1699*61e22dc2SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1700*61e22dc2SBarry Smith 
1701*61e22dc2SBarry Smith   PetscFunctionBegin;
1702*61e22dc2SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
1703*61e22dc2SBarry Smith 
1704*61e22dc2SBarry Smith   *B = 0;
1705*61e22dc2SBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView);
1706*61e22dc2SBarry Smith   PLogObjectCreate(C);
1707*61e22dc2SBarry Smith   C->data           = (void*)(c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c);
1708*61e22dc2SBarry Smith   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1709*61e22dc2SBarry Smith   C->factor         = A->factor;
1710*61e22dc2SBarry Smith   c->row            = 0;
1711*61e22dc2SBarry Smith   c->col            = 0;
1712*61e22dc2SBarry Smith   c->icol           = 0;
1713*61e22dc2SBarry Smith   c->saved_values   = 0;
1714*61e22dc2SBarry Smith   c->keepzeroedrows = a->keepzeroedrows;
1715*61e22dc2SBarry Smith   C->assembled      = PETSC_TRUE;
1716*61e22dc2SBarry Smith 
1717*61e22dc2SBarry Smith   c->m = C->m   = a->m;
1718*61e22dc2SBarry Smith   c->n = C->n   = a->n;
1719*61e22dc2SBarry Smith   C->M          = a->m;
1720*61e22dc2SBarry Smith   C->N          = a->n;
1721*61e22dc2SBarry Smith 
1722*61e22dc2SBarry Smith   c->bs         = a->bs;
1723*61e22dc2SBarry Smith   c->bs2        = a->bs2;
1724*61e22dc2SBarry Smith   c->mbs        = a->mbs;
1725*61e22dc2SBarry Smith   c->nbs        = a->nbs;
1726*61e22dc2SBarry Smith 
1727*61e22dc2SBarry Smith   c->imax       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax);
1728*61e22dc2SBarry Smith   c->ilen       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen);
1729*61e22dc2SBarry Smith   for (i=0; i<mbs; i++) {
1730*61e22dc2SBarry Smith     c->imax[i] = a->imax[i];
1731*61e22dc2SBarry Smith     c->ilen[i] = a->ilen[i];
1732*61e22dc2SBarry Smith   }
1733*61e22dc2SBarry Smith 
1734*61e22dc2SBarry Smith   /* allocate the matrix space */
1735*61e22dc2SBarry Smith   c->singlemalloc = PETSC_TRUE;
1736*61e22dc2SBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1737*61e22dc2SBarry Smith   c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a);
1738*61e22dc2SBarry Smith   c->j = (int*)(c->a + nz*bs2);
1739*61e22dc2SBarry Smith   c->i = c->j + nz;
1740*61e22dc2SBarry Smith   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1741*61e22dc2SBarry Smith   if (mbs > 0) {
1742*61e22dc2SBarry Smith     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
1743*61e22dc2SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1744*61e22dc2SBarry Smith       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
1745*61e22dc2SBarry Smith     } else {
1746*61e22dc2SBarry Smith       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
1747*61e22dc2SBarry Smith     }
1748*61e22dc2SBarry Smith   }
1749*61e22dc2SBarry Smith 
1750*61e22dc2SBarry Smith   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1751*61e22dc2SBarry Smith   c->sorted      = a->sorted;
1752*61e22dc2SBarry Smith   c->roworiented = a->roworiented;
1753*61e22dc2SBarry Smith   c->nonew       = a->nonew;
1754*61e22dc2SBarry Smith 
1755*61e22dc2SBarry Smith   if (a->diag) {
1756*61e22dc2SBarry Smith     c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag);
1757*61e22dc2SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1758*61e22dc2SBarry Smith     for (i=0; i<mbs; i++) {
1759*61e22dc2SBarry Smith       c->diag[i] = a->diag[i];
1760*61e22dc2SBarry Smith     }
1761*61e22dc2SBarry Smith   } else c->diag        = 0;
1762*61e22dc2SBarry Smith   c->nz                 = a->nz;
1763*61e22dc2SBarry Smith   c->maxnz              = a->maxnz;
1764*61e22dc2SBarry Smith   c->solve_work         = 0;
1765*61e22dc2SBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1766*61e22dc2SBarry Smith   c->mult_work          = 0;
1767*61e22dc2SBarry Smith   *B = C;
1768*61e22dc2SBarry Smith   ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
1769*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1770*61e22dc2SBarry Smith }
1771*61e22dc2SBarry Smith 
1772*61e22dc2SBarry Smith #undef __FUNC__
1773*61e22dc2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_SeqBAIJ"
1774*61e22dc2SBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
1775*61e22dc2SBarry Smith {
1776*61e22dc2SBarry Smith   Mat_SeqBAIJ  *a;
1777*61e22dc2SBarry Smith   Mat          B;
1778*61e22dc2SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1779*61e22dc2SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
1780*61e22dc2SBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1781*61e22dc2SBarry Smith   int          *masked,nmask,tmp,bs2,ishift;
1782*61e22dc2SBarry Smith   Scalar       *aa;
1783*61e22dc2SBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1784*61e22dc2SBarry Smith 
1785*61e22dc2SBarry Smith   PetscFunctionBegin;
1786*61e22dc2SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1787*61e22dc2SBarry Smith   bs2  = bs*bs;
1788*61e22dc2SBarry Smith 
1789*61e22dc2SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1790*61e22dc2SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
1791*61e22dc2SBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
1792*61e22dc2SBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1793*61e22dc2SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
1794*61e22dc2SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
1795*61e22dc2SBarry Smith 
1796*61e22dc2SBarry Smith   if (header[3] < 0) {
1797*61e22dc2SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1798*61e22dc2SBarry Smith   }
1799*61e22dc2SBarry Smith 
1800*61e22dc2SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
1801*61e22dc2SBarry Smith 
1802*61e22dc2SBarry Smith   /*
1803*61e22dc2SBarry Smith      This code adds extra rows to make sure the number of rows is
1804*61e22dc2SBarry Smith     divisible by the blocksize
1805*61e22dc2SBarry Smith   */
1806*61e22dc2SBarry Smith   mbs        = M/bs;
1807*61e22dc2SBarry Smith   extra_rows = bs - M + bs*(mbs);
1808*61e22dc2SBarry Smith   if (extra_rows == bs) extra_rows = 0;
1809*61e22dc2SBarry Smith   else                  mbs++;
1810*61e22dc2SBarry Smith   if (extra_rows) {
1811*61e22dc2SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
1812*61e22dc2SBarry Smith   }
1813*61e22dc2SBarry Smith 
1814*61e22dc2SBarry Smith   /* read in row lengths */
1815*61e22dc2SBarry Smith   rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
1816*61e22dc2SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1817*61e22dc2SBarry Smith   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
1818*61e22dc2SBarry Smith 
1819*61e22dc2SBarry Smith   /* read in column indices */
1820*61e22dc2SBarry Smith   jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj);
1821*61e22dc2SBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
1822*61e22dc2SBarry Smith   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
1823*61e22dc2SBarry Smith 
1824*61e22dc2SBarry Smith   /* loop over row lengths determining block row lengths */
1825*61e22dc2SBarry Smith   browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1826*61e22dc2SBarry Smith   ierr        = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
1827*61e22dc2SBarry Smith   mask        = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask);
1828*61e22dc2SBarry Smith   ierr        = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
1829*61e22dc2SBarry Smith   masked      = mask + mbs;
1830*61e22dc2SBarry Smith   rowcount    = 0; nzcount = 0;
1831*61e22dc2SBarry Smith   for (i=0; i<mbs; i++) {
1832*61e22dc2SBarry Smith     nmask = 0;
1833*61e22dc2SBarry Smith     for (j=0; j<bs; j++) {
1834*61e22dc2SBarry Smith       kmax = rowlengths[rowcount];
1835*61e22dc2SBarry Smith       for (k=0; k<kmax; k++) {
1836*61e22dc2SBarry Smith         tmp = jj[nzcount++]/bs;
1837*61e22dc2SBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1838*61e22dc2SBarry Smith       }
1839*61e22dc2SBarry Smith       rowcount++;
1840*61e22dc2SBarry Smith     }
1841*61e22dc2SBarry Smith     browlengths[i] += nmask;
1842*61e22dc2SBarry Smith     /* zero out the mask elements we set */
1843*61e22dc2SBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1844*61e22dc2SBarry Smith   }
1845*61e22dc2SBarry Smith 
1846*61e22dc2SBarry Smith   /* create our matrix */
1847*61e22dc2SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
1848*61e22dc2SBarry Smith   B = *A;
1849*61e22dc2SBarry Smith   a = (Mat_SeqBAIJ*)B->data;
1850*61e22dc2SBarry Smith 
1851*61e22dc2SBarry Smith   /* set matrix "i" values */
1852*61e22dc2SBarry Smith   a->i[0] = 0;
1853*61e22dc2SBarry Smith   for (i=1; i<= mbs; i++) {
1854*61e22dc2SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1855*61e22dc2SBarry Smith     a->ilen[i-1] = browlengths[i-1];
1856*61e22dc2SBarry Smith   }
1857*61e22dc2SBarry Smith   a->nz         = 0;
1858*61e22dc2SBarry Smith   for (i=0; i<mbs; i++) a->nz += browlengths[i];
1859*61e22dc2SBarry Smith 
1860*61e22dc2SBarry Smith   /* read in nonzero values */
1861*61e22dc2SBarry Smith   aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
1862*61e22dc2SBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
1863*61e22dc2SBarry Smith   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
1864*61e22dc2SBarry Smith 
1865*61e22dc2SBarry Smith   /* set "a" and "j" values into matrix */
1866*61e22dc2SBarry Smith   nzcount = 0; jcount = 0;
1867*61e22dc2SBarry Smith   for (i=0; i<mbs; i++) {
1868*61e22dc2SBarry Smith     nzcountb = nzcount;
1869*61e22dc2SBarry Smith     nmask    = 0;
1870*61e22dc2SBarry Smith     for (j=0; j<bs; j++) {
1871*61e22dc2SBarry Smith       kmax = rowlengths[i*bs+j];
1872*61e22dc2SBarry Smith       for (k=0; k<kmax; k++) {
1873*61e22dc2SBarry Smith         tmp = jj[nzcount++]/bs;
1874*61e22dc2SBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1875*61e22dc2SBarry Smith       }
1876*61e22dc2SBarry Smith       rowcount++;
1877*61e22dc2SBarry Smith     }
1878*61e22dc2SBarry Smith     /* sort the masked values */
1879*61e22dc2SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1880*61e22dc2SBarry Smith 
1881*61e22dc2SBarry Smith     /* set "j" values into matrix */
1882*61e22dc2SBarry Smith     maskcount = 1;
1883*61e22dc2SBarry Smith     for (j=0; j<nmask; j++) {
1884*61e22dc2SBarry Smith       a->j[jcount++]  = masked[j];
1885*61e22dc2SBarry Smith       mask[masked[j]] = maskcount++;
1886*61e22dc2SBarry Smith     }
1887*61e22dc2SBarry Smith     /* set "a" values into matrix */
1888*61e22dc2SBarry Smith     ishift = bs2*a->i[i];
1889*61e22dc2SBarry Smith     for (j=0; j<bs; j++) {
1890*61e22dc2SBarry Smith       kmax = rowlengths[i*bs+j];
1891*61e22dc2SBarry Smith       for (k=0; k<kmax; k++) {
1892*61e22dc2SBarry Smith         tmp       = jj[nzcountb]/bs ;
1893*61e22dc2SBarry Smith         block     = mask[tmp] - 1;
1894*61e22dc2SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1895*61e22dc2SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1896*61e22dc2SBarry Smith         a->a[idx] = (MatScalar)aa[nzcountb++];
1897*61e22dc2SBarry Smith       }
1898*61e22dc2SBarry Smith     }
1899*61e22dc2SBarry Smith     /* zero out the mask elements we set */
1900*61e22dc2SBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1901*61e22dc2SBarry Smith   }
1902*61e22dc2SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1903*61e22dc2SBarry Smith 
1904*61e22dc2SBarry Smith   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1905*61e22dc2SBarry Smith   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1906*61e22dc2SBarry Smith   ierr = PetscFree(aa);CHKERRQ(ierr);
1907*61e22dc2SBarry Smith   ierr = PetscFree(jj);CHKERRQ(ierr);
1908*61e22dc2SBarry Smith   ierr = PetscFree(mask);CHKERRQ(ierr);
1909*61e22dc2SBarry Smith 
1910*61e22dc2SBarry Smith   B->assembled = PETSC_TRUE;
1911*61e22dc2SBarry Smith 
1912*61e22dc2SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
1913*61e22dc2SBarry Smith   PetscFunctionReturn(0);
1914*61e22dc2SBarry Smith }
1915*61e22dc2SBarry Smith 
1916*61e22dc2SBarry Smith 
1917*61e22dc2SBarry Smith 
1918