xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 44a69424ff13a0ea173321b4754d03eaa842d22a)
1cb512458SBarry Smith #ifndef lint
2*44a69424SLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.17 1995/03/23 22:31:22 bsmith Exp curfman $";
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"
13289bc588SBarry Smith 
14289bc588SBarry Smith typedef struct {
15289bc588SBarry Smith   Scalar *v;
16289bc588SBarry Smith   int    roworiented;
17289bc588SBarry Smith   int    m,n,pad;
18289bc588SBarry Smith   int    *pivots;   /* pivots in LU factorization */
19*44a69424SLois Curfman McInnes } Mat_Dense;
20289bc588SBarry Smith 
21289bc588SBarry Smith 
22*44a69424SLois Curfman McInnes static int MatNz_Dense(Mat matin,int *nz)
23289bc588SBarry Smith {
24*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
25289bc588SBarry Smith   int    i,N = mat->m*mat->n,count = 0;
26289bc588SBarry Smith   Scalar *v = mat->v;
27289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
28289bc588SBarry Smith   *nz = count; return 0;
29289bc588SBarry Smith }
30*44a69424SLois Curfman McInnes static int MatMemory_Dense(Mat matin,int *mem)
31289bc588SBarry Smith {
32*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
33289bc588SBarry Smith   *mem = mat->m*mat->n*sizeof(Scalar); return 0;
34289bc588SBarry Smith }
35289bc588SBarry Smith 
36289bc588SBarry Smith /* ---------------------------------------------------------------*/
37289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
38289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
39*44a69424SLois Curfman McInnes static int MatLUFactor_Dense(Mat matin,IS row,IS col)
40289bc588SBarry Smith {
41*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
426abc6512SBarry Smith   int    info;
43289bc588SBarry Smith   if (!mat->pivots) {
44289bc588SBarry Smith     mat->pivots = (int *) MALLOC( mat->m*sizeof(int) );
45289bc588SBarry Smith     CHKPTR(mat->pivots);
46289bc588SBarry Smith   }
47289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
48289bc588SBarry Smith   if (info) SETERR(1,"Bad LU factorization");
49289bc588SBarry Smith   matin->factor = FACTOR_LU;
50289bc588SBarry Smith   return 0;
51289bc588SBarry Smith }
52*44a69424SLois Curfman McInnes static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact)
53289bc588SBarry Smith {
54289bc588SBarry Smith   int ierr;
556abc6512SBarry Smith   if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0);
56289bc588SBarry Smith   return 0;
57289bc588SBarry Smith }
58*44a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact)
59289bc588SBarry Smith {
6020563c6bSBarry Smith   return MatLUFactor(*fact,0,0);
61289bc588SBarry Smith }
62*44a69424SLois Curfman McInnes static int MatChFactorSymbolic_Dense(Mat matin,IS row,Mat *fact)
63289bc588SBarry Smith {
64289bc588SBarry Smith   int ierr;
656abc6512SBarry Smith   if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0);
66289bc588SBarry Smith   return 0;
67289bc588SBarry Smith }
68*44a69424SLois Curfman McInnes static int MatChFactorNumeric_Dense(Mat matin,Mat *fact)
69289bc588SBarry Smith {
7020563c6bSBarry Smith   return MatCholeskyFactor(*fact,0);
71289bc588SBarry Smith }
72*44a69424SLois Curfman McInnes static int MatChFactor_Dense(Mat matin,IS perm)
73289bc588SBarry Smith {
74*44a69424SLois Curfman McInnes   Mat_Dense    *mat = (Mat_Dense *) matin->data;
756abc6512SBarry Smith   int       info;
76289bc588SBarry Smith   if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;}
77289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
78289bc588SBarry Smith   if (info) SETERR(1,"Bad Cholesky factorization");
79289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
80289bc588SBarry Smith   return 0;
81289bc588SBarry Smith }
82289bc588SBarry Smith 
83*44a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy)
84289bc588SBarry Smith {
85*44a69424SLois 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);
89289bc588SBarry Smith   MEMCPY(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   }
989e25ed09SBarry Smith   else SETERR(1,"Matrix must be factored to solve");
99289bc588SBarry Smith   if (info) SETERR(1,"Bad solve");
100289bc588SBarry Smith   return 0;
101289bc588SBarry Smith }
102*44a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy)
103da3a660dSBarry Smith {
104*44a69424SLois 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);
108da3a660dSBarry Smith   MEMCPY(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   }
118da3a660dSBarry Smith   if (info) SETERR(1,"Bad solve");
119da3a660dSBarry Smith   return 0;
120da3a660dSBarry Smith }
121*44a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
122da3a660dSBarry Smith {
123*44a69424SLois 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) {
129da3a660dSBarry Smith     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
130da3a660dSBarry Smith     ierr = VecCopy(yy,tmp); CHKERR(ierr);
131da3a660dSBarry Smith   }
132da3a660dSBarry Smith   MEMCPY(y,x,mat->m*sizeof(Scalar));
133da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
134da3a660dSBarry Smith   if (mat->pivots) {
135da3a660dSBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
136da3a660dSBarry Smith               y, &mat->m, &info );
137da3a660dSBarry Smith   }
138da3a660dSBarry Smith   else {
139da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
140da3a660dSBarry Smith               y, &mat->m, &info );
141da3a660dSBarry Smith   }
142da3a660dSBarry Smith   if (info) SETERR(1,"Bad solve");
143da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
144da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
145da3a660dSBarry Smith   return 0;
146da3a660dSBarry Smith }
147*44a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy)
148da3a660dSBarry Smith {
149*44a69424SLois Curfman McInnes   Mat_Dense  *mat = (Mat_Dense *) matin->data;
1506abc6512SBarry Smith   int     one = 1, info,ierr;
1516abc6512SBarry Smith   Scalar  *x, *y, sone = 1.0;
152da3a660dSBarry Smith   Vec     tmp;
153da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
154da3a660dSBarry Smith   if (yy == zz) {
155da3a660dSBarry Smith     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
156da3a660dSBarry Smith     ierr = VecCopy(yy,tmp); CHKERR(ierr);
157da3a660dSBarry Smith   }
158da3a660dSBarry Smith   MEMCPY(y,x,mat->m*sizeof(Scalar));
159da3a660dSBarry Smith   /* assume if pivots exist then LU else Cholesky */
160da3a660dSBarry Smith   if (mat->pivots) {
161da3a660dSBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
162da3a660dSBarry Smith               y, &mat->m, &info );
163da3a660dSBarry Smith   }
164da3a660dSBarry Smith   else {
165da3a660dSBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
166da3a660dSBarry Smith               y, &mat->m, &info );
167da3a660dSBarry Smith   }
168da3a660dSBarry Smith   if (info) SETERR(1,"Bad solve");
169da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
170da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
171da3a660dSBarry Smith   return 0;
172da3a660dSBarry Smith }
173289bc588SBarry Smith /* ------------------------------------------------------------------*/
174*44a69424SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,int flag,double shift,
175289bc588SBarry Smith                        int its,Vec xx)
176289bc588SBarry Smith {
177*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
178289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
1796abc6512SBarry Smith   int    o = 1,ierr, m = mat->m, i;
180289bc588SBarry Smith 
181289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
182289bc588SBarry Smith     /* this is a hack fix, should have another version without
183289bc588SBarry Smith        the second BLdot */
1846abc6512SBarry Smith     if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0);
185289bc588SBarry Smith   }
186289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
187289bc588SBarry Smith   while (its--) {
188289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
189289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
190289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
191d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
192289bc588SBarry Smith       }
193289bc588SBarry Smith     }
194289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
195289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
196289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
197d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
198289bc588SBarry Smith       }
199289bc588SBarry Smith     }
200289bc588SBarry Smith   }
201289bc588SBarry Smith   return 0;
202289bc588SBarry Smith }
203289bc588SBarry Smith 
204289bc588SBarry Smith /* -----------------------------------------------------------------*/
205*44a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy)
206289bc588SBarry Smith {
207*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
208289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
209289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
210289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
211289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
212289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
213289bc588SBarry Smith   return 0;
214289bc588SBarry Smith }
215*44a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy)
216289bc588SBarry Smith {
217*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
218289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
219289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
220289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
221289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
222289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
223289bc588SBarry Smith   return 0;
224289bc588SBarry Smith }
225*44a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
226289bc588SBarry Smith {
227*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
228289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
2296abc6512SBarry Smith   int    _One=1; Scalar _DOne=1.0;
230289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
231289bc588SBarry Smith   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
232289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
233289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
234289bc588SBarry Smith   return 0;
235289bc588SBarry Smith }
236*44a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy)
237289bc588SBarry Smith {
238*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
239289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
240289bc588SBarry Smith   int    _One=1;
2416abc6512SBarry Smith   Scalar _DOne=1.0;
242289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
243289bc588SBarry Smith   VecGetArray(zz,&z);
244289bc588SBarry Smith   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
245289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
246289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
247289bc588SBarry Smith   return 0;
248289bc588SBarry Smith }
249289bc588SBarry Smith 
250289bc588SBarry Smith /* -----------------------------------------------------------------*/
251*44a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols,
252289bc588SBarry Smith                         Scalar **vals)
253289bc588SBarry Smith {
254*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
255289bc588SBarry Smith   Scalar *v;
256289bc588SBarry Smith   int    i;
257289bc588SBarry Smith   *ncols = mat->n;
258289bc588SBarry Smith   if (cols) {
259289bc588SBarry Smith     *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols);
260289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
261289bc588SBarry Smith   }
262289bc588SBarry Smith   if (vals) {
263289bc588SBarry Smith     *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals);
264289bc588SBarry Smith     v = mat->v + row;
265289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
266289bc588SBarry Smith   }
267289bc588SBarry Smith   return 0;
268289bc588SBarry Smith }
269*44a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols,
270289bc588SBarry Smith                             Scalar **vals)
271289bc588SBarry Smith {
272289bc588SBarry Smith   if (cols) { FREE(*cols); }
273289bc588SBarry Smith   if (vals) { FREE(*vals); }
274289bc588SBarry Smith   return 0;
275289bc588SBarry Smith }
276289bc588SBarry Smith /* ----------------------------------------------------------------*/
277*44a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n,
278e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
279289bc588SBarry Smith {
280*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
281289bc588SBarry Smith   int    i,j;
282d6dfbf8fSBarry Smith 
283289bc588SBarry Smith   if (!mat->roworiented) {
284e8d4e0b9SBarry Smith     if (addv == InsertValues) {
285289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
286d6dfbf8fSBarry Smith         if (indexn[j] < 0) {*v += m; continue;}
287289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
288d6dfbf8fSBarry Smith           if (indexm[i] < 0) {*v++; continue;}
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++ ) {
295d6dfbf8fSBarry Smith         if (indexn[j] < 0) {*v += m; continue;}
296289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
297d6dfbf8fSBarry Smith           if (indexm[i] < 0) {*v++; continue;}
298289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
299289bc588SBarry Smith         }
300289bc588SBarry Smith       }
301289bc588SBarry Smith     }
302e8d4e0b9SBarry Smith   }
303e8d4e0b9SBarry Smith   else {
304e8d4e0b9SBarry Smith     if (addv == InsertValues) {
305e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
306d6dfbf8fSBarry Smith         if (indexm[i] < 0) {*v += n; continue;}
307e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
308d6dfbf8fSBarry Smith           if (indexn[j] < 0) {*v++; continue;}
309e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
310e8d4e0b9SBarry Smith         }
311e8d4e0b9SBarry Smith       }
312e8d4e0b9SBarry Smith     }
313289bc588SBarry Smith     else {
314289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
315d6dfbf8fSBarry Smith         if (indexm[i] < 0) {*v += n; continue;}
316289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
317d6dfbf8fSBarry Smith           if (indexn[j] < 0) {*v++; continue;}
318289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
319289bc588SBarry Smith         }
320289bc588SBarry Smith       }
321289bc588SBarry Smith     }
322e8d4e0b9SBarry Smith   }
323289bc588SBarry Smith   return 0;
324289bc588SBarry Smith }
325e8d4e0b9SBarry Smith 
326289bc588SBarry Smith /* -----------------------------------------------------------------*/
327*44a69424SLois Curfman McInnes static int MatCopy_Dense(Mat matin,Mat *newmat)
328289bc588SBarry Smith {
329*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
330289bc588SBarry Smith   int ierr;
331289bc588SBarry Smith   Mat newi;
332*44a69424SLois Curfman McInnes   Mat_Dense *l;
3336abc6512SBarry Smith   if ((ierr = MatCreateSequentialDense(mat->m,mat->n,&newi))) SETERR(ierr,0);
334*44a69424SLois Curfman McInnes   l = (Mat_Dense *) newi->data;
335289bc588SBarry Smith   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
336289bc588SBarry Smith   *newmat = newi;
337289bc588SBarry Smith   return 0;
338289bc588SBarry Smith }
339da3a660dSBarry Smith #include "viewer.h"
340289bc588SBarry Smith 
341*44a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr)
342289bc588SBarry Smith {
343289bc588SBarry Smith   Mat         matin = (Mat) obj;
344*44a69424SLois Curfman McInnes   Mat_Dense      *mat = (Mat_Dense *) matin->data;
345289bc588SBarry Smith   Scalar      *v;
346289bc588SBarry Smith   int         i,j;
347da3a660dSBarry Smith   PetscObject ojb = (PetscObject) ptr;
348da3a660dSBarry Smith 
3496abc6512SBarry Smith   if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) {
350da3a660dSBarry Smith     return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v);
351da3a660dSBarry Smith   }
352da3a660dSBarry Smith   else {
353289bc588SBarry Smith     for ( i=0; i<mat->m; i++ ) {
354289bc588SBarry Smith       v = mat->v + i;
355289bc588SBarry Smith       for ( j=0; j<mat->n; j++ ) {
356289bc588SBarry Smith #if defined(PETSC_COMPLEX)
357da69df5fSBarry Smith         ViewerPrintf(ptr,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
358289bc588SBarry Smith #else
359da69df5fSBarry Smith         ViewerPrintf(ptr,"%6.4e ",*v); v += mat->m;
360289bc588SBarry Smith #endif
361289bc588SBarry Smith       }
362da69df5fSBarry Smith       ViewerPrintf(ptr,"\n");
363289bc588SBarry Smith     }
364da3a660dSBarry Smith   }
365289bc588SBarry Smith   return 0;
366289bc588SBarry Smith }
367289bc588SBarry Smith 
368289bc588SBarry Smith 
369*44a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj)
370289bc588SBarry Smith {
371289bc588SBarry Smith   Mat    mat = (Mat) obj;
372*44a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) mat->data;
373a5a9c739SBarry Smith #if defined(PETSC_LOG)
374a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
375a5a9c739SBarry Smith #endif
376289bc588SBarry Smith   if (l->pivots) FREE(l->pivots);
377289bc588SBarry Smith   FREE(l);
378a5a9c739SBarry Smith   PLogObjectDestroy(mat);
3799e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
380289bc588SBarry Smith   return 0;
381289bc588SBarry Smith }
382289bc588SBarry Smith 
383*44a69424SLois Curfman McInnes static int MatTrans_Dense(Mat matin)
384289bc588SBarry Smith {
385*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
386289bc588SBarry Smith   int    k,j;
387289bc588SBarry Smith   Scalar *v = mat->v, tmp;
388289bc588SBarry Smith   if (mat->m != mat->n) {
389289bc588SBarry Smith     SETERR(1,"Cannot transpose rectangular dense matrix");
390289bc588SBarry Smith   }
391289bc588SBarry Smith   for ( j=0; j<mat->m; j++ ) {
392289bc588SBarry Smith     for ( k=0; k<j; k++ ) {
393289bc588SBarry Smith       tmp = v[j + k*mat->n];
394289bc588SBarry Smith       v[j + k*mat->n] = v[k + j*mat->n];
395289bc588SBarry Smith       v[k + j*mat->n] = tmp;
396289bc588SBarry Smith     }
397289bc588SBarry Smith   }
398289bc588SBarry Smith   return 0;
399289bc588SBarry Smith }
400289bc588SBarry Smith 
401*44a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2)
402289bc588SBarry Smith {
403*44a69424SLois Curfman McInnes   Mat_Dense *mat1 = (Mat_Dense *) matin1->data;
404*44a69424SLois Curfman McInnes   Mat_Dense *mat2 = (Mat_Dense *) matin2->data;
405289bc588SBarry Smith   int    i;
406289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
407289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
408289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
409289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
410289bc588SBarry Smith     if (*v1 != *v2) return 0;
411289bc588SBarry Smith     v1++; v2++;
412289bc588SBarry Smith   }
413289bc588SBarry Smith   return 1;
414289bc588SBarry Smith }
415289bc588SBarry Smith 
416*44a69424SLois Curfman McInnes static int MatGetDiag_Dense(Mat matin,Vec v)
417289bc588SBarry Smith {
418*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
4196abc6512SBarry Smith   int    i, n;
4206abc6512SBarry Smith   Scalar *x;
421289bc588SBarry Smith   CHKTYPE(v,SEQVECTOR);
422289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
423289bc588SBarry Smith   if (n != mat->m) SETERR(1,"Nonconforming matrix and vector");
424289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
425289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
426289bc588SBarry Smith   }
427289bc588SBarry Smith   return 0;
428289bc588SBarry Smith }
429289bc588SBarry Smith 
430*44a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr)
431289bc588SBarry Smith {
432*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
433da3a660dSBarry Smith   Scalar *l,*r,x,*v;
434da3a660dSBarry Smith   int    i,j,m = mat->m, n = mat->n;
43528988994SBarry Smith   if (ll) {
436da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
437da3a660dSBarry Smith     if (m != mat->m) SETERR(1,"Left scaling vector wrong length");
438da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
439da3a660dSBarry Smith       x = l[i];
440da3a660dSBarry Smith       v = mat->v + i;
441da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
442da3a660dSBarry Smith     }
443da3a660dSBarry Smith   }
44428988994SBarry Smith   if (rr) {
445da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
446da3a660dSBarry Smith     if (n != mat->n) SETERR(1,"Right scaling vector wrong length");
447da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
448da3a660dSBarry Smith       x = r[i];
449da3a660dSBarry Smith       v = mat->v + i*m;
450da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
451da3a660dSBarry Smith     }
452da3a660dSBarry Smith   }
453289bc588SBarry Smith   return 0;
454289bc588SBarry Smith }
455289bc588SBarry Smith 
456da3a660dSBarry Smith 
457*44a69424SLois Curfman McInnes static int MatNorm_Dense(Mat matin,int type,double *norm)
458289bc588SBarry Smith {
459*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
460289bc588SBarry Smith   Scalar *v = mat->v;
461289bc588SBarry Smith   double sum = 0.0;
462289bc588SBarry Smith   int    i, j;
463289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
464289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
465289bc588SBarry Smith #if defined(PETSC_COMPLEX)
466289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
467289bc588SBarry Smith #else
468289bc588SBarry Smith       sum += (*v)*(*v); v++;
469289bc588SBarry Smith #endif
470289bc588SBarry Smith     }
471289bc588SBarry Smith     *norm = sqrt(sum);
472289bc588SBarry Smith   }
473289bc588SBarry Smith   else if (type == NORM_1) {
474289bc588SBarry Smith     *norm = 0.0;
475289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
476289bc588SBarry Smith       sum = 0.0;
477289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
478289bc588SBarry Smith #if defined(PETSC_COMPLEX)
479289bc588SBarry Smith         sum += abs(*v++);
480289bc588SBarry Smith #else
481289bc588SBarry Smith         sum += fabs(*v++);
482289bc588SBarry Smith #endif
483289bc588SBarry Smith       }
484289bc588SBarry Smith       if (sum > *norm) *norm = sum;
485289bc588SBarry Smith     }
486289bc588SBarry Smith   }
487289bc588SBarry Smith   else if (type == NORM_INFINITY) {
488289bc588SBarry Smith     *norm = 0.0;
489289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
490289bc588SBarry Smith       v = mat->v + j;
491289bc588SBarry Smith       sum = 0.0;
492289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
493289bc588SBarry Smith #if defined(PETSC_COMPLEX)
494289bc588SBarry Smith         sum += abs(*v); v += mat->m;
495289bc588SBarry Smith #else
496289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
497289bc588SBarry Smith #endif
498289bc588SBarry Smith       }
499289bc588SBarry Smith       if (sum > *norm) *norm = sum;
500289bc588SBarry Smith     }
501289bc588SBarry Smith   }
502289bc588SBarry Smith   else {
503289bc588SBarry Smith     SETERR(1,"No support for the two norm yet");
504289bc588SBarry Smith   }
505289bc588SBarry Smith   return 0;
506289bc588SBarry Smith }
507289bc588SBarry Smith 
508*44a69424SLois Curfman McInnes static int MatInsOpt_Dense(Mat aijin,int op)
509289bc588SBarry Smith {
510*44a69424SLois Curfman McInnes   Mat_Dense *aij = (Mat_Dense *) aijin->data;
511289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
512289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
513289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
514289bc588SBarry Smith   return 0;
515289bc588SBarry Smith }
516289bc588SBarry Smith 
517*44a69424SLois Curfman McInnes static int MatZero_Dense(Mat A)
5186f0a148fSBarry Smith {
519*44a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
5206f0a148fSBarry Smith   MEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
5216f0a148fSBarry Smith   return 0;
5226f0a148fSBarry Smith }
5236f0a148fSBarry Smith 
524*44a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag)
5256f0a148fSBarry Smith {
526*44a69424SLois Curfman McInnes   Mat_Dense *l = (Mat_Dense *) A->data;
5276abc6512SBarry Smith   int    n = l->n, i, j,ierr,N, *rows;
5286f0a148fSBarry Smith   Scalar *slot;
529260925b8SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERR(ierr);
530260925b8SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERR(ierr);
5316f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
5326f0a148fSBarry Smith     slot = l->v + rows[i];
5336f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
5346f0a148fSBarry Smith   }
5356f0a148fSBarry Smith   if (diag) {
5366f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
5376f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
538260925b8SBarry Smith       *slot = *diag;
5396f0a148fSBarry Smith     }
5406f0a148fSBarry Smith   }
541260925b8SBarry Smith   ISRestoreIndices(is,&rows);
5426f0a148fSBarry Smith   return 0;
5436f0a148fSBarry Smith }
544557bce09SLois Curfman McInnes 
545*44a69424SLois Curfman McInnes static int MatSize_Dense(Mat matin,int *m,int *n)
546557bce09SLois Curfman McInnes {
547*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
548557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
549557bce09SLois Curfman McInnes   return 0;
550557bce09SLois Curfman McInnes }
551557bce09SLois Curfman McInnes 
552*44a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array)
55364e87e97SBarry Smith {
554*44a69424SLois Curfman McInnes   Mat_Dense *mat = (Mat_Dense *) matin->data;
55564e87e97SBarry Smith   *array = mat->v;
55664e87e97SBarry Smith   return 0;
55764e87e97SBarry Smith }
558289bc588SBarry Smith /* -------------------------------------------------------------------*/
559*44a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense,
560*44a69424SLois Curfman McInnes        MatGetRow_Dense, MatRestoreRow_Dense,
561*44a69424SLois Curfman McInnes        MatMult_Dense, MatMultAdd_Dense,
562*44a69424SLois Curfman McInnes        MatMultTrans_Dense, MatMultTransAdd_Dense,
563*44a69424SLois Curfman McInnes        MatSolve_Dense,MatSolveAdd_Dense,
564*44a69424SLois Curfman McInnes        MatSolveTrans_Dense,MatSolveTransAdd_Dense,
565*44a69424SLois Curfman McInnes        MatLUFactor_Dense,MatChFactor_Dense,
566*44a69424SLois Curfman McInnes        MatRelax_Dense,
567*44a69424SLois Curfman McInnes        MatTrans_Dense,
568*44a69424SLois Curfman McInnes        MatNz_Dense,MatMemory_Dense,MatEqual_Dense,
569*44a69424SLois Curfman McInnes        MatCopy_Dense,
570*44a69424SLois Curfman McInnes        MatGetDiag_Dense,MatScale_Dense,MatNorm_Dense,
571289bc588SBarry Smith        0,0,
572*44a69424SLois Curfman McInnes        0, MatInsOpt_Dense,MatZero_Dense,MatZeroRows_Dense,0,
573*44a69424SLois Curfman McInnes        MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense,
574*44a69424SLois Curfman McInnes        MatChFactorSymbolic_Dense,MatChFactorNumeric_Dense,
575*44a69424SLois Curfman McInnes        MatSize_Dense,MatSize_Dense,0,
576*44a69424SLois Curfman McInnes        0,0,MatGetArray_Dense
577289bc588SBarry Smith };
57820563c6bSBarry Smith /*@
57920563c6bSBarry Smith     MatCreateSequentialDense - Creates a sequential dense matrix that
58020563c6bSBarry Smith         is stored in the usual Fortran 77 manner. Many of the matrix
58120563c6bSBarry Smith         operations use the BLAS and LAPACK routines.
582289bc588SBarry Smith 
58320563c6bSBarry Smith   Input Parameters:
58420563c6bSBarry Smith .   m, n - the number of rows and columns in the matrix.
58520563c6bSBarry Smith 
58620563c6bSBarry Smith   Output Parameter:
58720563c6bSBarry Smith .  newmat - the matrix created.
58820563c6bSBarry Smith 
58920563c6bSBarry Smith   Keywords: dense matrix, lapack, blas
59020563c6bSBarry Smith @*/
591289bc588SBarry Smith int MatCreateSequentialDense(int m,int n,Mat *newmat)
592289bc588SBarry Smith {
593*44a69424SLois Curfman McInnes   int       size = sizeof(Mat_Dense) + m*n*sizeof(Scalar);
594289bc588SBarry Smith   Mat mat;
595*44a69424SLois Curfman McInnes   Mat_Dense    *l;
596289bc588SBarry Smith   *newmat        = 0;
597c951ae0fSBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,MPI_COMM_SELF);
598a5a9c739SBarry Smith   PLogObjectCreate(mat);
599*44a69424SLois Curfman McInnes   l              = (Mat_Dense *) MALLOC(size); CHKPTR(l);
600289bc588SBarry Smith   mat->ops       = &MatOps;
601*44a69424SLois Curfman McInnes   mat->destroy   = MatDestroy_Dense;
602*44a69424SLois Curfman McInnes   mat->view      = MatView_Dense;
603289bc588SBarry Smith   mat->data      = (void *) l;
604289bc588SBarry Smith   mat->factor    = 0;
605289bc588SBarry Smith   mat->col       = 0;
606289bc588SBarry Smith   mat->row       = 0;
607d6dfbf8fSBarry Smith 
608289bc588SBarry Smith   l->m           = m;
609289bc588SBarry Smith   l->n           = n;
610289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
611289bc588SBarry Smith   l->pivots      = 0;
612289bc588SBarry Smith   l->roworiented = 1;
613d6dfbf8fSBarry Smith 
614289bc588SBarry Smith   MEMSET(l->v,0,m*n*sizeof(Scalar));
615289bc588SBarry Smith   *newmat = mat;
616289bc588SBarry Smith   return 0;
617289bc588SBarry Smith }
618289bc588SBarry Smith 
619*44a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat)
620289bc588SBarry Smith {
621*44a69424SLois Curfman McInnes   Mat_Dense *m = (Mat_Dense *) matin->data;
622289bc588SBarry Smith   return MatCreateSequentialDense(m->m,m->n,newmat);
623289bc588SBarry Smith }
624