xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 2b36279985d181bf1ae28d424577ae2c2fbebdf4)
1cb512458SBarry Smith #ifndef lint
2*2b362799SSatish Balay static char vcid[] = "$Id: dense.c,v 1.125 1997/03/13 16:33:41 curfman Exp balay $";
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 
125615d1e5SSatish Balay #undef __FUNC__
135615d1e5SSatish Balay #define __FUNC__ "MatAXPY_SeqDense"
141987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
151987afe7SBarry Smith {
161987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
171987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
181987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
191987afe7SBarry Smith   PLogFlops(2*N-1);
201987afe7SBarry Smith   return 0;
211987afe7SBarry Smith }
221987afe7SBarry Smith 
235615d1e5SSatish Balay #undef __FUNC__
245615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqDense"
258f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
26289bc588SBarry Smith {
274e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
284e220ebcSLois Curfman McInnes   int          i,N = a->m*a->n,count = 0;
294e220ebcSLois Curfman McInnes   Scalar       *v = a->v;
30289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
314e220ebcSLois Curfman McInnes 
324e220ebcSLois Curfman McInnes   info->rows_global       = (double)a->m;
334e220ebcSLois Curfman McInnes   info->columns_global    = (double)a->n;
344e220ebcSLois Curfman McInnes   info->rows_local        = (double)a->m;
354e220ebcSLois Curfman McInnes   info->columns_local     = (double)a->n;
364e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
374e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
384e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
394e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
404e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
414e220ebcSLois Curfman McInnes   info->mallocs           = 0;
424e220ebcSLois Curfman McInnes   info->memory            = A->mem;
434e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
444e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
454e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
464e220ebcSLois Curfman McInnes 
47fa9ec3f1SLois Curfman McInnes   return 0;
48289bc588SBarry Smith }
49289bc588SBarry Smith 
505615d1e5SSatish Balay #undef __FUNC__
515615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqDense"
528f6be9afSLois Curfman McInnes int MatScale_SeqDense(Scalar *alpha,Mat inA)
5380cd9d93SLois Curfman McInnes {
5480cd9d93SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
5580cd9d93SLois Curfman McInnes   int          one = 1, nz;
5680cd9d93SLois Curfman McInnes 
5780cd9d93SLois Curfman McInnes   nz = a->m*a->n;
5880cd9d93SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
5980cd9d93SLois Curfman McInnes   PLogFlops(nz);
6080cd9d93SLois Curfman McInnes   return 0;
6180cd9d93SLois Curfman McInnes }
6280cd9d93SLois Curfman McInnes 
63289bc588SBarry Smith /* ---------------------------------------------------------------*/
64289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
65289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
665615d1e5SSatish Balay #undef __FUNC__
675615d1e5SSatish Balay #define __FUNC__ "MatLUFactor_SeqDense"
688f6be9afSLois Curfman McInnes int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
69289bc588SBarry Smith {
70c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
716abc6512SBarry Smith   int          info;
7255659b69SBarry Smith 
73289bc588SBarry Smith   if (!mat->pivots) {
7445d48df9SBarry Smith     mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots);
75c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
76289bc588SBarry Smith   }
77289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
78e3372554SBarry Smith   if (info<0) SETERRQ(1,0,"Bad argument to LU factorization");
79e3372554SBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Bad LU factorization");
80c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
8155659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
82289bc588SBarry Smith   return 0;
83289bc588SBarry Smith }
846ee01492SSatish Balay 
855615d1e5SSatish Balay #undef __FUNC__
865615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqDense"
878f6be9afSLois Curfman McInnes int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
8802cad45dSBarry Smith {
8902cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
9002cad45dSBarry Smith   int          ierr;
9102cad45dSBarry Smith   Mat          newi;
9202cad45dSBarry Smith 
9302cad45dSBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
9402cad45dSBarry Smith   l = (Mat_SeqDense *) newi->data;
9502cad45dSBarry Smith   if (cpvalues == COPY_VALUES) {
9602cad45dSBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
9702cad45dSBarry Smith   }
9802cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
9902cad45dSBarry Smith   *newmat = newi;
10002cad45dSBarry Smith   return 0;
10102cad45dSBarry Smith }
10202cad45dSBarry Smith 
1035615d1e5SSatish Balay #undef __FUNC__
1045615d1e5SSatish Balay #define __FUNC__ "MatLUFactorSymbolic_SeqDense"
1058f6be9afSLois Curfman McInnes int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
106289bc588SBarry Smith {
10702cad45dSBarry Smith   return MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);
108289bc588SBarry Smith }
1096ee01492SSatish Balay 
1105615d1e5SSatish Balay #undef __FUNC__
1115615d1e5SSatish Balay #define __FUNC__ "MatLUFactorNumeric_SeqDense"
1128f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
113289bc588SBarry Smith {
11402cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data;
11502cad45dSBarry Smith   /* copy the numerical values */
11602cad45dSBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
11702cad45dSBarry Smith   (*fact)->factor = 0;
11849d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
119289bc588SBarry Smith }
1206ee01492SSatish Balay 
1215615d1e5SSatish Balay #undef __FUNC__
1225615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense"
1238f6be9afSLois Curfman McInnes int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
124289bc588SBarry Smith {
125a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
126289bc588SBarry Smith }
1276ee01492SSatish Balay 
1285615d1e5SSatish Balay #undef __FUNC__
1295615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense"
1308f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
131289bc588SBarry Smith {
13249d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
133289bc588SBarry Smith }
1346ee01492SSatish Balay 
1355615d1e5SSatish Balay #undef __FUNC__
1365615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense"
1378f6be9afSLois Curfman McInnes int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
138289bc588SBarry Smith {
139c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1406abc6512SBarry Smith   int           info;
14155659b69SBarry Smith 
142752f5567SLois Curfman McInnes   if (mat->pivots) {
1430452661fSBarry Smith     PetscFree(mat->pivots);
144c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
145752f5567SLois Curfman McInnes     mat->pivots = 0;
146752f5567SLois Curfman McInnes   }
147289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
148e3372554SBarry Smith   if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization");
149c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
15055659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
151289bc588SBarry Smith   return 0;
152289bc588SBarry Smith }
153289bc588SBarry Smith 
1545615d1e5SSatish Balay #undef __FUNC__
1555615d1e5SSatish Balay #define __FUNC__ "MatSolve_SeqDense"
1568f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
157289bc588SBarry Smith {
158c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
159a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
1606abc6512SBarry Smith   Scalar       *x, *y;
16167e560aaSBarry Smith 
162a57ad546SLois Curfman McInnes   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
163416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
164c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
16548d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
166289bc588SBarry Smith   }
167c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
16848d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
169289bc588SBarry Smith   }
170e3372554SBarry Smith   else SETERRQ(1,0,"Matrix must be factored to solve");
171e3372554SBarry Smith   if (info) SETERRQ(1,0,"MBad solve");
17255659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
173289bc588SBarry Smith   return 0;
174289bc588SBarry Smith }
1756ee01492SSatish Balay 
1765615d1e5SSatish Balay #undef __FUNC__
1775615d1e5SSatish Balay #define __FUNC__ "MatSolveTrans_SeqDense"
1788f6be9afSLois Curfman McInnes int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
179da3a660dSBarry Smith {
180c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1816abc6512SBarry Smith   int          one = 1, info;
1826abc6512SBarry Smith   Scalar       *x, *y;
18367e560aaSBarry Smith 
184da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
185416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
186752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
187da3a660dSBarry Smith   if (mat->pivots) {
18848d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
189da3a660dSBarry Smith   }
190da3a660dSBarry Smith   else {
19148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
192da3a660dSBarry Smith   }
193e3372554SBarry Smith   if (info) SETERRQ(1,0,"Bad solve");
19455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
195da3a660dSBarry Smith   return 0;
196da3a660dSBarry Smith }
1976ee01492SSatish Balay 
1985615d1e5SSatish Balay #undef __FUNC__
1995615d1e5SSatish Balay #define __FUNC__ "MatSolveAdd_SeqDense"
2008f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
201da3a660dSBarry Smith {
202c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2036abc6512SBarry Smith   int          one = 1, info,ierr;
2046abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
205da3a660dSBarry Smith   Vec          tmp = 0;
20667e560aaSBarry Smith 
207da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
208da3a660dSBarry Smith   if (yy == zz) {
20978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
210c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
21178b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
212da3a660dSBarry Smith   }
213416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
214752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
215da3a660dSBarry Smith   if (mat->pivots) {
21648d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
217da3a660dSBarry Smith   }
218da3a660dSBarry Smith   else {
21948d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
220da3a660dSBarry Smith   }
221e3372554SBarry Smith   if (info) SETERRQ(1,0,"Bad solve");
222da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
223da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
22455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
225da3a660dSBarry Smith   return 0;
226da3a660dSBarry Smith }
22767e560aaSBarry Smith 
2285615d1e5SSatish Balay #undef __FUNC__
2295615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense"
2308f6be9afSLois Curfman McInnes int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
231da3a660dSBarry Smith {
232c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
2336abc6512SBarry Smith   int           one = 1, info,ierr;
2346abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
235da3a660dSBarry Smith   Vec           tmp;
23667e560aaSBarry Smith 
237da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
238da3a660dSBarry Smith   if (yy == zz) {
23978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
240c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
24178b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
242da3a660dSBarry Smith   }
243416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
244752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
245da3a660dSBarry Smith   if (mat->pivots) {
24648d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
247da3a660dSBarry Smith   }
248da3a660dSBarry Smith   else {
24948d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
250da3a660dSBarry Smith   }
251e3372554SBarry Smith   if (info) SETERRQ(1,0,"Bad solve");
25290f02eecSBarry Smith   if (tmp) {
25390f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);  CHKERRQ(ierr);
25490f02eecSBarry Smith     ierr = VecDestroy(tmp); CHKERRQ(ierr);
25590f02eecSBarry Smith   }
25690f02eecSBarry Smith   else {
25790f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);
25890f02eecSBarry Smith   }
25955659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
260da3a660dSBarry Smith   return 0;
261da3a660dSBarry Smith }
262289bc588SBarry Smith /* ------------------------------------------------------------------*/
2635615d1e5SSatish Balay #undef __FUNC__
2645615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense"
2658f6be9afSLois Curfman McInnes int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
26620e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
267289bc588SBarry Smith {
268c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
269289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
2706abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
271289bc588SBarry Smith 
272289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
273289bc588SBarry Smith     /* this is a hack fix, should have another version without
274289bc588SBarry Smith        the second BLdot */
275bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
276289bc588SBarry Smith   }
277289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
278289bc588SBarry Smith   while (its--) {
279289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
280289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
281f1747703SBarry Smith #if defined(PETSC_COMPLEX)
282f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
283f1747703SBarry Smith            not happy about returning a double complex */
284f1747703SBarry Smith         int    _i;
285f1747703SBarry Smith         Scalar sum = b[i];
286f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
287f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
288f1747703SBarry Smith         }
289f1747703SBarry Smith         xt = sum;
290f1747703SBarry Smith #else
291289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
292f1747703SBarry Smith #endif
29355a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
294289bc588SBarry Smith       }
295289bc588SBarry Smith     }
296289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
297289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
298f1747703SBarry Smith #if defined(PETSC_COMPLEX)
299f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
300f1747703SBarry Smith            not happy about returning a double complex */
301f1747703SBarry Smith         int    _i;
302f1747703SBarry Smith         Scalar sum = b[i];
303f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
304f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
305f1747703SBarry Smith         }
306f1747703SBarry Smith         xt = sum;
307f1747703SBarry Smith #else
308289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
309f1747703SBarry Smith #endif
31055a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
311289bc588SBarry Smith       }
312289bc588SBarry Smith     }
313289bc588SBarry Smith   }
314289bc588SBarry Smith   return 0;
315289bc588SBarry Smith }
316289bc588SBarry Smith 
317289bc588SBarry Smith /* -----------------------------------------------------------------*/
3185615d1e5SSatish Balay #undef __FUNC__
3195615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqDense"
32044cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
321289bc588SBarry Smith {
322c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
323289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
324289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
325289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
32648d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
32755659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
328289bc588SBarry Smith   return 0;
329289bc588SBarry Smith }
3306ee01492SSatish Balay 
3315615d1e5SSatish Balay #undef __FUNC__
3325615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqDense"
33344cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
334289bc588SBarry Smith {
335c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
336289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
337289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
338289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
33948d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
34055659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
341289bc588SBarry Smith   return 0;
342289bc588SBarry Smith }
3436ee01492SSatish Balay 
3445615d1e5SSatish Balay #undef __FUNC__
3455615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqDense"
34644cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
347289bc588SBarry Smith {
348c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
349289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
3506abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
351289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
352416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
35348d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
35455659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
355289bc588SBarry Smith   return 0;
356289bc588SBarry Smith }
3576ee01492SSatish Balay 
3585615d1e5SSatish Balay #undef __FUNC__
3595615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense"
36044cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
361289bc588SBarry Smith {
362c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
363289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
364289bc588SBarry Smith   int          _One=1;
3656abc6512SBarry Smith   Scalar       _DOne=1.0;
36644cd7ae7SLois Curfman McInnes   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
367717eeb74SLois Curfman McInnes   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
36848d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
36955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
370289bc588SBarry Smith   return 0;
371289bc588SBarry Smith }
372289bc588SBarry Smith 
373289bc588SBarry Smith /* -----------------------------------------------------------------*/
3745615d1e5SSatish Balay #undef __FUNC__
3755615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense"
3768f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
377289bc588SBarry Smith {
378c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
379289bc588SBarry Smith   Scalar       *v;
380289bc588SBarry Smith   int          i;
38167e560aaSBarry Smith 
382289bc588SBarry Smith   *ncols = mat->n;
383289bc588SBarry Smith   if (cols) {
38445d48df9SBarry Smith     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols);
38573c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
386289bc588SBarry Smith   }
387289bc588SBarry Smith   if (vals) {
38845d48df9SBarry Smith     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals);
389289bc588SBarry Smith     v = mat->v + row;
39073c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
391289bc588SBarry Smith   }
392289bc588SBarry Smith   return 0;
393289bc588SBarry Smith }
3946ee01492SSatish Balay 
3955615d1e5SSatish Balay #undef __FUNC__
3965615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense"
3978f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
398289bc588SBarry Smith {
3990452661fSBarry Smith   if (cols) { PetscFree(*cols); }
4000452661fSBarry Smith   if (vals) { PetscFree(*vals); }
401289bc588SBarry Smith   return 0;
402289bc588SBarry Smith }
403289bc588SBarry Smith /* ----------------------------------------------------------------*/
4045615d1e5SSatish Balay #undef __FUNC__
4055615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense"
4068f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
407e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
408289bc588SBarry Smith {
409c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
410289bc588SBarry Smith   int          i,j;
411d6dfbf8fSBarry Smith 
412289bc588SBarry Smith   if (!mat->roworiented) {
413dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
414289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
41558804f6dSBarry Smith #if defined(PETSC_BOPT_g)
41658804f6dSBarry Smith         if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
41758804f6dSBarry Smith         if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
41858804f6dSBarry Smith #endif
419289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
42058804f6dSBarry Smith #if defined(PETSC_BOPT_g)
42158804f6dSBarry Smith           if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
42258804f6dSBarry Smith           if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
42358804f6dSBarry Smith #endif
424289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
425289bc588SBarry Smith         }
426289bc588SBarry Smith       }
427289bc588SBarry Smith     }
428289bc588SBarry Smith     else {
429289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
43058804f6dSBarry Smith #if defined(PETSC_BOPT_g)
43158804f6dSBarry Smith         if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
43258804f6dSBarry Smith         if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
43358804f6dSBarry Smith #endif
434289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
43558804f6dSBarry Smith #if defined(PETSC_BOPT_g)
43658804f6dSBarry Smith           if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
43758804f6dSBarry Smith           if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
43858804f6dSBarry Smith #endif
439289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
440289bc588SBarry Smith         }
441289bc588SBarry Smith       }
442289bc588SBarry Smith     }
443e8d4e0b9SBarry Smith   }
444e8d4e0b9SBarry Smith   else {
445dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
446e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
44758804f6dSBarry Smith #if defined(PETSC_BOPT_g)
44858804f6dSBarry Smith         if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
44958804f6dSBarry Smith         if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
45058804f6dSBarry Smith #endif
451e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
45258804f6dSBarry Smith #if defined(PETSC_BOPT_g)
45358804f6dSBarry Smith           if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
45458804f6dSBarry Smith           if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
45558804f6dSBarry Smith #endif
456e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
457e8d4e0b9SBarry Smith         }
458e8d4e0b9SBarry Smith       }
459e8d4e0b9SBarry Smith     }
460289bc588SBarry Smith     else {
461289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
46258804f6dSBarry Smith #if defined(PETSC_BOPT_g)
46358804f6dSBarry Smith         if (indexm[i] < 0) SETERRQ(1,0,"Negative row");
46458804f6dSBarry Smith         if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large");
46558804f6dSBarry Smith #endif
466289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
46758804f6dSBarry Smith #if defined(PETSC_BOPT_g)
46858804f6dSBarry Smith           if (indexn[j] < 0) SETERRQ(1,0,"Negative column");
46958804f6dSBarry Smith           if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large");
47058804f6dSBarry Smith #endif
471289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
472289bc588SBarry Smith         }
473289bc588SBarry Smith       }
474289bc588SBarry Smith     }
475e8d4e0b9SBarry Smith   }
476289bc588SBarry Smith   return 0;
477289bc588SBarry Smith }
478e8d4e0b9SBarry Smith 
4795615d1e5SSatish Balay #undef __FUNC__
4805615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense"
4818f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
482ae80bb75SLois Curfman McInnes {
483ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
484ae80bb75SLois Curfman McInnes   int          i, j;
485ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
486ae80bb75SLois Curfman McInnes 
487ae80bb75SLois Curfman McInnes   /* row-oriented output */
488ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
489ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
490ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
491ae80bb75SLois Curfman McInnes     }
492ae80bb75SLois Curfman McInnes   }
493ae80bb75SLois Curfman McInnes   return 0;
494ae80bb75SLois Curfman McInnes }
495ae80bb75SLois Curfman McInnes 
496289bc588SBarry Smith /* -----------------------------------------------------------------*/
497289bc588SBarry Smith 
49877c4ece6SBarry Smith #include "sys.h"
49988e32edaSLois Curfman McInnes 
5005615d1e5SSatish Balay #undef __FUNC__
5015615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense"
50219bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
50388e32edaSLois Curfman McInnes {
50488e32edaSLois Curfman McInnes   Mat_SeqDense *a;
50588e32edaSLois Curfman McInnes   Mat          B;
50688e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
50788e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
50888e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
50919bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
51088e32edaSLois Curfman McInnes 
51188e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
512e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"view must have one processor");
51390ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
51477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
515e3372554SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,0,"Not matrix object");
51688e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
51788e32edaSLois Curfman McInnes 
51890ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
51990ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
52090ace30eSBarry Smith     B = *A;
52190ace30eSBarry Smith     a = (Mat_SeqDense *) B->data;
52290ace30eSBarry Smith 
52390ace30eSBarry Smith     /* read in nonzero values */
52477c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
52590ace30eSBarry Smith 
5266d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5276d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
52890ace30eSBarry Smith     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
52990ace30eSBarry Smith   } else {
53088e32edaSLois Curfman McInnes     /* read row lengths */
53145d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
53277c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
53388e32edaSLois Curfman McInnes 
53488e32edaSLois Curfman McInnes     /* create our matrix */
535b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
53688e32edaSLois Curfman McInnes     B = *A;
53788e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
53888e32edaSLois Curfman McInnes     v = a->v;
53988e32edaSLois Curfman McInnes 
54088e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
54145d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
54277c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
54345d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
54477c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
54588e32edaSLois Curfman McInnes 
54688e32edaSLois Curfman McInnes     /* insert into matrix */
54788e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
54888e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
54988e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
55088e32edaSLois Curfman McInnes     }
5510452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
55288e32edaSLois Curfman McInnes 
5536d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5546d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
55590ace30eSBarry Smith   }
55688e32edaSLois Curfman McInnes   return 0;
55788e32edaSLois Curfman McInnes }
55888e32edaSLois Curfman McInnes 
559932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
56077c4ece6SBarry Smith #include "sys.h"
561932b0c3eSLois Curfman McInnes 
5625615d1e5SSatish Balay #undef __FUNC__
5635615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII"
564932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
565289bc588SBarry Smith {
566932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
567932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
568d636dbe3SBarry Smith   FILE         *fd;
569932b0c3eSLois Curfman McInnes   char         *outputname;
570932b0c3eSLois Curfman McInnes   Scalar       *v;
571932b0c3eSLois Curfman McInnes 
57290ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
573932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
57490ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
575639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5767ddc982cSLois Curfman McInnes     return 0;  /* do nothing for now */
577f72585f2SLois Curfman McInnes   }
57802cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
57980cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
58044cd7ae7SLois Curfman McInnes       v = a->v + i;
58180cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
58280cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
58380cd9d93SLois Curfman McInnes #if defined(PETSC_COMPLEX)
58480cd9d93SLois Curfman McInnes         if (real(*v) != 0.0 && imag(*v) != 0.0)
58580cd9d93SLois Curfman McInnes           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
58680cd9d93SLois Curfman McInnes         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
58780cd9d93SLois Curfman McInnes         v += a->m;
58880cd9d93SLois Curfman McInnes #else
58980cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
59080cd9d93SLois Curfman McInnes         v += a->m;
59180cd9d93SLois Curfman McInnes #endif
59280cd9d93SLois Curfman McInnes       }
59380cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
59480cd9d93SLois Curfman McInnes     }
59580cd9d93SLois Curfman McInnes   }
596f72585f2SLois Curfman McInnes   else {
59747989497SBarry Smith #if defined(PETSC_COMPLEX)
59847989497SBarry Smith     int allreal = 1;
59947989497SBarry Smith     /* determine if matrix has all real values */
60047989497SBarry Smith     v = a->v;
60147989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
60247989497SBarry Smith       if (imag(v[i])) { allreal = 0; break ;}
60347989497SBarry Smith     }
60447989497SBarry Smith #endif
605932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
606932b0c3eSLois Curfman McInnes       v = a->v + i;
607932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
608289bc588SBarry Smith #if defined(PETSC_COMPLEX)
60947989497SBarry Smith         if (allreal) {
61047989497SBarry Smith           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
61147989497SBarry Smith         } else {
612932b0c3eSLois Curfman McInnes           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
61347989497SBarry Smith         }
614289bc588SBarry Smith #else
615932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
616289bc588SBarry Smith #endif
617289bc588SBarry Smith       }
6188e37dc09SBarry Smith       fprintf(fd,"\n");
619289bc588SBarry Smith     }
620da3a660dSBarry Smith   }
621932b0c3eSLois Curfman McInnes   fflush(fd);
622289bc588SBarry Smith   return 0;
623289bc588SBarry Smith }
624289bc588SBarry Smith 
6255615d1e5SSatish Balay #undef __FUNC__
6265615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary"
627932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
628932b0c3eSLois Curfman McInnes {
629932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
630932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
63190ace30eSBarry Smith   int          format;
63290ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
633932b0c3eSLois Curfman McInnes 
63490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
63590ace30eSBarry Smith 
63690ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
63702cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
63890ace30eSBarry Smith     /* store the matrix as a dense matrix */
63990ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
64090ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
64190ace30eSBarry Smith     col_lens[1] = m;
64290ace30eSBarry Smith     col_lens[2] = n;
64390ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
64477c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
64590ace30eSBarry Smith     PetscFree(col_lens);
64690ace30eSBarry Smith 
64790ace30eSBarry Smith     /* write out matrix, by rows */
64845d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
64990ace30eSBarry Smith     v    = a->v;
65090ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
65190ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
65290ace30eSBarry Smith         vals[i + j*m] = *v++;
65390ace30eSBarry Smith       }
65490ace30eSBarry Smith     }
65577c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
65690ace30eSBarry Smith     PetscFree(vals);
65790ace30eSBarry Smith   } else {
6580452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
659932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
660932b0c3eSLois Curfman McInnes     col_lens[1] = m;
661932b0c3eSLois Curfman McInnes     col_lens[2] = n;
662932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
663932b0c3eSLois Curfman McInnes 
664932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
665932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
66677c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
667932b0c3eSLois Curfman McInnes 
668932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
669932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
670932b0c3eSLois Curfman McInnes     ict = 0;
671932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
672932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
673932b0c3eSLois Curfman McInnes     }
67477c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
6750452661fSBarry Smith     PetscFree(col_lens);
676932b0c3eSLois Curfman McInnes 
677932b0c3eSLois Curfman McInnes     /* store nonzero values */
67845d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
679932b0c3eSLois Curfman McInnes     ict = 0;
680932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
681932b0c3eSLois Curfman McInnes       v = a->v + i;
682932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
683932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
684932b0c3eSLois Curfman McInnes       }
685932b0c3eSLois Curfman McInnes     }
68677c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
6870452661fSBarry Smith     PetscFree(anonz);
68890ace30eSBarry Smith   }
689932b0c3eSLois Curfman McInnes   return 0;
690932b0c3eSLois Curfman McInnes }
691932b0c3eSLois Curfman McInnes 
6925615d1e5SSatish Balay #undef __FUNC__
6935615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense"
6948f6be9afSLois Curfman McInnes int MatView_SeqDense(PetscObject obj,Viewer viewer)
695932b0c3eSLois Curfman McInnes {
696932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
697932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
698bcd2baecSBarry Smith   ViewerType   vtype;
699bcd2baecSBarry Smith   int          ierr;
700932b0c3eSLois Curfman McInnes 
701bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
702bcd2baecSBarry Smith 
703bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
704932b0c3eSLois Curfman McInnes     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
705932b0c3eSLois Curfman McInnes   }
70619bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
707932b0c3eSLois Curfman McInnes     return MatView_SeqDense_ASCII(A,viewer);
708932b0c3eSLois Curfman McInnes   }
709bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
710932b0c3eSLois Curfman McInnes     return MatView_SeqDense_Binary(A,viewer);
711932b0c3eSLois Curfman McInnes   }
712932b0c3eSLois Curfman McInnes   return 0;
713932b0c3eSLois Curfman McInnes }
714289bc588SBarry Smith 
7155615d1e5SSatish Balay #undef __FUNC__
7165615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense"
7178f6be9afSLois Curfman McInnes int MatDestroy_SeqDense(PetscObject obj)
718289bc588SBarry Smith {
719289bc588SBarry Smith   Mat          mat = (Mat) obj;
720ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
72190f02eecSBarry Smith   int          ierr;
72290f02eecSBarry Smith 
723a5a9c739SBarry Smith #if defined(PETSC_LOG)
724a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
725a5a9c739SBarry Smith #endif
7260452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
7273439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
7280452661fSBarry Smith   PetscFree(l);
72990f02eecSBarry Smith   if (mat->mapping) {
73090f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
73190f02eecSBarry Smith   }
732a5a9c739SBarry Smith   PLogObjectDestroy(mat);
7330452661fSBarry Smith   PetscHeaderDestroy(mat);
734289bc588SBarry Smith   return 0;
735289bc588SBarry Smith }
736289bc588SBarry Smith 
7375615d1e5SSatish Balay #undef __FUNC__
7385615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense"
7398f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
740289bc588SBarry Smith {
741c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
742d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
743d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
74448b35521SBarry Smith 
745d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
7463638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
747e3372554SBarry Smith     if (m != n) SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place");
748d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
749289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
750d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
751d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
752d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
753289bc588SBarry Smith       }
754289bc588SBarry Smith     }
75548b35521SBarry Smith   }
756d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
757d3e5ee88SLois Curfman McInnes     int          ierr;
758d3e5ee88SLois Curfman McInnes     Mat          tmat;
759ec8511deSBarry Smith     Mat_SeqDense *tmatd;
760d3e5ee88SLois Curfman McInnes     Scalar       *v2;
761b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
762ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
7630de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
764d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
7650de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
766d3e5ee88SLois Curfman McInnes     }
7676d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
7686d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
769d3e5ee88SLois Curfman McInnes     *matout = tmat;
77048b35521SBarry Smith   }
771289bc588SBarry Smith   return 0;
772289bc588SBarry Smith }
773289bc588SBarry Smith 
7745615d1e5SSatish Balay #undef __FUNC__
7755615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense"
7768f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
777289bc588SBarry Smith {
778c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
779c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
780289bc588SBarry Smith   int          i;
781289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
7829ea5d5aeSSatish Balay 
783e3372554SBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Both matrices should be of type  MATSEQDENSE");
78477c4ece6SBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
78577c4ece6SBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
786289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
78777c4ece6SBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
788289bc588SBarry Smith     v1++; v2++;
789289bc588SBarry Smith   }
79077c4ece6SBarry Smith   *flg = PETSC_TRUE;
7919ea5d5aeSSatish Balay   return 0;
792289bc588SBarry Smith }
793289bc588SBarry Smith 
7945615d1e5SSatish Balay #undef __FUNC__
7955615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense"
7968f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
797289bc588SBarry Smith {
798c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
79944cd7ae7SLois Curfman McInnes   int          i, n, len;
80044cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
80144cd7ae7SLois Curfman McInnes 
80244cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
803289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
80444cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
805e3372554SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
80644cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
807289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
808289bc588SBarry Smith   }
809289bc588SBarry Smith   return 0;
810289bc588SBarry Smith }
811289bc588SBarry Smith 
8125615d1e5SSatish Balay #undef __FUNC__
8135615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense"
8148f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
815289bc588SBarry Smith {
816c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
817da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
818da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
81955659b69SBarry Smith 
82028988994SBarry Smith   if (ll) {
821da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
822e3372554SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
823da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
824da3a660dSBarry Smith       x = l[i];
825da3a660dSBarry Smith       v = mat->v + i;
826da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
827da3a660dSBarry Smith     }
82844cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
829da3a660dSBarry Smith   }
83028988994SBarry Smith   if (rr) {
831da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
832e3372554SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
833da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
834da3a660dSBarry Smith       x = r[i];
835da3a660dSBarry Smith       v = mat->v + i*m;
836da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
837da3a660dSBarry Smith     }
83844cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
839da3a660dSBarry Smith   }
840289bc588SBarry Smith   return 0;
841289bc588SBarry Smith }
842289bc588SBarry Smith 
8435615d1e5SSatish Balay #undef __FUNC__
8445615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense"
8458f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm)
846289bc588SBarry Smith {
847c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
848289bc588SBarry Smith   Scalar       *v = mat->v;
849289bc588SBarry Smith   double       sum = 0.0;
850289bc588SBarry Smith   int          i, j;
85155659b69SBarry Smith 
852289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
853289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
854289bc588SBarry Smith #if defined(PETSC_COMPLEX)
855289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
856289bc588SBarry Smith #else
857289bc588SBarry Smith       sum += (*v)*(*v); v++;
858289bc588SBarry Smith #endif
859289bc588SBarry Smith     }
860289bc588SBarry Smith     *norm = sqrt(sum);
86155659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
862289bc588SBarry Smith   }
863289bc588SBarry Smith   else if (type == NORM_1) {
864289bc588SBarry Smith     *norm = 0.0;
865289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
866289bc588SBarry Smith       sum = 0.0;
867289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
86833a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
869289bc588SBarry Smith       }
870289bc588SBarry Smith       if (sum > *norm) *norm = sum;
871289bc588SBarry Smith     }
87255659b69SBarry Smith     PLogFlops(mat->n*mat->m);
873289bc588SBarry Smith   }
874289bc588SBarry Smith   else if (type == NORM_INFINITY) {
875289bc588SBarry Smith     *norm = 0.0;
876289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
877289bc588SBarry Smith       v = mat->v + j;
878289bc588SBarry Smith       sum = 0.0;
879289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
880cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
881289bc588SBarry Smith       }
882289bc588SBarry Smith       if (sum > *norm) *norm = sum;
883289bc588SBarry Smith     }
88455659b69SBarry Smith     PLogFlops(mat->n*mat->m);
885289bc588SBarry Smith   }
886289bc588SBarry Smith   else {
887e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
888289bc588SBarry Smith   }
889289bc588SBarry Smith   return 0;
890289bc588SBarry Smith }
891289bc588SBarry Smith 
8925615d1e5SSatish Balay #undef __FUNC__
8935615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense"
8948f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
895289bc588SBarry Smith {
896c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
89767e560aaSBarry Smith 
8986d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
8996d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
9006d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
901219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
9026d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
903219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
9046d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
9056d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
9066d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
9076d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
908c2653b3dSLois Curfman McInnes            op == MAT_NEW_NONZERO_LOCATION_ERROR ||
9096d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
91090f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
911*2b362799SSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES)
91294a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
913c0bbcb79SLois Curfman McInnes   else
914e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
915289bc588SBarry Smith   return 0;
916289bc588SBarry Smith }
917289bc588SBarry Smith 
9185615d1e5SSatish Balay #undef __FUNC__
9195615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense"
9208f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
9216f0a148fSBarry Smith {
922ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
923cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
9246f0a148fSBarry Smith   return 0;
9256f0a148fSBarry Smith }
9266f0a148fSBarry Smith 
9275615d1e5SSatish Balay #undef __FUNC__
9285615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense"
9298f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
9304e220ebcSLois Curfman McInnes {
9314e220ebcSLois Curfman McInnes   *bs = 1;
9324e220ebcSLois Curfman McInnes   return 0;
9334e220ebcSLois Curfman McInnes }
9344e220ebcSLois Curfman McInnes 
9355615d1e5SSatish Balay #undef __FUNC__
9365615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense"
9378f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
9386f0a148fSBarry Smith {
939ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
9406abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
9416f0a148fSBarry Smith   Scalar       *slot;
94255659b69SBarry Smith 
94377c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
94478b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
9456f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
9466f0a148fSBarry Smith     slot = l->v + rows[i];
9476f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
9486f0a148fSBarry Smith   }
9496f0a148fSBarry Smith   if (diag) {
9506f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
9516f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
952260925b8SBarry Smith       *slot = *diag;
9536f0a148fSBarry Smith     }
9546f0a148fSBarry Smith   }
955260925b8SBarry Smith   ISRestoreIndices(is,&rows);
9566f0a148fSBarry Smith   return 0;
9576f0a148fSBarry Smith }
958557bce09SLois Curfman McInnes 
9595615d1e5SSatish Balay #undef __FUNC__
9605615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense"
9618f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n)
962557bce09SLois Curfman McInnes {
963c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
964557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
965557bce09SLois Curfman McInnes   return 0;
966557bce09SLois Curfman McInnes }
967557bce09SLois Curfman McInnes 
9685615d1e5SSatish Balay #undef __FUNC__
9695615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense"
9708f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
971d3e5ee88SLois Curfman McInnes {
972c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
973d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
974d3e5ee88SLois Curfman McInnes   return 0;
975d3e5ee88SLois Curfman McInnes }
976d3e5ee88SLois Curfman McInnes 
9775615d1e5SSatish Balay #undef __FUNC__
9785615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense"
9798f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array)
98064e87e97SBarry Smith {
981c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
98264e87e97SBarry Smith   *array = mat->v;
98364e87e97SBarry Smith   return 0;
98464e87e97SBarry Smith }
9850754003eSLois Curfman McInnes 
9865615d1e5SSatish Balay #undef __FUNC__
9875615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense"
9888f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array)
989ff14e315SSatish Balay {
990ff14e315SSatish Balay   return 0;
991ff14e315SSatish Balay }
9920754003eSLois Curfman McInnes 
9935615d1e5SSatish Balay #undef __FUNC__
9945615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense"
995cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
996cddf8d76SBarry Smith                                     Mat *submat)
9970754003eSLois Curfman McInnes {
998c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
9990754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
1000160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
10010754003eSLois Curfman McInnes   Scalar       *vwork, *val;
10020754003eSLois Curfman McInnes   Mat          newmat;
10030754003eSLois Curfman McInnes 
100478b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
100578b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
100678b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
100778b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
10080754003eSLois Curfman McInnes 
10090452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
10100452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
10110452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
1012cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
10130754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
10140754003eSLois Curfman McInnes 
10150754003eSLois Curfman McInnes   /* Create and fill new matrix */
1016b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
10170754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
10180754003eSLois Curfman McInnes     nznew = 0;
10190754003eSLois Curfman McInnes     val   = mat->v + irow[i];
10200754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
10210754003eSLois Curfman McInnes       if (smap[j]) {
10220754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
10230754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
10240754003eSLois Curfman McInnes       }
10250754003eSLois Curfman McInnes     }
1026dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
102778b31e54SBarry Smith            CHKERRQ(ierr);
10280754003eSLois Curfman McInnes   }
10296d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10306d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10310754003eSLois Curfman McInnes 
10320754003eSLois Curfman McInnes   /* Free work space */
10330452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
103478b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
103578b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10360754003eSLois Curfman McInnes   *submat = newmat;
10370754003eSLois Curfman McInnes   return 0;
10380754003eSLois Curfman McInnes }
10390754003eSLois Curfman McInnes 
10405615d1e5SSatish Balay #undef __FUNC__
10415615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense"
10428f6be9afSLois Curfman McInnes int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1043905e6a2fSBarry Smith                                     Mat **B)
1044905e6a2fSBarry Smith {
1045905e6a2fSBarry Smith   int ierr,i;
1046905e6a2fSBarry Smith 
1047905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1048905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1049905e6a2fSBarry Smith   }
1050905e6a2fSBarry Smith 
1051905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
1052905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1053905e6a2fSBarry Smith   }
1054905e6a2fSBarry Smith   return 0;
1055905e6a2fSBarry Smith }
1056905e6a2fSBarry Smith 
10575615d1e5SSatish Balay #undef __FUNC__
10585615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense"
10598f6be9afSLois Curfman McInnes int MatCopy_SeqDense(Mat A, Mat B)
10604b0e389bSBarry Smith {
10614b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
10624b0e389bSBarry Smith   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
1063e3372554SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
10644b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
10654b0e389bSBarry Smith   return 0;
10664b0e389bSBarry Smith }
10674b0e389bSBarry Smith 
1068289bc588SBarry Smith /* -------------------------------------------------------------------*/
1069ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
1070905e6a2fSBarry Smith        MatGetRow_SeqDense,
1071905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1072905e6a2fSBarry Smith        MatMult_SeqDense,
1073905e6a2fSBarry Smith        MatMultAdd_SeqDense,
1074905e6a2fSBarry Smith        MatMultTrans_SeqDense,
1075905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
1076905e6a2fSBarry Smith        MatSolve_SeqDense,
1077905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
1078905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
1079905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
1080905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1081905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1082ec8511deSBarry Smith        MatRelax_SeqDense,
1083ec8511deSBarry Smith        MatTranspose_SeqDense,
1084905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1085905e6a2fSBarry Smith        MatEqual_SeqDense,
1086905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1087905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1088905e6a2fSBarry Smith        MatNorm_SeqDense,
1089905e6a2fSBarry Smith        0,
1090905e6a2fSBarry Smith        0,
1091905e6a2fSBarry Smith        0,
1092905e6a2fSBarry Smith        MatSetOption_SeqDense,
1093905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1094905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1095905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1096905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1097905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1098905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1099905e6a2fSBarry Smith        MatGetSize_SeqDense,
1100905e6a2fSBarry Smith        MatGetSize_SeqDense,
1101905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1102905e6a2fSBarry Smith        0,
1103905e6a2fSBarry Smith        0,
1104905e6a2fSBarry Smith        MatGetArray_SeqDense,
1105905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
11063d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
1107905e6a2fSBarry Smith        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
11084b0e389bSBarry Smith        MatGetValues_SeqDense,
11094e220ebcSLois Curfman McInnes        MatCopy_SeqDense,0,MatScale_SeqDense,
11104e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_SeqDense};
111190ace30eSBarry Smith 
11125615d1e5SSatish Balay #undef __FUNC__
11135615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense"
11144b828684SBarry Smith /*@C
1115fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1116d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1117d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1118289bc588SBarry Smith 
111920563c6bSBarry Smith    Input Parameters:
11200c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
11210c775827SLois Curfman McInnes .  m - number of rows
112218f449edSLois Curfman McInnes .  n - number of columns
1123b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1124dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
112520563c6bSBarry Smith 
112620563c6bSBarry Smith    Output Parameter:
112744cd7ae7SLois Curfman McInnes .  A - the matrix
112820563c6bSBarry Smith 
112918f449edSLois Curfman McInnes   Notes:
113018f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
113118f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
1132b4fd4287SBarry Smith   set data=PETSC_NULL.
113318f449edSLois Curfman McInnes 
1134dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1135d65003e9SLois Curfman McInnes 
1136dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
113720563c6bSBarry Smith @*/
113844cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1139289bc588SBarry Smith {
114044cd7ae7SLois Curfman McInnes   Mat          B;
114144cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
11423b2fbd54SBarry Smith   int          ierr,flg,size;
11433b2fbd54SBarry Smith 
11443b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1145e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
114655659b69SBarry Smith 
114744cd7ae7SLois Curfman McInnes   *A            = 0;
114844cd7ae7SLois Curfman McInnes   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
114944cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
115044cd7ae7SLois Curfman McInnes   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
115144cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqDense));
115244cd7ae7SLois Curfman McInnes   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
115344cd7ae7SLois Curfman McInnes   B->destroy    = MatDestroy_SeqDense;
115444cd7ae7SLois Curfman McInnes   B->view       = MatView_SeqDense;
115544cd7ae7SLois Curfman McInnes   B->factor     = 0;
115690f02eecSBarry Smith   B->mapping    = 0;
115744cd7ae7SLois Curfman McInnes   PLogObjectMemory(B,sizeof(struct _Mat));
115844cd7ae7SLois Curfman McInnes   B->data       = (void *) b;
115918f449edSLois Curfman McInnes 
116044cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
116144cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
116244cd7ae7SLois Curfman McInnes   b->pivots       = 0;
116344cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1164b4fd4287SBarry Smith   if (data == PETSC_NULL) {
116544cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
116644cd7ae7SLois Curfman McInnes     PetscMemzero(b->v,m*n*sizeof(Scalar));
116744cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
116844cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
116918f449edSLois Curfman McInnes   }
11702dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
117144cd7ae7SLois Curfman McInnes     b->v = data;
117244cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
11732dd2b441SLois Curfman McInnes   }
117425cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
117544cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
11764e220ebcSLois Curfman McInnes 
117744cd7ae7SLois Curfman McInnes   *A = B;
1178289bc588SBarry Smith   return 0;
1179289bc588SBarry Smith }
1180289bc588SBarry Smith 
11813b2fbd54SBarry Smith 
11823b2fbd54SBarry Smith 
11833b2fbd54SBarry Smith 
11843b2fbd54SBarry Smith 
11853b2fbd54SBarry Smith 
1186