xref: /petsc/src/mat/impls/dense/seq/dense.c (revision f69a0ea3504314d029ee423e0de2950ece2dab86)
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"
12dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqDense(const PetscScalar *alpha,Mat X,Mat Y,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
1513f74950SBarry Smith   PetscInt     j;
164ce68768SBarry Smith   PetscBLASInt N = (PetscBLASInt)X->m*X->n,m=(PetscBLASInt)X->m,ldax = x->lda,lday=y->lda,one = 1;
17efee365bSSatish Balay   PetscErrorCode ierr;
183a40ed3dSBarry Smith 
193a40ed3dSBarry Smith   PetscFunctionBegin;
20a5ce6ee0Svictorle   if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
21a5ce6ee0Svictorle   if (ldax>m || lday>m) {
22a5ce6ee0Svictorle     for (j=0; j<X->n; j++) {
2371044d3cSBarry Smith       BLASaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one);
24a5ce6ee0Svictorle     }
25a5ce6ee0Svictorle   } else {
2671044d3cSBarry Smith     BLASaxpy_(&N,(PetscScalar*)alpha,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 {
364e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3713f74950SBarry Smith   PetscInt     i,N = A->m*A->n,count = 0;
3887828ca2SBarry Smith   PetscScalar  *v = a->v;
393a40ed3dSBarry Smith 
403a40ed3dSBarry Smith   PetscFunctionBegin;
41289bc588SBarry Smith   for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}
424e220ebcSLois Curfman McInnes 
43273d9f13SBarry Smith   info->rows_global       = (double)A->m;
44273d9f13SBarry Smith   info->columns_global    = (double)A->n;
45273d9f13SBarry Smith   info->rows_local        = (double)A->m;
46273d9f13SBarry Smith   info->columns_local     = (double)A->n;
474e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
484e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
494e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
504e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
514e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
524e220ebcSLois Curfman McInnes   info->mallocs           = 0;
534e220ebcSLois Curfman McInnes   info->memory            = A->mem;
544e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
554e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
564e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
574e220ebcSLois Curfman McInnes 
583a40ed3dSBarry Smith   PetscFunctionReturn(0);
59289bc588SBarry Smith }
60289bc588SBarry Smith 
614a2ae208SSatish Balay #undef __FUNCT__
624a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
63dfbe8321SBarry Smith PetscErrorCode MatScale_SeqDense(const PetscScalar *alpha,Mat A)
6480cd9d93SLois Curfman McInnes {
65273d9f13SBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
664ce68768SBarry Smith   PetscBLASInt one = 1,lda = a->lda,j,nz;
67efee365bSSatish Balay   PetscErrorCode ierr;
6880cd9d93SLois Curfman McInnes 
693a40ed3dSBarry Smith   PetscFunctionBegin;
70a5ce6ee0Svictorle   if (lda>A->m) {
714ce68768SBarry Smith     nz = (PetscBLASInt)A->m;
72a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
7371044d3cSBarry Smith       BLASscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one);
74a5ce6ee0Svictorle     }
75a5ce6ee0Svictorle   } else {
764ce68768SBarry Smith     nz = (PetscBLASInt)A->m*A->n;
7771044d3cSBarry Smith     BLASscal_(&nz,(PetscScalar*)alpha,a->v,&one);
78a5ce6ee0Svictorle   }
79efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
803a40ed3dSBarry Smith   PetscFunctionReturn(0);
8180cd9d93SLois Curfman McInnes }
8280cd9d93SLois Curfman McInnes 
83289bc588SBarry Smith /* ---------------------------------------------------------------*/
846831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
856831982aSBarry Smith    rather than put it in the Mat->row slot.*/
864a2ae208SSatish Balay #undef __FUNCT__
874a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
88dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
89289bc588SBarry Smith {
90a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
9181f479a6Svictorle   PetscFunctionBegin;
92a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
93a07ab388SBarry Smith #else
94c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
956849ba73SBarry Smith   PetscErrorCode ierr;
964ce68768SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->n,m = (PetscBLASInt)A->m,info;
9755659b69SBarry Smith 
983a40ed3dSBarry Smith   PetscFunctionBegin;
99289bc588SBarry Smith   if (!mat->pivots) {
1004ce68768SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
10152e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,A->m*sizeof(PetscBLASInt));CHKERRQ(ierr);
102289bc588SBarry Smith   }
103f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
104273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
10571044d3cSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
10629bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
10729bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
108efee365bSSatish Balay   ierr = PetscLogFlops((2*A->n*A->n*A->n)/3);CHKERRQ(ierr);
109a07ab388SBarry Smith #endif
1103a40ed3dSBarry Smith   PetscFunctionReturn(0);
111289bc588SBarry Smith }
1126ee01492SSatish Balay 
1134a2ae208SSatish Balay #undef __FUNCT__
1144a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
115dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11602cad45dSBarry Smith {
11702cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
1186849ba73SBarry Smith   PetscErrorCode ierr;
11913f74950SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
12002cad45dSBarry Smith   Mat            newi;
12102cad45dSBarry Smith 
1223a40ed3dSBarry Smith   PetscFunctionBegin;
123*f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&newi);CHKERRQ(ierr);
124*f69a0ea3SMatthew Knepley   ierr = MatSetSizes(newi,A->m,A->n,A->m,A->n);CHKERRQ(ierr);
1255c5985e7SKris Buschelman   ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr);
1265c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
1275609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
128a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
129a5ce6ee0Svictorle     if (lda>A->m) {
130a5ce6ee0Svictorle       m = A->m;
131a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
132a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
133a5ce6ee0Svictorle       }
134a5ce6ee0Svictorle     } else {
13587828ca2SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
13602cad45dSBarry Smith     }
137a5ce6ee0Svictorle   }
13802cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
13902cad45dSBarry Smith   *newmat = newi;
1403a40ed3dSBarry Smith   PetscFunctionReturn(0);
14102cad45dSBarry Smith }
14202cad45dSBarry Smith 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
145dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
146289bc588SBarry Smith {
147dfbe8321SBarry Smith   PetscErrorCode ierr;
1483a40ed3dSBarry Smith 
1493a40ed3dSBarry Smith   PetscFunctionBegin;
1505609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1513a40ed3dSBarry Smith   PetscFunctionReturn(0);
152289bc588SBarry Smith }
1536ee01492SSatish Balay 
1544a2ae208SSatish Balay #undef __FUNCT__
1554a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
156af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
157289bc588SBarry Smith {
15802cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1596849ba73SBarry Smith   PetscErrorCode ierr;
16013f74950SBarry Smith   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j;
1614482741eSBarry Smith   MatFactorInfo  info;
1623a40ed3dSBarry Smith 
1633a40ed3dSBarry Smith   PetscFunctionBegin;
16402cad45dSBarry Smith   /* copy the numerical values */
1651b807ce4Svictorle   if (lda1>m || lda2>m ) {
1661b807ce4Svictorle     for (j=0; j<n; j++) {
1671b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1681b807ce4Svictorle     }
1691b807ce4Svictorle   } else {
17087828ca2SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1711b807ce4Svictorle   }
17202cad45dSBarry Smith   (*fact)->factor = 0;
1734482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175289bc588SBarry Smith }
1766ee01492SSatish Balay 
1774a2ae208SSatish Balay #undef __FUNCT__
1784a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
179dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
180289bc588SBarry Smith {
181dfbe8321SBarry Smith   PetscErrorCode ierr;
1823a40ed3dSBarry Smith 
1833a40ed3dSBarry Smith   PetscFunctionBegin;
184ceb03754SKris Buschelman   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr);
1853a40ed3dSBarry Smith   PetscFunctionReturn(0);
186289bc588SBarry Smith }
1876ee01492SSatish Balay 
1884a2ae208SSatish Balay #undef __FUNCT__
1894a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
190dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
191289bc588SBarry Smith {
192a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
193a07ab388SBarry Smith   PetscFunctionBegin;
194a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
195a07ab388SBarry Smith #else
196c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1976849ba73SBarry Smith   PetscErrorCode ierr;
1984ce68768SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->n,info;
19955659b69SBarry Smith 
2003a40ed3dSBarry Smith   PetscFunctionBegin;
201752f5567SLois Curfman McInnes   if (mat->pivots) {
202606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
20352e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-A->m*sizeof(PetscInt));CHKERRQ(ierr);
204752f5567SLois Curfman McInnes     mat->pivots = 0;
205752f5567SLois Curfman McInnes   }
206273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
20771044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
20877431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
209c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
210efee365bSSatish Balay   ierr = PetscLogFlops((A->n*A->n*A->n)/3);CHKERRQ(ierr);
211a07ab388SBarry Smith #endif
2123a40ed3dSBarry Smith   PetscFunctionReturn(0);
213289bc588SBarry Smith }
214289bc588SBarry Smith 
2154a2ae208SSatish Balay #undef __FUNCT__
2160b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
217af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2180b4b3355SBarry Smith {
219dfbe8321SBarry Smith   PetscErrorCode ierr;
22015e8a5b3SHong Zhang   MatFactorInfo  info;
2210b4b3355SBarry Smith 
2220b4b3355SBarry Smith   PetscFunctionBegin;
22315e8a5b3SHong Zhang   info.fill = 1.0;
22415e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2250b4b3355SBarry Smith   PetscFunctionReturn(0);
2260b4b3355SBarry Smith }
2270b4b3355SBarry Smith 
2280b4b3355SBarry Smith #undef __FUNCT__
2294a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
230dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
231289bc588SBarry Smith {
232c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2336849ba73SBarry Smith   PetscErrorCode ierr;
2344ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, one = 1,info;
23587828ca2SBarry Smith   PetscScalar    *x,*y;
23667e560aaSBarry Smith 
2373a40ed3dSBarry Smith   PetscFunctionBegin;
238273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2391ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2401ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
24187828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
242c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
243ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
24480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
245ae7cfcebSSatish Balay #else
24671044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2474ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
248ae7cfcebSSatish Balay #endif
2497a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
250ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
25180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
252ae7cfcebSSatish Balay #else
25371044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2544ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
255ae7cfcebSSatish Balay #endif
256289bc588SBarry Smith   }
25729bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2591ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
260efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr);
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
262289bc588SBarry Smith }
2636ee01492SSatish Balay 
2644a2ae208SSatish Balay #undef __FUNCT__
2654a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
266dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
267da3a660dSBarry Smith {
268c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
269dfbe8321SBarry Smith   PetscErrorCode ierr;
2704ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt) A->m,one = 1,info;
27187828ca2SBarry Smith   PetscScalar    *x,*y;
27267e560aaSBarry Smith 
2733a40ed3dSBarry Smith   PetscFunctionBegin;
274273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2751ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2761ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
27787828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
278752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
279da3a660dSBarry Smith   if (mat->pivots) {
280ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
28180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
282ae7cfcebSSatish Balay #else
28371044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2844ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
285ae7cfcebSSatish Balay #endif
2867a97a34bSBarry Smith   } else {
287ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
28880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
289ae7cfcebSSatish Balay #else
29071044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2914ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
292ae7cfcebSSatish Balay #endif
293da3a660dSBarry Smith   }
2941ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2951ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
296efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr);
2973a40ed3dSBarry Smith   PetscFunctionReturn(0);
298da3a660dSBarry Smith }
2996ee01492SSatish Balay 
3004a2ae208SSatish Balay #undef __FUNCT__
3014a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
302dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
303da3a660dSBarry Smith {
304c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
305dfbe8321SBarry Smith   PetscErrorCode ierr;
3064ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
30787828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
308da3a660dSBarry Smith   Vec            tmp = 0;
30967e560aaSBarry Smith 
3103a40ed3dSBarry Smith   PetscFunctionBegin;
3111ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3121ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
313273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
314da3a660dSBarry Smith   if (yy == zz) {
31578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
31652e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
31778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
318da3a660dSBarry Smith   }
31987828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
320752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
321da3a660dSBarry Smith   if (mat->pivots) {
322ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
32380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
324ae7cfcebSSatish Balay #else
32571044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3262ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
327ae7cfcebSSatish Balay #endif
328a8c6a408SBarry Smith   } else {
329ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
33080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
331ae7cfcebSSatish Balay #else
33271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3332ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
334ae7cfcebSSatish Balay #endif
335da3a660dSBarry Smith   }
3362dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
3372dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
3381ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3391ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
340efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr);
3413a40ed3dSBarry Smith   PetscFunctionReturn(0);
342da3a660dSBarry Smith }
34367e560aaSBarry Smith 
3444a2ae208SSatish Balay #undef __FUNCT__
3454a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
346dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
347da3a660dSBarry Smith {
348c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3496849ba73SBarry Smith   PetscErrorCode ierr;
3504ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
35187828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
352da3a660dSBarry Smith   Vec            tmp;
35367e560aaSBarry Smith 
3543a40ed3dSBarry Smith   PetscFunctionBegin;
355273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3561ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3571ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
358da3a660dSBarry Smith   if (yy == zz) {
35978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
36052e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
36178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
362da3a660dSBarry Smith   }
36387828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
364752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
365da3a660dSBarry Smith   if (mat->pivots) {
366ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
36780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
368ae7cfcebSSatish Balay #else
36971044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3702ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
371ae7cfcebSSatish Balay #endif
3723a40ed3dSBarry Smith   } else {
373ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
37480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
375ae7cfcebSSatish Balay #else
37671044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3772ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
378ae7cfcebSSatish Balay #endif
379da3a660dSBarry Smith   }
38090f02eecSBarry Smith   if (tmp) {
3812dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
38290f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3833a40ed3dSBarry Smith   } else {
3842dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
38590f02eecSBarry Smith   }
3861ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3871ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
388efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr);
3893a40ed3dSBarry Smith   PetscFunctionReturn(0);
390da3a660dSBarry Smith }
391289bc588SBarry Smith /* ------------------------------------------------------------------*/
3924a2ae208SSatish Balay #undef __FUNCT__
3934a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
39413f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
395289bc588SBarry Smith {
396c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
39787828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
398dfbe8321SBarry Smith   PetscErrorCode ierr;
39913f74950SBarry Smith   PetscInt       m = A->m,i;
400aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4014ce68768SBarry Smith   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
402bc1b551cSSatish Balay #endif
403289bc588SBarry Smith 
4043a40ed3dSBarry Smith   PetscFunctionBegin;
405289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
40671044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4072dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
408289bc588SBarry Smith   }
4091ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4101ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
411b965ef7fSBarry Smith   its  = its*lits;
41277431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
413289bc588SBarry Smith   while (its--) {
414fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
415289bc588SBarry Smith       for (i=0; i<m; i++) {
416aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
417f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
418f1747703SBarry Smith            not happy about returning a double complex */
41913f74950SBarry Smith         PetscInt         _i;
42087828ca2SBarry Smith         PetscScalar sum = b[i];
421f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4223f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
423f1747703SBarry Smith         }
424f1747703SBarry Smith         xt = sum;
425f1747703SBarry Smith #else
42671044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
427f1747703SBarry Smith #endif
42855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
429289bc588SBarry Smith       }
430289bc588SBarry Smith     }
431fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
432289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
433aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
434f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
435f1747703SBarry Smith            not happy about returning a double complex */
43613f74950SBarry Smith         PetscInt         _i;
43787828ca2SBarry Smith         PetscScalar sum = b[i];
438f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4393f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
440f1747703SBarry Smith         }
441f1747703SBarry Smith         xt = sum;
442f1747703SBarry Smith #else
44371044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
444f1747703SBarry Smith #endif
44555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
446289bc588SBarry Smith       }
447289bc588SBarry Smith     }
448289bc588SBarry Smith   }
4491ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4501ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4513a40ed3dSBarry Smith   PetscFunctionReturn(0);
452289bc588SBarry Smith }
453289bc588SBarry Smith 
454289bc588SBarry Smith /* -----------------------------------------------------------------*/
4554a2ae208SSatish Balay #undef __FUNCT__
4564a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
457dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
458289bc588SBarry Smith {
459c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
46087828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
461dfbe8321SBarry Smith   PetscErrorCode ierr;
4624ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n,_One=1;
463ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
4643a40ed3dSBarry Smith 
4653a40ed3dSBarry Smith   PetscFunctionBegin;
466273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4671ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4681ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
46971044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
4701ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4711ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
472efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n - A->n);CHKERRQ(ierr);
4733a40ed3dSBarry Smith   PetscFunctionReturn(0);
474289bc588SBarry Smith }
4756ee01492SSatish Balay 
4764a2ae208SSatish Balay #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
478dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
479289bc588SBarry Smith {
480c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
48187828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
482dfbe8321SBarry Smith   PetscErrorCode ierr;
4834ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
4843a40ed3dSBarry Smith 
4853a40ed3dSBarry Smith   PetscFunctionBegin;
486273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4871ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4881ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
48971044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
4901ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4911ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
492efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n - A->m);CHKERRQ(ierr);
4933a40ed3dSBarry Smith   PetscFunctionReturn(0);
494289bc588SBarry Smith }
4956ee01492SSatish Balay 
4964a2ae208SSatish Balay #undef __FUNCT__
4974a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
498dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
499289bc588SBarry Smith {
500c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
50187828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
502dfbe8321SBarry Smith   PetscErrorCode ierr;
5034ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
5043a40ed3dSBarry Smith 
5053a40ed3dSBarry Smith   PetscFunctionBegin;
506273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
507600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5081ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5091ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
51071044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5111ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5121ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
513efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr);
5143a40ed3dSBarry Smith   PetscFunctionReturn(0);
515289bc588SBarry Smith }
5166ee01492SSatish Balay 
5174a2ae208SSatish Balay #undef __FUNCT__
5184a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
519dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
520289bc588SBarry Smith {
521c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
52287828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
523dfbe8321SBarry Smith   PetscErrorCode ierr;
5244ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
52587828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5263a40ed3dSBarry Smith 
5273a40ed3dSBarry Smith   PetscFunctionBegin;
528273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
529600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5301ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5311ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
53271044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5331ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5341ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
535efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr);
5363a40ed3dSBarry Smith   PetscFunctionReturn(0);
537289bc588SBarry Smith }
538289bc588SBarry Smith 
539289bc588SBarry Smith /* -----------------------------------------------------------------*/
5404a2ae208SSatish Balay #undef __FUNCT__
5414a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
54213f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
543289bc588SBarry Smith {
544c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
54587828ca2SBarry Smith   PetscScalar    *v;
5466849ba73SBarry Smith   PetscErrorCode ierr;
54713f74950SBarry Smith   PetscInt       i;
54867e560aaSBarry Smith 
5493a40ed3dSBarry Smith   PetscFunctionBegin;
550273d9f13SBarry Smith   *ncols = A->n;
551289bc588SBarry Smith   if (cols) {
55213f74950SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
553273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
554289bc588SBarry Smith   }
555289bc588SBarry Smith   if (vals) {
55687828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
557289bc588SBarry Smith     v    = mat->v + row;
5581b807ce4Svictorle     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;}
559289bc588SBarry Smith   }
5603a40ed3dSBarry Smith   PetscFunctionReturn(0);
561289bc588SBarry Smith }
5626ee01492SSatish Balay 
5634a2ae208SSatish Balay #undef __FUNCT__
5644a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
56513f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
566289bc588SBarry Smith {
567dfbe8321SBarry Smith   PetscErrorCode ierr;
568606d414cSSatish Balay   PetscFunctionBegin;
569606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
570606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5713a40ed3dSBarry Smith   PetscFunctionReturn(0);
572289bc588SBarry Smith }
573289bc588SBarry Smith /* ----------------------------------------------------------------*/
5744a2ae208SSatish Balay #undef __FUNCT__
5754a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
57613f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
577289bc588SBarry Smith {
578c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
57913f74950SBarry Smith   PetscInt     i,j,idx=0;
580d6dfbf8fSBarry Smith 
5813a40ed3dSBarry Smith   PetscFunctionBegin;
582289bc588SBarry Smith   if (!mat->roworiented) {
583dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
584289bc588SBarry Smith       for (j=0; j<n; j++) {
585cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5862515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
58777431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
58858804f6dSBarry Smith #endif
589289bc588SBarry Smith         for (i=0; i<m; i++) {
590cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5912515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
59277431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
59358804f6dSBarry Smith #endif
594cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
595289bc588SBarry Smith         }
596289bc588SBarry Smith       }
5973a40ed3dSBarry Smith     } else {
598289bc588SBarry Smith       for (j=0; j<n; j++) {
599cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6002515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60177431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
60258804f6dSBarry Smith #endif
603289bc588SBarry Smith         for (i=0; i<m; i++) {
604cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6052515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60677431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
60758804f6dSBarry Smith #endif
608cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
609289bc588SBarry Smith         }
610289bc588SBarry Smith       }
611289bc588SBarry Smith     }
6123a40ed3dSBarry Smith   } else {
613dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
614e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
615cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6162515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
61777431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
61858804f6dSBarry Smith #endif
619e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
620cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6212515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
62277431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
62358804f6dSBarry Smith #endif
624cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
625e8d4e0b9SBarry Smith         }
626e8d4e0b9SBarry Smith       }
6273a40ed3dSBarry Smith     } else {
628289bc588SBarry Smith       for (i=0; i<m; i++) {
629cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6302515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
63177431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
63258804f6dSBarry Smith #endif
633289bc588SBarry Smith         for (j=0; j<n; j++) {
634cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6352515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
63677431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
63758804f6dSBarry Smith #endif
638cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
639289bc588SBarry Smith         }
640289bc588SBarry Smith       }
641289bc588SBarry Smith     }
642e8d4e0b9SBarry Smith   }
6433a40ed3dSBarry Smith   PetscFunctionReturn(0);
644289bc588SBarry Smith }
645e8d4e0b9SBarry Smith 
6464a2ae208SSatish Balay #undef __FUNCT__
6474a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
64813f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
649ae80bb75SLois Curfman McInnes {
650ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
65113f74950SBarry Smith   PetscInt     i,j;
65287828ca2SBarry Smith   PetscScalar  *vpt = v;
653ae80bb75SLois Curfman McInnes 
6543a40ed3dSBarry Smith   PetscFunctionBegin;
655ae80bb75SLois Curfman McInnes   /* row-oriented output */
656ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
657ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6581b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
659ae80bb75SLois Curfman McInnes     }
660ae80bb75SLois Curfman McInnes   }
6613a40ed3dSBarry Smith   PetscFunctionReturn(0);
662ae80bb75SLois Curfman McInnes }
663ae80bb75SLois Curfman McInnes 
664289bc588SBarry Smith /* -----------------------------------------------------------------*/
665289bc588SBarry Smith 
666e090d566SSatish Balay #include "petscsys.h"
66788e32edaSLois Curfman McInnes 
6684a2ae208SSatish Balay #undef __FUNCT__
6694a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
670*f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A)
67188e32edaSLois Curfman McInnes {
67288e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
67388e32edaSLois Curfman McInnes   Mat            B;
6746849ba73SBarry Smith   PetscErrorCode ierr;
67513f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
67613f74950SBarry Smith   int            fd;
67713f74950SBarry Smith   PetscMPIInt    size;
67813f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
67987828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
68019bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
68188e32edaSLois Curfman McInnes 
6823a40ed3dSBarry Smith   PetscFunctionBegin;
683d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
68429bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
685b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6860752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
687552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
68888e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
68988e32edaSLois Curfman McInnes 
69090ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
691*f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
692*f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
693be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6945c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
69590ace30eSBarry Smith     B    = *A;
69690ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6974a41a4d3SSatish Balay     v    = a->v;
6984a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6994a41a4d3SSatish Balay        from row major to column major */
700b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
70190ace30eSBarry Smith     /* read in nonzero values */
7024a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7034a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7044a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7054a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7064a41a4d3SSatish Balay         *v++ =w[i*N+j];
7074a41a4d3SSatish Balay       }
7084a41a4d3SSatish Balay     }
709606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7106d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7116d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71290ace30eSBarry Smith   } else {
71388e32edaSLois Curfman McInnes     /* read row lengths */
71413f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7150752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
71688e32edaSLois Curfman McInnes 
71788e32edaSLois Curfman McInnes     /* create our matrix */
718*f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
719*f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
720be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7215c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
72288e32edaSLois Curfman McInnes     B = *A;
72388e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
72488e32edaSLois Curfman McInnes     v = a->v;
72588e32edaSLois Curfman McInnes 
72688e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
72713f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
728b0a32e0cSBarry Smith     cols = scols;
7290752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
73087828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
731b0a32e0cSBarry Smith     vals = svals;
7320752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
73388e32edaSLois Curfman McInnes 
73488e32edaSLois Curfman McInnes     /* insert into matrix */
73588e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
73688e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
73788e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
73888e32edaSLois Curfman McInnes     }
739606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
740606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
741606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
74288e32edaSLois Curfman McInnes 
7436d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7446d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74590ace30eSBarry Smith   }
7463a40ed3dSBarry Smith   PetscFunctionReturn(0);
74788e32edaSLois Curfman McInnes }
74888e32edaSLois Curfman McInnes 
749e090d566SSatish Balay #include "petscsys.h"
750932b0c3eSLois Curfman McInnes 
7514a2ae208SSatish Balay #undef __FUNCT__
7524a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
7536849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
754289bc588SBarry Smith {
755932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
756dfbe8321SBarry Smith   PetscErrorCode    ierr;
75713f74950SBarry Smith   PetscInt          i,j;
7582dcb1b2aSMatthew Knepley   const char        *name;
75987828ca2SBarry Smith   PetscScalar       *v;
760f3ef73ceSBarry Smith   PetscViewerFormat format;
761932b0c3eSLois Curfman McInnes 
7623a40ed3dSBarry Smith   PetscFunctionBegin;
763435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
764b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
765456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7663a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
767fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
768b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
769273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
77044cd7ae7SLois Curfman McInnes       v = a->v + i;
77177431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
772273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
773aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
774329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
77577431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
776329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
77777431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7786831982aSBarry Smith         }
77980cd9d93SLois Curfman McInnes #else
7806831982aSBarry Smith         if (*v) {
78177431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,*v);CHKERRQ(ierr);
7826831982aSBarry Smith         }
78380cd9d93SLois Curfman McInnes #endif
7841b807ce4Svictorle         v += a->lda;
78580cd9d93SLois Curfman McInnes       }
786b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
78780cd9d93SLois Curfman McInnes     }
788b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7893a40ed3dSBarry Smith   } else {
790b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
791aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
792ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
79347989497SBarry Smith     /* determine if matrix has all real values */
79447989497SBarry Smith     v = a->v;
795201f5f94SSatish Balay     for (i=0; i<A->m*A->n; i++) {
796ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
79747989497SBarry Smith     }
79847989497SBarry Smith #endif
799fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8003a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
80177431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->m,A->n);CHKERRQ(ierr);
80277431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->m,A->n);CHKERRQ(ierr);
803fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
804ffac6cdbSBarry Smith     }
805ffac6cdbSBarry Smith 
806273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
807932b0c3eSLois Curfman McInnes       v = a->v + i;
808273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
809aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
81047989497SBarry Smith         if (allreal) {
8111b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
81247989497SBarry Smith         } else {
8131b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
81447989497SBarry Smith         }
815289bc588SBarry Smith #else
8161b807ce4Svictorle         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
817289bc588SBarry Smith #endif
8181b807ce4Svictorle 	v += a->lda;
819289bc588SBarry Smith       }
820b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
821289bc588SBarry Smith     }
822fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
823b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
824ffac6cdbSBarry Smith     }
825b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
826da3a660dSBarry Smith   }
827b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8283a40ed3dSBarry Smith   PetscFunctionReturn(0);
829289bc588SBarry Smith }
830289bc588SBarry Smith 
8314a2ae208SSatish Balay #undef __FUNCT__
8324a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8336849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
834932b0c3eSLois Curfman McInnes {
835932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8366849ba73SBarry Smith   PetscErrorCode    ierr;
83713f74950SBarry Smith   int               fd;
83813f74950SBarry Smith   PetscInt          ict,j,n = A->n,m = A->m,i,*col_lens,nz = m*n;
83987828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
840f3ef73ceSBarry Smith   PetscViewerFormat format;
841932b0c3eSLois Curfman McInnes 
8423a40ed3dSBarry Smith   PetscFunctionBegin;
843b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
84490ace30eSBarry Smith 
845b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
846fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
84790ace30eSBarry Smith     /* store the matrix as a dense matrix */
84813f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
849552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
85090ace30eSBarry Smith     col_lens[1] = m;
85190ace30eSBarry Smith     col_lens[2] = n;
85290ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8536f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
854606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
85590ace30eSBarry Smith 
85690ace30eSBarry Smith     /* write out matrix, by rows */
85787828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
85890ace30eSBarry Smith     v    = a->v;
85990ace30eSBarry Smith     for (i=0; i<m; i++) {
86090ace30eSBarry Smith       for (j=0; j<n; j++) {
86190ace30eSBarry Smith         vals[i + j*m] = *v++;
86290ace30eSBarry Smith       }
86390ace30eSBarry Smith     }
8646f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
865606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
86690ace30eSBarry Smith   } else {
86713f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
868552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
869932b0c3eSLois Curfman McInnes     col_lens[1] = m;
870932b0c3eSLois Curfman McInnes     col_lens[2] = n;
871932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
872932b0c3eSLois Curfman McInnes 
873932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
874932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8756f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
876932b0c3eSLois Curfman McInnes 
877932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
878932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
879932b0c3eSLois Curfman McInnes     ict = 0;
880932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
881932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
882932b0c3eSLois Curfman McInnes     }
8836f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
884606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
885932b0c3eSLois Curfman McInnes 
886932b0c3eSLois Curfman McInnes     /* store nonzero values */
88787828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
888932b0c3eSLois Curfman McInnes     ict  = 0;
889932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
890932b0c3eSLois Curfman McInnes       v = a->v + i;
891932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
8921b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
893932b0c3eSLois Curfman McInnes       }
894932b0c3eSLois Curfman McInnes     }
8956f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
896606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
89790ace30eSBarry Smith   }
8983a40ed3dSBarry Smith   PetscFunctionReturn(0);
899932b0c3eSLois Curfman McInnes }
900932b0c3eSLois Curfman McInnes 
9014a2ae208SSatish Balay #undef __FUNCT__
9024a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
903dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
904f1af5d2fSBarry Smith {
905f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
906f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9076849ba73SBarry Smith   PetscErrorCode    ierr;
90813f74950SBarry Smith   PetscInt          m = A->m,n = A->n,color,i,j;
90987828ca2SBarry Smith   PetscScalar       *v = a->v;
910b0a32e0cSBarry Smith   PetscViewer       viewer;
911b0a32e0cSBarry Smith   PetscDraw         popup;
912329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
913f3ef73ceSBarry Smith   PetscViewerFormat format;
914f1af5d2fSBarry Smith 
915f1af5d2fSBarry Smith   PetscFunctionBegin;
916f1af5d2fSBarry Smith 
917f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
918b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
919b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
920f1af5d2fSBarry Smith 
921f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
922fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
923f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
924b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
925f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
926f1af5d2fSBarry Smith       x_l = j;
927f1af5d2fSBarry Smith       x_r = x_l + 1.0;
928f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
929f1af5d2fSBarry Smith         y_l = m - i - 1.0;
930f1af5d2fSBarry Smith         y_r = y_l + 1.0;
931f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
932329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
933b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
934329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
935b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
936f1af5d2fSBarry Smith         } else {
937f1af5d2fSBarry Smith           continue;
938f1af5d2fSBarry Smith         }
939f1af5d2fSBarry Smith #else
940f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
941b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
942f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
943b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
944f1af5d2fSBarry Smith         } else {
945f1af5d2fSBarry Smith           continue;
946f1af5d2fSBarry Smith         }
947f1af5d2fSBarry Smith #endif
948b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
949f1af5d2fSBarry Smith       }
950f1af5d2fSBarry Smith     }
951f1af5d2fSBarry Smith   } else {
952f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
953f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
954f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
955f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
956f1af5d2fSBarry Smith     }
957b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
958b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
959b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
960f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
961f1af5d2fSBarry Smith       x_l = j;
962f1af5d2fSBarry Smith       x_r = x_l + 1.0;
963f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
964f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
965f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
966b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
967b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
968f1af5d2fSBarry Smith       }
969f1af5d2fSBarry Smith     }
970f1af5d2fSBarry Smith   }
971f1af5d2fSBarry Smith   PetscFunctionReturn(0);
972f1af5d2fSBarry Smith }
973f1af5d2fSBarry Smith 
9744a2ae208SSatish Balay #undef __FUNCT__
9754a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
976dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
977f1af5d2fSBarry Smith {
978b0a32e0cSBarry Smith   PetscDraw      draw;
979f1af5d2fSBarry Smith   PetscTruth     isnull;
980329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
981dfbe8321SBarry Smith   PetscErrorCode ierr;
982f1af5d2fSBarry Smith 
983f1af5d2fSBarry Smith   PetscFunctionBegin;
984b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
985b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
986abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
987f1af5d2fSBarry Smith 
988f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
989273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
990f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
991b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
992b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
993f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
994f1af5d2fSBarry Smith   PetscFunctionReturn(0);
995f1af5d2fSBarry Smith }
996f1af5d2fSBarry Smith 
9974a2ae208SSatish Balay #undef __FUNCT__
9984a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
999dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1000932b0c3eSLois Curfman McInnes {
1001dfbe8321SBarry Smith   PetscErrorCode ierr;
100232077d6dSBarry Smith   PetscTruth     issocket,iascii,isbinary,isdraw;
1003932b0c3eSLois Curfman McInnes 
10043a40ed3dSBarry Smith   PetscFunctionBegin;
1005b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
100632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1007fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1008fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10090f5bd95cSBarry Smith 
1010c45a1595SBarry Smith   if (iascii) {
1011c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10124cf18111SSatish Balay #if defined(PETSC_USE_SOCKET_VIEWER)
1013c45a1595SBarry Smith   } else if (issocket) {
10144cf18111SSatish Balay     Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1015634064b4SBarry Smith     if (a->lda>A->m) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA");
1016b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
1017c45a1595SBarry Smith #endif
10180f5bd95cSBarry Smith   } else if (isbinary) {
10193a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1020f1af5d2fSBarry Smith   } else if (isdraw) {
1021f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10225cd90555SBarry Smith   } else {
1023958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1024932b0c3eSLois Curfman McInnes   }
10253a40ed3dSBarry Smith   PetscFunctionReturn(0);
1026932b0c3eSLois Curfman McInnes }
1027289bc588SBarry Smith 
10284a2ae208SSatish Balay #undef __FUNCT__
10294a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1030dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1031289bc588SBarry Smith {
1032ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1033dfbe8321SBarry Smith   PetscErrorCode ierr;
103490f02eecSBarry Smith 
10353a40ed3dSBarry Smith   PetscFunctionBegin;
1036aa482453SBarry Smith #if defined(PETSC_USE_LOG)
103777431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->m,mat->n);
1038a5a9c739SBarry Smith #endif
1039606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1040606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1041606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1042901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1044289bc588SBarry Smith }
1045289bc588SBarry Smith 
10464a2ae208SSatish Balay #undef __FUNCT__
10474a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1048dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1049289bc588SBarry Smith {
1050c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
10516849ba73SBarry Smith   PetscErrorCode ierr;
105213f74950SBarry Smith   PetscInt       k,j,m,n,M;
105387828ca2SBarry Smith   PetscScalar    *v,tmp;
105448b35521SBarry Smith 
10553a40ed3dSBarry Smith   PetscFunctionBegin;
10561b807ce4Svictorle   v = mat->v; m = A->m; M = mat->lda; n = A->n;
10577c922b88SBarry Smith   if (!matout) { /* in place transpose */
1058a5ce6ee0Svictorle     if (m != n) {
1059634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1060d64ed03dSBarry Smith     } else {
1061d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1062289bc588SBarry Smith         for (k=0; k<j; k++) {
10631b807ce4Svictorle           tmp = v[j + k*M];
10641b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10651b807ce4Svictorle           v[k + j*M] = tmp;
1066289bc588SBarry Smith         }
1067289bc588SBarry Smith       }
1068d64ed03dSBarry Smith     }
10693a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1070d3e5ee88SLois Curfman McInnes     Mat          tmat;
1071ec8511deSBarry Smith     Mat_SeqDense *tmatd;
107287828ca2SBarry Smith     PetscScalar  *v2;
1073ea709b57SSatish Balay 
1074*f69a0ea3SMatthew Knepley     ierr  = MatCreate(A->comm,&tmat);CHKERRQ(ierr);
1075*f69a0ea3SMatthew Knepley     ierr  = MatSetSizes(tmat,A->n,A->m,A->n,A->m);CHKERRQ(ierr);
10765c5985e7SKris Buschelman     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
10775c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1078ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10790de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1080d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10811b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1082d3e5ee88SLois Curfman McInnes     }
10836d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10846d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1085d3e5ee88SLois Curfman McInnes     *matout = tmat;
108648b35521SBarry Smith   }
10873a40ed3dSBarry Smith   PetscFunctionReturn(0);
1088289bc588SBarry Smith }
1089289bc588SBarry Smith 
10904a2ae208SSatish Balay #undef __FUNCT__
10914a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1092dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1093289bc588SBarry Smith {
1094c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1095c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
109613f74950SBarry Smith   PetscInt     i,j;
109787828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
10989ea5d5aeSSatish Balay 
10993a40ed3dSBarry Smith   PetscFunctionBegin;
1100273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1101273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11021b807ce4Svictorle   for (i=0; i<A1->m; i++) {
11031b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
11041b807ce4Svictorle     for (j=0; j<A1->n; j++) {
11053a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11061b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11071b807ce4Svictorle     }
1108289bc588SBarry Smith   }
110977c4ece6SBarry Smith   *flg = PETSC_TRUE;
11103a40ed3dSBarry Smith   PetscFunctionReturn(0);
1111289bc588SBarry Smith }
1112289bc588SBarry Smith 
11134a2ae208SSatish Balay #undef __FUNCT__
11144a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1115dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1116289bc588SBarry Smith {
1117c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1118dfbe8321SBarry Smith   PetscErrorCode ierr;
111913f74950SBarry Smith   PetscInt       i,n,len;
112087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
112144cd7ae7SLois Curfman McInnes 
11223a40ed3dSBarry Smith   PetscFunctionBegin;
11232dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11247a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11251ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1126273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1127273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
112844cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11291b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1130289bc588SBarry Smith   }
11311ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11323a40ed3dSBarry Smith   PetscFunctionReturn(0);
1133289bc588SBarry Smith }
1134289bc588SBarry Smith 
11354a2ae208SSatish Balay #undef __FUNCT__
11364a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1137dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1138289bc588SBarry Smith {
1139c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
114087828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1141dfbe8321SBarry Smith   PetscErrorCode ierr;
114213f74950SBarry Smith   PetscInt       i,j,m = A->m,n = A->n;
114355659b69SBarry Smith 
11443a40ed3dSBarry Smith   PetscFunctionBegin;
114528988994SBarry Smith   if (ll) {
11467a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11471ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1148273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1149da3a660dSBarry Smith     for (i=0; i<m; i++) {
1150da3a660dSBarry Smith       x = l[i];
1151da3a660dSBarry Smith       v = mat->v + i;
1152da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1153da3a660dSBarry Smith     }
11541ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1155efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1156da3a660dSBarry Smith   }
115728988994SBarry Smith   if (rr) {
11587a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11591ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1160273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1161da3a660dSBarry Smith     for (i=0; i<n; i++) {
1162da3a660dSBarry Smith       x = r[i];
1163da3a660dSBarry Smith       v = mat->v + i*m;
1164da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1165da3a660dSBarry Smith     }
11661ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1167efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1168da3a660dSBarry Smith   }
11693a40ed3dSBarry Smith   PetscFunctionReturn(0);
1170289bc588SBarry Smith }
1171289bc588SBarry Smith 
11724a2ae208SSatish Balay #undef __FUNCT__
11734a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1174dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1175289bc588SBarry Smith {
1176c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
117787828ca2SBarry Smith   PetscScalar  *v = mat->v;
1178329f5518SBarry Smith   PetscReal    sum = 0.0;
117913f74950SBarry Smith   PetscInt     lda=mat->lda,m=A->m,i,j;
1180efee365bSSatish Balay   PetscErrorCode ierr;
118155659b69SBarry Smith 
11823a40ed3dSBarry Smith   PetscFunctionBegin;
1183289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1184a5ce6ee0Svictorle     if (lda>m) {
1185a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
1186a5ce6ee0Svictorle 	v = mat->v+j*lda;
1187a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1188a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1189a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1190a5ce6ee0Svictorle #else
1191a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1192a5ce6ee0Svictorle #endif
1193a5ce6ee0Svictorle 	}
1194a5ce6ee0Svictorle       }
1195a5ce6ee0Svictorle     } else {
1196273d9f13SBarry Smith       for (i=0; i<A->n*A->m; i++) {
1197aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1198329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1199289bc588SBarry Smith #else
1200289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1201289bc588SBarry Smith #endif
1202289bc588SBarry Smith       }
1203a5ce6ee0Svictorle     }
1204064f8208SBarry Smith     *nrm = sqrt(sum);
1205efee365bSSatish Balay     ierr = PetscLogFlops(2*A->n*A->m);CHKERRQ(ierr);
12063a40ed3dSBarry Smith   } else if (type == NORM_1) {
1207064f8208SBarry Smith     *nrm = 0.0;
1208273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
12091b807ce4Svictorle       v = mat->v + j*mat->lda;
1210289bc588SBarry Smith       sum = 0.0;
1211273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
121233a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1213289bc588SBarry Smith       }
1214064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1215289bc588SBarry Smith     }
1216efee365bSSatish Balay     ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
12173a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1218064f8208SBarry Smith     *nrm = 0.0;
1219273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1220289bc588SBarry Smith       v = mat->v + j;
1221289bc588SBarry Smith       sum = 0.0;
1222273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
12231b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1224289bc588SBarry Smith       }
1225064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1226289bc588SBarry Smith     }
1227efee365bSSatish Balay     ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
12283a40ed3dSBarry Smith   } else {
122929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1230289bc588SBarry Smith   }
12313a40ed3dSBarry Smith   PetscFunctionReturn(0);
1232289bc588SBarry Smith }
1233289bc588SBarry Smith 
12344a2ae208SSatish Balay #undef __FUNCT__
12354a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1236dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op)
1237289bc588SBarry Smith {
1238c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
123963ba0a88SBarry Smith   PetscErrorCode ierr;
124067e560aaSBarry Smith 
12413a40ed3dSBarry Smith   PetscFunctionBegin;
1242b5a2b587SKris Buschelman   switch (op) {
1243b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1244b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1245b5a2b587SKris Buschelman     break;
1246b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1247b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1248b5a2b587SKris Buschelman     break;
1249b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1250b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1251b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1252b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1253b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1254b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1255b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1256b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1257b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1258b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1259b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
126063ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_SeqDense:Option ignored\n"));CHKERRQ(ierr);
1261b5a2b587SKris Buschelman     break;
126277e54ba9SKris Buschelman   case MAT_SYMMETRIC:
126377e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12649a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12659a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12669a4540c5SBarry Smith   case MAT_HERMITIAN:
12679a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12689a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12699a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
127077e54ba9SKris Buschelman     break;
1271b5a2b587SKris Buschelman   default:
127229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
12733a40ed3dSBarry Smith   }
12743a40ed3dSBarry Smith   PetscFunctionReturn(0);
1275289bc588SBarry Smith }
1276289bc588SBarry Smith 
12774a2ae208SSatish Balay #undef __FUNCT__
12784a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1279dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
12806f0a148fSBarry Smith {
1281ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12826849ba73SBarry Smith   PetscErrorCode ierr;
128313f74950SBarry Smith   PetscInt       lda=l->lda,m=A->m,j;
12843a40ed3dSBarry Smith 
12853a40ed3dSBarry Smith   PetscFunctionBegin;
1286a5ce6ee0Svictorle   if (lda>m) {
1287a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
1288a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1289a5ce6ee0Svictorle     }
1290a5ce6ee0Svictorle   } else {
129187828ca2SBarry Smith     ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1292a5ce6ee0Svictorle   }
12933a40ed3dSBarry Smith   PetscFunctionReturn(0);
12946f0a148fSBarry Smith }
12956f0a148fSBarry Smith 
12964a2ae208SSatish Balay #undef __FUNCT__
12974a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1298dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag)
12996f0a148fSBarry Smith {
1300ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
13016849ba73SBarry Smith   PetscErrorCode ierr;
130213f74950SBarry Smith   PetscInt       n = A->n,i,j,N,*rows;
130387828ca2SBarry Smith   PetscScalar    *slot;
130455659b69SBarry Smith 
13053a40ed3dSBarry Smith   PetscFunctionBegin;
1306e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
130778b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
13086f0a148fSBarry Smith   for (i=0; i<N; i++) {
13096f0a148fSBarry Smith     slot = l->v + rows[i];
13106f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13116f0a148fSBarry Smith   }
13126f0a148fSBarry Smith   if (diag) {
13136f0a148fSBarry Smith     for (i=0; i<N; i++) {
13146f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1315260925b8SBarry Smith       *slot = *diag;
13166f0a148fSBarry Smith     }
13176f0a148fSBarry Smith   }
1318260925b8SBarry Smith   ISRestoreIndices(is,&rows);
13193a40ed3dSBarry Smith   PetscFunctionReturn(0);
13206f0a148fSBarry Smith }
1321557bce09SLois Curfman McInnes 
13224a2ae208SSatish Balay #undef __FUNCT__
13234a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1324dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
132564e87e97SBarry Smith {
1326c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13273a40ed3dSBarry Smith 
13283a40ed3dSBarry Smith   PetscFunctionBegin;
132964e87e97SBarry Smith   *array = mat->v;
13303a40ed3dSBarry Smith   PetscFunctionReturn(0);
133164e87e97SBarry Smith }
13320754003eSLois Curfman McInnes 
13334a2ae208SSatish Balay #undef __FUNCT__
13344a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1335dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1336ff14e315SSatish Balay {
13373a40ed3dSBarry Smith   PetscFunctionBegin;
133809b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13393a40ed3dSBarry Smith   PetscFunctionReturn(0);
1340ff14e315SSatish Balay }
13410754003eSLois Curfman McInnes 
13424a2ae208SSatish Balay #undef __FUNCT__
13434a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
134413f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13450754003eSLois Curfman McInnes {
1346c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13476849ba73SBarry Smith   PetscErrorCode ierr;
134813f74950SBarry Smith   PetscInt       i,j,m = A->m,*irow,*icol,nrows,ncols;
134987828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13500754003eSLois Curfman McInnes   Mat            newmat;
13510754003eSLois Curfman McInnes 
13523a40ed3dSBarry Smith   PetscFunctionBegin;
135378b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
135478b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1355e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1356e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13570754003eSLois Curfman McInnes 
1358182d2002SSatish Balay   /* Check submatrixcall */
1359182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
136013f74950SBarry Smith     PetscInt n_cols,n_rows;
1361182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
136229bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1363182d2002SSatish Balay     newmat = *B;
1364182d2002SSatish Balay   } else {
13650754003eSLois Curfman McInnes     /* Create and fill new matrix */
1366*f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
1367*f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
13685c5985e7SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
13695c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1370182d2002SSatish Balay   }
1371182d2002SSatish Balay 
1372182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1373182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1374182d2002SSatish Balay 
1375182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1376182d2002SSatish Balay     av = v + m*icol[i];
1377182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1378182d2002SSatish Balay       *bv++ = av[irow[j]];
13790754003eSLois Curfman McInnes     }
13800754003eSLois Curfman McInnes   }
1381182d2002SSatish Balay 
1382182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13836d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13846d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13850754003eSLois Curfman McInnes 
13860754003eSLois Curfman McInnes   /* Free work space */
138778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
138878b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1389182d2002SSatish Balay   *B = newmat;
13903a40ed3dSBarry Smith   PetscFunctionReturn(0);
13910754003eSLois Curfman McInnes }
13920754003eSLois Curfman McInnes 
13934a2ae208SSatish Balay #undef __FUNCT__
13944a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
139513f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1396905e6a2fSBarry Smith {
13976849ba73SBarry Smith   PetscErrorCode ierr;
139813f74950SBarry Smith   PetscInt       i;
1399905e6a2fSBarry Smith 
14003a40ed3dSBarry Smith   PetscFunctionBegin;
1401905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1402b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1403905e6a2fSBarry Smith   }
1404905e6a2fSBarry Smith 
1405905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14066a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1407905e6a2fSBarry Smith   }
14083a40ed3dSBarry Smith   PetscFunctionReturn(0);
1409905e6a2fSBarry Smith }
1410905e6a2fSBarry Smith 
14114a2ae208SSatish Balay #undef __FUNCT__
14124a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1413dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14144b0e389bSBarry Smith {
14154b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14166849ba73SBarry Smith   PetscErrorCode ierr;
141713f74950SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j;
14183a40ed3dSBarry Smith 
14193a40ed3dSBarry Smith   PetscFunctionBegin;
142033f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
142133f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1422cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14233a40ed3dSBarry Smith     PetscFunctionReturn(0);
14243a40ed3dSBarry Smith   }
14250dbb7854Svictorle   if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1426a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14270dbb7854Svictorle     for (j=0; j<n; j++) {
1428a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1429a5ce6ee0Svictorle     }
1430a5ce6ee0Svictorle   } else {
143187828ca2SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1432a5ce6ee0Svictorle   }
1433273d9f13SBarry Smith   PetscFunctionReturn(0);
1434273d9f13SBarry Smith }
1435273d9f13SBarry Smith 
14364a2ae208SSatish Balay #undef __FUNCT__
14374a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1438dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1439273d9f13SBarry Smith {
1440dfbe8321SBarry Smith   PetscErrorCode ierr;
1441273d9f13SBarry Smith 
1442273d9f13SBarry Smith   PetscFunctionBegin;
1443273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14443a40ed3dSBarry Smith   PetscFunctionReturn(0);
14454b0e389bSBarry Smith }
14464b0e389bSBarry Smith 
1447289bc588SBarry Smith /* -------------------------------------------------------------------*/
1448a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1449905e6a2fSBarry Smith        MatGetRow_SeqDense,
1450905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1451905e6a2fSBarry Smith        MatMult_SeqDense,
145297304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
14537c922b88SBarry Smith        MatMultTranspose_SeqDense,
14547c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1455905e6a2fSBarry Smith        MatSolve_SeqDense,
1456905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
14577c922b88SBarry Smith        MatSolveTranspose_SeqDense,
145897304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1459905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1460905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1461ec8511deSBarry Smith        MatRelax_SeqDense,
1462ec8511deSBarry Smith        MatTranspose_SeqDense,
146397304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1464905e6a2fSBarry Smith        MatEqual_SeqDense,
1465905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1466905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1467905e6a2fSBarry Smith        MatNorm_SeqDense,
146897304618SKris Buschelman /*20*/ 0,
1469905e6a2fSBarry Smith        0,
1470905e6a2fSBarry Smith        0,
1471905e6a2fSBarry Smith        MatSetOption_SeqDense,
1472905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
147397304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1474905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1475905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1476905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1477905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
147897304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1479273d9f13SBarry Smith        0,
1480905e6a2fSBarry Smith        0,
1481905e6a2fSBarry Smith        MatGetArray_SeqDense,
1482905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
148397304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1484a5ae1ecdSBarry Smith        0,
1485a5ae1ecdSBarry Smith        0,
1486a5ae1ecdSBarry Smith        0,
1487a5ae1ecdSBarry Smith        0,
148897304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1489a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1490a5ae1ecdSBarry Smith        0,
14914b0e389bSBarry Smith        MatGetValues_SeqDense,
1492a5ae1ecdSBarry Smith        MatCopy_SeqDense,
149397304618SKris Buschelman /*45*/ 0,
1494a5ae1ecdSBarry Smith        MatScale_SeqDense,
1495a5ae1ecdSBarry Smith        0,
1496a5ae1ecdSBarry Smith        0,
1497a5ae1ecdSBarry Smith        0,
1498521d7252SBarry Smith /*50*/ 0,
1499a5ae1ecdSBarry Smith        0,
1500a5ae1ecdSBarry Smith        0,
1501a5ae1ecdSBarry Smith        0,
1502a5ae1ecdSBarry Smith        0,
150397304618SKris Buschelman /*55*/ 0,
1504a5ae1ecdSBarry Smith        0,
1505a5ae1ecdSBarry Smith        0,
1506a5ae1ecdSBarry Smith        0,
1507a5ae1ecdSBarry Smith        0,
150897304618SKris Buschelman /*60*/ 0,
1509e03a110bSBarry Smith        MatDestroy_SeqDense,
1510e03a110bSBarry Smith        MatView_SeqDense,
151197304618SKris Buschelman        MatGetPetscMaps_Petsc,
151297304618SKris Buschelman        0,
151397304618SKris Buschelman /*65*/ 0,
151497304618SKris Buschelman        0,
151597304618SKris Buschelman        0,
151697304618SKris Buschelman        0,
151797304618SKris Buschelman        0,
151897304618SKris Buschelman /*70*/ 0,
151997304618SKris Buschelman        0,
152097304618SKris Buschelman        0,
152197304618SKris Buschelman        0,
152297304618SKris Buschelman        0,
152397304618SKris Buschelman /*75*/ 0,
152497304618SKris Buschelman        0,
152597304618SKris Buschelman        0,
152697304618SKris Buschelman        0,
152797304618SKris Buschelman        0,
152897304618SKris Buschelman /*80*/ 0,
152997304618SKris Buschelman        0,
153097304618SKris Buschelman        0,
153197304618SKris Buschelman        0,
1532865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1533865e5f61SKris Buschelman        0,
1534865e5f61SKris Buschelman        0,
1535865e5f61SKris Buschelman        0,
1536865e5f61SKris Buschelman        0,
1537865e5f61SKris Buschelman        0,
1538865e5f61SKris Buschelman /*90*/ 0,
1539865e5f61SKris Buschelman        0,
1540865e5f61SKris Buschelman        0,
1541865e5f61SKris Buschelman        0,
1542865e5f61SKris Buschelman        0,
1543865e5f61SKris Buschelman /*95*/ 0,
1544865e5f61SKris Buschelman        0,
1545865e5f61SKris Buschelman        0,
1546865e5f61SKris Buschelman        0};
154790ace30eSBarry Smith 
15484a2ae208SSatish Balay #undef __FUNCT__
15494a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
15504b828684SBarry Smith /*@C
1551fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1552d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1553d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1554289bc588SBarry Smith 
1555db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1556db81eaa0SLois Curfman McInnes 
155720563c6bSBarry Smith    Input Parameters:
1558db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
15590c775827SLois Curfman McInnes .  m - number of rows
156018f449edSLois Curfman McInnes .  n - number of columns
1561db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1562dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
156320563c6bSBarry Smith 
156420563c6bSBarry Smith    Output Parameter:
156544cd7ae7SLois Curfman McInnes .  A - the matrix
156620563c6bSBarry Smith 
1567b259b22eSLois Curfman McInnes    Notes:
156818f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
156918f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1570b4fd4287SBarry Smith    set data=PETSC_NULL.
157118f449edSLois Curfman McInnes 
1572027ccd11SLois Curfman McInnes    Level: intermediate
1573027ccd11SLois Curfman McInnes 
1574dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1575d65003e9SLois Curfman McInnes 
1576db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
157720563c6bSBarry Smith @*/
1578be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1579289bc588SBarry Smith {
1580dfbe8321SBarry Smith   PetscErrorCode ierr;
15813b2fbd54SBarry Smith 
15823a40ed3dSBarry Smith   PetscFunctionBegin;
1583*f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1584*f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1585273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1586273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1587273d9f13SBarry Smith   PetscFunctionReturn(0);
1588273d9f13SBarry Smith }
1589273d9f13SBarry Smith 
15904a2ae208SSatish Balay #undef __FUNCT__
15914a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1592273d9f13SBarry Smith /*@C
1593273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1594273d9f13SBarry Smith 
1595273d9f13SBarry Smith    Collective on MPI_Comm
1596273d9f13SBarry Smith 
1597273d9f13SBarry Smith    Input Parameters:
1598273d9f13SBarry Smith +  A - the matrix
1599273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1600273d9f13SBarry Smith 
1601273d9f13SBarry Smith    Notes:
1602273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1603273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1604273d9f13SBarry Smith    set data=PETSC_NULL.
1605273d9f13SBarry Smith 
1606273d9f13SBarry Smith    Level: intermediate
1607273d9f13SBarry Smith 
1608273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1609273d9f13SBarry Smith 
1610273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1611273d9f13SBarry Smith @*/
1612be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1613273d9f13SBarry Smith {
1614dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1615a23d5eceSKris Buschelman 
1616a23d5eceSKris Buschelman   PetscFunctionBegin;
1617a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1618a23d5eceSKris Buschelman   if (f) {
1619a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1620a23d5eceSKris Buschelman   }
1621a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1622a23d5eceSKris Buschelman }
1623a23d5eceSKris Buschelman 
1624a23d5eceSKris Buschelman EXTERN_C_BEGIN
1625a23d5eceSKris Buschelman #undef __FUNCT__
1626a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1627be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1628a23d5eceSKris Buschelman {
1629273d9f13SBarry Smith   Mat_SeqDense   *b;
1630dfbe8321SBarry Smith   PetscErrorCode ierr;
1631273d9f13SBarry Smith 
1632273d9f13SBarry Smith   PetscFunctionBegin;
1633273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1634273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1635273d9f13SBarry Smith   if (!data) {
163687828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
163787828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1638273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
163952e6d16bSBarry Smith     ierr = PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));CHKERRQ(ierr);
1640273d9f13SBarry Smith   } else { /* user-allocated storage */
1641273d9f13SBarry Smith     b->v          = data;
1642273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1643273d9f13SBarry Smith   }
1644273d9f13SBarry Smith   PetscFunctionReturn(0);
1645273d9f13SBarry Smith }
1646a23d5eceSKris Buschelman EXTERN_C_END
1647273d9f13SBarry Smith 
16481b807ce4Svictorle #undef __FUNCT__
16491b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
16501b807ce4Svictorle /*@C
16511b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
16521b807ce4Svictorle 
16531b807ce4Svictorle   Input parameter:
16541b807ce4Svictorle + A - the matrix
16551b807ce4Svictorle - lda - the leading dimension
16561b807ce4Svictorle 
16571b807ce4Svictorle   Notes:
16581b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
16591b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
16601b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
16611b807ce4Svictorle 
16621b807ce4Svictorle   Level: intermediate
16631b807ce4Svictorle 
16641b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
16651b807ce4Svictorle 
16661b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation()
16671b807ce4Svictorle @*/
1668be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
16691b807ce4Svictorle {
16701b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
16711b807ce4Svictorle   PetscFunctionBegin;
167277431f27SBarry Smith   if (lda < B->m) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->m);
16731b807ce4Svictorle   b->lda = lda;
16741b807ce4Svictorle   PetscFunctionReturn(0);
16751b807ce4Svictorle }
16761b807ce4Svictorle 
16770bad9183SKris Buschelman /*MC
1678fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
16790bad9183SKris Buschelman 
16800bad9183SKris Buschelman    Options Database Keys:
16810bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
16820bad9183SKris Buschelman 
16830bad9183SKris Buschelman   Level: beginner
16840bad9183SKris Buschelman 
16850bad9183SKris Buschelman .seealso: MatCreateSeqDense
16860bad9183SKris Buschelman M*/
16870bad9183SKris Buschelman 
1688273d9f13SBarry Smith EXTERN_C_BEGIN
16894a2ae208SSatish Balay #undef __FUNCT__
16904a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1691be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
1692273d9f13SBarry Smith {
1693273d9f13SBarry Smith   Mat_SeqDense   *b;
1694dfbe8321SBarry Smith   PetscErrorCode ierr;
16957c334f02SBarry Smith   PetscMPIInt    size;
1696273d9f13SBarry Smith 
1697273d9f13SBarry Smith   PetscFunctionBegin;
1698273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
169929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
170055659b69SBarry Smith 
1701273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1702273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1703273d9f13SBarry Smith 
1704b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1705549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
170644cd7ae7SLois Curfman McInnes   B->factor       = 0;
170790f02eecSBarry Smith   B->mapping      = 0;
170852e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr);
170944cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
171018f449edSLois Curfman McInnes 
17118a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
17128a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1713a5ae1ecdSBarry Smith 
171444cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1715273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1716273d9f13SBarry Smith   b->v            = 0;
17171b807ce4Svictorle   b->lda          = B->m;
17184e220ebcSLois Curfman McInnes 
1719a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1720a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1721a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
17223a40ed3dSBarry Smith   PetscFunctionReturn(0);
1723289bc588SBarry Smith }
1724273d9f13SBarry Smith EXTERN_C_END
1725