xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 33a8263d2fdc82d6944056e76c023d5b91e5d5ea)
1cb512458SBarry Smith #ifndef lint
2*33a8263dSBarry Smith static char vcid[] = "$Id: dense.c,v 1.72 1995/11/01 23:17:41 bsmith Exp bsmith $";
3cb512458SBarry Smith #endif
4289bc588SBarry Smith 
517699dbbSLois Curfman McInnes #include "dense.h"
6b16a3bb1SBarry Smith #include "pinclude/plapack.h"
7b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
8289bc588SBarry Smith 
91987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
101987afe7SBarry Smith {
111987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
121987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
131987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
141987afe7SBarry Smith   PLogFlops(2*N-1);
151987afe7SBarry Smith   return 0;
161987afe7SBarry Smith }
171987afe7SBarry Smith 
18c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
19289bc588SBarry Smith {
20c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
21289bc588SBarry Smith   int          i,N = mat->m*mat->n,count = 0;
22289bc588SBarry Smith   Scalar       *v = mat->v;
23289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
24c0bbcb79SLois Curfman McInnes   *nz = count; *nzalloc = N; *mem = (int)A->mem;
25fa9ec3f1SLois Curfman McInnes   return 0;
26289bc588SBarry Smith }
27289bc588SBarry Smith 
28289bc588SBarry Smith /* ---------------------------------------------------------------*/
29289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
30289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
31c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
32289bc588SBarry Smith {
33c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
346abc6512SBarry Smith   int          info;
3555659b69SBarry Smith 
36289bc588SBarry Smith   if (!mat->pivots) {
370452661fSBarry Smith     mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
38c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
39289bc588SBarry Smith   }
40289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
41ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
42c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
4355659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
44289bc588SBarry Smith   return 0;
45289bc588SBarry Smith }
46c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
47289bc588SBarry Smith {
48289bc588SBarry Smith   int ierr;
49c0bbcb79SLois Curfman McInnes   ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr);
50289bc588SBarry Smith   return 0;
51289bc588SBarry Smith }
52c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
53289bc588SBarry Smith {
5449d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
55289bc588SBarry Smith }
56c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
57289bc588SBarry Smith {
58289bc588SBarry Smith   int ierr;
59c0bbcb79SLois Curfman McInnes   ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr);
60289bc588SBarry Smith   return 0;
61289bc588SBarry Smith }
62c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
63289bc588SBarry Smith {
6449d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
65289bc588SBarry Smith }
66c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
67289bc588SBarry Smith {
68c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
696abc6512SBarry Smith   int           info;
7055659b69SBarry Smith 
71752f5567SLois Curfman McInnes   if (mat->pivots) {
720452661fSBarry Smith     PetscFree(mat->pivots);
73c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
74752f5567SLois Curfman McInnes     mat->pivots = 0;
75752f5567SLois Curfman McInnes   }
76289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
77ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
78c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
7955659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
80289bc588SBarry Smith   return 0;
81289bc588SBarry Smith }
82289bc588SBarry Smith 
83c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
84289bc588SBarry Smith {
85c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
866abc6512SBarry Smith   int          one = 1, info;
876abc6512SBarry Smith   Scalar       *x, *y;
88289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
89416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
90c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
9148d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
92289bc588SBarry Smith   }
93c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
9448d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
95289bc588SBarry Smith   }
96ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
97ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
9855659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
99289bc588SBarry Smith   return 0;
100289bc588SBarry Smith }
101c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
102da3a660dSBarry Smith {
103c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1046abc6512SBarry Smith   int          one = 1, info;
1056abc6512SBarry Smith   Scalar       *x, *y;
106da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
107416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
108752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
109da3a660dSBarry Smith   if (mat->pivots) {
11048d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
111da3a660dSBarry Smith   }
112da3a660dSBarry Smith   else {
11348d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
114da3a660dSBarry Smith   }
115ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
11655659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
117da3a660dSBarry Smith   return 0;
118da3a660dSBarry Smith }
119c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
120da3a660dSBarry Smith {
121c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1226abc6512SBarry Smith   int          one = 1, info,ierr;
1236abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
124da3a660dSBarry Smith   Vec          tmp = 0;
125da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
126da3a660dSBarry Smith   if (yy == zz) {
12778b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
128c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
12978b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
130da3a660dSBarry Smith   }
131416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
132752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
133da3a660dSBarry Smith   if (mat->pivots) {
13448d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
135da3a660dSBarry Smith   }
136da3a660dSBarry Smith   else {
13748d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
138da3a660dSBarry Smith   }
139ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
140da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
141da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
14255659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
143da3a660dSBarry Smith   return 0;
144da3a660dSBarry Smith }
145c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
146da3a660dSBarry Smith {
147c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1486abc6512SBarry Smith   int           one = 1, info,ierr;
1496abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
150da3a660dSBarry Smith   Vec           tmp;
151da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
152da3a660dSBarry Smith   if (yy == zz) {
15378b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
154c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
15578b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
156da3a660dSBarry Smith   }
157416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
158752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
159da3a660dSBarry Smith   if (mat->pivots) {
16048d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
161da3a660dSBarry Smith   }
162da3a660dSBarry Smith   else {
16348d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
164da3a660dSBarry Smith   }
165ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
166da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
167da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
16855659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
169da3a660dSBarry Smith   return 0;
170da3a660dSBarry Smith }
171289bc588SBarry Smith /* ------------------------------------------------------------------*/
172c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
17320e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
174289bc588SBarry Smith {
175c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
176289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
1776abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
178289bc588SBarry Smith 
179289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
180289bc588SBarry Smith     /* this is a hack fix, should have another version without
181289bc588SBarry Smith        the second BLdot */
182bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
183289bc588SBarry Smith   }
184289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
185289bc588SBarry Smith   while (its--) {
186289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
187289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
188f1747703SBarry Smith #if defined(PETSC_COMPLEX)
189f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
190f1747703SBarry Smith            not happy about returning a double complex */
191f1747703SBarry Smith         int    _i;
192f1747703SBarry Smith         Scalar sum = b[i];
193f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
194f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
195f1747703SBarry Smith         }
196f1747703SBarry Smith         xt = sum;
197f1747703SBarry Smith #else
198289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
199f1747703SBarry Smith #endif
200d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
201289bc588SBarry Smith       }
202289bc588SBarry Smith     }
203289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
204289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
205f1747703SBarry Smith #if defined(PETSC_COMPLEX)
206f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
207f1747703SBarry Smith            not happy about returning a double complex */
208f1747703SBarry Smith         int    _i;
209f1747703SBarry Smith         Scalar sum = b[i];
210f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
211f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
212f1747703SBarry Smith         }
213f1747703SBarry Smith         xt = sum;
214f1747703SBarry Smith #else
215289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
216f1747703SBarry Smith #endif
217d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
218289bc588SBarry Smith       }
219289bc588SBarry Smith     }
220289bc588SBarry Smith   }
221289bc588SBarry Smith   return 0;
222289bc588SBarry Smith }
223289bc588SBarry Smith 
224289bc588SBarry Smith /* -----------------------------------------------------------------*/
225c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
226289bc588SBarry Smith {
227c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
228289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
229289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
230289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
23148d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
23255659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
233289bc588SBarry Smith   return 0;
234289bc588SBarry Smith }
235c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
236289bc588SBarry Smith {
237c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
238289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
239289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
240289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
24148d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
24255659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
243289bc588SBarry Smith   return 0;
244289bc588SBarry Smith }
245c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
246289bc588SBarry Smith {
247c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
248289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2496abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
250289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
251416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
25248d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
25355659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
254289bc588SBarry Smith   return 0;
255289bc588SBarry Smith }
256c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
257289bc588SBarry Smith {
258c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
259289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
260289bc588SBarry Smith   int          _One=1;
2616abc6512SBarry Smith   Scalar       _DOne=1.0;
262289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
263289bc588SBarry Smith   VecGetArray(zz,&z);
264416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
26548d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
26655659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
267289bc588SBarry Smith   return 0;
268289bc588SBarry Smith }
269289bc588SBarry Smith 
270289bc588SBarry Smith /* -----------------------------------------------------------------*/
271c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
272289bc588SBarry Smith {
273c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
274289bc588SBarry Smith   Scalar       *v;
275289bc588SBarry Smith   int          i;
276289bc588SBarry Smith   *ncols = mat->n;
277289bc588SBarry Smith   if (cols) {
2780452661fSBarry Smith     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
27973c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
280289bc588SBarry Smith   }
281289bc588SBarry Smith   if (vals) {
2820452661fSBarry Smith     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
283289bc588SBarry Smith     v = mat->v + row;
28473c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
285289bc588SBarry Smith   }
286289bc588SBarry Smith   return 0;
287289bc588SBarry Smith }
288c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
289289bc588SBarry Smith {
2900452661fSBarry Smith   if (cols) { PetscFree(*cols); }
2910452661fSBarry Smith   if (vals) { PetscFree(*vals); }
292289bc588SBarry Smith   return 0;
293289bc588SBarry Smith }
294289bc588SBarry Smith /* ----------------------------------------------------------------*/
295c0bbcb79SLois Curfman McInnes static int MatInsert_SeqDense(Mat A,int m,int *indexm,int n,
296e8d4e0b9SBarry Smith                                                  int *indexn,Scalar *v,InsertMode addv)
297289bc588SBarry Smith {
298c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
299289bc588SBarry Smith   int          i,j;
300d6dfbf8fSBarry Smith 
301289bc588SBarry Smith   if (!mat->roworiented) {
302dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
303289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
304289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
305289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
306289bc588SBarry Smith         }
307289bc588SBarry Smith       }
308289bc588SBarry Smith     }
309289bc588SBarry Smith     else {
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     }
316e8d4e0b9SBarry Smith   }
317e8d4e0b9SBarry Smith   else {
318dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
319e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
320e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
321e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
322e8d4e0b9SBarry Smith         }
323e8d4e0b9SBarry Smith       }
324e8d4e0b9SBarry Smith     }
325289bc588SBarry Smith     else {
326289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
327289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
328289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
329289bc588SBarry Smith         }
330289bc588SBarry Smith       }
331289bc588SBarry Smith     }
332e8d4e0b9SBarry Smith   }
333289bc588SBarry Smith   return 0;
334289bc588SBarry Smith }
335e8d4e0b9SBarry Smith 
336289bc588SBarry Smith /* -----------------------------------------------------------------*/
33755659b69SBarry Smith static int MatCopyPrivate_SeqDense(Mat A,Mat *newmat,int cpvalues)
338289bc588SBarry Smith {
339c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
340289bc588SBarry Smith   int          ierr;
341289bc588SBarry Smith   Mat          newi;
34248d91487SBarry Smith 
343c0bbcb79SLois Curfman McInnes   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,&newi); CHKERRQ(ierr);
344ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
34555659b69SBarry Smith   if (cpvalues == COPY_VALUES) {
346416022c9SBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
34755659b69SBarry Smith   }
348289bc588SBarry Smith   *newmat = newi;
349289bc588SBarry Smith   return 0;
350289bc588SBarry Smith }
351289bc588SBarry Smith 
35288e32edaSLois Curfman McInnes #include "sysio.h"
35388e32edaSLois Curfman McInnes 
35488e32edaSLois Curfman McInnes int MatLoad_SeqDense(Viewer bview,MatType type,Mat *A)
35588e32edaSLois Curfman McInnes {
35688e32edaSLois Curfman McInnes   Mat_SeqDense *a;
35788e32edaSLois Curfman McInnes   Mat          B;
35888e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
35988e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
36088e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
36188e32edaSLois Curfman McInnes   PetscObject  vobj = (PetscObject) bview;
36288e32edaSLois Curfman McInnes   MPI_Comm     comm = vobj->comm;
36388e32edaSLois Curfman McInnes 
36488e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
36588e32edaSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
36688e32edaSLois Curfman McInnes   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
36788e32edaSLois Curfman McInnes   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
36888e32edaSLois Curfman McInnes   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
36988e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
37088e32edaSLois Curfman McInnes 
37188e32edaSLois Curfman McInnes   /* read row lengths */
3720452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
37388e32edaSLois Curfman McInnes   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
37488e32edaSLois Curfman McInnes 
37588e32edaSLois Curfman McInnes   /* create our matrix */
37688e32edaSLois Curfman McInnes   ierr = MatCreateSeqDense(comm,M,N,A); CHKERRQ(ierr);
37788e32edaSLois Curfman McInnes   B = *A;
37888e32edaSLois Curfman McInnes   a = (Mat_SeqDense *) B->data;
37988e32edaSLois Curfman McInnes   v = a->v;
38088e32edaSLois Curfman McInnes 
38188e32edaSLois Curfman McInnes   /* read column indices and nonzeros */
3820452661fSBarry Smith   cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
38388e32edaSLois Curfman McInnes   ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr);
3840452661fSBarry Smith   vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
38588e32edaSLois Curfman McInnes   ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr);
38688e32edaSLois Curfman McInnes 
38788e32edaSLois Curfman McInnes   /* insert into matrix */
38888e32edaSLois Curfman McInnes   for ( i=0; i<M; i++ ) {
38988e32edaSLois Curfman McInnes     for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
39088e32edaSLois Curfman McInnes     svals += rowlengths[i]; scols += rowlengths[i];
39188e32edaSLois Curfman McInnes   }
3920452661fSBarry Smith   PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
39388e32edaSLois Curfman McInnes 
39488e32edaSLois Curfman McInnes   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
39588e32edaSLois Curfman McInnes   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
39688e32edaSLois Curfman McInnes   return 0;
39788e32edaSLois Curfman McInnes }
39888e32edaSLois Curfman McInnes 
399932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
400932b0c3eSLois Curfman McInnes #include "sysio.h"
401932b0c3eSLois Curfman McInnes 
402932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
403289bc588SBarry Smith {
404932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
405932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
406d636dbe3SBarry Smith   FILE         *fd;
407932b0c3eSLois Curfman McInnes   char         *outputname;
408932b0c3eSLois Curfman McInnes   Scalar       *v;
409932b0c3eSLois Curfman McInnes 
410932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
411932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
412932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetFormat_Private(viewer,&format);
413f72585f2SLois Curfman McInnes   if (format == FILE_FORMAT_INFO) {
414932b0c3eSLois Curfman McInnes     ;  /* do nothing for now */
415f72585f2SLois Curfman McInnes   }
416f72585f2SLois Curfman McInnes   else {
417932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
418932b0c3eSLois Curfman McInnes       v = a->v + i;
419932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
420289bc588SBarry Smith #if defined(PETSC_COMPLEX)
421932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
422289bc588SBarry Smith #else
423932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
424289bc588SBarry Smith #endif
425289bc588SBarry Smith       }
4268e37dc09SBarry Smith       fprintf(fd,"\n");
427289bc588SBarry Smith     }
428da3a660dSBarry Smith   }
429932b0c3eSLois Curfman McInnes   fflush(fd);
430289bc588SBarry Smith   return 0;
431289bc588SBarry Smith }
432289bc588SBarry Smith 
433932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
434932b0c3eSLois Curfman McInnes {
435932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
436932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
437932b0c3eSLois Curfman McInnes   Scalar       *v, *anonz;
438932b0c3eSLois Curfman McInnes 
439932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
4400452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
441932b0c3eSLois Curfman McInnes   col_lens[0] = MAT_COOKIE;
442932b0c3eSLois Curfman McInnes   col_lens[1] = m;
443932b0c3eSLois Curfman McInnes   col_lens[2] = n;
444932b0c3eSLois Curfman McInnes   col_lens[3] = nz;
445932b0c3eSLois Curfman McInnes 
446932b0c3eSLois Curfman McInnes   /* store lengths of each row and write (including header) to file */
447932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) col_lens[4+i] = n;
448932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr);
449932b0c3eSLois Curfman McInnes 
450932b0c3eSLois Curfman McInnes   /* Possibly should write in smaller increments, not whole matrix at once? */
451932b0c3eSLois Curfman McInnes  /* store column indices (zero start index) */
452932b0c3eSLois Curfman McInnes   ict = 0;
453932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
454932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) col_lens[ict++] = j;
455932b0c3eSLois Curfman McInnes   }
456932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr);
4570452661fSBarry Smith   PetscFree(col_lens);
458932b0c3eSLois Curfman McInnes 
459932b0c3eSLois Curfman McInnes   /* store nonzero values */
4600452661fSBarry Smith   anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
461932b0c3eSLois Curfman McInnes   ict = 0;
462932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
463932b0c3eSLois Curfman McInnes     v = a->v + i;
464932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) {
465932b0c3eSLois Curfman McInnes       anonz[ict++] = *v; v += a->m;
466932b0c3eSLois Curfman McInnes     }
467932b0c3eSLois Curfman McInnes   }
468932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr);
4690452661fSBarry Smith   PetscFree(anonz);
470932b0c3eSLois Curfman McInnes   return 0;
471932b0c3eSLois Curfman McInnes }
472932b0c3eSLois Curfman McInnes 
473932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
474932b0c3eSLois Curfman McInnes {
475932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
476932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
477932b0c3eSLois Curfman McInnes   PetscObject  vobj = (PetscObject) viewer;
478932b0c3eSLois Curfman McInnes 
479932b0c3eSLois Curfman McInnes   if (!viewer) {
480932b0c3eSLois Curfman McInnes     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
481932b0c3eSLois Curfman McInnes   }
482932b0c3eSLois Curfman McInnes   if (vobj->cookie == VIEWER_COOKIE) {
483932b0c3eSLois Curfman McInnes     if (vobj->type == MATLAB_VIEWER) {
484932b0c3eSLois Curfman McInnes       return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
485932b0c3eSLois Curfman McInnes     }
486932b0c3eSLois Curfman McInnes     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) {
487932b0c3eSLois Curfman McInnes       return MatView_SeqDense_ASCII(A,viewer);
488932b0c3eSLois Curfman McInnes     }
489932b0c3eSLois Curfman McInnes     else if (vobj->type == BINARY_FILE_VIEWER) {
490932b0c3eSLois Curfman McInnes       return MatView_SeqDense_Binary(A,viewer);
491932b0c3eSLois Curfman McInnes     }
492932b0c3eSLois Curfman McInnes   }
493932b0c3eSLois Curfman McInnes   return 0;
494932b0c3eSLois Curfman McInnes }
495289bc588SBarry Smith 
496ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
497289bc588SBarry Smith {
498289bc588SBarry Smith   Mat          mat = (Mat) obj;
499ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
500a5a9c739SBarry Smith #if defined(PETSC_LOG)
501a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
502a5a9c739SBarry Smith #endif
5030452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
5040452661fSBarry Smith   PetscFree(l);
505a5a9c739SBarry Smith   PLogObjectDestroy(mat);
5060452661fSBarry Smith   PetscHeaderDestroy(mat);
507289bc588SBarry Smith   return 0;
508289bc588SBarry Smith }
509289bc588SBarry Smith 
510c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
511289bc588SBarry Smith {
512c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
513d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
514d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
51548b35521SBarry Smith 
516d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
51748b35521SBarry Smith   if (!matout) { /* in place transpose */
518d9f96c7cSLois Curfman McInnes     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
519d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
520289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
521d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
522d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
523d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
524289bc588SBarry Smith       }
525289bc588SBarry Smith     }
52648b35521SBarry Smith   }
527d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
528d3e5ee88SLois Curfman McInnes     int          ierr;
529d3e5ee88SLois Curfman McInnes     Mat          tmat;
530ec8511deSBarry Smith     Mat_SeqDense *tmatd;
531d3e5ee88SLois Curfman McInnes     Scalar       *v2;
532c0bbcb79SLois Curfman McInnes     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr);
533ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
5340de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
535d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
5360de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
537d3e5ee88SLois Curfman McInnes     }
538d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
539d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
540d3e5ee88SLois Curfman McInnes     *matout = tmat;
54148b35521SBarry Smith   }
542289bc588SBarry Smith   return 0;
543289bc588SBarry Smith }
544289bc588SBarry Smith 
545c0bbcb79SLois Curfman McInnes static int MatEqual_SeqDense(Mat A1,Mat A2)
546289bc588SBarry Smith {
547c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
548c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
549289bc588SBarry Smith   int          i;
550289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
551289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
552289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
553289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
554289bc588SBarry Smith     if (*v1 != *v2) return 0;
555289bc588SBarry Smith     v1++; v2++;
556289bc588SBarry Smith   }
557289bc588SBarry Smith   return 1;
558289bc588SBarry Smith }
559289bc588SBarry Smith 
560c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
561289bc588SBarry Smith {
562c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
5636abc6512SBarry Smith   int          i, n;
5646abc6512SBarry Smith   Scalar       *x;
565289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
566ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
567289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
568289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
569289bc588SBarry Smith   }
570289bc588SBarry Smith   return 0;
571289bc588SBarry Smith }
572289bc588SBarry Smith 
573c0bbcb79SLois Curfman McInnes static int MatScale_SeqDense(Mat A,Vec ll,Vec rr)
574289bc588SBarry Smith {
575c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
576da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
577da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
57855659b69SBarry Smith 
57928988994SBarry Smith   if (ll) {
580da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
581ec8511deSBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size");
58255659b69SBarry Smith     PLogFlops(n*m);
583da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
584da3a660dSBarry Smith       x = l[i];
585da3a660dSBarry Smith       v = mat->v + i;
586da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
587da3a660dSBarry Smith     }
588da3a660dSBarry Smith   }
58928988994SBarry Smith   if (rr) {
590da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
591ec8511deSBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size");
59255659b69SBarry Smith     PLogFlops(n*m);
593da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
594da3a660dSBarry Smith       x = r[i];
595da3a660dSBarry Smith       v = mat->v + i*m;
596da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
597da3a660dSBarry Smith     }
598da3a660dSBarry Smith   }
599289bc588SBarry Smith   return 0;
600289bc588SBarry Smith }
601289bc588SBarry Smith 
602cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
603289bc588SBarry Smith {
604c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
605289bc588SBarry Smith   Scalar       *v = mat->v;
606289bc588SBarry Smith   double       sum = 0.0;
607289bc588SBarry Smith   int          i, j;
60855659b69SBarry Smith 
609289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
610289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
611289bc588SBarry Smith #if defined(PETSC_COMPLEX)
612289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
613289bc588SBarry Smith #else
614289bc588SBarry Smith       sum += (*v)*(*v); v++;
615289bc588SBarry Smith #endif
616289bc588SBarry Smith     }
617289bc588SBarry Smith     *norm = sqrt(sum);
61855659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
619289bc588SBarry Smith   }
620289bc588SBarry Smith   else if (type == NORM_1) {
621289bc588SBarry Smith     *norm = 0.0;
622289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
623289bc588SBarry Smith       sum = 0.0;
624289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
625*33a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
626289bc588SBarry Smith       }
627289bc588SBarry Smith       if (sum > *norm) *norm = sum;
628289bc588SBarry Smith     }
62955659b69SBarry Smith     PLogFlops(mat->n*mat->m);
630289bc588SBarry Smith   }
631289bc588SBarry Smith   else if (type == NORM_INFINITY) {
632289bc588SBarry Smith     *norm = 0.0;
633289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
634289bc588SBarry Smith       v = mat->v + j;
635289bc588SBarry Smith       sum = 0.0;
636289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
637cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
638289bc588SBarry Smith       }
639289bc588SBarry Smith       if (sum > *norm) *norm = sum;
640289bc588SBarry Smith     }
64155659b69SBarry Smith     PLogFlops(mat->n*mat->m);
642289bc588SBarry Smith   }
643289bc588SBarry Smith   else {
64448d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
645289bc588SBarry Smith   }
646289bc588SBarry Smith   return 0;
647289bc588SBarry Smith }
648289bc588SBarry Smith 
649c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
650289bc588SBarry Smith {
651c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
652289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
653289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
654c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
65588e32edaSLois Curfman McInnes            op == COLUMNS_SORTED ||
656c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
657c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
658c0bbcb79SLois Curfman McInnes            op == NO_NEW_NONZERO_LOCATIONS ||
659c0bbcb79SLois Curfman McInnes            op == YES_NEW_NONZERO_LOCATIONS ||
660c0bbcb79SLois Curfman McInnes            op == NO_NEW_DIAGONALS ||
661c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
662c0bbcb79SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n");
663c0bbcb79SLois Curfman McInnes   else
6644d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
665289bc588SBarry Smith   return 0;
666289bc588SBarry Smith }
667289bc588SBarry Smith 
668ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
6696f0a148fSBarry Smith {
670ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
671cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
6726f0a148fSBarry Smith   return 0;
6736f0a148fSBarry Smith }
6746f0a148fSBarry Smith 
675ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
6766f0a148fSBarry Smith {
677ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
6786abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
6796f0a148fSBarry Smith   Scalar       *slot;
68055659b69SBarry Smith 
68178b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
68278b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
6836f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
6846f0a148fSBarry Smith     slot = l->v + rows[i];
6856f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
6866f0a148fSBarry Smith   }
6876f0a148fSBarry Smith   if (diag) {
6886f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
6896f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
690260925b8SBarry Smith       *slot = *diag;
6916f0a148fSBarry Smith     }
6926f0a148fSBarry Smith   }
693260925b8SBarry Smith   ISRestoreIndices(is,&rows);
6946f0a148fSBarry Smith   return 0;
6956f0a148fSBarry Smith }
696557bce09SLois Curfman McInnes 
697c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
698557bce09SLois Curfman McInnes {
699c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
700557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
701557bce09SLois Curfman McInnes   return 0;
702557bce09SLois Curfman McInnes }
703557bce09SLois Curfman McInnes 
704c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
705d3e5ee88SLois Curfman McInnes {
706c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
707d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
708d3e5ee88SLois Curfman McInnes   return 0;
709d3e5ee88SLois Curfman McInnes }
710d3e5ee88SLois Curfman McInnes 
711c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
71264e87e97SBarry Smith {
713c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
71464e87e97SBarry Smith   *array = mat->v;
71564e87e97SBarry Smith   return 0;
71664e87e97SBarry Smith }
7170754003eSLois Curfman McInnes 
7180754003eSLois Curfman McInnes 
719c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
7200754003eSLois Curfman McInnes {
721ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
7220754003eSLois Curfman McInnes }
7230754003eSLois Curfman McInnes 
724cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
725cddf8d76SBarry Smith                                     Mat *submat)
7260754003eSLois Curfman McInnes {
727c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
7280754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
729160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
7300754003eSLois Curfman McInnes   Scalar       *vwork, *val;
7310754003eSLois Curfman McInnes   Mat          newmat;
7320754003eSLois Curfman McInnes 
73378b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
73478b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
73578b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
73678b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
7370754003eSLois Curfman McInnes 
7380452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
7390452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
7400452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
741cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
7420754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
7430754003eSLois Curfman McInnes 
7440754003eSLois Curfman McInnes   /* Create and fill new matrix */
745c0bbcb79SLois Curfman McInnes   ierr = MatCreateSeqDense(A->comm,nrows,ncols,&newmat);CHKERRQ(ierr);
7460754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
7470754003eSLois Curfman McInnes     nznew = 0;
7480754003eSLois Curfman McInnes     val   = mat->v + irow[i];
7490754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
7500754003eSLois Curfman McInnes       if (smap[j]) {
7510754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
7520754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
7530754003eSLois Curfman McInnes       }
7540754003eSLois Curfman McInnes     }
755dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
75678b31e54SBarry Smith            CHKERRQ(ierr);
7570754003eSLois Curfman McInnes   }
75878b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
75978b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
7600754003eSLois Curfman McInnes 
7610754003eSLois Curfman McInnes   /* Free work space */
7620452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
76378b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
76478b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
7650754003eSLois Curfman McInnes   *submat = newmat;
7660754003eSLois Curfman McInnes   return 0;
7670754003eSLois Curfman McInnes }
7680754003eSLois Curfman McInnes 
769289bc588SBarry Smith /* -------------------------------------------------------------------*/
770ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense,
771ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
772ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
773ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
774ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
775ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
776ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
777ec8511deSBarry Smith        MatRelax_SeqDense,
778ec8511deSBarry Smith        MatTranspose_SeqDense,
779ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
780ec8511deSBarry Smith        MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense,
781289bc588SBarry Smith        0,0,
782ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
783ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
784ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
785ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
786ec8511deSBarry Smith        0,0,MatGetArray_SeqDense,0,0,
787ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
7881987afe7SBarry Smith        MatCopyPrivate_SeqDense,0,0,0,0,
7891987afe7SBarry Smith        MatAXPY_SeqDense};
7900754003eSLois Curfman McInnes 
7914b828684SBarry Smith /*@C
792fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
793d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
794d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
795289bc588SBarry Smith 
79620563c6bSBarry Smith    Input Parameters:
7970c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
7980c775827SLois Curfman McInnes .  m - number of rows
7990c775827SLois Curfman McInnes .  n - number of column
80020563c6bSBarry Smith 
80120563c6bSBarry Smith    Output Parameter:
8020c775827SLois Curfman McInnes .  newmat - the matrix
80320563c6bSBarry Smith 
804dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
805d65003e9SLois Curfman McInnes 
806dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
80720563c6bSBarry Smith @*/
808fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat)
809289bc588SBarry Smith {
810ec8511deSBarry Smith   int          size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar);
811289bc588SBarry Smith   Mat          mat;
812ec8511deSBarry Smith   Mat_SeqDense *l;
81355659b69SBarry Smith 
814289bc588SBarry Smith   *newmat        = 0;
8150452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
816a5a9c739SBarry Smith   PLogObjectCreate(mat);
8170452661fSBarry Smith   l              = (Mat_SeqDense *) PetscMalloc(size); CHKPTRQ(l);
818416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
819ec8511deSBarry Smith   mat->destroy   = MatDestroy_SeqDense;
820ec8511deSBarry Smith   mat->view      = MatView_SeqDense;
821289bc588SBarry Smith   mat->data      = (void *) l;
822289bc588SBarry Smith   mat->factor    = 0;
823752f5567SLois Curfman McInnes   PLogObjectMemory(mat,sizeof(struct _Mat) + size);
824d6dfbf8fSBarry Smith 
825289bc588SBarry Smith   l->m           = m;
826289bc588SBarry Smith   l->n           = n;
827289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
828289bc588SBarry Smith   l->pivots      = 0;
829289bc588SBarry Smith   l->roworiented = 1;
830d6dfbf8fSBarry Smith 
831cddf8d76SBarry Smith   PetscMemzero(l->v,m*n*sizeof(Scalar));
832289bc588SBarry Smith   *newmat = mat;
833289bc588SBarry Smith   return 0;
834289bc588SBarry Smith }
835289bc588SBarry Smith 
836c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
837289bc588SBarry Smith {
838c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
839c0bbcb79SLois Curfman McInnes   return MatCreateSeqDense(A->comm,m->m,m->n,newmat);
840289bc588SBarry Smith }
841