xref: /petsc/src/mat/impls/dense/seq/dense.c (revision bcd2baec2099d24a3a9d1d6033b0e7e45ccc83b9)
1cb512458SBarry Smith #ifndef lint
2*bcd2baecSBarry Smith static char vcid[] = "$Id: dense.c,v 1.92 1996/03/06 21:47:29 balay Exp bsmith $";
3cb512458SBarry Smith #endif
467e560aaSBarry Smith /*
567e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
667e560aaSBarry Smith */
7289bc588SBarry Smith 
817699dbbSLois Curfman McInnes #include "dense.h"
9b16a3bb1SBarry Smith #include "pinclude/plapack.h"
10b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
11289bc588SBarry Smith 
121987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
151987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
161987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
171987afe7SBarry Smith   PLogFlops(2*N-1);
181987afe7SBarry Smith   return 0;
191987afe7SBarry Smith }
201987afe7SBarry Smith 
21c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
22289bc588SBarry Smith {
23c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
24289bc588SBarry Smith   int          i,N = mat->m*mat->n,count = 0;
25289bc588SBarry Smith   Scalar       *v = mat->v;
26289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
27*bcd2baecSBarry Smith   if (nz)      *nz      = count;
28*bcd2baecSBarry Smith   if (nzalloc) *nzalloc = N;
29*bcd2baecSBarry Smith   if (mem)     *mem     = (int)A->mem;
30fa9ec3f1SLois Curfman McInnes   return 0;
31289bc588SBarry Smith }
32289bc588SBarry Smith 
33289bc588SBarry Smith /* ---------------------------------------------------------------*/
34289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
35289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
36c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
37289bc588SBarry Smith {
38c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
396abc6512SBarry Smith   int          info;
4055659b69SBarry Smith 
41289bc588SBarry Smith   if (!mat->pivots) {
420452661fSBarry Smith     mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
43c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
44289bc588SBarry Smith   }
45289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
46ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
47c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
4855659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
49289bc588SBarry Smith   return 0;
50289bc588SBarry Smith }
51c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
52289bc588SBarry Smith {
53a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
54289bc588SBarry Smith }
55c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
56289bc588SBarry Smith {
5749d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
58289bc588SBarry Smith }
59c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
60289bc588SBarry Smith {
61a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
62289bc588SBarry Smith }
63c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
64289bc588SBarry Smith {
6549d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
66289bc588SBarry Smith }
67c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
68289bc588SBarry Smith {
69c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
706abc6512SBarry Smith   int           info;
7155659b69SBarry Smith 
72752f5567SLois Curfman McInnes   if (mat->pivots) {
730452661fSBarry Smith     PetscFree(mat->pivots);
74c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
75752f5567SLois Curfman McInnes     mat->pivots = 0;
76752f5567SLois Curfman McInnes   }
77289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
78ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
79c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
8055659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
81289bc588SBarry Smith   return 0;
82289bc588SBarry Smith }
83289bc588SBarry Smith 
84c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
85289bc588SBarry Smith {
86c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
87a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
886abc6512SBarry Smith   Scalar       *x, *y;
8967e560aaSBarry Smith 
90a57ad546SLois Curfman McInnes   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
91416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
92c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
9348d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
94289bc588SBarry Smith   }
95c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
9648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
97289bc588SBarry Smith   }
98ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
99ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
10055659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
101289bc588SBarry Smith   return 0;
102289bc588SBarry Smith }
103c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
104da3a660dSBarry Smith {
105c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1066abc6512SBarry Smith   int          one = 1, info;
1076abc6512SBarry Smith   Scalar       *x, *y;
10867e560aaSBarry Smith 
109da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
110416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
111752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
112da3a660dSBarry Smith   if (mat->pivots) {
11348d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
114da3a660dSBarry Smith   }
115da3a660dSBarry Smith   else {
11648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
117da3a660dSBarry Smith   }
118ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
11955659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
120da3a660dSBarry Smith   return 0;
121da3a660dSBarry Smith }
122c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
123da3a660dSBarry Smith {
124c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1256abc6512SBarry Smith   int          one = 1, info,ierr;
1266abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
127da3a660dSBarry Smith   Vec          tmp = 0;
12867e560aaSBarry Smith 
129da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
130da3a660dSBarry Smith   if (yy == zz) {
13178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
132c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
13378b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
134da3a660dSBarry Smith   }
135416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
136752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
137da3a660dSBarry Smith   if (mat->pivots) {
13848d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
139da3a660dSBarry Smith   }
140da3a660dSBarry Smith   else {
14148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
142da3a660dSBarry Smith   }
143ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
144da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
145da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
14655659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
147da3a660dSBarry Smith   return 0;
148da3a660dSBarry Smith }
14967e560aaSBarry Smith 
150c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
151da3a660dSBarry Smith {
152c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1536abc6512SBarry Smith   int           one = 1, info,ierr;
1546abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
155da3a660dSBarry Smith   Vec           tmp;
15667e560aaSBarry Smith 
157da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
158da3a660dSBarry Smith   if (yy == zz) {
15978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
160c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
16178b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
162da3a660dSBarry Smith   }
163416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
164752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
165da3a660dSBarry Smith   if (mat->pivots) {
16648d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
167da3a660dSBarry Smith   }
168da3a660dSBarry Smith   else {
16948d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
170da3a660dSBarry Smith   }
171ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
172da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
173da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
17455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
175da3a660dSBarry Smith   return 0;
176da3a660dSBarry Smith }
177289bc588SBarry Smith /* ------------------------------------------------------------------*/
178c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
17920e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
180289bc588SBarry Smith {
181c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
182289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
1836abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
184289bc588SBarry Smith 
185289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
186289bc588SBarry Smith     /* this is a hack fix, should have another version without
187289bc588SBarry Smith        the second BLdot */
188bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
189289bc588SBarry Smith   }
190289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
191289bc588SBarry Smith   while (its--) {
192289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
193289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
194f1747703SBarry Smith #if defined(PETSC_COMPLEX)
195f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
196f1747703SBarry Smith            not happy about returning a double complex */
197f1747703SBarry Smith         int    _i;
198f1747703SBarry Smith         Scalar sum = b[i];
199f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
200f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
201f1747703SBarry Smith         }
202f1747703SBarry Smith         xt = sum;
203f1747703SBarry Smith #else
204289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
205f1747703SBarry Smith #endif
206d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
207289bc588SBarry Smith       }
208289bc588SBarry Smith     }
209289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
210289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
211f1747703SBarry Smith #if defined(PETSC_COMPLEX)
212f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
213f1747703SBarry Smith            not happy about returning a double complex */
214f1747703SBarry Smith         int    _i;
215f1747703SBarry Smith         Scalar sum = b[i];
216f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
217f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
218f1747703SBarry Smith         }
219f1747703SBarry Smith         xt = sum;
220f1747703SBarry Smith #else
221289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
222f1747703SBarry Smith #endif
223d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
224289bc588SBarry Smith       }
225289bc588SBarry Smith     }
226289bc588SBarry Smith   }
227289bc588SBarry Smith   return 0;
228289bc588SBarry Smith }
229289bc588SBarry Smith 
230289bc588SBarry Smith /* -----------------------------------------------------------------*/
231c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
232289bc588SBarry Smith {
233c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
234289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
235289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
236289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
23748d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
23855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
239289bc588SBarry Smith   return 0;
240289bc588SBarry Smith }
241c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
242289bc588SBarry Smith {
243c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
244289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
245289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
246289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
24748d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
24855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
249289bc588SBarry Smith   return 0;
250289bc588SBarry Smith }
251c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
252289bc588SBarry Smith {
253c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
254289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2556abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
256289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
257416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
25848d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
25955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
260289bc588SBarry Smith   return 0;
261289bc588SBarry Smith }
262c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
263289bc588SBarry Smith {
264c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
265289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
266289bc588SBarry Smith   int          _One=1;
2676abc6512SBarry Smith   Scalar       _DOne=1.0;
268289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
269289bc588SBarry Smith   VecGetArray(zz,&z);
270416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
27148d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
27255659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
273289bc588SBarry Smith   return 0;
274289bc588SBarry Smith }
275289bc588SBarry Smith 
276289bc588SBarry Smith /* -----------------------------------------------------------------*/
277c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
278289bc588SBarry Smith {
279c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
280289bc588SBarry Smith   Scalar       *v;
281289bc588SBarry Smith   int          i;
28267e560aaSBarry Smith 
283289bc588SBarry Smith   *ncols = mat->n;
284289bc588SBarry Smith   if (cols) {
2850452661fSBarry Smith     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
28673c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
287289bc588SBarry Smith   }
288289bc588SBarry Smith   if (vals) {
2890452661fSBarry Smith     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
290289bc588SBarry Smith     v = mat->v + row;
29173c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
292289bc588SBarry Smith   }
293289bc588SBarry Smith   return 0;
294289bc588SBarry Smith }
295c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
296289bc588SBarry Smith {
2970452661fSBarry Smith   if (cols) { PetscFree(*cols); }
2980452661fSBarry Smith   if (vals) { PetscFree(*vals); }
299289bc588SBarry Smith   return 0;
300289bc588SBarry Smith }
301289bc588SBarry Smith /* ----------------------------------------------------------------*/
302ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
303e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
304289bc588SBarry Smith {
305c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
306289bc588SBarry Smith   int          i,j;
307d6dfbf8fSBarry Smith 
308289bc588SBarry Smith   if (!mat->roworiented) {
309dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
310289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
311289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
312289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
313289bc588SBarry Smith         }
314289bc588SBarry Smith       }
315289bc588SBarry Smith     }
316289bc588SBarry Smith     else {
317289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
318289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
319289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
320289bc588SBarry Smith         }
321289bc588SBarry Smith       }
322289bc588SBarry Smith     }
323e8d4e0b9SBarry Smith   }
324e8d4e0b9SBarry Smith   else {
325dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
326e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
327e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
328e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
329e8d4e0b9SBarry Smith         }
330e8d4e0b9SBarry Smith       }
331e8d4e0b9SBarry Smith     }
332289bc588SBarry Smith     else {
333289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
334289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
335289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
336289bc588SBarry Smith         }
337289bc588SBarry Smith       }
338289bc588SBarry Smith     }
339e8d4e0b9SBarry Smith   }
340289bc588SBarry Smith   return 0;
341289bc588SBarry Smith }
342e8d4e0b9SBarry Smith 
343ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
344ae80bb75SLois Curfman McInnes {
345ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
346ae80bb75SLois Curfman McInnes   int          i, j;
347ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
348ae80bb75SLois Curfman McInnes 
349ae80bb75SLois Curfman McInnes   /* row-oriented output */
350ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
351ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
352ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
353ae80bb75SLois Curfman McInnes     }
354ae80bb75SLois Curfman McInnes   }
355ae80bb75SLois Curfman McInnes   return 0;
356ae80bb75SLois Curfman McInnes }
357ae80bb75SLois Curfman McInnes 
358289bc588SBarry Smith /* -----------------------------------------------------------------*/
3593d1612f7SBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
360289bc588SBarry Smith {
361c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
362289bc588SBarry Smith   int          ierr;
363289bc588SBarry Smith   Mat          newi;
36448d91487SBarry Smith 
365b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
366ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
36755659b69SBarry Smith   if (cpvalues == COPY_VALUES) {
368416022c9SBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
36955659b69SBarry Smith   }
370227d817aSBarry Smith   newi->assembled = PETSC_TRUE;
371289bc588SBarry Smith   *newmat = newi;
372289bc588SBarry Smith   return 0;
373289bc588SBarry Smith }
374289bc588SBarry Smith 
37588e32edaSLois Curfman McInnes #include "sysio.h"
37688e32edaSLois Curfman McInnes 
37788e32edaSLois Curfman McInnes int MatLoad_SeqDense(Viewer bview,MatType type,Mat *A)
37888e32edaSLois Curfman McInnes {
37988e32edaSLois Curfman McInnes   Mat_SeqDense *a;
38088e32edaSLois Curfman McInnes   Mat          B;
38188e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
38288e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
38388e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
38488e32edaSLois Curfman McInnes   PetscObject  vobj = (PetscObject) bview;
38588e32edaSLois Curfman McInnes   MPI_Comm     comm = vobj->comm;
38688e32edaSLois Curfman McInnes 
38788e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
38888e32edaSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
38988e32edaSLois Curfman McInnes   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
39088e32edaSLois Curfman McInnes   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
39188e32edaSLois Curfman McInnes   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
39288e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
39388e32edaSLois Curfman McInnes 
39488e32edaSLois Curfman McInnes   /* read row lengths */
3950452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
39688e32edaSLois Curfman McInnes   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
39788e32edaSLois Curfman McInnes 
39888e32edaSLois Curfman McInnes   /* create our matrix */
399b4fd4287SBarry Smith   ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
40088e32edaSLois Curfman McInnes   B = *A;
40188e32edaSLois Curfman McInnes   a = (Mat_SeqDense *) B->data;
40288e32edaSLois Curfman McInnes   v = a->v;
40388e32edaSLois Curfman McInnes 
40488e32edaSLois Curfman McInnes   /* read column indices and nonzeros */
4050452661fSBarry Smith   cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
40688e32edaSLois Curfman McInnes   ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr);
4070452661fSBarry Smith   vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
40888e32edaSLois Curfman McInnes   ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr);
40988e32edaSLois Curfman McInnes 
41088e32edaSLois Curfman McInnes   /* insert into matrix */
41188e32edaSLois Curfman McInnes   for ( i=0; i<M; i++ ) {
41288e32edaSLois Curfman McInnes     for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
41388e32edaSLois Curfman McInnes     svals += rowlengths[i]; scols += rowlengths[i];
41488e32edaSLois Curfman McInnes   }
4150452661fSBarry Smith   PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
41688e32edaSLois Curfman McInnes 
41788e32edaSLois Curfman McInnes   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
41888e32edaSLois Curfman McInnes   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
41988e32edaSLois Curfman McInnes   return 0;
42088e32edaSLois Curfman McInnes }
42188e32edaSLois Curfman McInnes 
422932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
423932b0c3eSLois Curfman McInnes #include "sysio.h"
424932b0c3eSLois Curfman McInnes 
425932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
426289bc588SBarry Smith {
427932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
428932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
429d636dbe3SBarry Smith   FILE         *fd;
430932b0c3eSLois Curfman McInnes   char         *outputname;
431932b0c3eSLois Curfman McInnes   Scalar       *v;
432932b0c3eSLois Curfman McInnes 
433227d817aSBarry Smith   ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr);
434932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
435932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetFormat_Private(viewer,&format);
436f72585f2SLois Curfman McInnes   if (format == FILE_FORMAT_INFO) {
437932b0c3eSLois Curfman McInnes     ;  /* do nothing for now */
438f72585f2SLois Curfman McInnes   }
439f72585f2SLois Curfman McInnes   else {
440932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
441932b0c3eSLois Curfman McInnes       v = a->v + i;
442932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
443289bc588SBarry Smith #if defined(PETSC_COMPLEX)
444932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
445289bc588SBarry Smith #else
446932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
447289bc588SBarry Smith #endif
448289bc588SBarry Smith       }
4498e37dc09SBarry Smith       fprintf(fd,"\n");
450289bc588SBarry Smith     }
451da3a660dSBarry Smith   }
452932b0c3eSLois Curfman McInnes   fflush(fd);
453289bc588SBarry Smith   return 0;
454289bc588SBarry Smith }
455289bc588SBarry Smith 
456932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
457932b0c3eSLois Curfman McInnes {
458932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
459932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
460932b0c3eSLois Curfman McInnes   Scalar       *v, *anonz;
461932b0c3eSLois Curfman McInnes 
462932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
4630452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
464932b0c3eSLois Curfman McInnes   col_lens[0] = MAT_COOKIE;
465932b0c3eSLois Curfman McInnes   col_lens[1] = m;
466932b0c3eSLois Curfman McInnes   col_lens[2] = n;
467932b0c3eSLois Curfman McInnes   col_lens[3] = nz;
468932b0c3eSLois Curfman McInnes 
469932b0c3eSLois Curfman McInnes   /* store lengths of each row and write (including header) to file */
470932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) col_lens[4+i] = n;
471932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr);
472932b0c3eSLois Curfman McInnes 
473932b0c3eSLois Curfman McInnes   /* Possibly should write in smaller increments, not whole matrix at once? */
474932b0c3eSLois Curfman McInnes  /* store column indices (zero start index) */
475932b0c3eSLois Curfman McInnes   ict = 0;
476932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
477932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) col_lens[ict++] = j;
478932b0c3eSLois Curfman McInnes   }
479932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr);
4800452661fSBarry Smith   PetscFree(col_lens);
481932b0c3eSLois Curfman McInnes 
482932b0c3eSLois Curfman McInnes   /* store nonzero values */
4830452661fSBarry Smith   anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
484932b0c3eSLois Curfman McInnes   ict = 0;
485932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
486932b0c3eSLois Curfman McInnes     v = a->v + i;
487932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) {
488932b0c3eSLois Curfman McInnes       anonz[ict++] = *v; v += a->m;
489932b0c3eSLois Curfman McInnes     }
490932b0c3eSLois Curfman McInnes   }
491932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr);
4920452661fSBarry Smith   PetscFree(anonz);
493932b0c3eSLois Curfman McInnes   return 0;
494932b0c3eSLois Curfman McInnes }
495932b0c3eSLois Curfman McInnes 
496932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
497932b0c3eSLois Curfman McInnes {
498932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
499932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
500932b0c3eSLois Curfman McInnes   PetscObject  vobj = (PetscObject) viewer;
501*bcd2baecSBarry Smith   ViewerType   vtype;
502*bcd2baecSBarry Smith   int          ierr;
503932b0c3eSLois Curfman McInnes 
504932b0c3eSLois Curfman McInnes   if (!viewer) {
505932b0c3eSLois Curfman McInnes     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
506932b0c3eSLois Curfman McInnes   }
507*bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
508*bcd2baecSBarry Smith 
509*bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
510932b0c3eSLois Curfman McInnes     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
511932b0c3eSLois Curfman McInnes   }
512*bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) {
513932b0c3eSLois Curfman McInnes     return MatView_SeqDense_ASCII(A,viewer);
514932b0c3eSLois Curfman McInnes   }
515*bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
516932b0c3eSLois Curfman McInnes     return MatView_SeqDense_Binary(A,viewer);
517932b0c3eSLois Curfman McInnes   }
518932b0c3eSLois Curfman McInnes   return 0;
519932b0c3eSLois Curfman McInnes }
520289bc588SBarry Smith 
521ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
522289bc588SBarry Smith {
523289bc588SBarry Smith   Mat          mat = (Mat) obj;
524ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
525a5a9c739SBarry Smith #if defined(PETSC_LOG)
526a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
527a5a9c739SBarry Smith #endif
5280452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
5293439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
5300452661fSBarry Smith   PetscFree(l);
531a5a9c739SBarry Smith   PLogObjectDestroy(mat);
5320452661fSBarry Smith   PetscHeaderDestroy(mat);
533289bc588SBarry Smith   return 0;
534289bc588SBarry Smith }
535289bc588SBarry Smith 
536c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
537289bc588SBarry Smith {
538c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
539d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
540d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
54148b35521SBarry Smith 
542d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
5433638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
544d9f96c7cSLois Curfman McInnes     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
545d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
546289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
547d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
548d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
549d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
550289bc588SBarry Smith       }
551289bc588SBarry Smith     }
55248b35521SBarry Smith   }
553d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
554d3e5ee88SLois Curfman McInnes     int          ierr;
555d3e5ee88SLois Curfman McInnes     Mat          tmat;
556ec8511deSBarry Smith     Mat_SeqDense *tmatd;
557d3e5ee88SLois Curfman McInnes     Scalar       *v2;
558b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
559ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
5600de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
561d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
5620de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
563d3e5ee88SLois Curfman McInnes     }
564d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
565d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
566d3e5ee88SLois Curfman McInnes     *matout = tmat;
56748b35521SBarry Smith   }
568289bc588SBarry Smith   return 0;
569289bc588SBarry Smith }
570289bc588SBarry Smith 
5719ea5d5aeSSatish Balay static int MatEqual_SeqDense(Mat A1,Mat A2, int *flg)
572289bc588SBarry Smith {
573c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
574c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
575289bc588SBarry Smith   int          i;
576289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
5779ea5d5aeSSatish Balay 
5789ea5d5aeSSatish Balay   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
5799ea5d5aeSSatish Balay   if (mat1->m != mat2->m) { *flg = 0; return 0;}
5809ea5d5aeSSatish Balay   if (mat1->n != mat2->n) {*flg =0; return 0;}
581289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
5829ea5d5aeSSatish Balay     if (*v1 != *v2) {*flg =0; return 0;}
583289bc588SBarry Smith     v1++; v2++;
584289bc588SBarry Smith   }
5859ea5d5aeSSatish Balay   *flg = 1;
5869ea5d5aeSSatish Balay   return 0;
587289bc588SBarry Smith }
588289bc588SBarry Smith 
589c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
590289bc588SBarry Smith {
591c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
5926abc6512SBarry Smith   int          i, n;
5936abc6512SBarry Smith   Scalar       *x;
594289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
595ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
596289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
597289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
598289bc588SBarry Smith   }
599289bc588SBarry Smith   return 0;
600289bc588SBarry Smith }
601289bc588SBarry Smith 
602052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
603289bc588SBarry Smith {
604c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
605da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
606da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
60755659b69SBarry Smith 
60828988994SBarry Smith   if (ll) {
609da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
610052efed2SBarry Smith     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
61155659b69SBarry Smith     PLogFlops(n*m);
612da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
613da3a660dSBarry Smith       x = l[i];
614da3a660dSBarry Smith       v = mat->v + i;
615da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
616da3a660dSBarry Smith     }
617da3a660dSBarry Smith   }
61828988994SBarry Smith   if (rr) {
619da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
620052efed2SBarry Smith     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
62155659b69SBarry Smith     PLogFlops(n*m);
622da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
623da3a660dSBarry Smith       x = r[i];
624da3a660dSBarry Smith       v = mat->v + i*m;
625da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
626da3a660dSBarry Smith     }
627da3a660dSBarry Smith   }
628289bc588SBarry Smith   return 0;
629289bc588SBarry Smith }
630289bc588SBarry Smith 
631cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
632289bc588SBarry Smith {
633c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
634289bc588SBarry Smith   Scalar       *v = mat->v;
635289bc588SBarry Smith   double       sum = 0.0;
636289bc588SBarry Smith   int          i, j;
63755659b69SBarry Smith 
638289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
639289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
640289bc588SBarry Smith #if defined(PETSC_COMPLEX)
641289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
642289bc588SBarry Smith #else
643289bc588SBarry Smith       sum += (*v)*(*v); v++;
644289bc588SBarry Smith #endif
645289bc588SBarry Smith     }
646289bc588SBarry Smith     *norm = sqrt(sum);
64755659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
648289bc588SBarry Smith   }
649289bc588SBarry Smith   else if (type == NORM_1) {
650289bc588SBarry Smith     *norm = 0.0;
651289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
652289bc588SBarry Smith       sum = 0.0;
653289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
65433a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
655289bc588SBarry Smith       }
656289bc588SBarry Smith       if (sum > *norm) *norm = sum;
657289bc588SBarry Smith     }
65855659b69SBarry Smith     PLogFlops(mat->n*mat->m);
659289bc588SBarry Smith   }
660289bc588SBarry Smith   else if (type == NORM_INFINITY) {
661289bc588SBarry Smith     *norm = 0.0;
662289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
663289bc588SBarry Smith       v = mat->v + j;
664289bc588SBarry Smith       sum = 0.0;
665289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
666cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
667289bc588SBarry Smith       }
668289bc588SBarry Smith       if (sum > *norm) *norm = sum;
669289bc588SBarry Smith     }
67055659b69SBarry Smith     PLogFlops(mat->n*mat->m);
671289bc588SBarry Smith   }
672289bc588SBarry Smith   else {
67348d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
674289bc588SBarry Smith   }
675289bc588SBarry Smith   return 0;
676289bc588SBarry Smith }
677289bc588SBarry Smith 
678c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
679289bc588SBarry Smith {
680c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
68167e560aaSBarry Smith 
682289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
683289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
684c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
68588e32edaSLois Curfman McInnes            op == COLUMNS_SORTED ||
686c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
687c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
688c0bbcb79SLois Curfman McInnes            op == NO_NEW_NONZERO_LOCATIONS ||
689c0bbcb79SLois Curfman McInnes            op == YES_NEW_NONZERO_LOCATIONS ||
690c0bbcb79SLois Curfman McInnes            op == NO_NEW_DIAGONALS ||
691c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
692c0bbcb79SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n");
693c0bbcb79SLois Curfman McInnes   else
6944d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
695289bc588SBarry Smith   return 0;
696289bc588SBarry Smith }
697289bc588SBarry Smith 
698ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
6996f0a148fSBarry Smith {
700ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
701cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
7026f0a148fSBarry Smith   return 0;
7036f0a148fSBarry Smith }
7046f0a148fSBarry Smith 
705ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
7066f0a148fSBarry Smith {
707ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
7086abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
7096f0a148fSBarry Smith   Scalar       *slot;
71055659b69SBarry Smith 
71178b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
71278b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
7136f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
7146f0a148fSBarry Smith     slot = l->v + rows[i];
7156f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
7166f0a148fSBarry Smith   }
7176f0a148fSBarry Smith   if (diag) {
7186f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
7196f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
720260925b8SBarry Smith       *slot = *diag;
7216f0a148fSBarry Smith     }
7226f0a148fSBarry Smith   }
723260925b8SBarry Smith   ISRestoreIndices(is,&rows);
7246f0a148fSBarry Smith   return 0;
7256f0a148fSBarry Smith }
726557bce09SLois Curfman McInnes 
727c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
728557bce09SLois Curfman McInnes {
729c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
730557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
731557bce09SLois Curfman McInnes   return 0;
732557bce09SLois Curfman McInnes }
733557bce09SLois Curfman McInnes 
734c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
735d3e5ee88SLois Curfman McInnes {
736c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
737d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
738d3e5ee88SLois Curfman McInnes   return 0;
739d3e5ee88SLois Curfman McInnes }
740d3e5ee88SLois Curfman McInnes 
741c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
74264e87e97SBarry Smith {
743c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
74464e87e97SBarry Smith   *array = mat->v;
74564e87e97SBarry Smith   return 0;
74664e87e97SBarry Smith }
7470754003eSLois Curfman McInnes 
748ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
749ff14e315SSatish Balay {
750ff14e315SSatish Balay   return 0;
751ff14e315SSatish Balay }
7520754003eSLois Curfman McInnes 
753c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
7540754003eSLois Curfman McInnes {
755ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
7560754003eSLois Curfman McInnes }
7570754003eSLois Curfman McInnes 
758cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
759cddf8d76SBarry Smith                                     Mat *submat)
7600754003eSLois Curfman McInnes {
761c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
7620754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
763160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
7640754003eSLois Curfman McInnes   Scalar       *vwork, *val;
7650754003eSLois Curfman McInnes   Mat          newmat;
7660754003eSLois Curfman McInnes 
76778b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
76878b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
76978b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
77078b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
7710754003eSLois Curfman McInnes 
7720452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
7730452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
7740452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
775cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
7760754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
7770754003eSLois Curfman McInnes 
7780754003eSLois Curfman McInnes   /* Create and fill new matrix */
779b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
7800754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
7810754003eSLois Curfman McInnes     nznew = 0;
7820754003eSLois Curfman McInnes     val   = mat->v + irow[i];
7830754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
7840754003eSLois Curfman McInnes       if (smap[j]) {
7850754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
7860754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
7870754003eSLois Curfman McInnes       }
7880754003eSLois Curfman McInnes     }
789dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
79078b31e54SBarry Smith            CHKERRQ(ierr);
7910754003eSLois Curfman McInnes   }
79278b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
79378b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
7940754003eSLois Curfman McInnes 
7950754003eSLois Curfman McInnes   /* Free work space */
7960452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
79778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
79878b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
7990754003eSLois Curfman McInnes   *submat = newmat;
8000754003eSLois Curfman McInnes   return 0;
8010754003eSLois Curfman McInnes }
8020754003eSLois Curfman McInnes 
8034b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B)
8044b0e389bSBarry Smith {
8054b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
8064b0e389bSBarry Smith   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
8074b0e389bSBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
8084b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
8094b0e389bSBarry Smith   return 0;
8104b0e389bSBarry Smith }
8114b0e389bSBarry Smith 
812289bc588SBarry Smith /* -------------------------------------------------------------------*/
813ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
814ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
815ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
816ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
817ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
818ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
819ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
820ec8511deSBarry Smith        MatRelax_SeqDense,
821ec8511deSBarry Smith        MatTranspose_SeqDense,
822ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
823052efed2SBarry Smith        MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense,
824289bc588SBarry Smith        0,0,
825ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
826ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
827ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
828ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
829ff14e315SSatish Balay        0,0, MatGetArray_SeqDense, MatRestoreArray_SeqDense,0,
830ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
8313d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
832ae80bb75SLois Curfman McInnes        MatAXPY_SeqDense,0,0,
8334b0e389bSBarry Smith        MatGetValues_SeqDense,
8344b0e389bSBarry Smith        MatCopy_SeqDense};
8350754003eSLois Curfman McInnes 
8364b828684SBarry Smith /*@C
837fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
838d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
839d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
840289bc588SBarry Smith 
84120563c6bSBarry Smith    Input Parameters:
8420c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
8430c775827SLois Curfman McInnes .  m - number of rows
84418f449edSLois Curfman McInnes .  n - number of columns
845b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
846dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
84720563c6bSBarry Smith 
84820563c6bSBarry Smith    Output Parameter:
8490c775827SLois Curfman McInnes .  newmat - the matrix
85020563c6bSBarry Smith 
85118f449edSLois Curfman McInnes   Notes:
85218f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
85318f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
854b4fd4287SBarry Smith   set data=PETSC_NULL.
85518f449edSLois Curfman McInnes 
856dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
857d65003e9SLois Curfman McInnes 
858dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
85920563c6bSBarry Smith @*/
86018f449edSLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat)
861289bc588SBarry Smith {
862289bc588SBarry Smith   Mat          mat;
863ec8511deSBarry Smith   Mat_SeqDense *l;
86425cdf11fSBarry Smith   int          ierr,flg;
86555659b69SBarry Smith 
866289bc588SBarry Smith   *newmat        = 0;
86718f449edSLois Curfman McInnes 
8680452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
869a5a9c739SBarry Smith   PLogObjectCreate(mat);
8703439631bSBarry Smith   l               = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l);
871416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
872ec8511deSBarry Smith   mat->destroy    = MatDestroy_SeqDense;
873ec8511deSBarry Smith   mat->view       = MatView_SeqDense;
874289bc588SBarry Smith   mat->factor     = 0;
8753439631bSBarry Smith   PLogObjectMemory(mat,sizeof(struct _Mat));
87618f449edSLois Curfman McInnes   mat->data       = (void *) l;
877d6dfbf8fSBarry Smith 
878289bc588SBarry Smith   l->m            = m;
879289bc588SBarry Smith   l->n            = n;
880289bc588SBarry Smith   l->pivots       = 0;
881289bc588SBarry Smith   l->roworiented  = 1;
882b4fd4287SBarry Smith   if (data == PETSC_NULL) {
883ba7af9b1SBarry Smith     l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v);
884cddf8d76SBarry Smith     PetscMemzero(l->v,m*n*sizeof(Scalar));
8852dd2b441SLois Curfman McInnes     l->user_alloc = 0;
8863439631bSBarry Smith     PLogObjectMemory(mat,n*m*sizeof(Scalar));
88718f449edSLois Curfman McInnes   }
8882dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
8892dd2b441SLois Curfman McInnes     l->v = data;
8902dd2b441SLois Curfman McInnes     l->user_alloc = 1;
8912dd2b441SLois Curfman McInnes   }
89218f449edSLois Curfman McInnes 
89325cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
89425cdf11fSBarry Smith   if (flg) {
895682d7d0cSBarry Smith     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
896682d7d0cSBarry Smith   }
897289bc588SBarry Smith   *newmat = mat;
898289bc588SBarry Smith   return 0;
899289bc588SBarry Smith }
900289bc588SBarry Smith 
901c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
902289bc588SBarry Smith {
903c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
904b4fd4287SBarry Smith   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
905289bc588SBarry Smith }
906