xref: /petsc/src/mat/impls/dense/seq/dense.c (revision bda8bf91218abc8014db6f91a57cb182fae6533f)
1be1d678aSKris Buschelman 
267e560aaSBarry Smith /*
367e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
467e560aaSBarry Smith */
5289bc588SBarry Smith 
6c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h>
7c6db04a5SJed Brown #include <petscblaslapack.h>
8289bc588SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
11f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
121987afe7SBarry Smith {
131987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
14f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1513f74950SBarry Smith   PetscInt       j;
160805154bSBarry Smith   PetscBLASInt   N,m,ldax,lday,one = 1;
17efee365bSSatish Balay   PetscErrorCode ierr;
183a40ed3dSBarry Smith 
193a40ed3dSBarry Smith   PetscFunctionBegin;
20d0f46423SBarry Smith   N    = PetscBLASIntCast(X->rmap->n*X->cmap->n);
21d0f46423SBarry Smith   m    = PetscBLASIntCast(X->rmap->n);
220805154bSBarry Smith   ldax = PetscBLASIntCast(x->lda);
230805154bSBarry Smith   lday = PetscBLASIntCast(y->lda);
24a5ce6ee0Svictorle   if (ldax>m || lday>m) {
25d0f46423SBarry Smith     for (j=0; j<X->cmap->n; j++) {
26f4df32b1SMatthew Knepley       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
27a5ce6ee0Svictorle     }
28a5ce6ee0Svictorle   } else {
29f4df32b1SMatthew Knepley     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
30a5ce6ee0Svictorle   }
310450473dSBarry Smith   ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr);
323a40ed3dSBarry Smith   PetscFunctionReturn(0);
331987afe7SBarry Smith }
341987afe7SBarry Smith 
354a2ae208SSatish Balay #undef __FUNCT__
364a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
37dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
38289bc588SBarry Smith {
39d0f46423SBarry Smith   PetscInt     N = A->rmap->n*A->cmap->n;
403a40ed3dSBarry Smith 
413a40ed3dSBarry Smith   PetscFunctionBegin;
424e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
434e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
446de62eeeSBarry Smith   info->nz_used           = (double)N;
456de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
464e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
474e220ebcSLois Curfman McInnes   info->mallocs           = 0;
487adad957SLisandro Dalcin   info->memory            = ((PetscObject)A)->mem;
494e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
504e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
514e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
523a40ed3dSBarry Smith   PetscFunctionReturn(0);
53289bc588SBarry Smith }
54289bc588SBarry Smith 
554a2ae208SSatish Balay #undef __FUNCT__
564a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
57f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
5880cd9d93SLois Curfman McInnes {
59273d9f13SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
60f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
61efee365bSSatish Balay   PetscErrorCode ierr;
620805154bSBarry Smith   PetscBLASInt   one = 1,j,nz,lda = PetscBLASIntCast(a->lda);
6380cd9d93SLois Curfman McInnes 
643a40ed3dSBarry Smith   PetscFunctionBegin;
65d0f46423SBarry Smith   if (lda>A->rmap->n) {
66d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n);
67d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
68f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
69a5ce6ee0Svictorle     }
70a5ce6ee0Svictorle   } else {
71d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n*A->cmap->n);
72f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
73a5ce6ee0Svictorle   }
74efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
753a40ed3dSBarry Smith   PetscFunctionReturn(0);
7680cd9d93SLois Curfman McInnes }
7780cd9d93SLois Curfman McInnes 
781cbb95d3SBarry Smith #undef __FUNCT__
791cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense"
80ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool  *fl)
811cbb95d3SBarry Smith {
821cbb95d3SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
83d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,N;
841cbb95d3SBarry Smith   PetscScalar    *v = a->v;
851cbb95d3SBarry Smith 
861cbb95d3SBarry Smith   PetscFunctionBegin;
871cbb95d3SBarry Smith   *fl = PETSC_FALSE;
88d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0);
891cbb95d3SBarry Smith   N = a->lda;
901cbb95d3SBarry Smith 
911cbb95d3SBarry Smith   for (i=0; i<m; i++) {
921cbb95d3SBarry Smith     for (j=i+1; j<m; j++) {
931cbb95d3SBarry Smith       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
941cbb95d3SBarry Smith     }
951cbb95d3SBarry Smith   }
961cbb95d3SBarry Smith   *fl = PETSC_TRUE;
971cbb95d3SBarry Smith   PetscFunctionReturn(0);
981cbb95d3SBarry Smith }
991cbb95d3SBarry Smith 
100b24902e0SBarry Smith #undef __FUNCT__
101b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense"
102719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
103b24902e0SBarry Smith {
104b24902e0SBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
105b24902e0SBarry Smith   PetscErrorCode ierr;
106b24902e0SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
107b24902e0SBarry Smith 
108b24902e0SBarry Smith   PetscFunctionBegin;
109aa5ea44dSBarry Smith   ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr);
110aa5ea44dSBarry Smith   ierr = PetscLayoutReference(A->cmap,&newi->cmap);CHKERRQ(ierr);
111719d5645SBarry Smith   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
112b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
113b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
114d0f46423SBarry Smith     if (lda>A->rmap->n) {
115d0f46423SBarry Smith       m = A->rmap->n;
116d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
117b24902e0SBarry Smith 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
118b24902e0SBarry Smith       }
119b24902e0SBarry Smith     } else {
120d0f46423SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
121b24902e0SBarry Smith     }
122b24902e0SBarry Smith   }
123b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
124b24902e0SBarry Smith   PetscFunctionReturn(0);
125b24902e0SBarry Smith }
126b24902e0SBarry Smith 
1274a2ae208SSatish Balay #undef __FUNCT__
1284a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
129dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
13002cad45dSBarry Smith {
1316849ba73SBarry Smith   PetscErrorCode ierr;
13202cad45dSBarry Smith 
1333a40ed3dSBarry Smith   PetscFunctionBegin;
1345c9eb25fSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
135d0f46423SBarry Smith   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1365c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
137719d5645SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
138b24902e0SBarry Smith   PetscFunctionReturn(0);
139b24902e0SBarry Smith }
140b24902e0SBarry Smith 
1416ee01492SSatish Balay 
1420481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
143719d5645SBarry Smith 
1444a2ae208SSatish Balay #undef __FUNCT__
1454a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1460481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
147289bc588SBarry Smith {
1484482741eSBarry Smith   MatFactorInfo  info;
149a093e273SMatthew Knepley   PetscErrorCode ierr;
1503a40ed3dSBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
152c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
153719d5645SBarry Smith   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
1543a40ed3dSBarry Smith   PetscFunctionReturn(0);
155289bc588SBarry Smith }
1566ee01492SSatish Balay 
1570b4b3355SBarry Smith #undef __FUNCT__
1584a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
159dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
160289bc588SBarry Smith {
161c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1626849ba73SBarry Smith   PetscErrorCode ierr;
16387828ca2SBarry Smith   PetscScalar    *x,*y;
164d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
16567e560aaSBarry Smith 
1663a40ed3dSBarry Smith   PetscFunctionBegin;
1671ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1681ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
169d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
170d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
171ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
172e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
173ae7cfcebSSatish Balay #else
17471044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
175e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
176ae7cfcebSSatish Balay #endif
177d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_CHOLESKY){
178ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
179e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
180ae7cfcebSSatish Balay #else
18171044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
182e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
183ae7cfcebSSatish Balay #endif
184289bc588SBarry Smith   }
185e32f2f54SBarry Smith   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
1861ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1871ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
188dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
1893a40ed3dSBarry Smith   PetscFunctionReturn(0);
190289bc588SBarry Smith }
1916ee01492SSatish Balay 
1924a2ae208SSatish Balay #undef __FUNCT__
19385e2c93fSHong Zhang #define __FUNCT__ "MatMatSolve_SeqDense"
19485e2c93fSHong Zhang PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
19585e2c93fSHong Zhang {
19685e2c93fSHong Zhang   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
19785e2c93fSHong Zhang   PetscErrorCode ierr;
19885e2c93fSHong Zhang   PetscScalar    *b,*x;
199efb80c78SLisandro Dalcin   PetscInt       n;
20085e2c93fSHong Zhang   PetscBLASInt   nrhs,info,m=PetscBLASIntCast(A->rmap->n);
201*bda8bf91SBarry Smith   PetscBool      flg;
20285e2c93fSHong Zhang 
20385e2c93fSHong Zhang   PetscFunctionBegin;
204*bda8bf91SBarry Smith   ierr = PetscTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
205*bda8bf91SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
206*bda8bf91SBarry Smith   ierr = PetscTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
207*bda8bf91SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
208*bda8bf91SBarry Smith 
209efb80c78SLisandro Dalcin   ierr = MatGetSize(B,PETSC_NULL,&n);CHKERRQ(ierr);
210efb80c78SLisandro Dalcin   nrhs = PetscBLASIntCast(n);
21185e2c93fSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
21285e2c93fSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
21385e2c93fSHong Zhang 
21485e2c93fSHong Zhang   ierr = PetscMemcpy(x,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
21585e2c93fSHong Zhang 
21685e2c93fSHong Zhang   if (A->factortype == MAT_FACTOR_LU) {
21785e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS)
21885e2c93fSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
21985e2c93fSHong Zhang #else
22085e2c93fSHong Zhang     LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info);
22185e2c93fSHong Zhang     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
22285e2c93fSHong Zhang #endif
22385e2c93fSHong Zhang   } else if (A->factortype == MAT_FACTOR_CHOLESKY){
22485e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS)
22585e2c93fSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
22685e2c93fSHong Zhang #else
22785e2c93fSHong Zhang     LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info);
22885e2c93fSHong Zhang     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
22985e2c93fSHong Zhang #endif
23085e2c93fSHong Zhang   }
23185e2c93fSHong Zhang   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
23285e2c93fSHong Zhang 
23385e2c93fSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
23485e2c93fSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
23585e2c93fSHong Zhang   ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr);
23685e2c93fSHong Zhang   PetscFunctionReturn(0);
23785e2c93fSHong Zhang }
23885e2c93fSHong Zhang 
23985e2c93fSHong Zhang #undef __FUNCT__
2404a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
241dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
242da3a660dSBarry Smith {
243c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
244dfbe8321SBarry Smith   PetscErrorCode ierr;
24587828ca2SBarry Smith   PetscScalar    *x,*y;
246d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
24767e560aaSBarry Smith 
2483a40ed3dSBarry Smith   PetscFunctionBegin;
2491ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2501ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
251d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
252752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
253da3a660dSBarry Smith   if (mat->pivots) {
254ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
255e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
256ae7cfcebSSatish Balay #else
25771044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
258e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
259ae7cfcebSSatish Balay #endif
2607a97a34bSBarry Smith   } else {
261ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
262e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
263ae7cfcebSSatish Balay #else
26471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
265e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
266ae7cfcebSSatish Balay #endif
267da3a660dSBarry Smith   }
2681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2691ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
270dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
2713a40ed3dSBarry Smith   PetscFunctionReturn(0);
272da3a660dSBarry Smith }
2736ee01492SSatish Balay 
2744a2ae208SSatish Balay #undef __FUNCT__
2754a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
276dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
277da3a660dSBarry Smith {
278c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
279dfbe8321SBarry Smith   PetscErrorCode ierr;
28087828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
281da3a660dSBarry Smith   Vec            tmp = 0;
282d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
28367e560aaSBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
2851ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2861ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
287d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
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)
297e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
298ae7cfcebSSatish Balay #else
29971044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
300e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
301ae7cfcebSSatish Balay #endif
302a8c6a408SBarry Smith   } else {
303ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
304e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
305ae7cfcebSSatish Balay #else
30671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
307e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
308ae7cfcebSSatish Balay #endif
309da3a660dSBarry Smith   }
3106bf464f9SBarry Smith   if (tmp) {
3116bf464f9SBarry Smith     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
3126bf464f9SBarry Smith     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
3136bf464f9SBarry Smith   } else {
3146bf464f9SBarry Smith     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
3156bf464f9SBarry Smith   }
3161ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3171ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
318dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
320da3a660dSBarry Smith }
32167e560aaSBarry Smith 
3224a2ae208SSatish Balay #undef __FUNCT__
3234a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
324dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
325da3a660dSBarry Smith {
326c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3276849ba73SBarry Smith   PetscErrorCode ierr;
32887828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
329da3a660dSBarry Smith   Vec            tmp;
330d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
33167e560aaSBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
333d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
3341ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3351ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
336da3a660dSBarry Smith   if (yy == zz) {
33778b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
33852e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
33978b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
340da3a660dSBarry Smith   }
341d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
342752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
343da3a660dSBarry Smith   if (mat->pivots) {
344ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
345e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
346ae7cfcebSSatish Balay #else
34771044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
348e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
349ae7cfcebSSatish Balay #endif
3503a40ed3dSBarry Smith   } else {
351ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
352e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
353ae7cfcebSSatish Balay #else
35471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
355e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
356ae7cfcebSSatish Balay #endif
357da3a660dSBarry Smith   }
35890f02eecSBarry Smith   if (tmp) {
3592dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
3606bf464f9SBarry Smith     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
3613a40ed3dSBarry Smith   } else {
3622dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
36390f02eecSBarry Smith   }
3641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3651ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
366dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
3673a40ed3dSBarry Smith   PetscFunctionReturn(0);
368da3a660dSBarry Smith }
369db4efbfdSBarry Smith 
370db4efbfdSBarry Smith /* ---------------------------------------------------------------*/
371db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
372db4efbfdSBarry Smith    rather than put it in the Mat->row slot.*/
373db4efbfdSBarry Smith #undef __FUNCT__
374db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense"
3750481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
376db4efbfdSBarry Smith {
377db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
378db4efbfdSBarry Smith   PetscFunctionBegin;
379e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
380db4efbfdSBarry Smith #else
381db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
382db4efbfdSBarry Smith   PetscErrorCode ierr;
383db4efbfdSBarry Smith   PetscBLASInt   n,m,info;
384db4efbfdSBarry Smith 
385db4efbfdSBarry Smith   PetscFunctionBegin;
386db4efbfdSBarry Smith   n = PetscBLASIntCast(A->cmap->n);
387db4efbfdSBarry Smith   m = PetscBLASIntCast(A->rmap->n);
388db4efbfdSBarry Smith   if (!mat->pivots) {
389db4efbfdSBarry Smith     ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
390db4efbfdSBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
391db4efbfdSBarry Smith   }
392db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
393db4efbfdSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
394e32f2f54SBarry Smith   if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
395e32f2f54SBarry Smith   if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
396db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
397db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
398db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
399db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
400d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_LU;
401db4efbfdSBarry Smith 
402dc0b31edSSatish Balay   ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
403db4efbfdSBarry Smith #endif
404db4efbfdSBarry Smith   PetscFunctionReturn(0);
405db4efbfdSBarry Smith }
406db4efbfdSBarry Smith 
407db4efbfdSBarry Smith #undef __FUNCT__
408db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense"
4090481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
410db4efbfdSBarry Smith {
411db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
412db4efbfdSBarry Smith   PetscFunctionBegin;
413e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
414db4efbfdSBarry Smith #else
415db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
416db4efbfdSBarry Smith   PetscErrorCode ierr;
417db4efbfdSBarry Smith   PetscBLASInt   info,n = PetscBLASIntCast(A->cmap->n);
418db4efbfdSBarry Smith 
419db4efbfdSBarry Smith   PetscFunctionBegin;
420db4efbfdSBarry Smith   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
421db4efbfdSBarry Smith 
422db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
423db4efbfdSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
424e32f2f54SBarry Smith   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
425db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
426db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
427db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
428db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
429d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_CHOLESKY;
430dc0b31edSSatish Balay   ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr);
431db4efbfdSBarry Smith #endif
432db4efbfdSBarry Smith   PetscFunctionReturn(0);
433db4efbfdSBarry Smith }
434db4efbfdSBarry Smith 
435db4efbfdSBarry Smith 
436db4efbfdSBarry Smith #undef __FUNCT__
437db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
4380481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
439db4efbfdSBarry Smith {
440db4efbfdSBarry Smith   PetscErrorCode ierr;
441db4efbfdSBarry Smith   MatFactorInfo  info;
442db4efbfdSBarry Smith 
443db4efbfdSBarry Smith   PetscFunctionBegin;
444db4efbfdSBarry Smith   info.fill = 1.0;
445c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
446719d5645SBarry Smith   ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr);
447db4efbfdSBarry Smith   PetscFunctionReturn(0);
448db4efbfdSBarry Smith }
449db4efbfdSBarry Smith 
450db4efbfdSBarry Smith #undef __FUNCT__
451db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
4520481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
453db4efbfdSBarry Smith {
454db4efbfdSBarry Smith   PetscFunctionBegin;
455c3ef05f6SHong Zhang   fact->assembled                  = PETSC_TRUE;
456719d5645SBarry Smith   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
457db4efbfdSBarry Smith   PetscFunctionReturn(0);
458db4efbfdSBarry Smith }
459db4efbfdSBarry Smith 
460db4efbfdSBarry Smith #undef __FUNCT__
461db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
4620481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
463db4efbfdSBarry Smith {
464db4efbfdSBarry Smith   PetscFunctionBegin;
465c3ef05f6SHong Zhang   fact->assembled            = PETSC_TRUE;
466719d5645SBarry Smith   fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
467db4efbfdSBarry Smith   PetscFunctionReturn(0);
468db4efbfdSBarry Smith }
469db4efbfdSBarry Smith 
470bb5747d9SMatthew Knepley EXTERN_C_BEGIN
471db4efbfdSBarry Smith #undef __FUNCT__
472db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc"
473db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
474db4efbfdSBarry Smith {
475db4efbfdSBarry Smith   PetscErrorCode ierr;
476db4efbfdSBarry Smith 
477db4efbfdSBarry Smith   PetscFunctionBegin;
478db4efbfdSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
479db4efbfdSBarry Smith   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
480db4efbfdSBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
481db4efbfdSBarry Smith   if (ftype == MAT_FACTOR_LU){
482db4efbfdSBarry Smith     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
483db4efbfdSBarry Smith   } else {
484db4efbfdSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
485db4efbfdSBarry Smith   }
486d5f3da31SBarry Smith   (*fact)->factortype = ftype;
487db4efbfdSBarry Smith   PetscFunctionReturn(0);
488db4efbfdSBarry Smith }
489bb5747d9SMatthew Knepley EXTERN_C_END
490db4efbfdSBarry Smith 
491289bc588SBarry Smith /* ------------------------------------------------------------------*/
4924a2ae208SSatish Balay #undef __FUNCT__
49341f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense"
49441f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
495289bc588SBarry Smith {
496c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
49787828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
498dfbe8321SBarry Smith   PetscErrorCode ierr;
499d0f46423SBarry Smith   PetscInt       m = A->rmap->n,i;
500aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
5010805154bSBarry Smith   PetscBLASInt   o = 1,bm = PetscBLASIntCast(m);
502bc1b551cSSatish Balay #endif
503289bc588SBarry Smith 
5043a40ed3dSBarry Smith   PetscFunctionBegin;
505289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
50671044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
5072dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
508289bc588SBarry Smith   }
5091ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5101ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
511b965ef7fSBarry Smith   its  = its*lits;
512e32f2f54SBarry Smith   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
513289bc588SBarry Smith   while (its--) {
514fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
515289bc588SBarry Smith       for (i=0; i<m; i++) {
516aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
517f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
518f1747703SBarry Smith            not happy about returning a double complex */
51913f74950SBarry Smith         PetscInt    _i;
52087828ca2SBarry Smith         PetscScalar sum = b[i];
521f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
5223f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
523f1747703SBarry Smith         }
524f1747703SBarry Smith         xt = sum;
525f1747703SBarry Smith #else
52671044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
527f1747703SBarry Smith #endif
52855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
529289bc588SBarry Smith       }
530289bc588SBarry Smith     }
531fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
532289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
533aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
534f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
535f1747703SBarry Smith            not happy about returning a double complex */
53613f74950SBarry Smith         PetscInt    _i;
53787828ca2SBarry Smith         PetscScalar sum = b[i];
538f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
5393f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
540f1747703SBarry Smith         }
541f1747703SBarry Smith         xt = sum;
542f1747703SBarry Smith #else
54371044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
544f1747703SBarry Smith #endif
54555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
546289bc588SBarry Smith       }
547289bc588SBarry Smith     }
548289bc588SBarry Smith   }
5491ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
5501ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5513a40ed3dSBarry Smith   PetscFunctionReturn(0);
552289bc588SBarry Smith }
553289bc588SBarry Smith 
554289bc588SBarry Smith /* -----------------------------------------------------------------*/
5554a2ae208SSatish Balay #undef __FUNCT__
5564a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
557dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
558289bc588SBarry Smith {
559c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
56087828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
561dfbe8321SBarry Smith   PetscErrorCode ierr;
5620805154bSBarry Smith   PetscBLASInt   m, n,_One=1;
563ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
5643a40ed3dSBarry Smith 
5653a40ed3dSBarry Smith   PetscFunctionBegin;
566d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
567d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
568d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5691ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5701ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
57171044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
5721ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5731ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
574dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
5753a40ed3dSBarry Smith   PetscFunctionReturn(0);
576289bc588SBarry Smith }
577800995b7SMatthew Knepley 
5784a2ae208SSatish Balay #undef __FUNCT__
5794a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
580dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
581289bc588SBarry Smith {
582c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
58387828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
584dfbe8321SBarry Smith   PetscErrorCode ierr;
5850805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5863a40ed3dSBarry Smith 
5873a40ed3dSBarry Smith   PetscFunctionBegin;
588d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
589d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
590d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5911ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5921ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
59371044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5941ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5951ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
596dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
5973a40ed3dSBarry Smith   PetscFunctionReturn(0);
598289bc588SBarry Smith }
5996ee01492SSatish Balay 
6004a2ae208SSatish Balay #undef __FUNCT__
6014a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
602dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
603289bc588SBarry Smith {
604c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
60587828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
606dfbe8321SBarry Smith   PetscErrorCode ierr;
6070805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
6083a40ed3dSBarry Smith 
6093a40ed3dSBarry Smith   PetscFunctionBegin;
610d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
611d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
612d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
613600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
6141ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6151ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
61671044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
6171ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6181ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
619dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
6203a40ed3dSBarry Smith   PetscFunctionReturn(0);
621289bc588SBarry Smith }
6226ee01492SSatish Balay 
6234a2ae208SSatish Balay #undef __FUNCT__
6244a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
625dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
626289bc588SBarry Smith {
627c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
62887828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
629dfbe8321SBarry Smith   PetscErrorCode ierr;
6300805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
63187828ca2SBarry Smith   PetscScalar    _DOne=1.0;
6323a40ed3dSBarry Smith 
6333a40ed3dSBarry Smith   PetscFunctionBegin;
634d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
635d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
636d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
637600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
6381ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6391ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
64071044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
6411ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6421ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
643dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
6443a40ed3dSBarry Smith   PetscFunctionReturn(0);
645289bc588SBarry Smith }
646289bc588SBarry Smith 
647289bc588SBarry Smith /* -----------------------------------------------------------------*/
6484a2ae208SSatish Balay #undef __FUNCT__
6494a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
65013f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
651289bc588SBarry Smith {
652c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
65387828ca2SBarry Smith   PetscScalar    *v;
6546849ba73SBarry Smith   PetscErrorCode ierr;
65513f74950SBarry Smith   PetscInt       i;
65667e560aaSBarry Smith 
6573a40ed3dSBarry Smith   PetscFunctionBegin;
658d0f46423SBarry Smith   *ncols = A->cmap->n;
659289bc588SBarry Smith   if (cols) {
660d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
661d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
662289bc588SBarry Smith   }
663289bc588SBarry Smith   if (vals) {
664d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
665289bc588SBarry Smith     v    = mat->v + row;
666d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
667289bc588SBarry Smith   }
6683a40ed3dSBarry Smith   PetscFunctionReturn(0);
669289bc588SBarry Smith }
6706ee01492SSatish Balay 
6714a2ae208SSatish Balay #undef __FUNCT__
6724a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
67313f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
674289bc588SBarry Smith {
675dfbe8321SBarry Smith   PetscErrorCode ierr;
676606d414cSSatish Balay   PetscFunctionBegin;
677606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
678606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6793a40ed3dSBarry Smith   PetscFunctionReturn(0);
680289bc588SBarry Smith }
681289bc588SBarry Smith /* ----------------------------------------------------------------*/
6824a2ae208SSatish Balay #undef __FUNCT__
6834a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
68413f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
685289bc588SBarry Smith {
686c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
68713f74950SBarry Smith   PetscInt     i,j,idx=0;
688d6dfbf8fSBarry Smith 
6893a40ed3dSBarry Smith   PetscFunctionBegin;
69071fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
691289bc588SBarry Smith   if (!mat->roworiented) {
692dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
693289bc588SBarry Smith       for (j=0; j<n; j++) {
694cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6952515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
696e32f2f54SBarry Smith         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
69758804f6dSBarry Smith #endif
698289bc588SBarry Smith         for (i=0; i<m; i++) {
699cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
7002515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
701e32f2f54SBarry Smith           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
70258804f6dSBarry Smith #endif
703cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
704289bc588SBarry Smith         }
705289bc588SBarry Smith       }
7063a40ed3dSBarry Smith     } else {
707289bc588SBarry Smith       for (j=0; j<n; j++) {
708cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
7092515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
710e32f2f54SBarry Smith         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
71158804f6dSBarry Smith #endif
712289bc588SBarry Smith         for (i=0; i<m; i++) {
713cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
7142515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
715e32f2f54SBarry Smith           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
71658804f6dSBarry Smith #endif
717cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
718289bc588SBarry Smith         }
719289bc588SBarry Smith       }
720289bc588SBarry Smith     }
7213a40ed3dSBarry Smith   } else {
722dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
723e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
724cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
7252515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
726e32f2f54SBarry Smith         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
72758804f6dSBarry Smith #endif
728e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
729cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
7302515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
731e32f2f54SBarry Smith           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
73258804f6dSBarry Smith #endif
733cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
734e8d4e0b9SBarry Smith         }
735e8d4e0b9SBarry Smith       }
7363a40ed3dSBarry Smith     } else {
737289bc588SBarry Smith       for (i=0; i<m; i++) {
738cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
7392515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
740e32f2f54SBarry Smith         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
74158804f6dSBarry Smith #endif
742289bc588SBarry Smith         for (j=0; j<n; j++) {
743cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
7442515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
745e32f2f54SBarry Smith           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
74658804f6dSBarry Smith #endif
747cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
748289bc588SBarry Smith         }
749289bc588SBarry Smith       }
750289bc588SBarry Smith     }
751e8d4e0b9SBarry Smith   }
7523a40ed3dSBarry Smith   PetscFunctionReturn(0);
753289bc588SBarry Smith }
754e8d4e0b9SBarry Smith 
7554a2ae208SSatish Balay #undef __FUNCT__
7564a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
75713f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
758ae80bb75SLois Curfman McInnes {
759ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
76013f74950SBarry Smith   PetscInt     i,j;
761ae80bb75SLois Curfman McInnes 
7623a40ed3dSBarry Smith   PetscFunctionBegin;
763ae80bb75SLois Curfman McInnes   /* row-oriented output */
764ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
76597e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
766e32f2f54SBarry Smith     if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n);
767ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7686f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
769e32f2f54SBarry Smith       if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n);
77097e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
771ae80bb75SLois Curfman McInnes     }
772ae80bb75SLois Curfman McInnes   }
7733a40ed3dSBarry Smith   PetscFunctionReturn(0);
774ae80bb75SLois Curfman McInnes }
775ae80bb75SLois Curfman McInnes 
776289bc588SBarry Smith /* -----------------------------------------------------------------*/
777289bc588SBarry Smith 
7784a2ae208SSatish Balay #undef __FUNCT__
7795bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense"
780112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer)
781aabbc4fbSShri Abhyankar {
782aabbc4fbSShri Abhyankar   Mat_SeqDense   *a;
783aabbc4fbSShri Abhyankar   PetscErrorCode ierr;
784aabbc4fbSShri Abhyankar   PetscInt       *scols,i,j,nz,header[4];
785aabbc4fbSShri Abhyankar   int            fd;
786aabbc4fbSShri Abhyankar   PetscMPIInt    size;
787aabbc4fbSShri Abhyankar   PetscInt       *rowlengths = 0,M,N,*cols,grows,gcols;
788aabbc4fbSShri Abhyankar   PetscScalar    *vals,*svals,*v,*w;
789aabbc4fbSShri Abhyankar   MPI_Comm       comm = ((PetscObject)viewer)->comm;
790aabbc4fbSShri Abhyankar 
791aabbc4fbSShri Abhyankar   PetscFunctionBegin;
792aabbc4fbSShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
793aabbc4fbSShri Abhyankar   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
794aabbc4fbSShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
795aabbc4fbSShri Abhyankar   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
796aabbc4fbSShri Abhyankar   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
797aabbc4fbSShri Abhyankar   M = header[1]; N = header[2]; nz = header[3];
798aabbc4fbSShri Abhyankar 
799aabbc4fbSShri Abhyankar   /* set global size if not set already*/
800aabbc4fbSShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
801aabbc4fbSShri Abhyankar     ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr);
802aabbc4fbSShri Abhyankar   } else {
803aabbc4fbSShri Abhyankar     /* if sizes and type are already set, check if the vector global sizes are correct */
804aabbc4fbSShri Abhyankar     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
805aabbc4fbSShri Abhyankar     if (M != grows ||  N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols);
806aabbc4fbSShri Abhyankar   }
807aabbc4fbSShri Abhyankar   ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
808aabbc4fbSShri Abhyankar 
809aabbc4fbSShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
810aabbc4fbSShri Abhyankar     a    = (Mat_SeqDense*)newmat->data;
811aabbc4fbSShri Abhyankar     v    = a->v;
812aabbc4fbSShri Abhyankar     /* Allocate some temp space to read in the values and then flip them
813aabbc4fbSShri Abhyankar        from row major to column major */
814aabbc4fbSShri Abhyankar     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
815aabbc4fbSShri Abhyankar     /* read in nonzero values */
816aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
817aabbc4fbSShri Abhyankar     /* now flip the values and store them in the matrix*/
818aabbc4fbSShri Abhyankar     for (j=0; j<N; j++) {
819aabbc4fbSShri Abhyankar       for (i=0; i<M; i++) {
820aabbc4fbSShri Abhyankar         *v++ =w[i*N+j];
821aabbc4fbSShri Abhyankar       }
822aabbc4fbSShri Abhyankar     }
823aabbc4fbSShri Abhyankar     ierr = PetscFree(w);CHKERRQ(ierr);
824aabbc4fbSShri Abhyankar   } else {
825aabbc4fbSShri Abhyankar     /* read row lengths */
826aabbc4fbSShri Abhyankar     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
827aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
828aabbc4fbSShri Abhyankar 
829aabbc4fbSShri Abhyankar     a = (Mat_SeqDense*)newmat->data;
830aabbc4fbSShri Abhyankar     v = a->v;
831aabbc4fbSShri Abhyankar 
832aabbc4fbSShri Abhyankar     /* read column indices and nonzeros */
833aabbc4fbSShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
834aabbc4fbSShri Abhyankar     cols = scols;
835aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
836aabbc4fbSShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
837aabbc4fbSShri Abhyankar     vals = svals;
838aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
839aabbc4fbSShri Abhyankar 
840aabbc4fbSShri Abhyankar     /* insert into matrix */
841aabbc4fbSShri Abhyankar     for (i=0; i<M; i++) {
842aabbc4fbSShri Abhyankar       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
843aabbc4fbSShri Abhyankar       svals += rowlengths[i]; scols += rowlengths[i];
844aabbc4fbSShri Abhyankar     }
845aabbc4fbSShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
846aabbc4fbSShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
847aabbc4fbSShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
848aabbc4fbSShri Abhyankar   }
849aabbc4fbSShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
850aabbc4fbSShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
851aabbc4fbSShri Abhyankar 
852aabbc4fbSShri Abhyankar   PetscFunctionReturn(0);
853aabbc4fbSShri Abhyankar }
854aabbc4fbSShri Abhyankar 
855aabbc4fbSShri Abhyankar #undef __FUNCT__
8564a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8576849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
858289bc588SBarry Smith {
859932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
860dfbe8321SBarry Smith   PetscErrorCode    ierr;
86113f74950SBarry Smith   PetscInt          i,j;
8622dcb1b2aSMatthew Knepley   const char        *name;
86387828ca2SBarry Smith   PetscScalar       *v;
864f3ef73ceSBarry Smith   PetscViewerFormat format;
8655f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
866ace3abfcSBarry Smith   PetscBool         allreal = PETSC_TRUE;
8675f481a85SSatish Balay #endif
868932b0c3eSLois Curfman McInnes 
8693a40ed3dSBarry Smith   PetscFunctionBegin;
870b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
871456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8723a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
873fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
874d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
8757566de4bSShri Abhyankar     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
876d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
87744cd7ae7SLois Curfman McInnes       v = a->v + i;
87877431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
879d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
880aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
881329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
882a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
883329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
884a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8856831982aSBarry Smith         }
88680cd9d93SLois Curfman McInnes #else
8876831982aSBarry Smith         if (*v) {
888a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8896831982aSBarry Smith         }
89080cd9d93SLois Curfman McInnes #endif
8911b807ce4Svictorle         v += a->lda;
89280cd9d93SLois Curfman McInnes       }
893b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
89480cd9d93SLois Curfman McInnes     }
895d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
8963a40ed3dSBarry Smith   } else {
897d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
898aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
89947989497SBarry Smith     /* determine if matrix has all real values */
90047989497SBarry Smith     v = a->v;
901d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
902ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
90347989497SBarry Smith     }
90447989497SBarry Smith #endif
905fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
9063a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
907d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
908d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
909fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
9107566de4bSShri Abhyankar     } else {
9117566de4bSShri Abhyankar       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
912ffac6cdbSBarry Smith     }
913ffac6cdbSBarry Smith 
914d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
915932b0c3eSLois Curfman McInnes       v = a->v + i;
916d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
917aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
91847989497SBarry Smith         if (allreal) {
919f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
92047989497SBarry Smith         } else {
921f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
92247989497SBarry Smith         }
923289bc588SBarry Smith #else
924f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
925289bc588SBarry Smith #endif
9261b807ce4Svictorle 	v += a->lda;
927289bc588SBarry Smith       }
928b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
929289bc588SBarry Smith     }
930fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
931b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
932ffac6cdbSBarry Smith     }
933d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
934da3a660dSBarry Smith   }
935b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
9363a40ed3dSBarry Smith   PetscFunctionReturn(0);
937289bc588SBarry Smith }
938289bc588SBarry Smith 
9394a2ae208SSatish Balay #undef __FUNCT__
9404a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
9416849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
942932b0c3eSLois Curfman McInnes {
943932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9446849ba73SBarry Smith   PetscErrorCode    ierr;
94513f74950SBarry Smith   int               fd;
946d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
947f4403165SShri Abhyankar   PetscScalar       *v,*anonz,*vals;
948f4403165SShri Abhyankar   PetscViewerFormat format;
949932b0c3eSLois Curfman McInnes 
9503a40ed3dSBarry Smith   PetscFunctionBegin;
951b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
95290ace30eSBarry Smith 
953f4403165SShri Abhyankar   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
954f4403165SShri Abhyankar   if (format == PETSC_VIEWER_NATIVE) {
955f4403165SShri Abhyankar     /* store the matrix as a dense matrix */
956f4403165SShri Abhyankar     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
957f4403165SShri Abhyankar     col_lens[0] = MAT_FILE_CLASSID;
958f4403165SShri Abhyankar     col_lens[1] = m;
959f4403165SShri Abhyankar     col_lens[2] = n;
960f4403165SShri Abhyankar     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
961f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
962f4403165SShri Abhyankar     ierr = PetscFree(col_lens);CHKERRQ(ierr);
963f4403165SShri Abhyankar 
964f4403165SShri Abhyankar     /* write out matrix, by rows */
965f4403165SShri Abhyankar     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
966f4403165SShri Abhyankar     v    = a->v;
967f4403165SShri Abhyankar     for (j=0; j<n; j++) {
968f4403165SShri Abhyankar       for (i=0; i<m; i++) {
969f4403165SShri Abhyankar         vals[j + i*n] = *v++;
970f4403165SShri Abhyankar       }
971f4403165SShri Abhyankar     }
972f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
973f4403165SShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
974f4403165SShri Abhyankar   } else {
97513f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
9760700a824SBarry Smith     col_lens[0] = MAT_FILE_CLASSID;
977932b0c3eSLois Curfman McInnes     col_lens[1] = m;
978932b0c3eSLois Curfman McInnes     col_lens[2] = n;
979932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
980932b0c3eSLois Curfman McInnes 
981932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
982932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9836f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
984932b0c3eSLois Curfman McInnes 
985932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
986932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
987932b0c3eSLois Curfman McInnes     ict = 0;
988932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
989932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
990932b0c3eSLois Curfman McInnes     }
9916f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
992606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
993932b0c3eSLois Curfman McInnes 
994932b0c3eSLois Curfman McInnes     /* store nonzero values */
99587828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
996932b0c3eSLois Curfman McInnes     ict  = 0;
997932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
998932b0c3eSLois Curfman McInnes       v = a->v + i;
999932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
10001b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
1001932b0c3eSLois Curfman McInnes       }
1002932b0c3eSLois Curfman McInnes     }
10036f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
1004606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
1005f4403165SShri Abhyankar   }
10063a40ed3dSBarry Smith   PetscFunctionReturn(0);
1007932b0c3eSLois Curfman McInnes }
1008932b0c3eSLois Curfman McInnes 
10094a2ae208SSatish Balay #undef __FUNCT__
10104a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
1011dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1012f1af5d2fSBarry Smith {
1013f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
1014f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
10156849ba73SBarry Smith   PetscErrorCode    ierr;
1016d0f46423SBarry Smith   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
101787828ca2SBarry Smith   PetscScalar       *v = a->v;
1018b0a32e0cSBarry Smith   PetscViewer       viewer;
1019b0a32e0cSBarry Smith   PetscDraw         popup;
1020329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
1021f3ef73ceSBarry Smith   PetscViewerFormat format;
1022f1af5d2fSBarry Smith 
1023f1af5d2fSBarry Smith   PetscFunctionBegin;
1024f1af5d2fSBarry Smith 
1025f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
1026b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1027b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1028f1af5d2fSBarry Smith 
1029f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
1030fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1031f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
1032b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
1033f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1034f1af5d2fSBarry Smith       x_l = j;
1035f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1036f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1037f1af5d2fSBarry Smith         y_l = m - i - 1.0;
1038f1af5d2fSBarry Smith         y_r = y_l + 1.0;
1039f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
1040329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
1041b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1042329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
1043b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1044f1af5d2fSBarry Smith         } else {
1045f1af5d2fSBarry Smith           continue;
1046f1af5d2fSBarry Smith         }
1047f1af5d2fSBarry Smith #else
1048f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
1049b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1050f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1051b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1052f1af5d2fSBarry Smith         } else {
1053f1af5d2fSBarry Smith           continue;
1054f1af5d2fSBarry Smith         }
1055f1af5d2fSBarry Smith #endif
1056b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1057f1af5d2fSBarry Smith       }
1058f1af5d2fSBarry Smith     }
1059f1af5d2fSBarry Smith   } else {
1060f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1061f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1062f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1063f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1064f1af5d2fSBarry Smith     }
1065b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1066b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1067b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1068f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1069f1af5d2fSBarry Smith       x_l = j;
1070f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1071f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1072f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1073f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1074b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1075b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1076f1af5d2fSBarry Smith       }
1077f1af5d2fSBarry Smith     }
1078f1af5d2fSBarry Smith   }
1079f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1080f1af5d2fSBarry Smith }
1081f1af5d2fSBarry Smith 
10824a2ae208SSatish Balay #undef __FUNCT__
10834a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1084dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1085f1af5d2fSBarry Smith {
1086b0a32e0cSBarry Smith   PetscDraw      draw;
1087ace3abfcSBarry Smith   PetscBool      isnull;
1088329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1089dfbe8321SBarry Smith   PetscErrorCode ierr;
1090f1af5d2fSBarry Smith 
1091f1af5d2fSBarry Smith   PetscFunctionBegin;
1092b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1093b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1094abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1095f1af5d2fSBarry Smith 
1096f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1097d0f46423SBarry Smith   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1098f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1099b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1100b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1101f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1102f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1103f1af5d2fSBarry Smith }
1104f1af5d2fSBarry Smith 
11054a2ae208SSatish Balay #undef __FUNCT__
11064a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1107dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1108932b0c3eSLois Curfman McInnes {
1109dfbe8321SBarry Smith   PetscErrorCode ierr;
1110ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw;
1111932b0c3eSLois Curfman McInnes 
11123a40ed3dSBarry Smith   PetscFunctionBegin;
11132692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11142692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
11152692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
11160f5bd95cSBarry Smith 
1117c45a1595SBarry Smith   if (iascii) {
1118c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
11190f5bd95cSBarry Smith   } else if (isbinary) {
11203a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1121f1af5d2fSBarry Smith   } else if (isdraw) {
1122f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
11235cd90555SBarry Smith   } else {
1124e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1125932b0c3eSLois Curfman McInnes   }
11263a40ed3dSBarry Smith   PetscFunctionReturn(0);
1127932b0c3eSLois Curfman McInnes }
1128289bc588SBarry Smith 
11294a2ae208SSatish Balay #undef __FUNCT__
11304a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1131dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1132289bc588SBarry Smith {
1133ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1134dfbe8321SBarry Smith   PetscErrorCode ierr;
113590f02eecSBarry Smith 
11363a40ed3dSBarry Smith   PetscFunctionBegin;
1137aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1138d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1139a5a9c739SBarry Smith #endif
114005b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
11416857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1142bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
1143dbd8c25aSHong Zhang 
1144dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1145901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
11464ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11474ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11484ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11493a40ed3dSBarry Smith   PetscFunctionReturn(0);
1150289bc588SBarry Smith }
1151289bc588SBarry Smith 
11524a2ae208SSatish Balay #undef __FUNCT__
11534a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1154fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1155289bc588SBarry Smith {
1156c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11576849ba73SBarry Smith   PetscErrorCode ierr;
115813f74950SBarry Smith   PetscInt       k,j,m,n,M;
115987828ca2SBarry Smith   PetscScalar    *v,tmp;
116048b35521SBarry Smith 
11613a40ed3dSBarry Smith   PetscFunctionBegin;
1162d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1163e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1164e7e72b3dSBarry Smith     if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1165e7e72b3dSBarry Smith     else {
1166d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1167289bc588SBarry Smith         for (k=0; k<j; k++) {
11681b807ce4Svictorle           tmp = v[j + k*M];
11691b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11701b807ce4Svictorle           v[k + j*M] = tmp;
1171289bc588SBarry Smith         }
1172289bc588SBarry Smith       }
1173d64ed03dSBarry Smith     }
11743a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1175d3e5ee88SLois Curfman McInnes     Mat          tmat;
1176ec8511deSBarry Smith     Mat_SeqDense *tmatd;
117787828ca2SBarry Smith     PetscScalar  *v2;
1178ea709b57SSatish Balay 
1179fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11807adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1181d0f46423SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
11827adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11835c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1184fc4dec0aSBarry Smith     } else {
1185fc4dec0aSBarry Smith       tmat = *matout;
1186fc4dec0aSBarry Smith     }
1187ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11880de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1189d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11901b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1191d3e5ee88SLois Curfman McInnes     }
11926d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11936d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1194d3e5ee88SLois Curfman McInnes     *matout = tmat;
119548b35521SBarry Smith   }
11963a40ed3dSBarry Smith   PetscFunctionReturn(0);
1197289bc588SBarry Smith }
1198289bc588SBarry Smith 
11994a2ae208SSatish Balay #undef __FUNCT__
12004a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1201ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool  *flg)
1202289bc588SBarry Smith {
1203c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1204c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
120513f74950SBarry Smith   PetscInt     i,j;
120687828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
12079ea5d5aeSSatish Balay 
12083a40ed3dSBarry Smith   PetscFunctionBegin;
1209d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1210d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1211d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
12121b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1213d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
12143a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
12151b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
12161b807ce4Svictorle     }
1217289bc588SBarry Smith   }
121877c4ece6SBarry Smith   *flg = PETSC_TRUE;
12193a40ed3dSBarry Smith   PetscFunctionReturn(0);
1220289bc588SBarry Smith }
1221289bc588SBarry Smith 
12224a2ae208SSatish Balay #undef __FUNCT__
12234a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1224dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1225289bc588SBarry Smith {
1226c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1227dfbe8321SBarry Smith   PetscErrorCode ierr;
122813f74950SBarry Smith   PetscInt       i,n,len;
122987828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
123044cd7ae7SLois Curfman McInnes 
12313a40ed3dSBarry Smith   PetscFunctionBegin;
12322dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
12337a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
12341ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1235d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1236e32f2f54SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
123744cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
12381b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1239289bc588SBarry Smith   }
12401ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
12413a40ed3dSBarry Smith   PetscFunctionReturn(0);
1242289bc588SBarry Smith }
1243289bc588SBarry Smith 
12444a2ae208SSatish Balay #undef __FUNCT__
12454a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1246dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1247289bc588SBarry Smith {
1248c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
124987828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1250dfbe8321SBarry Smith   PetscErrorCode ierr;
1251d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
125255659b69SBarry Smith 
12533a40ed3dSBarry Smith   PetscFunctionBegin;
125428988994SBarry Smith   if (ll) {
12557a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12561ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1257e32f2f54SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1258da3a660dSBarry Smith     for (i=0; i<m; i++) {
1259da3a660dSBarry Smith       x = l[i];
1260da3a660dSBarry Smith       v = mat->v + i;
1261da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1262da3a660dSBarry Smith     }
12631ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1264efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1265da3a660dSBarry Smith   }
126628988994SBarry Smith   if (rr) {
12677a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12681ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1269e32f2f54SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1270da3a660dSBarry Smith     for (i=0; i<n; i++) {
1271da3a660dSBarry Smith       x = r[i];
1272da3a660dSBarry Smith       v = mat->v + i*m;
1273da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1274da3a660dSBarry Smith     }
12751ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1276efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1277da3a660dSBarry Smith   }
12783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1279289bc588SBarry Smith }
1280289bc588SBarry Smith 
12814a2ae208SSatish Balay #undef __FUNCT__
12824a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1283dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1284289bc588SBarry Smith {
1285c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
128687828ca2SBarry Smith   PetscScalar  *v = mat->v;
1287329f5518SBarry Smith   PetscReal    sum = 0.0;
1288d0f46423SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1289efee365bSSatish Balay   PetscErrorCode ierr;
129055659b69SBarry Smith 
12913a40ed3dSBarry Smith   PetscFunctionBegin;
1292289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1293a5ce6ee0Svictorle     if (lda>m) {
1294d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1295a5ce6ee0Svictorle 	v = mat->v+j*lda;
1296a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1297a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1298a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1299a5ce6ee0Svictorle #else
1300a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1301a5ce6ee0Svictorle #endif
1302a5ce6ee0Svictorle 	}
1303a5ce6ee0Svictorle       }
1304a5ce6ee0Svictorle     } else {
1305d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1306aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1307329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1308289bc588SBarry Smith #else
1309289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1310289bc588SBarry Smith #endif
1311289bc588SBarry Smith       }
1312a5ce6ee0Svictorle     }
13138f1a2a5eSBarry Smith     *nrm = PetscSqrtReal(sum);
1314dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
13153a40ed3dSBarry Smith   } else if (type == NORM_1) {
1316064f8208SBarry Smith     *nrm = 0.0;
1317d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
13181b807ce4Svictorle       v = mat->v + j*mat->lda;
1319289bc588SBarry Smith       sum = 0.0;
1320d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
132133a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1322289bc588SBarry Smith       }
1323064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1324289bc588SBarry Smith     }
1325d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
13263a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1327064f8208SBarry Smith     *nrm = 0.0;
1328d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1329289bc588SBarry Smith       v = mat->v + j;
1330289bc588SBarry Smith       sum = 0.0;
1331d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
13321b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1333289bc588SBarry Smith       }
1334064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1335289bc588SBarry Smith     }
1336d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1337e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
13383a40ed3dSBarry Smith   PetscFunctionReturn(0);
1339289bc588SBarry Smith }
1340289bc588SBarry Smith 
13414a2ae208SSatish Balay #undef __FUNCT__
13424a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1343ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool  flg)
1344289bc588SBarry Smith {
1345c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
134663ba0a88SBarry Smith   PetscErrorCode ierr;
134767e560aaSBarry Smith 
13483a40ed3dSBarry Smith   PetscFunctionBegin;
1349b5a2b587SKris Buschelman   switch (op) {
1350b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13514e0d8c25SBarry Smith     aij->roworiented = flg;
1352b5a2b587SKris Buschelman     break;
1353512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1354b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13553971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13564e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1357b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1358b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
135977e54ba9SKris Buschelman   case MAT_SYMMETRIC:
136077e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13619a4540c5SBarry Smith   case MAT_HERMITIAN:
13629a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1363600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1364290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
136577e54ba9SKris Buschelman     break;
1366b5a2b587SKris Buschelman   default:
1367e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13683a40ed3dSBarry Smith   }
13693a40ed3dSBarry Smith   PetscFunctionReturn(0);
1370289bc588SBarry Smith }
1371289bc588SBarry Smith 
13724a2ae208SSatish Balay #undef __FUNCT__
13734a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1374dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13756f0a148fSBarry Smith {
1376ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13776849ba73SBarry Smith   PetscErrorCode ierr;
1378d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13793a40ed3dSBarry Smith 
13803a40ed3dSBarry Smith   PetscFunctionBegin;
1381a5ce6ee0Svictorle   if (lda>m) {
1382d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1383a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1384a5ce6ee0Svictorle     }
1385a5ce6ee0Svictorle   } else {
1386d0f46423SBarry Smith     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1387a5ce6ee0Svictorle   }
13883a40ed3dSBarry Smith   PetscFunctionReturn(0);
13896f0a148fSBarry Smith }
13906f0a148fSBarry Smith 
13914a2ae208SSatish Balay #undef __FUNCT__
13924a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
13932b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
13946f0a148fSBarry Smith {
139597b48c8fSBarry Smith   PetscErrorCode    ierr;
1396ec8511deSBarry Smith   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
1397d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,j;
139897b48c8fSBarry Smith   PetscScalar       *slot,*bb;
139997b48c8fSBarry Smith   const PetscScalar *xx;
140055659b69SBarry Smith 
14013a40ed3dSBarry Smith   PetscFunctionBegin;
140297b48c8fSBarry Smith   /* fix right hand side if needed */
140397b48c8fSBarry Smith   if (x && b) {
140497b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
140597b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
140697b48c8fSBarry Smith     for (i=0; i<N; i++) {
140797b48c8fSBarry Smith       bb[rows[i]] = diag*xx[rows[i]];
140897b48c8fSBarry Smith     }
140997b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
141097b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
141197b48c8fSBarry Smith   }
141297b48c8fSBarry Smith 
14136f0a148fSBarry Smith   for (i=0; i<N; i++) {
14146f0a148fSBarry Smith     slot = l->v + rows[i];
14156f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
14166f0a148fSBarry Smith   }
1417f4df32b1SMatthew Knepley   if (diag != 0.0) {
14186f0a148fSBarry Smith     for (i=0; i<N; i++) {
14196f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1420f4df32b1SMatthew Knepley       *slot = diag;
14216f0a148fSBarry Smith     }
14226f0a148fSBarry Smith   }
14233a40ed3dSBarry Smith   PetscFunctionReturn(0);
14246f0a148fSBarry Smith }
1425557bce09SLois Curfman McInnes 
14264a2ae208SSatish Balay #undef __FUNCT__
14274a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1428dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
142964e87e97SBarry Smith {
1430c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
14313a40ed3dSBarry Smith 
14323a40ed3dSBarry Smith   PetscFunctionBegin;
1433e32f2f54SBarry Smith   if (mat->lda != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
143464e87e97SBarry Smith   *array = mat->v;
14353a40ed3dSBarry Smith   PetscFunctionReturn(0);
143664e87e97SBarry Smith }
14370754003eSLois Curfman McInnes 
14384a2ae208SSatish Balay #undef __FUNCT__
14394a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1440dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1441ff14e315SSatish Balay {
14423a40ed3dSBarry Smith   PetscFunctionBegin;
144309b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
14443a40ed3dSBarry Smith   PetscFunctionReturn(0);
1445ff14e315SSatish Balay }
14460754003eSLois Curfman McInnes 
14474a2ae208SSatish Balay #undef __FUNCT__
14484a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
144913f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
14500754003eSLois Curfman McInnes {
1451c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
14526849ba73SBarry Smith   PetscErrorCode ierr;
14535d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
14545d0c19d7SBarry Smith   const PetscInt *irow,*icol;
145587828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
14560754003eSLois Curfman McInnes   Mat            newmat;
14570754003eSLois Curfman McInnes 
14583a40ed3dSBarry Smith   PetscFunctionBegin;
145978b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
146078b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1461e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1462e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14630754003eSLois Curfman McInnes 
1464182d2002SSatish Balay   /* Check submatrixcall */
1465182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
146613f74950SBarry Smith     PetscInt n_cols,n_rows;
1467182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
146821a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
1469f746d493SDmitry Karpeev       /* resize the result matrix to match number of requested rows/columns */
1470c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
147121a2c019SBarry Smith     }
1472182d2002SSatish Balay     newmat = *B;
1473182d2002SSatish Balay   } else {
14740754003eSLois Curfman McInnes     /* Create and fill new matrix */
14757adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1476f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14777adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14785c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1479182d2002SSatish Balay   }
1480182d2002SSatish Balay 
1481182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1482182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1483182d2002SSatish Balay 
1484182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14856de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1486182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1487182d2002SSatish Balay       *bv++ = av[irow[j]];
14880754003eSLois Curfman McInnes     }
14890754003eSLois Curfman McInnes   }
1490182d2002SSatish Balay 
1491182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14926d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14936d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14940754003eSLois Curfman McInnes 
14950754003eSLois Curfman McInnes   /* Free work space */
149678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
149778b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1498182d2002SSatish Balay   *B = newmat;
14993a40ed3dSBarry Smith   PetscFunctionReturn(0);
15000754003eSLois Curfman McInnes }
15010754003eSLois Curfman McInnes 
15024a2ae208SSatish Balay #undef __FUNCT__
15034a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
150413f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1505905e6a2fSBarry Smith {
15066849ba73SBarry Smith   PetscErrorCode ierr;
150713f74950SBarry Smith   PetscInt       i;
1508905e6a2fSBarry Smith 
15093a40ed3dSBarry Smith   PetscFunctionBegin;
1510905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1511b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1512905e6a2fSBarry Smith   }
1513905e6a2fSBarry Smith 
1514905e6a2fSBarry Smith   for (i=0; i<n; i++) {
15156a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1516905e6a2fSBarry Smith   }
15173a40ed3dSBarry Smith   PetscFunctionReturn(0);
1518905e6a2fSBarry Smith }
1519905e6a2fSBarry Smith 
15204a2ae208SSatish Balay #undef __FUNCT__
1521c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1522c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1523c0aa2d19SHong Zhang {
1524c0aa2d19SHong Zhang   PetscFunctionBegin;
1525c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1526c0aa2d19SHong Zhang }
1527c0aa2d19SHong Zhang 
1528c0aa2d19SHong Zhang #undef __FUNCT__
1529c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1530c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1531c0aa2d19SHong Zhang {
1532c0aa2d19SHong Zhang   PetscFunctionBegin;
1533c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1534c0aa2d19SHong Zhang }
1535c0aa2d19SHong Zhang 
1536c0aa2d19SHong Zhang #undef __FUNCT__
15374a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1538dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
15394b0e389bSBarry Smith {
15404b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
15416849ba73SBarry Smith   PetscErrorCode ierr;
1542d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
15433a40ed3dSBarry Smith 
15443a40ed3dSBarry Smith   PetscFunctionBegin;
154533f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
154633f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1547cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
15483a40ed3dSBarry Smith     PetscFunctionReturn(0);
15493a40ed3dSBarry Smith   }
1550e32f2f54SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1551a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
15520dbb7854Svictorle     for (j=0; j<n; j++) {
1553a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1554a5ce6ee0Svictorle     }
1555a5ce6ee0Svictorle   } else {
1556d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1557a5ce6ee0Svictorle   }
1558273d9f13SBarry Smith   PetscFunctionReturn(0);
1559273d9f13SBarry Smith }
1560273d9f13SBarry Smith 
15614a2ae208SSatish Balay #undef __FUNCT__
15624a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1563dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1564273d9f13SBarry Smith {
1565dfbe8321SBarry Smith   PetscErrorCode ierr;
1566273d9f13SBarry Smith 
1567273d9f13SBarry Smith   PetscFunctionBegin;
1568273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15693a40ed3dSBarry Smith   PetscFunctionReturn(0);
15704b0e389bSBarry Smith }
15714b0e389bSBarry Smith 
1572284134d9SBarry Smith #undef __FUNCT__
1573284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1574284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1575284134d9SBarry Smith {
1576284134d9SBarry Smith   PetscFunctionBegin;
157721a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1578284134d9SBarry Smith   m = PetscMax(m,M);
1579284134d9SBarry Smith   n = PetscMax(n,N);
1580a868139aSShri Abhyankar 
158186d161a7SShri Abhyankar   /*  if (m > a->Mmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m);
158286d161a7SShri Abhyankar     if (n > a->Nmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n);
158386d161a7SShri Abhyankar   */
1584dc5cefdeSJed Brown   A->rmap->n = A->rmap->N = m;
1585d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
1586284134d9SBarry Smith   PetscFunctionReturn(0);
1587284134d9SBarry Smith }
1588170fe5c8SBarry Smith 
1589ba337c44SJed Brown #undef __FUNCT__
1590ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense"
1591ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A)
1592ba337c44SJed Brown {
1593ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1594ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1595ba337c44SJed Brown   PetscScalar    *aa = a->v;
1596ba337c44SJed Brown 
1597ba337c44SJed Brown   PetscFunctionBegin;
1598ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
1599ba337c44SJed Brown   PetscFunctionReturn(0);
1600ba337c44SJed Brown }
1601ba337c44SJed Brown 
1602ba337c44SJed Brown #undef __FUNCT__
1603ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense"
1604ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A)
1605ba337c44SJed Brown {
1606ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1607ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1608ba337c44SJed Brown   PetscScalar    *aa = a->v;
1609ba337c44SJed Brown 
1610ba337c44SJed Brown   PetscFunctionBegin;
1611ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1612ba337c44SJed Brown   PetscFunctionReturn(0);
1613ba337c44SJed Brown }
1614ba337c44SJed Brown 
1615ba337c44SJed Brown #undef __FUNCT__
1616ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense"
1617ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
1618ba337c44SJed Brown {
1619ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1620ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1621ba337c44SJed Brown   PetscScalar    *aa = a->v;
1622ba337c44SJed Brown 
1623ba337c44SJed Brown   PetscFunctionBegin;
1624ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1625ba337c44SJed Brown   PetscFunctionReturn(0);
1626ba337c44SJed Brown }
1627284134d9SBarry Smith 
1628a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1629a9fe9ddaSSatish Balay #undef __FUNCT__
1630a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1631a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1632a9fe9ddaSSatish Balay {
1633a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1634a9fe9ddaSSatish Balay 
1635a9fe9ddaSSatish Balay   PetscFunctionBegin;
1636a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1637a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1638a9fe9ddaSSatish Balay   }
1639a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1640a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1641a9fe9ddaSSatish Balay }
1642a9fe9ddaSSatish Balay 
1643a9fe9ddaSSatish Balay #undef __FUNCT__
1644a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1645a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1646a9fe9ddaSSatish Balay {
1647ee16a9a1SHong Zhang   PetscErrorCode ierr;
1648d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1649ee16a9a1SHong Zhang   Mat            Cmat;
1650a9fe9ddaSSatish Balay 
1651ee16a9a1SHong Zhang   PetscFunctionBegin;
1652e32f2f54SBarry Smith   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n);
1653ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1654ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1655ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1656ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1657ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1658ee16a9a1SHong Zhang   *C = Cmat;
1659ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1660ee16a9a1SHong Zhang }
1661a9fe9ddaSSatish Balay 
166298a3b096SSatish Balay #undef __FUNCT__
1663a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1664a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1665a9fe9ddaSSatish Balay {
1666a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1667a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1668a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16690805154bSBarry Smith   PetscBLASInt   m,n,k;
1670a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1671a9fe9ddaSSatish Balay 
1672a9fe9ddaSSatish Balay   PetscFunctionBegin;
1673d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1674d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1675d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1676a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1677a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1678a9fe9ddaSSatish Balay }
1679a9fe9ddaSSatish Balay 
1680a9fe9ddaSSatish Balay #undef __FUNCT__
1681a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1682a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1683a9fe9ddaSSatish Balay {
1684a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1685a9fe9ddaSSatish Balay 
1686a9fe9ddaSSatish Balay   PetscFunctionBegin;
1687a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1688a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1689a9fe9ddaSSatish Balay   }
1690a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1691a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1692a9fe9ddaSSatish Balay }
1693a9fe9ddaSSatish Balay 
1694a9fe9ddaSSatish Balay #undef __FUNCT__
1695a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1696a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1697a9fe9ddaSSatish Balay {
1698ee16a9a1SHong Zhang   PetscErrorCode ierr;
1699d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1700ee16a9a1SHong Zhang   Mat            Cmat;
1701a9fe9ddaSSatish Balay 
1702ee16a9a1SHong Zhang   PetscFunctionBegin;
1703e32f2f54SBarry Smith   if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n);
1704ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1705ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1706ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1707ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1708ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1709ee16a9a1SHong Zhang   *C = Cmat;
1710ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1711ee16a9a1SHong Zhang }
1712a9fe9ddaSSatish Balay 
1713a9fe9ddaSSatish Balay #undef __FUNCT__
1714a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1715a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1716a9fe9ddaSSatish Balay {
1717a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1718a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1719a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
17200805154bSBarry Smith   PetscBLASInt   m,n,k;
1721a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1722a9fe9ddaSSatish Balay 
1723a9fe9ddaSSatish Balay   PetscFunctionBegin;
1724d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1725d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1726d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
17272fbe02b9SBarry Smith   /*
17282fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
17292fbe02b9SBarry Smith   */
1730a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1731a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1732a9fe9ddaSSatish Balay }
1733985db425SBarry Smith 
1734985db425SBarry Smith #undef __FUNCT__
1735985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1736985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1737985db425SBarry Smith {
1738985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1739985db425SBarry Smith   PetscErrorCode ierr;
1740d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1741985db425SBarry Smith   PetscScalar    *x;
1742985db425SBarry Smith   MatScalar      *aa = a->v;
1743985db425SBarry Smith 
1744985db425SBarry Smith   PetscFunctionBegin;
1745e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1746985db425SBarry Smith 
1747985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1748985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1749985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1750e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1751985db425SBarry Smith   for (i=0; i<m; i++) {
1752985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1753985db425SBarry Smith     for (j=1; j<n; j++){
1754985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1755985db425SBarry Smith     }
1756985db425SBarry Smith   }
1757985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1758985db425SBarry Smith   PetscFunctionReturn(0);
1759985db425SBarry Smith }
1760985db425SBarry Smith 
1761985db425SBarry Smith #undef __FUNCT__
1762985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1763985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1764985db425SBarry Smith {
1765985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1766985db425SBarry Smith   PetscErrorCode ierr;
1767d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1768985db425SBarry Smith   PetscScalar    *x;
1769985db425SBarry Smith   PetscReal      atmp;
1770985db425SBarry Smith   MatScalar      *aa = a->v;
1771985db425SBarry Smith 
1772985db425SBarry Smith   PetscFunctionBegin;
1773e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1774985db425SBarry Smith 
1775985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1776985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1777985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1778e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1779985db425SBarry Smith   for (i=0; i<m; i++) {
17809189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1781985db425SBarry Smith     for (j=1; j<n; j++){
1782985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1783985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1784985db425SBarry Smith     }
1785985db425SBarry Smith   }
1786985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1787985db425SBarry Smith   PetscFunctionReturn(0);
1788985db425SBarry Smith }
1789985db425SBarry Smith 
1790985db425SBarry Smith #undef __FUNCT__
1791985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1792985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1793985db425SBarry Smith {
1794985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1795985db425SBarry Smith   PetscErrorCode ierr;
1796d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1797985db425SBarry Smith   PetscScalar    *x;
1798985db425SBarry Smith   MatScalar      *aa = a->v;
1799985db425SBarry Smith 
1800985db425SBarry Smith   PetscFunctionBegin;
1801e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1802985db425SBarry Smith 
1803985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1804985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1805985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1806e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1807985db425SBarry Smith   for (i=0; i<m; i++) {
1808985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1809985db425SBarry Smith     for (j=1; j<n; j++){
1810985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1811985db425SBarry Smith     }
1812985db425SBarry Smith   }
1813985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1814985db425SBarry Smith   PetscFunctionReturn(0);
1815985db425SBarry Smith }
1816985db425SBarry Smith 
18178d0534beSBarry Smith #undef __FUNCT__
18188d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
18198d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
18208d0534beSBarry Smith {
18218d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
18228d0534beSBarry Smith   PetscErrorCode ierr;
18238d0534beSBarry Smith   PetscScalar    *x;
18248d0534beSBarry Smith 
18258d0534beSBarry Smith   PetscFunctionBegin;
1826e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
18278d0534beSBarry Smith 
18288d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1829d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
18308d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
18318d0534beSBarry Smith   PetscFunctionReturn(0);
18328d0534beSBarry Smith }
18338d0534beSBarry Smith 
18340716a85fSBarry Smith 
18350716a85fSBarry Smith #undef __FUNCT__
18360716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense"
18370716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
18380716a85fSBarry Smith {
18390716a85fSBarry Smith   PetscErrorCode ierr;
18400716a85fSBarry Smith   PetscInt       i,j,m,n;
18410716a85fSBarry Smith   PetscScalar    *a;
18420716a85fSBarry Smith 
18430716a85fSBarry Smith   PetscFunctionBegin;
18440716a85fSBarry Smith   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
18450716a85fSBarry Smith   ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr);
18460716a85fSBarry Smith   ierr = MatGetArray(A,&a);CHKERRQ(ierr);
18470716a85fSBarry Smith   if (type == NORM_2) {
18480716a85fSBarry Smith     for (i=0; i<n; i++ ){
18490716a85fSBarry Smith       for (j=0; j<m; j++) {
18500716a85fSBarry Smith 	norms[i] += PetscAbsScalar(a[j]*a[j]);
18510716a85fSBarry Smith       }
18520716a85fSBarry Smith       a += m;
18530716a85fSBarry Smith     }
18540716a85fSBarry Smith   } else if (type == NORM_1) {
18550716a85fSBarry Smith     for (i=0; i<n; i++ ){
18560716a85fSBarry Smith       for (j=0; j<m; j++) {
18570716a85fSBarry Smith 	norms[i] += PetscAbsScalar(a[j]);
18580716a85fSBarry Smith       }
18590716a85fSBarry Smith       a += m;
18600716a85fSBarry Smith     }
18610716a85fSBarry Smith   } else if (type == NORM_INFINITY) {
18620716a85fSBarry Smith     for (i=0; i<n; i++ ){
18630716a85fSBarry Smith       for (j=0; j<m; j++) {
18640716a85fSBarry Smith 	norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
18650716a85fSBarry Smith       }
18660716a85fSBarry Smith       a += m;
18670716a85fSBarry Smith     }
18680716a85fSBarry Smith   } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType");
18690716a85fSBarry Smith   if (type == NORM_2) {
18708f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
18710716a85fSBarry Smith   }
18720716a85fSBarry Smith   PetscFunctionReturn(0);
18730716a85fSBarry Smith }
18740716a85fSBarry Smith 
1875289bc588SBarry Smith /* -------------------------------------------------------------------*/
1876a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1877905e6a2fSBarry Smith        MatGetRow_SeqDense,
1878905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1879905e6a2fSBarry Smith        MatMult_SeqDense,
188097304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
18817c922b88SBarry Smith        MatMultTranspose_SeqDense,
18827c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1883db4efbfdSBarry Smith        0,
1884db4efbfdSBarry Smith        0,
1885db4efbfdSBarry Smith        0,
1886db4efbfdSBarry Smith /*10*/ 0,
1887905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1888905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
188941f059aeSBarry Smith        MatSOR_SeqDense,
1890ec8511deSBarry Smith        MatTranspose_SeqDense,
189197304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1892905e6a2fSBarry Smith        MatEqual_SeqDense,
1893905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1894905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1895905e6a2fSBarry Smith        MatNorm_SeqDense,
1896c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1897c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1898905e6a2fSBarry Smith        MatSetOption_SeqDense,
1899905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1900d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense,
1901db4efbfdSBarry Smith        0,
1902db4efbfdSBarry Smith        0,
1903db4efbfdSBarry Smith        0,
1904db4efbfdSBarry Smith        0,
1905d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense,
1906273d9f13SBarry Smith        0,
1907905e6a2fSBarry Smith        0,
1908905e6a2fSBarry Smith        MatGetArray_SeqDense,
1909905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1910d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense,
1911a5ae1ecdSBarry Smith        0,
1912a5ae1ecdSBarry Smith        0,
1913a5ae1ecdSBarry Smith        0,
1914a5ae1ecdSBarry Smith        0,
1915d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense,
1916a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1917a5ae1ecdSBarry Smith        0,
19184b0e389bSBarry Smith        MatGetValues_SeqDense,
1919a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1920d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense,
1921a5ae1ecdSBarry Smith        MatScale_SeqDense,
1922a5ae1ecdSBarry Smith        0,
1923a5ae1ecdSBarry Smith        0,
1924a5ae1ecdSBarry Smith        0,
1925d519adbfSMatthew Knepley /*49*/ 0,
1926a5ae1ecdSBarry Smith        0,
1927a5ae1ecdSBarry Smith        0,
1928a5ae1ecdSBarry Smith        0,
1929a5ae1ecdSBarry Smith        0,
1930d519adbfSMatthew Knepley /*54*/ 0,
1931a5ae1ecdSBarry Smith        0,
1932a5ae1ecdSBarry Smith        0,
1933a5ae1ecdSBarry Smith        0,
1934a5ae1ecdSBarry Smith        0,
1935d519adbfSMatthew Knepley /*59*/ 0,
1936e03a110bSBarry Smith        MatDestroy_SeqDense,
1937e03a110bSBarry Smith        MatView_SeqDense,
1938357abbc8SBarry Smith        0,
193997304618SKris Buschelman        0,
1940d519adbfSMatthew Knepley /*64*/ 0,
194197304618SKris Buschelman        0,
194297304618SKris Buschelman        0,
194397304618SKris Buschelman        0,
194497304618SKris Buschelman        0,
1945d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense,
194697304618SKris Buschelman        0,
194797304618SKris Buschelman        0,
194897304618SKris Buschelman        0,
194997304618SKris Buschelman        0,
1950d519adbfSMatthew Knepley /*74*/ 0,
195197304618SKris Buschelman        0,
195297304618SKris Buschelman        0,
195397304618SKris Buschelman        0,
195497304618SKris Buschelman        0,
1955d519adbfSMatthew Knepley /*79*/ 0,
195697304618SKris Buschelman        0,
195797304618SKris Buschelman        0,
195897304618SKris Buschelman        0,
19595bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense,
1960865e5f61SKris Buschelman        0,
19611cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1962865e5f61SKris Buschelman        0,
1963865e5f61SKris Buschelman        0,
1964865e5f61SKris Buschelman        0,
1965d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense,
1966a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1967a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1968865e5f61SKris Buschelman        0,
1969865e5f61SKris Buschelman        0,
1970d519adbfSMatthew Knepley /*94*/ 0,
1971a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1972a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1973a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1974284134d9SBarry Smith        0,
1975d519adbfSMatthew Knepley /*99*/ 0,
1976284134d9SBarry Smith        0,
1977284134d9SBarry Smith        0,
1978ba337c44SJed Brown        MatConjugate_SeqDense,
1979985db425SBarry Smith        MatSetSizes_SeqDense,
1980ba337c44SJed Brown /*104*/0,
1981ba337c44SJed Brown        MatRealPart_SeqDense,
1982ba337c44SJed Brown        MatImaginaryPart_SeqDense,
1983985db425SBarry Smith        0,
1984985db425SBarry Smith        0,
198585e2c93fSHong Zhang /*109*/MatMatSolve_SeqDense,
1986985db425SBarry Smith        0,
19878d0534beSBarry Smith        MatGetRowMin_SeqDense,
1988aabbc4fbSShri Abhyankar        MatGetColumnVector_SeqDense,
1989aabbc4fbSShri Abhyankar        0,
1990aabbc4fbSShri Abhyankar /*114*/0,
1991aabbc4fbSShri Abhyankar        0,
1992aabbc4fbSShri Abhyankar        0,
1993aabbc4fbSShri Abhyankar        0,
1994aabbc4fbSShri Abhyankar        0,
1995aabbc4fbSShri Abhyankar /*119*/0,
1996aabbc4fbSShri Abhyankar        0,
1997aabbc4fbSShri Abhyankar        0,
19980716a85fSBarry Smith        0,
19990716a85fSBarry Smith        0,
20000716a85fSBarry Smith /*124*/0,
20010716a85fSBarry Smith        MatGetColumnNorms_SeqDense
2002985db425SBarry Smith };
200390ace30eSBarry Smith 
20044a2ae208SSatish Balay #undef __FUNCT__
20054a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
20064b828684SBarry Smith /*@C
2007fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
2008d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
2009d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
2010289bc588SBarry Smith 
2011db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
2012db81eaa0SLois Curfman McInnes 
201320563c6bSBarry Smith    Input Parameters:
2014db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
20150c775827SLois Curfman McInnes .  m - number of rows
201618f449edSLois Curfman McInnes .  n - number of columns
2017c0235b3cSMatthew Knepley -  data - optional location of matrix data in column major order.  Set data=PETSC_NULL for PETSc
2018dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
201920563c6bSBarry Smith 
202020563c6bSBarry Smith    Output Parameter:
202144cd7ae7SLois Curfman McInnes .  A - the matrix
202220563c6bSBarry Smith 
2023b259b22eSLois Curfman McInnes    Notes:
202418f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
202518f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
2026b4fd4287SBarry Smith    set data=PETSC_NULL.
202718f449edSLois Curfman McInnes 
2028027ccd11SLois Curfman McInnes    Level: intermediate
2029027ccd11SLois Curfman McInnes 
2030dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
2031d65003e9SLois Curfman McInnes 
2032db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
203320563c6bSBarry Smith @*/
20347087cfbeSBarry Smith PetscErrorCode  MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2035289bc588SBarry Smith {
2036dfbe8321SBarry Smith   PetscErrorCode ierr;
20373b2fbd54SBarry Smith 
20383a40ed3dSBarry Smith   PetscFunctionBegin;
2039f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2040f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
2041273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
2042273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
2043273d9f13SBarry Smith   PetscFunctionReturn(0);
2044273d9f13SBarry Smith }
2045273d9f13SBarry Smith 
20464a2ae208SSatish Balay #undef __FUNCT__
2047afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
2048273d9f13SBarry Smith /*@C
2049273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
2050273d9f13SBarry Smith 
2051273d9f13SBarry Smith    Collective on MPI_Comm
2052273d9f13SBarry Smith 
2053273d9f13SBarry Smith    Input Parameters:
2054273d9f13SBarry Smith +  A - the matrix
2055273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
2056273d9f13SBarry Smith 
2057273d9f13SBarry Smith    Notes:
2058273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
2059273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
2060284134d9SBarry Smith    need not call this routine.
2061273d9f13SBarry Smith 
2062273d9f13SBarry Smith    Level: intermediate
2063273d9f13SBarry Smith 
2064273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
2065273d9f13SBarry Smith 
2066867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
2067867c911aSBarry Smith 
2068273d9f13SBarry Smith @*/
20697087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2070273d9f13SBarry Smith {
20714ac538c5SBarry Smith   PetscErrorCode ierr;
2072a23d5eceSKris Buschelman 
2073a23d5eceSKris Buschelman   PetscFunctionBegin;
20744ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr);
2075a23d5eceSKris Buschelman   PetscFunctionReturn(0);
2076a23d5eceSKris Buschelman }
2077a23d5eceSKris Buschelman 
2078a23d5eceSKris Buschelman EXTERN_C_BEGIN
2079a23d5eceSKris Buschelman #undef __FUNCT__
2080afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
20817087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2082a23d5eceSKris Buschelman {
2083273d9f13SBarry Smith   Mat_SeqDense   *b;
2084dfbe8321SBarry Smith   PetscErrorCode ierr;
2085273d9f13SBarry Smith 
2086273d9f13SBarry Smith   PetscFunctionBegin;
2087273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
2088a868139aSShri Abhyankar 
208934ef9618SShri Abhyankar   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
209034ef9618SShri Abhyankar   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
209134ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
209234ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
209334ef9618SShri Abhyankar 
2094273d9f13SBarry Smith   b       = (Mat_SeqDense*)B->data;
209586d161a7SShri Abhyankar   b->Mmax = B->rmap->n;
209686d161a7SShri Abhyankar   b->Nmax = B->cmap->n;
209786d161a7SShri Abhyankar   if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n;
209886d161a7SShri Abhyankar 
20999e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
21009e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
21015afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
2102284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
2103284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
21049e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
2105273d9f13SBarry Smith   } else { /* user-allocated storage */
21069e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
2107273d9f13SBarry Smith     b->v          = data;
2108273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
2109273d9f13SBarry Smith   }
21100450473dSBarry Smith   B->assembled = PETSC_TRUE;
2111273d9f13SBarry Smith   PetscFunctionReturn(0);
2112273d9f13SBarry Smith }
2113a23d5eceSKris Buschelman EXTERN_C_END
2114273d9f13SBarry Smith 
21151b807ce4Svictorle #undef __FUNCT__
21161b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
21171b807ce4Svictorle /*@C
21181b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
21191b807ce4Svictorle 
21201b807ce4Svictorle   Input parameter:
21211b807ce4Svictorle + A - the matrix
21221b807ce4Svictorle - lda - the leading dimension
21231b807ce4Svictorle 
21241b807ce4Svictorle   Notes:
2125867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
21261b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
21271b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
21281b807ce4Svictorle 
21291b807ce4Svictorle   Level: intermediate
21301b807ce4Svictorle 
21311b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
21321b807ce4Svictorle 
2133284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2134867c911aSBarry Smith 
21351b807ce4Svictorle @*/
21367087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetLDA(Mat B,PetscInt lda)
21371b807ce4Svictorle {
21381b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
213921a2c019SBarry Smith 
21401b807ce4Svictorle   PetscFunctionBegin;
2141e32f2f54SBarry Smith   if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
21421b807ce4Svictorle   b->lda       = lda;
214321a2c019SBarry Smith   b->changelda = PETSC_FALSE;
214421a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
21451b807ce4Svictorle   PetscFunctionReturn(0);
21461b807ce4Svictorle }
21471b807ce4Svictorle 
21480bad9183SKris Buschelman /*MC
2149fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
21500bad9183SKris Buschelman 
21510bad9183SKris Buschelman    Options Database Keys:
21520bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
21530bad9183SKris Buschelman 
21540bad9183SKris Buschelman   Level: beginner
21550bad9183SKris Buschelman 
215689665df3SBarry Smith .seealso: MatCreateSeqDense()
215789665df3SBarry Smith 
21580bad9183SKris Buschelman M*/
21590bad9183SKris Buschelman 
2160273d9f13SBarry Smith EXTERN_C_BEGIN
21614a2ae208SSatish Balay #undef __FUNCT__
21624a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
21637087cfbeSBarry Smith PetscErrorCode  MatCreate_SeqDense(Mat B)
2164273d9f13SBarry Smith {
2165273d9f13SBarry Smith   Mat_SeqDense   *b;
2166dfbe8321SBarry Smith   PetscErrorCode ierr;
21677c334f02SBarry Smith   PetscMPIInt    size;
2168273d9f13SBarry Smith 
2169273d9f13SBarry Smith   PetscFunctionBegin;
21707adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
2171e32f2f54SBarry Smith   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
217255659b69SBarry Smith 
217338f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2174549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
217544cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
217618f449edSLois Curfman McInnes 
217744cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2178273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2179273d9f13SBarry Smith   b->v            = 0;
218021a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
21814e220ebcSLois Curfman McInnes 
2182b24902e0SBarry Smith 
2183ec1065edSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
2184b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2185b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2186a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2187a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2188a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
21894ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
21904ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
21914ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
21924ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
21934ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
21944ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
21954ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
21964ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
21974ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
219817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
21993a40ed3dSBarry Smith   PetscFunctionReturn(0);
2200289bc588SBarry Smith }
2201273d9f13SBarry Smith EXTERN_C_END
2202