xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
77c4f633dSBarry Smith #include "../src/mat/impls/dense/seq/dense.h"
8f3da1532SBarry Smith #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;
21d0f46423SBarry Smith   N    = PetscBLASIntCast(X->rmap->n*X->cmap->n);
22d0f46423SBarry Smith   m    = PetscBLASIntCast(X->rmap->n);
230805154bSBarry Smith   ldax = PetscBLASIntCast(x->lda);
240805154bSBarry Smith   lday = PetscBLASIntCast(y->lda);
25a5ce6ee0Svictorle   if (ldax>m || lday>m) {
26d0f46423SBarry 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   }
320450473dSBarry Smith   ierr = PetscLogFlops(PetscMax(2*N-1,0));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 {
40d0f46423SBarry Smith   PetscInt     N = A->rmap->n*A->cmap->n;
413a40ed3dSBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
434e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
444e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
456de62eeeSBarry Smith   info->nz_used           = (double)N;
466de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
474e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
484e220ebcSLois Curfman McInnes   info->mallocs           = 0;
497adad957SLisandro Dalcin   info->memory            = ((PetscObject)A)->mem;
504e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
514e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
524e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
54289bc588SBarry Smith }
55289bc588SBarry Smith 
564a2ae208SSatish Balay #undef __FUNCT__
574a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
58f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
5980cd9d93SLois Curfman McInnes {
60273d9f13SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
61f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
62efee365bSSatish Balay   PetscErrorCode ierr;
630805154bSBarry Smith   PetscBLASInt   one = 1,j,nz,lda = PetscBLASIntCast(a->lda);
6480cd9d93SLois Curfman McInnes 
653a40ed3dSBarry Smith   PetscFunctionBegin;
66d0f46423SBarry Smith   if (lda>A->rmap->n) {
67d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n);
68d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
69f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
72d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n*A->cmap->n);
73f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
74a5ce6ee0Svictorle   }
75efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
7780cd9d93SLois Curfman McInnes }
7880cd9d93SLois Curfman McInnes 
791cbb95d3SBarry Smith #undef __FUNCT__
801cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense"
811cbb95d3SBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscTruth *fl)
821cbb95d3SBarry Smith {
831cbb95d3SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
84d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,N;
851cbb95d3SBarry Smith   PetscScalar    *v = a->v;
861cbb95d3SBarry Smith 
871cbb95d3SBarry Smith   PetscFunctionBegin;
881cbb95d3SBarry Smith   *fl = PETSC_FALSE;
89d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0);
901cbb95d3SBarry Smith   N = a->lda;
911cbb95d3SBarry Smith 
921cbb95d3SBarry Smith   for (i=0; i<m; i++) {
931cbb95d3SBarry Smith     for (j=i+1; j<m; j++) {
941cbb95d3SBarry Smith       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
951cbb95d3SBarry Smith     }
961cbb95d3SBarry Smith   }
971cbb95d3SBarry Smith   *fl = PETSC_TRUE;
981cbb95d3SBarry Smith   PetscFunctionReturn(0);
991cbb95d3SBarry Smith }
1001cbb95d3SBarry Smith 
101b24902e0SBarry Smith #undef __FUNCT__
102b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense"
103719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
104b24902e0SBarry Smith {
105b24902e0SBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
106b24902e0SBarry Smith   PetscErrorCode ierr;
107b24902e0SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
108b24902e0SBarry Smith 
109b24902e0SBarry Smith   PetscFunctionBegin;
110719d5645SBarry Smith   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
111b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
112b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
113d0f46423SBarry Smith     if (lda>A->rmap->n) {
114d0f46423SBarry Smith       m = A->rmap->n;
115d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
116b24902e0SBarry Smith 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
117b24902e0SBarry Smith       }
118b24902e0SBarry Smith     } else {
119d0f46423SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
120b24902e0SBarry Smith     }
121b24902e0SBarry Smith   }
122b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
123b24902e0SBarry Smith   PetscFunctionReturn(0);
124b24902e0SBarry Smith }
125b24902e0SBarry Smith 
1264a2ae208SSatish Balay #undef __FUNCT__
1274a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
128dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12902cad45dSBarry Smith {
1306849ba73SBarry Smith   PetscErrorCode ierr;
13102cad45dSBarry Smith 
1323a40ed3dSBarry Smith   PetscFunctionBegin;
1335c9eb25fSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
134d0f46423SBarry Smith   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1355c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
136719d5645SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
137b24902e0SBarry Smith   PetscFunctionReturn(0);
138b24902e0SBarry Smith }
139b24902e0SBarry Smith 
1406ee01492SSatish Balay 
1410481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
142719d5645SBarry Smith 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1450481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
146289bc588SBarry Smith {
1474482741eSBarry Smith   MatFactorInfo  info;
148a093e273SMatthew Knepley   PetscErrorCode ierr;
1493a40ed3dSBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
151c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
152719d5645SBarry Smith   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154289bc588SBarry Smith }
1556ee01492SSatish Balay 
1560b4b3355SBarry Smith #undef __FUNCT__
1574a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
158dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
159289bc588SBarry Smith {
160c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1616849ba73SBarry Smith   PetscErrorCode ierr;
16287828ca2SBarry Smith   PetscScalar    *x,*y;
163d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
16467e560aaSBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
1661ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1671ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
168d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
169d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
170ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
17180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
172ae7cfcebSSatish Balay #else
17371044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
1744ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
175ae7cfcebSSatish Balay #endif
176d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_CHOLESKY){
177ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
17880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
179ae7cfcebSSatish Balay #else
18071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
1814ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
182ae7cfcebSSatish Balay #endif
183289bc588SBarry Smith   }
18429bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
1851ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1861ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
187dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
189289bc588SBarry Smith }
1906ee01492SSatish Balay 
1914a2ae208SSatish Balay #undef __FUNCT__
1924a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
193dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
194da3a660dSBarry Smith {
195c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
196dfbe8321SBarry Smith   PetscErrorCode ierr;
19787828ca2SBarry Smith   PetscScalar    *x,*y;
198d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
19967e560aaSBarry Smith 
2003a40ed3dSBarry Smith   PetscFunctionBegin;
2011ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2021ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
203d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
204752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
205da3a660dSBarry Smith   if (mat->pivots) {
206ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
20780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
208ae7cfcebSSatish Balay #else
20971044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2104ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
211ae7cfcebSSatish Balay #endif
2127a97a34bSBarry Smith   } else {
213ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
21480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
215ae7cfcebSSatish Balay #else
21671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2174ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
218ae7cfcebSSatish Balay #endif
219da3a660dSBarry Smith   }
2201ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2211ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
222dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
2233a40ed3dSBarry Smith   PetscFunctionReturn(0);
224da3a660dSBarry Smith }
2256ee01492SSatish Balay 
2264a2ae208SSatish Balay #undef __FUNCT__
2274a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
228dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
229da3a660dSBarry Smith {
230c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
231dfbe8321SBarry Smith   PetscErrorCode ierr;
23287828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
233da3a660dSBarry Smith   Vec            tmp = 0;
234d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
23567e560aaSBarry Smith 
2363a40ed3dSBarry Smith   PetscFunctionBegin;
2371ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2381ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
239d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
240da3a660dSBarry Smith   if (yy == zz) {
24178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
24252e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
24378b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
244da3a660dSBarry Smith   }
245d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
246752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
247da3a660dSBarry Smith   if (mat->pivots) {
248ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
24980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
250ae7cfcebSSatish Balay #else
25171044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2522ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
253ae7cfcebSSatish Balay #endif
254a8c6a408SBarry Smith   } else {
255ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
25680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
257ae7cfcebSSatish Balay #else
25871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2592ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
260ae7cfcebSSatish Balay #endif
261da3a660dSBarry Smith   }
2622dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
2632dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
2641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2651ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
266dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
2673a40ed3dSBarry Smith   PetscFunctionReturn(0);
268da3a660dSBarry Smith }
26967e560aaSBarry Smith 
2704a2ae208SSatish Balay #undef __FUNCT__
2714a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
272dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
273da3a660dSBarry Smith {
274c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2756849ba73SBarry Smith   PetscErrorCode ierr;
27687828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
277da3a660dSBarry Smith   Vec            tmp;
278d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
27967e560aaSBarry Smith 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
281d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
2821ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2831ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
284da3a660dSBarry Smith   if (yy == zz) {
28578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
28652e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
28778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
288da3a660dSBarry Smith   }
289d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
290752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
291da3a660dSBarry Smith   if (mat->pivots) {
292ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
29380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
294ae7cfcebSSatish Balay #else
29571044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2962ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
297ae7cfcebSSatish Balay #endif
2983a40ed3dSBarry Smith   } else {
299ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
30080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
301ae7cfcebSSatish Balay #else
30271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3032ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
304ae7cfcebSSatish Balay #endif
305da3a660dSBarry Smith   }
30690f02eecSBarry Smith   if (tmp) {
3072dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
30890f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3093a40ed3dSBarry Smith   } else {
3102dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
31190f02eecSBarry Smith   }
3121ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3131ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
314dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
316da3a660dSBarry Smith }
317db4efbfdSBarry Smith 
318db4efbfdSBarry Smith /* ---------------------------------------------------------------*/
319db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
320db4efbfdSBarry Smith    rather than put it in the Mat->row slot.*/
321db4efbfdSBarry Smith #undef __FUNCT__
322db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense"
3230481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
324db4efbfdSBarry Smith {
325db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
326db4efbfdSBarry Smith   PetscFunctionBegin;
327db4efbfdSBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
328db4efbfdSBarry Smith #else
329db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
330db4efbfdSBarry Smith   PetscErrorCode ierr;
331db4efbfdSBarry Smith   PetscBLASInt   n,m,info;
332db4efbfdSBarry Smith 
333db4efbfdSBarry Smith   PetscFunctionBegin;
334db4efbfdSBarry Smith   n = PetscBLASIntCast(A->cmap->n);
335db4efbfdSBarry Smith   m = PetscBLASIntCast(A->rmap->n);
336db4efbfdSBarry Smith   if (!mat->pivots) {
337db4efbfdSBarry Smith     ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
338db4efbfdSBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
339db4efbfdSBarry Smith   }
340db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
341db4efbfdSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
342db4efbfdSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
343db4efbfdSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
344db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
345db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
346db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
347db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
348d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_LU;
349db4efbfdSBarry Smith 
350dc0b31edSSatish Balay   ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
351db4efbfdSBarry Smith #endif
352db4efbfdSBarry Smith   PetscFunctionReturn(0);
353db4efbfdSBarry Smith }
354db4efbfdSBarry Smith 
355db4efbfdSBarry Smith #undef __FUNCT__
356db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense"
3570481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
358db4efbfdSBarry Smith {
359db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
360db4efbfdSBarry Smith   PetscFunctionBegin;
361db4efbfdSBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
362db4efbfdSBarry Smith #else
363db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
364db4efbfdSBarry Smith   PetscErrorCode ierr;
365db4efbfdSBarry Smith   PetscBLASInt   info,n = PetscBLASIntCast(A->cmap->n);
366db4efbfdSBarry Smith 
367db4efbfdSBarry Smith   PetscFunctionBegin;
368db4efbfdSBarry Smith   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
369db4efbfdSBarry Smith   mat->pivots = 0;
370db4efbfdSBarry Smith 
371db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
372db4efbfdSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
373db4efbfdSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
374db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
375db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
376db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
377db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
378d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_CHOLESKY;
379dc0b31edSSatish Balay   ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr);
380db4efbfdSBarry Smith #endif
381db4efbfdSBarry Smith   PetscFunctionReturn(0);
382db4efbfdSBarry Smith }
383db4efbfdSBarry Smith 
384db4efbfdSBarry Smith 
385db4efbfdSBarry Smith #undef __FUNCT__
386db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
3870481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
388db4efbfdSBarry Smith {
389db4efbfdSBarry Smith   PetscErrorCode ierr;
390db4efbfdSBarry Smith   MatFactorInfo  info;
391db4efbfdSBarry Smith 
392db4efbfdSBarry Smith   PetscFunctionBegin;
393db4efbfdSBarry Smith   info.fill = 1.0;
394c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
395719d5645SBarry Smith   ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr);
396db4efbfdSBarry Smith   PetscFunctionReturn(0);
397db4efbfdSBarry Smith }
398db4efbfdSBarry Smith 
399db4efbfdSBarry Smith #undef __FUNCT__
400db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
4010481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
402db4efbfdSBarry Smith {
403db4efbfdSBarry Smith   PetscFunctionBegin;
404c3ef05f6SHong Zhang   fact->assembled                  = PETSC_TRUE;
405719d5645SBarry Smith   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
406db4efbfdSBarry Smith   PetscFunctionReturn(0);
407db4efbfdSBarry Smith }
408db4efbfdSBarry Smith 
409db4efbfdSBarry Smith #undef __FUNCT__
410db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
4110481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
412db4efbfdSBarry Smith {
413db4efbfdSBarry Smith   PetscFunctionBegin;
414c3ef05f6SHong Zhang   fact->assembled            = PETSC_TRUE;
415719d5645SBarry Smith   fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
416db4efbfdSBarry Smith   PetscFunctionReturn(0);
417db4efbfdSBarry Smith }
418db4efbfdSBarry Smith 
419bb5747d9SMatthew Knepley EXTERN_C_BEGIN
420db4efbfdSBarry Smith #undef __FUNCT__
421db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc"
422db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
423db4efbfdSBarry Smith {
424db4efbfdSBarry Smith   PetscErrorCode ierr;
425db4efbfdSBarry Smith 
426db4efbfdSBarry Smith   PetscFunctionBegin;
427db4efbfdSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
428db4efbfdSBarry Smith   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
429db4efbfdSBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
430db4efbfdSBarry Smith   if (ftype == MAT_FACTOR_LU){
431db4efbfdSBarry Smith     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
432db4efbfdSBarry Smith   } else {
433db4efbfdSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
434db4efbfdSBarry Smith   }
435d5f3da31SBarry Smith   (*fact)->factortype = ftype;
436db4efbfdSBarry Smith   PetscFunctionReturn(0);
437db4efbfdSBarry Smith }
438bb5747d9SMatthew Knepley EXTERN_C_END
439db4efbfdSBarry Smith 
440289bc588SBarry Smith /* ------------------------------------------------------------------*/
4414a2ae208SSatish Balay #undef __FUNCT__
44241f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense"
44341f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
444289bc588SBarry Smith {
445c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
44687828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
447dfbe8321SBarry Smith   PetscErrorCode ierr;
448d0f46423SBarry Smith   PetscInt       m = A->rmap->n,i;
449aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4500805154bSBarry Smith   PetscBLASInt   o = 1,bm = PetscBLASIntCast(m);
451bc1b551cSSatish Balay #endif
452289bc588SBarry Smith 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
454289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
45571044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4562dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
457289bc588SBarry Smith   }
4581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4591ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
460b965ef7fSBarry Smith   its  = its*lits;
46177431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
462289bc588SBarry Smith   while (its--) {
463fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
464289bc588SBarry Smith       for (i=0; i<m; i++) {
465aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
466f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
467f1747703SBarry Smith            not happy about returning a double complex */
46813f74950SBarry Smith         PetscInt    _i;
46987828ca2SBarry Smith         PetscScalar sum = b[i];
470f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4713f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
472f1747703SBarry Smith         }
473f1747703SBarry Smith         xt = sum;
474f1747703SBarry Smith #else
47571044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
476f1747703SBarry Smith #endif
47755a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
478289bc588SBarry Smith       }
479289bc588SBarry Smith     }
480fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
481289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
483f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
484f1747703SBarry Smith            not happy about returning a double complex */
48513f74950SBarry Smith         PetscInt    _i;
48687828ca2SBarry Smith         PetscScalar sum = b[i];
487f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4883f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
489f1747703SBarry Smith         }
490f1747703SBarry Smith         xt = sum;
491f1747703SBarry Smith #else
49271044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
493f1747703SBarry Smith #endif
49455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
495289bc588SBarry Smith       }
496289bc588SBarry Smith     }
497289bc588SBarry Smith   }
4981ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4991ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5003a40ed3dSBarry Smith   PetscFunctionReturn(0);
501289bc588SBarry Smith }
502289bc588SBarry Smith 
503289bc588SBarry Smith /* -----------------------------------------------------------------*/
5044a2ae208SSatish Balay #undef __FUNCT__
5054a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
506dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
507289bc588SBarry Smith {
508c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
50987828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
510dfbe8321SBarry Smith   PetscErrorCode ierr;
5110805154bSBarry Smith   PetscBLASInt   m, n,_One=1;
512ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
5133a40ed3dSBarry Smith 
5143a40ed3dSBarry Smith   PetscFunctionBegin;
515d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
516d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
517d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5191ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
52071044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
5211ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5221ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
523dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
5243a40ed3dSBarry Smith   PetscFunctionReturn(0);
525289bc588SBarry Smith }
526800995b7SMatthew Knepley 
5274a2ae208SSatish Balay #undef __FUNCT__
5284a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
529dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
530289bc588SBarry Smith {
531c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
53287828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
533dfbe8321SBarry Smith   PetscErrorCode ierr;
5340805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5353a40ed3dSBarry Smith 
5363a40ed3dSBarry Smith   PetscFunctionBegin;
537d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
538d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
539d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5401ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5411ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
54271044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5431ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5441ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
545dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
5463a40ed3dSBarry Smith   PetscFunctionReturn(0);
547289bc588SBarry Smith }
5486ee01492SSatish Balay 
5494a2ae208SSatish Balay #undef __FUNCT__
5504a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
551dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
552289bc588SBarry Smith {
553c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
55487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
555dfbe8321SBarry Smith   PetscErrorCode ierr;
5560805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5573a40ed3dSBarry Smith 
5583a40ed3dSBarry Smith   PetscFunctionBegin;
559d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
560d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
561d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
562600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5631ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5641ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
56571044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5661ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5671ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
568dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5693a40ed3dSBarry Smith   PetscFunctionReturn(0);
570289bc588SBarry Smith }
5716ee01492SSatish Balay 
5724a2ae208SSatish Balay #undef __FUNCT__
5734a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
574dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
575289bc588SBarry Smith {
576c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
57787828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
578dfbe8321SBarry Smith   PetscErrorCode ierr;
5790805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
58087828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5813a40ed3dSBarry Smith 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
583d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
584d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
585d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
586600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5871ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5881ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
58971044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5901ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5911ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
592dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5933a40ed3dSBarry Smith   PetscFunctionReturn(0);
594289bc588SBarry Smith }
595289bc588SBarry Smith 
596289bc588SBarry Smith /* -----------------------------------------------------------------*/
5974a2ae208SSatish Balay #undef __FUNCT__
5984a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
59913f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
600289bc588SBarry Smith {
601c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
60287828ca2SBarry Smith   PetscScalar    *v;
6036849ba73SBarry Smith   PetscErrorCode ierr;
60413f74950SBarry Smith   PetscInt       i;
60567e560aaSBarry Smith 
6063a40ed3dSBarry Smith   PetscFunctionBegin;
607d0f46423SBarry Smith   *ncols = A->cmap->n;
608289bc588SBarry Smith   if (cols) {
609d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
610d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
611289bc588SBarry Smith   }
612289bc588SBarry Smith   if (vals) {
613d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
614289bc588SBarry Smith     v    = mat->v + row;
615d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
616289bc588SBarry Smith   }
6173a40ed3dSBarry Smith   PetscFunctionReturn(0);
618289bc588SBarry Smith }
6196ee01492SSatish Balay 
6204a2ae208SSatish Balay #undef __FUNCT__
6214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
62213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
623289bc588SBarry Smith {
624dfbe8321SBarry Smith   PetscErrorCode ierr;
625606d414cSSatish Balay   PetscFunctionBegin;
626606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
627606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6283a40ed3dSBarry Smith   PetscFunctionReturn(0);
629289bc588SBarry Smith }
630289bc588SBarry Smith /* ----------------------------------------------------------------*/
6314a2ae208SSatish Balay #undef __FUNCT__
6324a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
63313f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
634289bc588SBarry Smith {
635c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
63613f74950SBarry Smith   PetscInt     i,j,idx=0;
637d6dfbf8fSBarry Smith 
6383a40ed3dSBarry Smith   PetscFunctionBegin;
63971fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
640289bc588SBarry Smith   if (!mat->roworiented) {
641dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
642289bc588SBarry Smith       for (j=0; j<n; j++) {
643cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6442515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
645d0f46423SBarry 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);
64658804f6dSBarry Smith #endif
647289bc588SBarry Smith         for (i=0; i<m; i++) {
648cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6492515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
650d0f46423SBarry 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);
65158804f6dSBarry Smith #endif
652cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
653289bc588SBarry Smith         }
654289bc588SBarry Smith       }
6553a40ed3dSBarry Smith     } else {
656289bc588SBarry Smith       for (j=0; j<n; j++) {
657cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6582515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
659d0f46423SBarry 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);
66058804f6dSBarry Smith #endif
661289bc588SBarry Smith         for (i=0; i<m; i++) {
662cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6632515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
664d0f46423SBarry 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);
66558804f6dSBarry Smith #endif
666cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
667289bc588SBarry Smith         }
668289bc588SBarry Smith       }
669289bc588SBarry Smith     }
6703a40ed3dSBarry Smith   } else {
671dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
672e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
673cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6742515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
675d0f46423SBarry 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);
67658804f6dSBarry Smith #endif
677e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
678cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6792515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
680d0f46423SBarry 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);
68158804f6dSBarry Smith #endif
682cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
683e8d4e0b9SBarry Smith         }
684e8d4e0b9SBarry Smith       }
6853a40ed3dSBarry Smith     } else {
686289bc588SBarry Smith       for (i=0; i<m; i++) {
687cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6882515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
689d0f46423SBarry 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);
69058804f6dSBarry Smith #endif
691289bc588SBarry Smith         for (j=0; j<n; j++) {
692cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6932515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
694d0f46423SBarry 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);
69558804f6dSBarry Smith #endif
696cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
697289bc588SBarry Smith         }
698289bc588SBarry Smith       }
699289bc588SBarry Smith     }
700e8d4e0b9SBarry Smith   }
7013a40ed3dSBarry Smith   PetscFunctionReturn(0);
702289bc588SBarry Smith }
703e8d4e0b9SBarry Smith 
7044a2ae208SSatish Balay #undef __FUNCT__
7054a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
70613f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
707ae80bb75SLois Curfman McInnes {
708ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
70913f74950SBarry Smith   PetscInt     i,j;
710ae80bb75SLois Curfman McInnes 
7113a40ed3dSBarry Smith   PetscFunctionBegin;
712ae80bb75SLois Curfman McInnes   /* row-oriented output */
713ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
71497e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
715d0f46423SBarry 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);
716ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7176f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
718d0f46423SBarry 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);
71997e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
720ae80bb75SLois Curfman McInnes     }
721ae80bb75SLois Curfman McInnes   }
7223a40ed3dSBarry Smith   PetscFunctionReturn(0);
723ae80bb75SLois Curfman McInnes }
724ae80bb75SLois Curfman McInnes 
725289bc588SBarry Smith /* -----------------------------------------------------------------*/
726289bc588SBarry Smith 
7274a2ae208SSatish Balay #undef __FUNCT__
7284a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
729a313700dSBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, const MatType type,Mat *A)
73088e32edaSLois Curfman McInnes {
73188e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
73288e32edaSLois Curfman McInnes   Mat            B;
7336849ba73SBarry Smith   PetscErrorCode ierr;
73413f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
73513f74950SBarry Smith   int            fd;
73613f74950SBarry Smith   PetscMPIInt    size;
73713f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
73887828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
73919bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
74088e32edaSLois Curfman McInnes 
7413a40ed3dSBarry Smith   PetscFunctionBegin;
742d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
74329bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
744b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
7450752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
746*0700a824SBarry Smith   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
74788e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
74888e32edaSLois Curfman McInnes 
74990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
750f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
751f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
752be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7535c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
75490ace30eSBarry Smith     B    = *A;
75590ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
7564a41a4d3SSatish Balay     v    = a->v;
7574a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7584a41a4d3SSatish Balay        from row major to column major */
759b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
76090ace30eSBarry Smith     /* read in nonzero values */
7614a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7624a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7634a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7644a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7654a41a4d3SSatish Balay         *v++ =w[i*N+j];
7664a41a4d3SSatish Balay       }
7674a41a4d3SSatish Balay     }
768606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7696d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7706d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77190ace30eSBarry Smith   } else {
77288e32edaSLois Curfman McInnes     /* read row lengths */
77313f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7740752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
77588e32edaSLois Curfman McInnes 
77688e32edaSLois Curfman McInnes     /* create our matrix */
777f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
778f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
779be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7805c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
78188e32edaSLois Curfman McInnes     B = *A;
78288e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
78388e32edaSLois Curfman McInnes     v = a->v;
78488e32edaSLois Curfman McInnes 
78588e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
78613f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
787b0a32e0cSBarry Smith     cols = scols;
7880752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
78987828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
790b0a32e0cSBarry Smith     vals = svals;
7910752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
79288e32edaSLois Curfman McInnes 
79388e32edaSLois Curfman McInnes     /* insert into matrix */
79488e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
79588e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
79688e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
79788e32edaSLois Curfman McInnes     }
798606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
799606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
800606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
80188e32edaSLois Curfman McInnes 
8026d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8036d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80490ace30eSBarry Smith   }
8053a40ed3dSBarry Smith   PetscFunctionReturn(0);
80688e32edaSLois Curfman McInnes }
80788e32edaSLois 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);
829d0f46423SBarry 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);
832d0f46423SBarry 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;
854d0f46423SBarry 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);
860d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
861d0f46423SBarry 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 
865d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
866932b0c3eSLois Curfman McInnes       v = a->v + i;
867d0f46423SBarry 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;
897d0f46423SBarry 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);
9057973ad04SBarry Smith   if (format == PETSC_VIEWER_NATIVE) {
90690ace30eSBarry Smith     /* store the matrix as a dense matrix */
90713f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
908*0700a824SBarry Smith     col_lens[0] = MAT_FILE_CLASSID;
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);
927*0700a824SBarry Smith     col_lens[0] = MAT_FILE_CLASSID;
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;
967d0f46423SBarry 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);
1048d0f46423SBarry 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)
1089d0f46423SBarry 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;
1113d0f46423SBarry 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);
1133d0f46423SBarry 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;
1161d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1162d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1163d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
11641b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1165d0f46423SBarry 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);
1187d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1188d0f46423SBarry 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;
1203d0f46423SBarry 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);
1209d0f46423SBarry 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);
1221d0f46423SBarry 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;
1240d0f46423SBarry 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) {
1246d0f46423SBarry 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 {
1257d0f46423SBarry 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);
1266dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12673a40ed3dSBarry Smith   } else if (type == NORM_1) {
1268064f8208SBarry Smith     *nrm = 0.0;
1269d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
12701b807ce4Svictorle       v = mat->v + j*mat->lda;
1271289bc588SBarry Smith       sum = 0.0;
1272d0f46423SBarry 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     }
1277d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12783a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1279064f8208SBarry Smith     *nrm = 0.0;
1280d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1281289bc588SBarry Smith       v = mat->v + j;
1282289bc588SBarry Smith       sum = 0.0;
1283d0f46423SBarry 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     }
1288d0f46423SBarry 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;
1332d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13333a40ed3dSBarry Smith 
13343a40ed3dSBarry Smith   PetscFunctionBegin;
1335a5ce6ee0Svictorle   if (lda>m) {
1336d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1337a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1338a5ce6ee0Svictorle     }
1339a5ce6ee0Svictorle   } else {
1340d0f46423SBarry 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;
1350d0f46423SBarry 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;
1374d0f46423SBarry 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;
13945d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
13955d0c19d7SBarry Smith   const PetscInt *irow,*icol;
139687828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13970754003eSLois Curfman McInnes   Mat            newmat;
13980754003eSLois Curfman McInnes 
13993a40ed3dSBarry Smith   PetscFunctionBegin;
140078b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
140178b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1402e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1403e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14040754003eSLois Curfman McInnes 
1405182d2002SSatish Balay   /* Check submatrixcall */
1406182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
140713f74950SBarry Smith     PetscInt n_cols,n_rows;
1408182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
140921a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
141021a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
1411c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
141221a2c019SBarry Smith     }
1413182d2002SSatish Balay     newmat = *B;
1414182d2002SSatish Balay   } else {
14150754003eSLois Curfman McInnes     /* Create and fill new matrix */
14167adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1417f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14187adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14195c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1420182d2002SSatish Balay   }
1421182d2002SSatish Balay 
1422182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1423182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1424182d2002SSatish Balay 
1425182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14266de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1427182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1428182d2002SSatish Balay       *bv++ = av[irow[j]];
14290754003eSLois Curfman McInnes     }
14300754003eSLois Curfman McInnes   }
1431182d2002SSatish Balay 
1432182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14336d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14346d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14350754003eSLois Curfman McInnes 
14360754003eSLois Curfman McInnes   /* Free work space */
143778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
143878b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1439182d2002SSatish Balay   *B = newmat;
14403a40ed3dSBarry Smith   PetscFunctionReturn(0);
14410754003eSLois Curfman McInnes }
14420754003eSLois Curfman McInnes 
14434a2ae208SSatish Balay #undef __FUNCT__
14444a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
144513f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1446905e6a2fSBarry Smith {
14476849ba73SBarry Smith   PetscErrorCode ierr;
144813f74950SBarry Smith   PetscInt       i;
1449905e6a2fSBarry Smith 
14503a40ed3dSBarry Smith   PetscFunctionBegin;
1451905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1452b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1453905e6a2fSBarry Smith   }
1454905e6a2fSBarry Smith 
1455905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14566a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1457905e6a2fSBarry Smith   }
14583a40ed3dSBarry Smith   PetscFunctionReturn(0);
1459905e6a2fSBarry Smith }
1460905e6a2fSBarry Smith 
14614a2ae208SSatish Balay #undef __FUNCT__
1462c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1463c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1464c0aa2d19SHong Zhang {
1465c0aa2d19SHong Zhang   PetscFunctionBegin;
1466c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1467c0aa2d19SHong Zhang }
1468c0aa2d19SHong Zhang 
1469c0aa2d19SHong Zhang #undef __FUNCT__
1470c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1471c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1472c0aa2d19SHong Zhang {
1473c0aa2d19SHong Zhang   PetscFunctionBegin;
1474c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1475c0aa2d19SHong Zhang }
1476c0aa2d19SHong Zhang 
1477c0aa2d19SHong Zhang #undef __FUNCT__
14784a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1479dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14804b0e389bSBarry Smith {
14814b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14826849ba73SBarry Smith   PetscErrorCode ierr;
1483d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
14843a40ed3dSBarry Smith 
14853a40ed3dSBarry Smith   PetscFunctionBegin;
148633f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
148733f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1488cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14893a40ed3dSBarry Smith     PetscFunctionReturn(0);
14903a40ed3dSBarry Smith   }
1491d0f46423SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1492a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14930dbb7854Svictorle     for (j=0; j<n; j++) {
1494a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1495a5ce6ee0Svictorle     }
1496a5ce6ee0Svictorle   } else {
1497d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1498a5ce6ee0Svictorle   }
1499273d9f13SBarry Smith   PetscFunctionReturn(0);
1500273d9f13SBarry Smith }
1501273d9f13SBarry Smith 
15024a2ae208SSatish Balay #undef __FUNCT__
15034a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1504dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1505273d9f13SBarry Smith {
1506dfbe8321SBarry Smith   PetscErrorCode ierr;
1507273d9f13SBarry Smith 
1508273d9f13SBarry Smith   PetscFunctionBegin;
1509273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15103a40ed3dSBarry Smith   PetscFunctionReturn(0);
15114b0e389bSBarry Smith }
15124b0e389bSBarry Smith 
1513284134d9SBarry Smith #undef __FUNCT__
1514284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1515284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1516284134d9SBarry Smith {
1517284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
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);
1524dc5cefdeSJed Brown   A->rmap->n = A->rmap->N = m;
1525d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
152621a2c019SBarry Smith   if (a->changelda) a->lda = m;
1527284134d9SBarry Smith   PetscFunctionReturn(0);
1528284134d9SBarry Smith }
1529170fe5c8SBarry Smith 
1530ba337c44SJed Brown #undef __FUNCT__
1531ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense"
1532ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A)
1533ba337c44SJed Brown {
1534ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1535ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1536ba337c44SJed Brown   PetscScalar    *aa = a->v;
1537ba337c44SJed Brown 
1538ba337c44SJed Brown   PetscFunctionBegin;
1539ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
1540ba337c44SJed Brown   PetscFunctionReturn(0);
1541ba337c44SJed Brown }
1542ba337c44SJed Brown 
1543ba337c44SJed Brown #undef __FUNCT__
1544ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense"
1545ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A)
1546ba337c44SJed Brown {
1547ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1548ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1549ba337c44SJed Brown   PetscScalar    *aa = a->v;
1550ba337c44SJed Brown 
1551ba337c44SJed Brown   PetscFunctionBegin;
1552ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1553ba337c44SJed Brown   PetscFunctionReturn(0);
1554ba337c44SJed Brown }
1555ba337c44SJed Brown 
1556ba337c44SJed Brown #undef __FUNCT__
1557ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense"
1558ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
1559ba337c44SJed Brown {
1560ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1561ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1562ba337c44SJed Brown   PetscScalar    *aa = a->v;
1563ba337c44SJed Brown 
1564ba337c44SJed Brown   PetscFunctionBegin;
1565ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1566ba337c44SJed Brown   PetscFunctionReturn(0);
1567ba337c44SJed Brown }
1568284134d9SBarry Smith 
1569a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1570a9fe9ddaSSatish Balay #undef __FUNCT__
1571a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1572a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1573a9fe9ddaSSatish Balay {
1574a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1575a9fe9ddaSSatish Balay 
1576a9fe9ddaSSatish Balay   PetscFunctionBegin;
1577a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1578a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1579a9fe9ddaSSatish Balay   }
1580a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1581a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1582a9fe9ddaSSatish Balay }
1583a9fe9ddaSSatish Balay 
1584a9fe9ddaSSatish Balay #undef __FUNCT__
1585a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1586a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1587a9fe9ddaSSatish Balay {
1588ee16a9a1SHong Zhang   PetscErrorCode ierr;
1589d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1590ee16a9a1SHong Zhang   Mat            Cmat;
1591a9fe9ddaSSatish Balay 
1592ee16a9a1SHong Zhang   PetscFunctionBegin;
1593d0f46423SBarry 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);
1594ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1595ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1596ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1597ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1598ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1599ee16a9a1SHong Zhang   *C = Cmat;
1600ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1601ee16a9a1SHong Zhang }
1602a9fe9ddaSSatish Balay 
160398a3b096SSatish Balay #undef __FUNCT__
1604a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1605a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1606a9fe9ddaSSatish Balay {
1607a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1608a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1609a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16100805154bSBarry Smith   PetscBLASInt   m,n,k;
1611a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1612a9fe9ddaSSatish Balay 
1613a9fe9ddaSSatish Balay   PetscFunctionBegin;
1614d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1615d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1616d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1617a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1618a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1619a9fe9ddaSSatish Balay }
1620a9fe9ddaSSatish Balay 
1621a9fe9ddaSSatish Balay #undef __FUNCT__
1622a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1623a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1624a9fe9ddaSSatish Balay {
1625a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1626a9fe9ddaSSatish Balay 
1627a9fe9ddaSSatish Balay   PetscFunctionBegin;
1628a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1629a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1630a9fe9ddaSSatish Balay   }
1631a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1632a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1633a9fe9ddaSSatish Balay }
1634a9fe9ddaSSatish Balay 
1635a9fe9ddaSSatish Balay #undef __FUNCT__
1636a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1637a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1638a9fe9ddaSSatish Balay {
1639ee16a9a1SHong Zhang   PetscErrorCode ierr;
1640d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1641ee16a9a1SHong Zhang   Mat            Cmat;
1642a9fe9ddaSSatish Balay 
1643ee16a9a1SHong Zhang   PetscFunctionBegin;
1644d0f46423SBarry 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);
1645ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1646ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1647ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1648ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1649ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1650ee16a9a1SHong Zhang   *C = Cmat;
1651ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1652ee16a9a1SHong Zhang }
1653a9fe9ddaSSatish Balay 
1654a9fe9ddaSSatish Balay #undef __FUNCT__
1655a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1656a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1657a9fe9ddaSSatish Balay {
1658a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1659a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1660a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16610805154bSBarry Smith   PetscBLASInt   m,n,k;
1662a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1663a9fe9ddaSSatish Balay 
1664a9fe9ddaSSatish Balay   PetscFunctionBegin;
1665d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1666d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1667d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
16682fbe02b9SBarry Smith   /*
16692fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
16702fbe02b9SBarry Smith   */
1671a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1672a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1673a9fe9ddaSSatish Balay }
1674985db425SBarry Smith 
1675985db425SBarry Smith #undef __FUNCT__
1676985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1677985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1678985db425SBarry Smith {
1679985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1680985db425SBarry Smith   PetscErrorCode ierr;
1681d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1682985db425SBarry Smith   PetscScalar    *x;
1683985db425SBarry Smith   MatScalar      *aa = a->v;
1684985db425SBarry Smith 
1685985db425SBarry Smith   PetscFunctionBegin;
1686d5f3da31SBarry Smith   if (A->factortype) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1687985db425SBarry Smith 
1688985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1689985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1690985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1691d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1692985db425SBarry Smith   for (i=0; i<m; i++) {
1693985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1694985db425SBarry Smith     for (j=1; j<n; j++){
1695985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1696985db425SBarry Smith     }
1697985db425SBarry Smith   }
1698985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1699985db425SBarry Smith   PetscFunctionReturn(0);
1700985db425SBarry Smith }
1701985db425SBarry Smith 
1702985db425SBarry Smith #undef __FUNCT__
1703985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1704985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1705985db425SBarry Smith {
1706985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1707985db425SBarry Smith   PetscErrorCode ierr;
1708d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1709985db425SBarry Smith   PetscScalar    *x;
1710985db425SBarry Smith   PetscReal      atmp;
1711985db425SBarry Smith   MatScalar      *aa = a->v;
1712985db425SBarry Smith 
1713985db425SBarry Smith   PetscFunctionBegin;
1714d5f3da31SBarry Smith   if (A->factortype) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1715985db425SBarry Smith 
1716985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1717985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1718985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1719d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1720985db425SBarry Smith   for (i=0; i<m; i++) {
17219189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1722985db425SBarry Smith     for (j=1; j<n; j++){
1723985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1724985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1725985db425SBarry Smith     }
1726985db425SBarry Smith   }
1727985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1728985db425SBarry Smith   PetscFunctionReturn(0);
1729985db425SBarry Smith }
1730985db425SBarry Smith 
1731985db425SBarry Smith #undef __FUNCT__
1732985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1733985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1734985db425SBarry Smith {
1735985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1736985db425SBarry Smith   PetscErrorCode ierr;
1737d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1738985db425SBarry Smith   PetscScalar    *x;
1739985db425SBarry Smith   MatScalar      *aa = a->v;
1740985db425SBarry Smith 
1741985db425SBarry Smith   PetscFunctionBegin;
1742d5f3da31SBarry Smith   if (A->factortype) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1743985db425SBarry Smith 
1744985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1745985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1746985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1747d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1748985db425SBarry Smith   for (i=0; i<m; i++) {
1749985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1750985db425SBarry Smith     for (j=1; j<n; j++){
1751985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1752985db425SBarry Smith     }
1753985db425SBarry Smith   }
1754985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1755985db425SBarry Smith   PetscFunctionReturn(0);
1756985db425SBarry Smith }
1757985db425SBarry Smith 
17588d0534beSBarry Smith #undef __FUNCT__
17598d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
17608d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
17618d0534beSBarry Smith {
17628d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17638d0534beSBarry Smith   PetscErrorCode ierr;
17648d0534beSBarry Smith   PetscScalar    *x;
17658d0534beSBarry Smith 
17668d0534beSBarry Smith   PetscFunctionBegin;
1767d5f3da31SBarry Smith   if (A->factortype) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
17688d0534beSBarry Smith 
17698d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1770d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
17718d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
17728d0534beSBarry Smith   PetscFunctionReturn(0);
17738d0534beSBarry Smith }
17748d0534beSBarry Smith 
1775289bc588SBarry Smith /* -------------------------------------------------------------------*/
1776a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1777905e6a2fSBarry Smith        MatGetRow_SeqDense,
1778905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1779905e6a2fSBarry Smith        MatMult_SeqDense,
178097304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
17817c922b88SBarry Smith        MatMultTranspose_SeqDense,
17827c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1783db4efbfdSBarry Smith        0,
1784db4efbfdSBarry Smith        0,
1785db4efbfdSBarry Smith        0,
1786db4efbfdSBarry Smith /*10*/ 0,
1787905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1788905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
178941f059aeSBarry Smith        MatSOR_SeqDense,
1790ec8511deSBarry Smith        MatTranspose_SeqDense,
179197304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1792905e6a2fSBarry Smith        MatEqual_SeqDense,
1793905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1794905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1795905e6a2fSBarry Smith        MatNorm_SeqDense,
1796c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1797c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1798905e6a2fSBarry Smith        MatSetOption_SeqDense,
1799905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1800d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense,
1801db4efbfdSBarry Smith        0,
1802db4efbfdSBarry Smith        0,
1803db4efbfdSBarry Smith        0,
1804db4efbfdSBarry Smith        0,
1805d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense,
1806273d9f13SBarry Smith        0,
1807905e6a2fSBarry Smith        0,
1808905e6a2fSBarry Smith        MatGetArray_SeqDense,
1809905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1810d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense,
1811a5ae1ecdSBarry Smith        0,
1812a5ae1ecdSBarry Smith        0,
1813a5ae1ecdSBarry Smith        0,
1814a5ae1ecdSBarry Smith        0,
1815d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense,
1816a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1817a5ae1ecdSBarry Smith        0,
18184b0e389bSBarry Smith        MatGetValues_SeqDense,
1819a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1820d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense,
1821a5ae1ecdSBarry Smith        MatScale_SeqDense,
1822a5ae1ecdSBarry Smith        0,
1823a5ae1ecdSBarry Smith        0,
1824a5ae1ecdSBarry Smith        0,
1825d519adbfSMatthew Knepley /*49*/ 0,
1826a5ae1ecdSBarry Smith        0,
1827a5ae1ecdSBarry Smith        0,
1828a5ae1ecdSBarry Smith        0,
1829a5ae1ecdSBarry Smith        0,
1830d519adbfSMatthew Knepley /*54*/ 0,
1831a5ae1ecdSBarry Smith        0,
1832a5ae1ecdSBarry Smith        0,
1833a5ae1ecdSBarry Smith        0,
1834a5ae1ecdSBarry Smith        0,
1835d519adbfSMatthew Knepley /*59*/ 0,
1836e03a110bSBarry Smith        MatDestroy_SeqDense,
1837e03a110bSBarry Smith        MatView_SeqDense,
1838357abbc8SBarry Smith        0,
183997304618SKris Buschelman        0,
1840d519adbfSMatthew Knepley /*64*/ 0,
184197304618SKris Buschelman        0,
184297304618SKris Buschelman        0,
184397304618SKris Buschelman        0,
184497304618SKris Buschelman        0,
1845d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense,
184697304618SKris Buschelman        0,
184797304618SKris Buschelman        0,
184897304618SKris Buschelman        0,
184997304618SKris Buschelman        0,
1850d519adbfSMatthew Knepley /*74*/ 0,
185197304618SKris Buschelman        0,
185297304618SKris Buschelman        0,
185397304618SKris Buschelman        0,
185497304618SKris Buschelman        0,
1855d519adbfSMatthew Knepley /*79*/ 0,
185697304618SKris Buschelman        0,
185797304618SKris Buschelman        0,
185897304618SKris Buschelman        0,
1859d519adbfSMatthew Knepley /*83*/ MatLoad_SeqDense,
1860865e5f61SKris Buschelman        0,
18611cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1862865e5f61SKris Buschelman        0,
1863865e5f61SKris Buschelman        0,
1864865e5f61SKris Buschelman        0,
1865d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense,
1866a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1867a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1868865e5f61SKris Buschelman        0,
1869865e5f61SKris Buschelman        0,
1870d519adbfSMatthew Knepley /*94*/ 0,
1871a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1872a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1873a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1874284134d9SBarry Smith        0,
1875d519adbfSMatthew Knepley /*99*/ 0,
1876284134d9SBarry Smith        0,
1877284134d9SBarry Smith        0,
1878ba337c44SJed Brown        MatConjugate_SeqDense,
1879985db425SBarry Smith        MatSetSizes_SeqDense,
1880ba337c44SJed Brown /*104*/0,
1881ba337c44SJed Brown        MatRealPart_SeqDense,
1882ba337c44SJed Brown        MatImaginaryPart_SeqDense,
1883985db425SBarry Smith        0,
1884985db425SBarry Smith        0,
1885d519adbfSMatthew Knepley /*109*/0,
1886985db425SBarry Smith        0,
18878d0534beSBarry Smith        MatGetRowMin_SeqDense,
18888d0534beSBarry Smith        MatGetColumnVector_SeqDense
1889985db425SBarry Smith };
189090ace30eSBarry Smith 
18914a2ae208SSatish Balay #undef __FUNCT__
18924a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18934b828684SBarry Smith /*@C
1894fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1895d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1896d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1897289bc588SBarry Smith 
1898db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1899db81eaa0SLois Curfman McInnes 
190020563c6bSBarry Smith    Input Parameters:
1901db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
19020c775827SLois Curfman McInnes .  m - number of rows
190318f449edSLois Curfman McInnes .  n - number of columns
1904c0235b3cSMatthew Knepley -  data - optional location of matrix data in column major order.  Set data=PETSC_NULL for PETSc
1905dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
190620563c6bSBarry Smith 
190720563c6bSBarry Smith    Output Parameter:
190844cd7ae7SLois Curfman McInnes .  A - the matrix
190920563c6bSBarry Smith 
1910b259b22eSLois Curfman McInnes    Notes:
191118f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
191218f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1913b4fd4287SBarry Smith    set data=PETSC_NULL.
191418f449edSLois Curfman McInnes 
1915027ccd11SLois Curfman McInnes    Level: intermediate
1916027ccd11SLois Curfman McInnes 
1917dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1918d65003e9SLois Curfman McInnes 
1919db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
192020563c6bSBarry Smith @*/
1921be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1922289bc588SBarry Smith {
1923dfbe8321SBarry Smith   PetscErrorCode ierr;
19243b2fbd54SBarry Smith 
19253a40ed3dSBarry Smith   PetscFunctionBegin;
1926f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1927f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1928273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1929273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1930273d9f13SBarry Smith   PetscFunctionReturn(0);
1931273d9f13SBarry Smith }
1932273d9f13SBarry Smith 
19334a2ae208SSatish Balay #undef __FUNCT__
1934afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1935273d9f13SBarry Smith /*@C
1936273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1937273d9f13SBarry Smith 
1938273d9f13SBarry Smith    Collective on MPI_Comm
1939273d9f13SBarry Smith 
1940273d9f13SBarry Smith    Input Parameters:
1941273d9f13SBarry Smith +  A - the matrix
1942273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1943273d9f13SBarry Smith 
1944273d9f13SBarry Smith    Notes:
1945273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1946273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1947284134d9SBarry Smith    need not call this routine.
1948273d9f13SBarry Smith 
1949273d9f13SBarry Smith    Level: intermediate
1950273d9f13SBarry Smith 
1951273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1952273d9f13SBarry Smith 
1953867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
1954867c911aSBarry Smith 
1955273d9f13SBarry Smith @*/
1956be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1957273d9f13SBarry Smith {
1958dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1959a23d5eceSKris Buschelman 
1960a23d5eceSKris Buschelman   PetscFunctionBegin;
1961a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1962a23d5eceSKris Buschelman   if (f) {
1963a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1964a23d5eceSKris Buschelman   }
1965a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1966a23d5eceSKris Buschelman }
1967a23d5eceSKris Buschelman 
1968a23d5eceSKris Buschelman EXTERN_C_BEGIN
1969a23d5eceSKris Buschelman #undef __FUNCT__
1970afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1971be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1972a23d5eceSKris Buschelman {
1973273d9f13SBarry Smith   Mat_SeqDense   *b;
1974dfbe8321SBarry Smith   PetscErrorCode ierr;
1975273d9f13SBarry Smith 
1976273d9f13SBarry Smith   PetscFunctionBegin;
1977273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1978273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1979d0f46423SBarry Smith   if (b->lda <= 0) b->lda = B->rmap->n;
19809e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
19819e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
19825afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1983284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1984284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
19859e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
1986273d9f13SBarry Smith   } else { /* user-allocated storage */
19879e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1988273d9f13SBarry Smith     b->v          = data;
1989273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1990273d9f13SBarry Smith   }
19910450473dSBarry Smith   B->assembled = PETSC_TRUE;
1992273d9f13SBarry Smith   PetscFunctionReturn(0);
1993273d9f13SBarry Smith }
1994a23d5eceSKris Buschelman EXTERN_C_END
1995273d9f13SBarry Smith 
19961b807ce4Svictorle #undef __FUNCT__
19971b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
19981b807ce4Svictorle /*@C
19991b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
20001b807ce4Svictorle 
20011b807ce4Svictorle   Input parameter:
20021b807ce4Svictorle + A - the matrix
20031b807ce4Svictorle - lda - the leading dimension
20041b807ce4Svictorle 
20051b807ce4Svictorle   Notes:
2006867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
20071b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
20081b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
20091b807ce4Svictorle 
20101b807ce4Svictorle   Level: intermediate
20111b807ce4Svictorle 
20121b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
20131b807ce4Svictorle 
2014284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2015867c911aSBarry Smith 
20161b807ce4Svictorle @*/
2017be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
20181b807ce4Svictorle {
20191b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
202021a2c019SBarry Smith 
20211b807ce4Svictorle   PetscFunctionBegin;
2022d0f46423SBarry Smith   if (lda < B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
20231b807ce4Svictorle   b->lda       = lda;
202421a2c019SBarry Smith   b->changelda = PETSC_FALSE;
202521a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
20261b807ce4Svictorle   PetscFunctionReturn(0);
20271b807ce4Svictorle }
20281b807ce4Svictorle 
20290bad9183SKris Buschelman /*MC
2030fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
20310bad9183SKris Buschelman 
20320bad9183SKris Buschelman    Options Database Keys:
20330bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
20340bad9183SKris Buschelman 
20350bad9183SKris Buschelman   Level: beginner
20360bad9183SKris Buschelman 
203789665df3SBarry Smith .seealso: MatCreateSeqDense()
203889665df3SBarry Smith 
20390bad9183SKris Buschelman M*/
20400bad9183SKris Buschelman 
2041273d9f13SBarry Smith EXTERN_C_BEGIN
20424a2ae208SSatish Balay #undef __FUNCT__
20434a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
2044be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2045273d9f13SBarry Smith {
2046273d9f13SBarry Smith   Mat_SeqDense   *b;
2047dfbe8321SBarry Smith   PetscErrorCode ierr;
20487c334f02SBarry Smith   PetscMPIInt    size;
2049273d9f13SBarry Smith 
2050273d9f13SBarry Smith   PetscFunctionBegin;
20517adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
205229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
205355659b69SBarry Smith 
205426283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
205526283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
205626283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
205726283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
2058273d9f13SBarry Smith 
205938f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2060549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
206190f02eecSBarry Smith   B->mapping      = 0;
206244cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
206318f449edSLois Curfman McInnes 
2064a5ae1ecdSBarry Smith 
206544cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2066273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2067273d9f13SBarry Smith   b->v            = 0;
2068d0f46423SBarry Smith   b->lda          = B->rmap->n;
206921a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
2070d0f46423SBarry Smith   b->Mmax         = B->rmap->n;
2071d0f46423SBarry Smith   b->Nmax         = B->cmap->n;
20724e220ebcSLois Curfman McInnes 
2073b24902e0SBarry Smith 
2074ec1065edSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
2075b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2076b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2077a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2078a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2079a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
20804ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
20814ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
20824ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
20834ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
20844ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
20854ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
20864ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
20874ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
20884ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
208917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
20903a40ed3dSBarry Smith   PetscFunctionReturn(0);
2091289bc588SBarry Smith }
2092c0aa2d19SHong Zhang 
2093c0aa2d19SHong Zhang 
2094273d9f13SBarry Smith EXTERN_C_END
2095