xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 6805f65bba57675c6478fc5d0aad24746736e0fa)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
8b4c862e0SSatish Balay #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
17899cda47SBarry Smith   PetscBLASInt   N = (PetscBLASInt)X->rmap.n*X->cmap.n,m=(PetscBLASInt)X->rmap.n,ldax = x->lda,lday=y->lda,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21a5ce6ee0Svictorle   if (ldax>m || lday>m) {
22899cda47SBarry Smith     for (j=0; j<X->cmap.n; j++) {
23f4df32b1SMatthew Knepley       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
24a5ce6ee0Svictorle     }
25a5ce6ee0Svictorle   } else {
26f4df32b1SMatthew Knepley     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
27a5ce6ee0Svictorle   }
28efee365bSSatish Balay   ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr);
293a40ed3dSBarry Smith   PetscFunctionReturn(0);
301987afe7SBarry Smith }
311987afe7SBarry Smith 
324a2ae208SSatish Balay #undef __FUNCT__
334a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
34dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
35289bc588SBarry Smith {
36899cda47SBarry Smith   PetscInt     N = A->rmap.n*A->cmap.n;
373a40ed3dSBarry Smith 
383a40ed3dSBarry Smith   PetscFunctionBegin;
39899cda47SBarry Smith   info->rows_global       = (double)A->rmap.n;
40899cda47SBarry Smith   info->columns_global    = (double)A->cmap.n;
41899cda47SBarry Smith   info->rows_local        = (double)A->rmap.n;
42899cda47SBarry Smith   info->columns_local     = (double)A->cmap.n;
434e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
444e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
456de62eeeSBarry Smith   info->nz_used           = (double)N;
466de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
474e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
484e220ebcSLois Curfman McInnes   info->mallocs           = 0;
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;
66899cda47SBarry Smith   if (lda>A->rmap.n) {
67899cda47SBarry Smith     nz = (PetscBLASInt)A->rmap.n;
68899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
69f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
72899cda47SBarry Smith     nz = (PetscBLASInt)A->rmap.n*A->cmap.n;
73f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
74a5ce6ee0Svictorle   }
75efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
7780cd9d93SLois Curfman McInnes }
7880cd9d93SLois Curfman McInnes 
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;
92899cda47SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,m = (PetscBLASInt)A->rmap.n,info;
9355659b69SBarry Smith 
943a40ed3dSBarry Smith   PetscFunctionBegin;
95289bc588SBarry Smith   if (!mat->pivots) {
96899cda47SBarry Smith     ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
97899cda47SBarry Smith     ierr = PetscLogObjectMemory(A,A->rmap.n*sizeof(PetscBLASInt));CHKERRQ(ierr);
98289bc588SBarry Smith   }
99f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
100899cda47SBarry 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");
104899cda47SBarry 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);
120899cda47SBarry 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;
125899cda47SBarry Smith     if (lda>A->rmap.n) {
126899cda47SBarry Smith       m = A->rmap.n;
127899cda47SBarry 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 {
131899cda47SBarry 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;
156899cda47SBarry 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 {
166899cda47SBarry 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;
194899cda47SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->cmap.n,info;
19555659b69SBarry Smith 
1963a40ed3dSBarry Smith   PetscFunctionBegin;
197606d414cSSatish Balay   ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
198752f5567SLois Curfman McInnes   mat->pivots = 0;
19905b42c5fSBarry Smith 
200899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
20171044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
20277431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
203c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
204899cda47SBarry Smith   ierr = PetscLogFlops((A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr);
205a07ab388SBarry Smith #endif
2063a40ed3dSBarry Smith   PetscFunctionReturn(0);
207289bc588SBarry Smith }
208289bc588SBarry Smith 
2094a2ae208SSatish Balay #undef __FUNCT__
2100b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
211af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2120b4b3355SBarry Smith {
213dfbe8321SBarry Smith   PetscErrorCode ierr;
21415e8a5b3SHong Zhang   MatFactorInfo  info;
2150b4b3355SBarry Smith 
2160b4b3355SBarry Smith   PetscFunctionBegin;
21715e8a5b3SHong Zhang   info.fill = 1.0;
21815e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2190b4b3355SBarry Smith   PetscFunctionReturn(0);
2200b4b3355SBarry Smith }
2210b4b3355SBarry Smith 
2220b4b3355SBarry Smith #undef __FUNCT__
2234a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
224dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
225289bc588SBarry Smith {
226c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2276849ba73SBarry Smith   PetscErrorCode ierr;
228899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, one = 1,info;
22987828ca2SBarry Smith   PetscScalar    *x,*y;
23067e560aaSBarry Smith 
2313a40ed3dSBarry Smith   PetscFunctionBegin;
2321ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2331ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
234899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
235c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
236ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
23780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
238ae7cfcebSSatish Balay #else
23971044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2404ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
241ae7cfcebSSatish Balay #endif
2427a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
243ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
24480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
245ae7cfcebSSatish Balay #else
24671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2474ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
248ae7cfcebSSatish Balay #endif
249289bc588SBarry Smith   }
25029bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2511ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2521ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
253899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
2543a40ed3dSBarry Smith   PetscFunctionReturn(0);
255289bc588SBarry Smith }
2566ee01492SSatish Balay 
2574a2ae208SSatish Balay #undef __FUNCT__
2584a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
259dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
260da3a660dSBarry Smith {
261c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
262dfbe8321SBarry Smith   PetscErrorCode ierr;
263899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt) A->rmap.n,one = 1,info;
26487828ca2SBarry Smith   PetscScalar    *x,*y;
26567e560aaSBarry Smith 
2663a40ed3dSBarry Smith   PetscFunctionBegin;
2671ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2681ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
269899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
270752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
271da3a660dSBarry Smith   if (mat->pivots) {
272ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
27380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
274ae7cfcebSSatish Balay #else
27571044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2764ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
277ae7cfcebSSatish Balay #endif
2787a97a34bSBarry Smith   } else {
279ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
28080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
281ae7cfcebSSatish Balay #else
28271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2834ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
284ae7cfcebSSatish Balay #endif
285da3a660dSBarry Smith   }
2861ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2871ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
288899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
290da3a660dSBarry Smith }
2916ee01492SSatish Balay 
2924a2ae208SSatish Balay #undef __FUNCT__
2934a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
294dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
295da3a660dSBarry Smith {
296c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
297dfbe8321SBarry Smith   PetscErrorCode ierr;
298899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
29987828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
300da3a660dSBarry Smith   Vec            tmp = 0;
30167e560aaSBarry Smith 
3023a40ed3dSBarry Smith   PetscFunctionBegin;
3031ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3041ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
305899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
306da3a660dSBarry Smith   if (yy == zz) {
30778b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
30852e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
30978b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
310da3a660dSBarry Smith   }
311899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
312752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
313da3a660dSBarry Smith   if (mat->pivots) {
314ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
31580444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
316ae7cfcebSSatish Balay #else
31771044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3182ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
319ae7cfcebSSatish Balay #endif
320a8c6a408SBarry Smith   } else {
321ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
32280444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
323ae7cfcebSSatish Balay #else
32471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3252ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
326ae7cfcebSSatish Balay #endif
327da3a660dSBarry Smith   }
3282dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
3292dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
3301ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3311ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
332899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
3333a40ed3dSBarry Smith   PetscFunctionReturn(0);
334da3a660dSBarry Smith }
33567e560aaSBarry Smith 
3364a2ae208SSatish Balay #undef __FUNCT__
3374a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
338dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
339da3a660dSBarry Smith {
340c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3416849ba73SBarry Smith   PetscErrorCode ierr;
342899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n,one = 1,info;
34387828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
344da3a660dSBarry Smith   Vec            tmp;
34567e560aaSBarry Smith 
3463a40ed3dSBarry Smith   PetscFunctionBegin;
347899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
3481ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3491ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
350da3a660dSBarry Smith   if (yy == zz) {
35178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
35252e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
35378b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
354da3a660dSBarry Smith   }
355899cda47SBarry Smith   ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
356752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
357da3a660dSBarry Smith   if (mat->pivots) {
358ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
35980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
360ae7cfcebSSatish Balay #else
36171044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3622ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
363ae7cfcebSSatish Balay #endif
3643a40ed3dSBarry Smith   } else {
365ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
36680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
367ae7cfcebSSatish Balay #else
36871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3692ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
370ae7cfcebSSatish Balay #endif
371da3a660dSBarry Smith   }
37290f02eecSBarry Smith   if (tmp) {
3732dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
37490f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3753a40ed3dSBarry Smith   } else {
3762dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
37790f02eecSBarry Smith   }
3781ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3791ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
380899cda47SBarry Smith   ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr);
3813a40ed3dSBarry Smith   PetscFunctionReturn(0);
382da3a660dSBarry Smith }
383289bc588SBarry Smith /* ------------------------------------------------------------------*/
3844a2ae208SSatish Balay #undef __FUNCT__
3854a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
38613f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
387289bc588SBarry Smith {
388c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
38987828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
390dfbe8321SBarry Smith   PetscErrorCode ierr;
391899cda47SBarry Smith   PetscInt       m = A->rmap.n,i;
392aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3934ce68768SBarry Smith   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
394bc1b551cSSatish Balay #endif
395289bc588SBarry Smith 
3963a40ed3dSBarry Smith   PetscFunctionBegin;
397289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
39871044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
3992dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
400289bc588SBarry Smith   }
4011ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4021ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
403b965ef7fSBarry Smith   its  = its*lits;
40477431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
405289bc588SBarry Smith   while (its--) {
406fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
407289bc588SBarry Smith       for (i=0; i<m; i++) {
408aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
409f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
410f1747703SBarry Smith            not happy about returning a double complex */
41113f74950SBarry Smith         PetscInt         _i;
41287828ca2SBarry Smith         PetscScalar sum = b[i];
413f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4143f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
415f1747703SBarry Smith         }
416f1747703SBarry Smith         xt = sum;
417f1747703SBarry Smith #else
41871044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
419f1747703SBarry Smith #endif
42055a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
421289bc588SBarry Smith       }
422289bc588SBarry Smith     }
423fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
424289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
425aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
426f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
427f1747703SBarry Smith            not happy about returning a double complex */
42813f74950SBarry Smith         PetscInt         _i;
42987828ca2SBarry Smith         PetscScalar sum = b[i];
430f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4313f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
432f1747703SBarry Smith         }
433f1747703SBarry Smith         xt = sum;
434f1747703SBarry Smith #else
43571044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
436f1747703SBarry Smith #endif
43755a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
438289bc588SBarry Smith       }
439289bc588SBarry Smith     }
440289bc588SBarry Smith   }
4411ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4421ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4433a40ed3dSBarry Smith   PetscFunctionReturn(0);
444289bc588SBarry Smith }
445289bc588SBarry Smith 
446289bc588SBarry Smith /* -----------------------------------------------------------------*/
4474a2ae208SSatish Balay #undef __FUNCT__
4484a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
449dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
450289bc588SBarry Smith {
451c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
45287828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
453dfbe8321SBarry Smith   PetscErrorCode ierr;
454899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n,_One=1;
455ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
4563a40ed3dSBarry Smith 
4573a40ed3dSBarry Smith   PetscFunctionBegin;
458899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
4591ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4601ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
46171044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
4621ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4631ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
464899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr);
4653a40ed3dSBarry Smith   PetscFunctionReturn(0);
466289bc588SBarry Smith }
4676ee01492SSatish Balay 
4684a2ae208SSatish Balay #undef __FUNCT__
4694a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
470dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
471289bc588SBarry Smith {
472c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
47387828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
474dfbe8321SBarry Smith   PetscErrorCode ierr;
475899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
4763a40ed3dSBarry Smith 
4773a40ed3dSBarry Smith   PetscFunctionBegin;
478899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
4791ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4801ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
48171044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
4821ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4831ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
484899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->rmap.n);CHKERRQ(ierr);
4853a40ed3dSBarry Smith   PetscFunctionReturn(0);
486289bc588SBarry Smith }
4876ee01492SSatish Balay 
4884a2ae208SSatish Balay #undef __FUNCT__
4894a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
490dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
491289bc588SBarry Smith {
492c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
49387828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
494dfbe8321SBarry Smith   PetscErrorCode ierr;
495899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
4963a40ed3dSBarry Smith 
4973a40ed3dSBarry Smith   PetscFunctionBegin;
498899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
499600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5001ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5011ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
50271044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5031ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5041ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
505899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5063a40ed3dSBarry Smith   PetscFunctionReturn(0);
507289bc588SBarry Smith }
5086ee01492SSatish Balay 
5094a2ae208SSatish Balay #undef __FUNCT__
5104a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
511dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
512289bc588SBarry Smith {
513c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
51487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
515dfbe8321SBarry Smith   PetscErrorCode ierr;
516899cda47SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1;
51787828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5183a40ed3dSBarry Smith 
5193a40ed3dSBarry Smith   PetscFunctionBegin;
520899cda47SBarry Smith   if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0);
521600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5221ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5231ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
52471044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5251ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5261ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
527899cda47SBarry Smith   ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr);
5283a40ed3dSBarry Smith   PetscFunctionReturn(0);
529289bc588SBarry Smith }
530289bc588SBarry Smith 
531289bc588SBarry Smith /* -----------------------------------------------------------------*/
5324a2ae208SSatish Balay #undef __FUNCT__
5334a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
53413f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
535289bc588SBarry Smith {
536c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
53787828ca2SBarry Smith   PetscScalar    *v;
5386849ba73SBarry Smith   PetscErrorCode ierr;
53913f74950SBarry Smith   PetscInt       i;
54067e560aaSBarry Smith 
5413a40ed3dSBarry Smith   PetscFunctionBegin;
542899cda47SBarry Smith   *ncols = A->cmap.n;
543289bc588SBarry Smith   if (cols) {
544899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
545899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) (*cols)[i] = i;
546289bc588SBarry Smith   }
547289bc588SBarry Smith   if (vals) {
548899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
549289bc588SBarry Smith     v    = mat->v + row;
550899cda47SBarry Smith     for (i=0; i<A->cmap.n; i++) {(*vals)[i] = *v; v += mat->lda;}
551289bc588SBarry Smith   }
5523a40ed3dSBarry Smith   PetscFunctionReturn(0);
553289bc588SBarry Smith }
5546ee01492SSatish Balay 
5554a2ae208SSatish Balay #undef __FUNCT__
5564a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
55713f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
558289bc588SBarry Smith {
559dfbe8321SBarry Smith   PetscErrorCode ierr;
560606d414cSSatish Balay   PetscFunctionBegin;
561606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
562606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5633a40ed3dSBarry Smith   PetscFunctionReturn(0);
564289bc588SBarry Smith }
565289bc588SBarry Smith /* ----------------------------------------------------------------*/
5664a2ae208SSatish Balay #undef __FUNCT__
5674a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
56813f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
569289bc588SBarry Smith {
570c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
57113f74950SBarry Smith   PetscInt     i,j,idx=0;
572d6dfbf8fSBarry Smith 
5733a40ed3dSBarry Smith   PetscFunctionBegin;
574289bc588SBarry Smith   if (!mat->roworiented) {
575dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
576289bc588SBarry Smith       for (j=0; j<n; j++) {
577cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5782515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
579899cda47SBarry 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);
58058804f6dSBarry Smith #endif
581289bc588SBarry Smith         for (i=0; i<m; i++) {
582cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5832515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
584899cda47SBarry 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);
58558804f6dSBarry Smith #endif
586cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
587289bc588SBarry Smith         }
588289bc588SBarry Smith       }
5893a40ed3dSBarry Smith     } else {
590289bc588SBarry Smith       for (j=0; j<n; j++) {
591cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5922515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
593899cda47SBarry 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);
59458804f6dSBarry Smith #endif
595289bc588SBarry Smith         for (i=0; i<m; i++) {
596cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5972515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
598899cda47SBarry 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);
59958804f6dSBarry Smith #endif
600cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
601289bc588SBarry Smith         }
602289bc588SBarry Smith       }
603289bc588SBarry Smith     }
6043a40ed3dSBarry Smith   } else {
605dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
606e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
607cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6082515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
609899cda47SBarry 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);
61058804f6dSBarry Smith #endif
611e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
612cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6132515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
614899cda47SBarry 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);
61558804f6dSBarry Smith #endif
616cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
617e8d4e0b9SBarry Smith         }
618e8d4e0b9SBarry Smith       }
6193a40ed3dSBarry Smith     } else {
620289bc588SBarry Smith       for (i=0; i<m; i++) {
621cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6222515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
623899cda47SBarry 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);
62458804f6dSBarry Smith #endif
625289bc588SBarry Smith         for (j=0; j<n; j++) {
626cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6272515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
628899cda47SBarry 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);
62958804f6dSBarry Smith #endif
630cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
631289bc588SBarry Smith         }
632289bc588SBarry Smith       }
633289bc588SBarry Smith     }
634e8d4e0b9SBarry Smith   }
6353a40ed3dSBarry Smith   PetscFunctionReturn(0);
636289bc588SBarry Smith }
637e8d4e0b9SBarry Smith 
6384a2ae208SSatish Balay #undef __FUNCT__
6394a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
64013f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
641ae80bb75SLois Curfman McInnes {
642ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
64313f74950SBarry Smith   PetscInt     i,j;
64487828ca2SBarry Smith   PetscScalar  *vpt = v;
645ae80bb75SLois Curfman McInnes 
6463a40ed3dSBarry Smith   PetscFunctionBegin;
647ae80bb75SLois Curfman McInnes   /* row-oriented output */
648ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
649ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6501b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
651ae80bb75SLois Curfman McInnes     }
652ae80bb75SLois Curfman McInnes   }
6533a40ed3dSBarry Smith   PetscFunctionReturn(0);
654ae80bb75SLois Curfman McInnes }
655ae80bb75SLois Curfman McInnes 
656289bc588SBarry Smith /* -----------------------------------------------------------------*/
657289bc588SBarry Smith 
658e090d566SSatish Balay #include "petscsys.h"
65988e32edaSLois Curfman McInnes 
6604a2ae208SSatish Balay #undef __FUNCT__
6614a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
662f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A)
66388e32edaSLois Curfman McInnes {
66488e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
66588e32edaSLois Curfman McInnes   Mat            B;
6666849ba73SBarry Smith   PetscErrorCode ierr;
66713f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
66813f74950SBarry Smith   int            fd;
66913f74950SBarry Smith   PetscMPIInt    size;
67013f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
67187828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
67219bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
67388e32edaSLois Curfman McInnes 
6743a40ed3dSBarry Smith   PetscFunctionBegin;
675d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
67629bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
677b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6780752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
679552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
68088e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
68188e32edaSLois Curfman McInnes 
68290ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
683f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
684f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
685be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6865c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
68790ace30eSBarry Smith     B    = *A;
68890ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6894a41a4d3SSatish Balay     v    = a->v;
6904a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6914a41a4d3SSatish Balay        from row major to column major */
692b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
69390ace30eSBarry Smith     /* read in nonzero values */
6944a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6954a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6964a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6974a41a4d3SSatish Balay       for (i=0; i<M; i++) {
6984a41a4d3SSatish Balay         *v++ =w[i*N+j];
6994a41a4d3SSatish Balay       }
7004a41a4d3SSatish Balay     }
701606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7026d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7036d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70490ace30eSBarry Smith   } else {
70588e32edaSLois Curfman McInnes     /* read row lengths */
70613f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7070752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
70888e32edaSLois Curfman McInnes 
70988e32edaSLois Curfman McInnes     /* create our matrix */
710f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
711f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
712be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7135c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
71488e32edaSLois Curfman McInnes     B = *A;
71588e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
71688e32edaSLois Curfman McInnes     v = a->v;
71788e32edaSLois Curfman McInnes 
71888e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
71913f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
720b0a32e0cSBarry Smith     cols = scols;
7210752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
72287828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
723b0a32e0cSBarry Smith     vals = svals;
7240752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
72588e32edaSLois Curfman McInnes 
72688e32edaSLois Curfman McInnes     /* insert into matrix */
72788e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
72888e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
72988e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
73088e32edaSLois Curfman McInnes     }
731606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
732606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
733606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
73488e32edaSLois Curfman McInnes 
7356d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7366d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73790ace30eSBarry Smith   }
7383a40ed3dSBarry Smith   PetscFunctionReturn(0);
73988e32edaSLois Curfman McInnes }
74088e32edaSLois Curfman McInnes 
741e090d566SSatish Balay #include "petscsys.h"
742932b0c3eSLois Curfman McInnes 
7434a2ae208SSatish Balay #undef __FUNCT__
7444a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
7456849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
746289bc588SBarry Smith {
747932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
748dfbe8321SBarry Smith   PetscErrorCode    ierr;
74913f74950SBarry Smith   PetscInt          i,j;
7502dcb1b2aSMatthew Knepley   const char        *name;
75187828ca2SBarry Smith   PetscScalar       *v;
752f3ef73ceSBarry Smith   PetscViewerFormat format;
7535f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
7545f481a85SSatish Balay   PetscTruth allreal = PETSC_TRUE;
7555f481a85SSatish Balay #endif
756932b0c3eSLois Curfman McInnes 
7573a40ed3dSBarry Smith   PetscFunctionBegin;
758435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
759b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
760456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7613a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
762fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
763b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
764899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
76544cd7ae7SLois Curfman McInnes       v = a->v + i;
76677431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
767899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
768aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
769329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
770a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
771329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
772a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7736831982aSBarry Smith         }
77480cd9d93SLois Curfman McInnes #else
7756831982aSBarry Smith         if (*v) {
776a83599f4SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr);
7776831982aSBarry Smith         }
77880cd9d93SLois Curfman McInnes #endif
7791b807ce4Svictorle         v += a->lda;
78080cd9d93SLois Curfman McInnes       }
781b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
78280cd9d93SLois Curfman McInnes     }
783b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7843a40ed3dSBarry Smith   } else {
785b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
786aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
78747989497SBarry Smith     /* determine if matrix has all real values */
78847989497SBarry Smith     v = a->v;
789899cda47SBarry 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);
795899cda47SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr);
796899cda47SBarry 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 
800899cda47SBarry Smith     for (i=0; i<A->rmap.n; i++) {
801932b0c3eSLois Curfman McInnes       v = a->v + i;
802899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
803aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
80447989497SBarry Smith         if (allreal) {
805f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr);
80647989497SBarry Smith         } else {
807f32d5d43SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
80847989497SBarry Smith         }
809289bc588SBarry Smith #else
810f32d5d43SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*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;
832899cda47SBarry 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;
902899cda47SBarry 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);
983899cda47SBarry 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;
996*6805f65bSBarry Smith   PetscTruth     iascii,isbinary,isdraw;
997932b0c3eSLois Curfman McInnes 
9983a40ed3dSBarry Smith   PetscFunctionBegin;
99932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1000fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1001fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10020f5bd95cSBarry Smith 
1003c45a1595SBarry Smith   if (iascii) {
1004c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10050f5bd95cSBarry Smith   } else if (isbinary) {
10063a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1007f1af5d2fSBarry Smith   } else if (isdraw) {
1008f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10095cd90555SBarry Smith   } else {
1010958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1011932b0c3eSLois Curfman McInnes   }
10123a40ed3dSBarry Smith   PetscFunctionReturn(0);
1013932b0c3eSLois Curfman McInnes }
1014289bc588SBarry Smith 
10154a2ae208SSatish Balay #undef __FUNCT__
10164a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1017dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1018289bc588SBarry Smith {
1019ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1020dfbe8321SBarry Smith   PetscErrorCode ierr;
102190f02eecSBarry Smith 
10223a40ed3dSBarry Smith   PetscFunctionBegin;
1023aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1024899cda47SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n);
1025a5a9c739SBarry Smith #endif
102605b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
10276857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1028606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1029901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10304ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10314ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10324ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr);
10333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1034289bc588SBarry Smith }
1035289bc588SBarry Smith 
10364a2ae208SSatish Balay #undef __FUNCT__
10374a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1038dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1039289bc588SBarry Smith {
1040c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
10416849ba73SBarry Smith   PetscErrorCode ierr;
104213f74950SBarry Smith   PetscInt       k,j,m,n,M;
104387828ca2SBarry Smith   PetscScalar    *v,tmp;
104448b35521SBarry Smith 
10453a40ed3dSBarry Smith   PetscFunctionBegin;
1046899cda47SBarry Smith   v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n;
10477c922b88SBarry Smith   if (!matout) { /* in place transpose */
1048a5ce6ee0Svictorle     if (m != n) {
1049634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1050d64ed03dSBarry Smith     } else {
1051d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1052289bc588SBarry Smith         for (k=0; k<j; k++) {
10531b807ce4Svictorle           tmp = v[j + k*M];
10541b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10551b807ce4Svictorle           v[k + j*M] = tmp;
1056289bc588SBarry Smith         }
1057289bc588SBarry Smith       }
1058d64ed03dSBarry Smith     }
10593a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1060d3e5ee88SLois Curfman McInnes     Mat          tmat;
1061ec8511deSBarry Smith     Mat_SeqDense *tmatd;
106287828ca2SBarry Smith     PetscScalar  *v2;
1063ea709b57SSatish Balay 
1064f69a0ea3SMatthew Knepley     ierr  = MatCreate(A->comm,&tmat);CHKERRQ(ierr);
1065899cda47SBarry Smith     ierr  = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr);
10665c5985e7SKris Buschelman     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
10675c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1068ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10690de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1070d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10711b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1072d3e5ee88SLois Curfman McInnes     }
10736d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10746d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1075d3e5ee88SLois Curfman McInnes     *matout = tmat;
107648b35521SBarry Smith   }
10773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1078289bc588SBarry Smith }
1079289bc588SBarry Smith 
10804a2ae208SSatish Balay #undef __FUNCT__
10814a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1082dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1083289bc588SBarry Smith {
1084c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1085c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
108613f74950SBarry Smith   PetscInt     i,j;
108787828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
10889ea5d5aeSSatish Balay 
10893a40ed3dSBarry Smith   PetscFunctionBegin;
1090899cda47SBarry Smith   if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1091899cda47SBarry Smith   if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1092899cda47SBarry Smith   for (i=0; i<A1->rmap.n; i++) {
10931b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1094899cda47SBarry Smith     for (j=0; j<A1->cmap.n; j++) {
10953a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10961b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
10971b807ce4Svictorle     }
1098289bc588SBarry Smith   }
109977c4ece6SBarry Smith   *flg = PETSC_TRUE;
11003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1101289bc588SBarry Smith }
1102289bc588SBarry Smith 
11034a2ae208SSatish Balay #undef __FUNCT__
11044a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1105dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1106289bc588SBarry Smith {
1107c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1108dfbe8321SBarry Smith   PetscErrorCode ierr;
110913f74950SBarry Smith   PetscInt       i,n,len;
111087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
111144cd7ae7SLois Curfman McInnes 
11123a40ed3dSBarry Smith   PetscFunctionBegin;
11132dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11147a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11151ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1116899cda47SBarry Smith   len = PetscMin(A->rmap.n,A->cmap.n);
1117899cda47SBarry Smith   if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
111844cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11191b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1120289bc588SBarry Smith   }
11211ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11223a40ed3dSBarry Smith   PetscFunctionReturn(0);
1123289bc588SBarry Smith }
1124289bc588SBarry Smith 
11254a2ae208SSatish Balay #undef __FUNCT__
11264a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1127dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1128289bc588SBarry Smith {
1129c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
113087828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1131dfbe8321SBarry Smith   PetscErrorCode ierr;
1132899cda47SBarry Smith   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n;
113355659b69SBarry Smith 
11343a40ed3dSBarry Smith   PetscFunctionBegin;
113528988994SBarry Smith   if (ll) {
11367a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11371ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1138899cda47SBarry Smith     if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1139da3a660dSBarry Smith     for (i=0; i<m; i++) {
1140da3a660dSBarry Smith       x = l[i];
1141da3a660dSBarry Smith       v = mat->v + i;
1142da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1143da3a660dSBarry Smith     }
11441ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1145efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1146da3a660dSBarry Smith   }
114728988994SBarry Smith   if (rr) {
11487a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11491ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1150899cda47SBarry Smith     if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1151da3a660dSBarry Smith     for (i=0; i<n; i++) {
1152da3a660dSBarry Smith       x = r[i];
1153da3a660dSBarry Smith       v = mat->v + i*m;
1154da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1155da3a660dSBarry Smith     }
11561ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1157efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1158da3a660dSBarry Smith   }
11593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1160289bc588SBarry Smith }
1161289bc588SBarry Smith 
11624a2ae208SSatish Balay #undef __FUNCT__
11634a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1164dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1165289bc588SBarry Smith {
1166c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
116787828ca2SBarry Smith   PetscScalar  *v = mat->v;
1168329f5518SBarry Smith   PetscReal    sum = 0.0;
1169899cda47SBarry Smith   PetscInt     lda=mat->lda,m=A->rmap.n,i,j;
1170efee365bSSatish Balay   PetscErrorCode ierr;
117155659b69SBarry Smith 
11723a40ed3dSBarry Smith   PetscFunctionBegin;
1173289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1174a5ce6ee0Svictorle     if (lda>m) {
1175899cda47SBarry Smith       for (j=0; j<A->cmap.n; j++) {
1176a5ce6ee0Svictorle 	v = mat->v+j*lda;
1177a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1178a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1179a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1180a5ce6ee0Svictorle #else
1181a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1182a5ce6ee0Svictorle #endif
1183a5ce6ee0Svictorle 	}
1184a5ce6ee0Svictorle       }
1185a5ce6ee0Svictorle     } else {
1186899cda47SBarry Smith       for (i=0; i<A->cmap.n*A->rmap.n; i++) {
1187aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1188329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1189289bc588SBarry Smith #else
1190289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1191289bc588SBarry Smith #endif
1192289bc588SBarry Smith       }
1193a5ce6ee0Svictorle     }
1194064f8208SBarry Smith     *nrm = sqrt(sum);
1195899cda47SBarry Smith     ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr);
11963a40ed3dSBarry Smith   } else if (type == NORM_1) {
1197064f8208SBarry Smith     *nrm = 0.0;
1198899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
11991b807ce4Svictorle       v = mat->v + j*mat->lda;
1200289bc588SBarry Smith       sum = 0.0;
1201899cda47SBarry Smith       for (i=0; i<A->rmap.n; i++) {
120233a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1203289bc588SBarry Smith       }
1204064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1205289bc588SBarry Smith     }
1206899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12073a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1208064f8208SBarry Smith     *nrm = 0.0;
1209899cda47SBarry Smith     for (j=0; j<A->rmap.n; j++) {
1210289bc588SBarry Smith       v = mat->v + j;
1211289bc588SBarry Smith       sum = 0.0;
1212899cda47SBarry Smith       for (i=0; i<A->cmap.n; i++) {
12131b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1214289bc588SBarry Smith       }
1215064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1216289bc588SBarry Smith     }
1217899cda47SBarry Smith     ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
12183a40ed3dSBarry Smith   } else {
121929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1220289bc588SBarry Smith   }
12213a40ed3dSBarry Smith   PetscFunctionReturn(0);
1222289bc588SBarry Smith }
1223289bc588SBarry Smith 
12244a2ae208SSatish Balay #undef __FUNCT__
12254a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1226dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op)
1227289bc588SBarry Smith {
1228c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
122963ba0a88SBarry Smith   PetscErrorCode ierr;
123067e560aaSBarry Smith 
12313a40ed3dSBarry Smith   PetscFunctionBegin;
1232b5a2b587SKris Buschelman   switch (op) {
1233b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1234b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1235b5a2b587SKris Buschelman     break;
1236b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1237b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1238b5a2b587SKris Buschelman     break;
1239b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1240b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1241b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1242b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1243b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1244b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1245b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1246b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1247b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1248b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1249b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
1250ad86a440SBarry Smith     ierr = PetscInfo1(A,"Option %d ignored\n",op);CHKERRQ(ierr);
1251b5a2b587SKris Buschelman     break;
125277e54ba9SKris Buschelman   case MAT_SYMMETRIC:
125377e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12549a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12559a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12569a4540c5SBarry Smith   case MAT_HERMITIAN:
12579a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12589a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12599a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
126077e54ba9SKris Buschelman     break;
1261b5a2b587SKris Buschelman   default:
1262ad86a440SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
12633a40ed3dSBarry Smith   }
12643a40ed3dSBarry Smith   PetscFunctionReturn(0);
1265289bc588SBarry Smith }
1266289bc588SBarry Smith 
12674a2ae208SSatish Balay #undef __FUNCT__
12684a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1269dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
12706f0a148fSBarry Smith {
1271ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12726849ba73SBarry Smith   PetscErrorCode ierr;
1273899cda47SBarry Smith   PetscInt       lda=l->lda,m=A->rmap.n,j;
12743a40ed3dSBarry Smith 
12753a40ed3dSBarry Smith   PetscFunctionBegin;
1276a5ce6ee0Svictorle   if (lda>m) {
1277899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
1278a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1279a5ce6ee0Svictorle     }
1280a5ce6ee0Svictorle   } else {
1281899cda47SBarry Smith     ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1282a5ce6ee0Svictorle   }
12833a40ed3dSBarry Smith   PetscFunctionReturn(0);
12846f0a148fSBarry Smith }
12856f0a148fSBarry Smith 
12864a2ae208SSatish Balay #undef __FUNCT__
12874a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1288f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
12896f0a148fSBarry Smith {
1290ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1291899cda47SBarry Smith   PetscInt       n = A->cmap.n,i,j;
129287828ca2SBarry Smith   PetscScalar    *slot;
129355659b69SBarry Smith 
12943a40ed3dSBarry Smith   PetscFunctionBegin;
12956f0a148fSBarry Smith   for (i=0; i<N; i++) {
12966f0a148fSBarry Smith     slot = l->v + rows[i];
12976f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
12986f0a148fSBarry Smith   }
1299f4df32b1SMatthew Knepley   if (diag != 0.0) {
13006f0a148fSBarry Smith     for (i=0; i<N; i++) {
13016f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1302f4df32b1SMatthew Knepley       *slot = diag;
13036f0a148fSBarry Smith     }
13046f0a148fSBarry Smith   }
13053a40ed3dSBarry Smith   PetscFunctionReturn(0);
13066f0a148fSBarry Smith }
1307557bce09SLois Curfman McInnes 
13084a2ae208SSatish Balay #undef __FUNCT__
13094a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1310dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
131164e87e97SBarry Smith {
1312c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13133a40ed3dSBarry Smith 
13143a40ed3dSBarry Smith   PetscFunctionBegin;
1315899cda47SBarry Smith   if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
131664e87e97SBarry Smith   *array = mat->v;
13173a40ed3dSBarry Smith   PetscFunctionReturn(0);
131864e87e97SBarry Smith }
13190754003eSLois Curfman McInnes 
13204a2ae208SSatish Balay #undef __FUNCT__
13214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1322dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1323ff14e315SSatish Balay {
13243a40ed3dSBarry Smith   PetscFunctionBegin;
132509b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13263a40ed3dSBarry Smith   PetscFunctionReturn(0);
1327ff14e315SSatish Balay }
13280754003eSLois Curfman McInnes 
13294a2ae208SSatish Balay #undef __FUNCT__
13304a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
133113f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13320754003eSLois Curfman McInnes {
1333c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13346849ba73SBarry Smith   PetscErrorCode ierr;
133521a2c019SBarry Smith   PetscInt       i,j,*irow,*icol,nrows,ncols;
133687828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13370754003eSLois Curfman McInnes   Mat            newmat;
13380754003eSLois Curfman McInnes 
13393a40ed3dSBarry Smith   PetscFunctionBegin;
134078b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
134178b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1342e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1343e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13440754003eSLois Curfman McInnes 
1345182d2002SSatish Balay   /* Check submatrixcall */
1346182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
134713f74950SBarry Smith     PetscInt n_cols,n_rows;
1348182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
134921a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
135021a2c019SBarry Smith       /* resize the result result matrix to match number of requested rows/columns */
135121a2c019SBarry Smith       ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr);
135221a2c019SBarry Smith     }
1353182d2002SSatish Balay     newmat = *B;
1354182d2002SSatish Balay   } else {
13550754003eSLois Curfman McInnes     /* Create and fill new matrix */
1356f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
1357f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
13585c5985e7SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
13595c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1360182d2002SSatish Balay   }
1361182d2002SSatish Balay 
1362182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1363182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1364182d2002SSatish Balay 
1365182d2002SSatish Balay   for (i=0; i<ncols; i++) {
13666de62eeeSBarry Smith     av = v + mat->lda*icol[i];
1367182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1368182d2002SSatish Balay       *bv++ = av[irow[j]];
13690754003eSLois Curfman McInnes     }
13700754003eSLois Curfman McInnes   }
1371182d2002SSatish Balay 
1372182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13736d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13746d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13750754003eSLois Curfman McInnes 
13760754003eSLois Curfman McInnes   /* Free work space */
137778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
137878b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1379182d2002SSatish Balay   *B = newmat;
13803a40ed3dSBarry Smith   PetscFunctionReturn(0);
13810754003eSLois Curfman McInnes }
13820754003eSLois Curfman McInnes 
13834a2ae208SSatish Balay #undef __FUNCT__
13844a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
138513f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1386905e6a2fSBarry Smith {
13876849ba73SBarry Smith   PetscErrorCode ierr;
138813f74950SBarry Smith   PetscInt       i;
1389905e6a2fSBarry Smith 
13903a40ed3dSBarry Smith   PetscFunctionBegin;
1391905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1392b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1393905e6a2fSBarry Smith   }
1394905e6a2fSBarry Smith 
1395905e6a2fSBarry Smith   for (i=0; i<n; i++) {
13966a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1397905e6a2fSBarry Smith   }
13983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1399905e6a2fSBarry Smith }
1400905e6a2fSBarry Smith 
14014a2ae208SSatish Balay #undef __FUNCT__
1402c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense"
1403c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1404c0aa2d19SHong Zhang {
1405c0aa2d19SHong Zhang   PetscFunctionBegin;
1406c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1407c0aa2d19SHong Zhang }
1408c0aa2d19SHong Zhang 
1409c0aa2d19SHong Zhang #undef __FUNCT__
1410c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense"
1411c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1412c0aa2d19SHong Zhang {
1413c0aa2d19SHong Zhang   PetscFunctionBegin;
1414c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1415c0aa2d19SHong Zhang }
1416c0aa2d19SHong Zhang 
1417c0aa2d19SHong Zhang #undef __FUNCT__
14184a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1419dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14204b0e389bSBarry Smith {
14214b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14226849ba73SBarry Smith   PetscErrorCode ierr;
1423899cda47SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j;
14243a40ed3dSBarry Smith 
14253a40ed3dSBarry Smith   PetscFunctionBegin;
142633f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
142733f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1428cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14293a40ed3dSBarry Smith     PetscFunctionReturn(0);
14303a40ed3dSBarry Smith   }
1431899cda47SBarry Smith   if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1432a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14330dbb7854Svictorle     for (j=0; j<n; j++) {
1434a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1435a5ce6ee0Svictorle     }
1436a5ce6ee0Svictorle   } else {
1437899cda47SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr);
1438a5ce6ee0Svictorle   }
1439273d9f13SBarry Smith   PetscFunctionReturn(0);
1440273d9f13SBarry Smith }
1441273d9f13SBarry Smith 
14424a2ae208SSatish Balay #undef __FUNCT__
14434a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1444dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1445273d9f13SBarry Smith {
1446dfbe8321SBarry Smith   PetscErrorCode ierr;
1447273d9f13SBarry Smith 
1448273d9f13SBarry Smith   PetscFunctionBegin;
1449273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14503a40ed3dSBarry Smith   PetscFunctionReturn(0);
14514b0e389bSBarry Smith }
14524b0e389bSBarry Smith 
1453284134d9SBarry Smith #undef __FUNCT__
1454284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense"
1455284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
1456284134d9SBarry Smith {
1457284134d9SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
145821a2c019SBarry Smith   PetscErrorCode ierr;
1459284134d9SBarry Smith   PetscFunctionBegin;
146021a2c019SBarry Smith   /* this will not be called before lda, Mmax,  and Nmax have been set */
1461284134d9SBarry Smith   m = PetscMax(m,M);
1462284134d9SBarry Smith   n = PetscMax(n,N);
146321a2c019SBarry 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);
1464284134d9SBarry 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);
1465899cda47SBarry Smith   A->rmap.n = A->rmap.n = m;
1466899cda47SBarry Smith   A->cmap.n = A->cmap.N = n;
146721a2c019SBarry Smith   if (a->changelda) a->lda = m;
146821a2c019SBarry Smith   ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
1469284134d9SBarry Smith   PetscFunctionReturn(0);
1470284134d9SBarry Smith }
1471284134d9SBarry Smith 
1472a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
1473a9fe9ddaSSatish Balay #undef __FUNCT__
1474a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense"
1475a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1476a9fe9ddaSSatish Balay {
1477a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1478a9fe9ddaSSatish Balay 
1479a9fe9ddaSSatish Balay   PetscFunctionBegin;
1480a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1481a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1482a9fe9ddaSSatish Balay   }
1483a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1484a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1485a9fe9ddaSSatish Balay }
1486a9fe9ddaSSatish Balay 
1487a9fe9ddaSSatish Balay #undef __FUNCT__
1488a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense"
1489a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1490a9fe9ddaSSatish Balay {
1491ee16a9a1SHong Zhang   PetscErrorCode ierr;
1492899cda47SBarry Smith   PetscInt       m=A->rmap.n,n=B->cmap.n;
1493ee16a9a1SHong Zhang   Mat            Cmat;
1494a9fe9ddaSSatish Balay 
1495ee16a9a1SHong Zhang   PetscFunctionBegin;
1496899cda47SBarry 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);
1497ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1498ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1499ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1500ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1501ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1502ee16a9a1SHong Zhang   *C = Cmat;
1503ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1504ee16a9a1SHong Zhang }
1505a9fe9ddaSSatish Balay 
1506a9fe9ddaSSatish Balay #undef __FUNCT__
1507a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense"
1508a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1509a9fe9ddaSSatish Balay {
1510a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1511a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1512a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1513899cda47SBarry Smith   PetscBLASInt   m=(PetscBLASInt)A->rmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->cmap.n;
1514a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1515a9fe9ddaSSatish Balay 
1516a9fe9ddaSSatish Balay   PetscFunctionBegin;
1517a9fe9ddaSSatish Balay   BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1518a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1519a9fe9ddaSSatish Balay }
1520a9fe9ddaSSatish Balay 
1521a9fe9ddaSSatish Balay #undef __FUNCT__
1522a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense"
1523a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1524a9fe9ddaSSatish Balay {
1525a9fe9ddaSSatish Balay   PetscErrorCode ierr;
1526a9fe9ddaSSatish Balay 
1527a9fe9ddaSSatish Balay   PetscFunctionBegin;
1528a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX){
1529a9fe9ddaSSatish Balay     ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1530a9fe9ddaSSatish Balay   }
1531a9fe9ddaSSatish Balay   ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
1532a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1533a9fe9ddaSSatish Balay }
1534a9fe9ddaSSatish Balay 
1535a9fe9ddaSSatish Balay #undef __FUNCT__
1536a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense"
1537a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1538a9fe9ddaSSatish Balay {
1539ee16a9a1SHong Zhang   PetscErrorCode ierr;
1540899cda47SBarry Smith   PetscInt       m=A->cmap.n,n=B->cmap.n;
1541ee16a9a1SHong Zhang   Mat            Cmat;
1542a9fe9ddaSSatish Balay 
1543ee16a9a1SHong Zhang   PetscFunctionBegin;
1544899cda47SBarry 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);
1545ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
1546ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
1547ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
1548ee16a9a1SHong Zhang   ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
1549ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
1550ee16a9a1SHong Zhang   *C = Cmat;
1551ee16a9a1SHong Zhang   PetscFunctionReturn(0);
1552ee16a9a1SHong Zhang }
1553a9fe9ddaSSatish Balay 
1554a9fe9ddaSSatish Balay #undef __FUNCT__
1555a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense"
1556a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
1557a9fe9ddaSSatish Balay {
1558a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1559a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
1560a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
1561899cda47SBarry Smith   PetscBLASInt   m=(PetscBLASInt)A->cmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->rmap.n;
1562a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
1563a9fe9ddaSSatish Balay 
1564a9fe9ddaSSatish Balay   PetscFunctionBegin;
1565a9fe9ddaSSatish Balay   BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda);
1566a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
1567a9fe9ddaSSatish Balay }
1568289bc588SBarry Smith /* -------------------------------------------------------------------*/
1569a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1570905e6a2fSBarry Smith        MatGetRow_SeqDense,
1571905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1572905e6a2fSBarry Smith        MatMult_SeqDense,
157397304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
15747c922b88SBarry Smith        MatMultTranspose_SeqDense,
15757c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1576905e6a2fSBarry Smith        MatSolve_SeqDense,
1577905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
15787c922b88SBarry Smith        MatSolveTranspose_SeqDense,
157997304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1580905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1581905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1582ec8511deSBarry Smith        MatRelax_SeqDense,
1583ec8511deSBarry Smith        MatTranspose_SeqDense,
158497304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1585905e6a2fSBarry Smith        MatEqual_SeqDense,
1586905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1587905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1588905e6a2fSBarry Smith        MatNorm_SeqDense,
1589c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense,
1590c0aa2d19SHong Zhang        MatAssemblyEnd_SeqDense,
1591905e6a2fSBarry Smith        0,
1592905e6a2fSBarry Smith        MatSetOption_SeqDense,
1593905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
159497304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1595905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1596905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1597905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1598905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
159997304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1600273d9f13SBarry Smith        0,
1601905e6a2fSBarry Smith        0,
1602905e6a2fSBarry Smith        MatGetArray_SeqDense,
1603905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
160497304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1605a5ae1ecdSBarry Smith        0,
1606a5ae1ecdSBarry Smith        0,
1607a5ae1ecdSBarry Smith        0,
1608a5ae1ecdSBarry Smith        0,
160997304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1610a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1611a5ae1ecdSBarry Smith        0,
16124b0e389bSBarry Smith        MatGetValues_SeqDense,
1613a5ae1ecdSBarry Smith        MatCopy_SeqDense,
161497304618SKris Buschelman /*45*/ 0,
1615a5ae1ecdSBarry Smith        MatScale_SeqDense,
1616a5ae1ecdSBarry Smith        0,
1617a5ae1ecdSBarry Smith        0,
1618a5ae1ecdSBarry Smith        0,
1619521d7252SBarry Smith /*50*/ 0,
1620a5ae1ecdSBarry Smith        0,
1621a5ae1ecdSBarry Smith        0,
1622a5ae1ecdSBarry Smith        0,
1623a5ae1ecdSBarry Smith        0,
162497304618SKris Buschelman /*55*/ 0,
1625a5ae1ecdSBarry Smith        0,
1626a5ae1ecdSBarry Smith        0,
1627a5ae1ecdSBarry Smith        0,
1628a5ae1ecdSBarry Smith        0,
162997304618SKris Buschelman /*60*/ 0,
1630e03a110bSBarry Smith        MatDestroy_SeqDense,
1631e03a110bSBarry Smith        MatView_SeqDense,
1632357abbc8SBarry Smith        0,
163397304618SKris Buschelman        0,
163497304618SKris Buschelman /*65*/ 0,
163597304618SKris Buschelman        0,
163697304618SKris Buschelman        0,
163797304618SKris Buschelman        0,
163897304618SKris Buschelman        0,
163997304618SKris Buschelman /*70*/ 0,
164097304618SKris Buschelman        0,
164197304618SKris Buschelman        0,
164297304618SKris Buschelman        0,
164397304618SKris Buschelman        0,
164497304618SKris Buschelman /*75*/ 0,
164597304618SKris Buschelman        0,
164697304618SKris Buschelman        0,
164797304618SKris Buschelman        0,
164897304618SKris Buschelman        0,
164997304618SKris Buschelman /*80*/ 0,
165097304618SKris Buschelman        0,
165197304618SKris Buschelman        0,
165297304618SKris Buschelman        0,
1653865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1654865e5f61SKris Buschelman        0,
1655865e5f61SKris Buschelman        0,
1656865e5f61SKris Buschelman        0,
1657865e5f61SKris Buschelman        0,
1658865e5f61SKris Buschelman        0,
1659a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense,
1660a9fe9ddaSSatish Balay        MatMatMultSymbolic_SeqDense_SeqDense,
1661a9fe9ddaSSatish Balay        MatMatMultNumeric_SeqDense_SeqDense,
1662865e5f61SKris Buschelman        0,
1663865e5f61SKris Buschelman        0,
1664865e5f61SKris Buschelman /*95*/ 0,
1665a9fe9ddaSSatish Balay        MatMatMultTranspose_SeqDense_SeqDense,
1666a9fe9ddaSSatish Balay        MatMatMultTransposeSymbolic_SeqDense_SeqDense,
1667a9fe9ddaSSatish Balay        MatMatMultTransposeNumeric_SeqDense_SeqDense,
1668284134d9SBarry Smith        0,
1669284134d9SBarry Smith /*100*/0,
1670284134d9SBarry Smith        0,
1671284134d9SBarry Smith        0,
1672284134d9SBarry Smith        0,
1673284134d9SBarry Smith        MatSetSizes_SeqDense};
167490ace30eSBarry Smith 
16754a2ae208SSatish Balay #undef __FUNCT__
16764a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
16774b828684SBarry Smith /*@C
1678fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1679d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1680d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1681289bc588SBarry Smith 
1682db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1683db81eaa0SLois Curfman McInnes 
168420563c6bSBarry Smith    Input Parameters:
1685db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
16860c775827SLois Curfman McInnes .  m - number of rows
168718f449edSLois Curfman McInnes .  n - number of columns
1688db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1689dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
169020563c6bSBarry Smith 
169120563c6bSBarry Smith    Output Parameter:
169244cd7ae7SLois Curfman McInnes .  A - the matrix
169320563c6bSBarry Smith 
1694b259b22eSLois Curfman McInnes    Notes:
169518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
169618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1697b4fd4287SBarry Smith    set data=PETSC_NULL.
169818f449edSLois Curfman McInnes 
1699027ccd11SLois Curfman McInnes    Level: intermediate
1700027ccd11SLois Curfman McInnes 
1701dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1702d65003e9SLois Curfman McInnes 
1703db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
170420563c6bSBarry Smith @*/
1705be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1706289bc588SBarry Smith {
1707dfbe8321SBarry Smith   PetscErrorCode ierr;
17083b2fbd54SBarry Smith 
17093a40ed3dSBarry Smith   PetscFunctionBegin;
1710f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1711f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1712273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1713273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1714273d9f13SBarry Smith   PetscFunctionReturn(0);
1715273d9f13SBarry Smith }
1716273d9f13SBarry Smith 
17174a2ae208SSatish Balay #undef __FUNCT__
17184a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1719273d9f13SBarry Smith /*@C
1720273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1721273d9f13SBarry Smith 
1722273d9f13SBarry Smith    Collective on MPI_Comm
1723273d9f13SBarry Smith 
1724273d9f13SBarry Smith    Input Parameters:
1725273d9f13SBarry Smith +  A - the matrix
1726273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1727273d9f13SBarry Smith 
1728273d9f13SBarry Smith    Notes:
1729273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1730273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1731284134d9SBarry Smith    need not call this routine.
1732273d9f13SBarry Smith 
1733273d9f13SBarry Smith    Level: intermediate
1734273d9f13SBarry Smith 
1735273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1736273d9f13SBarry Smith 
1737273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1738273d9f13SBarry Smith @*/
1739be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1740273d9f13SBarry Smith {
1741dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1742a23d5eceSKris Buschelman 
1743a23d5eceSKris Buschelman   PetscFunctionBegin;
1744a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1745a23d5eceSKris Buschelman   if (f) {
1746a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1747a23d5eceSKris Buschelman   }
1748a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1749a23d5eceSKris Buschelman }
1750a23d5eceSKris Buschelman 
1751a23d5eceSKris Buschelman EXTERN_C_BEGIN
1752a23d5eceSKris Buschelman #undef __FUNCT__
1753a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1754be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1755a23d5eceSKris Buschelman {
1756273d9f13SBarry Smith   Mat_SeqDense   *b;
1757dfbe8321SBarry Smith   PetscErrorCode ierr;
1758273d9f13SBarry Smith 
1759273d9f13SBarry Smith   PetscFunctionBegin;
1760273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1761273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1762273d9f13SBarry Smith   if (!data) {
1763284134d9SBarry Smith     ierr          = PetscMalloc((b->lda*b->Nmax+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
1764284134d9SBarry Smith     ierr          = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1765273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
1766284134d9SBarry Smith     ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
1767273d9f13SBarry Smith   } else { /* user-allocated storage */
1768273d9f13SBarry Smith     b->v          = data;
1769273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1770273d9f13SBarry Smith   }
1771273d9f13SBarry Smith   PetscFunctionReturn(0);
1772273d9f13SBarry Smith }
1773a23d5eceSKris Buschelman EXTERN_C_END
1774273d9f13SBarry Smith 
17751b807ce4Svictorle #undef __FUNCT__
17761b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
17771b807ce4Svictorle /*@C
17781b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
17791b807ce4Svictorle 
17801b807ce4Svictorle   Input parameter:
17811b807ce4Svictorle + A - the matrix
17821b807ce4Svictorle - lda - the leading dimension
17831b807ce4Svictorle 
17841b807ce4Svictorle   Notes:
17851b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
17861b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
17871b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
17881b807ce4Svictorle 
17891b807ce4Svictorle   Level: intermediate
17901b807ce4Svictorle 
17911b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
17921b807ce4Svictorle 
1793284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
17941b807ce4Svictorle @*/
1795be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
17961b807ce4Svictorle {
17971b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
179821a2c019SBarry Smith 
17991b807ce4Svictorle   PetscFunctionBegin;
1800899cda47SBarry Smith   if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n);
18011b807ce4Svictorle   b->lda       = lda;
180221a2c019SBarry Smith   b->changelda = PETSC_FALSE;
180321a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
18041b807ce4Svictorle   PetscFunctionReturn(0);
18051b807ce4Svictorle }
18061b807ce4Svictorle 
18070bad9183SKris Buschelman /*MC
1808fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
18090bad9183SKris Buschelman 
18100bad9183SKris Buschelman    Options Database Keys:
18110bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
18120bad9183SKris Buschelman 
18130bad9183SKris Buschelman   Level: beginner
18140bad9183SKris Buschelman 
181589665df3SBarry Smith .seealso: MatCreateSeqDense()
181689665df3SBarry Smith 
18170bad9183SKris Buschelman M*/
18180bad9183SKris Buschelman 
1819273d9f13SBarry Smith EXTERN_C_BEGIN
18204a2ae208SSatish Balay #undef __FUNCT__
18214a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1822be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
1823273d9f13SBarry Smith {
1824273d9f13SBarry Smith   Mat_SeqDense   *b;
1825dfbe8321SBarry Smith   PetscErrorCode ierr;
18267c334f02SBarry Smith   PetscMPIInt    size;
1827273d9f13SBarry Smith 
1828273d9f13SBarry Smith   PetscFunctionBegin;
1829273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
183029bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
183155659b69SBarry Smith 
1832899cda47SBarry Smith   B->rmap.bs = B->cmap.bs = 1;
1833899cda47SBarry Smith   ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr);
1834899cda47SBarry Smith   ierr = PetscMapInitialize(B->comm,&B->cmap);CHKERRQ(ierr);
1835273d9f13SBarry Smith 
1836b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1837549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
183844cd7ae7SLois Curfman McInnes   B->factor       = 0;
183990f02eecSBarry Smith   B->mapping      = 0;
184052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr);
184144cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
184218f449edSLois Curfman McInnes 
1843a5ae1ecdSBarry Smith 
184444cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1845273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1846273d9f13SBarry Smith   b->v            = 0;
1847899cda47SBarry Smith   b->lda          = B->rmap.n;
184821a2c019SBarry Smith   b->changelda    = PETSC_FALSE;
1849899cda47SBarry Smith   b->Mmax         = B->rmap.n;
1850899cda47SBarry Smith   b->Nmax         = B->cmap.n;
18514e220ebcSLois Curfman McInnes 
1852a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1853a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1854a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
18554ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C",
18564ae313f4SHong Zhang                                      "MatMatMult_SeqAIJ_SeqDense",
18574ae313f4SHong Zhang                                       MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
18584ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",
18594ae313f4SHong Zhang                                      "MatMatMultSymbolic_SeqAIJ_SeqDense",
18604ae313f4SHong Zhang                                       MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
18614ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",
18624ae313f4SHong Zhang                                      "MatMatMultNumeric_SeqAIJ_SeqDense",
18634ae313f4SHong Zhang                                       MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
186417667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
18653a40ed3dSBarry Smith   PetscFunctionReturn(0);
1866289bc588SBarry Smith }
1867c0aa2d19SHong Zhang 
1868c0aa2d19SHong Zhang 
1869273d9f13SBarry Smith EXTERN_C_END
1870