xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 606d414c6e034e02e67059b83ebaefc3ebe99698)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*606d414cSSatish Balay static char vcid[] = "$Id: dense.c,v 1.171 1999/06/08 22:55:39 balay Exp balay $";
3cb512458SBarry Smith #endif
467e560aaSBarry Smith /*
567e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
667e560aaSBarry Smith */
7289bc588SBarry Smith 
870f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
9eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
10289bc588SBarry Smith 
115615d1e5SSatish Balay #undef __FUNC__
125615d1e5SSatish Balay #define __FUNC__ "MatAXPY_SeqDense"
131987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
141987afe7SBarry Smith {
151987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
161987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
173a40ed3dSBarry Smith 
183a40ed3dSBarry Smith   PetscFunctionBegin;
191987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
201987afe7SBarry Smith   PLogFlops(2*N-1);
213a40ed3dSBarry Smith   PetscFunctionReturn(0);
221987afe7SBarry Smith }
231987afe7SBarry Smith 
245615d1e5SSatish Balay #undef __FUNC__
255615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqDense"
268f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
27289bc588SBarry Smith {
284e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
294e220ebcSLois Curfman McInnes   int          i,N = a->m*a->n,count = 0;
304e220ebcSLois Curfman McInnes   Scalar       *v = a->v;
313a40ed3dSBarry Smith 
323a40ed3dSBarry Smith   PetscFunctionBegin;
33289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
344e220ebcSLois Curfman McInnes 
354e220ebcSLois Curfman McInnes   info->rows_global       = (double)a->m;
364e220ebcSLois Curfman McInnes   info->columns_global    = (double)a->n;
374e220ebcSLois Curfman McInnes   info->rows_local        = (double)a->m;
384e220ebcSLois Curfman McInnes   info->columns_local     = (double)a->n;
394e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
404e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
414e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
424e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
434e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
444e220ebcSLois Curfman McInnes   info->mallocs           = 0;
454e220ebcSLois Curfman McInnes   info->memory            = A->mem;
464e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
474e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
484e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
494e220ebcSLois Curfman McInnes 
503a40ed3dSBarry Smith   PetscFunctionReturn(0);
51289bc588SBarry Smith }
52289bc588SBarry Smith 
535615d1e5SSatish Balay #undef __FUNC__
545615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqDense"
558f6be9afSLois Curfman McInnes int MatScale_SeqDense(Scalar *alpha,Mat inA)
5680cd9d93SLois Curfman McInnes {
5780cd9d93SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
5880cd9d93SLois Curfman McInnes   int          one = 1, nz;
5980cd9d93SLois Curfman McInnes 
603a40ed3dSBarry Smith   PetscFunctionBegin;
6180cd9d93SLois Curfman McInnes   nz = a->m*a->n;
6280cd9d93SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
6380cd9d93SLois Curfman McInnes   PLogFlops(nz);
643a40ed3dSBarry Smith   PetscFunctionReturn(0);
6580cd9d93SLois Curfman McInnes }
6680cd9d93SLois Curfman McInnes 
67289bc588SBarry Smith /* ---------------------------------------------------------------*/
68289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
69289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
705615d1e5SSatish Balay #undef __FUNC__
715615d1e5SSatish Balay #define __FUNC__ "MatLUFactor_SeqDense"
728f6be9afSLois Curfman McInnes int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
73289bc588SBarry Smith {
74c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
756abc6512SBarry Smith   int          info;
7655659b69SBarry Smith 
773a40ed3dSBarry Smith   PetscFunctionBegin;
78289bc588SBarry Smith   if (!mat->pivots) {
7945d48df9SBarry Smith     mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots);
80c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
81289bc588SBarry Smith   }
827a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
83289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
84a8c6a408SBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,0,"Bad argument to LU factorization");
85e3372554SBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Bad LU factorization");
86c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
8755659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
883a40ed3dSBarry Smith   PetscFunctionReturn(0);
89289bc588SBarry Smith }
906ee01492SSatish Balay 
915615d1e5SSatish Balay #undef __FUNC__
925609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_SeqDense"
935609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
9402cad45dSBarry Smith {
9502cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
9602cad45dSBarry Smith   int          ierr;
9702cad45dSBarry Smith   Mat          newi;
9802cad45dSBarry Smith 
993a40ed3dSBarry Smith   PetscFunctionBegin;
10002cad45dSBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi);CHKERRQ(ierr);
10102cad45dSBarry Smith   l = (Mat_SeqDense *) newi->data;
1025609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
103549d3d68SSatish Balay     ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr);
10402cad45dSBarry Smith   }
10502cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
10602cad45dSBarry Smith   *newmat = newi;
1073a40ed3dSBarry Smith   PetscFunctionReturn(0);
10802cad45dSBarry Smith }
10902cad45dSBarry Smith 
1105615d1e5SSatish Balay #undef __FUNC__
1115615d1e5SSatish Balay #define __FUNC__ "MatLUFactorSymbolic_SeqDense"
1128f6be9afSLois Curfman McInnes int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
113289bc588SBarry Smith {
1143a40ed3dSBarry Smith   int ierr;
1153a40ed3dSBarry Smith 
1163a40ed3dSBarry Smith   PetscFunctionBegin;
1175609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1183a40ed3dSBarry Smith   PetscFunctionReturn(0);
119289bc588SBarry Smith }
1206ee01492SSatish Balay 
1215615d1e5SSatish Balay #undef __FUNC__
1225615d1e5SSatish Balay #define __FUNC__ "MatLUFactorNumeric_SeqDense"
1238f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
124289bc588SBarry Smith {
12502cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data;
1263a40ed3dSBarry Smith   int          ierr;
1273a40ed3dSBarry Smith 
1283a40ed3dSBarry Smith   PetscFunctionBegin;
12902cad45dSBarry Smith   /* copy the numerical values */
130549d3d68SSatish Balay   ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr);
13102cad45dSBarry Smith   (*fact)->factor = 0;
1323a40ed3dSBarry Smith   ierr = MatLUFactor(*fact,0,0,1.0);CHKERRQ(ierr);
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
134289bc588SBarry Smith }
1356ee01492SSatish Balay 
1365615d1e5SSatish Balay #undef __FUNC__
1375615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense"
1388f6be9afSLois Curfman McInnes int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
139289bc588SBarry Smith {
1403a40ed3dSBarry Smith   int ierr;
1413a40ed3dSBarry Smith 
1423a40ed3dSBarry Smith   PetscFunctionBegin;
1433a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1443a40ed3dSBarry Smith   PetscFunctionReturn(0);
145289bc588SBarry Smith }
1466ee01492SSatish Balay 
1475615d1e5SSatish Balay #undef __FUNC__
1485615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense"
1498f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
150289bc588SBarry Smith {
1513a40ed3dSBarry Smith   int ierr;
1523a40ed3dSBarry Smith 
1533a40ed3dSBarry Smith   PetscFunctionBegin;
1543a40ed3dSBarry Smith   ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr);
1553a40ed3dSBarry Smith   PetscFunctionReturn(0);
156289bc588SBarry Smith }
1576ee01492SSatish Balay 
1585615d1e5SSatish Balay #undef __FUNC__
1595615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense"
1608f6be9afSLois Curfman McInnes int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
161289bc588SBarry Smith {
162c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
163*606d414cSSatish Balay   int           info,ierr;
16455659b69SBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
166752f5567SLois Curfman McInnes   if (mat->pivots) {
167*606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
168c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
169752f5567SLois Curfman McInnes     mat->pivots = 0;
170752f5567SLois Curfman McInnes   }
1717a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
172289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
173e3372554SBarry Smith   if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization");
174c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
17555659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177289bc588SBarry Smith }
178289bc588SBarry Smith 
1795615d1e5SSatish Balay #undef __FUNC__
1805615d1e5SSatish Balay #define __FUNC__ "MatSolve_SeqDense"
1818f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
182289bc588SBarry Smith {
183c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
184a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
1856abc6512SBarry Smith   Scalar       *x, *y;
18667e560aaSBarry Smith 
1873a40ed3dSBarry Smith   PetscFunctionBegin;
1887a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
1897a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1907a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
191549d3d68SSatish Balay   ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr);
192c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
19348d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
1947a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
19548d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
196289bc588SBarry Smith   }
197a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Matrix must be factored to solve");
198a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"MBad solve");
1997a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2007a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
20155659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
203289bc588SBarry Smith }
2046ee01492SSatish Balay 
2055615d1e5SSatish Balay #undef __FUNC__
2065615d1e5SSatish Balay #define __FUNC__ "MatSolveTrans_SeqDense"
2078f6be9afSLois Curfman McInnes int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
208da3a660dSBarry Smith {
209c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2107a97a34bSBarry Smith   int          ierr,one = 1, info;
2116abc6512SBarry Smith   Scalar       *x, *y;
21267e560aaSBarry Smith 
2133a40ed3dSBarry Smith   PetscFunctionBegin;
2147a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
2157a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2167a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
217549d3d68SSatish Balay   ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr);
218752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
219da3a660dSBarry Smith   if (mat->pivots) {
22048d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
2217a97a34bSBarry Smith   } else {
22248d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
223da3a660dSBarry Smith   }
224a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
2257a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2267a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
22755659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2283a40ed3dSBarry Smith   PetscFunctionReturn(0);
229da3a660dSBarry Smith }
2306ee01492SSatish Balay 
2315615d1e5SSatish Balay #undef __FUNC__
2325615d1e5SSatish Balay #define __FUNC__ "MatSolveAdd_SeqDense"
2338f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
234da3a660dSBarry Smith {
235c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2367a97a34bSBarry Smith   int          ierr,one = 1, info;
2376abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
238da3a660dSBarry Smith   Vec          tmp = 0;
23967e560aaSBarry Smith 
2403a40ed3dSBarry Smith   PetscFunctionBegin;
2417a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2427a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
2437a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
244da3a660dSBarry Smith   if (yy == zz) {
24578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
246c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
24778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
248da3a660dSBarry Smith   }
249549d3d68SSatish Balay   ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr);
250752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
251da3a660dSBarry Smith   if (mat->pivots) {
25248d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
253a8c6a408SBarry Smith   } else {
25448d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
255da3a660dSBarry Smith   }
256a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
257a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
258a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
2597a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2607a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
26155659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2623a40ed3dSBarry Smith   PetscFunctionReturn(0);
263da3a660dSBarry Smith }
26467e560aaSBarry Smith 
2655615d1e5SSatish Balay #undef __FUNC__
2665615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense"
2678f6be9afSLois Curfman McInnes int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
268da3a660dSBarry Smith {
269c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
2706abc6512SBarry Smith   int           one = 1, info,ierr;
2716abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
272da3a660dSBarry Smith   Vec           tmp;
27367e560aaSBarry Smith 
2743a40ed3dSBarry Smith   PetscFunctionBegin;
2757a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
2767a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2777a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
278da3a660dSBarry Smith   if (yy == zz) {
27978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
280c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
28178b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
282da3a660dSBarry Smith   }
283549d3d68SSatish Balay   ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr);
284752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
285da3a660dSBarry Smith   if (mat->pivots) {
28648d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
2873a40ed3dSBarry Smith   } else {
28848d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
289da3a660dSBarry Smith   }
290a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
29190f02eecSBarry Smith   if (tmp) {
29290f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
29390f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
2943a40ed3dSBarry Smith   } else {
29590f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
29690f02eecSBarry Smith   }
2977a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2987a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
29955659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
3003a40ed3dSBarry Smith   PetscFunctionReturn(0);
301da3a660dSBarry Smith }
302289bc588SBarry Smith /* ------------------------------------------------------------------*/
3035615d1e5SSatish Balay #undef __FUNC__
3045615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense"
3058f6be9afSLois Curfman McInnes int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
30620e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
307289bc588SBarry Smith {
308c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
309289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
310bc1b551cSSatish Balay   int          ierr, m = mat->m, i;
311aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
312bc1b551cSSatish Balay   int          o = 1;
313bc1b551cSSatish Balay #endif
314289bc588SBarry Smith 
3153a40ed3dSBarry Smith   PetscFunctionBegin;
316289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
3173a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
318bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
319289bc588SBarry Smith   }
3207a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3217a97a34bSBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
322289bc588SBarry Smith   while (its--) {
323289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
324289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
325aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
326f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
327f1747703SBarry Smith            not happy about returning a double complex */
328f1747703SBarry Smith         int    _i;
329f1747703SBarry Smith         Scalar sum = b[i];
330f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
3313f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
332f1747703SBarry Smith         }
333f1747703SBarry Smith         xt = sum;
334f1747703SBarry Smith #else
335289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
336f1747703SBarry Smith #endif
33755a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
338289bc588SBarry Smith       }
339289bc588SBarry Smith     }
340289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
341289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
342aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
343f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
344f1747703SBarry Smith            not happy about returning a double complex */
345f1747703SBarry Smith         int    _i;
346f1747703SBarry Smith         Scalar sum = b[i];
347f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
3483f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
349f1747703SBarry Smith         }
350f1747703SBarry Smith         xt = sum;
351f1747703SBarry Smith #else
352289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
353f1747703SBarry Smith #endif
35455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
355289bc588SBarry Smith       }
356289bc588SBarry Smith     }
357289bc588SBarry Smith   }
358600d6d0bSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
3597a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3603a40ed3dSBarry Smith   PetscFunctionReturn(0);
361289bc588SBarry Smith }
362289bc588SBarry Smith 
363289bc588SBarry Smith /* -----------------------------------------------------------------*/
3645615d1e5SSatish Balay #undef __FUNC__
3655615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqDense"
36644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
367289bc588SBarry Smith {
368c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
369289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
3707a97a34bSBarry Smith   int          ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0;
3713a40ed3dSBarry Smith 
3723a40ed3dSBarry Smith   PetscFunctionBegin;
3737a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
3747a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3757a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
37648d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
3777a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3787a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
37955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
3803a40ed3dSBarry Smith   PetscFunctionReturn(0);
381289bc588SBarry Smith }
3826ee01492SSatish Balay 
3835615d1e5SSatish Balay #undef __FUNC__
3845615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqDense"
38544cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
386289bc588SBarry Smith {
387c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
388289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
3897a97a34bSBarry Smith   int          ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0;
3903a40ed3dSBarry Smith 
3913a40ed3dSBarry Smith   PetscFunctionBegin;
3927a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
3937a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3947a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
39548d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
3967a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3977a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
39855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
3993a40ed3dSBarry Smith   PetscFunctionReturn(0);
400289bc588SBarry Smith }
4016ee01492SSatish Balay 
4025615d1e5SSatish Balay #undef __FUNC__
4035615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqDense"
40444cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
405289bc588SBarry Smith {
406c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
407600d6d0bSBarry Smith   Scalar       *v = mat->v, *x, *y;
4087a97a34bSBarry Smith   int          ierr,_One=1; Scalar _DOne=1.0;
4093a40ed3dSBarry Smith 
4103a40ed3dSBarry Smith   PetscFunctionBegin;
4117a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
412600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4137a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4147a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
41548d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
4167a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4177a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
41855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
4193a40ed3dSBarry Smith   PetscFunctionReturn(0);
420289bc588SBarry Smith }
4216ee01492SSatish Balay 
4225615d1e5SSatish Balay #undef __FUNC__
4235615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense"
42444cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
425289bc588SBarry Smith {
426c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
427600d6d0bSBarry Smith   Scalar       *v = mat->v, *x, *y;
4287a97a34bSBarry Smith   int          ierr,_One=1;
4296abc6512SBarry Smith   Scalar       _DOne=1.0;
4303a40ed3dSBarry Smith 
4313a40ed3dSBarry Smith   PetscFunctionBegin;
4327a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
433600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4347a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4357a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
43648d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
4377a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4387a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
43955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
4403a40ed3dSBarry Smith   PetscFunctionReturn(0);
441289bc588SBarry Smith }
442289bc588SBarry Smith 
443289bc588SBarry Smith /* -----------------------------------------------------------------*/
4445615d1e5SSatish Balay #undef __FUNC__
4455615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense"
4468f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
447289bc588SBarry Smith {
448c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
449289bc588SBarry Smith   Scalar       *v;
450289bc588SBarry Smith   int          i;
45167e560aaSBarry Smith 
4523a40ed3dSBarry Smith   PetscFunctionBegin;
453289bc588SBarry Smith   *ncols = mat->n;
454289bc588SBarry Smith   if (cols) {
45545d48df9SBarry Smith     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int));CHKPTRQ(*cols);
45673c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
457289bc588SBarry Smith   }
458289bc588SBarry Smith   if (vals) {
45945d48df9SBarry Smith     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar));CHKPTRQ(*vals);
460289bc588SBarry Smith     v = mat->v + row;
46173c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
462289bc588SBarry Smith   }
4633a40ed3dSBarry Smith   PetscFunctionReturn(0);
464289bc588SBarry Smith }
4656ee01492SSatish Balay 
4665615d1e5SSatish Balay #undef __FUNC__
4675615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense"
4688f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
469289bc588SBarry Smith {
470*606d414cSSatish Balay   int ierr;
471*606d414cSSatish Balay   PetscFunctionBegin;
472*606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
473*606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
4743a40ed3dSBarry Smith   PetscFunctionReturn(0);
475289bc588SBarry Smith }
476289bc588SBarry Smith /* ----------------------------------------------------------------*/
4775615d1e5SSatish Balay #undef __FUNC__
4785615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense"
4798f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
480e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
481289bc588SBarry Smith {
482c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
483289bc588SBarry Smith   int          i,j;
484d6dfbf8fSBarry Smith 
4853a40ed3dSBarry Smith   PetscFunctionBegin;
486289bc588SBarry Smith   if (!mat->roworiented) {
487dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
488289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
4895ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
490aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
491a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
49258804f6dSBarry Smith #endif
493289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
4945ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
495aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
496a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
49758804f6dSBarry Smith #endif
498289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
499289bc588SBarry Smith         }
500289bc588SBarry Smith       }
5013a40ed3dSBarry Smith     } else {
502289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
5035ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
504aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
505a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
50658804f6dSBarry Smith #endif
507289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
5085ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
509aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
510a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
51158804f6dSBarry Smith #endif
512289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
513289bc588SBarry Smith         }
514289bc588SBarry Smith       }
515289bc588SBarry Smith     }
5163a40ed3dSBarry Smith   } else {
517dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
518e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
5195ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
520aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
521a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
52258804f6dSBarry Smith #endif
523e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
5245ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
525aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
526a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
52758804f6dSBarry Smith #endif
528e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
529e8d4e0b9SBarry Smith         }
530e8d4e0b9SBarry Smith       }
5313a40ed3dSBarry Smith     } else {
532289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
5335ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
534aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
535a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
53658804f6dSBarry Smith #endif
537289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
5385ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
539aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
540a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
54158804f6dSBarry Smith #endif
542289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
543289bc588SBarry Smith         }
544289bc588SBarry Smith       }
545289bc588SBarry Smith     }
546e8d4e0b9SBarry Smith   }
5473a40ed3dSBarry Smith   PetscFunctionReturn(0);
548289bc588SBarry Smith }
549e8d4e0b9SBarry Smith 
5505615d1e5SSatish Balay #undef __FUNC__
5515615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense"
5528f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
553ae80bb75SLois Curfman McInnes {
554ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
555ae80bb75SLois Curfman McInnes   int          i, j;
556ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
557ae80bb75SLois Curfman McInnes 
5583a40ed3dSBarry Smith   PetscFunctionBegin;
559ae80bb75SLois Curfman McInnes   /* row-oriented output */
560ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
561ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
562ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
563ae80bb75SLois Curfman McInnes     }
564ae80bb75SLois Curfman McInnes   }
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
566ae80bb75SLois Curfman McInnes }
567ae80bb75SLois Curfman McInnes 
568289bc588SBarry Smith /* -----------------------------------------------------------------*/
569289bc588SBarry Smith 
57077c4ece6SBarry Smith #include "sys.h"
57188e32edaSLois Curfman McInnes 
5725615d1e5SSatish Balay #undef __FUNC__
5735615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense"
57419bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
57588e32edaSLois Curfman McInnes {
57688e32edaSLois Curfman McInnes   Mat_SeqDense *a;
57788e32edaSLois Curfman McInnes   Mat          B;
57888e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
57988e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
5804a41a4d3SSatish Balay   Scalar       *vals, *svals, *v,*w;
58119bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
58288e32edaSLois Curfman McInnes 
5833a40ed3dSBarry Smith   PetscFunctionBegin;
584d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
585a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
58690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
5870752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
588a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object");
58988e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
59088e32edaSLois Curfman McInnes 
59190ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
59290ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
59390ace30eSBarry Smith     B    = *A;
59490ace30eSBarry Smith     a    = (Mat_SeqDense *) B->data;
5954a41a4d3SSatish Balay     v    = a->v;
5964a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
5974a41a4d3SSatish Balay        from row major to column major */
5984a41a4d3SSatish Balay     w = (Scalar *) PetscMalloc((M*N+1)*sizeof(Scalar));CHKPTRQ(w);
59990ace30eSBarry Smith     /* read in nonzero values */
6004a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6014a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6024a41a4d3SSatish Balay     for ( j=0; j<N; j++ ) {
6034a41a4d3SSatish Balay       for ( i=0; i<M; i++ ) {
6044a41a4d3SSatish Balay         *v++ =w[i*N+j];
6054a41a4d3SSatish Balay       }
6064a41a4d3SSatish Balay     }
607*606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
6086d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6096d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61090ace30eSBarry Smith   } else {
61188e32edaSLois Curfman McInnes     /* read row lengths */
61245d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) );CHKPTRQ(rowlengths);
6130752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
61488e32edaSLois Curfman McInnes 
61588e32edaSLois Curfman McInnes     /* create our matrix */
616b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr);
61788e32edaSLois Curfman McInnes     B = *A;
61888e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
61988e32edaSLois Curfman McInnes     v = a->v;
62088e32edaSLois Curfman McInnes 
62188e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
62245d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) );CHKPTRQ(cols);
6230752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
62445d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) );CHKPTRQ(vals);
6250752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
62688e32edaSLois Curfman McInnes 
62788e32edaSLois Curfman McInnes     /* insert into matrix */
62888e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
62988e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
63088e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
63188e32edaSLois Curfman McInnes     }
632*606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
633*606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
634*606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
63588e32edaSLois Curfman McInnes 
6366d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6376d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63890ace30eSBarry Smith   }
6393a40ed3dSBarry Smith   PetscFunctionReturn(0);
64088e32edaSLois Curfman McInnes }
64188e32edaSLois Curfman McInnes 
64277c4ece6SBarry Smith #include "sys.h"
643932b0c3eSLois Curfman McInnes 
6445615d1e5SSatish Balay #undef __FUNC__
6455615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII"
646932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
647289bc588SBarry Smith {
648932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
649932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
650d636dbe3SBarry Smith   FILE         *fd;
651932b0c3eSLois Curfman McInnes   char         *outputname;
652932b0c3eSLois Curfman McInnes   Scalar       *v;
653932b0c3eSLois Curfman McInnes 
6543a40ed3dSBarry Smith   PetscFunctionBegin;
65590ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr);
65677ed5343SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
657888f2ed8SSatish Balay   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
658639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
6593a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
660f72585f2SLois Curfman McInnes   }
66102cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
66280cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
66344cd7ae7SLois Curfman McInnes       v = a->v + i;
66480cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
66580cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
666aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
667e20fef11SSatish Balay         if (PetscReal(*v) != 0.0 && PetscImaginary(*v) != 0.0) fprintf(fd," %d %g + %g i",j,PetscReal(*v),PetscImaginary(*v));
6683f6de6efSSatish Balay         else if (PetscReal(*v)) fprintf(fd," %d %g ",j,PetscReal(*v));
66980cd9d93SLois Curfman McInnes #else
67080cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
67180cd9d93SLois Curfman McInnes #endif
672d64ed03dSBarry Smith         v += a->m;
67380cd9d93SLois Curfman McInnes       }
67480cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
67580cd9d93SLois Curfman McInnes     }
6763a40ed3dSBarry Smith   } else {
677aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
67847989497SBarry Smith     int allreal = 1;
67947989497SBarry Smith     /* determine if matrix has all real values */
68047989497SBarry Smith     v = a->v;
68147989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
682e20fef11SSatish Balay       if (PetscImaginary(v[i])) { allreal = 0; break ;}
68347989497SBarry Smith     }
68447989497SBarry Smith #endif
685932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
686932b0c3eSLois Curfman McInnes       v = a->v + i;
687932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
688aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
68947989497SBarry Smith         if (allreal) {
6903f6de6efSSatish Balay           fprintf(fd,"%6.4e ",PetscReal(*v)); v += a->m;
69147989497SBarry Smith         } else {
692e20fef11SSatish Balay           fprintf(fd,"%6.4e + %6.4e i ",PetscReal(*v),PetscImaginary(*v)); v += a->m;
69347989497SBarry Smith         }
694289bc588SBarry Smith #else
695932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
696289bc588SBarry Smith #endif
697289bc588SBarry Smith       }
6988e37dc09SBarry Smith       fprintf(fd,"\n");
699289bc588SBarry Smith     }
700da3a660dSBarry Smith   }
701932b0c3eSLois Curfman McInnes   fflush(fd);
7023a40ed3dSBarry Smith   PetscFunctionReturn(0);
703289bc588SBarry Smith }
704289bc588SBarry Smith 
7055615d1e5SSatish Balay #undef __FUNC__
7065615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary"
707932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
708932b0c3eSLois Curfman McInnes {
709932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
710932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
71190ace30eSBarry Smith   int          format;
71290ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
713932b0c3eSLois Curfman McInnes 
7143a40ed3dSBarry Smith   PetscFunctionBegin;
71590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
71690ace30eSBarry Smith 
71790ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
71802cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
71990ace30eSBarry Smith     /* store the matrix as a dense matrix */
72090ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) );CHKPTRQ(col_lens);
72190ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
72290ace30eSBarry Smith     col_lens[1] = m;
72390ace30eSBarry Smith     col_lens[2] = n;
72490ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
7250752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr);
726*606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
72790ace30eSBarry Smith 
72890ace30eSBarry Smith     /* write out matrix, by rows */
72945d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(vals);
73090ace30eSBarry Smith     v    = a->v;
73190ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
73290ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
73390ace30eSBarry Smith         vals[i + j*m] = *v++;
73490ace30eSBarry Smith       }
73590ace30eSBarry Smith     }
7360752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr);
737*606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
73890ace30eSBarry Smith   } else {
7390452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) );CHKPTRQ(col_lens);
740932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
741932b0c3eSLois Curfman McInnes     col_lens[1] = m;
742932b0c3eSLois Curfman McInnes     col_lens[2] = n;
743932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
744932b0c3eSLois Curfman McInnes 
745932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
746932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
7470752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr);
748932b0c3eSLois Curfman McInnes 
749932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
750932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
751932b0c3eSLois Curfman McInnes     ict = 0;
752932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
753932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
754932b0c3eSLois Curfman McInnes     }
7550752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr);
756*606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
757932b0c3eSLois Curfman McInnes 
758932b0c3eSLois Curfman McInnes     /* store nonzero values */
75945d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar));CHKPTRQ(anonz);
760932b0c3eSLois Curfman McInnes     ict = 0;
761932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
762932b0c3eSLois Curfman McInnes       v = a->v + i;
763932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
764932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
765932b0c3eSLois Curfman McInnes       }
766932b0c3eSLois Curfman McInnes     }
7670752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr);
768*606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
76990ace30eSBarry Smith   }
7703a40ed3dSBarry Smith   PetscFunctionReturn(0);
771932b0c3eSLois Curfman McInnes }
772932b0c3eSLois Curfman McInnes 
7735615d1e5SSatish Balay #undef __FUNC__
7745615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense"
775c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer)
776932b0c3eSLois Curfman McInnes {
777932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
778bcd2baecSBarry Smith   ViewerType   vtype;
779bcd2baecSBarry Smith   int          ierr;
780932b0c3eSLois Curfman McInnes 
7813a40ed3dSBarry Smith   PetscFunctionBegin;
782bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
783bcd2baecSBarry Smith 
7847b2a1423SBarry Smith   if (PetscTypeCompare(vtype,SOCKET_VIEWER)) {
7857b2a1423SBarry Smith     ierr = ViewerSocketPutScalar_Private(viewer,a->m,a->n,a->v);CHKERRQ(ierr);
7863f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
7873a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
7883f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) {
7893a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
7905cd90555SBarry Smith   } else {
7915cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
792932b0c3eSLois Curfman McInnes   }
7933a40ed3dSBarry Smith   PetscFunctionReturn(0);
794932b0c3eSLois Curfman McInnes }
795289bc588SBarry Smith 
7965615d1e5SSatish Balay #undef __FUNC__
7975615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense"
798c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
799289bc588SBarry Smith {
800ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
80190f02eecSBarry Smith   int          ierr;
80290f02eecSBarry Smith 
8033a40ed3dSBarry Smith   PetscFunctionBegin;
80494d884c6SBarry Smith   if (--mat->refct > 0) PetscFunctionReturn(0);
80594d884c6SBarry Smith 
80694d884c6SBarry Smith   if (mat->mapping) {
80794d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr);
80894d884c6SBarry Smith   }
80994d884c6SBarry Smith   if (mat->bmapping) {
81094d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr);
81194d884c6SBarry Smith   }
812aa482453SBarry Smith #if defined(PETSC_USE_LOG)
813c6e7dd08SBarry Smith   PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n);
814a5a9c739SBarry Smith #endif
815*606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
816*606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
817*606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
818a5ae1ecdSBarry Smith   if (mat->rmap) {
819a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
820a5ae1ecdSBarry Smith   }
821a5ae1ecdSBarry Smith   if (mat->cmap) {
822a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
82390f02eecSBarry Smith   }
824a5a9c739SBarry Smith   PLogObjectDestroy(mat);
8250452661fSBarry Smith   PetscHeaderDestroy(mat);
8263a40ed3dSBarry Smith   PetscFunctionReturn(0);
827289bc588SBarry Smith }
828289bc588SBarry Smith 
8295615d1e5SSatish Balay #undef __FUNC__
8305615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense"
8318f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
832289bc588SBarry Smith {
833c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
834549d3d68SSatish Balay   int          k, j, m, n, ierr;
835d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
83648b35521SBarry Smith 
8373a40ed3dSBarry Smith   PetscFunctionBegin;
838d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
8393638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
840d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
841d64ed03dSBarry Smith       Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w);
842d64ed03dSBarry Smith       for ( j=0; j<m; j++ ) {
843d64ed03dSBarry Smith         for ( k=0; k<n; k++ ) {
844d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
845d64ed03dSBarry Smith         }
846d64ed03dSBarry Smith       }
847549d3d68SSatish Balay       ierr = PetscMemcpy(v,w,m*n*sizeof(Scalar));CHKERRQ(ierr);
848*606d414cSSatish Balay       ierr = PetscFree(w);CHKERRQ(ierr);
849d64ed03dSBarry Smith     } else {
850d3e5ee88SLois Curfman McInnes       for ( j=0; j<m; j++ ) {
851289bc588SBarry Smith         for ( k=0; k<j; k++ ) {
852d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
853d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
854d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
855289bc588SBarry Smith         }
856289bc588SBarry Smith       }
857d64ed03dSBarry Smith     }
8583a40ed3dSBarry Smith   } else { /* out-of-place transpose */
859d3e5ee88SLois Curfman McInnes     Mat          tmat;
860ec8511deSBarry Smith     Mat_SeqDense *tmatd;
861d3e5ee88SLois Curfman McInnes     Scalar       *v2;
862b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat);CHKERRQ(ierr);
863ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
8640de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
865d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
8660de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
867d3e5ee88SLois Curfman McInnes     }
8686d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8696d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
870d3e5ee88SLois Curfman McInnes     *matout = tmat;
87148b35521SBarry Smith   }
8723a40ed3dSBarry Smith   PetscFunctionReturn(0);
873289bc588SBarry Smith }
874289bc588SBarry Smith 
8755615d1e5SSatish Balay #undef __FUNC__
8765615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense"
8778f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
878289bc588SBarry Smith {
879c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
880c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
881289bc588SBarry Smith   int          i;
882289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
8839ea5d5aeSSatish Balay 
8843a40ed3dSBarry Smith   PetscFunctionBegin;
8853a40ed3dSBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
8863a40ed3dSBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
8873a40ed3dSBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
888289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
8893a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
890289bc588SBarry Smith     v1++; v2++;
891289bc588SBarry Smith   }
89277c4ece6SBarry Smith   *flg = PETSC_TRUE;
8933a40ed3dSBarry Smith   PetscFunctionReturn(0);
894289bc588SBarry Smith }
895289bc588SBarry Smith 
8965615d1e5SSatish Balay #undef __FUNC__
8975615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense"
8988f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
899289bc588SBarry Smith {
900c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
9017a97a34bSBarry Smith   int          ierr,i, n, len;
90244cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
90344cd7ae7SLois Curfman McInnes 
9043a40ed3dSBarry Smith   PetscFunctionBegin;
9057a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
9067a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
907600d6d0bSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
90844cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
909e3372554SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
91044cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
911289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
912289bc588SBarry Smith   }
9137a97a34bSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
9143a40ed3dSBarry Smith   PetscFunctionReturn(0);
915289bc588SBarry Smith }
916289bc588SBarry Smith 
9175615d1e5SSatish Balay #undef __FUNC__
9185615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense"
9198f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
920289bc588SBarry Smith {
921c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
922da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
9237a97a34bSBarry Smith   int          ierr,i,j,m = mat->m, n = mat->n;
92455659b69SBarry Smith 
9253a40ed3dSBarry Smith   PetscFunctionBegin;
92628988994SBarry Smith   if (ll) {
9277a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
928600d6d0bSBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
929e3372554SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
930da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
931da3a660dSBarry Smith       x = l[i];
932da3a660dSBarry Smith       v = mat->v + i;
933da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
934da3a660dSBarry Smith     }
9357a97a34bSBarry Smith     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
93644cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
937da3a660dSBarry Smith   }
93828988994SBarry Smith   if (rr) {
9397a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
940600d6d0bSBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
941e3372554SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
942da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
943da3a660dSBarry Smith       x = r[i];
944da3a660dSBarry Smith       v = mat->v + i*m;
945da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
946da3a660dSBarry Smith     }
9477a97a34bSBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
94844cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
949da3a660dSBarry Smith   }
9503a40ed3dSBarry Smith   PetscFunctionReturn(0);
951289bc588SBarry Smith }
952289bc588SBarry Smith 
9535615d1e5SSatish Balay #undef __FUNC__
9545615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense"
9558f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm)
956289bc588SBarry Smith {
957c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
958289bc588SBarry Smith   Scalar       *v = mat->v;
959289bc588SBarry Smith   double       sum = 0.0;
960289bc588SBarry Smith   int          i, j;
96155659b69SBarry Smith 
9623a40ed3dSBarry Smith   PetscFunctionBegin;
963289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
964289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
965aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
9663f6de6efSSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
967289bc588SBarry Smith #else
968289bc588SBarry Smith       sum += (*v)*(*v); v++;
969289bc588SBarry Smith #endif
970289bc588SBarry Smith     }
971289bc588SBarry Smith     *norm = sqrt(sum);
97255659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
9733a40ed3dSBarry Smith   } else if (type == NORM_1) {
974289bc588SBarry Smith     *norm = 0.0;
975289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
976289bc588SBarry Smith       sum = 0.0;
977289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
97833a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
979289bc588SBarry Smith       }
980289bc588SBarry Smith       if (sum > *norm) *norm = sum;
981289bc588SBarry Smith     }
98255659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9833a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
984289bc588SBarry Smith     *norm = 0.0;
985289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
986289bc588SBarry Smith       v = mat->v + j;
987289bc588SBarry Smith       sum = 0.0;
988289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
989cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
990289bc588SBarry Smith       }
991289bc588SBarry Smith       if (sum > *norm) *norm = sum;
992289bc588SBarry Smith     }
99355659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9943a40ed3dSBarry Smith   } else {
995e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
996289bc588SBarry Smith   }
9973a40ed3dSBarry Smith   PetscFunctionReturn(0);
998289bc588SBarry Smith }
999289bc588SBarry Smith 
10005615d1e5SSatish Balay #undef __FUNC__
10015615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense"
10028f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1003289bc588SBarry Smith {
1004c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
100567e560aaSBarry Smith 
10063a40ed3dSBarry Smith   PetscFunctionBegin;
10076d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
10086d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
10096d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
1010219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
10116d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
1012219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
10136d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
10146d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
10156d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
10166d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
10174787f768SSatish Balay            op == MAT_NEW_NONZERO_LOCATION_ERR ||
10186d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
101990f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
1020b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
1021b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
1022981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
10233a40ed3dSBarry Smith   else {
10243a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
10253a40ed3dSBarry Smith   }
10263a40ed3dSBarry Smith   PetscFunctionReturn(0);
1027289bc588SBarry Smith }
1028289bc588SBarry Smith 
10295615d1e5SSatish Balay #undef __FUNC__
10305615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense"
10318f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
10326f0a148fSBarry Smith {
1033ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
1034549d3d68SSatish Balay   int          ierr;
10353a40ed3dSBarry Smith 
10363a40ed3dSBarry Smith   PetscFunctionBegin;
1037549d3d68SSatish Balay   ierr = PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));CHKERRQ(ierr);
10383a40ed3dSBarry Smith   PetscFunctionReturn(0);
10396f0a148fSBarry Smith }
10406f0a148fSBarry Smith 
10415615d1e5SSatish Balay #undef __FUNC__
10425615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense"
10438f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
10444e220ebcSLois Curfman McInnes {
10453a40ed3dSBarry Smith   PetscFunctionBegin;
10464e220ebcSLois Curfman McInnes   *bs = 1;
10473a40ed3dSBarry Smith   PetscFunctionReturn(0);
10484e220ebcSLois Curfman McInnes }
10494e220ebcSLois Curfman McInnes 
10505615d1e5SSatish Balay #undef __FUNC__
10515615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense"
10528f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
10536f0a148fSBarry Smith {
1054ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
10556abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
10566f0a148fSBarry Smith   Scalar       *slot;
105755659b69SBarry Smith 
10583a40ed3dSBarry Smith   PetscFunctionBegin;
105977c4ece6SBarry Smith   ierr = ISGetSize(is,&N);CHKERRQ(ierr);
106078b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
10616f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
10626f0a148fSBarry Smith     slot = l->v + rows[i];
10636f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
10646f0a148fSBarry Smith   }
10656f0a148fSBarry Smith   if (diag) {
10666f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
10676f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1068260925b8SBarry Smith       *slot = *diag;
10696f0a148fSBarry Smith     }
10706f0a148fSBarry Smith   }
1071260925b8SBarry Smith   ISRestoreIndices(is,&rows);
10723a40ed3dSBarry Smith   PetscFunctionReturn(0);
10736f0a148fSBarry Smith }
1074557bce09SLois Curfman McInnes 
10755615d1e5SSatish Balay #undef __FUNC__
10765615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense"
10778f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n)
1078557bce09SLois Curfman McInnes {
1079c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10803a40ed3dSBarry Smith 
10813a40ed3dSBarry Smith   PetscFunctionBegin;
1082557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
10833a40ed3dSBarry Smith   PetscFunctionReturn(0);
1084557bce09SLois Curfman McInnes }
1085557bce09SLois Curfman McInnes 
10865615d1e5SSatish Balay #undef __FUNC__
10875615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense"
10888f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1089d3e5ee88SLois Curfman McInnes {
1090c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10913a40ed3dSBarry Smith 
10923a40ed3dSBarry Smith   PetscFunctionBegin;
1093d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
10943a40ed3dSBarry Smith   PetscFunctionReturn(0);
1095d3e5ee88SLois Curfman McInnes }
1096d3e5ee88SLois Curfman McInnes 
10975615d1e5SSatish Balay #undef __FUNC__
10985615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense"
10998f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array)
110064e87e97SBarry Smith {
1101c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
11023a40ed3dSBarry Smith 
11033a40ed3dSBarry Smith   PetscFunctionBegin;
110464e87e97SBarry Smith   *array = mat->v;
11053a40ed3dSBarry Smith   PetscFunctionReturn(0);
110664e87e97SBarry Smith }
11070754003eSLois Curfman McInnes 
11085615d1e5SSatish Balay #undef __FUNC__
11095615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense"
11108f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1111ff14e315SSatish Balay {
11123a40ed3dSBarry Smith   PetscFunctionBegin;
111309b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
11143a40ed3dSBarry Smith   PetscFunctionReturn(0);
1115ff14e315SSatish Balay }
11160754003eSLois Curfman McInnes 
11175615d1e5SSatish Balay #undef __FUNC__
11185615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense"
11197b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
11200754003eSLois Curfman McInnes {
1121c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1122182d2002SSatish Balay   int          i, j, ierr, m = mat->m, *irow, *icol, nrows, ncols;
1123182d2002SSatish Balay   Scalar       *av, *bv, *v = mat->v;
11240754003eSLois Curfman McInnes   Mat          newmat;
11250754003eSLois Curfman McInnes 
11263a40ed3dSBarry Smith   PetscFunctionBegin;
112778b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
112878b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
112978b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
113078b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr);
11310754003eSLois Curfman McInnes 
1132182d2002SSatish Balay   /* Check submatrixcall */
1133182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1134182d2002SSatish Balay     int n_cols,n_rows;
1135182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
1136182d2002SSatish Balay     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1137182d2002SSatish Balay     newmat = *B;
1138182d2002SSatish Balay   } else {
11390754003eSLois Curfman McInnes     /* Create and fill new matrix */
1140b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
1141182d2002SSatish Balay   }
1142182d2002SSatish Balay 
1143182d2002SSatish Balay   /* Now extract the data pointers and do the copy, column at a time */
1144182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1145182d2002SSatish Balay 
1146182d2002SSatish Balay   for ( i=0; i<ncols; i++ ) {
1147182d2002SSatish Balay     av = v + m*icol[i];
1148182d2002SSatish Balay     for (j=0; j<nrows; j++ ) {
1149182d2002SSatish Balay       *bv++ = av[irow[j]];
11500754003eSLois Curfman McInnes     }
11510754003eSLois Curfman McInnes   }
1152182d2002SSatish Balay 
1153182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
11546d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11556d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11560754003eSLois Curfman McInnes 
11570754003eSLois Curfman McInnes   /* Free work space */
115878b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
115978b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1160182d2002SSatish Balay   *B = newmat;
11613a40ed3dSBarry Smith   PetscFunctionReturn(0);
11620754003eSLois Curfman McInnes }
11630754003eSLois Curfman McInnes 
11645615d1e5SSatish Balay #undef __FUNC__
11655615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense"
11667b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B)
1167905e6a2fSBarry Smith {
1168905e6a2fSBarry Smith   int ierr,i;
1169905e6a2fSBarry Smith 
11703a40ed3dSBarry Smith   PetscFunctionBegin;
1171905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1172905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) );CHKPTRQ(*B);
1173905e6a2fSBarry Smith   }
1174905e6a2fSBarry Smith 
1175905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
11766a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1177905e6a2fSBarry Smith   }
11783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1179905e6a2fSBarry Smith }
1180905e6a2fSBarry Smith 
11815615d1e5SSatish Balay #undef __FUNC__
11825615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense"
1183cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A, Mat B,MatStructure str)
11844b0e389bSBarry Smith {
11854b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
11863a40ed3dSBarry Smith   int          ierr;
11873a40ed3dSBarry Smith 
11883a40ed3dSBarry Smith   PetscFunctionBegin;
11893a40ed3dSBarry Smith   if (B->type != MATSEQDENSE) {
1190cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
11913a40ed3dSBarry Smith     PetscFunctionReturn(0);
11923a40ed3dSBarry Smith   }
1193e3372554SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
1194549d3d68SSatish Balay   ierr = PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));CHKERRQ(ierr);
11953a40ed3dSBarry Smith   PetscFunctionReturn(0);
11964b0e389bSBarry Smith }
11974b0e389bSBarry Smith 
1198289bc588SBarry Smith /* -------------------------------------------------------------------*/
1199a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1200905e6a2fSBarry Smith        MatGetRow_SeqDense,
1201905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1202905e6a2fSBarry Smith        MatMult_SeqDense,
1203905e6a2fSBarry Smith        MatMultAdd_SeqDense,
1204905e6a2fSBarry Smith        MatMultTrans_SeqDense,
1205905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
1206905e6a2fSBarry Smith        MatSolve_SeqDense,
1207905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
1208905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
1209905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
1210905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1211905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1212ec8511deSBarry Smith        MatRelax_SeqDense,
1213ec8511deSBarry Smith        MatTranspose_SeqDense,
1214905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1215905e6a2fSBarry Smith        MatEqual_SeqDense,
1216905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1217905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1218905e6a2fSBarry Smith        MatNorm_SeqDense,
1219905e6a2fSBarry Smith        0,
1220905e6a2fSBarry Smith        0,
1221905e6a2fSBarry Smith        0,
1222905e6a2fSBarry Smith        MatSetOption_SeqDense,
1223905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1224905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1225905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1226905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1227905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1228905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1229905e6a2fSBarry Smith        MatGetSize_SeqDense,
1230905e6a2fSBarry Smith        MatGetSize_SeqDense,
1231905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1232905e6a2fSBarry Smith        0,
1233905e6a2fSBarry Smith        0,
1234905e6a2fSBarry Smith        MatGetArray_SeqDense,
1235905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
12365609ef8eSBarry Smith        MatDuplicate_SeqDense,
1237a5ae1ecdSBarry Smith        0,
1238a5ae1ecdSBarry Smith        0,
1239a5ae1ecdSBarry Smith        0,
1240a5ae1ecdSBarry Smith        0,
1241a5ae1ecdSBarry Smith        MatAXPY_SeqDense,
1242a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1243a5ae1ecdSBarry Smith        0,
12444b0e389bSBarry Smith        MatGetValues_SeqDense,
1245a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1246a5ae1ecdSBarry Smith        0,
1247a5ae1ecdSBarry Smith        MatScale_SeqDense,
1248a5ae1ecdSBarry Smith        0,
1249a5ae1ecdSBarry Smith        0,
1250a5ae1ecdSBarry Smith        0,
1251a5ae1ecdSBarry Smith        MatGetBlockSize_SeqDense,
1252a5ae1ecdSBarry Smith        0,
1253a5ae1ecdSBarry Smith        0,
1254a5ae1ecdSBarry Smith        0,
1255a5ae1ecdSBarry Smith        0,
1256a5ae1ecdSBarry Smith        0,
1257a5ae1ecdSBarry Smith        0,
1258a5ae1ecdSBarry Smith        0,
1259a5ae1ecdSBarry Smith        0,
1260a5ae1ecdSBarry Smith        0,
1261a5ae1ecdSBarry Smith        0,
1262a5ae1ecdSBarry Smith        0,
1263a5ae1ecdSBarry Smith        0,
1264a5ae1ecdSBarry Smith        MatGetMaps_Petsc};
126590ace30eSBarry Smith 
12665615d1e5SSatish Balay #undef __FUNC__
12675615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense"
12684b828684SBarry Smith /*@C
1269fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1270d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1271d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1272289bc588SBarry Smith 
1273db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1274db81eaa0SLois Curfman McInnes 
127520563c6bSBarry Smith    Input Parameters:
1276db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
12770c775827SLois Curfman McInnes .  m - number of rows
127818f449edSLois Curfman McInnes .  n - number of columns
1279db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1280dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
128120563c6bSBarry Smith 
128220563c6bSBarry Smith    Output Parameter:
128344cd7ae7SLois Curfman McInnes .  A - the matrix
128420563c6bSBarry Smith 
1285b259b22eSLois Curfman McInnes    Notes:
128618f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
128718f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1288b4fd4287SBarry Smith    set data=PETSC_NULL.
128918f449edSLois Curfman McInnes 
1290027ccd11SLois Curfman McInnes    Level: intermediate
1291027ccd11SLois Curfman McInnes 
1292dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1293d65003e9SLois Curfman McInnes 
1294db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
129520563c6bSBarry Smith @*/
129644cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1297289bc588SBarry Smith {
129844cd7ae7SLois Curfman McInnes   Mat          B;
129944cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
13003b2fbd54SBarry Smith   int          ierr,flg,size;
13013b2fbd54SBarry Smith 
13023a40ed3dSBarry Smith   PetscFunctionBegin;
1303d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1304a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
130555659b69SBarry Smith 
130644cd7ae7SLois Curfman McInnes   *A            = 0;
13073f1db9ecSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,"Mat",comm,MatDestroy,MatView);
130844cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
130944cd7ae7SLois Curfman McInnes   b               = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense));CHKPTRQ(b);
1310549d3d68SSatish Balay   ierr            = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr);
1311549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1312c6e7dd08SBarry Smith   B->ops->destroy = MatDestroy_SeqDense;
1313c6e7dd08SBarry Smith   B->ops->view    = MatView_SeqDense;
131444cd7ae7SLois Curfman McInnes   B->factor       = 0;
131590f02eecSBarry Smith   B->mapping      = 0;
1316f09e8eb9SSatish Balay   PLogObjectMemory(B,sizeof(struct _p_Mat));
131744cd7ae7SLois Curfman McInnes   B->data         = (void *) b;
131818f449edSLois Curfman McInnes 
131944cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
132044cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
1321a5ae1ecdSBarry Smith 
1322488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
1323488ecbafSBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1324a5ae1ecdSBarry Smith 
132544cd7ae7SLois Curfman McInnes   b->pivots       = 0;
132644cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1327b4fd4287SBarry Smith   if (data == PETSC_NULL) {
132844cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(b->v);
1329549d3d68SSatish Balay     ierr = PetscMemzero(b->v,m*n*sizeof(Scalar));CHKERRQ(ierr);
133044cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
133144cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
1332273fdc76SBarry Smith   } else { /* user-allocated storage */
133344cd7ae7SLois Curfman McInnes     b->v = data;
133444cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
13352dd2b441SLois Curfman McInnes   }
133625cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
133744cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr);}
13384e220ebcSLois Curfman McInnes 
133944cd7ae7SLois Curfman McInnes   *A = B;
13403a40ed3dSBarry Smith   PetscFunctionReturn(0);
1341289bc588SBarry Smith }
1342289bc588SBarry Smith 
13433b2fbd54SBarry Smith 
13443b2fbd54SBarry Smith 
13453b2fbd54SBarry Smith 
13463b2fbd54SBarry Smith 
13473b2fbd54SBarry Smith 
1348