xref: /petsc/src/mat/impls/dense/seq/dense.c (revision d9f96c7ca71afc96fd7309d95e12b2b98fa1305e)
1cb512458SBarry Smith #ifndef lint
2*d9f96c7cSLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.68 1995/10/22 04:19:49 bsmith Exp curfman $";
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) {
3748d91487SBarry 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) {
72752f5567SLois Curfman McInnes     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) {
27878b31e54SBarry 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) {
28278b31e54SBarry 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 {
29078b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
29178b31e54SBarry 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 
352932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
353932b0c3eSLois Curfman McInnes #include "sysio.h"
354932b0c3eSLois Curfman McInnes 
355932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
356289bc588SBarry Smith {
357932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
358932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
359d636dbe3SBarry Smith   FILE         *fd;
360932b0c3eSLois Curfman McInnes   char         *outputname;
361932b0c3eSLois Curfman McInnes   Scalar       *v;
362932b0c3eSLois Curfman McInnes 
363932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
364932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
365932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetFormat_Private(viewer,&format);
366f72585f2SLois Curfman McInnes   if (format == FILE_FORMAT_INFO) {
367932b0c3eSLois Curfman McInnes     ;  /* do nothing for now */
368f72585f2SLois Curfman McInnes   }
369f72585f2SLois Curfman McInnes   else {
370932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
371932b0c3eSLois Curfman McInnes       v = a->v + i;
372932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
373289bc588SBarry Smith #if defined(PETSC_COMPLEX)
374932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
375289bc588SBarry Smith #else
376932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
377289bc588SBarry Smith #endif
378289bc588SBarry Smith       }
3798e37dc09SBarry Smith       fprintf(fd,"\n");
380289bc588SBarry Smith     }
381da3a660dSBarry Smith   }
382932b0c3eSLois Curfman McInnes   fflush(fd);
383289bc588SBarry Smith   return 0;
384289bc588SBarry Smith }
385289bc588SBarry Smith 
386932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
387932b0c3eSLois Curfman McInnes {
388932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
389932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
390932b0c3eSLois Curfman McInnes   Scalar       *v, *anonz;
391932b0c3eSLois Curfman McInnes 
392932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
393932b0c3eSLois Curfman McInnes   col_lens = (int *) PETSCMALLOC( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
394932b0c3eSLois Curfman McInnes   col_lens[0] = MAT_COOKIE;
395932b0c3eSLois Curfman McInnes   col_lens[1] = m;
396932b0c3eSLois Curfman McInnes   col_lens[2] = n;
397932b0c3eSLois Curfman McInnes   col_lens[3] = nz;
398932b0c3eSLois Curfman McInnes 
399932b0c3eSLois Curfman McInnes   /* store lengths of each row and write (including header) to file */
400932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) col_lens[4+i] = n;
401932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr);
402932b0c3eSLois Curfman McInnes 
403932b0c3eSLois Curfman McInnes   /* Possibly should write in smaller increments, not whole matrix at once? */
404932b0c3eSLois Curfman McInnes  /* store column indices (zero start index) */
405932b0c3eSLois Curfman McInnes   ict = 0;
406932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
407932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) col_lens[ict++] = j;
408932b0c3eSLois Curfman McInnes   }
409932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr);
410932b0c3eSLois Curfman McInnes   PETSCFREE(col_lens);
411932b0c3eSLois Curfman McInnes 
412932b0c3eSLois Curfman McInnes   /* store nonzero values */
413932b0c3eSLois Curfman McInnes   anonz = (Scalar *) PETSCMALLOC((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
414932b0c3eSLois Curfman McInnes   ict = 0;
415932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
416932b0c3eSLois Curfman McInnes     v = a->v + i;
417932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) {
418932b0c3eSLois Curfman McInnes       anonz[ict++] = *v; v += a->m;
419932b0c3eSLois Curfman McInnes     }
420932b0c3eSLois Curfman McInnes   }
421932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr);
422932b0c3eSLois Curfman McInnes   PETSCFREE(anonz);
423932b0c3eSLois Curfman McInnes   return 0;
424932b0c3eSLois Curfman McInnes }
425932b0c3eSLois Curfman McInnes 
426932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
427932b0c3eSLois Curfman McInnes {
428932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
429932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
430932b0c3eSLois Curfman McInnes   PetscObject  vobj = (PetscObject) viewer;
431932b0c3eSLois Curfman McInnes 
432932b0c3eSLois Curfman McInnes   if (!viewer) {
433932b0c3eSLois Curfman McInnes     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
434932b0c3eSLois Curfman McInnes   }
435932b0c3eSLois Curfman McInnes   if (vobj->cookie == VIEWER_COOKIE) {
436932b0c3eSLois Curfman McInnes     if (vobj->type == MATLAB_VIEWER) {
437932b0c3eSLois Curfman McInnes       return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
438932b0c3eSLois Curfman McInnes     }
439932b0c3eSLois Curfman McInnes     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) {
440932b0c3eSLois Curfman McInnes       return MatView_SeqDense_ASCII(A,viewer);
441932b0c3eSLois Curfman McInnes     }
442932b0c3eSLois Curfman McInnes     else if (vobj->type == BINARY_FILE_VIEWER) {
443932b0c3eSLois Curfman McInnes       return MatView_SeqDense_Binary(A,viewer);
444932b0c3eSLois Curfman McInnes     }
445932b0c3eSLois Curfman McInnes   }
446932b0c3eSLois Curfman McInnes   return 0;
447932b0c3eSLois Curfman McInnes }
448289bc588SBarry Smith 
449ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
450289bc588SBarry Smith {
451289bc588SBarry Smith   Mat          mat = (Mat) obj;
452ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
453a5a9c739SBarry Smith #if defined(PETSC_LOG)
454a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
455a5a9c739SBarry Smith #endif
45678b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
45778b31e54SBarry Smith   PETSCFREE(l);
458a5a9c739SBarry Smith   PLogObjectDestroy(mat);
4599e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
460289bc588SBarry Smith   return 0;
461289bc588SBarry Smith }
462289bc588SBarry Smith 
463c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
464289bc588SBarry Smith {
465c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
466d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
467d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
46848b35521SBarry Smith 
469d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
47048b35521SBarry Smith   if (!matout) { /* in place transpose */
471*d9f96c7cSLois Curfman McInnes     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
472d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
473289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
474d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
475d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
476d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
477289bc588SBarry Smith       }
478289bc588SBarry Smith     }
47948b35521SBarry Smith   }
480d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
481d3e5ee88SLois Curfman McInnes     int          ierr;
482d3e5ee88SLois Curfman McInnes     Mat          tmat;
483ec8511deSBarry Smith     Mat_SeqDense *tmatd;
484d3e5ee88SLois Curfman McInnes     Scalar       *v2;
485c0bbcb79SLois Curfman McInnes     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr);
486ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
4870de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
488d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
4890de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
490d3e5ee88SLois Curfman McInnes     }
491d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
492d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
493d3e5ee88SLois Curfman McInnes     *matout = tmat;
49448b35521SBarry Smith   }
495289bc588SBarry Smith   return 0;
496289bc588SBarry Smith }
497289bc588SBarry Smith 
498c0bbcb79SLois Curfman McInnes static int MatEqual_SeqDense(Mat A1,Mat A2)
499289bc588SBarry Smith {
500c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
501c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
502289bc588SBarry Smith   int          i;
503289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
504289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
505289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
506289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
507289bc588SBarry Smith     if (*v1 != *v2) return 0;
508289bc588SBarry Smith     v1++; v2++;
509289bc588SBarry Smith   }
510289bc588SBarry Smith   return 1;
511289bc588SBarry Smith }
512289bc588SBarry Smith 
513c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
514289bc588SBarry Smith {
515c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
5166abc6512SBarry Smith   int          i, n;
5176abc6512SBarry Smith   Scalar       *x;
518289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
519ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
520289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
521289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
522289bc588SBarry Smith   }
523289bc588SBarry Smith   return 0;
524289bc588SBarry Smith }
525289bc588SBarry Smith 
526c0bbcb79SLois Curfman McInnes static int MatScale_SeqDense(Mat A,Vec ll,Vec rr)
527289bc588SBarry Smith {
528c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
529da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
530da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
53155659b69SBarry Smith 
53228988994SBarry Smith   if (ll) {
533da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
534ec8511deSBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size");
53555659b69SBarry Smith     PLogFlops(n*m);
536da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
537da3a660dSBarry Smith       x = l[i];
538da3a660dSBarry Smith       v = mat->v + i;
539da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
540da3a660dSBarry Smith     }
541da3a660dSBarry Smith   }
54228988994SBarry Smith   if (rr) {
543da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
544ec8511deSBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size");
54555659b69SBarry Smith     PLogFlops(n*m);
546da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
547da3a660dSBarry Smith       x = r[i];
548da3a660dSBarry Smith       v = mat->v + i*m;
549da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
550da3a660dSBarry Smith     }
551da3a660dSBarry Smith   }
552289bc588SBarry Smith   return 0;
553289bc588SBarry Smith }
554289bc588SBarry Smith 
555c0bbcb79SLois Curfman McInnes static int MatNorm_SeqDense(Mat A,MatNormType type,double *norm)
556289bc588SBarry Smith {
557c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
558289bc588SBarry Smith   Scalar       *v = mat->v;
559289bc588SBarry Smith   double       sum = 0.0;
560289bc588SBarry Smith   int          i, j;
56155659b69SBarry Smith 
562289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
563289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
564289bc588SBarry Smith #if defined(PETSC_COMPLEX)
565289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
566289bc588SBarry Smith #else
567289bc588SBarry Smith       sum += (*v)*(*v); v++;
568289bc588SBarry Smith #endif
569289bc588SBarry Smith     }
570289bc588SBarry Smith     *norm = sqrt(sum);
57155659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
572289bc588SBarry Smith   }
573289bc588SBarry Smith   else if (type == NORM_1) {
574289bc588SBarry Smith     *norm = 0.0;
575289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
576289bc588SBarry Smith       sum = 0.0;
577289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
578289bc588SBarry Smith #if defined(PETSC_COMPLEX)
579289bc588SBarry Smith         sum += abs(*v++);
580289bc588SBarry Smith #else
581289bc588SBarry Smith         sum += fabs(*v++);
582289bc588SBarry Smith #endif
583289bc588SBarry Smith       }
584289bc588SBarry Smith       if (sum > *norm) *norm = sum;
585289bc588SBarry Smith     }
58655659b69SBarry Smith     PLogFlops(mat->n*mat->m);
587289bc588SBarry Smith   }
588289bc588SBarry Smith   else if (type == NORM_INFINITY) {
589289bc588SBarry Smith     *norm = 0.0;
590289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
591289bc588SBarry Smith       v = mat->v + j;
592289bc588SBarry Smith       sum = 0.0;
593289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
594289bc588SBarry Smith #if defined(PETSC_COMPLEX)
595289bc588SBarry Smith         sum += abs(*v); v += mat->m;
596289bc588SBarry Smith #else
597289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
598289bc588SBarry Smith #endif
599289bc588SBarry Smith       }
600289bc588SBarry Smith       if (sum > *norm) *norm = sum;
601289bc588SBarry Smith     }
60255659b69SBarry Smith     PLogFlops(mat->n*mat->m);
603289bc588SBarry Smith   }
604289bc588SBarry Smith   else {
60548d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
606289bc588SBarry Smith   }
607289bc588SBarry Smith   return 0;
608289bc588SBarry Smith }
609289bc588SBarry Smith 
610c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
611289bc588SBarry Smith {
612c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
613289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
614289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
615c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
616c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
617c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
618c0bbcb79SLois Curfman McInnes            op == NO_NEW_NONZERO_LOCATIONS ||
619c0bbcb79SLois Curfman McInnes            op == YES_NEW_NONZERO_LOCATIONS ||
620c0bbcb79SLois Curfman McInnes            op == NO_NEW_DIAGONALS ||
621c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
622c0bbcb79SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n");
623c0bbcb79SLois Curfman McInnes   else
6244d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
625289bc588SBarry Smith   return 0;
626289bc588SBarry Smith }
627289bc588SBarry Smith 
628ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
6296f0a148fSBarry Smith {
630ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
631416022c9SBarry Smith   PetscZero(l->v,l->m*l->n*sizeof(Scalar));
6326f0a148fSBarry Smith   return 0;
6336f0a148fSBarry Smith }
6346f0a148fSBarry Smith 
635ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
6366f0a148fSBarry Smith {
637ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
6386abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
6396f0a148fSBarry Smith   Scalar       *slot;
64055659b69SBarry Smith 
64178b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
64278b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
6436f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
6446f0a148fSBarry Smith     slot = l->v + rows[i];
6456f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
6466f0a148fSBarry Smith   }
6476f0a148fSBarry Smith   if (diag) {
6486f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
6496f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
650260925b8SBarry Smith       *slot = *diag;
6516f0a148fSBarry Smith     }
6526f0a148fSBarry Smith   }
653260925b8SBarry Smith   ISRestoreIndices(is,&rows);
6546f0a148fSBarry Smith   return 0;
6556f0a148fSBarry Smith }
656557bce09SLois Curfman McInnes 
657c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
658557bce09SLois Curfman McInnes {
659c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
660557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
661557bce09SLois Curfman McInnes   return 0;
662557bce09SLois Curfman McInnes }
663557bce09SLois Curfman McInnes 
664c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
665d3e5ee88SLois Curfman McInnes {
666c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
667d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
668d3e5ee88SLois Curfman McInnes   return 0;
669d3e5ee88SLois Curfman McInnes }
670d3e5ee88SLois Curfman McInnes 
671c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
67264e87e97SBarry Smith {
673c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
67464e87e97SBarry Smith   *array = mat->v;
67564e87e97SBarry Smith   return 0;
67664e87e97SBarry Smith }
6770754003eSLois Curfman McInnes 
6780754003eSLois Curfman McInnes 
679c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
6800754003eSLois Curfman McInnes {
681ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
6820754003eSLois Curfman McInnes }
6830754003eSLois Curfman McInnes 
684c0bbcb79SLois Curfman McInnes static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,Mat *submat)
6850754003eSLois Curfman McInnes {
686c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
6870754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
688160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
6890754003eSLois Curfman McInnes   Scalar       *vwork, *val;
6900754003eSLois Curfman McInnes   Mat          newmat;
6910754003eSLois Curfman McInnes 
69278b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
69378b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
69478b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
69578b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
6960754003eSLois Curfman McInnes 
69778b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
69878b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
69978b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
700416022c9SBarry Smith   PetscZero((char*)smap,oldcols*sizeof(int));
7010754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
7020754003eSLois Curfman McInnes 
7030754003eSLois Curfman McInnes   /* Create and fill new matrix */
704c0bbcb79SLois Curfman McInnes   ierr = MatCreateSeqDense(A->comm,nrows,ncols,&newmat);CHKERRQ(ierr);
7050754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
7060754003eSLois Curfman McInnes     nznew = 0;
7070754003eSLois Curfman McInnes     val   = mat->v + irow[i];
7080754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
7090754003eSLois Curfman McInnes       if (smap[j]) {
7100754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
7110754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
7120754003eSLois Curfman McInnes       }
7130754003eSLois Curfman McInnes     }
714dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
71578b31e54SBarry Smith            CHKERRQ(ierr);
7160754003eSLois Curfman McInnes   }
71778b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
71878b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
7190754003eSLois Curfman McInnes 
7200754003eSLois Curfman McInnes   /* Free work space */
72178b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
72278b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
72378b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
7240754003eSLois Curfman McInnes   *submat = newmat;
7250754003eSLois Curfman McInnes   return 0;
7260754003eSLois Curfman McInnes }
7270754003eSLois Curfman McInnes 
728289bc588SBarry Smith /* -------------------------------------------------------------------*/
729ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense,
730ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
731ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
732ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
733ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
734ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
735ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
736ec8511deSBarry Smith        MatRelax_SeqDense,
737ec8511deSBarry Smith        MatTranspose_SeqDense,
738ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
739ec8511deSBarry Smith        MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense,
740289bc588SBarry Smith        0,0,
741ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
742ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
743ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
744ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
745ec8511deSBarry Smith        0,0,MatGetArray_SeqDense,0,0,
746ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
7471987afe7SBarry Smith        MatCopyPrivate_SeqDense,0,0,0,0,
7481987afe7SBarry Smith        MatAXPY_SeqDense};
7490754003eSLois Curfman McInnes 
7504b828684SBarry Smith /*@C
751fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
752d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
753d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
754289bc588SBarry Smith 
75520563c6bSBarry Smith    Input Parameters:
7560c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
7570c775827SLois Curfman McInnes .  m - number of rows
7580c775827SLois Curfman McInnes .  n - number of column
75920563c6bSBarry Smith 
76020563c6bSBarry Smith    Output Parameter:
7610c775827SLois Curfman McInnes .  newmat - the matrix
76220563c6bSBarry Smith 
763dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
764d65003e9SLois Curfman McInnes 
765dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
76620563c6bSBarry Smith @*/
767fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat)
768289bc588SBarry Smith {
769ec8511deSBarry Smith   int          size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar);
770289bc588SBarry Smith   Mat          mat;
771ec8511deSBarry Smith   Mat_SeqDense *l;
77255659b69SBarry Smith 
773289bc588SBarry Smith   *newmat        = 0;
774ec8511deSBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
775a5a9c739SBarry Smith   PLogObjectCreate(mat);
776ec8511deSBarry Smith   l              = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l);
777416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
778ec8511deSBarry Smith   mat->destroy   = MatDestroy_SeqDense;
779ec8511deSBarry Smith   mat->view      = MatView_SeqDense;
780289bc588SBarry Smith   mat->data      = (void *) l;
781289bc588SBarry Smith   mat->factor    = 0;
782752f5567SLois Curfman McInnes   PLogObjectMemory(mat,sizeof(struct _Mat) + size);
783d6dfbf8fSBarry Smith 
784289bc588SBarry Smith   l->m           = m;
785289bc588SBarry Smith   l->n           = n;
786289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
787289bc588SBarry Smith   l->pivots      = 0;
788289bc588SBarry Smith   l->roworiented = 1;
789d6dfbf8fSBarry Smith 
790416022c9SBarry Smith   PetscZero(l->v,m*n*sizeof(Scalar));
791289bc588SBarry Smith   *newmat = mat;
792289bc588SBarry Smith   return 0;
793289bc588SBarry Smith }
794289bc588SBarry Smith 
795c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
796289bc588SBarry Smith {
797c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
798c0bbcb79SLois Curfman McInnes   return MatCreateSeqDense(A->comm,m->m,m->n,newmat);
799289bc588SBarry Smith }
800