xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 6abc6512c105a871402f018420c4c4316bc1e68d)
1289bc588SBarry Smith 
2289bc588SBarry Smith /*
3289bc588SBarry Smith     Standard Fortran style matrices
4289bc588SBarry Smith */
5289bc588SBarry Smith #include "ptscimpl.h"
6289bc588SBarry Smith #include "plapack.h"
7289bc588SBarry Smith #include "matimpl.h"
8289bc588SBarry Smith #include "math.h"
9289bc588SBarry Smith #include "vec/vecimpl.h"
10289bc588SBarry Smith 
11289bc588SBarry Smith typedef struct {
12289bc588SBarry Smith   Scalar *v;
13289bc588SBarry Smith   int    roworiented;
14289bc588SBarry Smith   int    m,n,pad;
15289bc588SBarry Smith   int    *pivots;   /* pivots in LU factorization */
16289bc588SBarry Smith } MatiSD;
17289bc588SBarry Smith 
18289bc588SBarry Smith 
19289bc588SBarry Smith static int MatiSDnz(Mat matin,int *nz)
20289bc588SBarry Smith {
21289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
22289bc588SBarry Smith   int    i,N = mat->m*mat->n,count = 0;
23289bc588SBarry Smith   Scalar *v = mat->v;
24289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
25289bc588SBarry Smith   *nz = count; return 0;
26289bc588SBarry Smith }
27289bc588SBarry Smith static int MatiSDmemory(Mat matin,int *mem)
28289bc588SBarry Smith {
29289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
30289bc588SBarry Smith   *mem = mat->m*mat->n*sizeof(Scalar); 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.*/
36289bc588SBarry Smith static int MatiSDlufactor(Mat matin,IS row,IS col)
37289bc588SBarry Smith {
38289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
39*6abc6512SBarry Smith   int    info;
40289bc588SBarry Smith   if (!mat->pivots) {
41289bc588SBarry Smith     mat->pivots = (int *) MALLOC( mat->m*sizeof(int) );
42289bc588SBarry Smith     CHKPTR(mat->pivots);
43289bc588SBarry Smith   }
44289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
45289bc588SBarry Smith   if (info) SETERR(1,"Bad LU factorization");
46289bc588SBarry Smith   matin->factor = FACTOR_LU;
47289bc588SBarry Smith   return 0;
48289bc588SBarry Smith }
49289bc588SBarry Smith static int MatiSDlufactorsymbolic(Mat matin,IS row,IS col,Mat *fact)
50289bc588SBarry Smith {
51289bc588SBarry Smith   int ierr;
52*6abc6512SBarry Smith   if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0);
53289bc588SBarry Smith   return 0;
54289bc588SBarry Smith }
5520563c6bSBarry Smith static int MatiSDlufactornumeric(Mat matin,Mat *fact)
56289bc588SBarry Smith {
5720563c6bSBarry Smith   return MatLUFactor(*fact,0,0);
58289bc588SBarry Smith }
59289bc588SBarry Smith static int MatiSDchfactorsymbolic(Mat matin,IS row,Mat *fact)
60289bc588SBarry Smith {
61289bc588SBarry Smith   int ierr;
62*6abc6512SBarry Smith   if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0);
63289bc588SBarry Smith   return 0;
64289bc588SBarry Smith }
6520563c6bSBarry Smith static int MatiSDchfactornumeric(Mat matin,Mat *fact)
66289bc588SBarry Smith {
6720563c6bSBarry Smith   return MatCholeskyFactor(*fact,0);
68289bc588SBarry Smith }
69289bc588SBarry Smith static int MatiSDchfactor(Mat matin,IS perm)
70289bc588SBarry Smith {
71289bc588SBarry Smith   MatiSD    *mat = (MatiSD *) matin->data;
72*6abc6512SBarry Smith   int       info;
73289bc588SBarry Smith   if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;}
74289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
75289bc588SBarry Smith   if (info) SETERR(1,"Bad Cholesky factorization");
76289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
77289bc588SBarry Smith   return 0;
78289bc588SBarry Smith }
79289bc588SBarry Smith 
80289bc588SBarry Smith static int MatiSDsolve(Mat matin,Vec xx,Vec yy)
81289bc588SBarry Smith {
82289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
83*6abc6512SBarry Smith   int    one = 1, info;
84*6abc6512SBarry Smith   Scalar *x, *y;
85289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
86289bc588SBarry Smith   MEMCPY(y,x,mat->m*sizeof(Scalar));
87289bc588SBarry Smith   /* assume if pivots exist then LU else Cholesky */
88289bc588SBarry Smith   if (mat->pivots) {
89289bc588SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
90289bc588SBarry Smith               y, &mat->m, &info );
91289bc588SBarry Smith   }
92289bc588SBarry Smith   else {
93289bc588SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
94289bc588SBarry Smith               y, &mat->m, &info );
95289bc588SBarry Smith   }
96289bc588SBarry Smith   if (info) SETERR(1,"Bad solve");
97289bc588SBarry Smith   return 0;
98289bc588SBarry Smith }
99289bc588SBarry Smith static int MatiSDsolvetrans(Mat matin,Vec xx,Vec yy)
100da3a660dSBarry Smith {
101da3a660dSBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
102*6abc6512SBarry Smith   int    one = 1, info;
103*6abc6512SBarry Smith   Scalar *x, *y;
104da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
105da3a660dSBarry Smith   MEMCPY(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   }
115da3a660dSBarry Smith   if (info) SETERR(1,"Bad solve");
116da3a660dSBarry Smith   return 0;
117da3a660dSBarry Smith }
118da3a660dSBarry Smith static int MatiSDsolveadd(Mat matin,Vec xx,Vec zz,Vec yy)
119da3a660dSBarry Smith {
120da3a660dSBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
121*6abc6512SBarry Smith   int    one = 1, info,ierr;
122*6abc6512SBarry 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) {
126da3a660dSBarry Smith     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
127da3a660dSBarry Smith     ierr = VecCopy(yy,tmp); CHKERR(ierr);
128da3a660dSBarry Smith   }
129da3a660dSBarry Smith   MEMCPY(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   }
139da3a660dSBarry Smith   if (info) SETERR(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 }
144da3a660dSBarry Smith static int MatiSDsolvetransadd(Mat matin,Vec xx,Vec zz, Vec yy)
145da3a660dSBarry Smith {
146da3a660dSBarry Smith   MatiSD  *mat = (MatiSD *) matin->data;
147*6abc6512SBarry Smith   int     one = 1, info,ierr;
148*6abc6512SBarry Smith   Scalar  *x, *y, sone = 1.0;
149da3a660dSBarry Smith   Vec     tmp;
150da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
151da3a660dSBarry Smith   if (yy == zz) {
152da3a660dSBarry Smith     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
153da3a660dSBarry Smith     ierr = VecCopy(yy,tmp); CHKERR(ierr);
154da3a660dSBarry Smith   }
155da3a660dSBarry Smith   MEMCPY(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   }
165da3a660dSBarry Smith   if (info) SETERR(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 /* ------------------------------------------------------------------*/
171d6dfbf8fSBarry Smith static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,double shift,
172289bc588SBarry Smith                        int its,Vec xx)
173289bc588SBarry Smith {
174289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
175289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
176*6abc6512SBarry 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*6abc6512SBarry Smith     if ((ierr = VecSet(&zero,xx))) SETERR(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 /* -----------------------------------------------------------------*/
202289bc588SBarry Smith static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy)
203289bc588SBarry Smith {
204289bc588SBarry Smith   MatiSD *mat = (MatiSD *) 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 }
212289bc588SBarry Smith static int MatiSDmult(Mat matin,Vec xx,Vec yy)
213289bc588SBarry Smith {
214289bc588SBarry Smith   MatiSD *mat = (MatiSD *) 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 }
222289bc588SBarry Smith static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy)
223289bc588SBarry Smith {
224289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
225289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
226*6abc6512SBarry Smith   int    _One=1; Scalar _DOne=1.0;
227289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
228289bc588SBarry Smith   if (zz != yy) MEMCPY(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 }
233289bc588SBarry Smith static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy)
234289bc588SBarry Smith {
235289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
236289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
237289bc588SBarry Smith   int    _One=1;
238*6abc6512SBarry Smith   Scalar _DOne=1.0;
239289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
240289bc588SBarry Smith   VecGetArray(zz,&z);
241289bc588SBarry Smith   if (zz != yy) MEMCPY(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 /* -----------------------------------------------------------------*/
248289bc588SBarry Smith static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols,
249289bc588SBarry Smith                         Scalar **vals)
250289bc588SBarry Smith {
251289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
252289bc588SBarry Smith   Scalar *v;
253289bc588SBarry Smith   int    i;
254289bc588SBarry Smith   *ncols = mat->n;
255289bc588SBarry Smith   if (cols) {
256289bc588SBarry Smith     *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols);
257289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
258289bc588SBarry Smith   }
259289bc588SBarry Smith   if (vals) {
260289bc588SBarry Smith     *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*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 }
266289bc588SBarry Smith static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols,
267289bc588SBarry Smith                             Scalar **vals)
268289bc588SBarry Smith {
269289bc588SBarry Smith   if (cols) { FREE(*cols); }
270289bc588SBarry Smith   if (vals) { FREE(*vals); }
271289bc588SBarry Smith   return 0;
272289bc588SBarry Smith }
273289bc588SBarry Smith /* ----------------------------------------------------------------*/
274e8d4e0b9SBarry Smith static int MatiSDinsert(Mat matin,int m,int *indexm,int n,
275e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
276289bc588SBarry Smith {
277289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
278289bc588SBarry Smith   int    i,j;
279d6dfbf8fSBarry Smith 
280289bc588SBarry Smith   if (!mat->roworiented) {
281e8d4e0b9SBarry Smith     if (addv == InsertValues) {
282289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
283d6dfbf8fSBarry Smith         if (indexn[j] < 0) {*v += m; continue;}
284289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
285d6dfbf8fSBarry Smith           if (indexm[i] < 0) {*v++; continue;}
286289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
287289bc588SBarry Smith         }
288289bc588SBarry Smith       }
289289bc588SBarry Smith     }
290289bc588SBarry Smith     else {
291289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
292d6dfbf8fSBarry Smith         if (indexn[j] < 0) {*v += m; continue;}
293289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
294d6dfbf8fSBarry Smith           if (indexm[i] < 0) {*v++; continue;}
295289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
296289bc588SBarry Smith         }
297289bc588SBarry Smith       }
298289bc588SBarry Smith     }
299e8d4e0b9SBarry Smith   }
300e8d4e0b9SBarry Smith   else {
301e8d4e0b9SBarry Smith     if (addv == InsertValues) {
302e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
303d6dfbf8fSBarry Smith         if (indexm[i] < 0) {*v += n; continue;}
304e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
305d6dfbf8fSBarry Smith           if (indexn[j] < 0) {*v++; continue;}
306e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
307e8d4e0b9SBarry Smith         }
308e8d4e0b9SBarry Smith       }
309e8d4e0b9SBarry Smith     }
310289bc588SBarry Smith     else {
311289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
312d6dfbf8fSBarry Smith         if (indexm[i] < 0) {*v += n; continue;}
313289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
314d6dfbf8fSBarry Smith           if (indexn[j] < 0) {*v++; continue;}
315289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
316289bc588SBarry Smith         }
317289bc588SBarry Smith       }
318289bc588SBarry Smith     }
319e8d4e0b9SBarry Smith   }
320289bc588SBarry Smith   return 0;
321289bc588SBarry Smith }
322e8d4e0b9SBarry Smith 
323289bc588SBarry Smith /* -----------------------------------------------------------------*/
324289bc588SBarry Smith static int MatiSDcopy(Mat matin,Mat *newmat)
325289bc588SBarry Smith {
326289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
327289bc588SBarry Smith   int ierr;
328289bc588SBarry Smith   Mat newi;
329289bc588SBarry Smith   MatiSD *l;
330*6abc6512SBarry Smith   if ((ierr = MatCreateSequentialDense(mat->m,mat->n,&newi))) SETERR(ierr,0);
331289bc588SBarry Smith   l = (MatiSD *) newi->data;
332289bc588SBarry Smith   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
333289bc588SBarry Smith   *newmat = newi;
334289bc588SBarry Smith   return 0;
335289bc588SBarry Smith }
336da3a660dSBarry Smith #include "viewer.h"
337289bc588SBarry Smith 
338289bc588SBarry Smith int MatiSDview(PetscObject obj,Viewer ptr)
339289bc588SBarry Smith {
340289bc588SBarry Smith   Mat         matin = (Mat) obj;
341289bc588SBarry Smith   MatiSD      *mat = (MatiSD *) matin->data;
342289bc588SBarry Smith   Scalar      *v;
343289bc588SBarry Smith   int         i,j;
344da3a660dSBarry Smith   PetscObject ojb = (PetscObject) ptr;
345da3a660dSBarry Smith 
346*6abc6512SBarry Smith   if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) {
347da3a660dSBarry Smith     return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v);
348da3a660dSBarry Smith   }
349da3a660dSBarry Smith   else {
350289bc588SBarry Smith     for ( i=0; i<mat->m; i++ ) {
351289bc588SBarry Smith       v = mat->v + i;
352289bc588SBarry Smith       for ( j=0; j<mat->n; j++ ) {
353289bc588SBarry Smith #if defined(PETSC_COMPLEX)
354289bc588SBarry Smith         printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
355289bc588SBarry Smith #else
356289bc588SBarry Smith         printf("%6.4e ",*v); v += mat->m;
357289bc588SBarry Smith #endif
358289bc588SBarry Smith       }
359289bc588SBarry Smith       printf("\n");
360289bc588SBarry Smith     }
361da3a660dSBarry Smith   }
362289bc588SBarry Smith   return 0;
363289bc588SBarry Smith }
364289bc588SBarry Smith 
365289bc588SBarry Smith 
366289bc588SBarry Smith static int MatiSDdestroy(PetscObject obj)
367289bc588SBarry Smith {
368289bc588SBarry Smith   Mat mat = (Mat) obj;
369289bc588SBarry Smith   MatiSD *l = (MatiSD *) mat->data;
370289bc588SBarry Smith   if (l->pivots) FREE(l->pivots);
371289bc588SBarry Smith   FREE(l);
372289bc588SBarry Smith   FREE(mat);
373289bc588SBarry Smith   return 0;
374289bc588SBarry Smith }
375289bc588SBarry Smith 
376289bc588SBarry Smith static int MatiSDtrans(Mat matin)
377289bc588SBarry Smith {
378289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
379289bc588SBarry Smith   int    k,j;
380289bc588SBarry Smith   Scalar *v = mat->v, tmp;
381289bc588SBarry Smith   if (mat->m != mat->n) {
382289bc588SBarry Smith     SETERR(1,"Cannot transpose rectangular dense matrix");
383289bc588SBarry Smith   }
384289bc588SBarry Smith   for ( j=0; j<mat->m; j++ ) {
385289bc588SBarry Smith     for ( k=0; k<j; k++ ) {
386289bc588SBarry Smith       tmp = v[j + k*mat->n];
387289bc588SBarry Smith       v[j + k*mat->n] = v[k + j*mat->n];
388289bc588SBarry Smith       v[k + j*mat->n] = tmp;
389289bc588SBarry Smith     }
390289bc588SBarry Smith   }
391289bc588SBarry Smith   return 0;
392289bc588SBarry Smith }
393289bc588SBarry Smith 
394289bc588SBarry Smith static int MatiSDequal(Mat matin1,Mat matin2)
395289bc588SBarry Smith {
396289bc588SBarry Smith   MatiSD *mat1 = (MatiSD *) matin1->data;
397289bc588SBarry Smith   MatiSD *mat2 = (MatiSD *) matin2->data;
398289bc588SBarry Smith   int    i;
399289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
400289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
401289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
402289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
403289bc588SBarry Smith     if (*v1 != *v2) return 0;
404289bc588SBarry Smith     v1++; v2++;
405289bc588SBarry Smith   }
406289bc588SBarry Smith   return 1;
407289bc588SBarry Smith }
408289bc588SBarry Smith 
409289bc588SBarry Smith static int MatiSDgetdiag(Mat matin,Vec v)
410289bc588SBarry Smith {
411289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
412*6abc6512SBarry Smith   int    i, n;
413*6abc6512SBarry Smith   Scalar *x;
414289bc588SBarry Smith   CHKTYPE(v,SEQVECTOR);
415289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
416289bc588SBarry Smith   if (n != mat->m) SETERR(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 
423da3a660dSBarry Smith static int MatiSDscale(Mat matin,Vec ll,Vec rr)
424289bc588SBarry Smith {
425da3a660dSBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
426da3a660dSBarry Smith   Scalar *l,*r,x,*v;
427da3a660dSBarry Smith   int    i,j,m = mat->m, n = mat->n;
428da3a660dSBarry Smith   if (l) {
429da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
430da3a660dSBarry Smith     if (m != mat->m) SETERR(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   }
437da3a660dSBarry Smith   if (l) {
438da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
439da3a660dSBarry Smith     if (n != mat->n) SETERR(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 
450289bc588SBarry Smith static int MatiSDnorm(Mat matin,int type,double *norm)
451289bc588SBarry Smith {
452289bc588SBarry Smith   MatiSD *mat = (MatiSD *) 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 {
496289bc588SBarry Smith     SETERR(1,"No support for the two norm yet");
497289bc588SBarry Smith   }
498289bc588SBarry Smith   return 0;
499289bc588SBarry Smith }
500289bc588SBarry Smith 
501289bc588SBarry Smith static int MatiDenseinsopt(Mat aijin,int op)
502289bc588SBarry Smith {
503289bc588SBarry Smith   MatiSD *aij = (MatiSD *) 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 
5106f0a148fSBarry Smith static int MatiZero(Mat A)
5116f0a148fSBarry Smith {
5126f0a148fSBarry Smith   MatiSD *l = (MatiSD *) A->data;
5136f0a148fSBarry Smith   MEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
5146f0a148fSBarry Smith   return 0;
5156f0a148fSBarry Smith }
5166f0a148fSBarry Smith 
517260925b8SBarry Smith static int MatiZerorows(Mat A,IS is,Scalar *diag)
5186f0a148fSBarry Smith {
5196f0a148fSBarry Smith   MatiSD *l = (MatiSD *) A->data;
520*6abc6512SBarry Smith   int    n = l->n, i, j,ierr,N, *rows;
5216f0a148fSBarry Smith   Scalar *slot;
522260925b8SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERR(ierr);
523260925b8SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERR(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 }
537289bc588SBarry Smith /* -------------------------------------------------------------------*/
538e8d4e0b9SBarry Smith static struct _MatOps MatOps = {MatiSDinsert,
539289bc588SBarry Smith        MatiSDgetrow, MatiSDrestorerow,
540289bc588SBarry Smith        MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd,
541da3a660dSBarry Smith        MatiSDsolve,MatiSDsolveadd,MatiSDsolvetrans,MatiSDsolvetransadd,
542289bc588SBarry Smith        MatiSDlufactor,MatiSDchfactor,
543289bc588SBarry Smith        MatiSDrelax,
544289bc588SBarry Smith        MatiSDtrans,
545289bc588SBarry Smith        MatiSDnz,MatiSDmemory,MatiSDequal,
546289bc588SBarry Smith        MatiSDcopy,
547289bc588SBarry Smith        MatiSDgetdiag,MatiSDscale,MatiSDnorm,
548289bc588SBarry Smith        0,0,
5496f0a148fSBarry Smith        0, MatiDenseinsopt,MatiZero,MatiZerorows,0,
550289bc588SBarry Smith        MatiSDlufactorsymbolic,MatiSDlufactornumeric,
551289bc588SBarry Smith        MatiSDchfactorsymbolic,MatiSDchfactornumeric
552289bc588SBarry Smith };
55320563c6bSBarry Smith /*@
55420563c6bSBarry Smith     MatCreateSequentialDense - Creates a sequential dense matrix that
55520563c6bSBarry Smith         is stored in the usual Fortran 77 manner. Many of the matrix
55620563c6bSBarry Smith         operations use the BLAS and LAPACK routines.
557289bc588SBarry Smith 
55820563c6bSBarry Smith   Input Parameters:
55920563c6bSBarry Smith .   m, n - the number of rows and columns in the matrix.
56020563c6bSBarry Smith 
56120563c6bSBarry Smith   Output Parameter:
56220563c6bSBarry Smith .  newmat - the matrix created.
56320563c6bSBarry Smith 
56420563c6bSBarry Smith   Keywords: dense matrix, lapack, blas
56520563c6bSBarry Smith @*/
566289bc588SBarry Smith int MatCreateSequentialDense(int m,int n,Mat *newmat)
567289bc588SBarry Smith {
568289bc588SBarry Smith   int       size = sizeof(MatiSD) + m*n*sizeof(Scalar);
569289bc588SBarry Smith   Mat mat;
570289bc588SBarry Smith   MatiSD    *l;
571289bc588SBarry Smith   *newmat        = 0;
572289bc588SBarry Smith   CREATEHEADER(mat,_Mat);
573289bc588SBarry Smith   l              = (MatiSD *) MALLOC(size); CHKPTR(l);
574289bc588SBarry Smith   mat->cookie    = MAT_COOKIE;
575289bc588SBarry Smith   mat->ops       = &MatOps;
576289bc588SBarry Smith   mat->destroy   = MatiSDdestroy;
577289bc588SBarry Smith   mat->view      = MatiSDview;
578289bc588SBarry Smith   mat->data      = (void *) l;
579289bc588SBarry Smith   mat->type      = MATDENSESEQ;
580289bc588SBarry Smith   mat->factor    = 0;
581289bc588SBarry Smith   mat->col       = 0;
582289bc588SBarry Smith   mat->row       = 0;
583d6dfbf8fSBarry Smith 
584289bc588SBarry Smith   l->m           = m;
585289bc588SBarry Smith   l->n           = n;
586289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
587289bc588SBarry Smith   l->pivots      = 0;
588289bc588SBarry Smith   l->roworiented = 1;
589d6dfbf8fSBarry Smith 
590289bc588SBarry Smith   MEMSET(l->v,0,m*n*sizeof(Scalar));
591289bc588SBarry Smith   *newmat = mat;
592289bc588SBarry Smith   return 0;
593289bc588SBarry Smith }
594289bc588SBarry Smith 
595289bc588SBarry Smith int MatiSDCreate(Mat matin,Mat *newmat)
596289bc588SBarry Smith {
597289bc588SBarry Smith   MatiSD *m = (MatiSD *) matin->data;
598289bc588SBarry Smith   return MatCreateSequentialDense(m->m,m->n,newmat);
599289bc588SBarry Smith }
600