xref: /petsc/src/mat/impls/dense/seq/dense.c (revision e78daa7df64465c5433a8e915d2a72f3bcad4f37)
1 #ifndef lint
2 static char vcid[] = "$Id: dense.c,v 1.111 1996/09/23 18:21:03 bsmith Exp bsmith $";
3 #endif
4 /*
5      Defines the basic matrix operations for sequential dense.
6 */
7 
8 #include "src/mat/impls/dense/seq/dense.h"
9 #include "pinclude/plapack.h"
10 #include "pinclude/pviewer.h"
11 
12 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
13 {
14   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
15   int          N = x->m*x->n, one = 1;
16   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
17   PLogFlops(2*N-1);
18   return 0;
19 }
20 
21 static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
22 {
23   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
24   int          i,N = a->m*a->n,count = 0;
25   Scalar       *v = a->v;
26   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
27 
28   info->rows_global       = (double)a->m;
29   info->columns_global    = (double)a->n;
30   info->rows_local        = (double)a->m;
31   info->columns_local     = (double)a->n;
32   info->block_size        = 1.0;
33   info->nz_allocated      = (double)N;
34   info->nz_used           = (double)count;
35   info->nz_unneeded       = (double)(N-count);
36   info->assemblies        = (double)A->num_ass;
37   info->mallocs           = 0;
38   info->memory            = A->mem;
39   info->fill_ratio_given  = 0;
40   info->fill_ratio_needed = 0;
41   info->factor_mallocs    = 0;
42 
43   return 0;
44 }
45 
46 static int MatScale_SeqDense(Scalar *alpha,Mat inA)
47 {
48   Mat_SeqDense *a = (Mat_SeqDense *) inA->data;
49   int          one = 1, nz;
50 
51   nz = a->m*a->n;
52   BLscal_( &nz, alpha, a->v, &one );
53   PLogFlops(nz);
54   return 0;
55 }
56 
57 /* ---------------------------------------------------------------*/
58 /* COMMENT: I have chosen to hide column permutation in the pivots,
59    rather than put it in the Mat->col slot.*/
60 static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
61 {
62   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
63   int          info;
64 
65   if (!mat->pivots) {
66     mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
67     PLogObjectMemory(A,mat->m*sizeof(int));
68   }
69   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
70   if (info<0) SETERRQ(1,"MatLUFactor_SeqDense:Bad argument to LU factorization");
71   if (info>0) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
72   A->factor = FACTOR_LU;
73   PLogFlops((2*mat->n*mat->n*mat->n)/3);
74   return 0;
75 }
76 static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
77 {
78   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
79   int          ierr;
80   Mat          newi;
81 
82   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
83   l = (Mat_SeqDense *) newi->data;
84   if (cpvalues == COPY_VALUES) {
85     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
86   }
87   newi->assembled = PETSC_TRUE;
88   *newmat = newi;
89   return 0;
90 }
91 
92 static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
93 {
94   return MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);
95 }
96 static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
97 {
98   Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data;
99   /* copy the numerical values */
100   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
101   (*fact)->factor = 0;
102   return MatLUFactor(*fact,0,0,1.0);
103 }
104 static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
105 {
106   return MatConvert(A,MATSAME,fact);
107 }
108 static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
109 {
110   return MatCholeskyFactor(*fact,0,1.0);
111 }
112 static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
113 {
114   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
115   int           info;
116 
117   if (mat->pivots) {
118     PetscFree(mat->pivots);
119     PLogObjectMemory(A,-mat->m*sizeof(int));
120     mat->pivots = 0;
121   }
122   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
123   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
124   A->factor = FACTOR_CHOLESKY;
125   PLogFlops((mat->n*mat->n*mat->n)/3);
126   return 0;
127 }
128 
129 static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
130 {
131   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
132   int          one = 1, info, ierr;
133   Scalar       *x, *y;
134 
135   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
136   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
137   if (A->factor == FACTOR_LU) {
138     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
139   }
140   else if (A->factor == FACTOR_CHOLESKY){
141     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
142   }
143   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
144   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
145   PLogFlops(mat->n*mat->n - mat->n);
146   return 0;
147 }
148 static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
149 {
150   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
151   int          one = 1, info;
152   Scalar       *x, *y;
153 
154   VecGetArray(xx,&x); VecGetArray(yy,&y);
155   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
156   /* assume if pivots exist then use LU; else Cholesky */
157   if (mat->pivots) {
158     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
159   }
160   else {
161     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
162   }
163   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
164   PLogFlops(mat->n*mat->n - mat->n);
165   return 0;
166 }
167 static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
168 {
169   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
170   int          one = 1, info,ierr;
171   Scalar       *x, *y, sone = 1.0;
172   Vec          tmp = 0;
173 
174   VecGetArray(xx,&x); VecGetArray(yy,&y);
175   if (yy == zz) {
176     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
177     PLogObjectParent(A,tmp);
178     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
179   }
180   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
181   /* assume if pivots exist then use LU; else Cholesky */
182   if (mat->pivots) {
183     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
184   }
185   else {
186     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
187   }
188   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
189   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
190   else VecAXPY(&sone,zz,yy);
191   PLogFlops(mat->n*mat->n - mat->n);
192   return 0;
193 }
194 
195 static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
196 {
197   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
198   int           one = 1, info,ierr;
199   Scalar        *x, *y, sone = 1.0;
200   Vec           tmp;
201 
202   VecGetArray(xx,&x); VecGetArray(yy,&y);
203   if (yy == zz) {
204     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
205     PLogObjectParent(A,tmp);
206     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
207   }
208   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
209   /* assume if pivots exist then use LU; else Cholesky */
210   if (mat->pivots) {
211     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
212   }
213   else {
214     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
215   }
216   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
217   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
218   else VecAXPY(&sone,zz,yy);
219   PLogFlops(mat->n*mat->n - mat->n);
220   return 0;
221 }
222 /* ------------------------------------------------------------------*/
223 static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
224                           double shift,int its,Vec xx)
225 {
226   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
227   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
228   int          o = 1,ierr, m = mat->m, i;
229 
230   if (flag & SOR_ZERO_INITIAL_GUESS) {
231     /* this is a hack fix, should have another version without
232        the second BLdot */
233     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
234   }
235   VecGetArray(xx,&x); VecGetArray(bb,&b);
236   while (its--) {
237     if (flag & SOR_FORWARD_SWEEP){
238       for ( i=0; i<m; i++ ) {
239 #if defined(PETSC_COMPLEX)
240         /* cannot use BLAS dot for complex because compiler/linker is
241            not happy about returning a double complex */
242         int    _i;
243         Scalar sum = b[i];
244         for ( _i=0; _i<m; _i++ ) {
245           sum -= conj(v[i+_i*m])*x[_i];
246         }
247         xt = sum;
248 #else
249         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
250 #endif
251         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
252       }
253     }
254     if (flag & SOR_BACKWARD_SWEEP) {
255       for ( i=m-1; i>=0; i-- ) {
256 #if defined(PETSC_COMPLEX)
257         /* cannot use BLAS dot for complex because compiler/linker is
258            not happy about returning a double complex */
259         int    _i;
260         Scalar sum = b[i];
261         for ( _i=0; _i<m; _i++ ) {
262           sum -= conj(v[i+_i*m])*x[_i];
263         }
264         xt = sum;
265 #else
266         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
267 #endif
268         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
269       }
270     }
271   }
272   return 0;
273 }
274 
275 /* -----------------------------------------------------------------*/
276 int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
277 {
278   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
279   Scalar       *v = mat->v, *x, *y;
280   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
281   VecGetArray(xx,&x), VecGetArray(yy,&y);
282   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
283   PLogFlops(2*mat->m*mat->n - mat->n);
284   return 0;
285 }
286 int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
287 {
288   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
289   Scalar       *v = mat->v, *x, *y;
290   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
291   VecGetArray(xx,&x); VecGetArray(yy,&y);
292   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
293   PLogFlops(2*mat->m*mat->n - mat->m);
294   return 0;
295 }
296 int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
297 {
298   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
299   Scalar       *v = mat->v, *x, *y, *z;
300   int          _One=1; Scalar _DOne=1.0;
301   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
302   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
303   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
304   PLogFlops(2*mat->m*mat->n);
305   return 0;
306 }
307 int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
308 {
309   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
310   Scalar       *v = mat->v, *x, *y, *z;
311   int          _One=1;
312   Scalar       _DOne=1.0;
313   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
314   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
315   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
316   PLogFlops(2*mat->m*mat->n);
317   return 0;
318 }
319 
320 /* -----------------------------------------------------------------*/
321 static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
322 {
323   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
324   Scalar       *v;
325   int          i;
326 
327   *ncols = mat->n;
328   if (cols) {
329     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
330     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
331   }
332   if (vals) {
333     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
334     v = mat->v + row;
335     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
336   }
337   return 0;
338 }
339 static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
340 {
341   if (cols) { PetscFree(*cols); }
342   if (vals) { PetscFree(*vals); }
343   return 0;
344 }
345 /* ----------------------------------------------------------------*/
346 static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
347                                     int *indexn,Scalar *v,InsertMode addv)
348 {
349   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
350   int          i,j;
351 
352   if (!mat->roworiented) {
353     if (addv == INSERT_VALUES) {
354       for ( j=0; j<n; j++ ) {
355         for ( i=0; i<m; i++ ) {
356           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
357         }
358       }
359     }
360     else {
361       for ( j=0; j<n; j++ ) {
362         for ( i=0; i<m; i++ ) {
363           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
364         }
365       }
366     }
367   }
368   else {
369     if (addv == INSERT_VALUES) {
370       for ( i=0; i<m; i++ ) {
371         for ( j=0; j<n; j++ ) {
372           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
373         }
374       }
375     }
376     else {
377       for ( i=0; i<m; i++ ) {
378         for ( j=0; j<n; j++ ) {
379           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
380         }
381       }
382     }
383   }
384   return 0;
385 }
386 
387 static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
388 {
389   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
390   int          i, j;
391   Scalar       *vpt = v;
392 
393   /* row-oriented output */
394   for ( i=0; i<m; i++ ) {
395     for ( j=0; j<n; j++ ) {
396       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
397     }
398   }
399   return 0;
400 }
401 
402 /* -----------------------------------------------------------------*/
403 
404 #include "sys.h"
405 
406 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
407 {
408   Mat_SeqDense *a;
409   Mat          B;
410   int          *scols, i, j, nz, ierr, fd, header[4], size;
411   int          *rowlengths = 0, M, N, *cols;
412   Scalar       *vals, *svals, *v;
413   MPI_Comm     comm = ((PetscObject)viewer)->comm;
414 
415   MPI_Comm_size(comm,&size);
416   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
417   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
418   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
419   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
420   M = header[1]; N = header[2]; nz = header[3];
421 
422   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
423     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
424     B = *A;
425     a = (Mat_SeqDense *) B->data;
426 
427     /* read in nonzero values */
428     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
429 
430     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
431     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
432     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
433   } else {
434     /* read row lengths */
435     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
436     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
437 
438     /* create our matrix */
439     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
440     B = *A;
441     a = (Mat_SeqDense *) B->data;
442     v = a->v;
443 
444     /* read column indices and nonzeros */
445     cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
446     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
447     vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
448     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
449 
450     /* insert into matrix */
451     for ( i=0; i<M; i++ ) {
452       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
453       svals += rowlengths[i]; scols += rowlengths[i];
454     }
455     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
456 
457     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
458     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
459   }
460   return 0;
461 }
462 
463 #include "pinclude/pviewer.h"
464 #include "sys.h"
465 
466 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
467 {
468   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
469   int          ierr, i, j, format;
470   FILE         *fd;
471   char         *outputname;
472   Scalar       *v;
473 
474   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
475   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
476   ierr = ViewerGetFormat(viewer,&format);
477   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_DETAILED) {
478     return 0;  /* do nothing for now */
479   }
480   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
481     for ( i=0; i<a->m; i++ ) {
482       v = a->v + i;
483       fprintf(fd,"row %d:",i);
484       for ( j=0; j<a->n; j++ ) {
485 #if defined(PETSC_COMPLEX)
486         if (real(*v) != 0.0 && imag(*v) != 0.0)
487           fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
488         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
489         v += a->m;
490 #else
491         if (*v) fprintf(fd," %d %g ",j,*v);
492         v += a->m;
493 #endif
494       }
495       fprintf(fd,"\n");
496     }
497   }
498   else {
499 #if defined(PETSC_COMPLEX)
500     int allreal = 1;
501     /* determine if matrix has all real values */
502     v = a->v;
503     for ( i=0; i<a->m*a->n; i++ ) {
504       if (imag(v[i])) { allreal = 0; break ;}
505     }
506 #endif
507     for ( i=0; i<a->m; i++ ) {
508       v = a->v + i;
509       for ( j=0; j<a->n; j++ ) {
510 #if defined(PETSC_COMPLEX)
511         if (allreal) {
512           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
513         } else {
514           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
515         }
516 #else
517         fprintf(fd,"%6.4e ",*v); v += a->m;
518 #endif
519       }
520       fprintf(fd,"\n");
521     }
522   }
523   fflush(fd);
524   return 0;
525 }
526 
527 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
528 {
529   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
530   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
531   int          format;
532   Scalar       *v, *anonz,*vals;
533 
534   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
535 
536   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
537   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
538     /* store the matrix as a dense matrix */
539     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
540     col_lens[0] = MAT_COOKIE;
541     col_lens[1] = m;
542     col_lens[2] = n;
543     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
544     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
545     PetscFree(col_lens);
546 
547     /* write out matrix, by rows */
548     vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals);
549     v    = a->v;
550     for ( i=0; i<m; i++ ) {
551       for ( j=0; j<n; j++ ) {
552         vals[i + j*m] = *v++;
553       }
554     }
555     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
556     PetscFree(vals);
557   } else {
558     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
559     col_lens[0] = MAT_COOKIE;
560     col_lens[1] = m;
561     col_lens[2] = n;
562     col_lens[3] = nz;
563 
564     /* store lengths of each row and write (including header) to file */
565     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
566     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
567 
568     /* Possibly should write in smaller increments, not whole matrix at once? */
569     /* store column indices (zero start index) */
570     ict = 0;
571     for ( i=0; i<m; i++ ) {
572       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
573     }
574     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
575     PetscFree(col_lens);
576 
577     /* store nonzero values */
578     anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
579     ict = 0;
580     for ( i=0; i<m; i++ ) {
581       v = a->v + i;
582       for ( j=0; j<n; j++ ) {
583         anonz[ict++] = *v; v += a->m;
584       }
585     }
586     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
587     PetscFree(anonz);
588   }
589   return 0;
590 }
591 
592 static int MatView_SeqDense(PetscObject obj,Viewer viewer)
593 {
594   Mat          A = (Mat) obj;
595   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
596   ViewerType   vtype;
597   int          ierr;
598 
599   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
600 
601   if (vtype == MATLAB_VIEWER) {
602     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
603   }
604   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
605     return MatView_SeqDense_ASCII(A,viewer);
606   }
607   else if (vtype == BINARY_FILE_VIEWER) {
608     return MatView_SeqDense_Binary(A,viewer);
609   }
610   return 0;
611 }
612 
613 static int MatDestroy_SeqDense(PetscObject obj)
614 {
615   Mat          mat = (Mat) obj;
616   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
617 #if defined(PETSC_LOG)
618   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
619 #endif
620   if (l->pivots) PetscFree(l->pivots);
621   if (!l->user_alloc) PetscFree(l->v);
622   PetscFree(l);
623   PLogObjectDestroy(mat);
624   PetscHeaderDestroy(mat);
625   return 0;
626 }
627 
628 static int MatTranspose_SeqDense(Mat A,Mat *matout)
629 {
630   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
631   int          k, j, m, n;
632   Scalar       *v, tmp;
633 
634   v = mat->v; m = mat->m; n = mat->n;
635   if (matout == PETSC_NULL) { /* in place transpose */
636     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
637     for ( j=0; j<m; j++ ) {
638       for ( k=0; k<j; k++ ) {
639         tmp = v[j + k*n];
640         v[j + k*n] = v[k + j*n];
641         v[k + j*n] = tmp;
642       }
643     }
644   }
645   else { /* out-of-place transpose */
646     int          ierr;
647     Mat          tmat;
648     Mat_SeqDense *tmatd;
649     Scalar       *v2;
650     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
651     tmatd = (Mat_SeqDense *) tmat->data;
652     v = mat->v; v2 = tmatd->v;
653     for ( j=0; j<n; j++ ) {
654       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
655     }
656     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
657     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
658     *matout = tmat;
659   }
660   return 0;
661 }
662 
663 static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
664 {
665   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
666   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
667   int          i;
668   Scalar       *v1 = mat1->v, *v2 = mat2->v;
669 
670   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
671   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
672   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
673   for ( i=0; i<mat1->m*mat1->n; i++ ) {
674     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
675     v1++; v2++;
676   }
677   *flg = PETSC_TRUE;
678   return 0;
679 }
680 
681 static int MatGetDiagonal_SeqDense(Mat A,Vec v)
682 {
683   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
684   int          i, n, len;
685   Scalar       *x, zero = 0.0;
686 
687   VecSet(&zero,v);
688   VecGetArray(v,&x); VecGetSize(v,&n);
689   len = PetscMin(mat->m,mat->n);
690   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
691   for ( i=0; i<len; i++ ) {
692     x[i] = mat->v[i*mat->m + i];
693   }
694   return 0;
695 }
696 
697 static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
698 {
699   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
700   Scalar       *l,*r,x,*v;
701   int          i,j,m = mat->m, n = mat->n;
702 
703   if (ll) {
704     VecGetArray(ll,&l); VecGetSize(ll,&m);
705     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
706     for ( i=0; i<m; i++ ) {
707       x = l[i];
708       v = mat->v + i;
709       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
710     }
711     PLogFlops(n*m);
712   }
713   if (rr) {
714     VecGetArray(rr,&r); VecGetSize(rr,&n);
715     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
716     for ( i=0; i<n; i++ ) {
717       x = r[i];
718       v = mat->v + i*m;
719       for ( j=0; j<m; j++ ) { (*v++) *= x;}
720     }
721     PLogFlops(n*m);
722   }
723   return 0;
724 }
725 
726 static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
727 {
728   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
729   Scalar       *v = mat->v;
730   double       sum = 0.0;
731   int          i, j;
732 
733   if (type == NORM_FROBENIUS) {
734     for (i=0; i<mat->n*mat->m; i++ ) {
735 #if defined(PETSC_COMPLEX)
736       sum += real(conj(*v)*(*v)); v++;
737 #else
738       sum += (*v)*(*v); v++;
739 #endif
740     }
741     *norm = sqrt(sum);
742     PLogFlops(2*mat->n*mat->m);
743   }
744   else if (type == NORM_1) {
745     *norm = 0.0;
746     for ( j=0; j<mat->n; j++ ) {
747       sum = 0.0;
748       for ( i=0; i<mat->m; i++ ) {
749         sum += PetscAbsScalar(*v);  v++;
750       }
751       if (sum > *norm) *norm = sum;
752     }
753     PLogFlops(mat->n*mat->m);
754   }
755   else if (type == NORM_INFINITY) {
756     *norm = 0.0;
757     for ( j=0; j<mat->m; j++ ) {
758       v = mat->v + j;
759       sum = 0.0;
760       for ( i=0; i<mat->n; i++ ) {
761         sum += PetscAbsScalar(*v); v += mat->m;
762       }
763       if (sum > *norm) *norm = sum;
764     }
765     PLogFlops(mat->n*mat->m);
766   }
767   else {
768     SETERRQ(1,"MatNorm_SeqDense:No two norm");
769   }
770   return 0;
771 }
772 
773 static int MatSetOption_SeqDense(Mat A,MatOption op)
774 {
775   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
776 
777   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
778   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
779   else if (op == MAT_ROWS_SORTED ||
780            op == MAT_COLUMNS_SORTED ||
781            op == MAT_SYMMETRIC ||
782            op == MAT_STRUCTURALLY_SYMMETRIC ||
783            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
784            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
785            op == MAT_NO_NEW_DIAGONALS ||
786            op == MAT_YES_NEW_DIAGONALS)
787     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
788   else
789     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
790   return 0;
791 }
792 
793 static int MatZeroEntries_SeqDense(Mat A)
794 {
795   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
796   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
797   return 0;
798 }
799 
800 static int MatGetBlockSize_SeqDense(Mat A,int *bs)
801 {
802   *bs = 1;
803   return 0;
804 }
805 
806 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
807 {
808   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
809   int          n = l->n, i, j,ierr,N, *rows;
810   Scalar       *slot;
811 
812   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
813   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
814   for ( i=0; i<N; i++ ) {
815     slot = l->v + rows[i];
816     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
817   }
818   if (diag) {
819     for ( i=0; i<N; i++ ) {
820       slot = l->v + (n+1)*rows[i];
821       *slot = *diag;
822     }
823   }
824   ISRestoreIndices(is,&rows);
825   return 0;
826 }
827 
828 static int MatGetSize_SeqDense(Mat A,int *m,int *n)
829 {
830   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
831   *m = mat->m; *n = mat->n;
832   return 0;
833 }
834 
835 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
836 {
837   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
838   *m = 0; *n = mat->m;
839   return 0;
840 }
841 
842 static int MatGetArray_SeqDense(Mat A,Scalar **array)
843 {
844   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
845   *array = mat->v;
846   return 0;
847 }
848 
849 static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
850 {
851   return 0;
852 }
853 
854 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
855                                     Mat *submat)
856 {
857   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
858   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
859   int          *irow, *icol, nrows, ncols, *cwork;
860   Scalar       *vwork, *val;
861   Mat          newmat;
862 
863   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
864   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
865   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
866   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
867 
868   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
869   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
870   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
871   PetscMemzero((char*)smap,oldcols*sizeof(int));
872   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
873 
874   /* Create and fill new matrix */
875   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
876   for (i=0; i<nrows; i++) {
877     nznew = 0;
878     val   = mat->v + irow[i];
879     for (j=0; j<oldcols; j++) {
880       if (smap[j]) {
881         cwork[nznew]   = smap[j] - 1;
882         vwork[nznew++] = val[j * mat->m];
883       }
884     }
885     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
886            CHKERRQ(ierr);
887   }
888   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
889   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
890 
891   /* Free work space */
892   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
893   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
894   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
895   *submat = newmat;
896   return 0;
897 }
898 
899 static int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
900                                     Mat **B)
901 {
902   int ierr,i;
903 
904   if (scall == MAT_INITIAL_MATRIX) {
905     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
906   }
907 
908   for ( i=0; i<n; i++ ) {
909     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
910   }
911   return 0;
912 }
913 
914 static int MatCopy_SeqDense(Mat A, Mat B)
915 {
916   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
917   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
918   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
919   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
920   return 0;
921 }
922 
923 /* -------------------------------------------------------------------*/
924 static struct _MatOps MatOps = {MatSetValues_SeqDense,
925        MatGetRow_SeqDense,
926        MatRestoreRow_SeqDense,
927        MatMult_SeqDense,
928        MatMultAdd_SeqDense,
929        MatMultTrans_SeqDense,
930        MatMultTransAdd_SeqDense,
931        MatSolve_SeqDense,
932        MatSolveAdd_SeqDense,
933        MatSolveTrans_SeqDense,
934        MatSolveTransAdd_SeqDense,
935        MatLUFactor_SeqDense,
936        MatCholeskyFactor_SeqDense,
937        MatRelax_SeqDense,
938        MatTranspose_SeqDense,
939        MatGetInfo_SeqDense,
940        MatEqual_SeqDense,
941        MatGetDiagonal_SeqDense,
942        MatDiagonalScale_SeqDense,
943        MatNorm_SeqDense,
944        0,
945        0,
946        0,
947        MatSetOption_SeqDense,
948        MatZeroEntries_SeqDense,
949        MatZeroRows_SeqDense,
950        MatLUFactorSymbolic_SeqDense,
951        MatLUFactorNumeric_SeqDense,
952        MatCholeskyFactorSymbolic_SeqDense,
953        MatCholeskyFactorNumeric_SeqDense,
954        MatGetSize_SeqDense,
955        MatGetSize_SeqDense,
956        MatGetOwnershipRange_SeqDense,
957        0,
958        0,
959        MatGetArray_SeqDense,
960        MatRestoreArray_SeqDense,
961        0,
962        MatConvertSameType_SeqDense,0,0,0,0,
963        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
964        MatGetValues_SeqDense,
965        MatCopy_SeqDense,0,MatScale_SeqDense,
966        0,0,0,MatGetBlockSize_SeqDense};
967 
968 /*@C
969    MatCreateSeqDense - Creates a sequential dense matrix that
970    is stored in column major order (the usual Fortran 77 manner). Many
971    of the matrix operations use the BLAS and LAPACK routines.
972 
973    Input Parameters:
974 .  comm - MPI communicator, set to MPI_COMM_SELF
975 .  m - number of rows
976 .  n - number of columns
977 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
978    to control all matrix memory allocation.
979 
980    Output Parameter:
981 .  A - the matrix
982 
983   Notes:
984   The data input variable is intended primarily for Fortran programmers
985   who wish to allocate their own matrix memory space.  Most users should
986   set data=PETSC_NULL.
987 
988 .keywords: dense, matrix, LAPACK, BLAS
989 
990 .seealso: MatCreate(), MatSetValues()
991 @*/
992 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
993 {
994   Mat          B;
995   Mat_SeqDense *b;
996   int          ierr,flg,size;
997 
998   MPI_Comm_size(comm,&size);
999   if (size > 1) SETERRQ(1,"MatCreateSeqDense:Comm must be of size 1");
1000 
1001   *A            = 0;
1002   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
1003   PLogObjectCreate(B);
1004   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
1005   PetscMemzero(b,sizeof(Mat_SeqDense));
1006   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1007   B->destroy    = MatDestroy_SeqDense;
1008   B->view       = MatView_SeqDense;
1009   B->factor     = 0;
1010   PLogObjectMemory(B,sizeof(struct _Mat));
1011   B->data       = (void *) b;
1012 
1013   b->m = m;  B->m = m; B->M = m;
1014   b->n = n;  B->n = n; B->N = n;
1015   b->pivots       = 0;
1016   b->roworiented  = 1;
1017   if (data == PETSC_NULL) {
1018     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
1019     PetscMemzero(b->v,m*n*sizeof(Scalar));
1020     b->user_alloc = 0;
1021     PLogObjectMemory(B,n*m*sizeof(Scalar));
1022   }
1023   else { /* user-allocated storage */
1024     b->v = data;
1025     b->user_alloc = 1;
1026   }
1027   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
1028   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
1029 
1030   *A = B;
1031   return 0;
1032 }
1033 
1034 int MatCreate_SeqDense(Mat A,Mat *newmat)
1035 {
1036   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
1037   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
1038 }
1039 
1040 
1041 
1042 
1043 
1044