xref: /petsc/src/mat/impls/dense/seq/dense.c (revision dc0b31ed242eb6637dca0468fbfa12899dc55c4b)
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);
187*dc0b31edSSatish 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);
222*dc0b31edSSatish 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);
266*dc0b31edSSatish 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);
314*dc0b31edSSatish 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 
350*dc0b31edSSatish 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;
379*dc0b31edSSatish 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 
419db4efbfdSBarry Smith #undef __FUNCT__
420db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc"
421db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
422db4efbfdSBarry Smith {
423db4efbfdSBarry Smith   PetscErrorCode ierr;
424db4efbfdSBarry Smith 
425db4efbfdSBarry Smith   PetscFunctionBegin;
426db4efbfdSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
427db4efbfdSBarry Smith   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
428db4efbfdSBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
429db4efbfdSBarry Smith   if (ftype == MAT_FACTOR_LU){
430db4efbfdSBarry Smith     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
431db4efbfdSBarry Smith   } else {
432db4efbfdSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
433db4efbfdSBarry Smith   }
434db4efbfdSBarry Smith   (*fact)->factor = ftype;
435db4efbfdSBarry Smith   PetscFunctionReturn(0);
436db4efbfdSBarry Smith }
437db4efbfdSBarry Smith 
438289bc588SBarry Smith /* ------------------------------------------------------------------*/
4394a2ae208SSatish Balay #undef __FUNCT__
4404a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
44113f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
442289bc588SBarry Smith {
443c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
44487828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
445dfbe8321SBarry Smith   PetscErrorCode ierr;
446d0f46423SBarry Smith   PetscInt       m = A->rmap->n,i;
447aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4480805154bSBarry Smith   PetscBLASInt   o = 1,bm = PetscBLASIntCast(m);
449bc1b551cSSatish Balay #endif
450289bc588SBarry Smith 
4513a40ed3dSBarry Smith   PetscFunctionBegin;
452289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
45371044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4542dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
455289bc588SBarry Smith   }
4561ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4571ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
458b965ef7fSBarry Smith   its  = its*lits;
45977431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
460289bc588SBarry Smith   while (its--) {
461fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
462289bc588SBarry Smith       for (i=0; i<m; i++) {
463aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
464f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
465f1747703SBarry Smith            not happy about returning a double complex */
46613f74950SBarry Smith         PetscInt    _i;
46787828ca2SBarry Smith         PetscScalar sum = b[i];
468f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4693f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
470f1747703SBarry Smith         }
471f1747703SBarry Smith         xt = sum;
472f1747703SBarry Smith #else
47371044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
474f1747703SBarry Smith #endif
47555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
476289bc588SBarry Smith       }
477289bc588SBarry Smith     }
478fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
479289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
480aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
481f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
482f1747703SBarry Smith            not happy about returning a double complex */
48313f74950SBarry Smith         PetscInt    _i;
48487828ca2SBarry Smith         PetscScalar sum = b[i];
485f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4863f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
487f1747703SBarry Smith         }
488f1747703SBarry Smith         xt = sum;
489f1747703SBarry Smith #else
49071044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
491f1747703SBarry Smith #endif
49255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
493289bc588SBarry Smith       }
494289bc588SBarry Smith     }
495289bc588SBarry Smith   }
4961ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4971ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4983a40ed3dSBarry Smith   PetscFunctionReturn(0);
499289bc588SBarry Smith }
500289bc588SBarry Smith 
501289bc588SBarry Smith /* -----------------------------------------------------------------*/
5024a2ae208SSatish Balay #undef __FUNCT__
5034a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
504dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
505289bc588SBarry Smith {
506c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
50787828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
508dfbe8321SBarry Smith   PetscErrorCode ierr;
5090805154bSBarry Smith   PetscBLASInt   m, n,_One=1;
510ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
5113a40ed3dSBarry Smith 
5123a40ed3dSBarry Smith   PetscFunctionBegin;
513d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
514d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
515d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5161ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5171ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
51871044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
5191ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5201ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
521*dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
5223a40ed3dSBarry Smith   PetscFunctionReturn(0);
523289bc588SBarry Smith }
5246ee01492SSatish Balay 
5254a2ae208SSatish Balay #undef __FUNCT__
5264a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
527dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
528289bc588SBarry Smith {
529c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
53087828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
531dfbe8321SBarry Smith   PetscErrorCode ierr;
5320805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5333a40ed3dSBarry Smith 
5343a40ed3dSBarry Smith   PetscFunctionBegin;
535d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
536d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
537d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5381ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5391ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
54071044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5411ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5421ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
543*dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
5443a40ed3dSBarry Smith   PetscFunctionReturn(0);
545289bc588SBarry Smith }
5466ee01492SSatish Balay 
5474a2ae208SSatish Balay #undef __FUNCT__
5484a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
549dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
550289bc588SBarry Smith {
551c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
55287828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
553dfbe8321SBarry Smith   PetscErrorCode ierr;
5540805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5553a40ed3dSBarry Smith 
5563a40ed3dSBarry Smith   PetscFunctionBegin;
557d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
558d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
559d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
560600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5611ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5621ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
56371044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5651ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
566*dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5673a40ed3dSBarry Smith   PetscFunctionReturn(0);
568289bc588SBarry Smith }
5696ee01492SSatish Balay 
5704a2ae208SSatish Balay #undef __FUNCT__
5714a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
572dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
573289bc588SBarry Smith {
574c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
57587828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
576dfbe8321SBarry Smith   PetscErrorCode ierr;
5770805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
57887828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5793a40ed3dSBarry Smith 
5803a40ed3dSBarry Smith   PetscFunctionBegin;
581d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
582d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
583d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
584600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5851ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5861ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
58771044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5881ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5891ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
590*dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5913a40ed3dSBarry Smith   PetscFunctionReturn(0);
592289bc588SBarry Smith }
593289bc588SBarry Smith 
594289bc588SBarry Smith /* -----------------------------------------------------------------*/
5954a2ae208SSatish Balay #undef __FUNCT__
5964a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
59713f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
598289bc588SBarry Smith {
599c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
60087828ca2SBarry Smith   PetscScalar    *v;
6016849ba73SBarry Smith   PetscErrorCode ierr;
60213f74950SBarry Smith   PetscInt       i;
60367e560aaSBarry Smith 
6043a40ed3dSBarry Smith   PetscFunctionBegin;
605d0f46423SBarry Smith   *ncols = A->cmap->n;
606289bc588SBarry Smith   if (cols) {
607d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
608d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
609289bc588SBarry Smith   }
610289bc588SBarry Smith   if (vals) {
611d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
612289bc588SBarry Smith     v    = mat->v + row;
613d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
614289bc588SBarry Smith   }
6153a40ed3dSBarry Smith   PetscFunctionReturn(0);
616289bc588SBarry Smith }
6176ee01492SSatish Balay 
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
62013f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
621289bc588SBarry Smith {
622dfbe8321SBarry Smith   PetscErrorCode ierr;
623606d414cSSatish Balay   PetscFunctionBegin;
624606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
625606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6263a40ed3dSBarry Smith   PetscFunctionReturn(0);
627289bc588SBarry Smith }
628289bc588SBarry Smith /* ----------------------------------------------------------------*/
6294a2ae208SSatish Balay #undef __FUNCT__
6304a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
63113f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
632289bc588SBarry Smith {
633c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
63413f74950SBarry Smith   PetscInt     i,j,idx=0;
635d6dfbf8fSBarry Smith 
6363a40ed3dSBarry Smith   PetscFunctionBegin;
637289bc588SBarry Smith   if (!mat->roworiented) {
638dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
639289bc588SBarry Smith       for (j=0; j<n; j++) {
640cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6412515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
642d0f46423SBarry 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);
64358804f6dSBarry Smith #endif
644289bc588SBarry Smith         for (i=0; i<m; i++) {
645cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6462515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
647d0f46423SBarry 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);
64858804f6dSBarry Smith #endif
649cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
650289bc588SBarry Smith         }
651289bc588SBarry Smith       }
6523a40ed3dSBarry Smith     } else {
653289bc588SBarry Smith       for (j=0; j<n; j++) {
654cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6552515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
656d0f46423SBarry 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);
65758804f6dSBarry Smith #endif
658289bc588SBarry Smith         for (i=0; i<m; i++) {
659cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6602515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
661d0f46423SBarry 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);
66258804f6dSBarry Smith #endif
663cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
664289bc588SBarry Smith         }
665289bc588SBarry Smith       }
666289bc588SBarry Smith     }
6673a40ed3dSBarry Smith   } else {
668dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
669e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
670cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6712515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
672d0f46423SBarry 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);
67358804f6dSBarry Smith #endif
674e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
675cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6762515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
677d0f46423SBarry 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);
67858804f6dSBarry Smith #endif
679cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
680e8d4e0b9SBarry Smith         }
681e8d4e0b9SBarry Smith       }
6823a40ed3dSBarry Smith     } else {
683289bc588SBarry Smith       for (i=0; i<m; i++) {
684cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6852515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
686d0f46423SBarry 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);
68758804f6dSBarry Smith #endif
688289bc588SBarry Smith         for (j=0; j<n; j++) {
689cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6902515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
691d0f46423SBarry 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);
69258804f6dSBarry Smith #endif
693cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
694289bc588SBarry Smith         }
695289bc588SBarry Smith       }
696289bc588SBarry Smith     }
697e8d4e0b9SBarry Smith   }
6983a40ed3dSBarry Smith   PetscFunctionReturn(0);
699289bc588SBarry Smith }
700e8d4e0b9SBarry Smith 
7014a2ae208SSatish Balay #undef __FUNCT__
7024a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
70313f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
704ae80bb75SLois Curfman McInnes {
705ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
70613f74950SBarry Smith   PetscInt     i,j;
707ae80bb75SLois Curfman McInnes 
7083a40ed3dSBarry Smith   PetscFunctionBegin;
709ae80bb75SLois Curfman McInnes   /* row-oriented output */
710ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
71197e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
712d0f46423SBarry 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);
713ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7146f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
715d0f46423SBarry 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);
71697e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
717ae80bb75SLois Curfman McInnes     }
718ae80bb75SLois Curfman McInnes   }
7193a40ed3dSBarry Smith   PetscFunctionReturn(0);
720ae80bb75SLois Curfman McInnes }
721ae80bb75SLois Curfman McInnes 
722289bc588SBarry Smith /* -----------------------------------------------------------------*/
723289bc588SBarry Smith 
724e090d566SSatish Balay #include "petscsys.h"
72588e32edaSLois Curfman McInnes 
7264a2ae208SSatish Balay #undef __FUNCT__
7274a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
728a313700dSBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, const MatType type,Mat *A)
72988e32edaSLois Curfman McInnes {
73088e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
73188e32edaSLois Curfman McInnes   Mat            B;
7326849ba73SBarry Smith   PetscErrorCode ierr;
73313f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
73413f74950SBarry Smith   int            fd;
73513f74950SBarry Smith   PetscMPIInt    size;
73613f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
73787828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
73819bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
73988e32edaSLois Curfman McInnes 
7403a40ed3dSBarry Smith   PetscFunctionBegin;
741d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
74229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
743b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
7440752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
745552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
74688e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
74788e32edaSLois Curfman McInnes 
74890ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
749f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
750f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
751be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7525c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
75390ace30eSBarry Smith     B    = *A;
75490ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
7554a41a4d3SSatish Balay     v    = a->v;
7564a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7574a41a4d3SSatish Balay        from row major to column major */
758b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
75990ace30eSBarry Smith     /* read in nonzero values */
7604a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7614a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7624a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7634a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7644a41a4d3SSatish Balay         *v++ =w[i*N+j];
7654a41a4d3SSatish Balay       }
7664a41a4d3SSatish Balay     }
767606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7686d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7696d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77090ace30eSBarry Smith   } else {
77188e32edaSLois Curfman McInnes     /* read row lengths */
77213f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7730752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
77488e32edaSLois Curfman McInnes 
77588e32edaSLois Curfman McInnes     /* create our matrix */
776f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
777f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
778be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7795c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
78088e32edaSLois Curfman McInnes     B = *A;
78188e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
78288e32edaSLois Curfman McInnes     v = a->v;
78388e32edaSLois Curfman McInnes 
78488e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
78513f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
786b0a32e0cSBarry Smith     cols = scols;
7870752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
78887828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
789b0a32e0cSBarry Smith     vals = svals;
7900752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
79188e32edaSLois Curfman McInnes 
79288e32edaSLois Curfman McInnes     /* insert into matrix */
79388e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
79488e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
79588e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
79688e32edaSLois Curfman McInnes     }
797606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
798606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
799606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
80088e32edaSLois Curfman McInnes 
8016d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8026d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80390ace30eSBarry Smith   }
8043a40ed3dSBarry Smith   PetscFunctionReturn(0);
80588e32edaSLois Curfman McInnes }
80688e32edaSLois Curfman McInnes 
807e090d566SSatish Balay #include "petscsys.h"
808932b0c3eSLois Curfman McInnes 
8094a2ae208SSatish Balay #undef __FUNCT__
8104a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8116849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
812289bc588SBarry Smith {
813932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
814dfbe8321SBarry Smith   PetscErrorCode    ierr;
81513f74950SBarry Smith   PetscInt          i,j;
8162dcb1b2aSMatthew Knepley   const char        *name;
81787828ca2SBarry Smith   PetscScalar       *v;
818f3ef73ceSBarry Smith   PetscViewerFormat format;
8195f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
8205f481a85SSatish Balay   PetscTruth        allreal = PETSC_TRUE;
8215f481a85SSatish Balay #endif
822932b0c3eSLois Curfman McInnes 
8233a40ed3dSBarry Smith   PetscFunctionBegin;
824435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
825b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
826456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8273a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
828fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
829b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
830d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
83144cd7ae7SLois Curfman McInnes       v = a->v + i;
83277431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
833d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
834aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
835329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
836a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
837329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
838a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8396831982aSBarry Smith         }
84080cd9d93SLois Curfman McInnes #else
8416831982aSBarry Smith         if (*v) {
842a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8436831982aSBarry Smith         }
84480cd9d93SLois Curfman McInnes #endif
8451b807ce4Svictorle         v += a->lda;
84680cd9d93SLois Curfman McInnes       }
847b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
84880cd9d93SLois Curfman McInnes     }
849b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
8503a40ed3dSBarry Smith   } else {
851b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
852aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
85347989497SBarry Smith     /* determine if matrix has all real values */
85447989497SBarry Smith     v = a->v;
855d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
856ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
85747989497SBarry Smith     }
85847989497SBarry Smith #endif
859fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8603a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
861d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
862d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
863fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
864ffac6cdbSBarry Smith     }
865ffac6cdbSBarry Smith 
866d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
867932b0c3eSLois Curfman McInnes       v = a->v + i;
868d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
869aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
87047989497SBarry Smith         if (allreal) {
871f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
87247989497SBarry Smith         } else {
873f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
87447989497SBarry Smith         }
875289bc588SBarry Smith #else
876f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
877289bc588SBarry Smith #endif
8781b807ce4Svictorle 	v += a->lda;
879289bc588SBarry Smith       }
880b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
881289bc588SBarry Smith     }
882fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
883b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
884ffac6cdbSBarry Smith     }
885b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
886da3a660dSBarry Smith   }
887b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8883a40ed3dSBarry Smith   PetscFunctionReturn(0);
889289bc588SBarry Smith }
890289bc588SBarry Smith 
8914a2ae208SSatish Balay #undef __FUNCT__
8924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8936849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
894932b0c3eSLois Curfman McInnes {
895932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8966849ba73SBarry Smith   PetscErrorCode    ierr;
89713f74950SBarry Smith   int               fd;
898d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
89987828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
900f3ef73ceSBarry Smith   PetscViewerFormat format;
901932b0c3eSLois Curfman McInnes 
9023a40ed3dSBarry Smith   PetscFunctionBegin;
903b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
90490ace30eSBarry Smith 
905b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
9067973ad04SBarry Smith   if (format == PETSC_VIEWER_NATIVE) {
90790ace30eSBarry Smith     /* store the matrix as a dense matrix */
90813f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
909552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
91090ace30eSBarry Smith     col_lens[1] = m;
91190ace30eSBarry Smith     col_lens[2] = n;
91290ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
9136f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
914606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
91590ace30eSBarry Smith 
91690ace30eSBarry Smith     /* write out matrix, by rows */
91787828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
91890ace30eSBarry Smith     v    = a->v;
91990ace30eSBarry Smith     for (j=0; j<n; j++) {
920578230a0SSatish Balay       for (i=0; i<m; i++) {
921578230a0SSatish Balay         vals[j + i*n] = *v++;
92290ace30eSBarry Smith       }
92390ace30eSBarry Smith     }
9246f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
925606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
92690ace30eSBarry Smith   } else {
92713f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
928552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
929932b0c3eSLois Curfman McInnes     col_lens[1] = m;
930932b0c3eSLois Curfman McInnes     col_lens[2] = n;
931932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
932932b0c3eSLois Curfman McInnes 
933932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
934932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9356f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
936932b0c3eSLois Curfman McInnes 
937932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
938932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
939932b0c3eSLois Curfman McInnes     ict = 0;
940932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
941932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
942932b0c3eSLois Curfman McInnes     }
9436f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
944606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
945932b0c3eSLois Curfman McInnes 
946932b0c3eSLois Curfman McInnes     /* store nonzero values */
94787828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
948932b0c3eSLois Curfman McInnes     ict  = 0;
949932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
950932b0c3eSLois Curfman McInnes       v = a->v + i;
951932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9521b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
953932b0c3eSLois Curfman McInnes       }
954932b0c3eSLois Curfman McInnes     }
9556f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
956606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
95790ace30eSBarry Smith   }
9583a40ed3dSBarry Smith   PetscFunctionReturn(0);
959932b0c3eSLois Curfman McInnes }
960932b0c3eSLois Curfman McInnes 
9614a2ae208SSatish Balay #undef __FUNCT__
9624a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
963dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
964f1af5d2fSBarry Smith {
965f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
966f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9676849ba73SBarry Smith   PetscErrorCode    ierr;
968d0f46423SBarry Smith   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
96987828ca2SBarry Smith   PetscScalar       *v = a->v;
970b0a32e0cSBarry Smith   PetscViewer       viewer;
971b0a32e0cSBarry Smith   PetscDraw         popup;
972329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
973f3ef73ceSBarry Smith   PetscViewerFormat format;
974f1af5d2fSBarry Smith 
975f1af5d2fSBarry Smith   PetscFunctionBegin;
976f1af5d2fSBarry Smith 
977f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
978b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
979b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
980f1af5d2fSBarry Smith 
981f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
982fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
983f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
984b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
985f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
986f1af5d2fSBarry Smith       x_l = j;
987f1af5d2fSBarry Smith       x_r = x_l + 1.0;
988f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
989f1af5d2fSBarry Smith         y_l = m - i - 1.0;
990f1af5d2fSBarry Smith         y_r = y_l + 1.0;
991f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
992329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
993b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
994329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
995b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
996f1af5d2fSBarry Smith         } else {
997f1af5d2fSBarry Smith           continue;
998f1af5d2fSBarry Smith         }
999f1af5d2fSBarry Smith #else
1000f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
1001b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1002f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1003b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1004f1af5d2fSBarry Smith         } else {
1005f1af5d2fSBarry Smith           continue;
1006f1af5d2fSBarry Smith         }
1007f1af5d2fSBarry Smith #endif
1008b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1009f1af5d2fSBarry Smith       }
1010f1af5d2fSBarry Smith     }
1011f1af5d2fSBarry Smith   } else {
1012f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1013f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1014f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1015f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1016f1af5d2fSBarry Smith     }
1017b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1018b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1019b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1020f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1021f1af5d2fSBarry Smith       x_l = j;
1022f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1023f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1024f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1025f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1026b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1027b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1028f1af5d2fSBarry Smith       }
1029f1af5d2fSBarry Smith     }
1030f1af5d2fSBarry Smith   }
1031f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1032f1af5d2fSBarry Smith }
1033f1af5d2fSBarry Smith 
10344a2ae208SSatish Balay #undef __FUNCT__
10354a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1036dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1037f1af5d2fSBarry Smith {
1038b0a32e0cSBarry Smith   PetscDraw      draw;
1039f1af5d2fSBarry Smith   PetscTruth     isnull;
1040329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1041dfbe8321SBarry Smith   PetscErrorCode ierr;
1042f1af5d2fSBarry Smith 
1043f1af5d2fSBarry Smith   PetscFunctionBegin;
1044b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1045b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1046abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1047f1af5d2fSBarry Smith 
1048f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1049d0f46423SBarry Smith   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1050f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1051b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1052b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1053f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1054f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1055f1af5d2fSBarry Smith }
1056f1af5d2fSBarry Smith 
10574a2ae208SSatish Balay #undef __FUNCT__
10584a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1059dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1060932b0c3eSLois Curfman McInnes {
1061dfbe8321SBarry Smith   PetscErrorCode ierr;
10626805f65bSBarry Smith   PetscTruth     iascii,isbinary,isdraw;
1063932b0c3eSLois Curfman McInnes 
10643a40ed3dSBarry Smith   PetscFunctionBegin;
106532077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1066fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1067fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10680f5bd95cSBarry Smith 
1069c45a1595SBarry Smith   if (iascii) {
1070c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10710f5bd95cSBarry Smith   } else if (isbinary) {
10723a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1073f1af5d2fSBarry Smith   } else if (isdraw) {
1074f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10755cd90555SBarry Smith   } else {
1076958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1077932b0c3eSLois Curfman McInnes   }
10783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1079932b0c3eSLois Curfman McInnes }
1080289bc588SBarry Smith 
10814a2ae208SSatish Balay #undef __FUNCT__
10824a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1083dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1084289bc588SBarry Smith {
1085ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1086dfbe8321SBarry Smith   PetscErrorCode ierr;
108790f02eecSBarry Smith 
10883a40ed3dSBarry Smith   PetscFunctionBegin;
1089aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1090d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1091a5a9c739SBarry Smith #endif
109205b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
10936857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1094606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1095dbd8c25aSHong Zhang 
1096dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1097901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10984ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10994ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11004ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11013a40ed3dSBarry Smith   PetscFunctionReturn(0);
1102289bc588SBarry Smith }
1103289bc588SBarry Smith 
11044a2ae208SSatish Balay #undef __FUNCT__
11054a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1106fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1107289bc588SBarry Smith {
1108c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11096849ba73SBarry Smith   PetscErrorCode ierr;
111013f74950SBarry Smith   PetscInt       k,j,m,n,M;
111187828ca2SBarry Smith   PetscScalar    *v,tmp;
111248b35521SBarry Smith 
11133a40ed3dSBarry Smith   PetscFunctionBegin;
1114d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1115e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1116a5ce6ee0Svictorle     if (m != n) {
1117634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1118d64ed03dSBarry Smith     } else {
1119d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1120289bc588SBarry Smith         for (k=0; k<j; k++) {
11211b807ce4Svictorle           tmp = v[j + k*M];
11221b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11231b807ce4Svictorle           v[k + j*M] = tmp;
1124289bc588SBarry Smith         }
1125289bc588SBarry Smith       }
1126d64ed03dSBarry Smith     }
11273a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1128d3e5ee88SLois Curfman McInnes     Mat          tmat;
1129ec8511deSBarry Smith     Mat_SeqDense *tmatd;
113087828ca2SBarry Smith     PetscScalar  *v2;
1131ea709b57SSatish Balay 
1132fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11337adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1134d0f46423SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
11357adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11365c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1137fc4dec0aSBarry Smith     } else {
1138fc4dec0aSBarry Smith       tmat = *matout;
1139fc4dec0aSBarry Smith     }
1140ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11410de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1142d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11431b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1144d3e5ee88SLois Curfman McInnes     }
11456d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11466d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1147d3e5ee88SLois Curfman McInnes     *matout = tmat;
114848b35521SBarry Smith   }
11493a40ed3dSBarry Smith   PetscFunctionReturn(0);
1150289bc588SBarry Smith }
1151289bc588SBarry Smith 
11524a2ae208SSatish Balay #undef __FUNCT__
11534a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1154dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1155289bc588SBarry Smith {
1156c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1157c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
115813f74950SBarry Smith   PetscInt     i,j;
115987828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11609ea5d5aeSSatish Balay 
11613a40ed3dSBarry Smith   PetscFunctionBegin;
1162d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1163d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1164d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
11651b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1166d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
11673a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11681b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11691b807ce4Svictorle     }
1170289bc588SBarry Smith   }
117177c4ece6SBarry Smith   *flg = PETSC_TRUE;
11723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1173289bc588SBarry Smith }
1174289bc588SBarry Smith 
11754a2ae208SSatish Balay #undef __FUNCT__
11764a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1177dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1178289bc588SBarry Smith {
1179c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1180dfbe8321SBarry Smith   PetscErrorCode ierr;
118113f74950SBarry Smith   PetscInt       i,n,len;
118287828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
118344cd7ae7SLois Curfman McInnes 
11843a40ed3dSBarry Smith   PetscFunctionBegin;
11852dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11867a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11871ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1188d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1189d0f46423SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
119044cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11911b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1192289bc588SBarry Smith   }
11931ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11943a40ed3dSBarry Smith   PetscFunctionReturn(0);
1195289bc588SBarry Smith }
1196289bc588SBarry Smith 
11974a2ae208SSatish Balay #undef __FUNCT__
11984a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1199dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1200289bc588SBarry Smith {
1201c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
120287828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1203dfbe8321SBarry Smith   PetscErrorCode ierr;
1204d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
120555659b69SBarry Smith 
12063a40ed3dSBarry Smith   PetscFunctionBegin;
120728988994SBarry Smith   if (ll) {
12087a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12091ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1210d0f46423SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1211da3a660dSBarry Smith     for (i=0; i<m; i++) {
1212da3a660dSBarry Smith       x = l[i];
1213da3a660dSBarry Smith       v = mat->v + i;
1214da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1215da3a660dSBarry Smith     }
12161ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1217efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1218da3a660dSBarry Smith   }
121928988994SBarry Smith   if (rr) {
12207a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12211ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1222d0f46423SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1223da3a660dSBarry Smith     for (i=0; i<n; i++) {
1224da3a660dSBarry Smith       x = r[i];
1225da3a660dSBarry Smith       v = mat->v + i*m;
1226da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1227da3a660dSBarry Smith     }
12281ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1229efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1230da3a660dSBarry Smith   }
12313a40ed3dSBarry Smith   PetscFunctionReturn(0);
1232289bc588SBarry Smith }
1233289bc588SBarry Smith 
12344a2ae208SSatish Balay #undef __FUNCT__
12354a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1236dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1237289bc588SBarry Smith {
1238c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
123987828ca2SBarry Smith   PetscScalar  *v = mat->v;
1240329f5518SBarry Smith   PetscReal    sum = 0.0;
1241d0f46423SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1242efee365bSSatish Balay   PetscErrorCode ierr;
124355659b69SBarry Smith 
12443a40ed3dSBarry Smith   PetscFunctionBegin;
1245289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1246a5ce6ee0Svictorle     if (lda>m) {
1247d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1248a5ce6ee0Svictorle 	v = mat->v+j*lda;
1249a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1250a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1251a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1252a5ce6ee0Svictorle #else
1253a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1254a5ce6ee0Svictorle #endif
1255a5ce6ee0Svictorle 	}
1256a5ce6ee0Svictorle       }
1257a5ce6ee0Svictorle     } else {
1258d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1259aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1260329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1261289bc588SBarry Smith #else
1262289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1263289bc588SBarry Smith #endif
1264289bc588SBarry Smith       }
1265a5ce6ee0Svictorle     }
1266064f8208SBarry Smith     *nrm = sqrt(sum);
1267*dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12683a40ed3dSBarry Smith   } else if (type == NORM_1) {
1269064f8208SBarry Smith     *nrm = 0.0;
1270d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
12711b807ce4Svictorle       v = mat->v + j*mat->lda;
1272289bc588SBarry Smith       sum = 0.0;
1273d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
127433a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1275289bc588SBarry Smith       }
1276064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1277289bc588SBarry Smith     }
1278d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12793a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1280064f8208SBarry Smith     *nrm = 0.0;
1281d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1282289bc588SBarry Smith       v = mat->v + j;
1283289bc588SBarry Smith       sum = 0.0;
1284d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
12851b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1286289bc588SBarry Smith       }
1287064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1288289bc588SBarry Smith     }
1289d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12903a40ed3dSBarry Smith   } else {
129129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1292289bc588SBarry Smith   }
12933a40ed3dSBarry Smith   PetscFunctionReturn(0);
1294289bc588SBarry Smith }
1295289bc588SBarry Smith 
12964a2ae208SSatish Balay #undef __FUNCT__
12974a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
12984e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg)
1299289bc588SBarry Smith {
1300c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
130163ba0a88SBarry Smith   PetscErrorCode ierr;
130267e560aaSBarry Smith 
13033a40ed3dSBarry Smith   PetscFunctionBegin;
1304b5a2b587SKris Buschelman   switch (op) {
1305b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13064e0d8c25SBarry Smith     aij->roworiented = flg;
1307b5a2b587SKris Buschelman     break;
1308512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1309b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13103971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13114e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1312b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1313b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
131477e54ba9SKris Buschelman   case MAT_SYMMETRIC:
131577e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13169a4540c5SBarry Smith   case MAT_HERMITIAN:
13179a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1318600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1319290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
132077e54ba9SKris Buschelman     break;
1321b5a2b587SKris Buschelman   default:
1322600fe468SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13233a40ed3dSBarry Smith   }
13243a40ed3dSBarry Smith   PetscFunctionReturn(0);
1325289bc588SBarry Smith }
1326289bc588SBarry Smith 
13274a2ae208SSatish Balay #undef __FUNCT__
13284a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1329dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13306f0a148fSBarry Smith {
1331ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13326849ba73SBarry Smith   PetscErrorCode ierr;
1333d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13343a40ed3dSBarry Smith 
13353a40ed3dSBarry Smith   PetscFunctionBegin;
1336a5ce6ee0Svictorle   if (lda>m) {
1337d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1338a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1339a5ce6ee0Svictorle     }
1340a5ce6ee0Svictorle   } else {
1341d0f46423SBarry Smith     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1342a5ce6ee0Svictorle   }
13433a40ed3dSBarry Smith   PetscFunctionReturn(0);
13446f0a148fSBarry Smith }
13456f0a148fSBarry Smith 
13464a2ae208SSatish Balay #undef __FUNCT__
13474a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1348f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13496f0a148fSBarry Smith {
1350ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1351d0f46423SBarry Smith   PetscInt       n = A->cmap->n,i,j;
135287828ca2SBarry Smith   PetscScalar    *slot;
135355659b69SBarry Smith 
13543a40ed3dSBarry Smith   PetscFunctionBegin;
13556f0a148fSBarry Smith   for (i=0; i<N; i++) {
13566f0a148fSBarry Smith     slot = l->v + rows[i];
13576f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13586f0a148fSBarry Smith   }
1359f4df32b1SMatthew Knepley   if (diag != 0.0) {
13606f0a148fSBarry Smith     for (i=0; i<N; i++) {
13616f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1362f4df32b1SMatthew Knepley       *slot = diag;
13636f0a148fSBarry Smith     }
13646f0a148fSBarry Smith   }
13653a40ed3dSBarry Smith   PetscFunctionReturn(0);
13666f0a148fSBarry Smith }
1367557bce09SLois Curfman McInnes 
13684a2ae208SSatish Balay #undef __FUNCT__
13694a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1370dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
137164e87e97SBarry Smith {
1372c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13733a40ed3dSBarry Smith 
13743a40ed3dSBarry Smith   PetscFunctionBegin;
1375d0f46423SBarry Smith   if (mat->lda != A->rmap->n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
137664e87e97SBarry Smith   *array = mat->v;
13773a40ed3dSBarry Smith   PetscFunctionReturn(0);
137864e87e97SBarry Smith }
13790754003eSLois Curfman McInnes 
13804a2ae208SSatish Balay #undef __FUNCT__
13814a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1382dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1383ff14e315SSatish Balay {
13843a40ed3dSBarry Smith   PetscFunctionBegin;
138509b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13863a40ed3dSBarry Smith   PetscFunctionReturn(0);
1387ff14e315SSatish Balay }
13880754003eSLois Curfman McInnes 
13894a2ae208SSatish Balay #undef __FUNCT__
13904a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
139113f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13920754003eSLois Curfman McInnes {
1393c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13946849ba73SBarry Smith   PetscErrorCode ierr;
13955d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
13965d0c19d7SBarry Smith   const PetscInt *irow,*icol;
139787828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13980754003eSLois Curfman McInnes   Mat            newmat;
13990754003eSLois Curfman McInnes 
14003a40ed3dSBarry Smith   PetscFunctionBegin;
140178b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
140278b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1403e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1404e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14050754003eSLois Curfman McInnes 
1406182d2002SSatish Balay   /* Check submatrixcall */
1407182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
140813f74950SBarry Smith     PetscInt n_cols,n_rows;
1409182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
141021a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
141121a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
1412c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
141321a2c019SBarry Smith     }
1414182d2002SSatish Balay     newmat = *B;
1415182d2002SSatish Balay   } else {
14160754003eSLois Curfman McInnes     /* Create and fill new matrix */
14177adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1418f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14197adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14205c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1421182d2002SSatish Balay   }
1422182d2002SSatish Balay 
1423182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1424182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1425182d2002SSatish Balay 
1426182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14276de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1428182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1429182d2002SSatish Balay       *bv++ = av[irow[j]];
14300754003eSLois Curfman McInnes     }
14310754003eSLois Curfman McInnes   }
1432182d2002SSatish Balay 
1433182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14346d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14356d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14360754003eSLois Curfman McInnes 
14370754003eSLois Curfman McInnes   /* Free work space */
143878b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
143978b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1440182d2002SSatish Balay   *B = newmat;
14413a40ed3dSBarry Smith   PetscFunctionReturn(0);
14420754003eSLois Curfman McInnes }
14430754003eSLois Curfman McInnes 
14444a2ae208SSatish Balay #undef __FUNCT__
14454a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
144613f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1447905e6a2fSBarry Smith {
14486849ba73SBarry Smith   PetscErrorCode ierr;
144913f74950SBarry Smith   PetscInt       i;
1450905e6a2fSBarry Smith 
14513a40ed3dSBarry Smith   PetscFunctionBegin;
1452905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1453b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1454905e6a2fSBarry Smith   }
1455905e6a2fSBarry Smith 
1456905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14576a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1458905e6a2fSBarry Smith   }
14593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1460905e6a2fSBarry Smith }
1461905e6a2fSBarry Smith 
14624a2ae208SSatish Balay #undef __FUNCT__
1463c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1464c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1465c0aa2d19SHong Zhang {
1466c0aa2d19SHong Zhang   PetscFunctionBegin;
1467c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1468c0aa2d19SHong Zhang }
1469c0aa2d19SHong Zhang 
1470c0aa2d19SHong Zhang #undef __FUNCT__
1471c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1472c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1473c0aa2d19SHong Zhang {
1474c0aa2d19SHong Zhang   PetscFunctionBegin;
1475c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1476c0aa2d19SHong Zhang }
1477c0aa2d19SHong Zhang 
1478c0aa2d19SHong Zhang #undef __FUNCT__
14794a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1480dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14814b0e389bSBarry Smith {
14824b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14836849ba73SBarry Smith   PetscErrorCode ierr;
1484d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
14853a40ed3dSBarry Smith 
14863a40ed3dSBarry Smith   PetscFunctionBegin;
148733f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
148833f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1489cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14903a40ed3dSBarry Smith     PetscFunctionReturn(0);
14913a40ed3dSBarry Smith   }
1492d0f46423SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1493a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14940dbb7854Svictorle     for (j=0; j<n; j++) {
1495a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1496a5ce6ee0Svictorle     }
1497a5ce6ee0Svictorle   } else {
1498d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1499a5ce6ee0Svictorle   }
1500273d9f13SBarry Smith   PetscFunctionReturn(0);
1501273d9f13SBarry Smith }
1502273d9f13SBarry Smith 
15034a2ae208SSatish Balay #undef __FUNCT__
15044a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1505dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1506273d9f13SBarry Smith {
1507dfbe8321SBarry Smith   PetscErrorCode ierr;
1508273d9f13SBarry Smith 
1509273d9f13SBarry Smith   PetscFunctionBegin;
1510273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15113a40ed3dSBarry Smith   PetscFunctionReturn(0);
15124b0e389bSBarry Smith }
15134b0e389bSBarry Smith 
1514284134d9SBarry Smith #undef __FUNCT__
1515284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1516284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1517284134d9SBarry Smith {
1518284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1519284134d9SBarry Smith   PetscFunctionBegin;
152021a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1521284134d9SBarry Smith   m = PetscMax(m,M);
1522284134d9SBarry Smith   n = PetscMax(n,N);
152321a2c019SBarry 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);
1524284134d9SBarry 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);
1525d0f46423SBarry Smith   A->rmap->n = A->rmap->n = m;
1526d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
152721a2c019SBarry Smith   if (a->changelda) a->lda = m;
1528284134d9SBarry Smith   PetscFunctionReturn(0);
1529284134d9SBarry Smith }
1530170fe5c8SBarry Smith 
1531284134d9SBarry Smith 
1532a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1533a9fe9ddaSSatish Balay #undef __FUNCT__
1534a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1535a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1536a9fe9ddaSSatish Balay {
1537a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1538a9fe9ddaSSatish Balay 
1539a9fe9ddaSSatish Balay   PetscFunctionBegin;
1540a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1541a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1542a9fe9ddaSSatish Balay   }
1543a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1544a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1545a9fe9ddaSSatish Balay }
1546a9fe9ddaSSatish Balay 
1547a9fe9ddaSSatish Balay #undef __FUNCT__
1548a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1549a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1550a9fe9ddaSSatish Balay {
1551ee16a9a1SHong Zhang   PetscErrorCode ierr;
1552d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1553ee16a9a1SHong Zhang   Mat            Cmat;
1554a9fe9ddaSSatish Balay 
1555ee16a9a1SHong Zhang   PetscFunctionBegin;
1556d0f46423SBarry Smith   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n);
1557ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1558ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1559ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1560ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1561ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1562ee16a9a1SHong Zhang   *C = Cmat;
1563ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1564ee16a9a1SHong Zhang }
1565a9fe9ddaSSatish Balay 
156698a3b096SSatish Balay #undef __FUNCT__
1567a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1568a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1569a9fe9ddaSSatish Balay {
1570a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1571a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1572a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
15730805154bSBarry Smith   PetscBLASInt   m,n,k;
1574a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1575a9fe9ddaSSatish Balay 
1576a9fe9ddaSSatish Balay   PetscFunctionBegin;
1577d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1578d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1579d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1580a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1581a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1582a9fe9ddaSSatish Balay }
1583a9fe9ddaSSatish Balay 
1584a9fe9ddaSSatish Balay #undef __FUNCT__
1585a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1586a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1587a9fe9ddaSSatish Balay {
1588a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1589a9fe9ddaSSatish Balay 
1590a9fe9ddaSSatish Balay   PetscFunctionBegin;
1591a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1592a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1593a9fe9ddaSSatish Balay   }
1594a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1595a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1596a9fe9ddaSSatish Balay }
1597a9fe9ddaSSatish Balay 
1598a9fe9ddaSSatish Balay #undef __FUNCT__
1599a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1600a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1601a9fe9ddaSSatish Balay {
1602ee16a9a1SHong Zhang   PetscErrorCode ierr;
1603d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1604ee16a9a1SHong Zhang   Mat            Cmat;
1605a9fe9ddaSSatish Balay 
1606ee16a9a1SHong Zhang   PetscFunctionBegin;
1607d0f46423SBarry Smith   if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n);
1608ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1609ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1610ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1611ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1612ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1613ee16a9a1SHong Zhang   *C = Cmat;
1614ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1615ee16a9a1SHong Zhang }
1616a9fe9ddaSSatish Balay 
1617a9fe9ddaSSatish Balay #undef __FUNCT__
1618a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1619a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1620a9fe9ddaSSatish Balay {
1621a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1622a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1623a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16240805154bSBarry Smith   PetscBLASInt   m,n,k;
1625a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1626a9fe9ddaSSatish Balay 
1627a9fe9ddaSSatish Balay   PetscFunctionBegin;
1628d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1629d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1630d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
16312fbe02b9SBarry Smith   /*
16322fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
16332fbe02b9SBarry Smith   */
1634a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1635a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1636a9fe9ddaSSatish Balay }
1637985db425SBarry Smith 
1638985db425SBarry Smith #undef __FUNCT__
1639985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1640985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1641985db425SBarry Smith {
1642985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1643985db425SBarry Smith   PetscErrorCode ierr;
1644d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1645985db425SBarry Smith   PetscScalar    *x;
1646985db425SBarry Smith   MatScalar      *aa = a->v;
1647985db425SBarry Smith 
1648985db425SBarry Smith   PetscFunctionBegin;
1649985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1650985db425SBarry Smith 
1651985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1652985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1653985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1654d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1655985db425SBarry Smith   for (i=0; i<m; i++) {
1656985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1657985db425SBarry Smith     for (j=1; j<n; j++){
1658985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1659985db425SBarry Smith     }
1660985db425SBarry Smith   }
1661985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1662985db425SBarry Smith   PetscFunctionReturn(0);
1663985db425SBarry Smith }
1664985db425SBarry Smith 
1665985db425SBarry Smith #undef __FUNCT__
1666985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1667985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1668985db425SBarry Smith {
1669985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1670985db425SBarry Smith   PetscErrorCode ierr;
1671d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1672985db425SBarry Smith   PetscScalar    *x;
1673985db425SBarry Smith   PetscReal      atmp;
1674985db425SBarry Smith   MatScalar      *aa = a->v;
1675985db425SBarry Smith 
1676985db425SBarry Smith   PetscFunctionBegin;
1677985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1678985db425SBarry Smith 
1679985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1680985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1681985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1682d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1683985db425SBarry Smith   for (i=0; i<m; i++) {
16849189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1685985db425SBarry Smith     for (j=1; j<n; j++){
1686985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1687985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1688985db425SBarry Smith     }
1689985db425SBarry Smith   }
1690985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1691985db425SBarry Smith   PetscFunctionReturn(0);
1692985db425SBarry Smith }
1693985db425SBarry Smith 
1694985db425SBarry Smith #undef __FUNCT__
1695985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1696985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1697985db425SBarry Smith {
1698985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1699985db425SBarry Smith   PetscErrorCode ierr;
1700d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1701985db425SBarry Smith   PetscScalar    *x;
1702985db425SBarry Smith   MatScalar      *aa = a->v;
1703985db425SBarry Smith 
1704985db425SBarry Smith   PetscFunctionBegin;
1705985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1706985db425SBarry Smith 
1707985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1708985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1709985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1710d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1711985db425SBarry Smith   for (i=0; i<m; i++) {
1712985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1713985db425SBarry Smith     for (j=1; j<n; j++){
1714985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1715985db425SBarry Smith     }
1716985db425SBarry Smith   }
1717985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1718985db425SBarry Smith   PetscFunctionReturn(0);
1719985db425SBarry Smith }
1720985db425SBarry Smith 
17218d0534beSBarry Smith #undef __FUNCT__
17228d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
17238d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
17248d0534beSBarry Smith {
17258d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17268d0534beSBarry Smith   PetscErrorCode ierr;
17278d0534beSBarry Smith   PetscScalar    *x;
17288d0534beSBarry Smith 
17298d0534beSBarry Smith   PetscFunctionBegin;
17308d0534beSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
17318d0534beSBarry Smith 
17328d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1733d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
17348d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
17358d0534beSBarry Smith   PetscFunctionReturn(0);
17368d0534beSBarry Smith }
17378d0534beSBarry Smith 
1738289bc588SBarry Smith /* -------------------------------------------------------------------*/
1739a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1740905e6a2fSBarry Smith        MatGetRow_SeqDense,
1741905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1742905e6a2fSBarry Smith        MatMult_SeqDense,
174397304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
17447c922b88SBarry Smith        MatMultTranspose_SeqDense,
17457c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1746db4efbfdSBarry Smith        0,
1747db4efbfdSBarry Smith        0,
1748db4efbfdSBarry Smith        0,
1749db4efbfdSBarry Smith /*10*/ 0,
1750905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1751905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1752ec8511deSBarry Smith        MatRelax_SeqDense,
1753ec8511deSBarry Smith        MatTranspose_SeqDense,
175497304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1755905e6a2fSBarry Smith        MatEqual_SeqDense,
1756905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1757905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1758905e6a2fSBarry Smith        MatNorm_SeqDense,
1759c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1760c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1761905e6a2fSBarry Smith        0,
1762905e6a2fSBarry Smith        MatSetOption_SeqDense,
1763905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
176497304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1765db4efbfdSBarry Smith        0,
1766db4efbfdSBarry Smith        0,
1767db4efbfdSBarry Smith        0,
1768db4efbfdSBarry Smith        0,
176997304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1770273d9f13SBarry Smith        0,
1771905e6a2fSBarry Smith        0,
1772905e6a2fSBarry Smith        MatGetArray_SeqDense,
1773905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
177497304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1775a5ae1ecdSBarry Smith        0,
1776a5ae1ecdSBarry Smith        0,
1777a5ae1ecdSBarry Smith        0,
1778a5ae1ecdSBarry Smith        0,
177997304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1780a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1781a5ae1ecdSBarry Smith        0,
17824b0e389bSBarry Smith        MatGetValues_SeqDense,
1783a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1784985db425SBarry Smith /*45*/ MatGetRowMax_SeqDense,
1785a5ae1ecdSBarry Smith        MatScale_SeqDense,
1786a5ae1ecdSBarry Smith        0,
1787a5ae1ecdSBarry Smith        0,
1788a5ae1ecdSBarry Smith        0,
1789521d7252SBarry Smith /*50*/ 0,
1790a5ae1ecdSBarry Smith        0,
1791a5ae1ecdSBarry Smith        0,
1792a5ae1ecdSBarry Smith        0,
1793a5ae1ecdSBarry Smith        0,
179497304618SKris Buschelman /*55*/ 0,
1795a5ae1ecdSBarry Smith        0,
1796a5ae1ecdSBarry Smith        0,
1797a5ae1ecdSBarry Smith        0,
1798a5ae1ecdSBarry Smith        0,
179997304618SKris Buschelman /*60*/ 0,
1800e03a110bSBarry Smith        MatDestroy_SeqDense,
1801e03a110bSBarry Smith        MatView_SeqDense,
1802357abbc8SBarry Smith        0,
180397304618SKris Buschelman        0,
180497304618SKris Buschelman /*65*/ 0,
180597304618SKris Buschelman        0,
180697304618SKris Buschelman        0,
180797304618SKris Buschelman        0,
180897304618SKris Buschelman        0,
1809985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqDense,
181097304618SKris Buschelman        0,
181197304618SKris Buschelman        0,
181297304618SKris Buschelman        0,
181397304618SKris Buschelman        0,
181497304618SKris Buschelman /*75*/ 0,
181597304618SKris Buschelman        0,
181697304618SKris Buschelman        0,
181797304618SKris Buschelman        0,
181897304618SKris Buschelman        0,
181997304618SKris Buschelman /*80*/ 0,
182097304618SKris Buschelman        0,
182197304618SKris Buschelman        0,
182297304618SKris Buschelman        0,
1823865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1824865e5f61SKris Buschelman        0,
18251cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1826865e5f61SKris Buschelman        0,
1827865e5f61SKris Buschelman        0,
1828865e5f61SKris Buschelman        0,
1829a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense,
1830a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1831a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1832865e5f61SKris Buschelman        0,
1833865e5f61SKris Buschelman        0,
1834865e5f61SKris Buschelman /*95*/ 0,
1835a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1836a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1837a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1838284134d9SBarry Smith        0,
1839284134d9SBarry Smith /*100*/0,
1840284134d9SBarry Smith        0,
1841284134d9SBarry Smith        0,
1842284134d9SBarry Smith        0,
1843985db425SBarry Smith        MatSetSizes_SeqDense,
1844985db425SBarry Smith        0,
1845985db425SBarry Smith        0,
1846985db425SBarry Smith        0,
1847985db425SBarry Smith        0,
1848985db425SBarry Smith        0,
1849985db425SBarry Smith /*110*/0,
1850985db425SBarry Smith        0,
18518d0534beSBarry Smith        MatGetRowMin_SeqDense,
18528d0534beSBarry Smith        MatGetColumnVector_SeqDense
1853985db425SBarry Smith };
185490ace30eSBarry Smith 
18554a2ae208SSatish Balay #undef __FUNCT__
18564a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18574b828684SBarry Smith /*@C
1858fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1859d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1860d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1861289bc588SBarry Smith 
1862db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1863db81eaa0SLois Curfman McInnes 
186420563c6bSBarry Smith    Input Parameters:
1865db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
18660c775827SLois Curfman McInnes .  m - number of rows
186718f449edSLois Curfman McInnes .  n - number of columns
1868db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1869dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
187020563c6bSBarry Smith 
187120563c6bSBarry Smith    Output Parameter:
187244cd7ae7SLois Curfman McInnes .  A - the matrix
187320563c6bSBarry Smith 
1874b259b22eSLois Curfman McInnes    Notes:
187518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
187618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1877b4fd4287SBarry Smith    set data=PETSC_NULL.
187818f449edSLois Curfman McInnes 
1879027ccd11SLois Curfman McInnes    Level: intermediate
1880027ccd11SLois Curfman McInnes 
1881dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1882d65003e9SLois Curfman McInnes 
1883db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
188420563c6bSBarry Smith @*/
1885be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1886289bc588SBarry Smith {
1887dfbe8321SBarry Smith   PetscErrorCode ierr;
18883b2fbd54SBarry Smith 
18893a40ed3dSBarry Smith   PetscFunctionBegin;
1890f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1891f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1892273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1893273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1894273d9f13SBarry Smith   PetscFunctionReturn(0);
1895273d9f13SBarry Smith }
1896273d9f13SBarry Smith 
18974a2ae208SSatish Balay #undef __FUNCT__
1898afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1899273d9f13SBarry Smith /*@C
1900273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1901273d9f13SBarry Smith 
1902273d9f13SBarry Smith    Collective on MPI_Comm
1903273d9f13SBarry Smith 
1904273d9f13SBarry Smith    Input Parameters:
1905273d9f13SBarry Smith +  A - the matrix
1906273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1907273d9f13SBarry Smith 
1908273d9f13SBarry Smith    Notes:
1909273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1910273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1911284134d9SBarry Smith    need not call this routine.
1912273d9f13SBarry Smith 
1913273d9f13SBarry Smith    Level: intermediate
1914273d9f13SBarry Smith 
1915273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1916273d9f13SBarry Smith 
1917867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
1918867c911aSBarry Smith 
1919273d9f13SBarry Smith @*/
1920be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1921273d9f13SBarry Smith {
1922dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1923a23d5eceSKris Buschelman 
1924a23d5eceSKris Buschelman   PetscFunctionBegin;
1925a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1926a23d5eceSKris Buschelman   if (f) {
1927a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1928a23d5eceSKris Buschelman   }
1929a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1930a23d5eceSKris Buschelman }
1931a23d5eceSKris Buschelman 
1932a23d5eceSKris Buschelman EXTERN_C_BEGIN
1933a23d5eceSKris Buschelman #undef __FUNCT__
1934afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1935be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1936a23d5eceSKris Buschelman {
1937273d9f13SBarry Smith   Mat_SeqDense   *b;
1938dfbe8321SBarry Smith   PetscErrorCode ierr;
1939273d9f13SBarry Smith 
1940273d9f13SBarry Smith   PetscFunctionBegin;
1941273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1942273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1943d0f46423SBarry Smith   if (b->lda <= 0) b->lda = B->rmap->n;
19449e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
19459e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
19465afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1947284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1948284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
19499e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
1950273d9f13SBarry Smith   } else { /* user-allocated storage */
19519e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1952273d9f13SBarry Smith     b->v          = data;
1953273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1954273d9f13SBarry Smith   }
19550450473dSBarry Smith   B->assembled = PETSC_TRUE;
1956273d9f13SBarry Smith   PetscFunctionReturn(0);
1957273d9f13SBarry Smith }
1958a23d5eceSKris Buschelman EXTERN_C_END
1959273d9f13SBarry Smith 
19601b807ce4Svictorle #undef __FUNCT__
19611b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
19621b807ce4Svictorle /*@C
19631b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
19641b807ce4Svictorle 
19651b807ce4Svictorle   Input parameter:
19661b807ce4Svictorle + A - the matrix
19671b807ce4Svictorle - lda - the leading dimension
19681b807ce4Svictorle 
19691b807ce4Svictorle   Notes:
1970867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
19711b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
19721b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
19731b807ce4Svictorle 
19741b807ce4Svictorle   Level: intermediate
19751b807ce4Svictorle 
19761b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
19771b807ce4Svictorle 
1978284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
1979867c911aSBarry Smith 
19801b807ce4Svictorle @*/
1981be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
19821b807ce4Svictorle {
19831b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
198421a2c019SBarry Smith 
19851b807ce4Svictorle   PetscFunctionBegin;
1986d0f46423SBarry Smith   if (lda < B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
19871b807ce4Svictorle   b->lda       = lda;
198821a2c019SBarry Smith   b->changelda = PETSC_FALSE;
198921a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
19901b807ce4Svictorle   PetscFunctionReturn(0);
19911b807ce4Svictorle }
19921b807ce4Svictorle 
19930bad9183SKris Buschelman /*MC
1994fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
19950bad9183SKris Buschelman 
19960bad9183SKris Buschelman    Options Database Keys:
19970bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
19980bad9183SKris Buschelman 
19990bad9183SKris Buschelman   Level: beginner
20000bad9183SKris Buschelman 
200189665df3SBarry Smith .seealso: MatCreateSeqDense()
200289665df3SBarry Smith 
20030bad9183SKris Buschelman M*/
20040bad9183SKris Buschelman 
2005273d9f13SBarry Smith EXTERN_C_BEGIN
20064a2ae208SSatish Balay #undef __FUNCT__
20074a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
2008be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2009273d9f13SBarry Smith {
2010273d9f13SBarry Smith   Mat_SeqDense   *b;
2011dfbe8321SBarry Smith   PetscErrorCode ierr;
20127c334f02SBarry Smith   PetscMPIInt    size;
2013273d9f13SBarry Smith 
2014273d9f13SBarry Smith   PetscFunctionBegin;
20157adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
201629bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
201755659b69SBarry Smith 
20187408324eSLisandro Dalcin   ierr = PetscMapSetBlockSize(B->rmap,1);CHKERRQ(ierr);
20197408324eSLisandro Dalcin   ierr = PetscMapSetBlockSize(B->cmap,1);CHKERRQ(ierr);
2020d0f46423SBarry Smith   ierr = PetscMapSetUp(B->rmap);CHKERRQ(ierr);
2021d0f46423SBarry Smith   ierr = PetscMapSetUp(B->cmap);CHKERRQ(ierr);
2022273d9f13SBarry Smith 
202338f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2024549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
202590f02eecSBarry Smith   B->mapping      = 0;
202644cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
202718f449edSLois Curfman McInnes 
2028a5ae1ecdSBarry Smith 
202944cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2030273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2031273d9f13SBarry Smith   b->v            = 0;
2032d0f46423SBarry Smith   b->lda          = B->rmap->n;
203321a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
2034d0f46423SBarry Smith   b->Mmax         = B->rmap->n;
2035d0f46423SBarry Smith   b->Nmax         = B->cmap->n;
20364e220ebcSLois Curfman McInnes 
2037b24902e0SBarry Smith 
2038b24902e0SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqdense_petsc_C",
2039b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2040b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2041a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2042a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2043a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
20444ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
20454ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
20464ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
20474ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
20484ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
20494ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
20504ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
20514ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
20524ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
205317667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
20543a40ed3dSBarry Smith   PetscFunctionReturn(0);
2055289bc588SBarry Smith }
2056c0aa2d19SHong Zhang 
2057c0aa2d19SHong Zhang 
2058273d9f13SBarry Smith EXTERN_C_END
2059