xref: /petsc/src/mat/impls/dense/seq/dense.c (revision e1311b9049e89cb3452dcd306fde571f4b440ff2)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dense.c,v 1.140 1998/04/03 21:48:26 bsmith Exp balay $";
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          ierr, m = mat->m, i;
296 #if !defined(USE_PETSC_COMPLEX)
297   int          o = 1;
298 #endif
299 
300   PetscFunctionBegin;
301   if (flag & SOR_ZERO_INITIAL_GUESS) {
302     /* this is a hack fix, should have another version without the second BLdot */
303     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
304   }
305   VecGetArray(xx,&x); VecGetArray(bb,&b);
306   while (its--) {
307     if (flag & SOR_FORWARD_SWEEP){
308       for ( i=0; i<m; i++ ) {
309 #if defined(USE_PETSC_COMPLEX)
310         /* cannot use BLAS dot for complex because compiler/linker is
311            not happy about returning a double complex */
312         int    _i;
313         Scalar sum = b[i];
314         for ( _i=0; _i<m; _i++ ) {
315           sum -= conj(v[i+_i*m])*x[_i];
316         }
317         xt = sum;
318 #else
319         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
320 #endif
321         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
322       }
323     }
324     if (flag & SOR_BACKWARD_SWEEP) {
325       for ( i=m-1; i>=0; i-- ) {
326 #if defined(USE_PETSC_COMPLEX)
327         /* cannot use BLAS dot for complex because compiler/linker is
328            not happy about returning a double complex */
329         int    _i;
330         Scalar sum = b[i];
331         for ( _i=0; _i<m; _i++ ) {
332           sum -= conj(v[i+_i*m])*x[_i];
333         }
334         xt = sum;
335 #else
336         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
337 #endif
338         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
339       }
340     }
341   }
342   PetscFunctionReturn(0);
343 }
344 
345 /* -----------------------------------------------------------------*/
346 #undef __FUNC__
347 #define __FUNC__ "MatMultTrans_SeqDense"
348 int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
349 {
350   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
351   Scalar       *v = mat->v, *x, *y;
352   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
353 
354   PetscFunctionBegin;
355   VecGetArray(xx,&x), VecGetArray(yy,&y);
356   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
357   PLogFlops(2*mat->m*mat->n - mat->n);
358   PetscFunctionReturn(0);
359 }
360 
361 #undef __FUNC__
362 #define __FUNC__ "MatMult_SeqDense"
363 int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
364 {
365   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
366   Scalar       *v = mat->v, *x, *y;
367   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
368 
369   PetscFunctionBegin;
370   VecGetArray(xx,&x); VecGetArray(yy,&y);
371   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
372   PLogFlops(2*mat->m*mat->n - mat->m);
373   PetscFunctionReturn(0);
374 }
375 
376 #undef __FUNC__
377 #define __FUNC__ "MatMultAdd_SeqDense"
378 int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
379 {
380   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
381   Scalar       *v = mat->v, *x, *y, *z;
382   int          _One=1; Scalar _DOne=1.0;
383 
384   PetscFunctionBegin;
385   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
386   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
387   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
388   PLogFlops(2*mat->m*mat->n);
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNC__
393 #define __FUNC__ "MatMultTransAdd_SeqDense"
394 int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
395 {
396   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
397   Scalar       *v = mat->v, *x, *y, *z;
398   int          _One=1;
399   Scalar       _DOne=1.0;
400 
401   PetscFunctionBegin;
402   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
403   if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar));
404   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
405   PLogFlops(2*mat->m*mat->n);
406   PetscFunctionReturn(0);
407 }
408 
409 /* -----------------------------------------------------------------*/
410 #undef __FUNC__
411 #define __FUNC__ "MatGetRow_SeqDense"
412 int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
413 {
414   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
415   Scalar       *v;
416   int          i;
417 
418   PetscFunctionBegin;
419   *ncols = mat->n;
420   if (cols) {
421     *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols);
422     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
423   }
424   if (vals) {
425     *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals);
426     v = mat->v + row;
427     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
428   }
429   PetscFunctionReturn(0);
430 }
431 
432 #undef __FUNC__
433 #define __FUNC__ "MatRestoreRow_SeqDense"
434 int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
435 {
436   if (cols) { PetscFree(*cols); }
437   if (vals) { PetscFree(*vals); }
438   PetscFunctionReturn(0);
439 }
440 /* ----------------------------------------------------------------*/
441 #undef __FUNC__
442 #define __FUNC__ "MatSetValues_SeqDense"
443 int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
444                                     int *indexn,Scalar *v,InsertMode addv)
445 {
446   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
447   int          i,j;
448 
449   PetscFunctionBegin;
450   if (!mat->roworiented) {
451     if (addv == INSERT_VALUES) {
452       for ( j=0; j<n; j++ ) {
453 #if defined(USE_PETSC_BOPT_g)
454         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
455         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
456 #endif
457         for ( i=0; i<m; i++ ) {
458 #if defined(USE_PETSC_BOPT_g)
459           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
460           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
461 #endif
462           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
463         }
464       }
465     } else {
466       for ( j=0; j<n; j++ ) {
467 #if defined(USE_PETSC_BOPT_g)
468         if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
469         if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
470 #endif
471         for ( i=0; i<m; i++ ) {
472 #if defined(USE_PETSC_BOPT_g)
473           if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
474           if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
475 #endif
476           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
477         }
478       }
479     }
480   } else {
481     if (addv == INSERT_VALUES) {
482       for ( i=0; i<m; i++ ) {
483 #if defined(USE_PETSC_BOPT_g)
484         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
485         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
486 #endif
487         for ( j=0; j<n; j++ ) {
488 #if defined(USE_PETSC_BOPT_g)
489           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
490           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
491 #endif
492           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
493         }
494       }
495     } else {
496       for ( i=0; i<m; i++ ) {
497 #if defined(USE_PETSC_BOPT_g)
498         if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
499         if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
500 #endif
501         for ( j=0; j<n; j++ ) {
502 #if defined(USE_PETSC_BOPT_g)
503           if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
504           if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
505 #endif
506           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
507         }
508       }
509     }
510   }
511   PetscFunctionReturn(0);
512 }
513 
514 #undef __FUNC__
515 #define __FUNC__ "MatGetValues_SeqDense"
516 int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
517 {
518   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
519   int          i, j;
520   Scalar       *vpt = v;
521 
522   PetscFunctionBegin;
523   /* row-oriented output */
524   for ( i=0; i<m; i++ ) {
525     for ( j=0; j<n; j++ ) {
526       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
527     }
528   }
529   PetscFunctionReturn(0);
530 }
531 
532 /* -----------------------------------------------------------------*/
533 
534 #include "sys.h"
535 
536 #undef __FUNC__
537 #define __FUNC__ "MatLoad_SeqDense"
538 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
539 {
540   Mat_SeqDense *a;
541   Mat          B;
542   int          *scols, i, j, nz, ierr, fd, header[4], size;
543   int          *rowlengths = 0, M, N, *cols;
544   Scalar       *vals, *svals, *v;
545   MPI_Comm     comm = ((PetscObject)viewer)->comm;
546 
547   PetscFunctionBegin;
548   MPI_Comm_size(comm,&size);
549   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
550   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
551   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
552   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object");
553   M = header[1]; N = header[2]; nz = header[3];
554 
555   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
556     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
557     B    = *A;
558     a    = (Mat_SeqDense *) B->data;
559 
560     /* read in nonzero values */
561     ierr = PetscBinaryRead(fd,a->v,M*N,PETSC_SCALAR); CHKERRQ(ierr);
562 
563     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
564     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
565     /* ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); */
566   } else {
567     /* read row lengths */
568     rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths);
569     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
570 
571     /* create our matrix */
572     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
573     B = *A;
574     a = (Mat_SeqDense *) B->data;
575     v = a->v;
576 
577     /* read column indices and nonzeros */
578     cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols);
579     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
580     vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals);
581     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
582 
583     /* insert into matrix */
584     for ( i=0; i<M; i++ ) {
585       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
586       svals += rowlengths[i]; scols += rowlengths[i];
587     }
588     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
589 
590     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
591     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
592   }
593   PetscFunctionReturn(0);
594 }
595 
596 #include "pinclude/pviewer.h"
597 #include "sys.h"
598 
599 #undef __FUNC__
600 #define __FUNC__ "MatView_SeqDense_ASCII"
601 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
602 {
603   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
604   int          ierr, i, j, format;
605   FILE         *fd;
606   char         *outputname;
607   Scalar       *v;
608 
609   PetscFunctionBegin;
610   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
611   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
612   ierr = ViewerGetFormat(viewer,&format);
613   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
614     PetscFunctionReturn(0);  /* do nothing for now */
615   }
616   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
617     for ( i=0; i<a->m; i++ ) {
618       v = a->v + i;
619       fprintf(fd,"row %d:",i);
620       for ( j=0; j<a->n; j++ ) {
621 #if defined(USE_PETSC_COMPLEX)
622         if (real(*v) != 0.0 && imag(*v) != 0.0) fprintf(fd," %d %g + %g i",j,real(*v),imag(*v));
623         else if (real(*v)) fprintf(fd," %d %g ",j,real(*v));
624 #else
625         if (*v) fprintf(fd," %d %g ",j,*v);
626 #endif
627         v += a->m;
628       }
629       fprintf(fd,"\n");
630     }
631   } else {
632 #if defined(USE_PETSC_COMPLEX)
633     int allreal = 1;
634     /* determine if matrix has all real values */
635     v = a->v;
636     for ( i=0; i<a->m*a->n; i++ ) {
637       if (imag(v[i])) { allreal = 0; break ;}
638     }
639 #endif
640     for ( i=0; i<a->m; i++ ) {
641       v = a->v + i;
642       for ( j=0; j<a->n; j++ ) {
643 #if defined(USE_PETSC_COMPLEX)
644         if (allreal) {
645           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
646         } else {
647           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
648         }
649 #else
650         fprintf(fd,"%6.4e ",*v); v += a->m;
651 #endif
652       }
653       fprintf(fd,"\n");
654     }
655   }
656   fflush(fd);
657   PetscFunctionReturn(0);
658 }
659 
660 #undef __FUNC__
661 #define __FUNC__ "MatView_SeqDense_Binary"
662 static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
663 {
664   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
665   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
666   int          format;
667   Scalar       *v, *anonz,*vals;
668 
669   PetscFunctionBegin;
670   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
671 
672   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
673   if (format == VIEWER_FORMAT_BINARY_NATIVE) {
674     /* store the matrix as a dense matrix */
675     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
676     col_lens[0] = MAT_COOKIE;
677     col_lens[1] = m;
678     col_lens[2] = n;
679     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
680     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr);
681     PetscFree(col_lens);
682 
683     /* write out matrix, by rows */
684     vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals);
685     v    = a->v;
686     for ( i=0; i<m; i++ ) {
687       for ( j=0; j<n; j++ ) {
688         vals[i + j*m] = *v++;
689       }
690     }
691     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr);
692     PetscFree(vals);
693   } else {
694     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
695     col_lens[0] = MAT_COOKIE;
696     col_lens[1] = m;
697     col_lens[2] = n;
698     col_lens[3] = nz;
699 
700     /* store lengths of each row and write (including header) to file */
701     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
702     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr);
703 
704     /* Possibly should write in smaller increments, not whole matrix at once? */
705     /* store column indices (zero start index) */
706     ict = 0;
707     for ( i=0; i<m; i++ ) {
708       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
709     }
710     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr);
711     PetscFree(col_lens);
712 
713     /* store nonzero values */
714     anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz);
715     ict = 0;
716     for ( i=0; i<m; i++ ) {
717       v = a->v + i;
718       for ( j=0; j<n; j++ ) {
719         anonz[ict++] = *v; v += a->m;
720       }
721     }
722     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr);
723     PetscFree(anonz);
724   }
725   PetscFunctionReturn(0);
726 }
727 
728 #undef __FUNC__
729 #define __FUNC__ "MatView_SeqDense"
730 int MatView_SeqDense(Mat A,Viewer viewer)
731 {
732   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
733   ViewerType   vtype;
734   int          ierr;
735 
736   PetscFunctionBegin;
737   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
738 
739   if (vtype == MATLAB_VIEWER) {
740     ierr = ViewerMatlabPutScalar_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr);
741   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
742     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
743   } else if (vtype == BINARY_FILE_VIEWER) {
744     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
745   }
746   PetscFunctionReturn(0);
747 }
748 
749 #undef __FUNC__
750 #define __FUNC__ "MatDestroy_SeqDense"
751 int MatDestroy_SeqDense(Mat mat)
752 {
753   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
754   int          ierr;
755 
756   PetscFunctionBegin;
757 #if defined(USE_PETSC_LOG)
758   PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n);
759 #endif
760   if (l->pivots) PetscFree(l->pivots);
761   if (!l->user_alloc) PetscFree(l->v);
762   PetscFree(l);
763   if (mat->mapping) {
764     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
765   }
766   PLogObjectDestroy(mat);
767   PetscHeaderDestroy(mat);
768   PetscFunctionReturn(0);
769 }
770 
771 #undef __FUNC__
772 #define __FUNC__ "MatTranspose_SeqDense"
773 int MatTranspose_SeqDense(Mat A,Mat *matout)
774 {
775   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
776   int          k, j, m, n;
777   Scalar       *v, tmp;
778 
779   PetscFunctionBegin;
780   v = mat->v; m = mat->m; n = mat->n;
781   if (matout == PETSC_NULL) { /* in place transpose */
782     if (m != n) { /* malloc temp to hold transpose */
783       Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w);
784       for ( j=0; j<m; j++ ) {
785         for ( k=0; k<n; k++ ) {
786           w[k + j*n] = v[j + k*m];
787         }
788       }
789       PetscMemcpy(v,w,m*n*sizeof(Scalar));
790       PetscFree(w);
791     } else {
792       for ( j=0; j<m; j++ ) {
793         for ( k=0; k<j; k++ ) {
794           tmp = v[j + k*n];
795           v[j + k*n] = v[k + j*n];
796           v[k + j*n] = tmp;
797         }
798       }
799     }
800   } else { /* out-of-place transpose */
801     int          ierr;
802     Mat          tmat;
803     Mat_SeqDense *tmatd;
804     Scalar       *v2;
805     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
806     tmatd = (Mat_SeqDense *) tmat->data;
807     v = mat->v; v2 = tmatd->v;
808     for ( j=0; j<n; j++ ) {
809       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
810     }
811     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
812     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
813     *matout = tmat;
814   }
815   PetscFunctionReturn(0);
816 }
817 
818 #undef __FUNC__
819 #define __FUNC__ "MatEqual_SeqDense"
820 int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
821 {
822   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
823   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
824   int          i;
825   Scalar       *v1 = mat1->v, *v2 = mat2->v;
826 
827   PetscFunctionBegin;
828   if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type  MATSEQDENSE");
829   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
830   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
831   for ( i=0; i<mat1->m*mat1->n; i++ ) {
832     if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
833     v1++; v2++;
834   }
835   *flg = PETSC_TRUE;
836   PetscFunctionReturn(0);
837 }
838 
839 #undef __FUNC__
840 #define __FUNC__ "MatGetDiagonal_SeqDense"
841 int MatGetDiagonal_SeqDense(Mat A,Vec v)
842 {
843   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
844   int          i, n, len;
845   Scalar       *x, zero = 0.0;
846 
847   PetscFunctionBegin;
848   VecSet(&zero,v);
849   VecGetArray(v,&x); VecGetSize(v,&n);
850   len = PetscMin(mat->m,mat->n);
851   if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
852   for ( i=0; i<len; i++ ) {
853     x[i] = mat->v[i*mat->m + i];
854   }
855   PetscFunctionReturn(0);
856 }
857 
858 #undef __FUNC__
859 #define __FUNC__ "MatDiagonalScale_SeqDense"
860 int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
861 {
862   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
863   Scalar       *l,*r,x,*v;
864   int          i,j,m = mat->m, n = mat->n;
865 
866   PetscFunctionBegin;
867   if (ll) {
868     VecGetArray(ll,&l); VecGetSize(ll,&m);
869     if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size");
870     for ( i=0; i<m; i++ ) {
871       x = l[i];
872       v = mat->v + i;
873       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
874     }
875     PLogFlops(n*m);
876   }
877   if (rr) {
878     VecGetArray(rr,&r); VecGetSize(rr,&n);
879     if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size");
880     for ( i=0; i<n; i++ ) {
881       x = r[i];
882       v = mat->v + i*m;
883       for ( j=0; j<m; j++ ) { (*v++) *= x;}
884     }
885     PLogFlops(n*m);
886   }
887   PetscFunctionReturn(0);
888 }
889 
890 #undef __FUNC__
891 #define __FUNC__ "MatNorm_SeqDense"
892 int MatNorm_SeqDense(Mat A,NormType type,double *norm)
893 {
894   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
895   Scalar       *v = mat->v;
896   double       sum = 0.0;
897   int          i, j;
898 
899   PetscFunctionBegin;
900   if (type == NORM_FROBENIUS) {
901     for (i=0; i<mat->n*mat->m; i++ ) {
902 #if defined(USE_PETSC_COMPLEX)
903       sum += real(conj(*v)*(*v)); v++;
904 #else
905       sum += (*v)*(*v); v++;
906 #endif
907     }
908     *norm = sqrt(sum);
909     PLogFlops(2*mat->n*mat->m);
910   } else if (type == NORM_1) {
911     *norm = 0.0;
912     for ( j=0; j<mat->n; j++ ) {
913       sum = 0.0;
914       for ( i=0; i<mat->m; i++ ) {
915         sum += PetscAbsScalar(*v);  v++;
916       }
917       if (sum > *norm) *norm = sum;
918     }
919     PLogFlops(mat->n*mat->m);
920   } else if (type == NORM_INFINITY) {
921     *norm = 0.0;
922     for ( j=0; j<mat->m; j++ ) {
923       v = mat->v + j;
924       sum = 0.0;
925       for ( i=0; i<mat->n; i++ ) {
926         sum += PetscAbsScalar(*v); v += mat->m;
927       }
928       if (sum > *norm) *norm = sum;
929     }
930     PLogFlops(mat->n*mat->m);
931   } else {
932     SETERRQ(PETSC_ERR_SUP,0,"No two norm");
933   }
934   PetscFunctionReturn(0);
935 }
936 
937 #undef __FUNC__
938 #define __FUNC__ "MatSetOption_SeqDense"
939 int MatSetOption_SeqDense(Mat A,MatOption op)
940 {
941   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
942 
943   PetscFunctionBegin;
944   if (op == MAT_ROW_ORIENTED)            aij->roworiented = 1;
945   else if (op == MAT_COLUMN_ORIENTED)    aij->roworiented = 0;
946   else if (op == MAT_ROWS_SORTED ||
947            op == MAT_ROWS_UNSORTED ||
948            op == MAT_COLUMNS_SORTED ||
949            op == MAT_COLUMNS_UNSORTED ||
950            op == MAT_SYMMETRIC ||
951            op == MAT_STRUCTURALLY_SYMMETRIC ||
952            op == MAT_NO_NEW_NONZERO_LOCATIONS ||
953            op == MAT_YES_NEW_NONZERO_LOCATIONS ||
954            op == MAT_NEW_NONZERO_LOCATION_ERROR ||
955            op == MAT_NO_NEW_DIAGONALS ||
956            op == MAT_YES_NEW_DIAGONALS ||
957            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
958            op == MAT_USE_HASH_TABLE)
959     PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
960   else {
961     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
962   }
963   PetscFunctionReturn(0);
964 }
965 
966 #undef __FUNC__
967 #define __FUNC__ "MatZeroEntries_SeqDense"
968 int MatZeroEntries_SeqDense(Mat A)
969 {
970   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
971 
972   PetscFunctionBegin;
973   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
974   PetscFunctionReturn(0);
975 }
976 
977 #undef __FUNC__
978 #define __FUNC__ "MatGetBlockSize_SeqDense"
979 int MatGetBlockSize_SeqDense(Mat A,int *bs)
980 {
981   PetscFunctionBegin;
982   *bs = 1;
983   PetscFunctionReturn(0);
984 }
985 
986 #undef __FUNC__
987 #define __FUNC__ "MatZeroRows_SeqDense"
988 int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
989 {
990   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
991   int          n = l->n, i, j,ierr,N, *rows;
992   Scalar       *slot;
993 
994   PetscFunctionBegin;
995   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
996   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
997   for ( i=0; i<N; i++ ) {
998     slot = l->v + rows[i];
999     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
1000   }
1001   if (diag) {
1002     for ( i=0; i<N; i++ ) {
1003       slot = l->v + (n+1)*rows[i];
1004       *slot = *diag;
1005     }
1006   }
1007   ISRestoreIndices(is,&rows);
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNC__
1012 #define __FUNC__ "MatGetSize_SeqDense"
1013 int MatGetSize_SeqDense(Mat A,int *m,int *n)
1014 {
1015   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1016 
1017   PetscFunctionBegin;
1018   *m = mat->m; *n = mat->n;
1019   PetscFunctionReturn(0);
1020 }
1021 
1022 #undef __FUNC__
1023 #define __FUNC__ "MatGetOwnershipRange_SeqDense"
1024 int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
1025 {
1026   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1027 
1028   PetscFunctionBegin;
1029   *m = 0; *n = mat->m;
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 #undef __FUNC__
1034 #define __FUNC__ "MatGetArray_SeqDense"
1035 int MatGetArray_SeqDense(Mat A,Scalar **array)
1036 {
1037   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1038 
1039   PetscFunctionBegin;
1040   *array = mat->v;
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 #undef __FUNC__
1045 #define __FUNC__ "MatRestoreArray_SeqDense"
1046 int MatRestoreArray_SeqDense(Mat A,Scalar **array)
1047 {
1048   PetscFunctionBegin;
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 #undef __FUNC__
1053 #define __FUNC__ "MatGetSubMatrix_SeqDense"
1054 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatGetSubMatrixCall scall,Mat *submat)
1055 {
1056   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1057   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
1058   int          *irow, *icol, nrows, ncols, *cwork;
1059   Scalar       *vwork, *val;
1060   Mat          newmat;
1061 
1062   PetscFunctionBegin;
1063   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
1064   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
1065   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
1066   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
1067 
1068   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
1069   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
1070   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
1071   PetscMemzero((char*)smap,oldcols*sizeof(int));
1072   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
1073 
1074   /* Create and fill new matrix */
1075   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
1076   for (i=0; i<nrows; i++) {
1077     nznew = 0;
1078     val   = mat->v + irow[i];
1079     for (j=0; j<oldcols; j++) {
1080       if (smap[j]) {
1081         cwork[nznew]   = smap[j] - 1;
1082         vwork[nznew++] = val[j * mat->m];
1083       }
1084     }
1085     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
1086   }
1087   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1088   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1089 
1090   /* Free work space */
1091   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
1092   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1093   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1094   *submat = newmat;
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 #undef __FUNC__
1099 #define __FUNC__ "MatGetSubMatrices_SeqDense"
1100 int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1101 {
1102   int ierr,i;
1103 
1104   PetscFunctionBegin;
1105   if (scall == MAT_INITIAL_MATRIX) {
1106     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1107   }
1108 
1109   for ( i=0; i<n; i++ ) {
1110     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1111   }
1112   PetscFunctionReturn(0);
1113 }
1114 
1115 #undef __FUNC__
1116 #define __FUNC__ "MatCopy_SeqDense"
1117 int MatCopy_SeqDense(Mat A, Mat B)
1118 {
1119   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
1120   int          ierr;
1121 
1122   PetscFunctionBegin;
1123   if (B->type != MATSEQDENSE) {
1124     ierr = MatCopy_Basic(A,B);CHKERRQ(ierr);
1125     PetscFunctionReturn(0);
1126   }
1127   if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)");
1128   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
1129   PetscFunctionReturn(0);
1130 }
1131 
1132 /* -------------------------------------------------------------------*/
1133 static struct _MatOps MatOps = {MatSetValues_SeqDense,
1134        MatGetRow_SeqDense,
1135        MatRestoreRow_SeqDense,
1136        MatMult_SeqDense,
1137        MatMultAdd_SeqDense,
1138        MatMultTrans_SeqDense,
1139        MatMultTransAdd_SeqDense,
1140        MatSolve_SeqDense,
1141        MatSolveAdd_SeqDense,
1142        MatSolveTrans_SeqDense,
1143        MatSolveTransAdd_SeqDense,
1144        MatLUFactor_SeqDense,
1145        MatCholeskyFactor_SeqDense,
1146        MatRelax_SeqDense,
1147        MatTranspose_SeqDense,
1148        MatGetInfo_SeqDense,
1149        MatEqual_SeqDense,
1150        MatGetDiagonal_SeqDense,
1151        MatDiagonalScale_SeqDense,
1152        MatNorm_SeqDense,
1153        0,
1154        0,
1155        0,
1156        MatSetOption_SeqDense,
1157        MatZeroEntries_SeqDense,
1158        MatZeroRows_SeqDense,
1159        MatLUFactorSymbolic_SeqDense,
1160        MatLUFactorNumeric_SeqDense,
1161        MatCholeskyFactorSymbolic_SeqDense,
1162        MatCholeskyFactorNumeric_SeqDense,
1163        MatGetSize_SeqDense,
1164        MatGetSize_SeqDense,
1165        MatGetOwnershipRange_SeqDense,
1166        0,
1167        0,
1168        MatGetArray_SeqDense,
1169        MatRestoreArray_SeqDense,
1170        MatConvertSameType_SeqDense,0,0,0,0,
1171        MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0,
1172        MatGetValues_SeqDense,
1173        MatCopy_SeqDense,0,MatScale_SeqDense,
1174        0,0,0,MatGetBlockSize_SeqDense};
1175 
1176 #undef __FUNC__
1177 #define __FUNC__ "MatCreateSeqDense"
1178 /*@C
1179    MatCreateSeqDense - Creates a sequential dense matrix that
1180    is stored in column major order (the usual Fortran 77 manner). Many
1181    of the matrix operations use the BLAS and LAPACK routines.
1182 
1183    Input Parameters:
1184 .  comm - MPI communicator, set to PETSC_COMM_SELF
1185 .  m - number of rows
1186 .  n - number of columns
1187 .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1188    to control all matrix memory allocation.
1189 
1190    Output Parameter:
1191 .  A - the matrix
1192 
1193   Notes:
1194   The data input variable is intended primarily for Fortran programmers
1195   who wish to allocate their own matrix memory space.  Most users should
1196   set data=PETSC_NULL.
1197 
1198 .keywords: dense, matrix, LAPACK, BLAS
1199 
1200 .seealso: MatCreate(), MatSetValues()
1201 @*/
1202 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A)
1203 {
1204   Mat          B;
1205   Mat_SeqDense *b;
1206   int          ierr,flg,size;
1207 
1208   PetscFunctionBegin;
1209   MPI_Comm_size(comm,&size);
1210   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1211 
1212   *A            = 0;
1213   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView);
1214   PLogObjectCreate(B);
1215   b             = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b);
1216   PetscMemzero(b,sizeof(Mat_SeqDense));
1217   PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps));
1218   B->ops->destroy    = MatDestroy_SeqDense;
1219   B->ops->view       = MatView_SeqDense;
1220   B->factor          = 0;
1221   B->mapping         = 0;
1222   PLogObjectMemory(B,sizeof(struct _p_Mat));
1223   B->data       = (void *) b;
1224 
1225   b->m = m;  B->m = m; B->M = m;
1226   b->n = n;  B->n = n; B->N = n;
1227   b->pivots       = 0;
1228   b->roworiented  = 1;
1229   if (data == PETSC_NULL) {
1230     b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v);
1231     PetscMemzero(b->v,m*n*sizeof(Scalar));
1232     b->user_alloc = 0;
1233     PLogObjectMemory(B,n*m*sizeof(Scalar));
1234   }
1235   else { /* user-allocated storage */
1236     b->v = data;
1237     b->user_alloc = 1;
1238   }
1239   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
1240   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);}
1241 
1242   *A = B;
1243   PetscFunctionReturn(0);
1244 }
1245 
1246 
1247 
1248 
1249 
1250 
1251