xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 899cda47322a0d0eb8e2428039961ef470104e3e)
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;
17*899cda47SBarry 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) {
22*899cda47SBarry 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 {
36*899cda47SBarry Smith   PetscInt     N = A->rmap.n*A->cmap.n;
373a40ed3dSBarry Smith 
383a40ed3dSBarry Smith   PetscFunctionBegin;
39*899cda47SBarry Smith   info->rows_global       = (double)A->rmap.n;
40*899cda47SBarry Smith   info->columns_global    = (double)A->cmap.n;
41*899cda47SBarry Smith   info->rows_local        = (double)A->rmap.n;
42*899cda47SBarry 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;
494e220ebcSLois Curfman McInnes   info->memory            = 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;
66*899cda47SBarry Smith   if (lda>A->rmap.n) {
67*899cda47SBarry Smith     nz = (PetscBLASInt)A->rmap.n;
68*899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
69f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
72*899cda47SBarry 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 
79289bc588SBarry Smith /* ---------------------------------------------------------------*/
806831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
816831982aSBarry Smith    rather than put it in the Mat->row slot.*/
824a2ae208SSatish Balay #undef __FUNCT__
834a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
84dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
85289bc588SBarry Smith {
86a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
8781f479a6Svictorle   PetscFunctionBegin;
88a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
89a07ab388SBarry Smith #else
90c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
916849ba73SBarry Smith   PetscErrorCode ierr;
92*899cda47SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,m = (PetscBLASInt)A->rmap.n,info;
9355659b69SBarry Smith 
943a40ed3dSBarry Smith   PetscFunctionBegin;
95289bc588SBarry Smith   if (!mat->pivots) {
96*899cda47SBarry Smith     ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
97*899cda47SBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap.n*sizeof(PetscBLASInt));CHKERRQ(ierr);
98289bc588SBarry Smith   }
99f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
100*899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
10171044d3cSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
10229bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
10329bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
104*899cda47SBarry Smith   ierr = PetscLogFlops((2*A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
105a07ab388SBarry Smith #endif
1063a40ed3dSBarry Smith   PetscFunctionReturn(0);
107289bc588SBarry Smith }
1086ee01492SSatish Balay 
1094a2ae208SSatish Balay #undef __FUNCT__
1104a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
111dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11202cad45dSBarry Smith {
11302cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
1146849ba73SBarry Smith   PetscErrorCode ierr;
11513f74950SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
11602cad45dSBarry Smith   Mat            newi;
11702cad45dSBarry Smith 
1183a40ed3dSBarry Smith   PetscFunctionBegin;
119f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&newi);CHKERRQ(ierr);
120*899cda47SBarry Smith   ierr = MatSetSizes(newi,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
1215c5985e7SKris Buschelman   ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr);
1225c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
1235609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
124a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
125*899cda47SBarry Smith     if (lda>A->rmap.n) {
126*899cda47SBarry Smith       m = A->rmap.n;
127*899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
128a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
129a5ce6ee0Svictorle       }
130a5ce6ee0Svictorle     } else {
131*899cda47SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
13202cad45dSBarry Smith     }
133a5ce6ee0Svictorle   }
13402cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
13502cad45dSBarry Smith   *newmat = newi;
1363a40ed3dSBarry Smith   PetscFunctionReturn(0);
13702cad45dSBarry Smith }
13802cad45dSBarry Smith 
1394a2ae208SSatish Balay #undef __FUNCT__
1404a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
141dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
142289bc588SBarry Smith {
143dfbe8321SBarry Smith   PetscErrorCode ierr;
1443a40ed3dSBarry Smith 
1453a40ed3dSBarry Smith   PetscFunctionBegin;
1465609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1473a40ed3dSBarry Smith   PetscFunctionReturn(0);
148289bc588SBarry Smith }
1496ee01492SSatish Balay 
1504a2ae208SSatish Balay #undef __FUNCT__
1514a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
152af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
153289bc588SBarry Smith {
15402cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1556849ba73SBarry Smith   PetscErrorCode ierr;
156*899cda47SBarry Smith   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->rmap.n,n=A->cmap.n, j;
1574482741eSBarry Smith   MatFactorInfo  info;
1583a40ed3dSBarry Smith 
1593a40ed3dSBarry Smith   PetscFunctionBegin;
16002cad45dSBarry Smith   /* copy the numerical values */
1611b807ce4Svictorle   if (lda1>m || lda2>m ) {
1621b807ce4Svictorle     for (j=0; j<n; j++) {
1631b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1641b807ce4Svictorle     }
1651b807ce4Svictorle   } else {
166*899cda47SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1671b807ce4Svictorle   }
16802cad45dSBarry Smith   (*fact)->factor = 0;
1694482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
1703a40ed3dSBarry Smith   PetscFunctionReturn(0);
171289bc588SBarry Smith }
1726ee01492SSatish Balay 
1734a2ae208SSatish Balay #undef __FUNCT__
1744a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
175dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
176289bc588SBarry Smith {
177dfbe8321SBarry Smith   PetscErrorCode ierr;
1783a40ed3dSBarry Smith 
1793a40ed3dSBarry Smith   PetscFunctionBegin;
180ceb03754SKris Buschelman   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr);
1813a40ed3dSBarry Smith   PetscFunctionReturn(0);
182289bc588SBarry Smith }
1836ee01492SSatish Balay 
1844a2ae208SSatish Balay #undef __FUNCT__
1854a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
186dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
187289bc588SBarry Smith {
188a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
189a07ab388SBarry Smith   PetscFunctionBegin;
190a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
191a07ab388SBarry Smith #else
192c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1936849ba73SBarry Smith   PetscErrorCode ierr;
194*899cda47SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,info;
19555659b69SBarry Smith 
1963a40ed3dSBarry Smith   PetscFunctionBegin;
197752f5567SLois Curfman McInnes   if (mat->pivots) {
198606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
199*899cda47SBarry Smith     ierr = PetscLogObjectMemory(A,-A->rmap.n*sizeof(PetscInt));CHKERRQ(ierr);
200752f5567SLois Curfman McInnes     mat->pivots = 0;
201752f5567SLois Curfman McInnes   }
202*899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
20371044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
20477431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
205c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
206*899cda47SBarry Smith   ierr = PetscLogFlops((A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
207a07ab388SBarry Smith #endif
2083a40ed3dSBarry Smith   PetscFunctionReturn(0);
209289bc588SBarry Smith }
210289bc588SBarry Smith 
2114a2ae208SSatish Balay #undef __FUNCT__
2120b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
213af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2140b4b3355SBarry Smith {
215dfbe8321SBarry Smith   PetscErrorCode ierr;
21615e8a5b3SHong Zhang   MatFactorInfo  info;
2170b4b3355SBarry Smith 
2180b4b3355SBarry Smith   PetscFunctionBegin;
21915e8a5b3SHong Zhang   info.fill = 1.0;
22015e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2210b4b3355SBarry Smith   PetscFunctionReturn(0);
2220b4b3355SBarry Smith }
2230b4b3355SBarry Smith 
2240b4b3355SBarry Smith #undef __FUNCT__
2254a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
226dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
227289bc588SBarry Smith {
228c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2296849ba73SBarry Smith   PetscErrorCode ierr;
230*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, one = 1,info;
23187828ca2SBarry Smith   PetscScalar    *x,*y;
23267e560aaSBarry Smith 
2333a40ed3dSBarry Smith   PetscFunctionBegin;
2341ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2351ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
236*899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
237c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
238ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
23980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
240ae7cfcebSSatish Balay #else
24171044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2424ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
243ae7cfcebSSatish Balay #endif
2447a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
245ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
24680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
247ae7cfcebSSatish Balay #else
24871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2494ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
250ae7cfcebSSatish Balay #endif
251289bc588SBarry Smith   }
25229bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2531ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2541ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
255*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
257289bc588SBarry Smith }
2586ee01492SSatish Balay 
2594a2ae208SSatish Balay #undef __FUNCT__
2604a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
261dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
262da3a660dSBarry Smith {
263c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
264dfbe8321SBarry Smith   PetscErrorCode ierr;
265*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt) A->rmap.n,one = 1,info;
26687828ca2SBarry Smith   PetscScalar    *x,*y;
26767e560aaSBarry Smith 
2683a40ed3dSBarry Smith   PetscFunctionBegin;
2691ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2701ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
271*899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
272752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
273da3a660dSBarry Smith   if (mat->pivots) {
274ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
27580444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
276ae7cfcebSSatish Balay #else
27771044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2784ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
279ae7cfcebSSatish Balay #endif
2807a97a34bSBarry Smith   } else {
281ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
28280444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
283ae7cfcebSSatish Balay #else
28471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2854ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
286ae7cfcebSSatish Balay #endif
287da3a660dSBarry Smith   }
2881ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2891ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
290*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
2913a40ed3dSBarry Smith   PetscFunctionReturn(0);
292da3a660dSBarry Smith }
2936ee01492SSatish Balay 
2944a2ae208SSatish Balay #undef __FUNCT__
2954a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
296dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
297da3a660dSBarry Smith {
298c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
299dfbe8321SBarry Smith   PetscErrorCode ierr;
300*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
30187828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
302da3a660dSBarry Smith   Vec            tmp = 0;
30367e560aaSBarry Smith 
3043a40ed3dSBarry Smith   PetscFunctionBegin;
3051ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3061ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
307*899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
308da3a660dSBarry Smith   if (yy == zz) {
30978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
31052e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
31178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
312da3a660dSBarry Smith   }
313*899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
314752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
315da3a660dSBarry Smith   if (mat->pivots) {
316ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
31780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
318ae7cfcebSSatish Balay #else
31971044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3202ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
321ae7cfcebSSatish Balay #endif
322a8c6a408SBarry Smith   } else {
323ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
32480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
325ae7cfcebSSatish Balay #else
32671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3272ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
328ae7cfcebSSatish Balay #endif
329da3a660dSBarry Smith   }
3302dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
3312dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
3321ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3331ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
334*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
3353a40ed3dSBarry Smith   PetscFunctionReturn(0);
336da3a660dSBarry Smith }
33767e560aaSBarry Smith 
3384a2ae208SSatish Balay #undef __FUNCT__
3394a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
340dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
341da3a660dSBarry Smith {
342c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3436849ba73SBarry Smith   PetscErrorCode ierr;
344*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
34587828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
346da3a660dSBarry Smith   Vec            tmp;
34767e560aaSBarry Smith 
3483a40ed3dSBarry Smith   PetscFunctionBegin;
349*899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
3501ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3511ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
352da3a660dSBarry Smith   if (yy == zz) {
35378b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
35452e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
35578b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
356da3a660dSBarry Smith   }
357*899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
358752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
359da3a660dSBarry Smith   if (mat->pivots) {
360ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
36180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
362ae7cfcebSSatish Balay #else
36371044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3642ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
365ae7cfcebSSatish Balay #endif
3663a40ed3dSBarry Smith   } else {
367ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
36880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
369ae7cfcebSSatish Balay #else
37071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3712ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
372ae7cfcebSSatish Balay #endif
373da3a660dSBarry Smith   }
37490f02eecSBarry Smith   if (tmp) {
3752dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
37690f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3773a40ed3dSBarry Smith   } else {
3782dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
37990f02eecSBarry Smith   }
3801ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3811ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
382*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
3833a40ed3dSBarry Smith   PetscFunctionReturn(0);
384da3a660dSBarry Smith }
385289bc588SBarry Smith /* ------------------------------------------------------------------*/
3864a2ae208SSatish Balay #undef __FUNCT__
3874a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
38813f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
389289bc588SBarry Smith {
390c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
39187828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
392dfbe8321SBarry Smith   PetscErrorCode ierr;
393*899cda47SBarry Smith   PetscInt       m = A->rmap.n,i;
394aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3954ce68768SBarry Smith   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
396bc1b551cSSatish Balay #endif
397289bc588SBarry Smith 
3983a40ed3dSBarry Smith   PetscFunctionBegin;
399289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
40071044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4012dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
402289bc588SBarry Smith   }
4031ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4041ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
405b965ef7fSBarry Smith   its  = its*lits;
40677431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
407289bc588SBarry Smith   while (its--) {
408fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
409289bc588SBarry Smith       for (i=0; i<m; i++) {
410aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
411f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
412f1747703SBarry Smith            not happy about returning a double complex */
41313f74950SBarry Smith         PetscInt         _i;
41487828ca2SBarry Smith         PetscScalar sum = b[i];
415f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4163f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
417f1747703SBarry Smith         }
418f1747703SBarry Smith         xt = sum;
419f1747703SBarry Smith #else
42071044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
421f1747703SBarry Smith #endif
42255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
423289bc588SBarry Smith       }
424289bc588SBarry Smith     }
425fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
426289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
427aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
428f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
429f1747703SBarry Smith            not happy about returning a double complex */
43013f74950SBarry Smith         PetscInt         _i;
43187828ca2SBarry Smith         PetscScalar sum = b[i];
432f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4333f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
434f1747703SBarry Smith         }
435f1747703SBarry Smith         xt = sum;
436f1747703SBarry Smith #else
43771044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
438f1747703SBarry Smith #endif
43955a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
440289bc588SBarry Smith       }
441289bc588SBarry Smith     }
442289bc588SBarry Smith   }
4431ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4441ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4453a40ed3dSBarry Smith   PetscFunctionReturn(0);
446289bc588SBarry Smith }
447289bc588SBarry Smith 
448289bc588SBarry Smith /* -----------------------------------------------------------------*/
4494a2ae208SSatish Balay #undef __FUNCT__
4504a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
451dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
452289bc588SBarry Smith {
453c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
45487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
455dfbe8321SBarry Smith   PetscErrorCode ierr;
456*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n,_One=1;
457ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
4583a40ed3dSBarry Smith 
4593a40ed3dSBarry Smith   PetscFunctionBegin;
460*899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
4611ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4621ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
46371044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
4641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4651ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
466*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
4673a40ed3dSBarry Smith   PetscFunctionReturn(0);
468289bc588SBarry Smith }
4696ee01492SSatish Balay 
4704a2ae208SSatish Balay #undef __FUNCT__
4714a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
472dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
473289bc588SBarry Smith {
474c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
47587828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
476dfbe8321SBarry Smith   PetscErrorCode ierr;
477*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
4783a40ed3dSBarry Smith 
4793a40ed3dSBarry Smith   PetscFunctionBegin;
480*899cda47SBarry 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_("N",&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);
486*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->rmap.n);CHKERRQ(ierr);
4873a40ed3dSBarry Smith   PetscFunctionReturn(0);
488289bc588SBarry Smith }
4896ee01492SSatish Balay 
4904a2ae208SSatish Balay #undef __FUNCT__
4914a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
492dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,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;
496dfbe8321SBarry Smith   PetscErrorCode ierr;
497*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
4983a40ed3dSBarry Smith 
4993a40ed3dSBarry Smith   PetscFunctionBegin;
500*899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
501600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5021ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5031ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
50471044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5051ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5061ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
507*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5083a40ed3dSBarry Smith   PetscFunctionReturn(0);
509289bc588SBarry Smith }
5106ee01492SSatish Balay 
5114a2ae208SSatish Balay #undef __FUNCT__
5124a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
513dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
514289bc588SBarry Smith {
515c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
51687828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
517dfbe8321SBarry Smith   PetscErrorCode ierr;
518*899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
51987828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5203a40ed3dSBarry Smith 
5213a40ed3dSBarry Smith   PetscFunctionBegin;
522*899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
523600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5241ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5251ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
52671044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5271ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5281ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
529*899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5303a40ed3dSBarry Smith   PetscFunctionReturn(0);
531289bc588SBarry Smith }
532289bc588SBarry Smith 
533289bc588SBarry Smith /* -----------------------------------------------------------------*/
5344a2ae208SSatish Balay #undef __FUNCT__
5354a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
53613f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
537289bc588SBarry Smith {
538c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
53987828ca2SBarry Smith   PetscScalar    *v;
5406849ba73SBarry Smith   PetscErrorCode ierr;
54113f74950SBarry Smith   PetscInt       i;
54267e560aaSBarry Smith 
5433a40ed3dSBarry Smith   PetscFunctionBegin;
544*899cda47SBarry Smith   *ncols = A->cmap.n;
545289bc588SBarry Smith   if (cols) {
546*899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
547*899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) (*cols)[i] = i;
548289bc588SBarry Smith   }
549289bc588SBarry Smith   if (vals) {
550*899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
551289bc588SBarry Smith     v    = mat->v + row;
552*899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) {(*vals)[i] = *v; v += mat->lda;}
553289bc588SBarry Smith   }
5543a40ed3dSBarry Smith   PetscFunctionReturn(0);
555289bc588SBarry Smith }
5566ee01492SSatish Balay 
5574a2ae208SSatish Balay #undef __FUNCT__
5584a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
55913f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
560289bc588SBarry Smith {
561dfbe8321SBarry Smith   PetscErrorCode ierr;
562606d414cSSatish Balay   PetscFunctionBegin;
563606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
564606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
566289bc588SBarry Smith }
567289bc588SBarry Smith /* ----------------------------------------------------------------*/
5684a2ae208SSatish Balay #undef __FUNCT__
5694a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
57013f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
571289bc588SBarry Smith {
572c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
57313f74950SBarry Smith   PetscInt     i,j,idx=0;
574d6dfbf8fSBarry Smith 
5753a40ed3dSBarry Smith   PetscFunctionBegin;
576289bc588SBarry Smith   if (!mat->roworiented) {
577dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
578289bc588SBarry Smith       for (j=0; j<n; j++) {
579cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5802515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
581*899cda47SBarry 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);
58258804f6dSBarry Smith #endif
583289bc588SBarry Smith         for (i=0; i<m; i++) {
584cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5852515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
586*899cda47SBarry 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);
58758804f6dSBarry Smith #endif
588cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
589289bc588SBarry Smith         }
590289bc588SBarry Smith       }
5913a40ed3dSBarry Smith     } else {
592289bc588SBarry Smith       for (j=0; j<n; j++) {
593cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5942515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
595*899cda47SBarry 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);
59658804f6dSBarry Smith #endif
597289bc588SBarry Smith         for (i=0; i<m; i++) {
598cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5992515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
600*899cda47SBarry 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);
60158804f6dSBarry Smith #endif
602cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
603289bc588SBarry Smith         }
604289bc588SBarry Smith       }
605289bc588SBarry Smith     }
6063a40ed3dSBarry Smith   } else {
607dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
608e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
609cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6102515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
611*899cda47SBarry 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);
61258804f6dSBarry Smith #endif
613e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
614cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6152515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
616*899cda47SBarry 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);
61758804f6dSBarry Smith #endif
618cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
619e8d4e0b9SBarry Smith         }
620e8d4e0b9SBarry Smith       }
6213a40ed3dSBarry Smith     } else {
622289bc588SBarry Smith       for (i=0; i<m; i++) {
623cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6242515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
625*899cda47SBarry 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);
62658804f6dSBarry Smith #endif
627289bc588SBarry Smith         for (j=0; j<n; j++) {
628cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6292515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
630*899cda47SBarry 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);
63158804f6dSBarry Smith #endif
632cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
633289bc588SBarry Smith         }
634289bc588SBarry Smith       }
635289bc588SBarry Smith     }
636e8d4e0b9SBarry Smith   }
6373a40ed3dSBarry Smith   PetscFunctionReturn(0);
638289bc588SBarry Smith }
639e8d4e0b9SBarry Smith 
6404a2ae208SSatish Balay #undef __FUNCT__
6414a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
64213f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
643ae80bb75SLois Curfman McInnes {
644ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
64513f74950SBarry Smith   PetscInt     i,j;
64687828ca2SBarry Smith   PetscScalar  *vpt = v;
647ae80bb75SLois Curfman McInnes 
6483a40ed3dSBarry Smith   PetscFunctionBegin;
649ae80bb75SLois Curfman McInnes   /* row-oriented output */
650ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
651ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6521b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
653ae80bb75SLois Curfman McInnes     }
654ae80bb75SLois Curfman McInnes   }
6553a40ed3dSBarry Smith   PetscFunctionReturn(0);
656ae80bb75SLois Curfman McInnes }
657ae80bb75SLois Curfman McInnes 
658289bc588SBarry Smith /* -----------------------------------------------------------------*/
659289bc588SBarry Smith 
660e090d566SSatish Balay #include "petscsys.h"
66188e32edaSLois Curfman McInnes 
6624a2ae208SSatish Balay #undef __FUNCT__
6634a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
664f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A)
66588e32edaSLois Curfman McInnes {
66688e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
66788e32edaSLois Curfman McInnes   Mat            B;
6686849ba73SBarry Smith   PetscErrorCode ierr;
66913f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
67013f74950SBarry Smith   int            fd;
67113f74950SBarry Smith   PetscMPIInt    size;
67213f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
67387828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
67419bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
67588e32edaSLois Curfman McInnes 
6763a40ed3dSBarry Smith   PetscFunctionBegin;
677d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
67829bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
679b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6800752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
681552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
68288e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
68388e32edaSLois Curfman McInnes 
68490ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
685f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
686f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
687be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6885c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
68990ace30eSBarry Smith     B    = *A;
69090ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6914a41a4d3SSatish Balay     v    = a->v;
6924a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6934a41a4d3SSatish Balay        from row major to column major */
694b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
69590ace30eSBarry Smith     /* read in nonzero values */
6964a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6974a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6984a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6994a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7004a41a4d3SSatish Balay         *v++ =w[i*N+j];
7014a41a4d3SSatish Balay       }
7024a41a4d3SSatish Balay     }
703606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7046d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7056d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70690ace30eSBarry Smith   } else {
70788e32edaSLois Curfman McInnes     /* read row lengths */
70813f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7090752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
71088e32edaSLois Curfman McInnes 
71188e32edaSLois Curfman McInnes     /* create our matrix */
712f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
713f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
714be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7155c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
71688e32edaSLois Curfman McInnes     B = *A;
71788e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
71888e32edaSLois Curfman McInnes     v = a->v;
71988e32edaSLois Curfman McInnes 
72088e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
72113f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
722b0a32e0cSBarry Smith     cols = scols;
7230752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
72487828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
725b0a32e0cSBarry Smith     vals = svals;
7260752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
72788e32edaSLois Curfman McInnes 
72888e32edaSLois Curfman McInnes     /* insert into matrix */
72988e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
73088e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
73188e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
73288e32edaSLois Curfman McInnes     }
733606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
734606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
735606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
73688e32edaSLois Curfman McInnes 
7376d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7386d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73990ace30eSBarry Smith   }
7403a40ed3dSBarry Smith   PetscFunctionReturn(0);
74188e32edaSLois Curfman McInnes }
74288e32edaSLois Curfman McInnes 
743e090d566SSatish Balay #include "petscsys.h"
744932b0c3eSLois Curfman McInnes 
7454a2ae208SSatish Balay #undef __FUNCT__
7464a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
7476849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
748289bc588SBarry Smith {
749932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
750dfbe8321SBarry Smith   PetscErrorCode    ierr;
75113f74950SBarry Smith   PetscInt          i,j;
7522dcb1b2aSMatthew Knepley   const char        *name;
75387828ca2SBarry Smith   PetscScalar       *v;
754f3ef73ceSBarry Smith   PetscViewerFormat format;
755932b0c3eSLois Curfman McInnes 
7563a40ed3dSBarry Smith   PetscFunctionBegin;
757435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
758b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
759456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7603a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
761fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
762b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
763*899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
76444cd7ae7SLois Curfman McInnes       v = a->v + i;
76577431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
766*899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
767aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
768329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
769a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
770329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
771a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7726831982aSBarry Smith         }
77380cd9d93SLois Curfman McInnes #else
7746831982aSBarry Smith         if (*v) {
775a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
7766831982aSBarry Smith         }
77780cd9d93SLois Curfman McInnes #endif
7781b807ce4Svictorle         v += a->lda;
77980cd9d93SLois Curfman McInnes       }
780b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
78180cd9d93SLois Curfman McInnes     }
782b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7833a40ed3dSBarry Smith   } else {
784b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
785aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
786ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
78747989497SBarry Smith     /* determine if matrix has all real values */
78847989497SBarry Smith     v = a->v;
789*899cda47SBarry Smith     for (i=0; i<A->rmap.n*A->cmap.n; i++) {
790ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
79147989497SBarry Smith     }
79247989497SBarry Smith #endif
793fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
7943a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
795*899cda47SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr);
796*899cda47SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap.n,A->cmap.n);CHKERRQ(ierr);
797fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
798ffac6cdbSBarry Smith     }
799ffac6cdbSBarry Smith 
800*899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
801932b0c3eSLois Curfman McInnes       v = a->v + i;
802*899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
803aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
80447989497SBarry Smith         if (allreal) {
8051b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
80647989497SBarry Smith         } else {
8071b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
80847989497SBarry Smith         }
809289bc588SBarry Smith #else
8101b807ce4Svictorle         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
811289bc588SBarry Smith #endif
8121b807ce4Svictorle 	v += a->lda;
813289bc588SBarry Smith       }
814b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
815289bc588SBarry Smith     }
816fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
817b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
818ffac6cdbSBarry Smith     }
819b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
820da3a660dSBarry Smith   }
821b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8223a40ed3dSBarry Smith   PetscFunctionReturn(0);
823289bc588SBarry Smith }
824289bc588SBarry Smith 
8254a2ae208SSatish Balay #undef __FUNCT__
8264a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8276849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
828932b0c3eSLois Curfman McInnes {
829932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8306849ba73SBarry Smith   PetscErrorCode    ierr;
83113f74950SBarry Smith   int               fd;
832*899cda47SBarry Smith   PetscInt          ict,j,n = A->cmap.n,m = A->rmap.n,i,*col_lens,nz = m*n;
83387828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
834f3ef73ceSBarry Smith   PetscViewerFormat format;
835932b0c3eSLois Curfman McInnes 
8363a40ed3dSBarry Smith   PetscFunctionBegin;
837b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
83890ace30eSBarry Smith 
839b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
840fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
84190ace30eSBarry Smith     /* store the matrix as a dense matrix */
84213f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
843552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
84490ace30eSBarry Smith     col_lens[1] = m;
84590ace30eSBarry Smith     col_lens[2] = n;
84690ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8476f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
848606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
84990ace30eSBarry Smith 
85090ace30eSBarry Smith     /* write out matrix, by rows */
85187828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
85290ace30eSBarry Smith     v    = a->v;
85390ace30eSBarry Smith     for (i=0; i<m; i++) {
85490ace30eSBarry Smith       for (j=0; j<n; j++) {
85590ace30eSBarry Smith         vals[i + j*m] = *v++;
85690ace30eSBarry Smith       }
85790ace30eSBarry Smith     }
8586f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
859606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
86090ace30eSBarry Smith   } else {
86113f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
862552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
863932b0c3eSLois Curfman McInnes     col_lens[1] = m;
864932b0c3eSLois Curfman McInnes     col_lens[2] = n;
865932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
866932b0c3eSLois Curfman McInnes 
867932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
868932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8696f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
870932b0c3eSLois Curfman McInnes 
871932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
872932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
873932b0c3eSLois Curfman McInnes     ict = 0;
874932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
875932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
876932b0c3eSLois Curfman McInnes     }
8776f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
878606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
879932b0c3eSLois Curfman McInnes 
880932b0c3eSLois Curfman McInnes     /* store nonzero values */
88187828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
882932b0c3eSLois Curfman McInnes     ict  = 0;
883932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
884932b0c3eSLois Curfman McInnes       v = a->v + i;
885932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
8861b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
887932b0c3eSLois Curfman McInnes       }
888932b0c3eSLois Curfman McInnes     }
8896f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
890606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
89190ace30eSBarry Smith   }
8923a40ed3dSBarry Smith   PetscFunctionReturn(0);
893932b0c3eSLois Curfman McInnes }
894932b0c3eSLois Curfman McInnes 
8954a2ae208SSatish Balay #undef __FUNCT__
8964a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
897dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
898f1af5d2fSBarry Smith {
899f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
900f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9016849ba73SBarry Smith   PetscErrorCode    ierr;
902*899cda47SBarry Smith   PetscInt          m = A->rmap.n,n = A->cmap.n,color,i,j;
90387828ca2SBarry Smith   PetscScalar       *v = a->v;
904b0a32e0cSBarry Smith   PetscViewer       viewer;
905b0a32e0cSBarry Smith   PetscDraw         popup;
906329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
907f3ef73ceSBarry Smith   PetscViewerFormat format;
908f1af5d2fSBarry Smith 
909f1af5d2fSBarry Smith   PetscFunctionBegin;
910f1af5d2fSBarry Smith 
911f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
912b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
913b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
914f1af5d2fSBarry Smith 
915f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
916fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
917f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
918b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
919f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
920f1af5d2fSBarry Smith       x_l = j;
921f1af5d2fSBarry Smith       x_r = x_l + 1.0;
922f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
923f1af5d2fSBarry Smith         y_l = m - i - 1.0;
924f1af5d2fSBarry Smith         y_r = y_l + 1.0;
925f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
926329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
927b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
928329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
929b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
930f1af5d2fSBarry Smith         } else {
931f1af5d2fSBarry Smith           continue;
932f1af5d2fSBarry Smith         }
933f1af5d2fSBarry Smith #else
934f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
935b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
936f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
937b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
938f1af5d2fSBarry Smith         } else {
939f1af5d2fSBarry Smith           continue;
940f1af5d2fSBarry Smith         }
941f1af5d2fSBarry Smith #endif
942b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
943f1af5d2fSBarry Smith       }
944f1af5d2fSBarry Smith     }
945f1af5d2fSBarry Smith   } else {
946f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
947f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
948f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
949f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
950f1af5d2fSBarry Smith     }
951b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
952b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
953b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
954f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
955f1af5d2fSBarry Smith       x_l = j;
956f1af5d2fSBarry Smith       x_r = x_l + 1.0;
957f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
958f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
959f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
960b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
961b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
962f1af5d2fSBarry Smith       }
963f1af5d2fSBarry Smith     }
964f1af5d2fSBarry Smith   }
965f1af5d2fSBarry Smith   PetscFunctionReturn(0);
966f1af5d2fSBarry Smith }
967f1af5d2fSBarry Smith 
9684a2ae208SSatish Balay #undef __FUNCT__
9694a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
970dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
971f1af5d2fSBarry Smith {
972b0a32e0cSBarry Smith   PetscDraw      draw;
973f1af5d2fSBarry Smith   PetscTruth     isnull;
974329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
975dfbe8321SBarry Smith   PetscErrorCode ierr;
976f1af5d2fSBarry Smith 
977f1af5d2fSBarry Smith   PetscFunctionBegin;
978b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
979b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
980abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
981f1af5d2fSBarry Smith 
982f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
983*899cda47SBarry Smith   xr  = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
984f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
985b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
986b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
987f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
988f1af5d2fSBarry Smith   PetscFunctionReturn(0);
989f1af5d2fSBarry Smith }
990f1af5d2fSBarry Smith 
9914a2ae208SSatish Balay #undef __FUNCT__
9924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
993dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
994932b0c3eSLois Curfman McInnes {
995dfbe8321SBarry Smith   PetscErrorCode ierr;
99632077d6dSBarry Smith   PetscTruth     issocket,iascii,isbinary,isdraw;
997932b0c3eSLois Curfman McInnes 
9983a40ed3dSBarry Smith   PetscFunctionBegin;
999b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
100032077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1001fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1002fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10030f5bd95cSBarry Smith 
1004c45a1595SBarry Smith   if (iascii) {
1005c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10064cf18111SSatish Balay #if defined(PETSC_USE_SOCKET_VIEWER)
1007c45a1595SBarry Smith   } else if (issocket) {
10084cf18111SSatish Balay     Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1009*899cda47SBarry Smith     if (a->lda>A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA");
1010*899cda47SBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->rmap.n,A->cmap.n,a->v);CHKERRQ(ierr);
1011c45a1595SBarry Smith #endif
10120f5bd95cSBarry Smith   } else if (isbinary) {
10133a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1014f1af5d2fSBarry Smith   } else if (isdraw) {
1015f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10165cd90555SBarry Smith   } else {
1017958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1018932b0c3eSLois Curfman McInnes   }
10193a40ed3dSBarry Smith   PetscFunctionReturn(0);
1020932b0c3eSLois Curfman McInnes }
1021289bc588SBarry Smith 
10224a2ae208SSatish Balay #undef __FUNCT__
10234a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1024dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1025289bc588SBarry Smith {
1026ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1027dfbe8321SBarry Smith   PetscErrorCode ierr;
102890f02eecSBarry Smith 
10293a40ed3dSBarry Smith   PetscFunctionBegin;
1030aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1031*899cda47SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n);
1032a5a9c739SBarry Smith #endif
1033606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1034606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1035606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1036901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10373a40ed3dSBarry Smith   PetscFunctionReturn(0);
1038289bc588SBarry Smith }
1039289bc588SBarry Smith 
10404a2ae208SSatish Balay #undef __FUNCT__
10414a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1042dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1043289bc588SBarry Smith {
1044c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
10456849ba73SBarry Smith   PetscErrorCode ierr;
104613f74950SBarry Smith   PetscInt       k,j,m,n,M;
104787828ca2SBarry Smith   PetscScalar    *v,tmp;
104848b35521SBarry Smith 
10493a40ed3dSBarry Smith   PetscFunctionBegin;
1050*899cda47SBarry Smith   v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n;
10517c922b88SBarry Smith   if (!matout) { /* in place transpose */
1052a5ce6ee0Svictorle     if (m != n) {
1053634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1054d64ed03dSBarry Smith     } else {
1055d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1056289bc588SBarry Smith         for (k=0; k<j; k++) {
10571b807ce4Svictorle           tmp = v[j + k*M];
10581b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10591b807ce4Svictorle           v[k + j*M] = tmp;
1060289bc588SBarry Smith         }
1061289bc588SBarry Smith       }
1062d64ed03dSBarry Smith     }
10633a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1064d3e5ee88SLois Curfman McInnes     Mat          tmat;
1065ec8511deSBarry Smith     Mat_SeqDense *tmatd;
106687828ca2SBarry Smith     PetscScalar  *v2;
1067ea709b57SSatish Balay 
1068f69a0ea3SMatthew Knepley     ierr  = MatCreate(A->comm,&tmat);CHKERRQ(ierr);
1069*899cda47SBarry Smith     ierr  = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr);
10705c5985e7SKris Buschelman     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
10715c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1072ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10730de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1074d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10751b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1076d3e5ee88SLois Curfman McInnes     }
10776d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10786d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1079d3e5ee88SLois Curfman McInnes     *matout = tmat;
108048b35521SBarry Smith   }
10813a40ed3dSBarry Smith   PetscFunctionReturn(0);
1082289bc588SBarry Smith }
1083289bc588SBarry Smith 
10844a2ae208SSatish Balay #undef __FUNCT__
10854a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1086dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1087289bc588SBarry Smith {
1088c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1089c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
109013f74950SBarry Smith   PetscInt     i,j;
109187828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
10929ea5d5aeSSatish Balay 
10933a40ed3dSBarry Smith   PetscFunctionBegin;
1094*899cda47SBarry Smith   if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1095*899cda47SBarry Smith   if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1096*899cda47SBarry Smith   for (i=0; i<A1->rmap.n; i++) {
10971b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1098*899cda47SBarry Smith     for (j=0; j<A1->cmap.n; j++) {
10993a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11001b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11011b807ce4Svictorle     }
1102289bc588SBarry Smith   }
110377c4ece6SBarry Smith   *flg = PETSC_TRUE;
11043a40ed3dSBarry Smith   PetscFunctionReturn(0);
1105289bc588SBarry Smith }
1106289bc588SBarry Smith 
11074a2ae208SSatish Balay #undef __FUNCT__
11084a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1109dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1110289bc588SBarry Smith {
1111c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1112dfbe8321SBarry Smith   PetscErrorCode ierr;
111313f74950SBarry Smith   PetscInt       i,n,len;
111487828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
111544cd7ae7SLois Curfman McInnes 
11163a40ed3dSBarry Smith   PetscFunctionBegin;
11172dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11187a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11191ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1120*899cda47SBarry Smith   len = PetscMin(A->rmap.n,A->cmap.n);
1121*899cda47SBarry Smith   if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
112244cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11231b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1124289bc588SBarry Smith   }
11251ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11263a40ed3dSBarry Smith   PetscFunctionReturn(0);
1127289bc588SBarry Smith }
1128289bc588SBarry Smith 
11294a2ae208SSatish Balay #undef __FUNCT__
11304a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1131dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1132289bc588SBarry Smith {
1133c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
113487828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1135dfbe8321SBarry Smith   PetscErrorCode ierr;
1136*899cda47SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n;
113755659b69SBarry Smith 
11383a40ed3dSBarry Smith   PetscFunctionBegin;
113928988994SBarry Smith   if (ll) {
11407a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11411ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1142*899cda47SBarry Smith     if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1143da3a660dSBarry Smith     for (i=0; i<m; i++) {
1144da3a660dSBarry Smith       x = l[i];
1145da3a660dSBarry Smith       v = mat->v + i;
1146da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1147da3a660dSBarry Smith     }
11481ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1149efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1150da3a660dSBarry Smith   }
115128988994SBarry Smith   if (rr) {
11527a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11531ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1154*899cda47SBarry Smith     if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1155da3a660dSBarry Smith     for (i=0; i<n; i++) {
1156da3a660dSBarry Smith       x = r[i];
1157da3a660dSBarry Smith       v = mat->v + i*m;
1158da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1159da3a660dSBarry Smith     }
11601ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1161efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1162da3a660dSBarry Smith   }
11633a40ed3dSBarry Smith   PetscFunctionReturn(0);
1164289bc588SBarry Smith }
1165289bc588SBarry Smith 
11664a2ae208SSatish Balay #undef __FUNCT__
11674a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1168dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1169289bc588SBarry Smith {
1170c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
117187828ca2SBarry Smith   PetscScalar  *v = mat->v;
1172329f5518SBarry Smith   PetscReal    sum = 0.0;
1173*899cda47SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap.n,i,j;
1174efee365bSSatish Balay   PetscErrorCode ierr;
117555659b69SBarry Smith 
11763a40ed3dSBarry Smith   PetscFunctionBegin;
1177289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1178a5ce6ee0Svictorle     if (lda>m) {
1179*899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
1180a5ce6ee0Svictorle 	v = mat->v+j*lda;
1181a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1182a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1183a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1184a5ce6ee0Svictorle #else
1185a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1186a5ce6ee0Svictorle #endif
1187a5ce6ee0Svictorle 	}
1188a5ce6ee0Svictorle       }
1189a5ce6ee0Svictorle     } else {
1190*899cda47SBarry Smith       for (i=0; i<A->cmap.n*A->rmap.n; i++) {
1191aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1192329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1193289bc588SBarry Smith #else
1194289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1195289bc588SBarry Smith #endif
1196289bc588SBarry Smith       }
1197a5ce6ee0Svictorle     }
1198064f8208SBarry Smith     *nrm = sqrt(sum);
1199*899cda47SBarry Smith     ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12003a40ed3dSBarry Smith   } else if (type == NORM_1) {
1201064f8208SBarry Smith     *nrm = 0.0;
1202*899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
12031b807ce4Svictorle       v = mat->v + j*mat->lda;
1204289bc588SBarry Smith       sum = 0.0;
1205*899cda47SBarry Smith       for (i=0; i<A->rmap.n; i++) {
120633a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1207289bc588SBarry Smith       }
1208064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1209289bc588SBarry Smith     }
1210*899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12113a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1212064f8208SBarry Smith     *nrm = 0.0;
1213*899cda47SBarry Smith     for (j=0; j<A->rmap.n; j++) {
1214289bc588SBarry Smith       v = mat->v + j;
1215289bc588SBarry Smith       sum = 0.0;
1216*899cda47SBarry Smith       for (i=0; i<A->cmap.n; i++) {
12171b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1218289bc588SBarry Smith       }
1219064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1220289bc588SBarry Smith     }
1221*899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12223a40ed3dSBarry Smith   } else {
122329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1224289bc588SBarry Smith   }
12253a40ed3dSBarry Smith   PetscFunctionReturn(0);
1226289bc588SBarry Smith }
1227289bc588SBarry Smith 
12284a2ae208SSatish Balay #undef __FUNCT__
12294a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1230dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op)
1231289bc588SBarry Smith {
1232c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
123363ba0a88SBarry Smith   PetscErrorCode ierr;
123467e560aaSBarry Smith 
12353a40ed3dSBarry Smith   PetscFunctionBegin;
1236b5a2b587SKris Buschelman   switch (op) {
1237b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1238b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1239b5a2b587SKris Buschelman     break;
1240b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1241b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1242b5a2b587SKris Buschelman     break;
1243b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1244b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1245b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1246b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1247b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1248b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1249b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1250b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1251b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1252b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1253b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
1254ae15b995SBarry Smith     ierr = PetscInfo(A,"Option ignored\n");CHKERRQ(ierr);
1255b5a2b587SKris Buschelman     break;
125677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
125777e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12589a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12599a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12609a4540c5SBarry Smith   case MAT_HERMITIAN:
12619a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12629a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12639a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
126477e54ba9SKris Buschelman     break;
1265b5a2b587SKris Buschelman   default:
126629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
12673a40ed3dSBarry Smith   }
12683a40ed3dSBarry Smith   PetscFunctionReturn(0);
1269289bc588SBarry Smith }
1270289bc588SBarry Smith 
12714a2ae208SSatish Balay #undef __FUNCT__
12724a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1273dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
12746f0a148fSBarry Smith {
1275ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12766849ba73SBarry Smith   PetscErrorCode ierr;
1277*899cda47SBarry Smith   PetscInt       lda=l->lda,m=A->rmap.n,j;
12783a40ed3dSBarry Smith 
12793a40ed3dSBarry Smith   PetscFunctionBegin;
1280a5ce6ee0Svictorle   if (lda>m) {
1281*899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
1282a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1283a5ce6ee0Svictorle     }
1284a5ce6ee0Svictorle   } else {
1285*899cda47SBarry Smith     ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1286a5ce6ee0Svictorle   }
12873a40ed3dSBarry Smith   PetscFunctionReturn(0);
12886f0a148fSBarry Smith }
12896f0a148fSBarry Smith 
12904a2ae208SSatish Balay #undef __FUNCT__
12914a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1292f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
12936f0a148fSBarry Smith {
1294ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1295*899cda47SBarry Smith   PetscInt       n = A->cmap.n,i,j;
129687828ca2SBarry Smith   PetscScalar    *slot;
129755659b69SBarry Smith 
12983a40ed3dSBarry Smith   PetscFunctionBegin;
12996f0a148fSBarry Smith   for (i=0; i<N; i++) {
13006f0a148fSBarry Smith     slot = l->v + rows[i];
13016f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13026f0a148fSBarry Smith   }
1303f4df32b1SMatthew Knepley   if (diag != 0.0) {
13046f0a148fSBarry Smith     for (i=0; i<N; i++) {
13056f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1306f4df32b1SMatthew Knepley       *slot = diag;
13076f0a148fSBarry Smith     }
13086f0a148fSBarry Smith   }
13093a40ed3dSBarry Smith   PetscFunctionReturn(0);
13106f0a148fSBarry Smith }
1311557bce09SLois Curfman McInnes 
13124a2ae208SSatish Balay #undef __FUNCT__
13134a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1314dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
131564e87e97SBarry Smith {
1316c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13173a40ed3dSBarry Smith 
13183a40ed3dSBarry Smith   PetscFunctionBegin;
1319*899cda47SBarry Smith   if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
132064e87e97SBarry Smith   *array = mat->v;
13213a40ed3dSBarry Smith   PetscFunctionReturn(0);
132264e87e97SBarry Smith }
13230754003eSLois Curfman McInnes 
13244a2ae208SSatish Balay #undef __FUNCT__
13254a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1326dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1327ff14e315SSatish Balay {
13283a40ed3dSBarry Smith   PetscFunctionBegin;
132909b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1331ff14e315SSatish Balay }
13320754003eSLois Curfman McInnes 
13334a2ae208SSatish Balay #undef __FUNCT__
13344a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
133513f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13360754003eSLois Curfman McInnes {
1337c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13386849ba73SBarry Smith   PetscErrorCode ierr;
133921a2c019SBarry Smith   PetscInt       i,j,*irow,*icol,nrows,ncols;
134087828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13410754003eSLois Curfman McInnes   Mat            newmat;
13420754003eSLois Curfman McInnes 
13433a40ed3dSBarry Smith   PetscFunctionBegin;
134478b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
134578b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1346e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1347e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13480754003eSLois Curfman McInnes 
1349182d2002SSatish Balay   /* Check submatrixcall */
1350182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
135113f74950SBarry Smith     PetscInt n_cols,n_rows;
1352182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
135321a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
135421a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
135521a2c019SBarry Smith       ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr);
135621a2c019SBarry Smith     }
1357182d2002SSatish Balay     newmat = *B;
1358182d2002SSatish Balay   } else {
13590754003eSLois Curfman McInnes     /* Create and fill new matrix */
1360f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
1361f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
13625c5985e7SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
13635c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1364182d2002SSatish Balay   }
1365182d2002SSatish Balay 
1366182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1367182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1368182d2002SSatish Balay 
1369182d2002SSatish Balay   for (i=0; i<ncols; i++) {
13706de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1371182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1372182d2002SSatish Balay       *bv++ = av[irow[j]];
13730754003eSLois Curfman McInnes     }
13740754003eSLois Curfman McInnes   }
1375182d2002SSatish Balay 
1376182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13776d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13786d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13790754003eSLois Curfman McInnes 
13800754003eSLois Curfman McInnes   /* Free work space */
138178b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
138278b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1383182d2002SSatish Balay   *B = newmat;
13843a40ed3dSBarry Smith   PetscFunctionReturn(0);
13850754003eSLois Curfman McInnes }
13860754003eSLois Curfman McInnes 
13874a2ae208SSatish Balay #undef __FUNCT__
13884a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
138913f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1390905e6a2fSBarry Smith {
13916849ba73SBarry Smith   PetscErrorCode ierr;
139213f74950SBarry Smith   PetscInt       i;
1393905e6a2fSBarry Smith 
13943a40ed3dSBarry Smith   PetscFunctionBegin;
1395905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1396b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1397905e6a2fSBarry Smith   }
1398905e6a2fSBarry Smith 
1399905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14006a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1401905e6a2fSBarry Smith   }
14023a40ed3dSBarry Smith   PetscFunctionReturn(0);
1403905e6a2fSBarry Smith }
1404905e6a2fSBarry Smith 
14054a2ae208SSatish Balay #undef __FUNCT__
1406c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1407c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1408c0aa2d19SHong Zhang {
1409c0aa2d19SHong Zhang   PetscFunctionBegin;
1410c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1411c0aa2d19SHong Zhang }
1412c0aa2d19SHong Zhang 
1413c0aa2d19SHong Zhang #undef __FUNCT__
1414c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1415c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1416c0aa2d19SHong Zhang {
1417c0aa2d19SHong Zhang   PetscFunctionBegin;
1418c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1419c0aa2d19SHong Zhang }
1420c0aa2d19SHong Zhang 
1421c0aa2d19SHong Zhang #undef __FUNCT__
14224a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1423dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14244b0e389bSBarry Smith {
14254b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14266849ba73SBarry Smith   PetscErrorCode ierr;
1427*899cda47SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j;
14283a40ed3dSBarry Smith 
14293a40ed3dSBarry Smith   PetscFunctionBegin;
143033f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
143133f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1432cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14333a40ed3dSBarry Smith     PetscFunctionReturn(0);
14343a40ed3dSBarry Smith   }
1435*899cda47SBarry Smith   if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1436a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14370dbb7854Svictorle     for (j=0; j<n; j++) {
1438a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1439a5ce6ee0Svictorle     }
1440a5ce6ee0Svictorle   } else {
1441*899cda47SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1442a5ce6ee0Svictorle   }
1443273d9f13SBarry Smith   PetscFunctionReturn(0);
1444273d9f13SBarry Smith }
1445273d9f13SBarry Smith 
14464a2ae208SSatish Balay #undef __FUNCT__
14474a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1448dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1449273d9f13SBarry Smith {
1450dfbe8321SBarry Smith   PetscErrorCode ierr;
1451273d9f13SBarry Smith 
1452273d9f13SBarry Smith   PetscFunctionBegin;
1453273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14543a40ed3dSBarry Smith   PetscFunctionReturn(0);
14554b0e389bSBarry Smith }
14564b0e389bSBarry Smith 
1457284134d9SBarry Smith #undef __FUNCT__
1458284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1459284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1460284134d9SBarry Smith {
1461284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
146221a2c019SBarry Smith   PetscErrorCode ierr;
1463284134d9SBarry Smith   PetscFunctionBegin;
146421a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1465284134d9SBarry Smith   m = PetscMax(m,M);
1466284134d9SBarry Smith   n = PetscMax(n,N);
146721a2c019SBarry 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);
1468284134d9SBarry 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);
1469*899cda47SBarry Smith   A->rmap.n = A->rmap.n = m;
1470*899cda47SBarry Smith   A->cmap.n = A->cmap.N = n;
147121a2c019SBarry Smith   if (a->changelda) a->lda = m;
147221a2c019SBarry Smith   ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
1473284134d9SBarry Smith   PetscFunctionReturn(0);
1474284134d9SBarry Smith }
1475284134d9SBarry Smith 
1476a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1477a9fe9ddaSSatish Balay #undef __FUNCT__
1478a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1479a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1480a9fe9ddaSSatish Balay {
1481a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1482a9fe9ddaSSatish Balay 
1483a9fe9ddaSSatish Balay   PetscFunctionBegin;
1484a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1485a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1486a9fe9ddaSSatish Balay   }
1487a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1488a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1489a9fe9ddaSSatish Balay }
1490a9fe9ddaSSatish Balay 
1491a9fe9ddaSSatish Balay #undef __FUNCT__
1492a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1493a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1494a9fe9ddaSSatish Balay {
1495ee16a9a1SHong Zhang   PetscErrorCode ierr;
1496*899cda47SBarry Smith   PetscInt       m=A->rmap.n,n=B->cmap.n;
1497ee16a9a1SHong Zhang   Mat            Cmat;
1498a9fe9ddaSSatish Balay 
1499ee16a9a1SHong Zhang   PetscFunctionBegin;
1500*899cda47SBarry 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);
1501ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1502ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1503ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1504ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1505ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1506ee16a9a1SHong Zhang   *C = Cmat;
1507ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1508ee16a9a1SHong Zhang }
1509a9fe9ddaSSatish Balay 
1510a9fe9ddaSSatish Balay #undef __FUNCT__
1511a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1512a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1513a9fe9ddaSSatish Balay {
1514a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1515a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1516a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1517*899cda47SBarry Smith   PetscBLASInt   m=(PetscBLASInt)A->rmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->cmap.n;
1518a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1519a9fe9ddaSSatish Balay 
1520a9fe9ddaSSatish Balay   PetscFunctionBegin;
1521a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1522a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1523a9fe9ddaSSatish Balay }
1524a9fe9ddaSSatish Balay 
1525a9fe9ddaSSatish Balay #undef __FUNCT__
1526a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1527a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1528a9fe9ddaSSatish Balay {
1529a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1530a9fe9ddaSSatish Balay 
1531a9fe9ddaSSatish Balay   PetscFunctionBegin;
1532a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1533a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1534a9fe9ddaSSatish Balay   }
1535a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1536a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1537a9fe9ddaSSatish Balay }
1538a9fe9ddaSSatish Balay 
1539a9fe9ddaSSatish Balay #undef __FUNCT__
1540a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1541a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1542a9fe9ddaSSatish Balay {
1543ee16a9a1SHong Zhang   PetscErrorCode ierr;
1544*899cda47SBarry Smith   PetscInt       m=A->cmap.n,n=B->cmap.n;
1545ee16a9a1SHong Zhang   Mat            Cmat;
1546a9fe9ddaSSatish Balay 
1547ee16a9a1SHong Zhang   PetscFunctionBegin;
1548*899cda47SBarry 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);
1549ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1550ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1551ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1552ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1553ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1554ee16a9a1SHong Zhang   *C = Cmat;
1555ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1556ee16a9a1SHong Zhang }
1557a9fe9ddaSSatish Balay 
1558a9fe9ddaSSatish Balay #undef __FUNCT__
1559a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1560a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1561a9fe9ddaSSatish Balay {
1562a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1563a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1564a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1565*899cda47SBarry Smith   PetscBLASInt   m=(PetscBLASInt)A->cmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->rmap.n;
1566a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1567a9fe9ddaSSatish Balay 
1568a9fe9ddaSSatish Balay   PetscFunctionBegin;
1569a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1570a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1571a9fe9ddaSSatish Balay }
1572289bc588SBarry Smith /* -------------------------------------------------------------------*/
1573a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1574905e6a2fSBarry Smith        MatGetRow_SeqDense,
1575905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1576905e6a2fSBarry Smith        MatMult_SeqDense,
157797304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
15787c922b88SBarry Smith        MatMultTranspose_SeqDense,
15797c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1580905e6a2fSBarry Smith        MatSolve_SeqDense,
1581905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
15827c922b88SBarry Smith        MatSolveTranspose_SeqDense,
158397304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1584905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1585905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1586ec8511deSBarry Smith        MatRelax_SeqDense,
1587ec8511deSBarry Smith        MatTranspose_SeqDense,
158897304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1589905e6a2fSBarry Smith        MatEqual_SeqDense,
1590905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1591905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1592905e6a2fSBarry Smith        MatNorm_SeqDense,
1593c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1594c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1595905e6a2fSBarry Smith        0,
1596905e6a2fSBarry Smith        MatSetOption_SeqDense,
1597905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
159897304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1599905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1600905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1601905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1602905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
160397304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1604273d9f13SBarry Smith        0,
1605905e6a2fSBarry Smith        0,
1606905e6a2fSBarry Smith        MatGetArray_SeqDense,
1607905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
160897304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1609a5ae1ecdSBarry Smith        0,
1610a5ae1ecdSBarry Smith        0,
1611a5ae1ecdSBarry Smith        0,
1612a5ae1ecdSBarry Smith        0,
161397304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1614a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1615a5ae1ecdSBarry Smith        0,
16164b0e389bSBarry Smith        MatGetValues_SeqDense,
1617a5ae1ecdSBarry Smith        MatCopy_SeqDense,
161897304618SKris Buschelman /*45*/ 0,
1619a5ae1ecdSBarry Smith        MatScale_SeqDense,
1620a5ae1ecdSBarry Smith        0,
1621a5ae1ecdSBarry Smith        0,
1622a5ae1ecdSBarry Smith        0,
1623521d7252SBarry Smith /*50*/ 0,
1624a5ae1ecdSBarry Smith        0,
1625a5ae1ecdSBarry Smith        0,
1626a5ae1ecdSBarry Smith        0,
1627a5ae1ecdSBarry Smith        0,
162897304618SKris Buschelman /*55*/ 0,
1629a5ae1ecdSBarry Smith        0,
1630a5ae1ecdSBarry Smith        0,
1631a5ae1ecdSBarry Smith        0,
1632a5ae1ecdSBarry Smith        0,
163397304618SKris Buschelman /*60*/ 0,
1634e03a110bSBarry Smith        MatDestroy_SeqDense,
1635e03a110bSBarry Smith        MatView_SeqDense,
1636357abbc8SBarry Smith        0,
163797304618SKris Buschelman        0,
163897304618SKris Buschelman /*65*/ 0,
163997304618SKris Buschelman        0,
164097304618SKris Buschelman        0,
164197304618SKris Buschelman        0,
164297304618SKris Buschelman        0,
164397304618SKris Buschelman /*70*/ 0,
164497304618SKris Buschelman        0,
164597304618SKris Buschelman        0,
164697304618SKris Buschelman        0,
164797304618SKris Buschelman        0,
164897304618SKris Buschelman /*75*/ 0,
164997304618SKris Buschelman        0,
165097304618SKris Buschelman        0,
165197304618SKris Buschelman        0,
165297304618SKris Buschelman        0,
165397304618SKris Buschelman /*80*/ 0,
165497304618SKris Buschelman        0,
165597304618SKris Buschelman        0,
165697304618SKris Buschelman        0,
1657865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1658865e5f61SKris Buschelman        0,
1659865e5f61SKris Buschelman        0,
1660865e5f61SKris Buschelman        0,
1661865e5f61SKris Buschelman        0,
1662865e5f61SKris Buschelman        0,
1663a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense,
1664a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1665a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1666865e5f61SKris Buschelman        0,
1667865e5f61SKris Buschelman        0,
1668865e5f61SKris Buschelman /*95*/ 0,
1669a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1670a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1671a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1672284134d9SBarry Smith        0,
1673284134d9SBarry Smith /*100*/0,
1674284134d9SBarry Smith        0,
1675284134d9SBarry Smith        0,
1676284134d9SBarry Smith        0,
1677284134d9SBarry Smith        MatSetSizes_SeqDense};
167890ace30eSBarry Smith 
16794a2ae208SSatish Balay #undef __FUNCT__
16804a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
16814b828684SBarry Smith /*@C
1682fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1683d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1684d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1685289bc588SBarry Smith 
1686db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1687db81eaa0SLois Curfman McInnes 
168820563c6bSBarry Smith    Input Parameters:
1689db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
16900c775827SLois Curfman McInnes .  m - number of rows
169118f449edSLois Curfman McInnes .  n - number of columns
1692db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1693dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
169420563c6bSBarry Smith 
169520563c6bSBarry Smith    Output Parameter:
169644cd7ae7SLois Curfman McInnes .  A - the matrix
169720563c6bSBarry Smith 
1698b259b22eSLois Curfman McInnes    Notes:
169918f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
170018f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1701b4fd4287SBarry Smith    set data=PETSC_NULL.
170218f449edSLois Curfman McInnes 
1703027ccd11SLois Curfman McInnes    Level: intermediate
1704027ccd11SLois Curfman McInnes 
1705dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1706d65003e9SLois Curfman McInnes 
1707db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
170820563c6bSBarry Smith @*/
1709be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1710289bc588SBarry Smith {
1711dfbe8321SBarry Smith   PetscErrorCode ierr;
17123b2fbd54SBarry Smith 
17133a40ed3dSBarry Smith   PetscFunctionBegin;
1714f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1715f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1716273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1717273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1718273d9f13SBarry Smith   PetscFunctionReturn(0);
1719273d9f13SBarry Smith }
1720273d9f13SBarry Smith 
17214a2ae208SSatish Balay #undef __FUNCT__
17224a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1723273d9f13SBarry Smith /*@C
1724273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1725273d9f13SBarry Smith 
1726273d9f13SBarry Smith    Collective on MPI_Comm
1727273d9f13SBarry Smith 
1728273d9f13SBarry Smith    Input Parameters:
1729273d9f13SBarry Smith +  A - the matrix
1730273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1731273d9f13SBarry Smith 
1732273d9f13SBarry Smith    Notes:
1733273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1734273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1735284134d9SBarry Smith    need not call this routine.
1736273d9f13SBarry Smith 
1737273d9f13SBarry Smith    Level: intermediate
1738273d9f13SBarry Smith 
1739273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1740273d9f13SBarry Smith 
1741273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1742273d9f13SBarry Smith @*/
1743be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1744273d9f13SBarry Smith {
1745dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1746a23d5eceSKris Buschelman 
1747a23d5eceSKris Buschelman   PetscFunctionBegin;
1748a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1749a23d5eceSKris Buschelman   if (f) {
1750a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1751a23d5eceSKris Buschelman   }
1752a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1753a23d5eceSKris Buschelman }
1754a23d5eceSKris Buschelman 
1755a23d5eceSKris Buschelman EXTERN_C_BEGIN
1756a23d5eceSKris Buschelman #undef __FUNCT__
1757a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1758be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1759a23d5eceSKris Buschelman {
1760273d9f13SBarry Smith   Mat_SeqDense   *b;
1761dfbe8321SBarry Smith   PetscErrorCode ierr;
1762273d9f13SBarry Smith 
1763273d9f13SBarry Smith   PetscFunctionBegin;
1764273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1765273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1766273d9f13SBarry Smith   if (!data) {
1767284134d9SBarry Smith     ierr          = PetscMalloc((b->lda*b->Nmax+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1768284134d9SBarry Smith     ierr          = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1769273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
1770284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1771273d9f13SBarry Smith   } else { /* user-allocated storage */
1772273d9f13SBarry Smith     b->v          = data;
1773273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1774273d9f13SBarry Smith   }
1775273d9f13SBarry Smith   PetscFunctionReturn(0);
1776273d9f13SBarry Smith }
1777a23d5eceSKris Buschelman EXTERN_C_END
1778273d9f13SBarry Smith 
17791b807ce4Svictorle #undef __FUNCT__
17801b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
17811b807ce4Svictorle /*@C
17821b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
17831b807ce4Svictorle 
17841b807ce4Svictorle   Input parameter:
17851b807ce4Svictorle + A - the matrix
17861b807ce4Svictorle - lda - the leading dimension
17871b807ce4Svictorle 
17881b807ce4Svictorle   Notes:
17891b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
17901b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
17911b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
17921b807ce4Svictorle 
17931b807ce4Svictorle   Level: intermediate
17941b807ce4Svictorle 
17951b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
17961b807ce4Svictorle 
1797284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
17981b807ce4Svictorle @*/
1799be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
18001b807ce4Svictorle {
18011b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
180221a2c019SBarry Smith 
18031b807ce4Svictorle   PetscFunctionBegin;
1804*899cda47SBarry Smith   if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n);
18051b807ce4Svictorle   b->lda       = lda;
180621a2c019SBarry Smith   b->changelda = PETSC_FALSE;
180721a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
18081b807ce4Svictorle   PetscFunctionReturn(0);
18091b807ce4Svictorle }
18101b807ce4Svictorle 
18110bad9183SKris Buschelman /*MC
1812fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
18130bad9183SKris Buschelman 
18140bad9183SKris Buschelman    Options Database Keys:
18150bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
18160bad9183SKris Buschelman 
18170bad9183SKris Buschelman   Level: beginner
18180bad9183SKris Buschelman 
181989665df3SBarry Smith .seealso: MatCreateSeqDense()
182089665df3SBarry Smith 
18210bad9183SKris Buschelman M*/
18220bad9183SKris Buschelman 
1823273d9f13SBarry Smith EXTERN_C_BEGIN
18244a2ae208SSatish Balay #undef __FUNCT__
18254a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1826be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
1827273d9f13SBarry Smith {
1828273d9f13SBarry Smith   Mat_SeqDense   *b;
1829dfbe8321SBarry Smith   PetscErrorCode ierr;
18307c334f02SBarry Smith   PetscMPIInt    size;
1831273d9f13SBarry Smith 
1832273d9f13SBarry Smith   PetscFunctionBegin;
1833273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
183429bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
183555659b69SBarry Smith 
1836*899cda47SBarry Smith   B->rmap.bs = B->cmap.bs = 1;
1837*899cda47SBarry Smith   ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr);
1838*899cda47SBarry Smith   ierr = PetscMapInitialize(B->comm,&B->cmap);CHKERRQ(ierr);
1839273d9f13SBarry Smith 
1840b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1841549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
184244cd7ae7SLois Curfman McInnes   B->factor       = 0;
184390f02eecSBarry Smith   B->mapping      = 0;
184452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr);
184544cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
184618f449edSLois Curfman McInnes 
1847a5ae1ecdSBarry Smith 
184844cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1849273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1850273d9f13SBarry Smith   b->v            = 0;
1851*899cda47SBarry Smith   b->lda          = B->rmap.n;
185221a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
1853*899cda47SBarry Smith   b->Mmax         = B->rmap.n;
1854*899cda47SBarry Smith   b->Nmax         = B->cmap.n;
18554e220ebcSLois Curfman McInnes 
1856a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1857a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1858a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
18593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1860289bc588SBarry Smith }
1861c0aa2d19SHong Zhang 
1862c0aa2d19SHong Zhang 
1863273d9f13SBarry Smith EXTERN_C_END
1864