xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 6a6a5d1dfdd5fcb6e31127fc1e7bbfb8f4c48c0a)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*6a6a5d1dSBarry Smith static char vcid[] = "$Id: dense.c,v 1.135 1998/01/06 20:10:09 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"
9b16a3bb1SBarry Smith #include "pinclude/plapack.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   }
83289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
84a8c6a408SBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,0,"Bad argument to LU factorization");
85e3372554SBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Bad LU factorization");
86c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
8755659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
883a40ed3dSBarry Smith   PetscFunctionReturn(0);
89289bc588SBarry Smith }
906ee01492SSatish Balay 
915615d1e5SSatish Balay #undef __FUNC__
925615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqDense"
938f6be9afSLois Curfman McInnes int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
9402cad45dSBarry Smith {
9502cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
9602cad45dSBarry Smith   int          ierr;
9702cad45dSBarry Smith   Mat          newi;
9802cad45dSBarry Smith 
993a40ed3dSBarry Smith   PetscFunctionBegin;
10002cad45dSBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
10102cad45dSBarry Smith   l = (Mat_SeqDense *) newi->data;
10202cad45dSBarry Smith   if (cpvalues == COPY_VALUES) {
10302cad45dSBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
10402cad45dSBarry Smith   }
10502cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
10602cad45dSBarry Smith   *newmat = newi;
1073a40ed3dSBarry Smith   PetscFunctionReturn(0);
10802cad45dSBarry Smith }
10902cad45dSBarry Smith 
1105615d1e5SSatish Balay #undef __FUNC__
1115615d1e5SSatish Balay #define __FUNC__ "MatLUFactorSymbolic_SeqDense"
1128f6be9afSLois Curfman McInnes int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
113289bc588SBarry Smith {
1143a40ed3dSBarry Smith   int ierr;
1153a40ed3dSBarry Smith 
1163a40ed3dSBarry Smith   PetscFunctionBegin;
1173a40ed3dSBarry Smith   ierr = MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);CHKERRQ(ierr);
1183a40ed3dSBarry Smith   PetscFunctionReturn(0);
119289bc588SBarry Smith }
1206ee01492SSatish Balay 
1215615d1e5SSatish Balay #undef __FUNC__
1225615d1e5SSatish Balay #define __FUNC__ "MatLUFactorNumeric_SeqDense"
1238f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
124289bc588SBarry Smith {
12502cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data;
1263a40ed3dSBarry Smith   int          ierr;
1273a40ed3dSBarry Smith 
1283a40ed3dSBarry Smith   PetscFunctionBegin;
12902cad45dSBarry Smith   /* copy the numerical values */
13002cad45dSBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
13102cad45dSBarry Smith   (*fact)->factor = 0;
1323a40ed3dSBarry Smith   ierr = MatLUFactor(*fact,0,0,1.0);CHKERRQ(ierr);
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
134289bc588SBarry Smith }
1356ee01492SSatish Balay 
1365615d1e5SSatish Balay #undef __FUNC__
1375615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense"
1388f6be9afSLois Curfman McInnes int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
139289bc588SBarry Smith {
1403a40ed3dSBarry Smith   int ierr;
1413a40ed3dSBarry Smith 
1423a40ed3dSBarry Smith   PetscFunctionBegin;
1433a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1443a40ed3dSBarry Smith   PetscFunctionReturn(0);
145289bc588SBarry Smith }
1466ee01492SSatish Balay 
1475615d1e5SSatish Balay #undef __FUNC__
1485615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense"
1498f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
150289bc588SBarry Smith {
1513a40ed3dSBarry Smith   int ierr;
1523a40ed3dSBarry Smith 
1533a40ed3dSBarry Smith   PetscFunctionBegin;
1543a40ed3dSBarry Smith   ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr);
1553a40ed3dSBarry Smith   PetscFunctionReturn(0);
156289bc588SBarry Smith }
1576ee01492SSatish Balay 
1585615d1e5SSatish Balay #undef __FUNC__
1595615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense"
1608f6be9afSLois Curfman McInnes int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
161289bc588SBarry Smith {
162c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1636abc6512SBarry Smith   int           info;
16455659b69SBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
166752f5567SLois Curfman McInnes   if (mat->pivots) {
1670452661fSBarry Smith     PetscFree(mat->pivots);
168c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
169752f5567SLois Curfman McInnes     mat->pivots = 0;
170752f5567SLois Curfman McInnes   }
171289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
172e3372554SBarry Smith   if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization");
173c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
17455659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
1753a40ed3dSBarry Smith   PetscFunctionReturn(0);
176289bc588SBarry Smith }
177289bc588SBarry Smith 
1785615d1e5SSatish Balay #undef __FUNC__
1795615d1e5SSatish Balay #define __FUNC__ "MatSolve_SeqDense"
1808f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
181289bc588SBarry Smith {
182c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
183a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
1846abc6512SBarry Smith   Scalar       *x, *y;
18567e560aaSBarry Smith 
1863a40ed3dSBarry Smith   PetscFunctionBegin;
187a57ad546SLois Curfman McInnes   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
188416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
189c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
19048d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
191289bc588SBarry Smith   }
192c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
19348d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
194289bc588SBarry Smith   }
195a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Matrix must be factored to solve");
196a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"MBad solve");
19755659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
199289bc588SBarry Smith }
2006ee01492SSatish Balay 
2015615d1e5SSatish Balay #undef __FUNC__
2025615d1e5SSatish Balay #define __FUNC__ "MatSolveTrans_SeqDense"
2038f6be9afSLois Curfman McInnes int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
204da3a660dSBarry Smith {
205c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2066abc6512SBarry Smith   int          one = 1, info;
2076abc6512SBarry Smith   Scalar       *x, *y;
20867e560aaSBarry Smith 
2093a40ed3dSBarry Smith   PetscFunctionBegin;
210da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
211416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
212752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
213da3a660dSBarry Smith   if (mat->pivots) {
21448d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
215da3a660dSBarry Smith   }
216da3a660dSBarry Smith   else {
21748d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
218da3a660dSBarry Smith   }
219a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
22055659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2213a40ed3dSBarry Smith   PetscFunctionReturn(0);
222da3a660dSBarry Smith }
2236ee01492SSatish Balay 
2245615d1e5SSatish Balay #undef __FUNC__
2255615d1e5SSatish Balay #define __FUNC__ "MatSolveAdd_SeqDense"
2268f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
227da3a660dSBarry Smith {
228c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2296abc6512SBarry Smith   int          one = 1, info,ierr;
2306abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
231da3a660dSBarry Smith   Vec          tmp = 0;
23267e560aaSBarry Smith 
2333a40ed3dSBarry Smith   PetscFunctionBegin;
234da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
235da3a660dSBarry Smith   if (yy == zz) {
23678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
237c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
23878b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
239da3a660dSBarry Smith   }
240416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
241752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
242da3a660dSBarry Smith   if (mat->pivots) {
24348d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
244a8c6a408SBarry Smith   } else {
24548d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
246da3a660dSBarry Smith   }
247a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
248a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy); CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
249a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);}
25055659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
252da3a660dSBarry Smith }
25367e560aaSBarry Smith 
2545615d1e5SSatish Balay #undef __FUNC__
2555615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense"
2568f6be9afSLois Curfman McInnes int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
257da3a660dSBarry Smith {
258c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
2596abc6512SBarry Smith   int           one = 1, info,ierr;
2606abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
261da3a660dSBarry Smith   Vec           tmp;
26267e560aaSBarry Smith 
2633a40ed3dSBarry Smith   PetscFunctionBegin;
264da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
265da3a660dSBarry Smith   if (yy == zz) {
26678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
267c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
26878b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
269da3a660dSBarry Smith   }
270416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
271752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
272da3a660dSBarry Smith   if (mat->pivots) {
27348d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
2743a40ed3dSBarry Smith   } else {
27548d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
276da3a660dSBarry Smith   }
277a8c6a408SBarry Smith   if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve");
27890f02eecSBarry Smith   if (tmp) {
27990f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);  CHKERRQ(ierr);
28090f02eecSBarry Smith     ierr = VecDestroy(tmp); CHKERRQ(ierr);
2813a40ed3dSBarry Smith   } else {
28290f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);
28390f02eecSBarry Smith   }
28455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
2853a40ed3dSBarry Smith   PetscFunctionReturn(0);
286da3a660dSBarry Smith }
287289bc588SBarry Smith /* ------------------------------------------------------------------*/
2885615d1e5SSatish Balay #undef __FUNC__
2895615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense"
2908f6be9afSLois Curfman McInnes int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
29120e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
292289bc588SBarry Smith {
293c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
294289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
2956abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
296289bc588SBarry Smith 
2973a40ed3dSBarry Smith   PetscFunctionBegin;
298289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
2993a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
300bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
301289bc588SBarry Smith   }
302289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
303289bc588SBarry Smith   while (its--) {
304289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
305289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
3063a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
307f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
308f1747703SBarry Smith            not happy about returning a double complex */
309f1747703SBarry Smith         int    _i;
310f1747703SBarry Smith         Scalar sum = b[i];
311f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
312f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
313f1747703SBarry Smith         }
314f1747703SBarry Smith         xt = sum;
315f1747703SBarry Smith #else
316289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
317f1747703SBarry Smith #endif
31855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
319289bc588SBarry Smith       }
320289bc588SBarry Smith     }
321289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
322289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
3233a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
324f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
325f1747703SBarry Smith            not happy about returning a double complex */
326f1747703SBarry Smith         int    _i;
327f1747703SBarry Smith         Scalar sum = b[i];
328f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
329f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
330f1747703SBarry Smith         }
331f1747703SBarry Smith         xt = sum;
332f1747703SBarry Smith #else
333289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
334f1747703SBarry Smith #endif
33555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
336289bc588SBarry Smith       }
337289bc588SBarry Smith     }
338289bc588SBarry Smith   }
3393a40ed3dSBarry Smith   PetscFunctionReturn(0);
340289bc588SBarry Smith }
341289bc588SBarry Smith 
342289bc588SBarry Smith /* -----------------------------------------------------------------*/
3435615d1e5SSatish Balay #undef __FUNC__
3445615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqDense"
34544cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
346289bc588SBarry Smith {
347c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
348289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
349289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
3503a40ed3dSBarry Smith 
3513a40ed3dSBarry Smith   PetscFunctionBegin;
352289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
35348d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
35455659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
3553a40ed3dSBarry Smith   PetscFunctionReturn(0);
356289bc588SBarry Smith }
3576ee01492SSatish Balay 
3585615d1e5SSatish Balay #undef __FUNC__
3595615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqDense"
36044cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
361289bc588SBarry Smith {
362c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
363289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
364289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
3653a40ed3dSBarry Smith 
3663a40ed3dSBarry Smith   PetscFunctionBegin;
367289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
36848d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
36955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
3703a40ed3dSBarry Smith   PetscFunctionReturn(0);
371289bc588SBarry Smith }
3726ee01492SSatish Balay 
3735615d1e5SSatish Balay #undef __FUNC__
3745615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqDense"
37544cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
376289bc588SBarry Smith {
377c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
378289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
3796abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
3803a40ed3dSBarry Smith 
3813a40ed3dSBarry Smith   PetscFunctionBegin;
382289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
383416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
38448d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
38555659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
3863a40ed3dSBarry Smith   PetscFunctionReturn(0);
387289bc588SBarry Smith }
3886ee01492SSatish Balay 
3895615d1e5SSatish Balay #undef __FUNC__
3905615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense"
39144cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
392289bc588SBarry Smith {
393c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
394289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
395289bc588SBarry Smith   int          _One=1;
3966abc6512SBarry Smith   Scalar       _DOne=1.0;
3973a40ed3dSBarry Smith 
3983a40ed3dSBarry Smith   PetscFunctionBegin;
39944cd7ae7SLois Curfman McInnes   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
400717eeb74SLois Curfman McInnes   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
40148d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
40255659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
4033a40ed3dSBarry Smith   PetscFunctionReturn(0);
404289bc588SBarry Smith }
405289bc588SBarry Smith 
406289bc588SBarry Smith /* -----------------------------------------------------------------*/
4075615d1e5SSatish Balay #undef __FUNC__
4085615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense"
4098f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
410289bc588SBarry Smith {
411c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
412289bc588SBarry Smith   Scalar       *v;
413289bc588SBarry Smith   int          i;
41467e560aaSBarry Smith 
4153a40ed3dSBarry Smith   PetscFunctionBegin;
416289bc588SBarry Smith   *ncols = mat->n;
417289bc588SBarry Smith   if (cols) {
41845d48df9SBarry Smith     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols);
41973c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
420289bc588SBarry Smith   }
421289bc588SBarry Smith   if (vals) {
42245d48df9SBarry Smith     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals);
423289bc588SBarry Smith     v = mat->v + row;
42473c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
425289bc588SBarry Smith   }
4263a40ed3dSBarry Smith   PetscFunctionReturn(0);
427289bc588SBarry Smith }
4286ee01492SSatish Balay 
4295615d1e5SSatish Balay #undef __FUNC__
4305615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense"
4318f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
432289bc588SBarry Smith {
4330452661fSBarry Smith   if (cols) { PetscFree(*cols); }
4340452661fSBarry Smith   if (vals) { PetscFree(*vals); }
4353a40ed3dSBarry Smith   PetscFunctionReturn(0);
436289bc588SBarry Smith }
437289bc588SBarry Smith /* ----------------------------------------------------------------*/
4385615d1e5SSatish Balay #undef __FUNC__
4395615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense"
4408f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
441e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
442289bc588SBarry Smith {
443c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
444289bc588SBarry Smith   int          i,j;
445d6dfbf8fSBarry Smith 
4463a40ed3dSBarry Smith   PetscFunctionBegin;
447289bc588SBarry Smith   if (!mat->roworiented) {
448dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
449289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
4503a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
451a8c6a408SBarry Smith         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
452a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
45358804f6dSBarry Smith #endif
454289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
4553a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
456a8c6a408SBarry Smith           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
457a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
45858804f6dSBarry Smith #endif
459289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
460289bc588SBarry Smith         }
461289bc588SBarry Smith       }
4623a40ed3dSBarry Smith     } else {
463289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
4643a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
465a8c6a408SBarry Smith         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
466a8c6a408SBarry Smith         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
46758804f6dSBarry Smith #endif
468289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
4693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
470a8c6a408SBarry Smith           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
471a8c6a408SBarry Smith           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
47258804f6dSBarry Smith #endif
473289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
474289bc588SBarry Smith         }
475289bc588SBarry Smith       }
476289bc588SBarry Smith     }
4773a40ed3dSBarry Smith   } else {
478dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
479e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
4803a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
481a8c6a408SBarry Smith         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
482a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
48358804f6dSBarry Smith #endif
484e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
4853a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
486a8c6a408SBarry Smith           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
487a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
48858804f6dSBarry Smith #endif
489e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
490e8d4e0b9SBarry Smith         }
491e8d4e0b9SBarry Smith       }
4923a40ed3dSBarry Smith     } else {
493289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
4943a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
495a8c6a408SBarry Smith         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
496a8c6a408SBarry Smith         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
49758804f6dSBarry Smith #endif
498289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
4993a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
500a8c6a408SBarry Smith           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
501a8c6a408SBarry Smith           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
50258804f6dSBarry Smith #endif
503289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
504289bc588SBarry Smith         }
505289bc588SBarry Smith       }
506289bc588SBarry Smith     }
507e8d4e0b9SBarry Smith   }
5083a40ed3dSBarry Smith   PetscFunctionReturn(0);
509289bc588SBarry Smith }
510e8d4e0b9SBarry Smith 
5115615d1e5SSatish Balay #undef __FUNC__
5125615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense"
5138f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
514ae80bb75SLois Curfman McInnes {
515ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
516ae80bb75SLois Curfman McInnes   int          i, j;
517ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
518ae80bb75SLois Curfman McInnes 
5193a40ed3dSBarry Smith   PetscFunctionBegin;
520ae80bb75SLois Curfman McInnes   /* row-oriented output */
521ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
522ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
523ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
524ae80bb75SLois Curfman McInnes     }
525ae80bb75SLois Curfman McInnes   }
5263a40ed3dSBarry Smith   PetscFunctionReturn(0);
527ae80bb75SLois Curfman McInnes }
528ae80bb75SLois Curfman McInnes 
529289bc588SBarry Smith /* -----------------------------------------------------------------*/
530289bc588SBarry Smith 
53177c4ece6SBarry Smith #include "sys.h"
53288e32edaSLois Curfman McInnes 
5335615d1e5SSatish Balay #undef __FUNC__
5345615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense"
53519bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
53688e32edaSLois Curfman McInnes {
53788e32edaSLois Curfman McInnes   Mat_SeqDense *a;
53888e32edaSLois Curfman McInnes   Mat          B;
53988e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
54088e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
54188e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
54219bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
54388e32edaSLois Curfman McInnes 
5443a40ed3dSBarry Smith   PetscFunctionBegin;
54588e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
546a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
54790ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
5480752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
549a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object");
55088e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
55188e32edaSLois Curfman McInnes 
55290ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
55390ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
55490ace30eSBarry Smith     B    = *A;
55590ace30eSBarry Smith     a    = (Mat_SeqDense *) B->data;
55690ace30eSBarry Smith 
55790ace30eSBarry Smith     /* read in nonzero values */
5580752156aSBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,PETSC_SCALAR); CHKERRQ(ierr);
55990ace30eSBarry Smith 
5606d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5616d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
562d64ed03dSBarry Smith     /* ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); */
56390ace30eSBarry Smith   } else {
56488e32edaSLois Curfman McInnes     /* read row lengths */
56545d48df9SBarry Smith     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
5660752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
56788e32edaSLois Curfman McInnes 
56888e32edaSLois Curfman McInnes     /* create our matrix */
569b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
57088e32edaSLois Curfman McInnes     B = *A;
57188e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
57288e32edaSLois Curfman McInnes     v = a->v;
57388e32edaSLois Curfman McInnes 
57488e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
57545d48df9SBarry Smith     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
5760752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
57745d48df9SBarry Smith     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
5780752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
57988e32edaSLois Curfman McInnes 
58088e32edaSLois Curfman McInnes     /* insert into matrix */
58188e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
58288e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
58388e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
58488e32edaSLois Curfman McInnes     }
5850452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
58688e32edaSLois Curfman McInnes 
5876d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5886d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
58990ace30eSBarry Smith   }
5903a40ed3dSBarry Smith   PetscFunctionReturn(0);
59188e32edaSLois Curfman McInnes }
59288e32edaSLois Curfman McInnes 
593932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
59477c4ece6SBarry Smith #include "sys.h"
595932b0c3eSLois Curfman McInnes 
5965615d1e5SSatish Balay #undef __FUNC__
5975615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII"
598932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
599289bc588SBarry Smith {
600932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
601932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
602d636dbe3SBarry Smith   FILE         *fd;
603932b0c3eSLois Curfman McInnes   char         *outputname;
604932b0c3eSLois Curfman McInnes   Scalar       *v;
605932b0c3eSLois Curfman McInnes 
6063a40ed3dSBarry Smith   PetscFunctionBegin;
60790ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
608932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
60990ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
610639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
6113a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
612f72585f2SLois Curfman McInnes   }
61302cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
61480cd9d93SLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
61544cd7ae7SLois Curfman McInnes       v = a->v + i;
61680cd9d93SLois Curfman McInnes       fprintf(fd,"row %d:",i);
61780cd9d93SLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6183a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
619d64ed03dSBarry Smith         if (real(*v) != 0.0 && imag(*v) != 0.0) fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
62080cd9d93SLois Curfman McInnes         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
62180cd9d93SLois Curfman McInnes #else
62280cd9d93SLois Curfman McInnes         if (*v) fprintf(fd," %d %g ",j,*v);
62380cd9d93SLois Curfman McInnes #endif
624d64ed03dSBarry Smith         v += a->m;
62580cd9d93SLois Curfman McInnes       }
62680cd9d93SLois Curfman McInnes       fprintf(fd,"\n");
62780cd9d93SLois Curfman McInnes     }
6283a40ed3dSBarry Smith   } else {
6293a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
63047989497SBarry Smith     int allreal = 1;
63147989497SBarry Smith     /* determine if matrix has all real values */
63247989497SBarry Smith     v = a->v;
63347989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
63447989497SBarry Smith       if (imag(v[i])) { allreal = 0; break ;}
63547989497SBarry Smith     }
63647989497SBarry Smith #endif
637932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
638932b0c3eSLois Curfman McInnes       v = a->v + i;
639932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
6403a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
64147989497SBarry Smith         if (allreal) {
64247989497SBarry Smith           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
64347989497SBarry Smith         } else {
644932b0c3eSLois Curfman McInnes           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
64547989497SBarry Smith         }
646289bc588SBarry Smith #else
647932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
648289bc588SBarry Smith #endif
649289bc588SBarry Smith       }
6508e37dc09SBarry Smith       fprintf(fd,"\n");
651289bc588SBarry Smith     }
652da3a660dSBarry Smith   }
653932b0c3eSLois Curfman McInnes   fflush(fd);
6543a40ed3dSBarry Smith   PetscFunctionReturn(0);
655289bc588SBarry Smith }
656289bc588SBarry Smith 
6575615d1e5SSatish Balay #undef __FUNC__
6585615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary"
659932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
660932b0c3eSLois Curfman McInnes {
661932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
662932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
66390ace30eSBarry Smith   int          format;
66490ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
665932b0c3eSLois Curfman McInnes 
6663a40ed3dSBarry Smith   PetscFunctionBegin;
66790ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
66890ace30eSBarry Smith 
66990ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
67002cad45dSBarry Smith   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
67190ace30eSBarry Smith     /* store the matrix as a dense matrix */
67290ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
67390ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
67490ace30eSBarry Smith     col_lens[1] = m;
67590ace30eSBarry Smith     col_lens[2] = n;
67690ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
6770752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr);
67890ace30eSBarry Smith     PetscFree(col_lens);
67990ace30eSBarry Smith 
68090ace30eSBarry Smith     /* write out matrix, by rows */
68145d48df9SBarry Smith     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
68290ace30eSBarry Smith     v    = a->v;
68390ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
68490ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
68590ace30eSBarry Smith         vals[i + j*m] = *v++;
68690ace30eSBarry Smith       }
68790ace30eSBarry Smith     }
6880752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr);
68990ace30eSBarry Smith     PetscFree(vals);
69090ace30eSBarry Smith   } else {
6910452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
692932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
693932b0c3eSLois Curfman McInnes     col_lens[1] = m;
694932b0c3eSLois Curfman McInnes     col_lens[2] = n;
695932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
696932b0c3eSLois Curfman McInnes 
697932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
698932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
6990752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr);
700932b0c3eSLois Curfman McInnes 
701932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
702932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
703932b0c3eSLois Curfman McInnes     ict = 0;
704932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
705932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
706932b0c3eSLois Curfman McInnes     }
7070752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr);
7080452661fSBarry Smith     PetscFree(col_lens);
709932b0c3eSLois Curfman McInnes 
710932b0c3eSLois Curfman McInnes     /* store nonzero values */
71145d48df9SBarry Smith     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
712932b0c3eSLois Curfman McInnes     ict = 0;
713932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
714932b0c3eSLois Curfman McInnes       v = a->v + i;
715932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
716932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
717932b0c3eSLois Curfman McInnes       }
718932b0c3eSLois Curfman McInnes     }
7190752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr);
7200452661fSBarry Smith     PetscFree(anonz);
72190ace30eSBarry Smith   }
7223a40ed3dSBarry Smith   PetscFunctionReturn(0);
723932b0c3eSLois Curfman McInnes }
724932b0c3eSLois Curfman McInnes 
7255615d1e5SSatish Balay #undef __FUNC__
7265615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense"
7278f6be9afSLois Curfman McInnes int MatView_SeqDense(PetscObject obj,Viewer viewer)
728932b0c3eSLois Curfman McInnes {
729932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
730932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
731bcd2baecSBarry Smith   ViewerType   vtype;
732bcd2baecSBarry Smith   int          ierr;
733932b0c3eSLois Curfman McInnes 
7343a40ed3dSBarry Smith   PetscFunctionBegin;
735bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
736bcd2baecSBarry Smith 
737bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
7383a40ed3dSBarry Smith     ierr = ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr);
7393a40ed3dSBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
7403a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
7413a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
7423a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
743932b0c3eSLois Curfman McInnes   }
7443a40ed3dSBarry Smith   PetscFunctionReturn(0);
745932b0c3eSLois Curfman McInnes }
746289bc588SBarry Smith 
7475615d1e5SSatish Balay #undef __FUNC__
7485615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense"
7498f6be9afSLois Curfman McInnes int MatDestroy_SeqDense(PetscObject obj)
750289bc588SBarry Smith {
751289bc588SBarry Smith   Mat          mat = (Mat) obj;
752ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
75390f02eecSBarry Smith   int          ierr;
75490f02eecSBarry Smith 
7553a40ed3dSBarry Smith   PetscFunctionBegin;
7563a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
757a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
758a5a9c739SBarry Smith #endif
7590452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
7603439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
7610452661fSBarry Smith   PetscFree(l);
76290f02eecSBarry Smith   if (mat->mapping) {
76390f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
76490f02eecSBarry Smith   }
765a5a9c739SBarry Smith   PLogObjectDestroy(mat);
7660452661fSBarry Smith   PetscHeaderDestroy(mat);
7673a40ed3dSBarry Smith   PetscFunctionReturn(0);
768289bc588SBarry Smith }
769289bc588SBarry Smith 
7705615d1e5SSatish Balay #undef __FUNC__
7715615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense"
7728f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
773289bc588SBarry Smith {
774c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
775d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
776d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
77748b35521SBarry Smith 
7783a40ed3dSBarry Smith   PetscFunctionBegin;
779d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
7803638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
781d64ed03dSBarry Smith     if (m != n) { /* malloc temp to hold transpose */
782d64ed03dSBarry Smith       Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w);
783d64ed03dSBarry Smith       for ( j=0; j<m; j++ ) {
784d64ed03dSBarry Smith         for ( k=0; k<n; k++ ) {
785d64ed03dSBarry Smith           w[k + j*n] = v[j + k*m];
786d64ed03dSBarry Smith         }
787d64ed03dSBarry Smith       }
788d64ed03dSBarry Smith       PetscMemcpy(v,w,m*n*sizeof(Scalar));
789d64ed03dSBarry Smith       PetscFree(w);
790d64ed03dSBarry Smith     } else {
791d3e5ee88SLois Curfman McInnes       for ( j=0; j<m; j++ ) {
792289bc588SBarry Smith         for ( k=0; k<j; k++ ) {
793d3e5ee88SLois Curfman McInnes           tmp = v[j + k*n];
794d3e5ee88SLois Curfman McInnes           v[j + k*n] = v[k + j*n];
795d3e5ee88SLois Curfman McInnes           v[k + j*n] = tmp;
796289bc588SBarry Smith         }
797289bc588SBarry Smith       }
798d64ed03dSBarry Smith     }
7993a40ed3dSBarry Smith   } else { /* out-of-place transpose */
800d3e5ee88SLois Curfman McInnes     int          ierr;
801d3e5ee88SLois Curfman McInnes     Mat          tmat;
802ec8511deSBarry Smith     Mat_SeqDense *tmatd;
803d3e5ee88SLois Curfman McInnes     Scalar       *v2;
804b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
805ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
8060de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
807d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
8080de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
809d3e5ee88SLois Curfman McInnes     }
8106d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8116d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
812d3e5ee88SLois Curfman McInnes     *matout = tmat;
81348b35521SBarry Smith   }
8143a40ed3dSBarry Smith   PetscFunctionReturn(0);
815289bc588SBarry Smith }
816289bc588SBarry Smith 
8175615d1e5SSatish Balay #undef __FUNC__
8185615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense"
8198f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
820289bc588SBarry Smith {
821c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
822c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
823289bc588SBarry Smith   int          i;
824289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
8259ea5d5aeSSatish Balay 
8263a40ed3dSBarry Smith   PetscFunctionBegin;
8273a40ed3dSBarry Smith   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
8283a40ed3dSBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
8293a40ed3dSBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
830289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
8313a40ed3dSBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
832289bc588SBarry Smith     v1++; v2++;
833289bc588SBarry Smith   }
83477c4ece6SBarry Smith   *flg = PETSC_TRUE;
8353a40ed3dSBarry Smith   PetscFunctionReturn(0);
836289bc588SBarry Smith }
837289bc588SBarry Smith 
8385615d1e5SSatish Balay #undef __FUNC__
8395615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense"
8408f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
841289bc588SBarry Smith {
842c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
84344cd7ae7SLois Curfman McInnes   int          i, n, len;
84444cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
84544cd7ae7SLois Curfman McInnes 
8463a40ed3dSBarry Smith   PetscFunctionBegin;
84744cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
848289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
84944cd7ae7SLois Curfman McInnes   len = PetscMin(mat->m,mat->n);
850e3372554SBarry Smith   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
85144cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
852289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
853289bc588SBarry Smith   }
8543a40ed3dSBarry Smith   PetscFunctionReturn(0);
855289bc588SBarry Smith }
856289bc588SBarry Smith 
8575615d1e5SSatish Balay #undef __FUNC__
8585615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense"
8598f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
860289bc588SBarry Smith {
861c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
862da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
863da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
86455659b69SBarry Smith 
8653a40ed3dSBarry Smith   PetscFunctionBegin;
86628988994SBarry Smith   if (ll) {
867da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
868e3372554SBarry Smith     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
869da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
870da3a660dSBarry Smith       x = l[i];
871da3a660dSBarry Smith       v = mat->v + i;
872da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
873da3a660dSBarry Smith     }
87444cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
875da3a660dSBarry Smith   }
87628988994SBarry Smith   if (rr) {
877da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
878e3372554SBarry Smith     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
879da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
880da3a660dSBarry Smith       x = r[i];
881da3a660dSBarry Smith       v = mat->v + i*m;
882da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
883da3a660dSBarry Smith     }
88444cd7ae7SLois Curfman McInnes     PLogFlops(n*m);
885da3a660dSBarry Smith   }
8863a40ed3dSBarry Smith   PetscFunctionReturn(0);
887289bc588SBarry Smith }
888289bc588SBarry Smith 
8895615d1e5SSatish Balay #undef __FUNC__
8905615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense"
8918f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm)
892289bc588SBarry Smith {
893c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
894289bc588SBarry Smith   Scalar       *v = mat->v;
895289bc588SBarry Smith   double       sum = 0.0;
896289bc588SBarry Smith   int          i, j;
89755659b69SBarry Smith 
8983a40ed3dSBarry Smith   PetscFunctionBegin;
899289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
900289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
9013a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
902289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
903289bc588SBarry Smith #else
904289bc588SBarry Smith       sum += (*v)*(*v); v++;
905289bc588SBarry Smith #endif
906289bc588SBarry Smith     }
907289bc588SBarry Smith     *norm = sqrt(sum);
90855659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
9093a40ed3dSBarry Smith   } else if (type == NORM_1) {
910289bc588SBarry Smith     *norm = 0.0;
911289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
912289bc588SBarry Smith       sum = 0.0;
913289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
91433a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
915289bc588SBarry Smith       }
916289bc588SBarry Smith       if (sum > *norm) *norm = sum;
917289bc588SBarry Smith     }
91855659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9193a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
920289bc588SBarry Smith     *norm = 0.0;
921289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
922289bc588SBarry Smith       v = mat->v + j;
923289bc588SBarry Smith       sum = 0.0;
924289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
925cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
926289bc588SBarry Smith       }
927289bc588SBarry Smith       if (sum > *norm) *norm = sum;
928289bc588SBarry Smith     }
92955659b69SBarry Smith     PLogFlops(mat->n*mat->m);
9303a40ed3dSBarry Smith   } else {
931e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
932289bc588SBarry Smith   }
9333a40ed3dSBarry Smith   PetscFunctionReturn(0);
934289bc588SBarry Smith }
935289bc588SBarry Smith 
9365615d1e5SSatish Balay #undef __FUNC__
9375615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense"
9388f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
939289bc588SBarry Smith {
940c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
94167e560aaSBarry Smith 
9423a40ed3dSBarry Smith   PetscFunctionBegin;
9436d4a8577SBarry Smith   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
9446d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
9456d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
946219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
9476d4a8577SBarry Smith            op == MAT_COLUMNS_SORTED ||
948219d9a1aSLois Curfman McInnes            op == MAT_COLUMNS_UNSORTED ||
9496d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
9506d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
9516d4a8577SBarry Smith            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
9526d4a8577SBarry Smith            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
953c2653b3dSLois Curfman McInnes            op == MAT_NEW_NONZERO_LOCATION_ERROR ||
9546d4a8577SBarry Smith            op == MAT_NO_NEW_DIAGONALS ||
95590f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
9562b362799SSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES)
957981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
9583a40ed3dSBarry Smith   else {
9593a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
9603a40ed3dSBarry Smith   }
9613a40ed3dSBarry Smith   PetscFunctionReturn(0);
962289bc588SBarry Smith }
963289bc588SBarry Smith 
9645615d1e5SSatish Balay #undef __FUNC__
9655615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense"
9668f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
9676f0a148fSBarry Smith {
968ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
9693a40ed3dSBarry Smith 
9703a40ed3dSBarry Smith   PetscFunctionBegin;
971cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
9723a40ed3dSBarry Smith   PetscFunctionReturn(0);
9736f0a148fSBarry Smith }
9746f0a148fSBarry Smith 
9755615d1e5SSatish Balay #undef __FUNC__
9765615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense"
9778f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
9784e220ebcSLois Curfman McInnes {
9793a40ed3dSBarry Smith   PetscFunctionBegin;
9804e220ebcSLois Curfman McInnes   *bs = 1;
9813a40ed3dSBarry Smith   PetscFunctionReturn(0);
9824e220ebcSLois Curfman McInnes }
9834e220ebcSLois Curfman McInnes 
9845615d1e5SSatish Balay #undef __FUNC__
9855615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense"
9868f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
9876f0a148fSBarry Smith {
988ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
9896abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
9906f0a148fSBarry Smith   Scalar       *slot;
99155659b69SBarry Smith 
9923a40ed3dSBarry Smith   PetscFunctionBegin;
99377c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
99478b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
9956f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
9966f0a148fSBarry Smith     slot = l->v + rows[i];
9976f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
9986f0a148fSBarry Smith   }
9996f0a148fSBarry Smith   if (diag) {
10006f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
10016f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1002260925b8SBarry Smith       *slot = *diag;
10036f0a148fSBarry Smith     }
10046f0a148fSBarry Smith   }
1005260925b8SBarry Smith   ISRestoreIndices(is,&rows);
10063a40ed3dSBarry Smith   PetscFunctionReturn(0);
10076f0a148fSBarry Smith }
1008557bce09SLois Curfman McInnes 
10095615d1e5SSatish Balay #undef __FUNC__
10105615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense"
10118f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n)
1012557bce09SLois Curfman McInnes {
1013c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10143a40ed3dSBarry Smith 
10153a40ed3dSBarry Smith   PetscFunctionBegin;
1016557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
10173a40ed3dSBarry Smith   PetscFunctionReturn(0);
1018557bce09SLois Curfman McInnes }
1019557bce09SLois Curfman McInnes 
10205615d1e5SSatish Balay #undef __FUNC__
10215615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense"
10228f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1023d3e5ee88SLois Curfman McInnes {
1024c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10253a40ed3dSBarry Smith 
10263a40ed3dSBarry Smith   PetscFunctionBegin;
1027d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
10283a40ed3dSBarry Smith   PetscFunctionReturn(0);
1029d3e5ee88SLois Curfman McInnes }
1030d3e5ee88SLois Curfman McInnes 
10315615d1e5SSatish Balay #undef __FUNC__
10325615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense"
10338f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array)
103464e87e97SBarry Smith {
1035c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10363a40ed3dSBarry Smith 
10373a40ed3dSBarry Smith   PetscFunctionBegin;
103864e87e97SBarry Smith   *array = mat->v;
10393a40ed3dSBarry Smith   PetscFunctionReturn(0);
104064e87e97SBarry Smith }
10410754003eSLois Curfman McInnes 
10425615d1e5SSatish Balay #undef __FUNC__
10435615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense"
10448f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1045ff14e315SSatish Balay {
10463a40ed3dSBarry Smith   PetscFunctionBegin;
10473a40ed3dSBarry Smith   PetscFunctionReturn(0);
1048ff14e315SSatish Balay }
10490754003eSLois Curfman McInnes 
10505615d1e5SSatish Balay #undef __FUNC__
10515615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense"
1052*6a6a5d1dSBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatGetSubMatrixCall scall,Mat *submat)
10530754003eSLois Curfman McInnes {
1054c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
10550754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
1056160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
10570754003eSLois Curfman McInnes   Scalar       *vwork, *val;
10580754003eSLois Curfman McInnes   Mat          newmat;
10590754003eSLois Curfman McInnes 
10603a40ed3dSBarry Smith   PetscFunctionBegin;
106178b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
106278b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
106378b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
106478b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
10650754003eSLois Curfman McInnes 
10660452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
10670452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
10680452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
1069cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
10700754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
10710754003eSLois Curfman McInnes 
10720754003eSLois Curfman McInnes   /* Create and fill new matrix */
1073b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
10740754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
10750754003eSLois Curfman McInnes     nznew = 0;
10760754003eSLois Curfman McInnes     val   = mat->v + irow[i];
10770754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
10780754003eSLois Curfman McInnes       if (smap[j]) {
10790754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
10800754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
10810754003eSLois Curfman McInnes       }
10820754003eSLois Curfman McInnes     }
10833a40ed3dSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
10840754003eSLois Curfman McInnes   }
10856d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10866d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10870754003eSLois Curfman McInnes 
10880754003eSLois Curfman McInnes   /* Free work space */
10890452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
109078b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
109178b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10920754003eSLois Curfman McInnes   *submat = newmat;
10933a40ed3dSBarry Smith   PetscFunctionReturn(0);
10940754003eSLois Curfman McInnes }
10950754003eSLois Curfman McInnes 
10965615d1e5SSatish Balay #undef __FUNC__
10975615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense"
1098*6a6a5d1dSBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1099905e6a2fSBarry Smith {
1100905e6a2fSBarry Smith   int ierr,i;
1101905e6a2fSBarry Smith 
11023a40ed3dSBarry Smith   PetscFunctionBegin;
1103905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1104905e6a2fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1105905e6a2fSBarry Smith   }
1106905e6a2fSBarry Smith 
1107905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) {
1108*6a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1109905e6a2fSBarry Smith   }
11103a40ed3dSBarry Smith   PetscFunctionReturn(0);
1111905e6a2fSBarry Smith }
1112905e6a2fSBarry Smith 
11135615d1e5SSatish Balay #undef __FUNC__
11145615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense"
11158f6be9afSLois Curfman McInnes int MatCopy_SeqDense(Mat A, Mat B)
11164b0e389bSBarry Smith {
11174b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
11183a40ed3dSBarry Smith   int          ierr;
11193a40ed3dSBarry Smith 
11203a40ed3dSBarry Smith   PetscFunctionBegin;
11213a40ed3dSBarry Smith   if (B->type != MATSEQDENSE) {
11223a40ed3dSBarry Smith     ierr = MatCopy_Basic(A,B);CHKERRQ(ierr);
11233a40ed3dSBarry Smith     PetscFunctionReturn(0);
11243a40ed3dSBarry Smith   }
1125e3372554SBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
11264b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
11273a40ed3dSBarry Smith   PetscFunctionReturn(0);
11284b0e389bSBarry Smith }
11294b0e389bSBarry Smith 
1130289bc588SBarry Smith /* -------------------------------------------------------------------*/
1131ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
1132905e6a2fSBarry Smith        MatGetRow_SeqDense,
1133905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1134905e6a2fSBarry Smith        MatMult_SeqDense,
1135905e6a2fSBarry Smith        MatMultAdd_SeqDense,
1136905e6a2fSBarry Smith        MatMultTrans_SeqDense,
1137905e6a2fSBarry Smith        MatMultTransAdd_SeqDense,
1138905e6a2fSBarry Smith        MatSolve_SeqDense,
1139905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
1140905e6a2fSBarry Smith        MatSolveTrans_SeqDense,
1141905e6a2fSBarry Smith        MatSolveTransAdd_SeqDense,
1142905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1143905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1144ec8511deSBarry Smith        MatRelax_SeqDense,
1145ec8511deSBarry Smith        MatTranspose_SeqDense,
1146905e6a2fSBarry Smith        MatGetInfo_SeqDense,
1147905e6a2fSBarry Smith        MatEqual_SeqDense,
1148905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1149905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1150905e6a2fSBarry Smith        MatNorm_SeqDense,
1151905e6a2fSBarry Smith        0,
1152905e6a2fSBarry Smith        0,
1153905e6a2fSBarry Smith        0,
1154905e6a2fSBarry Smith        MatSetOption_SeqDense,
1155905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
1156905e6a2fSBarry Smith        MatZeroRows_SeqDense,
1157905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1158905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1159905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1160905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
1161905e6a2fSBarry Smith        MatGetSize_SeqDense,
1162905e6a2fSBarry Smith        MatGetSize_SeqDense,
1163905e6a2fSBarry Smith        MatGetOwnershipRange_SeqDense,
1164905e6a2fSBarry Smith        0,
1165905e6a2fSBarry Smith        0,
1166905e6a2fSBarry Smith        MatGetArray_SeqDense,
1167905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
11683d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
1169905e6a2fSBarry Smith        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
11704b0e389bSBarry Smith        MatGetValues_SeqDense,
11714e220ebcSLois Curfman McInnes        MatCopy_SeqDense,0,MatScale_SeqDense,
11724e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_SeqDense};
117390ace30eSBarry Smith 
11745615d1e5SSatish Balay #undef __FUNC__
11755615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense"
11764b828684SBarry Smith /*@C
1177fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1178d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1179d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1180289bc588SBarry Smith 
118120563c6bSBarry Smith    Input Parameters:
1182029af93fSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
11830c775827SLois Curfman McInnes .  m - number of rows
118418f449edSLois Curfman McInnes .  n - number of columns
1185b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1186dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
118720563c6bSBarry Smith 
118820563c6bSBarry Smith    Output Parameter:
118944cd7ae7SLois Curfman McInnes .  A - the matrix
119020563c6bSBarry Smith 
119118f449edSLois Curfman McInnes   Notes:
119218f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
119318f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
1194b4fd4287SBarry Smith   set data=PETSC_NULL.
119518f449edSLois Curfman McInnes 
1196dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1197d65003e9SLois Curfman McInnes 
1198dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
119920563c6bSBarry Smith @*/
120044cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1201289bc588SBarry Smith {
120244cd7ae7SLois Curfman McInnes   Mat          B;
120344cd7ae7SLois Curfman McInnes   Mat_SeqDense *b;
12043b2fbd54SBarry Smith   int          ierr,flg,size;
12053b2fbd54SBarry Smith 
12063a40ed3dSBarry Smith   PetscFunctionBegin;
12073b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1208a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
120955659b69SBarry Smith 
121044cd7ae7SLois Curfman McInnes   *A            = 0;
1211d4bb536fSBarry Smith   PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView);
121244cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
121344cd7ae7SLois Curfman McInnes   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
121444cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqDense));
121544cd7ae7SLois Curfman McInnes   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
121644cd7ae7SLois Curfman McInnes   B->destroy    = MatDestroy_SeqDense;
121744cd7ae7SLois Curfman McInnes   B->view       = MatView_SeqDense;
121844cd7ae7SLois Curfman McInnes   B->factor     = 0;
121990f02eecSBarry Smith   B->mapping    = 0;
1220f09e8eb9SSatish Balay   PLogObjectMemory(B,sizeof(struct _p_Mat));
122144cd7ae7SLois Curfman McInnes   B->data       = (void *) b;
122218f449edSLois Curfman McInnes 
122344cd7ae7SLois Curfman McInnes   b->m = m;  B->m = m; B->M = m;
122444cd7ae7SLois Curfman McInnes   b->n = n;  B->n = n; B->N = n;
122544cd7ae7SLois Curfman McInnes   b->pivots       = 0;
122644cd7ae7SLois Curfman McInnes   b->roworiented  = 1;
1227b4fd4287SBarry Smith   if (data == PETSC_NULL) {
122844cd7ae7SLois Curfman McInnes     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
122944cd7ae7SLois Curfman McInnes     PetscMemzero(b->v,m*n*sizeof(Scalar));
123044cd7ae7SLois Curfman McInnes     b->user_alloc = 0;
123144cd7ae7SLois Curfman McInnes     PLogObjectMemory(B,n*m*sizeof(Scalar));
123218f449edSLois Curfman McInnes   }
12332dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
123444cd7ae7SLois Curfman McInnes     b->v = data;
123544cd7ae7SLois Curfman McInnes     b->user_alloc = 1;
12362dd2b441SLois Curfman McInnes   }
123725cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
123844cd7ae7SLois Curfman McInnes   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
12394e220ebcSLois Curfman McInnes 
124044cd7ae7SLois Curfman McInnes   *A = B;
12413a40ed3dSBarry Smith   PetscFunctionReturn(0);
1242289bc588SBarry Smith }
1243289bc588SBarry Smith 
12443b2fbd54SBarry Smith 
12453b2fbd54SBarry Smith 
12463b2fbd54SBarry Smith 
12473b2fbd54SBarry Smith 
12483b2fbd54SBarry Smith 
1249