xref: /petsc/src/mat/impls/dense/seq/dense.c (revision fb9695e52cb52bdd61c020aa3052f0ef9b41684c)
1*fb9695e5SSatish Balay /*$Id: dense.c,v 1.192 2001/01/15 21:45:30 bsmith Exp balay $*/
267e560aaSBarry Smith /*
367e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
467e560aaSBarry Smith */
5289bc588SBarry Smith 
670f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
7b4c862e0SSatish Balay #include "petscblaslapack.h"
8289bc588SBarry Smith 
95615d1e5SSatish Balay #undef __FUNC__
10b0a32e0cSBarry Smith #define __FUNC__ "MatAXPY_SeqDense"
111987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
121987afe7SBarry Smith {
131987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
14273d9f13SBarry Smith   int          N = X->m*X->n,one = 1;
153a40ed3dSBarry Smith 
163a40ed3dSBarry Smith   PetscFunctionBegin;
171987afe7SBarry Smith   BLaxpy_(&N,alpha,x->v,&one,y->v,&one);
18b0a32e0cSBarry Smith   PetscLogFlops(2*N-1);
193a40ed3dSBarry Smith   PetscFunctionReturn(0);
201987afe7SBarry Smith }
211987afe7SBarry Smith 
225615d1e5SSatish Balay #undef __FUNC__
23b0a32e0cSBarry Smith #define __FUNC__ "MatGetInfo_SeqDense"
248f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
25289bc588SBarry Smith {
264e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
27273d9f13SBarry Smith   int          i,N = A->m*A->n,count = 0;
284e220ebcSLois Curfman McInnes   Scalar       *v = a->v;
293a40ed3dSBarry Smith 
303a40ed3dSBarry Smith   PetscFunctionBegin;
31289bc588SBarry Smith   for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}
324e220ebcSLois Curfman McInnes 
33273d9f13SBarry Smith   info->rows_global       = (double)A->m;
34273d9f13SBarry Smith   info->columns_global    = (double)A->n;
35273d9f13SBarry Smith   info->rows_local        = (double)A->m;
36273d9f13SBarry Smith   info->columns_local     = (double)A->n;
374e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
384e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
394e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
404e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
414e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
424e220ebcSLois Curfman McInnes   info->mallocs           = 0;
434e220ebcSLois Curfman McInnes   info->memory            = A->mem;
444e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
454e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
464e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
474e220ebcSLois Curfman McInnes 
483a40ed3dSBarry Smith   PetscFunctionReturn(0);
49289bc588SBarry Smith }
50289bc588SBarry Smith 
515615d1e5SSatish Balay #undef __FUNC__
52b0a32e0cSBarry Smith #define __FUNC__ "MatScale_SeqDense"
53273d9f13SBarry Smith int MatScale_SeqDense(Scalar *alpha,Mat A)
5480cd9d93SLois Curfman McInnes {
55273d9f13SBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
5680cd9d93SLois Curfman McInnes   int          one = 1,nz;
5780cd9d93SLois Curfman McInnes 
583a40ed3dSBarry Smith   PetscFunctionBegin;
59273d9f13SBarry Smith   nz = A->m*A->n;
6080cd9d93SLois Curfman McInnes   BLscal_(&nz,alpha,a->v,&one);
61b0a32e0cSBarry Smith   PetscLogFlops(nz);
623a40ed3dSBarry Smith   PetscFunctionReturn(0);
6380cd9d93SLois Curfman McInnes }
6480cd9d93SLois Curfman McInnes 
65289bc588SBarry Smith /* ---------------------------------------------------------------*/
666831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
676831982aSBarry Smith    rather than put it in the Mat->row slot.*/
685615d1e5SSatish Balay #undef __FUNC__
69b0a32e0cSBarry Smith #define __FUNC__ "MatLUFactor_SeqDense"
70e03a110bSBarry Smith int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatLUInfo *minfo)
71289bc588SBarry Smith {
72c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
73b0a32e0cSBarry Smith   int          info,ierr;
7455659b69SBarry Smith 
753a40ed3dSBarry Smith   PetscFunctionBegin;
76289bc588SBarry Smith   if (!mat->pivots) {
77b0a32e0cSBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr);
78b0a32e0cSBarry Smith     PetscLogObjectMemory(A,A->m*sizeof(int));
79289bc588SBarry Smith   }
80f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
81273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
82273d9f13SBarry Smith   LAgetrf_(&A->m,&A->n,mat->v,&A->m,mat->pivots,&info);
8329bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
8429bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
85b0a32e0cSBarry Smith   PetscLogFlops((2*A->n*A->n*A->n)/3);
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
87289bc588SBarry Smith }
886ee01492SSatish Balay 
895615d1e5SSatish Balay #undef __FUNC__
90b0a32e0cSBarry Smith #define __FUNC__ "MatDuplicate_SeqDense"
915609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
9202cad45dSBarry Smith {
9302cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l;
9402cad45dSBarry Smith   int          ierr;
9502cad45dSBarry Smith   Mat          newi;
9602cad45dSBarry Smith 
973a40ed3dSBarry Smith   PetscFunctionBegin;
98273d9f13SBarry Smith   ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr);
9902cad45dSBarry Smith   l = (Mat_SeqDense*)newi->data;
1005609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
101273d9f13SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr);
10202cad45dSBarry Smith   }
10302cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
10402cad45dSBarry Smith   *newmat = newi;
1053a40ed3dSBarry Smith   PetscFunctionReturn(0);
10602cad45dSBarry Smith }
10702cad45dSBarry Smith 
1085615d1e5SSatish Balay #undef __FUNC__
109b0a32e0cSBarry Smith #define __FUNC__ "MatLUFactorSymbolic_SeqDense"
110e03a110bSBarry Smith int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatLUInfo *info,Mat *fact)
111289bc588SBarry Smith {
1123a40ed3dSBarry Smith   int ierr;
1133a40ed3dSBarry Smith 
1143a40ed3dSBarry Smith   PetscFunctionBegin;
1155609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1163a40ed3dSBarry Smith   PetscFunctionReturn(0);
117289bc588SBarry Smith }
1186ee01492SSatish Balay 
1195615d1e5SSatish Balay #undef __FUNC__
120b0a32e0cSBarry Smith #define __FUNC__ "MatLUFactorNumeric_SeqDense"
1218f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
122289bc588SBarry Smith {
12302cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1243a40ed3dSBarry Smith   int          ierr;
1253a40ed3dSBarry Smith 
1263a40ed3dSBarry Smith   PetscFunctionBegin;
12702cad45dSBarry Smith   /* copy the numerical values */
128273d9f13SBarry Smith   ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr);
12902cad45dSBarry Smith   (*fact)->factor = 0;
130e03a110bSBarry Smith   ierr = MatLUFactor(*fact,0,0,PETSC_NULL);CHKERRQ(ierr);
1313a40ed3dSBarry Smith   PetscFunctionReturn(0);
132289bc588SBarry Smith }
1336ee01492SSatish Balay 
1345615d1e5SSatish Balay #undef __FUNC__
135b0a32e0cSBarry Smith #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense"
136329f5518SBarry Smith int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,PetscReal f,Mat *fact)
137289bc588SBarry Smith {
1383a40ed3dSBarry Smith   int ierr;
1393a40ed3dSBarry Smith 
1403a40ed3dSBarry Smith   PetscFunctionBegin;
1413a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1423a40ed3dSBarry Smith   PetscFunctionReturn(0);
143289bc588SBarry Smith }
1446ee01492SSatish Balay 
1455615d1e5SSatish Balay #undef __FUNC__
146b0a32e0cSBarry Smith #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense"
1478f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
148289bc588SBarry Smith {
1493a40ed3dSBarry Smith   int ierr;
1503a40ed3dSBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
1523a40ed3dSBarry Smith   ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr);
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154289bc588SBarry Smith }
1556ee01492SSatish Balay 
1565615d1e5SSatish Balay #undef __FUNC__
1575615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense"
158329f5518SBarry Smith int MatCholeskyFactor_SeqDense(Mat A,IS perm,PetscReal f)
159289bc588SBarry Smith {
160c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
161606d414cSSatish Balay   int           info,ierr;
16255659b69SBarry Smith 
1633a40ed3dSBarry Smith   PetscFunctionBegin;
164752f5567SLois Curfman McInnes   if (mat->pivots) {
165606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
166b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-A->m*sizeof(int));
167752f5567SLois Curfman McInnes     mat->pivots = 0;
168752f5567SLois Curfman McInnes   }
169273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
170273d9f13SBarry Smith   LApotrf_("L",&A->n,mat->v,&A->m,&info);
17129bbc08cSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
172c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
173b0a32e0cSBarry Smith   PetscLogFlops((A->n*A->n*A->n)/3);
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175289bc588SBarry Smith }
176289bc588SBarry Smith 
1775615d1e5SSatish Balay #undef __FUNC__
178b0a32e0cSBarry Smith #define __FUNC__ "MatSolve_SeqDense"
1798f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
180289bc588SBarry Smith {
181c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
182a57ad546SLois Curfman McInnes   int          one = 1,info,ierr;
1836abc6512SBarry Smith   Scalar       *x,*y;
18467e560aaSBarry Smith 
1853a40ed3dSBarry Smith   PetscFunctionBegin;
186273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
1877a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1887a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
189273d9f13SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr);
190c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
191273d9f13SBarry Smith     LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
1927a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
193273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
194289bc588SBarry Smith   }
19529bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
19629bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
1977a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1987a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
199b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2003a40ed3dSBarry Smith   PetscFunctionReturn(0);
201289bc588SBarry Smith }
2026ee01492SSatish Balay 
2035615d1e5SSatish Balay #undef __FUNC__
204b0a32e0cSBarry Smith #define __FUNC__ "MatSolveTranspose_SeqDense"
2057c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
206da3a660dSBarry Smith {
207c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2087a97a34bSBarry Smith   int          ierr,one = 1,info;
2096abc6512SBarry Smith   Scalar       *x,*y;
21067e560aaSBarry Smith 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
212273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2137a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2147a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
215273d9f13SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr);
216752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
217da3a660dSBarry Smith   if (mat->pivots) {
218273d9f13SBarry Smith     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
2197a97a34bSBarry Smith   } else {
220273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
221da3a660dSBarry Smith   }
22229bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
2237a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2247a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
225b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2263a40ed3dSBarry Smith   PetscFunctionReturn(0);
227da3a660dSBarry Smith }
2286ee01492SSatish Balay 
2295615d1e5SSatish Balay #undef __FUNC__
230b0a32e0cSBarry Smith #define __FUNC__ "MatSolveAdd_SeqDense"
2318f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
232da3a660dSBarry Smith {
233c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2347a97a34bSBarry Smith   int          ierr,one = 1,info;
2356abc6512SBarry Smith   Scalar       *x,*y,sone = 1.0;
236da3a660dSBarry Smith   Vec          tmp = 0;
23767e560aaSBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
2397a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2407a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
241273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
242da3a660dSBarry Smith   if (yy == zz) {
24378b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
244b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
24578b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
246da3a660dSBarry Smith   }
247273d9f13SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr);
248752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
249da3a660dSBarry Smith   if (mat->pivots) {
250273d9f13SBarry Smith     LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
251a8c6a408SBarry Smith   } else {
252273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
253da3a660dSBarry Smith   }
25429bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
255a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
256a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
2577a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2587a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
259b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
261da3a660dSBarry Smith }
26267e560aaSBarry Smith 
2635615d1e5SSatish Balay #undef __FUNC__
264b0a32e0cSBarry Smith #define __FUNC__ "MatSolveTransposeAdd_SeqDense"
2657c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
266da3a660dSBarry Smith {
267c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
2686abc6512SBarry Smith   int           one = 1,info,ierr;
2696abc6512SBarry Smith   Scalar        *x,*y,sone = 1.0;
270da3a660dSBarry Smith   Vec           tmp;
27167e560aaSBarry Smith 
2723a40ed3dSBarry Smith   PetscFunctionBegin;
273273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2747a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2757a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
276da3a660dSBarry Smith   if (yy == zz) {
27778b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
278b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
27978b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
280da3a660dSBarry Smith   }
281273d9f13SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr);
282752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
283da3a660dSBarry Smith   if (mat->pivots) {
284273d9f13SBarry Smith     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
2853a40ed3dSBarry Smith   } else {
286273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
287da3a660dSBarry Smith   }
28829bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
28990f02eecSBarry Smith   if (tmp) {
29090f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
29190f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
2923a40ed3dSBarry Smith   } else {
29390f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
29490f02eecSBarry Smith   }
2957a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2967a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
297b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
2983a40ed3dSBarry Smith   PetscFunctionReturn(0);
299da3a660dSBarry Smith }
300289bc588SBarry Smith /* ------------------------------------------------------------------*/
3015615d1e5SSatish Balay #undef __FUNC__
302b0a32e0cSBarry Smith #define __FUNC__ "MatRelax_SeqDense"
303329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
304329f5518SBarry Smith                           PetscReal shift,int its,Vec xx)
305289bc588SBarry Smith {
306c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
307289bc588SBarry Smith   Scalar       *x,*b,*v = mat->v,zero = 0.0,xt;
308273d9f13SBarry Smith   int          ierr,m = A->m,i;
309aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
310bc1b551cSSatish Balay   int          o = 1;
311bc1b551cSSatish Balay #endif
312289bc588SBarry Smith 
3133a40ed3dSBarry Smith   PetscFunctionBegin;
314289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
3153a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
316bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
317289bc588SBarry Smith   }
3187a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3197a97a34bSBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
320289bc588SBarry Smith   while (its--) {
321289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
322289bc588SBarry Smith       for (i=0; i<m; i++) {
323aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
324f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
325f1747703SBarry Smith            not happy about returning a double complex */
326f1747703SBarry Smith         int    _i;
327f1747703SBarry Smith         Scalar sum = b[i];
328f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
3293f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
330f1747703SBarry Smith         }
331f1747703SBarry Smith         xt = sum;
332f1747703SBarry Smith #else
333289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
334f1747703SBarry Smith #endif
33555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
336289bc588SBarry Smith       }
337289bc588SBarry Smith     }
338289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
339289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
340aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
341f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
342f1747703SBarry Smith            not happy about returning a double complex */
343f1747703SBarry Smith         int    _i;
344f1747703SBarry Smith         Scalar sum = b[i];
345f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
3463f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
347f1747703SBarry Smith         }
348f1747703SBarry Smith         xt = sum;
349f1747703SBarry Smith #else
350289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
351f1747703SBarry Smith #endif
35255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
353289bc588SBarry Smith       }
354289bc588SBarry Smith     }
355289bc588SBarry Smith   }
356600d6d0bSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
3577a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3583a40ed3dSBarry Smith   PetscFunctionReturn(0);
359289bc588SBarry Smith }
360289bc588SBarry Smith 
361289bc588SBarry Smith /* -----------------------------------------------------------------*/
3625615d1e5SSatish Balay #undef __FUNC__
363b0a32e0cSBarry Smith #define __FUNC__ "MatMultTranspose_SeqDense"
3647c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
365289bc588SBarry Smith {
366c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
367289bc588SBarry Smith   Scalar       *v = mat->v,*x,*y;
3687a97a34bSBarry Smith   int          ierr,_One=1;Scalar _DOne=1.0,_DZero=0.0;
3693a40ed3dSBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
371273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3727a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3737a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
374273d9f13SBarry Smith   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
3757a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3767a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
377b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->n);
3783a40ed3dSBarry Smith   PetscFunctionReturn(0);
379289bc588SBarry Smith }
3806ee01492SSatish Balay 
3815615d1e5SSatish Balay #undef __FUNC__
382b0a32e0cSBarry Smith #define __FUNC__ "MatMult_SeqDense"
38344cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
384289bc588SBarry Smith {
385c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
386329f5518SBarry Smith   Scalar       *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
387329f5518SBarry Smith   int          ierr,_One=1;
3883a40ed3dSBarry Smith 
3893a40ed3dSBarry Smith   PetscFunctionBegin;
390273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3917a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3927a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
393273d9f13SBarry Smith   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
3947a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3957a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
396b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->m);
3973a40ed3dSBarry Smith   PetscFunctionReturn(0);
398289bc588SBarry Smith }
3996ee01492SSatish Balay 
4005615d1e5SSatish Balay #undef __FUNC__
401b0a32e0cSBarry Smith #define __FUNC__ "MatMultAdd_SeqDense"
40244cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
403289bc588SBarry Smith {
404c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
405329f5518SBarry Smith   Scalar       *v = mat->v,*x,*y,_DOne=1.0;
406329f5518SBarry Smith   int          ierr,_One=1;
4073a40ed3dSBarry Smith 
4083a40ed3dSBarry Smith   PetscFunctionBegin;
409273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
410600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4117a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4127a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
413273d9f13SBarry Smith   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
4147a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4157a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
416b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4173a40ed3dSBarry Smith   PetscFunctionReturn(0);
418289bc588SBarry Smith }
4196ee01492SSatish Balay 
4205615d1e5SSatish Balay #undef __FUNC__
421b0a32e0cSBarry Smith #define __FUNC__ "MatMultTransposeAdd_SeqDense"
4227c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
423289bc588SBarry Smith {
424c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
425600d6d0bSBarry Smith   Scalar       *v = mat->v,*x,*y;
4267a97a34bSBarry Smith   int          ierr,_One=1;
4276abc6512SBarry Smith   Scalar       _DOne=1.0;
4283a40ed3dSBarry Smith 
4293a40ed3dSBarry Smith   PetscFunctionBegin;
430273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
431600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4327a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4337a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
434273d9f13SBarry Smith   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
4357a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4367a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
437b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4383a40ed3dSBarry Smith   PetscFunctionReturn(0);
439289bc588SBarry Smith }
440289bc588SBarry Smith 
441289bc588SBarry Smith /* -----------------------------------------------------------------*/
4425615d1e5SSatish Balay #undef __FUNC__
443b0a32e0cSBarry Smith #define __FUNC__ "MatGetRow_SeqDense"
4448f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
445289bc588SBarry Smith {
446c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
447289bc588SBarry Smith   Scalar       *v;
448b0a32e0cSBarry Smith   int          i,ierr;
44967e560aaSBarry Smith 
4503a40ed3dSBarry Smith   PetscFunctionBegin;
451273d9f13SBarry Smith   *ncols = A->n;
452289bc588SBarry Smith   if (cols) {
453b0a32e0cSBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr);
454273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
455289bc588SBarry Smith   }
456289bc588SBarry Smith   if (vals) {
457b0a32e0cSBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(Scalar),vals);CHKERRQ(ierr);
458289bc588SBarry Smith     v    = mat->v + row;
459273d9f13SBarry Smith     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->m;}
460289bc588SBarry Smith   }
4613a40ed3dSBarry Smith   PetscFunctionReturn(0);
462289bc588SBarry Smith }
4636ee01492SSatish Balay 
4645615d1e5SSatish Balay #undef __FUNC__
465b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreRow_SeqDense"
4668f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
467289bc588SBarry Smith {
468606d414cSSatish Balay   int ierr;
469606d414cSSatish Balay   PetscFunctionBegin;
470606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
471606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
4723a40ed3dSBarry Smith   PetscFunctionReturn(0);
473289bc588SBarry Smith }
474289bc588SBarry Smith /* ----------------------------------------------------------------*/
4755615d1e5SSatish Balay #undef __FUNC__
476b0a32e0cSBarry Smith #define __FUNC__ "MatSetValues_SeqDense"
4778f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
478e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
479289bc588SBarry Smith {
480c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
481289bc588SBarry Smith   int          i,j;
482d6dfbf8fSBarry Smith 
4833a40ed3dSBarry Smith   PetscFunctionBegin;
484289bc588SBarry Smith   if (!mat->roworiented) {
485dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
486289bc588SBarry Smith       for (j=0; j<n; j++) {
4875ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
488aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
48929bbc08cSBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
49058804f6dSBarry Smith #endif
491289bc588SBarry Smith         for (i=0; i<m; i++) {
4925ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
493aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
49429bbc08cSBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
49558804f6dSBarry Smith #endif
496273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
497289bc588SBarry Smith         }
498289bc588SBarry Smith       }
4993a40ed3dSBarry Smith     } else {
500289bc588SBarry Smith       for (j=0; j<n; j++) {
5015ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
502aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
50329bbc08cSBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
50458804f6dSBarry Smith #endif
505289bc588SBarry Smith         for (i=0; i<m; i++) {
5065ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
507aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
50829bbc08cSBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
50958804f6dSBarry Smith #endif
510273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
511289bc588SBarry Smith         }
512289bc588SBarry Smith       }
513289bc588SBarry Smith     }
5143a40ed3dSBarry Smith   } else {
515dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
516e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
5175ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
518aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
51929bbc08cSBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
52058804f6dSBarry Smith #endif
521e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
5225ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
523aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
52429bbc08cSBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
52558804f6dSBarry Smith #endif
526273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
527e8d4e0b9SBarry Smith         }
528e8d4e0b9SBarry Smith       }
5293a40ed3dSBarry Smith     } else {
530289bc588SBarry Smith       for (i=0; i<m; i++) {
5315ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
532aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
53329bbc08cSBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
53458804f6dSBarry Smith #endif
535289bc588SBarry Smith         for (j=0; j<n; j++) {
5365ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
537aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
53829bbc08cSBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
53958804f6dSBarry Smith #endif
540273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
541289bc588SBarry Smith         }
542289bc588SBarry Smith       }
543289bc588SBarry Smith     }
544e8d4e0b9SBarry Smith   }
5453a40ed3dSBarry Smith   PetscFunctionReturn(0);
546289bc588SBarry Smith }
547e8d4e0b9SBarry Smith 
5485615d1e5SSatish Balay #undef __FUNC__
549b0a32e0cSBarry Smith #define __FUNC__ "MatGetValues_SeqDense"
5508f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
551ae80bb75SLois Curfman McInnes {
552ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
553ae80bb75SLois Curfman McInnes   int          i,j;
554ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
555ae80bb75SLois Curfman McInnes 
5563a40ed3dSBarry Smith   PetscFunctionBegin;
557ae80bb75SLois Curfman McInnes   /* row-oriented output */
558ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
559ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
560273d9f13SBarry Smith       *vpt++ = mat->v[indexn[j]*A->m + indexm[i]];
561ae80bb75SLois Curfman McInnes     }
562ae80bb75SLois Curfman McInnes   }
5633a40ed3dSBarry Smith   PetscFunctionReturn(0);
564ae80bb75SLois Curfman McInnes }
565ae80bb75SLois Curfman McInnes 
566289bc588SBarry Smith /* -----------------------------------------------------------------*/
567289bc588SBarry Smith 
568e090d566SSatish Balay #include "petscsys.h"
56988e32edaSLois Curfman McInnes 
570273d9f13SBarry Smith EXTERN_C_BEGIN
5715615d1e5SSatish Balay #undef __FUNC__
572b0a32e0cSBarry Smith #define __FUNC__ "MatLoad_SeqDense"
573b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A)
57488e32edaSLois Curfman McInnes {
57588e32edaSLois Curfman McInnes   Mat_SeqDense *a;
57688e32edaSLois Curfman McInnes   Mat          B;
57788e32edaSLois Curfman McInnes   int          *scols,i,j,nz,ierr,fd,header[4],size;
57888e32edaSLois Curfman McInnes   int          *rowlengths = 0,M,N,*cols;
5794a41a4d3SSatish Balay   Scalar       *vals,*svals,*v,*w;
58019bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
58188e32edaSLois Curfman McInnes 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
583d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
58429bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
585b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
5860752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
58729bbc08cSBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
58888e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
58988e32edaSLois Curfman McInnes 
59090ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
59190ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
59290ace30eSBarry Smith     B    = *A;
59390ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
5944a41a4d3SSatish Balay     v    = a->v;
5954a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
5964a41a4d3SSatish Balay        from row major to column major */
597b0a32e0cSBarry Smith     ierr = PetscMalloc((M*N+1)*sizeof(Scalar),&w);CHKERRQ(ierr);
59890ace30eSBarry Smith     /* read in nonzero values */
5994a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6004a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6014a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6024a41a4d3SSatish Balay       for (i=0; i<M; i++) {
6034a41a4d3SSatish Balay         *v++ =w[i*N+j];
6044a41a4d3SSatish Balay       }
6054a41a4d3SSatish Balay     }
606606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
6076d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6086d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60990ace30eSBarry Smith   } else {
61088e32edaSLois Curfman McInnes     /* read row lengths */
611b0a32e0cSBarry Smith     ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr);
6120752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
61388e32edaSLois Curfman McInnes 
61488e32edaSLois Curfman McInnes     /* create our matrix */
615b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
61688e32edaSLois Curfman McInnes     B = *A;
61788e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
61888e32edaSLois Curfman McInnes     v = a->v;
61988e32edaSLois Curfman McInnes 
62088e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
621b0a32e0cSBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr);
622b0a32e0cSBarry Smith     cols = scols;
6230752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
624b0a32e0cSBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(Scalar),&svals);CHKERRQ(ierr);
625b0a32e0cSBarry Smith     vals = svals;
6260752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
62788e32edaSLois Curfman McInnes 
62888e32edaSLois Curfman McInnes     /* insert into matrix */
62988e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
63088e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
63188e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
63288e32edaSLois Curfman McInnes     }
633606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
634606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
635606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
63688e32edaSLois Curfman McInnes 
6376d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6386d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63990ace30eSBarry Smith   }
6403a40ed3dSBarry Smith   PetscFunctionReturn(0);
64188e32edaSLois Curfman McInnes }
642273d9f13SBarry Smith EXTERN_C_END
64388e32edaSLois Curfman McInnes 
644e090d566SSatish Balay #include "petscsys.h"
645932b0c3eSLois Curfman McInnes 
6465615d1e5SSatish Balay #undef __FUNC__
647b0a32e0cSBarry Smith #define __FUNC__ "MatView_SeqDense_ASCII"
648b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
649289bc588SBarry Smith {
650932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
651*fb9695e5SSatish Balay   int          ierr,i,j;
652*fb9695e5SSatish Balay   char         *name;
653932b0c3eSLois Curfman McInnes   Scalar       *v;
654*fb9695e5SSatish Balay   PetscViewerFormatType  format;
655932b0c3eSLois Curfman McInnes 
6563a40ed3dSBarry Smith   PetscFunctionBegin;
657*fb9695e5SSatish Balay   ierr = PetscObjectGetName((PetscObject)viewer,&name);CHKERRQ(ierr);
658b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
659*fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
6603a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
661*fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
662b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
663273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
66444cd7ae7SLois Curfman McInnes       v = a->v + i;
665b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr);
666273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
667aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
668329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
669b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
670329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
671b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %d %g ",j,PetscRealPart(*v));CHKERRQ(ierr);
6726831982aSBarry Smith         }
67380cd9d93SLois Curfman McInnes #else
6746831982aSBarry Smith         if (*v) {
675b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %d %g ",j,*v);CHKERRQ(ierr);
6766831982aSBarry Smith         }
67780cd9d93SLois Curfman McInnes #endif
678273d9f13SBarry Smith         v += A->m;
67980cd9d93SLois Curfman McInnes       }
680b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
68180cd9d93SLois Curfman McInnes     }
682b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
6833a40ed3dSBarry Smith   } else {
684b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
685aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
686ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
68747989497SBarry Smith     /* determine if matrix has all real values */
68847989497SBarry Smith     v = a->v;
689273d9f13SBarry Smith     for (i=0; i<A->m*A->n; i++) {
690ffac6cdbSBarry Smith       if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
69147989497SBarry Smith     }
69247989497SBarry Smith #endif
693*fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
694*fb9695e5SSatish Balay       ierr = PetscObjectGetName((PetscObject)viewer,&name);CHKERRQ(ierr);
695b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr);
696*fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr);
697*fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
698ffac6cdbSBarry Smith     }
699ffac6cdbSBarry Smith 
700273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
701932b0c3eSLois Curfman McInnes       v = a->v + i;
702273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
703aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
70447989497SBarry Smith         if (allreal) {
705b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); v += A->m;
70647989497SBarry Smith         } else {
707b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); v += A->m;
70847989497SBarry Smith         }
709289bc588SBarry Smith #else
710b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += A->m;
711289bc588SBarry Smith #endif
712289bc588SBarry Smith       }
713b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
714289bc588SBarry Smith     }
715*fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
716b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
717ffac6cdbSBarry Smith     }
718b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
719da3a660dSBarry Smith   }
720b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7213a40ed3dSBarry Smith   PetscFunctionReturn(0);
722289bc588SBarry Smith }
723289bc588SBarry Smith 
7245615d1e5SSatish Balay #undef __FUNC__
725b0a32e0cSBarry Smith #define __FUNC__ "MatView_SeqDense_Binary"
726b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
727932b0c3eSLois Curfman McInnes {
728932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
729273d9f13SBarry Smith   int          ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
73090ace30eSBarry Smith   Scalar       *v,*anonz,*vals;
731*fb9695e5SSatish Balay   PetscViewerFormatType  format;
732932b0c3eSLois Curfman McInnes 
7333a40ed3dSBarry Smith   PetscFunctionBegin;
734b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
73590ace30eSBarry Smith 
736b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
737*fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
73890ace30eSBarry Smith     /* store the matrix as a dense matrix */
739b0a32e0cSBarry Smith     ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr);
74090ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
74190ace30eSBarry Smith     col_lens[1] = m;
74290ace30eSBarry Smith     col_lens[2] = n;
74390ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
7440752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr);
745606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
74690ace30eSBarry Smith 
74790ace30eSBarry Smith     /* write out matrix, by rows */
748b0a32e0cSBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(Scalar),&vals);CHKERRQ(ierr);
74990ace30eSBarry Smith     v    = a->v;
75090ace30eSBarry Smith     for (i=0; i<m; i++) {
75190ace30eSBarry Smith       for (j=0; j<n; j++) {
75290ace30eSBarry Smith         vals[i + j*m] = *v++;
75390ace30eSBarry Smith       }
75490ace30eSBarry Smith     }
7550752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr);
756606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
75790ace30eSBarry Smith   } else {
758b0a32e0cSBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr);
759932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
760932b0c3eSLois Curfman McInnes     col_lens[1] = m;
761932b0c3eSLois Curfman McInnes     col_lens[2] = n;
762932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
763932b0c3eSLois Curfman McInnes 
764932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
765932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
7660752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr);
767932b0c3eSLois Curfman McInnes 
768932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
769932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
770932b0c3eSLois Curfman McInnes     ict = 0;
771932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
772932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
773932b0c3eSLois Curfman McInnes     }
7740752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr);
775606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
776932b0c3eSLois Curfman McInnes 
777932b0c3eSLois Curfman McInnes     /* store nonzero values */
778b0a32e0cSBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(Scalar),&anonz);CHKERRQ(ierr);
779932b0c3eSLois Curfman McInnes     ict  = 0;
780932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
781932b0c3eSLois Curfman McInnes       v = a->v + i;
782932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
783273d9f13SBarry Smith         anonz[ict++] = *v; v += A->m;
784932b0c3eSLois Curfman McInnes       }
785932b0c3eSLois Curfman McInnes     }
7860752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr);
787606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
78890ace30eSBarry Smith   }
7893a40ed3dSBarry Smith   PetscFunctionReturn(0);
790932b0c3eSLois Curfman McInnes }
791932b0c3eSLois Curfman McInnes 
7925615d1e5SSatish Balay #undef __FUNC__
793b0a32e0cSBarry Smith #define __FUNC__ "MatView_SeqDense_Draw_Zoom"
794b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
795f1af5d2fSBarry Smith {
796f1af5d2fSBarry Smith   Mat           A = (Mat) Aa;
797f1af5d2fSBarry Smith   Mat_SeqDense  *a = (Mat_SeqDense*)A->data;
798*fb9695e5SSatish Balay   int           m = A->m,n = A->n,color,i,j,ierr;
799f1af5d2fSBarry Smith   Scalar        *v = a->v;
800b0a32e0cSBarry Smith   PetscViewer        viewer;
801b0a32e0cSBarry Smith   PetscDraw          popup;
802329f5518SBarry Smith   PetscReal          xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
803*fb9695e5SSatish Balay   PetscViewerFormatType  format;
804f1af5d2fSBarry Smith 
805f1af5d2fSBarry Smith   PetscFunctionBegin;
806f1af5d2fSBarry Smith 
807f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
808b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
809b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
810f1af5d2fSBarry Smith 
811f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
812*fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
813f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
814b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
815f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
816f1af5d2fSBarry Smith       x_l = j;
817f1af5d2fSBarry Smith       x_r = x_l + 1.0;
818f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
819f1af5d2fSBarry Smith         y_l = m - i - 1.0;
820f1af5d2fSBarry Smith         y_r = y_l + 1.0;
821f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
822329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
823b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
824329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
825b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
826f1af5d2fSBarry Smith         } else {
827f1af5d2fSBarry Smith           continue;
828f1af5d2fSBarry Smith         }
829f1af5d2fSBarry Smith #else
830f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
831b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
832f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
833b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
834f1af5d2fSBarry Smith         } else {
835f1af5d2fSBarry Smith           continue;
836f1af5d2fSBarry Smith         }
837f1af5d2fSBarry Smith #endif
838b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
839f1af5d2fSBarry Smith       }
840f1af5d2fSBarry Smith     }
841f1af5d2fSBarry Smith   } else {
842f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
843f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
844f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
845f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
846f1af5d2fSBarry Smith     }
847b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
848b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
849b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
850f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
851f1af5d2fSBarry Smith       x_l = j;
852f1af5d2fSBarry Smith       x_r = x_l + 1.0;
853f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
854f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
855f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
856b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
857b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
858f1af5d2fSBarry Smith       }
859f1af5d2fSBarry Smith     }
860f1af5d2fSBarry Smith   }
861f1af5d2fSBarry Smith   PetscFunctionReturn(0);
862f1af5d2fSBarry Smith }
863f1af5d2fSBarry Smith 
864f1af5d2fSBarry Smith #undef __FUNC__
865b0a32e0cSBarry Smith #define __FUNC__ "MatView_SeqDense_Draw"
866b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
867f1af5d2fSBarry Smith {
868b0a32e0cSBarry Smith   PetscDraw       draw;
869f1af5d2fSBarry Smith   PetscTruth isnull;
870329f5518SBarry Smith   PetscReal  xr,yr,xl,yl,h,w;
871f1af5d2fSBarry Smith   int        ierr;
872f1af5d2fSBarry Smith 
873f1af5d2fSBarry Smith   PetscFunctionBegin;
874b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
875b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
876f1af5d2fSBarry Smith   if (isnull == PETSC_TRUE) PetscFunctionReturn(0);
877f1af5d2fSBarry Smith 
878f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
879273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
880f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
881b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
882b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
883f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
884f1af5d2fSBarry Smith   PetscFunctionReturn(0);
885f1af5d2fSBarry Smith }
886f1af5d2fSBarry Smith 
887f1af5d2fSBarry Smith #undef __FUNC__
888b0a32e0cSBarry Smith #define __FUNC__ "MatView_SeqDense"
889b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer)
890932b0c3eSLois Curfman McInnes {
891932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
892bcd2baecSBarry Smith   int          ierr;
893f1af5d2fSBarry Smith   PetscTruth   issocket,isascii,isbinary,isdraw;
894932b0c3eSLois Curfman McInnes 
8953a40ed3dSBarry Smith   PetscFunctionBegin;
896b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
897b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
898*fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
899*fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
9000f5bd95cSBarry Smith 
9010f5bd95cSBarry Smith   if (issocket) {
902b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
9030f5bd95cSBarry Smith   } else if (isascii) {
9043a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
9050f5bd95cSBarry Smith   } else if (isbinary) {
9063a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
907f1af5d2fSBarry Smith   } else if (isdraw) {
908f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
9095cd90555SBarry Smith   } else {
91029bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
911932b0c3eSLois Curfman McInnes   }
9123a40ed3dSBarry Smith   PetscFunctionReturn(0);
913932b0c3eSLois Curfman McInnes }
914289bc588SBarry Smith 
9155615d1e5SSatish Balay #undef __FUNC__
916b0a32e0cSBarry Smith #define __FUNC__ "MatDestroy_SeqDense"
917c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
918289bc588SBarry Smith {
919ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
92090f02eecSBarry Smith   int          ierr;
92190f02eecSBarry Smith 
9223a40ed3dSBarry Smith   PetscFunctionBegin;
923aa482453SBarry Smith #if defined(PETSC_USE_LOG)
924b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
925a5a9c739SBarry Smith #endif
926606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
927606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
928606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
9293a40ed3dSBarry Smith   PetscFunctionReturn(0);
930289bc588SBarry Smith }
931289bc588SBarry Smith 
9325615d1e5SSatish Balay #undef __FUNC__
933b0a32e0cSBarry Smith #define __FUNC__ "MatTranspose_SeqDense"
9348f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
935289bc588SBarry Smith {
936c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
937549d3d68SSatish Balay   int          k,j,m,n,ierr;
938d3e5ee88SLois Curfman McInnes   Scalar       *v,tmp;
93948b35521SBarry Smith 
9403a40ed3dSBarry Smith   PetscFunctionBegin;
941273d9f13SBarry Smith   v = mat->v; m = A->m; n = A->n;
9427c922b88SBarry Smith   if (!matout) { /* in place transpose */
943d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
944b0a32e0cSBarry Smith       Scalar *w;
945b0a32e0cSBarry Smith       ierr = PetscMalloc((m+1)*(n+1)*sizeof(Scalar),&w);CHKERRQ(ierr);
946d64ed03dSBarry Smith       for (j=0; j<m; j++) {
947d64ed03dSBarry Smith         for (k=0; k<n; k++) {
948d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
949d64ed03dSBarry Smith         }
950d64ed03dSBarry Smith       }
951549d3d68SSatish Balay       ierr = PetscMemcpy(v,w,m*n*sizeof(Scalar));CHKERRQ(ierr);
952606d414cSSatish Balay       ierr = PetscFree(w);CHKERRQ(ierr);
953d64ed03dSBarry Smith     } else {
954d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
955289bc588SBarry Smith         for (k=0; k<j; k++) {
956d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
957d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
958d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
959289bc588SBarry Smith         }
960289bc588SBarry Smith       }
961d64ed03dSBarry Smith     }
9623a40ed3dSBarry Smith   } else { /* out-of-place transpose */
963d3e5ee88SLois Curfman McInnes     Mat          tmat;
964ec8511deSBarry Smith     Mat_SeqDense *tmatd;
965d3e5ee88SLois Curfman McInnes     Scalar       *v2;
966273d9f13SBarry Smith     ierr  = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr);
967ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
9680de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
969d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
9700de55854SLois Curfman McInnes       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m];
971d3e5ee88SLois Curfman McInnes     }
9726d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9736d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
974d3e5ee88SLois Curfman McInnes     *matout = tmat;
97548b35521SBarry Smith   }
9763a40ed3dSBarry Smith   PetscFunctionReturn(0);
977289bc588SBarry Smith }
978289bc588SBarry Smith 
9795615d1e5SSatish Balay #undef __FUNC__
980b0a32e0cSBarry Smith #define __FUNC__ "MatEqual_SeqDense"
9818f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
982289bc588SBarry Smith {
983c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
984c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
985273d9f13SBarry Smith   int          ierr,i;
986289bc588SBarry Smith   Scalar       *v1 = mat1->v,*v2 = mat2->v;
987273d9f13SBarry Smith   PetscTruth   flag;
9889ea5d5aeSSatish Balay 
9893a40ed3dSBarry Smith   PetscFunctionBegin;
990273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr);
991273d9f13SBarry Smith   if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type  MATSEQDENSE");
992273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
993273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
994273d9f13SBarry Smith   for (i=0; i<A1->m*A1->n; i++) {
9953a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
996289bc588SBarry Smith     v1++; v2++;
997289bc588SBarry Smith   }
99877c4ece6SBarry Smith   *flg = PETSC_TRUE;
9993a40ed3dSBarry Smith   PetscFunctionReturn(0);
1000289bc588SBarry Smith }
1001289bc588SBarry Smith 
10025615d1e5SSatish Balay #undef __FUNC__
1003b0a32e0cSBarry Smith #define __FUNC__ "MatGetDiagonal_SeqDense"
10048f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
1005289bc588SBarry Smith {
1006c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
10077a97a34bSBarry Smith   int          ierr,i,n,len;
100844cd7ae7SLois Curfman McInnes   Scalar       *x,zero = 0.0;
100944cd7ae7SLois Curfman McInnes 
10103a40ed3dSBarry Smith   PetscFunctionBegin;
10117a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
10127a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
1013600d6d0bSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1014273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1015273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
101644cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
1017273d9f13SBarry Smith     x[i] = mat->v[i*A->m + i];
1018289bc588SBarry Smith   }
10197a97a34bSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
10203a40ed3dSBarry Smith   PetscFunctionReturn(0);
1021289bc588SBarry Smith }
1022289bc588SBarry Smith 
10235615d1e5SSatish Balay #undef __FUNC__
1024b0a32e0cSBarry Smith #define __FUNC__ "MatDiagonalScale_SeqDense"
10258f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1026289bc588SBarry Smith {
1027c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1028da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
1029273d9f13SBarry Smith   int          ierr,i,j,m = A->m,n = A->n;
103055659b69SBarry Smith 
10313a40ed3dSBarry Smith   PetscFunctionBegin;
103228988994SBarry Smith   if (ll) {
10337a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1034600d6d0bSBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1035273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1036da3a660dSBarry Smith     for (i=0; i<m; i++) {
1037da3a660dSBarry Smith       x = l[i];
1038da3a660dSBarry Smith       v = mat->v + i;
1039da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1040da3a660dSBarry Smith     }
10417a97a34bSBarry Smith     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1042b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1043da3a660dSBarry Smith   }
104428988994SBarry Smith   if (rr) {
10457a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1046600d6d0bSBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1047273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1048da3a660dSBarry Smith     for (i=0; i<n; i++) {
1049da3a660dSBarry Smith       x = r[i];
1050da3a660dSBarry Smith       v = mat->v + i*m;
1051da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1052da3a660dSBarry Smith     }
10537a97a34bSBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1054b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1055da3a660dSBarry Smith   }
10563a40ed3dSBarry Smith   PetscFunctionReturn(0);
1057289bc588SBarry Smith }
1058289bc588SBarry Smith 
10595615d1e5SSatish Balay #undef __FUNC__
1060b0a32e0cSBarry Smith #define __FUNC__ "MatNorm_SeqDense"
1061329f5518SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *norm)
1062289bc588SBarry Smith {
1063c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1064289bc588SBarry Smith   Scalar       *v = mat->v;
1065329f5518SBarry Smith   PetscReal    sum = 0.0;
1066289bc588SBarry Smith   int          i,j;
106755659b69SBarry Smith 
10683a40ed3dSBarry Smith   PetscFunctionBegin;
1069289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1070273d9f13SBarry Smith     for (i=0; i<A->n*A->m; i++) {
1071aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1072329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1073289bc588SBarry Smith #else
1074289bc588SBarry Smith       sum += (*v)*(*v); v++;
1075289bc588SBarry Smith #endif
1076289bc588SBarry Smith     }
1077289bc588SBarry Smith     *norm = sqrt(sum);
1078b0a32e0cSBarry Smith     PetscLogFlops(2*A->n*A->m);
10793a40ed3dSBarry Smith   } else if (type == NORM_1) {
1080289bc588SBarry Smith     *norm = 0.0;
1081273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
1082289bc588SBarry Smith       sum = 0.0;
1083273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
108433a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1085289bc588SBarry Smith       }
1086289bc588SBarry Smith       if (sum > *norm) *norm = sum;
1087289bc588SBarry Smith     }
1088b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
10893a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1090289bc588SBarry Smith     *norm = 0.0;
1091273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1092289bc588SBarry Smith       v = mat->v + j;
1093289bc588SBarry Smith       sum = 0.0;
1094273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
1095273d9f13SBarry Smith         sum += PetscAbsScalar(*v); v += A->m;
1096289bc588SBarry Smith       }
1097289bc588SBarry Smith       if (sum > *norm) *norm = sum;
1098289bc588SBarry Smith     }
1099b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11003a40ed3dSBarry Smith   } else {
110129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1102289bc588SBarry Smith   }
11033a40ed3dSBarry Smith   PetscFunctionReturn(0);
1104289bc588SBarry Smith }
1105289bc588SBarry Smith 
11065615d1e5SSatish Balay #undef __FUNC__
1107b0a32e0cSBarry Smith #define __FUNC__ "MatSetOption_SeqDense"
11088f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1109289bc588SBarry Smith {
1110c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
111167e560aaSBarry Smith 
11123a40ed3dSBarry Smith   PetscFunctionBegin;
1113273d9f13SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = PETSC_TRUE;
1114273d9f13SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = PETSC_FALSE;
11156d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
1116219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
11176d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
1118219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
11196d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
11206d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
11216d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
11226d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
11234787f768SSatish Balay            op == MAT_NEW_NONZERO_LOCATION_ERR ||
11246d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
112590f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
1126b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
1127c38d4ed2SBarry Smith            op == MAT_USE_HASH_TABLE) {
1128b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
1129c38d4ed2SBarry Smith   } else {
113029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
11313a40ed3dSBarry Smith   }
11323a40ed3dSBarry Smith   PetscFunctionReturn(0);
1133289bc588SBarry Smith }
1134289bc588SBarry Smith 
11355615d1e5SSatish Balay #undef __FUNC__
1136b0a32e0cSBarry Smith #define __FUNC__ "MatZeroEntries_SeqDense"
11378f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
11386f0a148fSBarry Smith {
1139ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1140549d3d68SSatish Balay   int          ierr;
11413a40ed3dSBarry Smith 
11423a40ed3dSBarry Smith   PetscFunctionBegin;
1143273d9f13SBarry Smith   ierr = PetscMemzero(l->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr);
11443a40ed3dSBarry Smith   PetscFunctionReturn(0);
11456f0a148fSBarry Smith }
11466f0a148fSBarry Smith 
11475615d1e5SSatish Balay #undef __FUNC__
1148b0a32e0cSBarry Smith #define __FUNC__ "MatGetBlockSize_SeqDense"
11498f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
11504e220ebcSLois Curfman McInnes {
11513a40ed3dSBarry Smith   PetscFunctionBegin;
11524e220ebcSLois Curfman McInnes   *bs = 1;
11533a40ed3dSBarry Smith   PetscFunctionReturn(0);
11544e220ebcSLois Curfman McInnes }
11554e220ebcSLois Curfman McInnes 
11565615d1e5SSatish Balay #undef __FUNC__
1157b0a32e0cSBarry Smith #define __FUNC__ "MatZeroRows_SeqDense"
11588f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
11596f0a148fSBarry Smith {
1160ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1161273d9f13SBarry Smith   int          n = A->n,i,j,ierr,N,*rows;
11626f0a148fSBarry Smith   Scalar       *slot;
116355659b69SBarry Smith 
11643a40ed3dSBarry Smith   PetscFunctionBegin;
1165e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
116678b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
11676f0a148fSBarry Smith   for (i=0; i<N; i++) {
11686f0a148fSBarry Smith     slot = l->v + rows[i];
11696f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
11706f0a148fSBarry Smith   }
11716f0a148fSBarry Smith   if (diag) {
11726f0a148fSBarry Smith     for (i=0; i<N; i++) {
11736f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1174260925b8SBarry Smith       *slot = *diag;
11756f0a148fSBarry Smith     }
11766f0a148fSBarry Smith   }
1177260925b8SBarry Smith   ISRestoreIndices(is,&rows);
11783a40ed3dSBarry Smith   PetscFunctionReturn(0);
11796f0a148fSBarry Smith }
1180557bce09SLois Curfman McInnes 
11815615d1e5SSatish Balay #undef __FUNC__
1182b0a32e0cSBarry Smith #define __FUNC__ "MatGetOwnershipRange_SeqDense"
11838f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1184d3e5ee88SLois Curfman McInnes {
11853a40ed3dSBarry Smith   PetscFunctionBegin;
11864c49b128SBarry Smith   if (m) *m = 0;
1187273d9f13SBarry Smith   if (n) *n = A->m;
11883a40ed3dSBarry Smith   PetscFunctionReturn(0);
1189d3e5ee88SLois Curfman McInnes }
1190d3e5ee88SLois Curfman McInnes 
11915615d1e5SSatish Balay #undef __FUNC__
1192b0a32e0cSBarry Smith #define __FUNC__ "MatGetArray_SeqDense"
11938f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array)
119464e87e97SBarry Smith {
1195c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
11963a40ed3dSBarry Smith 
11973a40ed3dSBarry Smith   PetscFunctionBegin;
119864e87e97SBarry Smith   *array = mat->v;
11993a40ed3dSBarry Smith   PetscFunctionReturn(0);
120064e87e97SBarry Smith }
12010754003eSLois Curfman McInnes 
12025615d1e5SSatish Balay #undef __FUNC__
1203b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreArray_SeqDense"
12048f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1205ff14e315SSatish Balay {
12063a40ed3dSBarry Smith   PetscFunctionBegin;
120709b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
12083a40ed3dSBarry Smith   PetscFunctionReturn(0);
1209ff14e315SSatish Balay }
12100754003eSLois Curfman McInnes 
12115615d1e5SSatish Balay #undef __FUNC__
1212b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrix_SeqDense"
12137b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
12140754003eSLois Curfman McInnes {
1215c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1216273d9f13SBarry Smith   int          i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
1217182d2002SSatish Balay   Scalar       *av,*bv,*v = mat->v;
12180754003eSLois Curfman McInnes   Mat          newmat;
12190754003eSLois Curfman McInnes 
12203a40ed3dSBarry Smith   PetscFunctionBegin;
122178b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
122278b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1223e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1224e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
12250754003eSLois Curfman McInnes 
1226182d2002SSatish Balay   /* Check submatrixcall */
1227182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1228182d2002SSatish Balay     int n_cols,n_rows;
1229182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
123029bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1231182d2002SSatish Balay     newmat = *B;
1232182d2002SSatish Balay   } else {
12330754003eSLois Curfman McInnes     /* Create and fill new matrix */
1234b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
1235182d2002SSatish Balay   }
1236182d2002SSatish Balay 
1237182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1238182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1239182d2002SSatish Balay 
1240182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1241182d2002SSatish Balay     av = v + m*icol[i];
1242182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1243182d2002SSatish Balay       *bv++ = av[irow[j]];
12440754003eSLois Curfman McInnes     }
12450754003eSLois Curfman McInnes   }
1246182d2002SSatish Balay 
1247182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
12486d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12496d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12500754003eSLois Curfman McInnes 
12510754003eSLois Curfman McInnes   /* Free work space */
125278b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
125378b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1254182d2002SSatish Balay   *B = newmat;
12553a40ed3dSBarry Smith   PetscFunctionReturn(0);
12560754003eSLois Curfman McInnes }
12570754003eSLois Curfman McInnes 
12585615d1e5SSatish Balay #undef __FUNC__
1259b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrices_SeqDense"
12607b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1261905e6a2fSBarry Smith {
1262905e6a2fSBarry Smith   int ierr,i;
1263905e6a2fSBarry Smith 
12643a40ed3dSBarry Smith   PetscFunctionBegin;
1265905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1266b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1267905e6a2fSBarry Smith   }
1268905e6a2fSBarry Smith 
1269905e6a2fSBarry Smith   for (i=0; i<n; i++) {
12706a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1271905e6a2fSBarry Smith   }
12723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1273905e6a2fSBarry Smith }
1274905e6a2fSBarry Smith 
12755615d1e5SSatish Balay #undef __FUNC__
1276b0a32e0cSBarry Smith #define __FUNC__ "MatCopy_SeqDense"
1277cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
12784b0e389bSBarry Smith {
12794b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
12803a40ed3dSBarry Smith   int          ierr;
1281273d9f13SBarry Smith   PetscTruth   flag;
12823a40ed3dSBarry Smith 
12833a40ed3dSBarry Smith   PetscFunctionBegin;
1284273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr);
1285273d9f13SBarry Smith   if (!flag) {
1286cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
12873a40ed3dSBarry Smith     PetscFunctionReturn(0);
12883a40ed3dSBarry Smith   }
1289273d9f13SBarry Smith   if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1290273d9f13SBarry Smith   ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr);
1291273d9f13SBarry Smith   PetscFunctionReturn(0);
1292273d9f13SBarry Smith }
1293273d9f13SBarry Smith 
1294273d9f13SBarry Smith #undef __FUNC__
1295b0a32e0cSBarry Smith #define __FUNC__ "MatSetUpPreallocation_SeqDense"
1296273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A)
1297273d9f13SBarry Smith {
1298273d9f13SBarry Smith   int        ierr;
1299273d9f13SBarry Smith 
1300273d9f13SBarry Smith   PetscFunctionBegin;
1301273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
13023a40ed3dSBarry Smith   PetscFunctionReturn(0);
13034b0e389bSBarry Smith }
13044b0e389bSBarry Smith 
1305289bc588SBarry Smith /* -------------------------------------------------------------------*/
1306a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1307905e6a2fSBarry Smith        MatGetRow_SeqDense,
1308905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1309905e6a2fSBarry Smith        MatMult_SeqDense,
1310905e6a2fSBarry Smith        MatMultAdd_SeqDense,
13117c922b88SBarry Smith        MatMultTranspose_SeqDense,
13127c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1313905e6a2fSBarry Smith        MatSolve_SeqDense,
1314905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
13157c922b88SBarry Smith        MatSolveTranspose_SeqDense,
13167c922b88SBarry Smith        MatSolveTransposeAdd_SeqDense,
1317905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1318905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1319ec8511deSBarry Smith        MatRelax_SeqDense,
1320ec8511deSBarry Smith        MatTranspose_SeqDense,
1321905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1322905e6a2fSBarry Smith        MatEqual_SeqDense,
1323905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1324905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1325905e6a2fSBarry Smith        MatNorm_SeqDense,
1326905e6a2fSBarry Smith        0,
1327905e6a2fSBarry Smith        0,
1328905e6a2fSBarry Smith        0,
1329905e6a2fSBarry Smith        MatSetOption_SeqDense,
1330905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1331905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1332905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1333905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1334905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1335905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1336273d9f13SBarry Smith        MatSetUpPreallocation_SeqDense,
1337273d9f13SBarry Smith        0,
1338905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1339905e6a2fSBarry Smith        0,
1340905e6a2fSBarry Smith        0,
1341905e6a2fSBarry Smith        MatGetArray_SeqDense,
1342905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
13435609ef8eSBarry Smith        MatDuplicate_SeqDense,
1344a5ae1ecdSBarry Smith        0,
1345a5ae1ecdSBarry Smith        0,
1346a5ae1ecdSBarry Smith        0,
1347a5ae1ecdSBarry Smith        0,
1348a5ae1ecdSBarry Smith        MatAXPY_SeqDense,
1349a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1350a5ae1ecdSBarry Smith        0,
13514b0e389bSBarry Smith        MatGetValues_SeqDense,
1352a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1353a5ae1ecdSBarry Smith        0,
1354a5ae1ecdSBarry Smith        MatScale_SeqDense,
1355a5ae1ecdSBarry Smith        0,
1356a5ae1ecdSBarry Smith        0,
1357a5ae1ecdSBarry Smith        0,
1358a5ae1ecdSBarry Smith        MatGetBlockSize_SeqDense,
1359a5ae1ecdSBarry Smith        0,
1360a5ae1ecdSBarry Smith        0,
1361a5ae1ecdSBarry Smith        0,
1362a5ae1ecdSBarry Smith        0,
1363a5ae1ecdSBarry Smith        0,
1364a5ae1ecdSBarry Smith        0,
1365a5ae1ecdSBarry Smith        0,
1366a5ae1ecdSBarry Smith        0,
1367a5ae1ecdSBarry Smith        0,
1368a5ae1ecdSBarry Smith        0,
1369e03a110bSBarry Smith        MatDestroy_SeqDense,
1370e03a110bSBarry Smith        MatView_SeqDense,
1371a5ae1ecdSBarry Smith        MatGetMaps_Petsc};
137290ace30eSBarry Smith 
13735615d1e5SSatish Balay #undef __FUNC__
1374b0a32e0cSBarry Smith #define __FUNC__ "MatCreateSeqDense"
13754b828684SBarry Smith /*@C
1376fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1377d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1378d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1379289bc588SBarry Smith 
1380db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1381db81eaa0SLois Curfman McInnes 
138220563c6bSBarry Smith    Input Parameters:
1383db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
13840c775827SLois Curfman McInnes .  m - number of rows
138518f449edSLois Curfman McInnes .  n - number of columns
1386db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1387dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
138820563c6bSBarry Smith 
138920563c6bSBarry Smith    Output Parameter:
139044cd7ae7SLois Curfman McInnes .  A - the matrix
139120563c6bSBarry Smith 
1392b259b22eSLois Curfman McInnes    Notes:
139318f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
139418f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1395b4fd4287SBarry Smith    set data=PETSC_NULL.
139618f449edSLois Curfman McInnes 
1397027ccd11SLois Curfman McInnes    Level: intermediate
1398027ccd11SLois Curfman McInnes 
1399dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1400d65003e9SLois Curfman McInnes 
1401db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
140220563c6bSBarry Smith @*/
140344cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1404289bc588SBarry Smith {
1405273d9f13SBarry Smith   int ierr;
14063b2fbd54SBarry Smith 
14073a40ed3dSBarry Smith   PetscFunctionBegin;
1408273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1409273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1410273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1411273d9f13SBarry Smith   PetscFunctionReturn(0);
1412273d9f13SBarry Smith }
1413273d9f13SBarry Smith 
1414273d9f13SBarry Smith #undef __FUNC__
1415b0a32e0cSBarry Smith #define __FUNC__ "MatSeqDensePreallocation"
1416273d9f13SBarry Smith /*@C
1417273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1418273d9f13SBarry Smith 
1419273d9f13SBarry Smith    Collective on MPI_Comm
1420273d9f13SBarry Smith 
1421273d9f13SBarry Smith    Input Parameters:
1422273d9f13SBarry Smith +  A - the matrix
1423273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1424273d9f13SBarry Smith 
1425273d9f13SBarry Smith    Notes:
1426273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1427273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1428273d9f13SBarry Smith    set data=PETSC_NULL.
1429273d9f13SBarry Smith 
1430273d9f13SBarry Smith    Level: intermediate
1431273d9f13SBarry Smith 
1432273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1433273d9f13SBarry Smith 
1434273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1435273d9f13SBarry Smith @*/
1436273d9f13SBarry Smith int MatSeqDenseSetPreallocation(Mat B,Scalar *data)
1437273d9f13SBarry Smith {
1438273d9f13SBarry Smith   Mat_SeqDense *b;
1439273d9f13SBarry Smith   int          ierr;
1440273d9f13SBarry Smith   PetscTruth   flg2;
1441273d9f13SBarry Smith 
1442273d9f13SBarry Smith   PetscFunctionBegin;
1443273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);CHKERRQ(ierr);
1444273d9f13SBarry Smith   if (!flg2) PetscFunctionReturn(0);
1445273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1446273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1447273d9f13SBarry Smith   if (!data) {
1448b0a32e0cSBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(Scalar),&b->v);CHKERRQ(ierr);
1449273d9f13SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(Scalar));CHKERRQ(ierr);
1450273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
1451b0a32e0cSBarry Smith     PetscLogObjectMemory(B,B->n*B->m*sizeof(Scalar));
1452273d9f13SBarry Smith   } else { /* user-allocated storage */
1453273d9f13SBarry Smith     b->v          = data;
1454273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1455273d9f13SBarry Smith   }
1456273d9f13SBarry Smith   PetscFunctionReturn(0);
1457273d9f13SBarry Smith }
1458273d9f13SBarry Smith 
1459273d9f13SBarry Smith EXTERN_C_BEGIN
1460273d9f13SBarry Smith #undef __FUNC__
1461b0a32e0cSBarry Smith #define __FUNC__ "MatCreate_SeqDense"
1462273d9f13SBarry Smith int MatCreate_SeqDense(Mat B)
1463273d9f13SBarry Smith {
1464273d9f13SBarry Smith   Mat_SeqDense *b;
1465273d9f13SBarry Smith   int          ierr,size;
1466273d9f13SBarry Smith 
1467273d9f13SBarry Smith   PetscFunctionBegin;
1468273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
146929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
147055659b69SBarry Smith 
1471273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1472273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1473273d9f13SBarry Smith 
1474b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1475549d3d68SSatish Balay   ierr            = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr);
1476549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
147744cd7ae7SLois Curfman McInnes   B->factor       = 0;
147890f02eecSBarry Smith   B->mapping      = 0;
1479b0a32e0cSBarry Smith   PetscLogObjectMemory(B,sizeof(struct _p_Mat));
148044cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
148118f449edSLois Curfman McInnes 
1482273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
1483273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1484a5ae1ecdSBarry Smith 
148544cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1486273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1487273d9f13SBarry Smith   b->v            = 0;
14884e220ebcSLois Curfman McInnes 
14893a40ed3dSBarry Smith   PetscFunctionReturn(0);
1490289bc588SBarry Smith }
1491273d9f13SBarry Smith EXTERN_C_END
14923b2fbd54SBarry Smith 
14933b2fbd54SBarry Smith 
14943b2fbd54SBarry Smith 
14953b2fbd54SBarry Smith 
1496