xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 5c9eb25f753630cfd361293e05e7358a00a954ac)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
8b4c862e0SSatish Balay #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
170805154bSBarry Smith   PetscBLASInt   N,m,ldax,lday,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
210805154bSBarry Smith   N    = PetscBLASIntCast(X->rmap.n*X->cmap.n);
220805154bSBarry Smith   m    = PetscBLASIntCast(X->rmap.n);
230805154bSBarry Smith   ldax = PetscBLASIntCast(x->lda);
240805154bSBarry Smith   lday = PetscBLASIntCast(y->lda);
25a5ce6ee0Svictorle   if (ldax>m || lday>m) {
26899cda47SBarry Smith     for (j=0; j<X->cmap.n; j++) {
27f4df32b1SMatthew Knepley       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
28a5ce6ee0Svictorle     }
29a5ce6ee0Svictorle   } else {
30f4df32b1SMatthew Knepley     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
31a5ce6ee0Svictorle   }
32efee365bSSatish Balay   ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr);
333a40ed3dSBarry Smith   PetscFunctionReturn(0);
341987afe7SBarry Smith }
351987afe7SBarry Smith 
364a2ae208SSatish Balay #undef __FUNCT__
374a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
38dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
39289bc588SBarry Smith {
40899cda47SBarry Smith   PetscInt     N = A->rmap.n*A->cmap.n;
413a40ed3dSBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
43899cda47SBarry Smith   info->rows_global       = (double)A->rmap.n;
44899cda47SBarry Smith   info->columns_global    = (double)A->cmap.n;
45899cda47SBarry Smith   info->rows_local        = (double)A->rmap.n;
46899cda47SBarry Smith   info->columns_local     = (double)A->cmap.n;
474e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
484e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
496de62eeeSBarry Smith   info->nz_used           = (double)N;
506de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
514e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
524e220ebcSLois Curfman McInnes   info->mallocs           = 0;
537adad957SLisandro Dalcin   info->memory            = ((PetscObject)A)->mem;
544e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
554e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
564e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
573a40ed3dSBarry Smith   PetscFunctionReturn(0);
58289bc588SBarry Smith }
59289bc588SBarry Smith 
604a2ae208SSatish Balay #undef __FUNCT__
614a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
62f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
6380cd9d93SLois Curfman McInnes {
64273d9f13SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
65f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
66efee365bSSatish Balay   PetscErrorCode ierr;
670805154bSBarry Smith   PetscBLASInt   one = 1,j,nz,lda = PetscBLASIntCast(a->lda);
6880cd9d93SLois Curfman McInnes 
693a40ed3dSBarry Smith   PetscFunctionBegin;
70899cda47SBarry Smith   if (lda>A->rmap.n) {
710805154bSBarry Smith     nz = PetscBLASIntCast(A->rmap.n);
72899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
73f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
74a5ce6ee0Svictorle     }
75a5ce6ee0Svictorle   } else {
760805154bSBarry Smith     nz = PetscBLASIntCast(A->rmap.n*A->cmap.n);
77f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
78a5ce6ee0Svictorle   }
79efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
803a40ed3dSBarry Smith   PetscFunctionReturn(0);
8180cd9d93SLois Curfman McInnes }
8280cd9d93SLois Curfman McInnes 
831cbb95d3SBarry Smith #undef __FUNCT__
841cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense"
851cbb95d3SBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscTruth *fl)
861cbb95d3SBarry Smith {
871cbb95d3SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
881cbb95d3SBarry Smith   PetscInt       i,j,m = A->rmap.n,N;
891cbb95d3SBarry Smith   PetscScalar    *v = a->v;
901cbb95d3SBarry Smith 
911cbb95d3SBarry Smith   PetscFunctionBegin;
921cbb95d3SBarry Smith   *fl = PETSC_FALSE;
931cbb95d3SBarry Smith   if (A->rmap.n != A->cmap.n) PetscFunctionReturn(0);
941cbb95d3SBarry Smith   N = a->lda;
951cbb95d3SBarry Smith 
961cbb95d3SBarry Smith   for (i=0; i<m; i++) {
971cbb95d3SBarry Smith     for (j=i+1; j<m; j++) {
981cbb95d3SBarry Smith       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
991cbb95d3SBarry Smith     }
1001cbb95d3SBarry Smith   }
1011cbb95d3SBarry Smith   *fl = PETSC_TRUE;
1021cbb95d3SBarry Smith   PetscFunctionReturn(0);
1031cbb95d3SBarry Smith }
1041cbb95d3SBarry Smith 
105289bc588SBarry Smith /* ---------------------------------------------------------------*/
1066831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
1076831982aSBarry Smith    rather than put it in the Mat->row slot.*/
1084a2ae208SSatish Balay #undef __FUNCT__
1094a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
110dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
111289bc588SBarry Smith {
112a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
11381f479a6Svictorle   PetscFunctionBegin;
114a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
115a07ab388SBarry Smith #else
116c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1176849ba73SBarry Smith   PetscErrorCode ierr;
1180805154bSBarry Smith   PetscBLASInt   n,m,info;
11955659b69SBarry Smith 
1203a40ed3dSBarry Smith   PetscFunctionBegin;
1210805154bSBarry Smith   n = PetscBLASIntCast(A->cmap.n);
1220805154bSBarry Smith   m = PetscBLASIntCast(A->rmap.n);
123289bc588SBarry Smith   if (!mat->pivots) {
124899cda47SBarry Smith     ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
125899cda47SBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap.n*sizeof(PetscBLASInt));CHKERRQ(ierr);
126289bc588SBarry Smith   }
127*5c9eb25fSBarry Smith   A->factor = MAT_FACTOR_LU;
128899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
12971044d3cSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
13029bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
13129bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
132899cda47SBarry Smith   ierr = PetscLogFlops((2*A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
133a07ab388SBarry Smith #endif
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
135289bc588SBarry Smith }
1366ee01492SSatish Balay 
137b24902e0SBarry Smith 
138b24902e0SBarry Smith #undef __FUNCT__
139b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense"
140b24902e0SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
141b24902e0SBarry Smith {
142b24902e0SBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
143b24902e0SBarry Smith   PetscErrorCode ierr;
144b24902e0SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
145b24902e0SBarry Smith   Mat            newi = *newmat;
146b24902e0SBarry Smith 
147b24902e0SBarry Smith   PetscFunctionBegin;
148*5c9eb25fSBarry Smith   ierr = MatSeqDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
149b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
150b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
151b24902e0SBarry Smith     if (lda>A->rmap.n) {
152b24902e0SBarry Smith       m = A->rmap.n;
153b24902e0SBarry Smith       for (j=0; j<A->cmap.n; j++) {
154b24902e0SBarry Smith 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
155b24902e0SBarry Smith       }
156b24902e0SBarry Smith     } else {
157b24902e0SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
158b24902e0SBarry Smith     }
159b24902e0SBarry Smith   }
160b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
161b24902e0SBarry Smith   PetscFunctionReturn(0);
162b24902e0SBarry Smith }
163b24902e0SBarry Smith 
1644a2ae208SSatish Balay #undef __FUNCT__
1654a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
166dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
16702cad45dSBarry Smith {
16802cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
1696849ba73SBarry Smith   PetscErrorCode ierr;
17013f74950SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
17102cad45dSBarry Smith 
1723a40ed3dSBarry Smith   PetscFunctionBegin;
173*5c9eb25fSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
174*5c9eb25fSBarry Smith   ierr = MatSetSizes(*newmat,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
175*5c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
176b24902e0SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(A,cpvalues,newmat);CHKERRQ(ierr);
177b24902e0SBarry Smith   PetscFunctionReturn(0);
178b24902e0SBarry Smith }
179b24902e0SBarry Smith 
180b24902e0SBarry Smith #undef __FUNCT__
181b24902e0SBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
182b24902e0SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
183b24902e0SBarry Smith {
184b24902e0SBarry Smith   PetscErrorCode ierr;
185b24902e0SBarry Smith 
186b24902e0SBarry Smith   PetscFunctionBegin;
187*5c9eb25fSBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(A,MAT_COPY_VALUES,fact);CHKERRQ(ierr);
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
18902cad45dSBarry Smith }
19002cad45dSBarry Smith 
1914a2ae208SSatish Balay #undef __FUNCT__
1924a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
193dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
194289bc588SBarry Smith {
195dfbe8321SBarry Smith   PetscErrorCode ierr;
1963a40ed3dSBarry Smith 
1973a40ed3dSBarry Smith   PetscFunctionBegin;
198b24902e0SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
199b24902e0SBarry Smith   PetscFunctionReturn(0);
200b24902e0SBarry Smith }
201b24902e0SBarry Smith 
202b24902e0SBarry Smith #undef __FUNCT__
203b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc"
204*5c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
205b24902e0SBarry Smith {
206b24902e0SBarry Smith   PetscErrorCode ierr;
207b24902e0SBarry Smith 
208b24902e0SBarry Smith   PetscFunctionBegin;
209b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
210b24902e0SBarry Smith   ierr = MatSetSizes(*fact,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
211b24902e0SBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
2123a40ed3dSBarry Smith   PetscFunctionReturn(0);
213289bc588SBarry Smith }
2146ee01492SSatish Balay 
2154a2ae208SSatish Balay #undef __FUNCT__
2164a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
217af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
218289bc588SBarry Smith {
21902cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
2206849ba73SBarry Smith   PetscErrorCode ierr;
221899cda47SBarry Smith   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->rmap.n,n=A->cmap.n, j;
2224482741eSBarry Smith   MatFactorInfo  info;
2233a40ed3dSBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionBegin;
22502cad45dSBarry Smith   /* copy the numerical values */
2261b807ce4Svictorle   if (lda1>m || lda2>m ) {
2271b807ce4Svictorle     for (j=0; j<n; j++) {
2281b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
2291b807ce4Svictorle     }
2301b807ce4Svictorle   } else {
231899cda47SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
2321b807ce4Svictorle   }
233*5c9eb25fSBarry Smith   (*fact)->factor = MAT_FACTOR_LU;
2344482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
2353a40ed3dSBarry Smith   PetscFunctionReturn(0);
236289bc588SBarry Smith }
2376ee01492SSatish Balay 
2384a2ae208SSatish Balay #undef __FUNCT__
2394a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
240dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
241289bc588SBarry Smith {
242a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
243a07ab388SBarry Smith   PetscFunctionBegin;
244a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
245a07ab388SBarry Smith #else
246c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2476849ba73SBarry Smith   PetscErrorCode ierr;
2480805154bSBarry Smith   PetscBLASInt   info,n = PetscBLASIntCast(A->cmap.n);
24955659b69SBarry Smith 
2503a40ed3dSBarry Smith   PetscFunctionBegin;
251606d414cSSatish Balay   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
252752f5567SLois Curfman McInnes   mat->pivots = 0;
25305b42c5fSBarry Smith 
254899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
25571044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
25677431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
257*5c9eb25fSBarry Smith   A->factor = MAT_FACTOR_CHOLESKY;
258899cda47SBarry Smith   ierr = PetscLogFlops((A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
259a07ab388SBarry Smith #endif
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
261289bc588SBarry Smith }
262289bc588SBarry Smith 
2634a2ae208SSatish Balay #undef __FUNCT__
2640b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
265af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2660b4b3355SBarry Smith {
267dfbe8321SBarry Smith   PetscErrorCode ierr;
26815e8a5b3SHong Zhang   MatFactorInfo  info;
2690b4b3355SBarry Smith 
2700b4b3355SBarry Smith   PetscFunctionBegin;
27115e8a5b3SHong Zhang   info.fill = 1.0;
27215e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2730b4b3355SBarry Smith   PetscFunctionReturn(0);
2740b4b3355SBarry Smith }
2750b4b3355SBarry Smith 
2760b4b3355SBarry Smith #undef __FUNCT__
2774a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
278dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
279289bc588SBarry Smith {
280c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2816849ba73SBarry Smith   PetscErrorCode ierr;
28287828ca2SBarry Smith   PetscScalar    *x,*y;
2830805154bSBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap.n);
28467e560aaSBarry Smith 
2853a40ed3dSBarry Smith   PetscFunctionBegin;
2861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2871ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
288899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
289*5c9eb25fSBarry Smith   if (A->factor == MAT_FACTOR_LU) {
290ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
29180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
292ae7cfcebSSatish Balay #else
29371044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2944ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
295ae7cfcebSSatish Balay #endif
296*5c9eb25fSBarry Smith   } else if (A->factor == MAT_FACTOR_CHOLESKY){
297ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
29880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
299ae7cfcebSSatish Balay #else
30071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3014ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
302ae7cfcebSSatish Balay #endif
303289bc588SBarry Smith   }
30429bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
3051ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3061ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
307899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
309289bc588SBarry Smith }
3106ee01492SSatish Balay 
3114a2ae208SSatish Balay #undef __FUNCT__
3124a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
313dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
314da3a660dSBarry Smith {
315c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
316dfbe8321SBarry Smith   PetscErrorCode ierr;
31787828ca2SBarry Smith   PetscScalar    *x,*y;
3180805154bSBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap.n);
31967e560aaSBarry Smith 
3203a40ed3dSBarry Smith   PetscFunctionBegin;
3211ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3221ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
323899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
324752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
325da3a660dSBarry Smith   if (mat->pivots) {
326ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
32780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
328ae7cfcebSSatish Balay #else
32971044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3304ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
331ae7cfcebSSatish Balay #endif
3327a97a34bSBarry Smith   } else {
333ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
33480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
335ae7cfcebSSatish Balay #else
33671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3374ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
338ae7cfcebSSatish Balay #endif
339da3a660dSBarry Smith   }
3401ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3411ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
342899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
3433a40ed3dSBarry Smith   PetscFunctionReturn(0);
344da3a660dSBarry Smith }
3456ee01492SSatish Balay 
3464a2ae208SSatish Balay #undef __FUNCT__
3474a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
348dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
349da3a660dSBarry Smith {
350c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
351dfbe8321SBarry Smith   PetscErrorCode ierr;
35287828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
353da3a660dSBarry Smith   Vec            tmp = 0;
3540805154bSBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap.n);
35567e560aaSBarry Smith 
3563a40ed3dSBarry Smith   PetscFunctionBegin;
3571ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3581ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
359899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
360da3a660dSBarry Smith   if (yy == zz) {
36178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
36252e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
36378b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
364da3a660dSBarry Smith   }
365899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
366752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
367da3a660dSBarry Smith   if (mat->pivots) {
368ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
36980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
370ae7cfcebSSatish Balay #else
37171044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3722ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
373ae7cfcebSSatish Balay #endif
374a8c6a408SBarry Smith   } else {
375ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
37680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
377ae7cfcebSSatish Balay #else
37871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3792ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
380ae7cfcebSSatish Balay #endif
381da3a660dSBarry Smith   }
3822dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
3832dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
3841ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3851ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
386899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
388da3a660dSBarry Smith }
38967e560aaSBarry Smith 
3904a2ae208SSatish Balay #undef __FUNCT__
3914a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
392dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
393da3a660dSBarry Smith {
394c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3956849ba73SBarry Smith   PetscErrorCode ierr;
39687828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
397da3a660dSBarry Smith   Vec            tmp;
3980805154bSBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap.n);
39967e560aaSBarry Smith 
4003a40ed3dSBarry Smith   PetscFunctionBegin;
401899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
4021ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4031ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
404da3a660dSBarry Smith   if (yy == zz) {
40578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
40652e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
40778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
408da3a660dSBarry Smith   }
409899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
410752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
411da3a660dSBarry Smith   if (mat->pivots) {
412ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
41380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
414ae7cfcebSSatish Balay #else
41571044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
4162ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
417ae7cfcebSSatish Balay #endif
4183a40ed3dSBarry Smith   } else {
419ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
42080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
421ae7cfcebSSatish Balay #else
42271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
4232ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
424ae7cfcebSSatish Balay #endif
425da3a660dSBarry Smith   }
42690f02eecSBarry Smith   if (tmp) {
4272dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
42890f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
4293a40ed3dSBarry Smith   } else {
4302dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
43190f02eecSBarry Smith   }
4321ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4331ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
434899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
4353a40ed3dSBarry Smith   PetscFunctionReturn(0);
436da3a660dSBarry Smith }
437289bc588SBarry Smith /* ------------------------------------------------------------------*/
4384a2ae208SSatish Balay #undef __FUNCT__
4394a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
44013f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
441289bc588SBarry Smith {
442c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
44387828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
444dfbe8321SBarry Smith   PetscErrorCode ierr;
445899cda47SBarry Smith   PetscInt       m = A->rmap.n,i;
446aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4470805154bSBarry Smith   PetscBLASInt   o = 1,bm = PetscBLASIntCast(m);
448bc1b551cSSatish Balay #endif
449289bc588SBarry Smith 
4503a40ed3dSBarry Smith   PetscFunctionBegin;
451289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
45271044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4532dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
454289bc588SBarry Smith   }
4551ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4561ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
457b965ef7fSBarry Smith   its  = its*lits;
45877431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
459289bc588SBarry Smith   while (its--) {
460fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
461289bc588SBarry Smith       for (i=0; i<m; i++) {
462aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
463f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
464f1747703SBarry Smith            not happy about returning a double complex */
46513f74950SBarry Smith         PetscInt    _i;
46687828ca2SBarry Smith         PetscScalar sum = b[i];
467f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4683f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
469f1747703SBarry Smith         }
470f1747703SBarry Smith         xt = sum;
471f1747703SBarry Smith #else
47271044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
473f1747703SBarry Smith #endif
47455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
475289bc588SBarry Smith       }
476289bc588SBarry Smith     }
477fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
478289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
479aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
480f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
481f1747703SBarry Smith            not happy about returning a double complex */
48213f74950SBarry Smith         PetscInt    _i;
48387828ca2SBarry Smith         PetscScalar sum = b[i];
484f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4853f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
486f1747703SBarry Smith         }
487f1747703SBarry Smith         xt = sum;
488f1747703SBarry Smith #else
48971044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
490f1747703SBarry Smith #endif
49155a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
492289bc588SBarry Smith       }
493289bc588SBarry Smith     }
494289bc588SBarry Smith   }
4951ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4961ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4973a40ed3dSBarry Smith   PetscFunctionReturn(0);
498289bc588SBarry Smith }
499289bc588SBarry Smith 
500289bc588SBarry Smith /* -----------------------------------------------------------------*/
5014a2ae208SSatish Balay #undef __FUNCT__
5024a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
503dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
504289bc588SBarry Smith {
505c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
50687828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
507dfbe8321SBarry Smith   PetscErrorCode ierr;
5080805154bSBarry Smith   PetscBLASInt   m, n,_One=1;
509ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
5103a40ed3dSBarry Smith 
5113a40ed3dSBarry Smith   PetscFunctionBegin;
5120805154bSBarry Smith   m = PetscBLASIntCast(A->rmap.n);
5130805154bSBarry Smith   n = PetscBLASIntCast(A->cmap.n);
514899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
5151ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5161ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
51771044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
5181ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5191ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
520899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
5213a40ed3dSBarry Smith   PetscFunctionReturn(0);
522289bc588SBarry Smith }
5236ee01492SSatish Balay 
5244a2ae208SSatish Balay #undef __FUNCT__
5254a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
526dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
527289bc588SBarry Smith {
528c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
52987828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
530dfbe8321SBarry Smith   PetscErrorCode ierr;
5310805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5323a40ed3dSBarry Smith 
5333a40ed3dSBarry Smith   PetscFunctionBegin;
5340805154bSBarry Smith   m = PetscBLASIntCast(A->rmap.n);
5350805154bSBarry Smith   n = PetscBLASIntCast(A->cmap.n);
536899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
5371ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5381ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
53971044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5401ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5411ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
542899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->rmap.n);CHKERRQ(ierr);
5433a40ed3dSBarry Smith   PetscFunctionReturn(0);
544289bc588SBarry Smith }
5456ee01492SSatish Balay 
5464a2ae208SSatish Balay #undef __FUNCT__
5474a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
548dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
549289bc588SBarry Smith {
550c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
55187828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
552dfbe8321SBarry Smith   PetscErrorCode ierr;
5530805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5543a40ed3dSBarry Smith 
5553a40ed3dSBarry Smith   PetscFunctionBegin;
5560805154bSBarry Smith   m = PetscBLASIntCast(A->rmap.n);
5570805154bSBarry Smith   n = PetscBLASIntCast(A->cmap.n);
558899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
559600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5601ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5611ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
56271044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5631ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5641ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
565899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5663a40ed3dSBarry Smith   PetscFunctionReturn(0);
567289bc588SBarry Smith }
5686ee01492SSatish Balay 
5694a2ae208SSatish Balay #undef __FUNCT__
5704a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
571dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
572289bc588SBarry Smith {
573c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
57487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
575dfbe8321SBarry Smith   PetscErrorCode ierr;
5760805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
57787828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5783a40ed3dSBarry Smith 
5793a40ed3dSBarry Smith   PetscFunctionBegin;
5800805154bSBarry Smith   m = PetscBLASIntCast(A->rmap.n);
5810805154bSBarry Smith   n = PetscBLASIntCast(A->cmap.n);
582899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
583600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5841ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5851ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
58671044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5871ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5881ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
589899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5903a40ed3dSBarry Smith   PetscFunctionReturn(0);
591289bc588SBarry Smith }
592289bc588SBarry Smith 
593289bc588SBarry Smith /* -----------------------------------------------------------------*/
5944a2ae208SSatish Balay #undef __FUNCT__
5954a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
59613f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
597289bc588SBarry Smith {
598c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
59987828ca2SBarry Smith   PetscScalar    *v;
6006849ba73SBarry Smith   PetscErrorCode ierr;
60113f74950SBarry Smith   PetscInt       i;
60267e560aaSBarry Smith 
6033a40ed3dSBarry Smith   PetscFunctionBegin;
604899cda47SBarry Smith   *ncols = A->cmap.n;
605289bc588SBarry Smith   if (cols) {
606899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
607899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) (*cols)[i] = i;
608289bc588SBarry Smith   }
609289bc588SBarry Smith   if (vals) {
610899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
611289bc588SBarry Smith     v    = mat->v + row;
612899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) {(*vals)[i] = *v; v += mat->lda;}
613289bc588SBarry Smith   }
6143a40ed3dSBarry Smith   PetscFunctionReturn(0);
615289bc588SBarry Smith }
6166ee01492SSatish Balay 
6174a2ae208SSatish Balay #undef __FUNCT__
6184a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
61913f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
620289bc588SBarry Smith {
621dfbe8321SBarry Smith   PetscErrorCode ierr;
622606d414cSSatish Balay   PetscFunctionBegin;
623606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
624606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6253a40ed3dSBarry Smith   PetscFunctionReturn(0);
626289bc588SBarry Smith }
627289bc588SBarry Smith /* ----------------------------------------------------------------*/
6284a2ae208SSatish Balay #undef __FUNCT__
6294a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
63013f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
631289bc588SBarry Smith {
632c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
63313f74950SBarry Smith   PetscInt     i,j,idx=0;
634d6dfbf8fSBarry Smith 
6353a40ed3dSBarry Smith   PetscFunctionBegin;
636289bc588SBarry Smith   if (!mat->roworiented) {
637dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
638289bc588SBarry Smith       for (j=0; j<n; j++) {
639cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6402515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
641899cda47SBarry Smith         if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
64258804f6dSBarry Smith #endif
643289bc588SBarry Smith         for (i=0; i<m; i++) {
644cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6452515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
646899cda47SBarry Smith           if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
64758804f6dSBarry Smith #endif
648cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
649289bc588SBarry Smith         }
650289bc588SBarry Smith       }
6513a40ed3dSBarry Smith     } else {
652289bc588SBarry Smith       for (j=0; j<n; j++) {
653cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6542515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
655899cda47SBarry Smith         if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
65658804f6dSBarry Smith #endif
657289bc588SBarry Smith         for (i=0; i<m; i++) {
658cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6592515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
660899cda47SBarry Smith           if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
66158804f6dSBarry Smith #endif
662cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
663289bc588SBarry Smith         }
664289bc588SBarry Smith       }
665289bc588SBarry Smith     }
6663a40ed3dSBarry Smith   } else {
667dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
668e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
669cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6702515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
671899cda47SBarry Smith         if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
67258804f6dSBarry Smith #endif
673e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
674cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6752515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
676899cda47SBarry Smith           if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
67758804f6dSBarry Smith #endif
678cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
679e8d4e0b9SBarry Smith         }
680e8d4e0b9SBarry Smith       }
6813a40ed3dSBarry Smith     } else {
682289bc588SBarry Smith       for (i=0; i<m; i++) {
683cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6842515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
685899cda47SBarry Smith         if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1);
68658804f6dSBarry Smith #endif
687289bc588SBarry Smith         for (j=0; j<n; j++) {
688cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6892515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
690899cda47SBarry Smith           if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1);
69158804f6dSBarry Smith #endif
692cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
693289bc588SBarry Smith         }
694289bc588SBarry Smith       }
695289bc588SBarry Smith     }
696e8d4e0b9SBarry Smith   }
6973a40ed3dSBarry Smith   PetscFunctionReturn(0);
698289bc588SBarry Smith }
699e8d4e0b9SBarry Smith 
7004a2ae208SSatish Balay #undef __FUNCT__
7014a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
70213f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
703ae80bb75SLois Curfman McInnes {
704ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
70513f74950SBarry Smith   PetscInt     i,j;
706ae80bb75SLois Curfman McInnes 
7073a40ed3dSBarry Smith   PetscFunctionBegin;
708ae80bb75SLois Curfman McInnes   /* row-oriented output */
709ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
71097e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
711f12f8aa1SBarry Smith     if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap.n);
712ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7136f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
7146f31f424SBarry Smith       if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap.n);
71597e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
716ae80bb75SLois Curfman McInnes     }
717ae80bb75SLois Curfman McInnes   }
7183a40ed3dSBarry Smith   PetscFunctionReturn(0);
719ae80bb75SLois Curfman McInnes }
720ae80bb75SLois Curfman McInnes 
721289bc588SBarry Smith /* -----------------------------------------------------------------*/
722289bc588SBarry Smith 
723e090d566SSatish Balay #include "petscsys.h"
72488e32edaSLois Curfman McInnes 
7254a2ae208SSatish Balay #undef __FUNCT__
7264a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
727f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A)
72888e32edaSLois Curfman McInnes {
72988e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
73088e32edaSLois Curfman McInnes   Mat            B;
7316849ba73SBarry Smith   PetscErrorCode ierr;
73213f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
73313f74950SBarry Smith   int            fd;
73413f74950SBarry Smith   PetscMPIInt    size;
73513f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
73687828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
73719bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
73888e32edaSLois Curfman McInnes 
7393a40ed3dSBarry Smith   PetscFunctionBegin;
740d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
74129bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
742b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
7430752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
744552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
74588e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
74688e32edaSLois Curfman McInnes 
74790ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
748f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
749f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
750be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7515c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
75290ace30eSBarry Smith     B    = *A;
75390ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
7544a41a4d3SSatish Balay     v    = a->v;
7554a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7564a41a4d3SSatish Balay        from row major to column major */
757b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
75890ace30eSBarry Smith     /* read in nonzero values */
7594a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7604a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7614a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7624a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7634a41a4d3SSatish Balay         *v++ =w[i*N+j];
7644a41a4d3SSatish Balay       }
7654a41a4d3SSatish Balay     }
766606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7676d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7686d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76990ace30eSBarry Smith   } else {
77088e32edaSLois Curfman McInnes     /* read row lengths */
77113f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7720752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
77388e32edaSLois Curfman McInnes 
77488e32edaSLois Curfman McInnes     /* create our matrix */
775f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
776f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
777be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7785c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
77988e32edaSLois Curfman McInnes     B = *A;
78088e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
78188e32edaSLois Curfman McInnes     v = a->v;
78288e32edaSLois Curfman McInnes 
78388e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
78413f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
785b0a32e0cSBarry Smith     cols = scols;
7860752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
78787828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
788b0a32e0cSBarry Smith     vals = svals;
7890752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
79088e32edaSLois Curfman McInnes 
79188e32edaSLois Curfman McInnes     /* insert into matrix */
79288e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
79388e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
79488e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
79588e32edaSLois Curfman McInnes     }
796606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
797606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
798606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
79988e32edaSLois Curfman McInnes 
8006d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8016d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80290ace30eSBarry Smith   }
8033a40ed3dSBarry Smith   PetscFunctionReturn(0);
80488e32edaSLois Curfman McInnes }
80588e32edaSLois Curfman McInnes 
806e090d566SSatish Balay #include "petscsys.h"
807932b0c3eSLois Curfman McInnes 
8084a2ae208SSatish Balay #undef __FUNCT__
8094a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8106849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
811289bc588SBarry Smith {
812932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
813dfbe8321SBarry Smith   PetscErrorCode    ierr;
81413f74950SBarry Smith   PetscInt          i,j;
8152dcb1b2aSMatthew Knepley   const char        *name;
81687828ca2SBarry Smith   PetscScalar       *v;
817f3ef73ceSBarry Smith   PetscViewerFormat format;
8185f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
8195f481a85SSatish Balay   PetscTruth        allreal = PETSC_TRUE;
8205f481a85SSatish Balay #endif
821932b0c3eSLois Curfman McInnes 
8223a40ed3dSBarry Smith   PetscFunctionBegin;
823435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
824b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
825456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8263a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
827fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
828b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
829899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
83044cd7ae7SLois Curfman McInnes       v = a->v + i;
83177431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
832899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
833aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
834329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
835a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
836329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
837a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8386831982aSBarry Smith         }
83980cd9d93SLois Curfman McInnes #else
8406831982aSBarry Smith         if (*v) {
841a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8426831982aSBarry Smith         }
84380cd9d93SLois Curfman McInnes #endif
8441b807ce4Svictorle         v += a->lda;
84580cd9d93SLois Curfman McInnes       }
846b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
84780cd9d93SLois Curfman McInnes     }
848b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
8493a40ed3dSBarry Smith   } else {
850b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
851aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
85247989497SBarry Smith     /* determine if matrix has all real values */
85347989497SBarry Smith     v = a->v;
854899cda47SBarry Smith     for (i=0; i<A->rmap.n*A->cmap.n; i++) {
855ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
85647989497SBarry Smith     }
85747989497SBarry Smith #endif
858fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8593a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
860899cda47SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr);
861899cda47SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
862fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
863ffac6cdbSBarry Smith     }
864ffac6cdbSBarry Smith 
865899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
866932b0c3eSLois Curfman McInnes       v = a->v + i;
867899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
868aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
86947989497SBarry Smith         if (allreal) {
870f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
87147989497SBarry Smith         } else {
872f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
87347989497SBarry Smith         }
874289bc588SBarry Smith #else
875f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
876289bc588SBarry Smith #endif
8771b807ce4Svictorle 	v += a->lda;
878289bc588SBarry Smith       }
879b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
880289bc588SBarry Smith     }
881fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
882b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
883ffac6cdbSBarry Smith     }
884b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
885da3a660dSBarry Smith   }
886b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8873a40ed3dSBarry Smith   PetscFunctionReturn(0);
888289bc588SBarry Smith }
889289bc588SBarry Smith 
8904a2ae208SSatish Balay #undef __FUNCT__
8914a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8926849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
893932b0c3eSLois Curfman McInnes {
894932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8956849ba73SBarry Smith   PetscErrorCode    ierr;
89613f74950SBarry Smith   int               fd;
897899cda47SBarry Smith   PetscInt          ict,j,n = A->cmap.n,m = A->rmap.n,i,*col_lens,nz = m*n;
89887828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
899f3ef73ceSBarry Smith   PetscViewerFormat format;
900932b0c3eSLois Curfman McInnes 
9013a40ed3dSBarry Smith   PetscFunctionBegin;
902b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
90390ace30eSBarry Smith 
904b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
905fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
90690ace30eSBarry Smith     /* store the matrix as a dense matrix */
90713f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
908552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
90990ace30eSBarry Smith     col_lens[1] = m;
91090ace30eSBarry Smith     col_lens[2] = n;
91190ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
9126f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
913606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
91490ace30eSBarry Smith 
91590ace30eSBarry Smith     /* write out matrix, by rows */
91687828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
91790ace30eSBarry Smith     v    = a->v;
91890ace30eSBarry Smith     for (j=0; j<n; j++) {
919578230a0SSatish Balay       for (i=0; i<m; i++) {
920578230a0SSatish Balay         vals[j + i*n] = *v++;
92190ace30eSBarry Smith       }
92290ace30eSBarry Smith     }
9236f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
924606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
92590ace30eSBarry Smith   } else {
92613f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
927552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
928932b0c3eSLois Curfman McInnes     col_lens[1] = m;
929932b0c3eSLois Curfman McInnes     col_lens[2] = n;
930932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
931932b0c3eSLois Curfman McInnes 
932932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
933932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9346f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
935932b0c3eSLois Curfman McInnes 
936932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
937932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
938932b0c3eSLois Curfman McInnes     ict = 0;
939932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
940932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
941932b0c3eSLois Curfman McInnes     }
9426f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
943606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
944932b0c3eSLois Curfman McInnes 
945932b0c3eSLois Curfman McInnes     /* store nonzero values */
94687828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
947932b0c3eSLois Curfman McInnes     ict  = 0;
948932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
949932b0c3eSLois Curfman McInnes       v = a->v + i;
950932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9511b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
952932b0c3eSLois Curfman McInnes       }
953932b0c3eSLois Curfman McInnes     }
9546f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
955606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
95690ace30eSBarry Smith   }
9573a40ed3dSBarry Smith   PetscFunctionReturn(0);
958932b0c3eSLois Curfman McInnes }
959932b0c3eSLois Curfman McInnes 
9604a2ae208SSatish Balay #undef __FUNCT__
9614a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
962dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
963f1af5d2fSBarry Smith {
964f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
965f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9666849ba73SBarry Smith   PetscErrorCode    ierr;
967899cda47SBarry Smith   PetscInt          m = A->rmap.n,n = A->cmap.n,color,i,j;
96887828ca2SBarry Smith   PetscScalar       *v = a->v;
969b0a32e0cSBarry Smith   PetscViewer       viewer;
970b0a32e0cSBarry Smith   PetscDraw         popup;
971329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
972f3ef73ceSBarry Smith   PetscViewerFormat format;
973f1af5d2fSBarry Smith 
974f1af5d2fSBarry Smith   PetscFunctionBegin;
975f1af5d2fSBarry Smith 
976f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
977b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
978b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
979f1af5d2fSBarry Smith 
980f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
981fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
982f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
983b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
984f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
985f1af5d2fSBarry Smith       x_l = j;
986f1af5d2fSBarry Smith       x_r = x_l + 1.0;
987f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
988f1af5d2fSBarry Smith         y_l = m - i - 1.0;
989f1af5d2fSBarry Smith         y_r = y_l + 1.0;
990f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
991329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
992b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
993329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
994b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
995f1af5d2fSBarry Smith         } else {
996f1af5d2fSBarry Smith           continue;
997f1af5d2fSBarry Smith         }
998f1af5d2fSBarry Smith #else
999f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
1000b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1001f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1002b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1003f1af5d2fSBarry Smith         } else {
1004f1af5d2fSBarry Smith           continue;
1005f1af5d2fSBarry Smith         }
1006f1af5d2fSBarry Smith #endif
1007b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1008f1af5d2fSBarry Smith       }
1009f1af5d2fSBarry Smith     }
1010f1af5d2fSBarry Smith   } else {
1011f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1012f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1013f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1014f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1015f1af5d2fSBarry Smith     }
1016b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1017b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1018b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1019f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1020f1af5d2fSBarry Smith       x_l = j;
1021f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1022f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1023f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1024f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1025b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1026b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1027f1af5d2fSBarry Smith       }
1028f1af5d2fSBarry Smith     }
1029f1af5d2fSBarry Smith   }
1030f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1031f1af5d2fSBarry Smith }
1032f1af5d2fSBarry Smith 
10334a2ae208SSatish Balay #undef __FUNCT__
10344a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1035dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1036f1af5d2fSBarry Smith {
1037b0a32e0cSBarry Smith   PetscDraw      draw;
1038f1af5d2fSBarry Smith   PetscTruth     isnull;
1039329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1040dfbe8321SBarry Smith   PetscErrorCode ierr;
1041f1af5d2fSBarry Smith 
1042f1af5d2fSBarry Smith   PetscFunctionBegin;
1043b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1044b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1045abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1046f1af5d2fSBarry Smith 
1047f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1048899cda47SBarry Smith   xr  = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
1049f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1050b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1051b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1052f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1053f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1054f1af5d2fSBarry Smith }
1055f1af5d2fSBarry Smith 
10564a2ae208SSatish Balay #undef __FUNCT__
10574a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1058dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1059932b0c3eSLois Curfman McInnes {
1060dfbe8321SBarry Smith   PetscErrorCode ierr;
10616805f65bSBarry Smith   PetscTruth     iascii,isbinary,isdraw;
1062932b0c3eSLois Curfman McInnes 
10633a40ed3dSBarry Smith   PetscFunctionBegin;
106432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1065fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1066fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10670f5bd95cSBarry Smith 
1068c45a1595SBarry Smith   if (iascii) {
1069c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10700f5bd95cSBarry Smith   } else if (isbinary) {
10713a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1072f1af5d2fSBarry Smith   } else if (isdraw) {
1073f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10745cd90555SBarry Smith   } else {
1075958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1076932b0c3eSLois Curfman McInnes   }
10773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1078932b0c3eSLois Curfman McInnes }
1079289bc588SBarry Smith 
10804a2ae208SSatish Balay #undef __FUNCT__
10814a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1082dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1083289bc588SBarry Smith {
1084ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1085dfbe8321SBarry Smith   PetscErrorCode ierr;
108690f02eecSBarry Smith 
10873a40ed3dSBarry Smith   PetscFunctionBegin;
1088aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1089899cda47SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n);
1090a5a9c739SBarry Smith #endif
109105b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
10926857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1093606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1094dbd8c25aSHong Zhang 
1095dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1096901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10974ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10984ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10994ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1101289bc588SBarry Smith }
1102289bc588SBarry Smith 
11034a2ae208SSatish Balay #undef __FUNCT__
11044a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1105fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1106289bc588SBarry Smith {
1107c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11086849ba73SBarry Smith   PetscErrorCode ierr;
110913f74950SBarry Smith   PetscInt       k,j,m,n,M;
111087828ca2SBarry Smith   PetscScalar    *v,tmp;
111148b35521SBarry Smith 
11123a40ed3dSBarry Smith   PetscFunctionBegin;
1113899cda47SBarry Smith   v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n;
1114e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1115a5ce6ee0Svictorle     if (m != n) {
1116634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1117d64ed03dSBarry Smith     } else {
1118d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1119289bc588SBarry Smith         for (k=0; k<j; k++) {
11201b807ce4Svictorle           tmp = v[j + k*M];
11211b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11221b807ce4Svictorle           v[k + j*M] = tmp;
1123289bc588SBarry Smith         }
1124289bc588SBarry Smith       }
1125d64ed03dSBarry Smith     }
11263a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1127d3e5ee88SLois Curfman McInnes     Mat          tmat;
1128ec8511deSBarry Smith     Mat_SeqDense *tmatd;
112987828ca2SBarry Smith     PetscScalar  *v2;
1130ea709b57SSatish Balay 
1131fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11327adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1133899cda47SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr);
11347adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11355c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1136fc4dec0aSBarry Smith     } else {
1137fc4dec0aSBarry Smith       tmat = *matout;
1138fc4dec0aSBarry Smith     }
1139ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11400de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1141d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11421b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1143d3e5ee88SLois Curfman McInnes     }
11446d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11456d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1146d3e5ee88SLois Curfman McInnes     *matout = tmat;
114748b35521SBarry Smith   }
11483a40ed3dSBarry Smith   PetscFunctionReturn(0);
1149289bc588SBarry Smith }
1150289bc588SBarry Smith 
11514a2ae208SSatish Balay #undef __FUNCT__
11524a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1153dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1154289bc588SBarry Smith {
1155c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1156c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
115713f74950SBarry Smith   PetscInt     i,j;
115887828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11599ea5d5aeSSatish Balay 
11603a40ed3dSBarry Smith   PetscFunctionBegin;
1161899cda47SBarry Smith   if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1162899cda47SBarry Smith   if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1163899cda47SBarry Smith   for (i=0; i<A1->rmap.n; i++) {
11641b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1165899cda47SBarry Smith     for (j=0; j<A1->cmap.n; j++) {
11663a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11671b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11681b807ce4Svictorle     }
1169289bc588SBarry Smith   }
117077c4ece6SBarry Smith   *flg = PETSC_TRUE;
11713a40ed3dSBarry Smith   PetscFunctionReturn(0);
1172289bc588SBarry Smith }
1173289bc588SBarry Smith 
11744a2ae208SSatish Balay #undef __FUNCT__
11754a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1176dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1177289bc588SBarry Smith {
1178c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1179dfbe8321SBarry Smith   PetscErrorCode ierr;
118013f74950SBarry Smith   PetscInt       i,n,len;
118187828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
118244cd7ae7SLois Curfman McInnes 
11833a40ed3dSBarry Smith   PetscFunctionBegin;
11842dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11857a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11861ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1187899cda47SBarry Smith   len = PetscMin(A->rmap.n,A->cmap.n);
1188899cda47SBarry Smith   if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
118944cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11901b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1191289bc588SBarry Smith   }
11921ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11933a40ed3dSBarry Smith   PetscFunctionReturn(0);
1194289bc588SBarry Smith }
1195289bc588SBarry Smith 
11964a2ae208SSatish Balay #undef __FUNCT__
11974a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1198dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1199289bc588SBarry Smith {
1200c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
120187828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1202dfbe8321SBarry Smith   PetscErrorCode ierr;
1203899cda47SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n;
120455659b69SBarry Smith 
12053a40ed3dSBarry Smith   PetscFunctionBegin;
120628988994SBarry Smith   if (ll) {
12077a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12081ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1209899cda47SBarry Smith     if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1210da3a660dSBarry Smith     for (i=0; i<m; i++) {
1211da3a660dSBarry Smith       x = l[i];
1212da3a660dSBarry Smith       v = mat->v + i;
1213da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1214da3a660dSBarry Smith     }
12151ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1216efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1217da3a660dSBarry Smith   }
121828988994SBarry Smith   if (rr) {
12197a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12201ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1221899cda47SBarry Smith     if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1222da3a660dSBarry Smith     for (i=0; i<n; i++) {
1223da3a660dSBarry Smith       x = r[i];
1224da3a660dSBarry Smith       v = mat->v + i*m;
1225da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1226da3a660dSBarry Smith     }
12271ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1228efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1229da3a660dSBarry Smith   }
12303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1231289bc588SBarry Smith }
1232289bc588SBarry Smith 
12334a2ae208SSatish Balay #undef __FUNCT__
12344a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1235dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1236289bc588SBarry Smith {
1237c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
123887828ca2SBarry Smith   PetscScalar  *v = mat->v;
1239329f5518SBarry Smith   PetscReal    sum = 0.0;
1240899cda47SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap.n,i,j;
1241efee365bSSatish Balay   PetscErrorCode ierr;
124255659b69SBarry Smith 
12433a40ed3dSBarry Smith   PetscFunctionBegin;
1244289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1245a5ce6ee0Svictorle     if (lda>m) {
1246899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
1247a5ce6ee0Svictorle 	v = mat->v+j*lda;
1248a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1249a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1250a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1251a5ce6ee0Svictorle #else
1252a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1253a5ce6ee0Svictorle #endif
1254a5ce6ee0Svictorle 	}
1255a5ce6ee0Svictorle       }
1256a5ce6ee0Svictorle     } else {
1257899cda47SBarry Smith       for (i=0; i<A->cmap.n*A->rmap.n; i++) {
1258aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1259329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1260289bc588SBarry Smith #else
1261289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1262289bc588SBarry Smith #endif
1263289bc588SBarry Smith       }
1264a5ce6ee0Svictorle     }
1265064f8208SBarry Smith     *nrm = sqrt(sum);
1266899cda47SBarry Smith     ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12673a40ed3dSBarry Smith   } else if (type == NORM_1) {
1268064f8208SBarry Smith     *nrm = 0.0;
1269899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
12701b807ce4Svictorle       v = mat->v + j*mat->lda;
1271289bc588SBarry Smith       sum = 0.0;
1272899cda47SBarry Smith       for (i=0; i<A->rmap.n; i++) {
127333a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1274289bc588SBarry Smith       }
1275064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1276289bc588SBarry Smith     }
1277899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12783a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1279064f8208SBarry Smith     *nrm = 0.0;
1280899cda47SBarry Smith     for (j=0; j<A->rmap.n; j++) {
1281289bc588SBarry Smith       v = mat->v + j;
1282289bc588SBarry Smith       sum = 0.0;
1283899cda47SBarry Smith       for (i=0; i<A->cmap.n; i++) {
12841b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1285289bc588SBarry Smith       }
1286064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1287289bc588SBarry Smith     }
1288899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12893a40ed3dSBarry Smith   } else {
129029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1291289bc588SBarry Smith   }
12923a40ed3dSBarry Smith   PetscFunctionReturn(0);
1293289bc588SBarry Smith }
1294289bc588SBarry Smith 
12954a2ae208SSatish Balay #undef __FUNCT__
12964a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
12974e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg)
1298289bc588SBarry Smith {
1299c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
130063ba0a88SBarry Smith   PetscErrorCode ierr;
130167e560aaSBarry Smith 
13023a40ed3dSBarry Smith   PetscFunctionBegin;
1303b5a2b587SKris Buschelman   switch (op) {
1304b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13054e0d8c25SBarry Smith     aij->roworiented = flg;
1306b5a2b587SKris Buschelman     break;
1307512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1308b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13093971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13104e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1311b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1312b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
131377e54ba9SKris Buschelman   case MAT_SYMMETRIC:
131477e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13159a4540c5SBarry Smith   case MAT_HERMITIAN:
13169a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1317600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1318290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
131977e54ba9SKris Buschelman     break;
1320b5a2b587SKris Buschelman   default:
1321600fe468SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13223a40ed3dSBarry Smith   }
13233a40ed3dSBarry Smith   PetscFunctionReturn(0);
1324289bc588SBarry Smith }
1325289bc588SBarry Smith 
13264a2ae208SSatish Balay #undef __FUNCT__
13274a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1328dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13296f0a148fSBarry Smith {
1330ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13316849ba73SBarry Smith   PetscErrorCode ierr;
1332899cda47SBarry Smith   PetscInt       lda=l->lda,m=A->rmap.n,j;
13333a40ed3dSBarry Smith 
13343a40ed3dSBarry Smith   PetscFunctionBegin;
1335a5ce6ee0Svictorle   if (lda>m) {
1336899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
1337a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1338a5ce6ee0Svictorle     }
1339a5ce6ee0Svictorle   } else {
1340899cda47SBarry Smith     ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1341a5ce6ee0Svictorle   }
13423a40ed3dSBarry Smith   PetscFunctionReturn(0);
13436f0a148fSBarry Smith }
13446f0a148fSBarry Smith 
13454a2ae208SSatish Balay #undef __FUNCT__
13464a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1347f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13486f0a148fSBarry Smith {
1349ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1350899cda47SBarry Smith   PetscInt       n = A->cmap.n,i,j;
135187828ca2SBarry Smith   PetscScalar    *slot;
135255659b69SBarry Smith 
13533a40ed3dSBarry Smith   PetscFunctionBegin;
13546f0a148fSBarry Smith   for (i=0; i<N; i++) {
13556f0a148fSBarry Smith     slot = l->v + rows[i];
13566f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13576f0a148fSBarry Smith   }
1358f4df32b1SMatthew Knepley   if (diag != 0.0) {
13596f0a148fSBarry Smith     for (i=0; i<N; i++) {
13606f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1361f4df32b1SMatthew Knepley       *slot = diag;
13626f0a148fSBarry Smith     }
13636f0a148fSBarry Smith   }
13643a40ed3dSBarry Smith   PetscFunctionReturn(0);
13656f0a148fSBarry Smith }
1366557bce09SLois Curfman McInnes 
13674a2ae208SSatish Balay #undef __FUNCT__
13684a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1369dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
137064e87e97SBarry Smith {
1371c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13723a40ed3dSBarry Smith 
13733a40ed3dSBarry Smith   PetscFunctionBegin;
1374899cda47SBarry Smith   if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
137564e87e97SBarry Smith   *array = mat->v;
13763a40ed3dSBarry Smith   PetscFunctionReturn(0);
137764e87e97SBarry Smith }
13780754003eSLois Curfman McInnes 
13794a2ae208SSatish Balay #undef __FUNCT__
13804a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1381dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1382ff14e315SSatish Balay {
13833a40ed3dSBarry Smith   PetscFunctionBegin;
138409b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13853a40ed3dSBarry Smith   PetscFunctionReturn(0);
1386ff14e315SSatish Balay }
13870754003eSLois Curfman McInnes 
13884a2ae208SSatish Balay #undef __FUNCT__
13894a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
139013f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13910754003eSLois Curfman McInnes {
1392c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13936849ba73SBarry Smith   PetscErrorCode ierr;
139421a2c019SBarry Smith   PetscInt       i,j,*irow,*icol,nrows,ncols;
139587828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13960754003eSLois Curfman McInnes   Mat            newmat;
13970754003eSLois Curfman McInnes 
13983a40ed3dSBarry Smith   PetscFunctionBegin;
139978b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
140078b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1401e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1402e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14030754003eSLois Curfman McInnes 
1404182d2002SSatish Balay   /* Check submatrixcall */
1405182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
140613f74950SBarry Smith     PetscInt n_cols,n_rows;
1407182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
140821a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
140921a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
141021a2c019SBarry Smith       ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr);
141121a2c019SBarry Smith     }
1412182d2002SSatish Balay     newmat = *B;
1413182d2002SSatish Balay   } else {
14140754003eSLois Curfman McInnes     /* Create and fill new matrix */
14157adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1416f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14177adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14185c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1419182d2002SSatish Balay   }
1420182d2002SSatish Balay 
1421182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1422182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1423182d2002SSatish Balay 
1424182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14256de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1426182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1427182d2002SSatish Balay       *bv++ = av[irow[j]];
14280754003eSLois Curfman McInnes     }
14290754003eSLois Curfman McInnes   }
1430182d2002SSatish Balay 
1431182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14326d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14336d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14340754003eSLois Curfman McInnes 
14350754003eSLois Curfman McInnes   /* Free work space */
143678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
143778b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1438182d2002SSatish Balay   *B = newmat;
14393a40ed3dSBarry Smith   PetscFunctionReturn(0);
14400754003eSLois Curfman McInnes }
14410754003eSLois Curfman McInnes 
14424a2ae208SSatish Balay #undef __FUNCT__
14434a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
144413f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1445905e6a2fSBarry Smith {
14466849ba73SBarry Smith   PetscErrorCode ierr;
144713f74950SBarry Smith   PetscInt       i;
1448905e6a2fSBarry Smith 
14493a40ed3dSBarry Smith   PetscFunctionBegin;
1450905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1451b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1452905e6a2fSBarry Smith   }
1453905e6a2fSBarry Smith 
1454905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14556a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1456905e6a2fSBarry Smith   }
14573a40ed3dSBarry Smith   PetscFunctionReturn(0);
1458905e6a2fSBarry Smith }
1459905e6a2fSBarry Smith 
14604a2ae208SSatish Balay #undef __FUNCT__
1461c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1462c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1463c0aa2d19SHong Zhang {
1464c0aa2d19SHong Zhang   PetscFunctionBegin;
1465c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1466c0aa2d19SHong Zhang }
1467c0aa2d19SHong Zhang 
1468c0aa2d19SHong Zhang #undef __FUNCT__
1469c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1470c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1471c0aa2d19SHong Zhang {
1472c0aa2d19SHong Zhang   PetscFunctionBegin;
1473c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1474c0aa2d19SHong Zhang }
1475c0aa2d19SHong Zhang 
1476c0aa2d19SHong Zhang #undef __FUNCT__
14774a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1478dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14794b0e389bSBarry Smith {
14804b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14816849ba73SBarry Smith   PetscErrorCode ierr;
1482899cda47SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j;
14833a40ed3dSBarry Smith 
14843a40ed3dSBarry Smith   PetscFunctionBegin;
148533f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
148633f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1487cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14883a40ed3dSBarry Smith     PetscFunctionReturn(0);
14893a40ed3dSBarry Smith   }
1490899cda47SBarry Smith   if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1491a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14920dbb7854Svictorle     for (j=0; j<n; j++) {
1493a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1494a5ce6ee0Svictorle     }
1495a5ce6ee0Svictorle   } else {
1496899cda47SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1497a5ce6ee0Svictorle   }
1498273d9f13SBarry Smith   PetscFunctionReturn(0);
1499273d9f13SBarry Smith }
1500273d9f13SBarry Smith 
15014a2ae208SSatish Balay #undef __FUNCT__
15024a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1503dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1504273d9f13SBarry Smith {
1505dfbe8321SBarry Smith   PetscErrorCode ierr;
1506273d9f13SBarry Smith 
1507273d9f13SBarry Smith   PetscFunctionBegin;
1508273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15093a40ed3dSBarry Smith   PetscFunctionReturn(0);
15104b0e389bSBarry Smith }
15114b0e389bSBarry Smith 
1512284134d9SBarry Smith #undef __FUNCT__
1513284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1514284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1515284134d9SBarry Smith {
1516284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
151721a2c019SBarry Smith   PetscErrorCode ierr;
1518284134d9SBarry Smith   PetscFunctionBegin;
151921a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1520284134d9SBarry Smith   m = PetscMax(m,M);
1521284134d9SBarry Smith   n = PetscMax(n,N);
152221a2c019SBarry Smith   if (m > a->Mmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m);
1523284134d9SBarry Smith   if (n > a->Nmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n);
1524899cda47SBarry Smith   A->rmap.n = A->rmap.n = m;
1525899cda47SBarry Smith   A->cmap.n = A->cmap.N = n;
152621a2c019SBarry Smith   if (a->changelda) a->lda = m;
152721a2c019SBarry Smith   ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
1528284134d9SBarry Smith   PetscFunctionReturn(0);
1529284134d9SBarry Smith }
1530170fe5c8SBarry Smith 
1531284134d9SBarry Smith 
1532a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1533a9fe9ddaSSatish Balay #undef __FUNCT__
1534a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1535a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1536a9fe9ddaSSatish Balay {
1537a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1538a9fe9ddaSSatish Balay 
1539a9fe9ddaSSatish Balay   PetscFunctionBegin;
1540a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1541a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1542a9fe9ddaSSatish Balay   }
1543a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1544a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1545a9fe9ddaSSatish Balay }
1546a9fe9ddaSSatish Balay 
1547a9fe9ddaSSatish Balay #undef __FUNCT__
1548a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1549a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1550a9fe9ddaSSatish Balay {
1551ee16a9a1SHong Zhang   PetscErrorCode ierr;
1552899cda47SBarry Smith   PetscInt       m=A->rmap.n,n=B->cmap.n;
1553ee16a9a1SHong Zhang   Mat            Cmat;
1554a9fe9ddaSSatish Balay 
1555ee16a9a1SHong Zhang   PetscFunctionBegin;
1556899cda47SBarry Smith   if (A->cmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap.n %d != B->rmap.n %d\n",A->cmap.n,B->rmap.n);
1557ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1558ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1559ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1560ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1561ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1562ee16a9a1SHong Zhang   *C = Cmat;
1563ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1564ee16a9a1SHong Zhang }
1565a9fe9ddaSSatish Balay 
156698a3b096SSatish Balay #undef __FUNCT__
1567a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1568a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1569a9fe9ddaSSatish Balay {
1570a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1571a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1572a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
15730805154bSBarry Smith   PetscBLASInt   m,n,k;
1574a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1575a9fe9ddaSSatish Balay 
1576a9fe9ddaSSatish Balay   PetscFunctionBegin;
15770805154bSBarry Smith   m = PetscBLASIntCast(A->rmap.n);
15780805154bSBarry Smith   n = PetscBLASIntCast(B->cmap.n);
15790805154bSBarry Smith   k = PetscBLASIntCast(A->cmap.n);
1580a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1581a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1582a9fe9ddaSSatish Balay }
1583a9fe9ddaSSatish Balay 
1584a9fe9ddaSSatish Balay #undef __FUNCT__
1585a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1586a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1587a9fe9ddaSSatish Balay {
1588a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1589a9fe9ddaSSatish Balay 
1590a9fe9ddaSSatish Balay   PetscFunctionBegin;
1591a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1592a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1593a9fe9ddaSSatish Balay   }
1594a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1595a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1596a9fe9ddaSSatish Balay }
1597a9fe9ddaSSatish Balay 
1598a9fe9ddaSSatish Balay #undef __FUNCT__
1599a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1600a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1601a9fe9ddaSSatish Balay {
1602ee16a9a1SHong Zhang   PetscErrorCode ierr;
1603899cda47SBarry Smith   PetscInt       m=A->cmap.n,n=B->cmap.n;
1604ee16a9a1SHong Zhang   Mat            Cmat;
1605a9fe9ddaSSatish Balay 
1606ee16a9a1SHong Zhang   PetscFunctionBegin;
1607899cda47SBarry Smith   if (A->rmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->rmap.n %d != B->rmap.n %d\n",A->rmap.n,B->rmap.n);
1608ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1609ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1610ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1611ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1612ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1613ee16a9a1SHong Zhang   *C = Cmat;
1614ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1615ee16a9a1SHong Zhang }
1616a9fe9ddaSSatish Balay 
1617a9fe9ddaSSatish Balay #undef __FUNCT__
1618a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1619a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1620a9fe9ddaSSatish Balay {
1621a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1622a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1623a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16240805154bSBarry Smith   PetscBLASInt   m,n,k;
1625a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1626a9fe9ddaSSatish Balay 
1627a9fe9ddaSSatish Balay   PetscFunctionBegin;
16282fbe02b9SBarry Smith   m = PetscBLASIntCast(A->cmap.n);
16290805154bSBarry Smith   n = PetscBLASIntCast(B->cmap.n);
16302fbe02b9SBarry Smith   k = PetscBLASIntCast(A->rmap.n);
16312fbe02b9SBarry Smith   /*
16322fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
16332fbe02b9SBarry Smith   */
1634a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1635a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1636a9fe9ddaSSatish Balay }
1637985db425SBarry Smith 
1638985db425SBarry Smith #undef __FUNCT__
1639985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1640985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1641985db425SBarry Smith {
1642985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1643985db425SBarry Smith   PetscErrorCode ierr;
1644985db425SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n,p;
1645985db425SBarry Smith   PetscScalar    *x;
1646985db425SBarry Smith   MatScalar      *aa = a->v;
1647985db425SBarry Smith 
1648985db425SBarry Smith   PetscFunctionBegin;
1649985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1650985db425SBarry Smith 
1651985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1652985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1653985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1654985db425SBarry Smith   if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1655985db425SBarry Smith   for (i=0; i<m; i++) {
1656985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1657985db425SBarry Smith     for (j=1; j<n; j++){
1658985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1659985db425SBarry Smith     }
1660985db425SBarry Smith   }
1661985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1662985db425SBarry Smith   PetscFunctionReturn(0);
1663985db425SBarry Smith }
1664985db425SBarry Smith 
1665985db425SBarry Smith #undef __FUNCT__
1666985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1667985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1668985db425SBarry Smith {
1669985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1670985db425SBarry Smith   PetscErrorCode ierr;
1671985db425SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n,p;
1672985db425SBarry Smith   PetscScalar    *x;
1673985db425SBarry Smith   PetscReal      atmp;
1674985db425SBarry Smith   MatScalar      *aa = a->v;
1675985db425SBarry Smith 
1676985db425SBarry Smith   PetscFunctionBegin;
1677985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1678985db425SBarry Smith 
1679985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1680985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1681985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1682985db425SBarry Smith   if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1683985db425SBarry Smith   for (i=0; i<m; i++) {
16849189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1685985db425SBarry Smith     for (j=1; j<n; j++){
1686985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1687985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1688985db425SBarry Smith     }
1689985db425SBarry Smith   }
1690985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1691985db425SBarry Smith   PetscFunctionReturn(0);
1692985db425SBarry Smith }
1693985db425SBarry Smith 
1694985db425SBarry Smith #undef __FUNCT__
1695985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1696985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1697985db425SBarry Smith {
1698985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1699985db425SBarry Smith   PetscErrorCode ierr;
1700985db425SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n,p;
1701985db425SBarry Smith   PetscScalar    *x;
1702985db425SBarry Smith   MatScalar      *aa = a->v;
1703985db425SBarry Smith 
1704985db425SBarry Smith   PetscFunctionBegin;
1705985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1706985db425SBarry Smith 
1707985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1708985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1709985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1710985db425SBarry Smith   if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1711985db425SBarry Smith   for (i=0; i<m; i++) {
1712985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1713985db425SBarry Smith     for (j=1; j<n; j++){
1714985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1715985db425SBarry Smith     }
1716985db425SBarry Smith   }
1717985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1718985db425SBarry Smith   PetscFunctionReturn(0);
1719985db425SBarry Smith }
1720985db425SBarry Smith 
17218d0534beSBarry Smith #undef __FUNCT__
17228d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
17238d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
17248d0534beSBarry Smith {
17258d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17268d0534beSBarry Smith   PetscErrorCode ierr;
17278d0534beSBarry Smith   PetscScalar    *x;
17288d0534beSBarry Smith 
17298d0534beSBarry Smith   PetscFunctionBegin;
17308d0534beSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
17318d0534beSBarry Smith 
17328d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
17338d0534beSBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
17348d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
17358d0534beSBarry Smith   PetscFunctionReturn(0);
17368d0534beSBarry Smith }
17378d0534beSBarry Smith 
1738289bc588SBarry Smith /* -------------------------------------------------------------------*/
1739a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1740905e6a2fSBarry Smith        MatGetRow_SeqDense,
1741905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1742905e6a2fSBarry Smith        MatMult_SeqDense,
174397304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
17447c922b88SBarry Smith        MatMultTranspose_SeqDense,
17457c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1746905e6a2fSBarry Smith        MatSolve_SeqDense,
1747905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
17487c922b88SBarry Smith        MatSolveTranspose_SeqDense,
174997304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1750905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1751905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1752ec8511deSBarry Smith        MatRelax_SeqDense,
1753ec8511deSBarry Smith        MatTranspose_SeqDense,
175497304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1755905e6a2fSBarry Smith        MatEqual_SeqDense,
1756905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1757905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1758905e6a2fSBarry Smith        MatNorm_SeqDense,
1759c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1760c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1761905e6a2fSBarry Smith        0,
1762905e6a2fSBarry Smith        MatSetOption_SeqDense,
1763905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
176497304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1765905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1766905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1767905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1768905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
176997304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1770273d9f13SBarry Smith        0,
1771905e6a2fSBarry Smith        0,
1772905e6a2fSBarry Smith        MatGetArray_SeqDense,
1773905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
177497304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1775a5ae1ecdSBarry Smith        0,
1776a5ae1ecdSBarry Smith        0,
1777a5ae1ecdSBarry Smith        0,
1778a5ae1ecdSBarry Smith        0,
177997304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1780a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1781a5ae1ecdSBarry Smith        0,
17824b0e389bSBarry Smith        MatGetValues_SeqDense,
1783a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1784985db425SBarry Smith /*45*/ MatGetRowMax_SeqDense,
1785a5ae1ecdSBarry Smith        MatScale_SeqDense,
1786a5ae1ecdSBarry Smith        0,
1787a5ae1ecdSBarry Smith        0,
1788a5ae1ecdSBarry Smith        0,
1789521d7252SBarry Smith /*50*/ 0,
1790a5ae1ecdSBarry Smith        0,
1791a5ae1ecdSBarry Smith        0,
1792a5ae1ecdSBarry Smith        0,
1793a5ae1ecdSBarry Smith        0,
179497304618SKris Buschelman /*55*/ 0,
1795a5ae1ecdSBarry Smith        0,
1796a5ae1ecdSBarry Smith        0,
1797a5ae1ecdSBarry Smith        0,
1798a5ae1ecdSBarry Smith        0,
179997304618SKris Buschelman /*60*/ 0,
1800e03a110bSBarry Smith        MatDestroy_SeqDense,
1801e03a110bSBarry Smith        MatView_SeqDense,
1802357abbc8SBarry Smith        0,
180397304618SKris Buschelman        0,
180497304618SKris Buschelman /*65*/ 0,
180597304618SKris Buschelman        0,
180697304618SKris Buschelman        0,
180797304618SKris Buschelman        0,
180897304618SKris Buschelman        0,
1809985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqDense,
181097304618SKris Buschelman        0,
181197304618SKris Buschelman        0,
181297304618SKris Buschelman        0,
181397304618SKris Buschelman        0,
181497304618SKris Buschelman /*75*/ 0,
181597304618SKris Buschelman        0,
181697304618SKris Buschelman        0,
181797304618SKris Buschelman        0,
181897304618SKris Buschelman        0,
181997304618SKris Buschelman /*80*/ 0,
182097304618SKris Buschelman        0,
182197304618SKris Buschelman        0,
182297304618SKris Buschelman        0,
1823865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1824865e5f61SKris Buschelman        0,
18251cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1826865e5f61SKris Buschelman        0,
1827865e5f61SKris Buschelman        0,
1828865e5f61SKris Buschelman        0,
1829a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense,
1830a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1831a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1832865e5f61SKris Buschelman        0,
1833865e5f61SKris Buschelman        0,
1834865e5f61SKris Buschelman /*95*/ 0,
1835a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1836a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1837a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1838284134d9SBarry Smith        0,
1839284134d9SBarry Smith /*100*/0,
1840284134d9SBarry Smith        0,
1841284134d9SBarry Smith        0,
1842284134d9SBarry Smith        0,
1843985db425SBarry Smith        MatSetSizes_SeqDense,
1844985db425SBarry Smith        0,
1845985db425SBarry Smith        0,
1846985db425SBarry Smith        0,
1847985db425SBarry Smith        0,
1848985db425SBarry Smith        0,
1849985db425SBarry Smith /*110*/0,
1850985db425SBarry Smith        0,
18518d0534beSBarry Smith        MatGetRowMin_SeqDense,
18528d0534beSBarry Smith        MatGetColumnVector_SeqDense
1853985db425SBarry Smith };
185490ace30eSBarry Smith 
18554a2ae208SSatish Balay #undef __FUNCT__
18564a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18574b828684SBarry Smith /*@C
1858fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1859d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1860d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1861289bc588SBarry Smith 
1862db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1863db81eaa0SLois Curfman McInnes 
186420563c6bSBarry Smith    Input Parameters:
1865db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
18660c775827SLois Curfman McInnes .  m - number of rows
186718f449edSLois Curfman McInnes .  n - number of columns
1868db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1869dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
187020563c6bSBarry Smith 
187120563c6bSBarry Smith    Output Parameter:
187244cd7ae7SLois Curfman McInnes .  A - the matrix
187320563c6bSBarry Smith 
1874b259b22eSLois Curfman McInnes    Notes:
187518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
187618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1877b4fd4287SBarry Smith    set data=PETSC_NULL.
187818f449edSLois Curfman McInnes 
1879027ccd11SLois Curfman McInnes    Level: intermediate
1880027ccd11SLois Curfman McInnes 
1881dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1882d65003e9SLois Curfman McInnes 
1883db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
188420563c6bSBarry Smith @*/
1885be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1886289bc588SBarry Smith {
1887dfbe8321SBarry Smith   PetscErrorCode ierr;
18883b2fbd54SBarry Smith 
18893a40ed3dSBarry Smith   PetscFunctionBegin;
1890f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1891f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1892273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1893273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1894273d9f13SBarry Smith   PetscFunctionReturn(0);
1895273d9f13SBarry Smith }
1896273d9f13SBarry Smith 
18974a2ae208SSatish Balay #undef __FUNCT__
1898afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1899273d9f13SBarry Smith /*@C
1900273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1901273d9f13SBarry Smith 
1902273d9f13SBarry Smith    Collective on MPI_Comm
1903273d9f13SBarry Smith 
1904273d9f13SBarry Smith    Input Parameters:
1905273d9f13SBarry Smith +  A - the matrix
1906273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1907273d9f13SBarry Smith 
1908273d9f13SBarry Smith    Notes:
1909273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1910273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1911284134d9SBarry Smith    need not call this routine.
1912273d9f13SBarry Smith 
1913273d9f13SBarry Smith    Level: intermediate
1914273d9f13SBarry Smith 
1915273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1916273d9f13SBarry Smith 
1917273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1918273d9f13SBarry Smith @*/
1919be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1920273d9f13SBarry Smith {
1921dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1922a23d5eceSKris Buschelman 
1923a23d5eceSKris Buschelman   PetscFunctionBegin;
1924a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1925a23d5eceSKris Buschelman   if (f) {
1926a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1927a23d5eceSKris Buschelman   }
1928a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1929a23d5eceSKris Buschelman }
1930a23d5eceSKris Buschelman 
1931a23d5eceSKris Buschelman EXTERN_C_BEGIN
1932a23d5eceSKris Buschelman #undef __FUNCT__
1933afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1934be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1935a23d5eceSKris Buschelman {
1936273d9f13SBarry Smith   Mat_SeqDense   *b;
1937dfbe8321SBarry Smith   PetscErrorCode ierr;
1938273d9f13SBarry Smith 
1939273d9f13SBarry Smith   PetscFunctionBegin;
1940273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1941273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
19425afd5e0cSBarry Smith   if (b->lda <= 0) b->lda = B->rmap.n;
19439e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
19449e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
19455afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1946284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1947284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
19489e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
1949273d9f13SBarry Smith   } else { /* user-allocated storage */
19509e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1951273d9f13SBarry Smith     b->v          = data;
1952273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1953273d9f13SBarry Smith   }
1954273d9f13SBarry Smith   PetscFunctionReturn(0);
1955273d9f13SBarry Smith }
1956a23d5eceSKris Buschelman EXTERN_C_END
1957273d9f13SBarry Smith 
19581b807ce4Svictorle #undef __FUNCT__
19591b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
19601b807ce4Svictorle /*@C
19611b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
19621b807ce4Svictorle 
19631b807ce4Svictorle   Input parameter:
19641b807ce4Svictorle + A - the matrix
19651b807ce4Svictorle - lda - the leading dimension
19661b807ce4Svictorle 
19671b807ce4Svictorle   Notes:
19681b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
19691b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
19701b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
19711b807ce4Svictorle 
19721b807ce4Svictorle   Level: intermediate
19731b807ce4Svictorle 
19741b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
19751b807ce4Svictorle 
1976284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
19771b807ce4Svictorle @*/
1978be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
19791b807ce4Svictorle {
19801b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
198121a2c019SBarry Smith 
19821b807ce4Svictorle   PetscFunctionBegin;
1983899cda47SBarry Smith   if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n);
19841b807ce4Svictorle   b->lda       = lda;
198521a2c019SBarry Smith   b->changelda = PETSC_FALSE;
198621a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
19871b807ce4Svictorle   PetscFunctionReturn(0);
19881b807ce4Svictorle }
19891b807ce4Svictorle 
19900bad9183SKris Buschelman /*MC
1991fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
19920bad9183SKris Buschelman 
19930bad9183SKris Buschelman    Options Database Keys:
19940bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
19950bad9183SKris Buschelman 
19960bad9183SKris Buschelman   Level: beginner
19970bad9183SKris Buschelman 
199889665df3SBarry Smith .seealso: MatCreateSeqDense()
199989665df3SBarry Smith 
20000bad9183SKris Buschelman M*/
20010bad9183SKris Buschelman 
2002273d9f13SBarry Smith EXTERN_C_BEGIN
20034a2ae208SSatish Balay #undef __FUNCT__
20044a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
2005be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2006273d9f13SBarry Smith {
2007273d9f13SBarry Smith   Mat_SeqDense   *b;
2008dfbe8321SBarry Smith   PetscErrorCode ierr;
20097c334f02SBarry Smith   PetscMPIInt    size;
2010273d9f13SBarry Smith 
2011273d9f13SBarry Smith   PetscFunctionBegin;
20127adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
201329bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
201455659b69SBarry Smith 
2015899cda47SBarry Smith   B->rmap.bs = B->cmap.bs = 1;
20166148ca0dSBarry Smith   ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr);
20176148ca0dSBarry Smith   ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr);
2018273d9f13SBarry Smith 
201938f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2020549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
202190f02eecSBarry Smith   B->mapping      = 0;
202244cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
202318f449edSLois Curfman McInnes 
2024a5ae1ecdSBarry Smith 
202544cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2026273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2027273d9f13SBarry Smith   b->v            = 0;
2028899cda47SBarry Smith   b->lda          = B->rmap.n;
202921a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
2030899cda47SBarry Smith   b->Mmax         = B->rmap.n;
2031899cda47SBarry Smith   b->Nmax         = B->cmap.n;
20324e220ebcSLois Curfman McInnes 
2033b24902e0SBarry Smith 
2034b24902e0SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqdense_petsc_C",
2035b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2036b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2037a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2038a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2039a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
20404ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
20414ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
20424ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
20434ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
20444ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
20454ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
20464ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
20474ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
20484ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
204917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
20503a40ed3dSBarry Smith   PetscFunctionReturn(0);
2051289bc588SBarry Smith }
2052c0aa2d19SHong Zhang 
2053c0aa2d19SHong Zhang 
2054273d9f13SBarry Smith EXTERN_C_END
2055