xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 6abc6512c105a871402f018420c4c4316bc1e68d)
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    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       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    one = 1, info;
84   Scalar *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 {
101   MatiSD *mat = (MatiSD *) matin->data;
102   int    one = 1, info;
103   Scalar *x, *y;
104   VecGetArray(xx,&x); VecGetArray(yy,&y);
105   MEMCPY(y,x,mat->m*sizeof(Scalar));
106   /* assume if pivots exist then LU else Cholesky */
107   if (mat->pivots) {
108     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
109               y, &mat->m, &info );
110   }
111   else {
112     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
113               y, &mat->m, &info );
114   }
115   if (info) SETERR(1,"Bad solve");
116   return 0;
117 }
118 static int MatiSDsolveadd(Mat matin,Vec xx,Vec zz,Vec yy)
119 {
120   MatiSD *mat = (MatiSD *) matin->data;
121   int    one = 1, info,ierr;
122   Scalar *x, *y, sone = 1.0;
123   Vec    tmp = 0;
124   VecGetArray(xx,&x); VecGetArray(yy,&y);
125   if (yy == zz) {
126     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
127     ierr = VecCopy(yy,tmp); CHKERR(ierr);
128   }
129   MEMCPY(y,x,mat->m*sizeof(Scalar));
130   /* assume if pivots exist then LU else Cholesky */
131   if (mat->pivots) {
132     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,
133               y, &mat->m, &info );
134   }
135   else {
136     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
137               y, &mat->m, &info );
138   }
139   if (info) SETERR(1,"Bad solve");
140   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
141   else VecAXPY(&sone,zz,yy);
142   return 0;
143 }
144 static int MatiSDsolvetransadd(Mat matin,Vec xx,Vec zz, Vec yy)
145 {
146   MatiSD  *mat = (MatiSD *) matin->data;
147   int     one = 1, info,ierr;
148   Scalar  *x, *y, sone = 1.0;
149   Vec     tmp;
150   VecGetArray(xx,&x); VecGetArray(yy,&y);
151   if (yy == zz) {
152     ierr = VecCreate(yy,&tmp); CHKERR(ierr);
153     ierr = VecCopy(yy,tmp); CHKERR(ierr);
154   }
155   MEMCPY(y,x,mat->m*sizeof(Scalar));
156   /* assume if pivots exist then LU else Cholesky */
157   if (mat->pivots) {
158     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,
159               y, &mat->m, &info );
160   }
161   else {
162     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,
163               y, &mat->m, &info );
164   }
165   if (info) SETERR(1,"Bad solve");
166   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
167   else VecAXPY(&sone,zz,yy);
168   return 0;
169 }
170 /* ------------------------------------------------------------------*/
171 static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,double shift,
172                        int its,Vec xx)
173 {
174   MatiSD *mat = (MatiSD *) matin->data;
175   Scalar *x, *b, *v = mat->v, zero = 0.0, xt;
176   int    o = 1,ierr, m = mat->m, i;
177 
178   if (flag & SOR_ZERO_INITIAL_GUESS) {
179     /* this is a hack fix, should have another version without
180        the second BLdot */
181     if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0);
182   }
183   VecGetArray(xx,&x); VecGetArray(bb,&b);
184   while (its--) {
185     if (flag & SOR_FORWARD_SWEEP){
186       for ( i=0; i<m; i++ ) {
187         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
188         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
189       }
190     }
191     if (flag & SOR_BACKWARD_SWEEP) {
192       for ( i=m-1; i>=0; i-- ) {
193         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
194         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
195       }
196     }
197   }
198   return 0;
199 }
200 
201 /* -----------------------------------------------------------------*/
202 static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy)
203 {
204   MatiSD *mat = (MatiSD *) matin->data;
205   Scalar *v = mat->v, *x, *y;
206   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
207   VecGetArray(xx,&x), VecGetArray(yy,&y);
208   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
209          x, &_One, &_DZero, y, &_One );
210   return 0;
211 }
212 static int MatiSDmult(Mat matin,Vec xx,Vec yy)
213 {
214   MatiSD *mat = (MatiSD *) matin->data;
215   Scalar *v = mat->v, *x, *y;
216   int _One=1;Scalar _DOne=1.0, _DZero=0.0;
217   VecGetArray(xx,&x); VecGetArray(yy,&y);
218   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
219          x, &_One, &_DZero, y, &_One );
220   return 0;
221 }
222 static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy)
223 {
224   MatiSD *mat = (MatiSD *) matin->data;
225   Scalar *v = mat->v, *x, *y, *z;
226   int    _One=1; Scalar _DOne=1.0;
227   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
228   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
229   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
230          x, &_One, &_DOne, y, &_One );
231   return 0;
232 }
233 static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy)
234 {
235   MatiSD *mat = (MatiSD *) matin->data;
236   Scalar *v = mat->v, *x, *y, *z;
237   int    _One=1;
238   Scalar _DOne=1.0;
239   VecGetArray(xx,&x); VecGetArray(yy,&y);
240   VecGetArray(zz,&z);
241   if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar));
242   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),
243          x, &_One, &_DOne, y, &_One );
244   return 0;
245 }
246 
247 /* -----------------------------------------------------------------*/
248 static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols,
249                         Scalar **vals)
250 {
251   MatiSD *mat = (MatiSD *) matin->data;
252   Scalar *v;
253   int    i;
254   *ncols = mat->n;
255   if (cols) {
256     *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols);
257     for ( i=0; i<mat->n; i++ ) *cols[i] = i;
258   }
259   if (vals) {
260     *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals);
261     v = mat->v + row;
262     for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;}
263   }
264   return 0;
265 }
266 static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols,
267                             Scalar **vals)
268 {
269   if (cols) { FREE(*cols); }
270   if (vals) { FREE(*vals); }
271   return 0;
272 }
273 /* ----------------------------------------------------------------*/
274 static int MatiSDinsert(Mat matin,int m,int *indexm,int n,
275                         int *indexn,Scalar *v,InsertMode addv)
276 {
277   MatiSD *mat = (MatiSD *) matin->data;
278   int    i,j;
279 
280   if (!mat->roworiented) {
281     if (addv == InsertValues) {
282       for ( j=0; j<n; j++ ) {
283         if (indexn[j] < 0) {*v += m; continue;}
284         for ( i=0; i<m; i++ ) {
285           if (indexm[i] < 0) {*v++; continue;}
286           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
287         }
288       }
289     }
290     else {
291       for ( j=0; j<n; j++ ) {
292         if (indexn[j] < 0) {*v += m; continue;}
293         for ( i=0; i<m; i++ ) {
294           if (indexm[i] < 0) {*v++; continue;}
295           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
296         }
297       }
298     }
299   }
300   else {
301     if (addv == InsertValues) {
302       for ( i=0; i<m; i++ ) {
303         if (indexm[i] < 0) {*v += n; continue;}
304         for ( j=0; j<n; j++ ) {
305           if (indexn[j] < 0) {*v++; continue;}
306           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
307         }
308       }
309     }
310     else {
311       for ( i=0; i<m; i++ ) {
312         if (indexm[i] < 0) {*v += n; continue;}
313         for ( j=0; j<n; j++ ) {
314           if (indexn[j] < 0) {*v++; continue;}
315           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
316         }
317       }
318     }
319   }
320   return 0;
321 }
322 
323 /* -----------------------------------------------------------------*/
324 static int MatiSDcopy(Mat matin,Mat *newmat)
325 {
326   MatiSD *mat = (MatiSD *) matin->data;
327   int ierr;
328   Mat newi;
329   MatiSD *l;
330   if ((ierr = MatCreateSequentialDense(mat->m,mat->n,&newi))) SETERR(ierr,0);
331   l = (MatiSD *) newi->data;
332   MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
333   *newmat = newi;
334   return 0;
335 }
336 #include "viewer.h"
337 
338 int MatiSDview(PetscObject obj,Viewer ptr)
339 {
340   Mat         matin = (Mat) obj;
341   MatiSD      *mat = (MatiSD *) matin->data;
342   Scalar      *v;
343   int         i,j;
344   PetscObject ojb = (PetscObject) ptr;
345 
346   if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) {
347     return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v);
348   }
349   else {
350     for ( i=0; i<mat->m; i++ ) {
351       v = mat->v + i;
352       for ( j=0; j<mat->n; j++ ) {
353 #if defined(PETSC_COMPLEX)
354         printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m;
355 #else
356         printf("%6.4e ",*v); v += mat->m;
357 #endif
358       }
359       printf("\n");
360     }
361   }
362   return 0;
363 }
364 
365 
366 static int MatiSDdestroy(PetscObject obj)
367 {
368   Mat mat = (Mat) obj;
369   MatiSD *l = (MatiSD *) mat->data;
370   if (l->pivots) FREE(l->pivots);
371   FREE(l);
372   FREE(mat);
373   return 0;
374 }
375 
376 static int MatiSDtrans(Mat matin)
377 {
378   MatiSD *mat = (MatiSD *) matin->data;
379   int    k,j;
380   Scalar *v = mat->v, tmp;
381   if (mat->m != mat->n) {
382     SETERR(1,"Cannot transpose rectangular dense matrix");
383   }
384   for ( j=0; j<mat->m; j++ ) {
385     for ( k=0; k<j; k++ ) {
386       tmp = v[j + k*mat->n];
387       v[j + k*mat->n] = v[k + j*mat->n];
388       v[k + j*mat->n] = tmp;
389     }
390   }
391   return 0;
392 }
393 
394 static int MatiSDequal(Mat matin1,Mat matin2)
395 {
396   MatiSD *mat1 = (MatiSD *) matin1->data;
397   MatiSD *mat2 = (MatiSD *) matin2->data;
398   int    i;
399   Scalar *v1 = mat1->v, *v2 = mat2->v;
400   if (mat1->m != mat2->m) return 0;
401   if (mat1->n != mat2->n) return 0;
402   for ( i=0; i<mat1->m*mat1->n; i++ ) {
403     if (*v1 != *v2) return 0;
404     v1++; v2++;
405   }
406   return 1;
407 }
408 
409 static int MatiSDgetdiag(Mat matin,Vec v)
410 {
411   MatiSD *mat = (MatiSD *) matin->data;
412   int    i, n;
413   Scalar *x;
414   CHKTYPE(v,SEQVECTOR);
415   VecGetArray(v,&x); VecGetSize(v,&n);
416   if (n != mat->m) SETERR(1,"Nonconforming matrix and vector");
417   for ( i=0; i<mat->m; i++ ) {
418     x[i] = mat->v[i*mat->m + i];
419   }
420   return 0;
421 }
422 
423 static int MatiSDscale(Mat matin,Vec ll,Vec rr)
424 {
425   MatiSD *mat = (MatiSD *) matin->data;
426   Scalar *l,*r,x,*v;
427   int    i,j,m = mat->m, n = mat->n;
428   if (l) {
429     VecGetArray(ll,&l); VecGetSize(ll,&m);
430     if (m != mat->m) SETERR(1,"Left scaling vector wrong length");
431     for ( i=0; i<m; i++ ) {
432       x = l[i];
433       v = mat->v + i;
434       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
435     }
436   }
437   if (l) {
438     VecGetArray(rr,&r); VecGetSize(rr,&n);
439     if (n != mat->n) SETERR(1,"Right scaling vector wrong length");
440     for ( i=0; i<n; i++ ) {
441       x = r[i];
442       v = mat->v + i*m;
443       for ( j=0; j<m; j++ ) { (*v++) *= x;}
444     }
445   }
446   return 0;
447 }
448 
449 
450 static int MatiSDnorm(Mat matin,int type,double *norm)
451 {
452   MatiSD *mat = (MatiSD *) matin->data;
453   Scalar *v = mat->v;
454   double sum = 0.0;
455   int    i, j;
456   if (type == NORM_FROBENIUS) {
457     for (i=0; i<mat->n*mat->m; i++ ) {
458 #if defined(PETSC_COMPLEX)
459       sum += real(conj(*v)*(*v)); v++;
460 #else
461       sum += (*v)*(*v); v++;
462 #endif
463     }
464     *norm = sqrt(sum);
465   }
466   else if (type == NORM_1) {
467     *norm = 0.0;
468     for ( j=0; j<mat->n; j++ ) {
469       sum = 0.0;
470       for ( i=0; i<mat->m; i++ ) {
471 #if defined(PETSC_COMPLEX)
472         sum += abs(*v++);
473 #else
474         sum += fabs(*v++);
475 #endif
476       }
477       if (sum > *norm) *norm = sum;
478     }
479   }
480   else if (type == NORM_INFINITY) {
481     *norm = 0.0;
482     for ( j=0; j<mat->m; j++ ) {
483       v = mat->v + j;
484       sum = 0.0;
485       for ( i=0; i<mat->n; i++ ) {
486 #if defined(PETSC_COMPLEX)
487         sum += abs(*v); v += mat->m;
488 #else
489         sum += fabs(*v); v += mat->m;
490 #endif
491       }
492       if (sum > *norm) *norm = sum;
493     }
494   }
495   else {
496     SETERR(1,"No support for the two norm yet");
497   }
498   return 0;
499 }
500 
501 static int MatiDenseinsopt(Mat aijin,int op)
502 {
503   MatiSD *aij = (MatiSD *) aijin->data;
504   if (op == ROW_ORIENTED)            aij->roworiented = 1;
505   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
506   /* doesn't care about sorted rows or columns */
507   return 0;
508 }
509 
510 static int MatiZero(Mat A)
511 {
512   MatiSD *l = (MatiSD *) A->data;
513   MEMSET(l->v,0,l->m*l->n*sizeof(Scalar));
514   return 0;
515 }
516 
517 static int MatiZerorows(Mat A,IS is,Scalar *diag)
518 {
519   MatiSD *l = (MatiSD *) A->data;
520   int    n = l->n, i, j,ierr,N, *rows;
521   Scalar *slot;
522   ierr = ISGetLocalSize(is,&N); CHKERR(ierr);
523   ierr = ISGetIndices(is,&rows); CHKERR(ierr);
524   for ( i=0; i<N; i++ ) {
525     slot = l->v + rows[i];
526     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
527   }
528   if (diag) {
529     for ( i=0; i<N; i++ ) {
530       slot = l->v + (n+1)*rows[i];
531       *slot = *diag;
532     }
533   }
534   ISRestoreIndices(is,&rows);
535   return 0;
536 }
537 /* -------------------------------------------------------------------*/
538 static struct _MatOps MatOps = {MatiSDinsert,
539        MatiSDgetrow, MatiSDrestorerow,
540        MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd,
541        MatiSDsolve,MatiSDsolveadd,MatiSDsolvetrans,MatiSDsolvetransadd,
542        MatiSDlufactor,MatiSDchfactor,
543        MatiSDrelax,
544        MatiSDtrans,
545        MatiSDnz,MatiSDmemory,MatiSDequal,
546        MatiSDcopy,
547        MatiSDgetdiag,MatiSDscale,MatiSDnorm,
548        0,0,
549        0, MatiDenseinsopt,MatiZero,MatiZerorows,0,
550        MatiSDlufactorsymbolic,MatiSDlufactornumeric,
551        MatiSDchfactorsymbolic,MatiSDchfactornumeric
552 };
553 /*@
554     MatCreateSequentialDense - Creates a sequential dense matrix that
555         is stored in the usual Fortran 77 manner. Many of the matrix
556         operations use the BLAS and LAPACK routines.
557 
558   Input Parameters:
559 .   m, n - the number of rows and columns in the matrix.
560 
561   Output Parameter:
562 .  newmat - the matrix created.
563 
564   Keywords: dense matrix, lapack, blas
565 @*/
566 int MatCreateSequentialDense(int m,int n,Mat *newmat)
567 {
568   int       size = sizeof(MatiSD) + m*n*sizeof(Scalar);
569   Mat mat;
570   MatiSD    *l;
571   *newmat        = 0;
572   CREATEHEADER(mat,_Mat);
573   l              = (MatiSD *) MALLOC(size); CHKPTR(l);
574   mat->cookie    = MAT_COOKIE;
575   mat->ops       = &MatOps;
576   mat->destroy   = MatiSDdestroy;
577   mat->view      = MatiSDview;
578   mat->data      = (void *) l;
579   mat->type      = MATDENSESEQ;
580   mat->factor    = 0;
581   mat->col       = 0;
582   mat->row       = 0;
583 
584   l->m           = m;
585   l->n           = n;
586   l->v           = (Scalar *) (l + 1);
587   l->pivots      = 0;
588   l->roworiented = 1;
589 
590   MEMSET(l->v,0,m*n*sizeof(Scalar));
591   *newmat = mat;
592   return 0;
593 }
594 
595 int MatiSDCreate(Mat matin,Mat *newmat)
596 {
597   MatiSD *m = (MatiSD *) matin->data;
598   return MatCreateSequentialDense(m->m,m->n,newmat);
599 }
600