xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 456192e25527e6dfd94efe19f31bdabe85ed0893)
173f4d377SMatthew Knepley /*$Id: dense.c,v 1.208 2001/09/07 20:09:16 bsmith Exp $*/
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"
11607cd303SBarry Smith int MatAXPY_SeqDense(PetscScalar *alpha,Mat X,Mat Y,MatStructure str)
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);
8629bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
8729bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
882ffef20aSBarry Smith #endif
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);
16729bbc08cSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
1682ffef20aSBarry Smith #endif
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);
2032ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
204ae7cfcebSSatish Balay #endif
2057a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
206ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
207ae7cfcebSSatish Balay     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
208ae7cfcebSSatish Balay #else
209273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
2102ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
211ae7cfcebSSatish Balay #endif
212289bc588SBarry Smith   }
21329bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2147a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2157a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
216b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2173a40ed3dSBarry Smith   PetscFunctionReturn(0);
218289bc588SBarry Smith }
2196ee01492SSatish Balay 
2204a2ae208SSatish Balay #undef __FUNCT__
2214a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
2227c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
223da3a660dSBarry Smith {
224c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2257a97a34bSBarry Smith   int          ierr,one = 1,info;
22687828ca2SBarry Smith   PetscScalar  *x,*y;
22767e560aaSBarry Smith 
2283a40ed3dSBarry Smith   PetscFunctionBegin;
229273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2307a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2317a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
23287828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
233752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
234da3a660dSBarry Smith   if (mat->pivots) {
235ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
236ae7cfcebSSatish Balay     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
237ae7cfcebSSatish Balay #else
238273d9f13SBarry Smith     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
2392ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
240ae7cfcebSSatish Balay #endif
2417a97a34bSBarry Smith   } else {
242ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
243ae7cfcebSSatish Balay     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
244ae7cfcebSSatish Balay #else
245273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
2462ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
247ae7cfcebSSatish Balay #endif
248da3a660dSBarry Smith   }
2497a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2507a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
251b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
253da3a660dSBarry Smith }
2546ee01492SSatish Balay 
2554a2ae208SSatish Balay #undef __FUNCT__
2564a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
2578f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
258da3a660dSBarry Smith {
259c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2607a97a34bSBarry Smith   int          ierr,one = 1,info;
26187828ca2SBarry Smith   PetscScalar  *x,*y,sone = 1.0;
262da3a660dSBarry Smith   Vec          tmp = 0;
26367e560aaSBarry Smith 
2643a40ed3dSBarry Smith   PetscFunctionBegin;
2657a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2667a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
267273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
268da3a660dSBarry Smith   if (yy == zz) {
26978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
270b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
27178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
272da3a660dSBarry Smith   }
27387828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
274752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
275da3a660dSBarry Smith   if (mat->pivots) {
276ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
277ae7cfcebSSatish Balay     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
278ae7cfcebSSatish Balay #else
279273d9f13SBarry Smith     LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
2802ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
281ae7cfcebSSatish Balay #endif
282a8c6a408SBarry Smith   } else {
283ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
284ae7cfcebSSatish Balay     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
285ae7cfcebSSatish Balay #else
286273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
2872ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
288ae7cfcebSSatish Balay #endif
289da3a660dSBarry Smith   }
290a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
291a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
2927a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2937a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
294b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
2953a40ed3dSBarry Smith   PetscFunctionReturn(0);
296da3a660dSBarry Smith }
29767e560aaSBarry Smith 
2984a2ae208SSatish Balay #undef __FUNCT__
2994a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
3007c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
301da3a660dSBarry Smith {
302c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
3036abc6512SBarry Smith   int           one = 1,info,ierr;
30487828ca2SBarry Smith   PetscScalar   *x,*y,sone = 1.0;
305da3a660dSBarry Smith   Vec           tmp;
30667e560aaSBarry Smith 
3073a40ed3dSBarry Smith   PetscFunctionBegin;
308273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3097a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3107a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
311da3a660dSBarry Smith   if (yy == zz) {
31278b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
313b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
31478b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
315da3a660dSBarry Smith   }
31687828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
317752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
318da3a660dSBarry Smith   if (mat->pivots) {
319ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
320ae7cfcebSSatish Balay     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
321ae7cfcebSSatish Balay #else
322273d9f13SBarry Smith     LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
3232ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
324ae7cfcebSSatish Balay #endif
3253a40ed3dSBarry Smith   } else {
326ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
327ae7cfcebSSatish Balay     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
328ae7cfcebSSatish Balay #else
329273d9f13SBarry Smith     LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
3302ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
331ae7cfcebSSatish Balay #endif
332da3a660dSBarry Smith   }
33390f02eecSBarry Smith   if (tmp) {
33490f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
33590f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3363a40ed3dSBarry Smith   } else {
33790f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
33890f02eecSBarry Smith   }
3397a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3407a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
341b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3423a40ed3dSBarry Smith   PetscFunctionReturn(0);
343da3a660dSBarry Smith }
344289bc588SBarry Smith /* ------------------------------------------------------------------*/
3454a2ae208SSatish Balay #undef __FUNCT__
3464a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
347329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
348c14dc6b6SHong Zhang                           PetscReal shift,int its,int lits,Vec xx)
349289bc588SBarry Smith {
350c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
35187828ca2SBarry Smith   PetscScalar  *x,*b,*v = mat->v,zero = 0.0,xt;
352273d9f13SBarry Smith   int          ierr,m = A->m,i;
353aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
354bc1b551cSSatish Balay   int          o = 1;
355bc1b551cSSatish Balay #endif
356289bc588SBarry Smith 
3573a40ed3dSBarry Smith   PetscFunctionBegin;
358289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
3593a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
360bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
361289bc588SBarry Smith   }
3627a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3637a97a34bSBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
364b965ef7fSBarry Smith   its  = its*lits;
36591723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
366289bc588SBarry Smith   while (its--) {
367289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
368289bc588SBarry Smith       for (i=0; i<m; i++) {
369aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
370f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
371f1747703SBarry Smith            not happy about returning a double complex */
372f1747703SBarry Smith         int         _i;
37387828ca2SBarry Smith         PetscScalar sum = b[i];
374f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
3753f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
376f1747703SBarry Smith         }
377f1747703SBarry Smith         xt = sum;
378f1747703SBarry Smith #else
379289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
380f1747703SBarry Smith #endif
38155a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
382289bc588SBarry Smith       }
383289bc588SBarry Smith     }
384289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
385289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
386aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
387f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
388f1747703SBarry Smith            not happy about returning a double complex */
389f1747703SBarry Smith         int         _i;
39087828ca2SBarry Smith         PetscScalar sum = b[i];
391f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
3923f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
393f1747703SBarry Smith         }
394f1747703SBarry Smith         xt = sum;
395f1747703SBarry Smith #else
396289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
397f1747703SBarry Smith #endif
39855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
399289bc588SBarry Smith       }
400289bc588SBarry Smith     }
401289bc588SBarry Smith   }
402600d6d0bSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4037a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4043a40ed3dSBarry Smith   PetscFunctionReturn(0);
405289bc588SBarry Smith }
406289bc588SBarry Smith 
407289bc588SBarry Smith /* -----------------------------------------------------------------*/
4084a2ae208SSatish Balay #undef __FUNCT__
4094a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
4107c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
411289bc588SBarry Smith {
412c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
41387828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
414ea709b57SSatish Balay   int          ierr,_One=1;
415ea709b57SSatish Balay   PetscScalar  _DOne=1.0,_DZero=0.0;
4163a40ed3dSBarry Smith 
4173a40ed3dSBarry Smith   PetscFunctionBegin;
418273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4197a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4207a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
421273d9f13SBarry Smith   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
4227a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4237a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
424b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->n);
4253a40ed3dSBarry Smith   PetscFunctionReturn(0);
426289bc588SBarry Smith }
4276ee01492SSatish Balay 
4284a2ae208SSatish Balay #undef __FUNCT__
4294a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
43044cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
431289bc588SBarry Smith {
432c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
43387828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
434329f5518SBarry Smith   int          ierr,_One=1;
4353a40ed3dSBarry Smith 
4363a40ed3dSBarry Smith   PetscFunctionBegin;
437273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4387a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4397a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
440273d9f13SBarry Smith   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
4417a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4427a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
443b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->m);
4443a40ed3dSBarry Smith   PetscFunctionReturn(0);
445289bc588SBarry Smith }
4466ee01492SSatish Balay 
4474a2ae208SSatish Balay #undef __FUNCT__
4484a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
44944cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
450289bc588SBarry Smith {
451c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
45287828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0;
453329f5518SBarry Smith   int          ierr,_One=1;
4543a40ed3dSBarry Smith 
4553a40ed3dSBarry Smith   PetscFunctionBegin;
456273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
457600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4587a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4597a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
460273d9f13SBarry Smith   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
4617a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4627a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
463b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4643a40ed3dSBarry Smith   PetscFunctionReturn(0);
465289bc588SBarry Smith }
4666ee01492SSatish Balay 
4674a2ae208SSatish Balay #undef __FUNCT__
4684a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
4697c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
470289bc588SBarry Smith {
471c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
47287828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
4737a97a34bSBarry Smith   int          ierr,_One=1;
47487828ca2SBarry Smith   PetscScalar  _DOne=1.0;
4753a40ed3dSBarry Smith 
4763a40ed3dSBarry Smith   PetscFunctionBegin;
477273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
478600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4797a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4807a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
481273d9f13SBarry Smith   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
4827a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4837a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
484b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4853a40ed3dSBarry Smith   PetscFunctionReturn(0);
486289bc588SBarry Smith }
487289bc588SBarry Smith 
488289bc588SBarry Smith /* -----------------------------------------------------------------*/
4894a2ae208SSatish Balay #undef __FUNCT__
4904a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
49187828ca2SBarry Smith int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
492289bc588SBarry Smith {
493c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
49487828ca2SBarry Smith   PetscScalar  *v;
495b0a32e0cSBarry Smith   int          i,ierr;
49667e560aaSBarry Smith 
4973a40ed3dSBarry Smith   PetscFunctionBegin;
498273d9f13SBarry Smith   *ncols = A->n;
499289bc588SBarry Smith   if (cols) {
500b0a32e0cSBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr);
501273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
502289bc588SBarry Smith   }
503289bc588SBarry Smith   if (vals) {
50487828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
505289bc588SBarry Smith     v    = mat->v + row;
506273d9f13SBarry Smith     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->m;}
507289bc588SBarry Smith   }
5083a40ed3dSBarry Smith   PetscFunctionReturn(0);
509289bc588SBarry Smith }
5106ee01492SSatish Balay 
5114a2ae208SSatish Balay #undef __FUNCT__
5124a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
51387828ca2SBarry Smith int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
514289bc588SBarry Smith {
515606d414cSSatish Balay   int ierr;
516606d414cSSatish Balay   PetscFunctionBegin;
517606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
518606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5193a40ed3dSBarry Smith   PetscFunctionReturn(0);
520289bc588SBarry Smith }
521289bc588SBarry Smith /* ----------------------------------------------------------------*/
5224a2ae208SSatish Balay #undef __FUNCT__
5234a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
5248f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
52587828ca2SBarry Smith                                     int *indexn,PetscScalar *v,InsertMode addv)
526289bc588SBarry Smith {
527c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
528289bc588SBarry Smith   int          i,j;
529d6dfbf8fSBarry Smith 
5303a40ed3dSBarry Smith   PetscFunctionBegin;
531289bc588SBarry Smith   if (!mat->roworiented) {
532dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
533289bc588SBarry Smith       for (j=0; j<n; j++) {
5345ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
535aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
53629bbc08cSBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
53758804f6dSBarry Smith #endif
538289bc588SBarry Smith         for (i=0; i<m; i++) {
5395ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
540aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
54129bbc08cSBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
54258804f6dSBarry Smith #endif
543273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
544289bc588SBarry Smith         }
545289bc588SBarry Smith       }
5463a40ed3dSBarry Smith     } else {
547289bc588SBarry Smith       for (j=0; j<n; j++) {
5485ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
549aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
55029bbc08cSBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
55158804f6dSBarry Smith #endif
552289bc588SBarry Smith         for (i=0; i<m; i++) {
5535ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
554aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
55529bbc08cSBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
55658804f6dSBarry Smith #endif
557273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
558289bc588SBarry Smith         }
559289bc588SBarry Smith       }
560289bc588SBarry Smith     }
5613a40ed3dSBarry Smith   } else {
562dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
563e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
5645ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
565aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
56629bbc08cSBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
56758804f6dSBarry Smith #endif
568e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
5695ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
570aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
57129bbc08cSBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
57258804f6dSBarry Smith #endif
573273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] = *v++;
574e8d4e0b9SBarry Smith         }
575e8d4e0b9SBarry Smith       }
5763a40ed3dSBarry Smith     } else {
577289bc588SBarry Smith       for (i=0; i<m; i++) {
5785ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
579aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
58029bbc08cSBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
58158804f6dSBarry Smith #endif
582289bc588SBarry Smith         for (j=0; j<n; j++) {
5835ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
584aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
58529bbc08cSBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
58658804f6dSBarry Smith #endif
587273d9f13SBarry Smith           mat->v[indexn[j]*A->m + indexm[i]] += *v++;
588289bc588SBarry Smith         }
589289bc588SBarry Smith       }
590289bc588SBarry Smith     }
591e8d4e0b9SBarry Smith   }
5923a40ed3dSBarry Smith   PetscFunctionReturn(0);
593289bc588SBarry Smith }
594e8d4e0b9SBarry Smith 
5954a2ae208SSatish Balay #undef __FUNCT__
5964a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
59787828ca2SBarry Smith int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,PetscScalar *v)
598ae80bb75SLois Curfman McInnes {
599ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
600ae80bb75SLois Curfman McInnes   int          i,j;
60187828ca2SBarry Smith   PetscScalar  *vpt = v;
602ae80bb75SLois Curfman McInnes 
6033a40ed3dSBarry Smith   PetscFunctionBegin;
604ae80bb75SLois Curfman McInnes   /* row-oriented output */
605ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
606ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
607273d9f13SBarry Smith       *vpt++ = mat->v[indexn[j]*A->m + indexm[i]];
608ae80bb75SLois Curfman McInnes     }
609ae80bb75SLois Curfman McInnes   }
6103a40ed3dSBarry Smith   PetscFunctionReturn(0);
611ae80bb75SLois Curfman McInnes }
612ae80bb75SLois Curfman McInnes 
613289bc588SBarry Smith /* -----------------------------------------------------------------*/
614289bc588SBarry Smith 
615e090d566SSatish Balay #include "petscsys.h"
61688e32edaSLois Curfman McInnes 
617273d9f13SBarry Smith EXTERN_C_BEGIN
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
620b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A)
62188e32edaSLois Curfman McInnes {
62288e32edaSLois Curfman McInnes   Mat_SeqDense *a;
62388e32edaSLois Curfman McInnes   Mat          B;
62488e32edaSLois Curfman McInnes   int          *scols,i,j,nz,ierr,fd,header[4],size;
62588e32edaSLois Curfman McInnes   int          *rowlengths = 0,M,N,*cols;
62687828ca2SBarry Smith   PetscScalar  *vals,*svals,*v,*w;
62719bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
62888e32edaSLois Curfman McInnes 
6293a40ed3dSBarry Smith   PetscFunctionBegin;
630d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
63129bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
632b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6330752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
634552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
63588e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
63688e32edaSLois Curfman McInnes 
63790ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
63890ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
63990ace30eSBarry Smith     B    = *A;
64090ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6414a41a4d3SSatish Balay     v    = a->v;
6424a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6434a41a4d3SSatish Balay        from row major to column major */
64487828ca2SBarry Smith     ierr = PetscMalloc((M*N+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
64590ace30eSBarry Smith     /* read in nonzero values */
6464a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6474a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6484a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6494a41a4d3SSatish Balay       for (i=0; i<M; i++) {
6504a41a4d3SSatish Balay         *v++ =w[i*N+j];
6514a41a4d3SSatish Balay       }
6524a41a4d3SSatish Balay     }
653606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
6546d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6556d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65690ace30eSBarry Smith   } else {
65788e32edaSLois Curfman McInnes     /* read row lengths */
658b0a32e0cSBarry Smith     ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr);
6590752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
66088e32edaSLois Curfman McInnes 
66188e32edaSLois Curfman McInnes     /* create our matrix */
662b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
66388e32edaSLois Curfman McInnes     B = *A;
66488e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
66588e32edaSLois Curfman McInnes     v = a->v;
66688e32edaSLois Curfman McInnes 
66788e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
668b0a32e0cSBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr);
669b0a32e0cSBarry Smith     cols = scols;
6700752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
67187828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
672b0a32e0cSBarry Smith     vals = svals;
6730752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
67488e32edaSLois Curfman McInnes 
67588e32edaSLois Curfman McInnes     /* insert into matrix */
67688e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
67788e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
67888e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
67988e32edaSLois Curfman McInnes     }
680606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
681606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
682606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
68388e32edaSLois Curfman McInnes 
6846d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6856d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68690ace30eSBarry Smith   }
6873a40ed3dSBarry Smith   PetscFunctionReturn(0);
68888e32edaSLois Curfman McInnes }
689273d9f13SBarry Smith EXTERN_C_END
69088e32edaSLois Curfman McInnes 
691e090d566SSatish Balay #include "petscsys.h"
692932b0c3eSLois Curfman McInnes 
6934a2ae208SSatish Balay #undef __FUNCT__
6944a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
695b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
696289bc588SBarry Smith {
697932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
698fb9695e5SSatish Balay   int               ierr,i,j;
699fb9695e5SSatish Balay   char              *name;
70087828ca2SBarry Smith   PetscScalar       *v;
701f3ef73ceSBarry Smith   PetscViewerFormat format;
702932b0c3eSLois Curfman McInnes 
7033a40ed3dSBarry Smith   PetscFunctionBegin;
704435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
705b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
706*456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7073a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
708fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
709b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
710273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
71144cd7ae7SLois Curfman McInnes       v = a->v + i;
712b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr);
713273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
714aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
715329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
71662b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
717329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
71862b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7196831982aSBarry Smith         }
72080cd9d93SLois Curfman McInnes #else
7216831982aSBarry Smith         if (*v) {
72262b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr);
7236831982aSBarry Smith         }
72480cd9d93SLois Curfman McInnes #endif
725273d9f13SBarry Smith         v += A->m;
72680cd9d93SLois Curfman McInnes       }
727b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
72880cd9d93SLois Curfman McInnes     }
729b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7303a40ed3dSBarry Smith   } else {
731b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
732aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
733ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
73447989497SBarry Smith     /* determine if matrix has all real values */
73547989497SBarry Smith     v = a->v;
736273d9f13SBarry Smith     for (i=0; i<A->m*A->n; i++) {
737ffac6cdbSBarry Smith       if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
73847989497SBarry Smith     }
73947989497SBarry Smith #endif
740fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
7413a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
742b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr);
743fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr);
744fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
745ffac6cdbSBarry Smith     }
746ffac6cdbSBarry Smith 
747273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
748932b0c3eSLois Curfman McInnes       v = a->v + i;
749273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
750aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
75147989497SBarry Smith         if (allreal) {
752b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); v += A->m;
75347989497SBarry Smith         } else {
754b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); v += A->m;
75547989497SBarry Smith         }
756289bc588SBarry Smith #else
757b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += A->m;
758289bc588SBarry Smith #endif
759289bc588SBarry Smith       }
760b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
761289bc588SBarry Smith     }
762fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
763b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
764ffac6cdbSBarry Smith     }
765b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
766da3a660dSBarry Smith   }
767b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7683a40ed3dSBarry Smith   PetscFunctionReturn(0);
769289bc588SBarry Smith }
770289bc588SBarry Smith 
7714a2ae208SSatish Balay #undef __FUNCT__
7724a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
773b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
774932b0c3eSLois Curfman McInnes {
775932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
776273d9f13SBarry Smith   int               ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
77787828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
778f3ef73ceSBarry Smith   PetscViewerFormat format;
779932b0c3eSLois Curfman McInnes 
7803a40ed3dSBarry Smith   PetscFunctionBegin;
781b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
78290ace30eSBarry Smith 
783b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
784fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
78590ace30eSBarry Smith     /* store the matrix as a dense matrix */
786b0a32e0cSBarry Smith     ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr);
787552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
78890ace30eSBarry Smith     col_lens[1] = m;
78990ace30eSBarry Smith     col_lens[2] = n;
79090ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
7910752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr);
792606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
79390ace30eSBarry Smith 
79490ace30eSBarry Smith     /* write out matrix, by rows */
79587828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
79690ace30eSBarry Smith     v    = a->v;
79790ace30eSBarry Smith     for (i=0; i<m; i++) {
79890ace30eSBarry Smith       for (j=0; j<n; j++) {
79990ace30eSBarry Smith         vals[i + j*m] = *v++;
80090ace30eSBarry Smith       }
80190ace30eSBarry Smith     }
8020752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr);
803606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
80490ace30eSBarry Smith   } else {
805b0a32e0cSBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr);
806552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
807932b0c3eSLois Curfman McInnes     col_lens[1] = m;
808932b0c3eSLois Curfman McInnes     col_lens[2] = n;
809932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
810932b0c3eSLois Curfman McInnes 
811932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
812932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8130752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr);
814932b0c3eSLois Curfman McInnes 
815932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
816932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
817932b0c3eSLois Curfman McInnes     ict = 0;
818932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
819932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
820932b0c3eSLois Curfman McInnes     }
8210752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr);
822606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
823932b0c3eSLois Curfman McInnes 
824932b0c3eSLois Curfman McInnes     /* store nonzero values */
82587828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
826932b0c3eSLois Curfman McInnes     ict  = 0;
827932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
828932b0c3eSLois Curfman McInnes       v = a->v + i;
829932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
830273d9f13SBarry Smith         anonz[ict++] = *v; v += A->m;
831932b0c3eSLois Curfman McInnes       }
832932b0c3eSLois Curfman McInnes     }
8330752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr);
834606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
83590ace30eSBarry Smith   }
8363a40ed3dSBarry Smith   PetscFunctionReturn(0);
837932b0c3eSLois Curfman McInnes }
838932b0c3eSLois Curfman McInnes 
8394a2ae208SSatish Balay #undef __FUNCT__
8404a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
841b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
842f1af5d2fSBarry Smith {
843f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
844f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
845fb9695e5SSatish Balay   int               m = A->m,n = A->n,color,i,j,ierr;
84687828ca2SBarry Smith   PetscScalar       *v = a->v;
847b0a32e0cSBarry Smith   PetscViewer       viewer;
848b0a32e0cSBarry Smith   PetscDraw         popup;
849329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
850f3ef73ceSBarry Smith   PetscViewerFormat format;
851f1af5d2fSBarry Smith 
852f1af5d2fSBarry Smith   PetscFunctionBegin;
853f1af5d2fSBarry Smith 
854f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
855b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
856b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
857f1af5d2fSBarry Smith 
858f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
859fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
860f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
861b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
862f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
863f1af5d2fSBarry Smith       x_l = j;
864f1af5d2fSBarry Smith       x_r = x_l + 1.0;
865f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
866f1af5d2fSBarry Smith         y_l = m - i - 1.0;
867f1af5d2fSBarry Smith         y_r = y_l + 1.0;
868f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
869329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
870b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
871329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
872b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
873f1af5d2fSBarry Smith         } else {
874f1af5d2fSBarry Smith           continue;
875f1af5d2fSBarry Smith         }
876f1af5d2fSBarry Smith #else
877f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
878b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
879f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
880b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
881f1af5d2fSBarry Smith         } else {
882f1af5d2fSBarry Smith           continue;
883f1af5d2fSBarry Smith         }
884f1af5d2fSBarry Smith #endif
885b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
886f1af5d2fSBarry Smith       }
887f1af5d2fSBarry Smith     }
888f1af5d2fSBarry Smith   } else {
889f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
890f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
891f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
892f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
893f1af5d2fSBarry Smith     }
894b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
895b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
896b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
897f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
898f1af5d2fSBarry Smith       x_l = j;
899f1af5d2fSBarry Smith       x_r = x_l + 1.0;
900f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
901f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
902f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
903b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
904b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
905f1af5d2fSBarry Smith       }
906f1af5d2fSBarry Smith     }
907f1af5d2fSBarry Smith   }
908f1af5d2fSBarry Smith   PetscFunctionReturn(0);
909f1af5d2fSBarry Smith }
910f1af5d2fSBarry Smith 
9114a2ae208SSatish Balay #undef __FUNCT__
9124a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
913b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
914f1af5d2fSBarry Smith {
915b0a32e0cSBarry Smith   PetscDraw  draw;
916f1af5d2fSBarry Smith   PetscTruth isnull;
917329f5518SBarry Smith   PetscReal  xr,yr,xl,yl,h,w;
918f1af5d2fSBarry Smith   int        ierr;
919f1af5d2fSBarry Smith 
920f1af5d2fSBarry Smith   PetscFunctionBegin;
921b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
922b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
923f1af5d2fSBarry Smith   if (isnull == PETSC_TRUE) PetscFunctionReturn(0);
924f1af5d2fSBarry Smith 
925f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
926273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
927f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
928b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
929b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
930f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
931f1af5d2fSBarry Smith   PetscFunctionReturn(0);
932f1af5d2fSBarry Smith }
933f1af5d2fSBarry Smith 
9344a2ae208SSatish Balay #undef __FUNCT__
9354a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
936b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer)
937932b0c3eSLois Curfman McInnes {
938932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
939bcd2baecSBarry Smith   int          ierr;
940f1af5d2fSBarry Smith   PetscTruth   issocket,isascii,isbinary,isdraw;
941932b0c3eSLois Curfman McInnes 
9423a40ed3dSBarry Smith   PetscFunctionBegin;
943b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
944b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
945fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
946fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
9470f5bd95cSBarry Smith 
9480f5bd95cSBarry Smith   if (issocket) {
949b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
9500f5bd95cSBarry Smith   } else if (isascii) {
9513a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
9520f5bd95cSBarry Smith   } else if (isbinary) {
9533a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
954f1af5d2fSBarry Smith   } else if (isdraw) {
955f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
9565cd90555SBarry Smith   } else {
95729bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
958932b0c3eSLois Curfman McInnes   }
9593a40ed3dSBarry Smith   PetscFunctionReturn(0);
960932b0c3eSLois Curfman McInnes }
961289bc588SBarry Smith 
9624a2ae208SSatish Balay #undef __FUNCT__
9634a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
964c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
965289bc588SBarry Smith {
966ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
96790f02eecSBarry Smith   int          ierr;
96890f02eecSBarry Smith 
9693a40ed3dSBarry Smith   PetscFunctionBegin;
970aa482453SBarry Smith #if defined(PETSC_USE_LOG)
971b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
972a5a9c739SBarry Smith #endif
973606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
974606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
975606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
9763a40ed3dSBarry Smith   PetscFunctionReturn(0);
977289bc588SBarry Smith }
978289bc588SBarry Smith 
9794a2ae208SSatish Balay #undef __FUNCT__
9804a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
9818f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
982289bc588SBarry Smith {
983c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
984549d3d68SSatish Balay   int          k,j,m,n,ierr;
98587828ca2SBarry Smith   PetscScalar  *v,tmp;
98648b35521SBarry Smith 
9873a40ed3dSBarry Smith   PetscFunctionBegin;
988273d9f13SBarry Smith   v = mat->v; m = A->m; n = A->n;
9897c922b88SBarry Smith   if (!matout) { /* in place transpose */
990d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
99187828ca2SBarry Smith       PetscScalar *w;
99287828ca2SBarry Smith       ierr = PetscMalloc((m+1)*(n+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
993d64ed03dSBarry Smith       for (j=0; j<m; j++) {
994d64ed03dSBarry Smith         for (k=0; k<n; k++) {
995d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
996d64ed03dSBarry Smith         }
997d64ed03dSBarry Smith       }
99887828ca2SBarry Smith       ierr = PetscMemcpy(v,w,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
999606d414cSSatish Balay       ierr = PetscFree(w);CHKERRQ(ierr);
1000d64ed03dSBarry Smith     } else {
1001d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1002289bc588SBarry Smith         for (k=0; k<j; k++) {
1003d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
1004d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
1005d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
1006289bc588SBarry Smith         }
1007289bc588SBarry Smith       }
1008d64ed03dSBarry Smith     }
10093a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1010d3e5ee88SLois Curfman McInnes     Mat          tmat;
1011ec8511deSBarry Smith     Mat_SeqDense *tmatd;
101287828ca2SBarry Smith     PetscScalar  *v2;
1013ea709b57SSatish Balay 
1014273d9f13SBarry Smith     ierr  = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr);
1015ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10160de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1017d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10180de55854SLois Curfman McInnes       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m];
1019d3e5ee88SLois Curfman McInnes     }
10206d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10216d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1022d3e5ee88SLois Curfman McInnes     *matout = tmat;
102348b35521SBarry Smith   }
10243a40ed3dSBarry Smith   PetscFunctionReturn(0);
1025289bc588SBarry Smith }
1026289bc588SBarry Smith 
10274a2ae208SSatish Balay #undef __FUNCT__
10284a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
10298f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1030289bc588SBarry Smith {
1031c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1032c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1033273d9f13SBarry Smith   int          ierr,i;
103487828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
1035273d9f13SBarry Smith   PetscTruth   flag;
10369ea5d5aeSSatish Balay 
10373a40ed3dSBarry Smith   PetscFunctionBegin;
1038273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr);
1039273d9f13SBarry Smith   if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type  MATSEQDENSE");
1040273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1041273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1042273d9f13SBarry Smith   for (i=0; i<A1->m*A1->n; i++) {
10433a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1044289bc588SBarry Smith     v1++; v2++;
1045289bc588SBarry Smith   }
104677c4ece6SBarry Smith   *flg = PETSC_TRUE;
10473a40ed3dSBarry Smith   PetscFunctionReturn(0);
1048289bc588SBarry Smith }
1049289bc588SBarry Smith 
10504a2ae208SSatish Balay #undef __FUNCT__
10514a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
10528f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
1053289bc588SBarry Smith {
1054c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
10557a97a34bSBarry Smith   int          ierr,i,n,len;
105687828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
105744cd7ae7SLois Curfman McInnes 
10583a40ed3dSBarry Smith   PetscFunctionBegin;
10597a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
10607a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
1061600d6d0bSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1062273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1063273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
106444cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
1065273d9f13SBarry Smith     x[i] = mat->v[i*A->m + i];
1066289bc588SBarry Smith   }
10677a97a34bSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
10683a40ed3dSBarry Smith   PetscFunctionReturn(0);
1069289bc588SBarry Smith }
1070289bc588SBarry Smith 
10714a2ae208SSatish Balay #undef __FUNCT__
10724a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
10738f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1074289bc588SBarry Smith {
1075c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
107687828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
1077273d9f13SBarry Smith   int          ierr,i,j,m = A->m,n = A->n;
107855659b69SBarry Smith 
10793a40ed3dSBarry Smith   PetscFunctionBegin;
108028988994SBarry Smith   if (ll) {
10817a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1082600d6d0bSBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1083273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1084da3a660dSBarry Smith     for (i=0; i<m; i++) {
1085da3a660dSBarry Smith       x = l[i];
1086da3a660dSBarry Smith       v = mat->v + i;
1087da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1088da3a660dSBarry Smith     }
10897a97a34bSBarry Smith     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1090b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1091da3a660dSBarry Smith   }
109228988994SBarry Smith   if (rr) {
10937a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1094600d6d0bSBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1095273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1096da3a660dSBarry Smith     for (i=0; i<n; i++) {
1097da3a660dSBarry Smith       x = r[i];
1098da3a660dSBarry Smith       v = mat->v + i*m;
1099da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1100da3a660dSBarry Smith     }
11017a97a34bSBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1102b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1103da3a660dSBarry Smith   }
11043a40ed3dSBarry Smith   PetscFunctionReturn(0);
1105289bc588SBarry Smith }
1106289bc588SBarry Smith 
11074a2ae208SSatish Balay #undef __FUNCT__
11084a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1109064f8208SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1110289bc588SBarry Smith {
1111c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
111287828ca2SBarry Smith   PetscScalar  *v = mat->v;
1113329f5518SBarry Smith   PetscReal    sum = 0.0;
1114289bc588SBarry Smith   int          i,j;
111555659b69SBarry Smith 
11163a40ed3dSBarry Smith   PetscFunctionBegin;
1117289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1118273d9f13SBarry Smith     for (i=0; i<A->n*A->m; i++) {
1119aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1120329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1121289bc588SBarry Smith #else
1122289bc588SBarry Smith       sum += (*v)*(*v); v++;
1123289bc588SBarry Smith #endif
1124289bc588SBarry Smith     }
1125064f8208SBarry Smith     *nrm = sqrt(sum);
1126b0a32e0cSBarry Smith     PetscLogFlops(2*A->n*A->m);
11273a40ed3dSBarry Smith   } else if (type == NORM_1) {
1128064f8208SBarry Smith     *nrm = 0.0;
1129273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
1130289bc588SBarry Smith       sum = 0.0;
1131273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
113233a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1133289bc588SBarry Smith       }
1134064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1135289bc588SBarry Smith     }
1136b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11373a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1138064f8208SBarry Smith     *nrm = 0.0;
1139273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1140289bc588SBarry Smith       v = mat->v + j;
1141289bc588SBarry Smith       sum = 0.0;
1142273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
1143273d9f13SBarry Smith         sum += PetscAbsScalar(*v); v += A->m;
1144289bc588SBarry Smith       }
1145064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1146289bc588SBarry Smith     }
1147b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11483a40ed3dSBarry Smith   } else {
114929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1150289bc588SBarry Smith   }
11513a40ed3dSBarry Smith   PetscFunctionReturn(0);
1152289bc588SBarry Smith }
1153289bc588SBarry Smith 
11544a2ae208SSatish Balay #undef __FUNCT__
11554a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
11568f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1157289bc588SBarry Smith {
1158c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
115967e560aaSBarry Smith 
11603a40ed3dSBarry Smith   PetscFunctionBegin;
1161b5a2b587SKris Buschelman   switch (op) {
1162b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1163b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1164b5a2b587SKris Buschelman     break;
1165b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1166b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1167b5a2b587SKris Buschelman     break;
1168b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1169b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1170b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1171b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1172b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1173b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1174b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1175b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1176b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1177b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1178b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
1179d03495bdSKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
1180b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
1181b5a2b587SKris Buschelman     break;
1182b5a2b587SKris Buschelman   default:
118329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
11843a40ed3dSBarry Smith   }
11853a40ed3dSBarry Smith   PetscFunctionReturn(0);
1186289bc588SBarry Smith }
1187289bc588SBarry Smith 
11884a2ae208SSatish Balay #undef __FUNCT__
11894a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
11908f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
11916f0a148fSBarry Smith {
1192ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1193549d3d68SSatish Balay   int          ierr;
11943a40ed3dSBarry Smith 
11953a40ed3dSBarry Smith   PetscFunctionBegin;
119687828ca2SBarry Smith   ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
11973a40ed3dSBarry Smith   PetscFunctionReturn(0);
11986f0a148fSBarry Smith }
11996f0a148fSBarry Smith 
12004a2ae208SSatish Balay #undef __FUNCT__
12014a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense"
12028f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
12034e220ebcSLois Curfman McInnes {
12043a40ed3dSBarry Smith   PetscFunctionBegin;
12054e220ebcSLois Curfman McInnes   *bs = 1;
12063a40ed3dSBarry Smith   PetscFunctionReturn(0);
12074e220ebcSLois Curfman McInnes }
12084e220ebcSLois Curfman McInnes 
12094a2ae208SSatish Balay #undef __FUNCT__
12104a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
121187828ca2SBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,PetscScalar *diag)
12126f0a148fSBarry Smith {
1213ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1214273d9f13SBarry Smith   int          n = A->n,i,j,ierr,N,*rows;
121587828ca2SBarry Smith   PetscScalar  *slot;
121655659b69SBarry Smith 
12173a40ed3dSBarry Smith   PetscFunctionBegin;
1218e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
121978b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
12206f0a148fSBarry Smith   for (i=0; i<N; i++) {
12216f0a148fSBarry Smith     slot = l->v + rows[i];
12226f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
12236f0a148fSBarry Smith   }
12246f0a148fSBarry Smith   if (diag) {
12256f0a148fSBarry Smith     for (i=0; i<N; i++) {
12266f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1227260925b8SBarry Smith       *slot = *diag;
12286f0a148fSBarry Smith     }
12296f0a148fSBarry Smith   }
1230260925b8SBarry Smith   ISRestoreIndices(is,&rows);
12313a40ed3dSBarry Smith   PetscFunctionReturn(0);
12326f0a148fSBarry Smith }
1233557bce09SLois Curfman McInnes 
12344a2ae208SSatish Balay #undef __FUNCT__
12354a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
123687828ca2SBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar **array)
123764e87e97SBarry Smith {
1238c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
12393a40ed3dSBarry Smith 
12403a40ed3dSBarry Smith   PetscFunctionBegin;
124164e87e97SBarry Smith   *array = mat->v;
12423a40ed3dSBarry Smith   PetscFunctionReturn(0);
124364e87e97SBarry Smith }
12440754003eSLois Curfman McInnes 
12454a2ae208SSatish Balay #undef __FUNCT__
12464a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
124787828ca2SBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar **array)
1248ff14e315SSatish Balay {
12493a40ed3dSBarry Smith   PetscFunctionBegin;
125009b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
12513a40ed3dSBarry Smith   PetscFunctionReturn(0);
1252ff14e315SSatish Balay }
12530754003eSLois Curfman McInnes 
12544a2ae208SSatish Balay #undef __FUNCT__
12554a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
12567b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
12570754003eSLois Curfman McInnes {
1258c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1259273d9f13SBarry Smith   int          i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
126087828ca2SBarry Smith   PetscScalar  *av,*bv,*v = mat->v;
12610754003eSLois Curfman McInnes   Mat          newmat;
12620754003eSLois Curfman McInnes 
12633a40ed3dSBarry Smith   PetscFunctionBegin;
126478b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
126578b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1266e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1267e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
12680754003eSLois Curfman McInnes 
1269182d2002SSatish Balay   /* Check submatrixcall */
1270182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1271182d2002SSatish Balay     int n_cols,n_rows;
1272182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
127329bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1274182d2002SSatish Balay     newmat = *B;
1275182d2002SSatish Balay   } else {
12760754003eSLois Curfman McInnes     /* Create and fill new matrix */
1277b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
1278182d2002SSatish Balay   }
1279182d2002SSatish Balay 
1280182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1281182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1282182d2002SSatish Balay 
1283182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1284182d2002SSatish Balay     av = v + m*icol[i];
1285182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1286182d2002SSatish Balay       *bv++ = av[irow[j]];
12870754003eSLois Curfman McInnes     }
12880754003eSLois Curfman McInnes   }
1289182d2002SSatish Balay 
1290182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
12916d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12926d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12930754003eSLois Curfman McInnes 
12940754003eSLois Curfman McInnes   /* Free work space */
129578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
129678b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1297182d2002SSatish Balay   *B = newmat;
12983a40ed3dSBarry Smith   PetscFunctionReturn(0);
12990754003eSLois Curfman McInnes }
13000754003eSLois Curfman McInnes 
13014a2ae208SSatish Balay #undef __FUNCT__
13024a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
13037b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1304905e6a2fSBarry Smith {
1305905e6a2fSBarry Smith   int ierr,i;
1306905e6a2fSBarry Smith 
13073a40ed3dSBarry Smith   PetscFunctionBegin;
1308905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1309b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1310905e6a2fSBarry Smith   }
1311905e6a2fSBarry Smith 
1312905e6a2fSBarry Smith   for (i=0; i<n; i++) {
13136a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1314905e6a2fSBarry Smith   }
13153a40ed3dSBarry Smith   PetscFunctionReturn(0);
1316905e6a2fSBarry Smith }
1317905e6a2fSBarry Smith 
13184a2ae208SSatish Balay #undef __FUNCT__
13194a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1320cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
13214b0e389bSBarry Smith {
13224b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
13233a40ed3dSBarry Smith   int          ierr;
1324273d9f13SBarry Smith   PetscTruth   flag;
13253a40ed3dSBarry Smith 
13263a40ed3dSBarry Smith   PetscFunctionBegin;
1327273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr);
1328273d9f13SBarry Smith   if (!flag) {
1329cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
13303a40ed3dSBarry Smith     PetscFunctionReturn(0);
13313a40ed3dSBarry Smith   }
1332273d9f13SBarry Smith   if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
133387828ca2SBarry Smith   ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1334273d9f13SBarry Smith   PetscFunctionReturn(0);
1335273d9f13SBarry Smith }
1336273d9f13SBarry Smith 
13374a2ae208SSatish Balay #undef __FUNCT__
13384a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1339273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A)
1340273d9f13SBarry Smith {
1341273d9f13SBarry Smith   int        ierr;
1342273d9f13SBarry Smith 
1343273d9f13SBarry Smith   PetscFunctionBegin;
1344273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
13453a40ed3dSBarry Smith   PetscFunctionReturn(0);
13464b0e389bSBarry Smith }
13474b0e389bSBarry Smith 
1348289bc588SBarry Smith /* -------------------------------------------------------------------*/
1349a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1350905e6a2fSBarry Smith        MatGetRow_SeqDense,
1351905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1352905e6a2fSBarry Smith        MatMult_SeqDense,
1353905e6a2fSBarry Smith        MatMultAdd_SeqDense,
13547c922b88SBarry Smith        MatMultTranspose_SeqDense,
13557c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1356905e6a2fSBarry Smith        MatSolve_SeqDense,
1357905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
13587c922b88SBarry Smith        MatSolveTranspose_SeqDense,
13597c922b88SBarry Smith        MatSolveTransposeAdd_SeqDense,
1360905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1361905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1362ec8511deSBarry Smith        MatRelax_SeqDense,
1363ec8511deSBarry Smith        MatTranspose_SeqDense,
1364905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1365905e6a2fSBarry Smith        MatEqual_SeqDense,
1366905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1367905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1368905e6a2fSBarry Smith        MatNorm_SeqDense,
1369905e6a2fSBarry Smith        0,
1370905e6a2fSBarry Smith        0,
1371905e6a2fSBarry Smith        0,
1372905e6a2fSBarry Smith        MatSetOption_SeqDense,
1373905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1374905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1375905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1376905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1377905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1378905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1379273d9f13SBarry Smith        MatSetUpPreallocation_SeqDense,
1380273d9f13SBarry Smith        0,
1381905e6a2fSBarry Smith        0,
1382905e6a2fSBarry Smith        MatGetArray_SeqDense,
1383905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
13845609ef8eSBarry Smith        MatDuplicate_SeqDense,
1385a5ae1ecdSBarry Smith        0,
1386a5ae1ecdSBarry Smith        0,
1387a5ae1ecdSBarry Smith        0,
1388a5ae1ecdSBarry Smith        0,
1389a5ae1ecdSBarry Smith        MatAXPY_SeqDense,
1390a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1391a5ae1ecdSBarry Smith        0,
13924b0e389bSBarry Smith        MatGetValues_SeqDense,
1393a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1394a5ae1ecdSBarry Smith        0,
1395a5ae1ecdSBarry Smith        MatScale_SeqDense,
1396a5ae1ecdSBarry Smith        0,
1397a5ae1ecdSBarry Smith        0,
1398a5ae1ecdSBarry Smith        0,
1399a5ae1ecdSBarry Smith        MatGetBlockSize_SeqDense,
1400a5ae1ecdSBarry Smith        0,
1401a5ae1ecdSBarry Smith        0,
1402a5ae1ecdSBarry Smith        0,
1403a5ae1ecdSBarry Smith        0,
1404a5ae1ecdSBarry Smith        0,
1405a5ae1ecdSBarry Smith        0,
1406a5ae1ecdSBarry Smith        0,
1407a5ae1ecdSBarry Smith        0,
1408a5ae1ecdSBarry Smith        0,
1409a5ae1ecdSBarry Smith        0,
1410e03a110bSBarry Smith        MatDestroy_SeqDense,
1411e03a110bSBarry Smith        MatView_SeqDense,
14128a124369SBarry Smith        MatGetPetscMaps_Petsc};
141390ace30eSBarry Smith 
14144a2ae208SSatish Balay #undef __FUNCT__
14154a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
14164b828684SBarry Smith /*@C
1417fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1418d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1419d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1420289bc588SBarry Smith 
1421db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1422db81eaa0SLois Curfman McInnes 
142320563c6bSBarry Smith    Input Parameters:
1424db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
14250c775827SLois Curfman McInnes .  m - number of rows
142618f449edSLois Curfman McInnes .  n - number of columns
1427db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1428dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
142920563c6bSBarry Smith 
143020563c6bSBarry Smith    Output Parameter:
143144cd7ae7SLois Curfman McInnes .  A - the matrix
143220563c6bSBarry Smith 
1433b259b22eSLois Curfman McInnes    Notes:
143418f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
143518f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1436b4fd4287SBarry Smith    set data=PETSC_NULL.
143718f449edSLois Curfman McInnes 
1438027ccd11SLois Curfman McInnes    Level: intermediate
1439027ccd11SLois Curfman McInnes 
1440dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1441d65003e9SLois Curfman McInnes 
1442db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
144320563c6bSBarry Smith @*/
144487828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A)
1445289bc588SBarry Smith {
1446273d9f13SBarry Smith   int ierr;
14473b2fbd54SBarry Smith 
14483a40ed3dSBarry Smith   PetscFunctionBegin;
1449273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1450273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1451273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1452273d9f13SBarry Smith   PetscFunctionReturn(0);
1453273d9f13SBarry Smith }
1454273d9f13SBarry Smith 
14554a2ae208SSatish Balay #undef __FUNCT__
14564a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1457273d9f13SBarry Smith /*@C
1458273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1459273d9f13SBarry Smith 
1460273d9f13SBarry Smith    Collective on MPI_Comm
1461273d9f13SBarry Smith 
1462273d9f13SBarry Smith    Input Parameters:
1463273d9f13SBarry Smith +  A - the matrix
1464273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1465273d9f13SBarry Smith 
1466273d9f13SBarry Smith    Notes:
1467273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1468273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1469273d9f13SBarry Smith    set data=PETSC_NULL.
1470273d9f13SBarry Smith 
1471273d9f13SBarry Smith    Level: intermediate
1472273d9f13SBarry Smith 
1473273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1474273d9f13SBarry Smith 
1475273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1476273d9f13SBarry Smith @*/
147787828ca2SBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar *data)
1478273d9f13SBarry Smith {
1479273d9f13SBarry Smith   Mat_SeqDense *b;
1480273d9f13SBarry Smith   int          ierr;
1481273d9f13SBarry Smith   PetscTruth   flg2;
1482273d9f13SBarry Smith 
1483273d9f13SBarry Smith   PetscFunctionBegin;
1484273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);CHKERRQ(ierr);
1485273d9f13SBarry Smith   if (!flg2) PetscFunctionReturn(0);
1486273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1487273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1488273d9f13SBarry Smith   if (!data) {
148987828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
149087828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1491273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
149287828ca2SBarry Smith     PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));
1493273d9f13SBarry Smith   } else { /* user-allocated storage */
1494273d9f13SBarry Smith     b->v          = data;
1495273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1496273d9f13SBarry Smith   }
1497273d9f13SBarry Smith   PetscFunctionReturn(0);
1498273d9f13SBarry Smith }
1499273d9f13SBarry Smith 
1500273d9f13SBarry Smith EXTERN_C_BEGIN
15014a2ae208SSatish Balay #undef __FUNCT__
15024a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1503273d9f13SBarry Smith int MatCreate_SeqDense(Mat B)
1504273d9f13SBarry Smith {
1505273d9f13SBarry Smith   Mat_SeqDense *b;
1506273d9f13SBarry Smith   int          ierr,size;
1507273d9f13SBarry Smith 
1508273d9f13SBarry Smith   PetscFunctionBegin;
1509273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
151029bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
151155659b69SBarry Smith 
1512273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1513273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1514273d9f13SBarry Smith 
1515b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1516549d3d68SSatish Balay   ierr            = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr);
1517549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
151844cd7ae7SLois Curfman McInnes   B->factor       = 0;
151990f02eecSBarry Smith   B->mapping      = 0;
1520b0a32e0cSBarry Smith   PetscLogObjectMemory(B,sizeof(struct _p_Mat));
152144cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
152218f449edSLois Curfman McInnes 
15238a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
15248a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1525a5ae1ecdSBarry Smith 
152644cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1527273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1528273d9f13SBarry Smith   b->v            = 0;
15294e220ebcSLois Curfman McInnes 
15303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1531289bc588SBarry Smith }
1532273d9f13SBarry Smith EXTERN_C_END
15333b2fbd54SBarry Smith 
15343b2fbd54SBarry Smith 
15353b2fbd54SBarry Smith 
15363b2fbd54SBarry Smith 
1537