xref: /petsc/src/mat/impls/dense/seq/dense.c (revision d636dbe3a9c413b447e52fc0f7f06f5e4018d190)
1cb512458SBarry Smith #ifndef lint
2*d636dbe3SBarry Smith static char vcid[] = "$Id: dense.c,v 1.53 1995/09/04 17:24: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 */
2044a69424SLois Curfman McInnes } Mat_Dense;
21289bc588SBarry Smith 
2220e1a0d4SLois Curfman McInnes static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz,
2320e1a0d4SLois Curfman McInnes                             int *nzalloc,int *mem)
24289bc588SBarry Smith {
2544a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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.*/
3649d8b64dSBarry Smith static int MatLUFactor_Dense(Mat matin,IS row,IS col,double f)
37289bc588SBarry Smith {
3844a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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);
46bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatLUFactor_Dense: Bad LU factorization");
47289bc588SBarry Smith   matin->factor = FACTOR_LU;
48289bc588SBarry Smith   return 0;
49289bc588SBarry Smith }
5049d8b64dSBarry Smith static int MatLUFactorSymbolic_Dense(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 }
5744a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact)
58289bc588SBarry Smith {
5949d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
60289bc588SBarry Smith }
6149d8b64dSBarry Smith static int MatCholeskyFactorSymbolic_Dense(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 }
6746fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact)
68289bc588SBarry Smith {
6949d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
70289bc588SBarry Smith }
7149d8b64dSBarry Smith static int MatCholeskyFactor_Dense(Mat matin,IS perm,double f)
72289bc588SBarry Smith {
7344a69424SLois Curfman McInnes   Mat_Dense    *mat = (Mat_Dense *) 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);
81bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_Dense: Bad factorization");
82289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
83289bc588SBarry Smith   return 0;
84289bc588SBarry Smith }
85289bc588SBarry Smith 
8644a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy)
87289bc588SBarry Smith {
8844a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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   }
101bbb6d6a8SBarry Smith   else SETERRQ(1,"MatSolve_Dense: Matrix must be factored to solve");
102bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolve_Dense: Bad solve");
103289bc588SBarry Smith   return 0;
104289bc588SBarry Smith }
10544a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy)
106da3a660dSBarry Smith {
10744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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   }
121bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_Dense: Bad solve");
122da3a660dSBarry Smith   return 0;
123da3a660dSBarry Smith }
12444a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
125da3a660dSBarry Smith {
12644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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   }
146bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_Dense: 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 }
15144a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy)
152da3a660dSBarry Smith {
15344a69424SLois Curfman McInnes   Mat_Dense  *mat = (Mat_Dense *) 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   }
173bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_Dense: 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 /* ------------------------------------------------------------------*/
17920e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag,
18020e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
181289bc588SBarry Smith {
18244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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 /* -----------------------------------------------------------------*/
23244a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy)
233289bc588SBarry Smith {
23444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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 }
24244a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy)
243289bc588SBarry Smith {
24444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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 }
25244a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
253289bc588SBarry Smith {
25444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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 }
26344a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
264289bc588SBarry Smith {
26544a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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 /* -----------------------------------------------------------------*/
27844a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols,
279289bc588SBarry Smith                         Scalar **vals)
280289bc588SBarry Smith {
28144a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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 }
29644a69424SLois Curfman McInnes static int MatRestoreRow_Dense(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 /* ----------------------------------------------------------------*/
30444a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n,
305e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
306289bc588SBarry Smith {
30744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) 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 /* -----------------------------------------------------------------*/
346ff7509f2SBarry Smith static int MatCopyPrivate_Dense(Mat matin,Mat *newmat)
347289bc588SBarry Smith {
34844a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
349289bc588SBarry Smith   int ierr;
350289bc588SBarry Smith   Mat newi;
35144a69424SLois Curfman McInnes   Mat_Dense *l;
352bbb6d6a8SBarry Smith   ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi);
353bbb6d6a8SBarry Smith   CHKERRQ(ierr);
35444a69424SLois Curfman McInnes   l = (Mat_Dense *) 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 
36144a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr)
362289bc588SBarry Smith {
363289bc588SBarry Smith   Mat         matin = (Mat) obj;
36444a69424SLois Curfman McInnes   Mat_Dense   *mat = (Mat_Dense *) matin->data;
365289bc588SBarry Smith   Scalar      *v;
366*d636dbe3SBarry 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 {
376*d636dbe3SBarry Smith     FILE *fd;
377*d636dbe3SBarry Smith     int format;
378*d636dbe3SBarry Smith     ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr);
379*d636dbe3SBarry 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 
40144a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj)
402289bc588SBarry Smith {
403289bc588SBarry Smith   Mat    mat = (Mat) obj;
40444a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) 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 
41548b35521SBarry Smith static int MatTranspose_Dense(Mat matin,Mat *matout)
416289bc588SBarry Smith {
41744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
418289bc588SBarry Smith   int    k,j;
419289bc588SBarry Smith   Scalar *v = mat->v, tmp;
42048b35521SBarry Smith 
42148b35521SBarry Smith   if (!matout) { /* in place transpose */
422289bc588SBarry Smith     if (mat->m != mat->n) {
42348b35521SBarry Smith       SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix");
424289bc588SBarry Smith     }
425289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
426289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
427289bc588SBarry Smith         tmp = v[j + k*mat->n];
428289bc588SBarry Smith         v[j + k*mat->n] = v[k + j*mat->n];
429289bc588SBarry Smith         v[k + j*mat->n] = tmp;
430289bc588SBarry Smith       }
431289bc588SBarry Smith     }
43248b35521SBarry Smith   }
43348b35521SBarry Smith   else {
43448b35521SBarry Smith     SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet");
43548b35521SBarry Smith   }
436289bc588SBarry Smith   return 0;
437289bc588SBarry Smith }
438289bc588SBarry Smith 
43944a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2)
440289bc588SBarry Smith {
44144a69424SLois Curfman McInnes   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
44244a69424SLois Curfman McInnes   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
443289bc588SBarry Smith   int    i;
444289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
445289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
446289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
447289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
448289bc588SBarry Smith     if (*v1 != *v2) return 0;
449289bc588SBarry Smith     v1++; v2++;
450289bc588SBarry Smith   }
451289bc588SBarry Smith   return 1;
452289bc588SBarry Smith }
453289bc588SBarry Smith 
45446fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v)
455289bc588SBarry Smith {
45644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
4576abc6512SBarry Smith   int       i, n;
4586abc6512SBarry Smith   Scalar    *x;
459289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
460bbb6d6a8SBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_Dense:Nonconforming mat and vec");
461289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
462289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
463289bc588SBarry Smith   }
464289bc588SBarry Smith   return 0;
465289bc588SBarry Smith }
466289bc588SBarry Smith 
46744a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
468289bc588SBarry Smith {
46944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
470da3a660dSBarry Smith   Scalar *l,*r,x,*v;
471da3a660dSBarry Smith   int    i,j,m = mat->m, n = mat->n;
47228988994SBarry Smith   if (ll) {
473da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
474bbb6d6a8SBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_Dense:Left scaling vec wrong size");
475da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
476da3a660dSBarry Smith       x = l[i];
477da3a660dSBarry Smith       v = mat->v + i;
478da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
479da3a660dSBarry Smith     }
480da3a660dSBarry Smith   }
48128988994SBarry Smith   if (rr) {
482da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
483bbb6d6a8SBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_Dense:Right scaling vec wrong size");
484da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
485da3a660dSBarry Smith       x = r[i];
486da3a660dSBarry Smith       v = mat->v + i*m;
487da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
488da3a660dSBarry Smith     }
489da3a660dSBarry Smith   }
490289bc588SBarry Smith   return 0;
491289bc588SBarry Smith }
492289bc588SBarry Smith 
4937c17b2cfSLois Curfman McInnes static int MatNorm_Dense(Mat matin,MatNormType type,double *norm)
494289bc588SBarry Smith {
49544a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
496289bc588SBarry Smith   Scalar *v = mat->v;
497289bc588SBarry Smith   double sum = 0.0;
498289bc588SBarry Smith   int    i, j;
499289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
500289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
501289bc588SBarry Smith #if defined(PETSC_COMPLEX)
502289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
503289bc588SBarry Smith #else
504289bc588SBarry Smith       sum += (*v)*(*v); v++;
505289bc588SBarry Smith #endif
506289bc588SBarry Smith     }
507289bc588SBarry Smith     *norm = sqrt(sum);
508289bc588SBarry Smith   }
509289bc588SBarry Smith   else if (type == NORM_1) {
510289bc588SBarry Smith     *norm = 0.0;
511289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
512289bc588SBarry Smith       sum = 0.0;
513289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
514289bc588SBarry Smith #if defined(PETSC_COMPLEX)
515289bc588SBarry Smith         sum += abs(*v++);
516289bc588SBarry Smith #else
517289bc588SBarry Smith         sum += fabs(*v++);
518289bc588SBarry Smith #endif
519289bc588SBarry Smith       }
520289bc588SBarry Smith       if (sum > *norm) *norm = sum;
521289bc588SBarry Smith     }
522289bc588SBarry Smith   }
523289bc588SBarry Smith   else if (type == NORM_INFINITY) {
524289bc588SBarry Smith     *norm = 0.0;
525289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
526289bc588SBarry Smith       v = mat->v + j;
527289bc588SBarry Smith       sum = 0.0;
528289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
529289bc588SBarry Smith #if defined(PETSC_COMPLEX)
530289bc588SBarry Smith         sum += abs(*v); v += mat->m;
531289bc588SBarry Smith #else
532289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
533289bc588SBarry Smith #endif
534289bc588SBarry Smith       }
535289bc588SBarry Smith       if (sum > *norm) *norm = sum;
536289bc588SBarry Smith     }
537289bc588SBarry Smith   }
538289bc588SBarry Smith   else {
539bbb6d6a8SBarry Smith     SETERRQ(1,"MatNorm_Dense:No two norm yet");
540289bc588SBarry Smith   }
541289bc588SBarry Smith   return 0;
542289bc588SBarry Smith }
543289bc588SBarry Smith 
54420e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op)
545289bc588SBarry Smith {
54644a69424SLois Curfman McInnes   Mat_Dense *aij = (Mat_Dense *) aijin->data;
547289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
548289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
549289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
550289bc588SBarry Smith   return 0;
551289bc588SBarry Smith }
552289bc588SBarry Smith 
55344a69424SLois Curfman McInnes static int MatZero_Dense(Mat A)
5546f0a148fSBarry Smith {
55544a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
55678b31e54SBarry Smith   PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
5576f0a148fSBarry Smith   return 0;
5586f0a148fSBarry Smith }
5596f0a148fSBarry Smith 
56044a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag)
5616f0a148fSBarry Smith {
56244a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
5636abc6512SBarry Smith   int    n = l->n, i, j,ierr,N, *rows;
5646f0a148fSBarry Smith   Scalar *slot;
56578b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
56678b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5676f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5686f0a148fSBarry Smith     slot = l->v + rows[i];
5696f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5706f0a148fSBarry Smith   }
5716f0a148fSBarry Smith   if (diag) {
5726f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5736f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
574260925b8SBarry Smith       *slot = *diag;
5756f0a148fSBarry Smith     }
5766f0a148fSBarry Smith   }
577260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5786f0a148fSBarry Smith   return 0;
5796f0a148fSBarry Smith }
580557bce09SLois Curfman McInnes 
581fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n)
582557bce09SLois Curfman McInnes {
58344a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
584557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
585557bce09SLois Curfman McInnes   return 0;
586557bce09SLois Curfman McInnes }
587557bce09SLois Curfman McInnes 
58844a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array)
58964e87e97SBarry Smith {
59044a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
59164e87e97SBarry Smith   *array = mat->v;
59264e87e97SBarry Smith   return 0;
59364e87e97SBarry Smith }
5940754003eSLois Curfman McInnes 
5950754003eSLois Curfman McInnes 
5960754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol)
5970754003eSLois Curfman McInnes {
598bbb6d6a8SBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_Dense: not done");
5990754003eSLois Curfman McInnes }
6000754003eSLois Curfman McInnes 
6010754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat)
6020754003eSLois Curfman McInnes {
6030754003eSLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
6040754003eSLois Curfman McInnes   int     nznew, *smap, i, j, ierr, oldcols = mat->n;
605160018dcSBarry Smith   int     *irow, *icol, nrows, ncols, *cwork;
6060754003eSLois Curfman McInnes   Scalar  *vwork, *val;
6070754003eSLois Curfman McInnes   Mat     newmat;
6080754003eSLois Curfman McInnes 
60978b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
61078b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
61178b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
61278b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
6130754003eSLois Curfman McInnes 
61478b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
61578b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
61678b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
61778b31e54SBarry Smith   PETSCMEMSET((char*)smap,0,oldcols*sizeof(int));
6180754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
6190754003eSLois Curfman McInnes 
6200754003eSLois Curfman McInnes   /* Create and fill new matrix */
6210754003eSLois Curfman McInnes   ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat);
62278b31e54SBarry Smith          CHKERRQ(ierr);
6230754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
6240754003eSLois Curfman McInnes     nznew = 0;
6250754003eSLois Curfman McInnes     val   = mat->v + irow[i];
6260754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
6270754003eSLois Curfman McInnes       if (smap[j]) {
6280754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
6290754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
6300754003eSLois Curfman McInnes       }
6310754003eSLois Curfman McInnes     }
632edae2e7dSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES);
63378b31e54SBarry Smith            CHKERRQ(ierr);
6340754003eSLois Curfman McInnes   }
63578b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
63678b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
6370754003eSLois Curfman McInnes 
6380754003eSLois Curfman McInnes   /* Free work space */
63978b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
64078b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
64178b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
6420754003eSLois Curfman McInnes   *submat = newmat;
6430754003eSLois Curfman McInnes   return 0;
6440754003eSLois Curfman McInnes }
6450754003eSLois Curfman McInnes 
646289bc588SBarry Smith /* -------------------------------------------------------------------*/
64744a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense,
64844a69424SLois Curfman McInnes        MatGetRow_Dense, MatRestoreRow_Dense,
64944a69424SLois Curfman McInnes        MatMult_Dense, MatMultAdd_Dense,
65044a69424SLois Curfman McInnes        MatMultTrans_Dense, MatMultTransAdd_Dense,
65144a69424SLois Curfman McInnes        MatSolve_Dense,MatSolveAdd_Dense,
65244a69424SLois Curfman McInnes        MatSolveTrans_Dense,MatSolveTransAdd_Dense,
65346fc4a8fSLois Curfman McInnes        MatLUFactor_Dense,MatCholeskyFactor_Dense,
65444a69424SLois Curfman McInnes        MatRelax_Dense,
65548b35521SBarry Smith        MatTranspose_Dense,
656fa9ec3f1SLois Curfman McInnes        MatGetInfo_Dense,MatEqual_Dense,
65746fc4a8fSLois Curfman McInnes        MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense,
658289bc588SBarry Smith        0,0,
659681d1451SLois Curfman McInnes        0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0,
66044a69424SLois Curfman McInnes        MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense,
66146fc4a8fSLois Curfman McInnes        MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense,
662fa9ec3f1SLois Curfman McInnes        MatGetSize_Dense,MatGetSize_Dense,0,
66307a0e7adSLois Curfman McInnes        0,0,MatGetArray_Dense,0,0,
66407a0e7adSLois Curfman McInnes        MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense,
665ff7509f2SBarry Smith        MatCopyPrivate_Dense};
6660754003eSLois Curfman McInnes 
6674b828684SBarry Smith /*@C
66820563c6bSBarry Smith    MatCreateSequentialDense - Creates a sequential dense matrix that
669d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
670d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
671289bc588SBarry Smith 
67220563c6bSBarry Smith    Input Parameters:
6730c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
6740c775827SLois Curfman McInnes .  m - number of rows
6750c775827SLois Curfman McInnes .  n - number of column
67620563c6bSBarry Smith 
67720563c6bSBarry Smith    Output Parameter:
6780c775827SLois Curfman McInnes .  newmat - the matrix
67920563c6bSBarry Smith 
680dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
681d65003e9SLois Curfman McInnes 
682dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
68320563c6bSBarry Smith @*/
6846b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat)
685289bc588SBarry Smith {
68644a69424SLois Curfman McInnes   int       size = sizeof(Mat_Dense) + m*n*sizeof(Scalar);
687289bc588SBarry Smith   Mat mat;
68844a69424SLois Curfman McInnes   Mat_Dense    *l;
689289bc588SBarry Smith   *newmat        = 0;
6906b5873e3SBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm);
691a5a9c739SBarry Smith   PLogObjectCreate(mat);
69278b31e54SBarry Smith   l              = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l);
693289bc588SBarry Smith   mat->ops       = &MatOps;
69444a69424SLois Curfman McInnes   mat->destroy   = MatDestroy_Dense;
69544a69424SLois Curfman McInnes   mat->view      = MatView_Dense;
696289bc588SBarry Smith   mat->data      = (void *) l;
697289bc588SBarry Smith   mat->factor    = 0;
698752f5567SLois Curfman McInnes   PLogObjectMemory(mat,sizeof(struct _Mat) + size);
699d6dfbf8fSBarry Smith 
700289bc588SBarry Smith   l->m           = m;
701289bc588SBarry Smith   l->n           = n;
702289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
703289bc588SBarry Smith   l->pivots      = 0;
704289bc588SBarry Smith   l->roworiented = 1;
705d6dfbf8fSBarry Smith 
70678b31e54SBarry Smith   PETSCMEMSET(l->v,0,m*n*sizeof(Scalar));
707289bc588SBarry Smith   *newmat = mat;
708289bc588SBarry Smith   return 0;
709289bc588SBarry Smith }
710289bc588SBarry Smith 
71144a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat)
712289bc588SBarry Smith {
71344a69424SLois Curfman McInnes   Mat_Dense *m = (Mat_Dense *) matin->data;
7146b5873e3SBarry Smith   return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat);
715289bc588SBarry Smith }
716