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