xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 4787f7683dc5ea85e35b4269f7afceca0dae908f)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*4787f768SSatish Balay static char vcid[] = "$Id: dense.c,v 1.164 1999/02/03 03:18:16 curfman 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) {
10302cad45dSBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
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 */
13002cad45dSBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
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;
1636abc6512SBarry Smith   int           info;
16455659b69SBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
166752f5567SLois Curfman McInnes   if (mat->pivots) {
1670452661fSBarry Smith     PetscFree(mat->pivots);
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);
191416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
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);
217416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
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   }
249416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
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   }
283416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
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;
311bc1b551cSSatish Balay #if !defined(USE_PETSC_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++ ) {
3253a40ed3dSBarry Smith #if defined(USE_PETSC_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-- ) {
3423a40ed3dSBarry Smith #if defined(USE_PETSC_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 {
4700452661fSBarry Smith   if (cols) { PetscFree(*cols); }
4710452661fSBarry Smith   if (vals) { PetscFree(*vals); }
4723a40ed3dSBarry Smith   PetscFunctionReturn(0);
473289bc588SBarry Smith }
474289bc588SBarry Smith /* ----------------------------------------------------------------*/
4755615d1e5SSatish Balay #undef __FUNC__
4765615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense"
4778f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
478e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
479289bc588SBarry Smith {
480c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
481289bc588SBarry Smith   int          i,j;
482d6dfbf8fSBarry Smith 
4833a40ed3dSBarry Smith   PetscFunctionBegin;
484289bc588SBarry Smith   if (!mat->roworiented) {
485dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
486289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
4875ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
4883a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
489a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
49058804f6dSBarry Smith #endif
491289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
4925ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
4933a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
494a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
49558804f6dSBarry Smith #endif
496289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
497289bc588SBarry Smith         }
498289bc588SBarry Smith       }
4993a40ed3dSBarry Smith     } else {
500289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
5015ef9f2a5SBarry Smith         if (indexn[j] < 0) {v += m; continue;}
5023a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
503a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
50458804f6dSBarry Smith #endif
505289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
5065ef9f2a5SBarry Smith           if (indexm[i] < 0) {v++; continue;}
5073a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
508a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
50958804f6dSBarry Smith #endif
510289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
511289bc588SBarry Smith         }
512289bc588SBarry Smith       }
513289bc588SBarry Smith     }
5143a40ed3dSBarry Smith   } else {
515dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
516e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
5175ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
5183a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
519a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
52058804f6dSBarry Smith #endif
521e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
5225ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
5233a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
524a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
52558804f6dSBarry Smith #endif
526e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
527e8d4e0b9SBarry Smith         }
528e8d4e0b9SBarry Smith       }
5293a40ed3dSBarry Smith     } else {
530289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
5315ef9f2a5SBarry Smith         if (indexm[i] < 0) { v += n; continue;}
5323a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
533a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
53458804f6dSBarry Smith #endif
535289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
5365ef9f2a5SBarry Smith           if (indexn[j] < 0) { v++; continue;}
5373a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
538a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
53958804f6dSBarry Smith #endif
540289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
541289bc588SBarry Smith         }
542289bc588SBarry Smith       }
543289bc588SBarry Smith     }
544e8d4e0b9SBarry Smith   }
5453a40ed3dSBarry Smith   PetscFunctionReturn(0);
546289bc588SBarry Smith }
547e8d4e0b9SBarry Smith 
5485615d1e5SSatish Balay #undef __FUNC__
5495615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense"
5508f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
551ae80bb75SLois Curfman McInnes {
552ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
553ae80bb75SLois Curfman McInnes   int          i, j;
554ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
555ae80bb75SLois Curfman McInnes 
5563a40ed3dSBarry Smith   PetscFunctionBegin;
557ae80bb75SLois Curfman McInnes   /* row-oriented output */
558ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
559ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
560ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
561ae80bb75SLois Curfman McInnes     }
562ae80bb75SLois Curfman McInnes   }
5633a40ed3dSBarry Smith   PetscFunctionReturn(0);
564ae80bb75SLois Curfman McInnes }
565ae80bb75SLois Curfman McInnes 
566289bc588SBarry Smith /* -----------------------------------------------------------------*/
567289bc588SBarry Smith 
56877c4ece6SBarry Smith #include "sys.h"
56988e32edaSLois Curfman McInnes 
5705615d1e5SSatish Balay #undef __FUNC__
5715615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense"
57219bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
57388e32edaSLois Curfman McInnes {
57488e32edaSLois Curfman McInnes   Mat_SeqDense *a;
57588e32edaSLois Curfman McInnes   Mat          B;
57688e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
57788e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
5784a41a4d3SSatish Balay   Scalar       *vals, *svals, *v,*w;
57919bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
58088e32edaSLois Curfman McInnes 
5813a40ed3dSBarry Smith   PetscFunctionBegin;
58288e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
583a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
58490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
5850752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
586a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object");
58788e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
58888e32edaSLois Curfman McInnes 
58990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
59090ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
59190ace30eSBarry Smith     B    = *A;
59290ace30eSBarry Smith     a    = (Mat_SeqDense *) B->data;
5934a41a4d3SSatish Balay     v    = a->v;
5944a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
5954a41a4d3SSatish Balay        from row major to column major */
5964a41a4d3SSatish Balay     w = (Scalar *) PetscMalloc((M*N+1)*sizeof(Scalar));CHKPTRQ(w);
59790ace30eSBarry Smith     /* read in nonzero values */
5984a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR); CHKERRQ(ierr);
5994a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6004a41a4d3SSatish Balay     for ( j=0; j<N; j++ ) {
6014a41a4d3SSatish Balay       for ( i=0; i<M; i++ ) {
6024a41a4d3SSatish Balay         *v++ =w[i*N+j];
6034a41a4d3SSatish Balay       }
6044a41a4d3SSatish Balay     }
6054a41a4d3SSatish Balay     PetscFree(w);
6066d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6076d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
60890ace30eSBarry Smith   } else {
60988e32edaSLois Curfman McInnes     /* read row lengths */
61045d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
6110752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
61288e32edaSLois Curfman McInnes 
61388e32edaSLois Curfman McInnes     /* create our matrix */
614b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
61588e32edaSLois Curfman McInnes     B = *A;
61688e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
61788e32edaSLois Curfman McInnes     v = a->v;
61888e32edaSLois Curfman McInnes 
61988e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
62045d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
6210752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
62245d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
6230752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
62488e32edaSLois Curfman McInnes 
62588e32edaSLois Curfman McInnes     /* insert into matrix */
62688e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
62788e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
62888e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
62988e32edaSLois Curfman McInnes     }
6300452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
63188e32edaSLois Curfman McInnes 
6326d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6336d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
63490ace30eSBarry Smith   }
6353a40ed3dSBarry Smith   PetscFunctionReturn(0);
63688e32edaSLois Curfman McInnes }
63788e32edaSLois Curfman McInnes 
63877c4ece6SBarry Smith #include "sys.h"
639932b0c3eSLois Curfman McInnes 
6405615d1e5SSatish Balay #undef __FUNC__
6415615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII"
642932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
643289bc588SBarry Smith {
644932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
645932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
646d636dbe3SBarry Smith   FILE         *fd;
647932b0c3eSLois Curfman McInnes   char         *outputname;
648932b0c3eSLois Curfman McInnes   Scalar       *v;
649932b0c3eSLois Curfman McInnes 
6503a40ed3dSBarry Smith   PetscFunctionBegin;
65190ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
65277ed5343SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname); CHKERRQ(ierr);
65390ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
654639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
6553a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
656f72585f2SLois Curfman McInnes   }
65702cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
65880cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
65944cd7ae7SLois Curfman McInnes       v = a->v + i;
66080cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
66180cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6623a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
663e20fef11SSatish Balay         if (PetscReal(*v) != 0.0 && PetscImaginary(*v) != 0.0) fprintf(fd," %d %g + %g i",j,PetscReal(*v),PetscImaginary(*v));
6643f6de6efSSatish Balay         else if (PetscReal(*v)) fprintf(fd," %d %g ",j,PetscReal(*v));
66580cd9d93SLois Curfman McInnes #else
66680cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
66780cd9d93SLois Curfman McInnes #endif
668d64ed03dSBarry Smith         v += a->m;
66980cd9d93SLois Curfman McInnes       }
67080cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
67180cd9d93SLois Curfman McInnes     }
6723a40ed3dSBarry Smith   } else {
6733a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
67447989497SBarry Smith     int allreal = 1;
67547989497SBarry Smith     /* determine if matrix has all real values */
67647989497SBarry Smith     v = a->v;
67747989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
678e20fef11SSatish Balay       if (PetscImaginary(v[i])) { allreal = 0; break ;}
67947989497SBarry Smith     }
68047989497SBarry Smith #endif
681932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
682932b0c3eSLois Curfman McInnes       v = a->v + i;
683932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6843a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
68547989497SBarry Smith         if (allreal) {
6863f6de6efSSatish Balay           fprintf(fd,"%6.4e ",PetscReal(*v)); v += a->m;
68747989497SBarry Smith         } else {
688e20fef11SSatish Balay           fprintf(fd,"%6.4e + %6.4e i ",PetscReal(*v),PetscImaginary(*v)); v += a->m;
68947989497SBarry Smith         }
690289bc588SBarry Smith #else
691932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
692289bc588SBarry Smith #endif
693289bc588SBarry Smith       }
6948e37dc09SBarry Smith       fprintf(fd,"\n");
695289bc588SBarry Smith     }
696da3a660dSBarry Smith   }
697932b0c3eSLois Curfman McInnes   fflush(fd);
6983a40ed3dSBarry Smith   PetscFunctionReturn(0);
699289bc588SBarry Smith }
700289bc588SBarry Smith 
7015615d1e5SSatish Balay #undef __FUNC__
7025615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary"
703932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
704932b0c3eSLois Curfman McInnes {
705932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
706932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
70790ace30eSBarry Smith   int          format;
70890ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
709932b0c3eSLois Curfman McInnes 
7103a40ed3dSBarry Smith   PetscFunctionBegin;
71190ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
71290ace30eSBarry Smith 
71390ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
71402cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
71590ace30eSBarry Smith     /* store the matrix as a dense matrix */
71690ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
71790ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
71890ace30eSBarry Smith     col_lens[1] = m;
71990ace30eSBarry Smith     col_lens[2] = n;
72090ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
7210752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr);
72290ace30eSBarry Smith     PetscFree(col_lens);
72390ace30eSBarry Smith 
72490ace30eSBarry Smith     /* write out matrix, by rows */
72545d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
72690ace30eSBarry Smith     v    = a->v;
72790ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
72890ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
72990ace30eSBarry Smith         vals[i + j*m] = *v++;
73090ace30eSBarry Smith       }
73190ace30eSBarry Smith     }
7320752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr);
73390ace30eSBarry Smith     PetscFree(vals);
73490ace30eSBarry Smith   } else {
7350452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
736932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
737932b0c3eSLois Curfman McInnes     col_lens[1] = m;
738932b0c3eSLois Curfman McInnes     col_lens[2] = n;
739932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
740932b0c3eSLois Curfman McInnes 
741932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
742932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
7430752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr);
744932b0c3eSLois Curfman McInnes 
745932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
746932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
747932b0c3eSLois Curfman McInnes     ict = 0;
748932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
749932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
750932b0c3eSLois Curfman McInnes     }
7510752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr);
7520452661fSBarry Smith     PetscFree(col_lens);
753932b0c3eSLois Curfman McInnes 
754932b0c3eSLois Curfman McInnes     /* store nonzero values */
75545d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
756932b0c3eSLois Curfman McInnes     ict = 0;
757932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
758932b0c3eSLois Curfman McInnes       v = a->v + i;
759932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
760932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
761932b0c3eSLois Curfman McInnes       }
762932b0c3eSLois Curfman McInnes     }
7630752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr);
7640452661fSBarry Smith     PetscFree(anonz);
76590ace30eSBarry Smith   }
7663a40ed3dSBarry Smith   PetscFunctionReturn(0);
767932b0c3eSLois Curfman McInnes }
768932b0c3eSLois Curfman McInnes 
7695615d1e5SSatish Balay #undef __FUNC__
7705615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense"
771c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer)
772932b0c3eSLois Curfman McInnes {
773932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
774bcd2baecSBarry Smith   ViewerType   vtype;
775bcd2baecSBarry Smith   int          ierr;
776932b0c3eSLois Curfman McInnes 
7773a40ed3dSBarry Smith   PetscFunctionBegin;
778bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
779bcd2baecSBarry Smith 
7807b2a1423SBarry Smith   if (PetscTypeCompare(vtype,SOCKET_VIEWER)) {
7817b2a1423SBarry Smith     ierr = ViewerSocketPutScalar_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr);
7823f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
7833a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
7843f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) {
7853a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
7865cd90555SBarry Smith   } else {
7875cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
788932b0c3eSLois Curfman McInnes   }
7893a40ed3dSBarry Smith   PetscFunctionReturn(0);
790932b0c3eSLois Curfman McInnes }
791289bc588SBarry Smith 
7925615d1e5SSatish Balay #undef __FUNC__
7935615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense"
794c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
795289bc588SBarry Smith {
796ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
79790f02eecSBarry Smith   int          ierr;
79890f02eecSBarry Smith 
7993a40ed3dSBarry Smith   PetscFunctionBegin;
80094d884c6SBarry Smith   if (--mat->refct > 0) PetscFunctionReturn(0);
80194d884c6SBarry Smith 
80294d884c6SBarry Smith   if (mat->mapping) {
80394d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
80494d884c6SBarry Smith   }
80594d884c6SBarry Smith   if (mat->bmapping) {
80694d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr);
80794d884c6SBarry Smith   }
8083a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
809c6e7dd08SBarry Smith   PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n);
810a5a9c739SBarry Smith #endif
8110452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
8123439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
8130452661fSBarry Smith   PetscFree(l);
814a5ae1ecdSBarry Smith   if (mat->rmap) {
815a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
816a5ae1ecdSBarry Smith   }
817a5ae1ecdSBarry Smith   if (mat->cmap) {
818a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
81990f02eecSBarry Smith   }
820a5a9c739SBarry Smith   PLogObjectDestroy(mat);
8210452661fSBarry Smith   PetscHeaderDestroy(mat);
8223a40ed3dSBarry Smith   PetscFunctionReturn(0);
823289bc588SBarry Smith }
824289bc588SBarry Smith 
8255615d1e5SSatish Balay #undef __FUNC__
8265615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense"
8278f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
828289bc588SBarry Smith {
829c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
830d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
831d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
83248b35521SBarry Smith 
8333a40ed3dSBarry Smith   PetscFunctionBegin;
834d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
8353638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
836d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
837d64ed03dSBarry Smith       Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w);
838d64ed03dSBarry Smith       for ( j=0; j<m; j++ ) {
839d64ed03dSBarry Smith         for ( k=0; k<n; k++ ) {
840d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
841d64ed03dSBarry Smith         }
842d64ed03dSBarry Smith       }
843d64ed03dSBarry Smith       PetscMemcpy(v,w,m*n*sizeof(Scalar));
844d64ed03dSBarry Smith       PetscFree(w);
845d64ed03dSBarry Smith     } else {
846d3e5ee88SLois Curfman McInnes       for ( j=0; j<m; j++ ) {
847289bc588SBarry Smith         for ( k=0; k<j; k++ ) {
848d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
849d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
850d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
851289bc588SBarry Smith         }
852289bc588SBarry Smith       }
853d64ed03dSBarry Smith     }
8543a40ed3dSBarry Smith   } else { /* out-of-place transpose */
855d3e5ee88SLois Curfman McInnes     int          ierr;
856d3e5ee88SLois Curfman McInnes     Mat          tmat;
857ec8511deSBarry Smith     Mat_SeqDense *tmatd;
858d3e5ee88SLois Curfman McInnes     Scalar       *v2;
859b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
860ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
8610de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
862d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
8630de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
864d3e5ee88SLois Curfman McInnes     }
8656d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8666d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
867d3e5ee88SLois Curfman McInnes     *matout = tmat;
86848b35521SBarry Smith   }
8693a40ed3dSBarry Smith   PetscFunctionReturn(0);
870289bc588SBarry Smith }
871289bc588SBarry Smith 
8725615d1e5SSatish Balay #undef __FUNC__
8735615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense"
8748f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
875289bc588SBarry Smith {
876c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
877c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
878289bc588SBarry Smith   int          i;
879289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
8809ea5d5aeSSatish Balay 
8813a40ed3dSBarry Smith   PetscFunctionBegin;
8823a40ed3dSBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
8833a40ed3dSBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
8843a40ed3dSBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
885289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
8863a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
887289bc588SBarry Smith     v1++; v2++;
888289bc588SBarry Smith   }
88977c4ece6SBarry Smith   *flg = PETSC_TRUE;
8903a40ed3dSBarry Smith   PetscFunctionReturn(0);
891289bc588SBarry Smith }
892289bc588SBarry Smith 
8935615d1e5SSatish Balay #undef __FUNC__
8945615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense"
8958f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
896289bc588SBarry Smith {
897c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
8987a97a34bSBarry Smith   int          ierr,i, n, len;
89944cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
90044cd7ae7SLois Curfman McInnes 
9013a40ed3dSBarry Smith   PetscFunctionBegin;
9027a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
9037a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
904600d6d0bSBarry Smith   ierr = VecGetArray(v,&x); CHKERRQ(ierr);
90544cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
906e3372554SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
90744cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
908289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
909289bc588SBarry Smith   }
9107a97a34bSBarry Smith   ierr = VecRestoreArray(v,&x); CHKERRQ(ierr);
9113a40ed3dSBarry Smith   PetscFunctionReturn(0);
912289bc588SBarry Smith }
913289bc588SBarry Smith 
9145615d1e5SSatish Balay #undef __FUNC__
9155615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense"
9168f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
917289bc588SBarry Smith {
918c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
919da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
9207a97a34bSBarry Smith   int          ierr,i,j,m = mat->m, n = mat->n;
92155659b69SBarry Smith 
9223a40ed3dSBarry Smith   PetscFunctionBegin;
92328988994SBarry Smith   if (ll) {
9247a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
925600d6d0bSBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
926e3372554SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
927da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
928da3a660dSBarry Smith       x = l[i];
929da3a660dSBarry Smith       v = mat->v + i;
930da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
931da3a660dSBarry Smith     }
9327a97a34bSBarry Smith     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
93344cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
934da3a660dSBarry Smith   }
93528988994SBarry Smith   if (rr) {
9367a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
937600d6d0bSBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
938e3372554SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
939da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
940da3a660dSBarry Smith       x = r[i];
941da3a660dSBarry Smith       v = mat->v + i*m;
942da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
943da3a660dSBarry Smith     }
9447a97a34bSBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
94544cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
946da3a660dSBarry Smith   }
9473a40ed3dSBarry Smith   PetscFunctionReturn(0);
948289bc588SBarry Smith }
949289bc588SBarry Smith 
9505615d1e5SSatish Balay #undef __FUNC__
9515615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense"
9528f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm)
953289bc588SBarry Smith {
954c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
955289bc588SBarry Smith   Scalar       *v = mat->v;
956289bc588SBarry Smith   double       sum = 0.0;
957289bc588SBarry Smith   int          i, j;
95855659b69SBarry Smith 
9593a40ed3dSBarry Smith   PetscFunctionBegin;
960289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
961289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
9623a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
9633f6de6efSSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
964289bc588SBarry Smith #else
965289bc588SBarry Smith       sum += (*v)*(*v); v++;
966289bc588SBarry Smith #endif
967289bc588SBarry Smith     }
968289bc588SBarry Smith     *norm = sqrt(sum);
96955659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
9703a40ed3dSBarry Smith   } else if (type == NORM_1) {
971289bc588SBarry Smith     *norm = 0.0;
972289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
973289bc588SBarry Smith       sum = 0.0;
974289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
97533a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
976289bc588SBarry Smith       }
977289bc588SBarry Smith       if (sum > *norm) *norm = sum;
978289bc588SBarry Smith     }
97955659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9803a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
981289bc588SBarry Smith     *norm = 0.0;
982289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
983289bc588SBarry Smith       v = mat->v + j;
984289bc588SBarry Smith       sum = 0.0;
985289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
986cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
987289bc588SBarry Smith       }
988289bc588SBarry Smith       if (sum > *norm) *norm = sum;
989289bc588SBarry Smith     }
99055659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9913a40ed3dSBarry Smith   } else {
992e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
993289bc588SBarry Smith   }
9943a40ed3dSBarry Smith   PetscFunctionReturn(0);
995289bc588SBarry Smith }
996289bc588SBarry Smith 
9975615d1e5SSatish Balay #undef __FUNC__
9985615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense"
9998f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1000289bc588SBarry Smith {
1001c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
100267e560aaSBarry Smith 
10033a40ed3dSBarry Smith   PetscFunctionBegin;
10046d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
10056d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
10066d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
1007219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
10086d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
1009219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
10106d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
10116d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
10126d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
10136d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
1014*4787f768SSatish Balay            op == MAT_NEW_NONZERO_LOCATION_ERR ||
10156d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
101690f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
1017b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
1018b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
1019981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
10203a40ed3dSBarry Smith   else {
10213a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
10223a40ed3dSBarry Smith   }
10233a40ed3dSBarry Smith   PetscFunctionReturn(0);
1024289bc588SBarry Smith }
1025289bc588SBarry Smith 
10265615d1e5SSatish Balay #undef __FUNC__
10275615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense"
10288f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
10296f0a148fSBarry Smith {
1030ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
10313a40ed3dSBarry Smith 
10323a40ed3dSBarry Smith   PetscFunctionBegin;
1033cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
10343a40ed3dSBarry Smith   PetscFunctionReturn(0);
10356f0a148fSBarry Smith }
10366f0a148fSBarry Smith 
10375615d1e5SSatish Balay #undef __FUNC__
10385615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense"
10398f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
10404e220ebcSLois Curfman McInnes {
10413a40ed3dSBarry Smith   PetscFunctionBegin;
10424e220ebcSLois Curfman McInnes   *bs = 1;
10433a40ed3dSBarry Smith   PetscFunctionReturn(0);
10444e220ebcSLois Curfman McInnes }
10454e220ebcSLois Curfman McInnes 
10465615d1e5SSatish Balay #undef __FUNC__
10475615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense"
10488f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
10496f0a148fSBarry Smith {
1050ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
10516abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
10526f0a148fSBarry Smith   Scalar       *slot;
105355659b69SBarry Smith 
10543a40ed3dSBarry Smith   PetscFunctionBegin;
105577c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
105678b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
10576f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
10586f0a148fSBarry Smith     slot = l->v + rows[i];
10596f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
10606f0a148fSBarry Smith   }
10616f0a148fSBarry Smith   if (diag) {
10626f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
10636f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1064260925b8SBarry Smith       *slot = *diag;
10656f0a148fSBarry Smith     }
10666f0a148fSBarry Smith   }
1067260925b8SBarry Smith   ISRestoreIndices(is,&rows);
10683a40ed3dSBarry Smith   PetscFunctionReturn(0);
10696f0a148fSBarry Smith }
1070557bce09SLois Curfman McInnes 
10715615d1e5SSatish Balay #undef __FUNC__
10725615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense"
10738f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n)
1074557bce09SLois Curfman McInnes {
1075c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10763a40ed3dSBarry Smith 
10773a40ed3dSBarry Smith   PetscFunctionBegin;
1078557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
10793a40ed3dSBarry Smith   PetscFunctionReturn(0);
1080557bce09SLois Curfman McInnes }
1081557bce09SLois Curfman McInnes 
10825615d1e5SSatish Balay #undef __FUNC__
10835615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense"
10848f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1085d3e5ee88SLois Curfman McInnes {
1086c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10873a40ed3dSBarry Smith 
10883a40ed3dSBarry Smith   PetscFunctionBegin;
1089d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
10903a40ed3dSBarry Smith   PetscFunctionReturn(0);
1091d3e5ee88SLois Curfman McInnes }
1092d3e5ee88SLois Curfman McInnes 
10935615d1e5SSatish Balay #undef __FUNC__
10945615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense"
10958f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array)
109664e87e97SBarry Smith {
1097c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10983a40ed3dSBarry Smith 
10993a40ed3dSBarry Smith   PetscFunctionBegin;
110064e87e97SBarry Smith   *array = mat->v;
11013a40ed3dSBarry Smith   PetscFunctionReturn(0);
110264e87e97SBarry Smith }
11030754003eSLois Curfman McInnes 
11045615d1e5SSatish Balay #undef __FUNC__
11055615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense"
11068f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1107ff14e315SSatish Balay {
11083a40ed3dSBarry Smith   PetscFunctionBegin;
11093a40ed3dSBarry Smith   PetscFunctionReturn(0);
1110ff14e315SSatish Balay }
11110754003eSLois Curfman McInnes 
11125615d1e5SSatish Balay #undef __FUNC__
11135615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense"
11147b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
11150754003eSLois Curfman McInnes {
1116c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1117182d2002SSatish Balay   int          i, j, ierr, m = mat->m, *irow, *icol, nrows, ncols;
1118182d2002SSatish Balay   Scalar       *av, *bv, *v = mat->v;
11190754003eSLois Curfman McInnes   Mat          newmat;
11200754003eSLois Curfman McInnes 
11213a40ed3dSBarry Smith   PetscFunctionBegin;
112278b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
112378b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
112478b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
112578b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
11260754003eSLois Curfman McInnes 
1127182d2002SSatish Balay   /* Check submatrixcall */
1128182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1129182d2002SSatish Balay     int n_cols,n_rows;
1130182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1131182d2002SSatish Balay     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1132182d2002SSatish Balay     newmat = *B;
1133182d2002SSatish Balay   } else {
11340754003eSLois Curfman McInnes     /* Create and fill new matrix */
1135b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
1136182d2002SSatish Balay   }
1137182d2002SSatish Balay 
1138182d2002SSatish Balay   /* Now extract the data pointers and do the copy, column at a time */
1139182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1140182d2002SSatish Balay 
1141182d2002SSatish Balay   for ( i=0; i<ncols; i++ ) {
1142182d2002SSatish Balay     av = v + m*icol[i];
1143182d2002SSatish Balay     for (j=0; j<nrows; j++ ) {
1144182d2002SSatish Balay       *bv++ = av[irow[j]];
11450754003eSLois Curfman McInnes     }
11460754003eSLois Curfman McInnes   }
1147182d2002SSatish Balay 
1148182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
11496d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11506d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11510754003eSLois Curfman McInnes 
11520754003eSLois Curfman McInnes   /* Free work space */
115378b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
115478b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1155182d2002SSatish Balay   *B = newmat;
11563a40ed3dSBarry Smith   PetscFunctionReturn(0);
11570754003eSLois Curfman McInnes }
11580754003eSLois Curfman McInnes 
11595615d1e5SSatish Balay #undef __FUNC__
11605615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense"
11617b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B)
1162905e6a2fSBarry Smith {
1163905e6a2fSBarry Smith   int ierr,i;
1164905e6a2fSBarry Smith 
11653a40ed3dSBarry Smith   PetscFunctionBegin;
1166905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1167905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1168905e6a2fSBarry Smith   }
1169905e6a2fSBarry Smith 
1170905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
11716a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1172905e6a2fSBarry Smith   }
11733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1174905e6a2fSBarry Smith }
1175905e6a2fSBarry Smith 
11765615d1e5SSatish Balay #undef __FUNC__
11775615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense"
1178cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A, Mat B,MatStructure str)
11794b0e389bSBarry Smith {
11804b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
11813a40ed3dSBarry Smith   int          ierr;
11823a40ed3dSBarry Smith 
11833a40ed3dSBarry Smith   PetscFunctionBegin;
11843a40ed3dSBarry Smith   if (B->type != MATSEQDENSE) {
1185cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
11863a40ed3dSBarry Smith     PetscFunctionReturn(0);
11873a40ed3dSBarry Smith   }
1188e3372554SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
11894b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
11903a40ed3dSBarry Smith   PetscFunctionReturn(0);
11914b0e389bSBarry Smith }
11924b0e389bSBarry Smith 
1193289bc588SBarry Smith /* -------------------------------------------------------------------*/
1194a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1195905e6a2fSBarry Smith        MatGetRow_SeqDense,
1196905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1197905e6a2fSBarry Smith        MatMult_SeqDense,
1198905e6a2fSBarry Smith        MatMultAdd_SeqDense,
1199905e6a2fSBarry Smith        MatMultTrans_SeqDense,
1200905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
1201905e6a2fSBarry Smith        MatSolve_SeqDense,
1202905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
1203905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
1204905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
1205905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1206905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1207ec8511deSBarry Smith        MatRelax_SeqDense,
1208ec8511deSBarry Smith        MatTranspose_SeqDense,
1209905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1210905e6a2fSBarry Smith        MatEqual_SeqDense,
1211905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1212905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1213905e6a2fSBarry Smith        MatNorm_SeqDense,
1214905e6a2fSBarry Smith        0,
1215905e6a2fSBarry Smith        0,
1216905e6a2fSBarry Smith        0,
1217905e6a2fSBarry Smith        MatSetOption_SeqDense,
1218905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1219905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1220905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1221905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1222905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1223905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1224905e6a2fSBarry Smith        MatGetSize_SeqDense,
1225905e6a2fSBarry Smith        MatGetSize_SeqDense,
1226905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1227905e6a2fSBarry Smith        0,
1228905e6a2fSBarry Smith        0,
1229905e6a2fSBarry Smith        MatGetArray_SeqDense,
1230905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
12315609ef8eSBarry Smith        MatDuplicate_SeqDense,
1232a5ae1ecdSBarry Smith        0,
1233a5ae1ecdSBarry Smith        0,
1234a5ae1ecdSBarry Smith        0,
1235a5ae1ecdSBarry Smith        0,
1236a5ae1ecdSBarry Smith        MatAXPY_SeqDense,
1237a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1238a5ae1ecdSBarry Smith        0,
12394b0e389bSBarry Smith        MatGetValues_SeqDense,
1240a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1241a5ae1ecdSBarry Smith        0,
1242a5ae1ecdSBarry Smith        MatScale_SeqDense,
1243a5ae1ecdSBarry Smith        0,
1244a5ae1ecdSBarry Smith        0,
1245a5ae1ecdSBarry Smith        0,
1246a5ae1ecdSBarry Smith        MatGetBlockSize_SeqDense,
1247a5ae1ecdSBarry Smith        0,
1248a5ae1ecdSBarry Smith        0,
1249a5ae1ecdSBarry Smith        0,
1250a5ae1ecdSBarry Smith        0,
1251a5ae1ecdSBarry Smith        0,
1252a5ae1ecdSBarry Smith        0,
1253a5ae1ecdSBarry Smith        0,
1254a5ae1ecdSBarry Smith        0,
1255a5ae1ecdSBarry Smith        0,
1256a5ae1ecdSBarry Smith        0,
1257a5ae1ecdSBarry Smith        0,
1258a5ae1ecdSBarry Smith        0,
1259a5ae1ecdSBarry Smith        MatGetMaps_Petsc};
126090ace30eSBarry Smith 
12615615d1e5SSatish Balay #undef __FUNC__
12625615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense"
12634b828684SBarry Smith /*@C
1264fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1265d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1266d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1267289bc588SBarry Smith 
1268db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1269db81eaa0SLois Curfman McInnes 
127020563c6bSBarry Smith    Input Parameters:
1271db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
12720c775827SLois Curfman McInnes .  m - number of rows
127318f449edSLois Curfman McInnes .  n - number of columns
1274db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1275dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
127620563c6bSBarry Smith 
127720563c6bSBarry Smith    Output Parameter:
127844cd7ae7SLois Curfman McInnes .  A - the matrix
127920563c6bSBarry Smith 
1280b259b22eSLois Curfman McInnes    Notes:
128118f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
128218f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1283b4fd4287SBarry Smith    set data=PETSC_NULL.
128418f449edSLois Curfman McInnes 
1285027ccd11SLois Curfman McInnes    Level: intermediate
1286027ccd11SLois Curfman McInnes 
1287dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1288d65003e9SLois Curfman McInnes 
1289db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
129020563c6bSBarry Smith @*/
129144cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1292289bc588SBarry Smith {
129344cd7ae7SLois Curfman McInnes   Mat          B;
129444cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
12953b2fbd54SBarry Smith   int          ierr,flg,size;
12963b2fbd54SBarry Smith 
12973a40ed3dSBarry Smith   PetscFunctionBegin;
12983b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1299a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
130055659b69SBarry Smith 
130144cd7ae7SLois Curfman McInnes   *A            = 0;
13023f1db9ecSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,"Mat",comm,MatDestroy,MatView);
130344cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
130444cd7ae7SLois Curfman McInnes   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
130544cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqDense));
1306a5ae1ecdSBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1307c6e7dd08SBarry Smith   B->ops->destroy    = MatDestroy_SeqDense;
1308c6e7dd08SBarry Smith   B->ops->view       = MatView_SeqDense;
130944cd7ae7SLois Curfman McInnes   B->factor          = 0;
131090f02eecSBarry Smith   B->mapping         = 0;
1311f09e8eb9SSatish Balay   PLogObjectMemory(B,sizeof(struct _p_Mat));
131244cd7ae7SLois Curfman McInnes   B->data            = (void *) b;
131318f449edSLois Curfman McInnes 
131444cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
131544cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
1316a5ae1ecdSBarry Smith 
1317488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
1318488ecbafSBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1319a5ae1ecdSBarry Smith 
132044cd7ae7SLois Curfman McInnes   b->pivots       = 0;
132144cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1322b4fd4287SBarry Smith   if (data == PETSC_NULL) {
132344cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
132444cd7ae7SLois Curfman McInnes     PetscMemzero(b->v,m*n*sizeof(Scalar));
132544cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
132644cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
1327273fdc76SBarry Smith   } else { /* user-allocated storage */
132844cd7ae7SLois Curfman McInnes     b->v = data;
132944cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
13302dd2b441SLois Curfman McInnes   }
133125cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
133244cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
13334e220ebcSLois Curfman McInnes 
133444cd7ae7SLois Curfman McInnes   *A = B;
13353a40ed3dSBarry Smith   PetscFunctionReturn(0);
1336289bc588SBarry Smith }
1337289bc588SBarry Smith 
13383b2fbd54SBarry Smith 
13393b2fbd54SBarry Smith 
13403b2fbd54SBarry Smith 
13413b2fbd54SBarry Smith 
13423b2fbd54SBarry Smith 
1343