xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 7566de4b946d0c7b7a74549da1150676a63056f9)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
77c4f633dSBarry Smith #include "../src/mat/impls/dense/seq/dense.h"
8f3da1532SBarry Smith #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
170805154bSBarry Smith   PetscBLASInt   N,m,ldax,lday,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21d0f46423SBarry Smith   N    = PetscBLASIntCast(X->rmap->n*X->cmap->n);
22d0f46423SBarry Smith   m    = PetscBLASIntCast(X->rmap->n);
230805154bSBarry Smith   ldax = PetscBLASIntCast(x->lda);
240805154bSBarry Smith   lday = PetscBLASIntCast(y->lda);
25a5ce6ee0Svictorle   if (ldax>m || lday>m) {
26d0f46423SBarry Smith     for (j=0; j<X->cmap->n; j++) {
27f4df32b1SMatthew Knepley       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
28a5ce6ee0Svictorle     }
29a5ce6ee0Svictorle   } else {
30f4df32b1SMatthew Knepley     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
31a5ce6ee0Svictorle   }
320450473dSBarry Smith   ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr);
333a40ed3dSBarry Smith   PetscFunctionReturn(0);
341987afe7SBarry Smith }
351987afe7SBarry Smith 
364a2ae208SSatish Balay #undef __FUNCT__
374a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
38dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
39289bc588SBarry Smith {
40d0f46423SBarry Smith   PetscInt     N = A->rmap->n*A->cmap->n;
413a40ed3dSBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
434e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
444e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
456de62eeeSBarry Smith   info->nz_used           = (double)N;
466de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
474e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
484e220ebcSLois Curfman McInnes   info->mallocs           = 0;
497adad957SLisandro Dalcin   info->memory            = ((PetscObject)A)->mem;
504e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
514e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
524e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
54289bc588SBarry Smith }
55289bc588SBarry Smith 
564a2ae208SSatish Balay #undef __FUNCT__
574a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
58f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
5980cd9d93SLois Curfman McInnes {
60273d9f13SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
61f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
62efee365bSSatish Balay   PetscErrorCode ierr;
630805154bSBarry Smith   PetscBLASInt   one = 1,j,nz,lda = PetscBLASIntCast(a->lda);
6480cd9d93SLois Curfman McInnes 
653a40ed3dSBarry Smith   PetscFunctionBegin;
66d0f46423SBarry Smith   if (lda>A->rmap->n) {
67d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n);
68d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
69f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
72d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n*A->cmap->n);
73f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
74a5ce6ee0Svictorle   }
75efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
7780cd9d93SLois Curfman McInnes }
7880cd9d93SLois Curfman McInnes 
791cbb95d3SBarry Smith #undef __FUNCT__
801cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense"
81ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool  *fl)
821cbb95d3SBarry Smith {
831cbb95d3SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
84d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,N;
851cbb95d3SBarry Smith   PetscScalar    *v = a->v;
861cbb95d3SBarry Smith 
871cbb95d3SBarry Smith   PetscFunctionBegin;
881cbb95d3SBarry Smith   *fl = PETSC_FALSE;
89d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0);
901cbb95d3SBarry Smith   N = a->lda;
911cbb95d3SBarry Smith 
921cbb95d3SBarry Smith   for (i=0; i<m; i++) {
931cbb95d3SBarry Smith     for (j=i+1; j<m; j++) {
941cbb95d3SBarry Smith       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
951cbb95d3SBarry Smith     }
961cbb95d3SBarry Smith   }
971cbb95d3SBarry Smith   *fl = PETSC_TRUE;
981cbb95d3SBarry Smith   PetscFunctionReturn(0);
991cbb95d3SBarry Smith }
1001cbb95d3SBarry Smith 
101b24902e0SBarry Smith #undef __FUNCT__
102b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense"
103719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
104b24902e0SBarry Smith {
105b24902e0SBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
106b24902e0SBarry Smith   PetscErrorCode ierr;
107b24902e0SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
108b24902e0SBarry Smith 
109b24902e0SBarry Smith   PetscFunctionBegin;
110719d5645SBarry Smith   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
111b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
112b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
113d0f46423SBarry Smith     if (lda>A->rmap->n) {
114d0f46423SBarry Smith       m = A->rmap->n;
115d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
116b24902e0SBarry Smith 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
117b24902e0SBarry Smith       }
118b24902e0SBarry Smith     } else {
119d0f46423SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
120b24902e0SBarry Smith     }
121b24902e0SBarry Smith   }
122b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
123b24902e0SBarry Smith   PetscFunctionReturn(0);
124b24902e0SBarry Smith }
125b24902e0SBarry Smith 
1264a2ae208SSatish Balay #undef __FUNCT__
1274a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
128dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12902cad45dSBarry Smith {
1306849ba73SBarry Smith   PetscErrorCode ierr;
13102cad45dSBarry Smith 
1323a40ed3dSBarry Smith   PetscFunctionBegin;
1335c9eb25fSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
134d0f46423SBarry Smith   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1355c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
136719d5645SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
137b24902e0SBarry Smith   PetscFunctionReturn(0);
138b24902e0SBarry Smith }
139b24902e0SBarry Smith 
1406ee01492SSatish Balay 
1410481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
142719d5645SBarry Smith 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1450481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
146289bc588SBarry Smith {
1474482741eSBarry Smith   MatFactorInfo  info;
148a093e273SMatthew Knepley   PetscErrorCode ierr;
1493a40ed3dSBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
151c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
152719d5645SBarry Smith   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154289bc588SBarry Smith }
1556ee01492SSatish Balay 
1560b4b3355SBarry Smith #undef __FUNCT__
1574a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
158dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
159289bc588SBarry Smith {
160c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1616849ba73SBarry Smith   PetscErrorCode ierr;
16287828ca2SBarry Smith   PetscScalar    *x,*y;
163d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
16467e560aaSBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
1661ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1671ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
168d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
169d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
170ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
171e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
172ae7cfcebSSatish Balay #else
17371044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
174e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
175ae7cfcebSSatish Balay #endif
176d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_CHOLESKY){
177ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
178e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
179ae7cfcebSSatish Balay #else
18071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
181e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
182ae7cfcebSSatish Balay #endif
183289bc588SBarry Smith   }
184e32f2f54SBarry Smith   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
1851ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1861ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
187dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
189289bc588SBarry Smith }
1906ee01492SSatish Balay 
1914a2ae208SSatish Balay #undef __FUNCT__
1924a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
193dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
194da3a660dSBarry Smith {
195c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
196dfbe8321SBarry Smith   PetscErrorCode ierr;
19787828ca2SBarry Smith   PetscScalar    *x,*y;
198d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
19967e560aaSBarry Smith 
2003a40ed3dSBarry Smith   PetscFunctionBegin;
2011ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2021ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
203d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
204752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
205da3a660dSBarry Smith   if (mat->pivots) {
206ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
207e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
208ae7cfcebSSatish Balay #else
20971044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
210e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
211ae7cfcebSSatish Balay #endif
2127a97a34bSBarry Smith   } else {
213ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
214e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
215ae7cfcebSSatish Balay #else
21671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
217e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
218ae7cfcebSSatish Balay #endif
219da3a660dSBarry Smith   }
2201ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2211ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
222dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
2233a40ed3dSBarry Smith   PetscFunctionReturn(0);
224da3a660dSBarry Smith }
2256ee01492SSatish Balay 
2264a2ae208SSatish Balay #undef __FUNCT__
2274a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
228dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
229da3a660dSBarry Smith {
230c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
231dfbe8321SBarry Smith   PetscErrorCode ierr;
23287828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
233da3a660dSBarry Smith   Vec            tmp = 0;
234d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
23567e560aaSBarry Smith 
2363a40ed3dSBarry Smith   PetscFunctionBegin;
2371ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2381ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
239d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
240da3a660dSBarry Smith   if (yy == zz) {
24178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
24252e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
24378b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
244da3a660dSBarry Smith   }
245d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
246752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
247da3a660dSBarry Smith   if (mat->pivots) {
248ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
249e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
250ae7cfcebSSatish Balay #else
25171044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
252e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
253ae7cfcebSSatish Balay #endif
254a8c6a408SBarry Smith   } else {
255ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
256e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
257ae7cfcebSSatish Balay #else
25871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
259e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
260ae7cfcebSSatish Balay #endif
261da3a660dSBarry Smith   }
2622dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
2632dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
2641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2651ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
266dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
2673a40ed3dSBarry Smith   PetscFunctionReturn(0);
268da3a660dSBarry Smith }
26967e560aaSBarry Smith 
2704a2ae208SSatish Balay #undef __FUNCT__
2714a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
272dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
273da3a660dSBarry Smith {
274c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2756849ba73SBarry Smith   PetscErrorCode ierr;
27687828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
277da3a660dSBarry Smith   Vec            tmp;
278d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
27967e560aaSBarry Smith 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
281d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
2821ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2831ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
284da3a660dSBarry Smith   if (yy == zz) {
28578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
28652e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
28778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
288da3a660dSBarry Smith   }
289d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
290752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
291da3a660dSBarry Smith   if (mat->pivots) {
292ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
293e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
294ae7cfcebSSatish Balay #else
29571044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
296e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
297ae7cfcebSSatish Balay #endif
2983a40ed3dSBarry Smith   } else {
299ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
300e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
301ae7cfcebSSatish Balay #else
30271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
303e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve");
304ae7cfcebSSatish Balay #endif
305da3a660dSBarry Smith   }
30690f02eecSBarry Smith   if (tmp) {
3072dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
30890f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3093a40ed3dSBarry Smith   } else {
3102dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
31190f02eecSBarry Smith   }
3121ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3131ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
314dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr);
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
316da3a660dSBarry Smith }
317db4efbfdSBarry Smith 
318db4efbfdSBarry Smith /* ---------------------------------------------------------------*/
319db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
320db4efbfdSBarry Smith    rather than put it in the Mat->row slot.*/
321db4efbfdSBarry Smith #undef __FUNCT__
322db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense"
3230481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
324db4efbfdSBarry Smith {
325db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
326db4efbfdSBarry Smith   PetscFunctionBegin;
327e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
328db4efbfdSBarry Smith #else
329db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
330db4efbfdSBarry Smith   PetscErrorCode ierr;
331db4efbfdSBarry Smith   PetscBLASInt   n,m,info;
332db4efbfdSBarry Smith 
333db4efbfdSBarry Smith   PetscFunctionBegin;
334db4efbfdSBarry Smith   n = PetscBLASIntCast(A->cmap->n);
335db4efbfdSBarry Smith   m = PetscBLASIntCast(A->rmap->n);
336db4efbfdSBarry Smith   if (!mat->pivots) {
337db4efbfdSBarry Smith     ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
338db4efbfdSBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
339db4efbfdSBarry Smith   }
340db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
341db4efbfdSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
342e32f2f54SBarry Smith   if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
343e32f2f54SBarry Smith   if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
344db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
345db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
346db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
347db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
348d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_LU;
349db4efbfdSBarry Smith 
350dc0b31edSSatish Balay   ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
351db4efbfdSBarry Smith #endif
352db4efbfdSBarry Smith   PetscFunctionReturn(0);
353db4efbfdSBarry Smith }
354db4efbfdSBarry Smith 
355db4efbfdSBarry Smith #undef __FUNCT__
356db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense"
3570481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
358db4efbfdSBarry Smith {
359db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
360db4efbfdSBarry Smith   PetscFunctionBegin;
361e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
362db4efbfdSBarry Smith #else
363db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
364db4efbfdSBarry Smith   PetscErrorCode ierr;
365db4efbfdSBarry Smith   PetscBLASInt   info,n = PetscBLASIntCast(A->cmap->n);
366db4efbfdSBarry Smith 
367db4efbfdSBarry Smith   PetscFunctionBegin;
368db4efbfdSBarry Smith   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
369db4efbfdSBarry Smith   mat->pivots = 0;
370db4efbfdSBarry Smith 
371db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
372db4efbfdSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
373e32f2f54SBarry Smith   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
374db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
375db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
376db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
377db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
378d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_CHOLESKY;
379dc0b31edSSatish Balay   ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr);
380db4efbfdSBarry Smith #endif
381db4efbfdSBarry Smith   PetscFunctionReturn(0);
382db4efbfdSBarry Smith }
383db4efbfdSBarry Smith 
384db4efbfdSBarry Smith 
385db4efbfdSBarry Smith #undef __FUNCT__
386db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
3870481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
388db4efbfdSBarry Smith {
389db4efbfdSBarry Smith   PetscErrorCode ierr;
390db4efbfdSBarry Smith   MatFactorInfo  info;
391db4efbfdSBarry Smith 
392db4efbfdSBarry Smith   PetscFunctionBegin;
393db4efbfdSBarry Smith   info.fill = 1.0;
394c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
395719d5645SBarry Smith   ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr);
396db4efbfdSBarry Smith   PetscFunctionReturn(0);
397db4efbfdSBarry Smith }
398db4efbfdSBarry Smith 
399db4efbfdSBarry Smith #undef __FUNCT__
400db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
4010481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
402db4efbfdSBarry Smith {
403db4efbfdSBarry Smith   PetscFunctionBegin;
404c3ef05f6SHong Zhang   fact->assembled                  = PETSC_TRUE;
405719d5645SBarry Smith   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
406db4efbfdSBarry Smith   PetscFunctionReturn(0);
407db4efbfdSBarry Smith }
408db4efbfdSBarry Smith 
409db4efbfdSBarry Smith #undef __FUNCT__
410db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
4110481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
412db4efbfdSBarry Smith {
413db4efbfdSBarry Smith   PetscFunctionBegin;
414c3ef05f6SHong Zhang   fact->assembled            = PETSC_TRUE;
415719d5645SBarry Smith   fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
416db4efbfdSBarry Smith   PetscFunctionReturn(0);
417db4efbfdSBarry Smith }
418db4efbfdSBarry Smith 
419bb5747d9SMatthew Knepley EXTERN_C_BEGIN
420db4efbfdSBarry Smith #undef __FUNCT__
421db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc"
422db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
423db4efbfdSBarry Smith {
424db4efbfdSBarry Smith   PetscErrorCode ierr;
425db4efbfdSBarry Smith 
426db4efbfdSBarry Smith   PetscFunctionBegin;
427db4efbfdSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
428db4efbfdSBarry Smith   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
429db4efbfdSBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
430db4efbfdSBarry Smith   if (ftype == MAT_FACTOR_LU){
431db4efbfdSBarry Smith     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
432db4efbfdSBarry Smith   } else {
433db4efbfdSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
434db4efbfdSBarry Smith   }
435d5f3da31SBarry Smith   (*fact)->factortype = ftype;
436db4efbfdSBarry Smith   PetscFunctionReturn(0);
437db4efbfdSBarry Smith }
438bb5747d9SMatthew Knepley EXTERN_C_END
439db4efbfdSBarry Smith 
440289bc588SBarry Smith /* ------------------------------------------------------------------*/
4414a2ae208SSatish Balay #undef __FUNCT__
44241f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense"
44341f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
444289bc588SBarry Smith {
445c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
44687828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
447dfbe8321SBarry Smith   PetscErrorCode ierr;
448d0f46423SBarry Smith   PetscInt       m = A->rmap->n,i;
449aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4500805154bSBarry Smith   PetscBLASInt   o = 1,bm = PetscBLASIntCast(m);
451bc1b551cSSatish Balay #endif
452289bc588SBarry Smith 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
454289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
45571044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4562dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
457289bc588SBarry Smith   }
4581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4591ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
460b965ef7fSBarry Smith   its  = its*lits;
461e32f2f54SBarry 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);
462289bc588SBarry Smith   while (its--) {
463fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
464289bc588SBarry Smith       for (i=0; i<m; i++) {
465aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
466f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
467f1747703SBarry Smith            not happy about returning a double complex */
46813f74950SBarry Smith         PetscInt    _i;
46987828ca2SBarry Smith         PetscScalar sum = b[i];
470f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4713f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
472f1747703SBarry Smith         }
473f1747703SBarry Smith         xt = sum;
474f1747703SBarry Smith #else
47571044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
476f1747703SBarry Smith #endif
47755a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
478289bc588SBarry Smith       }
479289bc588SBarry Smith     }
480fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
481289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
483f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
484f1747703SBarry Smith            not happy about returning a double complex */
48513f74950SBarry Smith         PetscInt    _i;
48687828ca2SBarry Smith         PetscScalar sum = b[i];
487f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4883f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
489f1747703SBarry Smith         }
490f1747703SBarry Smith         xt = sum;
491f1747703SBarry Smith #else
49271044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
493f1747703SBarry Smith #endif
49455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
495289bc588SBarry Smith       }
496289bc588SBarry Smith     }
497289bc588SBarry Smith   }
4981ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4991ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5003a40ed3dSBarry Smith   PetscFunctionReturn(0);
501289bc588SBarry Smith }
502289bc588SBarry Smith 
503289bc588SBarry Smith /* -----------------------------------------------------------------*/
5044a2ae208SSatish Balay #undef __FUNCT__
5054a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
506dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
507289bc588SBarry Smith {
508c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
50987828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
510dfbe8321SBarry Smith   PetscErrorCode ierr;
5110805154bSBarry Smith   PetscBLASInt   m, n,_One=1;
512ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
5133a40ed3dSBarry Smith 
5143a40ed3dSBarry Smith   PetscFunctionBegin;
515d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
516d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
517d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5191ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
52071044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
5211ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5221ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
523dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
5243a40ed3dSBarry Smith   PetscFunctionReturn(0);
525289bc588SBarry Smith }
526800995b7SMatthew Knepley 
5274a2ae208SSatish Balay #undef __FUNCT__
5284a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
529dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
530289bc588SBarry Smith {
531c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
53287828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
533dfbe8321SBarry Smith   PetscErrorCode ierr;
5340805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5353a40ed3dSBarry Smith 
5363a40ed3dSBarry Smith   PetscFunctionBegin;
537d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
538d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
539d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5401ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5411ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
54271044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5431ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5441ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
545dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
5463a40ed3dSBarry Smith   PetscFunctionReturn(0);
547289bc588SBarry Smith }
5486ee01492SSatish Balay 
5494a2ae208SSatish Balay #undef __FUNCT__
5504a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
551dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
552289bc588SBarry Smith {
553c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
55487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
555dfbe8321SBarry Smith   PetscErrorCode ierr;
5560805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5573a40ed3dSBarry Smith 
5583a40ed3dSBarry Smith   PetscFunctionBegin;
559d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
560d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
561d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
562600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5631ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5641ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
56571044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5661ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5671ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
568dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5693a40ed3dSBarry Smith   PetscFunctionReturn(0);
570289bc588SBarry Smith }
5716ee01492SSatish Balay 
5724a2ae208SSatish Balay #undef __FUNCT__
5734a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
574dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
575289bc588SBarry Smith {
576c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
57787828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
578dfbe8321SBarry Smith   PetscErrorCode ierr;
5790805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
58087828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5813a40ed3dSBarry Smith 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
583d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
584d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
585d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
586600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5871ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5881ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
58971044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5901ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5911ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
592dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5933a40ed3dSBarry Smith   PetscFunctionReturn(0);
594289bc588SBarry Smith }
595289bc588SBarry Smith 
596289bc588SBarry Smith /* -----------------------------------------------------------------*/
5974a2ae208SSatish Balay #undef __FUNCT__
5984a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
59913f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
600289bc588SBarry Smith {
601c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
60287828ca2SBarry Smith   PetscScalar    *v;
6036849ba73SBarry Smith   PetscErrorCode ierr;
60413f74950SBarry Smith   PetscInt       i;
60567e560aaSBarry Smith 
6063a40ed3dSBarry Smith   PetscFunctionBegin;
607d0f46423SBarry Smith   *ncols = A->cmap->n;
608289bc588SBarry Smith   if (cols) {
609d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
610d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
611289bc588SBarry Smith   }
612289bc588SBarry Smith   if (vals) {
613d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
614289bc588SBarry Smith     v    = mat->v + row;
615d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
616289bc588SBarry Smith   }
6173a40ed3dSBarry Smith   PetscFunctionReturn(0);
618289bc588SBarry Smith }
6196ee01492SSatish Balay 
6204a2ae208SSatish Balay #undef __FUNCT__
6214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
62213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
623289bc588SBarry Smith {
624dfbe8321SBarry Smith   PetscErrorCode ierr;
625606d414cSSatish Balay   PetscFunctionBegin;
626606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
627606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6283a40ed3dSBarry Smith   PetscFunctionReturn(0);
629289bc588SBarry Smith }
630289bc588SBarry Smith /* ----------------------------------------------------------------*/
6314a2ae208SSatish Balay #undef __FUNCT__
6324a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
63313f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
634289bc588SBarry Smith {
635c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
63613f74950SBarry Smith   PetscInt     i,j,idx=0;
637d6dfbf8fSBarry Smith 
6383a40ed3dSBarry Smith   PetscFunctionBegin;
63971fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
640289bc588SBarry Smith   if (!mat->roworiented) {
641dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
642289bc588SBarry Smith       for (j=0; j<n; j++) {
643cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6442515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
645e32f2f54SBarry 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);
64658804f6dSBarry Smith #endif
647289bc588SBarry Smith         for (i=0; i<m; i++) {
648cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6492515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
650e32f2f54SBarry 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);
65158804f6dSBarry Smith #endif
652cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
653289bc588SBarry Smith         }
654289bc588SBarry Smith       }
6553a40ed3dSBarry Smith     } else {
656289bc588SBarry Smith       for (j=0; j<n; j++) {
657cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6582515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
659e32f2f54SBarry 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);
66058804f6dSBarry Smith #endif
661289bc588SBarry Smith         for (i=0; i<m; i++) {
662cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6632515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
664e32f2f54SBarry 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);
66558804f6dSBarry Smith #endif
666cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
667289bc588SBarry Smith         }
668289bc588SBarry Smith       }
669289bc588SBarry Smith     }
6703a40ed3dSBarry Smith   } else {
671dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
672e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
673cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6742515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
675e32f2f54SBarry 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);
67658804f6dSBarry Smith #endif
677e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
678cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6792515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
680e32f2f54SBarry 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);
68158804f6dSBarry Smith #endif
682cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
683e8d4e0b9SBarry Smith         }
684e8d4e0b9SBarry Smith       }
6853a40ed3dSBarry Smith     } else {
686289bc588SBarry Smith       for (i=0; i<m; i++) {
687cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6882515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
689e32f2f54SBarry 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);
69058804f6dSBarry Smith #endif
691289bc588SBarry Smith         for (j=0; j<n; j++) {
692cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6932515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
694e32f2f54SBarry 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);
69558804f6dSBarry Smith #endif
696cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
697289bc588SBarry Smith         }
698289bc588SBarry Smith       }
699289bc588SBarry Smith     }
700e8d4e0b9SBarry Smith   }
7013a40ed3dSBarry Smith   PetscFunctionReturn(0);
702289bc588SBarry Smith }
703e8d4e0b9SBarry Smith 
7044a2ae208SSatish Balay #undef __FUNCT__
7054a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
70613f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
707ae80bb75SLois Curfman McInnes {
708ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
70913f74950SBarry Smith   PetscInt     i,j;
710ae80bb75SLois Curfman McInnes 
7113a40ed3dSBarry Smith   PetscFunctionBegin;
712ae80bb75SLois Curfman McInnes   /* row-oriented output */
713ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
71497e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
715e32f2f54SBarry 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);
716ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7176f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
718e32f2f54SBarry 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);
71997e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
720ae80bb75SLois Curfman McInnes     }
721ae80bb75SLois Curfman McInnes   }
7223a40ed3dSBarry Smith   PetscFunctionReturn(0);
723ae80bb75SLois Curfman McInnes }
724ae80bb75SLois Curfman McInnes 
725289bc588SBarry Smith /* -----------------------------------------------------------------*/
726289bc588SBarry Smith 
7274a2ae208SSatish Balay #undef __FUNCT__
7285bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense"
729112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer)
730aabbc4fbSShri Abhyankar {
731aabbc4fbSShri Abhyankar   Mat_SeqDense   *a;
732aabbc4fbSShri Abhyankar   PetscErrorCode ierr;
733aabbc4fbSShri Abhyankar   PetscInt       *scols,i,j,nz,header[4];
734aabbc4fbSShri Abhyankar   int            fd;
735aabbc4fbSShri Abhyankar   PetscMPIInt    size;
736aabbc4fbSShri Abhyankar   PetscInt       *rowlengths = 0,M,N,*cols,grows,gcols;
737aabbc4fbSShri Abhyankar   PetscScalar    *vals,*svals,*v,*w;
738aabbc4fbSShri Abhyankar   MPI_Comm       comm = ((PetscObject)viewer)->comm;
739aabbc4fbSShri Abhyankar 
740aabbc4fbSShri Abhyankar   PetscFunctionBegin;
741aabbc4fbSShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
742aabbc4fbSShri Abhyankar   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
743aabbc4fbSShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
744aabbc4fbSShri Abhyankar   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
745aabbc4fbSShri Abhyankar   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
746aabbc4fbSShri Abhyankar   M = header[1]; N = header[2]; nz = header[3];
747aabbc4fbSShri Abhyankar 
748aabbc4fbSShri Abhyankar   /* set global size if not set already*/
749aabbc4fbSShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
750aabbc4fbSShri Abhyankar     ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr);
751aabbc4fbSShri Abhyankar   } else {
752aabbc4fbSShri Abhyankar     /* if sizes and type are already set, check if the vector global sizes are correct */
753aabbc4fbSShri Abhyankar     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
754aabbc4fbSShri 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);
755aabbc4fbSShri Abhyankar   }
756aabbc4fbSShri Abhyankar   ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
757aabbc4fbSShri Abhyankar 
758aabbc4fbSShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
759aabbc4fbSShri Abhyankar     a    = (Mat_SeqDense*)newmat->data;
760aabbc4fbSShri Abhyankar     v    = a->v;
761aabbc4fbSShri Abhyankar     /* Allocate some temp space to read in the values and then flip them
762aabbc4fbSShri Abhyankar        from row major to column major */
763aabbc4fbSShri Abhyankar     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
764aabbc4fbSShri Abhyankar     /* read in nonzero values */
765aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
766aabbc4fbSShri Abhyankar     /* now flip the values and store them in the matrix*/
767aabbc4fbSShri Abhyankar     for (j=0; j<N; j++) {
768aabbc4fbSShri Abhyankar       for (i=0; i<M; i++) {
769aabbc4fbSShri Abhyankar         *v++ =w[i*N+j];
770aabbc4fbSShri Abhyankar       }
771aabbc4fbSShri Abhyankar     }
772aabbc4fbSShri Abhyankar     ierr = PetscFree(w);CHKERRQ(ierr);
773aabbc4fbSShri Abhyankar   } else {
774aabbc4fbSShri Abhyankar     /* read row lengths */
775aabbc4fbSShri Abhyankar     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
776aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
777aabbc4fbSShri Abhyankar 
778aabbc4fbSShri Abhyankar     a = (Mat_SeqDense*)newmat->data;
779aabbc4fbSShri Abhyankar     v = a->v;
780aabbc4fbSShri Abhyankar 
781aabbc4fbSShri Abhyankar     /* read column indices and nonzeros */
782aabbc4fbSShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
783aabbc4fbSShri Abhyankar     cols = scols;
784aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
785aabbc4fbSShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
786aabbc4fbSShri Abhyankar     vals = svals;
787aabbc4fbSShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
788aabbc4fbSShri Abhyankar 
789aabbc4fbSShri Abhyankar     /* insert into matrix */
790aabbc4fbSShri Abhyankar     for (i=0; i<M; i++) {
791aabbc4fbSShri Abhyankar       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
792aabbc4fbSShri Abhyankar       svals += rowlengths[i]; scols += rowlengths[i];
793aabbc4fbSShri Abhyankar     }
794aabbc4fbSShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
795aabbc4fbSShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
796aabbc4fbSShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
797aabbc4fbSShri Abhyankar   }
798aabbc4fbSShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
799aabbc4fbSShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
800aabbc4fbSShri Abhyankar 
801aabbc4fbSShri Abhyankar   PetscFunctionReturn(0);
802aabbc4fbSShri Abhyankar }
803aabbc4fbSShri Abhyankar 
804aabbc4fbSShri Abhyankar #undef __FUNCT__
8054a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8066849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
807289bc588SBarry Smith {
808932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
809dfbe8321SBarry Smith   PetscErrorCode    ierr;
81013f74950SBarry Smith   PetscInt          i,j;
8112dcb1b2aSMatthew Knepley   const char        *name;
81287828ca2SBarry Smith   PetscScalar       *v;
813f3ef73ceSBarry Smith   PetscViewerFormat format;
8145f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
815ace3abfcSBarry Smith   PetscBool         allreal = PETSC_TRUE;
8165f481a85SSatish Balay #endif
817932b0c3eSLois Curfman McInnes 
8183a40ed3dSBarry Smith   PetscFunctionBegin;
819435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
820b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
821456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8223a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
823fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
824b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
825*7566de4bSShri Abhyankar     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
826d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
82744cd7ae7SLois Curfman McInnes       v = a->v + i;
82877431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
829d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
830aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
831329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
832a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
833329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
834a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8356831982aSBarry Smith         }
83680cd9d93SLois Curfman McInnes #else
8376831982aSBarry Smith         if (*v) {
838a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8396831982aSBarry Smith         }
84080cd9d93SLois Curfman McInnes #endif
8411b807ce4Svictorle         v += a->lda;
84280cd9d93SLois Curfman McInnes       }
843b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
84480cd9d93SLois Curfman McInnes     }
845b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
8463a40ed3dSBarry Smith   } else {
847b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
848aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
84947989497SBarry Smith     /* determine if matrix has all real values */
85047989497SBarry Smith     v = a->v;
851d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
852ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
85347989497SBarry Smith     }
85447989497SBarry Smith #endif
855fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8563a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
857d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
858d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
859fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
860*7566de4bSShri Abhyankar     } else {
861*7566de4bSShri Abhyankar       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
862ffac6cdbSBarry Smith     }
863ffac6cdbSBarry Smith 
864d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
865932b0c3eSLois Curfman McInnes       v = a->v + i;
866d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
867aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
86847989497SBarry Smith         if (allreal) {
869f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
87047989497SBarry Smith         } else {
871f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
87247989497SBarry Smith         }
873289bc588SBarry Smith #else
874f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
875289bc588SBarry Smith #endif
8761b807ce4Svictorle 	v += a->lda;
877289bc588SBarry Smith       }
878b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
879289bc588SBarry Smith     }
880fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
881b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
882ffac6cdbSBarry Smith     }
883b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
884da3a660dSBarry Smith   }
885b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8863a40ed3dSBarry Smith   PetscFunctionReturn(0);
887289bc588SBarry Smith }
888289bc588SBarry Smith 
8894a2ae208SSatish Balay #undef __FUNCT__
8904a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8916849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
892932b0c3eSLois Curfman McInnes {
893932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8946849ba73SBarry Smith   PetscErrorCode    ierr;
89513f74950SBarry Smith   int               fd;
896d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
897f4403165SShri Abhyankar   PetscScalar       *v,*anonz,*vals;
898f4403165SShri Abhyankar   PetscViewerFormat format;
899932b0c3eSLois Curfman McInnes 
9003a40ed3dSBarry Smith   PetscFunctionBegin;
901b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
90290ace30eSBarry Smith 
903f4403165SShri Abhyankar   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
904f4403165SShri Abhyankar   if (format == PETSC_VIEWER_NATIVE) {
905f4403165SShri Abhyankar     /* store the matrix as a dense matrix */
906f4403165SShri Abhyankar     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
907f4403165SShri Abhyankar     col_lens[0] = MAT_FILE_CLASSID;
908f4403165SShri Abhyankar     col_lens[1] = m;
909f4403165SShri Abhyankar     col_lens[2] = n;
910f4403165SShri Abhyankar     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
911f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
912f4403165SShri Abhyankar     ierr = PetscFree(col_lens);CHKERRQ(ierr);
913f4403165SShri Abhyankar 
914f4403165SShri Abhyankar     /* write out matrix, by rows */
915f4403165SShri Abhyankar     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
916f4403165SShri Abhyankar     v    = a->v;
917f4403165SShri Abhyankar     for (j=0; j<n; j++) {
918f4403165SShri Abhyankar       for (i=0; i<m; i++) {
919f4403165SShri Abhyankar         vals[j + i*n] = *v++;
920f4403165SShri Abhyankar       }
921f4403165SShri Abhyankar     }
922f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
923f4403165SShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
924f4403165SShri Abhyankar   } else {
92513f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
9260700a824SBarry Smith     col_lens[0] = MAT_FILE_CLASSID;
927932b0c3eSLois Curfman McInnes     col_lens[1] = m;
928932b0c3eSLois Curfman McInnes     col_lens[2] = n;
929932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
930932b0c3eSLois Curfman McInnes 
931932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
932932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9336f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
934932b0c3eSLois Curfman McInnes 
935932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
936932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
937932b0c3eSLois Curfman McInnes     ict = 0;
938932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
939932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
940932b0c3eSLois Curfman McInnes     }
9416f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
942606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
943932b0c3eSLois Curfman McInnes 
944932b0c3eSLois Curfman McInnes     /* store nonzero values */
94587828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
946932b0c3eSLois Curfman McInnes     ict  = 0;
947932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
948932b0c3eSLois Curfman McInnes       v = a->v + i;
949932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9501b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
951932b0c3eSLois Curfman McInnes       }
952932b0c3eSLois Curfman McInnes     }
9536f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
954606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
955f4403165SShri Abhyankar   }
9563a40ed3dSBarry Smith   PetscFunctionReturn(0);
957932b0c3eSLois Curfman McInnes }
958932b0c3eSLois Curfman McInnes 
9594a2ae208SSatish Balay #undef __FUNCT__
9604a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
961dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
962f1af5d2fSBarry Smith {
963f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
964f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9656849ba73SBarry Smith   PetscErrorCode    ierr;
966d0f46423SBarry Smith   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
96787828ca2SBarry Smith   PetscScalar       *v = a->v;
968b0a32e0cSBarry Smith   PetscViewer       viewer;
969b0a32e0cSBarry Smith   PetscDraw         popup;
970329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
971f3ef73ceSBarry Smith   PetscViewerFormat format;
972f1af5d2fSBarry Smith 
973f1af5d2fSBarry Smith   PetscFunctionBegin;
974f1af5d2fSBarry Smith 
975f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
976b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
977b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
978f1af5d2fSBarry Smith 
979f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
980fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
981f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
982b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
983f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
984f1af5d2fSBarry Smith       x_l = j;
985f1af5d2fSBarry Smith       x_r = x_l + 1.0;
986f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
987f1af5d2fSBarry Smith         y_l = m - i - 1.0;
988f1af5d2fSBarry Smith         y_r = y_l + 1.0;
989f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
990329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
991b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
992329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
993b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
994f1af5d2fSBarry Smith         } else {
995f1af5d2fSBarry Smith           continue;
996f1af5d2fSBarry Smith         }
997f1af5d2fSBarry Smith #else
998f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
999b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1000f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1001b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1002f1af5d2fSBarry Smith         } else {
1003f1af5d2fSBarry Smith           continue;
1004f1af5d2fSBarry Smith         }
1005f1af5d2fSBarry Smith #endif
1006b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1007f1af5d2fSBarry Smith       }
1008f1af5d2fSBarry Smith     }
1009f1af5d2fSBarry Smith   } else {
1010f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1011f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1012f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1013f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1014f1af5d2fSBarry Smith     }
1015b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1016b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1017b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1018f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1019f1af5d2fSBarry Smith       x_l = j;
1020f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1021f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1022f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1023f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1024b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1025b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1026f1af5d2fSBarry Smith       }
1027f1af5d2fSBarry Smith     }
1028f1af5d2fSBarry Smith   }
1029f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1030f1af5d2fSBarry Smith }
1031f1af5d2fSBarry Smith 
10324a2ae208SSatish Balay #undef __FUNCT__
10334a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1034dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1035f1af5d2fSBarry Smith {
1036b0a32e0cSBarry Smith   PetscDraw      draw;
1037ace3abfcSBarry Smith   PetscBool      isnull;
1038329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1039dfbe8321SBarry Smith   PetscErrorCode ierr;
1040f1af5d2fSBarry Smith 
1041f1af5d2fSBarry Smith   PetscFunctionBegin;
1042b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1043b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1044abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1045f1af5d2fSBarry Smith 
1046f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1047d0f46423SBarry Smith   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1048f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1049b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1050b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1051f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1052f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1053f1af5d2fSBarry Smith }
1054f1af5d2fSBarry Smith 
10554a2ae208SSatish Balay #undef __FUNCT__
10564a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1057dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1058932b0c3eSLois Curfman McInnes {
1059dfbe8321SBarry Smith   PetscErrorCode ierr;
1060ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw;
1061932b0c3eSLois Curfman McInnes 
10623a40ed3dSBarry Smith   PetscFunctionBegin;
10632692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10642692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
10652692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
10660f5bd95cSBarry Smith 
1067c45a1595SBarry Smith   if (iascii) {
1068c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10690f5bd95cSBarry Smith   } else if (isbinary) {
10703a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1071f1af5d2fSBarry Smith   } else if (isdraw) {
1072f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10735cd90555SBarry Smith   } else {
1074e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1075932b0c3eSLois Curfman McInnes   }
10763a40ed3dSBarry Smith   PetscFunctionReturn(0);
1077932b0c3eSLois Curfman McInnes }
1078289bc588SBarry Smith 
10794a2ae208SSatish Balay #undef __FUNCT__
10804a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1081dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1082289bc588SBarry Smith {
1083ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1084dfbe8321SBarry Smith   PetscErrorCode ierr;
108590f02eecSBarry Smith 
10863a40ed3dSBarry Smith   PetscFunctionBegin;
1087aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1088d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1089a5a9c739SBarry Smith #endif
109005b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
10916857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1092606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1093dbd8c25aSHong Zhang 
1094dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1095901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10964ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10974ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10984ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10993a40ed3dSBarry Smith   PetscFunctionReturn(0);
1100289bc588SBarry Smith }
1101289bc588SBarry Smith 
11024a2ae208SSatish Balay #undef __FUNCT__
11034a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1104fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1105289bc588SBarry Smith {
1106c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11076849ba73SBarry Smith   PetscErrorCode ierr;
110813f74950SBarry Smith   PetscInt       k,j,m,n,M;
110987828ca2SBarry Smith   PetscScalar    *v,tmp;
111048b35521SBarry Smith 
11113a40ed3dSBarry Smith   PetscFunctionBegin;
1112d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1113e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1114e7e72b3dSBarry Smith     if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1115e7e72b3dSBarry Smith     else {
1116d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1117289bc588SBarry Smith         for (k=0; k<j; k++) {
11181b807ce4Svictorle           tmp = v[j + k*M];
11191b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11201b807ce4Svictorle           v[k + j*M] = tmp;
1121289bc588SBarry Smith         }
1122289bc588SBarry Smith       }
1123d64ed03dSBarry Smith     }
11243a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1125d3e5ee88SLois Curfman McInnes     Mat          tmat;
1126ec8511deSBarry Smith     Mat_SeqDense *tmatd;
112787828ca2SBarry Smith     PetscScalar  *v2;
1128ea709b57SSatish Balay 
1129fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11307adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1131d0f46423SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
11327adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11335c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1134fc4dec0aSBarry Smith     } else {
1135fc4dec0aSBarry Smith       tmat = *matout;
1136fc4dec0aSBarry Smith     }
1137ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11380de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1139d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11401b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1141d3e5ee88SLois Curfman McInnes     }
11426d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11436d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1144d3e5ee88SLois Curfman McInnes     *matout = tmat;
114548b35521SBarry Smith   }
11463a40ed3dSBarry Smith   PetscFunctionReturn(0);
1147289bc588SBarry Smith }
1148289bc588SBarry Smith 
11494a2ae208SSatish Balay #undef __FUNCT__
11504a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1151ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool  *flg)
1152289bc588SBarry Smith {
1153c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1154c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
115513f74950SBarry Smith   PetscInt     i,j;
115687828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11579ea5d5aeSSatish Balay 
11583a40ed3dSBarry Smith   PetscFunctionBegin;
1159d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1160d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1161d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
11621b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1163d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
11643a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11651b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11661b807ce4Svictorle     }
1167289bc588SBarry Smith   }
116877c4ece6SBarry Smith   *flg = PETSC_TRUE;
11693a40ed3dSBarry Smith   PetscFunctionReturn(0);
1170289bc588SBarry Smith }
1171289bc588SBarry Smith 
11724a2ae208SSatish Balay #undef __FUNCT__
11734a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1174dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1175289bc588SBarry Smith {
1176c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1177dfbe8321SBarry Smith   PetscErrorCode ierr;
117813f74950SBarry Smith   PetscInt       i,n,len;
117987828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
118044cd7ae7SLois Curfman McInnes 
11813a40ed3dSBarry Smith   PetscFunctionBegin;
11822dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11837a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11841ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1185d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1186e32f2f54SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
118744cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11881b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1189289bc588SBarry Smith   }
11901ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11913a40ed3dSBarry Smith   PetscFunctionReturn(0);
1192289bc588SBarry Smith }
1193289bc588SBarry Smith 
11944a2ae208SSatish Balay #undef __FUNCT__
11954a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1196dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1197289bc588SBarry Smith {
1198c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
119987828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1200dfbe8321SBarry Smith   PetscErrorCode ierr;
1201d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
120255659b69SBarry Smith 
12033a40ed3dSBarry Smith   PetscFunctionBegin;
120428988994SBarry Smith   if (ll) {
12057a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12061ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1207e32f2f54SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1208da3a660dSBarry Smith     for (i=0; i<m; i++) {
1209da3a660dSBarry Smith       x = l[i];
1210da3a660dSBarry Smith       v = mat->v + i;
1211da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1212da3a660dSBarry Smith     }
12131ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1214efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1215da3a660dSBarry Smith   }
121628988994SBarry Smith   if (rr) {
12177a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12181ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1219e32f2f54SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1220da3a660dSBarry Smith     for (i=0; i<n; i++) {
1221da3a660dSBarry Smith       x = r[i];
1222da3a660dSBarry Smith       v = mat->v + i*m;
1223da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1224da3a660dSBarry Smith     }
12251ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1226efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1227da3a660dSBarry Smith   }
12283a40ed3dSBarry Smith   PetscFunctionReturn(0);
1229289bc588SBarry Smith }
1230289bc588SBarry Smith 
12314a2ae208SSatish Balay #undef __FUNCT__
12324a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1233dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1234289bc588SBarry Smith {
1235c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
123687828ca2SBarry Smith   PetscScalar  *v = mat->v;
1237329f5518SBarry Smith   PetscReal    sum = 0.0;
1238d0f46423SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1239efee365bSSatish Balay   PetscErrorCode ierr;
124055659b69SBarry Smith 
12413a40ed3dSBarry Smith   PetscFunctionBegin;
1242289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1243a5ce6ee0Svictorle     if (lda>m) {
1244d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1245a5ce6ee0Svictorle 	v = mat->v+j*lda;
1246a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1247a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1248a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1249a5ce6ee0Svictorle #else
1250a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1251a5ce6ee0Svictorle #endif
1252a5ce6ee0Svictorle 	}
1253a5ce6ee0Svictorle       }
1254a5ce6ee0Svictorle     } else {
1255d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1256aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1257329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1258289bc588SBarry Smith #else
1259289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1260289bc588SBarry Smith #endif
1261289bc588SBarry Smith       }
1262a5ce6ee0Svictorle     }
1263064f8208SBarry Smith     *nrm = sqrt(sum);
1264dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12653a40ed3dSBarry Smith   } else if (type == NORM_1) {
1266064f8208SBarry Smith     *nrm = 0.0;
1267d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
12681b807ce4Svictorle       v = mat->v + j*mat->lda;
1269289bc588SBarry Smith       sum = 0.0;
1270d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
127133a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1272289bc588SBarry Smith       }
1273064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1274289bc588SBarry Smith     }
1275d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12763a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1277064f8208SBarry Smith     *nrm = 0.0;
1278d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1279289bc588SBarry Smith       v = mat->v + j;
1280289bc588SBarry Smith       sum = 0.0;
1281d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
12821b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1283289bc588SBarry Smith       }
1284064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1285289bc588SBarry Smith     }
1286d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1287e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
12883a40ed3dSBarry Smith   PetscFunctionReturn(0);
1289289bc588SBarry Smith }
1290289bc588SBarry Smith 
12914a2ae208SSatish Balay #undef __FUNCT__
12924a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1293ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool  flg)
1294289bc588SBarry Smith {
1295c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
129663ba0a88SBarry Smith   PetscErrorCode ierr;
129767e560aaSBarry Smith 
12983a40ed3dSBarry Smith   PetscFunctionBegin;
1299b5a2b587SKris Buschelman   switch (op) {
1300b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13014e0d8c25SBarry Smith     aij->roworiented = flg;
1302b5a2b587SKris Buschelman     break;
1303512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1304b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13053971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13064e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1307b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1308b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
130977e54ba9SKris Buschelman   case MAT_SYMMETRIC:
131077e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13119a4540c5SBarry Smith   case MAT_HERMITIAN:
13129a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1313600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1314290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
131577e54ba9SKris Buschelman     break;
1316b5a2b587SKris Buschelman   default:
1317e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13183a40ed3dSBarry Smith   }
13193a40ed3dSBarry Smith   PetscFunctionReturn(0);
1320289bc588SBarry Smith }
1321289bc588SBarry Smith 
13224a2ae208SSatish Balay #undef __FUNCT__
13234a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1324dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13256f0a148fSBarry Smith {
1326ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13276849ba73SBarry Smith   PetscErrorCode ierr;
1328d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13293a40ed3dSBarry Smith 
13303a40ed3dSBarry Smith   PetscFunctionBegin;
1331a5ce6ee0Svictorle   if (lda>m) {
1332d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1333a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1334a5ce6ee0Svictorle     }
1335a5ce6ee0Svictorle   } else {
1336d0f46423SBarry Smith     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1337a5ce6ee0Svictorle   }
13383a40ed3dSBarry Smith   PetscFunctionReturn(0);
13396f0a148fSBarry Smith }
13406f0a148fSBarry Smith 
13414a2ae208SSatish Balay #undef __FUNCT__
13424a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1343f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13446f0a148fSBarry Smith {
1345ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1346d0f46423SBarry Smith   PetscInt       n = A->cmap->n,i,j;
134787828ca2SBarry Smith   PetscScalar    *slot;
134855659b69SBarry Smith 
13493a40ed3dSBarry Smith   PetscFunctionBegin;
13506f0a148fSBarry Smith   for (i=0; i<N; i++) {
13516f0a148fSBarry Smith     slot = l->v + rows[i];
13526f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13536f0a148fSBarry Smith   }
1354f4df32b1SMatthew Knepley   if (diag != 0.0) {
13556f0a148fSBarry Smith     for (i=0; i<N; i++) {
13566f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1357f4df32b1SMatthew Knepley       *slot = diag;
13586f0a148fSBarry Smith     }
13596f0a148fSBarry Smith   }
13603a40ed3dSBarry Smith   PetscFunctionReturn(0);
13616f0a148fSBarry Smith }
1362557bce09SLois Curfman McInnes 
13634a2ae208SSatish Balay #undef __FUNCT__
13644a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1365dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
136664e87e97SBarry Smith {
1367c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13683a40ed3dSBarry Smith 
13693a40ed3dSBarry Smith   PetscFunctionBegin;
1370e32f2f54SBarry 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");
137164e87e97SBarry Smith   *array = mat->v;
13723a40ed3dSBarry Smith   PetscFunctionReturn(0);
137364e87e97SBarry Smith }
13740754003eSLois Curfman McInnes 
13754a2ae208SSatish Balay #undef __FUNCT__
13764a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1377dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1378ff14e315SSatish Balay {
13793a40ed3dSBarry Smith   PetscFunctionBegin;
138009b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13813a40ed3dSBarry Smith   PetscFunctionReturn(0);
1382ff14e315SSatish Balay }
13830754003eSLois Curfman McInnes 
13844a2ae208SSatish Balay #undef __FUNCT__
13854a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
138613f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13870754003eSLois Curfman McInnes {
1388c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13896849ba73SBarry Smith   PetscErrorCode ierr;
13905d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
13915d0c19d7SBarry Smith   const PetscInt *irow,*icol;
139287828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13930754003eSLois Curfman McInnes   Mat            newmat;
13940754003eSLois Curfman McInnes 
13953a40ed3dSBarry Smith   PetscFunctionBegin;
139678b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
139778b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1398e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1399e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14000754003eSLois Curfman McInnes 
1401182d2002SSatish Balay   /* Check submatrixcall */
1402182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
140313f74950SBarry Smith     PetscInt n_cols,n_rows;
1404182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
140521a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
140621a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
1407c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
140821a2c019SBarry Smith     }
1409182d2002SSatish Balay     newmat = *B;
1410182d2002SSatish Balay   } else {
14110754003eSLois Curfman McInnes     /* Create and fill new matrix */
14127adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1413f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14147adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14155c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1416182d2002SSatish Balay   }
1417182d2002SSatish Balay 
1418182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1419182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1420182d2002SSatish Balay 
1421182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14226de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1423182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1424182d2002SSatish Balay       *bv++ = av[irow[j]];
14250754003eSLois Curfman McInnes     }
14260754003eSLois Curfman McInnes   }
1427182d2002SSatish Balay 
1428182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14296d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14306d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14310754003eSLois Curfman McInnes 
14320754003eSLois Curfman McInnes   /* Free work space */
143378b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
143478b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1435182d2002SSatish Balay   *B = newmat;
14363a40ed3dSBarry Smith   PetscFunctionReturn(0);
14370754003eSLois Curfman McInnes }
14380754003eSLois Curfman McInnes 
14394a2ae208SSatish Balay #undef __FUNCT__
14404a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
144113f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1442905e6a2fSBarry Smith {
14436849ba73SBarry Smith   PetscErrorCode ierr;
144413f74950SBarry Smith   PetscInt       i;
1445905e6a2fSBarry Smith 
14463a40ed3dSBarry Smith   PetscFunctionBegin;
1447905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1448b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1449905e6a2fSBarry Smith   }
1450905e6a2fSBarry Smith 
1451905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14526a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1453905e6a2fSBarry Smith   }
14543a40ed3dSBarry Smith   PetscFunctionReturn(0);
1455905e6a2fSBarry Smith }
1456905e6a2fSBarry Smith 
14574a2ae208SSatish Balay #undef __FUNCT__
1458c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1459c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1460c0aa2d19SHong Zhang {
1461c0aa2d19SHong Zhang   PetscFunctionBegin;
1462c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1463c0aa2d19SHong Zhang }
1464c0aa2d19SHong Zhang 
1465c0aa2d19SHong Zhang #undef __FUNCT__
1466c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1467c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1468c0aa2d19SHong Zhang {
1469c0aa2d19SHong Zhang   PetscFunctionBegin;
1470c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1471c0aa2d19SHong Zhang }
1472c0aa2d19SHong Zhang 
1473c0aa2d19SHong Zhang #undef __FUNCT__
14744a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1475dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14764b0e389bSBarry Smith {
14774b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14786849ba73SBarry Smith   PetscErrorCode ierr;
1479d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
14803a40ed3dSBarry Smith 
14813a40ed3dSBarry Smith   PetscFunctionBegin;
148233f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
148333f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1484cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14853a40ed3dSBarry Smith     PetscFunctionReturn(0);
14863a40ed3dSBarry Smith   }
1487e32f2f54SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1488a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14890dbb7854Svictorle     for (j=0; j<n; j++) {
1490a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1491a5ce6ee0Svictorle     }
1492a5ce6ee0Svictorle   } else {
1493d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1494a5ce6ee0Svictorle   }
1495273d9f13SBarry Smith   PetscFunctionReturn(0);
1496273d9f13SBarry Smith }
1497273d9f13SBarry Smith 
14984a2ae208SSatish Balay #undef __FUNCT__
14994a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1500dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1501273d9f13SBarry Smith {
1502dfbe8321SBarry Smith   PetscErrorCode ierr;
1503273d9f13SBarry Smith 
1504273d9f13SBarry Smith   PetscFunctionBegin;
1505273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15063a40ed3dSBarry Smith   PetscFunctionReturn(0);
15074b0e389bSBarry Smith }
15084b0e389bSBarry Smith 
1509284134d9SBarry Smith #undef __FUNCT__
1510284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1511284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1512284134d9SBarry Smith {
1513284134d9SBarry Smith   PetscFunctionBegin;
151421a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1515284134d9SBarry Smith   m = PetscMax(m,M);
1516284134d9SBarry Smith   n = PetscMax(n,N);
1517a868139aSShri Abhyankar 
151886d161a7SShri 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);
151986d161a7SShri 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);
152086d161a7SShri Abhyankar   */
1521dc5cefdeSJed Brown   A->rmap->n = A->rmap->N = m;
1522d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
1523284134d9SBarry Smith   PetscFunctionReturn(0);
1524284134d9SBarry Smith }
1525170fe5c8SBarry Smith 
1526ba337c44SJed Brown #undef __FUNCT__
1527ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense"
1528ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A)
1529ba337c44SJed Brown {
1530ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1531ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1532ba337c44SJed Brown   PetscScalar    *aa = a->v;
1533ba337c44SJed Brown 
1534ba337c44SJed Brown   PetscFunctionBegin;
1535ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
1536ba337c44SJed Brown   PetscFunctionReturn(0);
1537ba337c44SJed Brown }
1538ba337c44SJed Brown 
1539ba337c44SJed Brown #undef __FUNCT__
1540ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense"
1541ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A)
1542ba337c44SJed Brown {
1543ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1544ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1545ba337c44SJed Brown   PetscScalar    *aa = a->v;
1546ba337c44SJed Brown 
1547ba337c44SJed Brown   PetscFunctionBegin;
1548ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1549ba337c44SJed Brown   PetscFunctionReturn(0);
1550ba337c44SJed Brown }
1551ba337c44SJed Brown 
1552ba337c44SJed Brown #undef __FUNCT__
1553ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense"
1554ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
1555ba337c44SJed Brown {
1556ba337c44SJed Brown   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1557ba337c44SJed Brown   PetscInt       i,nz = A->rmap->n*A->cmap->n;
1558ba337c44SJed Brown   PetscScalar    *aa = a->v;
1559ba337c44SJed Brown 
1560ba337c44SJed Brown   PetscFunctionBegin;
1561ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1562ba337c44SJed Brown   PetscFunctionReturn(0);
1563ba337c44SJed Brown }
1564284134d9SBarry Smith 
1565a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1566a9fe9ddaSSatish Balay #undef __FUNCT__
1567a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1568a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1569a9fe9ddaSSatish Balay {
1570a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1571a9fe9ddaSSatish Balay 
1572a9fe9ddaSSatish Balay   PetscFunctionBegin;
1573a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1574a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1575a9fe9ddaSSatish Balay   }
1576a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1577a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1578a9fe9ddaSSatish Balay }
1579a9fe9ddaSSatish Balay 
1580a9fe9ddaSSatish Balay #undef __FUNCT__
1581a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1582a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1583a9fe9ddaSSatish Balay {
1584ee16a9a1SHong Zhang   PetscErrorCode ierr;
1585d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1586ee16a9a1SHong Zhang   Mat            Cmat;
1587a9fe9ddaSSatish Balay 
1588ee16a9a1SHong Zhang   PetscFunctionBegin;
1589e32f2f54SBarry 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);
1590ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1591ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1592ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1593ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1594ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1595ee16a9a1SHong Zhang   *C = Cmat;
1596ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1597ee16a9a1SHong Zhang }
1598a9fe9ddaSSatish Balay 
159998a3b096SSatish Balay #undef __FUNCT__
1600a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1601a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1602a9fe9ddaSSatish Balay {
1603a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1604a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1605a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16060805154bSBarry Smith   PetscBLASInt   m,n,k;
1607a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1608a9fe9ddaSSatish Balay 
1609a9fe9ddaSSatish Balay   PetscFunctionBegin;
1610d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1611d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1612d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1613a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1614a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1615a9fe9ddaSSatish Balay }
1616a9fe9ddaSSatish Balay 
1617a9fe9ddaSSatish Balay #undef __FUNCT__
1618a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1619a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1620a9fe9ddaSSatish Balay {
1621a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1622a9fe9ddaSSatish Balay 
1623a9fe9ddaSSatish Balay   PetscFunctionBegin;
1624a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1625a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1626a9fe9ddaSSatish Balay   }
1627a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1628a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1629a9fe9ddaSSatish Balay }
1630a9fe9ddaSSatish Balay 
1631a9fe9ddaSSatish Balay #undef __FUNCT__
1632a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1633a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1634a9fe9ddaSSatish Balay {
1635ee16a9a1SHong Zhang   PetscErrorCode ierr;
1636d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1637ee16a9a1SHong Zhang   Mat            Cmat;
1638a9fe9ddaSSatish Balay 
1639ee16a9a1SHong Zhang   PetscFunctionBegin;
1640e32f2f54SBarry 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);
1641ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1642ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1643ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1644ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1645ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1646ee16a9a1SHong Zhang   *C = Cmat;
1647ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1648ee16a9a1SHong Zhang }
1649a9fe9ddaSSatish Balay 
1650a9fe9ddaSSatish Balay #undef __FUNCT__
1651a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1652a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1653a9fe9ddaSSatish Balay {
1654a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1655a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1656a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16570805154bSBarry Smith   PetscBLASInt   m,n,k;
1658a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1659a9fe9ddaSSatish Balay 
1660a9fe9ddaSSatish Balay   PetscFunctionBegin;
1661d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1662d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1663d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
16642fbe02b9SBarry Smith   /*
16652fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
16662fbe02b9SBarry Smith   */
1667a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1668a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1669a9fe9ddaSSatish Balay }
1670985db425SBarry Smith 
1671985db425SBarry Smith #undef __FUNCT__
1672985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1673985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1674985db425SBarry Smith {
1675985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1676985db425SBarry Smith   PetscErrorCode ierr;
1677d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1678985db425SBarry Smith   PetscScalar    *x;
1679985db425SBarry Smith   MatScalar      *aa = a->v;
1680985db425SBarry Smith 
1681985db425SBarry Smith   PetscFunctionBegin;
1682e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1683985db425SBarry Smith 
1684985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1685985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1686985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1687e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1688985db425SBarry Smith   for (i=0; i<m; i++) {
1689985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1690985db425SBarry Smith     for (j=1; j<n; j++){
1691985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1692985db425SBarry Smith     }
1693985db425SBarry Smith   }
1694985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1695985db425SBarry Smith   PetscFunctionReturn(0);
1696985db425SBarry Smith }
1697985db425SBarry Smith 
1698985db425SBarry Smith #undef __FUNCT__
1699985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1700985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1701985db425SBarry Smith {
1702985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1703985db425SBarry Smith   PetscErrorCode ierr;
1704d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1705985db425SBarry Smith   PetscScalar    *x;
1706985db425SBarry Smith   PetscReal      atmp;
1707985db425SBarry Smith   MatScalar      *aa = a->v;
1708985db425SBarry Smith 
1709985db425SBarry Smith   PetscFunctionBegin;
1710e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1711985db425SBarry Smith 
1712985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1713985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1714985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1715e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1716985db425SBarry Smith   for (i=0; i<m; i++) {
17179189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1718985db425SBarry Smith     for (j=1; j<n; j++){
1719985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1720985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1721985db425SBarry Smith     }
1722985db425SBarry Smith   }
1723985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1724985db425SBarry Smith   PetscFunctionReturn(0);
1725985db425SBarry Smith }
1726985db425SBarry Smith 
1727985db425SBarry Smith #undef __FUNCT__
1728985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1729985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1730985db425SBarry Smith {
1731985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1732985db425SBarry Smith   PetscErrorCode ierr;
1733d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1734985db425SBarry Smith   PetscScalar    *x;
1735985db425SBarry Smith   MatScalar      *aa = a->v;
1736985db425SBarry Smith 
1737985db425SBarry Smith   PetscFunctionBegin;
1738e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1739985db425SBarry Smith 
1740985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1741985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1742985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1743e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1744985db425SBarry Smith   for (i=0; i<m; i++) {
1745985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1746985db425SBarry Smith     for (j=1; j<n; j++){
1747985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1748985db425SBarry Smith     }
1749985db425SBarry Smith   }
1750985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1751985db425SBarry Smith   PetscFunctionReturn(0);
1752985db425SBarry Smith }
1753985db425SBarry Smith 
17548d0534beSBarry Smith #undef __FUNCT__
17558d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
17568d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
17578d0534beSBarry Smith {
17588d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17598d0534beSBarry Smith   PetscErrorCode ierr;
17608d0534beSBarry Smith   PetscScalar    *x;
17618d0534beSBarry Smith 
17628d0534beSBarry Smith   PetscFunctionBegin;
1763e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
17648d0534beSBarry Smith 
17658d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1766d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
17678d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
17688d0534beSBarry Smith   PetscFunctionReturn(0);
17698d0534beSBarry Smith }
17708d0534beSBarry Smith 
1771289bc588SBarry Smith /* -------------------------------------------------------------------*/
1772a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1773905e6a2fSBarry Smith        MatGetRow_SeqDense,
1774905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1775905e6a2fSBarry Smith        MatMult_SeqDense,
177697304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
17777c922b88SBarry Smith        MatMultTranspose_SeqDense,
17787c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1779db4efbfdSBarry Smith        0,
1780db4efbfdSBarry Smith        0,
1781db4efbfdSBarry Smith        0,
1782db4efbfdSBarry Smith /*10*/ 0,
1783905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1784905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
178541f059aeSBarry Smith        MatSOR_SeqDense,
1786ec8511deSBarry Smith        MatTranspose_SeqDense,
178797304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1788905e6a2fSBarry Smith        MatEqual_SeqDense,
1789905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1790905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1791905e6a2fSBarry Smith        MatNorm_SeqDense,
1792c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1793c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1794905e6a2fSBarry Smith        MatSetOption_SeqDense,
1795905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1796d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense,
1797db4efbfdSBarry Smith        0,
1798db4efbfdSBarry Smith        0,
1799db4efbfdSBarry Smith        0,
1800db4efbfdSBarry Smith        0,
1801d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense,
1802273d9f13SBarry Smith        0,
1803905e6a2fSBarry Smith        0,
1804905e6a2fSBarry Smith        MatGetArray_SeqDense,
1805905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1806d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense,
1807a5ae1ecdSBarry Smith        0,
1808a5ae1ecdSBarry Smith        0,
1809a5ae1ecdSBarry Smith        0,
1810a5ae1ecdSBarry Smith        0,
1811d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense,
1812a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1813a5ae1ecdSBarry Smith        0,
18144b0e389bSBarry Smith        MatGetValues_SeqDense,
1815a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1816d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense,
1817a5ae1ecdSBarry Smith        MatScale_SeqDense,
1818a5ae1ecdSBarry Smith        0,
1819a5ae1ecdSBarry Smith        0,
1820a5ae1ecdSBarry Smith        0,
1821d519adbfSMatthew Knepley /*49*/ 0,
1822a5ae1ecdSBarry Smith        0,
1823a5ae1ecdSBarry Smith        0,
1824a5ae1ecdSBarry Smith        0,
1825a5ae1ecdSBarry Smith        0,
1826d519adbfSMatthew Knepley /*54*/ 0,
1827a5ae1ecdSBarry Smith        0,
1828a5ae1ecdSBarry Smith        0,
1829a5ae1ecdSBarry Smith        0,
1830a5ae1ecdSBarry Smith        0,
1831d519adbfSMatthew Knepley /*59*/ 0,
1832e03a110bSBarry Smith        MatDestroy_SeqDense,
1833e03a110bSBarry Smith        MatView_SeqDense,
1834357abbc8SBarry Smith        0,
183597304618SKris Buschelman        0,
1836d519adbfSMatthew Knepley /*64*/ 0,
183797304618SKris Buschelman        0,
183897304618SKris Buschelman        0,
183997304618SKris Buschelman        0,
184097304618SKris Buschelman        0,
1841d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense,
184297304618SKris Buschelman        0,
184397304618SKris Buschelman        0,
184497304618SKris Buschelman        0,
184597304618SKris Buschelman        0,
1846d519adbfSMatthew Knepley /*74*/ 0,
184797304618SKris Buschelman        0,
184897304618SKris Buschelman        0,
184997304618SKris Buschelman        0,
185097304618SKris Buschelman        0,
1851d519adbfSMatthew Knepley /*79*/ 0,
185297304618SKris Buschelman        0,
185397304618SKris Buschelman        0,
185497304618SKris Buschelman        0,
18555bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense,
1856865e5f61SKris Buschelman        0,
18571cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1858865e5f61SKris Buschelman        0,
1859865e5f61SKris Buschelman        0,
1860865e5f61SKris Buschelman        0,
1861d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense,
1862a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1863a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1864865e5f61SKris Buschelman        0,
1865865e5f61SKris Buschelman        0,
1866d519adbfSMatthew Knepley /*94*/ 0,
1867a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1868a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1869a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1870284134d9SBarry Smith        0,
1871d519adbfSMatthew Knepley /*99*/ 0,
1872284134d9SBarry Smith        0,
1873284134d9SBarry Smith        0,
1874ba337c44SJed Brown        MatConjugate_SeqDense,
1875985db425SBarry Smith        MatSetSizes_SeqDense,
1876ba337c44SJed Brown /*104*/0,
1877ba337c44SJed Brown        MatRealPart_SeqDense,
1878ba337c44SJed Brown        MatImaginaryPart_SeqDense,
1879985db425SBarry Smith        0,
1880985db425SBarry Smith        0,
1881d519adbfSMatthew Knepley /*109*/0,
1882985db425SBarry Smith        0,
18838d0534beSBarry Smith        MatGetRowMin_SeqDense,
1884aabbc4fbSShri Abhyankar        MatGetColumnVector_SeqDense,
1885aabbc4fbSShri Abhyankar        0,
1886aabbc4fbSShri Abhyankar /*114*/0,
1887aabbc4fbSShri Abhyankar        0,
1888aabbc4fbSShri Abhyankar        0,
1889aabbc4fbSShri Abhyankar        0,
1890aabbc4fbSShri Abhyankar        0,
1891aabbc4fbSShri Abhyankar /*119*/0,
1892aabbc4fbSShri Abhyankar        0,
1893aabbc4fbSShri Abhyankar        0,
18945bba2384SShri Abhyankar        0
1895985db425SBarry Smith };
189690ace30eSBarry Smith 
18974a2ae208SSatish Balay #undef __FUNCT__
18984a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18994b828684SBarry Smith /*@C
1900fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1901d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1902d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1903289bc588SBarry Smith 
1904db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1905db81eaa0SLois Curfman McInnes 
190620563c6bSBarry Smith    Input Parameters:
1907db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
19080c775827SLois Curfman McInnes .  m - number of rows
190918f449edSLois Curfman McInnes .  n - number of columns
1910c0235b3cSMatthew Knepley -  data - optional location of matrix data in column major order.  Set data=PETSC_NULL for PETSc
1911dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
191220563c6bSBarry Smith 
191320563c6bSBarry Smith    Output Parameter:
191444cd7ae7SLois Curfman McInnes .  A - the matrix
191520563c6bSBarry Smith 
1916b259b22eSLois Curfman McInnes    Notes:
191718f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
191818f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1919b4fd4287SBarry Smith    set data=PETSC_NULL.
192018f449edSLois Curfman McInnes 
1921027ccd11SLois Curfman McInnes    Level: intermediate
1922027ccd11SLois Curfman McInnes 
1923dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1924d65003e9SLois Curfman McInnes 
1925db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
192620563c6bSBarry Smith @*/
1927be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1928289bc588SBarry Smith {
1929dfbe8321SBarry Smith   PetscErrorCode ierr;
19303b2fbd54SBarry Smith 
19313a40ed3dSBarry Smith   PetscFunctionBegin;
1932f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1933f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1934273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1935273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1936273d9f13SBarry Smith   PetscFunctionReturn(0);
1937273d9f13SBarry Smith }
1938273d9f13SBarry Smith 
19394a2ae208SSatish Balay #undef __FUNCT__
1940afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1941273d9f13SBarry Smith /*@C
1942273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1943273d9f13SBarry Smith 
1944273d9f13SBarry Smith    Collective on MPI_Comm
1945273d9f13SBarry Smith 
1946273d9f13SBarry Smith    Input Parameters:
1947273d9f13SBarry Smith +  A - the matrix
1948273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1949273d9f13SBarry Smith 
1950273d9f13SBarry Smith    Notes:
1951273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1952273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1953284134d9SBarry Smith    need not call this routine.
1954273d9f13SBarry Smith 
1955273d9f13SBarry Smith    Level: intermediate
1956273d9f13SBarry Smith 
1957273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1958273d9f13SBarry Smith 
1959867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
1960867c911aSBarry Smith 
1961273d9f13SBarry Smith @*/
1962be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1963273d9f13SBarry Smith {
1964dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1965a23d5eceSKris Buschelman 
1966a23d5eceSKris Buschelman   PetscFunctionBegin;
1967a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1968a23d5eceSKris Buschelman   if (f) {
1969a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1970a23d5eceSKris Buschelman   }
1971a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1972a23d5eceSKris Buschelman }
1973a23d5eceSKris Buschelman 
1974a23d5eceSKris Buschelman EXTERN_C_BEGIN
1975a23d5eceSKris Buschelman #undef __FUNCT__
1976afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1977be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1978a23d5eceSKris Buschelman {
1979273d9f13SBarry Smith   Mat_SeqDense   *b;
1980dfbe8321SBarry Smith   PetscErrorCode ierr;
1981273d9f13SBarry Smith 
1982273d9f13SBarry Smith   PetscFunctionBegin;
1983273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1984a868139aSShri Abhyankar 
198534ef9618SShri Abhyankar   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
198634ef9618SShri Abhyankar   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
198734ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
198834ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
198934ef9618SShri Abhyankar 
1990273d9f13SBarry Smith   b       = (Mat_SeqDense*)B->data;
199186d161a7SShri Abhyankar   b->Mmax = B->rmap->n;
199286d161a7SShri Abhyankar   b->Nmax = B->cmap->n;
199386d161a7SShri Abhyankar   if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n;
199486d161a7SShri Abhyankar 
19959e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
19969e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
19975afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1998284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1999284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
20009e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
2001273d9f13SBarry Smith   } else { /* user-allocated storage */
20029e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
2003273d9f13SBarry Smith     b->v          = data;
2004273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
2005273d9f13SBarry Smith   }
20060450473dSBarry Smith   B->assembled = PETSC_TRUE;
2007273d9f13SBarry Smith   PetscFunctionReturn(0);
2008273d9f13SBarry Smith }
2009a23d5eceSKris Buschelman EXTERN_C_END
2010273d9f13SBarry Smith 
20111b807ce4Svictorle #undef __FUNCT__
20121b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
20131b807ce4Svictorle /*@C
20141b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
20151b807ce4Svictorle 
20161b807ce4Svictorle   Input parameter:
20171b807ce4Svictorle + A - the matrix
20181b807ce4Svictorle - lda - the leading dimension
20191b807ce4Svictorle 
20201b807ce4Svictorle   Notes:
2021867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
20221b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
20231b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
20241b807ce4Svictorle 
20251b807ce4Svictorle   Level: intermediate
20261b807ce4Svictorle 
20271b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
20281b807ce4Svictorle 
2029284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2030867c911aSBarry Smith 
20311b807ce4Svictorle @*/
2032be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
20331b807ce4Svictorle {
20341b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
203521a2c019SBarry Smith 
20361b807ce4Svictorle   PetscFunctionBegin;
2037e32f2f54SBarry 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);
20381b807ce4Svictorle   b->lda       = lda;
203921a2c019SBarry Smith   b->changelda = PETSC_FALSE;
204021a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
20411b807ce4Svictorle   PetscFunctionReturn(0);
20421b807ce4Svictorle }
20431b807ce4Svictorle 
20440bad9183SKris Buschelman /*MC
2045fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
20460bad9183SKris Buschelman 
20470bad9183SKris Buschelman    Options Database Keys:
20480bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
20490bad9183SKris Buschelman 
20500bad9183SKris Buschelman   Level: beginner
20510bad9183SKris Buschelman 
205289665df3SBarry Smith .seealso: MatCreateSeqDense()
205389665df3SBarry Smith 
20540bad9183SKris Buschelman M*/
20550bad9183SKris Buschelman 
2056273d9f13SBarry Smith EXTERN_C_BEGIN
20574a2ae208SSatish Balay #undef __FUNCT__
20584a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
2059be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2060273d9f13SBarry Smith {
2061273d9f13SBarry Smith   Mat_SeqDense   *b;
2062dfbe8321SBarry Smith   PetscErrorCode ierr;
20637c334f02SBarry Smith   PetscMPIInt    size;
2064273d9f13SBarry Smith 
2065273d9f13SBarry Smith   PetscFunctionBegin;
20667adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
2067e32f2f54SBarry Smith   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
206855659b69SBarry Smith 
206938f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2070549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
207190f02eecSBarry Smith   B->mapping      = 0;
207244cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
207318f449edSLois Curfman McInnes 
207444cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2075273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2076273d9f13SBarry Smith   b->v            = 0;
207721a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
20784e220ebcSLois Curfman McInnes 
2079b24902e0SBarry Smith 
2080ec1065edSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
2081b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2082b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2083a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2084a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2085a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
20864ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
20874ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
20884ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
20894ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
20904ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
20914ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
20924ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
20934ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
20944ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
209517667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
20963a40ed3dSBarry Smith   PetscFunctionReturn(0);
2097289bc588SBarry Smith }
2098273d9f13SBarry Smith EXTERN_C_END
2099