xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 45d48df93cc74c82ddf652f3b224ea65f4ac6a2b)
1cb512458SBarry Smith #ifndef lint
2*45d48df9SBarry Smith static char vcid[] = "$Id: dense.c,v 1.115 1996/11/29 22:20:34 curfman Exp bsmith $";
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 
214e220ebcSLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
22289bc588SBarry Smith {
234e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
244e220ebcSLois Curfman McInnes   int          i,N = a->m*a->n,count = 0;
254e220ebcSLois Curfman McInnes   Scalar       *v = a->v;
26289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
274e220ebcSLois Curfman McInnes 
284e220ebcSLois Curfman McInnes   info->rows_global       = (double)a->m;
294e220ebcSLois Curfman McInnes   info->columns_global    = (double)a->n;
304e220ebcSLois Curfman McInnes   info->rows_local        = (double)a->m;
314e220ebcSLois Curfman McInnes   info->columns_local     = (double)a->n;
324e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
334e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
344e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
354e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
364e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
374e220ebcSLois Curfman McInnes   info->mallocs           = 0;
384e220ebcSLois Curfman McInnes   info->memory            = A->mem;
394e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
404e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
414e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
424e220ebcSLois 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) {
66*45d48df9SBarry Smith     mat->pivots = (int *) PetscMalloc((mat->m+1)*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");
71*45d48df9SBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"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 }
7602cad45dSBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
7702cad45dSBarry Smith {
7802cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
7902cad45dSBarry Smith   int          ierr;
8002cad45dSBarry Smith   Mat          newi;
8102cad45dSBarry Smith 
8202cad45dSBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
8302cad45dSBarry Smith   l = (Mat_SeqDense *) newi->data;
8402cad45dSBarry Smith   if (cpvalues == COPY_VALUES) {
8502cad45dSBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
8602cad45dSBarry Smith   }
8702cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
8802cad45dSBarry Smith   *newmat = newi;
8902cad45dSBarry Smith   return 0;
9002cad45dSBarry Smith }
9102cad45dSBarry Smith 
92c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
93289bc588SBarry Smith {
9402cad45dSBarry Smith   return MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);
95289bc588SBarry Smith }
96c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
97289bc588SBarry Smith {
9802cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data;
9902cad45dSBarry Smith   /* copy the numerical values */
10002cad45dSBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
10102cad45dSBarry Smith   (*fact)->factor = 0;
10249d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
103289bc588SBarry Smith }
104c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
105289bc588SBarry Smith {
106a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
107289bc588SBarry Smith }
108c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
109289bc588SBarry Smith {
11049d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
111289bc588SBarry Smith }
112c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
113289bc588SBarry Smith {
114c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1156abc6512SBarry Smith   int           info;
11655659b69SBarry Smith 
117752f5567SLois Curfman McInnes   if (mat->pivots) {
1180452661fSBarry Smith     PetscFree(mat->pivots);
119c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
120752f5567SLois Curfman McInnes     mat->pivots = 0;
121752f5567SLois Curfman McInnes   }
122289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
123*45d48df9SBarry Smith   if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"MatCholeskyFactor_SeqDense:Bad factorization");
124c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
12555659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
126289bc588SBarry Smith   return 0;
127289bc588SBarry Smith }
128289bc588SBarry Smith 
129c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
130289bc588SBarry Smith {
131c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
132a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
1336abc6512SBarry Smith   Scalar       *x, *y;
13467e560aaSBarry Smith 
135a57ad546SLois Curfman McInnes   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
136416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
137c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
13848d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
139289bc588SBarry Smith   }
140c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
14148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
142289bc588SBarry Smith   }
143ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
144ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
14555659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
146289bc588SBarry Smith   return 0;
147289bc588SBarry Smith }
148c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
149da3a660dSBarry Smith {
150c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1516abc6512SBarry Smith   int          one = 1, info;
1526abc6512SBarry Smith   Scalar       *x, *y;
15367e560aaSBarry Smith 
154da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
155416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
156752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
157da3a660dSBarry Smith   if (mat->pivots) {
15848d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
159da3a660dSBarry Smith   }
160da3a660dSBarry Smith   else {
16148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
162da3a660dSBarry Smith   }
163ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
16455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
165da3a660dSBarry Smith   return 0;
166da3a660dSBarry Smith }
167c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
168da3a660dSBarry Smith {
169c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1706abc6512SBarry Smith   int          one = 1, info,ierr;
1716abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
172da3a660dSBarry Smith   Vec          tmp = 0;
17367e560aaSBarry Smith 
174da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
175da3a660dSBarry Smith   if (yy == zz) {
17678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
177c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
17878b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
179da3a660dSBarry Smith   }
180416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
181752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
182da3a660dSBarry Smith   if (mat->pivots) {
18348d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
184da3a660dSBarry Smith   }
185da3a660dSBarry Smith   else {
18648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
187da3a660dSBarry Smith   }
188ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
189da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
190da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
19155659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
192da3a660dSBarry Smith   return 0;
193da3a660dSBarry Smith }
19467e560aaSBarry Smith 
195c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
196da3a660dSBarry Smith {
197c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1986abc6512SBarry Smith   int           one = 1, info,ierr;
1996abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
200da3a660dSBarry Smith   Vec           tmp;
20167e560aaSBarry Smith 
202da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
203da3a660dSBarry Smith   if (yy == zz) {
20478b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
205c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
20678b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
207da3a660dSBarry Smith   }
208416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
209752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
210da3a660dSBarry Smith   if (mat->pivots) {
21148d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
212da3a660dSBarry Smith   }
213da3a660dSBarry Smith   else {
21448d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
215da3a660dSBarry Smith   }
216ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
21790f02eecSBarry Smith   if (tmp) {
21890f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);  CHKERRQ(ierr);
21990f02eecSBarry Smith     ierr = VecDestroy(tmp); CHKERRQ(ierr);
22090f02eecSBarry Smith   }
22190f02eecSBarry Smith   else {
22290f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);
22390f02eecSBarry Smith   }
22455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
225da3a660dSBarry Smith   return 0;
226da3a660dSBarry Smith }
227289bc588SBarry Smith /* ------------------------------------------------------------------*/
228c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
22920e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
230289bc588SBarry Smith {
231c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
232289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
2336abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
234289bc588SBarry Smith 
235289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
236289bc588SBarry Smith     /* this is a hack fix, should have another version without
237289bc588SBarry Smith        the second BLdot */
238bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
239289bc588SBarry Smith   }
240289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
241289bc588SBarry Smith   while (its--) {
242289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
243289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
244f1747703SBarry Smith #if defined(PETSC_COMPLEX)
245f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
246f1747703SBarry Smith            not happy about returning a double complex */
247f1747703SBarry Smith         int    _i;
248f1747703SBarry Smith         Scalar sum = b[i];
249f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
250f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
251f1747703SBarry Smith         }
252f1747703SBarry Smith         xt = sum;
253f1747703SBarry Smith #else
254289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
255f1747703SBarry Smith #endif
25655a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
257289bc588SBarry Smith       }
258289bc588SBarry Smith     }
259289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
260289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
261f1747703SBarry Smith #if defined(PETSC_COMPLEX)
262f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
263f1747703SBarry Smith            not happy about returning a double complex */
264f1747703SBarry Smith         int    _i;
265f1747703SBarry Smith         Scalar sum = b[i];
266f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
267f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
268f1747703SBarry Smith         }
269f1747703SBarry Smith         xt = sum;
270f1747703SBarry Smith #else
271289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
272f1747703SBarry Smith #endif
27355a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
274289bc588SBarry Smith       }
275289bc588SBarry Smith     }
276289bc588SBarry Smith   }
277289bc588SBarry Smith   return 0;
278289bc588SBarry Smith }
279289bc588SBarry Smith 
280289bc588SBarry Smith /* -----------------------------------------------------------------*/
28144cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
282289bc588SBarry Smith {
283c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
284289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
285289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
286289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
28748d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
28855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
289289bc588SBarry Smith   return 0;
290289bc588SBarry Smith }
29144cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
292289bc588SBarry Smith {
293c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
294289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
295289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
296289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
29748d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
29855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
299289bc588SBarry Smith   return 0;
300289bc588SBarry Smith }
30144cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
302289bc588SBarry Smith {
303c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
304289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
3056abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
306289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
307416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
30848d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
30955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
310289bc588SBarry Smith   return 0;
311289bc588SBarry Smith }
31244cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
313289bc588SBarry Smith {
314c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
315289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
316289bc588SBarry Smith   int          _One=1;
3176abc6512SBarry Smith   Scalar       _DOne=1.0;
31844cd7ae7SLois Curfman McInnes   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
319717eeb74SLois Curfman McInnes   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
32048d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
32155659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
322289bc588SBarry Smith   return 0;
323289bc588SBarry Smith }
324289bc588SBarry Smith 
325289bc588SBarry Smith /* -----------------------------------------------------------------*/
326c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
327289bc588SBarry Smith {
328c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
329289bc588SBarry Smith   Scalar       *v;
330289bc588SBarry Smith   int          i;
33167e560aaSBarry Smith 
332289bc588SBarry Smith   *ncols = mat->n;
333289bc588SBarry Smith   if (cols) {
334*45d48df9SBarry Smith     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols);
33573c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
336289bc588SBarry Smith   }
337289bc588SBarry Smith   if (vals) {
338*45d48df9SBarry Smith     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals);
339289bc588SBarry Smith     v = mat->v + row;
34073c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
341289bc588SBarry Smith   }
342289bc588SBarry Smith   return 0;
343289bc588SBarry Smith }
344c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
345289bc588SBarry Smith {
3460452661fSBarry Smith   if (cols) { PetscFree(*cols); }
3470452661fSBarry Smith   if (vals) { PetscFree(*vals); }
348289bc588SBarry Smith   return 0;
349289bc588SBarry Smith }
350289bc588SBarry Smith /* ----------------------------------------------------------------*/
351ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
352e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
353289bc588SBarry Smith {
354c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
355289bc588SBarry Smith   int          i,j;
356d6dfbf8fSBarry Smith 
357289bc588SBarry Smith   if (!mat->roworiented) {
358dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
359289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
360289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
361289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
362289bc588SBarry Smith         }
363289bc588SBarry Smith       }
364289bc588SBarry Smith     }
365289bc588SBarry Smith     else {
366289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
367289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
368289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
369289bc588SBarry Smith         }
370289bc588SBarry Smith       }
371289bc588SBarry Smith     }
372e8d4e0b9SBarry Smith   }
373e8d4e0b9SBarry Smith   else {
374dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
375e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
376e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
377e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
378e8d4e0b9SBarry Smith         }
379e8d4e0b9SBarry Smith       }
380e8d4e0b9SBarry Smith     }
381289bc588SBarry Smith     else {
382289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
383289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
384289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
385289bc588SBarry Smith         }
386289bc588SBarry Smith       }
387289bc588SBarry Smith     }
388e8d4e0b9SBarry Smith   }
389289bc588SBarry Smith   return 0;
390289bc588SBarry Smith }
391e8d4e0b9SBarry Smith 
392ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
393ae80bb75SLois Curfman McInnes {
394ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
395ae80bb75SLois Curfman McInnes   int          i, j;
396ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
397ae80bb75SLois Curfman McInnes 
398ae80bb75SLois Curfman McInnes   /* row-oriented output */
399ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
400ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
401ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
402ae80bb75SLois Curfman McInnes     }
403ae80bb75SLois Curfman McInnes   }
404ae80bb75SLois Curfman McInnes   return 0;
405ae80bb75SLois Curfman McInnes }
406ae80bb75SLois Curfman McInnes 
407289bc588SBarry Smith /* -----------------------------------------------------------------*/
408289bc588SBarry Smith 
40977c4ece6SBarry Smith #include "sys.h"
41088e32edaSLois Curfman McInnes 
41119bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
41288e32edaSLois Curfman McInnes {
41388e32edaSLois Curfman McInnes   Mat_SeqDense *a;
41488e32edaSLois Curfman McInnes   Mat          B;
41588e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
41688e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
41788e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
41819bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
41988e32edaSLois Curfman McInnes 
42088e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
42188e32edaSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
42290ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
42377c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
42488e32edaSLois Curfman McInnes   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
42588e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
42688e32edaSLois Curfman McInnes 
42790ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
42890ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
42990ace30eSBarry Smith     B = *A;
43090ace30eSBarry Smith     a = (Mat_SeqDense *) B->data;
43190ace30eSBarry Smith 
43290ace30eSBarry Smith     /* read in nonzero values */
43377c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
43490ace30eSBarry Smith 
4356d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4366d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
43790ace30eSBarry Smith     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
43890ace30eSBarry Smith   } else {
43988e32edaSLois Curfman McInnes     /* read row lengths */
440*45d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
44177c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
44288e32edaSLois Curfman McInnes 
44388e32edaSLois Curfman McInnes     /* create our matrix */
444b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
44588e32edaSLois Curfman McInnes     B = *A;
44688e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
44788e32edaSLois Curfman McInnes     v = a->v;
44888e32edaSLois Curfman McInnes 
44988e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
450*45d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
45177c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
452*45d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
45377c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
45488e32edaSLois Curfman McInnes 
45588e32edaSLois Curfman McInnes     /* insert into matrix */
45688e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
45788e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
45888e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
45988e32edaSLois Curfman McInnes     }
4600452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
46188e32edaSLois Curfman McInnes 
4626d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4636d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
46490ace30eSBarry Smith   }
46588e32edaSLois Curfman McInnes   return 0;
46688e32edaSLois Curfman McInnes }
46788e32edaSLois Curfman McInnes 
468932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
46977c4ece6SBarry Smith #include "sys.h"
470932b0c3eSLois Curfman McInnes 
471932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
472289bc588SBarry Smith {
473932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
474932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
475d636dbe3SBarry Smith   FILE         *fd;
476932b0c3eSLois Curfman McInnes   char         *outputname;
477932b0c3eSLois Curfman McInnes   Scalar       *v;
478932b0c3eSLois Curfman McInnes 
47990ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
480932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
48190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
482639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
4837ddc982cSLois Curfman McInnes     return 0;  /* do nothing for now */
484f72585f2SLois Curfman McInnes   }
48502cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
48680cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
48744cd7ae7SLois Curfman McInnes       v = a->v + i;
48880cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
48980cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
49080cd9d93SLois Curfman McInnes #if defined(PETSC_COMPLEX)
49180cd9d93SLois Curfman McInnes         if (real(*v) != 0.0 && imag(*v) != 0.0)
49280cd9d93SLois Curfman McInnes           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
49380cd9d93SLois Curfman McInnes         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
49480cd9d93SLois Curfman McInnes         v += a->m;
49580cd9d93SLois Curfman McInnes #else
49680cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
49780cd9d93SLois Curfman McInnes         v += a->m;
49880cd9d93SLois Curfman McInnes #endif
49980cd9d93SLois Curfman McInnes       }
50080cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
50180cd9d93SLois Curfman McInnes     }
50280cd9d93SLois Curfman McInnes   }
503f72585f2SLois Curfman McInnes   else {
50447989497SBarry Smith #if defined(PETSC_COMPLEX)
50547989497SBarry Smith     int allreal = 1;
50647989497SBarry Smith     /* determine if matrix has all real values */
50747989497SBarry Smith     v = a->v;
50847989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
50947989497SBarry Smith       if (imag(v[i])) { allreal = 0; break ;}
51047989497SBarry Smith     }
51147989497SBarry Smith #endif
512932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
513932b0c3eSLois Curfman McInnes       v = a->v + i;
514932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
515289bc588SBarry Smith #if defined(PETSC_COMPLEX)
51647989497SBarry Smith         if (allreal) {
51747989497SBarry Smith           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
51847989497SBarry Smith         } else {
519932b0c3eSLois Curfman McInnes           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
52047989497SBarry Smith         }
521289bc588SBarry Smith #else
522932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
523289bc588SBarry Smith #endif
524289bc588SBarry Smith       }
5258e37dc09SBarry Smith       fprintf(fd,"\n");
526289bc588SBarry Smith     }
527da3a660dSBarry Smith   }
528932b0c3eSLois Curfman McInnes   fflush(fd);
529289bc588SBarry Smith   return 0;
530289bc588SBarry Smith }
531289bc588SBarry Smith 
532932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
533932b0c3eSLois Curfman McInnes {
534932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
535932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
53690ace30eSBarry Smith   int          format;
53790ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
538932b0c3eSLois Curfman McInnes 
53990ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
54090ace30eSBarry Smith 
54190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
54202cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
54390ace30eSBarry Smith     /* store the matrix as a dense matrix */
54490ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
54590ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
54690ace30eSBarry Smith     col_lens[1] = m;
54790ace30eSBarry Smith     col_lens[2] = n;
54890ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
54977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
55090ace30eSBarry Smith     PetscFree(col_lens);
55190ace30eSBarry Smith 
55290ace30eSBarry Smith     /* write out matrix, by rows */
553*45d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
55490ace30eSBarry Smith     v    = a->v;
55590ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
55690ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
55790ace30eSBarry Smith         vals[i + j*m] = *v++;
55890ace30eSBarry Smith       }
55990ace30eSBarry Smith     }
56077c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
56190ace30eSBarry Smith     PetscFree(vals);
56290ace30eSBarry Smith   } else {
5630452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
564932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
565932b0c3eSLois Curfman McInnes     col_lens[1] = m;
566932b0c3eSLois Curfman McInnes     col_lens[2] = n;
567932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
568932b0c3eSLois Curfman McInnes 
569932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
570932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
57177c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
572932b0c3eSLois Curfman McInnes 
573932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
574932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
575932b0c3eSLois Curfman McInnes     ict = 0;
576932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
577932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
578932b0c3eSLois Curfman McInnes     }
57977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
5800452661fSBarry Smith     PetscFree(col_lens);
581932b0c3eSLois Curfman McInnes 
582932b0c3eSLois Curfman McInnes     /* store nonzero values */
583*45d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
584932b0c3eSLois Curfman McInnes     ict = 0;
585932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
586932b0c3eSLois Curfman McInnes       v = a->v + i;
587932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
588932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
589932b0c3eSLois Curfman McInnes       }
590932b0c3eSLois Curfman McInnes     }
59177c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
5920452661fSBarry Smith     PetscFree(anonz);
59390ace30eSBarry Smith   }
594932b0c3eSLois Curfman McInnes   return 0;
595932b0c3eSLois Curfman McInnes }
596932b0c3eSLois Curfman McInnes 
597932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
598932b0c3eSLois Curfman McInnes {
599932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
600932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
601bcd2baecSBarry Smith   ViewerType   vtype;
602bcd2baecSBarry Smith   int          ierr;
603932b0c3eSLois Curfman McInnes 
604bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
605bcd2baecSBarry Smith 
606bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
607932b0c3eSLois Curfman McInnes     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
608932b0c3eSLois Curfman McInnes   }
60919bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
610932b0c3eSLois Curfman McInnes     return MatView_SeqDense_ASCII(A,viewer);
611932b0c3eSLois Curfman McInnes   }
612bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
613932b0c3eSLois Curfman McInnes     return MatView_SeqDense_Binary(A,viewer);
614932b0c3eSLois Curfman McInnes   }
615932b0c3eSLois Curfman McInnes   return 0;
616932b0c3eSLois Curfman McInnes }
617289bc588SBarry Smith 
618ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
619289bc588SBarry Smith {
620289bc588SBarry Smith   Mat          mat = (Mat) obj;
621ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
62290f02eecSBarry Smith   int          ierr;
62390f02eecSBarry Smith 
624a5a9c739SBarry Smith #if defined(PETSC_LOG)
625a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
626a5a9c739SBarry Smith #endif
6270452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
6283439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
6290452661fSBarry Smith   PetscFree(l);
63090f02eecSBarry Smith   if (mat->mapping) {
63190f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
63290f02eecSBarry Smith   }
633a5a9c739SBarry Smith   PLogObjectDestroy(mat);
6340452661fSBarry Smith   PetscHeaderDestroy(mat);
635289bc588SBarry Smith   return 0;
636289bc588SBarry Smith }
637289bc588SBarry Smith 
638c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
639289bc588SBarry Smith {
640c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
641d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
642d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
64348b35521SBarry Smith 
644d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
6453638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
646*45d48df9SBarry Smith     if (m != n) SETERRQ(PETSC_ERR_SUP,"MatTranspose_SeqDense:Supports square matrix only in-place");
647d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
648289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
649d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
650d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
651d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
652289bc588SBarry Smith       }
653289bc588SBarry Smith     }
65448b35521SBarry Smith   }
655d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
656d3e5ee88SLois Curfman McInnes     int          ierr;
657d3e5ee88SLois Curfman McInnes     Mat          tmat;
658ec8511deSBarry Smith     Mat_SeqDense *tmatd;
659d3e5ee88SLois Curfman McInnes     Scalar       *v2;
660b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
661ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
6620de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
663d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
6640de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
665d3e5ee88SLois Curfman McInnes     }
6666d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6676d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
668d3e5ee88SLois Curfman McInnes     *matout = tmat;
66948b35521SBarry Smith   }
670289bc588SBarry Smith   return 0;
671289bc588SBarry Smith }
672289bc588SBarry Smith 
67377c4ece6SBarry Smith static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
674289bc588SBarry Smith {
675c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
676c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
677289bc588SBarry Smith   int          i;
678289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
6799ea5d5aeSSatish Balay 
680*45d48df9SBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
68177c4ece6SBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
68277c4ece6SBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
683289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
68477c4ece6SBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
685289bc588SBarry Smith     v1++; v2++;
686289bc588SBarry Smith   }
68777c4ece6SBarry Smith   *flg = PETSC_TRUE;
6889ea5d5aeSSatish Balay   return 0;
689289bc588SBarry Smith }
690289bc588SBarry Smith 
691c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
692289bc588SBarry Smith {
693c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
69444cd7ae7SLois Curfman McInnes   int          i, n, len;
69544cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
69644cd7ae7SLois Curfman McInnes 
69744cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
698289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
69944cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
700*45d48df9SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_SIZ,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
70144cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
702289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
703289bc588SBarry Smith   }
704289bc588SBarry Smith   return 0;
705289bc588SBarry Smith }
706289bc588SBarry Smith 
707052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
708289bc588SBarry Smith {
709c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
710da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
711da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
71255659b69SBarry Smith 
71328988994SBarry Smith   if (ll) {
714da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
715*45d48df9SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_SIZ,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
716da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
717da3a660dSBarry Smith       x = l[i];
718da3a660dSBarry Smith       v = mat->v + i;
719da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
720da3a660dSBarry Smith     }
72144cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
722da3a660dSBarry Smith   }
72328988994SBarry Smith   if (rr) {
724da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
725*45d48df9SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_SIZ,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
726da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
727da3a660dSBarry Smith       x = r[i];
728da3a660dSBarry Smith       v = mat->v + i*m;
729da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
730da3a660dSBarry Smith     }
73144cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
732da3a660dSBarry Smith   }
733289bc588SBarry Smith   return 0;
734289bc588SBarry Smith }
735289bc588SBarry Smith 
736cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
737289bc588SBarry Smith {
738c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
739289bc588SBarry Smith   Scalar       *v = mat->v;
740289bc588SBarry Smith   double       sum = 0.0;
741289bc588SBarry Smith   int          i, j;
74255659b69SBarry Smith 
743289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
744289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
745289bc588SBarry Smith #if defined(PETSC_COMPLEX)
746289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
747289bc588SBarry Smith #else
748289bc588SBarry Smith       sum += (*v)*(*v); v++;
749289bc588SBarry Smith #endif
750289bc588SBarry Smith     }
751289bc588SBarry Smith     *norm = sqrt(sum);
75255659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
753289bc588SBarry Smith   }
754289bc588SBarry Smith   else if (type == NORM_1) {
755289bc588SBarry Smith     *norm = 0.0;
756289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
757289bc588SBarry Smith       sum = 0.0;
758289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
75933a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
760289bc588SBarry Smith       }
761289bc588SBarry Smith       if (sum > *norm) *norm = sum;
762289bc588SBarry Smith     }
76355659b69SBarry Smith     PLogFlops(mat->n*mat->m);
764289bc588SBarry Smith   }
765289bc588SBarry Smith   else if (type == NORM_INFINITY) {
766289bc588SBarry Smith     *norm = 0.0;
767289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
768289bc588SBarry Smith       v = mat->v + j;
769289bc588SBarry Smith       sum = 0.0;
770289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
771cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
772289bc588SBarry Smith       }
773289bc588SBarry Smith       if (sum > *norm) *norm = sum;
774289bc588SBarry Smith     }
77555659b69SBarry Smith     PLogFlops(mat->n*mat->m);
776289bc588SBarry Smith   }
777289bc588SBarry Smith   else {
778*45d48df9SBarry Smith     SETERRQ(PETSC_ERR_SUP,"MatNorm_SeqDense:No two norm");
779289bc588SBarry Smith   }
780289bc588SBarry Smith   return 0;
781289bc588SBarry Smith }
782289bc588SBarry Smith 
783c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
784289bc588SBarry Smith {
785c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
78667e560aaSBarry Smith 
7876d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
7886d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
7896d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
790219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7916d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
792219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
7936d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7946d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
7956d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
7966d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
7976d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
79890f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
79990f02eecSBarry Smith            op == MAT_IGNORE_OFF_PROCESSOR_ENTRIES)
80094a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
801c0bbcb79SLois Curfman McInnes   else
8024d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
803289bc588SBarry Smith   return 0;
804289bc588SBarry Smith }
805289bc588SBarry Smith 
806ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
8076f0a148fSBarry Smith {
808ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
809cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
8106f0a148fSBarry Smith   return 0;
8116f0a148fSBarry Smith }
8126f0a148fSBarry Smith 
8134e220ebcSLois Curfman McInnes static int MatGetBlockSize_SeqDense(Mat A,int *bs)
8144e220ebcSLois Curfman McInnes {
8154e220ebcSLois Curfman McInnes   *bs = 1;
8164e220ebcSLois Curfman McInnes   return 0;
8174e220ebcSLois Curfman McInnes }
8184e220ebcSLois Curfman McInnes 
819ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
8206f0a148fSBarry Smith {
821ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
8226abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
8236f0a148fSBarry Smith   Scalar       *slot;
82455659b69SBarry Smith 
82577c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
82678b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
8276f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
8286f0a148fSBarry Smith     slot = l->v + rows[i];
8296f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
8306f0a148fSBarry Smith   }
8316f0a148fSBarry Smith   if (diag) {
8326f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
8336f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
834260925b8SBarry Smith       *slot = *diag;
8356f0a148fSBarry Smith     }
8366f0a148fSBarry Smith   }
837260925b8SBarry Smith   ISRestoreIndices(is,&rows);
8386f0a148fSBarry Smith   return 0;
8396f0a148fSBarry Smith }
840557bce09SLois Curfman McInnes 
841c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
842557bce09SLois Curfman McInnes {
843c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
844557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
845557bce09SLois Curfman McInnes   return 0;
846557bce09SLois Curfman McInnes }
847557bce09SLois Curfman McInnes 
848c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
849d3e5ee88SLois Curfman McInnes {
850c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
851d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
852d3e5ee88SLois Curfman McInnes   return 0;
853d3e5ee88SLois Curfman McInnes }
854d3e5ee88SLois Curfman McInnes 
855c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
85664e87e97SBarry Smith {
857c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
85864e87e97SBarry Smith   *array = mat->v;
85964e87e97SBarry Smith   return 0;
86064e87e97SBarry Smith }
8610754003eSLois Curfman McInnes 
862ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
863ff14e315SSatish Balay {
864ff14e315SSatish Balay   return 0;
865ff14e315SSatish Balay }
8660754003eSLois Curfman McInnes 
867cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
868cddf8d76SBarry Smith                                     Mat *submat)
8690754003eSLois Curfman McInnes {
870c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
8710754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
872160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
8730754003eSLois Curfman McInnes   Scalar       *vwork, *val;
8740754003eSLois Curfman McInnes   Mat          newmat;
8750754003eSLois Curfman McInnes 
87678b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
87778b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
87878b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
87978b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
8800754003eSLois Curfman McInnes 
8810452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
8820452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
8830452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
884cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
8850754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
8860754003eSLois Curfman McInnes 
8870754003eSLois Curfman McInnes   /* Create and fill new matrix */
888b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
8890754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
8900754003eSLois Curfman McInnes     nznew = 0;
8910754003eSLois Curfman McInnes     val   = mat->v + irow[i];
8920754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
8930754003eSLois Curfman McInnes       if (smap[j]) {
8940754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
8950754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
8960754003eSLois Curfman McInnes       }
8970754003eSLois Curfman McInnes     }
898dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
89978b31e54SBarry Smith            CHKERRQ(ierr);
9000754003eSLois Curfman McInnes   }
9016d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9026d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9030754003eSLois Curfman McInnes 
9040754003eSLois Curfman McInnes   /* Free work space */
9050452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
90678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
90778b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
9080754003eSLois Curfman McInnes   *submat = newmat;
9090754003eSLois Curfman McInnes   return 0;
9100754003eSLois Curfman McInnes }
9110754003eSLois Curfman McInnes 
912905e6a2fSBarry Smith static int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
913905e6a2fSBarry Smith                                     Mat **B)
914905e6a2fSBarry Smith {
915905e6a2fSBarry Smith   int ierr,i;
916905e6a2fSBarry Smith 
917905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
918905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
919905e6a2fSBarry Smith   }
920905e6a2fSBarry Smith 
921905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
922905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
923905e6a2fSBarry Smith   }
924905e6a2fSBarry Smith   return 0;
925905e6a2fSBarry Smith }
926905e6a2fSBarry Smith 
9274b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B)
9284b0e389bSBarry Smith {
9294b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
9304b0e389bSBarry Smith   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
931*45d48df9SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_SIZ,"MatCopy_SeqDense:size(B) != size(A)");
9324b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
9334b0e389bSBarry Smith   return 0;
9344b0e389bSBarry Smith }
9354b0e389bSBarry Smith 
936289bc588SBarry Smith /* -------------------------------------------------------------------*/
937ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
938905e6a2fSBarry Smith        MatGetRow_SeqDense,
939905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
940905e6a2fSBarry Smith        MatMult_SeqDense,
941905e6a2fSBarry Smith        MatMultAdd_SeqDense,
942905e6a2fSBarry Smith        MatMultTrans_SeqDense,
943905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
944905e6a2fSBarry Smith        MatSolve_SeqDense,
945905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
946905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
947905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
948905e6a2fSBarry Smith        MatLUFactor_SeqDense,
949905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
950ec8511deSBarry Smith        MatRelax_SeqDense,
951ec8511deSBarry Smith        MatTranspose_SeqDense,
952905e6a2fSBarry Smith        MatGetInfo_SeqDense,
953905e6a2fSBarry Smith        MatEqual_SeqDense,
954905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
955905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
956905e6a2fSBarry Smith        MatNorm_SeqDense,
957905e6a2fSBarry Smith        0,
958905e6a2fSBarry Smith        0,
959905e6a2fSBarry Smith        0,
960905e6a2fSBarry Smith        MatSetOption_SeqDense,
961905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
962905e6a2fSBarry Smith        MatZeroRows_SeqDense,
963905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
964905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
965905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
966905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
967905e6a2fSBarry Smith        MatGetSize_SeqDense,
968905e6a2fSBarry Smith        MatGetSize_SeqDense,
969905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
970905e6a2fSBarry Smith        0,
971905e6a2fSBarry Smith        0,
972905e6a2fSBarry Smith        MatGetArray_SeqDense,
973905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
974905e6a2fSBarry Smith        0,
9753d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
976905e6a2fSBarry Smith        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
9774b0e389bSBarry Smith        MatGetValues_SeqDense,
9784e220ebcSLois Curfman McInnes        MatCopy_SeqDense,0,MatScale_SeqDense,
9794e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_SeqDense};
98090ace30eSBarry Smith 
9814b828684SBarry Smith /*@C
982fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
983d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
984d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
985289bc588SBarry Smith 
98620563c6bSBarry Smith    Input Parameters:
9870c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
9880c775827SLois Curfman McInnes .  m - number of rows
98918f449edSLois Curfman McInnes .  n - number of columns
990b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
991dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
99220563c6bSBarry Smith 
99320563c6bSBarry Smith    Output Parameter:
99444cd7ae7SLois Curfman McInnes .  A - the matrix
99520563c6bSBarry Smith 
99618f449edSLois Curfman McInnes   Notes:
99718f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
99818f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
999b4fd4287SBarry Smith   set data=PETSC_NULL.
100018f449edSLois Curfman McInnes 
1001dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1002d65003e9SLois Curfman McInnes 
1003dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
100420563c6bSBarry Smith @*/
100544cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1006289bc588SBarry Smith {
100744cd7ae7SLois Curfman McInnes   Mat          B;
100844cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
10093b2fbd54SBarry Smith   int          ierr,flg,size;
10103b2fbd54SBarry Smith 
10113b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
10123b2fbd54SBarry Smith   if (size > 1) SETERRQ(1,"MatCreateSeqDense:Comm must be of size 1");
101355659b69SBarry Smith 
101444cd7ae7SLois Curfman McInnes   *A            = 0;
101544cd7ae7SLois Curfman McInnes   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
101644cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
101744cd7ae7SLois Curfman McInnes   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
101844cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqDense));
101944cd7ae7SLois Curfman McInnes   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
102044cd7ae7SLois Curfman McInnes   B->destroy    = MatDestroy_SeqDense;
102144cd7ae7SLois Curfman McInnes   B->view       = MatView_SeqDense;
102244cd7ae7SLois Curfman McInnes   B->factor     = 0;
102390f02eecSBarry Smith   B->mapping    = 0;
102444cd7ae7SLois Curfman McInnes   PLogObjectMemory(B,sizeof(struct _Mat));
102544cd7ae7SLois Curfman McInnes   B->data       = (void *) b;
102618f449edSLois Curfman McInnes 
102744cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
102844cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
102944cd7ae7SLois Curfman McInnes   b->pivots       = 0;
103044cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1031b4fd4287SBarry Smith   if (data == PETSC_NULL) {
103244cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
103344cd7ae7SLois Curfman McInnes     PetscMemzero(b->v,m*n*sizeof(Scalar));
103444cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
103544cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
103618f449edSLois Curfman McInnes   }
10372dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
103844cd7ae7SLois Curfman McInnes     b->v = data;
103944cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
10402dd2b441SLois Curfman McInnes   }
104125cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
104244cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
10434e220ebcSLois Curfman McInnes 
104444cd7ae7SLois Curfman McInnes   *A = B;
1045289bc588SBarry Smith   return 0;
1046289bc588SBarry Smith }
1047289bc588SBarry Smith 
1048c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
1049289bc588SBarry Smith {
1050c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
1051b4fd4287SBarry Smith   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
1052289bc588SBarry Smith }
10533b2fbd54SBarry Smith 
10543b2fbd54SBarry Smith 
10553b2fbd54SBarry Smith 
10563b2fbd54SBarry Smith 
10573b2fbd54SBarry Smith 
1058