xref: /petsc/src/mat/impls/dense/seq/dense.c (revision a07ab388b5b8675ffe5a01dc1bf9e23f3bee111b)
173f4d377SMatthew Knepley /*$Id: dense.c,v 1.208 2001/09/07 20:09:16 bsmith Exp $*/
267e560aaSBarry Smith /*
367e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
467e560aaSBarry Smith */
5289bc588SBarry Smith 
670f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
7b4c862e0SSatish Balay #include "petscblaslapack.h"
8289bc588SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
11268466fbSBarry Smith int MatAXPY_SeqDense(const PetscScalar *alpha,Mat X,Mat Y,MatStructure str)
121987afe7SBarry Smith {
131987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
14a5ce6ee0Svictorle   int          N = X->m*X->n,m=X->m,ldax=x->lda,lday=y->lda, j,one = 1;
153a40ed3dSBarry Smith 
163a40ed3dSBarry Smith   PetscFunctionBegin;
17a5ce6ee0Svictorle   if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
18a5ce6ee0Svictorle   if (ldax>m || lday>m) {
19a5ce6ee0Svictorle     for (j=0; j<X->n; j++) {
20268466fbSBarry Smith       BLaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one);
21a5ce6ee0Svictorle     }
22a5ce6ee0Svictorle   } else {
23268466fbSBarry Smith     BLaxpy_(&N,(PetscScalar*)alpha,x->v,&one,y->v,&one);
24a5ce6ee0Svictorle   }
25b0a32e0cSBarry Smith   PetscLogFlops(2*N-1);
263a40ed3dSBarry Smith   PetscFunctionReturn(0);
271987afe7SBarry Smith }
281987afe7SBarry Smith 
294a2ae208SSatish Balay #undef __FUNCT__
304a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
318f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
32289bc588SBarry Smith {
334e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
34273d9f13SBarry Smith   int          i,N = A->m*A->n,count = 0;
3587828ca2SBarry Smith   PetscScalar  *v = a->v;
363a40ed3dSBarry Smith 
373a40ed3dSBarry Smith   PetscFunctionBegin;
38289bc588SBarry Smith   for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}
394e220ebcSLois Curfman McInnes 
40273d9f13SBarry Smith   info->rows_global       = (double)A->m;
41273d9f13SBarry Smith   info->columns_global    = (double)A->n;
42273d9f13SBarry Smith   info->rows_local        = (double)A->m;
43273d9f13SBarry Smith   info->columns_local     = (double)A->n;
444e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
454e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
464e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
474e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
484e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
494e220ebcSLois Curfman McInnes   info->mallocs           = 0;
504e220ebcSLois Curfman McInnes   info->memory            = A->mem;
514e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
524e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
534e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
544e220ebcSLois Curfman McInnes 
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56289bc588SBarry Smith }
57289bc588SBarry Smith 
584a2ae208SSatish Balay #undef __FUNCT__
594a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
60268466fbSBarry Smith int MatScale_SeqDense(const PetscScalar *alpha,Mat A)
6180cd9d93SLois Curfman McInnes {
62273d9f13SBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
63a5ce6ee0Svictorle   int          one = 1,lda = a->lda,j,nz;
6480cd9d93SLois Curfman McInnes 
653a40ed3dSBarry Smith   PetscFunctionBegin;
66a5ce6ee0Svictorle   if (lda>A->m) {
67a5ce6ee0Svictorle     nz = A->m;
68a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
69268466fbSBarry Smith       BLscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
72273d9f13SBarry Smith     nz = A->m*A->n;
73268466fbSBarry Smith     BLscal_(&nz,(PetscScalar*)alpha,a->v,&one);
74a5ce6ee0Svictorle   }
75b0a32e0cSBarry Smith   PetscLogFlops(nz);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
7780cd9d93SLois Curfman McInnes }
7880cd9d93SLois Curfman McInnes 
79289bc588SBarry Smith /* ---------------------------------------------------------------*/
806831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
816831982aSBarry Smith    rather than put it in the Mat->row slot.*/
824a2ae208SSatish Balay #undef __FUNCT__
834a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
84b380c88cSHong Zhang int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
85289bc588SBarry Smith {
86*a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
87*a07ab388SBarry Smith   PetscFuntionBegin;
88*a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
89*a07ab388SBarry Smith #else
90c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
91b0a32e0cSBarry Smith   int          info,ierr;
9255659b69SBarry Smith 
933a40ed3dSBarry Smith   PetscFunctionBegin;
94289bc588SBarry Smith   if (!mat->pivots) {
95b0a32e0cSBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr);
96b0a32e0cSBarry Smith     PetscLogObjectMemory(A,A->m*sizeof(int));
97289bc588SBarry Smith   }
98f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
99273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
1001b807ce4Svictorle   LAgetrf_(&A->m,&A->n,mat->v,&mat->lda,mat->pivots,&info);
10129bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
10229bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
103b0a32e0cSBarry Smith   PetscLogFlops((2*A->n*A->n*A->n)/3);
104*a07ab388SBarry Smith #endif
1053a40ed3dSBarry Smith   PetscFunctionReturn(0);
106289bc588SBarry Smith }
1076ee01492SSatish Balay 
1084a2ae208SSatish Balay #undef __FUNCT__
1094a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
1105609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11102cad45dSBarry Smith {
11202cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l;
113a5ce6ee0Svictorle   int          lda = mat->lda,j,m,ierr;
11402cad45dSBarry Smith   Mat          newi;
11502cad45dSBarry Smith 
1163a40ed3dSBarry Smith   PetscFunctionBegin;
117273d9f13SBarry Smith   ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr);
1185609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
119a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
120a5ce6ee0Svictorle     if (lda>A->m) {
121a5ce6ee0Svictorle       m = A->m;
122a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
123a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
124a5ce6ee0Svictorle       }
125a5ce6ee0Svictorle     } else {
12687828ca2SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
12702cad45dSBarry Smith     }
128a5ce6ee0Svictorle   }
12902cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
13002cad45dSBarry Smith   *newmat = newi;
1313a40ed3dSBarry Smith   PetscFunctionReturn(0);
13202cad45dSBarry Smith }
13302cad45dSBarry Smith 
1344a2ae208SSatish Balay #undef __FUNCT__
1354a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
136b380c88cSHong Zhang int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
137289bc588SBarry Smith {
1383a40ed3dSBarry Smith   int ierr;
1393a40ed3dSBarry Smith 
1403a40ed3dSBarry Smith   PetscFunctionBegin;
1415609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1423a40ed3dSBarry Smith   PetscFunctionReturn(0);
143289bc588SBarry Smith }
1446ee01492SSatish Balay 
1454a2ae208SSatish Balay #undef __FUNCT__
1464a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1478f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
148289bc588SBarry Smith {
14902cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1501b807ce4Svictorle   int          lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j,ierr;
1513a40ed3dSBarry Smith 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
15302cad45dSBarry Smith   /* copy the numerical values */
1541b807ce4Svictorle   if (lda1>m || lda2>m ) {
1551b807ce4Svictorle     for (j=0; j<n; j++) {
1561b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar)); CHKERRQ(ierr);
1571b807ce4Svictorle     }
1581b807ce4Svictorle   } else {
15987828ca2SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1601b807ce4Svictorle   }
16102cad45dSBarry Smith   (*fact)->factor = 0;
162e03a110bSBarry Smith   ierr = MatLUFactor(*fact,0,0,PETSC_NULL);CHKERRQ(ierr);
1633a40ed3dSBarry Smith   PetscFunctionReturn(0);
164289bc588SBarry Smith }
1656ee01492SSatish Balay 
1664a2ae208SSatish Balay #undef __FUNCT__
1674a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
16815e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
169289bc588SBarry Smith {
1703a40ed3dSBarry Smith   int ierr;
1713a40ed3dSBarry Smith 
1723a40ed3dSBarry Smith   PetscFunctionBegin;
1733a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175289bc588SBarry Smith }
1766ee01492SSatish Balay 
1774a2ae208SSatish Balay #undef __FUNCT__
1784a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
17915e8a5b3SHong Zhang int MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
180289bc588SBarry Smith {
181*a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
182*a07ab388SBarry Smith   PetscFunctionBegin;
183*a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
184*a07ab388SBarry Smith #else
185c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
186606d414cSSatish Balay   int           info,ierr;
18755659b69SBarry Smith 
1883a40ed3dSBarry Smith   PetscFunctionBegin;
189752f5567SLois Curfman McInnes   if (mat->pivots) {
190606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
191b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-A->m*sizeof(int));
192752f5567SLois Curfman McInnes     mat->pivots = 0;
193752f5567SLois Curfman McInnes   }
194273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
1951b807ce4Svictorle   LApotrf_("L",&A->n,mat->v,&mat->lda,&info);
19629bbc08cSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
197c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
198b0a32e0cSBarry Smith   PetscLogFlops((A->n*A->n*A->n)/3);
199*a07ab388SBarry Smith #endif
2003a40ed3dSBarry Smith   PetscFunctionReturn(0);
201289bc588SBarry Smith }
202289bc588SBarry Smith 
2034a2ae208SSatish Balay #undef __FUNCT__
2040b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
2050b4b3355SBarry Smith int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
2060b4b3355SBarry Smith {
2070b4b3355SBarry Smith   int ierr;
20815e8a5b3SHong Zhang   MatFactorInfo info;
2090b4b3355SBarry Smith 
2100b4b3355SBarry Smith   PetscFunctionBegin;
21115e8a5b3SHong Zhang   info.fill = 1.0;
21215e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2130b4b3355SBarry Smith   PetscFunctionReturn(0);
2140b4b3355SBarry Smith }
2150b4b3355SBarry Smith 
2160b4b3355SBarry Smith #undef __FUNCT__
2174a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
2188f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
219289bc588SBarry Smith {
220c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
221a57ad546SLois Curfman McInnes   int          one = 1,info,ierr;
22287828ca2SBarry Smith   PetscScalar  *x,*y;
22367e560aaSBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionBegin;
225273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
226b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
227b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
22887828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
229c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
230ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
23180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
232ae7cfcebSSatish Balay #else
2331b807ce4Svictorle     LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
2342ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
235ae7cfcebSSatish Balay #endif
2367a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
237ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
23880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
239ae7cfcebSSatish Balay #else
2401b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
2412ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
242ae7cfcebSSatish Balay #endif
243289bc588SBarry Smith   }
24429bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
245b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
246b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
247b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
249289bc588SBarry Smith }
2506ee01492SSatish Balay 
2514a2ae208SSatish Balay #undef __FUNCT__
2524a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
2537c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
254da3a660dSBarry Smith {
255c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2567a97a34bSBarry Smith   int          ierr,one = 1,info;
25787828ca2SBarry Smith   PetscScalar  *x,*y;
25867e560aaSBarry Smith 
2593a40ed3dSBarry Smith   PetscFunctionBegin;
260273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
261b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
262b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
26387828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
264752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
265da3a660dSBarry Smith   if (mat->pivots) {
266ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
26780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
268ae7cfcebSSatish Balay #else
2691b807ce4Svictorle     LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
2702ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
271ae7cfcebSSatish Balay #endif
2727a97a34bSBarry Smith   } else {
273ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
27480444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
275ae7cfcebSSatish Balay #else
2761b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
2772ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
278ae7cfcebSSatish Balay #endif
279da3a660dSBarry Smith   }
280b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
281b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
282b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2833a40ed3dSBarry Smith   PetscFunctionReturn(0);
284da3a660dSBarry Smith }
2856ee01492SSatish Balay 
2864a2ae208SSatish Balay #undef __FUNCT__
2874a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
2888f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
289da3a660dSBarry Smith {
290c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2917a97a34bSBarry Smith   int          ierr,one = 1,info;
29287828ca2SBarry Smith   PetscScalar  *x,*y,sone = 1.0;
293da3a660dSBarry Smith   Vec          tmp = 0;
29467e560aaSBarry Smith 
2953a40ed3dSBarry Smith   PetscFunctionBegin;
296b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
297b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
298273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
299da3a660dSBarry Smith   if (yy == zz) {
30078b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
301b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
30278b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
303da3a660dSBarry Smith   }
30487828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
305752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
306da3a660dSBarry Smith   if (mat->pivots) {
307ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
30880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
309ae7cfcebSSatish Balay #else
3101b807ce4Svictorle     LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
3112ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
312ae7cfcebSSatish Balay #endif
313a8c6a408SBarry Smith   } else {
314ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
31580444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
316ae7cfcebSSatish Balay #else
3171b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
3182ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
319ae7cfcebSSatish Balay #endif
320da3a660dSBarry Smith   }
321a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
322a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
323b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
324b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
325b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3263a40ed3dSBarry Smith   PetscFunctionReturn(0);
327da3a660dSBarry Smith }
32867e560aaSBarry Smith 
3294a2ae208SSatish Balay #undef __FUNCT__
3304a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
3317c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
332da3a660dSBarry Smith {
333c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
3346abc6512SBarry Smith   int           one = 1,info,ierr;
33587828ca2SBarry Smith   PetscScalar   *x,*y,sone = 1.0;
336da3a660dSBarry Smith   Vec           tmp;
33767e560aaSBarry Smith 
3383a40ed3dSBarry Smith   PetscFunctionBegin;
339273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
340b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
341b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
342da3a660dSBarry Smith   if (yy == zz) {
34378b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
344b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
34578b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
346da3a660dSBarry Smith   }
34787828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
348752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
349da3a660dSBarry Smith   if (mat->pivots) {
350ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
35180444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
352ae7cfcebSSatish Balay #else
3531b807ce4Svictorle     LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
3542ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
355ae7cfcebSSatish Balay #endif
3563a40ed3dSBarry Smith   } else {
357ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
35880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
359ae7cfcebSSatish Balay #else
3601b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
3612ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
362ae7cfcebSSatish Balay #endif
363da3a660dSBarry Smith   }
36490f02eecSBarry Smith   if (tmp) {
36590f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
36690f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3673a40ed3dSBarry Smith   } else {
36890f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
36990f02eecSBarry Smith   }
370b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
371b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
372b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3733a40ed3dSBarry Smith   PetscFunctionReturn(0);
374da3a660dSBarry Smith }
375289bc588SBarry Smith /* ------------------------------------------------------------------*/
3764a2ae208SSatish Balay #undef __FUNCT__
3774a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
378329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
379c14dc6b6SHong Zhang                           PetscReal shift,int its,int lits,Vec xx)
380289bc588SBarry Smith {
381c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
38287828ca2SBarry Smith   PetscScalar  *x,*b,*v = mat->v,zero = 0.0,xt;
383273d9f13SBarry Smith   int          ierr,m = A->m,i;
384aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
385bc1b551cSSatish Balay   int          o = 1;
386bc1b551cSSatish Balay #endif
387289bc588SBarry Smith 
3883a40ed3dSBarry Smith   PetscFunctionBegin;
389289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
3903a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
391bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
392289bc588SBarry Smith   }
393b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
394b1d4fb26SBarry Smith   ierr = VecGetArrayFast(bb,&b);CHKERRQ(ierr);
395b965ef7fSBarry Smith   its  = its*lits;
39691723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
397289bc588SBarry Smith   while (its--) {
398289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
399289bc588SBarry Smith       for (i=0; i<m; i++) {
400aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
401f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
402f1747703SBarry Smith            not happy about returning a double complex */
403f1747703SBarry Smith         int         _i;
40487828ca2SBarry Smith         PetscScalar sum = b[i];
405f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4063f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
407f1747703SBarry Smith         }
408f1747703SBarry Smith         xt = sum;
409f1747703SBarry Smith #else
410289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
411f1747703SBarry Smith #endif
41255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
413289bc588SBarry Smith       }
414289bc588SBarry Smith     }
415289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
416289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
417aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
418f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
419f1747703SBarry Smith            not happy about returning a double complex */
420f1747703SBarry Smith         int         _i;
42187828ca2SBarry Smith         PetscScalar sum = b[i];
422f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4233f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
424f1747703SBarry Smith         }
425f1747703SBarry Smith         xt = sum;
426f1747703SBarry Smith #else
427289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
428f1747703SBarry Smith #endif
42955a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
430289bc588SBarry Smith       }
431289bc588SBarry Smith     }
432289bc588SBarry Smith   }
433b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(bb,&b);CHKERRQ(ierr);
434b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
4353a40ed3dSBarry Smith   PetscFunctionReturn(0);
436289bc588SBarry Smith }
437289bc588SBarry Smith 
438289bc588SBarry Smith /* -----------------------------------------------------------------*/
4394a2ae208SSatish Balay #undef __FUNCT__
4404a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
4417c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
442289bc588SBarry Smith {
443c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
44487828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
445ea709b57SSatish Balay   int          ierr,_One=1;
446ea709b57SSatish Balay   PetscScalar  _DOne=1.0,_DZero=0.0;
4473a40ed3dSBarry Smith 
4483a40ed3dSBarry Smith   PetscFunctionBegin;
449273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
450b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
451b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
4521b807ce4Svictorle   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
453b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
454b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
455b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->n);
4563a40ed3dSBarry Smith   PetscFunctionReturn(0);
457289bc588SBarry Smith }
4586ee01492SSatish Balay 
4594a2ae208SSatish Balay #undef __FUNCT__
4604a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
46144cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
462289bc588SBarry Smith {
463c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
46487828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
465329f5518SBarry Smith   int          ierr,_One=1;
4663a40ed3dSBarry Smith 
4673a40ed3dSBarry Smith   PetscFunctionBegin;
468273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
469b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
470b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
4711b807ce4Svictorle   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
472b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
473b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
474b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->m);
4753a40ed3dSBarry Smith   PetscFunctionReturn(0);
476289bc588SBarry Smith }
4776ee01492SSatish Balay 
4784a2ae208SSatish Balay #undef __FUNCT__
4794a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
48044cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
481289bc588SBarry Smith {
482c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
48387828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0;
484329f5518SBarry Smith   int          ierr,_One=1;
4853a40ed3dSBarry Smith 
4863a40ed3dSBarry Smith   PetscFunctionBegin;
487273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
488600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
489b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
490b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
4911b807ce4Svictorle   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
492b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
493b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
494b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4953a40ed3dSBarry Smith   PetscFunctionReturn(0);
496289bc588SBarry Smith }
4976ee01492SSatish Balay 
4984a2ae208SSatish Balay #undef __FUNCT__
4994a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
5007c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
501289bc588SBarry Smith {
502c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
50387828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
5047a97a34bSBarry Smith   int          ierr,_One=1;
50587828ca2SBarry Smith   PetscScalar  _DOne=1.0;
5063a40ed3dSBarry Smith 
5073a40ed3dSBarry Smith   PetscFunctionBegin;
508273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
509600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
510b1d4fb26SBarry Smith   ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr);
511b1d4fb26SBarry Smith   ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr);
5121b807ce4Svictorle   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
513b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr);
514b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr);
515b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
5163a40ed3dSBarry Smith   PetscFunctionReturn(0);
517289bc588SBarry Smith }
518289bc588SBarry Smith 
519289bc588SBarry Smith /* -----------------------------------------------------------------*/
5204a2ae208SSatish Balay #undef __FUNCT__
5214a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
52287828ca2SBarry Smith int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
523289bc588SBarry Smith {
524c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
52587828ca2SBarry Smith   PetscScalar  *v;
526b0a32e0cSBarry Smith   int          i,ierr;
52767e560aaSBarry Smith 
5283a40ed3dSBarry Smith   PetscFunctionBegin;
529273d9f13SBarry Smith   *ncols = A->n;
530289bc588SBarry Smith   if (cols) {
531b0a32e0cSBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr);
532273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
533289bc588SBarry Smith   }
534289bc588SBarry Smith   if (vals) {
53587828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
536289bc588SBarry Smith     v    = mat->v + row;
5371b807ce4Svictorle     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;}
538289bc588SBarry Smith   }
5393a40ed3dSBarry Smith   PetscFunctionReturn(0);
540289bc588SBarry Smith }
5416ee01492SSatish Balay 
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
54487828ca2SBarry Smith int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
545289bc588SBarry Smith {
546606d414cSSatish Balay   int ierr;
547606d414cSSatish Balay   PetscFunctionBegin;
548606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
549606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5503a40ed3dSBarry Smith   PetscFunctionReturn(0);
551289bc588SBarry Smith }
552289bc588SBarry Smith /* ----------------------------------------------------------------*/
5534a2ae208SSatish Balay #undef __FUNCT__
5544a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
555f15d580aSBarry Smith int MatSetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],const PetscScalar v[],InsertMode addv)
556289bc588SBarry Smith {
557c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
558cddbea37SSatish Balay   int          i,j,idx=0;
559d6dfbf8fSBarry Smith 
5603a40ed3dSBarry Smith   PetscFunctionBegin;
561289bc588SBarry Smith   if (!mat->roworiented) {
562dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
563289bc588SBarry Smith       for (j=0; j<n; j++) {
564cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
565aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
566590ac198SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1);
56758804f6dSBarry Smith #endif
568289bc588SBarry Smith         for (i=0; i<m; i++) {
569cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
570aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
5715c8e7597SSatish Balay           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1);
57258804f6dSBarry Smith #endif
573cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
574289bc588SBarry Smith         }
575289bc588SBarry Smith       }
5763a40ed3dSBarry Smith     } else {
577289bc588SBarry Smith       for (j=0; j<n; j++) {
578cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
579aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
580590ac198SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1);
58158804f6dSBarry Smith #endif
582289bc588SBarry Smith         for (i=0; i<m; i++) {
583cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
584aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
5855c8e7597SSatish Balay           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1);
58658804f6dSBarry Smith #endif
587cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
588289bc588SBarry Smith         }
589289bc588SBarry Smith       }
590289bc588SBarry Smith     }
5913a40ed3dSBarry Smith   } else {
592dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
593e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
594cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
595aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
5965c8e7597SSatish Balay         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1);
59758804f6dSBarry Smith #endif
598e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
599cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
600aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
6015c8e7597SSatish Balay           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1);
60258804f6dSBarry Smith #endif
603cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
604e8d4e0b9SBarry Smith         }
605e8d4e0b9SBarry Smith       }
6063a40ed3dSBarry Smith     } else {
607289bc588SBarry Smith       for (i=0; i<m; i++) {
608cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
609aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
6105c8e7597SSatish Balay         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1);
61158804f6dSBarry Smith #endif
612289bc588SBarry Smith         for (j=0; j<n; j++) {
613cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
614aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
6155c8e7597SSatish Balay           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1);
61658804f6dSBarry Smith #endif
617cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
618289bc588SBarry Smith         }
619289bc588SBarry Smith       }
620289bc588SBarry Smith     }
621e8d4e0b9SBarry Smith   }
6223a40ed3dSBarry Smith   PetscFunctionReturn(0);
623289bc588SBarry Smith }
624e8d4e0b9SBarry Smith 
6254a2ae208SSatish Balay #undef __FUNCT__
6264a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
627f15d580aSBarry Smith int MatGetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],PetscScalar v[])
628ae80bb75SLois Curfman McInnes {
629ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
630ae80bb75SLois Curfman McInnes   int          i,j;
63187828ca2SBarry Smith   PetscScalar  *vpt = v;
632ae80bb75SLois Curfman McInnes 
6333a40ed3dSBarry Smith   PetscFunctionBegin;
634ae80bb75SLois Curfman McInnes   /* row-oriented output */
635ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
636ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6371b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
638ae80bb75SLois Curfman McInnes     }
639ae80bb75SLois Curfman McInnes   }
6403a40ed3dSBarry Smith   PetscFunctionReturn(0);
641ae80bb75SLois Curfman McInnes }
642ae80bb75SLois Curfman McInnes 
643289bc588SBarry Smith /* -----------------------------------------------------------------*/
644289bc588SBarry Smith 
645e090d566SSatish Balay #include "petscsys.h"
64688e32edaSLois Curfman McInnes 
6474a2ae208SSatish Balay #undef __FUNCT__
6484a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
6498e9aea5cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,const MatType type,Mat *A)
65088e32edaSLois Curfman McInnes {
65188e32edaSLois Curfman McInnes   Mat_SeqDense *a;
65288e32edaSLois Curfman McInnes   Mat          B;
65388e32edaSLois Curfman McInnes   int          *scols,i,j,nz,ierr,fd,header[4],size;
65488e32edaSLois Curfman McInnes   int          *rowlengths = 0,M,N,*cols;
65587828ca2SBarry Smith   PetscScalar  *vals,*svals,*v,*w;
65619bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
65788e32edaSLois Curfman McInnes 
6583a40ed3dSBarry Smith   PetscFunctionBegin;
659d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
66029bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
661b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6620752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
663552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
66488e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
66588e32edaSLois Curfman McInnes 
66690ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
66790ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
66890ace30eSBarry Smith     B    = *A;
66990ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6704a41a4d3SSatish Balay     v    = a->v;
6714a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6724a41a4d3SSatish Balay        from row major to column major */
673b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
67490ace30eSBarry Smith     /* read in nonzero values */
6754a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6764a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6774a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6784a41a4d3SSatish Balay       for (i=0; i<M; i++) {
6794a41a4d3SSatish Balay         *v++ =w[i*N+j];
6804a41a4d3SSatish Balay       }
6814a41a4d3SSatish Balay     }
682606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
6836d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6846d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68590ace30eSBarry Smith   } else {
68688e32edaSLois Curfman McInnes     /* read row lengths */
687b0a32e0cSBarry Smith     ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr);
6880752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
68988e32edaSLois Curfman McInnes 
69088e32edaSLois Curfman McInnes     /* create our matrix */
691b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
69288e32edaSLois Curfman McInnes     B = *A;
69388e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
69488e32edaSLois Curfman McInnes     v = a->v;
69588e32edaSLois Curfman McInnes 
69688e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
697b0a32e0cSBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr);
698b0a32e0cSBarry Smith     cols = scols;
6990752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
70087828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
701b0a32e0cSBarry Smith     vals = svals;
7020752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
70388e32edaSLois Curfman McInnes 
70488e32edaSLois Curfman McInnes     /* insert into matrix */
70588e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
70688e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
70788e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
70888e32edaSLois Curfman McInnes     }
709606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
710606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
711606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
71288e32edaSLois Curfman McInnes 
7136d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7146d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71590ace30eSBarry Smith   }
7163a40ed3dSBarry Smith   PetscFunctionReturn(0);
71788e32edaSLois Curfman McInnes }
71888e32edaSLois Curfman McInnes 
719e090d566SSatish Balay #include "petscsys.h"
720932b0c3eSLois Curfman McInnes 
7214a2ae208SSatish Balay #undef __FUNCT__
7224a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
723b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
724289bc588SBarry Smith {
725932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
726fb9695e5SSatish Balay   int               ierr,i,j;
727fb9695e5SSatish Balay   char              *name;
72887828ca2SBarry Smith   PetscScalar       *v;
729f3ef73ceSBarry Smith   PetscViewerFormat format;
730932b0c3eSLois Curfman McInnes 
7313a40ed3dSBarry Smith   PetscFunctionBegin;
732435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
733b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
734456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7353a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
736fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
737b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
738273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
73944cd7ae7SLois Curfman McInnes       v = a->v + i;
740b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr);
741273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
742aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
743329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
74462b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
745329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
74662b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7476831982aSBarry Smith         }
74880cd9d93SLois Curfman McInnes #else
7496831982aSBarry Smith         if (*v) {
75062b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr);
7516831982aSBarry Smith         }
75280cd9d93SLois Curfman McInnes #endif
7531b807ce4Svictorle         v += a->lda;
75480cd9d93SLois Curfman McInnes       }
755b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
75680cd9d93SLois Curfman McInnes     }
757b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7583a40ed3dSBarry Smith   } else {
759b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
760aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
761ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
76247989497SBarry Smith     /* determine if matrix has all real values */
76347989497SBarry Smith     v = a->v;
764201f5f94SSatish Balay     for (i=0; i<A->m*A->n; i++) {
765ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
76647989497SBarry Smith     }
76747989497SBarry Smith #endif
768fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
7693a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
770b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr);
771fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr);
772fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
773ffac6cdbSBarry Smith     }
774ffac6cdbSBarry Smith 
775273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
776932b0c3eSLois Curfman McInnes       v = a->v + i;
777273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
778aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
77947989497SBarry Smith         if (allreal) {
7801b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
78147989497SBarry Smith         } else {
7821b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
78347989497SBarry Smith         }
784289bc588SBarry Smith #else
7851b807ce4Svictorle         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
786289bc588SBarry Smith #endif
7871b807ce4Svictorle 	v += a->lda;
788289bc588SBarry Smith       }
789b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
790289bc588SBarry Smith     }
791fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
792b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
793ffac6cdbSBarry Smith     }
794b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
795da3a660dSBarry Smith   }
796b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7973a40ed3dSBarry Smith   PetscFunctionReturn(0);
798289bc588SBarry Smith }
799289bc588SBarry Smith 
8004a2ae208SSatish Balay #undef __FUNCT__
8014a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
802b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
803932b0c3eSLois Curfman McInnes {
804932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
805273d9f13SBarry Smith   int               ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
80687828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
807f3ef73ceSBarry Smith   PetscViewerFormat format;
808932b0c3eSLois Curfman McInnes 
8093a40ed3dSBarry Smith   PetscFunctionBegin;
810b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
81190ace30eSBarry Smith 
812b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
813fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
81490ace30eSBarry Smith     /* store the matrix as a dense matrix */
815b0a32e0cSBarry Smith     ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr);
816552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
81790ace30eSBarry Smith     col_lens[1] = m;
81890ace30eSBarry Smith     col_lens[2] = n;
81990ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8200752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr);
821606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
82290ace30eSBarry Smith 
82390ace30eSBarry Smith     /* write out matrix, by rows */
82487828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
82590ace30eSBarry Smith     v    = a->v;
82690ace30eSBarry Smith     for (i=0; i<m; i++) {
82790ace30eSBarry Smith       for (j=0; j<n; j++) {
82890ace30eSBarry Smith         vals[i + j*m] = *v++;
82990ace30eSBarry Smith       }
83090ace30eSBarry Smith     }
8310752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr);
832606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
83390ace30eSBarry Smith   } else {
834b0a32e0cSBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr);
835552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
836932b0c3eSLois Curfman McInnes     col_lens[1] = m;
837932b0c3eSLois Curfman McInnes     col_lens[2] = n;
838932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
839932b0c3eSLois Curfman McInnes 
840932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
841932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8420752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr);
843932b0c3eSLois Curfman McInnes 
844932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
845932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
846932b0c3eSLois Curfman McInnes     ict = 0;
847932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
848932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
849932b0c3eSLois Curfman McInnes     }
8500752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr);
851606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
852932b0c3eSLois Curfman McInnes 
853932b0c3eSLois Curfman McInnes     /* store nonzero values */
85487828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
855932b0c3eSLois Curfman McInnes     ict  = 0;
856932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
857932b0c3eSLois Curfman McInnes       v = a->v + i;
858932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
8591b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
860932b0c3eSLois Curfman McInnes       }
861932b0c3eSLois Curfman McInnes     }
8620752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr);
863606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
86490ace30eSBarry Smith   }
8653a40ed3dSBarry Smith   PetscFunctionReturn(0);
866932b0c3eSLois Curfman McInnes }
867932b0c3eSLois Curfman McInnes 
8684a2ae208SSatish Balay #undef __FUNCT__
8694a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
870b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
871f1af5d2fSBarry Smith {
872f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
873f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
874fb9695e5SSatish Balay   int               m = A->m,n = A->n,color,i,j,ierr;
87587828ca2SBarry Smith   PetscScalar       *v = a->v;
876b0a32e0cSBarry Smith   PetscViewer       viewer;
877b0a32e0cSBarry Smith   PetscDraw         popup;
878329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
879f3ef73ceSBarry Smith   PetscViewerFormat format;
880f1af5d2fSBarry Smith 
881f1af5d2fSBarry Smith   PetscFunctionBegin;
882f1af5d2fSBarry Smith 
883f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
884b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
885b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
886f1af5d2fSBarry Smith 
887f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
888fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
889f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
890b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
891f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
892f1af5d2fSBarry Smith       x_l = j;
893f1af5d2fSBarry Smith       x_r = x_l + 1.0;
894f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
895f1af5d2fSBarry Smith         y_l = m - i - 1.0;
896f1af5d2fSBarry Smith         y_r = y_l + 1.0;
897f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
898329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
899b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
900329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
901b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
902f1af5d2fSBarry Smith         } else {
903f1af5d2fSBarry Smith           continue;
904f1af5d2fSBarry Smith         }
905f1af5d2fSBarry Smith #else
906f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
907b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
908f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
909b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
910f1af5d2fSBarry Smith         } else {
911f1af5d2fSBarry Smith           continue;
912f1af5d2fSBarry Smith         }
913f1af5d2fSBarry Smith #endif
914b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
915f1af5d2fSBarry Smith       }
916f1af5d2fSBarry Smith     }
917f1af5d2fSBarry Smith   } else {
918f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
919f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
920f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
921f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
922f1af5d2fSBarry Smith     }
923b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
924b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
925b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
926f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
927f1af5d2fSBarry Smith       x_l = j;
928f1af5d2fSBarry Smith       x_r = x_l + 1.0;
929f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
930f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
931f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
932b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
933b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
934f1af5d2fSBarry Smith       }
935f1af5d2fSBarry Smith     }
936f1af5d2fSBarry Smith   }
937f1af5d2fSBarry Smith   PetscFunctionReturn(0);
938f1af5d2fSBarry Smith }
939f1af5d2fSBarry Smith 
9404a2ae208SSatish Balay #undef __FUNCT__
9414a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
942b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
943f1af5d2fSBarry Smith {
944b0a32e0cSBarry Smith   PetscDraw  draw;
945f1af5d2fSBarry Smith   PetscTruth isnull;
946329f5518SBarry Smith   PetscReal  xr,yr,xl,yl,h,w;
947f1af5d2fSBarry Smith   int        ierr;
948f1af5d2fSBarry Smith 
949f1af5d2fSBarry Smith   PetscFunctionBegin;
950b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
951b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
952f1af5d2fSBarry Smith   if (isnull == PETSC_TRUE) PetscFunctionReturn(0);
953f1af5d2fSBarry Smith 
954f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
955273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
956f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
957b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
958b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
959f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
960f1af5d2fSBarry Smith   PetscFunctionReturn(0);
961f1af5d2fSBarry Smith }
962f1af5d2fSBarry Smith 
9634a2ae208SSatish Balay #undef __FUNCT__
9644a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
965b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer)
966932b0c3eSLois Curfman McInnes {
967932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
968bcd2baecSBarry Smith   int          ierr;
969f1af5d2fSBarry Smith   PetscTruth   issocket,isascii,isbinary,isdraw;
970932b0c3eSLois Curfman McInnes 
9713a40ed3dSBarry Smith   PetscFunctionBegin;
972b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
973b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
974fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
975fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
9760f5bd95cSBarry Smith 
9770f5bd95cSBarry Smith   if (issocket) {
9781b807ce4Svictorle     if (a->lda>A->m) SETERRQ(1,"Case can not handle LDA");
979b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
9800f5bd95cSBarry Smith   } else if (isascii) {
9813a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
9820f5bd95cSBarry Smith   } else if (isbinary) {
9833a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
984f1af5d2fSBarry Smith   } else if (isdraw) {
985f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
9865cd90555SBarry Smith   } else {
98729bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
988932b0c3eSLois Curfman McInnes   }
9893a40ed3dSBarry Smith   PetscFunctionReturn(0);
990932b0c3eSLois Curfman McInnes }
991289bc588SBarry Smith 
9924a2ae208SSatish Balay #undef __FUNCT__
9934a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
994c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
995289bc588SBarry Smith {
996ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
99790f02eecSBarry Smith   int          ierr;
99890f02eecSBarry Smith 
9993a40ed3dSBarry Smith   PetscFunctionBegin;
1000aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1001b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
1002a5a9c739SBarry Smith #endif
1003606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1004606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1005606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
10063a40ed3dSBarry Smith   PetscFunctionReturn(0);
1007289bc588SBarry Smith }
1008289bc588SBarry Smith 
10094a2ae208SSatish Balay #undef __FUNCT__
10104a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
10118f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
1012289bc588SBarry Smith {
1013c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
10141b807ce4Svictorle   int          k,j,m,n,M,ierr;
101587828ca2SBarry Smith   PetscScalar  *v,tmp;
101648b35521SBarry Smith 
10173a40ed3dSBarry Smith   PetscFunctionBegin;
10181b807ce4Svictorle   v = mat->v; m = A->m; M = mat->lda; n = A->n;
10197c922b88SBarry Smith   if (!matout) { /* in place transpose */
1020a5ce6ee0Svictorle     if (m != n) {
1021a5ce6ee0Svictorle       SETERRQ(1,"Can not transpose non-square matrix in place");
1022d64ed03dSBarry Smith     } else {
1023d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1024289bc588SBarry Smith         for (k=0; k<j; k++) {
10251b807ce4Svictorle           tmp = v[j + k*M];
10261b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10271b807ce4Svictorle           v[k + j*M] = tmp;
1028289bc588SBarry Smith         }
1029289bc588SBarry Smith       }
1030d64ed03dSBarry Smith     }
10313a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1032d3e5ee88SLois Curfman McInnes     Mat          tmat;
1033ec8511deSBarry Smith     Mat_SeqDense *tmatd;
103487828ca2SBarry Smith     PetscScalar  *v2;
1035ea709b57SSatish Balay 
1036273d9f13SBarry Smith     ierr  = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr);
1037ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10380de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1039d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10401b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1041d3e5ee88SLois Curfman McInnes     }
10426d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10436d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1044d3e5ee88SLois Curfman McInnes     *matout = tmat;
104548b35521SBarry Smith   }
10463a40ed3dSBarry Smith   PetscFunctionReturn(0);
1047289bc588SBarry Smith }
1048289bc588SBarry Smith 
10494a2ae208SSatish Balay #undef __FUNCT__
10504a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
10518f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1052289bc588SBarry Smith {
1053c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1054c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1055a43ee2ecSKris Buschelman   int          i,j;
105687828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
10579ea5d5aeSSatish Balay 
10583a40ed3dSBarry Smith   PetscFunctionBegin;
1059273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1060273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10611b807ce4Svictorle   for (i=0; i<A1->m; i++) {
10621b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
10631b807ce4Svictorle     for (j=0; j<A1->n; j++) {
10643a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10651b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
10661b807ce4Svictorle     }
1067289bc588SBarry Smith   }
106877c4ece6SBarry Smith   *flg = PETSC_TRUE;
10693a40ed3dSBarry Smith   PetscFunctionReturn(0);
1070289bc588SBarry Smith }
1071289bc588SBarry Smith 
10724a2ae208SSatish Balay #undef __FUNCT__
10734a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
10748f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
1075289bc588SBarry Smith {
1076c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
10777a97a34bSBarry Smith   int          ierr,i,n,len;
107887828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
107944cd7ae7SLois Curfman McInnes 
10803a40ed3dSBarry Smith   PetscFunctionBegin;
10817a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
10827a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
1083b1d4fb26SBarry Smith   ierr = VecGetArrayFast(v,&x);CHKERRQ(ierr);
1084273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1085273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
108644cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
10871b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1088289bc588SBarry Smith   }
1089b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(v,&x);CHKERRQ(ierr);
10903a40ed3dSBarry Smith   PetscFunctionReturn(0);
1091289bc588SBarry Smith }
1092289bc588SBarry Smith 
10934a2ae208SSatish Balay #undef __FUNCT__
10944a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
10958f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1096289bc588SBarry Smith {
1097c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
109887828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
1099273d9f13SBarry Smith   int          ierr,i,j,m = A->m,n = A->n;
110055659b69SBarry Smith 
11013a40ed3dSBarry Smith   PetscFunctionBegin;
110228988994SBarry Smith   if (ll) {
11037a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1104b1d4fb26SBarry Smith     ierr = VecGetArrayFast(ll,&l);CHKERRQ(ierr);
1105273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1106da3a660dSBarry Smith     for (i=0; i<m; i++) {
1107da3a660dSBarry Smith       x = l[i];
1108da3a660dSBarry Smith       v = mat->v + i;
1109da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1110da3a660dSBarry Smith     }
1111b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(ll,&l);CHKERRQ(ierr);
1112b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1113da3a660dSBarry Smith   }
111428988994SBarry Smith   if (rr) {
11157a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1116b1d4fb26SBarry Smith     ierr = VecGetArrayFast(rr,&r);CHKERRQ(ierr);
1117273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1118da3a660dSBarry Smith     for (i=0; i<n; i++) {
1119da3a660dSBarry Smith       x = r[i];
1120da3a660dSBarry Smith       v = mat->v + i*m;
1121da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1122da3a660dSBarry Smith     }
1123b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(rr,&r);CHKERRQ(ierr);
1124b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1125da3a660dSBarry Smith   }
11263a40ed3dSBarry Smith   PetscFunctionReturn(0);
1127289bc588SBarry Smith }
1128289bc588SBarry Smith 
11294a2ae208SSatish Balay #undef __FUNCT__
11304a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1131064f8208SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1132289bc588SBarry Smith {
1133c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
113487828ca2SBarry Smith   PetscScalar  *v = mat->v;
1135329f5518SBarry Smith   PetscReal    sum = 0.0;
1136a5ce6ee0Svictorle   int          lda=mat->lda,m=A->m,i,j;
113755659b69SBarry Smith 
11383a40ed3dSBarry Smith   PetscFunctionBegin;
1139289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1140a5ce6ee0Svictorle     if (lda>m) {
1141a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
1142a5ce6ee0Svictorle 	v = mat->v+j*lda;
1143a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1144a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1145a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1146a5ce6ee0Svictorle #else
1147a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1148a5ce6ee0Svictorle #endif
1149a5ce6ee0Svictorle 	}
1150a5ce6ee0Svictorle       }
1151a5ce6ee0Svictorle     } else {
1152273d9f13SBarry Smith       for (i=0; i<A->n*A->m; i++) {
1153aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1154329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1155289bc588SBarry Smith #else
1156289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1157289bc588SBarry Smith #endif
1158289bc588SBarry Smith       }
1159a5ce6ee0Svictorle     }
1160064f8208SBarry Smith     *nrm = sqrt(sum);
1161b0a32e0cSBarry Smith     PetscLogFlops(2*A->n*A->m);
11623a40ed3dSBarry Smith   } else if (type == NORM_1) {
1163064f8208SBarry Smith     *nrm = 0.0;
1164273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
11651b807ce4Svictorle       v = mat->v + j*mat->lda;
1166289bc588SBarry Smith       sum = 0.0;
1167273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
116833a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1169289bc588SBarry Smith       }
1170064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1171289bc588SBarry Smith     }
1172b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11733a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1174064f8208SBarry Smith     *nrm = 0.0;
1175273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1176289bc588SBarry Smith       v = mat->v + j;
1177289bc588SBarry Smith       sum = 0.0;
1178273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
11791b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1180289bc588SBarry Smith       }
1181064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1182289bc588SBarry Smith     }
1183b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11843a40ed3dSBarry Smith   } else {
118529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1186289bc588SBarry Smith   }
11873a40ed3dSBarry Smith   PetscFunctionReturn(0);
1188289bc588SBarry Smith }
1189289bc588SBarry Smith 
11904a2ae208SSatish Balay #undef __FUNCT__
11914a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
11928f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1193289bc588SBarry Smith {
1194c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
119567e560aaSBarry Smith 
11963a40ed3dSBarry Smith   PetscFunctionBegin;
1197b5a2b587SKris Buschelman   switch (op) {
1198b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1199b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1200b5a2b587SKris Buschelman     break;
1201b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1202b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1203b5a2b587SKris Buschelman     break;
1204b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1205b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1206b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1207b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1208b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1209b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1210b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1211b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1212b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1213b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1214b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
1215b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
1216b5a2b587SKris Buschelman     break;
121777e54ba9SKris Buschelman   case MAT_SYMMETRIC:
121877e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12199a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12209a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12219a4540c5SBarry Smith   case MAT_HERMITIAN:
12229a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12239a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12249a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
122577e54ba9SKris Buschelman     break;
1226b5a2b587SKris Buschelman   default:
122729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
12283a40ed3dSBarry Smith   }
12293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1230289bc588SBarry Smith }
1231289bc588SBarry Smith 
12324a2ae208SSatish Balay #undef __FUNCT__
12334a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
12348f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
12356f0a148fSBarry Smith {
1236ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1237a5ce6ee0Svictorle   int          lda=l->lda,m=A->m,j, ierr;
12383a40ed3dSBarry Smith 
12393a40ed3dSBarry Smith   PetscFunctionBegin;
1240a5ce6ee0Svictorle   if (lda>m) {
1241a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
1242a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar)); CHKERRQ(ierr);
1243a5ce6ee0Svictorle     }
1244a5ce6ee0Svictorle   } else {
124587828ca2SBarry Smith     ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1246a5ce6ee0Svictorle   }
12473a40ed3dSBarry Smith   PetscFunctionReturn(0);
12486f0a148fSBarry Smith }
12496f0a148fSBarry Smith 
12504a2ae208SSatish Balay #undef __FUNCT__
12514a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense"
12528f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
12534e220ebcSLois Curfman McInnes {
12543a40ed3dSBarry Smith   PetscFunctionBegin;
12554e220ebcSLois Curfman McInnes   *bs = 1;
12563a40ed3dSBarry Smith   PetscFunctionReturn(0);
12574e220ebcSLois Curfman McInnes }
12584e220ebcSLois Curfman McInnes 
12594a2ae208SSatish Balay #undef __FUNCT__
12604a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1261268466fbSBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag)
12626f0a148fSBarry Smith {
1263ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1264273d9f13SBarry Smith   int          n = A->n,i,j,ierr,N,*rows;
126587828ca2SBarry Smith   PetscScalar  *slot;
126655659b69SBarry Smith 
12673a40ed3dSBarry Smith   PetscFunctionBegin;
1268e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
126978b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
12706f0a148fSBarry Smith   for (i=0; i<N; i++) {
12716f0a148fSBarry Smith     slot = l->v + rows[i];
12726f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
12736f0a148fSBarry Smith   }
12746f0a148fSBarry Smith   if (diag) {
12756f0a148fSBarry Smith     for (i=0; i<N; i++) {
12766f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1277260925b8SBarry Smith       *slot = *diag;
12786f0a148fSBarry Smith     }
12796f0a148fSBarry Smith   }
1280260925b8SBarry Smith   ISRestoreIndices(is,&rows);
12813a40ed3dSBarry Smith   PetscFunctionReturn(0);
12826f0a148fSBarry Smith }
1283557bce09SLois Curfman McInnes 
12844a2ae208SSatish Balay #undef __FUNCT__
12854a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
12864e7234bfSBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar *array[])
128764e87e97SBarry Smith {
1288c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
12893a40ed3dSBarry Smith 
12903a40ed3dSBarry Smith   PetscFunctionBegin;
129164e87e97SBarry Smith   *array = mat->v;
12923a40ed3dSBarry Smith   PetscFunctionReturn(0);
129364e87e97SBarry Smith }
12940754003eSLois Curfman McInnes 
12954a2ae208SSatish Balay #undef __FUNCT__
12964a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
12974e7234bfSBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1298ff14e315SSatish Balay {
12993a40ed3dSBarry Smith   PetscFunctionBegin;
130009b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13013a40ed3dSBarry Smith   PetscFunctionReturn(0);
1302ff14e315SSatish Balay }
13030754003eSLois Curfman McInnes 
13044a2ae208SSatish Balay #undef __FUNCT__
13054a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
13067b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
13070754003eSLois Curfman McInnes {
1308c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1309273d9f13SBarry Smith   int          i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
131087828ca2SBarry Smith   PetscScalar  *av,*bv,*v = mat->v;
13110754003eSLois Curfman McInnes   Mat          newmat;
13120754003eSLois Curfman McInnes 
13133a40ed3dSBarry Smith   PetscFunctionBegin;
131478b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
131578b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1316e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1317e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13180754003eSLois Curfman McInnes 
1319182d2002SSatish Balay   /* Check submatrixcall */
1320182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1321182d2002SSatish Balay     int n_cols,n_rows;
1322182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
132329bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1324182d2002SSatish Balay     newmat = *B;
1325182d2002SSatish Balay   } else {
13260754003eSLois Curfman McInnes     /* Create and fill new matrix */
1327b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
1328182d2002SSatish Balay   }
1329182d2002SSatish Balay 
1330182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1331182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1332182d2002SSatish Balay 
1333182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1334182d2002SSatish Balay     av = v + m*icol[i];
1335182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1336182d2002SSatish Balay       *bv++ = av[irow[j]];
13370754003eSLois Curfman McInnes     }
13380754003eSLois Curfman McInnes   }
1339182d2002SSatish Balay 
1340182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13416d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13426d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13430754003eSLois Curfman McInnes 
13440754003eSLois Curfman McInnes   /* Free work space */
134578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
134678b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1347182d2002SSatish Balay   *B = newmat;
13483a40ed3dSBarry Smith   PetscFunctionReturn(0);
13490754003eSLois Curfman McInnes }
13500754003eSLois Curfman McInnes 
13514a2ae208SSatish Balay #undef __FUNCT__
13524a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
1353268466fbSBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1354905e6a2fSBarry Smith {
1355905e6a2fSBarry Smith   int ierr,i;
1356905e6a2fSBarry Smith 
13573a40ed3dSBarry Smith   PetscFunctionBegin;
1358905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1359b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1360905e6a2fSBarry Smith   }
1361905e6a2fSBarry Smith 
1362905e6a2fSBarry Smith   for (i=0; i<n; i++) {
13636a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1364905e6a2fSBarry Smith   }
13653a40ed3dSBarry Smith   PetscFunctionReturn(0);
1366905e6a2fSBarry Smith }
1367905e6a2fSBarry Smith 
13684a2ae208SSatish Balay #undef __FUNCT__
13694a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1370cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
13714b0e389bSBarry Smith {
13724b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1373a5ce6ee0Svictorle   int          lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j,ierr;
13743a40ed3dSBarry Smith 
13753a40ed3dSBarry Smith   PetscFunctionBegin;
137633f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
137733f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1378cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
13793a40ed3dSBarry Smith     PetscFunctionReturn(0);
13803a40ed3dSBarry Smith   }
13810dbb7854Svictorle   if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1382a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
13830dbb7854Svictorle     for (j=0; j<n; j++) {
1384a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1385a5ce6ee0Svictorle     }
1386a5ce6ee0Svictorle   } else {
138787828ca2SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1388a5ce6ee0Svictorle   }
1389273d9f13SBarry Smith   PetscFunctionReturn(0);
1390273d9f13SBarry Smith }
1391273d9f13SBarry Smith 
13924a2ae208SSatish Balay #undef __FUNCT__
13934a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1394273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A)
1395273d9f13SBarry Smith {
1396273d9f13SBarry Smith   int        ierr;
1397273d9f13SBarry Smith 
1398273d9f13SBarry Smith   PetscFunctionBegin;
1399273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14003a40ed3dSBarry Smith   PetscFunctionReturn(0);
14014b0e389bSBarry Smith }
14024b0e389bSBarry Smith 
1403289bc588SBarry Smith /* -------------------------------------------------------------------*/
1404a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1405905e6a2fSBarry Smith        MatGetRow_SeqDense,
1406905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1407905e6a2fSBarry Smith        MatMult_SeqDense,
140897304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
14097c922b88SBarry Smith        MatMultTranspose_SeqDense,
14107c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1411905e6a2fSBarry Smith        MatSolve_SeqDense,
1412905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
14137c922b88SBarry Smith        MatSolveTranspose_SeqDense,
141497304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1415905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1416905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1417ec8511deSBarry Smith        MatRelax_SeqDense,
1418ec8511deSBarry Smith        MatTranspose_SeqDense,
141997304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1420905e6a2fSBarry Smith        MatEqual_SeqDense,
1421905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1422905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1423905e6a2fSBarry Smith        MatNorm_SeqDense,
142497304618SKris Buschelman /*20*/ 0,
1425905e6a2fSBarry Smith        0,
1426905e6a2fSBarry Smith        0,
1427905e6a2fSBarry Smith        MatSetOption_SeqDense,
1428905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
142997304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1430905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1431905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1432905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1433905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
143497304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1435273d9f13SBarry Smith        0,
1436905e6a2fSBarry Smith        0,
1437905e6a2fSBarry Smith        MatGetArray_SeqDense,
1438905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
143997304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1440a5ae1ecdSBarry Smith        0,
1441a5ae1ecdSBarry Smith        0,
1442a5ae1ecdSBarry Smith        0,
1443a5ae1ecdSBarry Smith        0,
144497304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1445a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1446a5ae1ecdSBarry Smith        0,
14474b0e389bSBarry Smith        MatGetValues_SeqDense,
1448a5ae1ecdSBarry Smith        MatCopy_SeqDense,
144997304618SKris Buschelman /*45*/ 0,
1450a5ae1ecdSBarry Smith        MatScale_SeqDense,
1451a5ae1ecdSBarry Smith        0,
1452a5ae1ecdSBarry Smith        0,
1453a5ae1ecdSBarry Smith        0,
145497304618SKris Buschelman /*50*/ MatGetBlockSize_SeqDense,
1455a5ae1ecdSBarry Smith        0,
1456a5ae1ecdSBarry Smith        0,
1457a5ae1ecdSBarry Smith        0,
1458a5ae1ecdSBarry Smith        0,
145997304618SKris Buschelman /*55*/ 0,
1460a5ae1ecdSBarry Smith        0,
1461a5ae1ecdSBarry Smith        0,
1462a5ae1ecdSBarry Smith        0,
1463a5ae1ecdSBarry Smith        0,
146497304618SKris Buschelman /*60*/ 0,
1465e03a110bSBarry Smith        MatDestroy_SeqDense,
1466e03a110bSBarry Smith        MatView_SeqDense,
146797304618SKris Buschelman        MatGetPetscMaps_Petsc,
146897304618SKris Buschelman        0,
146997304618SKris Buschelman /*65*/ 0,
147097304618SKris Buschelman        0,
147197304618SKris Buschelman        0,
147297304618SKris Buschelman        0,
147397304618SKris Buschelman        0,
147497304618SKris Buschelman /*70*/ 0,
147597304618SKris Buschelman        0,
147697304618SKris Buschelman        0,
147797304618SKris Buschelman        0,
147897304618SKris Buschelman        0,
147997304618SKris Buschelman /*75*/ 0,
148097304618SKris Buschelman        0,
148197304618SKris Buschelman        0,
148297304618SKris Buschelman        0,
148397304618SKris Buschelman        0,
148497304618SKris Buschelman /*80*/ 0,
148597304618SKris Buschelman        0,
148697304618SKris Buschelman        0,
148797304618SKris Buschelman        0,
148897304618SKris Buschelman /*85*/ MatLoad_SeqDense};
148990ace30eSBarry Smith 
14904a2ae208SSatish Balay #undef __FUNCT__
14914a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
14924b828684SBarry Smith /*@C
1493fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1494d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1495d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1496289bc588SBarry Smith 
1497db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1498db81eaa0SLois Curfman McInnes 
149920563c6bSBarry Smith    Input Parameters:
1500db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
15010c775827SLois Curfman McInnes .  m - number of rows
150218f449edSLois Curfman McInnes .  n - number of columns
1503db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1504dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
150520563c6bSBarry Smith 
150620563c6bSBarry Smith    Output Parameter:
150744cd7ae7SLois Curfman McInnes .  A - the matrix
150820563c6bSBarry Smith 
1509b259b22eSLois Curfman McInnes    Notes:
151018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
151118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1512b4fd4287SBarry Smith    set data=PETSC_NULL.
151318f449edSLois Curfman McInnes 
1514027ccd11SLois Curfman McInnes    Level: intermediate
1515027ccd11SLois Curfman McInnes 
1516dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1517d65003e9SLois Curfman McInnes 
1518db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
151920563c6bSBarry Smith @*/
152087828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A)
1521289bc588SBarry Smith {
1522273d9f13SBarry Smith   int ierr;
15233b2fbd54SBarry Smith 
15243a40ed3dSBarry Smith   PetscFunctionBegin;
1525273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1526273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1527273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1528273d9f13SBarry Smith   PetscFunctionReturn(0);
1529273d9f13SBarry Smith }
1530273d9f13SBarry Smith 
15314a2ae208SSatish Balay #undef __FUNCT__
15324a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1533273d9f13SBarry Smith /*@C
1534273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1535273d9f13SBarry Smith 
1536273d9f13SBarry Smith    Collective on MPI_Comm
1537273d9f13SBarry Smith 
1538273d9f13SBarry Smith    Input Parameters:
1539273d9f13SBarry Smith +  A - the matrix
1540273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1541273d9f13SBarry Smith 
1542273d9f13SBarry Smith    Notes:
1543273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1544273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1545273d9f13SBarry Smith    set data=PETSC_NULL.
1546273d9f13SBarry Smith 
1547273d9f13SBarry Smith    Level: intermediate
1548273d9f13SBarry Smith 
1549273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1550273d9f13SBarry Smith 
1551273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1552273d9f13SBarry Smith @*/
1553ca01db9bSBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1554273d9f13SBarry Smith {
1555ca01db9bSBarry Smith   int ierr,(*f)(Mat,PetscScalar[]);
1556a23d5eceSKris Buschelman 
1557a23d5eceSKris Buschelman   PetscFunctionBegin;
1558a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1559a23d5eceSKris Buschelman   if (f) {
1560a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1561a23d5eceSKris Buschelman   }
1562a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1563a23d5eceSKris Buschelman }
1564a23d5eceSKris Buschelman 
1565a23d5eceSKris Buschelman EXTERN_C_BEGIN
1566a23d5eceSKris Buschelman #undef __FUNCT__
1567a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1568a23d5eceSKris Buschelman int MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1569a23d5eceSKris Buschelman {
1570273d9f13SBarry Smith   Mat_SeqDense *b;
1571273d9f13SBarry Smith   int          ierr;
1572273d9f13SBarry Smith 
1573273d9f13SBarry Smith   PetscFunctionBegin;
1574273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1575273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1576273d9f13SBarry Smith   if (!data) {
157787828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
157887828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1579273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
158087828ca2SBarry Smith     PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));
1581273d9f13SBarry Smith   } else { /* user-allocated storage */
1582273d9f13SBarry Smith     b->v          = data;
1583273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1584273d9f13SBarry Smith   }
1585273d9f13SBarry Smith   PetscFunctionReturn(0);
1586273d9f13SBarry Smith }
1587a23d5eceSKris Buschelman EXTERN_C_END
1588273d9f13SBarry Smith 
15891b807ce4Svictorle #undef __FUNCT__
15901b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
15911b807ce4Svictorle /*@C
15921b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
15931b807ce4Svictorle 
15941b807ce4Svictorle   Input parameter:
15951b807ce4Svictorle + A - the matrix
15961b807ce4Svictorle - lda - the leading dimension
15971b807ce4Svictorle 
15981b807ce4Svictorle   Notes:
15991b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
16001b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
16011b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
16021b807ce4Svictorle 
16031b807ce4Svictorle   Level: intermediate
16041b807ce4Svictorle 
16051b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
16061b807ce4Svictorle 
16071b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation()
16081b807ce4Svictorle @*/
16091b807ce4Svictorle int MatSeqDenseSetLDA(Mat B,int lda)
16101b807ce4Svictorle {
16111b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
16121b807ce4Svictorle   PetscFunctionBegin;
16131b807ce4Svictorle   if (lda<B->m) SETERRQ(1,"LDA must be at least matrix i dimension");
16141b807ce4Svictorle   b->lda = lda;
16151b807ce4Svictorle   PetscFunctionReturn(0);
16161b807ce4Svictorle }
16171b807ce4Svictorle 
16180bad9183SKris Buschelman /*MC
1619fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
16200bad9183SKris Buschelman 
16210bad9183SKris Buschelman    Options Database Keys:
16220bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
16230bad9183SKris Buschelman 
16240bad9183SKris Buschelman   Level: beginner
16250bad9183SKris Buschelman 
16260bad9183SKris Buschelman .seealso: MatCreateSeqDense
16270bad9183SKris Buschelman M*/
16280bad9183SKris Buschelman 
1629273d9f13SBarry Smith EXTERN_C_BEGIN
16304a2ae208SSatish Balay #undef __FUNCT__
16314a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1632273d9f13SBarry Smith int MatCreate_SeqDense(Mat B)
1633273d9f13SBarry Smith {
1634273d9f13SBarry Smith   Mat_SeqDense *b;
1635273d9f13SBarry Smith   int          ierr,size;
1636273d9f13SBarry Smith 
1637273d9f13SBarry Smith   PetscFunctionBegin;
1638273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
163929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
164055659b69SBarry Smith 
1641273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1642273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1643273d9f13SBarry Smith 
1644b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1645549d3d68SSatish Balay   ierr            = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr);
1646549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
164744cd7ae7SLois Curfman McInnes   B->factor       = 0;
164890f02eecSBarry Smith   B->mapping      = 0;
1649b0a32e0cSBarry Smith   PetscLogObjectMemory(B,sizeof(struct _p_Mat));
165044cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
165118f449edSLois Curfman McInnes 
16528a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
16538a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1654a5ae1ecdSBarry Smith 
165544cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1656273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1657273d9f13SBarry Smith   b->v            = 0;
16581b807ce4Svictorle   b->lda          = B->m;
16594e220ebcSLois Curfman McInnes 
1660a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1661a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1662a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
16633a40ed3dSBarry Smith   PetscFunctionReturn(0);
1664289bc588SBarry Smith }
1665273d9f13SBarry Smith EXTERN_C_END
1666