xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 1987afe70420aa833a01760c1264c6071036a247)
1cb512458SBarry Smith #ifndef lint
2*1987afe7SBarry Smith static char vcid[] = "$Id: dense.c,v 1.62 1995/10/10 16:14:22 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 */
20ec8511deSBarry Smith } Mat_SeqDense;
21289bc588SBarry Smith 
22*1987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
23*1987afe7SBarry Smith {
24*1987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
25*1987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
26*1987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
27*1987afe7SBarry Smith   PLogFlops(2*N-1);
28*1987afe7SBarry Smith   return 0;
29*1987afe7SBarry Smith }
30*1987afe7SBarry Smith 
31*1987afe7SBarry Smith 
3248d91487SBarry Smith static int MatGetInfo_SeqDense(Mat matin,MatInfoType flag,int *nz,int *nzalloc,int *mem)
33289bc588SBarry Smith {
34ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
35289bc588SBarry Smith   int          i,N = mat->m*mat->n,count = 0;
36289bc588SBarry Smith   Scalar       *v = mat->v;
37289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
38752f5567SLois Curfman McInnes   *nz = count; *nzalloc = N; *mem = (int)matin->mem;
39fa9ec3f1SLois Curfman McInnes   return 0;
40289bc588SBarry Smith }
41289bc588SBarry Smith 
42289bc588SBarry Smith /* ---------------------------------------------------------------*/
43289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
44289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
45ec8511deSBarry Smith static int MatLUFactor_SeqDense(Mat matin,IS row,IS col,double f)
46289bc588SBarry Smith {
47ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
486abc6512SBarry Smith   int          info;
49289bc588SBarry Smith   if (!mat->pivots) {
5048d91487SBarry Smith     mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
51752f5567SLois Curfman McInnes     PLogObjectMemory(matin,mat->m*sizeof(int));
52289bc588SBarry Smith   }
53289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
54ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
55289bc588SBarry Smith   matin->factor = FACTOR_LU;
56289bc588SBarry Smith   return 0;
57289bc588SBarry Smith }
5848d91487SBarry Smith static int MatLUFactorSymbolic_SeqDense(Mat matin,IS row,IS col,double f,Mat *fact)
59289bc588SBarry Smith {
60289bc588SBarry Smith   int ierr;
61bbb6d6a8SBarry Smith   ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr);
62289bc588SBarry Smith   return 0;
63289bc588SBarry Smith }
64ec8511deSBarry Smith static int MatLUFactorNumeric_SeqDense(Mat matin,Mat *fact)
65289bc588SBarry Smith {
6649d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
67289bc588SBarry Smith }
68ec8511deSBarry Smith static int MatCholeskyFactorSymbolic_SeqDense(Mat matin,IS row,double f,Mat *fact)
69289bc588SBarry Smith {
70289bc588SBarry Smith   int ierr;
71bbb6d6a8SBarry Smith   ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr);
72289bc588SBarry Smith   return 0;
73289bc588SBarry Smith }
74ec8511deSBarry Smith static int MatCholeskyFactorNumeric_SeqDense(Mat matin,Mat *fact)
75289bc588SBarry Smith {
7649d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
77289bc588SBarry Smith }
78ec8511deSBarry Smith static int MatCholeskyFactor_SeqDense(Mat matin,IS perm,double f)
79289bc588SBarry Smith {
80ec8511deSBarry Smith   Mat_SeqDense  *mat = (Mat_SeqDense *) matin->data;
816abc6512SBarry Smith   int           info;
82752f5567SLois Curfman McInnes   if (mat->pivots) {
83752f5567SLois Curfman McInnes     PETSCFREE(mat->pivots);
84752f5567SLois Curfman McInnes     PLogObjectMemory(matin,-mat->m*sizeof(int));
85752f5567SLois Curfman McInnes     mat->pivots = 0;
86752f5567SLois Curfman McInnes   }
87289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
88ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
89289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
90289bc588SBarry Smith   return 0;
91289bc588SBarry Smith }
92289bc588SBarry Smith 
93ec8511deSBarry Smith static int MatSolve_SeqDense(Mat matin,Vec xx,Vec yy)
94289bc588SBarry Smith {
95ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
966abc6512SBarry Smith   int          one = 1, info;
976abc6512SBarry Smith   Scalar       *x, *y;
98289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
99416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
1009e25ed09SBarry Smith   if (matin->factor == FACTOR_LU) {
10148d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
102289bc588SBarry Smith   }
1039e25ed09SBarry Smith   else if (matin->factor == FACTOR_CHOLESKY){
10448d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
105289bc588SBarry Smith   }
106ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
107ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
108289bc588SBarry Smith   return 0;
109289bc588SBarry Smith }
110ec8511deSBarry Smith static int MatSolveTrans_SeqDense(Mat matin,Vec xx,Vec yy)
111da3a660dSBarry Smith {
112ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
1136abc6512SBarry Smith   int          one = 1, info;
1146abc6512SBarry Smith   Scalar       *x, *y;
115da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
116416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
117752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
118da3a660dSBarry Smith   if (mat->pivots) {
11948d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
120da3a660dSBarry Smith   }
121da3a660dSBarry Smith   else {
12248d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
123da3a660dSBarry Smith   }
124ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
125da3a660dSBarry Smith   return 0;
126da3a660dSBarry Smith }
127ec8511deSBarry Smith static int MatSolveAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy)
128da3a660dSBarry Smith {
129ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
1306abc6512SBarry Smith   int          one = 1, info,ierr;
1316abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
132da3a660dSBarry Smith   Vec          tmp = 0;
133da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
134da3a660dSBarry Smith   if (yy == zz) {
13578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
136464493b3SBarry Smith     PLogObjectParent(matin,tmp);
13778b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
138da3a660dSBarry Smith   }
139416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
140752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
141da3a660dSBarry Smith   if (mat->pivots) {
14248d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
143da3a660dSBarry Smith   }
144da3a660dSBarry Smith   else {
14548d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
146da3a660dSBarry Smith   }
147ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
148da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
149da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
150da3a660dSBarry Smith   return 0;
151da3a660dSBarry Smith }
152ec8511deSBarry Smith static int MatSolveTransAdd_SeqDense(Mat matin,Vec xx,Vec zz, Vec yy)
153da3a660dSBarry Smith {
154ec8511deSBarry Smith   Mat_SeqDense  *mat = (Mat_SeqDense *) matin->data;
1556abc6512SBarry Smith   int           one = 1, info,ierr;
1566abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
157da3a660dSBarry Smith   Vec           tmp;
158da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
159da3a660dSBarry Smith   if (yy == zz) {
16078b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
161464493b3SBarry Smith     PLogObjectParent(matin,tmp);
16278b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
163da3a660dSBarry Smith   }
164416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
165752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
166da3a660dSBarry Smith   if (mat->pivots) {
16748d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
168da3a660dSBarry Smith   }
169da3a660dSBarry Smith   else {
17048d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
171da3a660dSBarry Smith   }
172ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
173da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
174da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
175da3a660dSBarry Smith   return 0;
176da3a660dSBarry Smith }
177289bc588SBarry Smith /* ------------------------------------------------------------------*/
178ec8511deSBarry Smith static int MatRelax_SeqDense(Mat matin,Vec bb,double omega,MatSORType flag,
17920e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
180289bc588SBarry Smith {
181ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
182289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
1836abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
184289bc588SBarry Smith 
185289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
186289bc588SBarry Smith     /* this is a hack fix, should have another version without
187289bc588SBarry Smith        the second BLdot */
188bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
189289bc588SBarry Smith   }
190289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
191289bc588SBarry Smith   while (its--) {
192289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
193289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
194f1747703SBarry Smith #if defined(PETSC_COMPLEX)
195f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
196f1747703SBarry Smith            not happy about returning a double complex */
197f1747703SBarry Smith         int    _i;
198f1747703SBarry Smith         Scalar sum = b[i];
199f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
200f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
201f1747703SBarry Smith         }
202f1747703SBarry Smith         xt = sum;
203f1747703SBarry Smith #else
204289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
205f1747703SBarry Smith #endif
206d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
207289bc588SBarry Smith       }
208289bc588SBarry Smith     }
209289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
210289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
211f1747703SBarry Smith #if defined(PETSC_COMPLEX)
212f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
213f1747703SBarry Smith            not happy about returning a double complex */
214f1747703SBarry Smith         int    _i;
215f1747703SBarry Smith         Scalar sum = b[i];
216f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
217f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
218f1747703SBarry Smith         }
219f1747703SBarry Smith         xt = sum;
220f1747703SBarry Smith #else
221289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
222f1747703SBarry Smith #endif
223d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
224289bc588SBarry Smith       }
225289bc588SBarry Smith     }
226289bc588SBarry Smith   }
227289bc588SBarry Smith   return 0;
228289bc588SBarry Smith }
229289bc588SBarry Smith 
230289bc588SBarry Smith /* -----------------------------------------------------------------*/
231ec8511deSBarry Smith static int MatMultTrans_SeqDense(Mat matin,Vec xx,Vec yy)
232289bc588SBarry Smith {
233ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
234289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
235289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
236289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
23748d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
238289bc588SBarry Smith   return 0;
239289bc588SBarry Smith }
240ec8511deSBarry Smith static int MatMult_SeqDense(Mat matin,Vec xx,Vec yy)
241289bc588SBarry Smith {
242ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
243289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
244289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
245289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
24648d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
247289bc588SBarry Smith   return 0;
248289bc588SBarry Smith }
249ec8511deSBarry Smith static int MatMultAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy)
250289bc588SBarry Smith {
251ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
252289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2536abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
254289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
255416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
25648d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
257289bc588SBarry Smith   return 0;
258289bc588SBarry Smith }
259ec8511deSBarry Smith static int MatMultTransAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy)
260289bc588SBarry Smith {
261ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
262289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
263289bc588SBarry Smith   int          _One=1;
2646abc6512SBarry Smith   Scalar       _DOne=1.0;
265289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
266289bc588SBarry Smith   VecGetArray(zz,&z);
267416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
26848d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
269289bc588SBarry Smith   return 0;
270289bc588SBarry Smith }
271289bc588SBarry Smith 
272289bc588SBarry Smith /* -----------------------------------------------------------------*/
27348d91487SBarry Smith static int MatGetRow_SeqDense(Mat matin,int row,int *ncols,int **cols,Scalar **vals)
274289bc588SBarry Smith {
275ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
276289bc588SBarry Smith   Scalar       *v;
277289bc588SBarry Smith   int          i;
278289bc588SBarry Smith   *ncols = mat->n;
279289bc588SBarry Smith   if (cols) {
28078b31e54SBarry Smith     *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols);
28173c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
282289bc588SBarry Smith   }
283289bc588SBarry Smith   if (vals) {
28478b31e54SBarry Smith     *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
285289bc588SBarry Smith     v = mat->v + row;
28673c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
287289bc588SBarry Smith   }
288289bc588SBarry Smith   return 0;
289289bc588SBarry Smith }
29048d91487SBarry Smith static int MatRestoreRow_SeqDense(Mat matin,int row,int *ncols,int **cols,Scalar **vals)
291289bc588SBarry Smith {
29278b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
29378b31e54SBarry Smith   if (vals) { PETSCFREE(*vals); }
294289bc588SBarry Smith   return 0;
295289bc588SBarry Smith }
296289bc588SBarry Smith /* ----------------------------------------------------------------*/
297ec8511deSBarry Smith static int MatInsert_SeqDense(Mat matin,int m,int *indexm,int n,
298e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
299289bc588SBarry Smith {
300ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
301289bc588SBarry Smith   int          i,j;
302d6dfbf8fSBarry Smith 
303289bc588SBarry Smith   if (!mat->roworiented) {
304dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
305289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
306289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
307289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
308289bc588SBarry Smith         }
309289bc588SBarry Smith       }
310289bc588SBarry Smith     }
311289bc588SBarry Smith     else {
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     }
318e8d4e0b9SBarry Smith   }
319e8d4e0b9SBarry Smith   else {
320dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
321e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
322e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
323e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
324e8d4e0b9SBarry Smith         }
325e8d4e0b9SBarry Smith       }
326e8d4e0b9SBarry Smith     }
327289bc588SBarry Smith     else {
328289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
329289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
330289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
331289bc588SBarry Smith         }
332289bc588SBarry Smith       }
333289bc588SBarry Smith     }
334e8d4e0b9SBarry Smith   }
335289bc588SBarry Smith   return 0;
336289bc588SBarry Smith }
337e8d4e0b9SBarry Smith 
338289bc588SBarry Smith /* -----------------------------------------------------------------*/
339ec8511deSBarry Smith static int MatCopyPrivate_SeqDense(Mat matin,Mat *newmat)
340289bc588SBarry Smith {
34148d91487SBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data,*l;
342289bc588SBarry Smith   int          ierr;
343289bc588SBarry Smith   Mat          newi;
34448d91487SBarry Smith 
34548d91487SBarry Smith   ierr = MatCreateSeqDense(matin->comm,mat->m,mat->n,&newi);CHKERRQ(ierr);
346ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
347416022c9SBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
348289bc588SBarry Smith   *newmat = newi;
349289bc588SBarry Smith   return 0;
350289bc588SBarry Smith }
351da3a660dSBarry Smith #include "viewer.h"
352289bc588SBarry Smith 
353ec8511deSBarry Smith int MatView_SeqDense(PetscObject obj,Viewer ptr)
354289bc588SBarry Smith {
355289bc588SBarry Smith   Mat           matin = (Mat) obj;
356ec8511deSBarry Smith   Mat_SeqDense  *mat = (Mat_SeqDense *) matin->data;
357289bc588SBarry Smith   Scalar        *v;
358d636dbe3SBarry Smith   int           i,j,ierr;
35923242f5aSBarry Smith   PetscObject   vobj = (PetscObject) ptr;
360da3a660dSBarry Smith 
36123242f5aSBarry Smith   if (ptr == 0) {
36221b0d8fbSLois Curfman McInnes     ptr = STDOUT_VIEWER_SELF; vobj = (PetscObject) ptr;
36323242f5aSBarry Smith   }
36423242f5aSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
3656f75cfa4SBarry Smith     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
366da3a660dSBarry Smith   }
367da3a660dSBarry Smith   else {
368d636dbe3SBarry Smith     FILE *fd;
369d636dbe3SBarry Smith     int format;
370d636dbe3SBarry Smith     ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr);
371d636dbe3SBarry Smith     ierr = ViewerFileGetFormat_Private(ptr,&format); CHKERRQ(ierr);
372f72585f2SLois Curfman McInnes     if (format == FILE_FORMAT_INFO) {
373f72585f2SLois Curfman McInnes       /* do nothing for now */
374f72585f2SLois Curfman McInnes     }
375f72585f2SLois Curfman McInnes     else {
376289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
377289bc588SBarry Smith         v = mat->v + i;
378289bc588SBarry Smith         for ( j=0; j<mat->n; j++ ) {
379289bc588SBarry Smith #if defined(PETSC_COMPLEX)
3808e37dc09SBarry Smith           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
381289bc588SBarry Smith #else
3828e37dc09SBarry Smith           fprintf(fd,"%6.4e ",*v); v += mat->m;
383289bc588SBarry Smith #endif
384289bc588SBarry Smith         }
3858e37dc09SBarry Smith         fprintf(fd,"\n");
386289bc588SBarry Smith       }
387da3a660dSBarry Smith     }
388f72585f2SLois Curfman McInnes   }
389289bc588SBarry Smith   return 0;
390289bc588SBarry Smith }
391289bc588SBarry Smith 
392289bc588SBarry Smith 
393ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
394289bc588SBarry Smith {
395289bc588SBarry Smith   Mat          mat = (Mat) obj;
396ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
397a5a9c739SBarry Smith #if defined(PETSC_LOG)
398a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
399a5a9c739SBarry Smith #endif
40078b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
40178b31e54SBarry Smith   PETSCFREE(l);
402a5a9c739SBarry Smith   PLogObjectDestroy(mat);
4039e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
404289bc588SBarry Smith   return 0;
405289bc588SBarry Smith }
406289bc588SBarry Smith 
407ec8511deSBarry Smith static int MatTranspose_SeqDense(Mat matin,Mat *matout)
408289bc588SBarry Smith {
409ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
410d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
411d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
41248b35521SBarry Smith 
413d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
41448b35521SBarry Smith   if (!matout) { /* in place transpose */
41548d91487SBarry Smith     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Not for rectangular matrix in place");
416d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
417289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
418d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
419d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
420d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
421289bc588SBarry Smith       }
422289bc588SBarry Smith     }
42348b35521SBarry Smith   }
424d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
425d3e5ee88SLois Curfman McInnes     int ierr;
426d3e5ee88SLois Curfman McInnes     Mat tmat;
427ec8511deSBarry Smith     Mat_SeqDense *tmatd;
428d3e5ee88SLois Curfman McInnes     Scalar *v2;
429fafbff53SBarry Smith     ierr = MatCreateSeqDense(matin->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr);
430ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
4310de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
432d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
4330de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
434d3e5ee88SLois Curfman McInnes     }
435d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
436d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
437d3e5ee88SLois Curfman McInnes     *matout = tmat;
43848b35521SBarry Smith   }
439289bc588SBarry Smith   return 0;
440289bc588SBarry Smith }
441289bc588SBarry Smith 
442ec8511deSBarry Smith static int MatEqual_SeqDense(Mat matin1,Mat matin2)
443289bc588SBarry Smith {
444ec8511deSBarry Smith   Mat_SeqDense *mat1 = (Mat_SeqDense *) matin1->data;
445ec8511deSBarry Smith   Mat_SeqDense *mat2 = (Mat_SeqDense *) matin2->data;
446289bc588SBarry Smith   int          i;
447289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
448289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
449289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
450289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
451289bc588SBarry Smith     if (*v1 != *v2) return 0;
452289bc588SBarry Smith     v1++; v2++;
453289bc588SBarry Smith   }
454289bc588SBarry Smith   return 1;
455289bc588SBarry Smith }
456289bc588SBarry Smith 
457ec8511deSBarry Smith static int MatGetDiagonal_SeqDense(Mat matin,Vec v)
458289bc588SBarry Smith {
459ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
4606abc6512SBarry Smith   int          i, n;
4616abc6512SBarry Smith   Scalar       *x;
462289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
463ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
464289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
465289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
466289bc588SBarry Smith   }
467289bc588SBarry Smith   return 0;
468289bc588SBarry Smith }
469289bc588SBarry Smith 
470ec8511deSBarry Smith static int MatScale_SeqDense(Mat matin,Vec ll,Vec rr)
471289bc588SBarry Smith {
472ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
473da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
474da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
47528988994SBarry Smith   if (ll) {
476da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
477ec8511deSBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size");
478da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
479da3a660dSBarry Smith       x = l[i];
480da3a660dSBarry Smith       v = mat->v + i;
481da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
482da3a660dSBarry Smith     }
483da3a660dSBarry Smith   }
48428988994SBarry Smith   if (rr) {
485da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
486ec8511deSBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size");
487da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
488da3a660dSBarry Smith       x = r[i];
489da3a660dSBarry Smith       v = mat->v + i*m;
490da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
491da3a660dSBarry Smith     }
492da3a660dSBarry Smith   }
493289bc588SBarry Smith   return 0;
494289bc588SBarry Smith }
495289bc588SBarry Smith 
496ec8511deSBarry Smith static int MatNorm_SeqDense(Mat matin,MatNormType type,double *norm)
497289bc588SBarry Smith {
498ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
499289bc588SBarry Smith   Scalar       *v = mat->v;
500289bc588SBarry Smith   double       sum = 0.0;
501289bc588SBarry Smith   int          i, j;
502289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
503289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
504289bc588SBarry Smith #if defined(PETSC_COMPLEX)
505289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
506289bc588SBarry Smith #else
507289bc588SBarry Smith       sum += (*v)*(*v); v++;
508289bc588SBarry Smith #endif
509289bc588SBarry Smith     }
510289bc588SBarry Smith     *norm = sqrt(sum);
511289bc588SBarry Smith   }
512289bc588SBarry Smith   else if (type == NORM_1) {
513289bc588SBarry Smith     *norm = 0.0;
514289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
515289bc588SBarry Smith       sum = 0.0;
516289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
517289bc588SBarry Smith #if defined(PETSC_COMPLEX)
518289bc588SBarry Smith         sum += abs(*v++);
519289bc588SBarry Smith #else
520289bc588SBarry Smith         sum += fabs(*v++);
521289bc588SBarry Smith #endif
522289bc588SBarry Smith       }
523289bc588SBarry Smith       if (sum > *norm) *norm = sum;
524289bc588SBarry Smith     }
525289bc588SBarry Smith   }
526289bc588SBarry Smith   else if (type == NORM_INFINITY) {
527289bc588SBarry Smith     *norm = 0.0;
528289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
529289bc588SBarry Smith       v = mat->v + j;
530289bc588SBarry Smith       sum = 0.0;
531289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
532289bc588SBarry Smith #if defined(PETSC_COMPLEX)
533289bc588SBarry Smith         sum += abs(*v); v += mat->m;
534289bc588SBarry Smith #else
535289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
536289bc588SBarry Smith #endif
537289bc588SBarry Smith       }
538289bc588SBarry Smith       if (sum > *norm) *norm = sum;
539289bc588SBarry Smith     }
540289bc588SBarry Smith   }
541289bc588SBarry Smith   else {
54248d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
543289bc588SBarry Smith   }
544289bc588SBarry Smith   return 0;
545289bc588SBarry Smith }
546289bc588SBarry Smith 
547ec8511deSBarry Smith static int MatSetOption_SeqDense(Mat aijin,MatOption op)
548289bc588SBarry Smith {
549ec8511deSBarry Smith   Mat_SeqDense *aij = (Mat_SeqDense *) aijin->data;
550289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
551289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
552289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
553289bc588SBarry Smith   return 0;
554289bc588SBarry Smith }
555289bc588SBarry Smith 
556ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
5576f0a148fSBarry Smith {
558ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
559416022c9SBarry Smith   PetscZero(l->v,l->m*l->n*sizeof(Scalar));
5606f0a148fSBarry Smith   return 0;
5616f0a148fSBarry Smith }
5626f0a148fSBarry Smith 
563ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
5646f0a148fSBarry Smith {
565ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
5666abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
5676f0a148fSBarry Smith   Scalar       *slot;
56878b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
56978b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5706f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5716f0a148fSBarry Smith     slot = l->v + rows[i];
5726f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5736f0a148fSBarry Smith   }
5746f0a148fSBarry Smith   if (diag) {
5756f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5766f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
577260925b8SBarry Smith       *slot = *diag;
5786f0a148fSBarry Smith     }
5796f0a148fSBarry Smith   }
580260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5816f0a148fSBarry Smith   return 0;
5826f0a148fSBarry Smith }
583557bce09SLois Curfman McInnes 
584ec8511deSBarry Smith static int MatGetSize_SeqDense(Mat matin,int *m,int *n)
585557bce09SLois Curfman McInnes {
586ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
587557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
588557bce09SLois Curfman McInnes   return 0;
589557bce09SLois Curfman McInnes }
590557bce09SLois Curfman McInnes 
591ec8511deSBarry Smith static int MatGetOwnershipRange_SeqDense(Mat matin,int *m,int *n)
592d3e5ee88SLois Curfman McInnes {
593ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
594d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
595d3e5ee88SLois Curfman McInnes   return 0;
596d3e5ee88SLois Curfman McInnes }
597d3e5ee88SLois Curfman McInnes 
598ec8511deSBarry Smith static int MatGetArray_SeqDense(Mat matin,Scalar **array)
59964e87e97SBarry Smith {
600ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
60164e87e97SBarry Smith   *array = mat->v;
60264e87e97SBarry Smith   return 0;
60364e87e97SBarry Smith }
6040754003eSLois Curfman McInnes 
6050754003eSLois Curfman McInnes 
606ec8511deSBarry Smith static int MatGetSubMatrixInPlace_SeqDense(Mat matin,IS isrow,IS iscol)
6070754003eSLois Curfman McInnes {
608ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
6090754003eSLois Curfman McInnes }
6100754003eSLois Curfman McInnes 
611ec8511deSBarry Smith static int MatGetSubMatrix_SeqDense(Mat matin,IS isrow,IS iscol,Mat *submat)
6120754003eSLois Curfman McInnes {
613ec8511deSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense *) matin->data;
6140754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
615160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
6160754003eSLois Curfman McInnes   Scalar       *vwork, *val;
6170754003eSLois Curfman McInnes   Mat          newmat;
6180754003eSLois Curfman McInnes 
61978b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
62078b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
62178b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
62278b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
6230754003eSLois Curfman McInnes 
62478b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
62578b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
62678b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
627416022c9SBarry Smith   PetscZero((char*)smap,oldcols*sizeof(int));
6280754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
6290754003eSLois Curfman McInnes 
6300754003eSLois Curfman McInnes   /* Create and fill new matrix */
63148d91487SBarry Smith   ierr = MatCreateSeqDense(matin->comm,nrows,ncols,&newmat);CHKERRQ(ierr);
6320754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
6330754003eSLois Curfman McInnes     nznew = 0;
6340754003eSLois Curfman McInnes     val   = mat->v + irow[i];
6350754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
6360754003eSLois Curfman McInnes       if (smap[j]) {
6370754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
6380754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
6390754003eSLois Curfman McInnes       }
6400754003eSLois Curfman McInnes     }
641dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
64278b31e54SBarry Smith            CHKERRQ(ierr);
6430754003eSLois Curfman McInnes   }
64478b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
64578b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
6460754003eSLois Curfman McInnes 
6470754003eSLois Curfman McInnes   /* Free work space */
64878b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
64978b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
65078b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
6510754003eSLois Curfman McInnes   *submat = newmat;
6520754003eSLois Curfman McInnes   return 0;
6530754003eSLois Curfman McInnes }
6540754003eSLois Curfman McInnes 
655289bc588SBarry Smith /* -------------------------------------------------------------------*/
656ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense,
657ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
658ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
659ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
660ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
661ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
662ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
663ec8511deSBarry Smith        MatRelax_SeqDense,
664ec8511deSBarry Smith        MatTranspose_SeqDense,
665ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
666ec8511deSBarry Smith        MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense,
667289bc588SBarry Smith        0,0,
668ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
669ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
670ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
671ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
672ec8511deSBarry Smith        0,0,MatGetArray_SeqDense,0,0,
673ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
674*1987afe7SBarry Smith        MatCopyPrivate_SeqDense,0,0,0,0,
675*1987afe7SBarry Smith        MatAXPY_SeqDense};
6760754003eSLois Curfman McInnes 
6774b828684SBarry Smith /*@C
678fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
679d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
680d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
681289bc588SBarry Smith 
68220563c6bSBarry Smith    Input Parameters:
6830c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
6840c775827SLois Curfman McInnes .  m - number of rows
6850c775827SLois Curfman McInnes .  n - number of column
68620563c6bSBarry Smith 
68720563c6bSBarry Smith    Output Parameter:
6880c775827SLois Curfman McInnes .  newmat - the matrix
68920563c6bSBarry Smith 
690dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
691d65003e9SLois Curfman McInnes 
692dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
69320563c6bSBarry Smith @*/
694fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat)
695289bc588SBarry Smith {
696ec8511deSBarry Smith   int          size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar);
697289bc588SBarry Smith   Mat          mat;
698ec8511deSBarry Smith   Mat_SeqDense *l;
699289bc588SBarry Smith   *newmat        = 0;
700ec8511deSBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
701a5a9c739SBarry Smith   PLogObjectCreate(mat);
702ec8511deSBarry Smith   l              = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l);
703416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
704ec8511deSBarry Smith   mat->destroy   = MatDestroy_SeqDense;
705ec8511deSBarry Smith   mat->view      = MatView_SeqDense;
706289bc588SBarry Smith   mat->data      = (void *) l;
707289bc588SBarry Smith   mat->factor    = 0;
708752f5567SLois Curfman McInnes   PLogObjectMemory(mat,sizeof(struct _Mat) + size);
709d6dfbf8fSBarry Smith 
710289bc588SBarry Smith   l->m           = m;
711289bc588SBarry Smith   l->n           = n;
712289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
713289bc588SBarry Smith   l->pivots      = 0;
714289bc588SBarry Smith   l->roworiented = 1;
715d6dfbf8fSBarry Smith 
716416022c9SBarry Smith   PetscZero(l->v,m*n*sizeof(Scalar));
717289bc588SBarry Smith   *newmat = mat;
718289bc588SBarry Smith   return 0;
719289bc588SBarry Smith }
720289bc588SBarry Smith 
721ec8511deSBarry Smith int MatCreate_SeqDense(Mat matin,Mat *newmat)
722289bc588SBarry Smith {
723ec8511deSBarry Smith   Mat_SeqDense *m = (Mat_SeqDense *) matin->data;
724fafbff53SBarry Smith   return MatCreateSeqDense(matin->comm,m->m,m->n,newmat);
725289bc588SBarry Smith }
726