xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 488ecbaffa467bb032d31d7eb20bc6d0ef6d986c)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*488ecbafSBarry Smith static char vcid[] = "$Id: dense.c,v 1.153 1998/07/14 02:49:23 bsmith Exp bsmith $";
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"
10b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
11289bc588SBarry Smith 
125615d1e5SSatish Balay #undef __FUNC__
135615d1e5SSatish Balay #define __FUNC__ "MatAXPY_SeqDense"
141987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
151987afe7SBarry Smith {
161987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
171987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
183a40ed3dSBarry Smith 
193a40ed3dSBarry Smith   PetscFunctionBegin;
201987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
211987afe7SBarry Smith   PLogFlops(2*N-1);
223a40ed3dSBarry Smith   PetscFunctionReturn(0);
231987afe7SBarry Smith }
241987afe7SBarry Smith 
255615d1e5SSatish Balay #undef __FUNC__
265615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqDense"
278f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
28289bc588SBarry Smith {
294e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
304e220ebcSLois Curfman McInnes   int          i,N = a->m*a->n,count = 0;
314e220ebcSLois Curfman McInnes   Scalar       *v = a->v;
323a40ed3dSBarry Smith 
333a40ed3dSBarry Smith   PetscFunctionBegin;
34289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
354e220ebcSLois Curfman McInnes 
364e220ebcSLois Curfman McInnes   info->rows_global       = (double)a->m;
374e220ebcSLois Curfman McInnes   info->columns_global    = (double)a->n;
384e220ebcSLois Curfman McInnes   info->rows_local        = (double)a->m;
394e220ebcSLois Curfman McInnes   info->columns_local     = (double)a->n;
404e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
414e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
424e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
434e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
444e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
454e220ebcSLois Curfman McInnes   info->mallocs           = 0;
464e220ebcSLois Curfman McInnes   info->memory            = A->mem;
474e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
484e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
494e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
504e220ebcSLois Curfman McInnes 
513a40ed3dSBarry Smith   PetscFunctionReturn(0);
52289bc588SBarry Smith }
53289bc588SBarry Smith 
545615d1e5SSatish Balay #undef __FUNC__
555615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqDense"
568f6be9afSLois Curfman McInnes int MatScale_SeqDense(Scalar *alpha,Mat inA)
5780cd9d93SLois Curfman McInnes {
5880cd9d93SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
5980cd9d93SLois Curfman McInnes   int          one = 1, nz;
6080cd9d93SLois Curfman McInnes 
613a40ed3dSBarry Smith   PetscFunctionBegin;
6280cd9d93SLois Curfman McInnes   nz = a->m*a->n;
6380cd9d93SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
6480cd9d93SLois Curfman McInnes   PLogFlops(nz);
653a40ed3dSBarry Smith   PetscFunctionReturn(0);
6680cd9d93SLois Curfman McInnes }
6780cd9d93SLois Curfman McInnes 
68289bc588SBarry Smith /* ---------------------------------------------------------------*/
69289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
70289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
715615d1e5SSatish Balay #undef __FUNC__
725615d1e5SSatish Balay #define __FUNC__ "MatLUFactor_SeqDense"
738f6be9afSLois Curfman McInnes int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
74289bc588SBarry Smith {
75c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
766abc6512SBarry Smith   int          info;
7755659b69SBarry Smith 
783a40ed3dSBarry Smith   PetscFunctionBegin;
79289bc588SBarry Smith   if (!mat->pivots) {
8045d48df9SBarry Smith     mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots);
81c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
82289bc588SBarry Smith   }
837a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
84289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
85a8c6a408SBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,0,"Bad argument to LU factorization");
86e3372554SBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Bad LU factorization");
87c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
8855659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
893a40ed3dSBarry Smith   PetscFunctionReturn(0);
90289bc588SBarry Smith }
916ee01492SSatish Balay 
925615d1e5SSatish Balay #undef __FUNC__
935615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqDense"
948f6be9afSLois Curfman McInnes int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
9502cad45dSBarry Smith {
9602cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
9702cad45dSBarry Smith   int          ierr;
9802cad45dSBarry Smith   Mat          newi;
9902cad45dSBarry Smith 
1003a40ed3dSBarry Smith   PetscFunctionBegin;
10102cad45dSBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
10202cad45dSBarry Smith   l = (Mat_SeqDense *) newi->data;
10302cad45dSBarry Smith   if (cpvalues == COPY_VALUES) {
10402cad45dSBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
10502cad45dSBarry Smith   }
10602cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
10702cad45dSBarry Smith   *newmat = newi;
1083a40ed3dSBarry Smith   PetscFunctionReturn(0);
10902cad45dSBarry Smith }
11002cad45dSBarry Smith 
1115615d1e5SSatish Balay #undef __FUNC__
1125615d1e5SSatish Balay #define __FUNC__ "MatLUFactorSymbolic_SeqDense"
1138f6be9afSLois Curfman McInnes int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
114289bc588SBarry Smith {
1153a40ed3dSBarry Smith   int ierr;
1163a40ed3dSBarry Smith 
1173a40ed3dSBarry Smith   PetscFunctionBegin;
1183a40ed3dSBarry Smith   ierr = MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);CHKERRQ(ierr);
1193a40ed3dSBarry Smith   PetscFunctionReturn(0);
120289bc588SBarry Smith }
1216ee01492SSatish Balay 
1225615d1e5SSatish Balay #undef __FUNC__
1235615d1e5SSatish Balay #define __FUNC__ "MatLUFactorNumeric_SeqDense"
1248f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
125289bc588SBarry Smith {
12602cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data;
1273a40ed3dSBarry Smith   int          ierr;
1283a40ed3dSBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionBegin;
13002cad45dSBarry Smith   /* copy the numerical values */
13102cad45dSBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
13202cad45dSBarry Smith   (*fact)->factor = 0;
1333a40ed3dSBarry Smith   ierr = MatLUFactor(*fact,0,0,1.0);CHKERRQ(ierr);
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
135289bc588SBarry Smith }
1366ee01492SSatish Balay 
1375615d1e5SSatish Balay #undef __FUNC__
1385615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense"
1398f6be9afSLois Curfman McInnes int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
140289bc588SBarry Smith {
1413a40ed3dSBarry Smith   int ierr;
1423a40ed3dSBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
1443a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1453a40ed3dSBarry Smith   PetscFunctionReturn(0);
146289bc588SBarry Smith }
1476ee01492SSatish Balay 
1485615d1e5SSatish Balay #undef __FUNC__
1495615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense"
1508f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
151289bc588SBarry Smith {
1523a40ed3dSBarry Smith   int ierr;
1533a40ed3dSBarry Smith 
1543a40ed3dSBarry Smith   PetscFunctionBegin;
1553a40ed3dSBarry Smith   ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr);
1563a40ed3dSBarry Smith   PetscFunctionReturn(0);
157289bc588SBarry Smith }
1586ee01492SSatish Balay 
1595615d1e5SSatish Balay #undef __FUNC__
1605615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense"
1618f6be9afSLois Curfman McInnes int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
162289bc588SBarry Smith {
163c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1646abc6512SBarry Smith   int           info;
16555659b69SBarry Smith 
1663a40ed3dSBarry Smith   PetscFunctionBegin;
167752f5567SLois Curfman McInnes   if (mat->pivots) {
1680452661fSBarry Smith     PetscFree(mat->pivots);
169c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
170752f5567SLois Curfman McInnes     mat->pivots = 0;
171752f5567SLois Curfman McInnes   }
1727a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
173289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
174e3372554SBarry Smith   if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization");
175c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
17655659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
1773a40ed3dSBarry Smith   PetscFunctionReturn(0);
178289bc588SBarry Smith }
179289bc588SBarry Smith 
1805615d1e5SSatish Balay #undef __FUNC__
1815615d1e5SSatish Balay #define __FUNC__ "MatSolve_SeqDense"
1828f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
183289bc588SBarry Smith {
184c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
185a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
1866abc6512SBarry Smith   Scalar       *x, *y;
18767e560aaSBarry Smith 
1883a40ed3dSBarry Smith   PetscFunctionBegin;
1897a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
1907a97a34bSBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
1917a97a34bSBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
192416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
193c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
19448d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
1957a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
19648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
197289bc588SBarry Smith   }
198a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Matrix must be factored to solve");
199a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"MBad solve");
2007a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
2017a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
20255659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2033a40ed3dSBarry Smith   PetscFunctionReturn(0);
204289bc588SBarry Smith }
2056ee01492SSatish Balay 
2065615d1e5SSatish Balay #undef __FUNC__
2075615d1e5SSatish Balay #define __FUNC__ "MatSolveTrans_SeqDense"
2088f6be9afSLois Curfman McInnes int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
209da3a660dSBarry Smith {
210c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2117a97a34bSBarry Smith   int          ierr,one = 1, info;
2126abc6512SBarry Smith   Scalar       *x, *y;
21367e560aaSBarry Smith 
2143a40ed3dSBarry Smith   PetscFunctionBegin;
2157a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
2167a97a34bSBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
2177a97a34bSBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
218416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
219752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
220da3a660dSBarry Smith   if (mat->pivots) {
22148d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
2227a97a34bSBarry Smith   } else {
22348d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
224da3a660dSBarry Smith   }
225a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
2267a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
2277a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
22855659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2293a40ed3dSBarry Smith   PetscFunctionReturn(0);
230da3a660dSBarry Smith }
2316ee01492SSatish Balay 
2325615d1e5SSatish Balay #undef __FUNC__
2335615d1e5SSatish Balay #define __FUNC__ "MatSolveAdd_SeqDense"
2348f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
235da3a660dSBarry Smith {
236c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2377a97a34bSBarry Smith   int          ierr,one = 1, info;
2386abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
239da3a660dSBarry Smith   Vec          tmp = 0;
24067e560aaSBarry Smith 
2413a40ed3dSBarry Smith   PetscFunctionBegin;
2427a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);  CHKERRQ(ierr);
2437a97a34bSBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
2447a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
245da3a660dSBarry Smith   if (yy == zz) {
24678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
247c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
24878b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
249da3a660dSBarry Smith   }
250416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
251752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
252da3a660dSBarry Smith   if (mat->pivots) {
25348d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
254a8c6a408SBarry Smith   } else {
25548d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
256da3a660dSBarry Smith   }
257a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
258a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy); CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
259a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);}
2607a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
2617a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
26255659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2633a40ed3dSBarry Smith   PetscFunctionReturn(0);
264da3a660dSBarry Smith }
26567e560aaSBarry Smith 
2665615d1e5SSatish Balay #undef __FUNC__
2675615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense"
2688f6be9afSLois Curfman McInnes int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
269da3a660dSBarry Smith {
270c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
2716abc6512SBarry Smith   int           one = 1, info,ierr;
2726abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
273da3a660dSBarry Smith   Vec           tmp;
27467e560aaSBarry Smith 
2753a40ed3dSBarry Smith   PetscFunctionBegin;
2767a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
2777a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2787a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
279da3a660dSBarry Smith   if (yy == zz) {
28078b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
281c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
28278b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
283da3a660dSBarry Smith   }
284416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
285752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
286da3a660dSBarry Smith   if (mat->pivots) {
28748d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
2883a40ed3dSBarry Smith   } else {
28948d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
290da3a660dSBarry Smith   }
291a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
29290f02eecSBarry Smith   if (tmp) {
29390f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);  CHKERRQ(ierr);
29490f02eecSBarry Smith     ierr = VecDestroy(tmp); CHKERRQ(ierr);
2953a40ed3dSBarry Smith   } else {
29690f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);
29790f02eecSBarry Smith   }
2987a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
2997a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
30055659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
3013a40ed3dSBarry Smith   PetscFunctionReturn(0);
302da3a660dSBarry Smith }
303289bc588SBarry Smith /* ------------------------------------------------------------------*/
3045615d1e5SSatish Balay #undef __FUNC__
3055615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense"
3068f6be9afSLois Curfman McInnes int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
30720e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
308289bc588SBarry Smith {
309c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
310289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
311bc1b551cSSatish Balay   int          ierr, m = mat->m, i;
312bc1b551cSSatish Balay #if !defined(USE_PETSC_COMPLEX)
313bc1b551cSSatish Balay   int          o = 1;
314bc1b551cSSatish Balay #endif
315289bc588SBarry Smith 
3163a40ed3dSBarry Smith   PetscFunctionBegin;
317289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
3183a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
319bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
320289bc588SBarry Smith   }
3217a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3227a97a34bSBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
323289bc588SBarry Smith   while (its--) {
324289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
325289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
3263a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
327f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
328f1747703SBarry Smith            not happy about returning a double complex */
329f1747703SBarry Smith         int    _i;
330f1747703SBarry Smith         Scalar sum = b[i];
331f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
3323f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
333f1747703SBarry Smith         }
334f1747703SBarry Smith         xt = sum;
335f1747703SBarry Smith #else
336289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
337f1747703SBarry Smith #endif
33855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
339289bc588SBarry Smith       }
340289bc588SBarry Smith     }
341289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
342289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
3433a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
344f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
345f1747703SBarry Smith            not happy about returning a double complex */
346f1747703SBarry Smith         int    _i;
347f1747703SBarry Smith         Scalar sum = b[i];
348f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
3493f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
350f1747703SBarry Smith         }
351f1747703SBarry Smith         xt = sum;
352f1747703SBarry Smith #else
353289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
354f1747703SBarry Smith #endif
35555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
356289bc588SBarry Smith       }
357289bc588SBarry Smith     }
358289bc588SBarry Smith   }
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;
407289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
4087a97a34bSBarry Smith   int          ierr,_One=1; Scalar _DOne=1.0;
4093a40ed3dSBarry Smith 
4103a40ed3dSBarry Smith   PetscFunctionBegin;
4117a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
4127a97a34bSBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
4137a97a34bSBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
4147a97a34bSBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
415416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
41648d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
4177a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
4187a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
4197a97a34bSBarry Smith   ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr);
42055659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
4213a40ed3dSBarry Smith   PetscFunctionReturn(0);
422289bc588SBarry Smith }
4236ee01492SSatish Balay 
4245615d1e5SSatish Balay #undef __FUNC__
4255615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense"
42644cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
427289bc588SBarry Smith {
428c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
429289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
4307a97a34bSBarry Smith   int          ierr,_One=1;
4316abc6512SBarry Smith   Scalar       _DOne=1.0;
4323a40ed3dSBarry Smith 
4333a40ed3dSBarry Smith   PetscFunctionBegin;
4347a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
4357a97a34bSBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
4367a97a34bSBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
4377a97a34bSBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
438717eeb74SLois Curfman McInnes   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
43948d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
4407a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
4417a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
4427a97a34bSBarry Smith   ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr);
44355659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
4443a40ed3dSBarry Smith   PetscFunctionReturn(0);
445289bc588SBarry Smith }
446289bc588SBarry Smith 
447289bc588SBarry Smith /* -----------------------------------------------------------------*/
4485615d1e5SSatish Balay #undef __FUNC__
4495615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense"
4508f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
451289bc588SBarry Smith {
452c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
453289bc588SBarry Smith   Scalar       *v;
454289bc588SBarry Smith   int          i;
45567e560aaSBarry Smith 
4563a40ed3dSBarry Smith   PetscFunctionBegin;
457289bc588SBarry Smith   *ncols = mat->n;
458289bc588SBarry Smith   if (cols) {
45945d48df9SBarry Smith     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols);
46073c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
461289bc588SBarry Smith   }
462289bc588SBarry Smith   if (vals) {
46345d48df9SBarry Smith     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals);
464289bc588SBarry Smith     v = mat->v + row;
46573c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
466289bc588SBarry Smith   }
4673a40ed3dSBarry Smith   PetscFunctionReturn(0);
468289bc588SBarry Smith }
4696ee01492SSatish Balay 
4705615d1e5SSatish Balay #undef __FUNC__
4715615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense"
4728f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
473289bc588SBarry Smith {
4740452661fSBarry Smith   if (cols) { PetscFree(*cols); }
4750452661fSBarry Smith   if (vals) { PetscFree(*vals); }
4763a40ed3dSBarry Smith   PetscFunctionReturn(0);
477289bc588SBarry Smith }
478289bc588SBarry Smith /* ----------------------------------------------------------------*/
4795615d1e5SSatish Balay #undef __FUNC__
4805615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense"
4818f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
482e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
483289bc588SBarry Smith {
484c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
485289bc588SBarry Smith   int          i,j;
486d6dfbf8fSBarry Smith 
4873a40ed3dSBarry Smith   PetscFunctionBegin;
488289bc588SBarry Smith   if (!mat->roworiented) {
489dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
490289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
4913a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
492a8c6a408SBarry Smith         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
493a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
49458804f6dSBarry Smith #endif
495289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
4963a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
497a8c6a408SBarry Smith           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
498a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
49958804f6dSBarry Smith #endif
500289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
501289bc588SBarry Smith         }
502289bc588SBarry Smith       }
5033a40ed3dSBarry Smith     } else {
504289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
5053a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
506a8c6a408SBarry Smith         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
507a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
50858804f6dSBarry Smith #endif
509289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
5103a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
511a8c6a408SBarry Smith           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
512a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
51358804f6dSBarry Smith #endif
514289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
515289bc588SBarry Smith         }
516289bc588SBarry Smith       }
517289bc588SBarry Smith     }
5183a40ed3dSBarry Smith   } else {
519dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
520e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
5213a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
522a8c6a408SBarry Smith         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
523a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
52458804f6dSBarry Smith #endif
525e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
5263a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
527a8c6a408SBarry Smith           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
528a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
52958804f6dSBarry Smith #endif
530e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
531e8d4e0b9SBarry Smith         }
532e8d4e0b9SBarry Smith       }
5333a40ed3dSBarry Smith     } else {
534289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
5353a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
536a8c6a408SBarry Smith         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
537a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
53858804f6dSBarry Smith #endif
539289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
5403a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
541a8c6a408SBarry Smith           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
542a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
54358804f6dSBarry Smith #endif
544289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
545289bc588SBarry Smith         }
546289bc588SBarry Smith       }
547289bc588SBarry Smith     }
548e8d4e0b9SBarry Smith   }
5493a40ed3dSBarry Smith   PetscFunctionReturn(0);
550289bc588SBarry Smith }
551e8d4e0b9SBarry Smith 
5525615d1e5SSatish Balay #undef __FUNC__
5535615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense"
5548f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
555ae80bb75SLois Curfman McInnes {
556ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
557ae80bb75SLois Curfman McInnes   int          i, j;
558ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
559ae80bb75SLois Curfman McInnes 
5603a40ed3dSBarry Smith   PetscFunctionBegin;
561ae80bb75SLois Curfman McInnes   /* row-oriented output */
562ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
563ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
564ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
565ae80bb75SLois Curfman McInnes     }
566ae80bb75SLois Curfman McInnes   }
5673a40ed3dSBarry Smith   PetscFunctionReturn(0);
568ae80bb75SLois Curfman McInnes }
569ae80bb75SLois Curfman McInnes 
570289bc588SBarry Smith /* -----------------------------------------------------------------*/
571289bc588SBarry Smith 
57277c4ece6SBarry Smith #include "sys.h"
57388e32edaSLois Curfman McInnes 
5745615d1e5SSatish Balay #undef __FUNC__
5755615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense"
57619bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
57788e32edaSLois Curfman McInnes {
57888e32edaSLois Curfman McInnes   Mat_SeqDense *a;
57988e32edaSLois Curfman McInnes   Mat          B;
58088e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
58188e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
58288e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
58319bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
58488e32edaSLois Curfman McInnes 
5853a40ed3dSBarry Smith   PetscFunctionBegin;
58688e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
587a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
58890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
5890752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
590a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object");
59188e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
59288e32edaSLois Curfman McInnes 
59390ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
59490ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
59590ace30eSBarry Smith     B    = *A;
59690ace30eSBarry Smith     a    = (Mat_SeqDense *) B->data;
59790ace30eSBarry Smith 
59890ace30eSBarry Smith     /* read in nonzero values */
5990752156aSBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,PETSC_SCALAR); CHKERRQ(ierr);
60090ace30eSBarry Smith 
6016d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6026d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
603d64ed03dSBarry Smith     /* ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); */
60490ace30eSBarry Smith   } else {
60588e32edaSLois Curfman McInnes     /* read row lengths */
60645d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
6070752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
60888e32edaSLois Curfman McInnes 
60988e32edaSLois Curfman McInnes     /* create our matrix */
610b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
61188e32edaSLois Curfman McInnes     B = *A;
61288e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
61388e32edaSLois Curfman McInnes     v = a->v;
61488e32edaSLois Curfman McInnes 
61588e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
61645d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
6170752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
61845d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
6190752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
62088e32edaSLois Curfman McInnes 
62188e32edaSLois Curfman McInnes     /* insert into matrix */
62288e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
62388e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
62488e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
62588e32edaSLois Curfman McInnes     }
6260452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
62788e32edaSLois Curfman McInnes 
6286d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6296d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
63090ace30eSBarry Smith   }
6313a40ed3dSBarry Smith   PetscFunctionReturn(0);
63288e32edaSLois Curfman McInnes }
63388e32edaSLois Curfman McInnes 
634932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
63577c4ece6SBarry Smith #include "sys.h"
636932b0c3eSLois Curfman McInnes 
6375615d1e5SSatish Balay #undef __FUNC__
6385615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII"
639932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
640289bc588SBarry Smith {
641932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
642932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
643d636dbe3SBarry Smith   FILE         *fd;
644932b0c3eSLois Curfman McInnes   char         *outputname;
645932b0c3eSLois Curfman McInnes   Scalar       *v;
646932b0c3eSLois Curfman McInnes 
6473a40ed3dSBarry Smith   PetscFunctionBegin;
64890ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
649932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
65090ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
651639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
6523a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
653f72585f2SLois Curfman McInnes   }
65402cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
65580cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
65644cd7ae7SLois Curfman McInnes       v = a->v + i;
65780cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
65880cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
660e20fef11SSatish Balay         if (PetscReal(*v) != 0.0 && PetscImaginary(*v) != 0.0) fprintf(fd," %d %g + %g i",j,PetscReal(*v),PetscImaginary(*v));
6613f6de6efSSatish Balay         else if (PetscReal(*v)) fprintf(fd," %d %g ",j,PetscReal(*v));
66280cd9d93SLois Curfman McInnes #else
66380cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
66480cd9d93SLois Curfman McInnes #endif
665d64ed03dSBarry Smith         v += a->m;
66680cd9d93SLois Curfman McInnes       }
66780cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
66880cd9d93SLois Curfman McInnes     }
6693a40ed3dSBarry Smith   } else {
6703a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
67147989497SBarry Smith     int allreal = 1;
67247989497SBarry Smith     /* determine if matrix has all real values */
67347989497SBarry Smith     v = a->v;
67447989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
675e20fef11SSatish Balay       if (PetscImaginary(v[i])) { allreal = 0; break ;}
67647989497SBarry Smith     }
67747989497SBarry Smith #endif
678932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
679932b0c3eSLois Curfman McInnes       v = a->v + i;
680932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6813a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
68247989497SBarry Smith         if (allreal) {
6833f6de6efSSatish Balay           fprintf(fd,"%6.4e ",PetscReal(*v)); v += a->m;
68447989497SBarry Smith         } else {
685e20fef11SSatish Balay           fprintf(fd,"%6.4e + %6.4e i ",PetscReal(*v),PetscImaginary(*v)); v += a->m;
68647989497SBarry Smith         }
687289bc588SBarry Smith #else
688932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
689289bc588SBarry Smith #endif
690289bc588SBarry Smith       }
6918e37dc09SBarry Smith       fprintf(fd,"\n");
692289bc588SBarry Smith     }
693da3a660dSBarry Smith   }
694932b0c3eSLois Curfman McInnes   fflush(fd);
6953a40ed3dSBarry Smith   PetscFunctionReturn(0);
696289bc588SBarry Smith }
697289bc588SBarry Smith 
6985615d1e5SSatish Balay #undef __FUNC__
6995615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary"
700932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
701932b0c3eSLois Curfman McInnes {
702932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
703932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
70490ace30eSBarry Smith   int          format;
70590ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
706932b0c3eSLois Curfman McInnes 
7073a40ed3dSBarry Smith   PetscFunctionBegin;
70890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
70990ace30eSBarry Smith 
71090ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
71102cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
71290ace30eSBarry Smith     /* store the matrix as a dense matrix */
71390ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
71490ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
71590ace30eSBarry Smith     col_lens[1] = m;
71690ace30eSBarry Smith     col_lens[2] = n;
71790ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
7180752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr);
71990ace30eSBarry Smith     PetscFree(col_lens);
72090ace30eSBarry Smith 
72190ace30eSBarry Smith     /* write out matrix, by rows */
72245d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
72390ace30eSBarry Smith     v    = a->v;
72490ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
72590ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
72690ace30eSBarry Smith         vals[i + j*m] = *v++;
72790ace30eSBarry Smith       }
72890ace30eSBarry Smith     }
7290752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr);
73090ace30eSBarry Smith     PetscFree(vals);
73190ace30eSBarry Smith   } else {
7320452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
733932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
734932b0c3eSLois Curfman McInnes     col_lens[1] = m;
735932b0c3eSLois Curfman McInnes     col_lens[2] = n;
736932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
737932b0c3eSLois Curfman McInnes 
738932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
739932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
7400752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr);
741932b0c3eSLois Curfman McInnes 
742932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
743932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
744932b0c3eSLois Curfman McInnes     ict = 0;
745932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
746932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
747932b0c3eSLois Curfman McInnes     }
7480752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr);
7490452661fSBarry Smith     PetscFree(col_lens);
750932b0c3eSLois Curfman McInnes 
751932b0c3eSLois Curfman McInnes     /* store nonzero values */
75245d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
753932b0c3eSLois Curfman McInnes     ict = 0;
754932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
755932b0c3eSLois Curfman McInnes       v = a->v + i;
756932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
757932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
758932b0c3eSLois Curfman McInnes       }
759932b0c3eSLois Curfman McInnes     }
7600752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr);
7610452661fSBarry Smith     PetscFree(anonz);
76290ace30eSBarry Smith   }
7633a40ed3dSBarry Smith   PetscFunctionReturn(0);
764932b0c3eSLois Curfman McInnes }
765932b0c3eSLois Curfman McInnes 
7665615d1e5SSatish Balay #undef __FUNC__
7675615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense"
768c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer)
769932b0c3eSLois Curfman McInnes {
770932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
771bcd2baecSBarry Smith   ViewerType   vtype;
772bcd2baecSBarry Smith   int          ierr;
773932b0c3eSLois Curfman McInnes 
7743a40ed3dSBarry Smith   PetscFunctionBegin;
775bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
776bcd2baecSBarry Smith 
777bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
7781a0f753aSSatish Balay     ierr = ViewerMatlabPutScalar_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr);
7793a40ed3dSBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
7803a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
7813a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
7823a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
7835cd90555SBarry Smith   } else {
7845cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
785932b0c3eSLois Curfman McInnes   }
7863a40ed3dSBarry Smith   PetscFunctionReturn(0);
787932b0c3eSLois Curfman McInnes }
788289bc588SBarry Smith 
7895615d1e5SSatish Balay #undef __FUNC__
7905615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense"
791c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
792289bc588SBarry Smith {
793ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
79490f02eecSBarry Smith   int          ierr;
79590f02eecSBarry Smith 
7963a40ed3dSBarry Smith   PetscFunctionBegin;
79794d884c6SBarry Smith   if (--mat->refct > 0) PetscFunctionReturn(0);
79894d884c6SBarry Smith 
79994d884c6SBarry Smith   if (mat->mapping) {
80094d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
80194d884c6SBarry Smith   }
80294d884c6SBarry Smith   if (mat->bmapping) {
80394d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr);
80494d884c6SBarry Smith   }
8053a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
806c6e7dd08SBarry Smith   PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n);
807a5a9c739SBarry Smith #endif
8080452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
8093439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
8100452661fSBarry Smith   PetscFree(l);
811a5ae1ecdSBarry Smith   if (mat->rmap) {
812a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
813a5ae1ecdSBarry Smith   }
814a5ae1ecdSBarry Smith   if (mat->cmap) {
815a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
81690f02eecSBarry Smith   }
817a5a9c739SBarry Smith   PLogObjectDestroy(mat);
8180452661fSBarry Smith   PetscHeaderDestroy(mat);
8193a40ed3dSBarry Smith   PetscFunctionReturn(0);
820289bc588SBarry Smith }
821289bc588SBarry Smith 
8225615d1e5SSatish Balay #undef __FUNC__
8235615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense"
8248f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
825289bc588SBarry Smith {
826c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
827d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
828d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
82948b35521SBarry Smith 
8303a40ed3dSBarry Smith   PetscFunctionBegin;
831d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
8323638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
833d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
834d64ed03dSBarry Smith       Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w);
835d64ed03dSBarry Smith       for ( j=0; j<m; j++ ) {
836d64ed03dSBarry Smith         for ( k=0; k<n; k++ ) {
837d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
838d64ed03dSBarry Smith         }
839d64ed03dSBarry Smith       }
840d64ed03dSBarry Smith       PetscMemcpy(v,w,m*n*sizeof(Scalar));
841d64ed03dSBarry Smith       PetscFree(w);
842d64ed03dSBarry Smith     } else {
843d3e5ee88SLois Curfman McInnes       for ( j=0; j<m; j++ ) {
844289bc588SBarry Smith         for ( k=0; k<j; k++ ) {
845d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
846d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
847d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
848289bc588SBarry Smith         }
849289bc588SBarry Smith       }
850d64ed03dSBarry Smith     }
8513a40ed3dSBarry Smith   } else { /* out-of-place transpose */
852d3e5ee88SLois Curfman McInnes     int          ierr;
853d3e5ee88SLois Curfman McInnes     Mat          tmat;
854ec8511deSBarry Smith     Mat_SeqDense *tmatd;
855d3e5ee88SLois Curfman McInnes     Scalar       *v2;
856b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
857ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
8580de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
859d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
8600de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
861d3e5ee88SLois Curfman McInnes     }
8626d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8636d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
864d3e5ee88SLois Curfman McInnes     *matout = tmat;
86548b35521SBarry Smith   }
8663a40ed3dSBarry Smith   PetscFunctionReturn(0);
867289bc588SBarry Smith }
868289bc588SBarry Smith 
8695615d1e5SSatish Balay #undef __FUNC__
8705615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense"
8718f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
872289bc588SBarry Smith {
873c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
874c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
875289bc588SBarry Smith   int          i;
876289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
8779ea5d5aeSSatish Balay 
8783a40ed3dSBarry Smith   PetscFunctionBegin;
8793a40ed3dSBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
8803a40ed3dSBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
8813a40ed3dSBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
882289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
8833a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
884289bc588SBarry Smith     v1++; v2++;
885289bc588SBarry Smith   }
88677c4ece6SBarry Smith   *flg = PETSC_TRUE;
8873a40ed3dSBarry Smith   PetscFunctionReturn(0);
888289bc588SBarry Smith }
889289bc588SBarry Smith 
8905615d1e5SSatish Balay #undef __FUNC__
8915615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense"
8928f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
893289bc588SBarry Smith {
894c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
8957a97a34bSBarry Smith   int          ierr,i, n, len;
89644cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
89744cd7ae7SLois Curfman McInnes 
8983a40ed3dSBarry Smith   PetscFunctionBegin;
8997a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
9007a97a34bSBarry Smith   ierr = VecGetArray(v,&x); CHKERRQ(ierr);
9017a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
90244cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
903e3372554SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
90444cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
905289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
906289bc588SBarry Smith   }
9077a97a34bSBarry Smith   ierr = VecRestoreArray(v,&x); CHKERRQ(ierr);
9083a40ed3dSBarry Smith   PetscFunctionReturn(0);
909289bc588SBarry Smith }
910289bc588SBarry Smith 
9115615d1e5SSatish Balay #undef __FUNC__
9125615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense"
9138f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
914289bc588SBarry Smith {
915c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
916da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
9177a97a34bSBarry Smith   int          ierr,i,j,m = mat->m, n = mat->n;
91855659b69SBarry Smith 
9193a40ed3dSBarry Smith   PetscFunctionBegin;
92028988994SBarry Smith   if (ll) {
9217a97a34bSBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
9227a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
923e3372554SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
924da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
925da3a660dSBarry Smith       x = l[i];
926da3a660dSBarry Smith       v = mat->v + i;
927da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
928da3a660dSBarry Smith     }
9297a97a34bSBarry Smith     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
93044cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
931da3a660dSBarry Smith   }
93228988994SBarry Smith   if (rr) {
9337a97a34bSBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
9347a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
935e3372554SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
936da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
937da3a660dSBarry Smith       x = r[i];
938da3a660dSBarry Smith       v = mat->v + i*m;
939da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
940da3a660dSBarry Smith     }
9417a97a34bSBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
94244cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
943da3a660dSBarry Smith   }
9443a40ed3dSBarry Smith   PetscFunctionReturn(0);
945289bc588SBarry Smith }
946289bc588SBarry Smith 
9475615d1e5SSatish Balay #undef __FUNC__
9485615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense"
9498f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm)
950289bc588SBarry Smith {
951c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
952289bc588SBarry Smith   Scalar       *v = mat->v;
953289bc588SBarry Smith   double       sum = 0.0;
954289bc588SBarry Smith   int          i, j;
95555659b69SBarry Smith 
9563a40ed3dSBarry Smith   PetscFunctionBegin;
957289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
958289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
9593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
9603f6de6efSSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
961289bc588SBarry Smith #else
962289bc588SBarry Smith       sum += (*v)*(*v); v++;
963289bc588SBarry Smith #endif
964289bc588SBarry Smith     }
965289bc588SBarry Smith     *norm = sqrt(sum);
96655659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
9673a40ed3dSBarry Smith   } else if (type == NORM_1) {
968289bc588SBarry Smith     *norm = 0.0;
969289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
970289bc588SBarry Smith       sum = 0.0;
971289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
97233a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
973289bc588SBarry Smith       }
974289bc588SBarry Smith       if (sum > *norm) *norm = sum;
975289bc588SBarry Smith     }
97655659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9773a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
978289bc588SBarry Smith     *norm = 0.0;
979289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
980289bc588SBarry Smith       v = mat->v + j;
981289bc588SBarry Smith       sum = 0.0;
982289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
983cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
984289bc588SBarry Smith       }
985289bc588SBarry Smith       if (sum > *norm) *norm = sum;
986289bc588SBarry Smith     }
98755659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9883a40ed3dSBarry Smith   } else {
989e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
990289bc588SBarry Smith   }
9913a40ed3dSBarry Smith   PetscFunctionReturn(0);
992289bc588SBarry Smith }
993289bc588SBarry Smith 
9945615d1e5SSatish Balay #undef __FUNC__
9955615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense"
9968f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
997289bc588SBarry Smith {
998c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
99967e560aaSBarry Smith 
10003a40ed3dSBarry Smith   PetscFunctionBegin;
10016d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
10026d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
10036d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
1004219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
10056d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
1006219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
10076d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
10086d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
10096d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
10106d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
1011c2653b3dSLois Curfman McInnes            op == MAT_NEW_NONZERO_LOCATION_ERROR ||
10126d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
101390f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
1014b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
1015b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
1016981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
10173a40ed3dSBarry Smith   else {
10183a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
10193a40ed3dSBarry Smith   }
10203a40ed3dSBarry Smith   PetscFunctionReturn(0);
1021289bc588SBarry Smith }
1022289bc588SBarry Smith 
10235615d1e5SSatish Balay #undef __FUNC__
10245615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense"
10258f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
10266f0a148fSBarry Smith {
1027ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
10283a40ed3dSBarry Smith 
10293a40ed3dSBarry Smith   PetscFunctionBegin;
1030cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
10313a40ed3dSBarry Smith   PetscFunctionReturn(0);
10326f0a148fSBarry Smith }
10336f0a148fSBarry Smith 
10345615d1e5SSatish Balay #undef __FUNC__
10355615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense"
10368f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
10374e220ebcSLois Curfman McInnes {
10383a40ed3dSBarry Smith   PetscFunctionBegin;
10394e220ebcSLois Curfman McInnes   *bs = 1;
10403a40ed3dSBarry Smith   PetscFunctionReturn(0);
10414e220ebcSLois Curfman McInnes }
10424e220ebcSLois Curfman McInnes 
10435615d1e5SSatish Balay #undef __FUNC__
10445615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense"
10458f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
10466f0a148fSBarry Smith {
1047ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
10486abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
10496f0a148fSBarry Smith   Scalar       *slot;
105055659b69SBarry Smith 
10513a40ed3dSBarry Smith   PetscFunctionBegin;
105277c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
105378b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
10546f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
10556f0a148fSBarry Smith     slot = l->v + rows[i];
10566f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
10576f0a148fSBarry Smith   }
10586f0a148fSBarry Smith   if (diag) {
10596f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
10606f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1061260925b8SBarry Smith       *slot = *diag;
10626f0a148fSBarry Smith     }
10636f0a148fSBarry Smith   }
1064260925b8SBarry Smith   ISRestoreIndices(is,&rows);
10653a40ed3dSBarry Smith   PetscFunctionReturn(0);
10666f0a148fSBarry Smith }
1067557bce09SLois Curfman McInnes 
10685615d1e5SSatish Balay #undef __FUNC__
10695615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense"
10708f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n)
1071557bce09SLois Curfman McInnes {
1072c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10733a40ed3dSBarry Smith 
10743a40ed3dSBarry Smith   PetscFunctionBegin;
1075557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
10763a40ed3dSBarry Smith   PetscFunctionReturn(0);
1077557bce09SLois Curfman McInnes }
1078557bce09SLois Curfman McInnes 
10795615d1e5SSatish Balay #undef __FUNC__
10805615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense"
10818f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1082d3e5ee88SLois Curfman McInnes {
1083c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10843a40ed3dSBarry Smith 
10853a40ed3dSBarry Smith   PetscFunctionBegin;
1086d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
10873a40ed3dSBarry Smith   PetscFunctionReturn(0);
1088d3e5ee88SLois Curfman McInnes }
1089d3e5ee88SLois Curfman McInnes 
10905615d1e5SSatish Balay #undef __FUNC__
10915615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense"
10928f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array)
109364e87e97SBarry Smith {
1094c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10953a40ed3dSBarry Smith 
10963a40ed3dSBarry Smith   PetscFunctionBegin;
109764e87e97SBarry Smith   *array = mat->v;
10983a40ed3dSBarry Smith   PetscFunctionReturn(0);
109964e87e97SBarry Smith }
11000754003eSLois Curfman McInnes 
11015615d1e5SSatish Balay #undef __FUNC__
11025615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense"
11038f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1104ff14e315SSatish Balay {
11053a40ed3dSBarry Smith   PetscFunctionBegin;
11063a40ed3dSBarry Smith   PetscFunctionReturn(0);
1107ff14e315SSatish Balay }
11080754003eSLois Curfman McInnes 
11095615d1e5SSatish Balay #undef __FUNC__
11105615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense"
11116a6a5d1dSBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatGetSubMatrixCall scall,Mat *submat)
11120754003eSLois Curfman McInnes {
1113c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
11140754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
1115160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
11160754003eSLois Curfman McInnes   Scalar       *vwork, *val;
11170754003eSLois Curfman McInnes   Mat          newmat;
11180754003eSLois Curfman McInnes 
11193a40ed3dSBarry Smith   PetscFunctionBegin;
112078b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
112178b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
112278b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
112378b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
11240754003eSLois Curfman McInnes 
11250452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
11260452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
11270452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
1128cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
11290754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
11300754003eSLois Curfman McInnes 
11310754003eSLois Curfman McInnes   /* Create and fill new matrix */
1132b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
11330754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
11340754003eSLois Curfman McInnes     nznew = 0;
11350754003eSLois Curfman McInnes     val   = mat->v + irow[i];
11360754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
11370754003eSLois Curfman McInnes       if (smap[j]) {
11380754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
11390754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
11400754003eSLois Curfman McInnes       }
11410754003eSLois Curfman McInnes     }
11423a40ed3dSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
11430754003eSLois Curfman McInnes   }
11446d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11456d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11460754003eSLois Curfman McInnes 
11470754003eSLois Curfman McInnes   /* Free work space */
11480452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
114978b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
115078b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
11510754003eSLois Curfman McInnes   *submat = newmat;
11523a40ed3dSBarry Smith   PetscFunctionReturn(0);
11530754003eSLois Curfman McInnes }
11540754003eSLois Curfman McInnes 
11555615d1e5SSatish Balay #undef __FUNC__
11565615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense"
11576a6a5d1dSBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1158905e6a2fSBarry Smith {
1159905e6a2fSBarry Smith   int ierr,i;
1160905e6a2fSBarry Smith 
11613a40ed3dSBarry Smith   PetscFunctionBegin;
1162905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1163905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1164905e6a2fSBarry Smith   }
1165905e6a2fSBarry Smith 
1166905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
11676a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1168905e6a2fSBarry Smith   }
11693a40ed3dSBarry Smith   PetscFunctionReturn(0);
1170905e6a2fSBarry Smith }
1171905e6a2fSBarry Smith 
11725615d1e5SSatish Balay #undef __FUNC__
11735615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense"
11748f6be9afSLois Curfman McInnes int MatCopy_SeqDense(Mat A, Mat B)
11754b0e389bSBarry Smith {
11764b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
11773a40ed3dSBarry Smith   int          ierr;
11783a40ed3dSBarry Smith 
11793a40ed3dSBarry Smith   PetscFunctionBegin;
11803a40ed3dSBarry Smith   if (B->type != MATSEQDENSE) {
11813a40ed3dSBarry Smith     ierr = MatCopy_Basic(A,B);CHKERRQ(ierr);
11823a40ed3dSBarry Smith     PetscFunctionReturn(0);
11833a40ed3dSBarry Smith   }
1184e3372554SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
11854b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
11863a40ed3dSBarry Smith   PetscFunctionReturn(0);
11874b0e389bSBarry Smith }
11884b0e389bSBarry Smith 
1189289bc588SBarry Smith /* -------------------------------------------------------------------*/
1190a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1191905e6a2fSBarry Smith        MatGetRow_SeqDense,
1192905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1193905e6a2fSBarry Smith        MatMult_SeqDense,
1194905e6a2fSBarry Smith        MatMultAdd_SeqDense,
1195905e6a2fSBarry Smith        MatMultTrans_SeqDense,
1196905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
1197905e6a2fSBarry Smith        MatSolve_SeqDense,
1198905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
1199905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
1200905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
1201905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1202905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1203ec8511deSBarry Smith        MatRelax_SeqDense,
1204ec8511deSBarry Smith        MatTranspose_SeqDense,
1205905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1206905e6a2fSBarry Smith        MatEqual_SeqDense,
1207905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1208905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1209905e6a2fSBarry Smith        MatNorm_SeqDense,
1210905e6a2fSBarry Smith        0,
1211905e6a2fSBarry Smith        0,
1212905e6a2fSBarry Smith        0,
1213905e6a2fSBarry Smith        MatSetOption_SeqDense,
1214905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1215905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1216905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1217905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1218905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1219905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1220905e6a2fSBarry Smith        MatGetSize_SeqDense,
1221905e6a2fSBarry Smith        MatGetSize_SeqDense,
1222905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1223905e6a2fSBarry Smith        0,
1224905e6a2fSBarry Smith        0,
1225905e6a2fSBarry Smith        MatGetArray_SeqDense,
1226905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
1227a5ae1ecdSBarry Smith        MatConvertSameType_SeqDense,
1228a5ae1ecdSBarry Smith        0,
1229a5ae1ecdSBarry Smith        0,
1230a5ae1ecdSBarry Smith        0,
1231a5ae1ecdSBarry Smith        0,
1232a5ae1ecdSBarry Smith        MatAXPY_SeqDense,
1233a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1234a5ae1ecdSBarry Smith        0,
12354b0e389bSBarry Smith        MatGetValues_SeqDense,
1236a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1237a5ae1ecdSBarry Smith        0,
1238a5ae1ecdSBarry Smith        MatScale_SeqDense,
1239a5ae1ecdSBarry Smith        0,
1240a5ae1ecdSBarry Smith        0,
1241a5ae1ecdSBarry Smith        0,
1242a5ae1ecdSBarry Smith        MatGetBlockSize_SeqDense,
1243a5ae1ecdSBarry Smith        0,
1244a5ae1ecdSBarry Smith        0,
1245a5ae1ecdSBarry Smith        0,
1246a5ae1ecdSBarry Smith        0,
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        MatGetMaps_Petsc};
125690ace30eSBarry Smith 
12575615d1e5SSatish Balay #undef __FUNC__
12585615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense"
12594b828684SBarry Smith /*@C
1260fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1261d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1262d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1263289bc588SBarry Smith 
1264db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1265db81eaa0SLois Curfman McInnes 
126620563c6bSBarry Smith    Input Parameters:
1267db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
12680c775827SLois Curfman McInnes .  m - number of rows
126918f449edSLois Curfman McInnes .  n - number of columns
1270db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1271dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
127220563c6bSBarry Smith 
127320563c6bSBarry Smith    Output Parameter:
127444cd7ae7SLois Curfman McInnes .  A - the matrix
127520563c6bSBarry Smith 
1276b259b22eSLois Curfman McInnes    Notes:
127718f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
127818f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1279b4fd4287SBarry Smith    set data=PETSC_NULL.
128018f449edSLois Curfman McInnes 
1281dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1282d65003e9SLois Curfman McInnes 
1283db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
128420563c6bSBarry Smith @*/
128544cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1286289bc588SBarry Smith {
128744cd7ae7SLois Curfman McInnes   Mat          B;
128844cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
12893b2fbd54SBarry Smith   int          ierr,flg,size;
12903b2fbd54SBarry Smith 
12913a40ed3dSBarry Smith   PetscFunctionBegin;
12923b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1293a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
129455659b69SBarry Smith 
129544cd7ae7SLois Curfman McInnes   *A            = 0;
1296f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView);
129744cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
129844cd7ae7SLois Curfman McInnes   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
129944cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqDense));
1300a5ae1ecdSBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1301c6e7dd08SBarry Smith   B->ops->destroy    = MatDestroy_SeqDense;
1302c6e7dd08SBarry Smith   B->ops->view       = MatView_SeqDense;
130344cd7ae7SLois Curfman McInnes   B->factor          = 0;
130490f02eecSBarry Smith   B->mapping         = 0;
1305f09e8eb9SSatish Balay   PLogObjectMemory(B,sizeof(struct _p_Mat));
130644cd7ae7SLois Curfman McInnes   B->data            = (void *) b;
130718f449edSLois Curfman McInnes 
130844cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
130944cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
1310a5ae1ecdSBarry Smith 
1311*488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
1312*488ecbafSBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1313a5ae1ecdSBarry Smith 
131444cd7ae7SLois Curfman McInnes   b->pivots       = 0;
131544cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1316b4fd4287SBarry Smith   if (data == PETSC_NULL) {
131744cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
131844cd7ae7SLois Curfman McInnes     PetscMemzero(b->v,m*n*sizeof(Scalar));
131944cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
132044cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
1321273fdc76SBarry Smith   } else { /* user-allocated storage */
132244cd7ae7SLois Curfman McInnes     b->v = data;
132344cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
13242dd2b441SLois Curfman McInnes   }
132525cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
132644cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
13274e220ebcSLois Curfman McInnes 
132844cd7ae7SLois Curfman McInnes   *A = B;
13293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1330289bc588SBarry Smith }
1331289bc588SBarry Smith 
13323b2fbd54SBarry Smith 
13333b2fbd54SBarry Smith 
13343b2fbd54SBarry Smith 
13353b2fbd54SBarry Smith 
13363b2fbd54SBarry Smith 
1337