xref: /petsc/src/mat/impls/dense/seq/dense.c (revision c0235b3c8bee801ea2f32fe0769d7268ab2d9066)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
77c4f633dSBarry Smith #include "../src/mat/impls/dense/seq/dense.h"
8b4c862e0SSatish Balay #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
170805154bSBarry Smith   PetscBLASInt   N,m,ldax,lday,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21d0f46423SBarry Smith   N    = PetscBLASIntCast(X->rmap->n*X->cmap->n);
22d0f46423SBarry Smith   m    = PetscBLASIntCast(X->rmap->n);
230805154bSBarry Smith   ldax = PetscBLASIntCast(x->lda);
240805154bSBarry Smith   lday = PetscBLASIntCast(y->lda);
25a5ce6ee0Svictorle   if (ldax>m || lday>m) {
26d0f46423SBarry Smith     for (j=0; j<X->cmap->n; j++) {
27f4df32b1SMatthew Knepley       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
28a5ce6ee0Svictorle     }
29a5ce6ee0Svictorle   } else {
30f4df32b1SMatthew Knepley     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
31a5ce6ee0Svictorle   }
320450473dSBarry Smith   ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr);
333a40ed3dSBarry Smith   PetscFunctionReturn(0);
341987afe7SBarry Smith }
351987afe7SBarry Smith 
364a2ae208SSatish Balay #undef __FUNCT__
374a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
38dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
39289bc588SBarry Smith {
40d0f46423SBarry Smith   PetscInt     N = A->rmap->n*A->cmap->n;
413a40ed3dSBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
434e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
444e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
456de62eeeSBarry Smith   info->nz_used           = (double)N;
466de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
474e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
484e220ebcSLois Curfman McInnes   info->mallocs           = 0;
497adad957SLisandro Dalcin   info->memory            = ((PetscObject)A)->mem;
504e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
514e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
524e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
54289bc588SBarry Smith }
55289bc588SBarry Smith 
564a2ae208SSatish Balay #undef __FUNCT__
574a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
58f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
5980cd9d93SLois Curfman McInnes {
60273d9f13SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
61f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
62efee365bSSatish Balay   PetscErrorCode ierr;
630805154bSBarry Smith   PetscBLASInt   one = 1,j,nz,lda = PetscBLASIntCast(a->lda);
6480cd9d93SLois Curfman McInnes 
653a40ed3dSBarry Smith   PetscFunctionBegin;
66d0f46423SBarry Smith   if (lda>A->rmap->n) {
67d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n);
68d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
69f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
72d0f46423SBarry Smith     nz = PetscBLASIntCast(A->rmap->n*A->cmap->n);
73f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
74a5ce6ee0Svictorle   }
75efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
7780cd9d93SLois Curfman McInnes }
7880cd9d93SLois Curfman McInnes 
791cbb95d3SBarry Smith #undef __FUNCT__
801cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense"
811cbb95d3SBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscTruth *fl)
821cbb95d3SBarry Smith {
831cbb95d3SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
84d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,N;
851cbb95d3SBarry Smith   PetscScalar    *v = a->v;
861cbb95d3SBarry Smith 
871cbb95d3SBarry Smith   PetscFunctionBegin;
881cbb95d3SBarry Smith   *fl = PETSC_FALSE;
89d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0);
901cbb95d3SBarry Smith   N = a->lda;
911cbb95d3SBarry Smith 
921cbb95d3SBarry Smith   for (i=0; i<m; i++) {
931cbb95d3SBarry Smith     for (j=i+1; j<m; j++) {
941cbb95d3SBarry Smith       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
951cbb95d3SBarry Smith     }
961cbb95d3SBarry Smith   }
971cbb95d3SBarry Smith   *fl = PETSC_TRUE;
981cbb95d3SBarry Smith   PetscFunctionReturn(0);
991cbb95d3SBarry Smith }
1001cbb95d3SBarry Smith 
101b24902e0SBarry Smith #undef __FUNCT__
102b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense"
103719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
104b24902e0SBarry Smith {
105b24902e0SBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
106b24902e0SBarry Smith   PetscErrorCode ierr;
107b24902e0SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
108b24902e0SBarry Smith 
109b24902e0SBarry Smith   PetscFunctionBegin;
110719d5645SBarry Smith   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
111b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
112b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
113d0f46423SBarry Smith     if (lda>A->rmap->n) {
114d0f46423SBarry Smith       m = A->rmap->n;
115d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
116b24902e0SBarry Smith 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
117b24902e0SBarry Smith       }
118b24902e0SBarry Smith     } else {
119d0f46423SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
120b24902e0SBarry Smith     }
121b24902e0SBarry Smith   }
122b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
123b24902e0SBarry Smith   PetscFunctionReturn(0);
124b24902e0SBarry Smith }
125b24902e0SBarry Smith 
1264a2ae208SSatish Balay #undef __FUNCT__
1274a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
128dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12902cad45dSBarry Smith {
1306849ba73SBarry Smith   PetscErrorCode ierr;
13102cad45dSBarry Smith 
1323a40ed3dSBarry Smith   PetscFunctionBegin;
1335c9eb25fSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr);
134d0f46423SBarry Smith   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1355c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
136719d5645SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
137b24902e0SBarry Smith   PetscFunctionReturn(0);
138b24902e0SBarry Smith }
139b24902e0SBarry Smith 
1406ee01492SSatish Balay 
1410481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
142719d5645SBarry Smith 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1450481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
146289bc588SBarry Smith {
1474482741eSBarry Smith   MatFactorInfo  info;
148a093e273SMatthew Knepley   PetscErrorCode ierr;
1493a40ed3dSBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
151c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
152719d5645SBarry Smith   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154289bc588SBarry Smith }
1556ee01492SSatish Balay 
1560b4b3355SBarry Smith #undef __FUNCT__
1574a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
158dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
159289bc588SBarry Smith {
160c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1616849ba73SBarry Smith   PetscErrorCode ierr;
16287828ca2SBarry Smith   PetscScalar    *x,*y;
163d0f46423SBarry Smith   PetscBLASInt   one = 1,info,m = PetscBLASIntCast(A->rmap->n);
16467e560aaSBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
1661ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1671ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
168d0f46423SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1695c9eb25fSBarry Smith   if (A->factor == MAT_FACTOR_LU) {
170ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
17180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
172ae7cfcebSSatish Balay #else
17371044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
1744ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
175ae7cfcebSSatish Balay #endif
1765c9eb25fSBarry Smith   } else if (A->factor == MAT_FACTOR_CHOLESKY){
177ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
17880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
179ae7cfcebSSatish Balay #else
18071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
1814ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
182ae7cfcebSSatish Balay #endif
183289bc588SBarry Smith   }
18429bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
1851ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1861ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
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)
20780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
208ae7cfcebSSatish Balay #else
20971044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2104ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
211ae7cfcebSSatish Balay #endif
2127a97a34bSBarry Smith   } else {
213ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
21480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
215ae7cfcebSSatish Balay #else
21671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2174ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
218ae7cfcebSSatish Balay #endif
219da3a660dSBarry Smith   }
2201ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2211ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
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)
24980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
250ae7cfcebSSatish Balay #else
25171044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2522ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
253ae7cfcebSSatish Balay #endif
254a8c6a408SBarry Smith   } else {
255ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
25680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
257ae7cfcebSSatish Balay #else
25871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2592ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
260ae7cfcebSSatish Balay #endif
261da3a660dSBarry Smith   }
2622dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
2632dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
2641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2651ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
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)
29380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
294ae7cfcebSSatish Balay #else
29571044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2962ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
297ae7cfcebSSatish Balay #endif
2983a40ed3dSBarry Smith   } else {
299ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
30080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
301ae7cfcebSSatish Balay #else
30271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3032ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
304ae7cfcebSSatish Balay #endif
305da3a660dSBarry Smith   }
30690f02eecSBarry Smith   if (tmp) {
3072dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
30890f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3093a40ed3dSBarry Smith   } else {
3102dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
31190f02eecSBarry Smith   }
3121ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3131ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
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;
327db4efbfdSBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
328db4efbfdSBarry Smith #else
329db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
330db4efbfdSBarry Smith   PetscErrorCode ierr;
331db4efbfdSBarry Smith   PetscBLASInt   n,m,info;
332db4efbfdSBarry Smith 
333db4efbfdSBarry Smith   PetscFunctionBegin;
334db4efbfdSBarry Smith   n = PetscBLASIntCast(A->cmap->n);
335db4efbfdSBarry Smith   m = PetscBLASIntCast(A->rmap->n);
336db4efbfdSBarry Smith   if (!mat->pivots) {
337db4efbfdSBarry Smith     ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
338db4efbfdSBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
339db4efbfdSBarry Smith   }
340db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
341db4efbfdSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
342db4efbfdSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
343db4efbfdSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
344db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
345db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
346db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
347db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
348db4efbfdSBarry Smith   A->factor = MAT_FACTOR_LU;
349db4efbfdSBarry Smith 
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;
361db4efbfdSBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
362db4efbfdSBarry Smith #else
363db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
364db4efbfdSBarry Smith   PetscErrorCode ierr;
365db4efbfdSBarry Smith   PetscBLASInt   info,n = PetscBLASIntCast(A->cmap->n);
366db4efbfdSBarry Smith 
367db4efbfdSBarry Smith   PetscFunctionBegin;
368db4efbfdSBarry Smith   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
369db4efbfdSBarry Smith   mat->pivots = 0;
370db4efbfdSBarry Smith 
371db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
372db4efbfdSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
373db4efbfdSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
374db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
375db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
376db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
377db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
378db4efbfdSBarry Smith   A->factor = MAT_FACTOR_CHOLESKY;
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   }
435db4efbfdSBarry Smith   (*fact)->factor = 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;
46177431f27SBarry Smith   if (its <= 0) SETERRQ2(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;
534800995b7SMatthew Knepley   PetscInt       i, j;
5350805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5363a40ed3dSBarry Smith 
5373a40ed3dSBarry Smith   PetscFunctionBegin;
538d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
539d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
540d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
5411ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5421ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
54361de5317SMatthew Knepley #if 0
54471044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
54561de5317SMatthew Knepley   {
54661de5317SMatthew Knepley     PetscScalar *yCheck = new PetscScalar[m];
54761de5317SMatthew Knepley     for(int i = 0; i < m; ++i) {
54861de5317SMatthew Knepley       yCheck[i] = 0.0;
54961de5317SMatthew Knepley       for(int j = 0; j < n; ++j) {
55061de5317SMatthew Knepley         yCheck[i] += v[i*n+j]*x[j];
55161de5317SMatthew Knepley       }
55261de5317SMatthew Knepley       if (fabs(y[i] - yCheck[i]) > 1.0e-8) {
55361de5317SMatthew Knepley         SETERRQ2(PETSC_ERR_LIB, "DGEMV call is fucked up: %g should be %g", y[i], yCheck[i]);
55461de5317SMatthew Knepley       }
55561de5317SMatthew Knepley     }
55661de5317SMatthew Knepley     delete [] y;
55761de5317SMatthew Knepley   }
55861de5317SMatthew Knepley #else
559800995b7SMatthew Knepley   for(i = 0; i < m; ++i) {
56061de5317SMatthew Knepley     y[i] = 0.0;
561800995b7SMatthew Knepley     for(j = 0; j < n; ++j) {
56261de5317SMatthew Knepley       y[i] += v[i*n+j]*x[j];
56361de5317SMatthew Knepley     }
56461de5317SMatthew Knepley   }
56561de5317SMatthew Knepley #endif
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 - A->rmap->n);CHKERRQ(ierr);
5693a40ed3dSBarry Smith   PetscFunctionReturn(0);
570289bc588SBarry Smith }
5716ee01492SSatish Balay 
5724a2ae208SSatish Balay #undef __FUNCT__
5734a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
574dfbe8321SBarry Smith PetscErrorCode MatMultAdd_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,_DOne=1.0;
578dfbe8321SBarry Smith   PetscErrorCode ierr;
5790805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
5803a40ed3dSBarry Smith 
5813a40ed3dSBarry Smith   PetscFunctionBegin;
582d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
583d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
584d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
585600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5871ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
58871044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5891ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5901ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
591dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
5923a40ed3dSBarry Smith   PetscFunctionReturn(0);
593289bc588SBarry Smith }
5946ee01492SSatish Balay 
5954a2ae208SSatish Balay #undef __FUNCT__
5964a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
597dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
598289bc588SBarry Smith {
599c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
60087828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
601dfbe8321SBarry Smith   PetscErrorCode ierr;
6020805154bSBarry Smith   PetscBLASInt   m, n, _One=1;
60387828ca2SBarry Smith   PetscScalar    _DOne=1.0;
6043a40ed3dSBarry Smith 
6053a40ed3dSBarry Smith   PetscFunctionBegin;
606d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
607d0f46423SBarry Smith   n = PetscBLASIntCast(A->cmap->n);
608d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
609600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
6101ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6111ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
61271044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
6131ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6141ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
615dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
6163a40ed3dSBarry Smith   PetscFunctionReturn(0);
617289bc588SBarry Smith }
618289bc588SBarry Smith 
619289bc588SBarry Smith /* -----------------------------------------------------------------*/
6204a2ae208SSatish Balay #undef __FUNCT__
6214a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
62213f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
623289bc588SBarry Smith {
624c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
62587828ca2SBarry Smith   PetscScalar    *v;
6266849ba73SBarry Smith   PetscErrorCode ierr;
62713f74950SBarry Smith   PetscInt       i;
62867e560aaSBarry Smith 
6293a40ed3dSBarry Smith   PetscFunctionBegin;
630d0f46423SBarry Smith   *ncols = A->cmap->n;
631289bc588SBarry Smith   if (cols) {
632d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
633d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
634289bc588SBarry Smith   }
635289bc588SBarry Smith   if (vals) {
636d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
637289bc588SBarry Smith     v    = mat->v + row;
638d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
639289bc588SBarry Smith   }
6403a40ed3dSBarry Smith   PetscFunctionReturn(0);
641289bc588SBarry Smith }
6426ee01492SSatish Balay 
6434a2ae208SSatish Balay #undef __FUNCT__
6444a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
64513f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
646289bc588SBarry Smith {
647dfbe8321SBarry Smith   PetscErrorCode ierr;
648606d414cSSatish Balay   PetscFunctionBegin;
649606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
650606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
6513a40ed3dSBarry Smith   PetscFunctionReturn(0);
652289bc588SBarry Smith }
653289bc588SBarry Smith /* ----------------------------------------------------------------*/
6544a2ae208SSatish Balay #undef __FUNCT__
6554a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
65613f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
657289bc588SBarry Smith {
658c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
65913f74950SBarry Smith   PetscInt     i,j,idx=0;
660d6dfbf8fSBarry Smith 
6613a40ed3dSBarry Smith   PetscFunctionBegin;
66271fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
663289bc588SBarry Smith   if (!mat->roworiented) {
664dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
665289bc588SBarry Smith       for (j=0; j<n; j++) {
666cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6672515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
668d0f46423SBarry Smith         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
66958804f6dSBarry Smith #endif
670289bc588SBarry Smith         for (i=0; i<m; i++) {
671cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6722515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
673d0f46423SBarry Smith           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
67458804f6dSBarry Smith #endif
675cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
676289bc588SBarry Smith         }
677289bc588SBarry Smith       }
6783a40ed3dSBarry Smith     } else {
679289bc588SBarry Smith       for (j=0; j<n; j++) {
680cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6812515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
682d0f46423SBarry Smith         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
68358804f6dSBarry Smith #endif
684289bc588SBarry Smith         for (i=0; i<m; i++) {
685cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6862515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
687d0f46423SBarry Smith           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
68858804f6dSBarry Smith #endif
689cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
690289bc588SBarry Smith         }
691289bc588SBarry Smith       }
692289bc588SBarry Smith     }
6933a40ed3dSBarry Smith   } else {
694dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
695e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
696cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6972515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
698d0f46423SBarry Smith         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
69958804f6dSBarry Smith #endif
700e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
701cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
7022515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
703d0f46423SBarry Smith           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
70458804f6dSBarry Smith #endif
705cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
706e8d4e0b9SBarry Smith         }
707e8d4e0b9SBarry Smith       }
7083a40ed3dSBarry Smith     } else {
709289bc588SBarry Smith       for (i=0; i<m; i++) {
710cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
7112515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
712d0f46423SBarry Smith         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
71358804f6dSBarry Smith #endif
714289bc588SBarry Smith         for (j=0; j<n; j++) {
715cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
7162515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
717d0f46423SBarry Smith           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
71858804f6dSBarry Smith #endif
719cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
720289bc588SBarry Smith         }
721289bc588SBarry Smith       }
722289bc588SBarry Smith     }
723e8d4e0b9SBarry Smith   }
7243a40ed3dSBarry Smith   PetscFunctionReturn(0);
725289bc588SBarry Smith }
726e8d4e0b9SBarry Smith 
7274a2ae208SSatish Balay #undef __FUNCT__
7284a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
72913f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
730ae80bb75SLois Curfman McInnes {
731ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
73213f74950SBarry Smith   PetscInt     i,j;
733ae80bb75SLois Curfman McInnes 
7343a40ed3dSBarry Smith   PetscFunctionBegin;
735ae80bb75SLois Curfman McInnes   /* row-oriented output */
736ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
73797e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
738d0f46423SBarry Smith     if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n);
739ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
7406f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
741d0f46423SBarry Smith       if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n);
74297e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
743ae80bb75SLois Curfman McInnes     }
744ae80bb75SLois Curfman McInnes   }
7453a40ed3dSBarry Smith   PetscFunctionReturn(0);
746ae80bb75SLois Curfman McInnes }
747ae80bb75SLois Curfman McInnes 
748289bc588SBarry Smith /* -----------------------------------------------------------------*/
749289bc588SBarry Smith 
750e090d566SSatish Balay #include "petscsys.h"
75188e32edaSLois Curfman McInnes 
7524a2ae208SSatish Balay #undef __FUNCT__
7534a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
754a313700dSBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, const MatType type,Mat *A)
75588e32edaSLois Curfman McInnes {
75688e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
75788e32edaSLois Curfman McInnes   Mat            B;
7586849ba73SBarry Smith   PetscErrorCode ierr;
75913f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
76013f74950SBarry Smith   int            fd;
76113f74950SBarry Smith   PetscMPIInt    size;
76213f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
76387828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
76419bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
76588e32edaSLois Curfman McInnes 
7663a40ed3dSBarry Smith   PetscFunctionBegin;
767d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
76829bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
769b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
7700752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
771552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
77288e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
77388e32edaSLois Curfman McInnes 
77490ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
775f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
776f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
777be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7785c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
77990ace30eSBarry Smith     B    = *A;
78090ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
7814a41a4d3SSatish Balay     v    = a->v;
7824a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7834a41a4d3SSatish Balay        from row major to column major */
784b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
78590ace30eSBarry Smith     /* read in nonzero values */
7864a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7874a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7884a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7894a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7904a41a4d3SSatish Balay         *v++ =w[i*N+j];
7914a41a4d3SSatish Balay       }
7924a41a4d3SSatish Balay     }
793606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7946d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7956d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79690ace30eSBarry Smith   } else {
79788e32edaSLois Curfman McInnes     /* read row lengths */
79813f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7990752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
80088e32edaSLois Curfman McInnes 
80188e32edaSLois Curfman McInnes     /* create our matrix */
802f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
803f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
804be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
8055c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
80688e32edaSLois Curfman McInnes     B = *A;
80788e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
80888e32edaSLois Curfman McInnes     v = a->v;
80988e32edaSLois Curfman McInnes 
81088e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
81113f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
812b0a32e0cSBarry Smith     cols = scols;
8130752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
81487828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
815b0a32e0cSBarry Smith     vals = svals;
8160752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
81788e32edaSLois Curfman McInnes 
81888e32edaSLois Curfman McInnes     /* insert into matrix */
81988e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
82088e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
82188e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
82288e32edaSLois Curfman McInnes     }
823606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
824606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
825606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
82688e32edaSLois Curfman McInnes 
8276d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8286d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
82990ace30eSBarry Smith   }
8303a40ed3dSBarry Smith   PetscFunctionReturn(0);
83188e32edaSLois Curfman McInnes }
83288e32edaSLois Curfman McInnes 
833e090d566SSatish Balay #include "petscsys.h"
834932b0c3eSLois Curfman McInnes 
8354a2ae208SSatish Balay #undef __FUNCT__
8364a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
8376849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
838289bc588SBarry Smith {
839932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
840dfbe8321SBarry Smith   PetscErrorCode    ierr;
84113f74950SBarry Smith   PetscInt          i,j;
8422dcb1b2aSMatthew Knepley   const char        *name;
84387828ca2SBarry Smith   PetscScalar       *v;
844f3ef73ceSBarry Smith   PetscViewerFormat format;
8455f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
8465f481a85SSatish Balay   PetscTruth        allreal = PETSC_TRUE;
8475f481a85SSatish Balay #endif
848932b0c3eSLois Curfman McInnes 
8493a40ed3dSBarry Smith   PetscFunctionBegin;
850435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
851b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
852456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
8533a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
854fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
855b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
856d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
85744cd7ae7SLois Curfman McInnes       v = a->v + i;
85877431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
859d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
860aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
861329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
862a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
863329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
864a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
8656831982aSBarry Smith         }
86680cd9d93SLois Curfman McInnes #else
8676831982aSBarry Smith         if (*v) {
868a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8696831982aSBarry Smith         }
87080cd9d93SLois Curfman McInnes #endif
8711b807ce4Svictorle         v += a->lda;
87280cd9d93SLois Curfman McInnes       }
873b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
87480cd9d93SLois Curfman McInnes     }
875b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
8763a40ed3dSBarry Smith   } else {
877b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
878aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
87947989497SBarry Smith     /* determine if matrix has all real values */
88047989497SBarry Smith     v = a->v;
881d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
882ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
88347989497SBarry Smith     }
88447989497SBarry Smith #endif
885fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8863a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
887d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
888d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
889fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
890ffac6cdbSBarry Smith     }
891ffac6cdbSBarry Smith 
892d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
893932b0c3eSLois Curfman McInnes       v = a->v + i;
894d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
895aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
89647989497SBarry Smith         if (allreal) {
897f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
89847989497SBarry Smith         } else {
899f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
90047989497SBarry Smith         }
901289bc588SBarry Smith #else
902f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
903289bc588SBarry Smith #endif
9041b807ce4Svictorle 	v += a->lda;
905289bc588SBarry Smith       }
906b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
907289bc588SBarry Smith     }
908fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
909b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
910ffac6cdbSBarry Smith     }
911b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
912da3a660dSBarry Smith   }
913b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
9143a40ed3dSBarry Smith   PetscFunctionReturn(0);
915289bc588SBarry Smith }
916289bc588SBarry Smith 
9174a2ae208SSatish Balay #undef __FUNCT__
9184a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
9196849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
920932b0c3eSLois Curfman McInnes {
921932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9226849ba73SBarry Smith   PetscErrorCode    ierr;
92313f74950SBarry Smith   int               fd;
924d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
92587828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
926f3ef73ceSBarry Smith   PetscViewerFormat format;
927932b0c3eSLois Curfman McInnes 
9283a40ed3dSBarry Smith   PetscFunctionBegin;
929b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
93090ace30eSBarry Smith 
931b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
9327973ad04SBarry Smith   if (format == PETSC_VIEWER_NATIVE) {
93390ace30eSBarry Smith     /* store the matrix as a dense matrix */
93413f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
935552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
93690ace30eSBarry Smith     col_lens[1] = m;
93790ace30eSBarry Smith     col_lens[2] = n;
93890ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
9396f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
940606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
94190ace30eSBarry Smith 
94290ace30eSBarry Smith     /* write out matrix, by rows */
94387828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
94490ace30eSBarry Smith     v    = a->v;
94590ace30eSBarry Smith     for (j=0; j<n; j++) {
946578230a0SSatish Balay       for (i=0; i<m; i++) {
947578230a0SSatish Balay         vals[j + i*n] = *v++;
94890ace30eSBarry Smith       }
94990ace30eSBarry Smith     }
9506f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
951606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
95290ace30eSBarry Smith   } else {
95313f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
954552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
955932b0c3eSLois Curfman McInnes     col_lens[1] = m;
956932b0c3eSLois Curfman McInnes     col_lens[2] = n;
957932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
958932b0c3eSLois Curfman McInnes 
959932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
960932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
9616f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
962932b0c3eSLois Curfman McInnes 
963932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
964932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
965932b0c3eSLois Curfman McInnes     ict = 0;
966932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
967932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
968932b0c3eSLois Curfman McInnes     }
9696f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
970606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
971932b0c3eSLois Curfman McInnes 
972932b0c3eSLois Curfman McInnes     /* store nonzero values */
97387828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
974932b0c3eSLois Curfman McInnes     ict  = 0;
975932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
976932b0c3eSLois Curfman McInnes       v = a->v + i;
977932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9781b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
979932b0c3eSLois Curfman McInnes       }
980932b0c3eSLois Curfman McInnes     }
9816f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
982606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
98390ace30eSBarry Smith   }
9843a40ed3dSBarry Smith   PetscFunctionReturn(0);
985932b0c3eSLois Curfman McInnes }
986932b0c3eSLois Curfman McInnes 
9874a2ae208SSatish Balay #undef __FUNCT__
9884a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
989dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
990f1af5d2fSBarry Smith {
991f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
992f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9936849ba73SBarry Smith   PetscErrorCode    ierr;
994d0f46423SBarry Smith   PetscInt          m = A->rmap->n,n = A->cmap->n,color,i,j;
99587828ca2SBarry Smith   PetscScalar       *v = a->v;
996b0a32e0cSBarry Smith   PetscViewer       viewer;
997b0a32e0cSBarry Smith   PetscDraw         popup;
998329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
999f3ef73ceSBarry Smith   PetscViewerFormat format;
1000f1af5d2fSBarry Smith 
1001f1af5d2fSBarry Smith   PetscFunctionBegin;
1002f1af5d2fSBarry Smith 
1003f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
1004b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1005b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1006f1af5d2fSBarry Smith 
1007f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
1008fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1009f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
1010b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
1011f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1012f1af5d2fSBarry Smith       x_l = j;
1013f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1014f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1015f1af5d2fSBarry Smith         y_l = m - i - 1.0;
1016f1af5d2fSBarry Smith         y_r = y_l + 1.0;
1017f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
1018329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
1019b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1020329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
1021b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1022f1af5d2fSBarry Smith         } else {
1023f1af5d2fSBarry Smith           continue;
1024f1af5d2fSBarry Smith         }
1025f1af5d2fSBarry Smith #else
1026f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
1027b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1028f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
1029b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1030f1af5d2fSBarry Smith         } else {
1031f1af5d2fSBarry Smith           continue;
1032f1af5d2fSBarry Smith         }
1033f1af5d2fSBarry Smith #endif
1034b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1035f1af5d2fSBarry Smith       }
1036f1af5d2fSBarry Smith     }
1037f1af5d2fSBarry Smith   } else {
1038f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1039f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1040f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
1041f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1042f1af5d2fSBarry Smith     }
1043b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
1044b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
1045b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
1046f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
1047f1af5d2fSBarry Smith       x_l = j;
1048f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1049f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
1050f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
1051f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
1052b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
1053b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1054f1af5d2fSBarry Smith       }
1055f1af5d2fSBarry Smith     }
1056f1af5d2fSBarry Smith   }
1057f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1058f1af5d2fSBarry Smith }
1059f1af5d2fSBarry Smith 
10604a2ae208SSatish Balay #undef __FUNCT__
10614a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
1062dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1063f1af5d2fSBarry Smith {
1064b0a32e0cSBarry Smith   PetscDraw      draw;
1065f1af5d2fSBarry Smith   PetscTruth     isnull;
1066329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1067dfbe8321SBarry Smith   PetscErrorCode ierr;
1068f1af5d2fSBarry Smith 
1069f1af5d2fSBarry Smith   PetscFunctionBegin;
1070b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1071b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1072abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1073f1af5d2fSBarry Smith 
1074f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1075d0f46423SBarry Smith   xr  = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1076f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1077b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1078b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1079f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1080f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1081f1af5d2fSBarry Smith }
1082f1af5d2fSBarry Smith 
10834a2ae208SSatish Balay #undef __FUNCT__
10844a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1085dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1086932b0c3eSLois Curfman McInnes {
1087dfbe8321SBarry Smith   PetscErrorCode ierr;
10886805f65bSBarry Smith   PetscTruth     iascii,isbinary,isdraw;
1089932b0c3eSLois Curfman McInnes 
10903a40ed3dSBarry Smith   PetscFunctionBegin;
109132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1092fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1093fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10940f5bd95cSBarry Smith 
1095c45a1595SBarry Smith   if (iascii) {
1096c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10970f5bd95cSBarry Smith   } else if (isbinary) {
10983a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1099f1af5d2fSBarry Smith   } else if (isdraw) {
1100f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
11015cd90555SBarry Smith   } else {
1102958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1103932b0c3eSLois Curfman McInnes   }
11043a40ed3dSBarry Smith   PetscFunctionReturn(0);
1105932b0c3eSLois Curfman McInnes }
1106289bc588SBarry Smith 
11074a2ae208SSatish Balay #undef __FUNCT__
11084a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1109dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1110289bc588SBarry Smith {
1111ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1112dfbe8321SBarry Smith   PetscErrorCode ierr;
111390f02eecSBarry Smith 
11143a40ed3dSBarry Smith   PetscFunctionBegin;
1115aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1116d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1117a5a9c739SBarry Smith #endif
111805b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
11196857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1120606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1121dbd8c25aSHong Zhang 
1122dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1123901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
11244ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11254ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11264ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
11273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1128289bc588SBarry Smith }
1129289bc588SBarry Smith 
11304a2ae208SSatish Balay #undef __FUNCT__
11314a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1132fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1133289bc588SBarry Smith {
1134c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
11356849ba73SBarry Smith   PetscErrorCode ierr;
113613f74950SBarry Smith   PetscInt       k,j,m,n,M;
113787828ca2SBarry Smith   PetscScalar    *v,tmp;
113848b35521SBarry Smith 
11393a40ed3dSBarry Smith   PetscFunctionBegin;
1140d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1141e9695a30SBarry Smith   if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */
1142a5ce6ee0Svictorle     if (m != n) {
1143634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1144d64ed03dSBarry Smith     } else {
1145d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1146289bc588SBarry Smith         for (k=0; k<j; k++) {
11471b807ce4Svictorle           tmp = v[j + k*M];
11481b807ce4Svictorle           v[j + k*M] = v[k + j*M];
11491b807ce4Svictorle           v[k + j*M] = tmp;
1150289bc588SBarry Smith         }
1151289bc588SBarry Smith       }
1152d64ed03dSBarry Smith     }
11533a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1154d3e5ee88SLois Curfman McInnes     Mat          tmat;
1155ec8511deSBarry Smith     Mat_SeqDense *tmatd;
115687828ca2SBarry Smith     PetscScalar  *v2;
1157ea709b57SSatish Balay 
1158fc4dec0aSBarry Smith     if (reuse == MAT_INITIAL_MATRIX) {
11597adad957SLisandro Dalcin       ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1160d0f46423SBarry Smith       ierr  = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
11617adad957SLisandro Dalcin       ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
11625c5985e7SKris Buschelman       ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1163fc4dec0aSBarry Smith     } else {
1164fc4dec0aSBarry Smith       tmat = *matout;
1165fc4dec0aSBarry Smith     }
1166ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
11670de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1168d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
11691b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1170d3e5ee88SLois Curfman McInnes     }
11716d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11726d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1173d3e5ee88SLois Curfman McInnes     *matout = tmat;
117448b35521SBarry Smith   }
11753a40ed3dSBarry Smith   PetscFunctionReturn(0);
1176289bc588SBarry Smith }
1177289bc588SBarry Smith 
11784a2ae208SSatish Balay #undef __FUNCT__
11794a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1180dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1181289bc588SBarry Smith {
1182c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1183c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
118413f74950SBarry Smith   PetscInt     i,j;
118587828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11869ea5d5aeSSatish Balay 
11873a40ed3dSBarry Smith   PetscFunctionBegin;
1188d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1189d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1190d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
11911b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1192d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
11933a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11941b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11951b807ce4Svictorle     }
1196289bc588SBarry Smith   }
119777c4ece6SBarry Smith   *flg = PETSC_TRUE;
11983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1199289bc588SBarry Smith }
1200289bc588SBarry Smith 
12014a2ae208SSatish Balay #undef __FUNCT__
12024a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1203dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1204289bc588SBarry Smith {
1205c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1206dfbe8321SBarry Smith   PetscErrorCode ierr;
120713f74950SBarry Smith   PetscInt       i,n,len;
120887828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
120944cd7ae7SLois Curfman McInnes 
12103a40ed3dSBarry Smith   PetscFunctionBegin;
12112dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
12127a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
12131ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1214d0f46423SBarry Smith   len = PetscMin(A->rmap->n,A->cmap->n);
1215d0f46423SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
121644cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
12171b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1218289bc588SBarry Smith   }
12191ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
12203a40ed3dSBarry Smith   PetscFunctionReturn(0);
1221289bc588SBarry Smith }
1222289bc588SBarry Smith 
12234a2ae208SSatish Balay #undef __FUNCT__
12244a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1225dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1226289bc588SBarry Smith {
1227c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
122887828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1229dfbe8321SBarry Smith   PetscErrorCode ierr;
1230d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n;
123155659b69SBarry Smith 
12323a40ed3dSBarry Smith   PetscFunctionBegin;
123328988994SBarry Smith   if (ll) {
12347a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
12351ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1236d0f46423SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1237da3a660dSBarry Smith     for (i=0; i<m; i++) {
1238da3a660dSBarry Smith       x = l[i];
1239da3a660dSBarry Smith       v = mat->v + i;
1240da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1241da3a660dSBarry Smith     }
12421ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1243efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1244da3a660dSBarry Smith   }
124528988994SBarry Smith   if (rr) {
12467a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
12471ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1248d0f46423SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1249da3a660dSBarry Smith     for (i=0; i<n; i++) {
1250da3a660dSBarry Smith       x = r[i];
1251da3a660dSBarry Smith       v = mat->v + i*m;
1252da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1253da3a660dSBarry Smith     }
12541ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1255efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1256da3a660dSBarry Smith   }
12573a40ed3dSBarry Smith   PetscFunctionReturn(0);
1258289bc588SBarry Smith }
1259289bc588SBarry Smith 
12604a2ae208SSatish Balay #undef __FUNCT__
12614a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1262dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1263289bc588SBarry Smith {
1264c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
126587828ca2SBarry Smith   PetscScalar  *v = mat->v;
1266329f5518SBarry Smith   PetscReal    sum = 0.0;
1267d0f46423SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap->n,i,j;
1268efee365bSSatish Balay   PetscErrorCode ierr;
126955659b69SBarry Smith 
12703a40ed3dSBarry Smith   PetscFunctionBegin;
1271289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1272a5ce6ee0Svictorle     if (lda>m) {
1273d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1274a5ce6ee0Svictorle 	v = mat->v+j*lda;
1275a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1276a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1277a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1278a5ce6ee0Svictorle #else
1279a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1280a5ce6ee0Svictorle #endif
1281a5ce6ee0Svictorle 	}
1282a5ce6ee0Svictorle       }
1283a5ce6ee0Svictorle     } else {
1284d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1285aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1286329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1287289bc588SBarry Smith #else
1288289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1289289bc588SBarry Smith #endif
1290289bc588SBarry Smith       }
1291a5ce6ee0Svictorle     }
1292064f8208SBarry Smith     *nrm = sqrt(sum);
1293dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
12943a40ed3dSBarry Smith   } else if (type == NORM_1) {
1295064f8208SBarry Smith     *nrm = 0.0;
1296d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
12971b807ce4Svictorle       v = mat->v + j*mat->lda;
1298289bc588SBarry Smith       sum = 0.0;
1299d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
130033a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1301289bc588SBarry Smith       }
1302064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1303289bc588SBarry Smith     }
1304d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
13053a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1306064f8208SBarry Smith     *nrm = 0.0;
1307d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1308289bc588SBarry Smith       v = mat->v + j;
1309289bc588SBarry Smith       sum = 0.0;
1310d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
13111b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1312289bc588SBarry Smith       }
1313064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1314289bc588SBarry Smith     }
1315d0f46423SBarry Smith     ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
13163a40ed3dSBarry Smith   } else {
131729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1318289bc588SBarry Smith   }
13193a40ed3dSBarry Smith   PetscFunctionReturn(0);
1320289bc588SBarry Smith }
1321289bc588SBarry Smith 
13224a2ae208SSatish Balay #undef __FUNCT__
13234a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
13244e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg)
1325289bc588SBarry Smith {
1326c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
132763ba0a88SBarry Smith   PetscErrorCode ierr;
132867e560aaSBarry Smith 
13293a40ed3dSBarry Smith   PetscFunctionBegin;
1330b5a2b587SKris Buschelman   switch (op) {
1331b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
13324e0d8c25SBarry Smith     aij->roworiented = flg;
1333b5a2b587SKris Buschelman     break;
1334512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1335b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
13363971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
13374e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1338b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1339b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
134077e54ba9SKris Buschelman   case MAT_SYMMETRIC:
134177e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
13429a4540c5SBarry Smith   case MAT_HERMITIAN:
13439a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1344600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
1345290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
134677e54ba9SKris Buschelman     break;
1347b5a2b587SKris Buschelman   default:
1348600fe468SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13493a40ed3dSBarry Smith   }
13503a40ed3dSBarry Smith   PetscFunctionReturn(0);
1351289bc588SBarry Smith }
1352289bc588SBarry Smith 
13534a2ae208SSatish Balay #undef __FUNCT__
13544a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1355dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
13566f0a148fSBarry Smith {
1357ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13586849ba73SBarry Smith   PetscErrorCode ierr;
1359d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
13603a40ed3dSBarry Smith 
13613a40ed3dSBarry Smith   PetscFunctionBegin;
1362a5ce6ee0Svictorle   if (lda>m) {
1363d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1364a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1365a5ce6ee0Svictorle     }
1366a5ce6ee0Svictorle   } else {
1367d0f46423SBarry Smith     ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1368a5ce6ee0Svictorle   }
13693a40ed3dSBarry Smith   PetscFunctionReturn(0);
13706f0a148fSBarry Smith }
13716f0a148fSBarry Smith 
13724a2ae208SSatish Balay #undef __FUNCT__
13734a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1374f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13756f0a148fSBarry Smith {
1376ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1377d0f46423SBarry Smith   PetscInt       n = A->cmap->n,i,j;
137887828ca2SBarry Smith   PetscScalar    *slot;
137955659b69SBarry Smith 
13803a40ed3dSBarry Smith   PetscFunctionBegin;
13816f0a148fSBarry Smith   for (i=0; i<N; i++) {
13826f0a148fSBarry Smith     slot = l->v + rows[i];
13836f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13846f0a148fSBarry Smith   }
1385f4df32b1SMatthew Knepley   if (diag != 0.0) {
13866f0a148fSBarry Smith     for (i=0; i<N; i++) {
13876f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1388f4df32b1SMatthew Knepley       *slot = diag;
13896f0a148fSBarry Smith     }
13906f0a148fSBarry Smith   }
13913a40ed3dSBarry Smith   PetscFunctionReturn(0);
13926f0a148fSBarry Smith }
1393557bce09SLois Curfman McInnes 
13944a2ae208SSatish Balay #undef __FUNCT__
13954a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1396dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
139764e87e97SBarry Smith {
1398c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13993a40ed3dSBarry Smith 
14003a40ed3dSBarry Smith   PetscFunctionBegin;
1401d0f46423SBarry Smith   if (mat->lda != A->rmap->n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
140264e87e97SBarry Smith   *array = mat->v;
14033a40ed3dSBarry Smith   PetscFunctionReturn(0);
140464e87e97SBarry Smith }
14050754003eSLois Curfman McInnes 
14064a2ae208SSatish Balay #undef __FUNCT__
14074a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1408dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1409ff14e315SSatish Balay {
14103a40ed3dSBarry Smith   PetscFunctionBegin;
141109b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
14123a40ed3dSBarry Smith   PetscFunctionReturn(0);
1413ff14e315SSatish Balay }
14140754003eSLois Curfman McInnes 
14154a2ae208SSatish Balay #undef __FUNCT__
14164a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
141713f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
14180754003eSLois Curfman McInnes {
1419c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
14206849ba73SBarry Smith   PetscErrorCode ierr;
14215d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
14225d0c19d7SBarry Smith   const PetscInt *irow,*icol;
142387828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
14240754003eSLois Curfman McInnes   Mat            newmat;
14250754003eSLois Curfman McInnes 
14263a40ed3dSBarry Smith   PetscFunctionBegin;
142778b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
142878b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1429e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1430e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
14310754003eSLois Curfman McInnes 
1432182d2002SSatish Balay   /* Check submatrixcall */
1433182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
143413f74950SBarry Smith     PetscInt n_cols,n_rows;
1435182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
143621a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
143721a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
1438c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
143921a2c019SBarry Smith     }
1440182d2002SSatish Balay     newmat = *B;
1441182d2002SSatish Balay   } else {
14420754003eSLois Curfman McInnes     /* Create and fill new matrix */
14437adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1444f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
14457adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
14465c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1447182d2002SSatish Balay   }
1448182d2002SSatish Balay 
1449182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1450182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1451182d2002SSatish Balay 
1452182d2002SSatish Balay   for (i=0; i<ncols; i++) {
14536de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1454182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1455182d2002SSatish Balay       *bv++ = av[irow[j]];
14560754003eSLois Curfman McInnes     }
14570754003eSLois Curfman McInnes   }
1458182d2002SSatish Balay 
1459182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
14606d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14616d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14620754003eSLois Curfman McInnes 
14630754003eSLois Curfman McInnes   /* Free work space */
146478b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
146578b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1466182d2002SSatish Balay   *B = newmat;
14673a40ed3dSBarry Smith   PetscFunctionReturn(0);
14680754003eSLois Curfman McInnes }
14690754003eSLois Curfman McInnes 
14704a2ae208SSatish Balay #undef __FUNCT__
14714a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
147213f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1473905e6a2fSBarry Smith {
14746849ba73SBarry Smith   PetscErrorCode ierr;
147513f74950SBarry Smith   PetscInt       i;
1476905e6a2fSBarry Smith 
14773a40ed3dSBarry Smith   PetscFunctionBegin;
1478905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1479b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1480905e6a2fSBarry Smith   }
1481905e6a2fSBarry Smith 
1482905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14836a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1484905e6a2fSBarry Smith   }
14853a40ed3dSBarry Smith   PetscFunctionReturn(0);
1486905e6a2fSBarry Smith }
1487905e6a2fSBarry Smith 
14884a2ae208SSatish Balay #undef __FUNCT__
1489c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1490c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1491c0aa2d19SHong Zhang {
1492c0aa2d19SHong Zhang   PetscFunctionBegin;
1493c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1494c0aa2d19SHong Zhang }
1495c0aa2d19SHong Zhang 
1496c0aa2d19SHong Zhang #undef __FUNCT__
1497c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1498c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1499c0aa2d19SHong Zhang {
1500c0aa2d19SHong Zhang   PetscFunctionBegin;
1501c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1502c0aa2d19SHong Zhang }
1503c0aa2d19SHong Zhang 
1504c0aa2d19SHong Zhang #undef __FUNCT__
15054a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1506dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
15074b0e389bSBarry Smith {
15084b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
15096849ba73SBarry Smith   PetscErrorCode ierr;
1510d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
15113a40ed3dSBarry Smith 
15123a40ed3dSBarry Smith   PetscFunctionBegin;
151333f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
151433f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1515cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
15163a40ed3dSBarry Smith     PetscFunctionReturn(0);
15173a40ed3dSBarry Smith   }
1518d0f46423SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1519a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
15200dbb7854Svictorle     for (j=0; j<n; j++) {
1521a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1522a5ce6ee0Svictorle     }
1523a5ce6ee0Svictorle   } else {
1524d0f46423SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
1525a5ce6ee0Svictorle   }
1526273d9f13SBarry Smith   PetscFunctionReturn(0);
1527273d9f13SBarry Smith }
1528273d9f13SBarry Smith 
15294a2ae208SSatish Balay #undef __FUNCT__
15304a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1531dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1532273d9f13SBarry Smith {
1533dfbe8321SBarry Smith   PetscErrorCode ierr;
1534273d9f13SBarry Smith 
1535273d9f13SBarry Smith   PetscFunctionBegin;
1536273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
15373a40ed3dSBarry Smith   PetscFunctionReturn(0);
15384b0e389bSBarry Smith }
15394b0e389bSBarry Smith 
1540284134d9SBarry Smith #undef __FUNCT__
1541284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1542284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1543284134d9SBarry Smith {
1544284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1545284134d9SBarry Smith   PetscFunctionBegin;
154621a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1547284134d9SBarry Smith   m = PetscMax(m,M);
1548284134d9SBarry Smith   n = PetscMax(n,N);
154921a2c019SBarry Smith   if (m > a->Mmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m);
1550284134d9SBarry Smith   if (n > a->Nmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n);
1551d0f46423SBarry Smith   A->rmap->n = A->rmap->n = m;
1552d0f46423SBarry Smith   A->cmap->n = A->cmap->N = n;
155321a2c019SBarry Smith   if (a->changelda) a->lda = m;
1554284134d9SBarry Smith   PetscFunctionReturn(0);
1555284134d9SBarry Smith }
1556170fe5c8SBarry Smith 
1557284134d9SBarry Smith 
1558a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1559a9fe9ddaSSatish Balay #undef __FUNCT__
1560a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1561a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1562a9fe9ddaSSatish Balay {
1563a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1564a9fe9ddaSSatish Balay 
1565a9fe9ddaSSatish Balay   PetscFunctionBegin;
1566a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1567a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1568a9fe9ddaSSatish Balay   }
1569a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1570a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1571a9fe9ddaSSatish Balay }
1572a9fe9ddaSSatish Balay 
1573a9fe9ddaSSatish Balay #undef __FUNCT__
1574a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1575a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1576a9fe9ddaSSatish Balay {
1577ee16a9a1SHong Zhang   PetscErrorCode ierr;
1578d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
1579ee16a9a1SHong Zhang   Mat            Cmat;
1580a9fe9ddaSSatish Balay 
1581ee16a9a1SHong Zhang   PetscFunctionBegin;
1582d0f46423SBarry Smith   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n);
1583ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1584ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1585ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1586ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1587ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1588ee16a9a1SHong Zhang   *C = Cmat;
1589ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1590ee16a9a1SHong Zhang }
1591a9fe9ddaSSatish Balay 
159298a3b096SSatish Balay #undef __FUNCT__
1593a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1594a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1595a9fe9ddaSSatish Balay {
1596a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1597a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1598a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
15990805154bSBarry Smith   PetscBLASInt   m,n,k;
1600a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1601a9fe9ddaSSatish Balay 
1602a9fe9ddaSSatish Balay   PetscFunctionBegin;
1603d0f46423SBarry Smith   m = PetscBLASIntCast(A->rmap->n);
1604d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1605d0f46423SBarry Smith   k = PetscBLASIntCast(A->cmap->n);
1606a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1607a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1608a9fe9ddaSSatish Balay }
1609a9fe9ddaSSatish Balay 
1610a9fe9ddaSSatish Balay #undef __FUNCT__
1611a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1612a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1613a9fe9ddaSSatish Balay {
1614a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1615a9fe9ddaSSatish Balay 
1616a9fe9ddaSSatish Balay   PetscFunctionBegin;
1617a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1618a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1619a9fe9ddaSSatish Balay   }
1620a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1621a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1622a9fe9ddaSSatish Balay }
1623a9fe9ddaSSatish Balay 
1624a9fe9ddaSSatish Balay #undef __FUNCT__
1625a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1626a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1627a9fe9ddaSSatish Balay {
1628ee16a9a1SHong Zhang   PetscErrorCode ierr;
1629d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
1630ee16a9a1SHong Zhang   Mat            Cmat;
1631a9fe9ddaSSatish Balay 
1632ee16a9a1SHong Zhang   PetscFunctionBegin;
1633d0f46423SBarry Smith   if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n);
1634ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1635ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1636ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1637ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1638ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1639ee16a9a1SHong Zhang   *C = Cmat;
1640ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1641ee16a9a1SHong Zhang }
1642a9fe9ddaSSatish Balay 
1643a9fe9ddaSSatish Balay #undef __FUNCT__
1644a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1645a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1646a9fe9ddaSSatish Balay {
1647a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1648a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1649a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
16500805154bSBarry Smith   PetscBLASInt   m,n,k;
1651a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1652a9fe9ddaSSatish Balay 
1653a9fe9ddaSSatish Balay   PetscFunctionBegin;
1654d0f46423SBarry Smith   m = PetscBLASIntCast(A->cmap->n);
1655d0f46423SBarry Smith   n = PetscBLASIntCast(B->cmap->n);
1656d0f46423SBarry Smith   k = PetscBLASIntCast(A->rmap->n);
16572fbe02b9SBarry Smith   /*
16582fbe02b9SBarry Smith      Note the m and n arguments below are the number rows and columns of A', not A!
16592fbe02b9SBarry Smith   */
1660a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1661a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1662a9fe9ddaSSatish Balay }
1663985db425SBarry Smith 
1664985db425SBarry Smith #undef __FUNCT__
1665985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1666985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1667985db425SBarry Smith {
1668985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1669985db425SBarry Smith   PetscErrorCode ierr;
1670d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1671985db425SBarry Smith   PetscScalar    *x;
1672985db425SBarry Smith   MatScalar      *aa = a->v;
1673985db425SBarry Smith 
1674985db425SBarry Smith   PetscFunctionBegin;
1675985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1676985db425SBarry Smith 
1677985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1678985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1679985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1680d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1681985db425SBarry Smith   for (i=0; i<m; i++) {
1682985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1683985db425SBarry Smith     for (j=1; j<n; j++){
1684985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1685985db425SBarry Smith     }
1686985db425SBarry Smith   }
1687985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1688985db425SBarry Smith   PetscFunctionReturn(0);
1689985db425SBarry Smith }
1690985db425SBarry Smith 
1691985db425SBarry Smith #undef __FUNCT__
1692985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1693985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1694985db425SBarry Smith {
1695985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1696985db425SBarry Smith   PetscErrorCode ierr;
1697d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1698985db425SBarry Smith   PetscScalar    *x;
1699985db425SBarry Smith   PetscReal      atmp;
1700985db425SBarry Smith   MatScalar      *aa = a->v;
1701985db425SBarry Smith 
1702985db425SBarry Smith   PetscFunctionBegin;
1703985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1704985db425SBarry Smith 
1705985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1706985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1707985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1708d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1709985db425SBarry Smith   for (i=0; i<m; i++) {
17109189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
1711985db425SBarry Smith     for (j=1; j<n; j++){
1712985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1713985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1714985db425SBarry Smith     }
1715985db425SBarry Smith   }
1716985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1717985db425SBarry Smith   PetscFunctionReturn(0);
1718985db425SBarry Smith }
1719985db425SBarry Smith 
1720985db425SBarry Smith #undef __FUNCT__
1721985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1722985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1723985db425SBarry Smith {
1724985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1725985db425SBarry Smith   PetscErrorCode ierr;
1726d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
1727985db425SBarry Smith   PetscScalar    *x;
1728985db425SBarry Smith   MatScalar      *aa = a->v;
1729985db425SBarry Smith 
1730985db425SBarry Smith   PetscFunctionBegin;
1731985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1732985db425SBarry Smith 
1733985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1734985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1735985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1736d0f46423SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1737985db425SBarry Smith   for (i=0; i<m; i++) {
1738985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1739985db425SBarry Smith     for (j=1; j<n; j++){
1740985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1741985db425SBarry Smith     }
1742985db425SBarry Smith   }
1743985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1744985db425SBarry Smith   PetscFunctionReturn(0);
1745985db425SBarry Smith }
1746985db425SBarry Smith 
17478d0534beSBarry Smith #undef __FUNCT__
17488d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
17498d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
17508d0534beSBarry Smith {
17518d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
17528d0534beSBarry Smith   PetscErrorCode ierr;
17538d0534beSBarry Smith   PetscScalar    *x;
17548d0534beSBarry Smith 
17558d0534beSBarry Smith   PetscFunctionBegin;
17568d0534beSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
17578d0534beSBarry Smith 
17588d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1759d0f46423SBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
17608d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
17618d0534beSBarry Smith   PetscFunctionReturn(0);
17628d0534beSBarry Smith }
17638d0534beSBarry Smith 
1764289bc588SBarry Smith /* -------------------------------------------------------------------*/
1765a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1766905e6a2fSBarry Smith        MatGetRow_SeqDense,
1767905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1768905e6a2fSBarry Smith        MatMult_SeqDense,
176997304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
17707c922b88SBarry Smith        MatMultTranspose_SeqDense,
17717c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1772db4efbfdSBarry Smith        0,
1773db4efbfdSBarry Smith        0,
1774db4efbfdSBarry Smith        0,
1775db4efbfdSBarry Smith /*10*/ 0,
1776905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1777905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
177841f059aeSBarry Smith        MatSOR_SeqDense,
1779ec8511deSBarry Smith        MatTranspose_SeqDense,
178097304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1781905e6a2fSBarry Smith        MatEqual_SeqDense,
1782905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1783905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1784905e6a2fSBarry Smith        MatNorm_SeqDense,
1785c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1786c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1787905e6a2fSBarry Smith        MatSetOption_SeqDense,
1788905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1789d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense,
1790db4efbfdSBarry Smith        0,
1791db4efbfdSBarry Smith        0,
1792db4efbfdSBarry Smith        0,
1793db4efbfdSBarry Smith        0,
1794d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense,
1795273d9f13SBarry Smith        0,
1796905e6a2fSBarry Smith        0,
1797905e6a2fSBarry Smith        MatGetArray_SeqDense,
1798905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1799d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense,
1800a5ae1ecdSBarry Smith        0,
1801a5ae1ecdSBarry Smith        0,
1802a5ae1ecdSBarry Smith        0,
1803a5ae1ecdSBarry Smith        0,
1804d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense,
1805a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1806a5ae1ecdSBarry Smith        0,
18074b0e389bSBarry Smith        MatGetValues_SeqDense,
1808a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1809d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense,
1810a5ae1ecdSBarry Smith        MatScale_SeqDense,
1811a5ae1ecdSBarry Smith        0,
1812a5ae1ecdSBarry Smith        0,
1813a5ae1ecdSBarry Smith        0,
1814d519adbfSMatthew Knepley /*49*/ 0,
1815a5ae1ecdSBarry Smith        0,
1816a5ae1ecdSBarry Smith        0,
1817a5ae1ecdSBarry Smith        0,
1818a5ae1ecdSBarry Smith        0,
1819d519adbfSMatthew Knepley /*54*/ 0,
1820a5ae1ecdSBarry Smith        0,
1821a5ae1ecdSBarry Smith        0,
1822a5ae1ecdSBarry Smith        0,
1823a5ae1ecdSBarry Smith        0,
1824d519adbfSMatthew Knepley /*59*/ 0,
1825e03a110bSBarry Smith        MatDestroy_SeqDense,
1826e03a110bSBarry Smith        MatView_SeqDense,
1827357abbc8SBarry Smith        0,
182897304618SKris Buschelman        0,
1829d519adbfSMatthew Knepley /*64*/ 0,
183097304618SKris Buschelman        0,
183197304618SKris Buschelman        0,
183297304618SKris Buschelman        0,
183397304618SKris Buschelman        0,
1834d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense,
183597304618SKris Buschelman        0,
183697304618SKris Buschelman        0,
183797304618SKris Buschelman        0,
183897304618SKris Buschelman        0,
1839d519adbfSMatthew Knepley /*74*/ 0,
184097304618SKris Buschelman        0,
184197304618SKris Buschelman        0,
184297304618SKris Buschelman        0,
184397304618SKris Buschelman        0,
1844d519adbfSMatthew Knepley /*79*/ 0,
184597304618SKris Buschelman        0,
184697304618SKris Buschelman        0,
184797304618SKris Buschelman        0,
1848d519adbfSMatthew Knepley /*83*/ MatLoad_SeqDense,
1849865e5f61SKris Buschelman        0,
18501cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1851865e5f61SKris Buschelman        0,
1852865e5f61SKris Buschelman        0,
1853865e5f61SKris Buschelman        0,
1854d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense,
1855a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1856a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1857865e5f61SKris Buschelman        0,
1858865e5f61SKris Buschelman        0,
1859d519adbfSMatthew Knepley /*94*/ 0,
1860a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1861a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1862a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1863284134d9SBarry Smith        0,
1864d519adbfSMatthew Knepley /*99*/ 0,
1865284134d9SBarry Smith        0,
1866284134d9SBarry Smith        0,
1867284134d9SBarry Smith        0,
1868985db425SBarry Smith        MatSetSizes_SeqDense,
1869985db425SBarry Smith        0,
1870985db425SBarry Smith        0,
1871985db425SBarry Smith        0,
1872985db425SBarry Smith        0,
1873985db425SBarry Smith        0,
1874d519adbfSMatthew Knepley /*109*/0,
1875985db425SBarry Smith        0,
18768d0534beSBarry Smith        MatGetRowMin_SeqDense,
18778d0534beSBarry Smith        MatGetColumnVector_SeqDense
1878985db425SBarry Smith };
187990ace30eSBarry Smith 
18804a2ae208SSatish Balay #undef __FUNCT__
18814a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18824b828684SBarry Smith /*@C
1883fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1884d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1885d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1886289bc588SBarry Smith 
1887db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1888db81eaa0SLois Curfman McInnes 
188920563c6bSBarry Smith    Input Parameters:
1890db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
18910c775827SLois Curfman McInnes .  m - number of rows
189218f449edSLois Curfman McInnes .  n - number of columns
1893*c0235b3cSMatthew Knepley -  data - optional location of matrix data in column major order.  Set data=PETSC_NULL for PETSc
1894dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
189520563c6bSBarry Smith 
189620563c6bSBarry Smith    Output Parameter:
189744cd7ae7SLois Curfman McInnes .  A - the matrix
189820563c6bSBarry Smith 
1899b259b22eSLois Curfman McInnes    Notes:
190018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
190118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1902b4fd4287SBarry Smith    set data=PETSC_NULL.
190318f449edSLois Curfman McInnes 
1904027ccd11SLois Curfman McInnes    Level: intermediate
1905027ccd11SLois Curfman McInnes 
1906dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1907d65003e9SLois Curfman McInnes 
1908db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
190920563c6bSBarry Smith @*/
1910be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1911289bc588SBarry Smith {
1912dfbe8321SBarry Smith   PetscErrorCode ierr;
19133b2fbd54SBarry Smith 
19143a40ed3dSBarry Smith   PetscFunctionBegin;
1915f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1916f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1917273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1918273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1919273d9f13SBarry Smith   PetscFunctionReturn(0);
1920273d9f13SBarry Smith }
1921273d9f13SBarry Smith 
19224a2ae208SSatish Balay #undef __FUNCT__
1923afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1924273d9f13SBarry Smith /*@C
1925273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1926273d9f13SBarry Smith 
1927273d9f13SBarry Smith    Collective on MPI_Comm
1928273d9f13SBarry Smith 
1929273d9f13SBarry Smith    Input Parameters:
1930273d9f13SBarry Smith +  A - the matrix
1931273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1932273d9f13SBarry Smith 
1933273d9f13SBarry Smith    Notes:
1934273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1935273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1936284134d9SBarry Smith    need not call this routine.
1937273d9f13SBarry Smith 
1938273d9f13SBarry Smith    Level: intermediate
1939273d9f13SBarry Smith 
1940273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1941273d9f13SBarry Smith 
1942867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA()
1943867c911aSBarry Smith 
1944273d9f13SBarry Smith @*/
1945be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1946273d9f13SBarry Smith {
1947dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1948a23d5eceSKris Buschelman 
1949a23d5eceSKris Buschelman   PetscFunctionBegin;
1950a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1951a23d5eceSKris Buschelman   if (f) {
1952a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1953a23d5eceSKris Buschelman   }
1954a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1955a23d5eceSKris Buschelman }
1956a23d5eceSKris Buschelman 
1957a23d5eceSKris Buschelman EXTERN_C_BEGIN
1958a23d5eceSKris Buschelman #undef __FUNCT__
1959afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1960be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1961a23d5eceSKris Buschelman {
1962273d9f13SBarry Smith   Mat_SeqDense   *b;
1963dfbe8321SBarry Smith   PetscErrorCode ierr;
1964273d9f13SBarry Smith 
1965273d9f13SBarry Smith   PetscFunctionBegin;
1966273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1967273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1968d0f46423SBarry Smith   if (b->lda <= 0) b->lda = B->rmap->n;
19699e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
19709e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
19715afd5e0cSBarry Smith     ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1972284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1973284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
19749e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
1975273d9f13SBarry Smith   } else { /* user-allocated storage */
19769e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1977273d9f13SBarry Smith     b->v          = data;
1978273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1979273d9f13SBarry Smith   }
19800450473dSBarry Smith   B->assembled = PETSC_TRUE;
1981273d9f13SBarry Smith   PetscFunctionReturn(0);
1982273d9f13SBarry Smith }
1983a23d5eceSKris Buschelman EXTERN_C_END
1984273d9f13SBarry Smith 
19851b807ce4Svictorle #undef __FUNCT__
19861b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
19871b807ce4Svictorle /*@C
19881b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
19891b807ce4Svictorle 
19901b807ce4Svictorle   Input parameter:
19911b807ce4Svictorle + A - the matrix
19921b807ce4Svictorle - lda - the leading dimension
19931b807ce4Svictorle 
19941b807ce4Svictorle   Notes:
1995867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
19961b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
19971b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
19981b807ce4Svictorle 
19991b807ce4Svictorle   Level: intermediate
20001b807ce4Svictorle 
20011b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
20021b807ce4Svictorle 
2003284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2004867c911aSBarry Smith 
20051b807ce4Svictorle @*/
2006be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
20071b807ce4Svictorle {
20081b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
200921a2c019SBarry Smith 
20101b807ce4Svictorle   PetscFunctionBegin;
2011d0f46423SBarry Smith   if (lda < B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
20121b807ce4Svictorle   b->lda       = lda;
201321a2c019SBarry Smith   b->changelda = PETSC_FALSE;
201421a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
20151b807ce4Svictorle   PetscFunctionReturn(0);
20161b807ce4Svictorle }
20171b807ce4Svictorle 
20180bad9183SKris Buschelman /*MC
2019fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
20200bad9183SKris Buschelman 
20210bad9183SKris Buschelman    Options Database Keys:
20220bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
20230bad9183SKris Buschelman 
20240bad9183SKris Buschelman   Level: beginner
20250bad9183SKris Buschelman 
202689665df3SBarry Smith .seealso: MatCreateSeqDense()
202789665df3SBarry Smith 
20280bad9183SKris Buschelman M*/
20290bad9183SKris Buschelman 
2030273d9f13SBarry Smith EXTERN_C_BEGIN
20314a2ae208SSatish Balay #undef __FUNCT__
20324a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
2033be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
2034273d9f13SBarry Smith {
2035273d9f13SBarry Smith   Mat_SeqDense   *b;
2036dfbe8321SBarry Smith   PetscErrorCode ierr;
20377c334f02SBarry Smith   PetscMPIInt    size;
2038273d9f13SBarry Smith 
2039273d9f13SBarry Smith   PetscFunctionBegin;
20407adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
204129bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
204255659b69SBarry Smith 
204326283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
204426283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
204526283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
204626283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
2047273d9f13SBarry Smith 
204838f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
2049549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
205090f02eecSBarry Smith   B->mapping      = 0;
205144cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
205218f449edSLois Curfman McInnes 
2053a5ae1ecdSBarry Smith 
205444cd7ae7SLois Curfman McInnes   b->pivots       = 0;
2055273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
2056273d9f13SBarry Smith   b->v            = 0;
2057d0f46423SBarry Smith   b->lda          = B->rmap->n;
205821a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
2059d0f46423SBarry Smith   b->Mmax         = B->rmap->n;
2060d0f46423SBarry Smith   b->Nmax         = B->cmap->n;
20614e220ebcSLois Curfman McInnes 
2062b24902e0SBarry Smith 
2063ec1065edSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
2064b24902e0SBarry Smith                                      "MatGetFactor_seqdense_petsc",
2065b24902e0SBarry Smith                                       MatGetFactor_seqdense_petsc);CHKERRQ(ierr);
2066a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
2067a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
2068a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
20694ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
20704ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
20714ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
20724ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
20734ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
20744ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
20754ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
20764ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
20774ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
207817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
20793a40ed3dSBarry Smith   PetscFunctionReturn(0);
2080289bc588SBarry Smith }
2081c0aa2d19SHong Zhang 
2082c0aa2d19SHong Zhang 
2083273d9f13SBarry Smith EXTERN_C_END
2084