xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 464493b3f8ca544de0565ce40d8044a414a4b8ef)
15392566eSBarry Smith 
25392566eSBarry Smith 
3cb512458SBarry Smith #ifndef lint
4*464493b3SBarry Smith static char vcid[] = "$Id: dense.c,v 1.46 1995/08/07 21:59:24 bsmith Exp bsmith $";
5cb512458SBarry Smith #endif
6289bc588SBarry Smith 
7289bc588SBarry Smith /*
8289bc588SBarry Smith     Standard Fortran style matrices
9289bc588SBarry Smith */
1019b02663SBarry Smith #include "petsc.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.*/
3849d8b64dSBarry Smith static int MatLUFactor_Dense(Mat matin,IS row,IS col,double f)
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);
47bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatLUFactor_Dense:Bad LU factorization");
48289bc588SBarry Smith   matin->factor = FACTOR_LU;
49289bc588SBarry Smith   return 0;
50289bc588SBarry Smith }
5149d8b64dSBarry Smith static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,double f,
5249d8b64dSBarry Smith                                      Mat *fact)
53289bc588SBarry Smith {
54289bc588SBarry Smith   int ierr;
55bbb6d6a8SBarry Smith   ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr);
56289bc588SBarry Smith   return 0;
57289bc588SBarry Smith }
5844a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact)
59289bc588SBarry Smith {
6049d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
61289bc588SBarry Smith }
6249d8b64dSBarry Smith static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,double f,Mat *fact)
63289bc588SBarry Smith {
64289bc588SBarry Smith   int ierr;
65bbb6d6a8SBarry Smith   ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr);
66289bc588SBarry Smith   return 0;
67289bc588SBarry Smith }
6846fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact)
69289bc588SBarry Smith {
7049d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
71289bc588SBarry Smith }
7249d8b64dSBarry Smith static int MatCholeskyFactor_Dense(Mat matin,IS perm,double f)
73289bc588SBarry Smith {
7444a69424SLois Curfman McInnes   Mat_Dense    *mat = (Mat_Dense *) matin->data;
756abc6512SBarry Smith   int       info;
7678b31e54SBarry Smith   if (mat->pivots) {PETSCFREE(mat->pivots); mat->pivots = 0;}
77289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
78bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_Dense:Bad factorization");
79289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
80289bc588SBarry Smith   return 0;
81289bc588SBarry Smith }
82289bc588SBarry Smith 
8344a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy)
84289bc588SBarry Smith {
8544a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
866abc6512SBarry Smith   int    one = 1, info;
876abc6512SBarry Smith   Scalar *x, *y;
88289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
8978b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
909e25ed09SBarry Smith   if (matin->factor == FACTOR_LU) {
91289bc588SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
92289bc588SBarry Smith               y, &mat->m, &info );
93289bc588SBarry Smith   }
949e25ed09SBarry Smith   else if (matin->factor == FACTOR_CHOLESKY){
95289bc588SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
96289bc588SBarry Smith               y, &mat->m, &info );
97289bc588SBarry Smith   }
98bbb6d6a8SBarry Smith   else SETERRQ(1,"MatSolve_Dense:Matrix must be factored to solve");
99bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolve_Dense:Bad solve");
100289bc588SBarry Smith   return 0;
101289bc588SBarry Smith }
10244a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy)
103da3a660dSBarry Smith {
10444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
1056abc6512SBarry Smith   int    one = 1, info;
1066abc6512SBarry Smith   Scalar *x, *y;
107da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
10878b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
109da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
110da3a660dSBarry Smith   if (mat->pivots) {
111da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
112da3a660dSBarry Smith               y, &mat->m, &info );
113da3a660dSBarry Smith   }
114da3a660dSBarry Smith   else {
115da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
116da3a660dSBarry Smith               y, &mat->m, &info );
117da3a660dSBarry Smith   }
118bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_Dense:Bad solve");
119da3a660dSBarry Smith   return 0;
120da3a660dSBarry Smith }
12144a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
122da3a660dSBarry Smith {
12344a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
1246abc6512SBarry Smith   int    one = 1, info,ierr;
1256abc6512SBarry Smith   Scalar *x, *y, sone = 1.0;
126da3a660dSBarry Smith   Vec    tmp = 0;
127da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
128da3a660dSBarry Smith   if (yy == zz) {
12978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
130*464493b3SBarry Smith     PLogObjectParent(matin,tmp);
13178b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
132da3a660dSBarry Smith   }
13378b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
134da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
135da3a660dSBarry Smith   if (mat->pivots) {
136da3a660dSBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
137da3a660dSBarry Smith               y, &mat->m, &info );
138da3a660dSBarry Smith   }
139da3a660dSBarry Smith   else {
140da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
141da3a660dSBarry Smith               y, &mat->m, &info );
142da3a660dSBarry Smith   }
143bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_Dense:Bad solve");
144da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
145da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
146da3a660dSBarry Smith   return 0;
147da3a660dSBarry Smith }
14844a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy)
149da3a660dSBarry Smith {
15044a69424SLois Curfman McInnes   Mat_Dense  *mat = (Mat_Dense *) matin->data;
1516abc6512SBarry Smith   int     one = 1, info,ierr;
1526abc6512SBarry Smith   Scalar  *x, *y, sone = 1.0;
153da3a660dSBarry Smith   Vec     tmp;
154da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
155da3a660dSBarry Smith   if (yy == zz) {
15678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
157*464493b3SBarry Smith     PLogObjectParent(matin,tmp);
15878b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
159da3a660dSBarry Smith   }
16078b31e54SBarry Smith   PETSCMEMCPY(y,x,mat->m*sizeof(Scalar));
161da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
162da3a660dSBarry Smith   if (mat->pivots) {
163da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
164da3a660dSBarry Smith               y, &mat->m, &info );
165da3a660dSBarry Smith   }
166da3a660dSBarry Smith   else {
167da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
168da3a660dSBarry Smith               y, &mat->m, &info );
169da3a660dSBarry Smith   }
170bbb6d6a8SBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_Dense:Bad solve");
171da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
172da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
173da3a660dSBarry Smith   return 0;
174da3a660dSBarry Smith }
175289bc588SBarry Smith /* ------------------------------------------------------------------*/
17620e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag,
17720e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
178289bc588SBarry Smith {
17944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
180289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
1816abc6512SBarry Smith   int    o = 1,ierr, m = mat->m, i;
182289bc588SBarry Smith 
183289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
184289bc588SBarry Smith     /* this is a hack fix, should have another version without
185289bc588SBarry Smith        the second BLdot */
186bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
187289bc588SBarry Smith   }
188289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
189289bc588SBarry Smith   while (its--) {
190289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
191289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
192289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
193d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
194289bc588SBarry Smith       }
195289bc588SBarry Smith     }
196289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
197289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
198289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
199d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
200289bc588SBarry Smith       }
201289bc588SBarry Smith     }
202289bc588SBarry Smith   }
203289bc588SBarry Smith   return 0;
204289bc588SBarry Smith }
205289bc588SBarry Smith 
206289bc588SBarry Smith /* -----------------------------------------------------------------*/
20744a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy)
208289bc588SBarry Smith {
20944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
210289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
211289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
212289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
213289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
214289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
215289bc588SBarry Smith   return 0;
216289bc588SBarry Smith }
21744a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy)
218289bc588SBarry Smith {
21944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
220289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
221289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
222289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
223289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
224289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
225289bc588SBarry Smith   return 0;
226289bc588SBarry Smith }
22744a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
228289bc588SBarry Smith {
22944a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
230289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
2316abc6512SBarry Smith   int    _One=1; Scalar _DOne=1.0;
232289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
23378b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
234289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
235289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
236289bc588SBarry Smith   return 0;
237289bc588SBarry Smith }
23844a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
239289bc588SBarry Smith {
24044a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
241289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
242289bc588SBarry Smith   int    _One=1;
2436abc6512SBarry Smith   Scalar _DOne=1.0;
244289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
245289bc588SBarry Smith   VecGetArray(zz,&z);
24678b31e54SBarry Smith   if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar));
247289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
248289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
249289bc588SBarry Smith   return 0;
250289bc588SBarry Smith }
251289bc588SBarry Smith 
252289bc588SBarry Smith /* -----------------------------------------------------------------*/
25344a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols,
254289bc588SBarry Smith                         Scalar **vals)
255289bc588SBarry Smith {
25644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
257289bc588SBarry Smith   Scalar *v;
258289bc588SBarry Smith   int    i;
259289bc588SBarry Smith   *ncols = mat->n;
260289bc588SBarry Smith   if (cols) {
26178b31e54SBarry Smith     *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols);
262289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
263289bc588SBarry Smith   }
264289bc588SBarry Smith   if (vals) {
26578b31e54SBarry Smith     *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
266289bc588SBarry Smith     v = mat->v + row;
267289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
268289bc588SBarry Smith   }
269289bc588SBarry Smith   return 0;
270289bc588SBarry Smith }
27144a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols,
272289bc588SBarry Smith                             Scalar **vals)
273289bc588SBarry Smith {
27478b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
27578b31e54SBarry Smith   if (vals) { PETSCFREE(*vals); }
276289bc588SBarry Smith   return 0;
277289bc588SBarry Smith }
278289bc588SBarry Smith /* ----------------------------------------------------------------*/
27944a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n,
280e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
281289bc588SBarry Smith {
28244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
283289bc588SBarry Smith   int    i,j;
284d6dfbf8fSBarry Smith 
285289bc588SBarry Smith   if (!mat->roworiented) {
286edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
287289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
288289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
289289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
290289bc588SBarry Smith         }
291289bc588SBarry Smith       }
292289bc588SBarry Smith     }
293289bc588SBarry Smith     else {
294289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
295289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
296289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
297289bc588SBarry Smith         }
298289bc588SBarry Smith       }
299289bc588SBarry Smith     }
300e8d4e0b9SBarry Smith   }
301e8d4e0b9SBarry Smith   else {
302edae2e7dSBarry Smith     if (addv == INSERTVALUES) {
303e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
304e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
305e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
306e8d4e0b9SBarry Smith         }
307e8d4e0b9SBarry Smith       }
308e8d4e0b9SBarry Smith     }
309289bc588SBarry Smith     else {
310289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
311289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
312289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
313289bc588SBarry Smith         }
314289bc588SBarry Smith       }
315289bc588SBarry Smith     }
316e8d4e0b9SBarry Smith   }
317289bc588SBarry Smith   return 0;
318289bc588SBarry Smith }
319e8d4e0b9SBarry Smith 
320289bc588SBarry Smith /* -----------------------------------------------------------------*/
321ff7509f2SBarry Smith static int MatCopyPrivate_Dense(Mat matin,Mat *newmat)
322289bc588SBarry Smith {
32344a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
324289bc588SBarry Smith   int ierr;
325289bc588SBarry Smith   Mat newi;
32644a69424SLois Curfman McInnes   Mat_Dense *l;
327bbb6d6a8SBarry Smith   ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi);
328bbb6d6a8SBarry Smith   CHKERRQ(ierr);
32944a69424SLois Curfman McInnes   l = (Mat_Dense *) newi->data;
33078b31e54SBarry Smith   PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
331289bc588SBarry Smith   *newmat = newi;
332289bc588SBarry Smith   return 0;
333289bc588SBarry Smith }
334da3a660dSBarry Smith #include "viewer.h"
335289bc588SBarry Smith 
33644a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr)
337289bc588SBarry Smith {
338289bc588SBarry Smith   Mat         matin = (Mat) obj;
33944a69424SLois Curfman McInnes   Mat_Dense   *mat = (Mat_Dense *) matin->data;
340289bc588SBarry Smith   Scalar      *v;
341289bc588SBarry Smith   int         i,j;
34223242f5aSBarry Smith   PetscObject vobj = (PetscObject) ptr;
343da3a660dSBarry Smith 
34423242f5aSBarry Smith   if (ptr == 0) {
34523242f5aSBarry Smith     ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr;
34623242f5aSBarry Smith   }
34723242f5aSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) {
3486f75cfa4SBarry Smith     return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v);
349da3a660dSBarry Smith   }
350da3a660dSBarry Smith   else {
3514a0ce102SLois Curfman McInnes     FILE *fd = ViewerFileGetPointer_Private(ptr);
352289bc588SBarry Smith     for ( i=0; i<mat->m; i++ ) {
353289bc588SBarry Smith       v = mat->v + i;
354289bc588SBarry Smith       for ( j=0; j<mat->n; j++ ) {
355289bc588SBarry Smith #if defined(PETSC_COMPLEX)
3568e37dc09SBarry Smith         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
357289bc588SBarry Smith #else
3588e37dc09SBarry Smith         fprintf(fd,"%6.4e ",*v); v += mat->m;
359289bc588SBarry Smith #endif
360289bc588SBarry Smith       }
3618e37dc09SBarry Smith       fprintf(fd,"\n");
362289bc588SBarry Smith     }
363da3a660dSBarry Smith   }
364289bc588SBarry Smith   return 0;
365289bc588SBarry Smith }
366289bc588SBarry Smith 
367289bc588SBarry Smith 
36844a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj)
369289bc588SBarry Smith {
370289bc588SBarry Smith   Mat    mat = (Mat) obj;
37144a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) mat->data;
372a5a9c739SBarry Smith #if defined(PETSC_LOG)
373a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
374a5a9c739SBarry Smith #endif
37578b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
37678b31e54SBarry Smith   PETSCFREE(l);
377a5a9c739SBarry Smith   PLogObjectDestroy(mat);
3789e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
379289bc588SBarry Smith   return 0;
380289bc588SBarry Smith }
381289bc588SBarry Smith 
38248b35521SBarry Smith static int MatTranspose_Dense(Mat matin,Mat *matout)
383289bc588SBarry Smith {
38444a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
385289bc588SBarry Smith   int    k,j;
386289bc588SBarry Smith   Scalar *v = mat->v, tmp;
38748b35521SBarry Smith 
38848b35521SBarry Smith   if (!matout) { /* in place transpose */
389289bc588SBarry Smith     if (mat->m != mat->n) {
39048b35521SBarry Smith       SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix");
391289bc588SBarry Smith     }
392289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
393289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
394289bc588SBarry Smith         tmp = v[j + k*mat->n];
395289bc588SBarry Smith         v[j + k*mat->n] = v[k + j*mat->n];
396289bc588SBarry Smith         v[k + j*mat->n] = tmp;
397289bc588SBarry Smith       }
398289bc588SBarry Smith     }
39948b35521SBarry Smith   }
40048b35521SBarry Smith   else {
40148b35521SBarry Smith     SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet");
40248b35521SBarry Smith   }
403289bc588SBarry Smith   return 0;
404289bc588SBarry Smith }
405289bc588SBarry Smith 
40644a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2)
407289bc588SBarry Smith {
40844a69424SLois Curfman McInnes   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
40944a69424SLois Curfman McInnes   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
410289bc588SBarry Smith   int    i;
411289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
412289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
413289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
414289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
415289bc588SBarry Smith     if (*v1 != *v2) return 0;
416289bc588SBarry Smith     v1++; v2++;
417289bc588SBarry Smith   }
418289bc588SBarry Smith   return 1;
419289bc588SBarry Smith }
420289bc588SBarry Smith 
42146fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v)
422289bc588SBarry Smith {
42344a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
4246abc6512SBarry Smith   int       i, n;
4256abc6512SBarry Smith   Scalar    *x;
426289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
427bbb6d6a8SBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_Dense:Nonconforming mat and vec");
428289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
429289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
430289bc588SBarry Smith   }
431289bc588SBarry Smith   return 0;
432289bc588SBarry Smith }
433289bc588SBarry Smith 
43444a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
435289bc588SBarry Smith {
43644a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
437da3a660dSBarry Smith   Scalar *l,*r,x,*v;
438da3a660dSBarry Smith   int    i,j,m = mat->m, n = mat->n;
43928988994SBarry Smith   if (ll) {
440da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
441bbb6d6a8SBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_Dense:Left scaling vec wrong size");
442da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
443da3a660dSBarry Smith       x = l[i];
444da3a660dSBarry Smith       v = mat->v + i;
445da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
446da3a660dSBarry Smith     }
447da3a660dSBarry Smith   }
44828988994SBarry Smith   if (rr) {
449da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
450bbb6d6a8SBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_Dense:Right scaling vec wrong size");
451da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
452da3a660dSBarry Smith       x = r[i];
453da3a660dSBarry Smith       v = mat->v + i*m;
454da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
455da3a660dSBarry Smith     }
456da3a660dSBarry Smith   }
457289bc588SBarry Smith   return 0;
458289bc588SBarry Smith }
459289bc588SBarry Smith 
4607c17b2cfSLois Curfman McInnes static int MatNorm_Dense(Mat matin,MatNormType type,double *norm)
461289bc588SBarry Smith {
46244a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
463289bc588SBarry Smith   Scalar *v = mat->v;
464289bc588SBarry Smith   double sum = 0.0;
465289bc588SBarry Smith   int    i, j;
466289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
467289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
468289bc588SBarry Smith #if defined(PETSC_COMPLEX)
469289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
470289bc588SBarry Smith #else
471289bc588SBarry Smith       sum += (*v)*(*v); v++;
472289bc588SBarry Smith #endif
473289bc588SBarry Smith     }
474289bc588SBarry Smith     *norm = sqrt(sum);
475289bc588SBarry Smith   }
476289bc588SBarry Smith   else if (type == NORM_1) {
477289bc588SBarry Smith     *norm = 0.0;
478289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
479289bc588SBarry Smith       sum = 0.0;
480289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
481289bc588SBarry Smith #if defined(PETSC_COMPLEX)
482289bc588SBarry Smith         sum += abs(*v++);
483289bc588SBarry Smith #else
484289bc588SBarry Smith         sum += fabs(*v++);
485289bc588SBarry Smith #endif
486289bc588SBarry Smith       }
487289bc588SBarry Smith       if (sum > *norm) *norm = sum;
488289bc588SBarry Smith     }
489289bc588SBarry Smith   }
490289bc588SBarry Smith   else if (type == NORM_INFINITY) {
491289bc588SBarry Smith     *norm = 0.0;
492289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
493289bc588SBarry Smith       v = mat->v + j;
494289bc588SBarry Smith       sum = 0.0;
495289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
496289bc588SBarry Smith #if defined(PETSC_COMPLEX)
497289bc588SBarry Smith         sum += abs(*v); v += mat->m;
498289bc588SBarry Smith #else
499289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
500289bc588SBarry Smith #endif
501289bc588SBarry Smith       }
502289bc588SBarry Smith       if (sum > *norm) *norm = sum;
503289bc588SBarry Smith     }
504289bc588SBarry Smith   }
505289bc588SBarry Smith   else {
506bbb6d6a8SBarry Smith     SETERRQ(1,"MatNorm_Dense:No two norm yet");
507289bc588SBarry Smith   }
508289bc588SBarry Smith   return 0;
509289bc588SBarry Smith }
510289bc588SBarry Smith 
51120e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op)
512289bc588SBarry Smith {
51344a69424SLois Curfman McInnes   Mat_Dense *aij = (Mat_Dense *) aijin->data;
514289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
515289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
516289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
517289bc588SBarry Smith   return 0;
518289bc588SBarry Smith }
519289bc588SBarry Smith 
52044a69424SLois Curfman McInnes static int MatZero_Dense(Mat A)
5216f0a148fSBarry Smith {
52244a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
52378b31e54SBarry Smith   PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
5246f0a148fSBarry Smith   return 0;
5256f0a148fSBarry Smith }
5266f0a148fSBarry Smith 
52744a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag)
5286f0a148fSBarry Smith {
52944a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
5306abc6512SBarry Smith   int    n = l->n, i, j,ierr,N, *rows;
5316f0a148fSBarry Smith   Scalar *slot;
53278b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
53378b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5346f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5356f0a148fSBarry Smith     slot = l->v + rows[i];
5366f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5376f0a148fSBarry Smith   }
5386f0a148fSBarry Smith   if (diag) {
5396f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5406f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
541260925b8SBarry Smith       *slot = *diag;
5426f0a148fSBarry Smith     }
5436f0a148fSBarry Smith   }
544260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5456f0a148fSBarry Smith   return 0;
5466f0a148fSBarry Smith }
547557bce09SLois Curfman McInnes 
548fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n)
549557bce09SLois Curfman McInnes {
55044a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
551557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
552557bce09SLois Curfman McInnes   return 0;
553557bce09SLois Curfman McInnes }
554557bce09SLois Curfman McInnes 
55544a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array)
55664e87e97SBarry Smith {
55744a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
55864e87e97SBarry Smith   *array = mat->v;
55964e87e97SBarry Smith   return 0;
56064e87e97SBarry Smith }
5610754003eSLois Curfman McInnes 
5620754003eSLois Curfman McInnes 
5630754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol)
5640754003eSLois Curfman McInnes {
565bbb6d6a8SBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_Dense: not done");
5660754003eSLois Curfman McInnes }
5670754003eSLois Curfman McInnes 
5680754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat)
5690754003eSLois Curfman McInnes {
5700754003eSLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
5710754003eSLois Curfman McInnes   int     nznew, *smap, i, j, ierr, oldcols = mat->n;
572160018dcSBarry Smith   int     *irow, *icol, nrows, ncols, *cwork;
5730754003eSLois Curfman McInnes   Scalar  *vwork, *val;
5740754003eSLois Curfman McInnes   Mat     newmat;
5750754003eSLois Curfman McInnes 
57678b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
57778b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
57878b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
57978b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
5800754003eSLois Curfman McInnes 
58178b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
58278b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
58378b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
58478b31e54SBarry Smith   PETSCMEMSET((char*)smap,0,oldcols*sizeof(int));
5850754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
5860754003eSLois Curfman McInnes 
5870754003eSLois Curfman McInnes   /* Create and fill new matrix */
5880754003eSLois Curfman McInnes   ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat);
58978b31e54SBarry Smith          CHKERRQ(ierr);
5900754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
5910754003eSLois Curfman McInnes     nznew = 0;
5920754003eSLois Curfman McInnes     val   = mat->v + irow[i];
5930754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
5940754003eSLois Curfman McInnes       if (smap[j]) {
5950754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
5960754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
5970754003eSLois Curfman McInnes       }
5980754003eSLois Curfman McInnes     }
599edae2e7dSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES);
60078b31e54SBarry Smith            CHKERRQ(ierr);
6010754003eSLois Curfman McInnes   }
60278b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
60378b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
6040754003eSLois Curfman McInnes 
6050754003eSLois Curfman McInnes   /* Free work space */
60678b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
60778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
60878b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
6090754003eSLois Curfman McInnes   *submat = newmat;
6100754003eSLois Curfman McInnes   return 0;
6110754003eSLois Curfman McInnes }
6120754003eSLois Curfman McInnes 
613289bc588SBarry Smith /* -------------------------------------------------------------------*/
61444a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense,
61544a69424SLois Curfman McInnes        MatGetRow_Dense, MatRestoreRow_Dense,
61644a69424SLois Curfman McInnes        MatMult_Dense, MatMultAdd_Dense,
61744a69424SLois Curfman McInnes        MatMultTrans_Dense, MatMultTransAdd_Dense,
61844a69424SLois Curfman McInnes        MatSolve_Dense,MatSolveAdd_Dense,
61944a69424SLois Curfman McInnes        MatSolveTrans_Dense,MatSolveTransAdd_Dense,
62046fc4a8fSLois Curfman McInnes        MatLUFactor_Dense,MatCholeskyFactor_Dense,
62144a69424SLois Curfman McInnes        MatRelax_Dense,
62248b35521SBarry Smith        MatTranspose_Dense,
623fa9ec3f1SLois Curfman McInnes        MatGetInfo_Dense,MatEqual_Dense,
62446fc4a8fSLois Curfman McInnes        MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense,
625289bc588SBarry Smith        0,0,
626681d1451SLois Curfman McInnes        0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0,
62744a69424SLois Curfman McInnes        MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense,
62846fc4a8fSLois Curfman McInnes        MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense,
629fa9ec3f1SLois Curfman McInnes        MatGetSize_Dense,MatGetSize_Dense,0,
63007a0e7adSLois Curfman McInnes        0,0,MatGetArray_Dense,0,0,
63107a0e7adSLois Curfman McInnes        MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense,
632ff7509f2SBarry Smith        MatCopyPrivate_Dense};
6330754003eSLois Curfman McInnes 
63420563c6bSBarry Smith /*@
63520563c6bSBarry Smith    MatCreateSequentialDense - Creates a sequential dense matrix that
636d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
637d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
638289bc588SBarry Smith 
63920563c6bSBarry Smith    Input Parameters:
6400c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
6410c775827SLois Curfman McInnes .  m - number of rows
6420c775827SLois Curfman McInnes .  n - number of column
64320563c6bSBarry Smith 
64420563c6bSBarry Smith    Output Parameter:
6450c775827SLois Curfman McInnes .  newmat - the matrix
64620563c6bSBarry Smith 
647dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
648d65003e9SLois Curfman McInnes 
649dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
65020563c6bSBarry Smith @*/
6516b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat)
652289bc588SBarry Smith {
65344a69424SLois Curfman McInnes   int       size = sizeof(Mat_Dense) + m*n*sizeof(Scalar);
654289bc588SBarry Smith   Mat mat;
65544a69424SLois Curfman McInnes   Mat_Dense    *l;
656289bc588SBarry Smith   *newmat        = 0;
6576b5873e3SBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm);
658a5a9c739SBarry Smith   PLogObjectCreate(mat);
65978b31e54SBarry Smith   l              = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l);
660289bc588SBarry Smith   mat->ops       = &MatOps;
66144a69424SLois Curfman McInnes   mat->destroy   = MatDestroy_Dense;
66244a69424SLois Curfman McInnes   mat->view      = MatView_Dense;
663289bc588SBarry Smith   mat->data      = (void *) l;
664289bc588SBarry Smith   mat->factor    = 0;
665d6dfbf8fSBarry Smith 
666289bc588SBarry Smith   l->m           = m;
667289bc588SBarry Smith   l->n           = n;
668289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
669289bc588SBarry Smith   l->pivots      = 0;
670289bc588SBarry Smith   l->roworiented = 1;
671d6dfbf8fSBarry Smith 
67278b31e54SBarry Smith   PETSCMEMSET(l->v,0,m*n*sizeof(Scalar));
673289bc588SBarry Smith   *newmat = mat;
674289bc588SBarry Smith   return 0;
675289bc588SBarry Smith }
676289bc588SBarry Smith 
67744a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat)
678289bc588SBarry Smith {
67944a69424SLois Curfman McInnes   Mat_Dense *m = (Mat_Dense *) matin->data;
6806b5873e3SBarry Smith   return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat);
681289bc588SBarry Smith }
682