xref: /petsc/src/mat/impls/dense/seq/dense.c (revision ec8511de82ec98b283fd006bf3680df2e1c57a57)
1cb512458SBarry Smith #ifndef lint
2*ec8511deSBarry Smith static char vcid[] = "$Id: dense.c,v 1.57 1995/09/11 18:47:37 bsmith Exp bsmith $";
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 {
16289bc588SBarry Smith   Scalar *v;
17289bc588SBarry Smith   int    roworiented;
18289bc588SBarry Smith   int    m,n,pad;
19289bc588SBarry Smith   int    *pivots;   /* pivots in LU factorization */
20*ec8511deSBarry Smith } Mat_SeqDense;
21289bc588SBarry Smith 
22*ec8511deSBarry Smith static int MatGetInfo_SeqDense(Mat matin,MatInfoType flag,int *nz,
2320e1a0d4SLois Curfman McInnes                             int *nzalloc,int *mem)
24289bc588SBarry Smith {
25*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
26289bc588SBarry Smith   int    i,N = mat->m*mat->n,count = 0;
27289bc588SBarry Smith   Scalar *v = mat->v;
28289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
29752f5567SLois Curfman McInnes   *nz = count; *nzalloc = N; *mem = (int)matin->mem;
30fa9ec3f1SLois Curfman McInnes   return 0;
31289bc588SBarry Smith }
32289bc588SBarry Smith 
33289bc588SBarry Smith /* ---------------------------------------------------------------*/
34289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
35289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
36*ec8511deSBarry Smith static int MatLUFactor_SeqDense(Mat matin,IS row,IS col,double f)
37289bc588SBarry Smith {
38*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
396abc6512SBarry Smith   int    info;
40289bc588SBarry Smith   if (!mat->pivots) {
4178b31e54SBarry Smith     mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int));
4278b31e54SBarry Smith     CHKPTRQ(mat->pivots);
43752f5567SLois Curfman McInnes     PLogObjectMemory(matin,mat->m*sizeof(int));
44289bc588SBarry Smith   }
45289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
46*ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense: Bad LU factorization");
47289bc588SBarry Smith   matin->factor = FACTOR_LU;
48289bc588SBarry Smith   return 0;
49289bc588SBarry Smith }
50*ec8511deSBarry Smith static int MatLUFactorSymbolic_SeqDense(Mat matin,IS row,IS col,double f,
5149d8b64dSBarry Smith                                      Mat *fact)
52289bc588SBarry Smith {
53289bc588SBarry Smith   int ierr;
54bbb6d6a8SBarry Smith   ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr);
55289bc588SBarry Smith   return 0;
56289bc588SBarry Smith }
57*ec8511deSBarry Smith static int MatLUFactorNumeric_SeqDense(Mat matin,Mat *fact)
58289bc588SBarry Smith {
5949d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
60289bc588SBarry Smith }
61*ec8511deSBarry Smith static int MatCholeskyFactorSymbolic_SeqDense(Mat matin,IS row,double f,Mat *fact)
62289bc588SBarry Smith {
63289bc588SBarry Smith   int ierr;
64bbb6d6a8SBarry Smith   ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr);
65289bc588SBarry Smith   return 0;
66289bc588SBarry Smith }
67*ec8511deSBarry Smith static int MatCholeskyFactorNumeric_SeqDense(Mat matin,Mat *fact)
68289bc588SBarry Smith {
6949d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
70289bc588SBarry Smith }
71*ec8511deSBarry Smith static int MatCholeskyFactor_SeqDense(Mat matin,IS perm,double f)
72289bc588SBarry Smith {
73*ec8511deSBarry Smith   Mat_SeqDense    *mat = (Mat_SeqDense *) matin->data;
746abc6512SBarry Smith   int       info;
75752f5567SLois Curfman McInnes   if (mat->pivots) {
76752f5567SLois Curfman McInnes     PETSCFREE(mat->pivots);
77752f5567SLois Curfman McInnes     PLogObjectMemory(matin,-mat->m*sizeof(int));
78752f5567SLois Curfman McInnes     mat->pivots = 0;
79752f5567SLois Curfman McInnes   }
80289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
81*ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense: Bad factorization");
82289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
83289bc588SBarry Smith   return 0;
84289bc588SBarry Smith }
85289bc588SBarry Smith 
86*ec8511deSBarry Smith static int MatSolve_SeqDense(Mat matin,Vec xx,Vec yy)
87289bc588SBarry Smith {
88*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
896abc6512SBarry Smith   int    one = 1, info;
906abc6512SBarry Smith   Scalar *x, *y;
91289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
9278b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
939e25ed09SBarry Smith   if (matin->factor == FACTOR_LU) {
94289bc588SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
95289bc588SBarry Smith               y, &mat->m, &info );
96289bc588SBarry Smith   }
979e25ed09SBarry Smith   else if (matin->factor == FACTOR_CHOLESKY){
98289bc588SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
99289bc588SBarry Smith               y, &mat->m, &info );
100289bc588SBarry Smith   }
101*ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense: Matrix must be factored to solve");
102*ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense: Bad solve");
103289bc588SBarry Smith   return 0;
104289bc588SBarry Smith }
105*ec8511deSBarry Smith static int MatSolveTrans_SeqDense(Mat matin,Vec xx,Vec yy)
106da3a660dSBarry Smith {
107*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
1086abc6512SBarry Smith   int    one = 1, info;
1096abc6512SBarry Smith   Scalar *x, *y;
110da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
11178b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
112752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
113da3a660dSBarry Smith   if (mat->pivots) {
114da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
115da3a660dSBarry Smith               y, &mat->m, &info );
116da3a660dSBarry Smith   }
117da3a660dSBarry Smith   else {
118da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
119da3a660dSBarry Smith               y, &mat->m, &info );
120da3a660dSBarry Smith   }
121*ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense: Bad solve");
122da3a660dSBarry Smith   return 0;
123da3a660dSBarry Smith }
124*ec8511deSBarry Smith static int MatSolveAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy)
125da3a660dSBarry Smith {
126*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
1276abc6512SBarry Smith   int    one = 1, info,ierr;
1286abc6512SBarry Smith   Scalar *x, *y, sone = 1.0;
129da3a660dSBarry Smith   Vec    tmp = 0;
130da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
131da3a660dSBarry Smith   if (yy == zz) {
13278b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
133464493b3SBarry Smith     PLogObjectParent(matin,tmp);
13478b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
135da3a660dSBarry Smith   }
13678b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
137752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
138da3a660dSBarry Smith   if (mat->pivots) {
139da3a660dSBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
140da3a660dSBarry Smith               y, &mat->m, &info );
141da3a660dSBarry Smith   }
142da3a660dSBarry Smith   else {
143da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
144da3a660dSBarry Smith               y, &mat->m, &info );
145da3a660dSBarry Smith   }
146*ec8511deSBarry 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*ec8511deSBarry Smith static int MatSolveTransAdd_SeqDense(Mat matin,Vec xx,Vec zz, Vec yy)
152da3a660dSBarry Smith {
153*ec8511deSBarry Smith   Mat_SeqDense  *mat = (Mat_SeqDense *) matin->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);
160464493b3SBarry Smith     PLogObjectParent(matin,tmp);
16178b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
162da3a660dSBarry Smith   }
16378b31e54SBarry 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) {
166da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
167da3a660dSBarry Smith               y, &mat->m, &info );
168da3a660dSBarry Smith   }
169da3a660dSBarry Smith   else {
170da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
171da3a660dSBarry Smith               y, &mat->m, &info );
172da3a660dSBarry Smith   }
173*ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense: Bad solve");
174da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
175da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
176da3a660dSBarry Smith   return 0;
177da3a660dSBarry Smith }
178289bc588SBarry Smith /* ------------------------------------------------------------------*/
179*ec8511deSBarry Smith static int MatRelax_SeqDense(Mat matin,Vec bb,double omega,MatSORType flag,
18020e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
181289bc588SBarry Smith {
182*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
183289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
1846abc6512SBarry Smith   int    o = 1,ierr, m = mat->m, i;
185289bc588SBarry Smith 
186289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
187289bc588SBarry Smith     /* this is a hack fix, should have another version without
188289bc588SBarry Smith        the second BLdot */
189bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
190289bc588SBarry Smith   }
191289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
192289bc588SBarry Smith   while (its--) {
193289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
194289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
195f1747703SBarry Smith #if defined(PETSC_COMPLEX)
196f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
197f1747703SBarry Smith            not happy about returning a double complex */
198f1747703SBarry Smith         int    _i;
199f1747703SBarry Smith         Scalar sum = b[i];
200f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
201f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
202f1747703SBarry Smith         }
203f1747703SBarry Smith         xt = sum;
204f1747703SBarry Smith #else
205289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
206f1747703SBarry Smith #endif
207d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
208289bc588SBarry Smith       }
209289bc588SBarry Smith     }
210289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
211289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
212f1747703SBarry Smith #if defined(PETSC_COMPLEX)
213f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
214f1747703SBarry Smith            not happy about returning a double complex */
215f1747703SBarry Smith         int    _i;
216f1747703SBarry Smith         Scalar sum = b[i];
217f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
218f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
219f1747703SBarry Smith         }
220f1747703SBarry Smith         xt = sum;
221f1747703SBarry Smith #else
222289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
223f1747703SBarry Smith #endif
224d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
225289bc588SBarry Smith       }
226289bc588SBarry Smith     }
227289bc588SBarry Smith   }
228289bc588SBarry Smith   return 0;
229289bc588SBarry Smith }
230289bc588SBarry Smith 
231289bc588SBarry Smith /* -----------------------------------------------------------------*/
232*ec8511deSBarry Smith static int MatMultTrans_SeqDense(Mat matin,Vec xx,Vec yy)
233289bc588SBarry Smith {
234*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
235289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
236289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
237289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
238289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
239289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
240289bc588SBarry Smith   return 0;
241289bc588SBarry Smith }
242*ec8511deSBarry Smith static int MatMult_SeqDense(Mat matin,Vec xx,Vec yy)
243289bc588SBarry Smith {
244*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
245289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
246289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
247289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
248289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
249289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
250289bc588SBarry Smith   return 0;
251289bc588SBarry Smith }
252*ec8511deSBarry Smith static int MatMultAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy)
253289bc588SBarry Smith {
254*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
255289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
2566abc6512SBarry Smith   int    _One=1; Scalar _DOne=1.0;
257289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
25878b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
259289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
260289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
261289bc588SBarry Smith   return 0;
262289bc588SBarry Smith }
263*ec8511deSBarry Smith static int MatMultTransAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy)
264289bc588SBarry Smith {
265*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
266289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
267289bc588SBarry Smith   int    _One=1;
2686abc6512SBarry Smith   Scalar _DOne=1.0;
269289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
270289bc588SBarry Smith   VecGetArray(zz,&z);
27178b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
272289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
273289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
274289bc588SBarry Smith   return 0;
275289bc588SBarry Smith }
276289bc588SBarry Smith 
277289bc588SBarry Smith /* -----------------------------------------------------------------*/
278*ec8511deSBarry Smith static int MatGetRow_SeqDense(Mat matin,int row,int *ncols,int **cols,
279289bc588SBarry Smith                         Scalar **vals)
280289bc588SBarry Smith {
281*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
282289bc588SBarry Smith   Scalar *v;
283289bc588SBarry Smith   int    i;
284289bc588SBarry Smith   *ncols = mat->n;
285289bc588SBarry Smith   if (cols) {
28678b31e54SBarry Smith     *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols);
287289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
288289bc588SBarry Smith   }
289289bc588SBarry Smith   if (vals) {
29078b31e54SBarry Smith     *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
291289bc588SBarry Smith     v = mat->v + row;
292289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
293289bc588SBarry Smith   }
294289bc588SBarry Smith   return 0;
295289bc588SBarry Smith }
296*ec8511deSBarry Smith static int MatRestoreRow_SeqDense(Mat matin,int row,int *ncols,int **cols,
297289bc588SBarry Smith                             Scalar **vals)
298289bc588SBarry Smith {
29978b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
30078b31e54SBarry Smith   if (vals) { PETSCFREE(*vals); }
301289bc588SBarry Smith   return 0;
302289bc588SBarry Smith }
303289bc588SBarry Smith /* ----------------------------------------------------------------*/
304*ec8511deSBarry Smith static int MatInsert_SeqDense(Mat matin,int m,int *indexm,int n,
305e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
306289bc588SBarry Smith {
307*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
308289bc588SBarry Smith   int    i,j;
309d6dfbf8fSBarry Smith 
310289bc588SBarry Smith   if (!mat->roworiented) {
311edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
312289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
313289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
314289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
315289bc588SBarry Smith         }
316289bc588SBarry Smith       }
317289bc588SBarry Smith     }
318289bc588SBarry Smith     else {
319289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
320289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
321289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
322289bc588SBarry Smith         }
323289bc588SBarry Smith       }
324289bc588SBarry Smith     }
325e8d4e0b9SBarry Smith   }
326e8d4e0b9SBarry Smith   else {
327edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
328e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
329e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
330e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
331e8d4e0b9SBarry Smith         }
332e8d4e0b9SBarry Smith       }
333e8d4e0b9SBarry Smith     }
334289bc588SBarry Smith     else {
335289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
336289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
337289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
338289bc588SBarry Smith         }
339289bc588SBarry Smith       }
340289bc588SBarry Smith     }
341e8d4e0b9SBarry Smith   }
342289bc588SBarry Smith   return 0;
343289bc588SBarry Smith }
344e8d4e0b9SBarry Smith 
345289bc588SBarry Smith /* -----------------------------------------------------------------*/
346*ec8511deSBarry Smith static int MatCopyPrivate_SeqDense(Mat matin,Mat *newmat)
347289bc588SBarry Smith {
348*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
349289bc588SBarry Smith   int ierr;
350289bc588SBarry Smith   Mat newi;
351*ec8511deSBarry Smith   Mat_SeqDense *l;
352fafbff53SBarry Smith   ierr = MatCreateSeqDense(matin->comm,mat->m,mat->n,&newi);
353bbb6d6a8SBarry Smith   CHKERRQ(ierr);
354*ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
35578b31e54SBarry Smith   PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
356289bc588SBarry Smith   *newmat = newi;
357289bc588SBarry Smith   return 0;
358289bc588SBarry Smith }
359da3a660dSBarry Smith #include "viewer.h"
360289bc588SBarry Smith 
361*ec8511deSBarry Smith int MatView_SeqDense(PetscObject obj,Viewer ptr)
362289bc588SBarry Smith {
363289bc588SBarry Smith   Mat         matin = (Mat) obj;
364*ec8511deSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense *) matin->data;
365289bc588SBarry Smith   Scalar      *v;
366d636dbe3SBarry Smith   int         i,j,ierr;
36723242f5aSBarry Smith   PetscObject vobj = (PetscObject) ptr;
368da3a660dSBarry Smith 
36923242f5aSBarry Smith   if (ptr == 0) {
37021b0d8fbSLois Curfman McInnes     ptr = STDOUT_VIEWER_SELF; vobj = (PetscObject) ptr;
37123242f5aSBarry Smith   }
37223242f5aSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
3736f75cfa4SBarry Smith     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
374da3a660dSBarry Smith   }
375da3a660dSBarry Smith   else {
376d636dbe3SBarry Smith     FILE *fd;
377d636dbe3SBarry Smith     int format;
378d636dbe3SBarry Smith     ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr);
379d636dbe3SBarry Smith     ierr = ViewerFileGetFormat_Private(ptr,&format); CHKERRQ(ierr);
380f72585f2SLois Curfman McInnes     if (format == FILE_FORMAT_INFO) {
381f72585f2SLois Curfman McInnes       /* do nothing for now */
382f72585f2SLois Curfman McInnes     }
383f72585f2SLois Curfman McInnes     else {
384289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
385289bc588SBarry Smith         v = mat->v + i;
386289bc588SBarry Smith         for ( j=0; j<mat->n; j++ ) {
387289bc588SBarry Smith #if defined(PETSC_COMPLEX)
3888e37dc09SBarry Smith           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
389289bc588SBarry Smith #else
3908e37dc09SBarry Smith           fprintf(fd,"%6.4e ",*v); v += mat->m;
391289bc588SBarry Smith #endif
392289bc588SBarry Smith         }
3938e37dc09SBarry Smith         fprintf(fd,"\n");
394289bc588SBarry Smith       }
395da3a660dSBarry Smith     }
396f72585f2SLois Curfman McInnes   }
397289bc588SBarry Smith   return 0;
398289bc588SBarry Smith }
399289bc588SBarry Smith 
400289bc588SBarry Smith 
401*ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
402289bc588SBarry Smith {
403289bc588SBarry Smith   Mat    mat = (Mat) obj;
404*ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
405a5a9c739SBarry Smith #if defined(PETSC_LOG)
406a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
407a5a9c739SBarry Smith #endif
40878b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
40978b31e54SBarry Smith   PETSCFREE(l);
410a5a9c739SBarry Smith   PLogObjectDestroy(mat);
4119e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
412289bc588SBarry Smith   return 0;
413289bc588SBarry Smith }
414289bc588SBarry Smith 
415*ec8511deSBarry Smith static int MatTranspose_SeqDense(Mat matin,Mat *matout)
416289bc588SBarry Smith {
417*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
418d3e5ee88SLois Curfman McInnes   int       k, j, m, n;
419d3e5ee88SLois Curfman McInnes   Scalar    *v, tmp;
42048b35521SBarry Smith 
421d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
42248b35521SBarry Smith   if (!matout) { /* in place transpose */
423d3e5ee88SLois Curfman McInnes     if (m != n) SETERRQ(1,
424*ec8511deSBarry Smith       "MatTranspose_SeqDense: Cannot transpose rectangular matrix in place");
425d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
426289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
427d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
428d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
429d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
430289bc588SBarry Smith       }
431289bc588SBarry Smith     }
43248b35521SBarry Smith   }
433d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
434d3e5ee88SLois Curfman McInnes     int ierr;
435d3e5ee88SLois Curfman McInnes     Mat tmat;
436*ec8511deSBarry Smith     Mat_SeqDense *tmatd;
437d3e5ee88SLois Curfman McInnes     Scalar *v2;
438fafbff53SBarry Smith     ierr = MatCreateSeqDense(matin->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr);
439*ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
4400de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
441d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
4420de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
443d3e5ee88SLois Curfman McInnes     }
444d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
445d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
446d3e5ee88SLois Curfman McInnes     *matout = tmat;
44748b35521SBarry Smith   }
448289bc588SBarry Smith   return 0;
449289bc588SBarry Smith }
450289bc588SBarry Smith 
451*ec8511deSBarry Smith static int MatEqual_SeqDense(Mat matin1,Mat matin2)
452289bc588SBarry Smith {
453*ec8511deSBarry Smith   Mat_SeqDense *mat1 = (Mat_SeqDense *) matin1->data;
454*ec8511deSBarry Smith   Mat_SeqDense *mat2 = (Mat_SeqDense *) matin2->data;
455289bc588SBarry Smith   int    i;
456289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
457289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
458289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
459289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
460289bc588SBarry Smith     if (*v1 != *v2) return 0;
461289bc588SBarry Smith     v1++; v2++;
462289bc588SBarry Smith   }
463289bc588SBarry Smith   return 1;
464289bc588SBarry Smith }
465289bc588SBarry Smith 
466*ec8511deSBarry Smith static int MatGetDiagonal_SeqDense(Mat matin,Vec v)
467289bc588SBarry Smith {
468*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
4696abc6512SBarry Smith   int       i, n;
4706abc6512SBarry Smith   Scalar    *x;
471289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
472*ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
473289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
474289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
475289bc588SBarry Smith   }
476289bc588SBarry Smith   return 0;
477289bc588SBarry Smith }
478289bc588SBarry Smith 
479*ec8511deSBarry Smith static int MatScale_SeqDense(Mat matin,Vec ll,Vec rr)
480289bc588SBarry Smith {
481*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
482da3a660dSBarry Smith   Scalar *l,*r,x,*v;
483da3a660dSBarry Smith   int    i,j,m = mat->m, n = mat->n;
48428988994SBarry Smith   if (ll) {
485da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
486*ec8511deSBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size");
487da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
488da3a660dSBarry Smith       x = l[i];
489da3a660dSBarry Smith       v = mat->v + i;
490da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
491da3a660dSBarry Smith     }
492da3a660dSBarry Smith   }
49328988994SBarry Smith   if (rr) {
494da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
495*ec8511deSBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size");
496da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
497da3a660dSBarry Smith       x = r[i];
498da3a660dSBarry Smith       v = mat->v + i*m;
499da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
500da3a660dSBarry Smith     }
501da3a660dSBarry Smith   }
502289bc588SBarry Smith   return 0;
503289bc588SBarry Smith }
504289bc588SBarry Smith 
505*ec8511deSBarry Smith static int MatNorm_SeqDense(Mat matin,MatNormType type,double *norm)
506289bc588SBarry Smith {
507*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
508289bc588SBarry Smith   Scalar *v = mat->v;
509289bc588SBarry Smith   double sum = 0.0;
510289bc588SBarry Smith   int    i, j;
511289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
512289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
513289bc588SBarry Smith #if defined(PETSC_COMPLEX)
514289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
515289bc588SBarry Smith #else
516289bc588SBarry Smith       sum += (*v)*(*v); v++;
517289bc588SBarry Smith #endif
518289bc588SBarry Smith     }
519289bc588SBarry Smith     *norm = sqrt(sum);
520289bc588SBarry Smith   }
521289bc588SBarry Smith   else if (type == NORM_1) {
522289bc588SBarry Smith     *norm = 0.0;
523289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
524289bc588SBarry Smith       sum = 0.0;
525289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
526289bc588SBarry Smith #if defined(PETSC_COMPLEX)
527289bc588SBarry Smith         sum += abs(*v++);
528289bc588SBarry Smith #else
529289bc588SBarry Smith         sum += fabs(*v++);
530289bc588SBarry Smith #endif
531289bc588SBarry Smith       }
532289bc588SBarry Smith       if (sum > *norm) *norm = sum;
533289bc588SBarry Smith     }
534289bc588SBarry Smith   }
535289bc588SBarry Smith   else if (type == NORM_INFINITY) {
536289bc588SBarry Smith     *norm = 0.0;
537289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
538289bc588SBarry Smith       v = mat->v + j;
539289bc588SBarry Smith       sum = 0.0;
540289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
541289bc588SBarry Smith #if defined(PETSC_COMPLEX)
542289bc588SBarry Smith         sum += abs(*v); v += mat->m;
543289bc588SBarry Smith #else
544289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
545289bc588SBarry Smith #endif
546289bc588SBarry Smith       }
547289bc588SBarry Smith       if (sum > *norm) *norm = sum;
548289bc588SBarry Smith     }
549289bc588SBarry Smith   }
550289bc588SBarry Smith   else {
551*ec8511deSBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm yet");
552289bc588SBarry Smith   }
553289bc588SBarry Smith   return 0;
554289bc588SBarry Smith }
555289bc588SBarry Smith 
556*ec8511deSBarry Smith static int MatSetOption_SeqDense(Mat aijin,MatOption op)
557289bc588SBarry Smith {
558*ec8511deSBarry Smith   Mat_SeqDense *aij = (Mat_SeqDense *) aijin->data;
559289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
560289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
561289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
562289bc588SBarry Smith   return 0;
563289bc588SBarry Smith }
564289bc588SBarry Smith 
565*ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
5666f0a148fSBarry Smith {
567*ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
56878b31e54SBarry Smith   PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
5696f0a148fSBarry Smith   return 0;
5706f0a148fSBarry Smith }
5716f0a148fSBarry Smith 
572*ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
5736f0a148fSBarry Smith {
574*ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
5756abc6512SBarry Smith   int    n = l->n, i, j,ierr,N, *rows;
5766f0a148fSBarry Smith   Scalar *slot;
57778b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
57878b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5796f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5806f0a148fSBarry Smith     slot = l->v + rows[i];
5816f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5826f0a148fSBarry Smith   }
5836f0a148fSBarry Smith   if (diag) {
5846f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5856f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
586260925b8SBarry Smith       *slot = *diag;
5876f0a148fSBarry Smith     }
5886f0a148fSBarry Smith   }
589260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5906f0a148fSBarry Smith   return 0;
5916f0a148fSBarry Smith }
592557bce09SLois Curfman McInnes 
593*ec8511deSBarry Smith static int MatGetSize_SeqDense(Mat matin,int *m,int *n)
594557bce09SLois Curfman McInnes {
595*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
596557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
597557bce09SLois Curfman McInnes   return 0;
598557bce09SLois Curfman McInnes }
599557bce09SLois Curfman McInnes 
600*ec8511deSBarry Smith static int MatGetOwnershipRange_SeqDense(Mat matin,int *m,int *n)
601d3e5ee88SLois Curfman McInnes {
602*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
603d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
604d3e5ee88SLois Curfman McInnes   return 0;
605d3e5ee88SLois Curfman McInnes }
606d3e5ee88SLois Curfman McInnes 
607*ec8511deSBarry Smith static int MatGetArray_SeqDense(Mat matin,Scalar **array)
60864e87e97SBarry Smith {
609*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
61064e87e97SBarry Smith   *array = mat->v;
61164e87e97SBarry Smith   return 0;
61264e87e97SBarry Smith }
6130754003eSLois Curfman McInnes 
6140754003eSLois Curfman McInnes 
615*ec8511deSBarry Smith static int MatGetSubMatrixInPlace_SeqDense(Mat matin,IS isrow,IS iscol)
6160754003eSLois Curfman McInnes {
617*ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense: not done");
6180754003eSLois Curfman McInnes }
6190754003eSLois Curfman McInnes 
620*ec8511deSBarry Smith static int MatGetSubMatrix_SeqDense(Mat matin,IS isrow,IS iscol,Mat *submat)
6210754003eSLois Curfman McInnes {
622*ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
6230754003eSLois Curfman McInnes   int     nznew, *smap, i, j, ierr, oldcols = mat->n;
624160018dcSBarry Smith   int     *irow, *icol, nrows, ncols, *cwork;
6250754003eSLois Curfman McInnes   Scalar  *vwork, *val;
6260754003eSLois Curfman McInnes   Mat     newmat;
6270754003eSLois Curfman McInnes 
62878b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
62978b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
63078b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
63178b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
6320754003eSLois Curfman McInnes 
63378b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
63478b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
63578b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
63678b31e54SBarry Smith   PETSCMEMSET((char*)smap,0,oldcols*sizeof(int));
6370754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
6380754003eSLois Curfman McInnes 
6390754003eSLois Curfman McInnes   /* Create and fill new matrix */
640fafbff53SBarry Smith   ierr = MatCreateSeqDense(matin->comm,nrows,ncols,&newmat);
64178b31e54SBarry Smith          CHKERRQ(ierr);
6420754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
6430754003eSLois Curfman McInnes     nznew = 0;
6440754003eSLois Curfman McInnes     val   = mat->v + irow[i];
6450754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
6460754003eSLois Curfman McInnes       if (smap[j]) {
6470754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
6480754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
6490754003eSLois Curfman McInnes       }
6500754003eSLois Curfman McInnes     }
651edae2e7dSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES);
65278b31e54SBarry Smith            CHKERRQ(ierr);
6530754003eSLois Curfman McInnes   }
65478b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
65578b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
6560754003eSLois Curfman McInnes 
6570754003eSLois Curfman McInnes   /* Free work space */
65878b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
65978b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
66078b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
6610754003eSLois Curfman McInnes   *submat = newmat;
6620754003eSLois Curfman McInnes   return 0;
6630754003eSLois Curfman McInnes }
6640754003eSLois Curfman McInnes 
665289bc588SBarry Smith /* -------------------------------------------------------------------*/
666*ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense,
667*ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
668*ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
669*ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
670*ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
671*ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
672*ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
673*ec8511deSBarry Smith        MatRelax_SeqDense,
674*ec8511deSBarry Smith        MatTranspose_SeqDense,
675*ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
676*ec8511deSBarry Smith        MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense,
677289bc588SBarry Smith        0,0,
678*ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
679*ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
680*ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
681*ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
682*ec8511deSBarry Smith        0,0,MatGetArray_SeqDense,0,0,
683*ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
684*ec8511deSBarry Smith        MatCopyPrivate_SeqDense};
6850754003eSLois Curfman McInnes 
6864b828684SBarry Smith /*@C
687fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
688d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
689d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
690289bc588SBarry Smith 
69120563c6bSBarry Smith    Input Parameters:
6920c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
6930c775827SLois Curfman McInnes .  m - number of rows
6940c775827SLois Curfman McInnes .  n - number of column
69520563c6bSBarry Smith 
69620563c6bSBarry Smith    Output Parameter:
6970c775827SLois Curfman McInnes .  newmat - the matrix
69820563c6bSBarry Smith 
699dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
700d65003e9SLois Curfman McInnes 
701dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
70220563c6bSBarry Smith @*/
703fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat)
704289bc588SBarry Smith {
705*ec8511deSBarry Smith   int       size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar);
706289bc588SBarry Smith   Mat mat;
707*ec8511deSBarry Smith   Mat_SeqDense    *l;
708289bc588SBarry Smith   *newmat        = 0;
709*ec8511deSBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
710a5a9c739SBarry Smith   PLogObjectCreate(mat);
711*ec8511deSBarry Smith   l              = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l);
712289bc588SBarry Smith   mat->ops       = &MatOps;
713*ec8511deSBarry Smith   mat->destroy   = MatDestroy_SeqDense;
714*ec8511deSBarry Smith   mat->view      = MatView_SeqDense;
715289bc588SBarry Smith   mat->data      = (void *) l;
716289bc588SBarry Smith   mat->factor    = 0;
717752f5567SLois Curfman McInnes   PLogObjectMemory(mat,sizeof(struct _Mat) + size);
718d6dfbf8fSBarry Smith 
719289bc588SBarry Smith   l->m           = m;
720289bc588SBarry Smith   l->n           = n;
721289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
722289bc588SBarry Smith   l->pivots      = 0;
723289bc588SBarry Smith   l->roworiented = 1;
724d6dfbf8fSBarry Smith 
72578b31e54SBarry Smith   PETSCMEMSET(l->v,0,m*n*sizeof(Scalar));
726289bc588SBarry Smith   *newmat = mat;
727289bc588SBarry Smith   return 0;
728289bc588SBarry Smith }
729289bc588SBarry Smith 
730*ec8511deSBarry Smith int MatCreate_SeqDense(Mat matin,Mat *newmat)
731289bc588SBarry Smith {
732*ec8511deSBarry Smith   Mat_SeqDense *m = (Mat_SeqDense *) matin->data;
733fafbff53SBarry Smith   return MatCreateSeqDense(matin->comm,m->m,m->n,newmat);
734289bc588SBarry Smith }
735