xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 4e220ebcc634eb0c8ccd78bed049803c9194c4d4)
1cb512458SBarry Smith #ifndef lint
2*4e220ebcSLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.108 1996/08/08 14:42:39 bsmith Exp curfman $";
3cb512458SBarry Smith #endif
467e560aaSBarry Smith /*
567e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
667e560aaSBarry Smith */
7289bc588SBarry Smith 
870f55243SBarry Smith #include "src/mat/impls/dense/seq/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 
21*4e220ebcSLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
22289bc588SBarry Smith {
23*4e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
24*4e220ebcSLois Curfman McInnes   int          i,N = a->m*a->n,count = 0;
25*4e220ebcSLois Curfman McInnes   Scalar       *v = a->v;
26289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
27*4e220ebcSLois Curfman McInnes 
28*4e220ebcSLois Curfman McInnes   info->rows_global       = (double)a->m;
29*4e220ebcSLois Curfman McInnes   info->columns_global    = (double)a->n;
30*4e220ebcSLois Curfman McInnes   info->rows_local        = (double)a->m;
31*4e220ebcSLois Curfman McInnes   info->columns_local     = (double)a->n;
32*4e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
33*4e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
34*4e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
35*4e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
36*4e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
37*4e220ebcSLois Curfman McInnes   info->mallocs           = 0;
38*4e220ebcSLois Curfman McInnes   info->memory            = A->mem;
39*4e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
40*4e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
41*4e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
42*4e220ebcSLois Curfman McInnes 
43fa9ec3f1SLois Curfman McInnes   return 0;
44289bc588SBarry Smith }
45289bc588SBarry Smith 
4680cd9d93SLois Curfman McInnes static int MatScale_SeqDense(Scalar *alpha,Mat inA)
4780cd9d93SLois Curfman McInnes {
4880cd9d93SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
4980cd9d93SLois Curfman McInnes   int          one = 1, nz;
5080cd9d93SLois Curfman McInnes 
5180cd9d93SLois Curfman McInnes   nz = a->m*a->n;
5280cd9d93SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
5380cd9d93SLois Curfman McInnes   PLogFlops(nz);
5480cd9d93SLois Curfman McInnes   return 0;
5580cd9d93SLois Curfman McInnes }
5680cd9d93SLois Curfman McInnes 
57289bc588SBarry Smith /* ---------------------------------------------------------------*/
58289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
59289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
60c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
61289bc588SBarry Smith {
62c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
636abc6512SBarry Smith   int          info;
6455659b69SBarry Smith 
65289bc588SBarry Smith   if (!mat->pivots) {
660452661fSBarry Smith     mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
67c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
68289bc588SBarry Smith   }
69289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
703d60af36SSatish Balay   if (info<0) SETERRQ(1,"MatLUFactor_SeqDense:Bad argument to LU factorization");
713d60af36SSatish Balay   if (info>0) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
72c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
7355659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
74289bc588SBarry Smith   return 0;
75289bc588SBarry Smith }
76c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
77289bc588SBarry Smith {
78a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
79289bc588SBarry Smith }
80c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
81289bc588SBarry Smith {
8249d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
83289bc588SBarry Smith }
84c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
85289bc588SBarry Smith {
86a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
87289bc588SBarry Smith }
88c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
89289bc588SBarry Smith {
9049d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
91289bc588SBarry Smith }
92c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
93289bc588SBarry Smith {
94c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
956abc6512SBarry Smith   int           info;
9655659b69SBarry Smith 
97752f5567SLois Curfman McInnes   if (mat->pivots) {
980452661fSBarry Smith     PetscFree(mat->pivots);
99c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
100752f5567SLois Curfman McInnes     mat->pivots = 0;
101752f5567SLois Curfman McInnes   }
102289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
103ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
104c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
10555659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
106289bc588SBarry Smith   return 0;
107289bc588SBarry Smith }
108289bc588SBarry Smith 
109c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
110289bc588SBarry Smith {
111c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
112a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
1136abc6512SBarry Smith   Scalar       *x, *y;
11467e560aaSBarry Smith 
115a57ad546SLois Curfman McInnes   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
116416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
117c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
11848d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
119289bc588SBarry Smith   }
120c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
12148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
122289bc588SBarry Smith   }
123ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
124ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
12555659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
126289bc588SBarry Smith   return 0;
127289bc588SBarry Smith }
128c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
129da3a660dSBarry Smith {
130c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1316abc6512SBarry Smith   int          one = 1, info;
1326abc6512SBarry Smith   Scalar       *x, *y;
13367e560aaSBarry Smith 
134da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
135416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
136752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
137da3a660dSBarry Smith   if (mat->pivots) {
13848d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
139da3a660dSBarry Smith   }
140da3a660dSBarry Smith   else {
14148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
142da3a660dSBarry Smith   }
143ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
14455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
145da3a660dSBarry Smith   return 0;
146da3a660dSBarry Smith }
147c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
148da3a660dSBarry Smith {
149c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1506abc6512SBarry Smith   int          one = 1, info,ierr;
1516abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
152da3a660dSBarry Smith   Vec          tmp = 0;
15367e560aaSBarry Smith 
154da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
155da3a660dSBarry Smith   if (yy == zz) {
15678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
157c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
15878b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
159da3a660dSBarry Smith   }
160416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
161752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
162da3a660dSBarry Smith   if (mat->pivots) {
16348d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
164da3a660dSBarry Smith   }
165da3a660dSBarry Smith   else {
16648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
167da3a660dSBarry Smith   }
168ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
169da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
170da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
17155659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
172da3a660dSBarry Smith   return 0;
173da3a660dSBarry Smith }
17467e560aaSBarry Smith 
175c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
176da3a660dSBarry Smith {
177c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1786abc6512SBarry Smith   int           one = 1, info,ierr;
1796abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
180da3a660dSBarry Smith   Vec           tmp;
18167e560aaSBarry Smith 
182da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
183da3a660dSBarry Smith   if (yy == zz) {
18478b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
185c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
18678b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
187da3a660dSBarry Smith   }
188416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
189752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
190da3a660dSBarry Smith   if (mat->pivots) {
19148d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
192da3a660dSBarry Smith   }
193da3a660dSBarry Smith   else {
19448d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
195da3a660dSBarry Smith   }
196ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
197da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
198da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
19955659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
200da3a660dSBarry Smith   return 0;
201da3a660dSBarry Smith }
202289bc588SBarry Smith /* ------------------------------------------------------------------*/
203c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
20420e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
205289bc588SBarry Smith {
206c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
207289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
2086abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
209289bc588SBarry Smith 
210289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
211289bc588SBarry Smith     /* this is a hack fix, should have another version without
212289bc588SBarry Smith        the second BLdot */
213bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
214289bc588SBarry Smith   }
215289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
216289bc588SBarry Smith   while (its--) {
217289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
218289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
219f1747703SBarry Smith #if defined(PETSC_COMPLEX)
220f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
221f1747703SBarry Smith            not happy about returning a double complex */
222f1747703SBarry Smith         int    _i;
223f1747703SBarry Smith         Scalar sum = b[i];
224f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
225f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
226f1747703SBarry Smith         }
227f1747703SBarry Smith         xt = sum;
228f1747703SBarry Smith #else
229289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
230f1747703SBarry Smith #endif
231d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
232289bc588SBarry Smith       }
233289bc588SBarry Smith     }
234289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
235289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
236f1747703SBarry Smith #if defined(PETSC_COMPLEX)
237f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
238f1747703SBarry Smith            not happy about returning a double complex */
239f1747703SBarry Smith         int    _i;
240f1747703SBarry Smith         Scalar sum = b[i];
241f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
242f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
243f1747703SBarry Smith         }
244f1747703SBarry Smith         xt = sum;
245f1747703SBarry Smith #else
246289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
247f1747703SBarry Smith #endif
248d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
249289bc588SBarry Smith       }
250289bc588SBarry Smith     }
251289bc588SBarry Smith   }
252289bc588SBarry Smith   return 0;
253289bc588SBarry Smith }
254289bc588SBarry Smith 
255289bc588SBarry Smith /* -----------------------------------------------------------------*/
25644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
257289bc588SBarry Smith {
258c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
259289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
260289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
261289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
26248d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
26355659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
264289bc588SBarry Smith   return 0;
265289bc588SBarry Smith }
26644cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
267289bc588SBarry Smith {
268c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
269289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
270289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
271289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
27248d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
27355659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
274289bc588SBarry Smith   return 0;
275289bc588SBarry Smith }
27644cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
277289bc588SBarry Smith {
278c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
279289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2806abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
281289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
282416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
28348d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
28455659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
285289bc588SBarry Smith   return 0;
286289bc588SBarry Smith }
28744cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
288289bc588SBarry Smith {
289c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
290289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
291289bc588SBarry Smith   int          _One=1;
2926abc6512SBarry Smith   Scalar       _DOne=1.0;
29344cd7ae7SLois Curfman McInnes   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
294717eeb74SLois Curfman McInnes   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
29548d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
29655659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
297289bc588SBarry Smith   return 0;
298289bc588SBarry Smith }
299289bc588SBarry Smith 
300289bc588SBarry Smith /* -----------------------------------------------------------------*/
301c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
302289bc588SBarry Smith {
303c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
304289bc588SBarry Smith   Scalar       *v;
305289bc588SBarry Smith   int          i;
30667e560aaSBarry Smith 
307289bc588SBarry Smith   *ncols = mat->n;
308289bc588SBarry Smith   if (cols) {
3090452661fSBarry Smith     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
31073c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
311289bc588SBarry Smith   }
312289bc588SBarry Smith   if (vals) {
3130452661fSBarry Smith     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
314289bc588SBarry Smith     v = mat->v + row;
31573c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
316289bc588SBarry Smith   }
317289bc588SBarry Smith   return 0;
318289bc588SBarry Smith }
319c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
320289bc588SBarry Smith {
3210452661fSBarry Smith   if (cols) { PetscFree(*cols); }
3220452661fSBarry Smith   if (vals) { PetscFree(*vals); }
323289bc588SBarry Smith   return 0;
324289bc588SBarry Smith }
325289bc588SBarry Smith /* ----------------------------------------------------------------*/
326ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
327e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
328289bc588SBarry Smith {
329c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
330289bc588SBarry Smith   int          i,j;
331d6dfbf8fSBarry Smith 
332289bc588SBarry Smith   if (!mat->roworiented) {
333dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
334289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
335289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
336289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
337289bc588SBarry Smith         }
338289bc588SBarry Smith       }
339289bc588SBarry Smith     }
340289bc588SBarry Smith     else {
341289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
342289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
343289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
344289bc588SBarry Smith         }
345289bc588SBarry Smith       }
346289bc588SBarry Smith     }
347e8d4e0b9SBarry Smith   }
348e8d4e0b9SBarry Smith   else {
349dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
350e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
351e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
352e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
353e8d4e0b9SBarry Smith         }
354e8d4e0b9SBarry Smith       }
355e8d4e0b9SBarry Smith     }
356289bc588SBarry Smith     else {
357289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
358289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
359289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
360289bc588SBarry Smith         }
361289bc588SBarry Smith       }
362289bc588SBarry Smith     }
363e8d4e0b9SBarry Smith   }
364289bc588SBarry Smith   return 0;
365289bc588SBarry Smith }
366e8d4e0b9SBarry Smith 
367ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
368ae80bb75SLois Curfman McInnes {
369ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
370ae80bb75SLois Curfman McInnes   int          i, j;
371ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
372ae80bb75SLois Curfman McInnes 
373ae80bb75SLois Curfman McInnes   /* row-oriented output */
374ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
375ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
376ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
377ae80bb75SLois Curfman McInnes     }
378ae80bb75SLois Curfman McInnes   }
379ae80bb75SLois Curfman McInnes   return 0;
380ae80bb75SLois Curfman McInnes }
381ae80bb75SLois Curfman McInnes 
382289bc588SBarry Smith /* -----------------------------------------------------------------*/
3833d1612f7SBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
384289bc588SBarry Smith {
385c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
386289bc588SBarry Smith   int          ierr;
387289bc588SBarry Smith   Mat          newi;
38848d91487SBarry Smith 
389b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
390ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
39155659b69SBarry Smith   if (cpvalues == COPY_VALUES) {
392416022c9SBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
39355659b69SBarry Smith   }
394227d817aSBarry Smith   newi->assembled = PETSC_TRUE;
395289bc588SBarry Smith   *newmat = newi;
396289bc588SBarry Smith   return 0;
397289bc588SBarry Smith }
398289bc588SBarry Smith 
39977c4ece6SBarry Smith #include "sys.h"
40088e32edaSLois Curfman McInnes 
40119bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
40288e32edaSLois Curfman McInnes {
40388e32edaSLois Curfman McInnes   Mat_SeqDense *a;
40488e32edaSLois Curfman McInnes   Mat          B;
40588e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
40688e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
40788e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
40819bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
40988e32edaSLois Curfman McInnes 
41088e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
41188e32edaSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
41290ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
41377c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
41488e32edaSLois Curfman McInnes   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
41588e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
41688e32edaSLois Curfman McInnes 
41790ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
41890ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
41990ace30eSBarry Smith     B = *A;
42090ace30eSBarry Smith     a = (Mat_SeqDense *) B->data;
42190ace30eSBarry Smith 
42290ace30eSBarry Smith     /* read in nonzero values */
42377c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
42490ace30eSBarry Smith 
4256d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4266d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
42790ace30eSBarry Smith     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
42890ace30eSBarry Smith   } else {
42988e32edaSLois Curfman McInnes     /* read row lengths */
4300452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
43177c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
43288e32edaSLois Curfman McInnes 
43388e32edaSLois Curfman McInnes     /* create our matrix */
434b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
43588e32edaSLois Curfman McInnes     B = *A;
43688e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
43788e32edaSLois Curfman McInnes     v = a->v;
43888e32edaSLois Curfman McInnes 
43988e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
4400452661fSBarry Smith     cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
44177c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
4420452661fSBarry Smith     vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
44377c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
44488e32edaSLois Curfman McInnes 
44588e32edaSLois Curfman McInnes     /* insert into matrix */
44688e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
44788e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
44888e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
44988e32edaSLois Curfman McInnes     }
4500452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
45188e32edaSLois Curfman McInnes 
4526d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4536d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
45490ace30eSBarry Smith   }
45588e32edaSLois Curfman McInnes   return 0;
45688e32edaSLois Curfman McInnes }
45788e32edaSLois Curfman McInnes 
458932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
45977c4ece6SBarry Smith #include "sys.h"
460932b0c3eSLois Curfman McInnes 
461932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
462289bc588SBarry Smith {
463932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
464932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
465d636dbe3SBarry Smith   FILE         *fd;
466932b0c3eSLois Curfman McInnes   char         *outputname;
467932b0c3eSLois Curfman McInnes   Scalar       *v;
468932b0c3eSLois Curfman McInnes 
46990ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
470932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
47190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
4727ddc982cSLois Curfman McInnes   if (format == ASCII_FORMAT_INFO || format == ASCII_FORMAT_INFO_DETAILED) {
4737ddc982cSLois Curfman McInnes     return 0;  /* do nothing for now */
474f72585f2SLois Curfman McInnes   }
47544cd7ae7SLois Curfman McInnes   else if (format == ASCII_FORMAT_COMMON) {
47680cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
47744cd7ae7SLois Curfman McInnes       v = a->v + i;
47880cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
47980cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
48080cd9d93SLois Curfman McInnes #if defined(PETSC_COMPLEX)
48180cd9d93SLois Curfman McInnes         if (real(*v) != 0.0 && imag(*v) != 0.0)
48280cd9d93SLois Curfman McInnes           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
48380cd9d93SLois Curfman McInnes         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
48480cd9d93SLois Curfman McInnes         v += a->m;
48580cd9d93SLois Curfman McInnes #else
48680cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
48780cd9d93SLois Curfman McInnes         v += a->m;
48880cd9d93SLois Curfman McInnes #endif
48980cd9d93SLois Curfman McInnes       }
49080cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
49180cd9d93SLois Curfman McInnes     }
49280cd9d93SLois Curfman McInnes   }
493f72585f2SLois Curfman McInnes   else {
49447989497SBarry Smith #if defined(PETSC_COMPLEX)
49547989497SBarry Smith     int allreal = 1;
49647989497SBarry Smith     /* determine if matrix has all real values */
49747989497SBarry Smith     v = a->v;
49847989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
49947989497SBarry Smith       if (imag(v[i])) { allreal = 0; break ;}
50047989497SBarry Smith     }
50147989497SBarry Smith #endif
502932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
503932b0c3eSLois Curfman McInnes       v = a->v + i;
504932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
505289bc588SBarry Smith #if defined(PETSC_COMPLEX)
50647989497SBarry Smith         if (allreal) {
50747989497SBarry Smith           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
50847989497SBarry Smith         } else {
509932b0c3eSLois Curfman McInnes           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
51047989497SBarry Smith         }
511289bc588SBarry Smith #else
512932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
513289bc588SBarry Smith #endif
514289bc588SBarry Smith       }
5158e37dc09SBarry Smith       fprintf(fd,"\n");
516289bc588SBarry Smith     }
517da3a660dSBarry Smith   }
518932b0c3eSLois Curfman McInnes   fflush(fd);
519289bc588SBarry Smith   return 0;
520289bc588SBarry Smith }
521289bc588SBarry Smith 
522932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
523932b0c3eSLois Curfman McInnes {
524932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
525932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
52690ace30eSBarry Smith   int          format;
52790ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
528932b0c3eSLois Curfman McInnes 
52990ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
53090ace30eSBarry Smith 
53190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
53290ace30eSBarry Smith   if (format == BINARY_FORMAT_NATIVE) {
53390ace30eSBarry Smith     /* store the matrix as a dense matrix */
53490ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
53590ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
53690ace30eSBarry Smith     col_lens[1] = m;
53790ace30eSBarry Smith     col_lens[2] = n;
53890ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
53977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
54090ace30eSBarry Smith     PetscFree(col_lens);
54190ace30eSBarry Smith 
54290ace30eSBarry Smith     /* write out matrix, by rows */
54390ace30eSBarry Smith     vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals);
54490ace30eSBarry Smith     v    = a->v;
54590ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
54690ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
54790ace30eSBarry Smith         vals[i + j*m] = *v++;
54890ace30eSBarry Smith       }
54990ace30eSBarry Smith     }
55077c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
55190ace30eSBarry Smith     PetscFree(vals);
55290ace30eSBarry Smith   } else {
5530452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
554932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
555932b0c3eSLois Curfman McInnes     col_lens[1] = m;
556932b0c3eSLois Curfman McInnes     col_lens[2] = n;
557932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
558932b0c3eSLois Curfman McInnes 
559932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
560932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
56177c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
562932b0c3eSLois Curfman McInnes 
563932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
564932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
565932b0c3eSLois Curfman McInnes     ict = 0;
566932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
567932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
568932b0c3eSLois Curfman McInnes     }
56977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
5700452661fSBarry Smith     PetscFree(col_lens);
571932b0c3eSLois Curfman McInnes 
572932b0c3eSLois Curfman McInnes     /* store nonzero values */
5730452661fSBarry Smith     anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
574932b0c3eSLois Curfman McInnes     ict = 0;
575932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
576932b0c3eSLois Curfman McInnes       v = a->v + i;
577932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
578932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
579932b0c3eSLois Curfman McInnes       }
580932b0c3eSLois Curfman McInnes     }
58177c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
5820452661fSBarry Smith     PetscFree(anonz);
58390ace30eSBarry Smith   }
584932b0c3eSLois Curfman McInnes   return 0;
585932b0c3eSLois Curfman McInnes }
586932b0c3eSLois Curfman McInnes 
587932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
588932b0c3eSLois Curfman McInnes {
589932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
590932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
591bcd2baecSBarry Smith   ViewerType   vtype;
592bcd2baecSBarry Smith   int          ierr;
593932b0c3eSLois Curfman McInnes 
594bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
595bcd2baecSBarry Smith 
596bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
597932b0c3eSLois Curfman McInnes     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
598932b0c3eSLois Curfman McInnes   }
59919bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
600932b0c3eSLois Curfman McInnes     return MatView_SeqDense_ASCII(A,viewer);
601932b0c3eSLois Curfman McInnes   }
602bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
603932b0c3eSLois Curfman McInnes     return MatView_SeqDense_Binary(A,viewer);
604932b0c3eSLois Curfman McInnes   }
605932b0c3eSLois Curfman McInnes   return 0;
606932b0c3eSLois Curfman McInnes }
607289bc588SBarry Smith 
608ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
609289bc588SBarry Smith {
610289bc588SBarry Smith   Mat          mat = (Mat) obj;
611ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
612a5a9c739SBarry Smith #if defined(PETSC_LOG)
613a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
614a5a9c739SBarry Smith #endif
6150452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
6163439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
6170452661fSBarry Smith   PetscFree(l);
618a5a9c739SBarry Smith   PLogObjectDestroy(mat);
6190452661fSBarry Smith   PetscHeaderDestroy(mat);
620289bc588SBarry Smith   return 0;
621289bc588SBarry Smith }
622289bc588SBarry Smith 
623c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
624289bc588SBarry Smith {
625c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
626d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
627d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
62848b35521SBarry Smith 
629d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
6303638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
631d9f96c7cSLois Curfman McInnes     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
632d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
633289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
634d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
635d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
636d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
637289bc588SBarry Smith       }
638289bc588SBarry Smith     }
63948b35521SBarry Smith   }
640d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
641d3e5ee88SLois Curfman McInnes     int          ierr;
642d3e5ee88SLois Curfman McInnes     Mat          tmat;
643ec8511deSBarry Smith     Mat_SeqDense *tmatd;
644d3e5ee88SLois Curfman McInnes     Scalar       *v2;
645b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
646ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
6470de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
648d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
6490de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
650d3e5ee88SLois Curfman McInnes     }
6516d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6526d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
653d3e5ee88SLois Curfman McInnes     *matout = tmat;
65448b35521SBarry Smith   }
655289bc588SBarry Smith   return 0;
656289bc588SBarry Smith }
657289bc588SBarry Smith 
65877c4ece6SBarry Smith static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
659289bc588SBarry Smith {
660c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
661c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
662289bc588SBarry Smith   int          i;
663289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
6649ea5d5aeSSatish Balay 
6659ea5d5aeSSatish Balay   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
66677c4ece6SBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
66777c4ece6SBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
668289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
66977c4ece6SBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
670289bc588SBarry Smith     v1++; v2++;
671289bc588SBarry Smith   }
67277c4ece6SBarry Smith   *flg = PETSC_TRUE;
6739ea5d5aeSSatish Balay   return 0;
674289bc588SBarry Smith }
675289bc588SBarry Smith 
676c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
677289bc588SBarry Smith {
678c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
67944cd7ae7SLois Curfman McInnes   int          i, n, len;
68044cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
68144cd7ae7SLois Curfman McInnes 
68244cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
683289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
68444cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
685ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
68644cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
687289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
688289bc588SBarry Smith   }
689289bc588SBarry Smith   return 0;
690289bc588SBarry Smith }
691289bc588SBarry Smith 
692052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
693289bc588SBarry Smith {
694c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
695da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
696da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
69755659b69SBarry Smith 
69828988994SBarry Smith   if (ll) {
699da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
700052efed2SBarry Smith     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
701da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
702da3a660dSBarry Smith       x = l[i];
703da3a660dSBarry Smith       v = mat->v + i;
704da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
705da3a660dSBarry Smith     }
70644cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
707da3a660dSBarry Smith   }
70828988994SBarry Smith   if (rr) {
709da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
710052efed2SBarry Smith     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
711da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
712da3a660dSBarry Smith       x = r[i];
713da3a660dSBarry Smith       v = mat->v + i*m;
714da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
715da3a660dSBarry Smith     }
71644cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
717da3a660dSBarry Smith   }
718289bc588SBarry Smith   return 0;
719289bc588SBarry Smith }
720289bc588SBarry Smith 
721cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
722289bc588SBarry Smith {
723c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
724289bc588SBarry Smith   Scalar       *v = mat->v;
725289bc588SBarry Smith   double       sum = 0.0;
726289bc588SBarry Smith   int          i, j;
72755659b69SBarry Smith 
728289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
729289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
730289bc588SBarry Smith #if defined(PETSC_COMPLEX)
731289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
732289bc588SBarry Smith #else
733289bc588SBarry Smith       sum += (*v)*(*v); v++;
734289bc588SBarry Smith #endif
735289bc588SBarry Smith     }
736289bc588SBarry Smith     *norm = sqrt(sum);
73755659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
738289bc588SBarry Smith   }
739289bc588SBarry Smith   else if (type == NORM_1) {
740289bc588SBarry Smith     *norm = 0.0;
741289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
742289bc588SBarry Smith       sum = 0.0;
743289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
74433a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
745289bc588SBarry Smith       }
746289bc588SBarry Smith       if (sum > *norm) *norm = sum;
747289bc588SBarry Smith     }
74855659b69SBarry Smith     PLogFlops(mat->n*mat->m);
749289bc588SBarry Smith   }
750289bc588SBarry Smith   else if (type == NORM_INFINITY) {
751289bc588SBarry Smith     *norm = 0.0;
752289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
753289bc588SBarry Smith       v = mat->v + j;
754289bc588SBarry Smith       sum = 0.0;
755289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
756cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
757289bc588SBarry Smith       }
758289bc588SBarry Smith       if (sum > *norm) *norm = sum;
759289bc588SBarry Smith     }
76055659b69SBarry Smith     PLogFlops(mat->n*mat->m);
761289bc588SBarry Smith   }
762289bc588SBarry Smith   else {
76348d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
764289bc588SBarry Smith   }
765289bc588SBarry Smith   return 0;
766289bc588SBarry Smith }
767289bc588SBarry Smith 
768c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
769289bc588SBarry Smith {
770c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
77167e560aaSBarry Smith 
7726d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
7736d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
7746d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
7756d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
7766d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7776d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
7786d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
7796d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
7806d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
7816d4a8577SBarry Smith            op == MAT_YES_NEW_DIAGONALS)
78294a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
783c0bbcb79SLois Curfman McInnes   else
7844d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
785289bc588SBarry Smith   return 0;
786289bc588SBarry Smith }
787289bc588SBarry Smith 
788ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
7896f0a148fSBarry Smith {
790ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
791cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
7926f0a148fSBarry Smith   return 0;
7936f0a148fSBarry Smith }
7946f0a148fSBarry Smith 
795*4e220ebcSLois Curfman McInnes static int MatGetBlockSize_SeqDense(Mat A,int *bs)
796*4e220ebcSLois Curfman McInnes {
797*4e220ebcSLois Curfman McInnes   *bs = 1;
798*4e220ebcSLois Curfman McInnes   return 0;
799*4e220ebcSLois Curfman McInnes }
800*4e220ebcSLois Curfman McInnes 
801ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
8026f0a148fSBarry Smith {
803ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
8046abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
8056f0a148fSBarry Smith   Scalar       *slot;
80655659b69SBarry Smith 
80777c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
80878b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
8096f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
8106f0a148fSBarry Smith     slot = l->v + rows[i];
8116f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
8126f0a148fSBarry Smith   }
8136f0a148fSBarry Smith   if (diag) {
8146f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
8156f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
816260925b8SBarry Smith       *slot = *diag;
8176f0a148fSBarry Smith     }
8186f0a148fSBarry Smith   }
819260925b8SBarry Smith   ISRestoreIndices(is,&rows);
8206f0a148fSBarry Smith   return 0;
8216f0a148fSBarry Smith }
822557bce09SLois Curfman McInnes 
823c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
824557bce09SLois Curfman McInnes {
825c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
826557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
827557bce09SLois Curfman McInnes   return 0;
828557bce09SLois Curfman McInnes }
829557bce09SLois Curfman McInnes 
830c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
831d3e5ee88SLois Curfman McInnes {
832c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
833d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
834d3e5ee88SLois Curfman McInnes   return 0;
835d3e5ee88SLois Curfman McInnes }
836d3e5ee88SLois Curfman McInnes 
837c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
83864e87e97SBarry Smith {
839c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
84064e87e97SBarry Smith   *array = mat->v;
84164e87e97SBarry Smith   return 0;
84264e87e97SBarry Smith }
8430754003eSLois Curfman McInnes 
844ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
845ff14e315SSatish Balay {
846ff14e315SSatish Balay   return 0;
847ff14e315SSatish Balay }
8480754003eSLois Curfman McInnes 
849cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
850cddf8d76SBarry Smith                                     Mat *submat)
8510754003eSLois Curfman McInnes {
852c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
8530754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
854160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
8550754003eSLois Curfman McInnes   Scalar       *vwork, *val;
8560754003eSLois Curfman McInnes   Mat          newmat;
8570754003eSLois Curfman McInnes 
85878b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
85978b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
86078b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
86178b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
8620754003eSLois Curfman McInnes 
8630452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
8640452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
8650452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
866cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
8670754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
8680754003eSLois Curfman McInnes 
8690754003eSLois Curfman McInnes   /* Create and fill new matrix */
870b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
8710754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
8720754003eSLois Curfman McInnes     nznew = 0;
8730754003eSLois Curfman McInnes     val   = mat->v + irow[i];
8740754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
8750754003eSLois Curfman McInnes       if (smap[j]) {
8760754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
8770754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
8780754003eSLois Curfman McInnes       }
8790754003eSLois Curfman McInnes     }
880dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
88178b31e54SBarry Smith            CHKERRQ(ierr);
8820754003eSLois Curfman McInnes   }
8836d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8846d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8850754003eSLois Curfman McInnes 
8860754003eSLois Curfman McInnes   /* Free work space */
8870452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
88878b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
88978b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
8900754003eSLois Curfman McInnes   *submat = newmat;
8910754003eSLois Curfman McInnes   return 0;
8920754003eSLois Curfman McInnes }
8930754003eSLois Curfman McInnes 
894905e6a2fSBarry Smith static int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
895905e6a2fSBarry Smith                                     Mat **B)
896905e6a2fSBarry Smith {
897905e6a2fSBarry Smith   int ierr,i;
898905e6a2fSBarry Smith 
899905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
900905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
901905e6a2fSBarry Smith   }
902905e6a2fSBarry Smith 
903905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
904905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
905905e6a2fSBarry Smith   }
906905e6a2fSBarry Smith   return 0;
907905e6a2fSBarry Smith }
908905e6a2fSBarry Smith 
9094b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B)
9104b0e389bSBarry Smith {
9114b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
9124b0e389bSBarry Smith   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
9134b0e389bSBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
9144b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
9154b0e389bSBarry Smith   return 0;
9164b0e389bSBarry Smith }
9174b0e389bSBarry Smith 
918289bc588SBarry Smith /* -------------------------------------------------------------------*/
919ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
920905e6a2fSBarry Smith        MatGetRow_SeqDense,
921905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
922905e6a2fSBarry Smith        MatMult_SeqDense,
923905e6a2fSBarry Smith        MatMultAdd_SeqDense,
924905e6a2fSBarry Smith        MatMultTrans_SeqDense,
925905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
926905e6a2fSBarry Smith        MatSolve_SeqDense,
927905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
928905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
929905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
930905e6a2fSBarry Smith        MatLUFactor_SeqDense,
931905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
932ec8511deSBarry Smith        MatRelax_SeqDense,
933ec8511deSBarry Smith        MatTranspose_SeqDense,
934905e6a2fSBarry Smith        MatGetInfo_SeqDense,
935905e6a2fSBarry Smith        MatEqual_SeqDense,
936905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
937905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
938905e6a2fSBarry Smith        MatNorm_SeqDense,
939905e6a2fSBarry Smith        0,
940905e6a2fSBarry Smith        0,
941905e6a2fSBarry Smith        0,
942905e6a2fSBarry Smith        MatSetOption_SeqDense,
943905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
944905e6a2fSBarry Smith        MatZeroRows_SeqDense,
945905e6a2fSBarry Smith        0,
946905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
947905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
948905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
949905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
950905e6a2fSBarry Smith        MatGetSize_SeqDense,
951905e6a2fSBarry Smith        MatGetSize_SeqDense,
952905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
953905e6a2fSBarry Smith        0,
954905e6a2fSBarry Smith        0,
955905e6a2fSBarry Smith        MatGetArray_SeqDense,
956905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
957905e6a2fSBarry Smith        0,
9583d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
959905e6a2fSBarry Smith        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
9604b0e389bSBarry Smith        MatGetValues_SeqDense,
961*4e220ebcSLois Curfman McInnes        MatCopy_SeqDense,0,MatScale_SeqDense,
962*4e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_SeqDense};
96390ace30eSBarry Smith 
9644b828684SBarry Smith /*@C
965fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
966d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
967d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
968289bc588SBarry Smith 
96920563c6bSBarry Smith    Input Parameters:
9700c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
9710c775827SLois Curfman McInnes .  m - number of rows
97218f449edSLois Curfman McInnes .  n - number of columns
973b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
974dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
97520563c6bSBarry Smith 
97620563c6bSBarry Smith    Output Parameter:
97744cd7ae7SLois Curfman McInnes .  A - the matrix
97820563c6bSBarry Smith 
97918f449edSLois Curfman McInnes   Notes:
98018f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
98118f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
982b4fd4287SBarry Smith   set data=PETSC_NULL.
98318f449edSLois Curfman McInnes 
984dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
985d65003e9SLois Curfman McInnes 
986dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
98720563c6bSBarry Smith @*/
98844cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
989289bc588SBarry Smith {
99044cd7ae7SLois Curfman McInnes   Mat          B;
99144cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
99225cdf11fSBarry Smith   int          ierr,flg;
99355659b69SBarry Smith 
99444cd7ae7SLois Curfman McInnes   *A            = 0;
99544cd7ae7SLois Curfman McInnes   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
99644cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
99744cd7ae7SLois Curfman McInnes   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
99844cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqDense));
99944cd7ae7SLois Curfman McInnes   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
100044cd7ae7SLois Curfman McInnes   B->destroy    = MatDestroy_SeqDense;
100144cd7ae7SLois Curfman McInnes   B->view       = MatView_SeqDense;
100244cd7ae7SLois Curfman McInnes   B->factor     = 0;
100344cd7ae7SLois Curfman McInnes   PLogObjectMemory(B,sizeof(struct _Mat));
100444cd7ae7SLois Curfman McInnes   B->data       = (void *) b;
100518f449edSLois Curfman McInnes 
100644cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
100744cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
100844cd7ae7SLois Curfman McInnes   b->pivots       = 0;
100944cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1010b4fd4287SBarry Smith   if (data == PETSC_NULL) {
101144cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
101244cd7ae7SLois Curfman McInnes     PetscMemzero(b->v,m*n*sizeof(Scalar));
101344cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
101444cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
101518f449edSLois Curfman McInnes   }
10162dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
101744cd7ae7SLois Curfman McInnes     b->v = data;
101844cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
10192dd2b441SLois Curfman McInnes   }
102025cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
102144cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
1022*4e220ebcSLois Curfman McInnes 
102344cd7ae7SLois Curfman McInnes   *A = B;
1024289bc588SBarry Smith   return 0;
1025289bc588SBarry Smith }
1026289bc588SBarry Smith 
1027c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
1028289bc588SBarry Smith {
1029c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
1030b4fd4287SBarry Smith   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
1031289bc588SBarry Smith }
1032