xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 78b31e5415da8002fdf1e4c89af509906889cf86)
1cb512458SBarry Smith #ifndef lint
2*78b31e54SBarry Smith static char vcid[] = "$Id: dense.c,v 1.39 1995/05/29 00:02:17 bsmith Exp bsmith $";
3cb512458SBarry Smith #endif
4289bc588SBarry Smith 
5289bc588SBarry Smith /*
6289bc588SBarry Smith     Standard Fortran style matrices
7289bc588SBarry Smith */
8289bc588SBarry Smith #include "ptscimpl.h"
9289bc588SBarry Smith #include "plapack.h"
10289bc588SBarry Smith #include "matimpl.h"
11289bc588SBarry Smith #include "math.h"
12289bc588SBarry Smith #include "vec/vecimpl.h"
1320e1a0d4SLois Curfman McInnes #include "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++;}
29fa9ec3f1SLois Curfman McInnes   *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar);
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.*/
3644a69424SLois Curfman McInnes static int MatLUFactor_Dense(Mat matin,IS row,IS col)
37289bc588SBarry Smith {
3844a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
396abc6512SBarry Smith   int    info;
40289bc588SBarry Smith   if (!mat->pivots) {
41*78b31e54SBarry Smith     mat->pivots = (int *) PETSCMALLOC( mat->m*sizeof(int) );
42*78b31e54SBarry Smith     CHKPTRQ(mat->pivots);
43289bc588SBarry Smith   }
44289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
45*78b31e54SBarry Smith   if (info) SETERRQ(1,"Bad LU factorization");
46289bc588SBarry Smith   matin->factor = FACTOR_LU;
47289bc588SBarry Smith   return 0;
48289bc588SBarry Smith }
4944a69424SLois Curfman McInnes static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact)
50289bc588SBarry Smith {
51289bc588SBarry Smith   int ierr;
52*78b31e54SBarry Smith   if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0);
53289bc588SBarry Smith   return 0;
54289bc588SBarry Smith }
5544a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact)
56289bc588SBarry Smith {
5720563c6bSBarry Smith   return MatLUFactor(*fact,0,0);
58289bc588SBarry Smith }
5946fc4a8fSLois Curfman McInnes static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,Mat *fact)
60289bc588SBarry Smith {
61289bc588SBarry Smith   int ierr;
62*78b31e54SBarry Smith   if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0);
63289bc588SBarry Smith   return 0;
64289bc588SBarry Smith }
6546fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact)
66289bc588SBarry Smith {
6720563c6bSBarry Smith   return MatCholeskyFactor(*fact,0);
68289bc588SBarry Smith }
6946fc4a8fSLois Curfman McInnes static int MatCholeskyFactor_Dense(Mat matin,IS perm)
70289bc588SBarry Smith {
7144a69424SLois Curfman McInnes   Mat_Dense    *mat = (Mat_Dense *) matin->data;
726abc6512SBarry Smith   int       info;
73*78b31e54SBarry Smith   if (mat->pivots) {PETSCFREE(mat->pivots); mat->pivots = 0;}
74289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
75*78b31e54SBarry Smith   if (info) SETERRQ(1,"Bad Cholesky factorization");
76289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
77289bc588SBarry Smith   return 0;
78289bc588SBarry Smith }
79289bc588SBarry Smith 
8044a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy)
81289bc588SBarry Smith {
8244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
836abc6512SBarry Smith   int    one = 1, info;
846abc6512SBarry Smith   Scalar *x, *y;
85289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
86*78b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
879e25ed09SBarry Smith   if (matin->factor == FACTOR_LU) {
88289bc588SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
89289bc588SBarry Smith               y, &mat->m, &info );
90289bc588SBarry Smith   }
919e25ed09SBarry Smith   else if (matin->factor == FACTOR_CHOLESKY){
92289bc588SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
93289bc588SBarry Smith               y, &mat->m, &info );
94289bc588SBarry Smith   }
95*78b31e54SBarry Smith   else SETERRQ(1,"Matrix must be factored to solve");
96*78b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
97289bc588SBarry Smith   return 0;
98289bc588SBarry Smith }
9944a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy)
100da3a660dSBarry Smith {
10144a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
1026abc6512SBarry Smith   int    one = 1, info;
1036abc6512SBarry Smith   Scalar *x, *y;
104da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
105*78b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
106da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
107da3a660dSBarry Smith   if (mat->pivots) {
108da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
109da3a660dSBarry Smith               y, &mat->m, &info );
110da3a660dSBarry Smith   }
111da3a660dSBarry Smith   else {
112da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
113da3a660dSBarry Smith               y, &mat->m, &info );
114da3a660dSBarry Smith   }
115*78b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
116da3a660dSBarry Smith   return 0;
117da3a660dSBarry Smith }
11844a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
119da3a660dSBarry Smith {
12044a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
1216abc6512SBarry Smith   int    one = 1, info,ierr;
1226abc6512SBarry Smith   Scalar *x, *y, sone = 1.0;
123da3a660dSBarry Smith   Vec    tmp = 0;
124da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
125da3a660dSBarry Smith   if (yy == zz) {
126*78b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
127*78b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
128da3a660dSBarry Smith   }
129*78b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
130da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
131da3a660dSBarry Smith   if (mat->pivots) {
132da3a660dSBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
133da3a660dSBarry Smith               y, &mat->m, &info );
134da3a660dSBarry Smith   }
135da3a660dSBarry Smith   else {
136da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
137da3a660dSBarry Smith               y, &mat->m, &info );
138da3a660dSBarry Smith   }
139*78b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
140da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
141da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
142da3a660dSBarry Smith   return 0;
143da3a660dSBarry Smith }
14444a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy)
145da3a660dSBarry Smith {
14644a69424SLois Curfman McInnes   Mat_Dense  *mat = (Mat_Dense *) matin->data;
1476abc6512SBarry Smith   int     one = 1, info,ierr;
1486abc6512SBarry Smith   Scalar  *x, *y, sone = 1.0;
149da3a660dSBarry Smith   Vec     tmp;
150da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
151da3a660dSBarry Smith   if (yy == zz) {
152*78b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
153*78b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
154da3a660dSBarry Smith   }
155*78b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
156da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
157da3a660dSBarry Smith   if (mat->pivots) {
158da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
159da3a660dSBarry Smith               y, &mat->m, &info );
160da3a660dSBarry Smith   }
161da3a660dSBarry Smith   else {
162da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
163da3a660dSBarry Smith               y, &mat->m, &info );
164da3a660dSBarry Smith   }
165*78b31e54SBarry Smith   if (info) SETERRQ(1,"Bad solve");
166da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
167da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
168da3a660dSBarry Smith   return 0;
169da3a660dSBarry Smith }
170289bc588SBarry Smith /* ------------------------------------------------------------------*/
17120e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag,
17220e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
173289bc588SBarry Smith {
17444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
175289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
1766abc6512SBarry Smith   int    o = 1,ierr, m = mat->m, i;
177289bc588SBarry Smith 
178289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
179289bc588SBarry Smith     /* this is a hack fix, should have another version without
180289bc588SBarry Smith        the second BLdot */
181*78b31e54SBarry Smith     if ((ierr = VecSet(&zero,xx))) SETERRQ(ierr,0);
182289bc588SBarry Smith   }
183289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
184289bc588SBarry Smith   while (its--) {
185289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
186289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
187289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
188d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
189289bc588SBarry Smith       }
190289bc588SBarry Smith     }
191289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
192289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
193289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
194d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
195289bc588SBarry Smith       }
196289bc588SBarry Smith     }
197289bc588SBarry Smith   }
198289bc588SBarry Smith   return 0;
199289bc588SBarry Smith }
200289bc588SBarry Smith 
201289bc588SBarry Smith /* -----------------------------------------------------------------*/
20244a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy)
203289bc588SBarry Smith {
20444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
205289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
206289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
207289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
208289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
209289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
210289bc588SBarry Smith   return 0;
211289bc588SBarry Smith }
21244a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy)
213289bc588SBarry Smith {
21444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
215289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
216289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
217289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
218289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
219289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
220289bc588SBarry Smith   return 0;
221289bc588SBarry Smith }
22244a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
223289bc588SBarry Smith {
22444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
225289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
2266abc6512SBarry Smith   int    _One=1; Scalar _DOne=1.0;
227289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
228*78b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
229289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
230289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
231289bc588SBarry Smith   return 0;
232289bc588SBarry Smith }
23344a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
234289bc588SBarry Smith {
23544a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
236289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
237289bc588SBarry Smith   int    _One=1;
2386abc6512SBarry Smith   Scalar _DOne=1.0;
239289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
240289bc588SBarry Smith   VecGetArray(zz,&z);
241*78b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
242289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
243289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
244289bc588SBarry Smith   return 0;
245289bc588SBarry Smith }
246289bc588SBarry Smith 
247289bc588SBarry Smith /* -----------------------------------------------------------------*/
24844a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols,
249289bc588SBarry Smith                         Scalar **vals)
250289bc588SBarry Smith {
25144a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
252289bc588SBarry Smith   Scalar *v;
253289bc588SBarry Smith   int    i;
254289bc588SBarry Smith   *ncols = mat->n;
255289bc588SBarry Smith   if (cols) {
256*78b31e54SBarry Smith     *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols);
257289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
258289bc588SBarry Smith   }
259289bc588SBarry Smith   if (vals) {
260*78b31e54SBarry Smith     *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
261289bc588SBarry Smith     v = mat->v + row;
262289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
263289bc588SBarry Smith   }
264289bc588SBarry Smith   return 0;
265289bc588SBarry Smith }
26644a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols,
267289bc588SBarry Smith                             Scalar **vals)
268289bc588SBarry Smith {
269*78b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
270*78b31e54SBarry Smith   if (vals) { PETSCFREE(*vals); }
271289bc588SBarry Smith   return 0;
272289bc588SBarry Smith }
273289bc588SBarry Smith /* ----------------------------------------------------------------*/
27444a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n,
275e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
276289bc588SBarry Smith {
27744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
278289bc588SBarry Smith   int    i,j;
279d6dfbf8fSBarry Smith 
280289bc588SBarry Smith   if (!mat->roworiented) {
281edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
282289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
283289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
284289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
285289bc588SBarry Smith         }
286289bc588SBarry Smith       }
287289bc588SBarry Smith     }
288289bc588SBarry Smith     else {
289289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
290289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
291289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
292289bc588SBarry Smith         }
293289bc588SBarry Smith       }
294289bc588SBarry Smith     }
295e8d4e0b9SBarry Smith   }
296e8d4e0b9SBarry Smith   else {
297edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
298e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
299e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
300e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
301e8d4e0b9SBarry Smith         }
302e8d4e0b9SBarry Smith       }
303e8d4e0b9SBarry Smith     }
304289bc588SBarry Smith     else {
305289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
306289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
307289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
308289bc588SBarry Smith         }
309289bc588SBarry Smith       }
310289bc588SBarry Smith     }
311e8d4e0b9SBarry Smith   }
312289bc588SBarry Smith   return 0;
313289bc588SBarry Smith }
314e8d4e0b9SBarry Smith 
315289bc588SBarry Smith /* -----------------------------------------------------------------*/
316ff7509f2SBarry Smith static int MatCopyPrivate_Dense(Mat matin,Mat *newmat)
317289bc588SBarry Smith {
31844a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
319289bc588SBarry Smith   int ierr;
320289bc588SBarry Smith   Mat newi;
32144a69424SLois Curfman McInnes   Mat_Dense *l;
3226b5873e3SBarry Smith   if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi)))
323*78b31e54SBarry Smith                                                           SETERRQ(ierr,0);
32444a69424SLois Curfman McInnes   l = (Mat_Dense *) newi->data;
325*78b31e54SBarry Smith   PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
326289bc588SBarry Smith   *newmat = newi;
327289bc588SBarry Smith   return 0;
328289bc588SBarry Smith }
329da3a660dSBarry Smith #include "viewer.h"
330289bc588SBarry Smith 
33144a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr)
332289bc588SBarry Smith {
333289bc588SBarry Smith   Mat         matin = (Mat) obj;
33444a69424SLois Curfman McInnes   Mat_Dense   *mat = (Mat_Dense *) matin->data;
335289bc588SBarry Smith   Scalar      *v;
336289bc588SBarry Smith   int         i,j;
33723242f5aSBarry Smith   PetscObject vobj = (PetscObject) ptr;
338da3a660dSBarry Smith 
33923242f5aSBarry Smith   if (ptr == 0) {
34023242f5aSBarry Smith     ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr;
34123242f5aSBarry Smith   }
34223242f5aSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
3436f75cfa4SBarry Smith     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
344da3a660dSBarry Smith   }
345da3a660dSBarry Smith   else {
3464a0ce102SLois Curfman McInnes     FILE *fd = ViewerFileGetPointer_Private(ptr);
347289bc588SBarry Smith     for ( i=0; i<mat->m; i++ ) {
348289bc588SBarry Smith       v = mat->v + i;
349289bc588SBarry Smith       for ( j=0; j<mat->n; j++ ) {
350289bc588SBarry Smith #if defined(PETSC_COMPLEX)
3518e37dc09SBarry Smith         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
352289bc588SBarry Smith #else
3538e37dc09SBarry Smith         fprintf(fd,"%6.4e ",*v); v += mat->m;
354289bc588SBarry Smith #endif
355289bc588SBarry Smith       }
3568e37dc09SBarry Smith       fprintf(fd,"\n");
357289bc588SBarry Smith     }
358da3a660dSBarry Smith   }
359289bc588SBarry Smith   return 0;
360289bc588SBarry Smith }
361289bc588SBarry Smith 
362289bc588SBarry Smith 
36344a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj)
364289bc588SBarry Smith {
365289bc588SBarry Smith   Mat    mat = (Mat) obj;
36644a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) mat->data;
367a5a9c739SBarry Smith #if defined(PETSC_LOG)
368a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
369a5a9c739SBarry Smith #endif
370*78b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
371*78b31e54SBarry Smith   PETSCFREE(l);
372a5a9c739SBarry Smith   PLogObjectDestroy(mat);
3739e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
374289bc588SBarry Smith   return 0;
375289bc588SBarry Smith }
376289bc588SBarry Smith 
37744a69424SLois Curfman McInnes static int MatTrans_Dense(Mat matin)
378289bc588SBarry Smith {
37944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
380289bc588SBarry Smith   int    k,j;
381289bc588SBarry Smith   Scalar *v = mat->v, tmp;
382289bc588SBarry Smith   if (mat->m != mat->n) {
383*78b31e54SBarry Smith     SETERRQ(1,"Cannot transpose rectangular dense matrix");
384289bc588SBarry Smith   }
385289bc588SBarry Smith   for ( j=0; j<mat->m; j++ ) {
386289bc588SBarry Smith     for ( k=0; k<j; k++ ) {
387289bc588SBarry Smith       tmp = v[j + k*mat->n];
388289bc588SBarry Smith       v[j + k*mat->n] = v[k + j*mat->n];
389289bc588SBarry Smith       v[k + j*mat->n] = tmp;
390289bc588SBarry Smith     }
391289bc588SBarry Smith   }
392289bc588SBarry Smith   return 0;
393289bc588SBarry Smith }
394289bc588SBarry Smith 
39544a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2)
396289bc588SBarry Smith {
39744a69424SLois Curfman McInnes   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
39844a69424SLois Curfman McInnes   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
399289bc588SBarry Smith   int    i;
400289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
401289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
402289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
403289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
404289bc588SBarry Smith     if (*v1 != *v2) return 0;
405289bc588SBarry Smith     v1++; v2++;
406289bc588SBarry Smith   }
407289bc588SBarry Smith   return 1;
408289bc588SBarry Smith }
409289bc588SBarry Smith 
41046fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v)
411289bc588SBarry Smith {
41244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
4136abc6512SBarry Smith   int       i, n;
4146abc6512SBarry Smith   Scalar    *x;
415289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
416*78b31e54SBarry Smith   if (n != mat->m) SETERRQ(1,"Nonconforming matrix and vector");
417289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
418289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
419289bc588SBarry Smith   }
420289bc588SBarry Smith   return 0;
421289bc588SBarry Smith }
422289bc588SBarry Smith 
42344a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
424289bc588SBarry Smith {
42544a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
426da3a660dSBarry Smith   Scalar *l,*r,x,*v;
427da3a660dSBarry Smith   int    i,j,m = mat->m, n = mat->n;
42828988994SBarry Smith   if (ll) {
429da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
430*78b31e54SBarry Smith     if (m != mat->m) SETERRQ(1,"Left scaling vector wrong length");
431da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
432da3a660dSBarry Smith       x = l[i];
433da3a660dSBarry Smith       v = mat->v + i;
434da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
435da3a660dSBarry Smith     }
436da3a660dSBarry Smith   }
43728988994SBarry Smith   if (rr) {
438da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
439*78b31e54SBarry Smith     if (n != mat->n) SETERRQ(1,"Right scaling vector wrong length");
440da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
441da3a660dSBarry Smith       x = r[i];
442da3a660dSBarry Smith       v = mat->v + i*m;
443da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
444da3a660dSBarry Smith     }
445da3a660dSBarry Smith   }
446289bc588SBarry Smith   return 0;
447289bc588SBarry Smith }
448289bc588SBarry Smith 
449da3a660dSBarry Smith 
45044a69424SLois Curfman McInnes static int MatNorm_Dense(Mat matin,int type,double *norm)
451289bc588SBarry Smith {
45244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
453289bc588SBarry Smith   Scalar *v = mat->v;
454289bc588SBarry Smith   double sum = 0.0;
455289bc588SBarry Smith   int    i, j;
456289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
457289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
458289bc588SBarry Smith #if defined(PETSC_COMPLEX)
459289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
460289bc588SBarry Smith #else
461289bc588SBarry Smith       sum += (*v)*(*v); v++;
462289bc588SBarry Smith #endif
463289bc588SBarry Smith     }
464289bc588SBarry Smith     *norm = sqrt(sum);
465289bc588SBarry Smith   }
466289bc588SBarry Smith   else if (type == NORM_1) {
467289bc588SBarry Smith     *norm = 0.0;
468289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
469289bc588SBarry Smith       sum = 0.0;
470289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
471289bc588SBarry Smith #if defined(PETSC_COMPLEX)
472289bc588SBarry Smith         sum += abs(*v++);
473289bc588SBarry Smith #else
474289bc588SBarry Smith         sum += fabs(*v++);
475289bc588SBarry Smith #endif
476289bc588SBarry Smith       }
477289bc588SBarry Smith       if (sum > *norm) *norm = sum;
478289bc588SBarry Smith     }
479289bc588SBarry Smith   }
480289bc588SBarry Smith   else if (type == NORM_INFINITY) {
481289bc588SBarry Smith     *norm = 0.0;
482289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
483289bc588SBarry Smith       v = mat->v + j;
484289bc588SBarry Smith       sum = 0.0;
485289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
486289bc588SBarry Smith #if defined(PETSC_COMPLEX)
487289bc588SBarry Smith         sum += abs(*v); v += mat->m;
488289bc588SBarry Smith #else
489289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
490289bc588SBarry Smith #endif
491289bc588SBarry Smith       }
492289bc588SBarry Smith       if (sum > *norm) *norm = sum;
493289bc588SBarry Smith     }
494289bc588SBarry Smith   }
495289bc588SBarry Smith   else {
496*78b31e54SBarry Smith     SETERRQ(1,"No support for the two norm yet");
497289bc588SBarry Smith   }
498289bc588SBarry Smith   return 0;
499289bc588SBarry Smith }
500289bc588SBarry Smith 
50120e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op)
502289bc588SBarry Smith {
50344a69424SLois Curfman McInnes   Mat_Dense *aij = (Mat_Dense *) aijin->data;
504289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
505289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
506289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
507289bc588SBarry Smith   return 0;
508289bc588SBarry Smith }
509289bc588SBarry Smith 
51044a69424SLois Curfman McInnes static int MatZero_Dense(Mat A)
5116f0a148fSBarry Smith {
51244a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
513*78b31e54SBarry Smith   PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
5146f0a148fSBarry Smith   return 0;
5156f0a148fSBarry Smith }
5166f0a148fSBarry Smith 
51744a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag)
5186f0a148fSBarry Smith {
51944a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
5206abc6512SBarry Smith   int    n = l->n, i, j,ierr,N, *rows;
5216f0a148fSBarry Smith   Scalar *slot;
522*78b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
523*78b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5246f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5256f0a148fSBarry Smith     slot = l->v + rows[i];
5266f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5276f0a148fSBarry Smith   }
5286f0a148fSBarry Smith   if (diag) {
5296f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5306f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
531260925b8SBarry Smith       *slot = *diag;
5326f0a148fSBarry Smith     }
5336f0a148fSBarry Smith   }
534260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5356f0a148fSBarry Smith   return 0;
5366f0a148fSBarry Smith }
537557bce09SLois Curfman McInnes 
538fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n)
539557bce09SLois Curfman McInnes {
54044a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
541557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
542557bce09SLois Curfman McInnes   return 0;
543557bce09SLois Curfman McInnes }
544557bce09SLois Curfman McInnes 
54544a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array)
54664e87e97SBarry Smith {
54744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
54864e87e97SBarry Smith   *array = mat->v;
54964e87e97SBarry Smith   return 0;
55064e87e97SBarry Smith }
5510754003eSLois Curfman McInnes 
5520754003eSLois Curfman McInnes 
5530754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol)
5540754003eSLois Curfman McInnes {
555*78b31e54SBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_Dense not yet done.");
5560754003eSLois Curfman McInnes }
5570754003eSLois Curfman McInnes 
5580754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat)
5590754003eSLois Curfman McInnes {
5600754003eSLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
5610754003eSLois Curfman McInnes   int     nznew, *smap, i, j, ierr, oldcols = mat->n;
562160018dcSBarry Smith   int     *irow, *icol, nrows, ncols, *cwork;
5630754003eSLois Curfman McInnes   Scalar  *vwork, *val;
5640754003eSLois Curfman McInnes   Mat     newmat;
5650754003eSLois Curfman McInnes 
566*78b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
567*78b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
568*78b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
569*78b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
5700754003eSLois Curfman McInnes 
571*78b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
572*78b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
573*78b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
574*78b31e54SBarry Smith   PETSCMEMSET((char*)smap,0,oldcols*sizeof(int));
5750754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
5760754003eSLois Curfman McInnes 
5770754003eSLois Curfman McInnes   /* Create and fill new matrix */
5780754003eSLois Curfman McInnes   ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat);
579*78b31e54SBarry Smith          CHKERRQ(ierr);
5800754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
5810754003eSLois Curfman McInnes     nznew = 0;
5820754003eSLois Curfman McInnes     val   = mat->v + irow[i];
5830754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
5840754003eSLois Curfman McInnes       if (smap[j]) {
5850754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
5860754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
5870754003eSLois Curfman McInnes       }
5880754003eSLois Curfman McInnes     }
589edae2e7dSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES);
590*78b31e54SBarry Smith            CHKERRQ(ierr);
5910754003eSLois Curfman McInnes   }
592*78b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
593*78b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
5940754003eSLois Curfman McInnes 
5950754003eSLois Curfman McInnes   /* Free work space */
596*78b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
597*78b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
598*78b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
5990754003eSLois Curfman McInnes   *submat = newmat;
6000754003eSLois Curfman McInnes   return 0;
6010754003eSLois Curfman McInnes }
6020754003eSLois Curfman McInnes 
603289bc588SBarry Smith /* -------------------------------------------------------------------*/
60444a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense,
60544a69424SLois Curfman McInnes        MatGetRow_Dense, MatRestoreRow_Dense,
60644a69424SLois Curfman McInnes        MatMult_Dense, MatMultAdd_Dense,
60744a69424SLois Curfman McInnes        MatMultTrans_Dense, MatMultTransAdd_Dense,
60844a69424SLois Curfman McInnes        MatSolve_Dense,MatSolveAdd_Dense,
60944a69424SLois Curfman McInnes        MatSolveTrans_Dense,MatSolveTransAdd_Dense,
61046fc4a8fSLois Curfman McInnes        MatLUFactor_Dense,MatCholeskyFactor_Dense,
61144a69424SLois Curfman McInnes        MatRelax_Dense,
61244a69424SLois Curfman McInnes        MatTrans_Dense,
613fa9ec3f1SLois Curfman McInnes        MatGetInfo_Dense,MatEqual_Dense,
61446fc4a8fSLois Curfman McInnes        MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense,
615289bc588SBarry Smith        0,0,
616681d1451SLois Curfman McInnes        0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0,
61744a69424SLois Curfman McInnes        MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense,
61846fc4a8fSLois Curfman McInnes        MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense,
619fa9ec3f1SLois Curfman McInnes        MatGetSize_Dense,MatGetSize_Dense,0,
62007a0e7adSLois Curfman McInnes        0,0,MatGetArray_Dense,0,0,
62107a0e7adSLois Curfman McInnes        MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense,
622ff7509f2SBarry Smith        MatCopyPrivate_Dense};
6230754003eSLois Curfman McInnes 
62420563c6bSBarry Smith /*@
62520563c6bSBarry Smith    MatCreateSequentialDense - Creates a sequential dense matrix that
626d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
627d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
628289bc588SBarry Smith 
62920563c6bSBarry Smith    Input Parameters:
6300c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
6310c775827SLois Curfman McInnes .  m - number of rows
6320c775827SLois Curfman McInnes .  n - number of column
63320563c6bSBarry Smith 
63420563c6bSBarry Smith    Output Parameter:
6350c775827SLois Curfman McInnes .  newmat - the matrix
63620563c6bSBarry Smith 
637dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
638d65003e9SLois Curfman McInnes 
639dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
64020563c6bSBarry Smith @*/
6416b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat)
642289bc588SBarry Smith {
64344a69424SLois Curfman McInnes   int       size = sizeof(Mat_Dense) + m*n*sizeof(Scalar);
644289bc588SBarry Smith   Mat mat;
64544a69424SLois Curfman McInnes   Mat_Dense    *l;
646289bc588SBarry Smith   *newmat        = 0;
6476b5873e3SBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm);
648a5a9c739SBarry Smith   PLogObjectCreate(mat);
649*78b31e54SBarry Smith   l              = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l);
650289bc588SBarry Smith   mat->ops       = &MatOps;
65144a69424SLois Curfman McInnes   mat->destroy   = MatDestroy_Dense;
65244a69424SLois Curfman McInnes   mat->view      = MatView_Dense;
653289bc588SBarry Smith   mat->data      = (void *) l;
654289bc588SBarry Smith   mat->factor    = 0;
655d6dfbf8fSBarry Smith 
656289bc588SBarry Smith   l->m           = m;
657289bc588SBarry Smith   l->n           = n;
658289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
659289bc588SBarry Smith   l->pivots      = 0;
660289bc588SBarry Smith   l->roworiented = 1;
661d6dfbf8fSBarry Smith 
662*78b31e54SBarry Smith   PETSCMEMSET(l->v,0,m*n*sizeof(Scalar));
663289bc588SBarry Smith   *newmat = mat;
664289bc588SBarry Smith   return 0;
665289bc588SBarry Smith }
666289bc588SBarry Smith 
66744a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat)
668289bc588SBarry Smith {
66944a69424SLois Curfman McInnes   Mat_Dense *m = (Mat_Dense *) matin->data;
6706b5873e3SBarry Smith   return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat);
671289bc588SBarry Smith }
672