xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 9e8f95c47d204a4cf6c5ba493eaddfe5153ef98d)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
8b4c862e0SSatish Balay #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
17899cda47SBarry Smith   PetscBLASInt   N = (PetscBLASInt)X->rmap.n*X->cmap.n,m=(PetscBLASInt)X->rmap.n,ldax = x->lda,lday=y->lda,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21a5ce6ee0Svictorle   if (ldax>m || lday>m) {
22899cda47SBarry Smith     for (j=0; j<X->cmap.n; j++) {
23f4df32b1SMatthew Knepley       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
24a5ce6ee0Svictorle     }
25a5ce6ee0Svictorle   } else {
26f4df32b1SMatthew Knepley     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
27a5ce6ee0Svictorle   }
28efee365bSSatish Balay   ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr);
293a40ed3dSBarry Smith   PetscFunctionReturn(0);
301987afe7SBarry Smith }
311987afe7SBarry Smith 
324a2ae208SSatish Balay #undef __FUNCT__
334a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
34dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
35289bc588SBarry Smith {
36899cda47SBarry Smith   PetscInt     N = A->rmap.n*A->cmap.n;
373a40ed3dSBarry Smith 
383a40ed3dSBarry Smith   PetscFunctionBegin;
39899cda47SBarry Smith   info->rows_global       = (double)A->rmap.n;
40899cda47SBarry Smith   info->columns_global    = (double)A->cmap.n;
41899cda47SBarry Smith   info->rows_local        = (double)A->rmap.n;
42899cda47SBarry Smith   info->columns_local     = (double)A->cmap.n;
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;
624ce68768SBarry Smith   PetscBLASInt   one = 1,lda = a->lda,j,nz;
63efee365bSSatish Balay   PetscErrorCode ierr;
6480cd9d93SLois Curfman McInnes 
653a40ed3dSBarry Smith   PetscFunctionBegin;
66899cda47SBarry Smith   if (lda>A->rmap.n) {
67899cda47SBarry Smith     nz = (PetscBLASInt)A->rmap.n;
68899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
69f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
72899cda47SBarry Smith     nz = (PetscBLASInt)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;
841cbb95d3SBarry 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;
891cbb95d3SBarry 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 
101289bc588SBarry Smith /* ---------------------------------------------------------------*/
1026831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
1036831982aSBarry Smith    rather than put it in the Mat->row slot.*/
1044a2ae208SSatish Balay #undef __FUNCT__
1054a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
106dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
107289bc588SBarry Smith {
108a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
10981f479a6Svictorle   PetscFunctionBegin;
110a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
111a07ab388SBarry Smith #else
112c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1136849ba73SBarry Smith   PetscErrorCode ierr;
114899cda47SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,m = (PetscBLASInt)A->rmap.n,info;
11555659b69SBarry Smith 
1163a40ed3dSBarry Smith   PetscFunctionBegin;
117289bc588SBarry Smith   if (!mat->pivots) {
118899cda47SBarry Smith     ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
119899cda47SBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap.n*sizeof(PetscBLASInt));CHKERRQ(ierr);
120289bc588SBarry Smith   }
121f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
122899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
12371044d3cSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
12429bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
12529bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
126899cda47SBarry Smith   ierr = PetscLogFlops((2*A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
127a07ab388SBarry Smith #endif
1283a40ed3dSBarry Smith   PetscFunctionReturn(0);
129289bc588SBarry Smith }
1306ee01492SSatish Balay 
1314a2ae208SSatish Balay #undef __FUNCT__
1324a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
133dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
13402cad45dSBarry Smith {
13502cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
1366849ba73SBarry Smith   PetscErrorCode ierr;
13713f74950SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
13802cad45dSBarry Smith   Mat            newi;
13902cad45dSBarry Smith 
1403a40ed3dSBarry Smith   PetscFunctionBegin;
1417adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&newi);CHKERRQ(ierr);
142899cda47SBarry Smith   ierr = MatSetSizes(newi,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
1437adad957SLisandro Dalcin   ierr = MatSetType(newi,((PetscObject)A)->type_name);CHKERRQ(ierr);
1445c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
1455609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
146a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
147899cda47SBarry Smith     if (lda>A->rmap.n) {
148899cda47SBarry Smith       m = A->rmap.n;
149899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
150a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
151a5ce6ee0Svictorle       }
152a5ce6ee0Svictorle     } else {
153899cda47SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
15402cad45dSBarry Smith     }
155a5ce6ee0Svictorle   }
15602cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
15702cad45dSBarry Smith   *newmat = newi;
1583a40ed3dSBarry Smith   PetscFunctionReturn(0);
15902cad45dSBarry Smith }
16002cad45dSBarry Smith 
1614a2ae208SSatish Balay #undef __FUNCT__
1624a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
163dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
164289bc588SBarry Smith {
165dfbe8321SBarry Smith   PetscErrorCode ierr;
1663a40ed3dSBarry Smith 
1673a40ed3dSBarry Smith   PetscFunctionBegin;
1685609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1693a40ed3dSBarry Smith   PetscFunctionReturn(0);
170289bc588SBarry Smith }
1716ee01492SSatish Balay 
1724a2ae208SSatish Balay #undef __FUNCT__
1734a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
174af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
175289bc588SBarry Smith {
17602cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1776849ba73SBarry Smith   PetscErrorCode ierr;
178899cda47SBarry Smith   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->rmap.n,n=A->cmap.n, j;
1794482741eSBarry Smith   MatFactorInfo  info;
1803a40ed3dSBarry Smith 
1813a40ed3dSBarry Smith   PetscFunctionBegin;
18202cad45dSBarry Smith   /* copy the numerical values */
1831b807ce4Svictorle   if (lda1>m || lda2>m ) {
1841b807ce4Svictorle     for (j=0; j<n; j++) {
1851b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1861b807ce4Svictorle     }
1871b807ce4Svictorle   } else {
188899cda47SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1891b807ce4Svictorle   }
19002cad45dSBarry Smith   (*fact)->factor = 0;
1914482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
1923a40ed3dSBarry Smith   PetscFunctionReturn(0);
193289bc588SBarry Smith }
1946ee01492SSatish Balay 
1954a2ae208SSatish Balay #undef __FUNCT__
1964a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
197dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
198289bc588SBarry Smith {
199dfbe8321SBarry Smith   PetscErrorCode ierr;
2003a40ed3dSBarry Smith 
2013a40ed3dSBarry Smith   PetscFunctionBegin;
202ceb03754SKris Buschelman   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr);
2033a40ed3dSBarry Smith   PetscFunctionReturn(0);
204289bc588SBarry Smith }
2056ee01492SSatish Balay 
2064a2ae208SSatish Balay #undef __FUNCT__
2074a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
208dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
209289bc588SBarry Smith {
210a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
211a07ab388SBarry Smith   PetscFunctionBegin;
212a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
213a07ab388SBarry Smith #else
214c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2156849ba73SBarry Smith   PetscErrorCode ierr;
216899cda47SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,info;
21755659b69SBarry Smith 
2183a40ed3dSBarry Smith   PetscFunctionBegin;
219606d414cSSatish Balay   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
220752f5567SLois Curfman McInnes   mat->pivots = 0;
22105b42c5fSBarry Smith 
222899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
22371044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
22477431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
225c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
226899cda47SBarry Smith   ierr = PetscLogFlops((A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
227a07ab388SBarry Smith #endif
2283a40ed3dSBarry Smith   PetscFunctionReturn(0);
229289bc588SBarry Smith }
230289bc588SBarry Smith 
2314a2ae208SSatish Balay #undef __FUNCT__
2320b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
233af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2340b4b3355SBarry Smith {
235dfbe8321SBarry Smith   PetscErrorCode ierr;
23615e8a5b3SHong Zhang   MatFactorInfo  info;
2370b4b3355SBarry Smith 
2380b4b3355SBarry Smith   PetscFunctionBegin;
23915e8a5b3SHong Zhang   info.fill = 1.0;
24015e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2410b4b3355SBarry Smith   PetscFunctionReturn(0);
2420b4b3355SBarry Smith }
2430b4b3355SBarry Smith 
2440b4b3355SBarry Smith #undef __FUNCT__
2454a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
246dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
247289bc588SBarry Smith {
248c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2496849ba73SBarry Smith   PetscErrorCode ierr;
250899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, one = 1,info;
25187828ca2SBarry Smith   PetscScalar    *x,*y;
25267e560aaSBarry Smith 
2533a40ed3dSBarry Smith   PetscFunctionBegin;
2541ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2551ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
256899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
257c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
258ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
25980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
260ae7cfcebSSatish Balay #else
26171044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2624ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
263ae7cfcebSSatish Balay #endif
2647a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
265ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
26680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
267ae7cfcebSSatish Balay #else
26871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2694ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
270ae7cfcebSSatish Balay #endif
271289bc588SBarry Smith   }
27229bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2731ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2741ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
275899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
2763a40ed3dSBarry Smith   PetscFunctionReturn(0);
277289bc588SBarry Smith }
2786ee01492SSatish Balay 
2794a2ae208SSatish Balay #undef __FUNCT__
2804a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
281dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
282da3a660dSBarry Smith {
283c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
284dfbe8321SBarry Smith   PetscErrorCode ierr;
285899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt) A->rmap.n,one = 1,info;
28687828ca2SBarry Smith   PetscScalar    *x,*y;
28767e560aaSBarry Smith 
2883a40ed3dSBarry Smith   PetscFunctionBegin;
2891ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2901ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
291899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
292752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
293da3a660dSBarry Smith   if (mat->pivots) {
294ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
29580444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
296ae7cfcebSSatish Balay #else
29771044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2984ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
299ae7cfcebSSatish Balay #endif
3007a97a34bSBarry Smith   } else {
301ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
30280444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
303ae7cfcebSSatish Balay #else
30471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3054ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
306ae7cfcebSSatish Balay #endif
307da3a660dSBarry Smith   }
3081ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3091ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
310899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
312da3a660dSBarry Smith }
3136ee01492SSatish Balay 
3144a2ae208SSatish Balay #undef __FUNCT__
3154a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
316dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
317da3a660dSBarry Smith {
318c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
319dfbe8321SBarry Smith   PetscErrorCode ierr;
320899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
32187828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
322da3a660dSBarry Smith   Vec            tmp = 0;
32367e560aaSBarry Smith 
3243a40ed3dSBarry Smith   PetscFunctionBegin;
3251ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3261ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
327899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
328da3a660dSBarry Smith   if (yy == zz) {
32978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
33052e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
33178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
332da3a660dSBarry Smith   }
333899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
334752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
335da3a660dSBarry Smith   if (mat->pivots) {
336ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
33780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
338ae7cfcebSSatish Balay #else
33971044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3402ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
341ae7cfcebSSatish Balay #endif
342a8c6a408SBarry Smith   } else {
343ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
34480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
345ae7cfcebSSatish Balay #else
34671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3472ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
348ae7cfcebSSatish Balay #endif
349da3a660dSBarry Smith   }
3502dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
3512dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
3521ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3531ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
354899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
3553a40ed3dSBarry Smith   PetscFunctionReturn(0);
356da3a660dSBarry Smith }
35767e560aaSBarry Smith 
3584a2ae208SSatish Balay #undef __FUNCT__
3594a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
360dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
361da3a660dSBarry Smith {
362c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3636849ba73SBarry Smith   PetscErrorCode ierr;
364899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
36587828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
366da3a660dSBarry Smith   Vec            tmp;
36767e560aaSBarry Smith 
3683a40ed3dSBarry Smith   PetscFunctionBegin;
369899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
3701ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3711ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
372da3a660dSBarry Smith   if (yy == zz) {
37378b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
37452e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
37578b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
376da3a660dSBarry Smith   }
377899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
378752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
379da3a660dSBarry Smith   if (mat->pivots) {
380ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
38180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
382ae7cfcebSSatish Balay #else
38371044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3842ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
385ae7cfcebSSatish Balay #endif
3863a40ed3dSBarry Smith   } else {
387ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
38880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
389ae7cfcebSSatish Balay #else
39071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3912ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
392ae7cfcebSSatish Balay #endif
393da3a660dSBarry Smith   }
39490f02eecSBarry Smith   if (tmp) {
3952dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
39690f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3973a40ed3dSBarry Smith   } else {
3982dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
39990f02eecSBarry Smith   }
4001ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4011ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
402899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
4033a40ed3dSBarry Smith   PetscFunctionReturn(0);
404da3a660dSBarry Smith }
405289bc588SBarry Smith /* ------------------------------------------------------------------*/
4064a2ae208SSatish Balay #undef __FUNCT__
4074a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
40813f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
409289bc588SBarry Smith {
410c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
41187828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
412dfbe8321SBarry Smith   PetscErrorCode ierr;
413899cda47SBarry Smith   PetscInt       m = A->rmap.n,i;
414aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4154ce68768SBarry Smith   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
416bc1b551cSSatish Balay #endif
417289bc588SBarry Smith 
4183a40ed3dSBarry Smith   PetscFunctionBegin;
419289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
42071044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4212dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
422289bc588SBarry Smith   }
4231ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4241ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
425b965ef7fSBarry Smith   its  = its*lits;
42677431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
427289bc588SBarry Smith   while (its--) {
428fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
429289bc588SBarry Smith       for (i=0; i<m; i++) {
430aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
431f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
432f1747703SBarry Smith            not happy about returning a double complex */
43313f74950SBarry Smith         PetscInt         _i;
43487828ca2SBarry Smith         PetscScalar sum = b[i];
435f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4363f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
437f1747703SBarry Smith         }
438f1747703SBarry Smith         xt = sum;
439f1747703SBarry Smith #else
44071044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
441f1747703SBarry Smith #endif
44255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
443289bc588SBarry Smith       }
444289bc588SBarry Smith     }
445fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
446289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
447aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
448f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
449f1747703SBarry Smith            not happy about returning a double complex */
45013f74950SBarry Smith         PetscInt         _i;
45187828ca2SBarry Smith         PetscScalar sum = b[i];
452f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4533f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
454f1747703SBarry Smith         }
455f1747703SBarry Smith         xt = sum;
456f1747703SBarry Smith #else
45771044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
458f1747703SBarry Smith #endif
45955a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
460289bc588SBarry Smith       }
461289bc588SBarry Smith     }
462289bc588SBarry Smith   }
4631ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4653a40ed3dSBarry Smith   PetscFunctionReturn(0);
466289bc588SBarry Smith }
467289bc588SBarry Smith 
468289bc588SBarry Smith /* -----------------------------------------------------------------*/
4694a2ae208SSatish Balay #undef __FUNCT__
4704a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
471dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
472289bc588SBarry Smith {
473c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
47487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
475dfbe8321SBarry Smith   PetscErrorCode ierr;
476899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n,_One=1;
477ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
4783a40ed3dSBarry Smith 
4793a40ed3dSBarry Smith   PetscFunctionBegin;
480899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
4811ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4821ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
48371044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
4841ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4851ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
486899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
4873a40ed3dSBarry Smith   PetscFunctionReturn(0);
488289bc588SBarry Smith }
4896ee01492SSatish Balay 
4904a2ae208SSatish Balay #undef __FUNCT__
4914a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
492dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
493289bc588SBarry Smith {
494c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
49587828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
496dfbe8321SBarry Smith   PetscErrorCode ierr;
497899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
4983a40ed3dSBarry Smith 
4993a40ed3dSBarry Smith   PetscFunctionBegin;
500899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
5011ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5021ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
50371044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
5041ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5051ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
506899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->rmap.n);CHKERRQ(ierr);
5073a40ed3dSBarry Smith   PetscFunctionReturn(0);
508289bc588SBarry Smith }
5096ee01492SSatish Balay 
5104a2ae208SSatish Balay #undef __FUNCT__
5114a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
512dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
513289bc588SBarry Smith {
514c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
51587828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
516dfbe8321SBarry Smith   PetscErrorCode ierr;
517899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
5183a40ed3dSBarry Smith 
5193a40ed3dSBarry Smith   PetscFunctionBegin;
520899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
521600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5221ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5231ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
52471044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5251ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5261ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
527899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5283a40ed3dSBarry Smith   PetscFunctionReturn(0);
529289bc588SBarry Smith }
5306ee01492SSatish Balay 
5314a2ae208SSatish Balay #undef __FUNCT__
5324a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
533dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
534289bc588SBarry Smith {
535c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
53687828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
537dfbe8321SBarry Smith   PetscErrorCode ierr;
538899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
53987828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5403a40ed3dSBarry Smith 
5413a40ed3dSBarry Smith   PetscFunctionBegin;
542899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
543600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5441ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5451ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
54671044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5471ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5481ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
549899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5503a40ed3dSBarry Smith   PetscFunctionReturn(0);
551289bc588SBarry Smith }
552289bc588SBarry Smith 
553289bc588SBarry Smith /* -----------------------------------------------------------------*/
5544a2ae208SSatish Balay #undef __FUNCT__
5554a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
55613f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
557289bc588SBarry Smith {
558c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
55987828ca2SBarry Smith   PetscScalar    *v;
5606849ba73SBarry Smith   PetscErrorCode ierr;
56113f74950SBarry Smith   PetscInt       i;
56267e560aaSBarry Smith 
5633a40ed3dSBarry Smith   PetscFunctionBegin;
564899cda47SBarry Smith   *ncols = A->cmap.n;
565289bc588SBarry Smith   if (cols) {
566899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
567899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) (*cols)[i] = i;
568289bc588SBarry Smith   }
569289bc588SBarry Smith   if (vals) {
570899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
571289bc588SBarry Smith     v    = mat->v + row;
572899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) {(*vals)[i] = *v; v += mat->lda;}
573289bc588SBarry Smith   }
5743a40ed3dSBarry Smith   PetscFunctionReturn(0);
575289bc588SBarry Smith }
5766ee01492SSatish Balay 
5774a2ae208SSatish Balay #undef __FUNCT__
5784a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
57913f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
580289bc588SBarry Smith {
581dfbe8321SBarry Smith   PetscErrorCode ierr;
582606d414cSSatish Balay   PetscFunctionBegin;
583606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
584606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5853a40ed3dSBarry Smith   PetscFunctionReturn(0);
586289bc588SBarry Smith }
587289bc588SBarry Smith /* ----------------------------------------------------------------*/
5884a2ae208SSatish Balay #undef __FUNCT__
5894a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
59013f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
591289bc588SBarry Smith {
592c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
59313f74950SBarry Smith   PetscInt     i,j,idx=0;
594d6dfbf8fSBarry Smith 
5953a40ed3dSBarry Smith   PetscFunctionBegin;
596289bc588SBarry Smith   if (!mat->roworiented) {
597dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
598289bc588SBarry Smith       for (j=0; j<n; j++) {
599cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6002515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
601899cda47SBarry 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);
60258804f6dSBarry Smith #endif
603289bc588SBarry Smith         for (i=0; i<m; i++) {
604cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6052515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
606899cda47SBarry 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);
60758804f6dSBarry Smith #endif
608cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
609289bc588SBarry Smith         }
610289bc588SBarry Smith       }
6113a40ed3dSBarry Smith     } else {
612289bc588SBarry Smith       for (j=0; j<n; j++) {
613cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6142515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
615899cda47SBarry 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);
61658804f6dSBarry Smith #endif
617289bc588SBarry Smith         for (i=0; i<m; i++) {
618cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6192515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
620899cda47SBarry 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);
62158804f6dSBarry Smith #endif
622cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
623289bc588SBarry Smith         }
624289bc588SBarry Smith       }
625289bc588SBarry Smith     }
6263a40ed3dSBarry Smith   } else {
627dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
628e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
629cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6302515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
631899cda47SBarry 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);
63258804f6dSBarry Smith #endif
633e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
634cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6352515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
636899cda47SBarry 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);
63758804f6dSBarry Smith #endif
638cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
639e8d4e0b9SBarry Smith         }
640e8d4e0b9SBarry Smith       }
6413a40ed3dSBarry Smith     } else {
642289bc588SBarry Smith       for (i=0; i<m; i++) {
643cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6442515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
645899cda47SBarry 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);
64658804f6dSBarry Smith #endif
647289bc588SBarry Smith         for (j=0; j<n; j++) {
648cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6492515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
650899cda47SBarry 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);
65158804f6dSBarry Smith #endif
652cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
653289bc588SBarry Smith         }
654289bc588SBarry Smith       }
655289bc588SBarry Smith     }
656e8d4e0b9SBarry Smith   }
6573a40ed3dSBarry Smith   PetscFunctionReturn(0);
658289bc588SBarry Smith }
659e8d4e0b9SBarry Smith 
6604a2ae208SSatish Balay #undef __FUNCT__
6614a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
66213f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
663ae80bb75SLois Curfman McInnes {
664ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
66513f74950SBarry Smith   PetscInt     i,j;
666ae80bb75SLois Curfman McInnes 
6673a40ed3dSBarry Smith   PetscFunctionBegin;
668ae80bb75SLois Curfman McInnes   /* row-oriented output */
669ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
67097e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
671f12f8aa1SBarry 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);
672ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6736f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
6746f31f424SBarry 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);
67597e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
676ae80bb75SLois Curfman McInnes     }
677ae80bb75SLois Curfman McInnes   }
6783a40ed3dSBarry Smith   PetscFunctionReturn(0);
679ae80bb75SLois Curfman McInnes }
680ae80bb75SLois Curfman McInnes 
681289bc588SBarry Smith /* -----------------------------------------------------------------*/
682289bc588SBarry Smith 
683e090d566SSatish Balay #include "petscsys.h"
68488e32edaSLois Curfman McInnes 
6854a2ae208SSatish Balay #undef __FUNCT__
6864a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
687f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A)
68888e32edaSLois Curfman McInnes {
68988e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
69088e32edaSLois Curfman McInnes   Mat            B;
6916849ba73SBarry Smith   PetscErrorCode ierr;
69213f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
69313f74950SBarry Smith   int            fd;
69413f74950SBarry Smith   PetscMPIInt    size;
69513f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
69687828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
69719bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
69888e32edaSLois Curfman McInnes 
6993a40ed3dSBarry Smith   PetscFunctionBegin;
700d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
70129bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
702b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
7030752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
704552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
70588e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
70688e32edaSLois Curfman McInnes 
70790ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
708f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
709f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
710be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7115c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
71290ace30eSBarry Smith     B    = *A;
71390ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
7144a41a4d3SSatish Balay     v    = a->v;
7154a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7164a41a4d3SSatish Balay        from row major to column major */
717b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
71890ace30eSBarry Smith     /* read in nonzero values */
7194a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7204a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7214a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7224a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7234a41a4d3SSatish Balay         *v++ =w[i*N+j];
7244a41a4d3SSatish Balay       }
7254a41a4d3SSatish Balay     }
726606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7276d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7286d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72990ace30eSBarry Smith   } else {
73088e32edaSLois Curfman McInnes     /* read row lengths */
73113f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7320752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
73388e32edaSLois Curfman McInnes 
73488e32edaSLois Curfman McInnes     /* create our matrix */
735f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
736f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
737be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7385c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
73988e32edaSLois Curfman McInnes     B = *A;
74088e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
74188e32edaSLois Curfman McInnes     v = a->v;
74288e32edaSLois Curfman McInnes 
74388e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
74413f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
745b0a32e0cSBarry Smith     cols = scols;
7460752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
74787828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
748b0a32e0cSBarry Smith     vals = svals;
7490752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
75088e32edaSLois Curfman McInnes 
75188e32edaSLois Curfman McInnes     /* insert into matrix */
75288e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
75388e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
75488e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
75588e32edaSLois Curfman McInnes     }
756606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
757606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
758606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
75988e32edaSLois Curfman McInnes 
7606d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7616d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76290ace30eSBarry Smith   }
7633a40ed3dSBarry Smith   PetscFunctionReturn(0);
76488e32edaSLois Curfman McInnes }
76588e32edaSLois Curfman McInnes 
766e090d566SSatish Balay #include "petscsys.h"
767932b0c3eSLois Curfman McInnes 
7684a2ae208SSatish Balay #undef __FUNCT__
7694a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
7706849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
771289bc588SBarry Smith {
772932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
773dfbe8321SBarry Smith   PetscErrorCode    ierr;
77413f74950SBarry Smith   PetscInt          i,j;
7752dcb1b2aSMatthew Knepley   const char        *name;
77687828ca2SBarry Smith   PetscScalar       *v;
777f3ef73ceSBarry Smith   PetscViewerFormat format;
7785f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
7795f481a85SSatish Balay   PetscTruth allreal = PETSC_TRUE;
7805f481a85SSatish Balay #endif
781932b0c3eSLois Curfman McInnes 
7823a40ed3dSBarry Smith   PetscFunctionBegin;
783435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
784b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
785456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7863a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
787fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
788b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
789899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
79044cd7ae7SLois Curfman McInnes       v = a->v + i;
79177431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
792899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
793aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
794329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
795a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
796329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
797a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7986831982aSBarry Smith         }
79980cd9d93SLois Curfman McInnes #else
8006831982aSBarry Smith         if (*v) {
801a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
8026831982aSBarry Smith         }
80380cd9d93SLois Curfman McInnes #endif
8041b807ce4Svictorle         v += a->lda;
80580cd9d93SLois Curfman McInnes       }
806b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
80780cd9d93SLois Curfman McInnes     }
808b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
8093a40ed3dSBarry Smith   } else {
810b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
811aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
81247989497SBarry Smith     /* determine if matrix has all real values */
81347989497SBarry Smith     v = a->v;
814899cda47SBarry Smith     for (i=0; i<A->rmap.n*A->cmap.n; i++) {
815ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
81647989497SBarry Smith     }
81747989497SBarry Smith #endif
818fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8193a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
820899cda47SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr);
821899cda47SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
822fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
823ffac6cdbSBarry Smith     }
824ffac6cdbSBarry Smith 
825899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
826932b0c3eSLois Curfman McInnes       v = a->v + i;
827899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
828aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
82947989497SBarry Smith         if (allreal) {
830f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
83147989497SBarry Smith         } else {
832f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
83347989497SBarry Smith         }
834289bc588SBarry Smith #else
835f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr);
836289bc588SBarry Smith #endif
8371b807ce4Svictorle 	v += a->lda;
838289bc588SBarry Smith       }
839b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
840289bc588SBarry Smith     }
841fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
842b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
843ffac6cdbSBarry Smith     }
844b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
845da3a660dSBarry Smith   }
846b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8473a40ed3dSBarry Smith   PetscFunctionReturn(0);
848289bc588SBarry Smith }
849289bc588SBarry Smith 
8504a2ae208SSatish Balay #undef __FUNCT__
8514a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8526849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
853932b0c3eSLois Curfman McInnes {
854932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8556849ba73SBarry Smith   PetscErrorCode    ierr;
85613f74950SBarry Smith   int               fd;
857899cda47SBarry Smith   PetscInt          ict,j,n = A->cmap.n,m = A->rmap.n,i,*col_lens,nz = m*n;
85887828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
859f3ef73ceSBarry Smith   PetscViewerFormat format;
860932b0c3eSLois Curfman McInnes 
8613a40ed3dSBarry Smith   PetscFunctionBegin;
862b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
86390ace30eSBarry Smith 
864b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
865fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
86690ace30eSBarry Smith     /* store the matrix as a dense matrix */
86713f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
868552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
86990ace30eSBarry Smith     col_lens[1] = m;
87090ace30eSBarry Smith     col_lens[2] = n;
87190ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8726f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
873606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
87490ace30eSBarry Smith 
87590ace30eSBarry Smith     /* write out matrix, by rows */
87687828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
87790ace30eSBarry Smith     v    = a->v;
87890ace30eSBarry Smith     for (j=0; j<n; j++) {
879578230a0SSatish Balay       for (i=0; i<m; i++) {
880578230a0SSatish Balay         vals[j + i*n] = *v++;
88190ace30eSBarry Smith       }
88290ace30eSBarry Smith     }
8836f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
884606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
88590ace30eSBarry Smith   } else {
88613f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
887552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
888932b0c3eSLois Curfman McInnes     col_lens[1] = m;
889932b0c3eSLois Curfman McInnes     col_lens[2] = n;
890932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
891932b0c3eSLois Curfman McInnes 
892932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
893932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8946f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
895932b0c3eSLois Curfman McInnes 
896932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
897932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
898932b0c3eSLois Curfman McInnes     ict = 0;
899932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
900932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
901932b0c3eSLois Curfman McInnes     }
9026f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
903606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
904932b0c3eSLois Curfman McInnes 
905932b0c3eSLois Curfman McInnes     /* store nonzero values */
90687828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
907932b0c3eSLois Curfman McInnes     ict  = 0;
908932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
909932b0c3eSLois Curfman McInnes       v = a->v + i;
910932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
9111b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
912932b0c3eSLois Curfman McInnes       }
913932b0c3eSLois Curfman McInnes     }
9146f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
915606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
91690ace30eSBarry Smith   }
9173a40ed3dSBarry Smith   PetscFunctionReturn(0);
918932b0c3eSLois Curfman McInnes }
919932b0c3eSLois Curfman McInnes 
9204a2ae208SSatish Balay #undef __FUNCT__
9214a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
922dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
923f1af5d2fSBarry Smith {
924f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
925f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9266849ba73SBarry Smith   PetscErrorCode    ierr;
927899cda47SBarry Smith   PetscInt          m = A->rmap.n,n = A->cmap.n,color,i,j;
92887828ca2SBarry Smith   PetscScalar       *v = a->v;
929b0a32e0cSBarry Smith   PetscViewer       viewer;
930b0a32e0cSBarry Smith   PetscDraw         popup;
931329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
932f3ef73ceSBarry Smith   PetscViewerFormat format;
933f1af5d2fSBarry Smith 
934f1af5d2fSBarry Smith   PetscFunctionBegin;
935f1af5d2fSBarry Smith 
936f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
937b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
938b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
939f1af5d2fSBarry Smith 
940f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
941fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
942f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
943b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
944f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
945f1af5d2fSBarry Smith       x_l = j;
946f1af5d2fSBarry Smith       x_r = x_l + 1.0;
947f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
948f1af5d2fSBarry Smith         y_l = m - i - 1.0;
949f1af5d2fSBarry Smith         y_r = y_l + 1.0;
950f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
951329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
952b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
953329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
954b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
955f1af5d2fSBarry Smith         } else {
956f1af5d2fSBarry Smith           continue;
957f1af5d2fSBarry Smith         }
958f1af5d2fSBarry Smith #else
959f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
960b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
961f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
962b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
963f1af5d2fSBarry Smith         } else {
964f1af5d2fSBarry Smith           continue;
965f1af5d2fSBarry Smith         }
966f1af5d2fSBarry Smith #endif
967b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
968f1af5d2fSBarry Smith       }
969f1af5d2fSBarry Smith     }
970f1af5d2fSBarry Smith   } else {
971f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
972f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
973f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
974f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
975f1af5d2fSBarry Smith     }
976b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
977b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
978b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
979f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
980f1af5d2fSBarry Smith       x_l = j;
981f1af5d2fSBarry Smith       x_r = x_l + 1.0;
982f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
983f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
984f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
985b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
986b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
987f1af5d2fSBarry Smith       }
988f1af5d2fSBarry Smith     }
989f1af5d2fSBarry Smith   }
990f1af5d2fSBarry Smith   PetscFunctionReturn(0);
991f1af5d2fSBarry Smith }
992f1af5d2fSBarry Smith 
9934a2ae208SSatish Balay #undef __FUNCT__
9944a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
995dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
996f1af5d2fSBarry Smith {
997b0a32e0cSBarry Smith   PetscDraw      draw;
998f1af5d2fSBarry Smith   PetscTruth     isnull;
999329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1000dfbe8321SBarry Smith   PetscErrorCode ierr;
1001f1af5d2fSBarry Smith 
1002f1af5d2fSBarry Smith   PetscFunctionBegin;
1003b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1004b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1005abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1006f1af5d2fSBarry Smith 
1007f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1008899cda47SBarry Smith   xr  = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
1009f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
1010b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1011b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
1012f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
1013f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1014f1af5d2fSBarry Smith }
1015f1af5d2fSBarry Smith 
10164a2ae208SSatish Balay #undef __FUNCT__
10174a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1018dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1019932b0c3eSLois Curfman McInnes {
1020dfbe8321SBarry Smith   PetscErrorCode ierr;
10216805f65bSBarry Smith   PetscTruth     iascii,isbinary,isdraw;
1022932b0c3eSLois Curfman McInnes 
10233a40ed3dSBarry Smith   PetscFunctionBegin;
102432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1025fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1026fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10270f5bd95cSBarry Smith 
1028c45a1595SBarry Smith   if (iascii) {
1029c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10300f5bd95cSBarry Smith   } else if (isbinary) {
10313a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1032f1af5d2fSBarry Smith   } else if (isdraw) {
1033f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10345cd90555SBarry Smith   } else {
1035958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1036932b0c3eSLois Curfman McInnes   }
10373a40ed3dSBarry Smith   PetscFunctionReturn(0);
1038932b0c3eSLois Curfman McInnes }
1039289bc588SBarry Smith 
10404a2ae208SSatish Balay #undef __FUNCT__
10414a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1042dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1043289bc588SBarry Smith {
1044ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1045dfbe8321SBarry Smith   PetscErrorCode ierr;
104690f02eecSBarry Smith 
10473a40ed3dSBarry Smith   PetscFunctionBegin;
1048aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1049899cda47SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n);
1050a5a9c739SBarry Smith #endif
105105b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
10526857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1053606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1054dbd8c25aSHong Zhang 
1055dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
1056901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10574ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10584ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10594ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10603a40ed3dSBarry Smith   PetscFunctionReturn(0);
1061289bc588SBarry Smith }
1062289bc588SBarry Smith 
10634a2ae208SSatish Balay #undef __FUNCT__
10644a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1065dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1066289bc588SBarry Smith {
1067c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
10686849ba73SBarry Smith   PetscErrorCode ierr;
106913f74950SBarry Smith   PetscInt       k,j,m,n,M;
107087828ca2SBarry Smith   PetscScalar    *v,tmp;
107148b35521SBarry Smith 
10723a40ed3dSBarry Smith   PetscFunctionBegin;
1073899cda47SBarry Smith   v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n;
10747c922b88SBarry Smith   if (!matout) { /* in place transpose */
1075a5ce6ee0Svictorle     if (m != n) {
1076634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1077d64ed03dSBarry Smith     } else {
1078d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1079289bc588SBarry Smith         for (k=0; k<j; k++) {
10801b807ce4Svictorle           tmp = v[j + k*M];
10811b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10821b807ce4Svictorle           v[k + j*M] = tmp;
1083289bc588SBarry Smith         }
1084289bc588SBarry Smith       }
1085d64ed03dSBarry Smith     }
10863a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1087d3e5ee88SLois Curfman McInnes     Mat          tmat;
1088ec8511deSBarry Smith     Mat_SeqDense *tmatd;
108987828ca2SBarry Smith     PetscScalar  *v2;
1090ea709b57SSatish Balay 
10917adad957SLisandro Dalcin     ierr  = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr);
1092899cda47SBarry Smith     ierr  = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr);
10937adad957SLisandro Dalcin     ierr  = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
10945c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1095ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10960de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1097d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10981b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1099d3e5ee88SLois Curfman McInnes     }
11006d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11016d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1102d3e5ee88SLois Curfman McInnes     *matout = tmat;
110348b35521SBarry Smith   }
11043a40ed3dSBarry Smith   PetscFunctionReturn(0);
1105289bc588SBarry Smith }
1106289bc588SBarry Smith 
11074a2ae208SSatish Balay #undef __FUNCT__
11084a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1109dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1110289bc588SBarry Smith {
1111c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1112c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
111313f74950SBarry Smith   PetscInt     i,j;
111487828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11159ea5d5aeSSatish Balay 
11163a40ed3dSBarry Smith   PetscFunctionBegin;
1117899cda47SBarry Smith   if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1118899cda47SBarry Smith   if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1119899cda47SBarry Smith   for (i=0; i<A1->rmap.n; i++) {
11201b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1121899cda47SBarry Smith     for (j=0; j<A1->cmap.n; j++) {
11223a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11231b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11241b807ce4Svictorle     }
1125289bc588SBarry Smith   }
112677c4ece6SBarry Smith   *flg = PETSC_TRUE;
11273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1128289bc588SBarry Smith }
1129289bc588SBarry Smith 
11304a2ae208SSatish Balay #undef __FUNCT__
11314a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1132dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1133289bc588SBarry Smith {
1134c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1135dfbe8321SBarry Smith   PetscErrorCode ierr;
113613f74950SBarry Smith   PetscInt       i,n,len;
113787828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
113844cd7ae7SLois Curfman McInnes 
11393a40ed3dSBarry Smith   PetscFunctionBegin;
11402dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11417a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11421ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1143899cda47SBarry Smith   len = PetscMin(A->rmap.n,A->cmap.n);
1144899cda47SBarry Smith   if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
114544cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11461b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1147289bc588SBarry Smith   }
11481ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11493a40ed3dSBarry Smith   PetscFunctionReturn(0);
1150289bc588SBarry Smith }
1151289bc588SBarry Smith 
11524a2ae208SSatish Balay #undef __FUNCT__
11534a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1154dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1155289bc588SBarry Smith {
1156c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
115787828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1158dfbe8321SBarry Smith   PetscErrorCode ierr;
1159899cda47SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n;
116055659b69SBarry Smith 
11613a40ed3dSBarry Smith   PetscFunctionBegin;
116228988994SBarry Smith   if (ll) {
11637a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11641ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1165899cda47SBarry Smith     if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1166da3a660dSBarry Smith     for (i=0; i<m; i++) {
1167da3a660dSBarry Smith       x = l[i];
1168da3a660dSBarry Smith       v = mat->v + i;
1169da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1170da3a660dSBarry Smith     }
11711ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1172efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1173da3a660dSBarry Smith   }
117428988994SBarry Smith   if (rr) {
11757a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11761ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1177899cda47SBarry Smith     if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1178da3a660dSBarry Smith     for (i=0; i<n; i++) {
1179da3a660dSBarry Smith       x = r[i];
1180da3a660dSBarry Smith       v = mat->v + i*m;
1181da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1182da3a660dSBarry Smith     }
11831ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1184efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1185da3a660dSBarry Smith   }
11863a40ed3dSBarry Smith   PetscFunctionReturn(0);
1187289bc588SBarry Smith }
1188289bc588SBarry Smith 
11894a2ae208SSatish Balay #undef __FUNCT__
11904a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1191dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1192289bc588SBarry Smith {
1193c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
119487828ca2SBarry Smith   PetscScalar  *v = mat->v;
1195329f5518SBarry Smith   PetscReal    sum = 0.0;
1196899cda47SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap.n,i,j;
1197efee365bSSatish Balay   PetscErrorCode ierr;
119855659b69SBarry Smith 
11993a40ed3dSBarry Smith   PetscFunctionBegin;
1200289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1201a5ce6ee0Svictorle     if (lda>m) {
1202899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
1203a5ce6ee0Svictorle 	v = mat->v+j*lda;
1204a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1205a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1206a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1207a5ce6ee0Svictorle #else
1208a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1209a5ce6ee0Svictorle #endif
1210a5ce6ee0Svictorle 	}
1211a5ce6ee0Svictorle       }
1212a5ce6ee0Svictorle     } else {
1213899cda47SBarry Smith       for (i=0; i<A->cmap.n*A->rmap.n; i++) {
1214aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1215329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1216289bc588SBarry Smith #else
1217289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1218289bc588SBarry Smith #endif
1219289bc588SBarry Smith       }
1220a5ce6ee0Svictorle     }
1221064f8208SBarry Smith     *nrm = sqrt(sum);
1222899cda47SBarry Smith     ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12233a40ed3dSBarry Smith   } else if (type == NORM_1) {
1224064f8208SBarry Smith     *nrm = 0.0;
1225899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
12261b807ce4Svictorle       v = mat->v + j*mat->lda;
1227289bc588SBarry Smith       sum = 0.0;
1228899cda47SBarry Smith       for (i=0; i<A->rmap.n; i++) {
122933a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1230289bc588SBarry Smith       }
1231064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1232289bc588SBarry Smith     }
1233899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12343a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1235064f8208SBarry Smith     *nrm = 0.0;
1236899cda47SBarry Smith     for (j=0; j<A->rmap.n; j++) {
1237289bc588SBarry Smith       v = mat->v + j;
1238289bc588SBarry Smith       sum = 0.0;
1239899cda47SBarry Smith       for (i=0; i<A->cmap.n; i++) {
12401b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1241289bc588SBarry Smith       }
1242064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1243289bc588SBarry Smith     }
1244899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12453a40ed3dSBarry Smith   } else {
124629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1247289bc588SBarry Smith   }
12483a40ed3dSBarry Smith   PetscFunctionReturn(0);
1249289bc588SBarry Smith }
1250289bc588SBarry Smith 
12514a2ae208SSatish Balay #undef __FUNCT__
12524a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
12534e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg)
1254289bc588SBarry Smith {
1255c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
125663ba0a88SBarry Smith   PetscErrorCode ierr;
125767e560aaSBarry Smith 
12583a40ed3dSBarry Smith   PetscFunctionBegin;
1259b5a2b587SKris Buschelman   switch (op) {
1260b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
12614e0d8c25SBarry Smith     aij->roworiented = flg;
1262b5a2b587SKris Buschelman     break;
1263512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1264b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
12653971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
12664e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1267b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1268b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
126977e54ba9SKris Buschelman   case MAT_SYMMETRIC:
127077e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12719a4540c5SBarry Smith   case MAT_HERMITIAN:
12729a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1273290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
127477e54ba9SKris Buschelman     break;
1275b5a2b587SKris Buschelman   default:
1276ad86a440SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
12773a40ed3dSBarry Smith   }
12783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1279289bc588SBarry Smith }
1280289bc588SBarry Smith 
12814a2ae208SSatish Balay #undef __FUNCT__
12824a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1283dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
12846f0a148fSBarry Smith {
1285ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12866849ba73SBarry Smith   PetscErrorCode ierr;
1287899cda47SBarry Smith   PetscInt       lda=l->lda,m=A->rmap.n,j;
12883a40ed3dSBarry Smith 
12893a40ed3dSBarry Smith   PetscFunctionBegin;
1290a5ce6ee0Svictorle   if (lda>m) {
1291899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
1292a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1293a5ce6ee0Svictorle     }
1294a5ce6ee0Svictorle   } else {
1295899cda47SBarry Smith     ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1296a5ce6ee0Svictorle   }
12973a40ed3dSBarry Smith   PetscFunctionReturn(0);
12986f0a148fSBarry Smith }
12996f0a148fSBarry Smith 
13004a2ae208SSatish Balay #undef __FUNCT__
13014a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1302f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13036f0a148fSBarry Smith {
1304ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1305899cda47SBarry Smith   PetscInt       n = A->cmap.n,i,j;
130687828ca2SBarry Smith   PetscScalar    *slot;
130755659b69SBarry Smith 
13083a40ed3dSBarry Smith   PetscFunctionBegin;
13096f0a148fSBarry Smith   for (i=0; i<N; i++) {
13106f0a148fSBarry Smith     slot = l->v + rows[i];
13116f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13126f0a148fSBarry Smith   }
1313f4df32b1SMatthew Knepley   if (diag != 0.0) {
13146f0a148fSBarry Smith     for (i=0; i<N; i++) {
13156f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1316f4df32b1SMatthew Knepley       *slot = diag;
13176f0a148fSBarry Smith     }
13186f0a148fSBarry Smith   }
13193a40ed3dSBarry Smith   PetscFunctionReturn(0);
13206f0a148fSBarry Smith }
1321557bce09SLois Curfman McInnes 
13224a2ae208SSatish Balay #undef __FUNCT__
13234a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1324dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
132564e87e97SBarry Smith {
1326c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13273a40ed3dSBarry Smith 
13283a40ed3dSBarry Smith   PetscFunctionBegin;
1329899cda47SBarry Smith   if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
133064e87e97SBarry Smith   *array = mat->v;
13313a40ed3dSBarry Smith   PetscFunctionReturn(0);
133264e87e97SBarry Smith }
13330754003eSLois Curfman McInnes 
13344a2ae208SSatish Balay #undef __FUNCT__
13354a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1336dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1337ff14e315SSatish Balay {
13383a40ed3dSBarry Smith   PetscFunctionBegin;
133909b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13403a40ed3dSBarry Smith   PetscFunctionReturn(0);
1341ff14e315SSatish Balay }
13420754003eSLois Curfman McInnes 
13434a2ae208SSatish Balay #undef __FUNCT__
13444a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
134513f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13460754003eSLois Curfman McInnes {
1347c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13486849ba73SBarry Smith   PetscErrorCode ierr;
134921a2c019SBarry Smith   PetscInt       i,j,*irow,*icol,nrows,ncols;
135087828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13510754003eSLois Curfman McInnes   Mat            newmat;
13520754003eSLois Curfman McInnes 
13533a40ed3dSBarry Smith   PetscFunctionBegin;
135478b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
135578b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1356e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1357e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13580754003eSLois Curfman McInnes 
1359182d2002SSatish Balay   /* Check submatrixcall */
1360182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
136113f74950SBarry Smith     PetscInt n_cols,n_rows;
1362182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
136321a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
136421a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
136521a2c019SBarry Smith       ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr);
136621a2c019SBarry Smith     }
1367182d2002SSatish Balay     newmat = *B;
1368182d2002SSatish Balay   } else {
13690754003eSLois Curfman McInnes     /* Create and fill new matrix */
13707adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
1371f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
13727adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
13735c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1374182d2002SSatish Balay   }
1375182d2002SSatish Balay 
1376182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1377182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1378182d2002SSatish Balay 
1379182d2002SSatish Balay   for (i=0; i<ncols; i++) {
13806de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1381182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1382182d2002SSatish Balay       *bv++ = av[irow[j]];
13830754003eSLois Curfman McInnes     }
13840754003eSLois Curfman McInnes   }
1385182d2002SSatish Balay 
1386182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13876d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13886d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13890754003eSLois Curfman McInnes 
13900754003eSLois Curfman McInnes   /* Free work space */
139178b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
139278b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1393182d2002SSatish Balay   *B = newmat;
13943a40ed3dSBarry Smith   PetscFunctionReturn(0);
13950754003eSLois Curfman McInnes }
13960754003eSLois Curfman McInnes 
13974a2ae208SSatish Balay #undef __FUNCT__
13984a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
139913f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1400905e6a2fSBarry Smith {
14016849ba73SBarry Smith   PetscErrorCode ierr;
140213f74950SBarry Smith   PetscInt       i;
1403905e6a2fSBarry Smith 
14043a40ed3dSBarry Smith   PetscFunctionBegin;
1405905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1406b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1407905e6a2fSBarry Smith   }
1408905e6a2fSBarry Smith 
1409905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14106a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1411905e6a2fSBarry Smith   }
14123a40ed3dSBarry Smith   PetscFunctionReturn(0);
1413905e6a2fSBarry Smith }
1414905e6a2fSBarry Smith 
14154a2ae208SSatish Balay #undef __FUNCT__
1416c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1417c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1418c0aa2d19SHong Zhang {
1419c0aa2d19SHong Zhang   PetscFunctionBegin;
1420c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1421c0aa2d19SHong Zhang }
1422c0aa2d19SHong Zhang 
1423c0aa2d19SHong Zhang #undef __FUNCT__
1424c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1425c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1426c0aa2d19SHong Zhang {
1427c0aa2d19SHong Zhang   PetscFunctionBegin;
1428c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1429c0aa2d19SHong Zhang }
1430c0aa2d19SHong Zhang 
1431c0aa2d19SHong Zhang #undef __FUNCT__
14324a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1433dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14344b0e389bSBarry Smith {
14354b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14366849ba73SBarry Smith   PetscErrorCode ierr;
1437899cda47SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j;
14383a40ed3dSBarry Smith 
14393a40ed3dSBarry Smith   PetscFunctionBegin;
144033f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
144133f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1442cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14433a40ed3dSBarry Smith     PetscFunctionReturn(0);
14443a40ed3dSBarry Smith   }
1445899cda47SBarry Smith   if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1446a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14470dbb7854Svictorle     for (j=0; j<n; j++) {
1448a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1449a5ce6ee0Svictorle     }
1450a5ce6ee0Svictorle   } else {
1451899cda47SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1452a5ce6ee0Svictorle   }
1453273d9f13SBarry Smith   PetscFunctionReturn(0);
1454273d9f13SBarry Smith }
1455273d9f13SBarry Smith 
14564a2ae208SSatish Balay #undef __FUNCT__
14574a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1458dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1459273d9f13SBarry Smith {
1460dfbe8321SBarry Smith   PetscErrorCode ierr;
1461273d9f13SBarry Smith 
1462273d9f13SBarry Smith   PetscFunctionBegin;
1463273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14643a40ed3dSBarry Smith   PetscFunctionReturn(0);
14654b0e389bSBarry Smith }
14664b0e389bSBarry Smith 
1467284134d9SBarry Smith #undef __FUNCT__
1468284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1469284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1470284134d9SBarry Smith {
1471284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
147221a2c019SBarry Smith   PetscErrorCode ierr;
1473284134d9SBarry Smith   PetscFunctionBegin;
147421a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1475284134d9SBarry Smith   m = PetscMax(m,M);
1476284134d9SBarry Smith   n = PetscMax(n,N);
147721a2c019SBarry 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);
1478284134d9SBarry 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);
1479899cda47SBarry Smith   A->rmap.n = A->rmap.n = m;
1480899cda47SBarry Smith   A->cmap.n = A->cmap.N = n;
148121a2c019SBarry Smith   if (a->changelda) a->lda = m;
148221a2c019SBarry Smith   ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
1483284134d9SBarry Smith   PetscFunctionReturn(0);
1484284134d9SBarry Smith }
1485284134d9SBarry Smith 
1486a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1487a9fe9ddaSSatish Balay #undef __FUNCT__
1488a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1489a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1490a9fe9ddaSSatish Balay {
1491a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1492a9fe9ddaSSatish Balay 
1493a9fe9ddaSSatish Balay   PetscFunctionBegin;
1494a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1495a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1496a9fe9ddaSSatish Balay   }
1497a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1498a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1499a9fe9ddaSSatish Balay }
1500a9fe9ddaSSatish Balay 
1501a9fe9ddaSSatish Balay #undef __FUNCT__
1502a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1503a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1504a9fe9ddaSSatish Balay {
1505ee16a9a1SHong Zhang   PetscErrorCode ierr;
1506899cda47SBarry Smith   PetscInt       m=A->rmap.n,n=B->cmap.n;
1507ee16a9a1SHong Zhang   Mat            Cmat;
1508a9fe9ddaSSatish Balay 
1509ee16a9a1SHong Zhang   PetscFunctionBegin;
1510899cda47SBarry 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);
1511ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1512ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1513ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1514ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1515ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1516ee16a9a1SHong Zhang   *C = Cmat;
1517ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1518ee16a9a1SHong Zhang }
1519a9fe9ddaSSatish Balay 
1520a9fe9ddaSSatish Balay #undef __FUNCT__
1521a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1522a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1523a9fe9ddaSSatish Balay {
1524a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1525a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1526a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1527899cda47SBarry Smith   PetscBLASInt   m=(PetscBLASInt)A->rmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->cmap.n;
1528a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1529a9fe9ddaSSatish Balay 
1530a9fe9ddaSSatish Balay   PetscFunctionBegin;
1531a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1532a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1533a9fe9ddaSSatish Balay }
1534a9fe9ddaSSatish Balay 
1535a9fe9ddaSSatish Balay #undef __FUNCT__
1536a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1537a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1538a9fe9ddaSSatish Balay {
1539a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1540a9fe9ddaSSatish Balay 
1541a9fe9ddaSSatish Balay   PetscFunctionBegin;
1542a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1543a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1544a9fe9ddaSSatish Balay   }
1545a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1546a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1547a9fe9ddaSSatish Balay }
1548a9fe9ddaSSatish Balay 
1549a9fe9ddaSSatish Balay #undef __FUNCT__
1550a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1551a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1552a9fe9ddaSSatish Balay {
1553ee16a9a1SHong Zhang   PetscErrorCode ierr;
1554899cda47SBarry Smith   PetscInt       m=A->cmap.n,n=B->cmap.n;
1555ee16a9a1SHong Zhang   Mat            Cmat;
1556a9fe9ddaSSatish Balay 
1557ee16a9a1SHong Zhang   PetscFunctionBegin;
1558899cda47SBarry 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);
1559ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1560ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1561ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1562ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1563ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1564ee16a9a1SHong Zhang   *C = Cmat;
1565ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1566ee16a9a1SHong Zhang }
1567a9fe9ddaSSatish Balay 
1568a9fe9ddaSSatish Balay #undef __FUNCT__
1569a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1570a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1571a9fe9ddaSSatish Balay {
1572a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1573a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1574a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1575899cda47SBarry Smith   PetscBLASInt   m=(PetscBLASInt)A->cmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->rmap.n;
1576a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1577a9fe9ddaSSatish Balay 
1578a9fe9ddaSSatish Balay   PetscFunctionBegin;
1579a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1580a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1581a9fe9ddaSSatish Balay }
1582985db425SBarry Smith 
1583985db425SBarry Smith #undef __FUNCT__
1584985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense"
1585985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
1586985db425SBarry Smith {
1587985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1588985db425SBarry Smith   PetscErrorCode ierr;
1589985db425SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n,p;
1590985db425SBarry Smith   PetscScalar    *x;
1591985db425SBarry Smith   MatScalar      *aa = a->v;
1592985db425SBarry Smith 
1593985db425SBarry Smith   PetscFunctionBegin;
1594985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1595985db425SBarry Smith 
1596985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1597985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1598985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1599985db425SBarry Smith   if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1600985db425SBarry Smith   for (i=0; i<m; i++) {
1601985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1602985db425SBarry Smith     for (j=1; j<n; j++){
1603985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1604985db425SBarry Smith     }
1605985db425SBarry Smith   }
1606985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1607985db425SBarry Smith   PetscFunctionReturn(0);
1608985db425SBarry Smith }
1609985db425SBarry Smith 
1610985db425SBarry Smith #undef __FUNCT__
1611985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense"
1612985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
1613985db425SBarry Smith {
1614985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1615985db425SBarry Smith   PetscErrorCode ierr;
1616985db425SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n,p;
1617985db425SBarry Smith   PetscScalar    *x;
1618985db425SBarry Smith   PetscReal      atmp;
1619985db425SBarry Smith   MatScalar      *aa = a->v;
1620985db425SBarry Smith 
1621985db425SBarry Smith   PetscFunctionBegin;
1622985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1623985db425SBarry Smith 
1624985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1625985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1626985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1627985db425SBarry Smith   if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1628985db425SBarry Smith   for (i=0; i<m; i++) {
1629985db425SBarry Smith     x[i] = PetscAbsScalar(aa[i]); if (idx) idx[i] = 0;
1630985db425SBarry Smith     for (j=1; j<n; j++){
1631985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
1632985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
1633985db425SBarry Smith     }
1634985db425SBarry Smith   }
1635985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1636985db425SBarry Smith   PetscFunctionReturn(0);
1637985db425SBarry Smith }
1638985db425SBarry Smith 
1639985db425SBarry Smith #undef __FUNCT__
1640985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense"
1641985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
1642985db425SBarry Smith {
1643985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1644985db425SBarry Smith   PetscErrorCode ierr;
1645985db425SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n,p;
1646985db425SBarry Smith   PetscScalar    *x;
1647985db425SBarry Smith   MatScalar      *aa = a->v;
1648985db425SBarry Smith 
1649985db425SBarry Smith   PetscFunctionBegin;
1650985db425SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1651985db425SBarry Smith 
1652985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1653985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1654985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
1655985db425SBarry Smith   if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1656985db425SBarry Smith   for (i=0; i<m; i++) {
1657985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
1658985db425SBarry Smith     for (j=1; j<n; j++){
1659985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
1660985db425SBarry Smith     }
1661985db425SBarry Smith   }
1662985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1663985db425SBarry Smith   PetscFunctionReturn(0);
1664985db425SBarry Smith }
1665985db425SBarry Smith 
16668d0534beSBarry Smith #undef __FUNCT__
16678d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense"
16688d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
16698d0534beSBarry Smith {
16708d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
16718d0534beSBarry Smith   PetscErrorCode ierr;
16728d0534beSBarry Smith   PetscScalar    *x;
16738d0534beSBarry Smith 
16748d0534beSBarry Smith   PetscFunctionBegin;
16758d0534beSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
16768d0534beSBarry Smith 
16778d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
16788d0534beSBarry Smith   ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
16798d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
16808d0534beSBarry Smith   PetscFunctionReturn(0);
16818d0534beSBarry Smith }
16828d0534beSBarry Smith 
1683289bc588SBarry Smith /* -------------------------------------------------------------------*/
1684a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1685905e6a2fSBarry Smith        MatGetRow_SeqDense,
1686905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1687905e6a2fSBarry Smith        MatMult_SeqDense,
168897304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
16897c922b88SBarry Smith        MatMultTranspose_SeqDense,
16907c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1691905e6a2fSBarry Smith        MatSolve_SeqDense,
1692905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
16937c922b88SBarry Smith        MatSolveTranspose_SeqDense,
169497304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1695905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1696905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1697ec8511deSBarry Smith        MatRelax_SeqDense,
1698ec8511deSBarry Smith        MatTranspose_SeqDense,
169997304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1700905e6a2fSBarry Smith        MatEqual_SeqDense,
1701905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1702905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1703905e6a2fSBarry Smith        MatNorm_SeqDense,
1704c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1705c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1706905e6a2fSBarry Smith        0,
1707905e6a2fSBarry Smith        MatSetOption_SeqDense,
1708905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
170997304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1710905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1711905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1712905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1713905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
171497304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1715273d9f13SBarry Smith        0,
1716905e6a2fSBarry Smith        0,
1717905e6a2fSBarry Smith        MatGetArray_SeqDense,
1718905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
171997304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1720a5ae1ecdSBarry Smith        0,
1721a5ae1ecdSBarry Smith        0,
1722a5ae1ecdSBarry Smith        0,
1723a5ae1ecdSBarry Smith        0,
172497304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1725a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1726a5ae1ecdSBarry Smith        0,
17274b0e389bSBarry Smith        MatGetValues_SeqDense,
1728a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1729985db425SBarry Smith /*45*/ MatGetRowMax_SeqDense,
1730a5ae1ecdSBarry Smith        MatScale_SeqDense,
1731a5ae1ecdSBarry Smith        0,
1732a5ae1ecdSBarry Smith        0,
1733a5ae1ecdSBarry Smith        0,
1734521d7252SBarry Smith /*50*/ 0,
1735a5ae1ecdSBarry Smith        0,
1736a5ae1ecdSBarry Smith        0,
1737a5ae1ecdSBarry Smith        0,
1738a5ae1ecdSBarry Smith        0,
173997304618SKris Buschelman /*55*/ 0,
1740a5ae1ecdSBarry Smith        0,
1741a5ae1ecdSBarry Smith        0,
1742a5ae1ecdSBarry Smith        0,
1743a5ae1ecdSBarry Smith        0,
174497304618SKris Buschelman /*60*/ 0,
1745e03a110bSBarry Smith        MatDestroy_SeqDense,
1746e03a110bSBarry Smith        MatView_SeqDense,
1747357abbc8SBarry Smith        0,
174897304618SKris Buschelman        0,
174997304618SKris Buschelman /*65*/ 0,
175097304618SKris Buschelman        0,
175197304618SKris Buschelman        0,
175297304618SKris Buschelman        0,
175397304618SKris Buschelman        0,
1754985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqDense,
175597304618SKris Buschelman        0,
175697304618SKris Buschelman        0,
175797304618SKris Buschelman        0,
175897304618SKris Buschelman        0,
175997304618SKris Buschelman /*75*/ 0,
176097304618SKris Buschelman        0,
176197304618SKris Buschelman        0,
176297304618SKris Buschelman        0,
176397304618SKris Buschelman        0,
176497304618SKris Buschelman /*80*/ 0,
176597304618SKris Buschelman        0,
176697304618SKris Buschelman        0,
176797304618SKris Buschelman        0,
1768865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1769865e5f61SKris Buschelman        0,
17701cbb95d3SBarry Smith        MatIsHermitian_SeqDense,
1771865e5f61SKris Buschelman        0,
1772865e5f61SKris Buschelman        0,
1773865e5f61SKris Buschelman        0,
1774a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense,
1775a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1776a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1777865e5f61SKris Buschelman        0,
1778865e5f61SKris Buschelman        0,
1779865e5f61SKris Buschelman /*95*/ 0,
1780a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1781a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1782a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1783284134d9SBarry Smith        0,
1784284134d9SBarry Smith /*100*/0,
1785284134d9SBarry Smith        0,
1786284134d9SBarry Smith        0,
1787284134d9SBarry Smith        0,
1788985db425SBarry Smith        MatSetSizes_SeqDense,
1789985db425SBarry Smith        0,
1790985db425SBarry Smith        0,
1791985db425SBarry Smith        0,
1792985db425SBarry Smith        0,
1793985db425SBarry Smith        0,
1794985db425SBarry Smith /*110*/0,
1795985db425SBarry Smith        0,
17968d0534beSBarry Smith        MatGetRowMin_SeqDense,
17978d0534beSBarry Smith        MatGetColumnVector_SeqDense
1798985db425SBarry Smith };
179990ace30eSBarry Smith 
18004a2ae208SSatish Balay #undef __FUNCT__
18014a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
18024b828684SBarry Smith /*@C
1803fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1804d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1805d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1806289bc588SBarry Smith 
1807db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1808db81eaa0SLois Curfman McInnes 
180920563c6bSBarry Smith    Input Parameters:
1810db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
18110c775827SLois Curfman McInnes .  m - number of rows
181218f449edSLois Curfman McInnes .  n - number of columns
1813db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1814dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
181520563c6bSBarry Smith 
181620563c6bSBarry Smith    Output Parameter:
181744cd7ae7SLois Curfman McInnes .  A - the matrix
181820563c6bSBarry Smith 
1819b259b22eSLois Curfman McInnes    Notes:
182018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
182118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1822b4fd4287SBarry Smith    set data=PETSC_NULL.
182318f449edSLois Curfman McInnes 
1824027ccd11SLois Curfman McInnes    Level: intermediate
1825027ccd11SLois Curfman McInnes 
1826dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1827d65003e9SLois Curfman McInnes 
1828db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
182920563c6bSBarry Smith @*/
1830be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1831289bc588SBarry Smith {
1832dfbe8321SBarry Smith   PetscErrorCode ierr;
18333b2fbd54SBarry Smith 
18343a40ed3dSBarry Smith   PetscFunctionBegin;
1835f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1836f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1837273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1838273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1839273d9f13SBarry Smith   PetscFunctionReturn(0);
1840273d9f13SBarry Smith }
1841273d9f13SBarry Smith 
18424a2ae208SSatish Balay #undef __FUNCT__
1843afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation"
1844273d9f13SBarry Smith /*@C
1845273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1846273d9f13SBarry Smith 
1847273d9f13SBarry Smith    Collective on MPI_Comm
1848273d9f13SBarry Smith 
1849273d9f13SBarry Smith    Input Parameters:
1850273d9f13SBarry Smith +  A - the matrix
1851273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1852273d9f13SBarry Smith 
1853273d9f13SBarry Smith    Notes:
1854273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1855273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1856284134d9SBarry Smith    need not call this routine.
1857273d9f13SBarry Smith 
1858273d9f13SBarry Smith    Level: intermediate
1859273d9f13SBarry Smith 
1860273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1861273d9f13SBarry Smith 
1862273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1863273d9f13SBarry Smith @*/
1864be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1865273d9f13SBarry Smith {
1866dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1867a23d5eceSKris Buschelman 
1868a23d5eceSKris Buschelman   PetscFunctionBegin;
1869a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1870a23d5eceSKris Buschelman   if (f) {
1871a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1872a23d5eceSKris Buschelman   }
1873a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1874a23d5eceSKris Buschelman }
1875a23d5eceSKris Buschelman 
1876a23d5eceSKris Buschelman EXTERN_C_BEGIN
1877a23d5eceSKris Buschelman #undef __FUNCT__
1878afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense"
1879be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1880a23d5eceSKris Buschelman {
1881273d9f13SBarry Smith   Mat_SeqDense   *b;
1882dfbe8321SBarry Smith   PetscErrorCode ierr;
1883273d9f13SBarry Smith 
1884273d9f13SBarry Smith   PetscFunctionBegin;
1885273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1886273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1887*9e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
1888*9e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1889284134d9SBarry Smith     ierr = PetscMalloc((b->lda*b->Nmax+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1890284134d9SBarry Smith     ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1891284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1892*9e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
1893273d9f13SBarry Smith   } else { /* user-allocated storage */
1894*9e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
1895273d9f13SBarry Smith     b->v          = data;
1896273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1897273d9f13SBarry Smith   }
1898273d9f13SBarry Smith   PetscFunctionReturn(0);
1899273d9f13SBarry Smith }
1900a23d5eceSKris Buschelman EXTERN_C_END
1901273d9f13SBarry Smith 
19021b807ce4Svictorle #undef __FUNCT__
19031b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
19041b807ce4Svictorle /*@C
19051b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
19061b807ce4Svictorle 
19071b807ce4Svictorle   Input parameter:
19081b807ce4Svictorle + A - the matrix
19091b807ce4Svictorle - lda - the leading dimension
19101b807ce4Svictorle 
19111b807ce4Svictorle   Notes:
19121b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
19131b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
19141b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
19151b807ce4Svictorle 
19161b807ce4Svictorle   Level: intermediate
19171b807ce4Svictorle 
19181b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
19191b807ce4Svictorle 
1920284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
19211b807ce4Svictorle @*/
1922be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
19231b807ce4Svictorle {
19241b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
192521a2c019SBarry Smith 
19261b807ce4Svictorle   PetscFunctionBegin;
1927899cda47SBarry Smith   if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n);
19281b807ce4Svictorle   b->lda       = lda;
192921a2c019SBarry Smith   b->changelda = PETSC_FALSE;
193021a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
19311b807ce4Svictorle   PetscFunctionReturn(0);
19321b807ce4Svictorle }
19331b807ce4Svictorle 
19340bad9183SKris Buschelman /*MC
1935fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
19360bad9183SKris Buschelman 
19370bad9183SKris Buschelman    Options Database Keys:
19380bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
19390bad9183SKris Buschelman 
19400bad9183SKris Buschelman   Level: beginner
19410bad9183SKris Buschelman 
194289665df3SBarry Smith .seealso: MatCreateSeqDense()
194389665df3SBarry Smith 
19440bad9183SKris Buschelman M*/
19450bad9183SKris Buschelman 
1946273d9f13SBarry Smith EXTERN_C_BEGIN
19474a2ae208SSatish Balay #undef __FUNCT__
19484a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1949be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
1950273d9f13SBarry Smith {
1951273d9f13SBarry Smith   Mat_SeqDense   *b;
1952dfbe8321SBarry Smith   PetscErrorCode ierr;
19537c334f02SBarry Smith   PetscMPIInt    size;
1954273d9f13SBarry Smith 
1955273d9f13SBarry Smith   PetscFunctionBegin;
19567adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
195729bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
195855659b69SBarry Smith 
1959899cda47SBarry Smith   B->rmap.bs = B->cmap.bs = 1;
19606148ca0dSBarry Smith   ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr);
19616148ca0dSBarry Smith   ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr);
1962273d9f13SBarry Smith 
196338f2d2fdSLisandro Dalcin   ierr            = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr);
1964549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
196544cd7ae7SLois Curfman McInnes   B->factor       = 0;
196690f02eecSBarry Smith   B->mapping      = 0;
196744cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
196818f449edSLois Curfman McInnes 
1969a5ae1ecdSBarry Smith 
197044cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1971273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1972273d9f13SBarry Smith   b->v            = 0;
1973899cda47SBarry Smith   b->lda          = B->rmap.n;
197421a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
1975899cda47SBarry Smith   b->Mmax         = B->rmap.n;
1976899cda47SBarry Smith   b->Nmax         = B->cmap.n;
19774e220ebcSLois Curfman McInnes 
1978a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1979a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1980a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
19814ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
19824ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
19834ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
19844ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
19854ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
19864ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
19874ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
19884ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
19894ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
199017667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
19913a40ed3dSBarry Smith   PetscFunctionReturn(0);
1992289bc588SBarry Smith }
1993c0aa2d19SHong Zhang 
1994c0aa2d19SHong Zhang 
1995273d9f13SBarry Smith EXTERN_C_END
1996