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