xref: /petsc/src/mat/impls/dense/seq/dense.c (revision c0bbcb798f27810f78e937d23f19ea55c7db8561)
1cb512458SBarry Smith #ifndef lint
2*c0bbcb79SLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.63 1995/10/11 15:19:29 bsmith Exp curfman $";
3cb512458SBarry Smith #endif
4289bc588SBarry Smith 
5289bc588SBarry Smith /*
6289bc588SBarry Smith     Standard Fortran style matrices
7289bc588SBarry Smith */
819b02663SBarry Smith #include "petsc.h"
9b16a3bb1SBarry Smith #include "pinclude/plapack.h"
10289bc588SBarry Smith #include "matimpl.h"
11289bc588SBarry Smith #include "math.h"
12289bc588SBarry Smith #include "vec/vecimpl.h"
13b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
14289bc588SBarry Smith 
15289bc588SBarry Smith typedef struct {
16*c0bbcb79SLois Curfman McInnes   Scalar *v;                /* matrix elements */
17*c0bbcb79SLois Curfman McInnes   int    roworiented;       /* if true, row oriented input (default) */
18*c0bbcb79SLois Curfman McInnes   int    m,n,pad;           /* rows, columns, padding */
19289bc588SBarry Smith   int    *pivots;           /* pivots in LU factorization */
20ec8511deSBarry Smith } Mat_SeqDense;
21289bc588SBarry Smith 
221987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
231987afe7SBarry Smith {
241987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
251987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
261987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
271987afe7SBarry Smith   PLogFlops(2*N-1);
281987afe7SBarry Smith   return 0;
291987afe7SBarry Smith }
301987afe7SBarry Smith 
31*c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
32289bc588SBarry Smith {
33*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
34289bc588SBarry Smith   int          i,N = mat->m*mat->n,count = 0;
35289bc588SBarry Smith   Scalar       *v = mat->v;
36289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
37*c0bbcb79SLois Curfman McInnes   *nz = count; *nzalloc = N; *mem = (int)A->mem;
38fa9ec3f1SLois Curfman McInnes   return 0;
39289bc588SBarry Smith }
40289bc588SBarry Smith 
41289bc588SBarry Smith /* ---------------------------------------------------------------*/
42289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
43289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
44*c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
45289bc588SBarry Smith {
46*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
476abc6512SBarry Smith   int          info;
48289bc588SBarry Smith   if (!mat->pivots) {
4948d91487SBarry Smith     mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
50*c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
51289bc588SBarry Smith   }
52289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
53ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
54*c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
55289bc588SBarry Smith   return 0;
56289bc588SBarry Smith }
57*c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
58289bc588SBarry Smith {
59289bc588SBarry Smith   int ierr;
60*c0bbcb79SLois Curfman McInnes   ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr);
61289bc588SBarry Smith   return 0;
62289bc588SBarry Smith }
63*c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
64289bc588SBarry Smith {
6549d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
66289bc588SBarry Smith }
67*c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
68289bc588SBarry Smith {
69289bc588SBarry Smith   int ierr;
70*c0bbcb79SLois Curfman McInnes   ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr);
71289bc588SBarry Smith   return 0;
72289bc588SBarry Smith }
73*c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
74289bc588SBarry Smith {
7549d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
76289bc588SBarry Smith }
77*c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
78289bc588SBarry Smith {
79*c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
806abc6512SBarry Smith   int           info;
81752f5567SLois Curfman McInnes   if (mat->pivots) {
82752f5567SLois Curfman McInnes     PETSCFREE(mat->pivots);
83*c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
84752f5567SLois Curfman McInnes     mat->pivots = 0;
85752f5567SLois Curfman McInnes   }
86289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
87ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
88*c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
89289bc588SBarry Smith   return 0;
90289bc588SBarry Smith }
91289bc588SBarry Smith 
92*c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
93289bc588SBarry Smith {
94*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
956abc6512SBarry Smith   int          one = 1, info;
966abc6512SBarry Smith   Scalar       *x, *y;
97289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
98416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
99*c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
10048d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
101289bc588SBarry Smith   }
102*c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
10348d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
104289bc588SBarry Smith   }
105ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
106ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
107289bc588SBarry Smith   return 0;
108289bc588SBarry Smith }
109*c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
110da3a660dSBarry Smith {
111*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1126abc6512SBarry Smith   int          one = 1, info;
1136abc6512SBarry Smith   Scalar       *x, *y;
114da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
115416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
116752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
117da3a660dSBarry Smith   if (mat->pivots) {
11848d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
119da3a660dSBarry Smith   }
120da3a660dSBarry Smith   else {
12148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
122da3a660dSBarry Smith   }
123ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
124da3a660dSBarry Smith   return 0;
125da3a660dSBarry Smith }
126*c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
127da3a660dSBarry Smith {
128*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1296abc6512SBarry Smith   int          one = 1, info,ierr;
1306abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
131da3a660dSBarry Smith   Vec          tmp = 0;
132da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
133da3a660dSBarry Smith   if (yy == zz) {
13478b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
135*c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
13678b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
137da3a660dSBarry Smith   }
138416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
139752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
140da3a660dSBarry Smith   if (mat->pivots) {
14148d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
142da3a660dSBarry Smith   }
143da3a660dSBarry Smith   else {
14448d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
145da3a660dSBarry Smith   }
146ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
147da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
148da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
149da3a660dSBarry Smith   return 0;
150da3a660dSBarry Smith }
151*c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
152da3a660dSBarry Smith {
153*c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1546abc6512SBarry Smith   int           one = 1, info,ierr;
1556abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
156da3a660dSBarry Smith   Vec           tmp;
157da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
158da3a660dSBarry Smith   if (yy == zz) {
15978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
160*c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
16178b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
162da3a660dSBarry Smith   }
163416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
164752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
165da3a660dSBarry Smith   if (mat->pivots) {
16648d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
167da3a660dSBarry Smith   }
168da3a660dSBarry Smith   else {
16948d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
170da3a660dSBarry Smith   }
171ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
172da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
173da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
174da3a660dSBarry Smith   return 0;
175da3a660dSBarry Smith }
176289bc588SBarry Smith /* ------------------------------------------------------------------*/
177*c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
17820e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
179289bc588SBarry Smith {
180*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
181289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
1826abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
183289bc588SBarry Smith 
184289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
185289bc588SBarry Smith     /* this is a hack fix, should have another version without
186289bc588SBarry Smith        the second BLdot */
187bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
188289bc588SBarry Smith   }
189289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
190289bc588SBarry Smith   while (its--) {
191289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
192289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
193f1747703SBarry Smith #if defined(PETSC_COMPLEX)
194f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
195f1747703SBarry Smith            not happy about returning a double complex */
196f1747703SBarry Smith         int    _i;
197f1747703SBarry Smith         Scalar sum = b[i];
198f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
199f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
200f1747703SBarry Smith         }
201f1747703SBarry Smith         xt = sum;
202f1747703SBarry Smith #else
203289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
204f1747703SBarry Smith #endif
205d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
206289bc588SBarry Smith       }
207289bc588SBarry Smith     }
208289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
209289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
210f1747703SBarry Smith #if defined(PETSC_COMPLEX)
211f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
212f1747703SBarry Smith            not happy about returning a double complex */
213f1747703SBarry Smith         int    _i;
214f1747703SBarry Smith         Scalar sum = b[i];
215f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
216f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
217f1747703SBarry Smith         }
218f1747703SBarry Smith         xt = sum;
219f1747703SBarry Smith #else
220289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
221f1747703SBarry Smith #endif
222d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
223289bc588SBarry Smith       }
224289bc588SBarry Smith     }
225289bc588SBarry Smith   }
226289bc588SBarry Smith   return 0;
227289bc588SBarry Smith }
228289bc588SBarry Smith 
229289bc588SBarry Smith /* -----------------------------------------------------------------*/
230*c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
231289bc588SBarry Smith {
232*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
233289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
234289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
235289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
23648d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
237289bc588SBarry Smith   return 0;
238289bc588SBarry Smith }
239*c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
240289bc588SBarry Smith {
241*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
242289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
243289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
244289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
24548d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
246289bc588SBarry Smith   return 0;
247289bc588SBarry Smith }
248*c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
249289bc588SBarry Smith {
250*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
251289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2526abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
253289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
254416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
25548d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
256289bc588SBarry Smith   return 0;
257289bc588SBarry Smith }
258*c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
259289bc588SBarry Smith {
260*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
261289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
262289bc588SBarry Smith   int          _One=1;
2636abc6512SBarry Smith   Scalar       _DOne=1.0;
264289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
265289bc588SBarry Smith   VecGetArray(zz,&z);
266416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
26748d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
268289bc588SBarry Smith   return 0;
269289bc588SBarry Smith }
270289bc588SBarry Smith 
271289bc588SBarry Smith /* -----------------------------------------------------------------*/
272*c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
273289bc588SBarry Smith {
274*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
275289bc588SBarry Smith   Scalar       *v;
276289bc588SBarry Smith   int          i;
277289bc588SBarry Smith   *ncols = mat->n;
278289bc588SBarry Smith   if (cols) {
27978b31e54SBarry Smith     *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols);
28073c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
281289bc588SBarry Smith   }
282289bc588SBarry Smith   if (vals) {
28378b31e54SBarry Smith     *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
284289bc588SBarry Smith     v = mat->v + row;
28573c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
286289bc588SBarry Smith   }
287289bc588SBarry Smith   return 0;
288289bc588SBarry Smith }
289*c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
290289bc588SBarry Smith {
29178b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
29278b31e54SBarry Smith   if (vals) { PETSCFREE(*vals); }
293289bc588SBarry Smith   return 0;
294289bc588SBarry Smith }
295289bc588SBarry Smith /* ----------------------------------------------------------------*/
296*c0bbcb79SLois Curfman McInnes static int MatInsert_SeqDense(Mat A,int m,int *indexm,int n,
297e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
298289bc588SBarry Smith {
299*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
300289bc588SBarry Smith   int          i,j;
301d6dfbf8fSBarry Smith 
302289bc588SBarry Smith   if (!mat->roworiented) {
303dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
304289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
305289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
306289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
307289bc588SBarry Smith         }
308289bc588SBarry Smith       }
309289bc588SBarry Smith     }
310289bc588SBarry Smith     else {
311289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
312289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
313289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
314289bc588SBarry Smith         }
315289bc588SBarry Smith       }
316289bc588SBarry Smith     }
317e8d4e0b9SBarry Smith   }
318e8d4e0b9SBarry Smith   else {
319dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
320e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
321e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
322e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
323e8d4e0b9SBarry Smith         }
324e8d4e0b9SBarry Smith       }
325e8d4e0b9SBarry Smith     }
326289bc588SBarry Smith     else {
327289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
328289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
329289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
330289bc588SBarry Smith         }
331289bc588SBarry Smith       }
332289bc588SBarry Smith     }
333e8d4e0b9SBarry Smith   }
334289bc588SBarry Smith   return 0;
335289bc588SBarry Smith }
336e8d4e0b9SBarry Smith 
337289bc588SBarry Smith /* -----------------------------------------------------------------*/
338*c0bbcb79SLois Curfman McInnes static int MatCopyPrivate_SeqDense(Mat A,Mat *newmat)
339289bc588SBarry Smith {
340*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data,*l;
341289bc588SBarry Smith   int          ierr;
342289bc588SBarry Smith   Mat          newi;
34348d91487SBarry Smith 
344*c0bbcb79SLois Curfman McInnes   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,&newi);CHKERRQ(ierr);
345ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
346416022c9SBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
347289bc588SBarry Smith   *newmat = newi;
348289bc588SBarry Smith   return 0;
349289bc588SBarry Smith }
350da3a660dSBarry Smith #include "viewer.h"
351289bc588SBarry Smith 
352ec8511deSBarry Smith int MatView_SeqDense(PetscObject obj,Viewer ptr)
353289bc588SBarry Smith {
354*c0bbcb79SLois Curfman McInnes   Mat           A = (Mat) obj;
355*c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
356289bc588SBarry Smith   Scalar        *v;
357d636dbe3SBarry Smith   int           i,j,ierr;
35823242f5aSBarry Smith   PetscObject   vobj = (PetscObject) ptr;
359da3a660dSBarry Smith 
36023242f5aSBarry Smith   if (ptr == 0) {
36121b0d8fbSLois Curfman McInnes     ptr = STDOUT_VIEWER_SELF; vobj = (PetscObject) ptr;
36223242f5aSBarry Smith   }
36323242f5aSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
3646f75cfa4SBarry Smith     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
365da3a660dSBarry Smith   }
366da3a660dSBarry Smith   else {
367d636dbe3SBarry Smith     FILE *fd;
368d636dbe3SBarry Smith     int format;
369d636dbe3SBarry Smith     ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr);
370d636dbe3SBarry Smith     ierr = ViewerFileGetFormat_Private(ptr,&format); CHKERRQ(ierr);
371f72585f2SLois Curfman McInnes     if (format == FILE_FORMAT_INFO) {
372f72585f2SLois Curfman McInnes       /* do nothing for now */
373f72585f2SLois Curfman McInnes     }
374f72585f2SLois Curfman McInnes     else {
375289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
376289bc588SBarry Smith         v = mat->v + i;
377289bc588SBarry Smith         for ( j=0; j<mat->n; j++ ) {
378289bc588SBarry Smith #if defined(PETSC_COMPLEX)
3798e37dc09SBarry Smith           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
380289bc588SBarry Smith #else
3818e37dc09SBarry Smith           fprintf(fd,"%6.4e ",*v); v += mat->m;
382289bc588SBarry Smith #endif
383289bc588SBarry Smith         }
3848e37dc09SBarry Smith         fprintf(fd,"\n");
385289bc588SBarry Smith       }
386da3a660dSBarry Smith     }
387f72585f2SLois Curfman McInnes   }
388289bc588SBarry Smith   return 0;
389289bc588SBarry Smith }
390289bc588SBarry Smith 
391289bc588SBarry Smith 
392ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
393289bc588SBarry Smith {
394289bc588SBarry Smith   Mat          mat = (Mat) obj;
395ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
396a5a9c739SBarry Smith #if defined(PETSC_LOG)
397a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
398a5a9c739SBarry Smith #endif
39978b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
40078b31e54SBarry Smith   PETSCFREE(l);
401a5a9c739SBarry Smith   PLogObjectDestroy(mat);
4029e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
403289bc588SBarry Smith   return 0;
404289bc588SBarry Smith }
405289bc588SBarry Smith 
406*c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
407289bc588SBarry Smith {
408*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
409d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
410d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
41148b35521SBarry Smith 
412d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
41348b35521SBarry Smith   if (!matout) { /* in place transpose */
41448d91487SBarry Smith     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Not for rectangular matrix in place");
415d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
416289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
417d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
418d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
419d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
420289bc588SBarry Smith       }
421289bc588SBarry Smith     }
42248b35521SBarry Smith   }
423d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
424d3e5ee88SLois Curfman McInnes     int ierr;
425d3e5ee88SLois Curfman McInnes     Mat tmat;
426ec8511deSBarry Smith     Mat_SeqDense *tmatd;
427d3e5ee88SLois Curfman McInnes     Scalar *v2;
428*c0bbcb79SLois Curfman McInnes     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr);
429ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
4300de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
431d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
4320de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
433d3e5ee88SLois Curfman McInnes     }
434d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
435d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
436d3e5ee88SLois Curfman McInnes     *matout = tmat;
43748b35521SBarry Smith   }
438289bc588SBarry Smith   return 0;
439289bc588SBarry Smith }
440289bc588SBarry Smith 
441*c0bbcb79SLois Curfman McInnes static int MatEqual_SeqDense(Mat A1,Mat A2)
442289bc588SBarry Smith {
443*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
444*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
445289bc588SBarry Smith   int          i;
446289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
447289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
448289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
449289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
450289bc588SBarry Smith     if (*v1 != *v2) return 0;
451289bc588SBarry Smith     v1++; v2++;
452289bc588SBarry Smith   }
453289bc588SBarry Smith   return 1;
454289bc588SBarry Smith }
455289bc588SBarry Smith 
456*c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
457289bc588SBarry Smith {
458*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
4596abc6512SBarry Smith   int          i, n;
4606abc6512SBarry Smith   Scalar       *x;
461289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
462ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
463289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
464289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
465289bc588SBarry Smith   }
466289bc588SBarry Smith   return 0;
467289bc588SBarry Smith }
468289bc588SBarry Smith 
469*c0bbcb79SLois Curfman McInnes static int MatScale_SeqDense(Mat A,Vec ll,Vec rr)
470289bc588SBarry Smith {
471*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
472da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
473da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
47428988994SBarry Smith   if (ll) {
475da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
476ec8511deSBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size");
477da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
478da3a660dSBarry Smith       x = l[i];
479da3a660dSBarry Smith       v = mat->v + i;
480da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
481da3a660dSBarry Smith     }
482da3a660dSBarry Smith   }
48328988994SBarry Smith   if (rr) {
484da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
485ec8511deSBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size");
486da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
487da3a660dSBarry Smith       x = r[i];
488da3a660dSBarry Smith       v = mat->v + i*m;
489da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
490da3a660dSBarry Smith     }
491da3a660dSBarry Smith   }
492289bc588SBarry Smith   return 0;
493289bc588SBarry Smith }
494289bc588SBarry Smith 
495*c0bbcb79SLois Curfman McInnes static int MatNorm_SeqDense(Mat A,MatNormType type,double *norm)
496289bc588SBarry Smith {
497*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
498289bc588SBarry Smith   Scalar       *v = mat->v;
499289bc588SBarry Smith   double       sum = 0.0;
500289bc588SBarry Smith   int          i, j;
501289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
502289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
503289bc588SBarry Smith #if defined(PETSC_COMPLEX)
504289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
505289bc588SBarry Smith #else
506289bc588SBarry Smith       sum += (*v)*(*v); v++;
507289bc588SBarry Smith #endif
508289bc588SBarry Smith     }
509289bc588SBarry Smith     *norm = sqrt(sum);
510289bc588SBarry Smith   }
511289bc588SBarry Smith   else if (type == NORM_1) {
512289bc588SBarry Smith     *norm = 0.0;
513289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
514289bc588SBarry Smith       sum = 0.0;
515289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
516289bc588SBarry Smith #if defined(PETSC_COMPLEX)
517289bc588SBarry Smith         sum += abs(*v++);
518289bc588SBarry Smith #else
519289bc588SBarry Smith         sum += fabs(*v++);
520289bc588SBarry Smith #endif
521289bc588SBarry Smith       }
522289bc588SBarry Smith       if (sum > *norm) *norm = sum;
523289bc588SBarry Smith     }
524289bc588SBarry Smith   }
525289bc588SBarry Smith   else if (type == NORM_INFINITY) {
526289bc588SBarry Smith     *norm = 0.0;
527289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
528289bc588SBarry Smith       v = mat->v + j;
529289bc588SBarry Smith       sum = 0.0;
530289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
531289bc588SBarry Smith #if defined(PETSC_COMPLEX)
532289bc588SBarry Smith         sum += abs(*v); v += mat->m;
533289bc588SBarry Smith #else
534289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
535289bc588SBarry Smith #endif
536289bc588SBarry Smith       }
537289bc588SBarry Smith       if (sum > *norm) *norm = sum;
538289bc588SBarry Smith     }
539289bc588SBarry Smith   }
540289bc588SBarry Smith   else {
54148d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
542289bc588SBarry Smith   }
543289bc588SBarry Smith   return 0;
544289bc588SBarry Smith }
545289bc588SBarry Smith 
546*c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
547289bc588SBarry Smith {
548*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
549289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
550289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
551*c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
552*c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
553*c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
554*c0bbcb79SLois Curfman McInnes            op == NO_NEW_NONZERO_LOCATIONS ||
555*c0bbcb79SLois Curfman McInnes            op == YES_NEW_NONZERO_LOCATIONS ||
556*c0bbcb79SLois Curfman McInnes            op == NO_NEW_DIAGONALS ||
557*c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
558*c0bbcb79SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n");
559*c0bbcb79SLois Curfman McInnes   else
560*c0bbcb79SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:Option not supported");}
561289bc588SBarry Smith   return 0;
562289bc588SBarry Smith }
563289bc588SBarry Smith 
564ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
5656f0a148fSBarry Smith {
566ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
567416022c9SBarry Smith   PetscZero(l->v,l->m*l->n*sizeof(Scalar));
5686f0a148fSBarry Smith   return 0;
5696f0a148fSBarry Smith }
5706f0a148fSBarry Smith 
571ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
5726f0a148fSBarry Smith {
573ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
5746abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
5756f0a148fSBarry Smith   Scalar       *slot;
57678b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
57778b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5786f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5796f0a148fSBarry Smith     slot = l->v + rows[i];
5806f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5816f0a148fSBarry Smith   }
5826f0a148fSBarry Smith   if (diag) {
5836f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5846f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
585260925b8SBarry Smith       *slot = *diag;
5866f0a148fSBarry Smith     }
5876f0a148fSBarry Smith   }
588260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5896f0a148fSBarry Smith   return 0;
5906f0a148fSBarry Smith }
591557bce09SLois Curfman McInnes 
592*c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
593557bce09SLois Curfman McInnes {
594*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
595557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
596557bce09SLois Curfman McInnes   return 0;
597557bce09SLois Curfman McInnes }
598557bce09SLois Curfman McInnes 
599*c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
600d3e5ee88SLois Curfman McInnes {
601*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
602d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
603d3e5ee88SLois Curfman McInnes   return 0;
604d3e5ee88SLois Curfman McInnes }
605d3e5ee88SLois Curfman McInnes 
606*c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
60764e87e97SBarry Smith {
608*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
60964e87e97SBarry Smith   *array = mat->v;
61064e87e97SBarry Smith   return 0;
61164e87e97SBarry Smith }
6120754003eSLois Curfman McInnes 
6130754003eSLois Curfman McInnes 
614*c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
6150754003eSLois Curfman McInnes {
616ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
6170754003eSLois Curfman McInnes }
6180754003eSLois Curfman McInnes 
619*c0bbcb79SLois Curfman McInnes static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,Mat *submat)
6200754003eSLois Curfman McInnes {
621*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
6220754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
623160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
6240754003eSLois Curfman McInnes   Scalar       *vwork, *val;
6250754003eSLois Curfman McInnes   Mat          newmat;
6260754003eSLois Curfman McInnes 
62778b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
62878b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
62978b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
63078b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
6310754003eSLois Curfman McInnes 
63278b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
63378b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
63478b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
635416022c9SBarry Smith   PetscZero((char*)smap,oldcols*sizeof(int));
6360754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
6370754003eSLois Curfman McInnes 
6380754003eSLois Curfman McInnes   /* Create and fill new matrix */
639*c0bbcb79SLois Curfman McInnes   ierr = MatCreateSeqDense(A->comm,nrows,ncols,&newmat);CHKERRQ(ierr);
6400754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
6410754003eSLois Curfman McInnes     nznew = 0;
6420754003eSLois Curfman McInnes     val   = mat->v + irow[i];
6430754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
6440754003eSLois Curfman McInnes       if (smap[j]) {
6450754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
6460754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
6470754003eSLois Curfman McInnes       }
6480754003eSLois Curfman McInnes     }
649dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
65078b31e54SBarry Smith            CHKERRQ(ierr);
6510754003eSLois Curfman McInnes   }
65278b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
65378b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
6540754003eSLois Curfman McInnes 
6550754003eSLois Curfman McInnes   /* Free work space */
65678b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
65778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
65878b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
6590754003eSLois Curfman McInnes   *submat = newmat;
6600754003eSLois Curfman McInnes   return 0;
6610754003eSLois Curfman McInnes }
6620754003eSLois Curfman McInnes 
663289bc588SBarry Smith /* -------------------------------------------------------------------*/
664ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense,
665ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
666ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
667ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
668ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
669ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
670ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
671ec8511deSBarry Smith        MatRelax_SeqDense,
672ec8511deSBarry Smith        MatTranspose_SeqDense,
673ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
674ec8511deSBarry Smith        MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense,
675289bc588SBarry Smith        0,0,
676ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
677ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
678ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
679ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
680ec8511deSBarry Smith        0,0,MatGetArray_SeqDense,0,0,
681ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
6821987afe7SBarry Smith        MatCopyPrivate_SeqDense,0,0,0,0,
6831987afe7SBarry Smith        MatAXPY_SeqDense};
6840754003eSLois Curfman McInnes 
6854b828684SBarry Smith /*@C
686fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
687d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
688d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
689289bc588SBarry Smith 
69020563c6bSBarry Smith    Input Parameters:
6910c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
6920c775827SLois Curfman McInnes .  m - number of rows
6930c775827SLois Curfman McInnes .  n - number of column
69420563c6bSBarry Smith 
69520563c6bSBarry Smith    Output Parameter:
6960c775827SLois Curfman McInnes .  newmat - the matrix
69720563c6bSBarry Smith 
698dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
699d65003e9SLois Curfman McInnes 
700dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
70120563c6bSBarry Smith @*/
702fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat)
703289bc588SBarry Smith {
704ec8511deSBarry Smith   int          size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar);
705289bc588SBarry Smith   Mat          mat;
706ec8511deSBarry Smith   Mat_SeqDense *l;
707289bc588SBarry Smith   *newmat        = 0;
708ec8511deSBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
709a5a9c739SBarry Smith   PLogObjectCreate(mat);
710ec8511deSBarry Smith   l              = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l);
711416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
712ec8511deSBarry Smith   mat->destroy   = MatDestroy_SeqDense;
713ec8511deSBarry Smith   mat->view      = MatView_SeqDense;
714289bc588SBarry Smith   mat->data      = (void *) l;
715289bc588SBarry Smith   mat->factor    = 0;
716752f5567SLois Curfman McInnes   PLogObjectMemory(mat,sizeof(struct _Mat) + size);
717d6dfbf8fSBarry Smith 
718289bc588SBarry Smith   l->m           = m;
719289bc588SBarry Smith   l->n           = n;
720289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
721289bc588SBarry Smith   l->pivots      = 0;
722289bc588SBarry Smith   l->roworiented = 1;
723d6dfbf8fSBarry Smith 
724416022c9SBarry Smith   PetscZero(l->v,m*n*sizeof(Scalar));
725289bc588SBarry Smith   *newmat = mat;
726289bc588SBarry Smith   return 0;
727289bc588SBarry Smith }
728289bc588SBarry Smith 
729*c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
730289bc588SBarry Smith {
731*c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
732*c0bbcb79SLois Curfman McInnes   return MatCreateSeqDense(A->comm,m->m,m->n,newmat);
733289bc588SBarry Smith }
734