xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 5392566e307f4c743d2e3a1f9d4485cd364fbdad)
1*5392566eSBarry Smith 
2*5392566eSBarry Smith 
3cb512458SBarry Smith #ifndef lint
4*5392566eSBarry Smith static char vcid[] = "$Id: dense.c,v 1.40 1995/06/08 03:09:05 bsmith Exp bsmith $";
5cb512458SBarry Smith #endif
6289bc588SBarry Smith 
7289bc588SBarry Smith /*
8289bc588SBarry Smith     Standard Fortran style matrices
9289bc588SBarry Smith */
10289bc588SBarry Smith #include "ptscimpl.h"
11289bc588SBarry Smith #include "plapack.h"
12289bc588SBarry Smith #include "matimpl.h"
13289bc588SBarry Smith #include "math.h"
14289bc588SBarry Smith #include "vec/vecimpl.h"
1520e1a0d4SLois Curfman McInnes #include "pviewer.h"
16289bc588SBarry Smith 
17289bc588SBarry Smith typedef struct {
18289bc588SBarry Smith   Scalar *v;
19289bc588SBarry Smith   int    roworiented;
20289bc588SBarry Smith   int    m,n,pad;
21289bc588SBarry Smith   int    *pivots;   /* pivots in LU factorization */
2244a69424SLois Curfman McInnes } Mat_Dense;
23289bc588SBarry Smith 
2420e1a0d4SLois Curfman McInnes static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz,
2520e1a0d4SLois Curfman McInnes                             int *nzalloc,int *mem)
26289bc588SBarry Smith {
2744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
28289bc588SBarry Smith   int    i,N = mat->m*mat->n,count = 0;
29289bc588SBarry Smith   Scalar *v = mat->v;
30289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
31fa9ec3f1SLois Curfman McInnes   *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar);
32fa9ec3f1SLois Curfman McInnes   return 0;
33289bc588SBarry Smith }
34289bc588SBarry Smith 
35289bc588SBarry Smith /* ---------------------------------------------------------------*/
36289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
37289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
3844a69424SLois Curfman McInnes static int MatLUFactor_Dense(Mat matin,IS row,IS col)
39289bc588SBarry Smith {
4044a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
416abc6512SBarry Smith   int    info;
42289bc588SBarry Smith   if (!mat->pivots) {
4378b31e54SBarry Smith     mat->pivots = (int *) PETSCMALLOC( mat->m*sizeof(int) );
4478b31e54SBarry Smith     CHKPTRQ(mat->pivots);
45289bc588SBarry Smith   }
46289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
4778b31e54SBarry Smith   if (info) SETERRQ(1,"Bad LU factorization");
48289bc588SBarry Smith   matin->factor = FACTOR_LU;
49289bc588SBarry Smith   return 0;
50289bc588SBarry Smith }
5144a69424SLois Curfman McInnes static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact)
52289bc588SBarry Smith {
53289bc588SBarry Smith   int ierr;
5478b31e54SBarry Smith   if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0);
55289bc588SBarry Smith   return 0;
56289bc588SBarry Smith }
5744a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact)
58289bc588SBarry Smith {
5920563c6bSBarry Smith   return MatLUFactor(*fact,0,0);
60289bc588SBarry Smith }
6146fc4a8fSLois Curfman McInnes static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,Mat *fact)
62289bc588SBarry Smith {
63289bc588SBarry Smith   int ierr;
6478b31e54SBarry Smith   if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0);
65289bc588SBarry Smith   return 0;
66289bc588SBarry Smith }
6746fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact)
68289bc588SBarry Smith {
6920563c6bSBarry Smith   return MatCholeskyFactor(*fact,0);
70289bc588SBarry Smith }
7146fc4a8fSLois Curfman McInnes static int MatCholeskyFactor_Dense(Mat matin,IS perm)
72289bc588SBarry Smith {
7344a69424SLois Curfman McInnes   Mat_Dense    *mat = (Mat_Dense *) matin->data;
746abc6512SBarry Smith   int       info;
7578b31e54SBarry Smith   if (mat->pivots) {PETSCFREE(mat->pivots); mat->pivots = 0;}
76289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
7778b31e54SBarry Smith   if (info) SETERRQ(1,"Bad Cholesky factorization");
78289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
79289bc588SBarry Smith   return 0;
80289bc588SBarry Smith }
81289bc588SBarry Smith 
8244a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy)
83289bc588SBarry Smith {
8444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
856abc6512SBarry Smith   int    one = 1, info;
866abc6512SBarry Smith   Scalar *x, *y;
87289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
8878b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
899e25ed09SBarry Smith   if (matin->factor == FACTOR_LU) {
90289bc588SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
91289bc588SBarry Smith               y, &mat->m, &info );
92289bc588SBarry Smith   }
939e25ed09SBarry Smith   else if (matin->factor == FACTOR_CHOLESKY){
94289bc588SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
95289bc588SBarry Smith               y, &mat->m, &info );
96289bc588SBarry Smith   }
9778b31e54SBarry Smith   else SETERRQ(1,"Matrix must be factored to solve");
9878b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
99289bc588SBarry Smith   return 0;
100289bc588SBarry Smith }
10144a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy)
102da3a660dSBarry Smith {
10344a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
1046abc6512SBarry Smith   int    one = 1, info;
1056abc6512SBarry Smith   Scalar *x, *y;
106da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
10778b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
108da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
109da3a660dSBarry Smith   if (mat->pivots) {
110da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
111da3a660dSBarry Smith               y, &mat->m, &info );
112da3a660dSBarry Smith   }
113da3a660dSBarry Smith   else {
114da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
115da3a660dSBarry Smith               y, &mat->m, &info );
116da3a660dSBarry Smith   }
11778b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
118da3a660dSBarry Smith   return 0;
119da3a660dSBarry Smith }
12044a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
121da3a660dSBarry Smith {
12244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
1236abc6512SBarry Smith   int    one = 1, info,ierr;
1246abc6512SBarry Smith   Scalar *x, *y, sone = 1.0;
125da3a660dSBarry Smith   Vec    tmp = 0;
126da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
127da3a660dSBarry Smith   if (yy == zz) {
12878b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
12978b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
130da3a660dSBarry Smith   }
13178b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
132da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
133da3a660dSBarry Smith   if (mat->pivots) {
134da3a660dSBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
135da3a660dSBarry Smith               y, &mat->m, &info );
136da3a660dSBarry Smith   }
137da3a660dSBarry Smith   else {
138da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
139da3a660dSBarry Smith               y, &mat->m, &info );
140da3a660dSBarry Smith   }
14178b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
142da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
143da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
144da3a660dSBarry Smith   return 0;
145da3a660dSBarry Smith }
14644a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy)
147da3a660dSBarry Smith {
14844a69424SLois Curfman McInnes   Mat_Dense  *mat = (Mat_Dense *) matin->data;
1496abc6512SBarry Smith   int     one = 1, info,ierr;
1506abc6512SBarry Smith   Scalar  *x, *y, sone = 1.0;
151da3a660dSBarry Smith   Vec     tmp;
152da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
153da3a660dSBarry Smith   if (yy == zz) {
15478b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
15578b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
156da3a660dSBarry Smith   }
15778b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
158da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
159da3a660dSBarry Smith   if (mat->pivots) {
160da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
161da3a660dSBarry Smith               y, &mat->m, &info );
162da3a660dSBarry Smith   }
163da3a660dSBarry Smith   else {
164da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
165da3a660dSBarry Smith               y, &mat->m, &info );
166da3a660dSBarry Smith   }
16778b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
168da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
169da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
170da3a660dSBarry Smith   return 0;
171da3a660dSBarry Smith }
172289bc588SBarry Smith /* ------------------------------------------------------------------*/
17320e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag,
17420e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
175289bc588SBarry Smith {
17644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
177289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
1786abc6512SBarry Smith   int    o = 1,ierr, m = mat->m, i;
179289bc588SBarry Smith 
180289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
181289bc588SBarry Smith     /* this is a hack fix, should have another version without
182289bc588SBarry Smith        the second BLdot */
18378b31e54SBarry Smith     if ((ierr = VecSet(&zero,xx))) SETERRQ(ierr,0);
184289bc588SBarry Smith   }
185289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
186289bc588SBarry Smith   while (its--) {
187289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
188289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
189289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
190d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
191289bc588SBarry Smith       }
192289bc588SBarry Smith     }
193289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
194289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
195289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
196d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
197289bc588SBarry Smith       }
198289bc588SBarry Smith     }
199289bc588SBarry Smith   }
200289bc588SBarry Smith   return 0;
201289bc588SBarry Smith }
202289bc588SBarry Smith 
203289bc588SBarry Smith /* -----------------------------------------------------------------*/
20444a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy)
205289bc588SBarry Smith {
20644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
207289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
208289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
209289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
210289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
211289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
212289bc588SBarry Smith   return 0;
213289bc588SBarry Smith }
21444a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy)
215289bc588SBarry Smith {
21644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
217289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
218289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
219289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
220289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
221289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
222289bc588SBarry Smith   return 0;
223289bc588SBarry Smith }
22444a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
225289bc588SBarry Smith {
22644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
227289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
2286abc6512SBarry Smith   int    _One=1; Scalar _DOne=1.0;
229289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
23078b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
231289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
232289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
233289bc588SBarry Smith   return 0;
234289bc588SBarry Smith }
23544a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
236289bc588SBarry Smith {
23744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
238289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
239289bc588SBarry Smith   int    _One=1;
2406abc6512SBarry Smith   Scalar _DOne=1.0;
241289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
242289bc588SBarry Smith   VecGetArray(zz,&z);
24378b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
244289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
245289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
246289bc588SBarry Smith   return 0;
247289bc588SBarry Smith }
248289bc588SBarry Smith 
249289bc588SBarry Smith /* -----------------------------------------------------------------*/
25044a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols,
251289bc588SBarry Smith                         Scalar **vals)
252289bc588SBarry Smith {
25344a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
254289bc588SBarry Smith   Scalar *v;
255289bc588SBarry Smith   int    i;
256289bc588SBarry Smith   *ncols = mat->n;
257289bc588SBarry Smith   if (cols) {
25878b31e54SBarry Smith     *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols);
259289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
260289bc588SBarry Smith   }
261289bc588SBarry Smith   if (vals) {
26278b31e54SBarry Smith     *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
263289bc588SBarry Smith     v = mat->v + row;
264289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
265289bc588SBarry Smith   }
266289bc588SBarry Smith   return 0;
267289bc588SBarry Smith }
26844a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols,
269289bc588SBarry Smith                             Scalar **vals)
270289bc588SBarry Smith {
27178b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
27278b31e54SBarry Smith   if (vals) { PETSCFREE(*vals); }
273289bc588SBarry Smith   return 0;
274289bc588SBarry Smith }
275289bc588SBarry Smith /* ----------------------------------------------------------------*/
27644a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n,
277e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
278289bc588SBarry Smith {
27944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
280289bc588SBarry Smith   int    i,j;
281d6dfbf8fSBarry Smith 
282289bc588SBarry Smith   if (!mat->roworiented) {
283edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
284289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
285289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
286289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
287289bc588SBarry Smith         }
288289bc588SBarry Smith       }
289289bc588SBarry Smith     }
290289bc588SBarry Smith     else {
291289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
292289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
293289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
294289bc588SBarry Smith         }
295289bc588SBarry Smith       }
296289bc588SBarry Smith     }
297e8d4e0b9SBarry Smith   }
298e8d4e0b9SBarry Smith   else {
299edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
300e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
301e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
302e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
303e8d4e0b9SBarry Smith         }
304e8d4e0b9SBarry Smith       }
305e8d4e0b9SBarry Smith     }
306289bc588SBarry Smith     else {
307289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
308289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
309289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
310289bc588SBarry Smith         }
311289bc588SBarry Smith       }
312289bc588SBarry Smith     }
313e8d4e0b9SBarry Smith   }
314289bc588SBarry Smith   return 0;
315289bc588SBarry Smith }
316e8d4e0b9SBarry Smith 
317289bc588SBarry Smith /* -----------------------------------------------------------------*/
318ff7509f2SBarry Smith static int MatCopyPrivate_Dense(Mat matin,Mat *newmat)
319289bc588SBarry Smith {
32044a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
321289bc588SBarry Smith   int ierr;
322289bc588SBarry Smith   Mat newi;
32344a69424SLois Curfman McInnes   Mat_Dense *l;
3246b5873e3SBarry Smith   if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi)))
32578b31e54SBarry Smith                                                           SETERRQ(ierr,0);
32644a69424SLois Curfman McInnes   l = (Mat_Dense *) newi->data;
32778b31e54SBarry Smith   PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
328289bc588SBarry Smith   *newmat = newi;
329289bc588SBarry Smith   return 0;
330289bc588SBarry Smith }
331da3a660dSBarry Smith #include "viewer.h"
332289bc588SBarry Smith 
33344a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr)
334289bc588SBarry Smith {
335289bc588SBarry Smith   Mat         matin = (Mat) obj;
33644a69424SLois Curfman McInnes   Mat_Dense   *mat = (Mat_Dense *) matin->data;
337289bc588SBarry Smith   Scalar      *v;
338289bc588SBarry Smith   int         i,j;
33923242f5aSBarry Smith   PetscObject vobj = (PetscObject) ptr;
340da3a660dSBarry Smith 
34123242f5aSBarry Smith   if (ptr == 0) {
34223242f5aSBarry Smith     ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr;
34323242f5aSBarry Smith   }
34423242f5aSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
3456f75cfa4SBarry Smith     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
346da3a660dSBarry Smith   }
347da3a660dSBarry Smith   else {
3484a0ce102SLois Curfman McInnes     FILE *fd = ViewerFileGetPointer_Private(ptr);
349289bc588SBarry Smith     for ( i=0; i<mat->m; i++ ) {
350289bc588SBarry Smith       v = mat->v + i;
351289bc588SBarry Smith       for ( j=0; j<mat->n; j++ ) {
352289bc588SBarry Smith #if defined(PETSC_COMPLEX)
3538e37dc09SBarry Smith         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
354289bc588SBarry Smith #else
3558e37dc09SBarry Smith         fprintf(fd,"%6.4e ",*v); v += mat->m;
356289bc588SBarry Smith #endif
357289bc588SBarry Smith       }
3588e37dc09SBarry Smith       fprintf(fd,"\n");
359289bc588SBarry Smith     }
360da3a660dSBarry Smith   }
361289bc588SBarry Smith   return 0;
362289bc588SBarry Smith }
363289bc588SBarry Smith 
364289bc588SBarry Smith 
36544a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj)
366289bc588SBarry Smith {
367289bc588SBarry Smith   Mat    mat = (Mat) obj;
36844a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) mat->data;
369a5a9c739SBarry Smith #if defined(PETSC_LOG)
370a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
371a5a9c739SBarry Smith #endif
37278b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
37378b31e54SBarry Smith   PETSCFREE(l);
374a5a9c739SBarry Smith   PLogObjectDestroy(mat);
3759e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
376289bc588SBarry Smith   return 0;
377289bc588SBarry Smith }
378289bc588SBarry Smith 
37944a69424SLois Curfman McInnes static int MatTrans_Dense(Mat matin)
380289bc588SBarry Smith {
38144a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
382289bc588SBarry Smith   int    k,j;
383289bc588SBarry Smith   Scalar *v = mat->v, tmp;
384289bc588SBarry Smith   if (mat->m != mat->n) {
38578b31e54SBarry Smith     SETERRQ(1,"Cannot transpose rectangular dense matrix");
386289bc588SBarry Smith   }
387289bc588SBarry Smith   for ( j=0; j<mat->m; j++ ) {
388289bc588SBarry Smith     for ( k=0; k<j; k++ ) {
389289bc588SBarry Smith       tmp = v[j + k*mat->n];
390289bc588SBarry Smith       v[j + k*mat->n] = v[k + j*mat->n];
391289bc588SBarry Smith       v[k + j*mat->n] = tmp;
392289bc588SBarry Smith     }
393289bc588SBarry Smith   }
394289bc588SBarry Smith   return 0;
395289bc588SBarry Smith }
396289bc588SBarry Smith 
39744a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2)
398289bc588SBarry Smith {
39944a69424SLois Curfman McInnes   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
40044a69424SLois Curfman McInnes   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
401289bc588SBarry Smith   int    i;
402289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
403289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
404289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
405289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
406289bc588SBarry Smith     if (*v1 != *v2) return 0;
407289bc588SBarry Smith     v1++; v2++;
408289bc588SBarry Smith   }
409289bc588SBarry Smith   return 1;
410289bc588SBarry Smith }
411289bc588SBarry Smith 
41246fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v)
413289bc588SBarry Smith {
41444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
4156abc6512SBarry Smith   int       i, n;
4166abc6512SBarry Smith   Scalar    *x;
417289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
41878b31e54SBarry Smith   if (n != mat->m) SETERRQ(1,"Nonconforming matrix and vector");
419289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
420289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
421289bc588SBarry Smith   }
422289bc588SBarry Smith   return 0;
423289bc588SBarry Smith }
424289bc588SBarry Smith 
42544a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
426289bc588SBarry Smith {
42744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
428da3a660dSBarry Smith   Scalar *l,*r,x,*v;
429da3a660dSBarry Smith   int    i,j,m = mat->m, n = mat->n;
43028988994SBarry Smith   if (ll) {
431da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
43278b31e54SBarry Smith     if (m != mat->m) SETERRQ(1,"Left scaling vector wrong length");
433da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
434da3a660dSBarry Smith       x = l[i];
435da3a660dSBarry Smith       v = mat->v + i;
436da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
437da3a660dSBarry Smith     }
438da3a660dSBarry Smith   }
43928988994SBarry Smith   if (rr) {
440da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
44178b31e54SBarry Smith     if (n != mat->n) SETERRQ(1,"Right scaling vector wrong length");
442da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
443da3a660dSBarry Smith       x = r[i];
444da3a660dSBarry Smith       v = mat->v + i*m;
445da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
446da3a660dSBarry Smith     }
447da3a660dSBarry Smith   }
448289bc588SBarry Smith   return 0;
449289bc588SBarry Smith }
450289bc588SBarry Smith 
451da3a660dSBarry Smith 
45244a69424SLois Curfman McInnes static int MatNorm_Dense(Mat matin,int type,double *norm)
453289bc588SBarry Smith {
45444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
455289bc588SBarry Smith   Scalar *v = mat->v;
456289bc588SBarry Smith   double sum = 0.0;
457289bc588SBarry Smith   int    i, j;
458289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
459289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
460289bc588SBarry Smith #if defined(PETSC_COMPLEX)
461289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
462289bc588SBarry Smith #else
463289bc588SBarry Smith       sum += (*v)*(*v); v++;
464289bc588SBarry Smith #endif
465289bc588SBarry Smith     }
466289bc588SBarry Smith     *norm = sqrt(sum);
467289bc588SBarry Smith   }
468289bc588SBarry Smith   else if (type == NORM_1) {
469289bc588SBarry Smith     *norm = 0.0;
470289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
471289bc588SBarry Smith       sum = 0.0;
472289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
473289bc588SBarry Smith #if defined(PETSC_COMPLEX)
474289bc588SBarry Smith         sum += abs(*v++);
475289bc588SBarry Smith #else
476289bc588SBarry Smith         sum += fabs(*v++);
477289bc588SBarry Smith #endif
478289bc588SBarry Smith       }
479289bc588SBarry Smith       if (sum > *norm) *norm = sum;
480289bc588SBarry Smith     }
481289bc588SBarry Smith   }
482289bc588SBarry Smith   else if (type == NORM_INFINITY) {
483289bc588SBarry Smith     *norm = 0.0;
484289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
485289bc588SBarry Smith       v = mat->v + j;
486289bc588SBarry Smith       sum = 0.0;
487289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
488289bc588SBarry Smith #if defined(PETSC_COMPLEX)
489289bc588SBarry Smith         sum += abs(*v); v += mat->m;
490289bc588SBarry Smith #else
491289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
492289bc588SBarry Smith #endif
493289bc588SBarry Smith       }
494289bc588SBarry Smith       if (sum > *norm) *norm = sum;
495289bc588SBarry Smith     }
496289bc588SBarry Smith   }
497289bc588SBarry Smith   else {
49878b31e54SBarry Smith     SETERRQ(1,"No support for the two norm yet");
499289bc588SBarry Smith   }
500289bc588SBarry Smith   return 0;
501289bc588SBarry Smith }
502289bc588SBarry Smith 
50320e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op)
504289bc588SBarry Smith {
50544a69424SLois Curfman McInnes   Mat_Dense *aij = (Mat_Dense *) aijin->data;
506289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
507289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
508289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
509289bc588SBarry Smith   return 0;
510289bc588SBarry Smith }
511289bc588SBarry Smith 
51244a69424SLois Curfman McInnes static int MatZero_Dense(Mat A)
5136f0a148fSBarry Smith {
51444a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
51578b31e54SBarry Smith   PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
5166f0a148fSBarry Smith   return 0;
5176f0a148fSBarry Smith }
5186f0a148fSBarry Smith 
51944a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag)
5206f0a148fSBarry Smith {
52144a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
5226abc6512SBarry Smith   int    n = l->n, i, j,ierr,N, *rows;
5236f0a148fSBarry Smith   Scalar *slot;
52478b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
52578b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5266f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5276f0a148fSBarry Smith     slot = l->v + rows[i];
5286f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5296f0a148fSBarry Smith   }
5306f0a148fSBarry Smith   if (diag) {
5316f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5326f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
533260925b8SBarry Smith       *slot = *diag;
5346f0a148fSBarry Smith     }
5356f0a148fSBarry Smith   }
536260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5376f0a148fSBarry Smith   return 0;
5386f0a148fSBarry Smith }
539557bce09SLois Curfman McInnes 
540fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n)
541557bce09SLois Curfman McInnes {
54244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
543557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
544557bce09SLois Curfman McInnes   return 0;
545557bce09SLois Curfman McInnes }
546557bce09SLois Curfman McInnes 
54744a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array)
54864e87e97SBarry Smith {
54944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
55064e87e97SBarry Smith   *array = mat->v;
55164e87e97SBarry Smith   return 0;
55264e87e97SBarry Smith }
5530754003eSLois Curfman McInnes 
5540754003eSLois Curfman McInnes 
5550754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol)
5560754003eSLois Curfman McInnes {
55778b31e54SBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_Dense not yet done.");
5580754003eSLois Curfman McInnes }
5590754003eSLois Curfman McInnes 
5600754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat)
5610754003eSLois Curfman McInnes {
5620754003eSLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
5630754003eSLois Curfman McInnes   int     nznew, *smap, i, j, ierr, oldcols = mat->n;
564160018dcSBarry Smith   int     *irow, *icol, nrows, ncols, *cwork;
5650754003eSLois Curfman McInnes   Scalar  *vwork, *val;
5660754003eSLois Curfman McInnes   Mat     newmat;
5670754003eSLois Curfman McInnes 
56878b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
56978b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
57078b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
57178b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
5720754003eSLois Curfman McInnes 
57378b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
57478b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
57578b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
57678b31e54SBarry Smith   PETSCMEMSET((char*)smap,0,oldcols*sizeof(int));
5770754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
5780754003eSLois Curfman McInnes 
5790754003eSLois Curfman McInnes   /* Create and fill new matrix */
5800754003eSLois Curfman McInnes   ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat);
58178b31e54SBarry Smith          CHKERRQ(ierr);
5820754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
5830754003eSLois Curfman McInnes     nznew = 0;
5840754003eSLois Curfman McInnes     val   = mat->v + irow[i];
5850754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
5860754003eSLois Curfman McInnes       if (smap[j]) {
5870754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
5880754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
5890754003eSLois Curfman McInnes       }
5900754003eSLois Curfman McInnes     }
591edae2e7dSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES);
59278b31e54SBarry Smith            CHKERRQ(ierr);
5930754003eSLois Curfman McInnes   }
59478b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
59578b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
5960754003eSLois Curfman McInnes 
5970754003eSLois Curfman McInnes   /* Free work space */
59878b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
59978b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
60078b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
6010754003eSLois Curfman McInnes   *submat = newmat;
6020754003eSLois Curfman McInnes   return 0;
6030754003eSLois Curfman McInnes }
6040754003eSLois Curfman McInnes 
605289bc588SBarry Smith /* -------------------------------------------------------------------*/
60644a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense,
60744a69424SLois Curfman McInnes        MatGetRow_Dense, MatRestoreRow_Dense,
60844a69424SLois Curfman McInnes        MatMult_Dense, MatMultAdd_Dense,
60944a69424SLois Curfman McInnes        MatMultTrans_Dense, MatMultTransAdd_Dense,
61044a69424SLois Curfman McInnes        MatSolve_Dense,MatSolveAdd_Dense,
61144a69424SLois Curfman McInnes        MatSolveTrans_Dense,MatSolveTransAdd_Dense,
61246fc4a8fSLois Curfman McInnes        MatLUFactor_Dense,MatCholeskyFactor_Dense,
61344a69424SLois Curfman McInnes        MatRelax_Dense,
61444a69424SLois Curfman McInnes        MatTrans_Dense,
615fa9ec3f1SLois Curfman McInnes        MatGetInfo_Dense,MatEqual_Dense,
61646fc4a8fSLois Curfman McInnes        MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense,
617289bc588SBarry Smith        0,0,
618681d1451SLois Curfman McInnes        0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0,
61944a69424SLois Curfman McInnes        MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense,
62046fc4a8fSLois Curfman McInnes        MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense,
621fa9ec3f1SLois Curfman McInnes        MatGetSize_Dense,MatGetSize_Dense,0,
62207a0e7adSLois Curfman McInnes        0,0,MatGetArray_Dense,0,0,
62307a0e7adSLois Curfman McInnes        MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense,
624ff7509f2SBarry Smith        MatCopyPrivate_Dense};
6250754003eSLois Curfman McInnes 
62620563c6bSBarry Smith /*@
62720563c6bSBarry Smith    MatCreateSequentialDense - Creates a sequential dense matrix that
628d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
629d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
630289bc588SBarry Smith 
63120563c6bSBarry Smith    Input Parameters:
6320c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
6330c775827SLois Curfman McInnes .  m - number of rows
6340c775827SLois Curfman McInnes .  n - number of column
63520563c6bSBarry Smith 
63620563c6bSBarry Smith    Output Parameter:
6370c775827SLois Curfman McInnes .  newmat - the matrix
63820563c6bSBarry Smith 
639dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
640d65003e9SLois Curfman McInnes 
641dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
64220563c6bSBarry Smith @*/
6436b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat)
644289bc588SBarry Smith {
64544a69424SLois Curfman McInnes   int       size = sizeof(Mat_Dense) + m*n*sizeof(Scalar);
646289bc588SBarry Smith   Mat mat;
64744a69424SLois Curfman McInnes   Mat_Dense    *l;
648289bc588SBarry Smith   *newmat        = 0;
6496b5873e3SBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm);
650a5a9c739SBarry Smith   PLogObjectCreate(mat);
65178b31e54SBarry Smith   l              = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l);
652289bc588SBarry Smith   mat->ops       = &MatOps;
65344a69424SLois Curfman McInnes   mat->destroy   = MatDestroy_Dense;
65444a69424SLois Curfman McInnes   mat->view      = MatView_Dense;
655289bc588SBarry Smith   mat->data      = (void *) l;
656289bc588SBarry Smith   mat->factor    = 0;
657d6dfbf8fSBarry Smith 
658289bc588SBarry Smith   l->m           = m;
659289bc588SBarry Smith   l->n           = n;
660289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
661289bc588SBarry Smith   l->pivots      = 0;
662289bc588SBarry Smith   l->roworiented = 1;
663d6dfbf8fSBarry Smith 
66478b31e54SBarry Smith   PETSCMEMSET(l->v,0,m*n*sizeof(Scalar));
665289bc588SBarry Smith   *newmat = mat;
666289bc588SBarry Smith   return 0;
667289bc588SBarry Smith }
668289bc588SBarry Smith 
66944a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat)
670289bc588SBarry Smith {
67144a69424SLois Curfman McInnes   Mat_Dense *m = (Mat_Dense *) matin->data;
6726b5873e3SBarry Smith   return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat);
673289bc588SBarry Smith }
674