xref: /petsc/src/mat/impls/dense/seq/dense.c (revision dc5cefde2fdcdb92db1db0a399b233ecd778849b)
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"
8b4c862e0SSatish Balay #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
170805154bSBarry Smith   PetscBLASInt   N,m,ldax,lday,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
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);
1695c9eb25fSBarry Smith   if (A->factor == 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
1765c9eb25fSBarry Smith   } else if (A->factor == 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;
348db4efbfdSBarry Smith   A->factor = 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;
378db4efbfdSBarry Smith   A->factor = 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   }
435db4efbfdSBarry Smith   (*fact)->factor = 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 
727e090d566SSatish Balay #include "petscsys.h"
72888e32edaSLois Curfman McInnes 
7294a2ae208SSatish Balay #undef __FUNCT__
7304a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
731a313700dSBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, const MatType type,Mat *A)
73288e32edaSLois Curfman McInnes {
73388e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
73488e32edaSLois Curfman McInnes   Mat            B;
7356849ba73SBarry Smith   PetscErrorCode ierr;
73613f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
73713f74950SBarry Smith   int            fd;
73813f74950SBarry Smith   PetscMPIInt    size;
73913f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
74087828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
74119bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
74288e32edaSLois Curfman McInnes 
7433a40ed3dSBarry Smith   PetscFunctionBegin;
744d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
74529bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
746b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
7470752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
748552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
74988e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
75088e32edaSLois Curfman McInnes 
75190ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
752f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
753f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
754be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7555c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
75690ace30eSBarry Smith     B    = *A;
75790ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
7584a41a4d3SSatish Balay     v    = a->v;
7594a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7604a41a4d3SSatish Balay        from row major to column major */
761b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
76290ace30eSBarry Smith     /* read in nonzero values */
7634a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7644a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7654a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7664a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7674a41a4d3SSatish Balay         *v++ =w[i*N+j];
7684a41a4d3SSatish Balay       }
7694a41a4d3SSatish Balay     }
770606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7716d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7726d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77390ace30eSBarry Smith   } else {
77488e32edaSLois Curfman McInnes     /* read row lengths */
77513f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7760752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
77788e32edaSLois Curfman McInnes 
77888e32edaSLois Curfman McInnes     /* create our matrix */
779f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
780f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
781be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7825c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
78388e32edaSLois Curfman McInnes     B = *A;
78488e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
78588e32edaSLois Curfman McInnes     v = a->v;
78688e32edaSLois Curfman McInnes 
78788e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
78813f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
789b0a32e0cSBarry Smith     cols = scols;
7900752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
79187828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
792b0a32e0cSBarry Smith     vals = svals;
7930752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
79488e32edaSLois Curfman McInnes 
79588e32edaSLois Curfman McInnes     /* insert into matrix */
79688e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
79788e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
79888e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
79988e32edaSLois Curfman McInnes     }
800606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
801606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
802606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
80388e32edaSLois Curfman McInnes 
8046d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8056d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80690ace30eSBarry Smith   }
8073a40ed3dSBarry Smith   PetscFunctionReturn(0);
80888e32edaSLois Curfman McInnes }
80988e32edaSLois Curfman McInnes 
810e090d566SSatish Balay #include "petscsys.h"
811932b0c3eSLois Curfman McInnes 
8124a2ae208SSatish Balay #undef __FUNCT__
8134a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8146849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
815289bc588SBarry Smith {
816932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
817dfbe8321SBarry Smith   PetscErrorCode    ierr;
81813f74950SBarry Smith   PetscInt          i,j;
8192dcb1b2aSMatthew Knepley   const char        *name;
82087828ca2SBarry Smith   PetscScalar       *v;
821f3ef73ceSBarry Smith   PetscViewerFormat format;
8225f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
8235f481a85SSatish Balay   PetscTruth        allreal = PETSC_TRUE;
8245f481a85SSatish Balay #endif
825932b0c3eSLois Curfman McInnes 
8263a40ed3dSBarry Smith   PetscFunctionBegin;
827435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
828b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
829456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8303a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
831fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
832b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
833d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
83444cd7ae7SLois Curfman McInnes       v = a->v + i;
83577431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
836d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
837aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
838329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
839a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
840329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
841a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8426831982aSBarry Smith         }
84380cd9d93SLois Curfman McInnes #else
8446831982aSBarry Smith         if (*v) {
845a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8466831982aSBarry Smith         }
84780cd9d93SLois Curfman McInnes #endif
8481b807ce4Svictorle         v += a->lda;
84980cd9d93SLois Curfman McInnes       }
850b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
85180cd9d93SLois Curfman McInnes     }
852b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
8533a40ed3dSBarry Smith   } else {
854b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
855aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
85647989497SBarry Smith     /* determine if matrix has all real values */
85747989497SBarry Smith     v = a->v;
858d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
859ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
86047989497SBarry Smith     }
86147989497SBarry Smith #endif
862fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8633a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
864d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
865d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
866fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
867ffac6cdbSBarry Smith     }
868ffac6cdbSBarry Smith 
869d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
870932b0c3eSLois Curfman McInnes       v = a->v + i;
871d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
872aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
87347989497SBarry Smith         if (allreal) {
874f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
87547989497SBarry Smith         } else {
876f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
87747989497SBarry Smith         }
878289bc588SBarry Smith #else
879f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
880289bc588SBarry Smith #endif
8811b807ce4Svictorle 	v += a->lda;
882289bc588SBarry Smith       }
883b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
884289bc588SBarry Smith     }
885fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
886b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
887ffac6cdbSBarry Smith     }
888b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
889da3a660dSBarry Smith   }
890b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8913a40ed3dSBarry Smith   PetscFunctionReturn(0);
892289bc588SBarry Smith }
893289bc588SBarry Smith 
8944a2ae208SSatish Balay #undef __FUNCT__
8954a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8966849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
897932b0c3eSLois Curfman McInnes {
898932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8996849ba73SBarry Smith   PetscErrorCode    ierr;
90013f74950SBarry Smith   int               fd;
901d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
90287828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
903f3ef73ceSBarry Smith   PetscViewerFormat format;
904932b0c3eSLois Curfman McInnes 
9053a40ed3dSBarry Smith   PetscFunctionBegin;
906b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
90790ace30eSBarry Smith 
908b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
9097973ad04SBarry Smith   if (format == PETSC_VIEWER_NATIVE) {
91090ace30eSBarry Smith     /* store the matrix as a dense matrix */
91113f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
912552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
91390ace30eSBarry Smith     col_lens[1] = m;
91490ace30eSBarry Smith     col_lens[2] = n;
91590ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
9166f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
917606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
91890ace30eSBarry Smith 
91990ace30eSBarry Smith     /* write out matrix, by rows */
92087828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
92190ace30eSBarry Smith     v    = a->v;
92290ace30eSBarry Smith     for (j=0; j<n; j++) {
923578230a0SSatish Balay       for (i=0; i<m; i++) {
924578230a0SSatish Balay         vals[j + i*n] = *v++;
92590ace30eSBarry Smith       }
92690ace30eSBarry Smith     }
9276f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
928606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
92990ace30eSBarry Smith   } else {
93013f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
931552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
932932b0c3eSLois Curfman McInnes     col_lens[1] = m;
933932b0c3eSLois Curfman McInnes     col_lens[2] = n;
934932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
935932b0c3eSLois Curfman McInnes 
936932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
937932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9386f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
939932b0c3eSLois Curfman McInnes 
940932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
941932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
942932b0c3eSLois Curfman McInnes     ict = 0;
943932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
944932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
945932b0c3eSLois Curfman McInnes     }
9466f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
947606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
948932b0c3eSLois Curfman McInnes 
949932b0c3eSLois Curfman McInnes     /* store nonzero values */
95087828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
951932b0c3eSLois Curfman McInnes     ict  = 0;
952932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
953932b0c3eSLois Curfman McInnes       v = a->v + i;
954932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9551b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
956932b0c3eSLois Curfman McInnes       }
957932b0c3eSLois Curfman McInnes     }
9586f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
959606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
96090ace30eSBarry Smith   }
9613a40ed3dSBarry Smith   PetscFunctionReturn(0);
962932b0c3eSLois Curfman McInnes }
963932b0c3eSLois Curfman McInnes 
9644a2ae208SSatish Balay #undef __FUNCT__
9654a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
966dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
967f1af5d2fSBarry Smith {
968f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
969f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9706849ba73SBarry Smith   PetscErrorCode    ierr;
971d0f46423SBarry Smith   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
97287828ca2SBarry Smith   PetscScalar       *v = a->v;
973b0a32e0cSBarry Smith   PetscViewer       viewer;
974b0a32e0cSBarry Smith   PetscDraw         popup;
975329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
976f3ef73ceSBarry Smith   PetscViewerFormat format;
977f1af5d2fSBarry Smith 
978f1af5d2fSBarry Smith   PetscFunctionBegin;
979f1af5d2fSBarry Smith 
980f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
981b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
982b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
983f1af5d2fSBarry Smith 
984f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
985fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
986f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
987b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
988f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
989f1af5d2fSBarry Smith       x_l = j;
990f1af5d2fSBarry Smith       x_r = x_l + 1.0;
991f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
992f1af5d2fSBarry Smith         y_l = m - i - 1.0;
993f1af5d2fSBarry Smith         y_r = y_l + 1.0;
994f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
995329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
996b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
997329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
998b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
999f1af5d2fSBarry Smith         } else {
1000f1af5d2fSBarry Smith           continue;
1001f1af5d2fSBarry Smith         }
1002f1af5d2fSBarry Smith #else
1003f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
1004b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1005f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1006b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1007f1af5d2fSBarry Smith         } else {
1008f1af5d2fSBarry Smith           continue;
1009f1af5d2fSBarry Smith         }
1010f1af5d2fSBarry Smith #endif
1011b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1012f1af5d2fSBarry Smith       }
1013f1af5d2fSBarry Smith     }
1014f1af5d2fSBarry Smith   } else {
1015f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1016f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1017f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1018f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1019f1af5d2fSBarry Smith     }
1020b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1021b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1022b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1023f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1024f1af5d2fSBarry Smith       x_l = j;
1025f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1026f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1027f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1028f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1029b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1030b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1031f1af5d2fSBarry Smith       }
1032f1af5d2fSBarry Smith     }
1033f1af5d2fSBarry Smith   }
1034f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1035f1af5d2fSBarry Smith }
1036f1af5d2fSBarry Smith 
10374a2ae208SSatish Balay #undef __FUNCT__
10384a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1039dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1040f1af5d2fSBarry Smith {
1041b0a32e0cSBarry Smith   PetscDraw      draw;
1042f1af5d2fSBarry Smith   PetscTruth     isnull;
1043329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1044dfbe8321SBarry Smith   PetscErrorCode ierr;
1045f1af5d2fSBarry Smith 
1046f1af5d2fSBarry Smith   PetscFunctionBegin;
1047b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1048b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1049abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1050f1af5d2fSBarry Smith 
1051f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1052d0f46423SBarry Smith   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1053f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1054b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1055b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1056f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1057f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1058f1af5d2fSBarry Smith }
1059f1af5d2fSBarry Smith 
10604a2ae208SSatish Balay #undef __FUNCT__
10614a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1062dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1063932b0c3eSLois Curfman McInnes {
1064dfbe8321SBarry Smith   PetscErrorCode ierr;
10656805f65bSBarry Smith   PetscTruth     iascii,isbinary,isdraw;
1066932b0c3eSLois Curfman McInnes 
10673a40ed3dSBarry Smith   PetscFunctionBegin;
106832077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1069fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1070fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10710f5bd95cSBarry Smith 
1072c45a1595SBarry Smith   if (iascii) {
1073c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10740f5bd95cSBarry Smith   } else if (isbinary) {
10753a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1076f1af5d2fSBarry Smith   } else if (isdraw) {
1077f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10785cd90555SBarry Smith   } else {
1079958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1080932b0c3eSLois Curfman McInnes   }
10813a40ed3dSBarry Smith   PetscFunctionReturn(0);
1082932b0c3eSLois Curfman McInnes }
1083289bc588SBarry Smith 
10844a2ae208SSatish Balay #undef __FUNCT__
10854a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1086dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1087289bc588SBarry Smith {
1088ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1089dfbe8321SBarry Smith   PetscErrorCode ierr;
109090f02eecSBarry Smith 
10913a40ed3dSBarry Smith   PetscFunctionBegin;
1092aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1093d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1094a5a9c739SBarry Smith #endif
109505b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
10966857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1097606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1098dbd8c25aSHong Zhang 
1099dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1100901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
11014ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11024ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11034ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11043a40ed3dSBarry Smith   PetscFunctionReturn(0);
1105289bc588SBarry Smith }
1106289bc588SBarry Smith 
11074a2ae208SSatish Balay #undef __FUNCT__
11084a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1109fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1110289bc588SBarry Smith {
1111c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11126849ba73SBarry Smith   PetscErrorCode ierr;
111313f74950SBarry Smith   PetscInt       k,j,m,n,M;
111487828ca2SBarry Smith   PetscScalar    *v,tmp;
111548b35521SBarry Smith 
11163a40ed3dSBarry Smith   PetscFunctionBegin;
1117d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1118e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1119a5ce6ee0Svictorle     if (m != n) {
1120634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1121d64ed03dSBarry Smith     } else {
1122d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1123289bc588SBarry Smith         for (k=0; k<j; k++) {
11241b807ce4Svictorle           tmp = v[j + k*M];
11251b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11261b807ce4Svictorle           v[k + j*M] = tmp;
1127289bc588SBarry Smith         }
1128289bc588SBarry Smith       }
1129d64ed03dSBarry Smith     }
11303a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1131d3e5ee88SLois Curfman McInnes     Mat          tmat;
1132ec8511deSBarry Smith     Mat_SeqDense *tmatd;
113387828ca2SBarry Smith     PetscScalar  *v2;
1134ea709b57SSatish Balay 
1135fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11367adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1137d0f46423SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
11387adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11395c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1140fc4dec0aSBarry Smith     } else {
1141fc4dec0aSBarry Smith       tmat = *matout;
1142fc4dec0aSBarry Smith     }
1143ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11440de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1145d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11461b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1147d3e5ee88SLois Curfman McInnes     }
11486d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11496d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1150d3e5ee88SLois Curfman McInnes     *matout = tmat;
115148b35521SBarry Smith   }
11523a40ed3dSBarry Smith   PetscFunctionReturn(0);
1153289bc588SBarry Smith }
1154289bc588SBarry Smith 
11554a2ae208SSatish Balay #undef __FUNCT__
11564a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1157dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1158289bc588SBarry Smith {
1159c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1160c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
116113f74950SBarry Smith   PetscInt     i,j;
116287828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11639ea5d5aeSSatish Balay 
11643a40ed3dSBarry Smith   PetscFunctionBegin;
1165d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1166d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1167d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
11681b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1169d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
11703a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11711b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11721b807ce4Svictorle     }
1173289bc588SBarry Smith   }
117477c4ece6SBarry Smith   *flg = PETSC_TRUE;
11753a40ed3dSBarry Smith   PetscFunctionReturn(0);
1176289bc588SBarry Smith }
1177289bc588SBarry Smith 
11784a2ae208SSatish Balay #undef __FUNCT__
11794a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1180dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1181289bc588SBarry Smith {
1182c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1183dfbe8321SBarry Smith   PetscErrorCode ierr;
118413f74950SBarry Smith   PetscInt       i,n,len;
118587828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
118644cd7ae7SLois Curfman McInnes 
11873a40ed3dSBarry Smith   PetscFunctionBegin;
11882dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11897a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11901ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1191d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1192d0f46423SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
119344cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11941b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1195289bc588SBarry Smith   }
11961ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1198289bc588SBarry Smith }
1199289bc588SBarry Smith 
12004a2ae208SSatish Balay #undef __FUNCT__
12014a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1202dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1203289bc588SBarry Smith {
1204c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
120587828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1206dfbe8321SBarry Smith   PetscErrorCode ierr;
1207d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
120855659b69SBarry Smith 
12093a40ed3dSBarry Smith   PetscFunctionBegin;
121028988994SBarry Smith   if (ll) {
12117a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12121ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1213d0f46423SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1214da3a660dSBarry Smith     for (i=0; i<m; i++) {
1215da3a660dSBarry Smith       x = l[i];
1216da3a660dSBarry Smith       v = mat->v + i;
1217da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1218da3a660dSBarry Smith     }
12191ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1220efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1221da3a660dSBarry Smith   }
122228988994SBarry Smith   if (rr) {
12237a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12241ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1225d0f46423SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1226da3a660dSBarry Smith     for (i=0; i<n; i++) {
1227da3a660dSBarry Smith       x = r[i];
1228da3a660dSBarry Smith       v = mat->v + i*m;
1229da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1230da3a660dSBarry Smith     }
12311ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1232efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1233da3a660dSBarry Smith   }
12343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1235289bc588SBarry Smith }
1236289bc588SBarry Smith 
12374a2ae208SSatish Balay #undef __FUNCT__
12384a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1239dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1240289bc588SBarry Smith {
1241c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
124287828ca2SBarry Smith   PetscScalar  *v = mat->v;
1243329f5518SBarry Smith   PetscReal    sum = 0.0;
1244d0f46423SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1245efee365bSSatish Balay   PetscErrorCode ierr;
124655659b69SBarry Smith 
12473a40ed3dSBarry Smith   PetscFunctionBegin;
1248289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1249a5ce6ee0Svictorle     if (lda>m) {
1250d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1251a5ce6ee0Svictorle 	v = mat->v+j*lda;
1252a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1253a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1254a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1255a5ce6ee0Svictorle #else
1256a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1257a5ce6ee0Svictorle #endif
1258a5ce6ee0Svictorle 	}
1259a5ce6ee0Svictorle       }
1260a5ce6ee0Svictorle     } else {
1261d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1262aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1263329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1264289bc588SBarry Smith #else
1265289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1266289bc588SBarry Smith #endif
1267289bc588SBarry Smith       }
1268a5ce6ee0Svictorle     }
1269064f8208SBarry Smith     *nrm = sqrt(sum);
1270dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12713a40ed3dSBarry Smith   } else if (type == NORM_1) {
1272064f8208SBarry Smith     *nrm = 0.0;
1273d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
12741b807ce4Svictorle       v = mat->v + j*mat->lda;
1275289bc588SBarry Smith       sum = 0.0;
1276d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
127733a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1278289bc588SBarry Smith       }
1279064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1280289bc588SBarry Smith     }
1281d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12823a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1283064f8208SBarry Smith     *nrm = 0.0;
1284d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1285289bc588SBarry Smith       v = mat->v + j;
1286289bc588SBarry Smith       sum = 0.0;
1287d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
12881b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1289289bc588SBarry Smith       }
1290064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1291289bc588SBarry Smith     }
1292d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12933a40ed3dSBarry Smith   } else {
129429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1295289bc588SBarry Smith   }
12963a40ed3dSBarry Smith   PetscFunctionReturn(0);
1297289bc588SBarry Smith }
1298289bc588SBarry Smith 
12994a2ae208SSatish Balay #undef __FUNCT__
13004a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
13014e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg)
1302289bc588SBarry Smith {
1303c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
130463ba0a88SBarry Smith   PetscErrorCode ierr;
130567e560aaSBarry Smith 
13063a40ed3dSBarry Smith   PetscFunctionBegin;
1307b5a2b587SKris Buschelman   switch (op) {
1308b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13094e0d8c25SBarry Smith     aij->roworiented = flg;
1310b5a2b587SKris Buschelman     break;
1311512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1312b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13133971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13144e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1315b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1316b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
131777e54ba9SKris Buschelman   case MAT_SYMMETRIC:
131877e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13199a4540c5SBarry Smith   case MAT_HERMITIAN:
13209a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1321600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1322290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
132377e54ba9SKris Buschelman     break;
1324b5a2b587SKris Buschelman   default:
1325600fe468SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13263a40ed3dSBarry Smith   }
13273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1328289bc588SBarry Smith }
1329289bc588SBarry Smith 
13304a2ae208SSatish Balay #undef __FUNCT__
13314a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1332dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13336f0a148fSBarry Smith {
1334ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13356849ba73SBarry Smith   PetscErrorCode ierr;
1336d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13373a40ed3dSBarry Smith 
13383a40ed3dSBarry Smith   PetscFunctionBegin;
1339a5ce6ee0Svictorle   if (lda>m) {
1340d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1341a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1342a5ce6ee0Svictorle     }
1343a5ce6ee0Svictorle   } else {
1344d0f46423SBarry Smith     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1345a5ce6ee0Svictorle   }
13463a40ed3dSBarry Smith   PetscFunctionReturn(0);
13476f0a148fSBarry Smith }
13486f0a148fSBarry Smith 
13494a2ae208SSatish Balay #undef __FUNCT__
13504a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1351f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13526f0a148fSBarry Smith {
1353ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1354d0f46423SBarry Smith   PetscInt       n = A->cmap->n,i,j;
135587828ca2SBarry Smith   PetscScalar    *slot;
135655659b69SBarry Smith 
13573a40ed3dSBarry Smith   PetscFunctionBegin;
13586f0a148fSBarry Smith   for (i=0; i<N; i++) {
13596f0a148fSBarry Smith     slot = l->v + rows[i];
13606f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13616f0a148fSBarry Smith   }
1362f4df32b1SMatthew Knepley   if (diag != 0.0) {
13636f0a148fSBarry Smith     for (i=0; i<N; i++) {
13646f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1365f4df32b1SMatthew Knepley       *slot = diag;
13666f0a148fSBarry Smith     }
13676f0a148fSBarry Smith   }
13683a40ed3dSBarry Smith   PetscFunctionReturn(0);
13696f0a148fSBarry Smith }
1370557bce09SLois Curfman McInnes 
13714a2ae208SSatish Balay #undef __FUNCT__
13724a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1373dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
137464e87e97SBarry Smith {
1375c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13763a40ed3dSBarry Smith 
13773a40ed3dSBarry Smith   PetscFunctionBegin;
1378d0f46423SBarry Smith   if (mat->lda != A->rmap->n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
137964e87e97SBarry Smith   *array = mat->v;
13803a40ed3dSBarry Smith   PetscFunctionReturn(0);
138164e87e97SBarry Smith }
13820754003eSLois Curfman McInnes 
13834a2ae208SSatish Balay #undef __FUNCT__
13844a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1385dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1386ff14e315SSatish Balay {
13873a40ed3dSBarry Smith   PetscFunctionBegin;
138809b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13893a40ed3dSBarry Smith   PetscFunctionReturn(0);
1390ff14e315SSatish Balay }
13910754003eSLois Curfman McInnes 
13924a2ae208SSatish Balay #undef __FUNCT__
13934a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
139413f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13950754003eSLois Curfman McInnes {
1396c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13976849ba73SBarry Smith   PetscErrorCode ierr;
13985d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
13995d0c19d7SBarry Smith   const PetscInt *irow,*icol;
140087828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
14010754003eSLois Curfman McInnes   Mat            newmat;
14020754003eSLois Curfman McInnes 
14033a40ed3dSBarry Smith   PetscFunctionBegin;
140478b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
140578b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1406e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1407e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14080754003eSLois Curfman McInnes 
1409182d2002SSatish Balay   /* Check submatrixcall */
1410182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
141113f74950SBarry Smith     PetscInt n_cols,n_rows;
1412182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
141321a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
141421a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
1415c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
141621a2c019SBarry Smith     }
1417182d2002SSatish Balay     newmat = *B;
1418182d2002SSatish Balay   } else {
14190754003eSLois Curfman McInnes     /* Create and fill new matrix */
14207adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1421f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14227adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14235c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1424182d2002SSatish Balay   }
1425182d2002SSatish Balay 
1426182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1427182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1428182d2002SSatish Balay 
1429182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14306de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1431182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1432182d2002SSatish Balay       *bv++ = av[irow[j]];
14330754003eSLois Curfman McInnes     }
14340754003eSLois Curfman McInnes   }
1435182d2002SSatish Balay 
1436182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14376d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14386d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14390754003eSLois Curfman McInnes 
14400754003eSLois Curfman McInnes   /* Free work space */
144178b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
144278b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1443182d2002SSatish Balay   *B = newmat;
14443a40ed3dSBarry Smith   PetscFunctionReturn(0);
14450754003eSLois Curfman McInnes }
14460754003eSLois Curfman McInnes 
14474a2ae208SSatish Balay #undef __FUNCT__
14484a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
144913f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1450905e6a2fSBarry Smith {
14516849ba73SBarry Smith   PetscErrorCode ierr;
145213f74950SBarry Smith   PetscInt       i;
1453905e6a2fSBarry Smith 
14543a40ed3dSBarry Smith   PetscFunctionBegin;
1455905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1456b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1457905e6a2fSBarry Smith   }
1458905e6a2fSBarry Smith 
1459905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14606a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1461905e6a2fSBarry Smith   }
14623a40ed3dSBarry Smith   PetscFunctionReturn(0);
1463905e6a2fSBarry Smith }
1464905e6a2fSBarry Smith 
14654a2ae208SSatish Balay #undef __FUNCT__
1466c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1467c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1468c0aa2d19SHong Zhang {
1469c0aa2d19SHong Zhang   PetscFunctionBegin;
1470c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1471c0aa2d19SHong Zhang }
1472c0aa2d19SHong Zhang 
1473c0aa2d19SHong Zhang #undef __FUNCT__
1474c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1475c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1476c0aa2d19SHong Zhang {
1477c0aa2d19SHong Zhang   PetscFunctionBegin;
1478c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1479c0aa2d19SHong Zhang }
1480c0aa2d19SHong Zhang 
1481c0aa2d19SHong Zhang #undef __FUNCT__
14824a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1483dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14844b0e389bSBarry Smith {
14854b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14866849ba73SBarry Smith   PetscErrorCode ierr;
1487d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
14883a40ed3dSBarry Smith 
14893a40ed3dSBarry Smith   PetscFunctionBegin;
149033f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
149133f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1492cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14933a40ed3dSBarry Smith     PetscFunctionReturn(0);
14943a40ed3dSBarry Smith   }
1495d0f46423SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1496a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14970dbb7854Svictorle     for (j=0; j<n; j++) {
1498a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1499a5ce6ee0Svictorle     }
1500a5ce6ee0Svictorle   } else {
1501d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1502a5ce6ee0Svictorle   }
1503273d9f13SBarry Smith   PetscFunctionReturn(0);
1504273d9f13SBarry Smith }
1505273d9f13SBarry Smith 
15064a2ae208SSatish Balay #undef __FUNCT__
15074a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1508dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1509273d9f13SBarry Smith {
1510dfbe8321SBarry Smith   PetscErrorCode ierr;
1511273d9f13SBarry Smith 
1512273d9f13SBarry Smith   PetscFunctionBegin;
1513273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15143a40ed3dSBarry Smith   PetscFunctionReturn(0);
15154b0e389bSBarry Smith }
15164b0e389bSBarry Smith 
1517284134d9SBarry Smith #undef __FUNCT__
1518284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1519284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1520284134d9SBarry Smith {
1521284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1522284134d9SBarry Smith   PetscFunctionBegin;
152321a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1524284134d9SBarry Smith   m = PetscMax(m,M);
1525284134d9SBarry Smith   n = PetscMax(n,N);
152621a2c019SBarry 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);
1527284134d9SBarry 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);
1528*dc5cefdeSJed Brown   A->rmap->n = A->rmap->N = m;
1529d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
153021a2c019SBarry Smith   if (a->changelda) a->lda = m;
1531284134d9SBarry Smith   PetscFunctionReturn(0);
1532284134d9SBarry Smith }
1533170fe5c8SBarry Smith 
1534284134d9SBarry Smith 
1535a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1536a9fe9ddaSSatish Balay #undef __FUNCT__
1537a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1538a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1539a9fe9ddaSSatish Balay {
1540a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1541a9fe9ddaSSatish Balay 
1542a9fe9ddaSSatish Balay   PetscFunctionBegin;
1543a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1544a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1545a9fe9ddaSSatish Balay   }
1546a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1547a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1548a9fe9ddaSSatish Balay }
1549a9fe9ddaSSatish Balay 
1550a9fe9ddaSSatish Balay #undef __FUNCT__
1551a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1552a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1553a9fe9ddaSSatish Balay {
1554ee16a9a1SHong Zhang   PetscErrorCode ierr;
1555d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1556ee16a9a1SHong Zhang   Mat            Cmat;
1557a9fe9ddaSSatish Balay 
1558ee16a9a1SHong Zhang   PetscFunctionBegin;
1559d0f46423SBarry 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);
1560ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1561ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1562ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1563ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1564ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1565ee16a9a1SHong Zhang   *C = Cmat;
1566ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1567ee16a9a1SHong Zhang }
1568a9fe9ddaSSatish Balay 
156998a3b096SSatish Balay #undef __FUNCT__
1570a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1571a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1572a9fe9ddaSSatish Balay {
1573a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1574a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1575a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
15760805154bSBarry Smith   PetscBLASInt   m,n,k;
1577a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1578a9fe9ddaSSatish Balay 
1579a9fe9ddaSSatish Balay   PetscFunctionBegin;
1580d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1581d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1582d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1583a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1584a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1585a9fe9ddaSSatish Balay }
1586a9fe9ddaSSatish Balay 
1587a9fe9ddaSSatish Balay #undef __FUNCT__
1588a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1589a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1590a9fe9ddaSSatish Balay {
1591a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1592a9fe9ddaSSatish Balay 
1593a9fe9ddaSSatish Balay   PetscFunctionBegin;
1594a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1595a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1596a9fe9ddaSSatish Balay   }
1597a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1598a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1599a9fe9ddaSSatish Balay }
1600a9fe9ddaSSatish Balay 
1601a9fe9ddaSSatish Balay #undef __FUNCT__
1602a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1603a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1604a9fe9ddaSSatish Balay {
1605ee16a9a1SHong Zhang   PetscErrorCode ierr;
1606d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1607ee16a9a1SHong Zhang   Mat            Cmat;
1608a9fe9ddaSSatish Balay 
1609ee16a9a1SHong Zhang   PetscFunctionBegin;
1610d0f46423SBarry 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);
1611ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1612ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1613ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1614ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1615ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1616ee16a9a1SHong Zhang   *C = Cmat;
1617ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1618ee16a9a1SHong Zhang }
1619a9fe9ddaSSatish Balay 
1620a9fe9ddaSSatish Balay #undef __FUNCT__
1621a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1622a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1623a9fe9ddaSSatish Balay {
1624a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1625a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1626a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16270805154bSBarry Smith   PetscBLASInt   m,n,k;
1628a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1629a9fe9ddaSSatish Balay 
1630a9fe9ddaSSatish Balay   PetscFunctionBegin;
1631d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1632d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1633d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
16342fbe02b9SBarry Smith   /*
16352fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
16362fbe02b9SBarry Smith   */
1637a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1638a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1639a9fe9ddaSSatish Balay }
1640985db425SBarry Smith 
1641985db425SBarry Smith #undef __FUNCT__
1642985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1643985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1644985db425SBarry Smith {
1645985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1646985db425SBarry Smith   PetscErrorCode ierr;
1647d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1648985db425SBarry Smith   PetscScalar    *x;
1649985db425SBarry Smith   MatScalar      *aa = a->v;
1650985db425SBarry Smith 
1651985db425SBarry Smith   PetscFunctionBegin;
1652985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1653985db425SBarry Smith 
1654985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1655985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1656985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1657d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1658985db425SBarry Smith   for (i=0; i<m; i++) {
1659985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1660985db425SBarry Smith     for (j=1; j<n; j++){
1661985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1662985db425SBarry Smith     }
1663985db425SBarry Smith   }
1664985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1665985db425SBarry Smith   PetscFunctionReturn(0);
1666985db425SBarry Smith }
1667985db425SBarry Smith 
1668985db425SBarry Smith #undef __FUNCT__
1669985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1670985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1671985db425SBarry Smith {
1672985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1673985db425SBarry Smith   PetscErrorCode ierr;
1674d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1675985db425SBarry Smith   PetscScalar    *x;
1676985db425SBarry Smith   PetscReal      atmp;
1677985db425SBarry Smith   MatScalar      *aa = a->v;
1678985db425SBarry Smith 
1679985db425SBarry Smith   PetscFunctionBegin;
1680985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1681985db425SBarry Smith 
1682985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1683985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1684985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1685d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1686985db425SBarry Smith   for (i=0; i<m; i++) {
16879189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1688985db425SBarry Smith     for (j=1; j<n; j++){
1689985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1690985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1691985db425SBarry Smith     }
1692985db425SBarry Smith   }
1693985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1694985db425SBarry Smith   PetscFunctionReturn(0);
1695985db425SBarry Smith }
1696985db425SBarry Smith 
1697985db425SBarry Smith #undef __FUNCT__
1698985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1699985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1700985db425SBarry Smith {
1701985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1702985db425SBarry Smith   PetscErrorCode ierr;
1703d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1704985db425SBarry Smith   PetscScalar    *x;
1705985db425SBarry Smith   MatScalar      *aa = a->v;
1706985db425SBarry Smith 
1707985db425SBarry Smith   PetscFunctionBegin;
1708985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1709985db425SBarry Smith 
1710985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1711985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1712985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1713d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1714985db425SBarry Smith   for (i=0; i<m; i++) {
1715985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1716985db425SBarry Smith     for (j=1; j<n; j++){
1717985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1718985db425SBarry Smith     }
1719985db425SBarry Smith   }
1720985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1721985db425SBarry Smith   PetscFunctionReturn(0);
1722985db425SBarry Smith }
1723985db425SBarry Smith 
17248d0534beSBarry Smith #undef __FUNCT__
17258d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
17268d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
17278d0534beSBarry Smith {
17288d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17298d0534beSBarry Smith   PetscErrorCode ierr;
17308d0534beSBarry Smith   PetscScalar    *x;
17318d0534beSBarry Smith 
17328d0534beSBarry Smith   PetscFunctionBegin;
17338d0534beSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
17348d0534beSBarry Smith 
17358d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1736d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
17378d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
17388d0534beSBarry Smith   PetscFunctionReturn(0);
17398d0534beSBarry Smith }
17408d0534beSBarry Smith 
1741289bc588SBarry Smith /* -------------------------------------------------------------------*/
1742a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1743905e6a2fSBarry Smith        MatGetRow_SeqDense,
1744905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1745905e6a2fSBarry Smith        MatMult_SeqDense,
174697304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
17477c922b88SBarry Smith        MatMultTranspose_SeqDense,
17487c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1749db4efbfdSBarry Smith        0,
1750db4efbfdSBarry Smith        0,
1751db4efbfdSBarry Smith        0,
1752db4efbfdSBarry Smith /*10*/ 0,
1753905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1754905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
175541f059aeSBarry Smith        MatSOR_SeqDense,
1756ec8511deSBarry Smith        MatTranspose_SeqDense,
175797304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1758905e6a2fSBarry Smith        MatEqual_SeqDense,
1759905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1760905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1761905e6a2fSBarry Smith        MatNorm_SeqDense,
1762c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1763c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1764905e6a2fSBarry Smith        MatSetOption_SeqDense,
1765905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1766d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense,
1767db4efbfdSBarry Smith        0,
1768db4efbfdSBarry Smith        0,
1769db4efbfdSBarry Smith        0,
1770db4efbfdSBarry Smith        0,
1771d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense,
1772273d9f13SBarry Smith        0,
1773905e6a2fSBarry Smith        0,
1774905e6a2fSBarry Smith        MatGetArray_SeqDense,
1775905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1776d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense,
1777a5ae1ecdSBarry Smith        0,
1778a5ae1ecdSBarry Smith        0,
1779a5ae1ecdSBarry Smith        0,
1780a5ae1ecdSBarry Smith        0,
1781d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense,
1782a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1783a5ae1ecdSBarry Smith        0,
17844b0e389bSBarry Smith        MatGetValues_SeqDense,
1785a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1786d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense,
1787a5ae1ecdSBarry Smith        MatScale_SeqDense,
1788a5ae1ecdSBarry Smith        0,
1789a5ae1ecdSBarry Smith        0,
1790a5ae1ecdSBarry Smith        0,
1791d519adbfSMatthew Knepley /*49*/ 0,
1792a5ae1ecdSBarry Smith        0,
1793a5ae1ecdSBarry Smith        0,
1794a5ae1ecdSBarry Smith        0,
1795a5ae1ecdSBarry Smith        0,
1796d519adbfSMatthew Knepley /*54*/ 0,
1797a5ae1ecdSBarry Smith        0,
1798a5ae1ecdSBarry Smith        0,
1799a5ae1ecdSBarry Smith        0,
1800a5ae1ecdSBarry Smith        0,
1801d519adbfSMatthew Knepley /*59*/ 0,
1802e03a110bSBarry Smith        MatDestroy_SeqDense,
1803e03a110bSBarry Smith        MatView_SeqDense,
1804357abbc8SBarry Smith        0,
180597304618SKris Buschelman        0,
1806d519adbfSMatthew Knepley /*64*/ 0,
180797304618SKris Buschelman        0,
180897304618SKris Buschelman        0,
180997304618SKris Buschelman        0,
181097304618SKris Buschelman        0,
1811d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense,
181297304618SKris Buschelman        0,
181397304618SKris Buschelman        0,
181497304618SKris Buschelman        0,
181597304618SKris Buschelman        0,
1816d519adbfSMatthew Knepley /*74*/ 0,
181797304618SKris Buschelman        0,
181897304618SKris Buschelman        0,
181997304618SKris Buschelman        0,
182097304618SKris Buschelman        0,
1821d519adbfSMatthew Knepley /*79*/ 0,
182297304618SKris Buschelman        0,
182397304618SKris Buschelman        0,
182497304618SKris Buschelman        0,
1825d519adbfSMatthew Knepley /*83*/ MatLoad_SeqDense,
1826865e5f61SKris Buschelman        0,
18271cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1828865e5f61SKris Buschelman        0,
1829865e5f61SKris Buschelman        0,
1830865e5f61SKris Buschelman        0,
1831d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense,
1832a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1833a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1834865e5f61SKris Buschelman        0,
1835865e5f61SKris Buschelman        0,
1836d519adbfSMatthew Knepley /*94*/ 0,
1837a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1838a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1839a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1840284134d9SBarry Smith        0,
1841d519adbfSMatthew Knepley /*99*/ 0,
1842284134d9SBarry Smith        0,
1843284134d9SBarry Smith        0,
1844284134d9SBarry Smith        0,
1845985db425SBarry Smith        MatSetSizes_SeqDense,
1846985db425SBarry Smith        0,
1847985db425SBarry Smith        0,
1848985db425SBarry Smith        0,
1849985db425SBarry Smith        0,
1850985db425SBarry Smith        0,
1851d519adbfSMatthew Knepley /*109*/0,
1852985db425SBarry Smith        0,
18538d0534beSBarry Smith        MatGetRowMin_SeqDense,
18548d0534beSBarry Smith        MatGetColumnVector_SeqDense
1855985db425SBarry Smith };
185690ace30eSBarry Smith 
18574a2ae208SSatish Balay #undef __FUNCT__
18584a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18594b828684SBarry Smith /*@C
1860fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1861d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1862d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1863289bc588SBarry Smith 
1864db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1865db81eaa0SLois Curfman McInnes 
186620563c6bSBarry Smith    Input Parameters:
1867db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
18680c775827SLois Curfman McInnes .  m - number of rows
186918f449edSLois Curfman McInnes .  n - number of columns
1870c0235b3cSMatthew Knepley -  data - optional location of matrix data in column major order.  Set data=PETSC_NULL for PETSc
1871dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
187220563c6bSBarry Smith 
187320563c6bSBarry Smith    Output Parameter:
187444cd7ae7SLois Curfman McInnes .  A - the matrix
187520563c6bSBarry Smith 
1876b259b22eSLois Curfman McInnes    Notes:
187718f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
187818f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1879b4fd4287SBarry Smith    set data=PETSC_NULL.
188018f449edSLois Curfman McInnes 
1881027ccd11SLois Curfman McInnes    Level: intermediate
1882027ccd11SLois Curfman McInnes 
1883dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1884d65003e9SLois Curfman McInnes 
1885db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
188620563c6bSBarry Smith @*/
1887be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1888289bc588SBarry Smith {
1889dfbe8321SBarry Smith   PetscErrorCode ierr;
18903b2fbd54SBarry Smith 
18913a40ed3dSBarry Smith   PetscFunctionBegin;
1892f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1893f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1894273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1895273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1896273d9f13SBarry Smith   PetscFunctionReturn(0);
1897273d9f13SBarry Smith }
1898273d9f13SBarry Smith 
18994a2ae208SSatish Balay #undef __FUNCT__
1900afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1901273d9f13SBarry Smith /*@C
1902273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1903273d9f13SBarry Smith 
1904273d9f13SBarry Smith    Collective on MPI_Comm
1905273d9f13SBarry Smith 
1906273d9f13SBarry Smith    Input Parameters:
1907273d9f13SBarry Smith +  A - the matrix
1908273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1909273d9f13SBarry Smith 
1910273d9f13SBarry Smith    Notes:
1911273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1912273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1913284134d9SBarry Smith    need not call this routine.
1914273d9f13SBarry Smith 
1915273d9f13SBarry Smith    Level: intermediate
1916273d9f13SBarry Smith 
1917273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1918273d9f13SBarry Smith 
1919867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
1920867c911aSBarry Smith 
1921273d9f13SBarry Smith @*/
1922be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1923273d9f13SBarry Smith {
1924dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1925a23d5eceSKris Buschelman 
1926a23d5eceSKris Buschelman   PetscFunctionBegin;
1927a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1928a23d5eceSKris Buschelman   if (f) {
1929a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1930a23d5eceSKris Buschelman   }
1931a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1932a23d5eceSKris Buschelman }
1933a23d5eceSKris Buschelman 
1934a23d5eceSKris Buschelman EXTERN_C_BEGIN
1935a23d5eceSKris Buschelman #undef __FUNCT__
1936afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1937be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1938a23d5eceSKris Buschelman {
1939273d9f13SBarry Smith   Mat_SeqDense   *b;
1940dfbe8321SBarry Smith   PetscErrorCode ierr;
1941273d9f13SBarry Smith 
1942273d9f13SBarry Smith   PetscFunctionBegin;
1943273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1944273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1945d0f46423SBarry Smith   if (b->lda <= 0) b->lda = B->rmap->n;
19469e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
19479e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
19485afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1949284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1950284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
19519e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
1952273d9f13SBarry Smith   } else { /* user-allocated storage */
19539e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1954273d9f13SBarry Smith     b->v          = data;
1955273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1956273d9f13SBarry Smith   }
19570450473dSBarry Smith   B->assembled = PETSC_TRUE;
1958273d9f13SBarry Smith   PetscFunctionReturn(0);
1959273d9f13SBarry Smith }
1960a23d5eceSKris Buschelman EXTERN_C_END
1961273d9f13SBarry Smith 
19621b807ce4Svictorle #undef __FUNCT__
19631b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
19641b807ce4Svictorle /*@C
19651b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
19661b807ce4Svictorle 
19671b807ce4Svictorle   Input parameter:
19681b807ce4Svictorle + A - the matrix
19691b807ce4Svictorle - lda - the leading dimension
19701b807ce4Svictorle 
19711b807ce4Svictorle   Notes:
1972867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
19731b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
19741b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
19751b807ce4Svictorle 
19761b807ce4Svictorle   Level: intermediate
19771b807ce4Svictorle 
19781b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
19791b807ce4Svictorle 
1980284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
1981867c911aSBarry Smith 
19821b807ce4Svictorle @*/
1983be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
19841b807ce4Svictorle {
19851b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
198621a2c019SBarry Smith 
19871b807ce4Svictorle   PetscFunctionBegin;
1988d0f46423SBarry Smith   if (lda < B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
19891b807ce4Svictorle   b->lda       = lda;
199021a2c019SBarry Smith   b->changelda = PETSC_FALSE;
199121a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
19921b807ce4Svictorle   PetscFunctionReturn(0);
19931b807ce4Svictorle }
19941b807ce4Svictorle 
19950bad9183SKris Buschelman /*MC
1996fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
19970bad9183SKris Buschelman 
19980bad9183SKris Buschelman    Options Database Keys:
19990bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
20000bad9183SKris Buschelman 
20010bad9183SKris Buschelman   Level: beginner
20020bad9183SKris Buschelman 
200389665df3SBarry Smith .seealso: MatCreateSeqDense()
200489665df3SBarry Smith 
20050bad9183SKris Buschelman M*/
20060bad9183SKris Buschelman 
2007273d9f13SBarry Smith EXTERN_C_BEGIN
20084a2ae208SSatish Balay #undef __FUNCT__
20094a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
2010be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2011273d9f13SBarry Smith {
2012273d9f13SBarry Smith   Mat_SeqDense   *b;
2013dfbe8321SBarry Smith   PetscErrorCode ierr;
20147c334f02SBarry Smith   PetscMPIInt    size;
2015273d9f13SBarry Smith 
2016273d9f13SBarry Smith   PetscFunctionBegin;
20177adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
201829bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
201955659b69SBarry Smith 
202026283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
202126283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
202226283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
202326283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
2024273d9f13SBarry Smith 
202538f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2026549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
202790f02eecSBarry Smith   B->mapping      = 0;
202844cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
202918f449edSLois Curfman McInnes 
2030a5ae1ecdSBarry Smith 
203144cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2032273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2033273d9f13SBarry Smith   b->v            = 0;
2034d0f46423SBarry Smith   b->lda          = B->rmap->n;
203521a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
2036d0f46423SBarry Smith   b->Mmax         = B->rmap->n;
2037d0f46423SBarry Smith   b->Nmax         = B->cmap->n;
20384e220ebcSLois Curfman McInnes 
2039b24902e0SBarry Smith 
2040ec1065edSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
2041b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2042b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2043a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2044a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2045a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
20464ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
20474ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
20484ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
20494ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
20504ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
20514ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
20524ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
20534ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
20544ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
205517667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
20563a40ed3dSBarry Smith   PetscFunctionReturn(0);
2057289bc588SBarry Smith }
2058c0aa2d19SHong Zhang 
2059c0aa2d19SHong Zhang 
2060273d9f13SBarry Smith EXTERN_C_END
2061