xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 600d6d0b0d3d6155c6326a1d80903795834b2e01)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*600d6d0bSBarry Smith static char vcid[] = "$Id: dense.c,v 1.158 1998/10/01 21:32:47 balay 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__
935609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_SeqDense"
945609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
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;
1035609ef8eSBarry Smith   if (cpvalues == MAT_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;
1185609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);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   }
359*600d6d0bSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
3607a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
3613a40ed3dSBarry Smith   PetscFunctionReturn(0);
362289bc588SBarry Smith }
363289bc588SBarry Smith 
364289bc588SBarry Smith /* -----------------------------------------------------------------*/
3655615d1e5SSatish Balay #undef __FUNC__
3665615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqDense"
36744cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
368289bc588SBarry Smith {
369c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
370289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
3717a97a34bSBarry Smith   int          ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0;
3723a40ed3dSBarry Smith 
3733a40ed3dSBarry Smith   PetscFunctionBegin;
3747a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
3757a97a34bSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3767a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
37748d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
3787a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
3797a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
38055659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
3813a40ed3dSBarry Smith   PetscFunctionReturn(0);
382289bc588SBarry Smith }
3836ee01492SSatish Balay 
3845615d1e5SSatish Balay #undef __FUNC__
3855615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqDense"
38644cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
387289bc588SBarry Smith {
388c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
389289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
3907a97a34bSBarry Smith   int          ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0;
3913a40ed3dSBarry Smith 
3923a40ed3dSBarry Smith   PetscFunctionBegin;
3937a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
3947a97a34bSBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
3957a97a34bSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
39648d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
3977a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
3987a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
39955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
4003a40ed3dSBarry Smith   PetscFunctionReturn(0);
401289bc588SBarry Smith }
4026ee01492SSatish Balay 
4035615d1e5SSatish Balay #undef __FUNC__
4045615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqDense"
40544cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
406289bc588SBarry Smith {
407c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
408*600d6d0bSBarry Smith   Scalar       *v = mat->v, *x, *y;
4097a97a34bSBarry Smith   int          ierr,_One=1; Scalar _DOne=1.0;
4103a40ed3dSBarry Smith 
4113a40ed3dSBarry Smith   PetscFunctionBegin;
4127a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
413*600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4147a97a34bSBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
4157a97a34bSBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
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);
41955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
4203a40ed3dSBarry Smith   PetscFunctionReturn(0);
421289bc588SBarry Smith }
4226ee01492SSatish Balay 
4235615d1e5SSatish Balay #undef __FUNC__
4245615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense"
42544cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
426289bc588SBarry Smith {
427c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
428*600d6d0bSBarry Smith   Scalar       *v = mat->v, *x, *y;
4297a97a34bSBarry Smith   int          ierr,_One=1;
4306abc6512SBarry Smith   Scalar       _DOne=1.0;
4313a40ed3dSBarry Smith 
4323a40ed3dSBarry Smith   PetscFunctionBegin;
4337a97a34bSBarry Smith   if (!mat->m || !mat->n) PetscFunctionReturn(0);
434*600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4357a97a34bSBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
4367a97a34bSBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
43748d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
4387a97a34bSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
4397a97a34bSBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
44055659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
4413a40ed3dSBarry Smith   PetscFunctionReturn(0);
442289bc588SBarry Smith }
443289bc588SBarry Smith 
444289bc588SBarry Smith /* -----------------------------------------------------------------*/
4455615d1e5SSatish Balay #undef __FUNC__
4465615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense"
4478f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
448289bc588SBarry Smith {
449c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
450289bc588SBarry Smith   Scalar       *v;
451289bc588SBarry Smith   int          i;
45267e560aaSBarry Smith 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
454289bc588SBarry Smith   *ncols = mat->n;
455289bc588SBarry Smith   if (cols) {
45645d48df9SBarry Smith     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols);
45773c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
458289bc588SBarry Smith   }
459289bc588SBarry Smith   if (vals) {
46045d48df9SBarry Smith     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals);
461289bc588SBarry Smith     v = mat->v + row;
46273c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
463289bc588SBarry Smith   }
4643a40ed3dSBarry Smith   PetscFunctionReturn(0);
465289bc588SBarry Smith }
4666ee01492SSatish Balay 
4675615d1e5SSatish Balay #undef __FUNC__
4685615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense"
4698f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
470289bc588SBarry Smith {
4710452661fSBarry Smith   if (cols) { PetscFree(*cols); }
4720452661fSBarry Smith   if (vals) { PetscFree(*vals); }
4733a40ed3dSBarry Smith   PetscFunctionReturn(0);
474289bc588SBarry Smith }
475289bc588SBarry Smith /* ----------------------------------------------------------------*/
4765615d1e5SSatish Balay #undef __FUNC__
4775615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense"
4788f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
479e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
480289bc588SBarry Smith {
481c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
482289bc588SBarry Smith   int          i,j;
483d6dfbf8fSBarry Smith 
4843a40ed3dSBarry Smith   PetscFunctionBegin;
485289bc588SBarry Smith   if (!mat->roworiented) {
486dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
487289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
4883a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
489a8c6a408SBarry Smith         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
490a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
49158804f6dSBarry Smith #endif
492289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
4933a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
494a8c6a408SBarry Smith           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
495a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
49658804f6dSBarry Smith #endif
497289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
498289bc588SBarry Smith         }
499289bc588SBarry Smith       }
5003a40ed3dSBarry Smith     } else {
501289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
5023a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
503a8c6a408SBarry Smith         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
504a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
50558804f6dSBarry Smith #endif
506289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
5073a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
508a8c6a408SBarry Smith           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
509a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
51058804f6dSBarry Smith #endif
511289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
512289bc588SBarry Smith         }
513289bc588SBarry Smith       }
514289bc588SBarry Smith     }
5153a40ed3dSBarry Smith   } else {
516dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
517e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
5183a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
519a8c6a408SBarry Smith         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
520a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
52158804f6dSBarry Smith #endif
522e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
5233a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
524a8c6a408SBarry Smith           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
525a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
52658804f6dSBarry Smith #endif
527e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
528e8d4e0b9SBarry Smith         }
529e8d4e0b9SBarry Smith       }
5303a40ed3dSBarry Smith     } else {
531289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
5323a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
533a8c6a408SBarry Smith         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
534a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
53558804f6dSBarry Smith #endif
536289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
5373a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
538a8c6a408SBarry Smith           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
539a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
54058804f6dSBarry Smith #endif
541289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
542289bc588SBarry Smith         }
543289bc588SBarry Smith       }
544289bc588SBarry Smith     }
545e8d4e0b9SBarry Smith   }
5463a40ed3dSBarry Smith   PetscFunctionReturn(0);
547289bc588SBarry Smith }
548e8d4e0b9SBarry Smith 
5495615d1e5SSatish Balay #undef __FUNC__
5505615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense"
5518f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
552ae80bb75SLois Curfman McInnes {
553ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
554ae80bb75SLois Curfman McInnes   int          i, j;
555ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
556ae80bb75SLois Curfman McInnes 
5573a40ed3dSBarry Smith   PetscFunctionBegin;
558ae80bb75SLois Curfman McInnes   /* row-oriented output */
559ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
560ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
561ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
562ae80bb75SLois Curfman McInnes     }
563ae80bb75SLois Curfman McInnes   }
5643a40ed3dSBarry Smith   PetscFunctionReturn(0);
565ae80bb75SLois Curfman McInnes }
566ae80bb75SLois Curfman McInnes 
567289bc588SBarry Smith /* -----------------------------------------------------------------*/
568289bc588SBarry Smith 
56977c4ece6SBarry Smith #include "sys.h"
57088e32edaSLois Curfman McInnes 
5715615d1e5SSatish Balay #undef __FUNC__
5725615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense"
57319bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
57488e32edaSLois Curfman McInnes {
57588e32edaSLois Curfman McInnes   Mat_SeqDense *a;
57688e32edaSLois Curfman McInnes   Mat          B;
57788e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
57888e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
5794a41a4d3SSatish Balay   Scalar       *vals, *svals, *v,*w;
58019bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
58188e32edaSLois Curfman McInnes 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
58388e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
584a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
58590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
5860752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
587a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object");
58888e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
58988e32edaSLois Curfman McInnes 
59090ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
59190ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
59290ace30eSBarry Smith     B    = *A;
59390ace30eSBarry Smith     a    = (Mat_SeqDense *) B->data;
5944a41a4d3SSatish Balay     v    = a->v;
5954a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
5964a41a4d3SSatish Balay        from row major to column major */
5974a41a4d3SSatish Balay     w = (Scalar *) PetscMalloc((M*N+1)*sizeof(Scalar));CHKPTRQ(w);
59890ace30eSBarry Smith     /* read in nonzero values */
5994a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR); CHKERRQ(ierr);
6004a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6014a41a4d3SSatish Balay     for ( j=0; j<N; j++ ) {
6024a41a4d3SSatish Balay       for ( i=0; i<M; i++ ) {
6034a41a4d3SSatish Balay         *v++ =w[i*N+j];
6044a41a4d3SSatish Balay       }
6054a41a4d3SSatish Balay     }
6064a41a4d3SSatish Balay     PetscFree(w);
6076d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6086d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
60990ace30eSBarry Smith   } else {
61088e32edaSLois Curfman McInnes     /* read row lengths */
61145d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
6120752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
61388e32edaSLois Curfman McInnes 
61488e32edaSLois Curfman McInnes     /* create our matrix */
615b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
61688e32edaSLois Curfman McInnes     B = *A;
61788e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
61888e32edaSLois Curfman McInnes     v = a->v;
61988e32edaSLois Curfman McInnes 
62088e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
62145d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
6220752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
62345d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
6240752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
62588e32edaSLois Curfman McInnes 
62688e32edaSLois Curfman McInnes     /* insert into matrix */
62788e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
62888e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
62988e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
63088e32edaSLois Curfman McInnes     }
6310452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
63288e32edaSLois Curfman McInnes 
6336d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6346d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
63590ace30eSBarry Smith   }
6363a40ed3dSBarry Smith   PetscFunctionReturn(0);
63788e32edaSLois Curfman McInnes }
63888e32edaSLois Curfman McInnes 
639932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
64077c4ece6SBarry Smith #include "sys.h"
641932b0c3eSLois Curfman McInnes 
6425615d1e5SSatish Balay #undef __FUNC__
6435615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII"
644932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
645289bc588SBarry Smith {
646932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
647932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
648d636dbe3SBarry Smith   FILE         *fd;
649932b0c3eSLois Curfman McInnes   char         *outputname;
650932b0c3eSLois Curfman McInnes   Scalar       *v;
651932b0c3eSLois Curfman McInnes 
6523a40ed3dSBarry Smith   PetscFunctionBegin;
65390ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
654932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
65590ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
656639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
6573a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
658f72585f2SLois Curfman McInnes   }
65902cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
66080cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
66144cd7ae7SLois Curfman McInnes       v = a->v + i;
66280cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
66380cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6643a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
665e20fef11SSatish Balay         if (PetscReal(*v) != 0.0 && PetscImaginary(*v) != 0.0) fprintf(fd," %d %g + %g i",j,PetscReal(*v),PetscImaginary(*v));
6663f6de6efSSatish Balay         else if (PetscReal(*v)) fprintf(fd," %d %g ",j,PetscReal(*v));
66780cd9d93SLois Curfman McInnes #else
66880cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
66980cd9d93SLois Curfman McInnes #endif
670d64ed03dSBarry Smith         v += a->m;
67180cd9d93SLois Curfman McInnes       }
67280cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
67380cd9d93SLois Curfman McInnes     }
6743a40ed3dSBarry Smith   } else {
6753a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
67647989497SBarry Smith     int allreal = 1;
67747989497SBarry Smith     /* determine if matrix has all real values */
67847989497SBarry Smith     v = a->v;
67947989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
680e20fef11SSatish Balay       if (PetscImaginary(v[i])) { allreal = 0; break ;}
68147989497SBarry Smith     }
68247989497SBarry Smith #endif
683932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
684932b0c3eSLois Curfman McInnes       v = a->v + i;
685932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6863a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
68747989497SBarry Smith         if (allreal) {
6883f6de6efSSatish Balay           fprintf(fd,"%6.4e ",PetscReal(*v)); v += a->m;
68947989497SBarry Smith         } else {
690e20fef11SSatish Balay           fprintf(fd,"%6.4e + %6.4e i ",PetscReal(*v),PetscImaginary(*v)); v += a->m;
69147989497SBarry Smith         }
692289bc588SBarry Smith #else
693932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
694289bc588SBarry Smith #endif
695289bc588SBarry Smith       }
6968e37dc09SBarry Smith       fprintf(fd,"\n");
697289bc588SBarry Smith     }
698da3a660dSBarry Smith   }
699932b0c3eSLois Curfman McInnes   fflush(fd);
7003a40ed3dSBarry Smith   PetscFunctionReturn(0);
701289bc588SBarry Smith }
702289bc588SBarry Smith 
7035615d1e5SSatish Balay #undef __FUNC__
7045615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary"
705932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
706932b0c3eSLois Curfman McInnes {
707932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
708932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
70990ace30eSBarry Smith   int          format;
71090ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
711932b0c3eSLois Curfman McInnes 
7123a40ed3dSBarry Smith   PetscFunctionBegin;
71390ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
71490ace30eSBarry Smith 
71590ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
71602cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
71790ace30eSBarry Smith     /* store the matrix as a dense matrix */
71890ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
71990ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
72090ace30eSBarry Smith     col_lens[1] = m;
72190ace30eSBarry Smith     col_lens[2] = n;
72290ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
7230752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr);
72490ace30eSBarry Smith     PetscFree(col_lens);
72590ace30eSBarry Smith 
72690ace30eSBarry Smith     /* write out matrix, by rows */
72745d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
72890ace30eSBarry Smith     v    = a->v;
72990ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
73090ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
73190ace30eSBarry Smith         vals[i + j*m] = *v++;
73290ace30eSBarry Smith       }
73390ace30eSBarry Smith     }
7340752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr);
73590ace30eSBarry Smith     PetscFree(vals);
73690ace30eSBarry Smith   } else {
7370452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
738932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
739932b0c3eSLois Curfman McInnes     col_lens[1] = m;
740932b0c3eSLois Curfman McInnes     col_lens[2] = n;
741932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
742932b0c3eSLois Curfman McInnes 
743932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
744932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
7450752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr);
746932b0c3eSLois Curfman McInnes 
747932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
748932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
749932b0c3eSLois Curfman McInnes     ict = 0;
750932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
751932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
752932b0c3eSLois Curfman McInnes     }
7530752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr);
7540452661fSBarry Smith     PetscFree(col_lens);
755932b0c3eSLois Curfman McInnes 
756932b0c3eSLois Curfman McInnes     /* store nonzero values */
75745d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
758932b0c3eSLois Curfman McInnes     ict = 0;
759932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
760932b0c3eSLois Curfman McInnes       v = a->v + i;
761932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
762932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
763932b0c3eSLois Curfman McInnes       }
764932b0c3eSLois Curfman McInnes     }
7650752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr);
7660452661fSBarry Smith     PetscFree(anonz);
76790ace30eSBarry Smith   }
7683a40ed3dSBarry Smith   PetscFunctionReturn(0);
769932b0c3eSLois Curfman McInnes }
770932b0c3eSLois Curfman McInnes 
7715615d1e5SSatish Balay #undef __FUNC__
7725615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense"
773c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer)
774932b0c3eSLois Curfman McInnes {
775932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
776bcd2baecSBarry Smith   ViewerType   vtype;
777bcd2baecSBarry Smith   int          ierr;
778932b0c3eSLois Curfman McInnes 
7793a40ed3dSBarry Smith   PetscFunctionBegin;
780bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
781bcd2baecSBarry Smith 
782bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
7831a0f753aSSatish Balay     ierr = ViewerMatlabPutScalar_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr);
7843a40ed3dSBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
7853a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
7863a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
7873a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
7885cd90555SBarry Smith   } else {
7895cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
790932b0c3eSLois Curfman McInnes   }
7913a40ed3dSBarry Smith   PetscFunctionReturn(0);
792932b0c3eSLois Curfman McInnes }
793289bc588SBarry Smith 
7945615d1e5SSatish Balay #undef __FUNC__
7955615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense"
796c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
797289bc588SBarry Smith {
798ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
79990f02eecSBarry Smith   int          ierr;
80090f02eecSBarry Smith 
8013a40ed3dSBarry Smith   PetscFunctionBegin;
80294d884c6SBarry Smith   if (--mat->refct > 0) PetscFunctionReturn(0);
80394d884c6SBarry Smith 
80494d884c6SBarry Smith   if (mat->mapping) {
80594d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
80694d884c6SBarry Smith   }
80794d884c6SBarry Smith   if (mat->bmapping) {
80894d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr);
80994d884c6SBarry Smith   }
8103a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
811c6e7dd08SBarry Smith   PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n);
812a5a9c739SBarry Smith #endif
8130452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
8143439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
8150452661fSBarry Smith   PetscFree(l);
816a5ae1ecdSBarry Smith   if (mat->rmap) {
817a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
818a5ae1ecdSBarry Smith   }
819a5ae1ecdSBarry Smith   if (mat->cmap) {
820a5ae1ecdSBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
82190f02eecSBarry Smith   }
822a5a9c739SBarry Smith   PLogObjectDestroy(mat);
8230452661fSBarry Smith   PetscHeaderDestroy(mat);
8243a40ed3dSBarry Smith   PetscFunctionReturn(0);
825289bc588SBarry Smith }
826289bc588SBarry Smith 
8275615d1e5SSatish Balay #undef __FUNC__
8285615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense"
8298f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
830289bc588SBarry Smith {
831c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
832d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
833d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
83448b35521SBarry Smith 
8353a40ed3dSBarry Smith   PetscFunctionBegin;
836d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
8373638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
838d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
839d64ed03dSBarry Smith       Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w);
840d64ed03dSBarry Smith       for ( j=0; j<m; j++ ) {
841d64ed03dSBarry Smith         for ( k=0; k<n; k++ ) {
842d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
843d64ed03dSBarry Smith         }
844d64ed03dSBarry Smith       }
845d64ed03dSBarry Smith       PetscMemcpy(v,w,m*n*sizeof(Scalar));
846d64ed03dSBarry Smith       PetscFree(w);
847d64ed03dSBarry Smith     } else {
848d3e5ee88SLois Curfman McInnes       for ( j=0; j<m; j++ ) {
849289bc588SBarry Smith         for ( k=0; k<j; k++ ) {
850d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
851d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
852d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
853289bc588SBarry Smith         }
854289bc588SBarry Smith       }
855d64ed03dSBarry Smith     }
8563a40ed3dSBarry Smith   } else { /* out-of-place transpose */
857d3e5ee88SLois Curfman McInnes     int          ierr;
858d3e5ee88SLois Curfman McInnes     Mat          tmat;
859ec8511deSBarry Smith     Mat_SeqDense *tmatd;
860d3e5ee88SLois Curfman McInnes     Scalar       *v2;
861b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
862ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
8630de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
864d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
8650de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
866d3e5ee88SLois Curfman McInnes     }
8676d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8686d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
869d3e5ee88SLois Curfman McInnes     *matout = tmat;
87048b35521SBarry Smith   }
8713a40ed3dSBarry Smith   PetscFunctionReturn(0);
872289bc588SBarry Smith }
873289bc588SBarry Smith 
8745615d1e5SSatish Balay #undef __FUNC__
8755615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense"
8768f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
877289bc588SBarry Smith {
878c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
879c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
880289bc588SBarry Smith   int          i;
881289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
8829ea5d5aeSSatish Balay 
8833a40ed3dSBarry Smith   PetscFunctionBegin;
8843a40ed3dSBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
8853a40ed3dSBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
8863a40ed3dSBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
887289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
8883a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
889289bc588SBarry Smith     v1++; v2++;
890289bc588SBarry Smith   }
89177c4ece6SBarry Smith   *flg = PETSC_TRUE;
8923a40ed3dSBarry Smith   PetscFunctionReturn(0);
893289bc588SBarry Smith }
894289bc588SBarry Smith 
8955615d1e5SSatish Balay #undef __FUNC__
8965615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense"
8978f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
898289bc588SBarry Smith {
899c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
9007a97a34bSBarry Smith   int          ierr,i, n, len;
90144cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
90244cd7ae7SLois Curfman McInnes 
9033a40ed3dSBarry Smith   PetscFunctionBegin;
9047a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
9057a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
906*600d6d0bSBarry Smith   ierr = VecGetArray(v,&x); CHKERRQ(ierr);
90744cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
908e3372554SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
90944cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
910289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
911289bc588SBarry Smith   }
9127a97a34bSBarry Smith   ierr = VecRestoreArray(v,&x); CHKERRQ(ierr);
9133a40ed3dSBarry Smith   PetscFunctionReturn(0);
914289bc588SBarry Smith }
915289bc588SBarry Smith 
9165615d1e5SSatish Balay #undef __FUNC__
9175615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense"
9188f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
919289bc588SBarry Smith {
920c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
921da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
9227a97a34bSBarry Smith   int          ierr,i,j,m = mat->m, n = mat->n;
92355659b69SBarry Smith 
9243a40ed3dSBarry Smith   PetscFunctionBegin;
92528988994SBarry Smith   if (ll) {
9267a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
927*600d6d0bSBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
928e3372554SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
929da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
930da3a660dSBarry Smith       x = l[i];
931da3a660dSBarry Smith       v = mat->v + i;
932da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
933da3a660dSBarry Smith     }
9347a97a34bSBarry Smith     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
93544cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
936da3a660dSBarry Smith   }
93728988994SBarry Smith   if (rr) {
9387a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
939*600d6d0bSBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
940e3372554SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
941da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
942da3a660dSBarry Smith       x = r[i];
943da3a660dSBarry Smith       v = mat->v + i*m;
944da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
945da3a660dSBarry Smith     }
9467a97a34bSBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
94744cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
948da3a660dSBarry Smith   }
9493a40ed3dSBarry Smith   PetscFunctionReturn(0);
950289bc588SBarry Smith }
951289bc588SBarry Smith 
9525615d1e5SSatish Balay #undef __FUNC__
9535615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense"
9548f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm)
955289bc588SBarry Smith {
956c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
957289bc588SBarry Smith   Scalar       *v = mat->v;
958289bc588SBarry Smith   double       sum = 0.0;
959289bc588SBarry Smith   int          i, j;
96055659b69SBarry Smith 
9613a40ed3dSBarry Smith   PetscFunctionBegin;
962289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
963289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
9643a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
9653f6de6efSSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
966289bc588SBarry Smith #else
967289bc588SBarry Smith       sum += (*v)*(*v); v++;
968289bc588SBarry Smith #endif
969289bc588SBarry Smith     }
970289bc588SBarry Smith     *norm = sqrt(sum);
97155659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
9723a40ed3dSBarry Smith   } else if (type == NORM_1) {
973289bc588SBarry Smith     *norm = 0.0;
974289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
975289bc588SBarry Smith       sum = 0.0;
976289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
97733a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
978289bc588SBarry Smith       }
979289bc588SBarry Smith       if (sum > *norm) *norm = sum;
980289bc588SBarry Smith     }
98155659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9823a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
983289bc588SBarry Smith     *norm = 0.0;
984289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
985289bc588SBarry Smith       v = mat->v + j;
986289bc588SBarry Smith       sum = 0.0;
987289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
988cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
989289bc588SBarry Smith       }
990289bc588SBarry Smith       if (sum > *norm) *norm = sum;
991289bc588SBarry Smith     }
99255659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9933a40ed3dSBarry Smith   } else {
994e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
995289bc588SBarry Smith   }
9963a40ed3dSBarry Smith   PetscFunctionReturn(0);
997289bc588SBarry Smith }
998289bc588SBarry Smith 
9995615d1e5SSatish Balay #undef __FUNC__
10005615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense"
10018f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1002289bc588SBarry Smith {
1003c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
100467e560aaSBarry Smith 
10053a40ed3dSBarry Smith   PetscFunctionBegin;
10066d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
10076d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
10086d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
1009219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
10106d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
1011219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
10126d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
10136d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
10146d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
10156d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
1016c2653b3dSLois Curfman McInnes            op == MAT_NEW_NONZERO_LOCATION_ERROR ||
10176d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
101890f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
1019b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
1020b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
1021981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
10223a40ed3dSBarry Smith   else {
10233a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
10243a40ed3dSBarry Smith   }
10253a40ed3dSBarry Smith   PetscFunctionReturn(0);
1026289bc588SBarry Smith }
1027289bc588SBarry Smith 
10285615d1e5SSatish Balay #undef __FUNC__
10295615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense"
10308f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
10316f0a148fSBarry Smith {
1032ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
10333a40ed3dSBarry Smith 
10343a40ed3dSBarry Smith   PetscFunctionBegin;
1035cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
10363a40ed3dSBarry Smith   PetscFunctionReturn(0);
10376f0a148fSBarry Smith }
10386f0a148fSBarry Smith 
10395615d1e5SSatish Balay #undef __FUNC__
10405615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense"
10418f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
10424e220ebcSLois Curfman McInnes {
10433a40ed3dSBarry Smith   PetscFunctionBegin;
10444e220ebcSLois Curfman McInnes   *bs = 1;
10453a40ed3dSBarry Smith   PetscFunctionReturn(0);
10464e220ebcSLois Curfman McInnes }
10474e220ebcSLois Curfman McInnes 
10485615d1e5SSatish Balay #undef __FUNC__
10495615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense"
10508f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
10516f0a148fSBarry Smith {
1052ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
10536abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
10546f0a148fSBarry Smith   Scalar       *slot;
105555659b69SBarry Smith 
10563a40ed3dSBarry Smith   PetscFunctionBegin;
105777c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
105878b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
10596f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
10606f0a148fSBarry Smith     slot = l->v + rows[i];
10616f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
10626f0a148fSBarry Smith   }
10636f0a148fSBarry Smith   if (diag) {
10646f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
10656f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1066260925b8SBarry Smith       *slot = *diag;
10676f0a148fSBarry Smith     }
10686f0a148fSBarry Smith   }
1069260925b8SBarry Smith   ISRestoreIndices(is,&rows);
10703a40ed3dSBarry Smith   PetscFunctionReturn(0);
10716f0a148fSBarry Smith }
1072557bce09SLois Curfman McInnes 
10735615d1e5SSatish Balay #undef __FUNC__
10745615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense"
10758f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n)
1076557bce09SLois Curfman McInnes {
1077c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10783a40ed3dSBarry Smith 
10793a40ed3dSBarry Smith   PetscFunctionBegin;
1080557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
10813a40ed3dSBarry Smith   PetscFunctionReturn(0);
1082557bce09SLois Curfman McInnes }
1083557bce09SLois Curfman McInnes 
10845615d1e5SSatish Balay #undef __FUNC__
10855615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense"
10868f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1087d3e5ee88SLois Curfman McInnes {
1088c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10893a40ed3dSBarry Smith 
10903a40ed3dSBarry Smith   PetscFunctionBegin;
1091d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
10923a40ed3dSBarry Smith   PetscFunctionReturn(0);
1093d3e5ee88SLois Curfman McInnes }
1094d3e5ee88SLois Curfman McInnes 
10955615d1e5SSatish Balay #undef __FUNC__
10965615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense"
10978f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array)
109864e87e97SBarry Smith {
1099c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
11003a40ed3dSBarry Smith 
11013a40ed3dSBarry Smith   PetscFunctionBegin;
110264e87e97SBarry Smith   *array = mat->v;
11033a40ed3dSBarry Smith   PetscFunctionReturn(0);
110464e87e97SBarry Smith }
11050754003eSLois Curfman McInnes 
11065615d1e5SSatish Balay #undef __FUNC__
11075615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense"
11088f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1109ff14e315SSatish Balay {
11103a40ed3dSBarry Smith   PetscFunctionBegin;
11113a40ed3dSBarry Smith   PetscFunctionReturn(0);
1112ff14e315SSatish Balay }
11130754003eSLois Curfman McInnes 
11145615d1e5SSatish Balay #undef __FUNC__
11155615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense"
1116182d2002SSatish Balay static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatGetSubMatrixCall scall,Mat *B)
11170754003eSLois Curfman McInnes {
1118c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1119182d2002SSatish Balay   int          i, j, ierr, m = mat->m, *irow, *icol, nrows, ncols;
1120182d2002SSatish Balay   Scalar       *av, *bv, *v = mat->v;
11210754003eSLois Curfman McInnes   Mat          newmat;
11220754003eSLois Curfman McInnes 
11233a40ed3dSBarry Smith   PetscFunctionBegin;
112478b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
112578b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
112678b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
112778b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
11280754003eSLois Curfman McInnes 
1129182d2002SSatish Balay   /* Check submatrixcall */
1130182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1131182d2002SSatish Balay     int n_cols,n_rows;
1132182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1133182d2002SSatish Balay     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1134182d2002SSatish Balay     newmat = *B;
1135182d2002SSatish Balay   } else {
11360754003eSLois Curfman McInnes     /* Create and fill new matrix */
1137b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
1138182d2002SSatish Balay   }
1139182d2002SSatish Balay 
1140182d2002SSatish Balay   /* Now extract the data pointers and do the copy, column at a time */
1141182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1142182d2002SSatish Balay 
1143182d2002SSatish Balay   for ( i=0; i<ncols; i++ ) {
1144182d2002SSatish Balay     av = v + m*icol[i];
1145182d2002SSatish Balay     for (j=0; j<nrows; j++ ) {
1146182d2002SSatish Balay       *bv++ = av[irow[j]];
11470754003eSLois Curfman McInnes     }
11480754003eSLois Curfman McInnes   }
1149182d2002SSatish Balay 
1150182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
11516d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11526d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11530754003eSLois Curfman McInnes 
11540754003eSLois Curfman McInnes   /* Free work space */
115578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
115678b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1157182d2002SSatish Balay   *B = newmat;
11583a40ed3dSBarry Smith   PetscFunctionReturn(0);
11590754003eSLois Curfman McInnes }
11600754003eSLois Curfman McInnes 
11615615d1e5SSatish Balay #undef __FUNC__
11625615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense"
11636a6a5d1dSBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1164905e6a2fSBarry Smith {
1165905e6a2fSBarry Smith   int ierr,i;
1166905e6a2fSBarry Smith 
11673a40ed3dSBarry Smith   PetscFunctionBegin;
1168905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1169905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1170905e6a2fSBarry Smith   }
1171905e6a2fSBarry Smith 
1172905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
11736a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1174905e6a2fSBarry Smith   }
11753a40ed3dSBarry Smith   PetscFunctionReturn(0);
1176905e6a2fSBarry Smith }
1177905e6a2fSBarry Smith 
11785615d1e5SSatish Balay #undef __FUNC__
11795615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense"
1180cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A, Mat B,MatStructure str)
11814b0e389bSBarry Smith {
11824b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
11833a40ed3dSBarry Smith   int          ierr;
11843a40ed3dSBarry Smith 
11853a40ed3dSBarry Smith   PetscFunctionBegin;
11863a40ed3dSBarry Smith   if (B->type != MATSEQDENSE) {
1187cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
11883a40ed3dSBarry Smith     PetscFunctionReturn(0);
11893a40ed3dSBarry Smith   }
1190e3372554SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
11914b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
11923a40ed3dSBarry Smith   PetscFunctionReturn(0);
11934b0e389bSBarry Smith }
11944b0e389bSBarry Smith 
1195289bc588SBarry Smith /* -------------------------------------------------------------------*/
1196a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1197905e6a2fSBarry Smith        MatGetRow_SeqDense,
1198905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1199905e6a2fSBarry Smith        MatMult_SeqDense,
1200905e6a2fSBarry Smith        MatMultAdd_SeqDense,
1201905e6a2fSBarry Smith        MatMultTrans_SeqDense,
1202905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
1203905e6a2fSBarry Smith        MatSolve_SeqDense,
1204905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
1205905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
1206905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
1207905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1208905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1209ec8511deSBarry Smith        MatRelax_SeqDense,
1210ec8511deSBarry Smith        MatTranspose_SeqDense,
1211905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1212905e6a2fSBarry Smith        MatEqual_SeqDense,
1213905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1214905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1215905e6a2fSBarry Smith        MatNorm_SeqDense,
1216905e6a2fSBarry Smith        0,
1217905e6a2fSBarry Smith        0,
1218905e6a2fSBarry Smith        0,
1219905e6a2fSBarry Smith        MatSetOption_SeqDense,
1220905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1221905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1222905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1223905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1224905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1225905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1226905e6a2fSBarry Smith        MatGetSize_SeqDense,
1227905e6a2fSBarry Smith        MatGetSize_SeqDense,
1228905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1229905e6a2fSBarry Smith        0,
1230905e6a2fSBarry Smith        0,
1231905e6a2fSBarry Smith        MatGetArray_SeqDense,
1232905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
12335609ef8eSBarry Smith        MatDuplicate_SeqDense,
1234a5ae1ecdSBarry Smith        0,
1235a5ae1ecdSBarry Smith        0,
1236a5ae1ecdSBarry Smith        0,
1237a5ae1ecdSBarry Smith        0,
1238a5ae1ecdSBarry Smith        MatAXPY_SeqDense,
1239a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1240a5ae1ecdSBarry Smith        0,
12414b0e389bSBarry Smith        MatGetValues_SeqDense,
1242a5ae1ecdSBarry Smith        MatCopy_SeqDense,
1243a5ae1ecdSBarry Smith        0,
1244a5ae1ecdSBarry Smith        MatScale_SeqDense,
1245a5ae1ecdSBarry Smith        0,
1246a5ae1ecdSBarry Smith        0,
1247a5ae1ecdSBarry Smith        0,
1248a5ae1ecdSBarry Smith        MatGetBlockSize_SeqDense,
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        0,
1260a5ae1ecdSBarry Smith        0,
1261a5ae1ecdSBarry Smith        MatGetMaps_Petsc};
126290ace30eSBarry Smith 
12635615d1e5SSatish Balay #undef __FUNC__
12645615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense"
12654b828684SBarry Smith /*@C
1266fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1267d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1268d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1269289bc588SBarry Smith 
1270db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1271db81eaa0SLois Curfman McInnes 
127220563c6bSBarry Smith    Input Parameters:
1273db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
12740c775827SLois Curfman McInnes .  m - number of rows
127518f449edSLois Curfman McInnes .  n - number of columns
1276db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1277dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
127820563c6bSBarry Smith 
127920563c6bSBarry Smith    Output Parameter:
128044cd7ae7SLois Curfman McInnes .  A - the matrix
128120563c6bSBarry Smith 
1282b259b22eSLois Curfman McInnes    Notes:
128318f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
128418f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1285b4fd4287SBarry Smith    set data=PETSC_NULL.
128618f449edSLois 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;
1302f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,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