xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 289bc588618de054f36015f5c3d330cf718be7c1)
1*289bc588SBarry Smith 
2*289bc588SBarry Smith /*
3*289bc588SBarry Smith     Standard Fortran style matrices
4*289bc588SBarry Smith */
5*289bc588SBarry Smith #include "ptscimpl.h"
6*289bc588SBarry Smith #include "plapack.h"
7*289bc588SBarry Smith #include "matimpl.h"
8*289bc588SBarry Smith #include "math.h"
9*289bc588SBarry Smith #include "vec/vecimpl.h"
10*289bc588SBarry Smith 
11*289bc588SBarry Smith typedef struct {
12*289bc588SBarry Smith   Scalar *v;
13*289bc588SBarry Smith   int    roworiented;
14*289bc588SBarry Smith   int    m,n,pad;
15*289bc588SBarry Smith   int    *pivots;   /* pivots in LU factorization */
16*289bc588SBarry Smith } MatiSD;
17*289bc588SBarry Smith 
18*289bc588SBarry Smith 
19*289bc588SBarry Smith static int MatiSDnz(Mat matin,int *nz)
20*289bc588SBarry Smith {
21*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
22*289bc588SBarry Smith   int    i,N = mat->m*mat->n,count = 0;
23*289bc588SBarry Smith   Scalar *v = mat->v;
24*289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
25*289bc588SBarry Smith   *nz = count; return 0;
26*289bc588SBarry Smith }
27*289bc588SBarry Smith static int MatiSDmemory(Mat matin,int *mem)
28*289bc588SBarry Smith {
29*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
30*289bc588SBarry Smith   *mem = mat->m*mat->n*sizeof(Scalar); return 0;
31*289bc588SBarry Smith }
32*289bc588SBarry Smith 
33*289bc588SBarry Smith /* ---------------------------------------------------------------*/
34*289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
35*289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
36*289bc588SBarry Smith static int MatiSDlufactor(Mat matin,IS row,IS col)
37*289bc588SBarry Smith {
38*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
39*289bc588SBarry Smith   int    ierr, one  = 1, info;
40*289bc588SBarry Smith   if (!mat->pivots) {
41*289bc588SBarry Smith     mat->pivots = (int *) MALLOC( mat->m*sizeof(int) );
42*289bc588SBarry Smith     CHKPTR(mat->pivots);
43*289bc588SBarry Smith   }
44*289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
45*289bc588SBarry Smith   if (info) SETERR(1,"Bad LU factorization");
46*289bc588SBarry Smith   matin->factor = FACTOR_LU;
47*289bc588SBarry Smith   return 0;
48*289bc588SBarry Smith }
49*289bc588SBarry Smith static int MatiSDlufactorsymbolic(Mat matin,IS row,IS col,Mat *fact)
50*289bc588SBarry Smith {
51*289bc588SBarry Smith   int ierr;
52*289bc588SBarry Smith   if (ierr = MatCopy(matin,fact)) SETERR(ierr,0);
53*289bc588SBarry Smith   return 0;
54*289bc588SBarry Smith }
55*289bc588SBarry Smith static int MatiSDlufactornumeric(Mat matin,Mat fact)
56*289bc588SBarry Smith {
57*289bc588SBarry Smith   return MatLUFactor(fact,0,0);
58*289bc588SBarry Smith }
59*289bc588SBarry Smith static int MatiSDchfactorsymbolic(Mat matin,IS row,Mat *fact)
60*289bc588SBarry Smith {
61*289bc588SBarry Smith   int ierr;
62*289bc588SBarry Smith   if (ierr = MatCopy(matin,fact)) SETERR(ierr,0);
63*289bc588SBarry Smith   return 0;
64*289bc588SBarry Smith }
65*289bc588SBarry Smith static int MatiSDchfactornumeric(Mat matin,Mat fact)
66*289bc588SBarry Smith {
67*289bc588SBarry Smith   return MatCholeskyFactor(fact,0);
68*289bc588SBarry Smith }
69*289bc588SBarry Smith static int MatiSDchfactor(Mat matin,IS perm)
70*289bc588SBarry Smith {
71*289bc588SBarry Smith   MatiSD    *mat = (MatiSD *) matin->data;
72*289bc588SBarry Smith   int       ierr, one  = 1, info;
73*289bc588SBarry Smith   if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;}
74*289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
75*289bc588SBarry Smith   if (info) SETERR(1,"Bad Cholesky factorization");
76*289bc588SBarry Smith   matin->factor = FACTOR_CHOLESKY;
77*289bc588SBarry Smith   return 0;
78*289bc588SBarry Smith }
79*289bc588SBarry Smith 
80*289bc588SBarry Smith static int MatiSDsolve(Mat matin,Vec xx,Vec yy)
81*289bc588SBarry Smith {
82*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
83*289bc588SBarry Smith   int i,j, one = 1, info;
84*289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
85*289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
86*289bc588SBarry Smith   MEMCPY(y,x,mat->m*sizeof(Scalar));
87*289bc588SBarry Smith   /* assume if pivots exist then LU else Cholesky */
88*289bc588SBarry Smith   if (mat->pivots) {
89*289bc588SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
90*289bc588SBarry Smith               y, &mat->m, &info );
91*289bc588SBarry Smith   }
92*289bc588SBarry Smith   else {
93*289bc588SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
94*289bc588SBarry Smith               y, &mat->m, &info );
95*289bc588SBarry Smith   }
96*289bc588SBarry Smith   if (info) SETERR(1,"Bad solve");
97*289bc588SBarry Smith   return 0;
98*289bc588SBarry Smith }
99*289bc588SBarry Smith static int MatiSDsolveadd(Mat matin,Vec xx,Vec yy,Vec zz)
100*289bc588SBarry Smith {
101*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
102*289bc588SBarry Smith   int    i;
103*289bc588SBarry Smith   Scalar *y, *z;
104*289bc588SBarry Smith   VecGetArray(yy,&y); VecGetArray(zz,&z);
105*289bc588SBarry Smith   MatiSDsolve(matin,xx,zz);
106*289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) z[i] += y[i];
107*289bc588SBarry Smith   return 0;
108*289bc588SBarry Smith }
109*289bc588SBarry Smith static int MatiSDsolvetrans(Mat matin,Vec xx,Vec yy)
110*289bc588SBarry Smith {return 0;}
111*289bc588SBarry Smith 
112*289bc588SBarry Smith static int MatiSDsolvetransadd(Mat matin,Vec xx,Vec yy,Vec zz)
113*289bc588SBarry Smith {return 0;}
114*289bc588SBarry Smith 
115*289bc588SBarry Smith /* ------------------------------------------------------------------*/
116*289bc588SBarry Smith static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,IS is,
117*289bc588SBarry Smith                        int its,Vec xx)
118*289bc588SBarry Smith {
119*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
120*289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
121*289bc588SBarry Smith   int    o = 1, tmp,n = mat->n,ierr, m = mat->m, i, j;
122*289bc588SBarry Smith 
123*289bc588SBarry Smith   if (is) SETERR(1,"No support for ordering in relaxation");
124*289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
125*289bc588SBarry Smith     /* this is a hack fix, should have another version without
126*289bc588SBarry Smith        the second BLdot */
127*289bc588SBarry Smith     if (ierr = VecSet(&zero,xx)) SETERR(ierr,0);
128*289bc588SBarry Smith   }
129*289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
130*289bc588SBarry Smith   while (its--) {
131*289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
132*289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
133*289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
134*289bc588SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/v[i + i*m] + x[i]);
135*289bc588SBarry Smith       }
136*289bc588SBarry Smith     }
137*289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
138*289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
139*289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
140*289bc588SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/v[i + i*m] + x[i]);
141*289bc588SBarry Smith       }
142*289bc588SBarry Smith     }
143*289bc588SBarry Smith   }
144*289bc588SBarry Smith   return 0;
145*289bc588SBarry Smith }
146*289bc588SBarry Smith 
147*289bc588SBarry Smith /* -----------------------------------------------------------------*/
148*289bc588SBarry Smith static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy)
149*289bc588SBarry Smith {
150*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
151*289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
152*289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
153*289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
154*289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
155*289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
156*289bc588SBarry Smith   return 0;
157*289bc588SBarry Smith }
158*289bc588SBarry Smith static int MatiSDmult(Mat matin,Vec xx,Vec yy)
159*289bc588SBarry Smith {
160*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
161*289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
162*289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
163*289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
164*289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
165*289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
166*289bc588SBarry Smith   return 0;
167*289bc588SBarry Smith }
168*289bc588SBarry Smith static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy)
169*289bc588SBarry Smith {
170*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
171*289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
172*289bc588SBarry Smith   int    _One=1; Scalar _DOne=1.0, _DZero=0.0;
173*289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
174*289bc588SBarry Smith   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
175*289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
176*289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
177*289bc588SBarry Smith   return 0;
178*289bc588SBarry Smith }
179*289bc588SBarry Smith static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy)
180*289bc588SBarry Smith {
181*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
182*289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
183*289bc588SBarry Smith   int    _One=1;
184*289bc588SBarry Smith   Scalar _DOne=1.0, _DZero=0.0;
185*289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
186*289bc588SBarry Smith   VecGetArray(zz,&z);
187*289bc588SBarry Smith   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
188*289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
189*289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
190*289bc588SBarry Smith   return 0;
191*289bc588SBarry Smith }
192*289bc588SBarry Smith 
193*289bc588SBarry Smith /* -----------------------------------------------------------------*/
194*289bc588SBarry Smith static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols,
195*289bc588SBarry Smith                         Scalar **vals)
196*289bc588SBarry Smith {
197*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
198*289bc588SBarry Smith   Scalar *v;
199*289bc588SBarry Smith   int    i;
200*289bc588SBarry Smith   *ncols = mat->n;
201*289bc588SBarry Smith   if (cols) {
202*289bc588SBarry Smith     *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols);
203*289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
204*289bc588SBarry Smith   }
205*289bc588SBarry Smith   if (vals) {
206*289bc588SBarry Smith     *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals);
207*289bc588SBarry Smith     v = mat->v + row;
208*289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
209*289bc588SBarry Smith   }
210*289bc588SBarry Smith   return 0;
211*289bc588SBarry Smith }
212*289bc588SBarry Smith static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols,
213*289bc588SBarry Smith                             Scalar **vals)
214*289bc588SBarry Smith {
215*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
216*289bc588SBarry Smith   if (cols) { FREE(*cols); }
217*289bc588SBarry Smith   if (vals) { FREE(*vals); }
218*289bc588SBarry Smith   return 0;
219*289bc588SBarry Smith }
220*289bc588SBarry Smith /* ----------------------------------------------------------------*/
221*289bc588SBarry Smith static int MatiSDinsert(Mat matin,Scalar *v,int m,int *indexm,int n,
222*289bc588SBarry Smith                         int *indexn)
223*289bc588SBarry Smith {
224*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
225*289bc588SBarry Smith   int i,j;
226*289bc588SBarry Smith   if (!mat->roworiented) {
227*289bc588SBarry Smith     for ( j=0; j<n; j++ ) {
228*289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
229*289bc588SBarry Smith         mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
230*289bc588SBarry Smith       }
231*289bc588SBarry Smith     }
232*289bc588SBarry Smith   }
233*289bc588SBarry Smith   else {
234*289bc588SBarry Smith     for ( i=0; i<m; i++ ) {
235*289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
236*289bc588SBarry Smith         mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
237*289bc588SBarry Smith       }
238*289bc588SBarry Smith     }
239*289bc588SBarry Smith   }
240*289bc588SBarry Smith   return 0;
241*289bc588SBarry Smith }
242*289bc588SBarry Smith static int MatiSDinsertadd(Mat matin,Scalar *v,int m,int *indexm,
243*289bc588SBarry Smith                            int n,int *indexn)
244*289bc588SBarry Smith {
245*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
246*289bc588SBarry Smith   int i,j;
247*289bc588SBarry Smith   if (!mat->roworiented) {
248*289bc588SBarry Smith     for ( j=0; j<n; j++ ) {
249*289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
250*289bc588SBarry Smith         mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
251*289bc588SBarry Smith       }
252*289bc588SBarry Smith     }
253*289bc588SBarry Smith   }
254*289bc588SBarry Smith   else {
255*289bc588SBarry Smith     for ( i=0; i<m; i++ ) {
256*289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
257*289bc588SBarry Smith         mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
258*289bc588SBarry Smith       }
259*289bc588SBarry Smith     }
260*289bc588SBarry Smith   }
261*289bc588SBarry Smith   return 0;
262*289bc588SBarry Smith }
263*289bc588SBarry Smith /* -----------------------------------------------------------------*/
264*289bc588SBarry Smith static int MatiSDcopy(Mat matin,Mat *newmat)
265*289bc588SBarry Smith {
266*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
267*289bc588SBarry Smith   int ierr;
268*289bc588SBarry Smith   Mat newi;
269*289bc588SBarry Smith   MatiSD *l;
270*289bc588SBarry Smith   if (ierr = MatCreateSequentialDense(mat->m,mat->n,&newi)) SETERR(ierr,0);
271*289bc588SBarry Smith   l = (MatiSD *) newi->data;
272*289bc588SBarry Smith   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
273*289bc588SBarry Smith   *newmat = newi;
274*289bc588SBarry Smith   return 0;
275*289bc588SBarry Smith }
276*289bc588SBarry Smith 
277*289bc588SBarry Smith int MatiSDview(PetscObject obj,Viewer ptr)
278*289bc588SBarry Smith {
279*289bc588SBarry Smith   Mat    matin = (Mat) obj;
280*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
281*289bc588SBarry Smith   Scalar *v;
282*289bc588SBarry Smith   int i,j;
283*289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
284*289bc588SBarry Smith     v = mat->v + i;
285*289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
286*289bc588SBarry Smith #if defined(PETSC_COMPLEX)
287*289bc588SBarry Smith       printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
288*289bc588SBarry Smith #else
289*289bc588SBarry Smith       printf("%6.4e ",*v); v += mat->m;
290*289bc588SBarry Smith #endif
291*289bc588SBarry Smith     }
292*289bc588SBarry Smith     printf("\n");
293*289bc588SBarry Smith   }
294*289bc588SBarry Smith   return 0;
295*289bc588SBarry Smith }
296*289bc588SBarry Smith 
297*289bc588SBarry Smith 
298*289bc588SBarry Smith static int MatiSDdestroy(PetscObject obj)
299*289bc588SBarry Smith {
300*289bc588SBarry Smith   Mat mat = (Mat) obj;
301*289bc588SBarry Smith   MatiSD *l = (MatiSD *) mat->data;
302*289bc588SBarry Smith   if (l->pivots) FREE(l->pivots);
303*289bc588SBarry Smith   FREE(l);
304*289bc588SBarry Smith   FREE(mat);
305*289bc588SBarry Smith   return 0;
306*289bc588SBarry Smith }
307*289bc588SBarry Smith 
308*289bc588SBarry Smith static int MatiSDtrans(Mat matin)
309*289bc588SBarry Smith {
310*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
311*289bc588SBarry Smith   int    k,j;
312*289bc588SBarry Smith   Scalar *v = mat->v, tmp;
313*289bc588SBarry Smith   if (mat->m != mat->n) {
314*289bc588SBarry Smith     SETERR(1,"Cannot transpose rectangular dense matrix");
315*289bc588SBarry Smith   }
316*289bc588SBarry Smith   for ( j=0; j<mat->m; j++ ) {
317*289bc588SBarry Smith     for ( k=0; k<j; k++ ) {
318*289bc588SBarry Smith       tmp = v[j + k*mat->n];
319*289bc588SBarry Smith       v[j + k*mat->n] = v[k + j*mat->n];
320*289bc588SBarry Smith       v[k + j*mat->n] = tmp;
321*289bc588SBarry Smith     }
322*289bc588SBarry Smith   }
323*289bc588SBarry Smith   return 0;
324*289bc588SBarry Smith }
325*289bc588SBarry Smith 
326*289bc588SBarry Smith static int MatiSDequal(Mat matin1,Mat matin2)
327*289bc588SBarry Smith {
328*289bc588SBarry Smith   MatiSD *mat1 = (MatiSD *) matin1->data;
329*289bc588SBarry Smith   MatiSD *mat2 = (MatiSD *) matin2->data;
330*289bc588SBarry Smith   int    i;
331*289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
332*289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
333*289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
334*289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
335*289bc588SBarry Smith     if (*v1 != *v2) return 0;
336*289bc588SBarry Smith     v1++; v2++;
337*289bc588SBarry Smith   }
338*289bc588SBarry Smith   return 1;
339*289bc588SBarry Smith }
340*289bc588SBarry Smith 
341*289bc588SBarry Smith static int MatiSDgetdiag(Mat matin,Vec v)
342*289bc588SBarry Smith {
343*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
344*289bc588SBarry Smith   int    i,j, n;
345*289bc588SBarry Smith   Scalar *x, zero = 0.0;
346*289bc588SBarry Smith   CHKTYPE(v,SEQVECTOR);
347*289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
348*289bc588SBarry Smith   if (n != mat->m) SETERR(1,"Nonconforming matrix and vector");
349*289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
350*289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
351*289bc588SBarry Smith   }
352*289bc588SBarry Smith   return 0;
353*289bc588SBarry Smith }
354*289bc588SBarry Smith 
355*289bc588SBarry Smith static int MatiSDscale(Mat matin,Vec l,Vec r)
356*289bc588SBarry Smith {
357*289bc588SBarry Smith return 0;
358*289bc588SBarry Smith }
359*289bc588SBarry Smith 
360*289bc588SBarry Smith static int MatiSDnorm(Mat matin,int type,double *norm)
361*289bc588SBarry Smith {
362*289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
363*289bc588SBarry Smith   Scalar *v = mat->v;
364*289bc588SBarry Smith   double sum = 0.0;
365*289bc588SBarry Smith   int    i, j;
366*289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
367*289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
368*289bc588SBarry Smith #if defined(PETSC_COMPLEX)
369*289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
370*289bc588SBarry Smith #else
371*289bc588SBarry Smith       sum += (*v)*(*v); v++;
372*289bc588SBarry Smith #endif
373*289bc588SBarry Smith     }
374*289bc588SBarry Smith     *norm = sqrt(sum);
375*289bc588SBarry Smith   }
376*289bc588SBarry Smith   else if (type == NORM_1) {
377*289bc588SBarry Smith     *norm = 0.0;
378*289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
379*289bc588SBarry Smith       sum = 0.0;
380*289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
381*289bc588SBarry Smith #if defined(PETSC_COMPLEX)
382*289bc588SBarry Smith         sum += abs(*v++);
383*289bc588SBarry Smith #else
384*289bc588SBarry Smith         sum += fabs(*v++);
385*289bc588SBarry Smith #endif
386*289bc588SBarry Smith       }
387*289bc588SBarry Smith       if (sum > *norm) *norm = sum;
388*289bc588SBarry Smith     }
389*289bc588SBarry Smith   }
390*289bc588SBarry Smith   else if (type == NORM_INFINITY) {
391*289bc588SBarry Smith     *norm = 0.0;
392*289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
393*289bc588SBarry Smith       v = mat->v + j;
394*289bc588SBarry Smith       sum = 0.0;
395*289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
396*289bc588SBarry Smith #if defined(PETSC_COMPLEX)
397*289bc588SBarry Smith         sum += abs(*v); v += mat->m;
398*289bc588SBarry Smith #else
399*289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
400*289bc588SBarry Smith #endif
401*289bc588SBarry Smith       }
402*289bc588SBarry Smith       if (sum > *norm) *norm = sum;
403*289bc588SBarry Smith     }
404*289bc588SBarry Smith   }
405*289bc588SBarry Smith   else {
406*289bc588SBarry Smith     SETERR(1,"No support for the two norm yet");
407*289bc588SBarry Smith   }
408*289bc588SBarry Smith   return 0;
409*289bc588SBarry Smith }
410*289bc588SBarry Smith 
411*289bc588SBarry Smith static int MatiDenseinsopt(Mat aijin,int op)
412*289bc588SBarry Smith {
413*289bc588SBarry Smith   MatiSD *aij = (MatiSD *) aijin->data;
414*289bc588SBarry Smith   if (op == ROW_ORIENTED)         aij->roworiented = 1;
415*289bc588SBarry Smith   else if (op == COLUMN_ORIENTED) aij->roworiented = 0;
416*289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
417*289bc588SBarry Smith   return 0;
418*289bc588SBarry Smith }
419*289bc588SBarry Smith 
420*289bc588SBarry Smith /* -------------------------------------------------------------------*/
421*289bc588SBarry Smith static struct _MatOps MatOps = {MatiSDinsert,MatiSDinsert,
422*289bc588SBarry Smith        MatiSDgetrow, MatiSDrestorerow,
423*289bc588SBarry Smith        MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd,
424*289bc588SBarry Smith        MatiSDsolve,MatiSDsolveadd,MatiSDsolvetrans,MatiSDsolvetransadd,
425*289bc588SBarry Smith        MatiSDlufactor,MatiSDchfactor,
426*289bc588SBarry Smith        MatiSDrelax,
427*289bc588SBarry Smith        MatiSDtrans,
428*289bc588SBarry Smith        MatiSDnz,MatiSDmemory,MatiSDequal,
429*289bc588SBarry Smith        MatiSDcopy,
430*289bc588SBarry Smith        MatiSDgetdiag,MatiSDscale,MatiSDnorm,
431*289bc588SBarry Smith        0,0,
432*289bc588SBarry Smith        0, MatiDenseinsopt,0,0,0,
433*289bc588SBarry Smith        MatiSDlufactorsymbolic,MatiSDlufactornumeric,
434*289bc588SBarry Smith        MatiSDchfactorsymbolic,MatiSDchfactornumeric
435*289bc588SBarry Smith };
436*289bc588SBarry Smith 
437*289bc588SBarry Smith int MatCreateSequentialDense(int m,int n,Mat *newmat)
438*289bc588SBarry Smith {
439*289bc588SBarry Smith   int       size = sizeof(MatiSD) + m*n*sizeof(Scalar);
440*289bc588SBarry Smith   Mat mat;
441*289bc588SBarry Smith   MatiSD    *l;
442*289bc588SBarry Smith   *newmat        = 0;
443*289bc588SBarry Smith   CREATEHEADER(mat,_Mat);
444*289bc588SBarry Smith   l              = (MatiSD *) MALLOC(size); CHKPTR(l);
445*289bc588SBarry Smith   mat->cookie    = MAT_COOKIE;
446*289bc588SBarry Smith   mat->ops       = &MatOps;
447*289bc588SBarry Smith   mat->destroy   = MatiSDdestroy;
448*289bc588SBarry Smith   mat->view      = MatiSDview;
449*289bc588SBarry Smith   mat->data      = (void *) l;
450*289bc588SBarry Smith   mat->type      = MATDENSESEQ;
451*289bc588SBarry Smith   mat->factor    = 0;
452*289bc588SBarry Smith   mat->col       = 0;
453*289bc588SBarry Smith   mat->row       = 0;
454*289bc588SBarry Smith   l->m           = m;
455*289bc588SBarry Smith   l->n           = n;
456*289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
457*289bc588SBarry Smith   l->pivots      = 0;
458*289bc588SBarry Smith   l->roworiented = 1;
459*289bc588SBarry Smith   MEMSET(l->v,0,m*n*sizeof(Scalar));
460*289bc588SBarry Smith   *newmat = mat;
461*289bc588SBarry Smith   return 0;
462*289bc588SBarry Smith }
463*289bc588SBarry Smith 
464*289bc588SBarry Smith int MatiSDCreate(Mat matin,Mat *newmat)
465*289bc588SBarry Smith {
466*289bc588SBarry Smith   MatiSD *m = (MatiSD *) matin->data;
467*289bc588SBarry Smith   return MatCreateSequentialDense(m->m,m->n,newmat);
468*289bc588SBarry Smith }
469