xref: /petsc/src/mat/impls/dense/seq/dense.c (revision c3ef05f6a9524b64ba2cd8698cb1d10c96c5a834)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
8b4c862e0SSatish Balay #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
170805154bSBarry Smith   PetscBLASInt   N,m,ldax,lday,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
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   }
32efee365bSSatish Balay   ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr);
333a40ed3dSBarry Smith   PetscFunctionReturn(0);
341987afe7SBarry Smith }
351987afe7SBarry Smith 
364a2ae208SSatish Balay #undef __FUNCT__
374a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
38dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
39289bc588SBarry Smith {
40d0f46423SBarry Smith   PetscInt     N = A->rmap->n*A->cmap->n;
413a40ed3dSBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
43d0f46423SBarry Smith   info->rows_global       = (double)A->rmap->n;
44d0f46423SBarry Smith   info->columns_global    = (double)A->cmap->n;
45d0f46423SBarry Smith   info->rows_local        = (double)A->rmap->n;
46d0f46423SBarry Smith   info->columns_local     = (double)A->cmap->n;
474e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
484e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
496de62eeeSBarry Smith   info->nz_used           = (double)N;
506de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
514e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
524e220ebcSLois Curfman McInnes   info->mallocs           = 0;
537adad957SLisandro Dalcin   info->memory            = ((PetscObject)A)->mem;
544e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
554e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
564e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
573a40ed3dSBarry Smith   PetscFunctionReturn(0);
58289bc588SBarry Smith }
59289bc588SBarry Smith 
604a2ae208SSatish Balay #undef __FUNCT__
614a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
62f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
6380cd9d93SLois Curfman McInnes {
64273d9f13SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
65f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
66efee365bSSatish Balay   PetscErrorCode ierr;
670805154bSBarry Smith   PetscBLASInt   one = 1,j,nz,lda = PetscBLASIntCast(a->lda);
6880cd9d93SLois Curfman McInnes 
693a40ed3dSBarry Smith   PetscFunctionBegin;
70d0f46423SBarry Smith   if (lda>A->rmap->n) {
71d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n);
72d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
73f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
74a5ce6ee0Svictorle     }
75a5ce6ee0Svictorle   } else {
76d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n*A->cmap->n);
77f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
78a5ce6ee0Svictorle   }
79efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
803a40ed3dSBarry Smith   PetscFunctionReturn(0);
8180cd9d93SLois Curfman McInnes }
8280cd9d93SLois Curfman McInnes 
831cbb95d3SBarry Smith #undef __FUNCT__
841cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense"
851cbb95d3SBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscTruth *fl)
861cbb95d3SBarry Smith {
871cbb95d3SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
88d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,N;
891cbb95d3SBarry Smith   PetscScalar    *v = a->v;
901cbb95d3SBarry Smith 
911cbb95d3SBarry Smith   PetscFunctionBegin;
921cbb95d3SBarry Smith   *fl = PETSC_FALSE;
93d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0);
941cbb95d3SBarry Smith   N = a->lda;
951cbb95d3SBarry Smith 
961cbb95d3SBarry Smith   for (i=0; i<m; i++) {
971cbb95d3SBarry Smith     for (j=i+1; j<m; j++) {
981cbb95d3SBarry Smith       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
991cbb95d3SBarry Smith     }
1001cbb95d3SBarry Smith   }
1011cbb95d3SBarry Smith   *fl = PETSC_TRUE;
1021cbb95d3SBarry Smith   PetscFunctionReturn(0);
1031cbb95d3SBarry Smith }
1041cbb95d3SBarry Smith 
105b24902e0SBarry Smith #undef __FUNCT__
106b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense"
107719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
108b24902e0SBarry Smith {
109b24902e0SBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
110b24902e0SBarry Smith   PetscErrorCode ierr;
111b24902e0SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
112b24902e0SBarry Smith 
113b24902e0SBarry Smith   PetscFunctionBegin;
114719d5645SBarry Smith   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
115b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
116b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
117d0f46423SBarry Smith     if (lda>A->rmap->n) {
118d0f46423SBarry Smith       m = A->rmap->n;
119d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
120b24902e0SBarry Smith 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
121b24902e0SBarry Smith       }
122b24902e0SBarry Smith     } else {
123d0f46423SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
124b24902e0SBarry Smith     }
125b24902e0SBarry Smith   }
126b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
127b24902e0SBarry Smith   PetscFunctionReturn(0);
128b24902e0SBarry Smith }
129b24902e0SBarry Smith 
1304a2ae208SSatish Balay #undef __FUNCT__
1314a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
132dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
13302cad45dSBarry Smith {
1346849ba73SBarry Smith   PetscErrorCode ierr;
13502cad45dSBarry Smith 
1363a40ed3dSBarry Smith   PetscFunctionBegin;
1375c9eb25fSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
138d0f46423SBarry Smith   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1395c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
140719d5645SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
141b24902e0SBarry Smith   PetscFunctionReturn(0);
142b24902e0SBarry Smith }
143b24902e0SBarry Smith 
1446ee01492SSatish Balay 
1450481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
146719d5645SBarry Smith 
1474a2ae208SSatish Balay #undef __FUNCT__
1484a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1490481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
150289bc588SBarry Smith {
1514482741eSBarry Smith   MatFactorInfo  info;
152a093e273SMatthew Knepley   PetscErrorCode ierr;
1533a40ed3dSBarry Smith 
1543a40ed3dSBarry Smith   PetscFunctionBegin;
155*c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
156719d5645SBarry Smith   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
1573a40ed3dSBarry Smith   PetscFunctionReturn(0);
158289bc588SBarry Smith }
1596ee01492SSatish Balay 
1600b4b3355SBarry Smith #undef __FUNCT__
1614a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
162dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
163289bc588SBarry Smith {
164c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1656849ba73SBarry Smith   PetscErrorCode ierr;
16687828ca2SBarry Smith   PetscScalar    *x,*y;
167d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
16867e560aaSBarry Smith 
1693a40ed3dSBarry Smith   PetscFunctionBegin;
1701ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1711ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
172d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1735c9eb25fSBarry Smith   if (A->factor == MAT_FACTOR_LU) {
174ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
17580444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
176ae7cfcebSSatish Balay #else
17771044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
1784ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
179ae7cfcebSSatish Balay #endif
1805c9eb25fSBarry Smith   } else if (A->factor == MAT_FACTOR_CHOLESKY){
181ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
18280444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
183ae7cfcebSSatish Balay #else
18471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
1854ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
186ae7cfcebSSatish Balay #endif
187289bc588SBarry Smith   }
18829bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
1891ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1901ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
191d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
1923a40ed3dSBarry Smith   PetscFunctionReturn(0);
193289bc588SBarry Smith }
1946ee01492SSatish Balay 
1954a2ae208SSatish Balay #undef __FUNCT__
1964a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
197dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
198da3a660dSBarry Smith {
199c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
200dfbe8321SBarry Smith   PetscErrorCode ierr;
20187828ca2SBarry Smith   PetscScalar    *x,*y;
202d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
20367e560aaSBarry Smith 
2043a40ed3dSBarry Smith   PetscFunctionBegin;
2051ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2061ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
207d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
208752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
209da3a660dSBarry Smith   if (mat->pivots) {
210ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
21180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
212ae7cfcebSSatish Balay #else
21371044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2144ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
215ae7cfcebSSatish Balay #endif
2167a97a34bSBarry Smith   } else {
217ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
21880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
219ae7cfcebSSatish Balay #else
22071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2214ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
222ae7cfcebSSatish Balay #endif
223da3a660dSBarry Smith   }
2241ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2251ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
226d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
2273a40ed3dSBarry Smith   PetscFunctionReturn(0);
228da3a660dSBarry Smith }
2296ee01492SSatish Balay 
2304a2ae208SSatish Balay #undef __FUNCT__
2314a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
232dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
233da3a660dSBarry Smith {
234c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
235dfbe8321SBarry Smith   PetscErrorCode ierr;
23687828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
237da3a660dSBarry Smith   Vec            tmp = 0;
238d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
23967e560aaSBarry Smith 
2403a40ed3dSBarry Smith   PetscFunctionBegin;
2411ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2421ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
243d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
244da3a660dSBarry Smith   if (yy == zz) {
24578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
24652e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
24778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
248da3a660dSBarry Smith   }
249d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
250752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
251da3a660dSBarry Smith   if (mat->pivots) {
252ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
25380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
254ae7cfcebSSatish Balay #else
25571044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2562ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
257ae7cfcebSSatish Balay #endif
258a8c6a408SBarry Smith   } else {
259ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
26080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
261ae7cfcebSSatish Balay #else
26271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2632ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
264ae7cfcebSSatish Balay #endif
265da3a660dSBarry Smith   }
2662dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
2672dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
2681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2691ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
270d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
2713a40ed3dSBarry Smith   PetscFunctionReturn(0);
272da3a660dSBarry Smith }
27367e560aaSBarry Smith 
2744a2ae208SSatish Balay #undef __FUNCT__
2754a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
276dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
277da3a660dSBarry Smith {
278c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2796849ba73SBarry Smith   PetscErrorCode ierr;
28087828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
281da3a660dSBarry Smith   Vec            tmp;
282d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
28367e560aaSBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
285d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
2861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2871ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
288da3a660dSBarry Smith   if (yy == zz) {
28978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
29052e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
29178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
292da3a660dSBarry Smith   }
293d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
294752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
295da3a660dSBarry Smith   if (mat->pivots) {
296ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
29780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
298ae7cfcebSSatish Balay #else
29971044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3002ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
301ae7cfcebSSatish Balay #endif
3023a40ed3dSBarry Smith   } else {
303ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
30480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
305ae7cfcebSSatish Balay #else
30671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3072ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
308ae7cfcebSSatish Balay #endif
309da3a660dSBarry Smith   }
31090f02eecSBarry Smith   if (tmp) {
3112dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
31290f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3133a40ed3dSBarry Smith   } else {
3142dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
31590f02eecSBarry Smith   }
3161ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3171ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
318d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
320da3a660dSBarry Smith }
321db4efbfdSBarry Smith 
322db4efbfdSBarry Smith /* ---------------------------------------------------------------*/
323db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
324db4efbfdSBarry Smith    rather than put it in the Mat->row slot.*/
325db4efbfdSBarry Smith #undef __FUNCT__
326db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense"
3270481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
328db4efbfdSBarry Smith {
329db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
330db4efbfdSBarry Smith   PetscFunctionBegin;
331db4efbfdSBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
332db4efbfdSBarry Smith #else
333db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
334db4efbfdSBarry Smith   PetscErrorCode ierr;
335db4efbfdSBarry Smith   PetscBLASInt   n,m,info;
336db4efbfdSBarry Smith 
337db4efbfdSBarry Smith   PetscFunctionBegin;
338db4efbfdSBarry Smith   n = PetscBLASIntCast(A->cmap->n);
339db4efbfdSBarry Smith   m = PetscBLASIntCast(A->rmap->n);
340db4efbfdSBarry Smith   if (!mat->pivots) {
341db4efbfdSBarry Smith     ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
342db4efbfdSBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
343db4efbfdSBarry Smith   }
344db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
345db4efbfdSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
346db4efbfdSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
347db4efbfdSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
348db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
349db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
350db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
351db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
352db4efbfdSBarry Smith   A->factor = MAT_FACTOR_LU;
353db4efbfdSBarry Smith 
354db4efbfdSBarry Smith   ierr = PetscLogFlops((2*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
355db4efbfdSBarry Smith #endif
356db4efbfdSBarry Smith   PetscFunctionReturn(0);
357db4efbfdSBarry Smith }
358db4efbfdSBarry Smith 
359db4efbfdSBarry Smith #undef __FUNCT__
360db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense"
3610481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
362db4efbfdSBarry Smith {
363db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
364db4efbfdSBarry Smith   PetscFunctionBegin;
365db4efbfdSBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
366db4efbfdSBarry Smith #else
367db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
368db4efbfdSBarry Smith   PetscErrorCode ierr;
369db4efbfdSBarry Smith   PetscBLASInt   info,n = PetscBLASIntCast(A->cmap->n);
370db4efbfdSBarry Smith 
371db4efbfdSBarry Smith   PetscFunctionBegin;
372db4efbfdSBarry Smith   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
373db4efbfdSBarry Smith   mat->pivots = 0;
374db4efbfdSBarry Smith 
375db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
376db4efbfdSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
377db4efbfdSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
378db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
379db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
380db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
381db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
382db4efbfdSBarry Smith   A->factor = MAT_FACTOR_CHOLESKY;
383db4efbfdSBarry Smith   ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
384db4efbfdSBarry Smith #endif
385db4efbfdSBarry Smith   PetscFunctionReturn(0);
386db4efbfdSBarry Smith }
387db4efbfdSBarry Smith 
388db4efbfdSBarry Smith 
389db4efbfdSBarry Smith #undef __FUNCT__
390db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
3910481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
392db4efbfdSBarry Smith {
393db4efbfdSBarry Smith   PetscErrorCode ierr;
394db4efbfdSBarry Smith   MatFactorInfo  info;
395db4efbfdSBarry Smith 
396db4efbfdSBarry Smith   PetscFunctionBegin;
397db4efbfdSBarry Smith   info.fill = 1.0;
398*c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
399719d5645SBarry Smith   ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr);
400db4efbfdSBarry Smith   PetscFunctionReturn(0);
401db4efbfdSBarry Smith }
402db4efbfdSBarry Smith 
403db4efbfdSBarry Smith #undef __FUNCT__
404db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
4050481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
406db4efbfdSBarry Smith {
407db4efbfdSBarry Smith   PetscFunctionBegin;
408*c3ef05f6SHong Zhang   fact->assembled                  = PETSC_TRUE;
409719d5645SBarry Smith   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
410db4efbfdSBarry Smith   PetscFunctionReturn(0);
411db4efbfdSBarry Smith }
412db4efbfdSBarry Smith 
413db4efbfdSBarry Smith #undef __FUNCT__
414db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
4150481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
416db4efbfdSBarry Smith {
417db4efbfdSBarry Smith   PetscFunctionBegin;
418*c3ef05f6SHong Zhang   fact->assembled            = PETSC_TRUE;
419719d5645SBarry Smith   fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
420db4efbfdSBarry Smith   PetscFunctionReturn(0);
421db4efbfdSBarry Smith }
422db4efbfdSBarry Smith 
423db4efbfdSBarry Smith #undef __FUNCT__
424db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc"
425db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
426db4efbfdSBarry Smith {
427db4efbfdSBarry Smith   PetscErrorCode ierr;
428db4efbfdSBarry Smith 
429db4efbfdSBarry Smith   PetscFunctionBegin;
430db4efbfdSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
431db4efbfdSBarry Smith   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
432db4efbfdSBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
433db4efbfdSBarry Smith   if (ftype == MAT_FACTOR_LU){
434db4efbfdSBarry Smith     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
435db4efbfdSBarry Smith   } else {
436db4efbfdSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
437db4efbfdSBarry Smith   }
438db4efbfdSBarry Smith   (*fact)->factor = ftype;
439db4efbfdSBarry Smith   PetscFunctionReturn(0);
440db4efbfdSBarry Smith }
441db4efbfdSBarry Smith 
442289bc588SBarry Smith /* ------------------------------------------------------------------*/
4434a2ae208SSatish Balay #undef __FUNCT__
4444a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
44513f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
446289bc588SBarry Smith {
447c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
44887828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
449dfbe8321SBarry Smith   PetscErrorCode ierr;
450d0f46423SBarry Smith   PetscInt       m = A->rmap->n,i;
451aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4520805154bSBarry Smith   PetscBLASInt   o = 1,bm = PetscBLASIntCast(m);
453bc1b551cSSatish Balay #endif
454289bc588SBarry Smith 
4553a40ed3dSBarry Smith   PetscFunctionBegin;
456289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
45771044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4582dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
459289bc588SBarry Smith   }
4601ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4611ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
462b965ef7fSBarry Smith   its  = its*lits;
46377431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
464289bc588SBarry Smith   while (its--) {
465fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
466289bc588SBarry Smith       for (i=0; i<m; i++) {
467aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
468f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
469f1747703SBarry Smith            not happy about returning a double complex */
47013f74950SBarry Smith         PetscInt    _i;
47187828ca2SBarry Smith         PetscScalar sum = b[i];
472f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4733f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
474f1747703SBarry Smith         }
475f1747703SBarry Smith         xt = sum;
476f1747703SBarry Smith #else
47771044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
478f1747703SBarry Smith #endif
47955a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
480289bc588SBarry Smith       }
481289bc588SBarry Smith     }
482fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
483289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
484aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
485f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
486f1747703SBarry Smith            not happy about returning a double complex */
48713f74950SBarry Smith         PetscInt    _i;
48887828ca2SBarry Smith         PetscScalar sum = b[i];
489f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4903f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
491f1747703SBarry Smith         }
492f1747703SBarry Smith         xt = sum;
493f1747703SBarry Smith #else
49471044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
495f1747703SBarry Smith #endif
49655a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
497289bc588SBarry Smith       }
498289bc588SBarry Smith     }
499289bc588SBarry Smith   }
5001ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
5011ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5023a40ed3dSBarry Smith   PetscFunctionReturn(0);
503289bc588SBarry Smith }
504289bc588SBarry Smith 
505289bc588SBarry Smith /* -----------------------------------------------------------------*/
5064a2ae208SSatish Balay #undef __FUNCT__
5074a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
508dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
509289bc588SBarry Smith {
510c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
51187828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
512dfbe8321SBarry Smith   PetscErrorCode ierr;
5130805154bSBarry Smith   PetscBLASInt   m, n,_One=1;
514ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
5153a40ed3dSBarry Smith 
5163a40ed3dSBarry Smith   PetscFunctionBegin;
517d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
518d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
519d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5201ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5211ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
52271044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
5231ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5241ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
525d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
5263a40ed3dSBarry Smith   PetscFunctionReturn(0);
527289bc588SBarry Smith }
5286ee01492SSatish Balay 
5294a2ae208SSatish Balay #undef __FUNCT__
5304a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
531dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
532289bc588SBarry Smith {
533c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
53487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
535dfbe8321SBarry Smith   PetscErrorCode ierr;
5360805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5373a40ed3dSBarry Smith 
5383a40ed3dSBarry Smith   PetscFunctionBegin;
539d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
540d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
541d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5421ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5431ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
54471044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5451ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5461ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
547d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
5483a40ed3dSBarry Smith   PetscFunctionReturn(0);
549289bc588SBarry Smith }
5506ee01492SSatish Balay 
5514a2ae208SSatish Balay #undef __FUNCT__
5524a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
553dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
554289bc588SBarry Smith {
555c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
55687828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
557dfbe8321SBarry Smith   PetscErrorCode ierr;
5580805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5593a40ed3dSBarry Smith 
5603a40ed3dSBarry Smith   PetscFunctionBegin;
561d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
562d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
563d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
564600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5651ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5661ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
56771044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5691ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
570d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5713a40ed3dSBarry Smith   PetscFunctionReturn(0);
572289bc588SBarry Smith }
5736ee01492SSatish Balay 
5744a2ae208SSatish Balay #undef __FUNCT__
5754a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
576dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
577289bc588SBarry Smith {
578c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
57987828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
580dfbe8321SBarry Smith   PetscErrorCode ierr;
5810805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
58287828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5833a40ed3dSBarry Smith 
5843a40ed3dSBarry Smith   PetscFunctionBegin;
585d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
586d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
587d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
588600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5891ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5901ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
59171044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5921ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5931ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
594d0f46423SBarry Smith   ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5953a40ed3dSBarry Smith   PetscFunctionReturn(0);
596289bc588SBarry Smith }
597289bc588SBarry Smith 
598289bc588SBarry Smith /* -----------------------------------------------------------------*/
5994a2ae208SSatish Balay #undef __FUNCT__
6004a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
60113f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
602289bc588SBarry Smith {
603c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
60487828ca2SBarry Smith   PetscScalar    *v;
6056849ba73SBarry Smith   PetscErrorCode ierr;
60613f74950SBarry Smith   PetscInt       i;
60767e560aaSBarry Smith 
6083a40ed3dSBarry Smith   PetscFunctionBegin;
609d0f46423SBarry Smith   *ncols = A->cmap->n;
610289bc588SBarry Smith   if (cols) {
611d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
612d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
613289bc588SBarry Smith   }
614289bc588SBarry Smith   if (vals) {
615d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
616289bc588SBarry Smith     v    = mat->v + row;
617d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
618289bc588SBarry Smith   }
6193a40ed3dSBarry Smith   PetscFunctionReturn(0);
620289bc588SBarry Smith }
6216ee01492SSatish Balay 
6224a2ae208SSatish Balay #undef __FUNCT__
6234a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
62413f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
625289bc588SBarry Smith {
626dfbe8321SBarry Smith   PetscErrorCode ierr;
627606d414cSSatish Balay   PetscFunctionBegin;
628606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
629606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6303a40ed3dSBarry Smith   PetscFunctionReturn(0);
631289bc588SBarry Smith }
632289bc588SBarry Smith /* ----------------------------------------------------------------*/
6334a2ae208SSatish Balay #undef __FUNCT__
6344a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
63513f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
636289bc588SBarry Smith {
637c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
63813f74950SBarry Smith   PetscInt     i,j,idx=0;
639d6dfbf8fSBarry Smith 
6403a40ed3dSBarry Smith   PetscFunctionBegin;
641289bc588SBarry Smith   if (!mat->roworiented) {
642dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
643289bc588SBarry Smith       for (j=0; j<n; j++) {
644cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6452515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
646d0f46423SBarry 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);
64758804f6dSBarry Smith #endif
648289bc588SBarry Smith         for (i=0; i<m; i++) {
649cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6502515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
651d0f46423SBarry 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);
65258804f6dSBarry Smith #endif
653cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
654289bc588SBarry Smith         }
655289bc588SBarry Smith       }
6563a40ed3dSBarry Smith     } else {
657289bc588SBarry Smith       for (j=0; j<n; j++) {
658cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6592515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
660d0f46423SBarry 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);
66158804f6dSBarry Smith #endif
662289bc588SBarry Smith         for (i=0; i<m; i++) {
663cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6642515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
665d0f46423SBarry 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);
66658804f6dSBarry Smith #endif
667cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
668289bc588SBarry Smith         }
669289bc588SBarry Smith       }
670289bc588SBarry Smith     }
6713a40ed3dSBarry Smith   } else {
672dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
673e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
674cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6752515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
676d0f46423SBarry 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);
67758804f6dSBarry Smith #endif
678e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
679cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6802515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
681d0f46423SBarry 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);
68258804f6dSBarry Smith #endif
683cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
684e8d4e0b9SBarry Smith         }
685e8d4e0b9SBarry Smith       }
6863a40ed3dSBarry Smith     } else {
687289bc588SBarry Smith       for (i=0; i<m; i++) {
688cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6892515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
690d0f46423SBarry 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);
69158804f6dSBarry Smith #endif
692289bc588SBarry Smith         for (j=0; j<n; j++) {
693cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6942515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
695d0f46423SBarry 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);
69658804f6dSBarry Smith #endif
697cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
698289bc588SBarry Smith         }
699289bc588SBarry Smith       }
700289bc588SBarry Smith     }
701e8d4e0b9SBarry Smith   }
7023a40ed3dSBarry Smith   PetscFunctionReturn(0);
703289bc588SBarry Smith }
704e8d4e0b9SBarry Smith 
7054a2ae208SSatish Balay #undef __FUNCT__
7064a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
70713f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
708ae80bb75SLois Curfman McInnes {
709ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
71013f74950SBarry Smith   PetscInt     i,j;
711ae80bb75SLois Curfman McInnes 
7123a40ed3dSBarry Smith   PetscFunctionBegin;
713ae80bb75SLois Curfman McInnes   /* row-oriented output */
714ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
71597e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
716d0f46423SBarry 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);
717ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7186f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
719d0f46423SBarry 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);
72097e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
721ae80bb75SLois Curfman McInnes     }
722ae80bb75SLois Curfman McInnes   }
7233a40ed3dSBarry Smith   PetscFunctionReturn(0);
724ae80bb75SLois Curfman McInnes }
725ae80bb75SLois Curfman McInnes 
726289bc588SBarry Smith /* -----------------------------------------------------------------*/
727289bc588SBarry Smith 
728e090d566SSatish Balay #include "petscsys.h"
72988e32edaSLois Curfman McInnes 
7304a2ae208SSatish Balay #undef __FUNCT__
7314a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
732a313700dSBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, const MatType type,Mat *A)
73388e32edaSLois Curfman McInnes {
73488e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
73588e32edaSLois Curfman McInnes   Mat            B;
7366849ba73SBarry Smith   PetscErrorCode ierr;
73713f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
73813f74950SBarry Smith   int            fd;
73913f74950SBarry Smith   PetscMPIInt    size;
74013f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
74187828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
74219bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
74388e32edaSLois Curfman McInnes 
7443a40ed3dSBarry Smith   PetscFunctionBegin;
745d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
74629bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
747b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
7480752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
749552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
75088e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
75188e32edaSLois Curfman McInnes 
75290ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
753f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
754f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
755be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7565c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
75790ace30eSBarry Smith     B    = *A;
75890ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
7594a41a4d3SSatish Balay     v    = a->v;
7604a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7614a41a4d3SSatish Balay        from row major to column major */
762b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
76390ace30eSBarry Smith     /* read in nonzero values */
7644a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7654a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7664a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7674a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7684a41a4d3SSatish Balay         *v++ =w[i*N+j];
7694a41a4d3SSatish Balay       }
7704a41a4d3SSatish Balay     }
771606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7726d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7736d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77490ace30eSBarry Smith   } else {
77588e32edaSLois Curfman McInnes     /* read row lengths */
77613f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7770752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
77888e32edaSLois Curfman McInnes 
77988e32edaSLois Curfman McInnes     /* create our matrix */
780f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
781f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
782be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7835c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
78488e32edaSLois Curfman McInnes     B = *A;
78588e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
78688e32edaSLois Curfman McInnes     v = a->v;
78788e32edaSLois Curfman McInnes 
78888e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
78913f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
790b0a32e0cSBarry Smith     cols = scols;
7910752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
79287828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
793b0a32e0cSBarry Smith     vals = svals;
7940752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
79588e32edaSLois Curfman McInnes 
79688e32edaSLois Curfman McInnes     /* insert into matrix */
79788e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
79888e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
79988e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
80088e32edaSLois Curfman McInnes     }
801606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
802606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
803606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
80488e32edaSLois Curfman McInnes 
8056d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8066d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80790ace30eSBarry Smith   }
8083a40ed3dSBarry Smith   PetscFunctionReturn(0);
80988e32edaSLois Curfman McInnes }
81088e32edaSLois Curfman McInnes 
811e090d566SSatish Balay #include "petscsys.h"
812932b0c3eSLois Curfman McInnes 
8134a2ae208SSatish Balay #undef __FUNCT__
8144a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8156849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
816289bc588SBarry Smith {
817932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
818dfbe8321SBarry Smith   PetscErrorCode    ierr;
81913f74950SBarry Smith   PetscInt          i,j;
8202dcb1b2aSMatthew Knepley   const char        *name;
82187828ca2SBarry Smith   PetscScalar       *v;
822f3ef73ceSBarry Smith   PetscViewerFormat format;
8235f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
8245f481a85SSatish Balay   PetscTruth        allreal = PETSC_TRUE;
8255f481a85SSatish Balay #endif
826932b0c3eSLois Curfman McInnes 
8273a40ed3dSBarry Smith   PetscFunctionBegin;
828435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
829b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
830456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8313a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
832fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
833b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
834d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
83544cd7ae7SLois Curfman McInnes       v = a->v + i;
83677431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
837d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
838aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
839329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
840a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
841329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
842a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8436831982aSBarry Smith         }
84480cd9d93SLois Curfman McInnes #else
8456831982aSBarry Smith         if (*v) {
846a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8476831982aSBarry Smith         }
84880cd9d93SLois Curfman McInnes #endif
8491b807ce4Svictorle         v += a->lda;
85080cd9d93SLois Curfman McInnes       }
851b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
85280cd9d93SLois Curfman McInnes     }
853b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
8543a40ed3dSBarry Smith   } else {
855b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
856aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
85747989497SBarry Smith     /* determine if matrix has all real values */
85847989497SBarry Smith     v = a->v;
859d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
860ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
86147989497SBarry Smith     }
86247989497SBarry Smith #endif
863fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8643a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
865d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
866d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
867fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
868ffac6cdbSBarry Smith     }
869ffac6cdbSBarry Smith 
870d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
871932b0c3eSLois Curfman McInnes       v = a->v + i;
872d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
873aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
87447989497SBarry Smith         if (allreal) {
875f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
87647989497SBarry Smith         } else {
877f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
87847989497SBarry Smith         }
879289bc588SBarry Smith #else
880f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
881289bc588SBarry Smith #endif
8821b807ce4Svictorle 	v += a->lda;
883289bc588SBarry Smith       }
884b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
885289bc588SBarry Smith     }
886fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
887b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
888ffac6cdbSBarry Smith     }
889b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
890da3a660dSBarry Smith   }
891b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8923a40ed3dSBarry Smith   PetscFunctionReturn(0);
893289bc588SBarry Smith }
894289bc588SBarry Smith 
8954a2ae208SSatish Balay #undef __FUNCT__
8964a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8976849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
898932b0c3eSLois Curfman McInnes {
899932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9006849ba73SBarry Smith   PetscErrorCode    ierr;
90113f74950SBarry Smith   int               fd;
902d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
90387828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
904f3ef73ceSBarry Smith   PetscViewerFormat format;
905932b0c3eSLois Curfman McInnes 
9063a40ed3dSBarry Smith   PetscFunctionBegin;
907b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
90890ace30eSBarry Smith 
909b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
910fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
91190ace30eSBarry Smith     /* store the matrix as a dense matrix */
91213f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
913552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
91490ace30eSBarry Smith     col_lens[1] = m;
91590ace30eSBarry Smith     col_lens[2] = n;
91690ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
9176f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
918606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
91990ace30eSBarry Smith 
92090ace30eSBarry Smith     /* write out matrix, by rows */
92187828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
92290ace30eSBarry Smith     v    = a->v;
92390ace30eSBarry Smith     for (j=0; j<n; j++) {
924578230a0SSatish Balay       for (i=0; i<m; i++) {
925578230a0SSatish Balay         vals[j + i*n] = *v++;
92690ace30eSBarry Smith       }
92790ace30eSBarry Smith     }
9286f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
929606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
93090ace30eSBarry Smith   } else {
93113f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
932552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
933932b0c3eSLois Curfman McInnes     col_lens[1] = m;
934932b0c3eSLois Curfman McInnes     col_lens[2] = n;
935932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
936932b0c3eSLois Curfman McInnes 
937932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
938932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9396f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
940932b0c3eSLois Curfman McInnes 
941932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
942932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
943932b0c3eSLois Curfman McInnes     ict = 0;
944932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
945932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
946932b0c3eSLois Curfman McInnes     }
9476f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
948606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
949932b0c3eSLois Curfman McInnes 
950932b0c3eSLois Curfman McInnes     /* store nonzero values */
95187828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
952932b0c3eSLois Curfman McInnes     ict  = 0;
953932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
954932b0c3eSLois Curfman McInnes       v = a->v + i;
955932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9561b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
957932b0c3eSLois Curfman McInnes       }
958932b0c3eSLois Curfman McInnes     }
9596f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
960606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
96190ace30eSBarry Smith   }
9623a40ed3dSBarry Smith   PetscFunctionReturn(0);
963932b0c3eSLois Curfman McInnes }
964932b0c3eSLois Curfman McInnes 
9654a2ae208SSatish Balay #undef __FUNCT__
9664a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
967dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
968f1af5d2fSBarry Smith {
969f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
970f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9716849ba73SBarry Smith   PetscErrorCode    ierr;
972d0f46423SBarry Smith   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
97387828ca2SBarry Smith   PetscScalar       *v = a->v;
974b0a32e0cSBarry Smith   PetscViewer       viewer;
975b0a32e0cSBarry Smith   PetscDraw         popup;
976329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
977f3ef73ceSBarry Smith   PetscViewerFormat format;
978f1af5d2fSBarry Smith 
979f1af5d2fSBarry Smith   PetscFunctionBegin;
980f1af5d2fSBarry Smith 
981f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
982b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
983b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
984f1af5d2fSBarry Smith 
985f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
986fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
987f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
988b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
989f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
990f1af5d2fSBarry Smith       x_l = j;
991f1af5d2fSBarry Smith       x_r = x_l + 1.0;
992f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
993f1af5d2fSBarry Smith         y_l = m - i - 1.0;
994f1af5d2fSBarry Smith         y_r = y_l + 1.0;
995f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
996329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
997b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
998329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
999b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1000f1af5d2fSBarry Smith         } else {
1001f1af5d2fSBarry Smith           continue;
1002f1af5d2fSBarry Smith         }
1003f1af5d2fSBarry Smith #else
1004f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
1005b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1006f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1007b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1008f1af5d2fSBarry Smith         } else {
1009f1af5d2fSBarry Smith           continue;
1010f1af5d2fSBarry Smith         }
1011f1af5d2fSBarry Smith #endif
1012b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1013f1af5d2fSBarry Smith       }
1014f1af5d2fSBarry Smith     }
1015f1af5d2fSBarry Smith   } else {
1016f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1017f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1018f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1019f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1020f1af5d2fSBarry Smith     }
1021b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1022b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1023b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1024f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1025f1af5d2fSBarry Smith       x_l = j;
1026f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1027f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1028f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1029f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1030b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1031b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1032f1af5d2fSBarry Smith       }
1033f1af5d2fSBarry Smith     }
1034f1af5d2fSBarry Smith   }
1035f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1036f1af5d2fSBarry Smith }
1037f1af5d2fSBarry Smith 
10384a2ae208SSatish Balay #undef __FUNCT__
10394a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1040dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1041f1af5d2fSBarry Smith {
1042b0a32e0cSBarry Smith   PetscDraw      draw;
1043f1af5d2fSBarry Smith   PetscTruth     isnull;
1044329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1045dfbe8321SBarry Smith   PetscErrorCode ierr;
1046f1af5d2fSBarry Smith 
1047f1af5d2fSBarry Smith   PetscFunctionBegin;
1048b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1049b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1050abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1051f1af5d2fSBarry Smith 
1052f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1053d0f46423SBarry Smith   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1054f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1055b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1056b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1057f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1058f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1059f1af5d2fSBarry Smith }
1060f1af5d2fSBarry Smith 
10614a2ae208SSatish Balay #undef __FUNCT__
10624a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1063dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1064932b0c3eSLois Curfman McInnes {
1065dfbe8321SBarry Smith   PetscErrorCode ierr;
10666805f65bSBarry Smith   PetscTruth     iascii,isbinary,isdraw;
1067932b0c3eSLois Curfman McInnes 
10683a40ed3dSBarry Smith   PetscFunctionBegin;
106932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1070fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1071fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10720f5bd95cSBarry Smith 
1073c45a1595SBarry Smith   if (iascii) {
1074c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10750f5bd95cSBarry Smith   } else if (isbinary) {
10763a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1077f1af5d2fSBarry Smith   } else if (isdraw) {
1078f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10795cd90555SBarry Smith   } else {
1080958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1081932b0c3eSLois Curfman McInnes   }
10823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1083932b0c3eSLois Curfman McInnes }
1084289bc588SBarry Smith 
10854a2ae208SSatish Balay #undef __FUNCT__
10864a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1087dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1088289bc588SBarry Smith {
1089ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1090dfbe8321SBarry Smith   PetscErrorCode ierr;
109190f02eecSBarry Smith 
10923a40ed3dSBarry Smith   PetscFunctionBegin;
1093aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1094d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1095a5a9c739SBarry Smith #endif
109605b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
10976857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1098606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1099dbd8c25aSHong Zhang 
1100dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1101901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
11024ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11034ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11044ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11053a40ed3dSBarry Smith   PetscFunctionReturn(0);
1106289bc588SBarry Smith }
1107289bc588SBarry Smith 
11084a2ae208SSatish Balay #undef __FUNCT__
11094a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1110fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1111289bc588SBarry Smith {
1112c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11136849ba73SBarry Smith   PetscErrorCode ierr;
111413f74950SBarry Smith   PetscInt       k,j,m,n,M;
111587828ca2SBarry Smith   PetscScalar    *v,tmp;
111648b35521SBarry Smith 
11173a40ed3dSBarry Smith   PetscFunctionBegin;
1118d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1119e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1120a5ce6ee0Svictorle     if (m != n) {
1121634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1122d64ed03dSBarry Smith     } else {
1123d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1124289bc588SBarry Smith         for (k=0; k<j; k++) {
11251b807ce4Svictorle           tmp = v[j + k*M];
11261b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11271b807ce4Svictorle           v[k + j*M] = tmp;
1128289bc588SBarry Smith         }
1129289bc588SBarry Smith       }
1130d64ed03dSBarry Smith     }
11313a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1132d3e5ee88SLois Curfman McInnes     Mat          tmat;
1133ec8511deSBarry Smith     Mat_SeqDense *tmatd;
113487828ca2SBarry Smith     PetscScalar  *v2;
1135ea709b57SSatish Balay 
1136fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11377adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1138d0f46423SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
11397adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11405c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1141fc4dec0aSBarry Smith     } else {
1142fc4dec0aSBarry Smith       tmat = *matout;
1143fc4dec0aSBarry Smith     }
1144ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11450de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1146d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11471b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1148d3e5ee88SLois Curfman McInnes     }
11496d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11506d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1151d3e5ee88SLois Curfman McInnes     *matout = tmat;
115248b35521SBarry Smith   }
11533a40ed3dSBarry Smith   PetscFunctionReturn(0);
1154289bc588SBarry Smith }
1155289bc588SBarry Smith 
11564a2ae208SSatish Balay #undef __FUNCT__
11574a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1158dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1159289bc588SBarry Smith {
1160c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1161c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
116213f74950SBarry Smith   PetscInt     i,j;
116387828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11649ea5d5aeSSatish Balay 
11653a40ed3dSBarry Smith   PetscFunctionBegin;
1166d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1167d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1168d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
11691b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1170d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
11713a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11721b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11731b807ce4Svictorle     }
1174289bc588SBarry Smith   }
117577c4ece6SBarry Smith   *flg = PETSC_TRUE;
11763a40ed3dSBarry Smith   PetscFunctionReturn(0);
1177289bc588SBarry Smith }
1178289bc588SBarry Smith 
11794a2ae208SSatish Balay #undef __FUNCT__
11804a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1181dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1182289bc588SBarry Smith {
1183c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1184dfbe8321SBarry Smith   PetscErrorCode ierr;
118513f74950SBarry Smith   PetscInt       i,n,len;
118687828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
118744cd7ae7SLois Curfman McInnes 
11883a40ed3dSBarry Smith   PetscFunctionBegin;
11892dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11907a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11911ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1192d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1193d0f46423SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
119444cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11951b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1196289bc588SBarry Smith   }
11971ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1199289bc588SBarry Smith }
1200289bc588SBarry Smith 
12014a2ae208SSatish Balay #undef __FUNCT__
12024a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1203dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1204289bc588SBarry Smith {
1205c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
120687828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1207dfbe8321SBarry Smith   PetscErrorCode ierr;
1208d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
120955659b69SBarry Smith 
12103a40ed3dSBarry Smith   PetscFunctionBegin;
121128988994SBarry Smith   if (ll) {
12127a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12131ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1214d0f46423SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1215da3a660dSBarry Smith     for (i=0; i<m; i++) {
1216da3a660dSBarry Smith       x = l[i];
1217da3a660dSBarry Smith       v = mat->v + i;
1218da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1219da3a660dSBarry Smith     }
12201ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1221efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1222da3a660dSBarry Smith   }
122328988994SBarry Smith   if (rr) {
12247a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12251ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1226d0f46423SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1227da3a660dSBarry Smith     for (i=0; i<n; i++) {
1228da3a660dSBarry Smith       x = r[i];
1229da3a660dSBarry Smith       v = mat->v + i*m;
1230da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1231da3a660dSBarry Smith     }
12321ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1233efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1234da3a660dSBarry Smith   }
12353a40ed3dSBarry Smith   PetscFunctionReturn(0);
1236289bc588SBarry Smith }
1237289bc588SBarry Smith 
12384a2ae208SSatish Balay #undef __FUNCT__
12394a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1240dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1241289bc588SBarry Smith {
1242c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
124387828ca2SBarry Smith   PetscScalar  *v = mat->v;
1244329f5518SBarry Smith   PetscReal    sum = 0.0;
1245d0f46423SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1246efee365bSSatish Balay   PetscErrorCode ierr;
124755659b69SBarry Smith 
12483a40ed3dSBarry Smith   PetscFunctionBegin;
1249289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1250a5ce6ee0Svictorle     if (lda>m) {
1251d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1252a5ce6ee0Svictorle 	v = mat->v+j*lda;
1253a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1254a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1255a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1256a5ce6ee0Svictorle #else
1257a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1258a5ce6ee0Svictorle #endif
1259a5ce6ee0Svictorle 	}
1260a5ce6ee0Svictorle       }
1261a5ce6ee0Svictorle     } else {
1262d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1263aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1264329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1265289bc588SBarry Smith #else
1266289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1267289bc588SBarry Smith #endif
1268289bc588SBarry Smith       }
1269a5ce6ee0Svictorle     }
1270064f8208SBarry Smith     *nrm = sqrt(sum);
1271d0f46423SBarry Smith     ierr = PetscLogFlops(2*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12723a40ed3dSBarry Smith   } else if (type == NORM_1) {
1273064f8208SBarry Smith     *nrm = 0.0;
1274d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
12751b807ce4Svictorle       v = mat->v + j*mat->lda;
1276289bc588SBarry Smith       sum = 0.0;
1277d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
127833a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1279289bc588SBarry Smith       }
1280064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1281289bc588SBarry Smith     }
1282d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12833a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1284064f8208SBarry Smith     *nrm = 0.0;
1285d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1286289bc588SBarry Smith       v = mat->v + j;
1287289bc588SBarry Smith       sum = 0.0;
1288d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
12891b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1290289bc588SBarry Smith       }
1291064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1292289bc588SBarry Smith     }
1293d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12943a40ed3dSBarry Smith   } else {
129529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1296289bc588SBarry Smith   }
12973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1298289bc588SBarry Smith }
1299289bc588SBarry Smith 
13004a2ae208SSatish Balay #undef __FUNCT__
13014a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
13024e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg)
1303289bc588SBarry Smith {
1304c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
130563ba0a88SBarry Smith   PetscErrorCode ierr;
130667e560aaSBarry Smith 
13073a40ed3dSBarry Smith   PetscFunctionBegin;
1308b5a2b587SKris Buschelman   switch (op) {
1309b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13104e0d8c25SBarry Smith     aij->roworiented = flg;
1311b5a2b587SKris Buschelman     break;
1312512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1313b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13143971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13154e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1316b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1317b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
131877e54ba9SKris Buschelman   case MAT_SYMMETRIC:
131977e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13209a4540c5SBarry Smith   case MAT_HERMITIAN:
13219a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1322600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1323290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
132477e54ba9SKris Buschelman     break;
1325b5a2b587SKris Buschelman   default:
1326600fe468SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13273a40ed3dSBarry Smith   }
13283a40ed3dSBarry Smith   PetscFunctionReturn(0);
1329289bc588SBarry Smith }
1330289bc588SBarry Smith 
13314a2ae208SSatish Balay #undef __FUNCT__
13324a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1333dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13346f0a148fSBarry Smith {
1335ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13366849ba73SBarry Smith   PetscErrorCode ierr;
1337d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13383a40ed3dSBarry Smith 
13393a40ed3dSBarry Smith   PetscFunctionBegin;
1340a5ce6ee0Svictorle   if (lda>m) {
1341d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1342a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1343a5ce6ee0Svictorle     }
1344a5ce6ee0Svictorle   } else {
1345d0f46423SBarry Smith     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1346a5ce6ee0Svictorle   }
13473a40ed3dSBarry Smith   PetscFunctionReturn(0);
13486f0a148fSBarry Smith }
13496f0a148fSBarry Smith 
13504a2ae208SSatish Balay #undef __FUNCT__
13514a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1352f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13536f0a148fSBarry Smith {
1354ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1355d0f46423SBarry Smith   PetscInt       n = A->cmap->n,i,j;
135687828ca2SBarry Smith   PetscScalar    *slot;
135755659b69SBarry Smith 
13583a40ed3dSBarry Smith   PetscFunctionBegin;
13596f0a148fSBarry Smith   for (i=0; i<N; i++) {
13606f0a148fSBarry Smith     slot = l->v + rows[i];
13616f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13626f0a148fSBarry Smith   }
1363f4df32b1SMatthew Knepley   if (diag != 0.0) {
13646f0a148fSBarry Smith     for (i=0; i<N; i++) {
13656f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1366f4df32b1SMatthew Knepley       *slot = diag;
13676f0a148fSBarry Smith     }
13686f0a148fSBarry Smith   }
13693a40ed3dSBarry Smith   PetscFunctionReturn(0);
13706f0a148fSBarry Smith }
1371557bce09SLois Curfman McInnes 
13724a2ae208SSatish Balay #undef __FUNCT__
13734a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1374dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
137564e87e97SBarry Smith {
1376c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13773a40ed3dSBarry Smith 
13783a40ed3dSBarry Smith   PetscFunctionBegin;
1379d0f46423SBarry Smith   if (mat->lda != A->rmap->n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
138064e87e97SBarry Smith   *array = mat->v;
13813a40ed3dSBarry Smith   PetscFunctionReturn(0);
138264e87e97SBarry Smith }
13830754003eSLois Curfman McInnes 
13844a2ae208SSatish Balay #undef __FUNCT__
13854a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1386dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1387ff14e315SSatish Balay {
13883a40ed3dSBarry Smith   PetscFunctionBegin;
138909b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13903a40ed3dSBarry Smith   PetscFunctionReturn(0);
1391ff14e315SSatish Balay }
13920754003eSLois Curfman McInnes 
13934a2ae208SSatish Balay #undef __FUNCT__
13944a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
139513f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13960754003eSLois Curfman McInnes {
1397c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13986849ba73SBarry Smith   PetscErrorCode ierr;
13995d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
14005d0c19d7SBarry Smith   const PetscInt *irow,*icol;
140187828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
14020754003eSLois Curfman McInnes   Mat            newmat;
14030754003eSLois Curfman McInnes 
14043a40ed3dSBarry Smith   PetscFunctionBegin;
140578b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
140678b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1407e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1408e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14090754003eSLois Curfman McInnes 
1410182d2002SSatish Balay   /* Check submatrixcall */
1411182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
141213f74950SBarry Smith     PetscInt n_cols,n_rows;
1413182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
141421a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
141521a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
141621a2c019SBarry Smith       ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr);
141721a2c019SBarry Smith     }
1418182d2002SSatish Balay     newmat = *B;
1419182d2002SSatish Balay   } else {
14200754003eSLois Curfman McInnes     /* Create and fill new matrix */
14217adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1422f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14237adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14245c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1425182d2002SSatish Balay   }
1426182d2002SSatish Balay 
1427182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1428182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1429182d2002SSatish Balay 
1430182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14316de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1432182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1433182d2002SSatish Balay       *bv++ = av[irow[j]];
14340754003eSLois Curfman McInnes     }
14350754003eSLois Curfman McInnes   }
1436182d2002SSatish Balay 
1437182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14386d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14396d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14400754003eSLois Curfman McInnes 
14410754003eSLois Curfman McInnes   /* Free work space */
144278b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
144378b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1444182d2002SSatish Balay   *B = newmat;
14453a40ed3dSBarry Smith   PetscFunctionReturn(0);
14460754003eSLois Curfman McInnes }
14470754003eSLois Curfman McInnes 
14484a2ae208SSatish Balay #undef __FUNCT__
14494a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
145013f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1451905e6a2fSBarry Smith {
14526849ba73SBarry Smith   PetscErrorCode ierr;
145313f74950SBarry Smith   PetscInt       i;
1454905e6a2fSBarry Smith 
14553a40ed3dSBarry Smith   PetscFunctionBegin;
1456905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1457b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1458905e6a2fSBarry Smith   }
1459905e6a2fSBarry Smith 
1460905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14616a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1462905e6a2fSBarry Smith   }
14633a40ed3dSBarry Smith   PetscFunctionReturn(0);
1464905e6a2fSBarry Smith }
1465905e6a2fSBarry Smith 
14664a2ae208SSatish Balay #undef __FUNCT__
1467c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1468c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1469c0aa2d19SHong Zhang {
1470c0aa2d19SHong Zhang   PetscFunctionBegin;
1471c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1472c0aa2d19SHong Zhang }
1473c0aa2d19SHong Zhang 
1474c0aa2d19SHong Zhang #undef __FUNCT__
1475c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1476c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1477c0aa2d19SHong Zhang {
1478c0aa2d19SHong Zhang   PetscFunctionBegin;
1479c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1480c0aa2d19SHong Zhang }
1481c0aa2d19SHong Zhang 
1482c0aa2d19SHong Zhang #undef __FUNCT__
14834a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1484dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14854b0e389bSBarry Smith {
14864b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14876849ba73SBarry Smith   PetscErrorCode ierr;
1488d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
14893a40ed3dSBarry Smith 
14903a40ed3dSBarry Smith   PetscFunctionBegin;
149133f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
149233f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1493cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14943a40ed3dSBarry Smith     PetscFunctionReturn(0);
14953a40ed3dSBarry Smith   }
1496d0f46423SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1497a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14980dbb7854Svictorle     for (j=0; j<n; j++) {
1499a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1500a5ce6ee0Svictorle     }
1501a5ce6ee0Svictorle   } else {
1502d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1503a5ce6ee0Svictorle   }
1504273d9f13SBarry Smith   PetscFunctionReturn(0);
1505273d9f13SBarry Smith }
1506273d9f13SBarry Smith 
15074a2ae208SSatish Balay #undef __FUNCT__
15084a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1509dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1510273d9f13SBarry Smith {
1511dfbe8321SBarry Smith   PetscErrorCode ierr;
1512273d9f13SBarry Smith 
1513273d9f13SBarry Smith   PetscFunctionBegin;
1514273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15153a40ed3dSBarry Smith   PetscFunctionReturn(0);
15164b0e389bSBarry Smith }
15174b0e389bSBarry Smith 
1518284134d9SBarry Smith #undef __FUNCT__
1519284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1520284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1521284134d9SBarry Smith {
1522284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
152321a2c019SBarry Smith   PetscErrorCode ierr;
1524284134d9SBarry Smith   PetscFunctionBegin;
152521a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1526284134d9SBarry Smith   m = PetscMax(m,M);
1527284134d9SBarry Smith   n = PetscMax(n,N);
152821a2c019SBarry 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);
1529284134d9SBarry 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);
1530d0f46423SBarry Smith   A->rmap->n = A->rmap->n = m;
1531d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
153221a2c019SBarry Smith   if (a->changelda) a->lda = m;
153321a2c019SBarry Smith   ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
1534284134d9SBarry Smith   PetscFunctionReturn(0);
1535284134d9SBarry Smith }
1536170fe5c8SBarry Smith 
1537284134d9SBarry Smith 
1538a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1539a9fe9ddaSSatish Balay #undef __FUNCT__
1540a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1541a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1542a9fe9ddaSSatish Balay {
1543a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1544a9fe9ddaSSatish Balay 
1545a9fe9ddaSSatish Balay   PetscFunctionBegin;
1546a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1547a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1548a9fe9ddaSSatish Balay   }
1549a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1550a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1551a9fe9ddaSSatish Balay }
1552a9fe9ddaSSatish Balay 
1553a9fe9ddaSSatish Balay #undef __FUNCT__
1554a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1555a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1556a9fe9ddaSSatish Balay {
1557ee16a9a1SHong Zhang   PetscErrorCode ierr;
1558d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1559ee16a9a1SHong Zhang   Mat            Cmat;
1560a9fe9ddaSSatish Balay 
1561ee16a9a1SHong Zhang   PetscFunctionBegin;
1562d0f46423SBarry 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);
1563ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1564ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1565ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1566ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1567ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1568ee16a9a1SHong Zhang   *C = Cmat;
1569ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1570ee16a9a1SHong Zhang }
1571a9fe9ddaSSatish Balay 
157298a3b096SSatish Balay #undef __FUNCT__
1573a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1574a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1575a9fe9ddaSSatish Balay {
1576a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1577a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1578a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
15790805154bSBarry Smith   PetscBLASInt   m,n,k;
1580a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1581a9fe9ddaSSatish Balay 
1582a9fe9ddaSSatish Balay   PetscFunctionBegin;
1583d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1584d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1585d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1586a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1587a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1588a9fe9ddaSSatish Balay }
1589a9fe9ddaSSatish Balay 
1590a9fe9ddaSSatish Balay #undef __FUNCT__
1591a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1592a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1593a9fe9ddaSSatish Balay {
1594a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1595a9fe9ddaSSatish Balay 
1596a9fe9ddaSSatish Balay   PetscFunctionBegin;
1597a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1598a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1599a9fe9ddaSSatish Balay   }
1600a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1601a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1602a9fe9ddaSSatish Balay }
1603a9fe9ddaSSatish Balay 
1604a9fe9ddaSSatish Balay #undef __FUNCT__
1605a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1606a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1607a9fe9ddaSSatish Balay {
1608ee16a9a1SHong Zhang   PetscErrorCode ierr;
1609d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1610ee16a9a1SHong Zhang   Mat            Cmat;
1611a9fe9ddaSSatish Balay 
1612ee16a9a1SHong Zhang   PetscFunctionBegin;
1613d0f46423SBarry 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);
1614ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1615ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1616ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1617ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1618ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1619ee16a9a1SHong Zhang   *C = Cmat;
1620ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1621ee16a9a1SHong Zhang }
1622a9fe9ddaSSatish Balay 
1623a9fe9ddaSSatish Balay #undef __FUNCT__
1624a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1625a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1626a9fe9ddaSSatish Balay {
1627a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1628a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1629a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16300805154bSBarry Smith   PetscBLASInt   m,n,k;
1631a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1632a9fe9ddaSSatish Balay 
1633a9fe9ddaSSatish Balay   PetscFunctionBegin;
1634d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1635d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1636d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
16372fbe02b9SBarry Smith   /*
16382fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
16392fbe02b9SBarry Smith   */
1640a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1641a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1642a9fe9ddaSSatish Balay }
1643985db425SBarry Smith 
1644985db425SBarry Smith #undef __FUNCT__
1645985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1646985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1647985db425SBarry Smith {
1648985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1649985db425SBarry Smith   PetscErrorCode ierr;
1650d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1651985db425SBarry Smith   PetscScalar    *x;
1652985db425SBarry Smith   MatScalar      *aa = a->v;
1653985db425SBarry Smith 
1654985db425SBarry Smith   PetscFunctionBegin;
1655985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1656985db425SBarry Smith 
1657985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1658985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1659985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1660d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1661985db425SBarry Smith   for (i=0; i<m; i++) {
1662985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1663985db425SBarry Smith     for (j=1; j<n; j++){
1664985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1665985db425SBarry Smith     }
1666985db425SBarry Smith   }
1667985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1668985db425SBarry Smith   PetscFunctionReturn(0);
1669985db425SBarry Smith }
1670985db425SBarry Smith 
1671985db425SBarry Smith #undef __FUNCT__
1672985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1673985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1674985db425SBarry Smith {
1675985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1676985db425SBarry Smith   PetscErrorCode ierr;
1677d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1678985db425SBarry Smith   PetscScalar    *x;
1679985db425SBarry Smith   PetscReal      atmp;
1680985db425SBarry Smith   MatScalar      *aa = a->v;
1681985db425SBarry Smith 
1682985db425SBarry Smith   PetscFunctionBegin;
1683985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1684985db425SBarry Smith 
1685985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1686985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1687985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1688d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1689985db425SBarry Smith   for (i=0; i<m; i++) {
16909189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1691985db425SBarry Smith     for (j=1; j<n; j++){
1692985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1693985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1694985db425SBarry Smith     }
1695985db425SBarry Smith   }
1696985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1697985db425SBarry Smith   PetscFunctionReturn(0);
1698985db425SBarry Smith }
1699985db425SBarry Smith 
1700985db425SBarry Smith #undef __FUNCT__
1701985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1702985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1703985db425SBarry Smith {
1704985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1705985db425SBarry Smith   PetscErrorCode ierr;
1706d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1707985db425SBarry Smith   PetscScalar    *x;
1708985db425SBarry Smith   MatScalar      *aa = a->v;
1709985db425SBarry Smith 
1710985db425SBarry Smith   PetscFunctionBegin;
1711985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1712985db425SBarry Smith 
1713985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1714985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1715985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1716d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1717985db425SBarry Smith   for (i=0; i<m; i++) {
1718985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1719985db425SBarry Smith     for (j=1; j<n; j++){
1720985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1721985db425SBarry Smith     }
1722985db425SBarry Smith   }
1723985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1724985db425SBarry Smith   PetscFunctionReturn(0);
1725985db425SBarry Smith }
1726985db425SBarry Smith 
17278d0534beSBarry Smith #undef __FUNCT__
17288d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
17298d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
17308d0534beSBarry Smith {
17318d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17328d0534beSBarry Smith   PetscErrorCode ierr;
17338d0534beSBarry Smith   PetscScalar    *x;
17348d0534beSBarry Smith 
17358d0534beSBarry Smith   PetscFunctionBegin;
17368d0534beSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
17378d0534beSBarry Smith 
17388d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1739d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
17408d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
17418d0534beSBarry Smith   PetscFunctionReturn(0);
17428d0534beSBarry Smith }
17438d0534beSBarry Smith 
1744289bc588SBarry Smith /* -------------------------------------------------------------------*/
1745a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1746905e6a2fSBarry Smith        MatGetRow_SeqDense,
1747905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1748905e6a2fSBarry Smith        MatMult_SeqDense,
174997304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
17507c922b88SBarry Smith        MatMultTranspose_SeqDense,
17517c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1752db4efbfdSBarry Smith        0,
1753db4efbfdSBarry Smith        0,
1754db4efbfdSBarry Smith        0,
1755db4efbfdSBarry Smith /*10*/ 0,
1756905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1757905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1758ec8511deSBarry Smith        MatRelax_SeqDense,
1759ec8511deSBarry Smith        MatTranspose_SeqDense,
176097304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1761905e6a2fSBarry Smith        MatEqual_SeqDense,
1762905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1763905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1764905e6a2fSBarry Smith        MatNorm_SeqDense,
1765c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1766c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1767905e6a2fSBarry Smith        0,
1768905e6a2fSBarry Smith        MatSetOption_SeqDense,
1769905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
177097304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1771db4efbfdSBarry Smith        0,
1772db4efbfdSBarry Smith        0,
1773db4efbfdSBarry Smith        0,
1774db4efbfdSBarry Smith        0,
177597304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1776273d9f13SBarry Smith        0,
1777905e6a2fSBarry Smith        0,
1778905e6a2fSBarry Smith        MatGetArray_SeqDense,
1779905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
178097304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1781a5ae1ecdSBarry Smith        0,
1782a5ae1ecdSBarry Smith        0,
1783a5ae1ecdSBarry Smith        0,
1784a5ae1ecdSBarry Smith        0,
178597304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1786a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1787a5ae1ecdSBarry Smith        0,
17884b0e389bSBarry Smith        MatGetValues_SeqDense,
1789a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1790985db425SBarry Smith /*45*/ MatGetRowMax_SeqDense,
1791a5ae1ecdSBarry Smith        MatScale_SeqDense,
1792a5ae1ecdSBarry Smith        0,
1793a5ae1ecdSBarry Smith        0,
1794a5ae1ecdSBarry Smith        0,
1795521d7252SBarry Smith /*50*/ 0,
1796a5ae1ecdSBarry Smith        0,
1797a5ae1ecdSBarry Smith        0,
1798a5ae1ecdSBarry Smith        0,
1799a5ae1ecdSBarry Smith        0,
180097304618SKris Buschelman /*55*/ 0,
1801a5ae1ecdSBarry Smith        0,
1802a5ae1ecdSBarry Smith        0,
1803a5ae1ecdSBarry Smith        0,
1804a5ae1ecdSBarry Smith        0,
180597304618SKris Buschelman /*60*/ 0,
1806e03a110bSBarry Smith        MatDestroy_SeqDense,
1807e03a110bSBarry Smith        MatView_SeqDense,
1808357abbc8SBarry Smith        0,
180997304618SKris Buschelman        0,
181097304618SKris Buschelman /*65*/ 0,
181197304618SKris Buschelman        0,
181297304618SKris Buschelman        0,
181397304618SKris Buschelman        0,
181497304618SKris Buschelman        0,
1815985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqDense,
181697304618SKris Buschelman        0,
181797304618SKris Buschelman        0,
181897304618SKris Buschelman        0,
181997304618SKris Buschelman        0,
182097304618SKris Buschelman /*75*/ 0,
182197304618SKris Buschelman        0,
182297304618SKris Buschelman        0,
182397304618SKris Buschelman        0,
182497304618SKris Buschelman        0,
182597304618SKris Buschelman /*80*/ 0,
182697304618SKris Buschelman        0,
182797304618SKris Buschelman        0,
182897304618SKris Buschelman        0,
1829865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1830865e5f61SKris Buschelman        0,
18311cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1832865e5f61SKris Buschelman        0,
1833865e5f61SKris Buschelman        0,
1834865e5f61SKris Buschelman        0,
1835a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense,
1836a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1837a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1838865e5f61SKris Buschelman        0,
1839865e5f61SKris Buschelman        0,
1840865e5f61SKris Buschelman /*95*/ 0,
1841a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1842a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1843a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1844284134d9SBarry Smith        0,
1845284134d9SBarry Smith /*100*/0,
1846284134d9SBarry Smith        0,
1847284134d9SBarry Smith        0,
1848284134d9SBarry Smith        0,
1849985db425SBarry Smith        MatSetSizes_SeqDense,
1850985db425SBarry Smith        0,
1851985db425SBarry Smith        0,
1852985db425SBarry Smith        0,
1853985db425SBarry Smith        0,
1854985db425SBarry Smith        0,
1855985db425SBarry Smith /*110*/0,
1856985db425SBarry Smith        0,
18578d0534beSBarry Smith        MatGetRowMin_SeqDense,
18588d0534beSBarry Smith        MatGetColumnVector_SeqDense
1859985db425SBarry Smith };
186090ace30eSBarry Smith 
18614a2ae208SSatish Balay #undef __FUNCT__
18624a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18634b828684SBarry Smith /*@C
1864fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1865d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1866d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1867289bc588SBarry Smith 
1868db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1869db81eaa0SLois Curfman McInnes 
187020563c6bSBarry Smith    Input Parameters:
1871db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
18720c775827SLois Curfman McInnes .  m - number of rows
187318f449edSLois Curfman McInnes .  n - number of columns
1874db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1875dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
187620563c6bSBarry Smith 
187720563c6bSBarry Smith    Output Parameter:
187844cd7ae7SLois Curfman McInnes .  A - the matrix
187920563c6bSBarry Smith 
1880b259b22eSLois Curfman McInnes    Notes:
188118f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
188218f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1883b4fd4287SBarry Smith    set data=PETSC_NULL.
188418f449edSLois Curfman McInnes 
1885027ccd11SLois Curfman McInnes    Level: intermediate
1886027ccd11SLois Curfman McInnes 
1887dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1888d65003e9SLois Curfman McInnes 
1889db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
189020563c6bSBarry Smith @*/
1891be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1892289bc588SBarry Smith {
1893dfbe8321SBarry Smith   PetscErrorCode ierr;
18943b2fbd54SBarry Smith 
18953a40ed3dSBarry Smith   PetscFunctionBegin;
1896f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1897f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1898273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1899273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1900273d9f13SBarry Smith   PetscFunctionReturn(0);
1901273d9f13SBarry Smith }
1902273d9f13SBarry Smith 
19034a2ae208SSatish Balay #undef __FUNCT__
1904afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1905273d9f13SBarry Smith /*@C
1906273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1907273d9f13SBarry Smith 
1908273d9f13SBarry Smith    Collective on MPI_Comm
1909273d9f13SBarry Smith 
1910273d9f13SBarry Smith    Input Parameters:
1911273d9f13SBarry Smith +  A - the matrix
1912273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1913273d9f13SBarry Smith 
1914273d9f13SBarry Smith    Notes:
1915273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1916273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1917284134d9SBarry Smith    need not call this routine.
1918273d9f13SBarry Smith 
1919273d9f13SBarry Smith    Level: intermediate
1920273d9f13SBarry Smith 
1921273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1922273d9f13SBarry Smith 
1923273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1924273d9f13SBarry Smith @*/
1925be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1926273d9f13SBarry Smith {
1927dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1928a23d5eceSKris Buschelman 
1929a23d5eceSKris Buschelman   PetscFunctionBegin;
1930a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1931a23d5eceSKris Buschelman   if (f) {
1932a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1933a23d5eceSKris Buschelman   }
1934a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1935a23d5eceSKris Buschelman }
1936a23d5eceSKris Buschelman 
1937a23d5eceSKris Buschelman EXTERN_C_BEGIN
1938a23d5eceSKris Buschelman #undef __FUNCT__
1939afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1940be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1941a23d5eceSKris Buschelman {
1942273d9f13SBarry Smith   Mat_SeqDense   *b;
1943dfbe8321SBarry Smith   PetscErrorCode ierr;
1944273d9f13SBarry Smith 
1945273d9f13SBarry Smith   PetscFunctionBegin;
1946273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1947273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1948d0f46423SBarry Smith   if (b->lda <= 0) b->lda = B->rmap->n;
19499e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
19509e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
19515afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1952284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1953284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
19549e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
1955273d9f13SBarry Smith   } else { /* user-allocated storage */
19569e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1957273d9f13SBarry Smith     b->v          = data;
1958273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1959273d9f13SBarry Smith   }
1960273d9f13SBarry Smith   PetscFunctionReturn(0);
1961273d9f13SBarry Smith }
1962a23d5eceSKris Buschelman EXTERN_C_END
1963273d9f13SBarry Smith 
19641b807ce4Svictorle #undef __FUNCT__
19651b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
19661b807ce4Svictorle /*@C
19671b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
19681b807ce4Svictorle 
19691b807ce4Svictorle   Input parameter:
19701b807ce4Svictorle + A - the matrix
19711b807ce4Svictorle - lda - the leading dimension
19721b807ce4Svictorle 
19731b807ce4Svictorle   Notes:
19741b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
19751b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
19761b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
19771b807ce4Svictorle 
19781b807ce4Svictorle   Level: intermediate
19791b807ce4Svictorle 
19801b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
19811b807ce4Svictorle 
1982284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
19831b807ce4Svictorle @*/
1984be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
19851b807ce4Svictorle {
19861b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
198721a2c019SBarry Smith 
19881b807ce4Svictorle   PetscFunctionBegin;
1989d0f46423SBarry Smith   if (lda < B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
19901b807ce4Svictorle   b->lda       = lda;
199121a2c019SBarry Smith   b->changelda = PETSC_FALSE;
199221a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
19931b807ce4Svictorle   PetscFunctionReturn(0);
19941b807ce4Svictorle }
19951b807ce4Svictorle 
19960bad9183SKris Buschelman /*MC
1997fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
19980bad9183SKris Buschelman 
19990bad9183SKris Buschelman    Options Database Keys:
20000bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
20010bad9183SKris Buschelman 
20020bad9183SKris Buschelman   Level: beginner
20030bad9183SKris Buschelman 
200489665df3SBarry Smith .seealso: MatCreateSeqDense()
200589665df3SBarry Smith 
20060bad9183SKris Buschelman M*/
20070bad9183SKris Buschelman 
2008273d9f13SBarry Smith EXTERN_C_BEGIN
20094a2ae208SSatish Balay #undef __FUNCT__
20104a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
2011be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2012273d9f13SBarry Smith {
2013273d9f13SBarry Smith   Mat_SeqDense   *b;
2014dfbe8321SBarry Smith   PetscErrorCode ierr;
20157c334f02SBarry Smith   PetscMPIInt    size;
2016273d9f13SBarry Smith 
2017273d9f13SBarry Smith   PetscFunctionBegin;
20187adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
201929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
202055659b69SBarry Smith 
2021d0f46423SBarry Smith   B->rmap->bs = B->cmap->bs = 1;
2022d0f46423SBarry Smith   ierr = PetscMapSetUp(B->rmap);CHKERRQ(ierr);
2023d0f46423SBarry Smith   ierr = PetscMapSetUp(B->cmap);CHKERRQ(ierr);
2024273d9f13SBarry Smith 
202538f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2026549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
202790f02eecSBarry Smith   B->mapping      = 0;
202844cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
202918f449edSLois Curfman McInnes 
2030a5ae1ecdSBarry Smith 
203144cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2032273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2033273d9f13SBarry Smith   b->v            = 0;
2034d0f46423SBarry Smith   b->lda          = B->rmap->n;
203521a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
2036d0f46423SBarry Smith   b->Mmax         = B->rmap->n;
2037d0f46423SBarry Smith   b->Nmax         = B->cmap->n;
20384e220ebcSLois Curfman McInnes 
2039b24902e0SBarry Smith 
2040b24902e0SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqdense_petsc_C",
2041b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2042b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2043a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2044a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2045a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
20464ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
20474ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
20484ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
20494ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
20504ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
20514ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
20524ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
20534ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
20544ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
205517667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
20563a40ed3dSBarry Smith   PetscFunctionReturn(0);
2057289bc588SBarry Smith }
2058c0aa2d19SHong Zhang 
2059c0aa2d19SHong Zhang 
2060273d9f13SBarry Smith EXTERN_C_END
2061