xref: /petsc/src/mat/impls/dense/seq/dense.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris 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;
1235c5985e7SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&newi);CHKERRQ(ierr);
1245c5985e7SKris Buschelman   ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr);
1255c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
1265609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
127a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
128a5ce6ee0Svictorle     if (lda>A->m) {
129a5ce6ee0Svictorle       m = A->m;
130a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
131a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
132a5ce6ee0Svictorle       }
133a5ce6ee0Svictorle     } else {
13487828ca2SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
13502cad45dSBarry Smith     }
136a5ce6ee0Svictorle   }
13702cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
13802cad45dSBarry Smith   *newmat = newi;
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
14002cad45dSBarry Smith }
14102cad45dSBarry Smith 
1424a2ae208SSatish Balay #undef __FUNCT__
1434a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
144dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
145289bc588SBarry Smith {
146dfbe8321SBarry Smith   PetscErrorCode ierr;
1473a40ed3dSBarry Smith 
1483a40ed3dSBarry Smith   PetscFunctionBegin;
1495609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1503a40ed3dSBarry Smith   PetscFunctionReturn(0);
151289bc588SBarry Smith }
1526ee01492SSatish Balay 
1534a2ae208SSatish Balay #undef __FUNCT__
1544a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
155af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
156289bc588SBarry Smith {
15702cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1586849ba73SBarry Smith   PetscErrorCode ierr;
15913f74950SBarry Smith   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j;
1604482741eSBarry Smith   MatFactorInfo  info;
1613a40ed3dSBarry Smith 
1623a40ed3dSBarry Smith   PetscFunctionBegin;
16302cad45dSBarry Smith   /* copy the numerical values */
1641b807ce4Svictorle   if (lda1>m || lda2>m ) {
1651b807ce4Svictorle     for (j=0; j<n; j++) {
1661b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1671b807ce4Svictorle     }
1681b807ce4Svictorle   } else {
16987828ca2SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1701b807ce4Svictorle   }
17102cad45dSBarry Smith   (*fact)->factor = 0;
1724482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
1733a40ed3dSBarry Smith   PetscFunctionReturn(0);
174289bc588SBarry Smith }
1756ee01492SSatish Balay 
1764a2ae208SSatish Balay #undef __FUNCT__
1774a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
178dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
179289bc588SBarry Smith {
180dfbe8321SBarry Smith   PetscErrorCode ierr;
1813a40ed3dSBarry Smith 
1823a40ed3dSBarry Smith   PetscFunctionBegin;
183ceb03754SKris Buschelman   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr);
1843a40ed3dSBarry Smith   PetscFunctionReturn(0);
185289bc588SBarry Smith }
1866ee01492SSatish Balay 
1874a2ae208SSatish Balay #undef __FUNCT__
1884a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
189dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
190289bc588SBarry Smith {
191a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
192a07ab388SBarry Smith   PetscFunctionBegin;
193a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
194a07ab388SBarry Smith #else
195c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1966849ba73SBarry Smith   PetscErrorCode ierr;
1974ce68768SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->n,info;
19855659b69SBarry Smith 
1993a40ed3dSBarry Smith   PetscFunctionBegin;
200752f5567SLois Curfman McInnes   if (mat->pivots) {
201606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
20252e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-A->m*sizeof(PetscInt));CHKERRQ(ierr);
203752f5567SLois Curfman McInnes     mat->pivots = 0;
204752f5567SLois Curfman McInnes   }
205273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
20671044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
20777431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
208c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
209efee365bSSatish Balay   ierr = PetscLogFlops((A->n*A->n*A->n)/3);CHKERRQ(ierr);
210a07ab388SBarry Smith #endif
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
212289bc588SBarry Smith }
213289bc588SBarry Smith 
2144a2ae208SSatish Balay #undef __FUNCT__
2150b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
216af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2170b4b3355SBarry Smith {
218dfbe8321SBarry Smith   PetscErrorCode ierr;
21915e8a5b3SHong Zhang   MatFactorInfo  info;
2200b4b3355SBarry Smith 
2210b4b3355SBarry Smith   PetscFunctionBegin;
22215e8a5b3SHong Zhang   info.fill = 1.0;
22315e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2240b4b3355SBarry Smith   PetscFunctionReturn(0);
2250b4b3355SBarry Smith }
2260b4b3355SBarry Smith 
2270b4b3355SBarry Smith #undef __FUNCT__
2284a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
229dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
230289bc588SBarry Smith {
231c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2326849ba73SBarry Smith   PetscErrorCode ierr;
2334ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, one = 1,info;
23487828ca2SBarry Smith   PetscScalar    *x,*y;
23567e560aaSBarry Smith 
2363a40ed3dSBarry Smith   PetscFunctionBegin;
237273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2381ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2391ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
24087828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
241c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
242ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
24380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
244ae7cfcebSSatish Balay #else
24571044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2464ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
247ae7cfcebSSatish Balay #endif
2487a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
249ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
25080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
251ae7cfcebSSatish Balay #else
25271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2534ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
254ae7cfcebSSatish Balay #endif
255289bc588SBarry Smith   }
25629bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2571ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2581ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
259efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr);
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
261289bc588SBarry Smith }
2626ee01492SSatish Balay 
2634a2ae208SSatish Balay #undef __FUNCT__
2644a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
265dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
266da3a660dSBarry Smith {
267c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
268dfbe8321SBarry Smith   PetscErrorCode ierr;
2694ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt) A->m,one = 1,info;
27087828ca2SBarry Smith   PetscScalar    *x,*y;
27167e560aaSBarry Smith 
2723a40ed3dSBarry Smith   PetscFunctionBegin;
273273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2741ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2751ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
27687828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
277752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
278da3a660dSBarry Smith   if (mat->pivots) {
279ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
28080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
281ae7cfcebSSatish Balay #else
28271044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2834ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
284ae7cfcebSSatish Balay #endif
2857a97a34bSBarry Smith   } else {
286ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
28780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
288ae7cfcebSSatish Balay #else
28971044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2904ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
291ae7cfcebSSatish Balay #endif
292da3a660dSBarry Smith   }
2931ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2941ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
295efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr);
2963a40ed3dSBarry Smith   PetscFunctionReturn(0);
297da3a660dSBarry Smith }
2986ee01492SSatish Balay 
2994a2ae208SSatish Balay #undef __FUNCT__
3004a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
301dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
302da3a660dSBarry Smith {
303c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
304dfbe8321SBarry Smith   PetscErrorCode ierr;
3054ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
30687828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
307da3a660dSBarry Smith   Vec            tmp = 0;
30867e560aaSBarry Smith 
3093a40ed3dSBarry Smith   PetscFunctionBegin;
3101ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3111ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
312273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
313da3a660dSBarry Smith   if (yy == zz) {
31478b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
31552e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
31678b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
317da3a660dSBarry Smith   }
31887828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
319752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
320da3a660dSBarry Smith   if (mat->pivots) {
321ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
32280444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
323ae7cfcebSSatish Balay #else
32471044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3252ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
326ae7cfcebSSatish Balay #endif
327a8c6a408SBarry Smith   } else {
328ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
32980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
330ae7cfcebSSatish Balay #else
33171044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3322ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
333ae7cfcebSSatish Balay #endif
334da3a660dSBarry Smith   }
335a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
336a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
3371ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3381ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
339efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr);
3403a40ed3dSBarry Smith   PetscFunctionReturn(0);
341da3a660dSBarry Smith }
34267e560aaSBarry Smith 
3434a2ae208SSatish Balay #undef __FUNCT__
3444a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
345dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
346da3a660dSBarry Smith {
347c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3486849ba73SBarry Smith   PetscErrorCode ierr;
3494ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
35087828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
351da3a660dSBarry Smith   Vec            tmp;
35267e560aaSBarry Smith 
3533a40ed3dSBarry Smith   PetscFunctionBegin;
354273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3551ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3561ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
357da3a660dSBarry Smith   if (yy == zz) {
35878b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
35952e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
36078b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
361da3a660dSBarry Smith   }
36287828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
363752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
364da3a660dSBarry Smith   if (mat->pivots) {
365ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
36680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
367ae7cfcebSSatish Balay #else
36871044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3692ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
370ae7cfcebSSatish Balay #endif
3713a40ed3dSBarry Smith   } else {
372ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
37380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
374ae7cfcebSSatish Balay #else
37571044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3762ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
377ae7cfcebSSatish Balay #endif
378da3a660dSBarry Smith   }
37990f02eecSBarry Smith   if (tmp) {
38090f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
38190f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3823a40ed3dSBarry Smith   } else {
38390f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
38490f02eecSBarry Smith   }
3851ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3861ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
387efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr);
3883a40ed3dSBarry Smith   PetscFunctionReturn(0);
389da3a660dSBarry Smith }
390289bc588SBarry Smith /* ------------------------------------------------------------------*/
3914a2ae208SSatish Balay #undef __FUNCT__
3924a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
39313f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
394289bc588SBarry Smith {
395c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
39687828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
397dfbe8321SBarry Smith   PetscErrorCode ierr;
39813f74950SBarry Smith   PetscInt       m = A->m,i;
399aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4004ce68768SBarry Smith   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
401bc1b551cSSatish Balay #endif
402289bc588SBarry Smith 
4033a40ed3dSBarry Smith   PetscFunctionBegin;
404289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
40571044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
406bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
407289bc588SBarry Smith   }
4081ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4091ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
410b965ef7fSBarry Smith   its  = its*lits;
41177431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
412289bc588SBarry Smith   while (its--) {
413289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
414289bc588SBarry Smith       for (i=0; i<m; i++) {
415aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
416f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
417f1747703SBarry Smith            not happy about returning a double complex */
41813f74950SBarry Smith         PetscInt         _i;
41987828ca2SBarry Smith         PetscScalar sum = b[i];
420f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4213f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
422f1747703SBarry Smith         }
423f1747703SBarry Smith         xt = sum;
424f1747703SBarry Smith #else
42571044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
426f1747703SBarry Smith #endif
42755a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
428289bc588SBarry Smith       }
429289bc588SBarry Smith     }
430289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
431289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
432aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
433f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
434f1747703SBarry Smith            not happy about returning a double complex */
43513f74950SBarry Smith         PetscInt         _i;
43687828ca2SBarry Smith         PetscScalar sum = b[i];
437f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4383f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
439f1747703SBarry Smith         }
440f1747703SBarry Smith         xt = sum;
441f1747703SBarry Smith #else
44271044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
443f1747703SBarry Smith #endif
44455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
445289bc588SBarry Smith       }
446289bc588SBarry Smith     }
447289bc588SBarry Smith   }
4481ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4491ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4503a40ed3dSBarry Smith   PetscFunctionReturn(0);
451289bc588SBarry Smith }
452289bc588SBarry Smith 
453289bc588SBarry Smith /* -----------------------------------------------------------------*/
4544a2ae208SSatish Balay #undef __FUNCT__
4554a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
456dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
457289bc588SBarry Smith {
458c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
45987828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
460dfbe8321SBarry Smith   PetscErrorCode ierr;
4614ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n,_One=1;
462ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
4633a40ed3dSBarry Smith 
4643a40ed3dSBarry Smith   PetscFunctionBegin;
465273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4661ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4671ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
46871044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
4691ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4701ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
471efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n - A->n);CHKERRQ(ierr);
4723a40ed3dSBarry Smith   PetscFunctionReturn(0);
473289bc588SBarry Smith }
4746ee01492SSatish Balay 
4754a2ae208SSatish Balay #undef __FUNCT__
4764a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
477dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
478289bc588SBarry Smith {
479c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
48087828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
481dfbe8321SBarry Smith   PetscErrorCode ierr;
4824ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
4833a40ed3dSBarry Smith 
4843a40ed3dSBarry Smith   PetscFunctionBegin;
485273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4871ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
48871044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
4891ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4901ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
491efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n - A->m);CHKERRQ(ierr);
4923a40ed3dSBarry Smith   PetscFunctionReturn(0);
493289bc588SBarry Smith }
4946ee01492SSatish Balay 
4954a2ae208SSatish Balay #undef __FUNCT__
4964a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
497dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
498289bc588SBarry Smith {
499c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
50087828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
501dfbe8321SBarry Smith   PetscErrorCode ierr;
5024ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
5033a40ed3dSBarry Smith 
5043a40ed3dSBarry Smith   PetscFunctionBegin;
505273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
506600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5071ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5081ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
50971044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5101ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5111ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
512efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr);
5133a40ed3dSBarry Smith   PetscFunctionReturn(0);
514289bc588SBarry Smith }
5156ee01492SSatish Balay 
5164a2ae208SSatish Balay #undef __FUNCT__
5174a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
518dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
519289bc588SBarry Smith {
520c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
52187828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
522dfbe8321SBarry Smith   PetscErrorCode ierr;
5234ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
52487828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5253a40ed3dSBarry Smith 
5263a40ed3dSBarry Smith   PetscFunctionBegin;
527273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
528600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5291ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5301ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
53171044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5321ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5331ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
534efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr);
5353a40ed3dSBarry Smith   PetscFunctionReturn(0);
536289bc588SBarry Smith }
537289bc588SBarry Smith 
538289bc588SBarry Smith /* -----------------------------------------------------------------*/
5394a2ae208SSatish Balay #undef __FUNCT__
5404a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
54113f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
542289bc588SBarry Smith {
543c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
54487828ca2SBarry Smith   PetscScalar    *v;
5456849ba73SBarry Smith   PetscErrorCode ierr;
54613f74950SBarry Smith   PetscInt       i;
54767e560aaSBarry Smith 
5483a40ed3dSBarry Smith   PetscFunctionBegin;
549273d9f13SBarry Smith   *ncols = A->n;
550289bc588SBarry Smith   if (cols) {
55113f74950SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
552273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
553289bc588SBarry Smith   }
554289bc588SBarry Smith   if (vals) {
55587828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
556289bc588SBarry Smith     v    = mat->v + row;
5571b807ce4Svictorle     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;}
558289bc588SBarry Smith   }
5593a40ed3dSBarry Smith   PetscFunctionReturn(0);
560289bc588SBarry Smith }
5616ee01492SSatish Balay 
5624a2ae208SSatish Balay #undef __FUNCT__
5634a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
56413f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
565289bc588SBarry Smith {
566dfbe8321SBarry Smith   PetscErrorCode ierr;
567606d414cSSatish Balay   PetscFunctionBegin;
568606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
569606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5703a40ed3dSBarry Smith   PetscFunctionReturn(0);
571289bc588SBarry Smith }
572289bc588SBarry Smith /* ----------------------------------------------------------------*/
5734a2ae208SSatish Balay #undef __FUNCT__
5744a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
57513f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
576289bc588SBarry Smith {
577c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
57813f74950SBarry Smith   PetscInt     i,j,idx=0;
579d6dfbf8fSBarry Smith 
5803a40ed3dSBarry Smith   PetscFunctionBegin;
581289bc588SBarry Smith   if (!mat->roworiented) {
582dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
583289bc588SBarry Smith       for (j=0; j<n; j++) {
584cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5852515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
58677431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
58758804f6dSBarry Smith #endif
588289bc588SBarry Smith         for (i=0; i<m; i++) {
589cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5902515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
59177431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
59258804f6dSBarry Smith #endif
593cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
594289bc588SBarry Smith         }
595289bc588SBarry Smith       }
5963a40ed3dSBarry Smith     } else {
597289bc588SBarry Smith       for (j=0; j<n; j++) {
598cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5992515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60077431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
60158804f6dSBarry Smith #endif
602289bc588SBarry Smith         for (i=0; i<m; i++) {
603cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6042515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60577431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
60658804f6dSBarry Smith #endif
607cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
608289bc588SBarry Smith         }
609289bc588SBarry Smith       }
610289bc588SBarry Smith     }
6113a40ed3dSBarry Smith   } else {
612dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
613e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
614cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6152515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
61677431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
61758804f6dSBarry Smith #endif
618e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
619cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6202515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
62177431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
62258804f6dSBarry Smith #endif
623cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
624e8d4e0b9SBarry Smith         }
625e8d4e0b9SBarry Smith       }
6263a40ed3dSBarry Smith     } else {
627289bc588SBarry Smith       for (i=0; i<m; i++) {
628cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6292515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
63077431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
63158804f6dSBarry Smith #endif
632289bc588SBarry Smith         for (j=0; j<n; j++) {
633cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6342515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
63577431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
63658804f6dSBarry Smith #endif
637cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
638289bc588SBarry Smith         }
639289bc588SBarry Smith       }
640289bc588SBarry Smith     }
641e8d4e0b9SBarry Smith   }
6423a40ed3dSBarry Smith   PetscFunctionReturn(0);
643289bc588SBarry Smith }
644e8d4e0b9SBarry Smith 
6454a2ae208SSatish Balay #undef __FUNCT__
6464a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
64713f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
648ae80bb75SLois Curfman McInnes {
649ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
65013f74950SBarry Smith   PetscInt     i,j;
65187828ca2SBarry Smith   PetscScalar  *vpt = v;
652ae80bb75SLois Curfman McInnes 
6533a40ed3dSBarry Smith   PetscFunctionBegin;
654ae80bb75SLois Curfman McInnes   /* row-oriented output */
655ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
656ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6571b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
658ae80bb75SLois Curfman McInnes     }
659ae80bb75SLois Curfman McInnes   }
6603a40ed3dSBarry Smith   PetscFunctionReturn(0);
661ae80bb75SLois Curfman McInnes }
662ae80bb75SLois Curfman McInnes 
663289bc588SBarry Smith /* -----------------------------------------------------------------*/
664289bc588SBarry Smith 
665e090d566SSatish Balay #include "petscsys.h"
66688e32edaSLois Curfman McInnes 
6674a2ae208SSatish Balay #undef __FUNCT__
6684a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
669dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer,const MatType type,Mat *A)
67088e32edaSLois Curfman McInnes {
67188e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
67288e32edaSLois Curfman McInnes   Mat            B;
6736849ba73SBarry Smith   PetscErrorCode ierr;
67413f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
67513f74950SBarry Smith   int            fd;
67613f74950SBarry Smith   PetscMPIInt    size;
67713f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
67887828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
67919bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
68088e32edaSLois Curfman McInnes 
6813a40ed3dSBarry Smith   PetscFunctionBegin;
682d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
68329bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
684b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6850752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
686552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
68788e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
68888e32edaSLois Curfman McInnes 
68990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
6905c5985e7SKris Buschelman     ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr);
691be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6925c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
69390ace30eSBarry Smith     B    = *A;
69490ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6954a41a4d3SSatish Balay     v    = a->v;
6964a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6974a41a4d3SSatish Balay        from row major to column major */
698b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
69990ace30eSBarry Smith     /* read in nonzero values */
7004a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7014a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7024a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7034a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7044a41a4d3SSatish Balay         *v++ =w[i*N+j];
7054a41a4d3SSatish Balay       }
7064a41a4d3SSatish Balay     }
707606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7086d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7096d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71090ace30eSBarry Smith   } else {
71188e32edaSLois Curfman McInnes     /* read row lengths */
71213f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7130752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
71488e32edaSLois Curfman McInnes 
71588e32edaSLois Curfman McInnes     /* create our matrix */
7165c5985e7SKris Buschelman     ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr);
717be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7185c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
71988e32edaSLois Curfman McInnes     B = *A;
72088e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
72188e32edaSLois Curfman McInnes     v = a->v;
72288e32edaSLois Curfman McInnes 
72388e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
72413f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
725b0a32e0cSBarry Smith     cols = scols;
7260752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
72787828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
728b0a32e0cSBarry Smith     vals = svals;
7290752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
73088e32edaSLois Curfman McInnes 
73188e32edaSLois Curfman McInnes     /* insert into matrix */
73288e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
73388e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
73488e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
73588e32edaSLois Curfman McInnes     }
736606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
737606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
738606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
73988e32edaSLois Curfman McInnes 
7406d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7416d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74290ace30eSBarry Smith   }
7433a40ed3dSBarry Smith   PetscFunctionReturn(0);
74488e32edaSLois Curfman McInnes }
74588e32edaSLois Curfman McInnes 
746e090d566SSatish Balay #include "petscsys.h"
747932b0c3eSLois Curfman McInnes 
7484a2ae208SSatish Balay #undef __FUNCT__
7494a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
7506849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
751289bc588SBarry Smith {
752932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
753dfbe8321SBarry Smith   PetscErrorCode    ierr;
75413f74950SBarry Smith   PetscInt          i,j;
755fb9695e5SSatish Balay   char              *name;
75687828ca2SBarry Smith   PetscScalar       *v;
757f3ef73ceSBarry Smith   PetscViewerFormat format;
758932b0c3eSLois Curfman McInnes 
7593a40ed3dSBarry Smith   PetscFunctionBegin;
760435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
761b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
762456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7633a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
764fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
765b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
766273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
76744cd7ae7SLois Curfman McInnes       v = a->v + i;
76877431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
769273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
770aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
771329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
77277431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
773329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
77477431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7756831982aSBarry Smith         }
77680cd9d93SLois Curfman McInnes #else
7776831982aSBarry Smith         if (*v) {
77877431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,*v);CHKERRQ(ierr);
7796831982aSBarry Smith         }
78080cd9d93SLois Curfman McInnes #endif
7811b807ce4Svictorle         v += a->lda;
78280cd9d93SLois Curfman McInnes       }
783b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
78480cd9d93SLois Curfman McInnes     }
785b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7863a40ed3dSBarry Smith   } else {
787b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
788aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
789ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
79047989497SBarry Smith     /* determine if matrix has all real values */
79147989497SBarry Smith     v = a->v;
792201f5f94SSatish Balay     for (i=0; i<A->m*A->n; i++) {
793ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
79447989497SBarry Smith     }
79547989497SBarry Smith #endif
796fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
7973a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
79877431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->m,A->n);CHKERRQ(ierr);
79977431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->m,A->n);CHKERRQ(ierr);
800fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
801ffac6cdbSBarry Smith     }
802ffac6cdbSBarry Smith 
803273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
804932b0c3eSLois Curfman McInnes       v = a->v + i;
805273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
806aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
80747989497SBarry Smith         if (allreal) {
8081b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
80947989497SBarry Smith         } else {
8101b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
81147989497SBarry Smith         }
812289bc588SBarry Smith #else
8131b807ce4Svictorle         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
814289bc588SBarry Smith #endif
8151b807ce4Svictorle 	v += a->lda;
816289bc588SBarry Smith       }
817b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
818289bc588SBarry Smith     }
819fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
820b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
821ffac6cdbSBarry Smith     }
822b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
823da3a660dSBarry Smith   }
824b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8253a40ed3dSBarry Smith   PetscFunctionReturn(0);
826289bc588SBarry Smith }
827289bc588SBarry Smith 
8284a2ae208SSatish Balay #undef __FUNCT__
8294a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8306849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
831932b0c3eSLois Curfman McInnes {
832932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8336849ba73SBarry Smith   PetscErrorCode    ierr;
83413f74950SBarry Smith   int               fd;
83513f74950SBarry Smith   PetscInt          ict,j,n = A->n,m = A->m,i,*col_lens,nz = m*n;
83687828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
837f3ef73ceSBarry Smith   PetscViewerFormat format;
838932b0c3eSLois Curfman McInnes 
8393a40ed3dSBarry Smith   PetscFunctionBegin;
840b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
84190ace30eSBarry Smith 
842b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
843fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
84490ace30eSBarry Smith     /* store the matrix as a dense matrix */
84513f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
846552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
84790ace30eSBarry Smith     col_lens[1] = m;
84890ace30eSBarry Smith     col_lens[2] = n;
84990ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8506f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
851606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
85290ace30eSBarry Smith 
85390ace30eSBarry Smith     /* write out matrix, by rows */
85487828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
85590ace30eSBarry Smith     v    = a->v;
85690ace30eSBarry Smith     for (i=0; i<m; i++) {
85790ace30eSBarry Smith       for (j=0; j<n; j++) {
85890ace30eSBarry Smith         vals[i + j*m] = *v++;
85990ace30eSBarry Smith       }
86090ace30eSBarry Smith     }
8616f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
862606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
86390ace30eSBarry Smith   } else {
86413f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
865552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
866932b0c3eSLois Curfman McInnes     col_lens[1] = m;
867932b0c3eSLois Curfman McInnes     col_lens[2] = n;
868932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
869932b0c3eSLois Curfman McInnes 
870932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
871932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8726f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
873932b0c3eSLois Curfman McInnes 
874932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
875932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
876932b0c3eSLois Curfman McInnes     ict = 0;
877932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
878932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
879932b0c3eSLois Curfman McInnes     }
8806f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
881606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
882932b0c3eSLois Curfman McInnes 
883932b0c3eSLois Curfman McInnes     /* store nonzero values */
88487828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
885932b0c3eSLois Curfman McInnes     ict  = 0;
886932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
887932b0c3eSLois Curfman McInnes       v = a->v + i;
888932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
8891b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
890932b0c3eSLois Curfman McInnes       }
891932b0c3eSLois Curfman McInnes     }
8926f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
893606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
89490ace30eSBarry Smith   }
8953a40ed3dSBarry Smith   PetscFunctionReturn(0);
896932b0c3eSLois Curfman McInnes }
897932b0c3eSLois Curfman McInnes 
8984a2ae208SSatish Balay #undef __FUNCT__
8994a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
900dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
901f1af5d2fSBarry Smith {
902f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
903f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9046849ba73SBarry Smith   PetscErrorCode    ierr;
90513f74950SBarry Smith   PetscInt          m = A->m,n = A->n,color,i,j;
90687828ca2SBarry Smith   PetscScalar       *v = a->v;
907b0a32e0cSBarry Smith   PetscViewer       viewer;
908b0a32e0cSBarry Smith   PetscDraw         popup;
909329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
910f3ef73ceSBarry Smith   PetscViewerFormat format;
911f1af5d2fSBarry Smith 
912f1af5d2fSBarry Smith   PetscFunctionBegin;
913f1af5d2fSBarry Smith 
914f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
915b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
916b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
917f1af5d2fSBarry Smith 
918f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
919fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
920f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
921b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
922f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
923f1af5d2fSBarry Smith       x_l = j;
924f1af5d2fSBarry Smith       x_r = x_l + 1.0;
925f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
926f1af5d2fSBarry Smith         y_l = m - i - 1.0;
927f1af5d2fSBarry Smith         y_r = y_l + 1.0;
928f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
929329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
930b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
931329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
932b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
933f1af5d2fSBarry Smith         } else {
934f1af5d2fSBarry Smith           continue;
935f1af5d2fSBarry Smith         }
936f1af5d2fSBarry Smith #else
937f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
938b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
939f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
940b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
941f1af5d2fSBarry Smith         } else {
942f1af5d2fSBarry Smith           continue;
943f1af5d2fSBarry Smith         }
944f1af5d2fSBarry Smith #endif
945b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
946f1af5d2fSBarry Smith       }
947f1af5d2fSBarry Smith     }
948f1af5d2fSBarry Smith   } else {
949f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
950f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
951f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
952f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
953f1af5d2fSBarry Smith     }
954b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
955b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
956b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
957f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
958f1af5d2fSBarry Smith       x_l = j;
959f1af5d2fSBarry Smith       x_r = x_l + 1.0;
960f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
961f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
962f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
963b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
964b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
965f1af5d2fSBarry Smith       }
966f1af5d2fSBarry Smith     }
967f1af5d2fSBarry Smith   }
968f1af5d2fSBarry Smith   PetscFunctionReturn(0);
969f1af5d2fSBarry Smith }
970f1af5d2fSBarry Smith 
9714a2ae208SSatish Balay #undef __FUNCT__
9724a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
973dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
974f1af5d2fSBarry Smith {
975b0a32e0cSBarry Smith   PetscDraw      draw;
976f1af5d2fSBarry Smith   PetscTruth     isnull;
977329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
978dfbe8321SBarry Smith   PetscErrorCode ierr;
979f1af5d2fSBarry Smith 
980f1af5d2fSBarry Smith   PetscFunctionBegin;
981b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
982b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
983abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
984f1af5d2fSBarry Smith 
985f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
986273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
987f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
988b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
989b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
990f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
991f1af5d2fSBarry Smith   PetscFunctionReturn(0);
992f1af5d2fSBarry Smith }
993f1af5d2fSBarry Smith 
9944a2ae208SSatish Balay #undef __FUNCT__
9954a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
996dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
997932b0c3eSLois Curfman McInnes {
998932b0c3eSLois Curfman McInnes   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
999dfbe8321SBarry Smith   PetscErrorCode ierr;
100032077d6dSBarry Smith   PetscTruth     issocket,iascii,isbinary,isdraw;
1001932b0c3eSLois Curfman McInnes 
10023a40ed3dSBarry Smith   PetscFunctionBegin;
1003b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
100432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1005fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1006fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10070f5bd95cSBarry Smith 
10080f5bd95cSBarry Smith   if (issocket) {
1009634064b4SBarry Smith     if (a->lda>A->m) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA");
1010b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
101132077d6dSBarry Smith   } else if (iascii) {
10123a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10130f5bd95cSBarry Smith   } else if (isbinary) {
10143a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1015f1af5d2fSBarry Smith   } else if (isdraw) {
1016f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10175cd90555SBarry Smith   } else {
1018958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1019932b0c3eSLois Curfman McInnes   }
10203a40ed3dSBarry Smith   PetscFunctionReturn(0);
1021932b0c3eSLois Curfman McInnes }
1022289bc588SBarry Smith 
10234a2ae208SSatish Balay #undef __FUNCT__
10244a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1025dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1026289bc588SBarry Smith {
1027ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1028dfbe8321SBarry Smith   PetscErrorCode ierr;
102990f02eecSBarry Smith 
10303a40ed3dSBarry Smith   PetscFunctionBegin;
1031aa482453SBarry Smith #if defined(PETSC_USE_LOG)
103277431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->m,mat->n);
1033a5a9c739SBarry Smith #endif
1034606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1035606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1036606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1037901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10383a40ed3dSBarry Smith   PetscFunctionReturn(0);
1039289bc588SBarry Smith }
1040289bc588SBarry Smith 
10414a2ae208SSatish Balay #undef __FUNCT__
10424a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1043dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1044289bc588SBarry Smith {
1045c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
10466849ba73SBarry Smith   PetscErrorCode ierr;
104713f74950SBarry Smith   PetscInt       k,j,m,n,M;
104887828ca2SBarry Smith   PetscScalar    *v,tmp;
104948b35521SBarry Smith 
10503a40ed3dSBarry Smith   PetscFunctionBegin;
10511b807ce4Svictorle   v = mat->v; m = A->m; M = mat->lda; n = A->n;
10527c922b88SBarry Smith   if (!matout) { /* in place transpose */
1053a5ce6ee0Svictorle     if (m != n) {
1054634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1055d64ed03dSBarry Smith     } else {
1056d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1057289bc588SBarry Smith         for (k=0; k<j; k++) {
10581b807ce4Svictorle           tmp = v[j + k*M];
10591b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10601b807ce4Svictorle           v[k + j*M] = tmp;
1061289bc588SBarry Smith         }
1062289bc588SBarry Smith       }
1063d64ed03dSBarry Smith     }
10643a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1065d3e5ee88SLois Curfman McInnes     Mat          tmat;
1066ec8511deSBarry Smith     Mat_SeqDense *tmatd;
106787828ca2SBarry Smith     PetscScalar  *v2;
1068ea709b57SSatish Balay 
10695c5985e7SKris Buschelman     ierr  = MatCreate(A->comm,A->n,A->m,A->n,A->m,&tmat);CHKERRQ(ierr);
10705c5985e7SKris Buschelman     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
10715c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1072ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10730de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1074d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10751b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1076d3e5ee88SLois Curfman McInnes     }
10776d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10786d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1079d3e5ee88SLois Curfman McInnes     *matout = tmat;
108048b35521SBarry Smith   }
10813a40ed3dSBarry Smith   PetscFunctionReturn(0);
1082289bc588SBarry Smith }
1083289bc588SBarry Smith 
10844a2ae208SSatish Balay #undef __FUNCT__
10854a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1086dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1087289bc588SBarry Smith {
1088c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1089c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
109013f74950SBarry Smith   PetscInt     i,j;
109187828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
10929ea5d5aeSSatish Balay 
10933a40ed3dSBarry Smith   PetscFunctionBegin;
1094273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1095273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10961b807ce4Svictorle   for (i=0; i<A1->m; i++) {
10971b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
10981b807ce4Svictorle     for (j=0; j<A1->n; j++) {
10993a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11001b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11011b807ce4Svictorle     }
1102289bc588SBarry Smith   }
110377c4ece6SBarry Smith   *flg = PETSC_TRUE;
11043a40ed3dSBarry Smith   PetscFunctionReturn(0);
1105289bc588SBarry Smith }
1106289bc588SBarry Smith 
11074a2ae208SSatish Balay #undef __FUNCT__
11084a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1109dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1110289bc588SBarry Smith {
1111c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1112dfbe8321SBarry Smith   PetscErrorCode ierr;
111313f74950SBarry Smith   PetscInt       i,n,len;
111487828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
111544cd7ae7SLois Curfman McInnes 
11163a40ed3dSBarry Smith   PetscFunctionBegin;
11177a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
11187a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11191ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1120273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1121273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
112244cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11231b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1124289bc588SBarry Smith   }
11251ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11263a40ed3dSBarry Smith   PetscFunctionReturn(0);
1127289bc588SBarry Smith }
1128289bc588SBarry Smith 
11294a2ae208SSatish Balay #undef __FUNCT__
11304a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1131dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1132289bc588SBarry Smith {
1133c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
113487828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1135dfbe8321SBarry Smith   PetscErrorCode ierr;
113613f74950SBarry Smith   PetscInt       i,j,m = A->m,n = A->n;
113755659b69SBarry Smith 
11383a40ed3dSBarry Smith   PetscFunctionBegin;
113928988994SBarry Smith   if (ll) {
11407a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11411ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1142273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1143da3a660dSBarry Smith     for (i=0; i<m; i++) {
1144da3a660dSBarry Smith       x = l[i];
1145da3a660dSBarry Smith       v = mat->v + i;
1146da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1147da3a660dSBarry Smith     }
11481ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1149efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1150da3a660dSBarry Smith   }
115128988994SBarry Smith   if (rr) {
11527a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11531ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1154273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1155da3a660dSBarry Smith     for (i=0; i<n; i++) {
1156da3a660dSBarry Smith       x = r[i];
1157da3a660dSBarry Smith       v = mat->v + i*m;
1158da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1159da3a660dSBarry Smith     }
11601ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1161efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1162da3a660dSBarry Smith   }
11633a40ed3dSBarry Smith   PetscFunctionReturn(0);
1164289bc588SBarry Smith }
1165289bc588SBarry Smith 
11664a2ae208SSatish Balay #undef __FUNCT__
11674a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1168dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1169289bc588SBarry Smith {
1170c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
117187828ca2SBarry Smith   PetscScalar  *v = mat->v;
1172329f5518SBarry Smith   PetscReal    sum = 0.0;
117313f74950SBarry Smith   PetscInt     lda=mat->lda,m=A->m,i,j;
1174efee365bSSatish Balay   PetscErrorCode ierr;
117555659b69SBarry Smith 
11763a40ed3dSBarry Smith   PetscFunctionBegin;
1177289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1178a5ce6ee0Svictorle     if (lda>m) {
1179a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
1180a5ce6ee0Svictorle 	v = mat->v+j*lda;
1181a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1182a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1183a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1184a5ce6ee0Svictorle #else
1185a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1186a5ce6ee0Svictorle #endif
1187a5ce6ee0Svictorle 	}
1188a5ce6ee0Svictorle       }
1189a5ce6ee0Svictorle     } else {
1190273d9f13SBarry Smith       for (i=0; i<A->n*A->m; i++) {
1191aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1192329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1193289bc588SBarry Smith #else
1194289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1195289bc588SBarry Smith #endif
1196289bc588SBarry Smith       }
1197a5ce6ee0Svictorle     }
1198064f8208SBarry Smith     *nrm = sqrt(sum);
1199efee365bSSatish Balay     ierr = PetscLogFlops(2*A->n*A->m);CHKERRQ(ierr);
12003a40ed3dSBarry Smith   } else if (type == NORM_1) {
1201064f8208SBarry Smith     *nrm = 0.0;
1202273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
12031b807ce4Svictorle       v = mat->v + j*mat->lda;
1204289bc588SBarry Smith       sum = 0.0;
1205273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
120633a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1207289bc588SBarry Smith       }
1208064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1209289bc588SBarry Smith     }
1210efee365bSSatish Balay     ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
12113a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1212064f8208SBarry Smith     *nrm = 0.0;
1213273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1214289bc588SBarry Smith       v = mat->v + j;
1215289bc588SBarry Smith       sum = 0.0;
1216273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
12171b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1218289bc588SBarry Smith       }
1219064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1220289bc588SBarry Smith     }
1221efee365bSSatish Balay     ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
12223a40ed3dSBarry Smith   } else {
122329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1224289bc588SBarry Smith   }
12253a40ed3dSBarry Smith   PetscFunctionReturn(0);
1226289bc588SBarry Smith }
1227289bc588SBarry Smith 
12284a2ae208SSatish Balay #undef __FUNCT__
12294a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1230dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op)
1231289bc588SBarry Smith {
1232c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
123363ba0a88SBarry Smith   PetscErrorCode ierr;
123467e560aaSBarry Smith 
12353a40ed3dSBarry Smith   PetscFunctionBegin;
1236b5a2b587SKris Buschelman   switch (op) {
1237b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1238b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1239b5a2b587SKris Buschelman     break;
1240b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1241b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1242b5a2b587SKris Buschelman     break;
1243b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1244b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1245b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1246b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1247b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1248b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1249b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1250b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1251b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1252b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1253b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
125463ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_SeqDense:Option ignored\n"));CHKERRQ(ierr);
1255b5a2b587SKris Buschelman     break;
125677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
125777e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12589a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12599a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12609a4540c5SBarry Smith   case MAT_HERMITIAN:
12619a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12629a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12639a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
126477e54ba9SKris Buschelman     break;
1265b5a2b587SKris Buschelman   default:
126629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
12673a40ed3dSBarry Smith   }
12683a40ed3dSBarry Smith   PetscFunctionReturn(0);
1269289bc588SBarry Smith }
1270289bc588SBarry Smith 
12714a2ae208SSatish Balay #undef __FUNCT__
12724a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1273dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
12746f0a148fSBarry Smith {
1275ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12766849ba73SBarry Smith   PetscErrorCode ierr;
127713f74950SBarry Smith   PetscInt       lda=l->lda,m=A->m,j;
12783a40ed3dSBarry Smith 
12793a40ed3dSBarry Smith   PetscFunctionBegin;
1280a5ce6ee0Svictorle   if (lda>m) {
1281a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
1282a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1283a5ce6ee0Svictorle     }
1284a5ce6ee0Svictorle   } else {
128587828ca2SBarry Smith     ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1286a5ce6ee0Svictorle   }
12873a40ed3dSBarry Smith   PetscFunctionReturn(0);
12886f0a148fSBarry Smith }
12896f0a148fSBarry Smith 
12904a2ae208SSatish Balay #undef __FUNCT__
12914a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1292dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag)
12936f0a148fSBarry Smith {
1294ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12956849ba73SBarry Smith   PetscErrorCode ierr;
129613f74950SBarry Smith   PetscInt       n = A->n,i,j,N,*rows;
129787828ca2SBarry Smith   PetscScalar    *slot;
129855659b69SBarry Smith 
12993a40ed3dSBarry Smith   PetscFunctionBegin;
1300e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
130178b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
13026f0a148fSBarry Smith   for (i=0; i<N; i++) {
13036f0a148fSBarry Smith     slot = l->v + rows[i];
13046f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13056f0a148fSBarry Smith   }
13066f0a148fSBarry Smith   if (diag) {
13076f0a148fSBarry Smith     for (i=0; i<N; i++) {
13086f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1309260925b8SBarry Smith       *slot = *diag;
13106f0a148fSBarry Smith     }
13116f0a148fSBarry Smith   }
1312260925b8SBarry Smith   ISRestoreIndices(is,&rows);
13133a40ed3dSBarry Smith   PetscFunctionReturn(0);
13146f0a148fSBarry Smith }
1315557bce09SLois Curfman McInnes 
13164a2ae208SSatish Balay #undef __FUNCT__
13174a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1318dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
131964e87e97SBarry Smith {
1320c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13213a40ed3dSBarry Smith 
13223a40ed3dSBarry Smith   PetscFunctionBegin;
132364e87e97SBarry Smith   *array = mat->v;
13243a40ed3dSBarry Smith   PetscFunctionReturn(0);
132564e87e97SBarry Smith }
13260754003eSLois Curfman McInnes 
13274a2ae208SSatish Balay #undef __FUNCT__
13284a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1329dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1330ff14e315SSatish Balay {
13313a40ed3dSBarry Smith   PetscFunctionBegin;
133209b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1334ff14e315SSatish Balay }
13350754003eSLois Curfman McInnes 
13364a2ae208SSatish Balay #undef __FUNCT__
13374a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
133813f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13390754003eSLois Curfman McInnes {
1340c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13416849ba73SBarry Smith   PetscErrorCode ierr;
134213f74950SBarry Smith   PetscInt       i,j,m = A->m,*irow,*icol,nrows,ncols;
134387828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13440754003eSLois Curfman McInnes   Mat            newmat;
13450754003eSLois Curfman McInnes 
13463a40ed3dSBarry Smith   PetscFunctionBegin;
134778b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
134878b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1349e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1350e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13510754003eSLois Curfman McInnes 
1352182d2002SSatish Balay   /* Check submatrixcall */
1353182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
135413f74950SBarry Smith     PetscInt n_cols,n_rows;
1355182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
135629bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1357182d2002SSatish Balay     newmat = *B;
1358182d2002SSatish Balay   } else {
13590754003eSLois Curfman McInnes     /* Create and fill new matrix */
13605c5985e7SKris Buschelman     ierr = MatCreate(A->comm,nrows,ncols,nrows,ncols,&newmat);CHKERRQ(ierr);
13615c5985e7SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
13625c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1363182d2002SSatish Balay   }
1364182d2002SSatish Balay 
1365182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1366182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1367182d2002SSatish Balay 
1368182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1369182d2002SSatish Balay     av = v + m*icol[i];
1370182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1371182d2002SSatish Balay       *bv++ = av[irow[j]];
13720754003eSLois Curfman McInnes     }
13730754003eSLois Curfman McInnes   }
1374182d2002SSatish Balay 
1375182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13766d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13776d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13780754003eSLois Curfman McInnes 
13790754003eSLois Curfman McInnes   /* Free work space */
138078b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
138178b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1382182d2002SSatish Balay   *B = newmat;
13833a40ed3dSBarry Smith   PetscFunctionReturn(0);
13840754003eSLois Curfman McInnes }
13850754003eSLois Curfman McInnes 
13864a2ae208SSatish Balay #undef __FUNCT__
13874a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
138813f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1389905e6a2fSBarry Smith {
13906849ba73SBarry Smith   PetscErrorCode ierr;
139113f74950SBarry Smith   PetscInt       i;
1392905e6a2fSBarry Smith 
13933a40ed3dSBarry Smith   PetscFunctionBegin;
1394905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1395b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1396905e6a2fSBarry Smith   }
1397905e6a2fSBarry Smith 
1398905e6a2fSBarry Smith   for (i=0; i<n; i++) {
13996a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1400905e6a2fSBarry Smith   }
14013a40ed3dSBarry Smith   PetscFunctionReturn(0);
1402905e6a2fSBarry Smith }
1403905e6a2fSBarry Smith 
14044a2ae208SSatish Balay #undef __FUNCT__
14054a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1406dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14074b0e389bSBarry Smith {
14084b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14096849ba73SBarry Smith   PetscErrorCode ierr;
141013f74950SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j;
14113a40ed3dSBarry Smith 
14123a40ed3dSBarry Smith   PetscFunctionBegin;
141333f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
141433f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1415cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14163a40ed3dSBarry Smith     PetscFunctionReturn(0);
14173a40ed3dSBarry Smith   }
14180dbb7854Svictorle   if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1419a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14200dbb7854Svictorle     for (j=0; j<n; j++) {
1421a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1422a5ce6ee0Svictorle     }
1423a5ce6ee0Svictorle   } else {
142487828ca2SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1425a5ce6ee0Svictorle   }
1426273d9f13SBarry Smith   PetscFunctionReturn(0);
1427273d9f13SBarry Smith }
1428273d9f13SBarry Smith 
14294a2ae208SSatish Balay #undef __FUNCT__
14304a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1431dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1432273d9f13SBarry Smith {
1433dfbe8321SBarry Smith   PetscErrorCode ierr;
1434273d9f13SBarry Smith 
1435273d9f13SBarry Smith   PetscFunctionBegin;
1436273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14373a40ed3dSBarry Smith   PetscFunctionReturn(0);
14384b0e389bSBarry Smith }
14394b0e389bSBarry Smith 
1440289bc588SBarry Smith /* -------------------------------------------------------------------*/
1441a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1442905e6a2fSBarry Smith        MatGetRow_SeqDense,
1443905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1444905e6a2fSBarry Smith        MatMult_SeqDense,
144597304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
14467c922b88SBarry Smith        MatMultTranspose_SeqDense,
14477c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1448905e6a2fSBarry Smith        MatSolve_SeqDense,
1449905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
14507c922b88SBarry Smith        MatSolveTranspose_SeqDense,
145197304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1452905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1453905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1454ec8511deSBarry Smith        MatRelax_SeqDense,
1455ec8511deSBarry Smith        MatTranspose_SeqDense,
145697304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1457905e6a2fSBarry Smith        MatEqual_SeqDense,
1458905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1459905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1460905e6a2fSBarry Smith        MatNorm_SeqDense,
146197304618SKris Buschelman /*20*/ 0,
1462905e6a2fSBarry Smith        0,
1463905e6a2fSBarry Smith        0,
1464905e6a2fSBarry Smith        MatSetOption_SeqDense,
1465905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
146697304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1467905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1468905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1469905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1470905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
147197304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1472273d9f13SBarry Smith        0,
1473905e6a2fSBarry Smith        0,
1474905e6a2fSBarry Smith        MatGetArray_SeqDense,
1475905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
147697304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1477a5ae1ecdSBarry Smith        0,
1478a5ae1ecdSBarry Smith        0,
1479a5ae1ecdSBarry Smith        0,
1480a5ae1ecdSBarry Smith        0,
148197304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1482a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1483a5ae1ecdSBarry Smith        0,
14844b0e389bSBarry Smith        MatGetValues_SeqDense,
1485a5ae1ecdSBarry Smith        MatCopy_SeqDense,
148697304618SKris Buschelman /*45*/ 0,
1487a5ae1ecdSBarry Smith        MatScale_SeqDense,
1488a5ae1ecdSBarry Smith        0,
1489a5ae1ecdSBarry Smith        0,
1490a5ae1ecdSBarry Smith        0,
1491521d7252SBarry Smith /*50*/ 0,
1492a5ae1ecdSBarry Smith        0,
1493a5ae1ecdSBarry Smith        0,
1494a5ae1ecdSBarry Smith        0,
1495a5ae1ecdSBarry Smith        0,
149697304618SKris Buschelman /*55*/ 0,
1497a5ae1ecdSBarry Smith        0,
1498a5ae1ecdSBarry Smith        0,
1499a5ae1ecdSBarry Smith        0,
1500a5ae1ecdSBarry Smith        0,
150197304618SKris Buschelman /*60*/ 0,
1502e03a110bSBarry Smith        MatDestroy_SeqDense,
1503e03a110bSBarry Smith        MatView_SeqDense,
150497304618SKris Buschelman        MatGetPetscMaps_Petsc,
150597304618SKris Buschelman        0,
150697304618SKris Buschelman /*65*/ 0,
150797304618SKris Buschelman        0,
150897304618SKris Buschelman        0,
150997304618SKris Buschelman        0,
151097304618SKris Buschelman        0,
151197304618SKris Buschelman /*70*/ 0,
151297304618SKris Buschelman        0,
151397304618SKris Buschelman        0,
151497304618SKris Buschelman        0,
151597304618SKris Buschelman        0,
151697304618SKris Buschelman /*75*/ 0,
151797304618SKris Buschelman        0,
151897304618SKris Buschelman        0,
151997304618SKris Buschelman        0,
152097304618SKris Buschelman        0,
152197304618SKris Buschelman /*80*/ 0,
152297304618SKris Buschelman        0,
152397304618SKris Buschelman        0,
152497304618SKris Buschelman        0,
1525865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1526865e5f61SKris Buschelman        0,
1527865e5f61SKris Buschelman        0,
1528865e5f61SKris Buschelman        0,
1529865e5f61SKris Buschelman        0,
1530865e5f61SKris Buschelman        0,
1531865e5f61SKris Buschelman /*90*/ 0,
1532865e5f61SKris Buschelman        0,
1533865e5f61SKris Buschelman        0,
1534865e5f61SKris Buschelman        0,
1535865e5f61SKris Buschelman        0,
1536865e5f61SKris Buschelman /*95*/ 0,
1537865e5f61SKris Buschelman        0,
1538865e5f61SKris Buschelman        0,
1539865e5f61SKris Buschelman        0};
154090ace30eSBarry Smith 
15414a2ae208SSatish Balay #undef __FUNCT__
15424a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
15434b828684SBarry Smith /*@C
1544fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1545d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1546d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1547289bc588SBarry Smith 
1548db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1549db81eaa0SLois Curfman McInnes 
155020563c6bSBarry Smith    Input Parameters:
1551db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
15520c775827SLois Curfman McInnes .  m - number of rows
155318f449edSLois Curfman McInnes .  n - number of columns
1554db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1555dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
155620563c6bSBarry Smith 
155720563c6bSBarry Smith    Output Parameter:
155844cd7ae7SLois Curfman McInnes .  A - the matrix
155920563c6bSBarry Smith 
1560b259b22eSLois Curfman McInnes    Notes:
156118f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
156218f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1563b4fd4287SBarry Smith    set data=PETSC_NULL.
156418f449edSLois Curfman McInnes 
1565027ccd11SLois Curfman McInnes    Level: intermediate
1566027ccd11SLois Curfman McInnes 
1567dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1568d65003e9SLois Curfman McInnes 
1569db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
157020563c6bSBarry Smith @*/
1571*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1572289bc588SBarry Smith {
1573dfbe8321SBarry Smith   PetscErrorCode ierr;
15743b2fbd54SBarry Smith 
15753a40ed3dSBarry Smith   PetscFunctionBegin;
1576273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1577273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1578273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1579273d9f13SBarry Smith   PetscFunctionReturn(0);
1580273d9f13SBarry Smith }
1581273d9f13SBarry Smith 
15824a2ae208SSatish Balay #undef __FUNCT__
15834a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1584273d9f13SBarry Smith /*@C
1585273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1586273d9f13SBarry Smith 
1587273d9f13SBarry Smith    Collective on MPI_Comm
1588273d9f13SBarry Smith 
1589273d9f13SBarry Smith    Input Parameters:
1590273d9f13SBarry Smith +  A - the matrix
1591273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1592273d9f13SBarry Smith 
1593273d9f13SBarry Smith    Notes:
1594273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1595273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1596273d9f13SBarry Smith    set data=PETSC_NULL.
1597273d9f13SBarry Smith 
1598273d9f13SBarry Smith    Level: intermediate
1599273d9f13SBarry Smith 
1600273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1601273d9f13SBarry Smith 
1602273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1603273d9f13SBarry Smith @*/
1604*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1605273d9f13SBarry Smith {
1606dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1607a23d5eceSKris Buschelman 
1608a23d5eceSKris Buschelman   PetscFunctionBegin;
1609a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1610a23d5eceSKris Buschelman   if (f) {
1611a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1612a23d5eceSKris Buschelman   }
1613a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1614a23d5eceSKris Buschelman }
1615a23d5eceSKris Buschelman 
1616a23d5eceSKris Buschelman EXTERN_C_BEGIN
1617a23d5eceSKris Buschelman #undef __FUNCT__
1618a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1619*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1620a23d5eceSKris Buschelman {
1621273d9f13SBarry Smith   Mat_SeqDense   *b;
1622dfbe8321SBarry Smith   PetscErrorCode ierr;
1623273d9f13SBarry Smith 
1624273d9f13SBarry Smith   PetscFunctionBegin;
1625273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1626273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1627273d9f13SBarry Smith   if (!data) {
162887828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
162987828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1630273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
163152e6d16bSBarry Smith     ierr = PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));CHKERRQ(ierr);
1632273d9f13SBarry Smith   } else { /* user-allocated storage */
1633273d9f13SBarry Smith     b->v          = data;
1634273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1635273d9f13SBarry Smith   }
1636273d9f13SBarry Smith   PetscFunctionReturn(0);
1637273d9f13SBarry Smith }
1638a23d5eceSKris Buschelman EXTERN_C_END
1639273d9f13SBarry Smith 
16401b807ce4Svictorle #undef __FUNCT__
16411b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
16421b807ce4Svictorle /*@C
16431b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
16441b807ce4Svictorle 
16451b807ce4Svictorle   Input parameter:
16461b807ce4Svictorle + A - the matrix
16471b807ce4Svictorle - lda - the leading dimension
16481b807ce4Svictorle 
16491b807ce4Svictorle   Notes:
16501b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
16511b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
16521b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
16531b807ce4Svictorle 
16541b807ce4Svictorle   Level: intermediate
16551b807ce4Svictorle 
16561b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
16571b807ce4Svictorle 
16581b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation()
16591b807ce4Svictorle @*/
1660*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
16611b807ce4Svictorle {
16621b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
16631b807ce4Svictorle   PetscFunctionBegin;
166477431f27SBarry Smith   if (lda < B->m) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->m);
16651b807ce4Svictorle   b->lda = lda;
16661b807ce4Svictorle   PetscFunctionReturn(0);
16671b807ce4Svictorle }
16681b807ce4Svictorle 
16690bad9183SKris Buschelman /*MC
1670fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
16710bad9183SKris Buschelman 
16720bad9183SKris Buschelman    Options Database Keys:
16730bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
16740bad9183SKris Buschelman 
16750bad9183SKris Buschelman   Level: beginner
16760bad9183SKris Buschelman 
16770bad9183SKris Buschelman .seealso: MatCreateSeqDense
16780bad9183SKris Buschelman M*/
16790bad9183SKris Buschelman 
1680273d9f13SBarry Smith EXTERN_C_BEGIN
16814a2ae208SSatish Balay #undef __FUNCT__
16824a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1683*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
1684273d9f13SBarry Smith {
1685273d9f13SBarry Smith   Mat_SeqDense   *b;
1686dfbe8321SBarry Smith   PetscErrorCode ierr;
16877c334f02SBarry Smith   PetscMPIInt    size;
1688273d9f13SBarry Smith 
1689273d9f13SBarry Smith   PetscFunctionBegin;
1690273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
169129bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
169255659b69SBarry Smith 
1693273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1694273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1695273d9f13SBarry Smith 
1696b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1697549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
169844cd7ae7SLois Curfman McInnes   B->factor       = 0;
169990f02eecSBarry Smith   B->mapping      = 0;
170052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr);
170144cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
170218f449edSLois Curfman McInnes 
17038a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
17048a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1705a5ae1ecdSBarry Smith 
170644cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1707273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1708273d9f13SBarry Smith   b->v            = 0;
17091b807ce4Svictorle   b->lda          = B->m;
17104e220ebcSLois Curfman McInnes 
1711a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1712a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1713a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
17143a40ed3dSBarry Smith   PetscFunctionReturn(0);
1715289bc588SBarry Smith }
1716273d9f13SBarry Smith EXTERN_C_END
1717