xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 0716a85f2ec13166f862a3128e4c2e598c442603)
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;
109719d5645SBarry Smith   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
110b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
111b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
112d0f46423SBarry Smith     if (lda>A->rmap->n) {
113d0f46423SBarry Smith       m = A->rmap->n;
114d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
115b24902e0SBarry Smith 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
116b24902e0SBarry Smith       }
117b24902e0SBarry Smith     } else {
118d0f46423SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
119b24902e0SBarry Smith     }
120b24902e0SBarry Smith   }
121b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
122b24902e0SBarry Smith   PetscFunctionReturn(0);
123b24902e0SBarry Smith }
124b24902e0SBarry Smith 
1254a2ae208SSatish Balay #undef __FUNCT__
1264a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
127dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12802cad45dSBarry Smith {
1296849ba73SBarry Smith   PetscErrorCode ierr;
13002cad45dSBarry Smith 
1313a40ed3dSBarry Smith   PetscFunctionBegin;
1325c9eb25fSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
133d0f46423SBarry Smith   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1345c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
135719d5645SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
136b24902e0SBarry Smith   PetscFunctionReturn(0);
137b24902e0SBarry Smith }
138b24902e0SBarry Smith 
1396ee01492SSatish Balay 
1400481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
141719d5645SBarry Smith 
1424a2ae208SSatish Balay #undef __FUNCT__
1434a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1440481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
145289bc588SBarry Smith {
1464482741eSBarry Smith   MatFactorInfo  info;
147a093e273SMatthew Knepley   PetscErrorCode ierr;
1483a40ed3dSBarry Smith 
1493a40ed3dSBarry Smith   PetscFunctionBegin;
150c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
151719d5645SBarry Smith   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
1523a40ed3dSBarry Smith   PetscFunctionReturn(0);
153289bc588SBarry Smith }
1546ee01492SSatish Balay 
1550b4b3355SBarry Smith #undef __FUNCT__
1564a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
157dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
158289bc588SBarry Smith {
159c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1606849ba73SBarry Smith   PetscErrorCode ierr;
16187828ca2SBarry Smith   PetscScalar    *x,*y;
162d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
16367e560aaSBarry Smith 
1643a40ed3dSBarry Smith   PetscFunctionBegin;
1651ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1661ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
167d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
168d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
169ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
170e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
171ae7cfcebSSatish Balay #else
17271044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
173e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
174ae7cfcebSSatish Balay #endif
175d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_CHOLESKY){
176ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
177e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
178ae7cfcebSSatish Balay #else
17971044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
180e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
181ae7cfcebSSatish Balay #endif
182289bc588SBarry Smith   }
183e32f2f54SBarry Smith   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
1841ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1851ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
186dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
1873a40ed3dSBarry Smith   PetscFunctionReturn(0);
188289bc588SBarry Smith }
1896ee01492SSatish Balay 
1904a2ae208SSatish Balay #undef __FUNCT__
19185e2c93fSHong Zhang #define __FUNCT__ "MatMatSolve_SeqDense"
19285e2c93fSHong Zhang PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
19385e2c93fSHong Zhang {
19485e2c93fSHong Zhang   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
19585e2c93fSHong Zhang   PetscErrorCode ierr;
19685e2c93fSHong Zhang   PetscScalar    *b,*x;
19785e2c93fSHong Zhang   PetscBLASInt   nrhs,info,m=PetscBLASIntCast(A->rmap->n);
19885e2c93fSHong Zhang 
19985e2c93fSHong Zhang   PetscFunctionBegin;
20085e2c93fSHong Zhang   ierr = MatGetSize(B,PETSC_NULL,&nrhs);CHKERRQ(ierr);
20185e2c93fSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
20285e2c93fSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
20385e2c93fSHong Zhang 
20485e2c93fSHong Zhang   ierr = PetscMemcpy(x,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
20585e2c93fSHong Zhang 
20685e2c93fSHong Zhang   if (A->factortype == MAT_FACTOR_LU) {
20785e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS)
20885e2c93fSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
20985e2c93fSHong Zhang #else
21085e2c93fSHong Zhang     LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info);
21185e2c93fSHong Zhang     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
21285e2c93fSHong Zhang #endif
21385e2c93fSHong Zhang   } else if (A->factortype == MAT_FACTOR_CHOLESKY){
21485e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS)
21585e2c93fSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
21685e2c93fSHong Zhang #else
21785e2c93fSHong Zhang     LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info);
21885e2c93fSHong Zhang     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
21985e2c93fSHong Zhang #endif
22085e2c93fSHong Zhang   }
22185e2c93fSHong Zhang   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
22285e2c93fSHong Zhang 
22385e2c93fSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
22485e2c93fSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
22585e2c93fSHong Zhang   ierr = PetscLogFlops(nrhs*(2.0*m*m- m));CHKERRQ(ierr);
22685e2c93fSHong Zhang   PetscFunctionReturn(0);
22785e2c93fSHong Zhang }
22885e2c93fSHong Zhang 
22985e2c93fSHong Zhang #undef __FUNCT__
2304a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
231dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
232da3a660dSBarry Smith {
233c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
234dfbe8321SBarry Smith   PetscErrorCode ierr;
23587828ca2SBarry Smith   PetscScalar    *x,*y;
236d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
23767e560aaSBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
2391ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2401ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
241d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
242752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
243da3a660dSBarry Smith   if (mat->pivots) {
244ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
245e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
246ae7cfcebSSatish Balay #else
24771044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
248e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
249ae7cfcebSSatish Balay #endif
2507a97a34bSBarry Smith   } else {
251ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
252e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
253ae7cfcebSSatish Balay #else
25471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
255e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
256ae7cfcebSSatish Balay #endif
257da3a660dSBarry Smith   }
2581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2591ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
260dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
262da3a660dSBarry Smith }
2636ee01492SSatish Balay 
2644a2ae208SSatish Balay #undef __FUNCT__
2654a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
266dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
267da3a660dSBarry Smith {
268c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
269dfbe8321SBarry Smith   PetscErrorCode ierr;
27087828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
271da3a660dSBarry Smith   Vec            tmp = 0;
272d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
27367e560aaSBarry Smith 
2743a40ed3dSBarry Smith   PetscFunctionBegin;
2751ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2761ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
277d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
278da3a660dSBarry Smith   if (yy == zz) {
27978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
28052e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
28178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
282da3a660dSBarry Smith   }
283d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
284752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
285da3a660dSBarry Smith   if (mat->pivots) {
286ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
287e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
288ae7cfcebSSatish Balay #else
28971044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
290e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
291ae7cfcebSSatish Balay #endif
292a8c6a408SBarry Smith   } else {
293ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
294e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
295ae7cfcebSSatish Balay #else
29671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
297e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
298ae7cfcebSSatish Balay #endif
299da3a660dSBarry Smith   }
3006bf464f9SBarry Smith   if (tmp) {
3016bf464f9SBarry Smith     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
3026bf464f9SBarry Smith     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
3036bf464f9SBarry Smith   } else {
3046bf464f9SBarry Smith     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
3056bf464f9SBarry Smith   }
3061ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3071ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
308dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
3093a40ed3dSBarry Smith   PetscFunctionReturn(0);
310da3a660dSBarry Smith }
31167e560aaSBarry Smith 
3124a2ae208SSatish Balay #undef __FUNCT__
3134a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
314dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
315da3a660dSBarry Smith {
316c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3176849ba73SBarry Smith   PetscErrorCode ierr;
31887828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
319da3a660dSBarry Smith   Vec            tmp;
320d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
32167e560aaSBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
323d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
3241ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3251ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
326da3a660dSBarry Smith   if (yy == zz) {
32778b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
32852e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
32978b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
330da3a660dSBarry Smith   }
331d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
332752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
333da3a660dSBarry Smith   if (mat->pivots) {
334ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
335e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
336ae7cfcebSSatish Balay #else
33771044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
338e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
339ae7cfcebSSatish Balay #endif
3403a40ed3dSBarry Smith   } else {
341ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
342e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
343ae7cfcebSSatish Balay #else
34471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
345e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
346ae7cfcebSSatish Balay #endif
347da3a660dSBarry Smith   }
34890f02eecSBarry Smith   if (tmp) {
3492dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
3506bf464f9SBarry Smith     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
3513a40ed3dSBarry Smith   } else {
3522dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
35390f02eecSBarry Smith   }
3541ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3551ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
356dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
3573a40ed3dSBarry Smith   PetscFunctionReturn(0);
358da3a660dSBarry Smith }
359db4efbfdSBarry Smith 
360db4efbfdSBarry Smith /* ---------------------------------------------------------------*/
361db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
362db4efbfdSBarry Smith    rather than put it in the Mat->row slot.*/
363db4efbfdSBarry Smith #undef __FUNCT__
364db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense"
3650481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
366db4efbfdSBarry Smith {
367db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
368db4efbfdSBarry Smith   PetscFunctionBegin;
369e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
370db4efbfdSBarry Smith #else
371db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
372db4efbfdSBarry Smith   PetscErrorCode ierr;
373db4efbfdSBarry Smith   PetscBLASInt   n,m,info;
374db4efbfdSBarry Smith 
375db4efbfdSBarry Smith   PetscFunctionBegin;
376db4efbfdSBarry Smith   n = PetscBLASIntCast(A->cmap->n);
377db4efbfdSBarry Smith   m = PetscBLASIntCast(A->rmap->n);
378db4efbfdSBarry Smith   if (!mat->pivots) {
379db4efbfdSBarry Smith     ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
380db4efbfdSBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
381db4efbfdSBarry Smith   }
382db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
383db4efbfdSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
384e32f2f54SBarry Smith   if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
385e32f2f54SBarry Smith   if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
386db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
387db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
388db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
389db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
390d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_LU;
391db4efbfdSBarry Smith 
392dc0b31edSSatish Balay   ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
393db4efbfdSBarry Smith #endif
394db4efbfdSBarry Smith   PetscFunctionReturn(0);
395db4efbfdSBarry Smith }
396db4efbfdSBarry Smith 
397db4efbfdSBarry Smith #undef __FUNCT__
398db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense"
3990481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
400db4efbfdSBarry Smith {
401db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
402db4efbfdSBarry Smith   PetscFunctionBegin;
403e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
404db4efbfdSBarry Smith #else
405db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
406db4efbfdSBarry Smith   PetscErrorCode ierr;
407db4efbfdSBarry Smith   PetscBLASInt   info,n = PetscBLASIntCast(A->cmap->n);
408db4efbfdSBarry Smith 
409db4efbfdSBarry Smith   PetscFunctionBegin;
410db4efbfdSBarry Smith   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
411db4efbfdSBarry Smith 
412db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
413db4efbfdSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
414e32f2f54SBarry Smith   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
415db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
416db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
417db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
418db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
419d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_CHOLESKY;
420dc0b31edSSatish Balay   ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr);
421db4efbfdSBarry Smith #endif
422db4efbfdSBarry Smith   PetscFunctionReturn(0);
423db4efbfdSBarry Smith }
424db4efbfdSBarry Smith 
425db4efbfdSBarry Smith 
426db4efbfdSBarry Smith #undef __FUNCT__
427db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
4280481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
429db4efbfdSBarry Smith {
430db4efbfdSBarry Smith   PetscErrorCode ierr;
431db4efbfdSBarry Smith   MatFactorInfo  info;
432db4efbfdSBarry Smith 
433db4efbfdSBarry Smith   PetscFunctionBegin;
434db4efbfdSBarry Smith   info.fill = 1.0;
435c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
436719d5645SBarry Smith   ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr);
437db4efbfdSBarry Smith   PetscFunctionReturn(0);
438db4efbfdSBarry Smith }
439db4efbfdSBarry Smith 
440db4efbfdSBarry Smith #undef __FUNCT__
441db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
4420481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
443db4efbfdSBarry Smith {
444db4efbfdSBarry Smith   PetscFunctionBegin;
445c3ef05f6SHong Zhang   fact->assembled                  = PETSC_TRUE;
446719d5645SBarry Smith   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
447db4efbfdSBarry Smith   PetscFunctionReturn(0);
448db4efbfdSBarry Smith }
449db4efbfdSBarry Smith 
450db4efbfdSBarry Smith #undef __FUNCT__
451db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
4520481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
453db4efbfdSBarry Smith {
454db4efbfdSBarry Smith   PetscFunctionBegin;
455c3ef05f6SHong Zhang   fact->assembled            = PETSC_TRUE;
456719d5645SBarry Smith   fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
457db4efbfdSBarry Smith   PetscFunctionReturn(0);
458db4efbfdSBarry Smith }
459db4efbfdSBarry Smith 
460bb5747d9SMatthew Knepley EXTERN_C_BEGIN
461db4efbfdSBarry Smith #undef __FUNCT__
462db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc"
463db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
464db4efbfdSBarry Smith {
465db4efbfdSBarry Smith   PetscErrorCode ierr;
466db4efbfdSBarry Smith 
467db4efbfdSBarry Smith   PetscFunctionBegin;
468db4efbfdSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
469db4efbfdSBarry Smith   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
470db4efbfdSBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
471db4efbfdSBarry Smith   if (ftype == MAT_FACTOR_LU){
472db4efbfdSBarry Smith     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
473db4efbfdSBarry Smith   } else {
474db4efbfdSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
475db4efbfdSBarry Smith   }
476d5f3da31SBarry Smith   (*fact)->factortype = ftype;
477db4efbfdSBarry Smith   PetscFunctionReturn(0);
478db4efbfdSBarry Smith }
479bb5747d9SMatthew Knepley EXTERN_C_END
480db4efbfdSBarry Smith 
481289bc588SBarry Smith /* ------------------------------------------------------------------*/
4824a2ae208SSatish Balay #undef __FUNCT__
48341f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense"
48441f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
485289bc588SBarry Smith {
486c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
48787828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
488dfbe8321SBarry Smith   PetscErrorCode ierr;
489d0f46423SBarry Smith   PetscInt       m = A->rmap->n,i;
490aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4910805154bSBarry Smith   PetscBLASInt   o = 1,bm = PetscBLASIntCast(m);
492bc1b551cSSatish Balay #endif
493289bc588SBarry Smith 
4943a40ed3dSBarry Smith   PetscFunctionBegin;
495289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
49671044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4972dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
498289bc588SBarry Smith   }
4991ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5001ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
501b965ef7fSBarry Smith   its  = its*lits;
502e32f2f54SBarry 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);
503289bc588SBarry Smith   while (its--) {
504fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
505289bc588SBarry Smith       for (i=0; i<m; i++) {
506aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
507f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
508f1747703SBarry Smith            not happy about returning a double complex */
50913f74950SBarry Smith         PetscInt    _i;
51087828ca2SBarry Smith         PetscScalar sum = b[i];
511f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
5123f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
513f1747703SBarry Smith         }
514f1747703SBarry Smith         xt = sum;
515f1747703SBarry Smith #else
51671044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
517f1747703SBarry Smith #endif
51855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
519289bc588SBarry Smith       }
520289bc588SBarry Smith     }
521fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
522289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
523aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
524f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
525f1747703SBarry Smith            not happy about returning a double complex */
52613f74950SBarry Smith         PetscInt    _i;
52787828ca2SBarry Smith         PetscScalar sum = b[i];
528f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
5293f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
530f1747703SBarry Smith         }
531f1747703SBarry Smith         xt = sum;
532f1747703SBarry Smith #else
53371044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
534f1747703SBarry Smith #endif
53555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
536289bc588SBarry Smith       }
537289bc588SBarry Smith     }
538289bc588SBarry Smith   }
5391ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
5401ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5413a40ed3dSBarry Smith   PetscFunctionReturn(0);
542289bc588SBarry Smith }
543289bc588SBarry Smith 
544289bc588SBarry Smith /* -----------------------------------------------------------------*/
5454a2ae208SSatish Balay #undef __FUNCT__
5464a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
547dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
548289bc588SBarry Smith {
549c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
55087828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
551dfbe8321SBarry Smith   PetscErrorCode ierr;
5520805154bSBarry Smith   PetscBLASInt   m, n,_One=1;
553ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
5543a40ed3dSBarry Smith 
5553a40ed3dSBarry Smith   PetscFunctionBegin;
556d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
557d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
558d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5591ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5601ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
56171044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
5621ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5631ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
564dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
566289bc588SBarry Smith }
567800995b7SMatthew Knepley 
5684a2ae208SSatish Balay #undef __FUNCT__
5694a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
570dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
571289bc588SBarry Smith {
572c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
57387828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
574dfbe8321SBarry Smith   PetscErrorCode ierr;
5750805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5763a40ed3dSBarry Smith 
5773a40ed3dSBarry Smith   PetscFunctionBegin;
578d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
579d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
580d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5811ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5821ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
58371044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5841ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5851ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
586dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
5873a40ed3dSBarry Smith   PetscFunctionReturn(0);
588289bc588SBarry Smith }
5896ee01492SSatish Balay 
5904a2ae208SSatish Balay #undef __FUNCT__
5914a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
592dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
593289bc588SBarry Smith {
594c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
59587828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
596dfbe8321SBarry Smith   PetscErrorCode ierr;
5970805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5983a40ed3dSBarry Smith 
5993a40ed3dSBarry Smith   PetscFunctionBegin;
600d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
601d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
602d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
603600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
6041ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6051ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
60671044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
6071ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6081ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
609dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
6103a40ed3dSBarry Smith   PetscFunctionReturn(0);
611289bc588SBarry Smith }
6126ee01492SSatish Balay 
6134a2ae208SSatish Balay #undef __FUNCT__
6144a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
615dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
616289bc588SBarry Smith {
617c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
61887828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
619dfbe8321SBarry Smith   PetscErrorCode ierr;
6200805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
62187828ca2SBarry Smith   PetscScalar    _DOne=1.0;
6223a40ed3dSBarry Smith 
6233a40ed3dSBarry Smith   PetscFunctionBegin;
624d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
625d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
626d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
627600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
6281ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6291ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
63071044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
6311ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6321ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
633dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
6343a40ed3dSBarry Smith   PetscFunctionReturn(0);
635289bc588SBarry Smith }
636289bc588SBarry Smith 
637289bc588SBarry Smith /* -----------------------------------------------------------------*/
6384a2ae208SSatish Balay #undef __FUNCT__
6394a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
64013f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
641289bc588SBarry Smith {
642c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
64387828ca2SBarry Smith   PetscScalar    *v;
6446849ba73SBarry Smith   PetscErrorCode ierr;
64513f74950SBarry Smith   PetscInt       i;
64667e560aaSBarry Smith 
6473a40ed3dSBarry Smith   PetscFunctionBegin;
648d0f46423SBarry Smith   *ncols = A->cmap->n;
649289bc588SBarry Smith   if (cols) {
650d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
651d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
652289bc588SBarry Smith   }
653289bc588SBarry Smith   if (vals) {
654d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
655289bc588SBarry Smith     v    = mat->v + row;
656d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
657289bc588SBarry Smith   }
6583a40ed3dSBarry Smith   PetscFunctionReturn(0);
659289bc588SBarry Smith }
6606ee01492SSatish Balay 
6614a2ae208SSatish Balay #undef __FUNCT__
6624a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
66313f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
664289bc588SBarry Smith {
665dfbe8321SBarry Smith   PetscErrorCode ierr;
666606d414cSSatish Balay   PetscFunctionBegin;
667606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
668606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6693a40ed3dSBarry Smith   PetscFunctionReturn(0);
670289bc588SBarry Smith }
671289bc588SBarry Smith /* ----------------------------------------------------------------*/
6724a2ae208SSatish Balay #undef __FUNCT__
6734a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
67413f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
675289bc588SBarry Smith {
676c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
67713f74950SBarry Smith   PetscInt     i,j,idx=0;
678d6dfbf8fSBarry Smith 
6793a40ed3dSBarry Smith   PetscFunctionBegin;
68071fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
681289bc588SBarry Smith   if (!mat->roworiented) {
682dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
683289bc588SBarry Smith       for (j=0; j<n; j++) {
684cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6852515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
686e32f2f54SBarry 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);
68758804f6dSBarry Smith #endif
688289bc588SBarry Smith         for (i=0; i<m; i++) {
689cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6902515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
691e32f2f54SBarry 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);
69258804f6dSBarry Smith #endif
693cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
694289bc588SBarry Smith         }
695289bc588SBarry Smith       }
6963a40ed3dSBarry Smith     } else {
697289bc588SBarry Smith       for (j=0; j<n; j++) {
698cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6992515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
700e32f2f54SBarry 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);
70158804f6dSBarry Smith #endif
702289bc588SBarry Smith         for (i=0; i<m; i++) {
703cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
7042515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
705e32f2f54SBarry 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);
70658804f6dSBarry Smith #endif
707cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
708289bc588SBarry Smith         }
709289bc588SBarry Smith       }
710289bc588SBarry Smith     }
7113a40ed3dSBarry Smith   } else {
712dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
713e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
714cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
7152515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
716e32f2f54SBarry 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);
71758804f6dSBarry Smith #endif
718e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
719cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
7202515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
721e32f2f54SBarry 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);
72258804f6dSBarry Smith #endif
723cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
724e8d4e0b9SBarry Smith         }
725e8d4e0b9SBarry Smith       }
7263a40ed3dSBarry Smith     } else {
727289bc588SBarry Smith       for (i=0; i<m; i++) {
728cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
7292515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
730e32f2f54SBarry 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);
73158804f6dSBarry Smith #endif
732289bc588SBarry Smith         for (j=0; j<n; j++) {
733cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
7342515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
735e32f2f54SBarry 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);
73658804f6dSBarry Smith #endif
737cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
738289bc588SBarry Smith         }
739289bc588SBarry Smith       }
740289bc588SBarry Smith     }
741e8d4e0b9SBarry Smith   }
7423a40ed3dSBarry Smith   PetscFunctionReturn(0);
743289bc588SBarry Smith }
744e8d4e0b9SBarry Smith 
7454a2ae208SSatish Balay #undef __FUNCT__
7464a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
74713f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
748ae80bb75SLois Curfman McInnes {
749ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
75013f74950SBarry Smith   PetscInt     i,j;
751ae80bb75SLois Curfman McInnes 
7523a40ed3dSBarry Smith   PetscFunctionBegin;
753ae80bb75SLois Curfman McInnes   /* row-oriented output */
754ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
75597e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
756e32f2f54SBarry 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);
757ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7586f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
759e32f2f54SBarry 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);
76097e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
761ae80bb75SLois Curfman McInnes     }
762ae80bb75SLois Curfman McInnes   }
7633a40ed3dSBarry Smith   PetscFunctionReturn(0);
764ae80bb75SLois Curfman McInnes }
765ae80bb75SLois Curfman McInnes 
766289bc588SBarry Smith /* -----------------------------------------------------------------*/
767289bc588SBarry Smith 
7684a2ae208SSatish Balay #undef __FUNCT__
7695bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense"
770112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer)
771aabbc4fbSShri Abhyankar {
772aabbc4fbSShri Abhyankar   Mat_SeqDense   *a;
773aabbc4fbSShri Abhyankar   PetscErrorCode ierr;
774aabbc4fbSShri Abhyankar   PetscInt       *scols,i,j,nz,header[4];
775aabbc4fbSShri Abhyankar   int            fd;
776aabbc4fbSShri Abhyankar   PetscMPIInt    size;
777aabbc4fbSShri Abhyankar   PetscInt       *rowlengths = 0,M,N,*cols,grows,gcols;
778aabbc4fbSShri Abhyankar   PetscScalar    *vals,*svals,*v,*w;
779aabbc4fbSShri Abhyankar   MPI_Comm       comm = ((PetscObject)viewer)->comm;
780aabbc4fbSShri Abhyankar 
781aabbc4fbSShri Abhyankar   PetscFunctionBegin;
782aabbc4fbSShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
783aabbc4fbSShri Abhyankar   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
784aabbc4fbSShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
785aabbc4fbSShri Abhyankar   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
786aabbc4fbSShri Abhyankar   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
787aabbc4fbSShri Abhyankar   M = header[1]; N = header[2]; nz = header[3];
788aabbc4fbSShri Abhyankar 
789aabbc4fbSShri Abhyankar   /* set global size if not set already*/
790aabbc4fbSShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
791aabbc4fbSShri Abhyankar     ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr);
792aabbc4fbSShri Abhyankar   } else {
793aabbc4fbSShri Abhyankar     /* if sizes and type are already set, check if the vector global sizes are correct */
794aabbc4fbSShri Abhyankar     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
795aabbc4fbSShri 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);
796aabbc4fbSShri Abhyankar   }
797aabbc4fbSShri Abhyankar   ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
798aabbc4fbSShri Abhyankar 
799aabbc4fbSShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
800aabbc4fbSShri Abhyankar     a    = (Mat_SeqDense*)newmat->data;
801aabbc4fbSShri Abhyankar     v    = a->v;
802aabbc4fbSShri Abhyankar     /* Allocate some temp space to read in the values and then flip them
803aabbc4fbSShri Abhyankar        from row major to column major */
804aabbc4fbSShri Abhyankar     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
805aabbc4fbSShri Abhyankar     /* read in nonzero values */
806aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
807aabbc4fbSShri Abhyankar     /* now flip the values and store them in the matrix*/
808aabbc4fbSShri Abhyankar     for (j=0; j<N; j++) {
809aabbc4fbSShri Abhyankar       for (i=0; i<M; i++) {
810aabbc4fbSShri Abhyankar         *v++ =w[i*N+j];
811aabbc4fbSShri Abhyankar       }
812aabbc4fbSShri Abhyankar     }
813aabbc4fbSShri Abhyankar     ierr = PetscFree(w);CHKERRQ(ierr);
814aabbc4fbSShri Abhyankar   } else {
815aabbc4fbSShri Abhyankar     /* read row lengths */
816aabbc4fbSShri Abhyankar     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
817aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
818aabbc4fbSShri Abhyankar 
819aabbc4fbSShri Abhyankar     a = (Mat_SeqDense*)newmat->data;
820aabbc4fbSShri Abhyankar     v = a->v;
821aabbc4fbSShri Abhyankar 
822aabbc4fbSShri Abhyankar     /* read column indices and nonzeros */
823aabbc4fbSShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
824aabbc4fbSShri Abhyankar     cols = scols;
825aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
826aabbc4fbSShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
827aabbc4fbSShri Abhyankar     vals = svals;
828aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
829aabbc4fbSShri Abhyankar 
830aabbc4fbSShri Abhyankar     /* insert into matrix */
831aabbc4fbSShri Abhyankar     for (i=0; i<M; i++) {
832aabbc4fbSShri Abhyankar       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
833aabbc4fbSShri Abhyankar       svals += rowlengths[i]; scols += rowlengths[i];
834aabbc4fbSShri Abhyankar     }
835aabbc4fbSShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
836aabbc4fbSShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
837aabbc4fbSShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
838aabbc4fbSShri Abhyankar   }
839aabbc4fbSShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
840aabbc4fbSShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
841aabbc4fbSShri Abhyankar 
842aabbc4fbSShri Abhyankar   PetscFunctionReturn(0);
843aabbc4fbSShri Abhyankar }
844aabbc4fbSShri Abhyankar 
845aabbc4fbSShri Abhyankar #undef __FUNCT__
8464a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8476849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
848289bc588SBarry Smith {
849932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
850dfbe8321SBarry Smith   PetscErrorCode    ierr;
85113f74950SBarry Smith   PetscInt          i,j;
8522dcb1b2aSMatthew Knepley   const char        *name;
85387828ca2SBarry Smith   PetscScalar       *v;
854f3ef73ceSBarry Smith   PetscViewerFormat format;
8555f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
856ace3abfcSBarry Smith   PetscBool         allreal = PETSC_TRUE;
8575f481a85SSatish Balay #endif
858932b0c3eSLois Curfman McInnes 
8593a40ed3dSBarry Smith   PetscFunctionBegin;
860b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
861456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8623a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
863fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
864d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
8657566de4bSShri Abhyankar     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
866d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
86744cd7ae7SLois Curfman McInnes       v = a->v + i;
86877431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
869d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
870aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
871329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
872a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
873329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
874a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8756831982aSBarry Smith         }
87680cd9d93SLois Curfman McInnes #else
8776831982aSBarry Smith         if (*v) {
878a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8796831982aSBarry Smith         }
88080cd9d93SLois Curfman McInnes #endif
8811b807ce4Svictorle         v += a->lda;
88280cd9d93SLois Curfman McInnes       }
883b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
88480cd9d93SLois Curfman McInnes     }
885d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
8863a40ed3dSBarry Smith   } else {
887d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
888aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
88947989497SBarry Smith     /* determine if matrix has all real values */
89047989497SBarry Smith     v = a->v;
891d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
892ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
89347989497SBarry Smith     }
89447989497SBarry Smith #endif
895fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8963a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
897d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
898d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
899fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
9007566de4bSShri Abhyankar     } else {
9017566de4bSShri Abhyankar       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
902ffac6cdbSBarry Smith     }
903ffac6cdbSBarry Smith 
904d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
905932b0c3eSLois Curfman McInnes       v = a->v + i;
906d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
907aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
90847989497SBarry Smith         if (allreal) {
909f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
91047989497SBarry Smith         } else {
911f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
91247989497SBarry Smith         }
913289bc588SBarry Smith #else
914f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
915289bc588SBarry Smith #endif
9161b807ce4Svictorle 	v += a->lda;
917289bc588SBarry Smith       }
918b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
919289bc588SBarry Smith     }
920fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
921b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
922ffac6cdbSBarry Smith     }
923d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
924da3a660dSBarry Smith   }
925b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
9263a40ed3dSBarry Smith   PetscFunctionReturn(0);
927289bc588SBarry Smith }
928289bc588SBarry Smith 
9294a2ae208SSatish Balay #undef __FUNCT__
9304a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
9316849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
932932b0c3eSLois Curfman McInnes {
933932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9346849ba73SBarry Smith   PetscErrorCode    ierr;
93513f74950SBarry Smith   int               fd;
936d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
937f4403165SShri Abhyankar   PetscScalar       *v,*anonz,*vals;
938f4403165SShri Abhyankar   PetscViewerFormat format;
939932b0c3eSLois Curfman McInnes 
9403a40ed3dSBarry Smith   PetscFunctionBegin;
941b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
94290ace30eSBarry Smith 
943f4403165SShri Abhyankar   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
944f4403165SShri Abhyankar   if (format == PETSC_VIEWER_NATIVE) {
945f4403165SShri Abhyankar     /* store the matrix as a dense matrix */
946f4403165SShri Abhyankar     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
947f4403165SShri Abhyankar     col_lens[0] = MAT_FILE_CLASSID;
948f4403165SShri Abhyankar     col_lens[1] = m;
949f4403165SShri Abhyankar     col_lens[2] = n;
950f4403165SShri Abhyankar     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
951f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
952f4403165SShri Abhyankar     ierr = PetscFree(col_lens);CHKERRQ(ierr);
953f4403165SShri Abhyankar 
954f4403165SShri Abhyankar     /* write out matrix, by rows */
955f4403165SShri Abhyankar     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
956f4403165SShri Abhyankar     v    = a->v;
957f4403165SShri Abhyankar     for (j=0; j<n; j++) {
958f4403165SShri Abhyankar       for (i=0; i<m; i++) {
959f4403165SShri Abhyankar         vals[j + i*n] = *v++;
960f4403165SShri Abhyankar       }
961f4403165SShri Abhyankar     }
962f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
963f4403165SShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
964f4403165SShri Abhyankar   } else {
96513f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
9660700a824SBarry Smith     col_lens[0] = MAT_FILE_CLASSID;
967932b0c3eSLois Curfman McInnes     col_lens[1] = m;
968932b0c3eSLois Curfman McInnes     col_lens[2] = n;
969932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
970932b0c3eSLois Curfman McInnes 
971932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
972932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9736f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
974932b0c3eSLois Curfman McInnes 
975932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
976932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
977932b0c3eSLois Curfman McInnes     ict = 0;
978932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
979932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
980932b0c3eSLois Curfman McInnes     }
9816f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
982606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
983932b0c3eSLois Curfman McInnes 
984932b0c3eSLois Curfman McInnes     /* store nonzero values */
98587828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
986932b0c3eSLois Curfman McInnes     ict  = 0;
987932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
988932b0c3eSLois Curfman McInnes       v = a->v + i;
989932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9901b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
991932b0c3eSLois Curfman McInnes       }
992932b0c3eSLois Curfman McInnes     }
9936f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
994606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
995f4403165SShri Abhyankar   }
9963a40ed3dSBarry Smith   PetscFunctionReturn(0);
997932b0c3eSLois Curfman McInnes }
998932b0c3eSLois Curfman McInnes 
9994a2ae208SSatish Balay #undef __FUNCT__
10004a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
1001dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1002f1af5d2fSBarry Smith {
1003f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
1004f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
10056849ba73SBarry Smith   PetscErrorCode    ierr;
1006d0f46423SBarry Smith   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
100787828ca2SBarry Smith   PetscScalar       *v = a->v;
1008b0a32e0cSBarry Smith   PetscViewer       viewer;
1009b0a32e0cSBarry Smith   PetscDraw         popup;
1010329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
1011f3ef73ceSBarry Smith   PetscViewerFormat format;
1012f1af5d2fSBarry Smith 
1013f1af5d2fSBarry Smith   PetscFunctionBegin;
1014f1af5d2fSBarry Smith 
1015f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
1016b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1017b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1018f1af5d2fSBarry Smith 
1019f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
1020fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1021f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
1022b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
1023f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1024f1af5d2fSBarry Smith       x_l = j;
1025f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1026f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1027f1af5d2fSBarry Smith         y_l = m - i - 1.0;
1028f1af5d2fSBarry Smith         y_r = y_l + 1.0;
1029f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
1030329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
1031b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1032329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
1033b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1034f1af5d2fSBarry Smith         } else {
1035f1af5d2fSBarry Smith           continue;
1036f1af5d2fSBarry Smith         }
1037f1af5d2fSBarry Smith #else
1038f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
1039b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1040f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1041b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1042f1af5d2fSBarry Smith         } else {
1043f1af5d2fSBarry Smith           continue;
1044f1af5d2fSBarry Smith         }
1045f1af5d2fSBarry Smith #endif
1046b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1047f1af5d2fSBarry Smith       }
1048f1af5d2fSBarry Smith     }
1049f1af5d2fSBarry Smith   } else {
1050f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1051f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1052f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1053f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1054f1af5d2fSBarry Smith     }
1055b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1056b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1057b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1058f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1059f1af5d2fSBarry Smith       x_l = j;
1060f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1061f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1062f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1063f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1064b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1065b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1066f1af5d2fSBarry Smith       }
1067f1af5d2fSBarry Smith     }
1068f1af5d2fSBarry Smith   }
1069f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1070f1af5d2fSBarry Smith }
1071f1af5d2fSBarry Smith 
10724a2ae208SSatish Balay #undef __FUNCT__
10734a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1074dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1075f1af5d2fSBarry Smith {
1076b0a32e0cSBarry Smith   PetscDraw      draw;
1077ace3abfcSBarry Smith   PetscBool      isnull;
1078329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1079dfbe8321SBarry Smith   PetscErrorCode ierr;
1080f1af5d2fSBarry Smith 
1081f1af5d2fSBarry Smith   PetscFunctionBegin;
1082b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1083b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1084abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1085f1af5d2fSBarry Smith 
1086f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1087d0f46423SBarry Smith   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1088f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1089b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1090b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1091f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1092f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1093f1af5d2fSBarry Smith }
1094f1af5d2fSBarry Smith 
10954a2ae208SSatish Balay #undef __FUNCT__
10964a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1097dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1098932b0c3eSLois Curfman McInnes {
1099dfbe8321SBarry Smith   PetscErrorCode ierr;
1100ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw;
1101932b0c3eSLois Curfman McInnes 
11023a40ed3dSBarry Smith   PetscFunctionBegin;
11032692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11042692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
11052692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
11060f5bd95cSBarry Smith 
1107c45a1595SBarry Smith   if (iascii) {
1108c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
11090f5bd95cSBarry Smith   } else if (isbinary) {
11103a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1111f1af5d2fSBarry Smith   } else if (isdraw) {
1112f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
11135cd90555SBarry Smith   } else {
1114e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1115932b0c3eSLois Curfman McInnes   }
11163a40ed3dSBarry Smith   PetscFunctionReturn(0);
1117932b0c3eSLois Curfman McInnes }
1118289bc588SBarry Smith 
11194a2ae208SSatish Balay #undef __FUNCT__
11204a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1121dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1122289bc588SBarry Smith {
1123ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1124dfbe8321SBarry Smith   PetscErrorCode ierr;
112590f02eecSBarry Smith 
11263a40ed3dSBarry Smith   PetscFunctionBegin;
1127aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1128d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1129a5a9c739SBarry Smith #endif
113005b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
11316857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1132bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
1133dbd8c25aSHong Zhang 
1134dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1135901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
11364ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11374ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11384ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11393a40ed3dSBarry Smith   PetscFunctionReturn(0);
1140289bc588SBarry Smith }
1141289bc588SBarry Smith 
11424a2ae208SSatish Balay #undef __FUNCT__
11434a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1144fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1145289bc588SBarry Smith {
1146c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11476849ba73SBarry Smith   PetscErrorCode ierr;
114813f74950SBarry Smith   PetscInt       k,j,m,n,M;
114987828ca2SBarry Smith   PetscScalar    *v,tmp;
115048b35521SBarry Smith 
11513a40ed3dSBarry Smith   PetscFunctionBegin;
1152d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1153e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1154e7e72b3dSBarry Smith     if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1155e7e72b3dSBarry Smith     else {
1156d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1157289bc588SBarry Smith         for (k=0; k<j; k++) {
11581b807ce4Svictorle           tmp = v[j + k*M];
11591b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11601b807ce4Svictorle           v[k + j*M] = tmp;
1161289bc588SBarry Smith         }
1162289bc588SBarry Smith       }
1163d64ed03dSBarry Smith     }
11643a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1165d3e5ee88SLois Curfman McInnes     Mat          tmat;
1166ec8511deSBarry Smith     Mat_SeqDense *tmatd;
116787828ca2SBarry Smith     PetscScalar  *v2;
1168ea709b57SSatish Balay 
1169fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11707adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1171d0f46423SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
11727adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11735c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1174fc4dec0aSBarry Smith     } else {
1175fc4dec0aSBarry Smith       tmat = *matout;
1176fc4dec0aSBarry Smith     }
1177ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11780de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1179d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11801b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1181d3e5ee88SLois Curfman McInnes     }
11826d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11836d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1184d3e5ee88SLois Curfman McInnes     *matout = tmat;
118548b35521SBarry Smith   }
11863a40ed3dSBarry Smith   PetscFunctionReturn(0);
1187289bc588SBarry Smith }
1188289bc588SBarry Smith 
11894a2ae208SSatish Balay #undef __FUNCT__
11904a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1191ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool  *flg)
1192289bc588SBarry Smith {
1193c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1194c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
119513f74950SBarry Smith   PetscInt     i,j;
119687828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11979ea5d5aeSSatish Balay 
11983a40ed3dSBarry Smith   PetscFunctionBegin;
1199d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1200d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1201d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
12021b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1203d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
12043a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
12051b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
12061b807ce4Svictorle     }
1207289bc588SBarry Smith   }
120877c4ece6SBarry Smith   *flg = PETSC_TRUE;
12093a40ed3dSBarry Smith   PetscFunctionReturn(0);
1210289bc588SBarry Smith }
1211289bc588SBarry Smith 
12124a2ae208SSatish Balay #undef __FUNCT__
12134a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1214dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1215289bc588SBarry Smith {
1216c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1217dfbe8321SBarry Smith   PetscErrorCode ierr;
121813f74950SBarry Smith   PetscInt       i,n,len;
121987828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
122044cd7ae7SLois Curfman McInnes 
12213a40ed3dSBarry Smith   PetscFunctionBegin;
12222dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
12237a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
12241ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1225d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1226e32f2f54SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
122744cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
12281b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1229289bc588SBarry Smith   }
12301ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
12313a40ed3dSBarry Smith   PetscFunctionReturn(0);
1232289bc588SBarry Smith }
1233289bc588SBarry Smith 
12344a2ae208SSatish Balay #undef __FUNCT__
12354a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1236dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1237289bc588SBarry Smith {
1238c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
123987828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1240dfbe8321SBarry Smith   PetscErrorCode ierr;
1241d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
124255659b69SBarry Smith 
12433a40ed3dSBarry Smith   PetscFunctionBegin;
124428988994SBarry Smith   if (ll) {
12457a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12461ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1247e32f2f54SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1248da3a660dSBarry Smith     for (i=0; i<m; i++) {
1249da3a660dSBarry Smith       x = l[i];
1250da3a660dSBarry Smith       v = mat->v + i;
1251da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1252da3a660dSBarry Smith     }
12531ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1254efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1255da3a660dSBarry Smith   }
125628988994SBarry Smith   if (rr) {
12577a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12581ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1259e32f2f54SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1260da3a660dSBarry Smith     for (i=0; i<n; i++) {
1261da3a660dSBarry Smith       x = r[i];
1262da3a660dSBarry Smith       v = mat->v + i*m;
1263da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1264da3a660dSBarry Smith     }
12651ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1266efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1267da3a660dSBarry Smith   }
12683a40ed3dSBarry Smith   PetscFunctionReturn(0);
1269289bc588SBarry Smith }
1270289bc588SBarry Smith 
12714a2ae208SSatish Balay #undef __FUNCT__
12724a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1273dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1274289bc588SBarry Smith {
1275c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
127687828ca2SBarry Smith   PetscScalar  *v = mat->v;
1277329f5518SBarry Smith   PetscReal    sum = 0.0;
1278d0f46423SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1279efee365bSSatish Balay   PetscErrorCode ierr;
128055659b69SBarry Smith 
12813a40ed3dSBarry Smith   PetscFunctionBegin;
1282289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1283a5ce6ee0Svictorle     if (lda>m) {
1284d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1285a5ce6ee0Svictorle 	v = mat->v+j*lda;
1286a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1287a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1288a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1289a5ce6ee0Svictorle #else
1290a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1291a5ce6ee0Svictorle #endif
1292a5ce6ee0Svictorle 	}
1293a5ce6ee0Svictorle       }
1294a5ce6ee0Svictorle     } else {
1295d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1296aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1297329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1298289bc588SBarry Smith #else
1299289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1300289bc588SBarry Smith #endif
1301289bc588SBarry Smith       }
1302a5ce6ee0Svictorle     }
1303064f8208SBarry Smith     *nrm = sqrt(sum);
1304dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
13053a40ed3dSBarry Smith   } else if (type == NORM_1) {
1306064f8208SBarry Smith     *nrm = 0.0;
1307d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
13081b807ce4Svictorle       v = mat->v + j*mat->lda;
1309289bc588SBarry Smith       sum = 0.0;
1310d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
131133a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1312289bc588SBarry Smith       }
1313064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1314289bc588SBarry Smith     }
1315d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
13163a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1317064f8208SBarry Smith     *nrm = 0.0;
1318d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1319289bc588SBarry Smith       v = mat->v + j;
1320289bc588SBarry Smith       sum = 0.0;
1321d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
13221b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1323289bc588SBarry Smith       }
1324064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1325289bc588SBarry Smith     }
1326d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1327e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
13283a40ed3dSBarry Smith   PetscFunctionReturn(0);
1329289bc588SBarry Smith }
1330289bc588SBarry Smith 
13314a2ae208SSatish Balay #undef __FUNCT__
13324a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1333ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool  flg)
1334289bc588SBarry Smith {
1335c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
133663ba0a88SBarry Smith   PetscErrorCode ierr;
133767e560aaSBarry Smith 
13383a40ed3dSBarry Smith   PetscFunctionBegin;
1339b5a2b587SKris Buschelman   switch (op) {
1340b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13414e0d8c25SBarry Smith     aij->roworiented = flg;
1342b5a2b587SKris Buschelman     break;
1343512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1344b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13453971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13464e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1347b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1348b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
134977e54ba9SKris Buschelman   case MAT_SYMMETRIC:
135077e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13519a4540c5SBarry Smith   case MAT_HERMITIAN:
13529a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1353600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1354290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
135577e54ba9SKris Buschelman     break;
1356b5a2b587SKris Buschelman   default:
1357e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13583a40ed3dSBarry Smith   }
13593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1360289bc588SBarry Smith }
1361289bc588SBarry Smith 
13624a2ae208SSatish Balay #undef __FUNCT__
13634a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1364dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13656f0a148fSBarry Smith {
1366ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13676849ba73SBarry Smith   PetscErrorCode ierr;
1368d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13693a40ed3dSBarry Smith 
13703a40ed3dSBarry Smith   PetscFunctionBegin;
1371a5ce6ee0Svictorle   if (lda>m) {
1372d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1373a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1374a5ce6ee0Svictorle     }
1375a5ce6ee0Svictorle   } else {
1376d0f46423SBarry Smith     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1377a5ce6ee0Svictorle   }
13783a40ed3dSBarry Smith   PetscFunctionReturn(0);
13796f0a148fSBarry Smith }
13806f0a148fSBarry Smith 
13814a2ae208SSatish Balay #undef __FUNCT__
13824a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
13832b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
13846f0a148fSBarry Smith {
138597b48c8fSBarry Smith   PetscErrorCode    ierr;
1386ec8511deSBarry Smith   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
1387d0f46423SBarry Smith   PetscInt          n = A->cmap->n,i,j;
138897b48c8fSBarry Smith   PetscScalar       *slot,*bb;
138997b48c8fSBarry Smith   const PetscScalar *xx;
139055659b69SBarry Smith 
13913a40ed3dSBarry Smith   PetscFunctionBegin;
139297b48c8fSBarry Smith   /* fix right hand side if needed */
139397b48c8fSBarry Smith   if (x && b) {
139497b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
139597b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
139697b48c8fSBarry Smith     for (i=0; i<N; i++) {
139797b48c8fSBarry Smith       bb[rows[i]] = diag*xx[rows[i]];
139897b48c8fSBarry Smith     }
139997b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
140097b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
140197b48c8fSBarry Smith   }
140297b48c8fSBarry Smith 
14036f0a148fSBarry Smith   for (i=0; i<N; i++) {
14046f0a148fSBarry Smith     slot = l->v + rows[i];
14056f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
14066f0a148fSBarry Smith   }
1407f4df32b1SMatthew Knepley   if (diag != 0.0) {
14086f0a148fSBarry Smith     for (i=0; i<N; i++) {
14096f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1410f4df32b1SMatthew Knepley       *slot = diag;
14116f0a148fSBarry Smith     }
14126f0a148fSBarry Smith   }
14133a40ed3dSBarry Smith   PetscFunctionReturn(0);
14146f0a148fSBarry Smith }
1415557bce09SLois Curfman McInnes 
14164a2ae208SSatish Balay #undef __FUNCT__
14174a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1418dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
141964e87e97SBarry Smith {
1420c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
14213a40ed3dSBarry Smith 
14223a40ed3dSBarry Smith   PetscFunctionBegin;
1423e32f2f54SBarry 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");
142464e87e97SBarry Smith   *array = mat->v;
14253a40ed3dSBarry Smith   PetscFunctionReturn(0);
142664e87e97SBarry Smith }
14270754003eSLois Curfman McInnes 
14284a2ae208SSatish Balay #undef __FUNCT__
14294a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1430dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1431ff14e315SSatish Balay {
14323a40ed3dSBarry Smith   PetscFunctionBegin;
143309b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
14343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1435ff14e315SSatish Balay }
14360754003eSLois Curfman McInnes 
14374a2ae208SSatish Balay #undef __FUNCT__
14384a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
143913f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
14400754003eSLois Curfman McInnes {
1441c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
14426849ba73SBarry Smith   PetscErrorCode ierr;
14435d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
14445d0c19d7SBarry Smith   const PetscInt *irow,*icol;
144587828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
14460754003eSLois Curfman McInnes   Mat            newmat;
14470754003eSLois Curfman McInnes 
14483a40ed3dSBarry Smith   PetscFunctionBegin;
144978b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
145078b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1451e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1452e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14530754003eSLois Curfman McInnes 
1454182d2002SSatish Balay   /* Check submatrixcall */
1455182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
145613f74950SBarry Smith     PetscInt n_cols,n_rows;
1457182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
145821a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
1459f746d493SDmitry Karpeev       /* resize the result matrix to match number of requested rows/columns */
1460c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
146121a2c019SBarry Smith     }
1462182d2002SSatish Balay     newmat = *B;
1463182d2002SSatish Balay   } else {
14640754003eSLois Curfman McInnes     /* Create and fill new matrix */
14657adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1466f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14677adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14685c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1469182d2002SSatish Balay   }
1470182d2002SSatish Balay 
1471182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1472182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1473182d2002SSatish Balay 
1474182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14756de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1476182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1477182d2002SSatish Balay       *bv++ = av[irow[j]];
14780754003eSLois Curfman McInnes     }
14790754003eSLois Curfman McInnes   }
1480182d2002SSatish Balay 
1481182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14826d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14836d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14840754003eSLois Curfman McInnes 
14850754003eSLois Curfman McInnes   /* Free work space */
148678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
148778b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1488182d2002SSatish Balay   *B = newmat;
14893a40ed3dSBarry Smith   PetscFunctionReturn(0);
14900754003eSLois Curfman McInnes }
14910754003eSLois Curfman McInnes 
14924a2ae208SSatish Balay #undef __FUNCT__
14934a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
149413f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1495905e6a2fSBarry Smith {
14966849ba73SBarry Smith   PetscErrorCode ierr;
149713f74950SBarry Smith   PetscInt       i;
1498905e6a2fSBarry Smith 
14993a40ed3dSBarry Smith   PetscFunctionBegin;
1500905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1501b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1502905e6a2fSBarry Smith   }
1503905e6a2fSBarry Smith 
1504905e6a2fSBarry Smith   for (i=0; i<n; i++) {
15056a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1506905e6a2fSBarry Smith   }
15073a40ed3dSBarry Smith   PetscFunctionReturn(0);
1508905e6a2fSBarry Smith }
1509905e6a2fSBarry Smith 
15104a2ae208SSatish Balay #undef __FUNCT__
1511c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1512c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1513c0aa2d19SHong Zhang {
1514c0aa2d19SHong Zhang   PetscFunctionBegin;
1515c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1516c0aa2d19SHong Zhang }
1517c0aa2d19SHong Zhang 
1518c0aa2d19SHong Zhang #undef __FUNCT__
1519c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1520c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1521c0aa2d19SHong Zhang {
1522c0aa2d19SHong Zhang   PetscFunctionBegin;
1523c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1524c0aa2d19SHong Zhang }
1525c0aa2d19SHong Zhang 
1526c0aa2d19SHong Zhang #undef __FUNCT__
15274a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1528dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
15294b0e389bSBarry Smith {
15304b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
15316849ba73SBarry Smith   PetscErrorCode ierr;
1532d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
15333a40ed3dSBarry Smith 
15343a40ed3dSBarry Smith   PetscFunctionBegin;
153533f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
153633f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1537cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
15383a40ed3dSBarry Smith     PetscFunctionReturn(0);
15393a40ed3dSBarry Smith   }
1540e32f2f54SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1541a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
15420dbb7854Svictorle     for (j=0; j<n; j++) {
1543a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1544a5ce6ee0Svictorle     }
1545a5ce6ee0Svictorle   } else {
1546d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1547a5ce6ee0Svictorle   }
1548273d9f13SBarry Smith   PetscFunctionReturn(0);
1549273d9f13SBarry Smith }
1550273d9f13SBarry Smith 
15514a2ae208SSatish Balay #undef __FUNCT__
15524a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1553dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1554273d9f13SBarry Smith {
1555dfbe8321SBarry Smith   PetscErrorCode ierr;
1556273d9f13SBarry Smith 
1557273d9f13SBarry Smith   PetscFunctionBegin;
1558273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15593a40ed3dSBarry Smith   PetscFunctionReturn(0);
15604b0e389bSBarry Smith }
15614b0e389bSBarry Smith 
1562284134d9SBarry Smith #undef __FUNCT__
1563284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1564284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1565284134d9SBarry Smith {
1566284134d9SBarry Smith   PetscFunctionBegin;
156721a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1568284134d9SBarry Smith   m = PetscMax(m,M);
1569284134d9SBarry Smith   n = PetscMax(n,N);
1570a868139aSShri Abhyankar 
157186d161a7SShri 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);
157286d161a7SShri 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);
157386d161a7SShri Abhyankar   */
1574dc5cefdeSJed Brown   A->rmap->n = A->rmap->N = m;
1575d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
1576284134d9SBarry Smith   PetscFunctionReturn(0);
1577284134d9SBarry Smith }
1578170fe5c8SBarry Smith 
1579ba337c44SJed Brown #undef __FUNCT__
1580ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense"
1581ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A)
1582ba337c44SJed Brown {
1583ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1584ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1585ba337c44SJed Brown   PetscScalar    *aa = a->v;
1586ba337c44SJed Brown 
1587ba337c44SJed Brown   PetscFunctionBegin;
1588ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
1589ba337c44SJed Brown   PetscFunctionReturn(0);
1590ba337c44SJed Brown }
1591ba337c44SJed Brown 
1592ba337c44SJed Brown #undef __FUNCT__
1593ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense"
1594ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A)
1595ba337c44SJed Brown {
1596ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1597ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1598ba337c44SJed Brown   PetscScalar    *aa = a->v;
1599ba337c44SJed Brown 
1600ba337c44SJed Brown   PetscFunctionBegin;
1601ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1602ba337c44SJed Brown   PetscFunctionReturn(0);
1603ba337c44SJed Brown }
1604ba337c44SJed Brown 
1605ba337c44SJed Brown #undef __FUNCT__
1606ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense"
1607ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
1608ba337c44SJed Brown {
1609ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1610ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1611ba337c44SJed Brown   PetscScalar    *aa = a->v;
1612ba337c44SJed Brown 
1613ba337c44SJed Brown   PetscFunctionBegin;
1614ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1615ba337c44SJed Brown   PetscFunctionReturn(0);
1616ba337c44SJed Brown }
1617284134d9SBarry Smith 
1618a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1619a9fe9ddaSSatish Balay #undef __FUNCT__
1620a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1621a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1622a9fe9ddaSSatish Balay {
1623a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1624a9fe9ddaSSatish Balay 
1625a9fe9ddaSSatish Balay   PetscFunctionBegin;
1626a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1627a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1628a9fe9ddaSSatish Balay   }
1629a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1630a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1631a9fe9ddaSSatish Balay }
1632a9fe9ddaSSatish Balay 
1633a9fe9ddaSSatish Balay #undef __FUNCT__
1634a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1635a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1636a9fe9ddaSSatish Balay {
1637ee16a9a1SHong Zhang   PetscErrorCode ierr;
1638d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1639ee16a9a1SHong Zhang   Mat            Cmat;
1640a9fe9ddaSSatish Balay 
1641ee16a9a1SHong Zhang   PetscFunctionBegin;
1642e32f2f54SBarry 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);
1643ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1644ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1645ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1646ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1647ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1648ee16a9a1SHong Zhang   *C = Cmat;
1649ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1650ee16a9a1SHong Zhang }
1651a9fe9ddaSSatish Balay 
165298a3b096SSatish Balay #undef __FUNCT__
1653a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1654a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1655a9fe9ddaSSatish Balay {
1656a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1657a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1658a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16590805154bSBarry Smith   PetscBLASInt   m,n,k;
1660a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1661a9fe9ddaSSatish Balay 
1662a9fe9ddaSSatish Balay   PetscFunctionBegin;
1663d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1664d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1665d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1666a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1667a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1668a9fe9ddaSSatish Balay }
1669a9fe9ddaSSatish Balay 
1670a9fe9ddaSSatish Balay #undef __FUNCT__
1671a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1672a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1673a9fe9ddaSSatish Balay {
1674a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1675a9fe9ddaSSatish Balay 
1676a9fe9ddaSSatish Balay   PetscFunctionBegin;
1677a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1678a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1679a9fe9ddaSSatish Balay   }
1680a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1681a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1682a9fe9ddaSSatish Balay }
1683a9fe9ddaSSatish Balay 
1684a9fe9ddaSSatish Balay #undef __FUNCT__
1685a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1686a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1687a9fe9ddaSSatish Balay {
1688ee16a9a1SHong Zhang   PetscErrorCode ierr;
1689d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1690ee16a9a1SHong Zhang   Mat            Cmat;
1691a9fe9ddaSSatish Balay 
1692ee16a9a1SHong Zhang   PetscFunctionBegin;
1693e32f2f54SBarry 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);
1694ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1695ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1696ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1697ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1698ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1699ee16a9a1SHong Zhang   *C = Cmat;
1700ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1701ee16a9a1SHong Zhang }
1702a9fe9ddaSSatish Balay 
1703a9fe9ddaSSatish Balay #undef __FUNCT__
1704a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1705a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1706a9fe9ddaSSatish Balay {
1707a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1708a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1709a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
17100805154bSBarry Smith   PetscBLASInt   m,n,k;
1711a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1712a9fe9ddaSSatish Balay 
1713a9fe9ddaSSatish Balay   PetscFunctionBegin;
1714d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1715d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1716d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
17172fbe02b9SBarry Smith   /*
17182fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
17192fbe02b9SBarry Smith   */
1720a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1721a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1722a9fe9ddaSSatish Balay }
1723985db425SBarry Smith 
1724985db425SBarry Smith #undef __FUNCT__
1725985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1726985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1727985db425SBarry Smith {
1728985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1729985db425SBarry Smith   PetscErrorCode ierr;
1730d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1731985db425SBarry Smith   PetscScalar    *x;
1732985db425SBarry Smith   MatScalar      *aa = a->v;
1733985db425SBarry Smith 
1734985db425SBarry Smith   PetscFunctionBegin;
1735e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1736985db425SBarry Smith 
1737985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1738985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1739985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1740e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1741985db425SBarry Smith   for (i=0; i<m; i++) {
1742985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1743985db425SBarry Smith     for (j=1; j<n; j++){
1744985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1745985db425SBarry Smith     }
1746985db425SBarry Smith   }
1747985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1748985db425SBarry Smith   PetscFunctionReturn(0);
1749985db425SBarry Smith }
1750985db425SBarry Smith 
1751985db425SBarry Smith #undef __FUNCT__
1752985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1753985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1754985db425SBarry Smith {
1755985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1756985db425SBarry Smith   PetscErrorCode ierr;
1757d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1758985db425SBarry Smith   PetscScalar    *x;
1759985db425SBarry Smith   PetscReal      atmp;
1760985db425SBarry Smith   MatScalar      *aa = a->v;
1761985db425SBarry Smith 
1762985db425SBarry Smith   PetscFunctionBegin;
1763e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1764985db425SBarry Smith 
1765985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1766985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1767985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1768e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1769985db425SBarry Smith   for (i=0; i<m; i++) {
17709189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1771985db425SBarry Smith     for (j=1; j<n; j++){
1772985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1773985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1774985db425SBarry Smith     }
1775985db425SBarry Smith   }
1776985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1777985db425SBarry Smith   PetscFunctionReturn(0);
1778985db425SBarry Smith }
1779985db425SBarry Smith 
1780985db425SBarry Smith #undef __FUNCT__
1781985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1782985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1783985db425SBarry Smith {
1784985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1785985db425SBarry Smith   PetscErrorCode ierr;
1786d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1787985db425SBarry Smith   PetscScalar    *x;
1788985db425SBarry Smith   MatScalar      *aa = a->v;
1789985db425SBarry Smith 
1790985db425SBarry Smith   PetscFunctionBegin;
1791e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1792985db425SBarry Smith 
1793985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1794985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1795985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1796e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1797985db425SBarry Smith   for (i=0; i<m; i++) {
1798985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1799985db425SBarry Smith     for (j=1; j<n; j++){
1800985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1801985db425SBarry Smith     }
1802985db425SBarry Smith   }
1803985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1804985db425SBarry Smith   PetscFunctionReturn(0);
1805985db425SBarry Smith }
1806985db425SBarry Smith 
18078d0534beSBarry Smith #undef __FUNCT__
18088d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
18098d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
18108d0534beSBarry Smith {
18118d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
18128d0534beSBarry Smith   PetscErrorCode ierr;
18138d0534beSBarry Smith   PetscScalar    *x;
18148d0534beSBarry Smith 
18158d0534beSBarry Smith   PetscFunctionBegin;
1816e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
18178d0534beSBarry Smith 
18188d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1819d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
18208d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
18218d0534beSBarry Smith   PetscFunctionReturn(0);
18228d0534beSBarry Smith }
18238d0534beSBarry Smith 
1824*0716a85fSBarry Smith 
1825*0716a85fSBarry Smith #undef __FUNCT__
1826*0716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense"
1827*0716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
1828*0716a85fSBarry Smith {
1829*0716a85fSBarry Smith   PetscErrorCode ierr;
1830*0716a85fSBarry Smith   PetscInt       i,j,m,n;
1831*0716a85fSBarry Smith   PetscScalar    *a;
1832*0716a85fSBarry Smith 
1833*0716a85fSBarry Smith   PetscFunctionBegin;
1834*0716a85fSBarry Smith   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
1835*0716a85fSBarry Smith   ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr);
1836*0716a85fSBarry Smith   ierr = MatGetArray(A,&a);CHKERRQ(ierr);
1837*0716a85fSBarry Smith   if (type == NORM_2) {
1838*0716a85fSBarry Smith     for (i=0; i<n; i++ ){
1839*0716a85fSBarry Smith       for (j=0; j<m; j++) {
1840*0716a85fSBarry Smith 	norms[i] += PetscAbsScalar(a[j]*a[j]);
1841*0716a85fSBarry Smith       }
1842*0716a85fSBarry Smith       a += m;
1843*0716a85fSBarry Smith     }
1844*0716a85fSBarry Smith   } else if (type == NORM_1) {
1845*0716a85fSBarry Smith     for (i=0; i<n; i++ ){
1846*0716a85fSBarry Smith       for (j=0; j<m; j++) {
1847*0716a85fSBarry Smith 	norms[i] += PetscAbsScalar(a[j]);
1848*0716a85fSBarry Smith       }
1849*0716a85fSBarry Smith       a += m;
1850*0716a85fSBarry Smith     }
1851*0716a85fSBarry Smith   } else if (type == NORM_INFINITY) {
1852*0716a85fSBarry Smith     for (i=0; i<n; i++ ){
1853*0716a85fSBarry Smith       for (j=0; j<m; j++) {
1854*0716a85fSBarry Smith 	norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
1855*0716a85fSBarry Smith       }
1856*0716a85fSBarry Smith       a += m;
1857*0716a85fSBarry Smith     }
1858*0716a85fSBarry Smith   } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType");
1859*0716a85fSBarry Smith   if (type == NORM_2) {
1860*0716a85fSBarry Smith     for (i=0; i<n; i++) norms[i] = sqrt(norms[i]);
1861*0716a85fSBarry Smith   }
1862*0716a85fSBarry Smith   PetscFunctionReturn(0);
1863*0716a85fSBarry Smith }
1864*0716a85fSBarry Smith 
1865289bc588SBarry Smith /* -------------------------------------------------------------------*/
1866a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1867905e6a2fSBarry Smith        MatGetRow_SeqDense,
1868905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1869905e6a2fSBarry Smith        MatMult_SeqDense,
187097304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
18717c922b88SBarry Smith        MatMultTranspose_SeqDense,
18727c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1873db4efbfdSBarry Smith        0,
1874db4efbfdSBarry Smith        0,
1875db4efbfdSBarry Smith        0,
1876db4efbfdSBarry Smith /*10*/ 0,
1877905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1878905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
187941f059aeSBarry Smith        MatSOR_SeqDense,
1880ec8511deSBarry Smith        MatTranspose_SeqDense,
188197304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1882905e6a2fSBarry Smith        MatEqual_SeqDense,
1883905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1884905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1885905e6a2fSBarry Smith        MatNorm_SeqDense,
1886c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1887c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1888905e6a2fSBarry Smith        MatSetOption_SeqDense,
1889905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1890d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense,
1891db4efbfdSBarry Smith        0,
1892db4efbfdSBarry Smith        0,
1893db4efbfdSBarry Smith        0,
1894db4efbfdSBarry Smith        0,
1895d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense,
1896273d9f13SBarry Smith        0,
1897905e6a2fSBarry Smith        0,
1898905e6a2fSBarry Smith        MatGetArray_SeqDense,
1899905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1900d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense,
1901a5ae1ecdSBarry Smith        0,
1902a5ae1ecdSBarry Smith        0,
1903a5ae1ecdSBarry Smith        0,
1904a5ae1ecdSBarry Smith        0,
1905d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense,
1906a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1907a5ae1ecdSBarry Smith        0,
19084b0e389bSBarry Smith        MatGetValues_SeqDense,
1909a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1910d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense,
1911a5ae1ecdSBarry Smith        MatScale_SeqDense,
1912a5ae1ecdSBarry Smith        0,
1913a5ae1ecdSBarry Smith        0,
1914a5ae1ecdSBarry Smith        0,
1915d519adbfSMatthew Knepley /*49*/ 0,
1916a5ae1ecdSBarry Smith        0,
1917a5ae1ecdSBarry Smith        0,
1918a5ae1ecdSBarry Smith        0,
1919a5ae1ecdSBarry Smith        0,
1920d519adbfSMatthew Knepley /*54*/ 0,
1921a5ae1ecdSBarry Smith        0,
1922a5ae1ecdSBarry Smith        0,
1923a5ae1ecdSBarry Smith        0,
1924a5ae1ecdSBarry Smith        0,
1925d519adbfSMatthew Knepley /*59*/ 0,
1926e03a110bSBarry Smith        MatDestroy_SeqDense,
1927e03a110bSBarry Smith        MatView_SeqDense,
1928357abbc8SBarry Smith        0,
192997304618SKris Buschelman        0,
1930d519adbfSMatthew Knepley /*64*/ 0,
193197304618SKris Buschelman        0,
193297304618SKris Buschelman        0,
193397304618SKris Buschelman        0,
193497304618SKris Buschelman        0,
1935d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense,
193697304618SKris Buschelman        0,
193797304618SKris Buschelman        0,
193897304618SKris Buschelman        0,
193997304618SKris Buschelman        0,
1940d519adbfSMatthew Knepley /*74*/ 0,
194197304618SKris Buschelman        0,
194297304618SKris Buschelman        0,
194397304618SKris Buschelman        0,
194497304618SKris Buschelman        0,
1945d519adbfSMatthew Knepley /*79*/ 0,
194697304618SKris Buschelman        0,
194797304618SKris Buschelman        0,
194897304618SKris Buschelman        0,
19495bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense,
1950865e5f61SKris Buschelman        0,
19511cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1952865e5f61SKris Buschelman        0,
1953865e5f61SKris Buschelman        0,
1954865e5f61SKris Buschelman        0,
1955d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense,
1956a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1957a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1958865e5f61SKris Buschelman        0,
1959865e5f61SKris Buschelman        0,
1960d519adbfSMatthew Knepley /*94*/ 0,
1961a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1962a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1963a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1964284134d9SBarry Smith        0,
1965d519adbfSMatthew Knepley /*99*/ 0,
1966284134d9SBarry Smith        0,
1967284134d9SBarry Smith        0,
1968ba337c44SJed Brown        MatConjugate_SeqDense,
1969985db425SBarry Smith        MatSetSizes_SeqDense,
1970ba337c44SJed Brown /*104*/0,
1971ba337c44SJed Brown        MatRealPart_SeqDense,
1972ba337c44SJed Brown        MatImaginaryPart_SeqDense,
1973985db425SBarry Smith        0,
1974985db425SBarry Smith        0,
197585e2c93fSHong Zhang /*109*/MatMatSolve_SeqDense,
1976985db425SBarry Smith        0,
19778d0534beSBarry Smith        MatGetRowMin_SeqDense,
1978aabbc4fbSShri Abhyankar        MatGetColumnVector_SeqDense,
1979aabbc4fbSShri Abhyankar        0,
1980aabbc4fbSShri Abhyankar /*114*/0,
1981aabbc4fbSShri Abhyankar        0,
1982aabbc4fbSShri Abhyankar        0,
1983aabbc4fbSShri Abhyankar        0,
1984aabbc4fbSShri Abhyankar        0,
1985aabbc4fbSShri Abhyankar /*119*/0,
1986aabbc4fbSShri Abhyankar        0,
1987aabbc4fbSShri Abhyankar        0,
1988*0716a85fSBarry Smith        0,
1989*0716a85fSBarry Smith        0,
1990*0716a85fSBarry Smith /*124*/0,
1991*0716a85fSBarry Smith        MatGetColumnNorms_SeqDense
1992985db425SBarry Smith };
199390ace30eSBarry Smith 
19944a2ae208SSatish Balay #undef __FUNCT__
19954a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
19964b828684SBarry Smith /*@C
1997fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1998d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1999d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
2000289bc588SBarry Smith 
2001db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
2002db81eaa0SLois Curfman McInnes 
200320563c6bSBarry Smith    Input Parameters:
2004db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
20050c775827SLois Curfman McInnes .  m - number of rows
200618f449edSLois Curfman McInnes .  n - number of columns
2007c0235b3cSMatthew Knepley -  data - optional location of matrix data in column major order.  Set data=PETSC_NULL for PETSc
2008dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
200920563c6bSBarry Smith 
201020563c6bSBarry Smith    Output Parameter:
201144cd7ae7SLois Curfman McInnes .  A - the matrix
201220563c6bSBarry Smith 
2013b259b22eSLois Curfman McInnes    Notes:
201418f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
201518f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
2016b4fd4287SBarry Smith    set data=PETSC_NULL.
201718f449edSLois Curfman McInnes 
2018027ccd11SLois Curfman McInnes    Level: intermediate
2019027ccd11SLois Curfman McInnes 
2020dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
2021d65003e9SLois Curfman McInnes 
2022db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
202320563c6bSBarry Smith @*/
20247087cfbeSBarry Smith PetscErrorCode  MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2025289bc588SBarry Smith {
2026dfbe8321SBarry Smith   PetscErrorCode ierr;
20273b2fbd54SBarry Smith 
20283a40ed3dSBarry Smith   PetscFunctionBegin;
2029f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2030f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
2031273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
2032273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
2033273d9f13SBarry Smith   PetscFunctionReturn(0);
2034273d9f13SBarry Smith }
2035273d9f13SBarry Smith 
20364a2ae208SSatish Balay #undef __FUNCT__
2037afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
2038273d9f13SBarry Smith /*@C
2039273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
2040273d9f13SBarry Smith 
2041273d9f13SBarry Smith    Collective on MPI_Comm
2042273d9f13SBarry Smith 
2043273d9f13SBarry Smith    Input Parameters:
2044273d9f13SBarry Smith +  A - the matrix
2045273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
2046273d9f13SBarry Smith 
2047273d9f13SBarry Smith    Notes:
2048273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
2049273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
2050284134d9SBarry Smith    need not call this routine.
2051273d9f13SBarry Smith 
2052273d9f13SBarry Smith    Level: intermediate
2053273d9f13SBarry Smith 
2054273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
2055273d9f13SBarry Smith 
2056867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
2057867c911aSBarry Smith 
2058273d9f13SBarry Smith @*/
20597087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2060273d9f13SBarry Smith {
20614ac538c5SBarry Smith   PetscErrorCode ierr;
2062a23d5eceSKris Buschelman 
2063a23d5eceSKris Buschelman   PetscFunctionBegin;
20644ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr);
2065a23d5eceSKris Buschelman   PetscFunctionReturn(0);
2066a23d5eceSKris Buschelman }
2067a23d5eceSKris Buschelman 
2068a23d5eceSKris Buschelman EXTERN_C_BEGIN
2069a23d5eceSKris Buschelman #undef __FUNCT__
2070afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
20717087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2072a23d5eceSKris Buschelman {
2073273d9f13SBarry Smith   Mat_SeqDense   *b;
2074dfbe8321SBarry Smith   PetscErrorCode ierr;
2075273d9f13SBarry Smith 
2076273d9f13SBarry Smith   PetscFunctionBegin;
2077273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
2078a868139aSShri Abhyankar 
207934ef9618SShri Abhyankar   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
208034ef9618SShri Abhyankar   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
208134ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
208234ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
208334ef9618SShri Abhyankar 
2084273d9f13SBarry Smith   b       = (Mat_SeqDense*)B->data;
208586d161a7SShri Abhyankar   b->Mmax = B->rmap->n;
208686d161a7SShri Abhyankar   b->Nmax = B->cmap->n;
208786d161a7SShri Abhyankar   if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n;
208886d161a7SShri Abhyankar 
20899e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
20909e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
20915afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
2092284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
2093284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
20949e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
2095273d9f13SBarry Smith   } else { /* user-allocated storage */
20969e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
2097273d9f13SBarry Smith     b->v          = data;
2098273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
2099273d9f13SBarry Smith   }
21000450473dSBarry Smith   B->assembled = PETSC_TRUE;
2101273d9f13SBarry Smith   PetscFunctionReturn(0);
2102273d9f13SBarry Smith }
2103a23d5eceSKris Buschelman EXTERN_C_END
2104273d9f13SBarry Smith 
21051b807ce4Svictorle #undef __FUNCT__
21061b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
21071b807ce4Svictorle /*@C
21081b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
21091b807ce4Svictorle 
21101b807ce4Svictorle   Input parameter:
21111b807ce4Svictorle + A - the matrix
21121b807ce4Svictorle - lda - the leading dimension
21131b807ce4Svictorle 
21141b807ce4Svictorle   Notes:
2115867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
21161b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
21171b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
21181b807ce4Svictorle 
21191b807ce4Svictorle   Level: intermediate
21201b807ce4Svictorle 
21211b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
21221b807ce4Svictorle 
2123284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2124867c911aSBarry Smith 
21251b807ce4Svictorle @*/
21267087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetLDA(Mat B,PetscInt lda)
21271b807ce4Svictorle {
21281b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
212921a2c019SBarry Smith 
21301b807ce4Svictorle   PetscFunctionBegin;
2131e32f2f54SBarry 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);
21321b807ce4Svictorle   b->lda       = lda;
213321a2c019SBarry Smith   b->changelda = PETSC_FALSE;
213421a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
21351b807ce4Svictorle   PetscFunctionReturn(0);
21361b807ce4Svictorle }
21371b807ce4Svictorle 
21380bad9183SKris Buschelman /*MC
2139fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
21400bad9183SKris Buschelman 
21410bad9183SKris Buschelman    Options Database Keys:
21420bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
21430bad9183SKris Buschelman 
21440bad9183SKris Buschelman   Level: beginner
21450bad9183SKris Buschelman 
214689665df3SBarry Smith .seealso: MatCreateSeqDense()
214789665df3SBarry Smith 
21480bad9183SKris Buschelman M*/
21490bad9183SKris Buschelman 
2150273d9f13SBarry Smith EXTERN_C_BEGIN
21514a2ae208SSatish Balay #undef __FUNCT__
21524a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
21537087cfbeSBarry Smith PetscErrorCode  MatCreate_SeqDense(Mat B)
2154273d9f13SBarry Smith {
2155273d9f13SBarry Smith   Mat_SeqDense   *b;
2156dfbe8321SBarry Smith   PetscErrorCode ierr;
21577c334f02SBarry Smith   PetscMPIInt    size;
2158273d9f13SBarry Smith 
2159273d9f13SBarry Smith   PetscFunctionBegin;
21607adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
2161e32f2f54SBarry Smith   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
216255659b69SBarry Smith 
216338f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2164549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
216544cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
216618f449edSLois Curfman McInnes 
216744cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2168273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2169273d9f13SBarry Smith   b->v            = 0;
217021a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
21714e220ebcSLois Curfman McInnes 
2172b24902e0SBarry Smith 
2173ec1065edSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
2174b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2175b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2176a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2177a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2178a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
21794ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
21804ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
21814ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
21824ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
21834ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
21844ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
21854ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
21864ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
21874ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
218817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
21893a40ed3dSBarry Smith   PetscFunctionReturn(0);
2190289bc588SBarry Smith }
2191273d9f13SBarry Smith EXTERN_C_END
2192