xref: /petsc/src/mat/impls/dense/seq/dense.c (revision d6dfbf8f6aa26f18ad3ebf3c96fa582aa9714225)
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;
39289bc588SBarry Smith   int    ierr, one  = 1, 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;
52289bc588SBarry 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;
62289bc588SBarry 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;
72289bc588SBarry Smith   int       ierr, one  = 1, 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;
83289bc588SBarry Smith   int i,j, one = 1, info;
84289bc588SBarry Smith   Scalar *v = mat->v, *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)
100289bc588SBarry Smith {return 0;}
101289bc588SBarry Smith 
102289bc588SBarry Smith /* ------------------------------------------------------------------*/
103*d6dfbf8fSBarry Smith static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,double shift,
104289bc588SBarry Smith                        int its,Vec xx)
105289bc588SBarry Smith {
106289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
107289bc588SBarry Smith   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
108289bc588SBarry Smith   int    o = 1, tmp,n = mat->n,ierr, m = mat->m, i, j;
109289bc588SBarry Smith 
110289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
111289bc588SBarry Smith     /* this is a hack fix, should have another version without
112289bc588SBarry Smith        the second BLdot */
113289bc588SBarry Smith     if (ierr = VecSet(&zero,xx)) SETERR(ierr,0);
114289bc588SBarry Smith   }
115289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
116289bc588SBarry Smith   while (its--) {
117289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
118289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
119289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
120*d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
121289bc588SBarry Smith       }
122289bc588SBarry Smith     }
123289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
124289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
125289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
126*d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
127289bc588SBarry Smith       }
128289bc588SBarry Smith     }
129289bc588SBarry Smith   }
130289bc588SBarry Smith   return 0;
131289bc588SBarry Smith }
132289bc588SBarry Smith 
133289bc588SBarry Smith /* -----------------------------------------------------------------*/
134289bc588SBarry Smith static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy)
135289bc588SBarry Smith {
136289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
137289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
138289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
139289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
140289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
141289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
142289bc588SBarry Smith   return 0;
143289bc588SBarry Smith }
144289bc588SBarry Smith static int MatiSDmult(Mat matin,Vec xx,Vec yy)
145289bc588SBarry Smith {
146289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
147289bc588SBarry Smith   Scalar *v = mat->v, *x, *y;
148289bc588SBarry Smith   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
149289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
150289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
151289bc588SBarry Smith          x, &_One, &_DZero, y, &_One );
152289bc588SBarry Smith   return 0;
153289bc588SBarry Smith }
154289bc588SBarry Smith static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy)
155289bc588SBarry Smith {
156289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
157289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
158289bc588SBarry Smith   int    _One=1; Scalar _DOne=1.0, _DZero=0.0;
159289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
160289bc588SBarry Smith   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
161289bc588SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
162289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
163289bc588SBarry Smith   return 0;
164289bc588SBarry Smith }
165289bc588SBarry Smith static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy)
166289bc588SBarry Smith {
167289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
168289bc588SBarry Smith   Scalar *v = mat->v, *x, *y, *z;
169289bc588SBarry Smith   int    _One=1;
170289bc588SBarry Smith   Scalar _DOne=1.0, _DZero=0.0;
171289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
172289bc588SBarry Smith   VecGetArray(zz,&z);
173289bc588SBarry Smith   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
174289bc588SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
175289bc588SBarry Smith          x, &_One, &_DOne, y, &_One );
176289bc588SBarry Smith   return 0;
177289bc588SBarry Smith }
178289bc588SBarry Smith 
179289bc588SBarry Smith /* -----------------------------------------------------------------*/
180289bc588SBarry Smith static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols,
181289bc588SBarry Smith                         Scalar **vals)
182289bc588SBarry Smith {
183289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
184289bc588SBarry Smith   Scalar *v;
185289bc588SBarry Smith   int    i;
186289bc588SBarry Smith   *ncols = mat->n;
187289bc588SBarry Smith   if (cols) {
188289bc588SBarry Smith     *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols);
189289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
190289bc588SBarry Smith   }
191289bc588SBarry Smith   if (vals) {
192289bc588SBarry Smith     *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals);
193289bc588SBarry Smith     v = mat->v + row;
194289bc588SBarry Smith     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
195289bc588SBarry Smith   }
196289bc588SBarry Smith   return 0;
197289bc588SBarry Smith }
198289bc588SBarry Smith static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols,
199289bc588SBarry Smith                             Scalar **vals)
200289bc588SBarry Smith {
201289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
202289bc588SBarry Smith   if (cols) { FREE(*cols); }
203289bc588SBarry Smith   if (vals) { FREE(*vals); }
204289bc588SBarry Smith   return 0;
205289bc588SBarry Smith }
206289bc588SBarry Smith /* ----------------------------------------------------------------*/
207e8d4e0b9SBarry Smith static int MatiSDinsert(Mat matin,int m,int *indexm,int n,
208e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
209289bc588SBarry Smith {
210289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
211289bc588SBarry Smith   int    i,j;
212*d6dfbf8fSBarry Smith 
213289bc588SBarry Smith   if (!mat->roworiented) {
214e8d4e0b9SBarry Smith     if (addv == InsertValues) {
215289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
216*d6dfbf8fSBarry Smith         if (indexn[j] < 0) {*v += m; continue;}
217289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
218*d6dfbf8fSBarry Smith           if (indexm[i] < 0) {*v++; continue;}
219289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
220289bc588SBarry Smith         }
221289bc588SBarry Smith       }
222289bc588SBarry Smith     }
223289bc588SBarry Smith     else {
224289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
225*d6dfbf8fSBarry Smith         if (indexn[j] < 0) {*v += m; continue;}
226289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
227*d6dfbf8fSBarry Smith           if (indexm[i] < 0) {*v++; continue;}
228289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
229289bc588SBarry Smith         }
230289bc588SBarry Smith       }
231289bc588SBarry Smith     }
232e8d4e0b9SBarry Smith   }
233e8d4e0b9SBarry Smith   else {
234e8d4e0b9SBarry Smith     if (addv == InsertValues) {
235e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
236*d6dfbf8fSBarry Smith         if (indexm[i] < 0) {*v += n; continue;}
237e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
238*d6dfbf8fSBarry Smith           if (indexn[j] < 0) {*v++; continue;}
239e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
240e8d4e0b9SBarry Smith         }
241e8d4e0b9SBarry Smith       }
242e8d4e0b9SBarry Smith     }
243289bc588SBarry Smith     else {
244289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
245*d6dfbf8fSBarry Smith         if (indexm[i] < 0) {*v += n; continue;}
246289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
247*d6dfbf8fSBarry Smith           if (indexn[j] < 0) {*v++; continue;}
248289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
249289bc588SBarry Smith         }
250289bc588SBarry Smith       }
251289bc588SBarry Smith     }
252e8d4e0b9SBarry Smith   }
253289bc588SBarry Smith   return 0;
254289bc588SBarry Smith }
255e8d4e0b9SBarry Smith 
256289bc588SBarry Smith /* -----------------------------------------------------------------*/
257289bc588SBarry Smith static int MatiSDcopy(Mat matin,Mat *newmat)
258289bc588SBarry Smith {
259289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
260289bc588SBarry Smith   int ierr;
261289bc588SBarry Smith   Mat newi;
262289bc588SBarry Smith   MatiSD *l;
263289bc588SBarry Smith   if (ierr = MatCreateSequentialDense(mat->m,mat->n,&newi)) SETERR(ierr,0);
264289bc588SBarry Smith   l = (MatiSD *) newi->data;
265289bc588SBarry Smith   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
266289bc588SBarry Smith   *newmat = newi;
267289bc588SBarry Smith   return 0;
268289bc588SBarry Smith }
269289bc588SBarry Smith 
270289bc588SBarry Smith int MatiSDview(PetscObject obj,Viewer ptr)
271289bc588SBarry Smith {
272289bc588SBarry Smith   Mat    matin = (Mat) obj;
273289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
274289bc588SBarry Smith   Scalar *v;
275289bc588SBarry Smith   int i,j;
276289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
277289bc588SBarry Smith     v = mat->v + i;
278289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
279289bc588SBarry Smith #if defined(PETSC_COMPLEX)
280289bc588SBarry Smith       printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
281289bc588SBarry Smith #else
282289bc588SBarry Smith       printf("%6.4e ",*v); v += mat->m;
283289bc588SBarry Smith #endif
284289bc588SBarry Smith     }
285289bc588SBarry Smith     printf("\n");
286289bc588SBarry Smith   }
287289bc588SBarry Smith   return 0;
288289bc588SBarry Smith }
289289bc588SBarry Smith 
290289bc588SBarry Smith 
291289bc588SBarry Smith static int MatiSDdestroy(PetscObject obj)
292289bc588SBarry Smith {
293289bc588SBarry Smith   Mat mat = (Mat) obj;
294289bc588SBarry Smith   MatiSD *l = (MatiSD *) mat->data;
295289bc588SBarry Smith   if (l->pivots) FREE(l->pivots);
296289bc588SBarry Smith   FREE(l);
297289bc588SBarry Smith   FREE(mat);
298289bc588SBarry Smith   return 0;
299289bc588SBarry Smith }
300289bc588SBarry Smith 
301289bc588SBarry Smith static int MatiSDtrans(Mat matin)
302289bc588SBarry Smith {
303289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
304289bc588SBarry Smith   int    k,j;
305289bc588SBarry Smith   Scalar *v = mat->v, tmp;
306289bc588SBarry Smith   if (mat->m != mat->n) {
307289bc588SBarry Smith     SETERR(1,"Cannot transpose rectangular dense matrix");
308289bc588SBarry Smith   }
309289bc588SBarry Smith   for ( j=0; j<mat->m; j++ ) {
310289bc588SBarry Smith     for ( k=0; k<j; k++ ) {
311289bc588SBarry Smith       tmp = v[j + k*mat->n];
312289bc588SBarry Smith       v[j + k*mat->n] = v[k + j*mat->n];
313289bc588SBarry Smith       v[k + j*mat->n] = tmp;
314289bc588SBarry Smith     }
315289bc588SBarry Smith   }
316289bc588SBarry Smith   return 0;
317289bc588SBarry Smith }
318289bc588SBarry Smith 
319289bc588SBarry Smith static int MatiSDequal(Mat matin1,Mat matin2)
320289bc588SBarry Smith {
321289bc588SBarry Smith   MatiSD *mat1 = (MatiSD *) matin1->data;
322289bc588SBarry Smith   MatiSD *mat2 = (MatiSD *) matin2->data;
323289bc588SBarry Smith   int    i;
324289bc588SBarry Smith   Scalar *v1 = mat1->v, *v2 = mat2->v;
325289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
326289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
327289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
328289bc588SBarry Smith     if (*v1 != *v2) return 0;
329289bc588SBarry Smith     v1++; v2++;
330289bc588SBarry Smith   }
331289bc588SBarry Smith   return 1;
332289bc588SBarry Smith }
333289bc588SBarry Smith 
334289bc588SBarry Smith static int MatiSDgetdiag(Mat matin,Vec v)
335289bc588SBarry Smith {
336289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
337289bc588SBarry Smith   int    i,j, n;
338289bc588SBarry Smith   Scalar *x, zero = 0.0;
339289bc588SBarry Smith   CHKTYPE(v,SEQVECTOR);
340289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
341289bc588SBarry Smith   if (n != mat->m) SETERR(1,"Nonconforming matrix and vector");
342289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
343289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
344289bc588SBarry Smith   }
345289bc588SBarry Smith   return 0;
346289bc588SBarry Smith }
347289bc588SBarry Smith 
348289bc588SBarry Smith static int MatiSDscale(Mat matin,Vec l,Vec r)
349289bc588SBarry Smith {
350289bc588SBarry Smith return 0;
351289bc588SBarry Smith }
352289bc588SBarry Smith 
353289bc588SBarry Smith static int MatiSDnorm(Mat matin,int type,double *norm)
354289bc588SBarry Smith {
355289bc588SBarry Smith   MatiSD *mat = (MatiSD *) matin->data;
356289bc588SBarry Smith   Scalar *v = mat->v;
357289bc588SBarry Smith   double sum = 0.0;
358289bc588SBarry Smith   int    i, j;
359289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
360289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
361289bc588SBarry Smith #if defined(PETSC_COMPLEX)
362289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
363289bc588SBarry Smith #else
364289bc588SBarry Smith       sum += (*v)*(*v); v++;
365289bc588SBarry Smith #endif
366289bc588SBarry Smith     }
367289bc588SBarry Smith     *norm = sqrt(sum);
368289bc588SBarry Smith   }
369289bc588SBarry Smith   else if (type == NORM_1) {
370289bc588SBarry Smith     *norm = 0.0;
371289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
372289bc588SBarry Smith       sum = 0.0;
373289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
374289bc588SBarry Smith #if defined(PETSC_COMPLEX)
375289bc588SBarry Smith         sum += abs(*v++);
376289bc588SBarry Smith #else
377289bc588SBarry Smith         sum += fabs(*v++);
378289bc588SBarry Smith #endif
379289bc588SBarry Smith       }
380289bc588SBarry Smith       if (sum > *norm) *norm = sum;
381289bc588SBarry Smith     }
382289bc588SBarry Smith   }
383289bc588SBarry Smith   else if (type == NORM_INFINITY) {
384289bc588SBarry Smith     *norm = 0.0;
385289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
386289bc588SBarry Smith       v = mat->v + j;
387289bc588SBarry Smith       sum = 0.0;
388289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
389289bc588SBarry Smith #if defined(PETSC_COMPLEX)
390289bc588SBarry Smith         sum += abs(*v); v += mat->m;
391289bc588SBarry Smith #else
392289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
393289bc588SBarry Smith #endif
394289bc588SBarry Smith       }
395289bc588SBarry Smith       if (sum > *norm) *norm = sum;
396289bc588SBarry Smith     }
397289bc588SBarry Smith   }
398289bc588SBarry Smith   else {
399289bc588SBarry Smith     SETERR(1,"No support for the two norm yet");
400289bc588SBarry Smith   }
401289bc588SBarry Smith   return 0;
402289bc588SBarry Smith }
403289bc588SBarry Smith 
404289bc588SBarry Smith static int MatiDenseinsopt(Mat aijin,int op)
405289bc588SBarry Smith {
406289bc588SBarry Smith   MatiSD *aij = (MatiSD *) aijin->data;
407289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
408289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
409289bc588SBarry Smith   /* doesn't care about sorted rows or columns */
410289bc588SBarry Smith   return 0;
411289bc588SBarry Smith }
412289bc588SBarry Smith 
4136f0a148fSBarry Smith static int MatiZero(Mat A)
4146f0a148fSBarry Smith {
4156f0a148fSBarry Smith   MatiSD *l = (MatiSD *) A->data;
4166f0a148fSBarry Smith   MEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
4176f0a148fSBarry Smith   return 0;
4186f0a148fSBarry Smith }
4196f0a148fSBarry Smith 
420260925b8SBarry Smith static int MatiZerorows(Mat A,IS is,Scalar *diag)
4216f0a148fSBarry Smith {
4226f0a148fSBarry Smith   MatiSD *l = (MatiSD *) A->data;
423260925b8SBarry Smith   int     m = l->m, n = l->n, i, j,ierr,N, *rows;
4246f0a148fSBarry Smith   Scalar  *slot;
425260925b8SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERR(ierr);
426260925b8SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERR(ierr);
4276f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
4286f0a148fSBarry Smith     slot = l->v + rows[i];
4296f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
4306f0a148fSBarry Smith   }
4316f0a148fSBarry Smith   if (diag) {
4326f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
4336f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
434260925b8SBarry Smith       *slot = *diag;
4356f0a148fSBarry Smith     }
4366f0a148fSBarry Smith   }
437260925b8SBarry Smith   ISRestoreIndices(is,&rows);
4386f0a148fSBarry Smith   return 0;
4396f0a148fSBarry Smith }
440289bc588SBarry Smith /* -------------------------------------------------------------------*/
441e8d4e0b9SBarry Smith static struct _MatOps MatOps = {MatiSDinsert,
442289bc588SBarry Smith        MatiSDgetrow, MatiSDrestorerow,
443289bc588SBarry Smith        MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd,
4448c37ef55SBarry Smith        MatiSDsolve,0,MatiSDsolvetrans,0,
445289bc588SBarry Smith        MatiSDlufactor,MatiSDchfactor,
446289bc588SBarry Smith        MatiSDrelax,
447289bc588SBarry Smith        MatiSDtrans,
448289bc588SBarry Smith        MatiSDnz,MatiSDmemory,MatiSDequal,
449289bc588SBarry Smith        MatiSDcopy,
450289bc588SBarry Smith        MatiSDgetdiag,MatiSDscale,MatiSDnorm,
451289bc588SBarry Smith        0,0,
4526f0a148fSBarry Smith        0, MatiDenseinsopt,MatiZero,MatiZerorows,0,
453289bc588SBarry Smith        MatiSDlufactorsymbolic,MatiSDlufactornumeric,
454289bc588SBarry Smith        MatiSDchfactorsymbolic,MatiSDchfactornumeric
455289bc588SBarry Smith };
45620563c6bSBarry Smith /*@
45720563c6bSBarry Smith     MatCreateSequentialDense - Creates a sequential dense matrix that
45820563c6bSBarry Smith         is stored in the usual Fortran 77 manner. Many of the matrix
45920563c6bSBarry Smith         operations use the BLAS and LAPACK routines.
460289bc588SBarry Smith 
46120563c6bSBarry Smith   Input Parameters:
46220563c6bSBarry Smith .   m, n - the number of rows and columns in the matrix.
46320563c6bSBarry Smith 
46420563c6bSBarry Smith   Output Parameter:
46520563c6bSBarry Smith .  newmat - the matrix created.
46620563c6bSBarry Smith 
46720563c6bSBarry Smith   Keywords: dense matrix, lapack, blas
46820563c6bSBarry Smith @*/
469289bc588SBarry Smith int MatCreateSequentialDense(int m,int n,Mat *newmat)
470289bc588SBarry Smith {
471289bc588SBarry Smith   int       size = sizeof(MatiSD) + m*n*sizeof(Scalar);
472289bc588SBarry Smith   Mat mat;
473289bc588SBarry Smith   MatiSD    *l;
474289bc588SBarry Smith   *newmat        = 0;
475289bc588SBarry Smith   CREATEHEADER(mat,_Mat);
476289bc588SBarry Smith   l              = (MatiSD *) MALLOC(size); CHKPTR(l);
477289bc588SBarry Smith   mat->cookie    = MAT_COOKIE;
478289bc588SBarry Smith   mat->ops       = &MatOps;
479289bc588SBarry Smith   mat->destroy   = MatiSDdestroy;
480289bc588SBarry Smith   mat->view      = MatiSDview;
481289bc588SBarry Smith   mat->data      = (void *) l;
482289bc588SBarry Smith   mat->type      = MATDENSESEQ;
483289bc588SBarry Smith   mat->factor    = 0;
484289bc588SBarry Smith   mat->col       = 0;
485289bc588SBarry Smith   mat->row       = 0;
486*d6dfbf8fSBarry Smith   mat->outofrange= 0;
487*d6dfbf8fSBarry Smith   mat->Mlow      = 0;
488*d6dfbf8fSBarry Smith   mat->Mhigh     = m;
489*d6dfbf8fSBarry Smith   mat->Nlow      = 0;
490*d6dfbf8fSBarry Smith   mat->Nhigh     = n;
491*d6dfbf8fSBarry Smith 
492289bc588SBarry Smith   l->m           = m;
493289bc588SBarry Smith   l->n           = n;
494289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
495289bc588SBarry Smith   l->pivots      = 0;
496289bc588SBarry Smith   l->roworiented = 1;
497*d6dfbf8fSBarry Smith 
498289bc588SBarry Smith   MEMSET(l->v,0,m*n*sizeof(Scalar));
499289bc588SBarry Smith   *newmat = mat;
500289bc588SBarry Smith   return 0;
501289bc588SBarry Smith }
502289bc588SBarry Smith 
503289bc588SBarry Smith int MatiSDCreate(Mat matin,Mat *newmat)
504289bc588SBarry Smith {
505289bc588SBarry Smith   MatiSD *m = (MatiSD *) matin->data;
506289bc588SBarry Smith   return MatCreateSequentialDense(m->m,m->n,newmat);
507289bc588SBarry Smith }
508