xref: /petsc/src/mat/impls/dense/seq/dense.c (revision ea709b57ddb8cc304f11df2e9a39d6a8fdb28b53)
1*ea709b57SSatish Balay /*$Id: dense.c,v 1.205 2001/08/06 21:15:12 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 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
1187828ca2SBarry Smith int MatAXPY_SeqDense(PetscScalar *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 
224a2ae208SSatish Balay #undef __FUNCT__
234a2ae208SSatish Balay #define __FUNCT__ "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;
2887828ca2SBarry Smith   PetscScalar  *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 
514a2ae208SSatish Balay #undef __FUNCT__
524a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
5387828ca2SBarry Smith int MatScale_SeqDense(PetscScalar *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.*/
684a2ae208SSatish Balay #undef __FUNCT__
694a2ae208SSatish Balay #define __FUNCT__ "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);
82ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRF)
83ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable.");
84ae7cfcebSSatish Balay #else
85273d9f13SBarry Smith   LAgetrf_(&A->m,&A->n,mat->v,&A->m,mat->pivots,&info);
86ae7cfcebSSatish Balay #endif
8729bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
8829bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
89b0a32e0cSBarry Smith   PetscLogFlops((2*A->n*A->n*A->n)/3);
903a40ed3dSBarry Smith   PetscFunctionReturn(0);
91289bc588SBarry Smith }
926ee01492SSatish Balay 
934a2ae208SSatish Balay #undef __FUNCT__
944a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
955609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
9602cad45dSBarry Smith {
9702cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l;
9802cad45dSBarry Smith   int          ierr;
9902cad45dSBarry Smith   Mat          newi;
10002cad45dSBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionBegin;
102273d9f13SBarry Smith   ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr);
10302cad45dSBarry Smith   l = (Mat_SeqDense*)newi->data;
1045609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
10587828ca2SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
10602cad45dSBarry Smith   }
10702cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
10802cad45dSBarry Smith   *newmat = newi;
1093a40ed3dSBarry Smith   PetscFunctionReturn(0);
11002cad45dSBarry Smith }
11102cad45dSBarry Smith 
1124a2ae208SSatish Balay #undef __FUNCT__
1134a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
114e03a110bSBarry Smith int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatLUInfo *info,Mat *fact)
115289bc588SBarry Smith {
1163a40ed3dSBarry Smith   int ierr;
1173a40ed3dSBarry Smith 
1183a40ed3dSBarry Smith   PetscFunctionBegin;
1195609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1203a40ed3dSBarry Smith   PetscFunctionReturn(0);
121289bc588SBarry Smith }
1226ee01492SSatish Balay 
1234a2ae208SSatish Balay #undef __FUNCT__
1244a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1258f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
126289bc588SBarry Smith {
12702cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1283a40ed3dSBarry Smith   int          ierr;
1293a40ed3dSBarry Smith 
1303a40ed3dSBarry Smith   PetscFunctionBegin;
13102cad45dSBarry Smith   /* copy the numerical values */
13287828ca2SBarry Smith   ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
13302cad45dSBarry Smith   (*fact)->factor = 0;
134e03a110bSBarry Smith   ierr = MatLUFactor(*fact,0,0,PETSC_NULL);CHKERRQ(ierr);
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
136289bc588SBarry Smith }
1376ee01492SSatish Balay 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
140329f5518SBarry Smith int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,PetscReal f,Mat *fact)
141289bc588SBarry Smith {
1423a40ed3dSBarry Smith   int ierr;
1433a40ed3dSBarry Smith 
1443a40ed3dSBarry Smith   PetscFunctionBegin;
1453a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1463a40ed3dSBarry Smith   PetscFunctionReturn(0);
147289bc588SBarry Smith }
1486ee01492SSatish Balay 
1494a2ae208SSatish Balay #undef __FUNCT__
1504a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
151329f5518SBarry Smith int MatCholeskyFactor_SeqDense(Mat A,IS perm,PetscReal f)
152289bc588SBarry Smith {
153c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
154606d414cSSatish Balay   int           info,ierr;
15555659b69SBarry Smith 
1563a40ed3dSBarry Smith   PetscFunctionBegin;
157752f5567SLois Curfman McInnes   if (mat->pivots) {
158606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
159b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-A->m*sizeof(int));
160752f5567SLois Curfman McInnes     mat->pivots = 0;
161752f5567SLois Curfman McInnes   }
162273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
163ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRF)
164ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable.");
165ae7cfcebSSatish Balay #else
166273d9f13SBarry Smith   LApotrf_("L",&A->n,mat->v,&A->m,&info);
167ae7cfcebSSatish Balay #endif
16829bbc08cSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
169c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
170b0a32e0cSBarry Smith   PetscLogFlops((A->n*A->n*A->n)/3);
1713a40ed3dSBarry Smith   PetscFunctionReturn(0);
172289bc588SBarry Smith }
173289bc588SBarry Smith 
1744a2ae208SSatish Balay #undef __FUNCT__
1750b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
1760b4b3355SBarry Smith int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
1770b4b3355SBarry Smith {
1780b4b3355SBarry Smith   int ierr;
1790b4b3355SBarry Smith 
1800b4b3355SBarry Smith   PetscFunctionBegin;
1810b4b3355SBarry Smith   ierr = MatCholeskyFactor_SeqDense(*fact,0,1.0);CHKERRQ(ierr);
1820b4b3355SBarry Smith   PetscFunctionReturn(0);
1830b4b3355SBarry Smith }
1840b4b3355SBarry Smith 
1850b4b3355SBarry Smith #undef __FUNCT__
1864a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
1878f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
188289bc588SBarry Smith {
189c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
190a57ad546SLois Curfman McInnes   int          one = 1,info,ierr;
19187828ca2SBarry Smith   PetscScalar  *x,*y;
19267e560aaSBarry Smith 
1933a40ed3dSBarry Smith   PetscFunctionBegin;
194273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
1957a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1967a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
19787828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
198c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
199ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
200ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
201ae7cfcebSSatish Balay #else
202273d9f13SBarry Smith     LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
203ae7cfcebSSatish Balay #endif
2047a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
205ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
206ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
207ae7cfcebSSatish Balay #else
208273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
209ae7cfcebSSatish Balay #endif
210289bc588SBarry Smith   }
21129bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
21229bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
2137a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2147a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
215b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2163a40ed3dSBarry Smith   PetscFunctionReturn(0);
217289bc588SBarry Smith }
2186ee01492SSatish Balay 
2194a2ae208SSatish Balay #undef __FUNCT__
2204a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
2217c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
222da3a660dSBarry Smith {
223c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2247a97a34bSBarry Smith   int          ierr,one = 1,info;
22587828ca2SBarry Smith   PetscScalar  *x,*y;
22667e560aaSBarry Smith 
2273a40ed3dSBarry Smith   PetscFunctionBegin;
228273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2297a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2307a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
23187828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
232752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
233da3a660dSBarry Smith   if (mat->pivots) {
234ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
235ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
236ae7cfcebSSatish Balay #else
237273d9f13SBarry Smith     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
238ae7cfcebSSatish Balay #endif
2397a97a34bSBarry Smith   } else {
240ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
241ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
242ae7cfcebSSatish Balay #else
243273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
244ae7cfcebSSatish Balay #endif
245da3a660dSBarry Smith   }
24629bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
2477a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2487a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
249b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
251da3a660dSBarry Smith }
2526ee01492SSatish Balay 
2534a2ae208SSatish Balay #undef __FUNCT__
2544a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
2558f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
256da3a660dSBarry Smith {
257c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2587a97a34bSBarry Smith   int          ierr,one = 1,info;
25987828ca2SBarry Smith   PetscScalar  *x,*y,sone = 1.0;
260da3a660dSBarry Smith   Vec          tmp = 0;
26167e560aaSBarry Smith 
2623a40ed3dSBarry Smith   PetscFunctionBegin;
2637a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2647a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
265273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
266da3a660dSBarry Smith   if (yy == zz) {
26778b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
268b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
26978b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
270da3a660dSBarry Smith   }
27187828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
272752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
273da3a660dSBarry Smith   if (mat->pivots) {
274ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
275ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
276ae7cfcebSSatish Balay #else
277273d9f13SBarry Smith     LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
278ae7cfcebSSatish Balay #endif
279a8c6a408SBarry Smith   } else {
280ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
281ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
282ae7cfcebSSatish Balay #else
283273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
284ae7cfcebSSatish Balay #endif
285da3a660dSBarry Smith   }
28629bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
287a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
288a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
2897a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2907a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
291b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
2923a40ed3dSBarry Smith   PetscFunctionReturn(0);
293da3a660dSBarry Smith }
29467e560aaSBarry Smith 
2954a2ae208SSatish Balay #undef __FUNCT__
2964a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
2977c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
298da3a660dSBarry Smith {
299c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
3006abc6512SBarry Smith   int           one = 1,info,ierr;
30187828ca2SBarry Smith   PetscScalar   *x,*y,sone = 1.0;
302da3a660dSBarry Smith   Vec           tmp;
30367e560aaSBarry Smith 
3043a40ed3dSBarry Smith   PetscFunctionBegin;
305273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3067a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3077a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
308da3a660dSBarry Smith   if (yy == zz) {
30978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
310b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
31178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
312da3a660dSBarry Smith   }
31387828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
314752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
315da3a660dSBarry Smith   if (mat->pivots) {
316ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
317ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
318ae7cfcebSSatish Balay #else
319273d9f13SBarry Smith     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
320ae7cfcebSSatish Balay #endif
3213a40ed3dSBarry Smith   } else {
322ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
323ae7cfcebSSatish Balay   SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
324ae7cfcebSSatish Balay #else
325273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
326ae7cfcebSSatish Balay #endif
327da3a660dSBarry Smith   }
32829bbc08cSBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
32990f02eecSBarry Smith   if (tmp) {
33090f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
33190f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3323a40ed3dSBarry Smith   } else {
33390f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
33490f02eecSBarry Smith   }
3357a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3367a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
337b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3383a40ed3dSBarry Smith   PetscFunctionReturn(0);
339da3a660dSBarry Smith }
340289bc588SBarry Smith /* ------------------------------------------------------------------*/
3414a2ae208SSatish Balay #undef __FUNCT__
3424a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
343329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
344329f5518SBarry Smith                           PetscReal shift,int its,Vec xx)
345289bc588SBarry Smith {
346c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
34787828ca2SBarry Smith   PetscScalar  *x,*b,*v = mat->v,zero = 0.0,xt;
348273d9f13SBarry Smith   int          ierr,m = A->m,i;
349aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
350bc1b551cSSatish Balay   int          o = 1;
351bc1b551cSSatish Balay #endif
352289bc588SBarry Smith 
3533a40ed3dSBarry Smith   PetscFunctionBegin;
354289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
3553a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
356bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
357289bc588SBarry Smith   }
3587a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3597a97a34bSBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
360289bc588SBarry Smith   while (its--) {
361289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
362289bc588SBarry Smith       for (i=0; i<m; i++) {
363aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
364f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
365f1747703SBarry Smith            not happy about returning a double complex */
366f1747703SBarry Smith         int    _i;
36787828ca2SBarry Smith         PetscScalar sum = b[i];
368f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
3693f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
370f1747703SBarry Smith         }
371f1747703SBarry Smith         xt = sum;
372f1747703SBarry Smith #else
373289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
374f1747703SBarry Smith #endif
37555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
376289bc588SBarry Smith       }
377289bc588SBarry Smith     }
378289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
379289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
380aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
381f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
382f1747703SBarry Smith            not happy about returning a double complex */
383f1747703SBarry Smith         int    _i;
38487828ca2SBarry Smith         PetscScalar sum = b[i];
385f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
3863f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
387f1747703SBarry Smith         }
388f1747703SBarry Smith         xt = sum;
389f1747703SBarry Smith #else
390289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
391f1747703SBarry Smith #endif
39255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
393289bc588SBarry Smith       }
394289bc588SBarry Smith     }
395289bc588SBarry Smith   }
396600d6d0bSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
3977a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3983a40ed3dSBarry Smith   PetscFunctionReturn(0);
399289bc588SBarry Smith }
400289bc588SBarry Smith 
401289bc588SBarry Smith /* -----------------------------------------------------------------*/
4024a2ae208SSatish Balay #undef __FUNCT__
4034a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
4047c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
405289bc588SBarry Smith {
406c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
40787828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
408*ea709b57SSatish Balay   int          ierr,_One=1;
409*ea709b57SSatish Balay   PetscScalar  _DOne=1.0,_DZero=0.0;
4103a40ed3dSBarry Smith 
4113a40ed3dSBarry Smith   PetscFunctionBegin;
412273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4137a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4147a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
415273d9f13SBarry Smith   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
4167a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4177a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
418b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->n);
4193a40ed3dSBarry Smith   PetscFunctionReturn(0);
420289bc588SBarry Smith }
4216ee01492SSatish Balay 
4224a2ae208SSatish Balay #undef __FUNCT__
4234a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
42444cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
425289bc588SBarry Smith {
426c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
42787828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
428329f5518SBarry Smith   int          ierr,_One=1;
4293a40ed3dSBarry Smith 
4303a40ed3dSBarry Smith   PetscFunctionBegin;
431273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4327a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4337a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
434273d9f13SBarry Smith   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,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 - A->m);
4383a40ed3dSBarry Smith   PetscFunctionReturn(0);
439289bc588SBarry Smith }
4406ee01492SSatish Balay 
4414a2ae208SSatish Balay #undef __FUNCT__
4424a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
44344cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
444289bc588SBarry Smith {
445c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
44687828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0;
447329f5518SBarry Smith   int          ierr,_One=1;
4483a40ed3dSBarry Smith 
4493a40ed3dSBarry Smith   PetscFunctionBegin;
450273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
451600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4527a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4537a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
454273d9f13SBarry Smith   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
4557a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4567a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
457b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4583a40ed3dSBarry Smith   PetscFunctionReturn(0);
459289bc588SBarry Smith }
4606ee01492SSatish Balay 
4614a2ae208SSatish Balay #undef __FUNCT__
4624a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
4637c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
464289bc588SBarry Smith {
465c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
46687828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
4677a97a34bSBarry Smith   int          ierr,_One=1;
46887828ca2SBarry Smith   PetscScalar  _DOne=1.0;
4693a40ed3dSBarry Smith 
4703a40ed3dSBarry Smith   PetscFunctionBegin;
471273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
472600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4737a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4747a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
475273d9f13SBarry Smith   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
4767a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4777a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
478b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4793a40ed3dSBarry Smith   PetscFunctionReturn(0);
480289bc588SBarry Smith }
481289bc588SBarry Smith 
482289bc588SBarry Smith /* -----------------------------------------------------------------*/
4834a2ae208SSatish Balay #undef __FUNCT__
4844a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
48587828ca2SBarry Smith int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
486289bc588SBarry Smith {
487c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
48887828ca2SBarry Smith   PetscScalar  *v;
489b0a32e0cSBarry Smith   int          i,ierr;
49067e560aaSBarry Smith 
4913a40ed3dSBarry Smith   PetscFunctionBegin;
492273d9f13SBarry Smith   *ncols = A->n;
493289bc588SBarry Smith   if (cols) {
494b0a32e0cSBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr);
495273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
496289bc588SBarry Smith   }
497289bc588SBarry Smith   if (vals) {
49887828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
499289bc588SBarry Smith     v    = mat->v + row;
500273d9f13SBarry Smith     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->m;}
501289bc588SBarry Smith   }
5023a40ed3dSBarry Smith   PetscFunctionReturn(0);
503289bc588SBarry Smith }
5046ee01492SSatish Balay 
5054a2ae208SSatish Balay #undef __FUNCT__
5064a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
50787828ca2SBarry Smith int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
508289bc588SBarry Smith {
509606d414cSSatish Balay   int ierr;
510606d414cSSatish Balay   PetscFunctionBegin;
511606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
512606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5133a40ed3dSBarry Smith   PetscFunctionReturn(0);
514289bc588SBarry Smith }
515289bc588SBarry Smith /* ----------------------------------------------------------------*/
5164a2ae208SSatish Balay #undef __FUNCT__
5174a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
5188f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
51987828ca2SBarry Smith                                     int *indexn,PetscScalar *v,InsertMode addv)
520289bc588SBarry Smith {
521c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
522289bc588SBarry Smith   int          i,j;
523d6dfbf8fSBarry Smith 
5243a40ed3dSBarry Smith   PetscFunctionBegin;
525289bc588SBarry Smith   if (!mat->roworiented) {
526dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
527289bc588SBarry Smith       for (j=0; j<n; j++) {
5285ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
529aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
53029bbc08cSBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
53158804f6dSBarry Smith #endif
532289bc588SBarry Smith         for (i=0; i<m; i++) {
5335ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
534aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
53529bbc08cSBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
53658804f6dSBarry Smith #endif
537273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
538289bc588SBarry Smith         }
539289bc588SBarry Smith       }
5403a40ed3dSBarry Smith     } else {
541289bc588SBarry Smith       for (j=0; j<n; j++) {
5425ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
543aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
54429bbc08cSBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
54558804f6dSBarry Smith #endif
546289bc588SBarry Smith         for (i=0; i<m; i++) {
5475ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
548aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
54929bbc08cSBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
55058804f6dSBarry Smith #endif
551273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
552289bc588SBarry Smith         }
553289bc588SBarry Smith       }
554289bc588SBarry Smith     }
5553a40ed3dSBarry Smith   } else {
556dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
557e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
5585ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
559aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
56029bbc08cSBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
56158804f6dSBarry Smith #endif
562e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
5635ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
564aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
56529bbc08cSBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
56658804f6dSBarry Smith #endif
567273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
568e8d4e0b9SBarry Smith         }
569e8d4e0b9SBarry Smith       }
5703a40ed3dSBarry Smith     } else {
571289bc588SBarry Smith       for (i=0; i<m; i++) {
5725ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
573aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
57429bbc08cSBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
57558804f6dSBarry Smith #endif
576289bc588SBarry Smith         for (j=0; j<n; j++) {
5775ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
578aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
57929bbc08cSBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
58058804f6dSBarry Smith #endif
581273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
582289bc588SBarry Smith         }
583289bc588SBarry Smith       }
584289bc588SBarry Smith     }
585e8d4e0b9SBarry Smith   }
5863a40ed3dSBarry Smith   PetscFunctionReturn(0);
587289bc588SBarry Smith }
588e8d4e0b9SBarry Smith 
5894a2ae208SSatish Balay #undef __FUNCT__
5904a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
59187828ca2SBarry Smith int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,PetscScalar *v)
592ae80bb75SLois Curfman McInnes {
593ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
594ae80bb75SLois Curfman McInnes   int          i,j;
59587828ca2SBarry Smith   PetscScalar  *vpt = v;
596ae80bb75SLois Curfman McInnes 
5973a40ed3dSBarry Smith   PetscFunctionBegin;
598ae80bb75SLois Curfman McInnes   /* row-oriented output */
599ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
600ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
601273d9f13SBarry Smith       *vpt++ = mat->v[indexn[j]*A->m + indexm[i]];
602ae80bb75SLois Curfman McInnes     }
603ae80bb75SLois Curfman McInnes   }
6043a40ed3dSBarry Smith   PetscFunctionReturn(0);
605ae80bb75SLois Curfman McInnes }
606ae80bb75SLois Curfman McInnes 
607289bc588SBarry Smith /* -----------------------------------------------------------------*/
608289bc588SBarry Smith 
609e090d566SSatish Balay #include "petscsys.h"
61088e32edaSLois Curfman McInnes 
611273d9f13SBarry Smith EXTERN_C_BEGIN
6124a2ae208SSatish Balay #undef __FUNCT__
6134a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
614b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A)
61588e32edaSLois Curfman McInnes {
61688e32edaSLois Curfman McInnes   Mat_SeqDense *a;
61788e32edaSLois Curfman McInnes   Mat          B;
61888e32edaSLois Curfman McInnes   int          *scols,i,j,nz,ierr,fd,header[4],size;
61988e32edaSLois Curfman McInnes   int          *rowlengths = 0,M,N,*cols;
62087828ca2SBarry Smith   PetscScalar  *vals,*svals,*v,*w;
62119bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
62288e32edaSLois Curfman McInnes 
6233a40ed3dSBarry Smith   PetscFunctionBegin;
624d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
62529bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
626b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6270752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
62829bbc08cSBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
62988e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
63088e32edaSLois Curfman McInnes 
63190ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
63290ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
63390ace30eSBarry Smith     B    = *A;
63490ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6354a41a4d3SSatish Balay     v    = a->v;
6364a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6374a41a4d3SSatish Balay        from row major to column major */
63887828ca2SBarry Smith     ierr = PetscMalloc((M*N+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
63990ace30eSBarry Smith     /* read in nonzero values */
6404a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6414a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6424a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6434a41a4d3SSatish Balay       for (i=0; i<M; i++) {
6444a41a4d3SSatish Balay         *v++ =w[i*N+j];
6454a41a4d3SSatish Balay       }
6464a41a4d3SSatish Balay     }
647606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
6486d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6496d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65090ace30eSBarry Smith   } else {
65188e32edaSLois Curfman McInnes     /* read row lengths */
652b0a32e0cSBarry Smith     ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr);
6530752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
65488e32edaSLois Curfman McInnes 
65588e32edaSLois Curfman McInnes     /* create our matrix */
656b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
65788e32edaSLois Curfman McInnes     B = *A;
65888e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
65988e32edaSLois Curfman McInnes     v = a->v;
66088e32edaSLois Curfman McInnes 
66188e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
662b0a32e0cSBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr);
663b0a32e0cSBarry Smith     cols = scols;
6640752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
66587828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
666b0a32e0cSBarry Smith     vals = svals;
6670752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
66888e32edaSLois Curfman McInnes 
66988e32edaSLois Curfman McInnes     /* insert into matrix */
67088e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
67188e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
67288e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
67388e32edaSLois Curfman McInnes     }
674606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
675606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
676606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
67788e32edaSLois Curfman McInnes 
6786d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6796d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68090ace30eSBarry Smith   }
6813a40ed3dSBarry Smith   PetscFunctionReturn(0);
68288e32edaSLois Curfman McInnes }
683273d9f13SBarry Smith EXTERN_C_END
68488e32edaSLois Curfman McInnes 
685e090d566SSatish Balay #include "petscsys.h"
686932b0c3eSLois Curfman McInnes 
6874a2ae208SSatish Balay #undef __FUNCT__
6884a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
689b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
690289bc588SBarry Smith {
691932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
692fb9695e5SSatish Balay   int               ierr,i,j;
693fb9695e5SSatish Balay   char              *name;
69487828ca2SBarry Smith   PetscScalar       *v;
695f3ef73ceSBarry Smith   PetscViewerFormat format;
696932b0c3eSLois Curfman McInnes 
6973a40ed3dSBarry Smith   PetscFunctionBegin;
698435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
699b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
700fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
7013a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
702fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
703b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
704273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
70544cd7ae7SLois Curfman McInnes       v = a->v + i;
706b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr);
707273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
708aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
709329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
710b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
711329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
712b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %d %g ",j,PetscRealPart(*v));CHKERRQ(ierr);
7136831982aSBarry Smith         }
71480cd9d93SLois Curfman McInnes #else
7156831982aSBarry Smith         if (*v) {
716b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %d %g ",j,*v);CHKERRQ(ierr);
7176831982aSBarry Smith         }
71880cd9d93SLois Curfman McInnes #endif
719273d9f13SBarry Smith         v += A->m;
72080cd9d93SLois Curfman McInnes       }
721b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
72280cd9d93SLois Curfman McInnes     }
723b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7243a40ed3dSBarry Smith   } else {
725b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
726aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
727ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
72847989497SBarry Smith     /* determine if matrix has all real values */
72947989497SBarry Smith     v = a->v;
730273d9f13SBarry Smith     for (i=0; i<A->m*A->n; i++) {
731ffac6cdbSBarry Smith       if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
73247989497SBarry Smith     }
73347989497SBarry Smith #endif
734fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
7353a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
736b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr);
737fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr);
738fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
739ffac6cdbSBarry Smith     }
740ffac6cdbSBarry Smith 
741273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
742932b0c3eSLois Curfman McInnes       v = a->v + i;
743273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
744aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
74547989497SBarry Smith         if (allreal) {
746b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); v += A->m;
74747989497SBarry Smith         } else {
748b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); v += A->m;
74947989497SBarry Smith         }
750289bc588SBarry Smith #else
751b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += A->m;
752289bc588SBarry Smith #endif
753289bc588SBarry Smith       }
754b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
755289bc588SBarry Smith     }
756fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
757b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
758ffac6cdbSBarry Smith     }
759b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
760da3a660dSBarry Smith   }
761b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7623a40ed3dSBarry Smith   PetscFunctionReturn(0);
763289bc588SBarry Smith }
764289bc588SBarry Smith 
7654a2ae208SSatish Balay #undef __FUNCT__
7664a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
767b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
768932b0c3eSLois Curfman McInnes {
769932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
770273d9f13SBarry Smith   int               ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
77187828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
772f3ef73ceSBarry Smith   PetscViewerFormat format;
773932b0c3eSLois Curfman McInnes 
7743a40ed3dSBarry Smith   PetscFunctionBegin;
775b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
77690ace30eSBarry Smith 
777b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
778fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
77990ace30eSBarry Smith     /* store the matrix as a dense matrix */
780b0a32e0cSBarry Smith     ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr);
78190ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
78290ace30eSBarry Smith     col_lens[1] = m;
78390ace30eSBarry Smith     col_lens[2] = n;
78490ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
7850752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr);
786606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
78790ace30eSBarry Smith 
78890ace30eSBarry Smith     /* write out matrix, by rows */
78987828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
79090ace30eSBarry Smith     v    = a->v;
79190ace30eSBarry Smith     for (i=0; i<m; i++) {
79290ace30eSBarry Smith       for (j=0; j<n; j++) {
79390ace30eSBarry Smith         vals[i + j*m] = *v++;
79490ace30eSBarry Smith       }
79590ace30eSBarry Smith     }
7960752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr);
797606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
79890ace30eSBarry Smith   } else {
799b0a32e0cSBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr);
800932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
801932b0c3eSLois Curfman McInnes     col_lens[1] = m;
802932b0c3eSLois Curfman McInnes     col_lens[2] = n;
803932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
804932b0c3eSLois Curfman McInnes 
805932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
806932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8070752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr);
808932b0c3eSLois Curfman McInnes 
809932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
810932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
811932b0c3eSLois Curfman McInnes     ict = 0;
812932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
813932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
814932b0c3eSLois Curfman McInnes     }
8150752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr);
816606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
817932b0c3eSLois Curfman McInnes 
818932b0c3eSLois Curfman McInnes     /* store nonzero values */
81987828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
820932b0c3eSLois Curfman McInnes     ict  = 0;
821932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
822932b0c3eSLois Curfman McInnes       v = a->v + i;
823932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
824273d9f13SBarry Smith         anonz[ict++] = *v; v += A->m;
825932b0c3eSLois Curfman McInnes       }
826932b0c3eSLois Curfman McInnes     }
8270752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr);
828606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
82990ace30eSBarry Smith   }
8303a40ed3dSBarry Smith   PetscFunctionReturn(0);
831932b0c3eSLois Curfman McInnes }
832932b0c3eSLois Curfman McInnes 
8334a2ae208SSatish Balay #undef __FUNCT__
8344a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
835b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
836f1af5d2fSBarry Smith {
837f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
838f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
839fb9695e5SSatish Balay   int               m = A->m,n = A->n,color,i,j,ierr;
84087828ca2SBarry Smith   PetscScalar       *v = a->v;
841b0a32e0cSBarry Smith   PetscViewer       viewer;
842b0a32e0cSBarry Smith   PetscDraw         popup;
843329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
844f3ef73ceSBarry Smith   PetscViewerFormat format;
845f1af5d2fSBarry Smith 
846f1af5d2fSBarry Smith   PetscFunctionBegin;
847f1af5d2fSBarry Smith 
848f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
849b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
850b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
851f1af5d2fSBarry Smith 
852f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
853fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
854f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
855b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
856f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
857f1af5d2fSBarry Smith       x_l = j;
858f1af5d2fSBarry Smith       x_r = x_l + 1.0;
859f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
860f1af5d2fSBarry Smith         y_l = m - i - 1.0;
861f1af5d2fSBarry Smith         y_r = y_l + 1.0;
862f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
863329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
864b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
865329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
866b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
867f1af5d2fSBarry Smith         } else {
868f1af5d2fSBarry Smith           continue;
869f1af5d2fSBarry Smith         }
870f1af5d2fSBarry Smith #else
871f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
872b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
873f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
874b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
875f1af5d2fSBarry Smith         } else {
876f1af5d2fSBarry Smith           continue;
877f1af5d2fSBarry Smith         }
878f1af5d2fSBarry Smith #endif
879b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
880f1af5d2fSBarry Smith       }
881f1af5d2fSBarry Smith     }
882f1af5d2fSBarry Smith   } else {
883f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
884f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
885f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
886f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
887f1af5d2fSBarry Smith     }
888b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
889b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
890b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
891f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
892f1af5d2fSBarry Smith       x_l = j;
893f1af5d2fSBarry Smith       x_r = x_l + 1.0;
894f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
895f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
896f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
897b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
898b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
899f1af5d2fSBarry Smith       }
900f1af5d2fSBarry Smith     }
901f1af5d2fSBarry Smith   }
902f1af5d2fSBarry Smith   PetscFunctionReturn(0);
903f1af5d2fSBarry Smith }
904f1af5d2fSBarry Smith 
9054a2ae208SSatish Balay #undef __FUNCT__
9064a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
907b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
908f1af5d2fSBarry Smith {
909b0a32e0cSBarry Smith   PetscDraw  draw;
910f1af5d2fSBarry Smith   PetscTruth isnull;
911329f5518SBarry Smith   PetscReal  xr,yr,xl,yl,h,w;
912f1af5d2fSBarry Smith   int        ierr;
913f1af5d2fSBarry Smith 
914f1af5d2fSBarry Smith   PetscFunctionBegin;
915b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
916b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
917f1af5d2fSBarry Smith   if (isnull == PETSC_TRUE) PetscFunctionReturn(0);
918f1af5d2fSBarry Smith 
919f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
920273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
921f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
922b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
923b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
924f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
925f1af5d2fSBarry Smith   PetscFunctionReturn(0);
926f1af5d2fSBarry Smith }
927f1af5d2fSBarry Smith 
9284a2ae208SSatish Balay #undef __FUNCT__
9294a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
930b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer)
931932b0c3eSLois Curfman McInnes {
932932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
933bcd2baecSBarry Smith   int          ierr;
934f1af5d2fSBarry Smith   PetscTruth   issocket,isascii,isbinary,isdraw;
935932b0c3eSLois Curfman McInnes 
9363a40ed3dSBarry Smith   PetscFunctionBegin;
937b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
938b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
939fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
940fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
9410f5bd95cSBarry Smith 
9420f5bd95cSBarry Smith   if (issocket) {
943b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
9440f5bd95cSBarry Smith   } else if (isascii) {
9453a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
9460f5bd95cSBarry Smith   } else if (isbinary) {
9473a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
948f1af5d2fSBarry Smith   } else if (isdraw) {
949f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
9505cd90555SBarry Smith   } else {
95129bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
952932b0c3eSLois Curfman McInnes   }
9533a40ed3dSBarry Smith   PetscFunctionReturn(0);
954932b0c3eSLois Curfman McInnes }
955289bc588SBarry Smith 
9564a2ae208SSatish Balay #undef __FUNCT__
9574a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
958c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
959289bc588SBarry Smith {
960ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
96190f02eecSBarry Smith   int          ierr;
96290f02eecSBarry Smith 
9633a40ed3dSBarry Smith   PetscFunctionBegin;
964aa482453SBarry Smith #if defined(PETSC_USE_LOG)
965b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
966a5a9c739SBarry Smith #endif
967606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
968606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
969606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
9703a40ed3dSBarry Smith   PetscFunctionReturn(0);
971289bc588SBarry Smith }
972289bc588SBarry Smith 
9734a2ae208SSatish Balay #undef __FUNCT__
9744a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
9758f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
976289bc588SBarry Smith {
977c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
978549d3d68SSatish Balay   int          k,j,m,n,ierr;
97987828ca2SBarry Smith   PetscScalar  *v,tmp;
98048b35521SBarry Smith 
9813a40ed3dSBarry Smith   PetscFunctionBegin;
982273d9f13SBarry Smith   v = mat->v; m = A->m; n = A->n;
9837c922b88SBarry Smith   if (!matout) { /* in place transpose */
984d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
98587828ca2SBarry Smith       PetscScalar *w;
98687828ca2SBarry Smith       ierr = PetscMalloc((m+1)*(n+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
987d64ed03dSBarry Smith       for (j=0; j<m; j++) {
988d64ed03dSBarry Smith         for (k=0; k<n; k++) {
989d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
990d64ed03dSBarry Smith         }
991d64ed03dSBarry Smith       }
99287828ca2SBarry Smith       ierr = PetscMemcpy(v,w,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
993606d414cSSatish Balay       ierr = PetscFree(w);CHKERRQ(ierr);
994d64ed03dSBarry Smith     } else {
995d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
996289bc588SBarry Smith         for (k=0; k<j; k++) {
997d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
998d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
999d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
1000289bc588SBarry Smith         }
1001289bc588SBarry Smith       }
1002d64ed03dSBarry Smith     }
10033a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1004d3e5ee88SLois Curfman McInnes     Mat          tmat;
1005ec8511deSBarry Smith     Mat_SeqDense *tmatd;
100687828ca2SBarry Smith     PetscScalar  *v2;
1007*ea709b57SSatish Balay 
1008273d9f13SBarry Smith     ierr  = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr);
1009ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10100de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1011d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10120de55854SLois Curfman McInnes       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m];
1013d3e5ee88SLois Curfman McInnes     }
10146d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10156d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1016d3e5ee88SLois Curfman McInnes     *matout = tmat;
101748b35521SBarry Smith   }
10183a40ed3dSBarry Smith   PetscFunctionReturn(0);
1019289bc588SBarry Smith }
1020289bc588SBarry Smith 
10214a2ae208SSatish Balay #undef __FUNCT__
10224a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
10238f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1024289bc588SBarry Smith {
1025c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1026c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1027273d9f13SBarry Smith   int          ierr,i;
102887828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
1029273d9f13SBarry Smith   PetscTruth   flag;
10309ea5d5aeSSatish Balay 
10313a40ed3dSBarry Smith   PetscFunctionBegin;
1032273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr);
1033273d9f13SBarry Smith   if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type  MATSEQDENSE");
1034273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1035273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1036273d9f13SBarry Smith   for (i=0; i<A1->m*A1->n; i++) {
10373a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1038289bc588SBarry Smith     v1++; v2++;
1039289bc588SBarry Smith   }
104077c4ece6SBarry Smith   *flg = PETSC_TRUE;
10413a40ed3dSBarry Smith   PetscFunctionReturn(0);
1042289bc588SBarry Smith }
1043289bc588SBarry Smith 
10444a2ae208SSatish Balay #undef __FUNCT__
10454a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
10468f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
1047289bc588SBarry Smith {
1048c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
10497a97a34bSBarry Smith   int          ierr,i,n,len;
105087828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
105144cd7ae7SLois Curfman McInnes 
10523a40ed3dSBarry Smith   PetscFunctionBegin;
10537a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
10547a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
1055600d6d0bSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1056273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1057273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
105844cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
1059273d9f13SBarry Smith     x[i] = mat->v[i*A->m + i];
1060289bc588SBarry Smith   }
10617a97a34bSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
10623a40ed3dSBarry Smith   PetscFunctionReturn(0);
1063289bc588SBarry Smith }
1064289bc588SBarry Smith 
10654a2ae208SSatish Balay #undef __FUNCT__
10664a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
10678f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1068289bc588SBarry Smith {
1069c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
107087828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
1071273d9f13SBarry Smith   int          ierr,i,j,m = A->m,n = A->n;
107255659b69SBarry Smith 
10733a40ed3dSBarry Smith   PetscFunctionBegin;
107428988994SBarry Smith   if (ll) {
10757a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1076600d6d0bSBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1077273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1078da3a660dSBarry Smith     for (i=0; i<m; i++) {
1079da3a660dSBarry Smith       x = l[i];
1080da3a660dSBarry Smith       v = mat->v + i;
1081da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1082da3a660dSBarry Smith     }
10837a97a34bSBarry Smith     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1084b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1085da3a660dSBarry Smith   }
108628988994SBarry Smith   if (rr) {
10877a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1088600d6d0bSBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1089273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1090da3a660dSBarry Smith     for (i=0; i<n; i++) {
1091da3a660dSBarry Smith       x = r[i];
1092da3a660dSBarry Smith       v = mat->v + i*m;
1093da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1094da3a660dSBarry Smith     }
10957a97a34bSBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1096b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1097da3a660dSBarry Smith   }
10983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1099289bc588SBarry Smith }
1100289bc588SBarry Smith 
11014a2ae208SSatish Balay #undef __FUNCT__
11024a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1103329f5518SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *norm)
1104289bc588SBarry Smith {
1105c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
110687828ca2SBarry Smith   PetscScalar  *v = mat->v;
1107329f5518SBarry Smith   PetscReal    sum = 0.0;
1108289bc588SBarry Smith   int          i,j;
110955659b69SBarry Smith 
11103a40ed3dSBarry Smith   PetscFunctionBegin;
1111289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1112273d9f13SBarry Smith     for (i=0; i<A->n*A->m; i++) {
1113aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1114329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1115289bc588SBarry Smith #else
1116289bc588SBarry Smith       sum += (*v)*(*v); v++;
1117289bc588SBarry Smith #endif
1118289bc588SBarry Smith     }
1119289bc588SBarry Smith     *norm = sqrt(sum);
1120b0a32e0cSBarry Smith     PetscLogFlops(2*A->n*A->m);
11213a40ed3dSBarry Smith   } else if (type == NORM_1) {
1122289bc588SBarry Smith     *norm = 0.0;
1123273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
1124289bc588SBarry Smith       sum = 0.0;
1125273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
112633a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1127289bc588SBarry Smith       }
1128289bc588SBarry Smith       if (sum > *norm) *norm = sum;
1129289bc588SBarry Smith     }
1130b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11313a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1132289bc588SBarry Smith     *norm = 0.0;
1133273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1134289bc588SBarry Smith       v = mat->v + j;
1135289bc588SBarry Smith       sum = 0.0;
1136273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
1137273d9f13SBarry Smith         sum += PetscAbsScalar(*v); v += A->m;
1138289bc588SBarry Smith       }
1139289bc588SBarry Smith       if (sum > *norm) *norm = sum;
1140289bc588SBarry Smith     }
1141b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11423a40ed3dSBarry Smith   } else {
114329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1144289bc588SBarry Smith   }
11453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1146289bc588SBarry Smith }
1147289bc588SBarry Smith 
11484a2ae208SSatish Balay #undef __FUNCT__
11494a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
11508f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1151289bc588SBarry Smith {
1152c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
115367e560aaSBarry Smith 
11543a40ed3dSBarry Smith   PetscFunctionBegin;
1155b5a2b587SKris Buschelman   switch (op) {
1156b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1157b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1158b5a2b587SKris Buschelman     break;
1159b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1160b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1161b5a2b587SKris Buschelman     break;
1162b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1163b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1164b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1165b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1166b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1167b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1168b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1169b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1170b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1171b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1172b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
1173d03495bdSKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
1174b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
1175b5a2b587SKris Buschelman     break;
1176b5a2b587SKris Buschelman   default:
117729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
11783a40ed3dSBarry Smith   }
11793a40ed3dSBarry Smith   PetscFunctionReturn(0);
1180289bc588SBarry Smith }
1181289bc588SBarry Smith 
11824a2ae208SSatish Balay #undef __FUNCT__
11834a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
11848f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
11856f0a148fSBarry Smith {
1186ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1187549d3d68SSatish Balay   int          ierr;
11883a40ed3dSBarry Smith 
11893a40ed3dSBarry Smith   PetscFunctionBegin;
119087828ca2SBarry Smith   ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
11913a40ed3dSBarry Smith   PetscFunctionReturn(0);
11926f0a148fSBarry Smith }
11936f0a148fSBarry Smith 
11944a2ae208SSatish Balay #undef __FUNCT__
11954a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense"
11968f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
11974e220ebcSLois Curfman McInnes {
11983a40ed3dSBarry Smith   PetscFunctionBegin;
11994e220ebcSLois Curfman McInnes   *bs = 1;
12003a40ed3dSBarry Smith   PetscFunctionReturn(0);
12014e220ebcSLois Curfman McInnes }
12024e220ebcSLois Curfman McInnes 
12034a2ae208SSatish Balay #undef __FUNCT__
12044a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
120587828ca2SBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,PetscScalar *diag)
12066f0a148fSBarry Smith {
1207ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1208273d9f13SBarry Smith   int          n = A->n,i,j,ierr,N,*rows;
120987828ca2SBarry Smith   PetscScalar  *slot;
121055659b69SBarry Smith 
12113a40ed3dSBarry Smith   PetscFunctionBegin;
1212e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
121378b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
12146f0a148fSBarry Smith   for (i=0; i<N; i++) {
12156f0a148fSBarry Smith     slot = l->v + rows[i];
12166f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
12176f0a148fSBarry Smith   }
12186f0a148fSBarry Smith   if (diag) {
12196f0a148fSBarry Smith     for (i=0; i<N; i++) {
12206f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1221260925b8SBarry Smith       *slot = *diag;
12226f0a148fSBarry Smith     }
12236f0a148fSBarry Smith   }
1224260925b8SBarry Smith   ISRestoreIndices(is,&rows);
12253a40ed3dSBarry Smith   PetscFunctionReturn(0);
12266f0a148fSBarry Smith }
1227557bce09SLois Curfman McInnes 
12284a2ae208SSatish Balay #undef __FUNCT__
12294a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_SeqDense"
12308f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1231d3e5ee88SLois Curfman McInnes {
12323a40ed3dSBarry Smith   PetscFunctionBegin;
12334c49b128SBarry Smith   if (m) *m = 0;
1234273d9f13SBarry Smith   if (n) *n = A->m;
12353a40ed3dSBarry Smith   PetscFunctionReturn(0);
1236d3e5ee88SLois Curfman McInnes }
1237d3e5ee88SLois Curfman McInnes 
12384a2ae208SSatish Balay #undef __FUNCT__
12394a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
124087828ca2SBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar **array)
124164e87e97SBarry Smith {
1242c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
12433a40ed3dSBarry Smith 
12443a40ed3dSBarry Smith   PetscFunctionBegin;
124564e87e97SBarry Smith   *array = mat->v;
12463a40ed3dSBarry Smith   PetscFunctionReturn(0);
124764e87e97SBarry Smith }
12480754003eSLois Curfman McInnes 
12494a2ae208SSatish Balay #undef __FUNCT__
12504a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
125187828ca2SBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar **array)
1252ff14e315SSatish Balay {
12533a40ed3dSBarry Smith   PetscFunctionBegin;
125409b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
12553a40ed3dSBarry Smith   PetscFunctionReturn(0);
1256ff14e315SSatish Balay }
12570754003eSLois Curfman McInnes 
12584a2ae208SSatish Balay #undef __FUNCT__
12594a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
12607b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
12610754003eSLois Curfman McInnes {
1262c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1263273d9f13SBarry Smith   int          i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
126487828ca2SBarry Smith   PetscScalar  *av,*bv,*v = mat->v;
12650754003eSLois Curfman McInnes   Mat          newmat;
12660754003eSLois Curfman McInnes 
12673a40ed3dSBarry Smith   PetscFunctionBegin;
126878b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
126978b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1270e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1271e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
12720754003eSLois Curfman McInnes 
1273182d2002SSatish Balay   /* Check submatrixcall */
1274182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1275182d2002SSatish Balay     int n_cols,n_rows;
1276182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
127729bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1278182d2002SSatish Balay     newmat = *B;
1279182d2002SSatish Balay   } else {
12800754003eSLois Curfman McInnes     /* Create and fill new matrix */
1281b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
1282182d2002SSatish Balay   }
1283182d2002SSatish Balay 
1284182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1285182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1286182d2002SSatish Balay 
1287182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1288182d2002SSatish Balay     av = v + m*icol[i];
1289182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1290182d2002SSatish Balay       *bv++ = av[irow[j]];
12910754003eSLois Curfman McInnes     }
12920754003eSLois Curfman McInnes   }
1293182d2002SSatish Balay 
1294182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
12956d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12966d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12970754003eSLois Curfman McInnes 
12980754003eSLois Curfman McInnes   /* Free work space */
129978b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
130078b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1301182d2002SSatish Balay   *B = newmat;
13023a40ed3dSBarry Smith   PetscFunctionReturn(0);
13030754003eSLois Curfman McInnes }
13040754003eSLois Curfman McInnes 
13054a2ae208SSatish Balay #undef __FUNCT__
13064a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
13077b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1308905e6a2fSBarry Smith {
1309905e6a2fSBarry Smith   int ierr,i;
1310905e6a2fSBarry Smith 
13113a40ed3dSBarry Smith   PetscFunctionBegin;
1312905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1313b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1314905e6a2fSBarry Smith   }
1315905e6a2fSBarry Smith 
1316905e6a2fSBarry Smith   for (i=0; i<n; i++) {
13176a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1318905e6a2fSBarry Smith   }
13193a40ed3dSBarry Smith   PetscFunctionReturn(0);
1320905e6a2fSBarry Smith }
1321905e6a2fSBarry Smith 
13224a2ae208SSatish Balay #undef __FUNCT__
13234a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1324cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
13254b0e389bSBarry Smith {
13264b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
13273a40ed3dSBarry Smith   int          ierr;
1328273d9f13SBarry Smith   PetscTruth   flag;
13293a40ed3dSBarry Smith 
13303a40ed3dSBarry Smith   PetscFunctionBegin;
1331273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr);
1332273d9f13SBarry Smith   if (!flag) {
1333cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
13343a40ed3dSBarry Smith     PetscFunctionReturn(0);
13353a40ed3dSBarry Smith   }
1336273d9f13SBarry Smith   if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
133787828ca2SBarry Smith   ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1338273d9f13SBarry Smith   PetscFunctionReturn(0);
1339273d9f13SBarry Smith }
1340273d9f13SBarry Smith 
13414a2ae208SSatish Balay #undef __FUNCT__
13424a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1343273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A)
1344273d9f13SBarry Smith {
1345273d9f13SBarry Smith   int        ierr;
1346273d9f13SBarry Smith 
1347273d9f13SBarry Smith   PetscFunctionBegin;
1348273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
13493a40ed3dSBarry Smith   PetscFunctionReturn(0);
13504b0e389bSBarry Smith }
13514b0e389bSBarry Smith 
1352289bc588SBarry Smith /* -------------------------------------------------------------------*/
1353a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1354905e6a2fSBarry Smith        MatGetRow_SeqDense,
1355905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1356905e6a2fSBarry Smith        MatMult_SeqDense,
1357905e6a2fSBarry Smith        MatMultAdd_SeqDense,
13587c922b88SBarry Smith        MatMultTranspose_SeqDense,
13597c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1360905e6a2fSBarry Smith        MatSolve_SeqDense,
1361905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
13627c922b88SBarry Smith        MatSolveTranspose_SeqDense,
13637c922b88SBarry Smith        MatSolveTransposeAdd_SeqDense,
1364905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1365905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1366ec8511deSBarry Smith        MatRelax_SeqDense,
1367ec8511deSBarry Smith        MatTranspose_SeqDense,
1368905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1369905e6a2fSBarry Smith        MatEqual_SeqDense,
1370905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1371905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1372905e6a2fSBarry Smith        MatNorm_SeqDense,
1373905e6a2fSBarry Smith        0,
1374905e6a2fSBarry Smith        0,
1375905e6a2fSBarry Smith        0,
1376905e6a2fSBarry Smith        MatSetOption_SeqDense,
1377905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1378905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1379905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1380905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1381905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1382905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1383273d9f13SBarry Smith        MatSetUpPreallocation_SeqDense,
1384273d9f13SBarry Smith        0,
1385905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1386905e6a2fSBarry Smith        0,
1387905e6a2fSBarry Smith        0,
1388905e6a2fSBarry Smith        MatGetArray_SeqDense,
1389905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
13905609ef8eSBarry Smith        MatDuplicate_SeqDense,
1391a5ae1ecdSBarry Smith        0,
1392a5ae1ecdSBarry Smith        0,
1393a5ae1ecdSBarry Smith        0,
1394a5ae1ecdSBarry Smith        0,
1395a5ae1ecdSBarry Smith        MatAXPY_SeqDense,
1396a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1397a5ae1ecdSBarry Smith        0,
13984b0e389bSBarry Smith        MatGetValues_SeqDense,
1399a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1400a5ae1ecdSBarry Smith        0,
1401a5ae1ecdSBarry Smith        MatScale_SeqDense,
1402a5ae1ecdSBarry Smith        0,
1403a5ae1ecdSBarry Smith        0,
1404a5ae1ecdSBarry Smith        0,
1405a5ae1ecdSBarry Smith        MatGetBlockSize_SeqDense,
1406a5ae1ecdSBarry Smith        0,
1407a5ae1ecdSBarry Smith        0,
1408a5ae1ecdSBarry Smith        0,
1409a5ae1ecdSBarry Smith        0,
1410a5ae1ecdSBarry Smith        0,
1411a5ae1ecdSBarry Smith        0,
1412a5ae1ecdSBarry Smith        0,
1413a5ae1ecdSBarry Smith        0,
1414a5ae1ecdSBarry Smith        0,
1415a5ae1ecdSBarry Smith        0,
1416e03a110bSBarry Smith        MatDestroy_SeqDense,
1417e03a110bSBarry Smith        MatView_SeqDense,
14188a124369SBarry Smith        MatGetPetscMaps_Petsc};
141990ace30eSBarry Smith 
14204a2ae208SSatish Balay #undef __FUNCT__
14214a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
14224b828684SBarry Smith /*@C
1423fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1424d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1425d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1426289bc588SBarry Smith 
1427db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1428db81eaa0SLois Curfman McInnes 
142920563c6bSBarry Smith    Input Parameters:
1430db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
14310c775827SLois Curfman McInnes .  m - number of rows
143218f449edSLois Curfman McInnes .  n - number of columns
1433db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1434dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
143520563c6bSBarry Smith 
143620563c6bSBarry Smith    Output Parameter:
143744cd7ae7SLois Curfman McInnes .  A - the matrix
143820563c6bSBarry Smith 
1439b259b22eSLois Curfman McInnes    Notes:
144018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
144118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1442b4fd4287SBarry Smith    set data=PETSC_NULL.
144318f449edSLois Curfman McInnes 
1444027ccd11SLois Curfman McInnes    Level: intermediate
1445027ccd11SLois Curfman McInnes 
1446dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1447d65003e9SLois Curfman McInnes 
1448db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
144920563c6bSBarry Smith @*/
145087828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A)
1451289bc588SBarry Smith {
1452273d9f13SBarry Smith   int ierr;
14533b2fbd54SBarry Smith 
14543a40ed3dSBarry Smith   PetscFunctionBegin;
1455273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1456273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1457273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1458273d9f13SBarry Smith   PetscFunctionReturn(0);
1459273d9f13SBarry Smith }
1460273d9f13SBarry Smith 
14614a2ae208SSatish Balay #undef __FUNCT__
14624a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1463273d9f13SBarry Smith /*@C
1464273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1465273d9f13SBarry Smith 
1466273d9f13SBarry Smith    Collective on MPI_Comm
1467273d9f13SBarry Smith 
1468273d9f13SBarry Smith    Input Parameters:
1469273d9f13SBarry Smith +  A - the matrix
1470273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1471273d9f13SBarry Smith 
1472273d9f13SBarry Smith    Notes:
1473273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1474273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1475273d9f13SBarry Smith    set data=PETSC_NULL.
1476273d9f13SBarry Smith 
1477273d9f13SBarry Smith    Level: intermediate
1478273d9f13SBarry Smith 
1479273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1480273d9f13SBarry Smith 
1481273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1482273d9f13SBarry Smith @*/
148387828ca2SBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar *data)
1484273d9f13SBarry Smith {
1485273d9f13SBarry Smith   Mat_SeqDense *b;
1486273d9f13SBarry Smith   int          ierr;
1487273d9f13SBarry Smith   PetscTruth   flg2;
1488273d9f13SBarry Smith 
1489273d9f13SBarry Smith   PetscFunctionBegin;
1490273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);CHKERRQ(ierr);
1491273d9f13SBarry Smith   if (!flg2) PetscFunctionReturn(0);
1492273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1493273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1494273d9f13SBarry Smith   if (!data) {
149587828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
149687828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1497273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
149887828ca2SBarry Smith     PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));
1499273d9f13SBarry Smith   } else { /* user-allocated storage */
1500273d9f13SBarry Smith     b->v          = data;
1501273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1502273d9f13SBarry Smith   }
1503273d9f13SBarry Smith   PetscFunctionReturn(0);
1504273d9f13SBarry Smith }
1505273d9f13SBarry Smith 
1506273d9f13SBarry Smith EXTERN_C_BEGIN
15074a2ae208SSatish Balay #undef __FUNCT__
15084a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1509273d9f13SBarry Smith int MatCreate_SeqDense(Mat B)
1510273d9f13SBarry Smith {
1511273d9f13SBarry Smith   Mat_SeqDense *b;
1512273d9f13SBarry Smith   int          ierr,size;
1513273d9f13SBarry Smith 
1514273d9f13SBarry Smith   PetscFunctionBegin;
1515273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
151629bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
151755659b69SBarry Smith 
1518273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1519273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1520273d9f13SBarry Smith 
1521b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1522549d3d68SSatish Balay   ierr            = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr);
1523549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
152444cd7ae7SLois Curfman McInnes   B->factor       = 0;
152590f02eecSBarry Smith   B->mapping      = 0;
1526b0a32e0cSBarry Smith   PetscLogObjectMemory(B,sizeof(struct _p_Mat));
152744cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
152818f449edSLois Curfman McInnes 
15298a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
15308a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1531a5ae1ecdSBarry Smith 
153244cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1533273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1534273d9f13SBarry Smith   b->v            = 0;
15354e220ebcSLois Curfman McInnes 
15363a40ed3dSBarry Smith   PetscFunctionReturn(0);
1537289bc588SBarry Smith }
1538273d9f13SBarry Smith EXTERN_C_END
15393b2fbd54SBarry Smith 
15403b2fbd54SBarry Smith 
15413b2fbd54SBarry Smith 
15423b2fbd54SBarry Smith 
1543