xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 5615d1e584023db9367fb782d85b1b4ebbb8df18)
1cb512458SBarry Smith #ifndef lint
2*5615d1e5SSatish Balay static char vcid[] = "$Id: dense.c,v 1.120 1997/01/01 03:37:31 bsmith 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 
12*5615d1e5SSatish Balay #undef __FUNC__
13*5615d1e5SSatish 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 
23*5615d1e5SSatish Balay #undef __FUNC__
24*5615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqDense"
254e220ebcSLois Curfman McInnes static 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 
50*5615d1e5SSatish Balay #undef __FUNC__
51*5615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqDense"
5280cd9d93SLois Curfman McInnes static 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.*/
66*5615d1e5SSatish Balay #undef __FUNC__
67*5615d1e5SSatish Balay #define __FUNC__ "MatLUFactor_SeqDense"
68c0bbcb79SLois Curfman McInnes static 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 
85*5615d1e5SSatish Balay #undef __FUNC__
86*5615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqDense"
8702cad45dSBarry Smith static 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 
103*5615d1e5SSatish Balay #undef __FUNC__
104*5615d1e5SSatish Balay #define __FUNC__ "MatLUFactorSymbolic_SeqDense"
105c0bbcb79SLois Curfman McInnes static 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 
110*5615d1e5SSatish Balay #undef __FUNC__
111*5615d1e5SSatish Balay #define __FUNC__ "MatLUFactorNumeric_SeqDense"
112c0bbcb79SLois Curfman McInnes static 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 
121*5615d1e5SSatish Balay #undef __FUNC__
122*5615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense"
123c0bbcb79SLois Curfman McInnes static 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 
128*5615d1e5SSatish Balay #undef __FUNC__
129*5615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense"
130c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
131289bc588SBarry Smith {
13249d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
133289bc588SBarry Smith }
1346ee01492SSatish Balay 
135*5615d1e5SSatish Balay #undef __FUNC__
136*5615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense"
137c0bbcb79SLois Curfman McInnes static 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 
154*5615d1e5SSatish Balay #undef __FUNC__
155*5615d1e5SSatish Balay #define __FUNC__ "MatSolve_SeqDense"
156c0bbcb79SLois Curfman McInnes static 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 
176*5615d1e5SSatish Balay #undef __FUNC__
177*5615d1e5SSatish Balay #define __FUNC__ "MatSolveTrans_SeqDense"
178c0bbcb79SLois Curfman McInnes static 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 
198*5615d1e5SSatish Balay #undef __FUNC__
199*5615d1e5SSatish Balay #define __FUNC__ "MatSolveAdd_SeqDense"
200c0bbcb79SLois Curfman McInnes static 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 
228*5615d1e5SSatish Balay #undef __FUNC__
229*5615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense"
230c0bbcb79SLois Curfman McInnes static 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 /* ------------------------------------------------------------------*/
263*5615d1e5SSatish Balay #undef __FUNC__
264*5615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense"
265c0bbcb79SLois Curfman McInnes static 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 /* -----------------------------------------------------------------*/
318*5615d1e5SSatish Balay #undef __FUNC__
319*5615d1e5SSatish 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 
331*5615d1e5SSatish Balay #undef __FUNC__
332*5615d1e5SSatish 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 
344*5615d1e5SSatish Balay #undef __FUNC__
345*5615d1e5SSatish 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 
358*5615d1e5SSatish Balay #undef __FUNC__
359*5615d1e5SSatish 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 /* -----------------------------------------------------------------*/
374*5615d1e5SSatish Balay #undef __FUNC__
375*5615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense"
376c0bbcb79SLois Curfman McInnes static 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 
395*5615d1e5SSatish Balay #undef __FUNC__
396*5615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense"
397c0bbcb79SLois Curfman McInnes static 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 /* ----------------------------------------------------------------*/
404*5615d1e5SSatish Balay #undef __FUNC__
405*5615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense"
406ae80bb75SLois Curfman McInnes static 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++ ) {
415289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
416289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
417289bc588SBarry Smith         }
418289bc588SBarry Smith       }
419289bc588SBarry Smith     }
420289bc588SBarry Smith     else {
421289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
422289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
423289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
424289bc588SBarry Smith         }
425289bc588SBarry Smith       }
426289bc588SBarry Smith     }
427e8d4e0b9SBarry Smith   }
428e8d4e0b9SBarry Smith   else {
429dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
430e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
431e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
432e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
433e8d4e0b9SBarry Smith         }
434e8d4e0b9SBarry Smith       }
435e8d4e0b9SBarry Smith     }
436289bc588SBarry Smith     else {
437289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
438289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
439289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
440289bc588SBarry Smith         }
441289bc588SBarry Smith       }
442289bc588SBarry Smith     }
443e8d4e0b9SBarry Smith   }
444289bc588SBarry Smith   return 0;
445289bc588SBarry Smith }
446e8d4e0b9SBarry Smith 
447*5615d1e5SSatish Balay #undef __FUNC__
448*5615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense"
449ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
450ae80bb75SLois Curfman McInnes {
451ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
452ae80bb75SLois Curfman McInnes   int          i, j;
453ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
454ae80bb75SLois Curfman McInnes 
455ae80bb75SLois Curfman McInnes   /* row-oriented output */
456ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
457ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
458ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
459ae80bb75SLois Curfman McInnes     }
460ae80bb75SLois Curfman McInnes   }
461ae80bb75SLois Curfman McInnes   return 0;
462ae80bb75SLois Curfman McInnes }
463ae80bb75SLois Curfman McInnes 
464289bc588SBarry Smith /* -----------------------------------------------------------------*/
465289bc588SBarry Smith 
46677c4ece6SBarry Smith #include "sys.h"
46788e32edaSLois Curfman McInnes 
468*5615d1e5SSatish Balay #undef __FUNC__
469*5615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense"
47019bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
47188e32edaSLois Curfman McInnes {
47288e32edaSLois Curfman McInnes   Mat_SeqDense *a;
47388e32edaSLois Curfman McInnes   Mat          B;
47488e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
47588e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
47688e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
47719bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
47888e32edaSLois Curfman McInnes 
47988e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
480e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"view must have one processor");
48190ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
48277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
483e3372554SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,0,"Not matrix object");
48488e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
48588e32edaSLois Curfman McInnes 
48690ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
48790ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
48890ace30eSBarry Smith     B = *A;
48990ace30eSBarry Smith     a = (Mat_SeqDense *) B->data;
49090ace30eSBarry Smith 
49190ace30eSBarry Smith     /* read in nonzero values */
49277c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
49390ace30eSBarry Smith 
4946d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4956d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
49690ace30eSBarry Smith     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
49790ace30eSBarry Smith   } else {
49888e32edaSLois Curfman McInnes     /* read row lengths */
49945d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
50077c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
50188e32edaSLois Curfman McInnes 
50288e32edaSLois Curfman McInnes     /* create our matrix */
503b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
50488e32edaSLois Curfman McInnes     B = *A;
50588e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
50688e32edaSLois Curfman McInnes     v = a->v;
50788e32edaSLois Curfman McInnes 
50888e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
50945d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
51077c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
51145d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
51277c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
51388e32edaSLois Curfman McInnes 
51488e32edaSLois Curfman McInnes     /* insert into matrix */
51588e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
51688e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
51788e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
51888e32edaSLois Curfman McInnes     }
5190452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
52088e32edaSLois Curfman McInnes 
5216d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5226d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
52390ace30eSBarry Smith   }
52488e32edaSLois Curfman McInnes   return 0;
52588e32edaSLois Curfman McInnes }
52688e32edaSLois Curfman McInnes 
527932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
52877c4ece6SBarry Smith #include "sys.h"
529932b0c3eSLois Curfman McInnes 
530*5615d1e5SSatish Balay #undef __FUNC__
531*5615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII"
532932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
533289bc588SBarry Smith {
534932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
535932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
536d636dbe3SBarry Smith   FILE         *fd;
537932b0c3eSLois Curfman McInnes   char         *outputname;
538932b0c3eSLois Curfman McInnes   Scalar       *v;
539932b0c3eSLois Curfman McInnes 
54090ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
541932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
54290ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
543639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5447ddc982cSLois Curfman McInnes     return 0;  /* do nothing for now */
545f72585f2SLois Curfman McInnes   }
54602cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
54780cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
54844cd7ae7SLois Curfman McInnes       v = a->v + i;
54980cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
55080cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
55180cd9d93SLois Curfman McInnes #if defined(PETSC_COMPLEX)
55280cd9d93SLois Curfman McInnes         if (real(*v) != 0.0 && imag(*v) != 0.0)
55380cd9d93SLois Curfman McInnes           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
55480cd9d93SLois Curfman McInnes         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
55580cd9d93SLois Curfman McInnes         v += a->m;
55680cd9d93SLois Curfman McInnes #else
55780cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
55880cd9d93SLois Curfman McInnes         v += a->m;
55980cd9d93SLois Curfman McInnes #endif
56080cd9d93SLois Curfman McInnes       }
56180cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
56280cd9d93SLois Curfman McInnes     }
56380cd9d93SLois Curfman McInnes   }
564f72585f2SLois Curfman McInnes   else {
56547989497SBarry Smith #if defined(PETSC_COMPLEX)
56647989497SBarry Smith     int allreal = 1;
56747989497SBarry Smith     /* determine if matrix has all real values */
56847989497SBarry Smith     v = a->v;
56947989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
57047989497SBarry Smith       if (imag(v[i])) { allreal = 0; break ;}
57147989497SBarry Smith     }
57247989497SBarry Smith #endif
573932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
574932b0c3eSLois Curfman McInnes       v = a->v + i;
575932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
576289bc588SBarry Smith #if defined(PETSC_COMPLEX)
57747989497SBarry Smith         if (allreal) {
57847989497SBarry Smith           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
57947989497SBarry Smith         } else {
580932b0c3eSLois Curfman McInnes           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
58147989497SBarry Smith         }
582289bc588SBarry Smith #else
583932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
584289bc588SBarry Smith #endif
585289bc588SBarry Smith       }
5868e37dc09SBarry Smith       fprintf(fd,"\n");
587289bc588SBarry Smith     }
588da3a660dSBarry Smith   }
589932b0c3eSLois Curfman McInnes   fflush(fd);
590289bc588SBarry Smith   return 0;
591289bc588SBarry Smith }
592289bc588SBarry Smith 
593*5615d1e5SSatish Balay #undef __FUNC__
594*5615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary"
595932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
596932b0c3eSLois Curfman McInnes {
597932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
598932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
59990ace30eSBarry Smith   int          format;
60090ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
601932b0c3eSLois Curfman McInnes 
60290ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
60390ace30eSBarry Smith 
60490ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
60502cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
60690ace30eSBarry Smith     /* store the matrix as a dense matrix */
60790ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
60890ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
60990ace30eSBarry Smith     col_lens[1] = m;
61090ace30eSBarry Smith     col_lens[2] = n;
61190ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
61277c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
61390ace30eSBarry Smith     PetscFree(col_lens);
61490ace30eSBarry Smith 
61590ace30eSBarry Smith     /* write out matrix, by rows */
61645d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
61790ace30eSBarry Smith     v    = a->v;
61890ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
61990ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
62090ace30eSBarry Smith         vals[i + j*m] = *v++;
62190ace30eSBarry Smith       }
62290ace30eSBarry Smith     }
62377c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
62490ace30eSBarry Smith     PetscFree(vals);
62590ace30eSBarry Smith   } else {
6260452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
627932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
628932b0c3eSLois Curfman McInnes     col_lens[1] = m;
629932b0c3eSLois Curfman McInnes     col_lens[2] = n;
630932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
631932b0c3eSLois Curfman McInnes 
632932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
633932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
63477c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
635932b0c3eSLois Curfman McInnes 
636932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
637932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
638932b0c3eSLois Curfman McInnes     ict = 0;
639932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
640932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
641932b0c3eSLois Curfman McInnes     }
64277c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
6430452661fSBarry Smith     PetscFree(col_lens);
644932b0c3eSLois Curfman McInnes 
645932b0c3eSLois Curfman McInnes     /* store nonzero values */
64645d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
647932b0c3eSLois Curfman McInnes     ict = 0;
648932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
649932b0c3eSLois Curfman McInnes       v = a->v + i;
650932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
651932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
652932b0c3eSLois Curfman McInnes       }
653932b0c3eSLois Curfman McInnes     }
65477c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
6550452661fSBarry Smith     PetscFree(anonz);
65690ace30eSBarry Smith   }
657932b0c3eSLois Curfman McInnes   return 0;
658932b0c3eSLois Curfman McInnes }
659932b0c3eSLois Curfman McInnes 
660*5615d1e5SSatish Balay #undef __FUNC__
661*5615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense"
662932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
663932b0c3eSLois Curfman McInnes {
664932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
665932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
666bcd2baecSBarry Smith   ViewerType   vtype;
667bcd2baecSBarry Smith   int          ierr;
668932b0c3eSLois Curfman McInnes 
669bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
670bcd2baecSBarry Smith 
671bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
672932b0c3eSLois Curfman McInnes     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
673932b0c3eSLois Curfman McInnes   }
67419bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
675932b0c3eSLois Curfman McInnes     return MatView_SeqDense_ASCII(A,viewer);
676932b0c3eSLois Curfman McInnes   }
677bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
678932b0c3eSLois Curfman McInnes     return MatView_SeqDense_Binary(A,viewer);
679932b0c3eSLois Curfman McInnes   }
680932b0c3eSLois Curfman McInnes   return 0;
681932b0c3eSLois Curfman McInnes }
682289bc588SBarry Smith 
683*5615d1e5SSatish Balay #undef __FUNC__
684*5615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense"
685ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
686289bc588SBarry Smith {
687289bc588SBarry Smith   Mat          mat = (Mat) obj;
688ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
68990f02eecSBarry Smith   int          ierr;
69090f02eecSBarry Smith 
691a5a9c739SBarry Smith #if defined(PETSC_LOG)
692a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
693a5a9c739SBarry Smith #endif
6940452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
6953439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
6960452661fSBarry Smith   PetscFree(l);
69790f02eecSBarry Smith   if (mat->mapping) {
69890f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
69990f02eecSBarry Smith   }
700a5a9c739SBarry Smith   PLogObjectDestroy(mat);
7010452661fSBarry Smith   PetscHeaderDestroy(mat);
702289bc588SBarry Smith   return 0;
703289bc588SBarry Smith }
704289bc588SBarry Smith 
705*5615d1e5SSatish Balay #undef __FUNC__
706*5615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense"
707c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
708289bc588SBarry Smith {
709c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
710d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
711d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
71248b35521SBarry Smith 
713d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
7143638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
715e3372554SBarry Smith     if (m != n) SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place");
716d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
717289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
718d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
719d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
720d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
721289bc588SBarry Smith       }
722289bc588SBarry Smith     }
72348b35521SBarry Smith   }
724d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
725d3e5ee88SLois Curfman McInnes     int          ierr;
726d3e5ee88SLois Curfman McInnes     Mat          tmat;
727ec8511deSBarry Smith     Mat_SeqDense *tmatd;
728d3e5ee88SLois Curfman McInnes     Scalar       *v2;
729b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
730ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
7310de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
732d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
7330de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
734d3e5ee88SLois Curfman McInnes     }
7356d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
7366d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
737d3e5ee88SLois Curfman McInnes     *matout = tmat;
73848b35521SBarry Smith   }
739289bc588SBarry Smith   return 0;
740289bc588SBarry Smith }
741289bc588SBarry Smith 
742*5615d1e5SSatish Balay #undef __FUNC__
743*5615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense"
74477c4ece6SBarry Smith static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
745289bc588SBarry Smith {
746c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
747c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
748289bc588SBarry Smith   int          i;
749289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
7509ea5d5aeSSatish Balay 
751e3372554SBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Both matrices should be of type  MATSEQDENSE");
75277c4ece6SBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
75377c4ece6SBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
754289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
75577c4ece6SBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
756289bc588SBarry Smith     v1++; v2++;
757289bc588SBarry Smith   }
75877c4ece6SBarry Smith   *flg = PETSC_TRUE;
7599ea5d5aeSSatish Balay   return 0;
760289bc588SBarry Smith }
761289bc588SBarry Smith 
762*5615d1e5SSatish Balay #undef __FUNC__
763*5615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense"
764c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
765289bc588SBarry Smith {
766c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
76744cd7ae7SLois Curfman McInnes   int          i, n, len;
76844cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
76944cd7ae7SLois Curfman McInnes 
77044cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
771289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
77244cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
773e3372554SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
77444cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
775289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
776289bc588SBarry Smith   }
777289bc588SBarry Smith   return 0;
778289bc588SBarry Smith }
779289bc588SBarry Smith 
780*5615d1e5SSatish Balay #undef __FUNC__
781*5615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense"
782052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
783289bc588SBarry Smith {
784c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
785da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
786da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
78755659b69SBarry Smith 
78828988994SBarry Smith   if (ll) {
789da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
790e3372554SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
791da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
792da3a660dSBarry Smith       x = l[i];
793da3a660dSBarry Smith       v = mat->v + i;
794da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
795da3a660dSBarry Smith     }
79644cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
797da3a660dSBarry Smith   }
79828988994SBarry Smith   if (rr) {
799da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
800e3372554SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
801da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
802da3a660dSBarry Smith       x = r[i];
803da3a660dSBarry Smith       v = mat->v + i*m;
804da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
805da3a660dSBarry Smith     }
80644cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
807da3a660dSBarry Smith   }
808289bc588SBarry Smith   return 0;
809289bc588SBarry Smith }
810289bc588SBarry Smith 
811*5615d1e5SSatish Balay #undef __FUNC__
812*5615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense"
813cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
814289bc588SBarry Smith {
815c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
816289bc588SBarry Smith   Scalar       *v = mat->v;
817289bc588SBarry Smith   double       sum = 0.0;
818289bc588SBarry Smith   int          i, j;
81955659b69SBarry Smith 
820289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
821289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
822289bc588SBarry Smith #if defined(PETSC_COMPLEX)
823289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
824289bc588SBarry Smith #else
825289bc588SBarry Smith       sum += (*v)*(*v); v++;
826289bc588SBarry Smith #endif
827289bc588SBarry Smith     }
828289bc588SBarry Smith     *norm = sqrt(sum);
82955659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
830289bc588SBarry Smith   }
831289bc588SBarry Smith   else if (type == NORM_1) {
832289bc588SBarry Smith     *norm = 0.0;
833289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
834289bc588SBarry Smith       sum = 0.0;
835289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
83633a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
837289bc588SBarry Smith       }
838289bc588SBarry Smith       if (sum > *norm) *norm = sum;
839289bc588SBarry Smith     }
84055659b69SBarry Smith     PLogFlops(mat->n*mat->m);
841289bc588SBarry Smith   }
842289bc588SBarry Smith   else if (type == NORM_INFINITY) {
843289bc588SBarry Smith     *norm = 0.0;
844289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
845289bc588SBarry Smith       v = mat->v + j;
846289bc588SBarry Smith       sum = 0.0;
847289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
848cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
849289bc588SBarry Smith       }
850289bc588SBarry Smith       if (sum > *norm) *norm = sum;
851289bc588SBarry Smith     }
85255659b69SBarry Smith     PLogFlops(mat->n*mat->m);
853289bc588SBarry Smith   }
854289bc588SBarry Smith   else {
855e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
856289bc588SBarry Smith   }
857289bc588SBarry Smith   return 0;
858289bc588SBarry Smith }
859289bc588SBarry Smith 
860*5615d1e5SSatish Balay #undef __FUNC__
861*5615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense"
862c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
863289bc588SBarry Smith {
864c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
86567e560aaSBarry Smith 
8666d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
8676d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
8686d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
869219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
8706d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
871219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
8726d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
8736d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
8746d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
8756d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
8766d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
87790f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
87890f02eecSBarry Smith            op == MAT_IGNORE_OFF_PROCESSOR_ENTRIES)
87994a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
880c0bbcb79SLois Curfman McInnes   else
881e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
882289bc588SBarry Smith   return 0;
883289bc588SBarry Smith }
884289bc588SBarry Smith 
885*5615d1e5SSatish Balay #undef __FUNC__
886*5615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense"
887ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
8886f0a148fSBarry Smith {
889ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
890cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
8916f0a148fSBarry Smith   return 0;
8926f0a148fSBarry Smith }
8936f0a148fSBarry Smith 
894*5615d1e5SSatish Balay #undef __FUNC__
895*5615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense"
8964e220ebcSLois Curfman McInnes static int MatGetBlockSize_SeqDense(Mat A,int *bs)
8974e220ebcSLois Curfman McInnes {
8984e220ebcSLois Curfman McInnes   *bs = 1;
8994e220ebcSLois Curfman McInnes   return 0;
9004e220ebcSLois Curfman McInnes }
9014e220ebcSLois Curfman McInnes 
902*5615d1e5SSatish Balay #undef __FUNC__
903*5615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense"
904ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
9056f0a148fSBarry Smith {
906ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
9076abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
9086f0a148fSBarry Smith   Scalar       *slot;
90955659b69SBarry Smith 
91077c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
91178b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
9126f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
9136f0a148fSBarry Smith     slot = l->v + rows[i];
9146f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
9156f0a148fSBarry Smith   }
9166f0a148fSBarry Smith   if (diag) {
9176f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
9186f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
919260925b8SBarry Smith       *slot = *diag;
9206f0a148fSBarry Smith     }
9216f0a148fSBarry Smith   }
922260925b8SBarry Smith   ISRestoreIndices(is,&rows);
9236f0a148fSBarry Smith   return 0;
9246f0a148fSBarry Smith }
925557bce09SLois Curfman McInnes 
926*5615d1e5SSatish Balay #undef __FUNC__
927*5615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense"
928c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
929557bce09SLois Curfman McInnes {
930c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
931557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
932557bce09SLois Curfman McInnes   return 0;
933557bce09SLois Curfman McInnes }
934557bce09SLois Curfman McInnes 
935*5615d1e5SSatish Balay #undef __FUNC__
936*5615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense"
937c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
938d3e5ee88SLois Curfman McInnes {
939c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
940d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
941d3e5ee88SLois Curfman McInnes   return 0;
942d3e5ee88SLois Curfman McInnes }
943d3e5ee88SLois Curfman McInnes 
944*5615d1e5SSatish Balay #undef __FUNC__
945*5615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense"
946c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
94764e87e97SBarry Smith {
948c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
94964e87e97SBarry Smith   *array = mat->v;
95064e87e97SBarry Smith   return 0;
95164e87e97SBarry Smith }
9520754003eSLois Curfman McInnes 
953*5615d1e5SSatish Balay #undef __FUNC__
954*5615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense"
955ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
956ff14e315SSatish Balay {
957ff14e315SSatish Balay   return 0;
958ff14e315SSatish Balay }
9590754003eSLois Curfman McInnes 
960*5615d1e5SSatish Balay #undef __FUNC__
961*5615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense"
962cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
963cddf8d76SBarry Smith                                     Mat *submat)
9640754003eSLois Curfman McInnes {
965c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
9660754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
967160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
9680754003eSLois Curfman McInnes   Scalar       *vwork, *val;
9690754003eSLois Curfman McInnes   Mat          newmat;
9700754003eSLois Curfman McInnes 
97178b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
97278b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
97378b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
97478b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
9750754003eSLois Curfman McInnes 
9760452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
9770452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
9780452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
979cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
9800754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
9810754003eSLois Curfman McInnes 
9820754003eSLois Curfman McInnes   /* Create and fill new matrix */
983b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
9840754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
9850754003eSLois Curfman McInnes     nznew = 0;
9860754003eSLois Curfman McInnes     val   = mat->v + irow[i];
9870754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
9880754003eSLois Curfman McInnes       if (smap[j]) {
9890754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
9900754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
9910754003eSLois Curfman McInnes       }
9920754003eSLois Curfman McInnes     }
993dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
99478b31e54SBarry Smith            CHKERRQ(ierr);
9950754003eSLois Curfman McInnes   }
9966d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9976d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9980754003eSLois Curfman McInnes 
9990754003eSLois Curfman McInnes   /* Free work space */
10000452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
100178b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
100278b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10030754003eSLois Curfman McInnes   *submat = newmat;
10040754003eSLois Curfman McInnes   return 0;
10050754003eSLois Curfman McInnes }
10060754003eSLois Curfman McInnes 
1007*5615d1e5SSatish Balay #undef __FUNC__
1008*5615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense"
1009905e6a2fSBarry Smith static int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1010905e6a2fSBarry Smith                                     Mat **B)
1011905e6a2fSBarry Smith {
1012905e6a2fSBarry Smith   int ierr,i;
1013905e6a2fSBarry Smith 
1014905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1015905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1016905e6a2fSBarry Smith   }
1017905e6a2fSBarry Smith 
1018905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
1019905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1020905e6a2fSBarry Smith   }
1021905e6a2fSBarry Smith   return 0;
1022905e6a2fSBarry Smith }
1023905e6a2fSBarry Smith 
1024*5615d1e5SSatish Balay #undef __FUNC__
1025*5615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense"
10264b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B)
10274b0e389bSBarry Smith {
10284b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
10294b0e389bSBarry Smith   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
1030e3372554SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
10314b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
10324b0e389bSBarry Smith   return 0;
10334b0e389bSBarry Smith }
10344b0e389bSBarry Smith 
1035289bc588SBarry Smith /* -------------------------------------------------------------------*/
1036ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
1037905e6a2fSBarry Smith        MatGetRow_SeqDense,
1038905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1039905e6a2fSBarry Smith        MatMult_SeqDense,
1040905e6a2fSBarry Smith        MatMultAdd_SeqDense,
1041905e6a2fSBarry Smith        MatMultTrans_SeqDense,
1042905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
1043905e6a2fSBarry Smith        MatSolve_SeqDense,
1044905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
1045905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
1046905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
1047905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1048905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1049ec8511deSBarry Smith        MatRelax_SeqDense,
1050ec8511deSBarry Smith        MatTranspose_SeqDense,
1051905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1052905e6a2fSBarry Smith        MatEqual_SeqDense,
1053905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1054905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1055905e6a2fSBarry Smith        MatNorm_SeqDense,
1056905e6a2fSBarry Smith        0,
1057905e6a2fSBarry Smith        0,
1058905e6a2fSBarry Smith        0,
1059905e6a2fSBarry Smith        MatSetOption_SeqDense,
1060905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1061905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1062905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1063905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1064905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1065905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1066905e6a2fSBarry Smith        MatGetSize_SeqDense,
1067905e6a2fSBarry Smith        MatGetSize_SeqDense,
1068905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1069905e6a2fSBarry Smith        0,
1070905e6a2fSBarry Smith        0,
1071905e6a2fSBarry Smith        MatGetArray_SeqDense,
1072905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1073905e6a2fSBarry Smith        0,
10743d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
1075905e6a2fSBarry Smith        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
10764b0e389bSBarry Smith        MatGetValues_SeqDense,
10774e220ebcSLois Curfman McInnes        MatCopy_SeqDense,0,MatScale_SeqDense,
10784e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_SeqDense};
107990ace30eSBarry Smith 
1080*5615d1e5SSatish Balay #undef __FUNC__
1081*5615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense"
10824b828684SBarry Smith /*@C
1083fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1084d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1085d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1086289bc588SBarry Smith 
108720563c6bSBarry Smith    Input Parameters:
10880c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
10890c775827SLois Curfman McInnes .  m - number of rows
109018f449edSLois Curfman McInnes .  n - number of columns
1091b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1092dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
109320563c6bSBarry Smith 
109420563c6bSBarry Smith    Output Parameter:
109544cd7ae7SLois Curfman McInnes .  A - the matrix
109620563c6bSBarry Smith 
109718f449edSLois Curfman McInnes   Notes:
109818f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
109918f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
1100b4fd4287SBarry Smith   set data=PETSC_NULL.
110118f449edSLois Curfman McInnes 
1102dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1103d65003e9SLois Curfman McInnes 
1104dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
110520563c6bSBarry Smith @*/
110644cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1107289bc588SBarry Smith {
110844cd7ae7SLois Curfman McInnes   Mat          B;
110944cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
11103b2fbd54SBarry Smith   int          ierr,flg,size;
11113b2fbd54SBarry Smith 
11123b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1113e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
111455659b69SBarry Smith 
111544cd7ae7SLois Curfman McInnes   *A            = 0;
111644cd7ae7SLois Curfman McInnes   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
111744cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
111844cd7ae7SLois Curfman McInnes   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
111944cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqDense));
112044cd7ae7SLois Curfman McInnes   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
112144cd7ae7SLois Curfman McInnes   B->destroy    = MatDestroy_SeqDense;
112244cd7ae7SLois Curfman McInnes   B->view       = MatView_SeqDense;
112344cd7ae7SLois Curfman McInnes   B->factor     = 0;
112490f02eecSBarry Smith   B->mapping    = 0;
112544cd7ae7SLois Curfman McInnes   PLogObjectMemory(B,sizeof(struct _Mat));
112644cd7ae7SLois Curfman McInnes   B->data       = (void *) b;
112718f449edSLois Curfman McInnes 
112844cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
112944cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
113044cd7ae7SLois Curfman McInnes   b->pivots       = 0;
113144cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1132b4fd4287SBarry Smith   if (data == PETSC_NULL) {
113344cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
113444cd7ae7SLois Curfman McInnes     PetscMemzero(b->v,m*n*sizeof(Scalar));
113544cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
113644cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
113718f449edSLois Curfman McInnes   }
11382dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
113944cd7ae7SLois Curfman McInnes     b->v = data;
114044cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
11412dd2b441SLois Curfman McInnes   }
114225cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
114344cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
11444e220ebcSLois Curfman McInnes 
114544cd7ae7SLois Curfman McInnes   *A = B;
1146289bc588SBarry Smith   return 0;
1147289bc588SBarry Smith }
1148289bc588SBarry Smith 
1149*5615d1e5SSatish Balay #undef __FUNC__
1150*5615d1e5SSatish Balay #define __FUNC__ "MatCreate_SeqDense"
1151c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
1152289bc588SBarry Smith {
1153c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
1154b4fd4287SBarry Smith   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
1155289bc588SBarry Smith }
11563b2fbd54SBarry Smith 
11573b2fbd54SBarry Smith 
11583b2fbd54SBarry Smith 
11593b2fbd54SBarry Smith 
11603b2fbd54SBarry Smith 
1161