Actual source code: dense.c
1: /*$Id: dense.c,v 1.208 2001/09/07 20:09:16 bsmith Exp $*/
2: /*
3: Defines the basic matrix operations for sequential dense.
4: */
6: #include src/mat/impls/dense/seq/dense.h
7: #include petscblaslapack.h
9: int MatAXPY_SeqDense(PetscScalar *alpha,Mat X,Mat Y,MatStructure str)
10: {
11: Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
12: int N = X->m*X->n,one = 1;
15: BLaxpy_(&N,alpha,x->v,&one,y->v,&one);
16: PetscLogFlops(2*N-1);
17: return(0);
18: }
20: int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
21: {
22: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
23: int i,N = A->m*A->n,count = 0;
24: PetscScalar *v = a->v;
27: for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}
29: info->rows_global = (double)A->m;
30: info->columns_global = (double)A->n;
31: info->rows_local = (double)A->m;
32: info->columns_local = (double)A->n;
33: info->block_size = 1.0;
34: info->nz_allocated = (double)N;
35: info->nz_used = (double)count;
36: info->nz_unneeded = (double)(N-count);
37: info->assemblies = (double)A->num_ass;
38: info->mallocs = 0;
39: info->memory = A->mem;
40: info->fill_ratio_given = 0;
41: info->fill_ratio_needed = 0;
42: info->factor_mallocs = 0;
44: return(0);
45: }
47: int MatScale_SeqDense(PetscScalar *alpha,Mat A)
48: {
49: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
50: int one = 1,nz;
53: nz = A->m*A->n;
54: BLscal_(&nz,alpha,a->v,&one);
55: PetscLogFlops(nz);
56: return(0);
57: }
58:
59: /* ---------------------------------------------------------------*/
60: /* COMMENT: I have chosen to hide row permutation in the pivots,
61: rather than put it in the Mat->row slot.*/
62: int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatLUInfo *minfo)
63: {
64: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
65: int info,ierr;
68: if (!mat->pivots) {
69: PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);
70: PetscLogObjectMemory(A,A->m*sizeof(int));
71: }
72: A->factor = FACTOR_LU;
73: if (!A->m || !A->n) return(0);
74: #if defined(PETSC_MISSING_LAPACK_GETRF)
75: SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable.");
76: #else
77: LAgetrf_(&A->m,&A->n,mat->v,&A->m,mat->pivots,&info);
78: if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
79: if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
80: #endif
81: PetscLogFlops((2*A->n*A->n*A->n)/3);
82: return(0);
83: }
85: int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
86: {
87: Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l;
88: int ierr;
89: Mat newi;
92: MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);
93: l = (Mat_SeqDense*)newi->data;
94: if (cpvalues == MAT_COPY_VALUES) {
95: PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));
96: }
97: newi->assembled = PETSC_TRUE;
98: *newmat = newi;
99: return(0);
100: }
102: int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatLUInfo *info,Mat *fact)
103: {
107: MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);
108: return(0);
109: }
111: int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
112: {
113: Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
114: int ierr;
117: /* copy the numerical values */
118: PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));
119: (*fact)->factor = 0;
120: MatLUFactor(*fact,0,0,PETSC_NULL);
121: return(0);
122: }
124: int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,PetscReal f,Mat *fact)
125: {
129: MatConvert(A,MATSAME,fact);
130: return(0);
131: }
133: int MatCholeskyFactor_SeqDense(Mat A,IS perm,PetscReal f)
134: {
135: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
136: int info,ierr;
137:
139: if (mat->pivots) {
140: PetscFree(mat->pivots);
141: PetscLogObjectMemory(A,-A->m*sizeof(int));
142: mat->pivots = 0;
143: }
144: if (!A->m || !A->n) return(0);
145: #if defined(PETSC_MISSING_LAPACK_POTRF)
146: SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable.");
147: #else
148: LApotrf_("L",&A->n,mat->v,&A->m,&info);
149: if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
150: #endif
151: A->factor = FACTOR_CHOLESKY;
152: PetscLogFlops((A->n*A->n*A->n)/3);
153: return(0);
154: }
156: int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
157: {
161: MatCholeskyFactor_SeqDense(*fact,0,1.0);
162: return(0);
163: }
165: int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
166: {
167: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
168: int one = 1,info,ierr;
169: PetscScalar *x,*y;
170:
172: if (!A->m || !A->n) return(0);
173: VecGetArray(xx,&x);
174: VecGetArray(yy,&y);
175: PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
176: if (A->factor == FACTOR_LU) {
177: #if defined(PETSC_MISSING_LAPACK_GETRS)
178: SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
179: #else
180: LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
181: if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
182: #endif
183: } else if (A->factor == FACTOR_CHOLESKY){
184: #if defined(PETSC_MISSING_LAPACK_POTRS)
185: SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
186: #else
187: LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
188: if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
189: #endif
190: }
191: else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
192: VecRestoreArray(xx,&x);
193: VecRestoreArray(yy,&y);
194: PetscLogFlops(2*A->n*A->n - A->n);
195: return(0);
196: }
198: int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
199: {
200: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
201: int ierr,one = 1,info;
202: PetscScalar *x,*y;
203:
205: if (!A->m || !A->n) return(0);
206: VecGetArray(xx,&x);
207: VecGetArray(yy,&y);
208: PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
209: /* assume if pivots exist then use LU; else Cholesky */
210: if (mat->pivots) {
211: #if defined(PETSC_MISSING_LAPACK_GETRS)
212: SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
213: #else
214: LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
215: if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
216: #endif
217: } else {
218: #if defined(PETSC_MISSING_LAPACK_POTRS)
219: SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
220: #else
221: LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
222: if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
223: #endif
224: }
225: VecRestoreArray(xx,&x);
226: VecRestoreArray(yy,&y);
227: PetscLogFlops(2*A->n*A->n - A->n);
228: return(0);
229: }
231: int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
232: {
233: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
234: int ierr,one = 1,info;
235: PetscScalar *x,*y,sone = 1.0;
236: Vec tmp = 0;
237:
239: VecGetArray(xx,&x);
240: VecGetArray(yy,&y);
241: if (!A->m || !A->n) return(0);
242: if (yy == zz) {
243: VecDuplicate(yy,&tmp);
244: PetscLogObjectParent(A,tmp);
245: VecCopy(yy,tmp);
246: }
247: PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
248: /* assume if pivots exist then use LU; else Cholesky */
249: if (mat->pivots) {
250: #if defined(PETSC_MISSING_LAPACK_GETRS)
251: SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
252: #else
253: LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
254: if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
255: #endif
256: } else {
257: #if defined(PETSC_MISSING_LAPACK_POTRS)
258: SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
259: #else
260: LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
261: if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
262: #endif
263: }
264: if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
265: else {VecAXPY(&sone,zz,yy);}
266: VecRestoreArray(xx,&x);
267: VecRestoreArray(yy,&y);
268: PetscLogFlops(2*A->n*A->n);
269: return(0);
270: }
272: int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
273: {
274: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
275: int one = 1,info,ierr;
276: PetscScalar *x,*y,sone = 1.0;
277: Vec tmp;
278:
280: if (!A->m || !A->n) return(0);
281: VecGetArray(xx,&x);
282: VecGetArray(yy,&y);
283: if (yy == zz) {
284: VecDuplicate(yy,&tmp);
285: PetscLogObjectParent(A,tmp);
286: VecCopy(yy,tmp);
287: }
288: PetscMemcpy(y,x,A->m*sizeof(PetscScalar));
289: /* assume if pivots exist then use LU; else Cholesky */
290: if (mat->pivots) {
291: #if defined(PETSC_MISSING_LAPACK_GETRS)
292: SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable.");
293: #else
294: LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info);
295: if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
296: #endif
297: } else {
298: #if defined(PETSC_MISSING_LAPACK_POTRS)
299: SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable.");
300: #else
301: LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info);
302: if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
303: #endif
304: }
305: if (tmp) {
306: VecAXPY(&sone,tmp,yy);
307: VecDestroy(tmp);
308: } else {
309: VecAXPY(&sone,zz,yy);
310: }
311: VecRestoreArray(xx,&x);
312: VecRestoreArray(yy,&y);
313: PetscLogFlops(2*A->n*A->n);
314: return(0);
315: }
316: /* ------------------------------------------------------------------*/
317: int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
318: PetscReal shift,int its,int lits,Vec xx)
319: {
320: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
321: PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt;
322: int ierr,m = A->m,i;
323: #if !defined(PETSC_USE_COMPLEX)
324: int o = 1;
325: #endif
328: if (flag & SOR_ZERO_INITIAL_GUESS) {
329: /* this is a hack fix, should have another version without the second BLdot */
330: VecSet(&zero,xx);
331: }
332: VecGetArray(xx,&x);
333: VecGetArray(bb,&b);
334: its = its*lits;
335: if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
336: while (its--) {
337: if (flag & SOR_FORWARD_SWEEP){
338: for (i=0; i<m; i++) {
339: #if defined(PETSC_USE_COMPLEX)
340: /* cannot use BLAS dot for complex because compiler/linker is
341: not happy about returning a double complex */
342: int _i;
343: PetscScalar sum = b[i];
344: for (_i=0; _i<m; _i++) {
345: sum -= PetscConj(v[i+_i*m])*x[_i];
346: }
347: xt = sum;
348: #else
349: xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
350: #endif
351: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
352: }
353: }
354: if (flag & SOR_BACKWARD_SWEEP) {
355: for (i=m-1; i>=0; i--) {
356: #if defined(PETSC_USE_COMPLEX)
357: /* cannot use BLAS dot for complex because compiler/linker is
358: not happy about returning a double complex */
359: int _i;
360: PetscScalar sum = b[i];
361: for (_i=0; _i<m; _i++) {
362: sum -= PetscConj(v[i+_i*m])*x[_i];
363: }
364: xt = sum;
365: #else
366: xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
367: #endif
368: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
369: }
370: }
371: }
372: VecRestoreArray(bb,&b);
373: VecRestoreArray(xx,&x);
374: return(0);
375: }
377: /* -----------------------------------------------------------------*/
378: int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
379: {
380: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
381: PetscScalar *v = mat->v,*x,*y;
382: int ierr,_One=1;
383: PetscScalar _DOne=1.0,_DZero=0.0;
386: if (!A->m || !A->n) return(0);
387: VecGetArray(xx,&x);
388: VecGetArray(yy,&y);
389: LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
390: VecRestoreArray(xx,&x);
391: VecRestoreArray(yy,&y);
392: PetscLogFlops(2*A->m*A->n - A->n);
393: return(0);
394: }
396: int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
397: {
398: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
399: PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
400: int ierr,_One=1;
403: if (!A->m || !A->n) return(0);
404: VecGetArray(xx,&x);
405: VecGetArray(yy,&y);
406: LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One);
407: VecRestoreArray(xx,&x);
408: VecRestoreArray(yy,&y);
409: PetscLogFlops(2*A->m*A->n - A->m);
410: return(0);
411: }
413: int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
414: {
415: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
416: PetscScalar *v = mat->v,*x,*y,_DOne=1.0;
417: int ierr,_One=1;
420: if (!A->m || !A->n) return(0);
421: if (zz != yy) {VecCopy(zz,yy);}
422: VecGetArray(xx,&x);
423: VecGetArray(yy,&y);
424: LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
425: VecRestoreArray(xx,&x);
426: VecRestoreArray(yy,&y);
427: PetscLogFlops(2*A->m*A->n);
428: return(0);
429: }
431: int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
432: {
433: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
434: PetscScalar *v = mat->v,*x,*y;
435: int ierr,_One=1;
436: PetscScalar _DOne=1.0;
439: if (!A->m || !A->n) return(0);
440: if (zz != yy) {VecCopy(zz,yy);}
441: VecGetArray(xx,&x);
442: VecGetArray(yy,&y);
443: LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One);
444: VecRestoreArray(xx,&x);
445: VecRestoreArray(yy,&y);
446: PetscLogFlops(2*A->m*A->n);
447: return(0);
448: }
450: /* -----------------------------------------------------------------*/
451: int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
452: {
453: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
454: PetscScalar *v;
455: int i,ierr;
456:
458: *ncols = A->n;
459: if (cols) {
460: PetscMalloc((A->n+1)*sizeof(int),cols);
461: for (i=0; i<A->n; i++) (*cols)[i] = i;
462: }
463: if (vals) {
464: PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);
465: v = mat->v + row;
466: for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->m;}
467: }
468: return(0);
469: }
471: int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
472: {
475: if (cols) {PetscFree(*cols);}
476: if (vals) {PetscFree(*vals); }
477: return(0);
478: }
479: /* ----------------------------------------------------------------*/
480: int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
481: int *indexn,PetscScalar *v,InsertMode addv)
482: {
483: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
484: int i,j;
485:
487: if (!mat->roworiented) {
488: if (addv == INSERT_VALUES) {
489: for (j=0; j<n; j++) {
490: if (indexn[j] < 0) {v += m; continue;}
491: #if defined(PETSC_USE_BOPT_g)
492: if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
493: #endif
494: for (i=0; i<m; i++) {
495: if (indexm[i] < 0) {v++; continue;}
496: #if defined(PETSC_USE_BOPT_g)
497: if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
498: #endif
499: mat->v[indexn[j]*A->m + indexm[i]] = *v++;
500: }
501: }
502: } else {
503: for (j=0; j<n; j++) {
504: if (indexn[j] < 0) {v += m; continue;}
505: #if defined(PETSC_USE_BOPT_g)
506: if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
507: #endif
508: for (i=0; i<m; i++) {
509: if (indexm[i] < 0) {v++; continue;}
510: #if defined(PETSC_USE_BOPT_g)
511: if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
512: #endif
513: mat->v[indexn[j]*A->m + indexm[i]] += *v++;
514: }
515: }
516: }
517: } else {
518: if (addv == INSERT_VALUES) {
519: for (i=0; i<m; i++) {
520: if (indexm[i] < 0) { v += n; continue;}
521: #if defined(PETSC_USE_BOPT_g)
522: if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
523: #endif
524: for (j=0; j<n; j++) {
525: if (indexn[j] < 0) { v++; continue;}
526: #if defined(PETSC_USE_BOPT_g)
527: if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
528: #endif
529: mat->v[indexn[j]*A->m + indexm[i]] = *v++;
530: }
531: }
532: } else {
533: for (i=0; i<m; i++) {
534: if (indexm[i] < 0) { v += n; continue;}
535: #if defined(PETSC_USE_BOPT_g)
536: if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
537: #endif
538: for (j=0; j<n; j++) {
539: if (indexn[j] < 0) { v++; continue;}
540: #if defined(PETSC_USE_BOPT_g)
541: if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
542: #endif
543: mat->v[indexn[j]*A->m + indexm[i]] += *v++;
544: }
545: }
546: }
547: }
548: return(0);
549: }
551: int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,PetscScalar *v)
552: {
553: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
554: int i,j;
555: PetscScalar *vpt = v;
558: /* row-oriented output */
559: for (i=0; i<m; i++) {
560: for (j=0; j<n; j++) {
561: *vpt++ = mat->v[indexn[j]*A->m + indexm[i]];
562: }
563: }
564: return(0);
565: }
567: /* -----------------------------------------------------------------*/
569: #include petscsys.h
571: EXTERN_C_BEGIN
572: int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A)
573: {
574: Mat_SeqDense *a;
575: Mat B;
576: int *scols,i,j,nz,ierr,fd,header[4],size;
577: int *rowlengths = 0,M,N,*cols;
578: PetscScalar *vals,*svals,*v,*w;
579: MPI_Comm comm = ((PetscObject)viewer)->comm;
582: MPI_Comm_size(comm,&size);
583: if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
584: PetscViewerBinaryGetDescriptor(viewer,&fd);
585: PetscBinaryRead(fd,header,4,PETSC_INT);
586: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
587: M = header[1]; N = header[2]; nz = header[3];
589: if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
590: MatCreateSeqDense(comm,M,N,PETSC_NULL,A);
591: B = *A;
592: a = (Mat_SeqDense*)B->data;
593: v = a->v;
594: /* Allocate some temp space to read in the values and then flip them
595: from row major to column major */
596: PetscMalloc((M*N+1)*sizeof(PetscScalar),&w);
597: /* read in nonzero values */
598: PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);
599: /* now flip the values and store them in the matrix*/
600: for (j=0; j<N; j++) {
601: for (i=0; i<M; i++) {
602: *v++ =w[i*N+j];
603: }
604: }
605: PetscFree(w);
606: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
607: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
608: } else {
609: /* read row lengths */
610: PetscMalloc((M+1)*sizeof(int),&rowlengths);
611: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
613: /* create our matrix */
614: MatCreateSeqDense(comm,M,N,PETSC_NULL,A);
615: B = *A;
616: a = (Mat_SeqDense*)B->data;
617: v = a->v;
619: /* read column indices and nonzeros */
620: PetscMalloc((nz+1)*sizeof(int),&scols);
621: cols = scols;
622: PetscBinaryRead(fd,cols,nz,PETSC_INT);
623: PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);
624: vals = svals;
625: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
627: /* insert into matrix */
628: for (i=0; i<M; i++) {
629: for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
630: svals += rowlengths[i]; scols += rowlengths[i];
631: }
632: PetscFree(vals);
633: PetscFree(cols);
634: PetscFree(rowlengths);
636: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
637: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
638: }
639: return(0);
640: }
641: EXTERN_C_END
643: #include petscsys.h
645: static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
646: {
647: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
648: int ierr,i,j;
649: char *name;
650: PetscScalar *v;
651: PetscViewerFormat format;
654: PetscObjectGetName((PetscObject)A,&name);
655: PetscViewerGetFormat(viewer,&format);
656: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
657: return(0); /* do nothing for now */
658: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
659: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
660: for (i=0; i<A->m; i++) {
661: v = a->v + i;
662: PetscViewerASCIIPrintf(viewer,"row %d:",i);
663: for (j=0; j<A->n; j++) {
664: #if defined(PETSC_USE_COMPLEX)
665: if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
666: PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));
667: } else if (PetscRealPart(*v)) {
668: PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));
669: }
670: #else
671: if (*v) {
672: PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);
673: }
674: #endif
675: v += A->m;
676: }
677: PetscViewerASCIIPrintf(viewer,"n");
678: }
679: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
680: } else {
681: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
682: #if defined(PETSC_USE_COMPLEX)
683: PetscTruth allreal = PETSC_TRUE;
684: /* determine if matrix has all real values */
685: v = a->v;
686: for (i=0; i<A->m*A->n; i++) {
687: if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
688: }
689: #endif
690: if (format == PETSC_VIEWER_ASCII_MATLAB) {
691: PetscObjectGetName((PetscObject)A,&name);
692: PetscViewerASCIIPrintf(viewer,"%% Size = %d %d n",A->m,A->n);
693: PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);n",name,A->m,A->n);
694: PetscViewerASCIIPrintf(viewer,"%s = [n",name);
695: }
697: for (i=0; i<A->m; i++) {
698: v = a->v + i;
699: for (j=0; j<A->n; j++) {
700: #if defined(PETSC_USE_COMPLEX)
701: if (allreal) {
702: PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v)); v += A->m;
703: } else {
704: PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v)); v += A->m;
705: }
706: #else
707: PetscViewerASCIIPrintf(viewer,"%6.4e ",*v); v += A->m;
708: #endif
709: }
710: PetscViewerASCIIPrintf(viewer,"n");
711: }
712: if (format == PETSC_VIEWER_ASCII_MATLAB) {
713: PetscViewerASCIIPrintf(viewer,"];n");
714: }
715: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
716: }
717: PetscViewerFlush(viewer);
718: return(0);
719: }
721: static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
722: {
723: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
724: int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
725: PetscScalar *v,*anonz,*vals;
726: PetscViewerFormat format;
727:
729: PetscViewerBinaryGetDescriptor(viewer,&fd);
731: PetscViewerGetFormat(viewer,&format);
732: if (format == PETSC_VIEWER_BINARY_NATIVE) {
733: /* store the matrix as a dense matrix */
734: PetscMalloc(4*sizeof(int),&col_lens);
735: col_lens[0] = MAT_FILE_COOKIE;
736: col_lens[1] = m;
737: col_lens[2] = n;
738: col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
739: PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);
740: PetscFree(col_lens);
742: /* write out matrix, by rows */
743: PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);
744: v = a->v;
745: for (i=0; i<m; i++) {
746: for (j=0; j<n; j++) {
747: vals[i + j*m] = *v++;
748: }
749: }
750: PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);
751: PetscFree(vals);
752: } else {
753: PetscMalloc((4+nz)*sizeof(int),&col_lens);
754: col_lens[0] = MAT_FILE_COOKIE;
755: col_lens[1] = m;
756: col_lens[2] = n;
757: col_lens[3] = nz;
759: /* store lengths of each row and write (including header) to file */
760: for (i=0; i<m; i++) col_lens[4+i] = n;
761: PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);
763: /* Possibly should write in smaller increments, not whole matrix at once? */
764: /* store column indices (zero start index) */
765: ict = 0;
766: for (i=0; i<m; i++) {
767: for (j=0; j<n; j++) col_lens[ict++] = j;
768: }
769: PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);
770: PetscFree(col_lens);
772: /* store nonzero values */
773: PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);
774: ict = 0;
775: for (i=0; i<m; i++) {
776: v = a->v + i;
777: for (j=0; j<n; j++) {
778: anonz[ict++] = *v; v += A->m;
779: }
780: }
781: PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);
782: PetscFree(anonz);
783: }
784: return(0);
785: }
787: int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
788: {
789: Mat A = (Mat) Aa;
790: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
791: int m = A->m,n = A->n,color,i,j,ierr;
792: PetscScalar *v = a->v;
793: PetscViewer viewer;
794: PetscDraw popup;
795: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
796: PetscViewerFormat format;
800: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
801: PetscViewerGetFormat(viewer,&format);
802: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
804: /* Loop over matrix elements drawing boxes */
805: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
806: /* Blue for negative and Red for positive */
807: color = PETSC_DRAW_BLUE;
808: for(j = 0; j < n; j++) {
809: x_l = j;
810: x_r = x_l + 1.0;
811: for(i = 0; i < m; i++) {
812: y_l = m - i - 1.0;
813: y_r = y_l + 1.0;
814: #if defined(PETSC_USE_COMPLEX)
815: if (PetscRealPart(v[j*m+i]) > 0.) {
816: color = PETSC_DRAW_RED;
817: } else if (PetscRealPart(v[j*m+i]) < 0.) {
818: color = PETSC_DRAW_BLUE;
819: } else {
820: continue;
821: }
822: #else
823: if (v[j*m+i] > 0.) {
824: color = PETSC_DRAW_RED;
825: } else if (v[j*m+i] < 0.) {
826: color = PETSC_DRAW_BLUE;
827: } else {
828: continue;
829: }
830: #endif
831: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
832: }
833: }
834: } else {
835: /* use contour shading to indicate magnitude of values */
836: /* first determine max of all nonzero values */
837: for(i = 0; i < m*n; i++) {
838: if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
839: }
840: scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
841: ierr = PetscDrawGetPopup(draw,&popup);
842: if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);}
843: for(j = 0; j < n; j++) {
844: x_l = j;
845: x_r = x_l + 1.0;
846: for(i = 0; i < m; i++) {
847: y_l = m - i - 1.0;
848: y_r = y_l + 1.0;
849: color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
850: ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
851: }
852: }
853: }
854: return(0);
855: }
857: int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
858: {
859: PetscDraw draw;
860: PetscTruth isnull;
861: PetscReal xr,yr,xl,yl,h,w;
862: int ierr;
865: PetscViewerDrawGetDraw(viewer,0,&draw);
866: PetscDrawIsNull(draw,&isnull);
867: if (isnull == PETSC_TRUE) return(0);
869: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
870: xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
871: xr += w; yr += h; xl = -w; yl = -h;
872: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
873: PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
874: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
875: return(0);
876: }
878: int MatView_SeqDense(Mat A,PetscViewer viewer)
879: {
880: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
881: int ierr;
882: PetscTruth issocket,isascii,isbinary,isdraw;
885: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
886: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
887: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
888: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
890: if (issocket) {
891: PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);
892: } else if (isascii) {
893: MatView_SeqDense_ASCII(A,viewer);
894: } else if (isbinary) {
895: MatView_SeqDense_Binary(A,viewer);
896: } else if (isdraw) {
897: MatView_SeqDense_Draw(A,viewer);
898: } else {
899: SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
900: }
901: return(0);
902: }
904: int MatDestroy_SeqDense(Mat mat)
905: {
906: Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
907: int ierr;
910: #if defined(PETSC_USE_LOG)
911: PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
912: #endif
913: if (l->pivots) {PetscFree(l->pivots);}
914: if (!l->user_alloc) {PetscFree(l->v);}
915: PetscFree(l);
916: return(0);
917: }
919: int MatTranspose_SeqDense(Mat A,Mat *matout)
920: {
921: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
922: int k,j,m,n,ierr;
923: PetscScalar *v,tmp;
926: v = mat->v; m = A->m; n = A->n;
927: if (!matout) { /* in place transpose */
928: if (m != n) { /* malloc temp to hold transpose */
929: PetscScalar *w;
930: PetscMalloc((m+1)*(n+1)*sizeof(PetscScalar),&w);
931: for (j=0; j<m; j++) {
932: for (k=0; k<n; k++) {
933: w[k + j*n] = v[j + k*m];
934: }
935: }
936: PetscMemcpy(v,w,m*n*sizeof(PetscScalar));
937: PetscFree(w);
938: } else {
939: for (j=0; j<m; j++) {
940: for (k=0; k<j; k++) {
941: tmp = v[j + k*n];
942: v[j + k*n] = v[k + j*n];
943: v[k + j*n] = tmp;
944: }
945: }
946: }
947: } else { /* out-of-place transpose */
948: Mat tmat;
949: Mat_SeqDense *tmatd;
950: PetscScalar *v2;
952: ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);
953: tmatd = (Mat_SeqDense*)tmat->data;
954: v = mat->v; v2 = tmatd->v;
955: for (j=0; j<n; j++) {
956: for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m];
957: }
958: MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
959: MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
960: *matout = tmat;
961: }
962: return(0);
963: }
965: int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
966: {
967: Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
968: Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
969: int ierr,i;
970: PetscScalar *v1 = mat1->v,*v2 = mat2->v;
971: PetscTruth flag;
974: PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);
975: if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type MATSEQDENSE");
976: if (A1->m != A2->m) {*flg = PETSC_FALSE; return(0);}
977: if (A1->n != A2->n) {*flg = PETSC_FALSE; return(0);}
978: for (i=0; i<A1->m*A1->n; i++) {
979: if (*v1 != *v2) {*flg = PETSC_FALSE; return(0);}
980: v1++; v2++;
981: }
982: *flg = PETSC_TRUE;
983: return(0);
984: }
986: int MatGetDiagonal_SeqDense(Mat A,Vec v)
987: {
988: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
989: int ierr,i,n,len;
990: PetscScalar *x,zero = 0.0;
993: VecSet(&zero,v);
994: VecGetSize(v,&n);
995: VecGetArray(v,&x);
996: len = PetscMin(A->m,A->n);
997: if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
998: for (i=0; i<len; i++) {
999: x[i] = mat->v[i*A->m + i];
1000: }
1001: VecRestoreArray(v,&x);
1002: return(0);
1003: }
1005: int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1006: {
1007: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1008: PetscScalar *l,*r,x,*v;
1009: int ierr,i,j,m = A->m,n = A->n;
1012: if (ll) {
1013: VecGetSize(ll,&m);
1014: VecGetArray(ll,&l);
1015: if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1016: for (i=0; i<m; i++) {
1017: x = l[i];
1018: v = mat->v + i;
1019: for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1020: }
1021: VecRestoreArray(ll,&l);
1022: PetscLogFlops(n*m);
1023: }
1024: if (rr) {
1025: VecGetSize(rr,&n);
1026: VecGetArray(rr,&r);
1027: if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1028: for (i=0; i<n; i++) {
1029: x = r[i];
1030: v = mat->v + i*m;
1031: for (j=0; j<m; j++) { (*v++) *= x;}
1032: }
1033: VecRestoreArray(rr,&r);
1034: PetscLogFlops(n*m);
1035: }
1036: return(0);
1037: }
1039: int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1040: {
1041: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1042: PetscScalar *v = mat->v;
1043: PetscReal sum = 0.0;
1044: int i,j;
1047: if (type == NORM_FROBENIUS) {
1048: for (i=0; i<A->n*A->m; i++) {
1049: #if defined(PETSC_USE_COMPLEX)
1050: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1051: #else
1052: sum += (*v)*(*v); v++;
1053: #endif
1054: }
1055: *nrm = sqrt(sum);
1056: PetscLogFlops(2*A->n*A->m);
1057: } else if (type == NORM_1) {
1058: *nrm = 0.0;
1059: for (j=0; j<A->n; j++) {
1060: sum = 0.0;
1061: for (i=0; i<A->m; i++) {
1062: sum += PetscAbsScalar(*v); v++;
1063: }
1064: if (sum > *nrm) *nrm = sum;
1065: }
1066: PetscLogFlops(A->n*A->m);
1067: } else if (type == NORM_INFINITY) {
1068: *nrm = 0.0;
1069: for (j=0; j<A->m; j++) {
1070: v = mat->v + j;
1071: sum = 0.0;
1072: for (i=0; i<A->n; i++) {
1073: sum += PetscAbsScalar(*v); v += A->m;
1074: }
1075: if (sum > *nrm) *nrm = sum;
1076: }
1077: PetscLogFlops(A->n*A->m);
1078: } else {
1079: SETERRQ(PETSC_ERR_SUP,"No two norm");
1080: }
1081: return(0);
1082: }
1084: int MatSetOption_SeqDense(Mat A,MatOption op)
1085: {
1086: Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
1087:
1089: switch (op) {
1090: case MAT_ROW_ORIENTED:
1091: aij->roworiented = PETSC_TRUE;
1092: break;
1093: case MAT_COLUMN_ORIENTED:
1094: aij->roworiented = PETSC_FALSE;
1095: break;
1096: case MAT_ROWS_SORTED:
1097: case MAT_ROWS_UNSORTED:
1098: case MAT_COLUMNS_SORTED:
1099: case MAT_COLUMNS_UNSORTED:
1100: case MAT_NO_NEW_NONZERO_LOCATIONS:
1101: case MAT_YES_NEW_NONZERO_LOCATIONS:
1102: case MAT_NEW_NONZERO_LOCATION_ERR:
1103: case MAT_NO_NEW_DIAGONALS:
1104: case MAT_YES_NEW_DIAGONALS:
1105: case MAT_IGNORE_OFF_PROC_ENTRIES:
1106: case MAT_USE_HASH_TABLE:
1107: case MAT_USE_SINGLE_PRECISION_SOLVES:
1108: PetscLogInfo(A,"MatSetOption_SeqDense:Option ignoredn");
1109: break;
1110: default:
1111: SETERRQ(PETSC_ERR_SUP,"unknown option");
1112: }
1113: return(0);
1114: }
1116: int MatZeroEntries_SeqDense(Mat A)
1117: {
1118: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1119: int ierr;
1122: PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));
1123: return(0);
1124: }
1126: int MatGetBlockSize_SeqDense(Mat A,int *bs)
1127: {
1129: *bs = 1;
1130: return(0);
1131: }
1133: int MatZeroRows_SeqDense(Mat A,IS is,PetscScalar *diag)
1134: {
1135: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1136: int n = A->n,i,j,ierr,N,*rows;
1137: PetscScalar *slot;
1140: ISGetLocalSize(is,&N);
1141: ISGetIndices(is,&rows);
1142: for (i=0; i<N; i++) {
1143: slot = l->v + rows[i];
1144: for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
1145: }
1146: if (diag) {
1147: for (i=0; i<N; i++) {
1148: slot = l->v + (n+1)*rows[i];
1149: *slot = *diag;
1150: }
1151: }
1152: ISRestoreIndices(is,&rows);
1153: return(0);
1154: }
1156: int MatGetArray_SeqDense(Mat A,PetscScalar **array)
1157: {
1158: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1161: *array = mat->v;
1162: return(0);
1163: }
1165: int MatRestoreArray_SeqDense(Mat A,PetscScalar **array)
1166: {
1168: *array = 0; /* user cannot accidently use the array later */
1169: return(0);
1170: }
1172: static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
1173: {
1174: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1175: int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
1176: PetscScalar *av,*bv,*v = mat->v;
1177: Mat newmat;
1180: ISGetIndices(isrow,&irow);
1181: ISGetIndices(iscol,&icol);
1182: ISGetLocalSize(isrow,&nrows);
1183: ISGetLocalSize(iscol,&ncols);
1184:
1185: /* Check submatrixcall */
1186: if (scall == MAT_REUSE_MATRIX) {
1187: int n_cols,n_rows;
1188: MatGetSize(*B,&n_rows,&n_cols);
1189: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1190: newmat = *B;
1191: } else {
1192: /* Create and fill new matrix */
1193: MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);
1194: }
1196: /* Now extract the data pointers and do the copy,column at a time */
1197: bv = ((Mat_SeqDense*)newmat->data)->v;
1198:
1199: for (i=0; i<ncols; i++) {
1200: av = v + m*icol[i];
1201: for (j=0; j<nrows; j++) {
1202: *bv++ = av[irow[j]];
1203: }
1204: }
1206: /* Assemble the matrices so that the correct flags are set */
1207: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1208: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1210: /* Free work space */
1211: ISRestoreIndices(isrow,&irow);
1212: ISRestoreIndices(iscol,&icol);
1213: *B = newmat;
1214: return(0);
1215: }
1217: int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1218: {
1219: int ierr,i;
1222: if (scall == MAT_INITIAL_MATRIX) {
1223: PetscMalloc((n+1)*sizeof(Mat),B);
1224: }
1226: for (i=0; i<n; i++) {
1227: MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1228: }
1229: return(0);
1230: }
1232: int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
1233: {
1234: Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1235: int ierr;
1236: PetscTruth flag;
1239: PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);
1240: if (!flag) {
1241: MatCopy_Basic(A,B,str);
1242: return(0);
1243: }
1244: if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1245: PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));
1246: return(0);
1247: }
1249: int MatSetUpPreallocation_SeqDense(Mat A)
1250: {
1251: int ierr;
1254: MatSeqDenseSetPreallocation(A,0);
1255: return(0);
1256: }
1258: /* -------------------------------------------------------------------*/
1259: static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1260: MatGetRow_SeqDense,
1261: MatRestoreRow_SeqDense,
1262: MatMult_SeqDense,
1263: MatMultAdd_SeqDense,
1264: MatMultTranspose_SeqDense,
1265: MatMultTransposeAdd_SeqDense,
1266: MatSolve_SeqDense,
1267: MatSolveAdd_SeqDense,
1268: MatSolveTranspose_SeqDense,
1269: MatSolveTransposeAdd_SeqDense,
1270: MatLUFactor_SeqDense,
1271: MatCholeskyFactor_SeqDense,
1272: MatRelax_SeqDense,
1273: MatTranspose_SeqDense,
1274: MatGetInfo_SeqDense,
1275: MatEqual_SeqDense,
1276: MatGetDiagonal_SeqDense,
1277: MatDiagonalScale_SeqDense,
1278: MatNorm_SeqDense,
1279: 0,
1280: 0,
1281: 0,
1282: MatSetOption_SeqDense,
1283: MatZeroEntries_SeqDense,
1284: MatZeroRows_SeqDense,
1285: MatLUFactorSymbolic_SeqDense,
1286: MatLUFactorNumeric_SeqDense,
1287: MatCholeskyFactorSymbolic_SeqDense,
1288: MatCholeskyFactorNumeric_SeqDense,
1289: MatSetUpPreallocation_SeqDense,
1290: 0,
1291: 0,
1292: MatGetArray_SeqDense,
1293: MatRestoreArray_SeqDense,
1294: MatDuplicate_SeqDense,
1295: 0,
1296: 0,
1297: 0,
1298: 0,
1299: MatAXPY_SeqDense,
1300: MatGetSubMatrices_SeqDense,
1301: 0,
1302: MatGetValues_SeqDense,
1303: MatCopy_SeqDense,
1304: 0,
1305: MatScale_SeqDense,
1306: 0,
1307: 0,
1308: 0,
1309: MatGetBlockSize_SeqDense,
1310: 0,
1311: 0,
1312: 0,
1313: 0,
1314: 0,
1315: 0,
1316: 0,
1317: 0,
1318: 0,
1319: 0,
1320: MatDestroy_SeqDense,
1321: MatView_SeqDense,
1322: MatGetPetscMaps_Petsc};
1324: /*@C
1325: MatCreateSeqDense - Creates a sequential dense matrix that
1326: is stored in column major order (the usual Fortran 77 manner). Many
1327: of the matrix operations use the BLAS and LAPACK routines.
1329: Collective on MPI_Comm
1331: Input Parameters:
1332: + comm - MPI communicator, set to PETSC_COMM_SELF
1333: . m - number of rows
1334: . n - number of columns
1335: - data - optional location of matrix data. Set data=PETSC_NULL for PETSc
1336: to control all matrix memory allocation.
1338: Output Parameter:
1339: . A - the matrix
1341: Notes:
1342: The data input variable is intended primarily for Fortran programmers
1343: who wish to allocate their own matrix memory space. Most users should
1344: set data=PETSC_NULL.
1346: Level: intermediate
1348: .keywords: dense, matrix, LAPACK, BLAS
1350: .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1351: @*/
1352: int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A)
1353: {
1357: MatCreate(comm,m,n,m,n,A);
1358: MatSetType(*A,MATSEQDENSE);
1359: MatSeqDenseSetPreallocation(*A,data);
1360: return(0);
1361: }
1363: /*@C
1364: MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1366: Collective on MPI_Comm
1368: Input Parameters:
1369: + A - the matrix
1370: - data - the array (or PETSC_NULL)
1372: Notes:
1373: The data input variable is intended primarily for Fortran programmers
1374: who wish to allocate their own matrix memory space. Most users should
1375: set data=PETSC_NULL.
1377: Level: intermediate
1379: .keywords: dense, matrix, LAPACK, BLAS
1381: .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1382: @*/
1383: int MatSeqDenseSetPreallocation(Mat B,PetscScalar *data)
1384: {
1385: Mat_SeqDense *b;
1386: int ierr;
1387: PetscTruth flg2;
1390: PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);
1391: if (!flg2) return(0);
1392: B->preallocated = PETSC_TRUE;
1393: b = (Mat_SeqDense*)B->data;
1394: if (!data) {
1395: ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);
1396: ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));
1397: b->user_alloc = PETSC_FALSE;
1398: PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));
1399: } else { /* user-allocated storage */
1400: b->v = data;
1401: b->user_alloc = PETSC_TRUE;
1402: }
1403: return(0);
1404: }
1406: EXTERN_C_BEGIN
1407: int MatCreate_SeqDense(Mat B)
1408: {
1409: Mat_SeqDense *b;
1410: int ierr,size;
1413: MPI_Comm_size(B->comm,&size);
1414: if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1416: B->m = B->M = PetscMax(B->m,B->M);
1417: B->n = B->N = PetscMax(B->n,B->N);
1419: ierr = PetscNew(Mat_SeqDense,&b);
1420: ierr = PetscMemzero(b,sizeof(Mat_SeqDense));
1421: ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1422: B->factor = 0;
1423: B->mapping = 0;
1424: PetscLogObjectMemory(B,sizeof(struct _p_Mat));
1425: B->data = (void*)b;
1427: PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
1428: PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);
1430: b->pivots = 0;
1431: b->roworiented = PETSC_TRUE;
1432: b->v = 0;
1434: return(0);
1435: }
1436: EXTERN_C_END