xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 80cd9d9306eb90d00e6fd99dae46147660087c21)
1cb512458SBarry Smith #ifndef lint
2*80cd9d93SLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.100 1996/04/05 19:41:25 curfman Exp curfman $";
3cb512458SBarry Smith #endif
467e560aaSBarry Smith /*
567e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
667e560aaSBarry Smith */
7289bc588SBarry Smith 
817699dbbSLois Curfman McInnes #include "dense.h"
9b16a3bb1SBarry Smith #include "pinclude/plapack.h"
10b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
11289bc588SBarry Smith 
121987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
151987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
161987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
171987afe7SBarry Smith   PLogFlops(2*N-1);
181987afe7SBarry Smith   return 0;
191987afe7SBarry Smith }
201987afe7SBarry Smith 
21c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
22289bc588SBarry Smith {
23c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
24289bc588SBarry Smith   int          i,N = mat->m*mat->n,count = 0;
25289bc588SBarry Smith   Scalar       *v = mat->v;
26289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
27bcd2baecSBarry Smith   if (nz)      *nz      = count;
28bcd2baecSBarry Smith   if (nzalloc) *nzalloc = N;
29bcd2baecSBarry Smith   if (mem)     *mem     = (int)A->mem;
30fa9ec3f1SLois Curfman McInnes   return 0;
31289bc588SBarry Smith }
32289bc588SBarry Smith 
33*80cd9d93SLois Curfman McInnes static int MatScale_SeqDense(Scalar *alpha,Mat inA)
34*80cd9d93SLois Curfman McInnes {
35*80cd9d93SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
36*80cd9d93SLois Curfman McInnes   int          one = 1, nz;
37*80cd9d93SLois Curfman McInnes 
38*80cd9d93SLois Curfman McInnes   nz = a->m*a->n;
39*80cd9d93SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
40*80cd9d93SLois Curfman McInnes   PLogFlops(nz);
41*80cd9d93SLois Curfman McInnes   return 0;
42*80cd9d93SLois Curfman McInnes }
43*80cd9d93SLois Curfman McInnes 
44289bc588SBarry Smith /* ---------------------------------------------------------------*/
45289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
46289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
47c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
48289bc588SBarry Smith {
49c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
506abc6512SBarry Smith   int          info;
5155659b69SBarry Smith 
52289bc588SBarry Smith   if (!mat->pivots) {
530452661fSBarry Smith     mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
54c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
55289bc588SBarry Smith   }
56289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
57ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
58c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
5955659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
60289bc588SBarry Smith   return 0;
61289bc588SBarry Smith }
62c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
63289bc588SBarry Smith {
64a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
65289bc588SBarry Smith }
66c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
67289bc588SBarry Smith {
6849d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
69289bc588SBarry Smith }
70c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
71289bc588SBarry Smith {
72a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
73289bc588SBarry Smith }
74c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
75289bc588SBarry Smith {
7649d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
77289bc588SBarry Smith }
78c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
79289bc588SBarry Smith {
80c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
816abc6512SBarry Smith   int           info;
8255659b69SBarry Smith 
83752f5567SLois Curfman McInnes   if (mat->pivots) {
840452661fSBarry Smith     PetscFree(mat->pivots);
85c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
86752f5567SLois Curfman McInnes     mat->pivots = 0;
87752f5567SLois Curfman McInnes   }
88289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
89ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
90c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
9155659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
92289bc588SBarry Smith   return 0;
93289bc588SBarry Smith }
94289bc588SBarry Smith 
95c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
96289bc588SBarry Smith {
97c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
98a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
996abc6512SBarry Smith   Scalar       *x, *y;
10067e560aaSBarry Smith 
101a57ad546SLois Curfman McInnes   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
102416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
103c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
10448d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
105289bc588SBarry Smith   }
106c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
10748d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
108289bc588SBarry Smith   }
109ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
110ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
11155659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
112289bc588SBarry Smith   return 0;
113289bc588SBarry Smith }
114c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
115da3a660dSBarry Smith {
116c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1176abc6512SBarry Smith   int          one = 1, info;
1186abc6512SBarry Smith   Scalar       *x, *y;
11967e560aaSBarry Smith 
120da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
121416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
122752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
123da3a660dSBarry Smith   if (mat->pivots) {
12448d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
125da3a660dSBarry Smith   }
126da3a660dSBarry Smith   else {
12748d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
128da3a660dSBarry Smith   }
129ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
13055659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
131da3a660dSBarry Smith   return 0;
132da3a660dSBarry Smith }
133c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
134da3a660dSBarry Smith {
135c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1366abc6512SBarry Smith   int          one = 1, info,ierr;
1376abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
138da3a660dSBarry Smith   Vec          tmp = 0;
13967e560aaSBarry Smith 
140da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
141da3a660dSBarry Smith   if (yy == zz) {
14278b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
143c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
14478b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
145da3a660dSBarry Smith   }
146416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
147752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
148da3a660dSBarry Smith   if (mat->pivots) {
14948d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
150da3a660dSBarry Smith   }
151da3a660dSBarry Smith   else {
15248d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
153da3a660dSBarry Smith   }
154ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
155da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
156da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
15755659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
158da3a660dSBarry Smith   return 0;
159da3a660dSBarry Smith }
16067e560aaSBarry Smith 
161c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
162da3a660dSBarry Smith {
163c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1646abc6512SBarry Smith   int           one = 1, info,ierr;
1656abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
166da3a660dSBarry Smith   Vec           tmp;
16767e560aaSBarry Smith 
168da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
169da3a660dSBarry Smith   if (yy == zz) {
17078b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
171c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
17278b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
173da3a660dSBarry Smith   }
174416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
175752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
176da3a660dSBarry Smith   if (mat->pivots) {
17748d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
178da3a660dSBarry Smith   }
179da3a660dSBarry Smith   else {
18048d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
181da3a660dSBarry Smith   }
182ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
183da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
184da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
18555659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
186da3a660dSBarry Smith   return 0;
187da3a660dSBarry Smith }
188289bc588SBarry Smith /* ------------------------------------------------------------------*/
189c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
19020e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
191289bc588SBarry Smith {
192c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
193289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
1946abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
195289bc588SBarry Smith 
196289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
197289bc588SBarry Smith     /* this is a hack fix, should have another version without
198289bc588SBarry Smith        the second BLdot */
199bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
200289bc588SBarry Smith   }
201289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
202289bc588SBarry Smith   while (its--) {
203289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
204289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
205f1747703SBarry Smith #if defined(PETSC_COMPLEX)
206f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
207f1747703SBarry Smith            not happy about returning a double complex */
208f1747703SBarry Smith         int    _i;
209f1747703SBarry Smith         Scalar sum = b[i];
210f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
211f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
212f1747703SBarry Smith         }
213f1747703SBarry Smith         xt = sum;
214f1747703SBarry Smith #else
215289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
216f1747703SBarry Smith #endif
217d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
218289bc588SBarry Smith       }
219289bc588SBarry Smith     }
220289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
221289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
222f1747703SBarry Smith #if defined(PETSC_COMPLEX)
223f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
224f1747703SBarry Smith            not happy about returning a double complex */
225f1747703SBarry Smith         int    _i;
226f1747703SBarry Smith         Scalar sum = b[i];
227f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
228f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
229f1747703SBarry Smith         }
230f1747703SBarry Smith         xt = sum;
231f1747703SBarry Smith #else
232289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
233f1747703SBarry Smith #endif
234d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
235289bc588SBarry Smith       }
236289bc588SBarry Smith     }
237289bc588SBarry Smith   }
238289bc588SBarry Smith   return 0;
239289bc588SBarry Smith }
240289bc588SBarry Smith 
241289bc588SBarry Smith /* -----------------------------------------------------------------*/
242c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
243289bc588SBarry Smith {
244c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
245289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
246289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
247289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
24848d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
24955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
250289bc588SBarry Smith   return 0;
251289bc588SBarry Smith }
252c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
253289bc588SBarry Smith {
254c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
255289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
256289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
257289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
25848d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
25955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
260289bc588SBarry Smith   return 0;
261289bc588SBarry Smith }
262c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
263289bc588SBarry Smith {
264c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
265289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2666abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
267289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
268416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
26948d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
27055659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
271289bc588SBarry Smith   return 0;
272289bc588SBarry Smith }
273c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
274289bc588SBarry Smith {
275c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
276289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
277289bc588SBarry Smith   int          _One=1;
2786abc6512SBarry Smith   Scalar       _DOne=1.0;
279289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
280289bc588SBarry Smith   VecGetArray(zz,&z);
281717eeb74SLois Curfman McInnes   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
28248d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
28355659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
284289bc588SBarry Smith   return 0;
285289bc588SBarry Smith }
286289bc588SBarry Smith 
287289bc588SBarry Smith /* -----------------------------------------------------------------*/
288c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
289289bc588SBarry Smith {
290c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
291289bc588SBarry Smith   Scalar       *v;
292289bc588SBarry Smith   int          i;
29367e560aaSBarry Smith 
294289bc588SBarry Smith   *ncols = mat->n;
295289bc588SBarry Smith   if (cols) {
2960452661fSBarry Smith     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
29773c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
298289bc588SBarry Smith   }
299289bc588SBarry Smith   if (vals) {
3000452661fSBarry Smith     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
301289bc588SBarry Smith     v = mat->v + row;
30273c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
303289bc588SBarry Smith   }
304289bc588SBarry Smith   return 0;
305289bc588SBarry Smith }
306c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
307289bc588SBarry Smith {
3080452661fSBarry Smith   if (cols) { PetscFree(*cols); }
3090452661fSBarry Smith   if (vals) { PetscFree(*vals); }
310289bc588SBarry Smith   return 0;
311289bc588SBarry Smith }
312289bc588SBarry Smith /* ----------------------------------------------------------------*/
313ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
314e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
315289bc588SBarry Smith {
316c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
317289bc588SBarry Smith   int          i,j;
318d6dfbf8fSBarry Smith 
319289bc588SBarry Smith   if (!mat->roworiented) {
320dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
321289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
322289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
323289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
324289bc588SBarry Smith         }
325289bc588SBarry Smith       }
326289bc588SBarry Smith     }
327289bc588SBarry Smith     else {
328289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
329289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
330289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
331289bc588SBarry Smith         }
332289bc588SBarry Smith       }
333289bc588SBarry Smith     }
334e8d4e0b9SBarry Smith   }
335e8d4e0b9SBarry Smith   else {
336dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
337e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
338e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
339e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
340e8d4e0b9SBarry Smith         }
341e8d4e0b9SBarry Smith       }
342e8d4e0b9SBarry Smith     }
343289bc588SBarry Smith     else {
344289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
345289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
346289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
347289bc588SBarry Smith         }
348289bc588SBarry Smith       }
349289bc588SBarry Smith     }
350e8d4e0b9SBarry Smith   }
351289bc588SBarry Smith   return 0;
352289bc588SBarry Smith }
353e8d4e0b9SBarry Smith 
354ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
355ae80bb75SLois Curfman McInnes {
356ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
357ae80bb75SLois Curfman McInnes   int          i, j;
358ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
359ae80bb75SLois Curfman McInnes 
360ae80bb75SLois Curfman McInnes   /* row-oriented output */
361ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
362ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
363ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
364ae80bb75SLois Curfman McInnes     }
365ae80bb75SLois Curfman McInnes   }
366ae80bb75SLois Curfman McInnes   return 0;
367ae80bb75SLois Curfman McInnes }
368ae80bb75SLois Curfman McInnes 
369289bc588SBarry Smith /* -----------------------------------------------------------------*/
3703d1612f7SBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
371289bc588SBarry Smith {
372c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
373289bc588SBarry Smith   int          ierr;
374289bc588SBarry Smith   Mat          newi;
37548d91487SBarry Smith 
376b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
377ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
37855659b69SBarry Smith   if (cpvalues == COPY_VALUES) {
379416022c9SBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
38055659b69SBarry Smith   }
381227d817aSBarry Smith   newi->assembled = PETSC_TRUE;
382289bc588SBarry Smith   *newmat = newi;
383289bc588SBarry Smith   return 0;
384289bc588SBarry Smith }
385289bc588SBarry Smith 
38677c4ece6SBarry Smith #include "sys.h"
38788e32edaSLois Curfman McInnes 
38819bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
38988e32edaSLois Curfman McInnes {
39088e32edaSLois Curfman McInnes   Mat_SeqDense *a;
39188e32edaSLois Curfman McInnes   Mat          B;
39288e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
39388e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
39488e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
39519bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
39688e32edaSLois Curfman McInnes 
39788e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
39888e32edaSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
39990ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
40077c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
40188e32edaSLois Curfman McInnes   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
40288e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
40388e32edaSLois Curfman McInnes 
40490ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
40590ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
40690ace30eSBarry Smith     B = *A;
40790ace30eSBarry Smith     a = (Mat_SeqDense *) B->data;
40890ace30eSBarry Smith 
40990ace30eSBarry Smith     /* read in nonzero values */
41077c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
41190ace30eSBarry Smith 
41290ace30eSBarry Smith     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
41390ace30eSBarry Smith     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
41490ace30eSBarry Smith     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
41590ace30eSBarry Smith   } else {
41688e32edaSLois Curfman McInnes     /* read row lengths */
4170452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
41877c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
41988e32edaSLois Curfman McInnes 
42088e32edaSLois Curfman McInnes     /* create our matrix */
421b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
42288e32edaSLois Curfman McInnes     B = *A;
42388e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
42488e32edaSLois Curfman McInnes     v = a->v;
42588e32edaSLois Curfman McInnes 
42688e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
4270452661fSBarry Smith     cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
42877c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
4290452661fSBarry Smith     vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
43077c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
43188e32edaSLois Curfman McInnes 
43288e32edaSLois Curfman McInnes     /* insert into matrix */
43388e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
43488e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
43588e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
43688e32edaSLois Curfman McInnes     }
4370452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
43888e32edaSLois Curfman McInnes 
43988e32edaSLois Curfman McInnes     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
44088e32edaSLois Curfman McInnes     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
44190ace30eSBarry Smith   }
44288e32edaSLois Curfman McInnes   return 0;
44388e32edaSLois Curfman McInnes }
44488e32edaSLois Curfman McInnes 
445932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
44677c4ece6SBarry Smith #include "sys.h"
447932b0c3eSLois Curfman McInnes 
448932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
449289bc588SBarry Smith {
450932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
451932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
452d636dbe3SBarry Smith   FILE         *fd;
453932b0c3eSLois Curfman McInnes   char         *outputname;
454932b0c3eSLois Curfman McInnes   Scalar       *v;
455932b0c3eSLois Curfman McInnes 
45690ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
457932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
45890ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
45990ace30eSBarry Smith   if (format == ASCII_FORMAT_INFO) {
460932b0c3eSLois Curfman McInnes     ;  /* do nothing for now */
461f72585f2SLois Curfman McInnes   }
462*80cd9d93SLois Curfman McInnes   /* We retain dense format as default; use ASCII_FORMAT_IMPL to get
463*80cd9d93SLois Curfman McInnes      sparse format output ... maybe should switch this? */
464*80cd9d93SLois Curfman McInnes   else if (format == ASCII_FORMAT_IMPL) {
465*80cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
466*80cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
467*80cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
468*80cd9d93SLois Curfman McInnes #if defined(PETSC_COMPLEX)
469*80cd9d93SLois Curfman McInnes         if (real(*v) != 0.0 && imag(*v) != 0.0)
470*80cd9d93SLois Curfman McInnes           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
471*80cd9d93SLois Curfman McInnes         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
472*80cd9d93SLois Curfman McInnes         v += a->m;
473*80cd9d93SLois Curfman McInnes #else
474*80cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
475*80cd9d93SLois Curfman McInnes         v += a->m;
476*80cd9d93SLois Curfman McInnes #endif
477*80cd9d93SLois Curfman McInnes       }
478*80cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
479*80cd9d93SLois Curfman McInnes     }
480*80cd9d93SLois Curfman McInnes   }
481f72585f2SLois Curfman McInnes   else {
48247989497SBarry Smith #if defined(PETSC_COMPLEX)
48347989497SBarry Smith     int allreal = 1;
48447989497SBarry Smith     /* determine if matrix has all real values */
48547989497SBarry Smith     v = a->v;
48647989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
48747989497SBarry Smith       if (imag(v[i])) { allreal = 0; break ;}
48847989497SBarry Smith     }
48947989497SBarry Smith #endif
490932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
491932b0c3eSLois Curfman McInnes       v = a->v + i;
492932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
493289bc588SBarry Smith #if defined(PETSC_COMPLEX)
49447989497SBarry Smith         if (allreal) {
49547989497SBarry Smith           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
49647989497SBarry Smith         } else {
497932b0c3eSLois Curfman McInnes           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
49847989497SBarry Smith         }
499289bc588SBarry Smith #else
500932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
501289bc588SBarry Smith #endif
502289bc588SBarry Smith       }
5038e37dc09SBarry Smith       fprintf(fd,"\n");
504289bc588SBarry Smith     }
505da3a660dSBarry Smith   }
506932b0c3eSLois Curfman McInnes   fflush(fd);
507289bc588SBarry Smith   return 0;
508289bc588SBarry Smith }
509289bc588SBarry Smith 
510932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
511932b0c3eSLois Curfman McInnes {
512932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
513932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
51490ace30eSBarry Smith   int          format;
51590ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
516932b0c3eSLois Curfman McInnes 
51790ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
51890ace30eSBarry Smith 
51990ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
52090ace30eSBarry Smith   if (format == BINARY_FORMAT_NATIVE) {
52190ace30eSBarry Smith     /* store the matrix as a dense matrix */
52290ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
52390ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
52490ace30eSBarry Smith     col_lens[1] = m;
52590ace30eSBarry Smith     col_lens[2] = n;
52690ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
52777c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
52890ace30eSBarry Smith     PetscFree(col_lens);
52990ace30eSBarry Smith 
53090ace30eSBarry Smith     /* write out matrix, by rows */
53190ace30eSBarry Smith     vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals);
53290ace30eSBarry Smith     v    = a->v;
53390ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
53490ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
53590ace30eSBarry Smith         vals[i + j*m] = *v++;
53690ace30eSBarry Smith       }
53790ace30eSBarry Smith     }
53877c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
53990ace30eSBarry Smith     PetscFree(vals);
54090ace30eSBarry Smith   } else {
5410452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
542932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
543932b0c3eSLois Curfman McInnes     col_lens[1] = m;
544932b0c3eSLois Curfman McInnes     col_lens[2] = n;
545932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
546932b0c3eSLois Curfman McInnes 
547932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
548932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
54977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
550932b0c3eSLois Curfman McInnes 
551932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
552932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
553932b0c3eSLois Curfman McInnes     ict = 0;
554932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
555932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
556932b0c3eSLois Curfman McInnes     }
55777c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
5580452661fSBarry Smith     PetscFree(col_lens);
559932b0c3eSLois Curfman McInnes 
560932b0c3eSLois Curfman McInnes     /* store nonzero values */
5610452661fSBarry Smith     anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
562932b0c3eSLois Curfman McInnes     ict = 0;
563932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
564932b0c3eSLois Curfman McInnes       v = a->v + i;
565932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
566932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
567932b0c3eSLois Curfman McInnes       }
568932b0c3eSLois Curfman McInnes     }
56977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
5700452661fSBarry Smith     PetscFree(anonz);
57190ace30eSBarry Smith   }
572932b0c3eSLois Curfman McInnes   return 0;
573932b0c3eSLois Curfman McInnes }
574932b0c3eSLois Curfman McInnes 
575932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
576932b0c3eSLois Curfman McInnes {
577932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
578932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
579bcd2baecSBarry Smith   ViewerType   vtype;
580bcd2baecSBarry Smith   int          ierr;
581932b0c3eSLois Curfman McInnes 
582932b0c3eSLois Curfman McInnes   if (!viewer) {
58319bcc07fSBarry Smith     viewer = STDOUT_VIEWER_SELF;
584932b0c3eSLois Curfman McInnes   }
585bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
586bcd2baecSBarry Smith 
587bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
588932b0c3eSLois Curfman McInnes     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
589932b0c3eSLois Curfman McInnes   }
59019bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
591932b0c3eSLois Curfman McInnes     return MatView_SeqDense_ASCII(A,viewer);
592932b0c3eSLois Curfman McInnes   }
593bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
594932b0c3eSLois Curfman McInnes     return MatView_SeqDense_Binary(A,viewer);
595932b0c3eSLois Curfman McInnes   }
596932b0c3eSLois Curfman McInnes   return 0;
597932b0c3eSLois Curfman McInnes }
598289bc588SBarry Smith 
599ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
600289bc588SBarry Smith {
601289bc588SBarry Smith   Mat          mat = (Mat) obj;
602ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
603a5a9c739SBarry Smith #if defined(PETSC_LOG)
604a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
605a5a9c739SBarry Smith #endif
6060452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
6073439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
6080452661fSBarry Smith   PetscFree(l);
609a5a9c739SBarry Smith   PLogObjectDestroy(mat);
6100452661fSBarry Smith   PetscHeaderDestroy(mat);
611289bc588SBarry Smith   return 0;
612289bc588SBarry Smith }
613289bc588SBarry Smith 
614c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
615289bc588SBarry Smith {
616c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
617d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
618d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
61948b35521SBarry Smith 
620d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
6213638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
622d9f96c7cSLois Curfman McInnes     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
623d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
624289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
625d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
626d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
627d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
628289bc588SBarry Smith       }
629289bc588SBarry Smith     }
63048b35521SBarry Smith   }
631d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
632d3e5ee88SLois Curfman McInnes     int          ierr;
633d3e5ee88SLois Curfman McInnes     Mat          tmat;
634ec8511deSBarry Smith     Mat_SeqDense *tmatd;
635d3e5ee88SLois Curfman McInnes     Scalar       *v2;
636b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
637ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
6380de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
639d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
6400de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
641d3e5ee88SLois Curfman McInnes     }
642d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
643d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
644d3e5ee88SLois Curfman McInnes     *matout = tmat;
64548b35521SBarry Smith   }
646289bc588SBarry Smith   return 0;
647289bc588SBarry Smith }
648289bc588SBarry Smith 
64977c4ece6SBarry Smith static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
650289bc588SBarry Smith {
651c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
652c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
653289bc588SBarry Smith   int          i;
654289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
6559ea5d5aeSSatish Balay 
6569ea5d5aeSSatish Balay   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
65777c4ece6SBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
65877c4ece6SBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
659289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
66077c4ece6SBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
661289bc588SBarry Smith     v1++; v2++;
662289bc588SBarry Smith   }
66377c4ece6SBarry Smith   *flg = PETSC_TRUE;
6649ea5d5aeSSatish Balay   return 0;
665289bc588SBarry Smith }
666289bc588SBarry Smith 
667c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
668289bc588SBarry Smith {
669c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
6706abc6512SBarry Smith   int          i, n;
6716abc6512SBarry Smith   Scalar       *x;
672289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
673ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
674289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
675289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
676289bc588SBarry Smith   }
677289bc588SBarry Smith   return 0;
678289bc588SBarry Smith }
679289bc588SBarry Smith 
680052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
681289bc588SBarry Smith {
682c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
683da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
684da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
68555659b69SBarry Smith 
68628988994SBarry Smith   if (ll) {
687da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
688052efed2SBarry Smith     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
68955659b69SBarry Smith     PLogFlops(n*m);
690da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
691da3a660dSBarry Smith       x = l[i];
692da3a660dSBarry Smith       v = mat->v + i;
693da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
694da3a660dSBarry Smith     }
695da3a660dSBarry Smith   }
69628988994SBarry Smith   if (rr) {
697da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
698052efed2SBarry Smith     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
69955659b69SBarry Smith     PLogFlops(n*m);
700da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
701da3a660dSBarry Smith       x = r[i];
702da3a660dSBarry Smith       v = mat->v + i*m;
703da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
704da3a660dSBarry Smith     }
705da3a660dSBarry Smith   }
706289bc588SBarry Smith   return 0;
707289bc588SBarry Smith }
708289bc588SBarry Smith 
709cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
710289bc588SBarry Smith {
711c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
712289bc588SBarry Smith   Scalar       *v = mat->v;
713289bc588SBarry Smith   double       sum = 0.0;
714289bc588SBarry Smith   int          i, j;
71555659b69SBarry Smith 
716289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
717289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
718289bc588SBarry Smith #if defined(PETSC_COMPLEX)
719289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
720289bc588SBarry Smith #else
721289bc588SBarry Smith       sum += (*v)*(*v); v++;
722289bc588SBarry Smith #endif
723289bc588SBarry Smith     }
724289bc588SBarry Smith     *norm = sqrt(sum);
72555659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
726289bc588SBarry Smith   }
727289bc588SBarry Smith   else if (type == NORM_1) {
728289bc588SBarry Smith     *norm = 0.0;
729289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
730289bc588SBarry Smith       sum = 0.0;
731289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
73233a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
733289bc588SBarry Smith       }
734289bc588SBarry Smith       if (sum > *norm) *norm = sum;
735289bc588SBarry Smith     }
73655659b69SBarry Smith     PLogFlops(mat->n*mat->m);
737289bc588SBarry Smith   }
738289bc588SBarry Smith   else if (type == NORM_INFINITY) {
739289bc588SBarry Smith     *norm = 0.0;
740289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
741289bc588SBarry Smith       v = mat->v + j;
742289bc588SBarry Smith       sum = 0.0;
743289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
744cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
745289bc588SBarry Smith       }
746289bc588SBarry Smith       if (sum > *norm) *norm = sum;
747289bc588SBarry Smith     }
74855659b69SBarry Smith     PLogFlops(mat->n*mat->m);
749289bc588SBarry Smith   }
750289bc588SBarry Smith   else {
75148d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
752289bc588SBarry Smith   }
753289bc588SBarry Smith   return 0;
754289bc588SBarry Smith }
755289bc588SBarry Smith 
756c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
757289bc588SBarry Smith {
758c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
75967e560aaSBarry Smith 
760289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
761289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
762c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
76388e32edaSLois Curfman McInnes            op == COLUMNS_SORTED ||
764c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
765c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
766c0bbcb79SLois Curfman McInnes            op == NO_NEW_NONZERO_LOCATIONS ||
767c0bbcb79SLois Curfman McInnes            op == YES_NEW_NONZERO_LOCATIONS ||
768c0bbcb79SLois Curfman McInnes            op == NO_NEW_DIAGONALS ||
769c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
77094a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
771c0bbcb79SLois Curfman McInnes   else
7724d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
773289bc588SBarry Smith   return 0;
774289bc588SBarry Smith }
775289bc588SBarry Smith 
776ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
7776f0a148fSBarry Smith {
778ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
779cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
7806f0a148fSBarry Smith   return 0;
7816f0a148fSBarry Smith }
7826f0a148fSBarry Smith 
783ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
7846f0a148fSBarry Smith {
785ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
7866abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
7876f0a148fSBarry Smith   Scalar       *slot;
78855659b69SBarry Smith 
78977c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
79078b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
7916f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
7926f0a148fSBarry Smith     slot = l->v + rows[i];
7936f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
7946f0a148fSBarry Smith   }
7956f0a148fSBarry Smith   if (diag) {
7966f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
7976f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
798260925b8SBarry Smith       *slot = *diag;
7996f0a148fSBarry Smith     }
8006f0a148fSBarry Smith   }
801260925b8SBarry Smith   ISRestoreIndices(is,&rows);
8026f0a148fSBarry Smith   return 0;
8036f0a148fSBarry Smith }
804557bce09SLois Curfman McInnes 
805c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
806557bce09SLois Curfman McInnes {
807c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
808557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
809557bce09SLois Curfman McInnes   return 0;
810557bce09SLois Curfman McInnes }
811557bce09SLois Curfman McInnes 
812c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
813d3e5ee88SLois Curfman McInnes {
814c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
815d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
816d3e5ee88SLois Curfman McInnes   return 0;
817d3e5ee88SLois Curfman McInnes }
818d3e5ee88SLois Curfman McInnes 
819c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
82064e87e97SBarry Smith {
821c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
82264e87e97SBarry Smith   *array = mat->v;
82364e87e97SBarry Smith   return 0;
82464e87e97SBarry Smith }
8250754003eSLois Curfman McInnes 
826ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
827ff14e315SSatish Balay {
828ff14e315SSatish Balay   return 0;
829ff14e315SSatish Balay }
8300754003eSLois Curfman McInnes 
831c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
8320754003eSLois Curfman McInnes {
833ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
8340754003eSLois Curfman McInnes }
8350754003eSLois Curfman McInnes 
836cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
837cddf8d76SBarry Smith                                     Mat *submat)
8380754003eSLois Curfman McInnes {
839c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
8400754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
841160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
8420754003eSLois Curfman McInnes   Scalar       *vwork, *val;
8430754003eSLois Curfman McInnes   Mat          newmat;
8440754003eSLois Curfman McInnes 
84578b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
84678b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
84778b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
84878b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
8490754003eSLois Curfman McInnes 
8500452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
8510452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
8520452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
853cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
8540754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
8550754003eSLois Curfman McInnes 
8560754003eSLois Curfman McInnes   /* Create and fill new matrix */
857b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
8580754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
8590754003eSLois Curfman McInnes     nznew = 0;
8600754003eSLois Curfman McInnes     val   = mat->v + irow[i];
8610754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
8620754003eSLois Curfman McInnes       if (smap[j]) {
8630754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
8640754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
8650754003eSLois Curfman McInnes       }
8660754003eSLois Curfman McInnes     }
867dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
86878b31e54SBarry Smith            CHKERRQ(ierr);
8690754003eSLois Curfman McInnes   }
87078b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
87178b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
8720754003eSLois Curfman McInnes 
8730754003eSLois Curfman McInnes   /* Free work space */
8740452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
87578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
87678b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
8770754003eSLois Curfman McInnes   *submat = newmat;
8780754003eSLois Curfman McInnes   return 0;
8790754003eSLois Curfman McInnes }
8800754003eSLois Curfman McInnes 
8814b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B)
8824b0e389bSBarry Smith {
8834b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
8844b0e389bSBarry Smith   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
8854b0e389bSBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
8864b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
8874b0e389bSBarry Smith   return 0;
8884b0e389bSBarry Smith }
8894b0e389bSBarry Smith 
890289bc588SBarry Smith /* -------------------------------------------------------------------*/
891ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
892ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
893ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
894ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
895ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
896ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
897ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
898ec8511deSBarry Smith        MatRelax_SeqDense,
899ec8511deSBarry Smith        MatTranspose_SeqDense,
900ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
901052efed2SBarry Smith        MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense,
902289bc588SBarry Smith        0,0,
903ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
904ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
905ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
906ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
907ff14e315SSatish Balay        0,0, MatGetArray_SeqDense, MatRestoreArray_SeqDense,0,
908ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
9093d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
910ae80bb75SLois Curfman McInnes        MatAXPY_SeqDense,0,0,
9114b0e389bSBarry Smith        MatGetValues_SeqDense,
912*80cd9d93SLois Curfman McInnes        MatCopy_SeqDense,0,MatScale_SeqDense};
9130754003eSLois Curfman McInnes 
91490ace30eSBarry Smith 
9154b828684SBarry Smith /*@C
916fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
917d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
918d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
919289bc588SBarry Smith 
92020563c6bSBarry Smith    Input Parameters:
9210c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
9220c775827SLois Curfman McInnes .  m - number of rows
92318f449edSLois Curfman McInnes .  n - number of columns
924b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
925dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
92620563c6bSBarry Smith 
92720563c6bSBarry Smith    Output Parameter:
9280c775827SLois Curfman McInnes .  newmat - the matrix
92920563c6bSBarry Smith 
93018f449edSLois Curfman McInnes   Notes:
93118f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
93218f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
933b4fd4287SBarry Smith   set data=PETSC_NULL.
93418f449edSLois Curfman McInnes 
935dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
936d65003e9SLois Curfman McInnes 
937dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
93820563c6bSBarry Smith @*/
93918f449edSLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat)
940289bc588SBarry Smith {
941289bc588SBarry Smith   Mat          mat;
942ec8511deSBarry Smith   Mat_SeqDense *l;
94325cdf11fSBarry Smith   int          ierr,flg;
94455659b69SBarry Smith 
945289bc588SBarry Smith   *newmat        = 0;
94618f449edSLois Curfman McInnes 
9470452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
948a5a9c739SBarry Smith   PLogObjectCreate(mat);
9493439631bSBarry Smith   l               = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l);
950416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
951ec8511deSBarry Smith   mat->destroy    = MatDestroy_SeqDense;
952ec8511deSBarry Smith   mat->view       = MatView_SeqDense;
953289bc588SBarry Smith   mat->factor     = 0;
9543439631bSBarry Smith   PLogObjectMemory(mat,sizeof(struct _Mat));
95518f449edSLois Curfman McInnes   mat->data       = (void *) l;
956d6dfbf8fSBarry Smith 
957289bc588SBarry Smith   l->m            = m;
958289bc588SBarry Smith   l->n            = n;
959289bc588SBarry Smith   l->pivots       = 0;
960289bc588SBarry Smith   l->roworiented  = 1;
961b4fd4287SBarry Smith   if (data == PETSC_NULL) {
962ba7af9b1SBarry Smith     l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v);
963cddf8d76SBarry Smith     PetscMemzero(l->v,m*n*sizeof(Scalar));
9642dd2b441SLois Curfman McInnes     l->user_alloc = 0;
9653439631bSBarry Smith     PLogObjectMemory(mat,n*m*sizeof(Scalar));
96618f449edSLois Curfman McInnes   }
9672dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
9682dd2b441SLois Curfman McInnes     l->v = data;
9692dd2b441SLois Curfman McInnes     l->user_alloc = 1;
9702dd2b441SLois Curfman McInnes   }
97118f449edSLois Curfman McInnes 
97225cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
97325cdf11fSBarry Smith   if (flg) {
974682d7d0cSBarry Smith     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
975682d7d0cSBarry Smith   }
976289bc588SBarry Smith   *newmat = mat;
977289bc588SBarry Smith   return 0;
978289bc588SBarry Smith }
979289bc588SBarry Smith 
980c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
981289bc588SBarry Smith {
982c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
983b4fd4287SBarry Smith   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
984289bc588SBarry Smith }
985