Actual source code: dense.c
petsc-3.13.0 2020-03-29
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: #include <../src/mat/impls/aij/seq/aij.h>
11: PetscErrorCode MatSeqDenseSymmetrize_Private(Mat A, PetscBool hermitian)
12: {
13: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
14: PetscInt j, k, n = A->rmap->n;
15: PetscScalar *v;
19: if (A->rmap->n != A->cmap->n) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot symmetrize a rectangular matrix");
20: MatDenseGetArray(A,&v);
21: if (!hermitian) {
22: for (k=0;k<n;k++) {
23: for (j=k;j<n;j++) {
24: v[j*mat->lda + k] = v[k*mat->lda + j];
25: }
26: }
27: } else {
28: for (k=0;k<n;k++) {
29: for (j=k;j<n;j++) {
30: v[j*mat->lda + k] = PetscConj(v[k*mat->lda + j]);
31: }
32: }
33: }
34: MatDenseRestoreArray(A,&v);
35: return(0);
36: }
38: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat A)
39: {
40: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
42: PetscBLASInt info,n;
45: if (!A->rmap->n || !A->cmap->n) return(0);
46: PetscBLASIntCast(A->cmap->n,&n);
47: if (A->factortype == MAT_FACTOR_LU) {
48: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
49: if (!mat->fwork) {
50: mat->lfwork = n;
51: PetscMalloc1(mat->lfwork,&mat->fwork);
52: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
53: }
54: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
55: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
56: PetscFPTrapPop();
57: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
58: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
59: if (A->spd) {
60: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
61: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&n,mat->v,&mat->lda,&info));
62: PetscFPTrapPop();
63: MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
64: #if defined(PETSC_USE_COMPLEX)
65: } else if (A->hermitian) {
66: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
67: if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
68: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
69: PetscStackCallBLAS("LAPACKhetri",LAPACKhetri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
70: PetscFPTrapPop();
71: MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
72: #endif
73: } else { /* symmetric case */
74: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
75: if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
76: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
77: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
78: PetscFPTrapPop();
79: MatSeqDenseSymmetrize_Private(A,PETSC_FALSE);
80: }
81: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad Inversion: zero pivot in row %D",(PetscInt)info-1);
82: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
83: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
85: A->ops->solve = NULL;
86: A->ops->matsolve = NULL;
87: A->ops->solvetranspose = NULL;
88: A->ops->matsolvetranspose = NULL;
89: A->ops->solveadd = NULL;
90: A->ops->solvetransposeadd = NULL;
91: A->factortype = MAT_FACTOR_NONE;
92: PetscFree(A->solvertype);
93: return(0);
94: }
96: PetscErrorCode MatZeroRowsColumns_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
97: {
98: PetscErrorCode ierr;
99: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
100: PetscInt m = l->lda, n = A->cmap->n,r = A->rmap->n, i,j;
101: PetscScalar *slot,*bb,*v;
102: const PetscScalar *xx;
105: #if defined(PETSC_USE_DEBUG)
106: for (i=0; i<N; i++) {
107: if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
108: if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n);
109: if (rows[i] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Col %D requested to be zeroed greater than or equal number of cols %D",rows[i],A->cmap->n);
110: }
111: #endif
112: if (!N) return(0);
114: /* fix right hand side if needed */
115: if (x && b) {
116: Vec xt;
118: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
119: VecDuplicate(x,&xt);
120: VecCopy(x,xt);
121: VecScale(xt,-1.0);
122: MatMultAdd(A,xt,b,b);
123: VecDestroy(&xt);
124: VecGetArrayRead(x,&xx);
125: VecGetArray(b,&bb);
126: for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
127: VecRestoreArrayRead(x,&xx);
128: VecRestoreArray(b,&bb);
129: }
131: MatDenseGetArray(A,&v);
132: for (i=0; i<N; i++) {
133: slot = v + rows[i]*m;
134: PetscArrayzero(slot,r);
135: }
136: for (i=0; i<N; i++) {
137: slot = v + rows[i];
138: for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
139: }
140: if (diag != 0.0) {
141: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
142: for (i=0; i<N; i++) {
143: slot = v + (m+1)*rows[i];
144: *slot = diag;
145: }
146: }
147: MatDenseRestoreArray(A,&v);
148: return(0);
149: }
151: PetscErrorCode MatPtAPNumeric_SeqDense_SeqDense(Mat A,Mat P,Mat C)
152: {
153: Mat_SeqDense *c = (Mat_SeqDense*)(C->data);
157: if (c->ptapwork) {
158: (*C->ops->matmultnumeric)(A,P,c->ptapwork);
159: (*C->ops->transposematmultnumeric)(P,c->ptapwork,C);
160: } else SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Must call MatPtAPSymbolic_SeqDense_SeqDense() first");
161: return(0);
162: }
164: PetscErrorCode MatPtAPSymbolic_SeqDense_SeqDense(Mat A,Mat P,PetscReal fill,Mat C)
165: {
166: Mat_SeqDense *c;
167: PetscBool flg;
171: PetscObjectTypeCompare((PetscObject)P,((PetscObject)A)->type_name,&flg);
172: MatSetSizes(C,P->cmap->n,P->cmap->n,P->cmap->N,P->cmap->N);
173: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
174: MatSeqDenseSetPreallocation(C,NULL);
175: c = (Mat_SeqDense*)C->data;
176: MatCreate(PetscObjectComm((PetscObject)A),&c->ptapwork);
177: MatSetSizes(c->ptapwork,A->rmap->n,P->cmap->n,A->rmap->N,P->cmap->N);
178: MatSetType(c->ptapwork,flg ? ((PetscObject)A)->type_name : MATDENSE);
179: MatSeqDenseSetPreallocation(c->ptapwork,NULL);
180: return(0);
181: }
183: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
184: {
185: Mat B = NULL;
186: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
187: Mat_SeqDense *b;
189: PetscInt *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->N,i;
190: MatScalar *av=a->a;
191: PetscBool isseqdense;
194: if (reuse == MAT_REUSE_MATRIX) {
195: PetscObjectTypeCompare((PetscObject)*newmat,MATSEQDENSE,&isseqdense);
196: if (!isseqdense) SETERRQ1(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type_name);
197: }
198: if (reuse != MAT_REUSE_MATRIX) {
199: MatCreate(PetscObjectComm((PetscObject)A),&B);
200: MatSetSizes(B,m,n,m,n);
201: MatSetType(B,MATSEQDENSE);
202: MatSeqDenseSetPreallocation(B,NULL);
203: b = (Mat_SeqDense*)(B->data);
204: } else {
205: b = (Mat_SeqDense*)((*newmat)->data);
206: PetscArrayzero(b->v,m*n);
207: }
208: for (i=0; i<m; i++) {
209: PetscInt j;
210: for (j=0;j<ai[1]-ai[0];j++) {
211: b->v[*aj*m+i] = *av;
212: aj++;
213: av++;
214: }
215: ai++;
216: }
218: if (reuse == MAT_INPLACE_MATRIX) {
219: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
220: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
221: MatHeaderReplace(A,&B);
222: } else {
223: if (B) *newmat = B;
224: MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
225: MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
226: }
227: return(0);
228: }
230: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
231: {
232: Mat B;
233: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
235: PetscInt i, j;
236: PetscInt *rows, *nnz;
237: MatScalar *aa = a->v, *vals;
240: MatCreate(PetscObjectComm((PetscObject)A),&B);
241: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
242: MatSetType(B,MATSEQAIJ);
243: PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals);
244: for (j=0; j<A->cmap->n; j++) {
245: for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) ++nnz[i];
246: aa += a->lda;
247: }
248: MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz);
249: aa = a->v;
250: for (j=0; j<A->cmap->n; j++) {
251: PetscInt numRows = 0;
252: for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) {rows[numRows] = i; vals[numRows++] = aa[i];}
253: MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES);
254: aa += a->lda;
255: }
256: PetscFree3(rows,nnz,vals);
257: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
258: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
260: if (reuse == MAT_INPLACE_MATRIX) {
261: MatHeaderReplace(A,&B);
262: } else {
263: *newmat = B;
264: }
265: return(0);
266: }
268: PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
269: {
270: Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
271: const PetscScalar *xv;
272: PetscScalar *yv;
273: PetscBLASInt N,m,ldax,lday,one = 1;
274: PetscErrorCode ierr;
277: MatDenseGetArrayRead(X,&xv);
278: MatDenseGetArray(Y,&yv);
279: PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);
280: PetscBLASIntCast(X->rmap->n,&m);
281: PetscBLASIntCast(x->lda,&ldax);
282: PetscBLASIntCast(y->lda,&lday);
283: if (ldax>m || lday>m) {
284: PetscInt j;
286: for (j=0; j<X->cmap->n; j++) {
287: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&alpha,xv+j*ldax,&one,yv+j*lday,&one));
288: }
289: } else {
290: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&alpha,xv,&one,yv,&one));
291: }
292: MatDenseRestoreArrayRead(X,&xv);
293: MatDenseRestoreArray(Y,&yv);
294: PetscLogFlops(PetscMax(2*N-1,0));
295: return(0);
296: }
298: static PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
299: {
300: PetscLogDouble N = A->rmap->n*A->cmap->n;
303: info->block_size = 1.0;
304: info->nz_allocated = N;
305: info->nz_used = N;
306: info->nz_unneeded = 0;
307: info->assemblies = A->num_ass;
308: info->mallocs = 0;
309: info->memory = ((PetscObject)A)->mem;
310: info->fill_ratio_given = 0;
311: info->fill_ratio_needed = 0;
312: info->factor_mallocs = 0;
313: return(0);
314: }
316: static PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
317: {
318: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
319: PetscScalar *v;
321: PetscBLASInt one = 1,j,nz,lda;
324: MatDenseGetArray(A,&v);
325: PetscBLASIntCast(a->lda,&lda);
326: if (lda>A->rmap->n) {
327: PetscBLASIntCast(A->rmap->n,&nz);
328: for (j=0; j<A->cmap->n; j++) {
329: PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v+j*lda,&one));
330: }
331: } else {
332: PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);
333: PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v,&one));
334: }
335: PetscLogFlops(nz);
336: MatDenseRestoreArray(A,&v);
337: return(0);
338: }
340: static PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl)
341: {
342: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
343: PetscInt i,j,m = A->rmap->n,N = a->lda;
344: const PetscScalar *v;
345: PetscErrorCode ierr;
348: *fl = PETSC_FALSE;
349: if (A->rmap->n != A->cmap->n) return(0);
350: MatDenseGetArrayRead(A,&v);
351: for (i=0; i<m; i++) {
352: for (j=i; j<m; j++) {
353: if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) return(0);
354: }
355: }
356: MatDenseRestoreArrayRead(A,&v);
357: *fl = PETSC_TRUE;
358: return(0);
359: }
361: PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
362: {
363: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
365: PetscInt lda = (PetscInt)mat->lda,j,m;
368: PetscLayoutReference(A->rmap,&newi->rmap);
369: PetscLayoutReference(A->cmap,&newi->cmap);
370: MatSeqDenseSetPreallocation(newi,NULL);
371: if (cpvalues == MAT_COPY_VALUES) {
372: const PetscScalar *av;
373: PetscScalar *v;
375: MatDenseGetArrayRead(A,&av);
376: MatDenseGetArray(newi,&v);
377: if (lda>A->rmap->n) {
378: m = A->rmap->n;
379: for (j=0; j<A->cmap->n; j++) {
380: PetscArraycpy(v+j*m,av+j*lda,m);
381: }
382: } else {
383: PetscArraycpy(v,av,A->rmap->n*A->cmap->n);
384: }
385: MatDenseRestoreArray(newi,&v);
386: MatDenseRestoreArrayRead(A,&av);
387: }
388: return(0);
389: }
391: PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
392: {
396: MatCreate(PetscObjectComm((PetscObject)A),newmat);
397: MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
398: MatSetType(*newmat,((PetscObject)A)->type_name);
399: MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);
400: return(0);
401: }
403: static PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
404: {
405: MatFactorInfo info;
409: MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
410: (*fact->ops->lufactor)(fact,0,0,&info);
411: return(0);
412: }
414: static PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
415: {
416: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
417: PetscErrorCode ierr;
418: const PetscScalar *x;
419: PetscScalar *y;
420: PetscBLASInt one = 1,info,m;
423: PetscBLASIntCast(A->rmap->n,&m);
424: VecGetArrayRead(xx,&x);
425: VecGetArray(yy,&y);
426: PetscArraycpy(y,x,A->rmap->n);
427: if (A->factortype == MAT_FACTOR_LU) {
428: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
429: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
430: PetscFPTrapPop();
431: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
432: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
433: if (A->spd) {
434: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
435: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
436: PetscFPTrapPop();
437: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
438: #if defined(PETSC_USE_COMPLEX)
439: } else if (A->hermitian) {
440: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
441: PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
442: PetscFPTrapPop();
443: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
444: #endif
445: } else { /* symmetric case */
446: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
447: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
448: PetscFPTrapPop();
449: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
450: }
451: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
452: VecRestoreArrayRead(xx,&x);
453: VecRestoreArray(yy,&y);
454: PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
455: return(0);
456: }
458: static PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
459: {
460: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
461: PetscErrorCode ierr;
462: const PetscScalar *b;
463: PetscScalar *x;
464: PetscInt n;
465: PetscBLASInt nrhs,info,m;
468: PetscBLASIntCast(A->rmap->n,&m);
469: MatGetSize(B,NULL,&n);
470: PetscBLASIntCast(n,&nrhs);
471: MatDenseGetArrayRead(B,&b);
472: MatDenseGetArray(X,&x);
474: PetscArraycpy(x,b,m*nrhs);
476: if (A->factortype == MAT_FACTOR_LU) {
477: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
478: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
479: PetscFPTrapPop();
480: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
481: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
482: if (A->spd) {
483: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
484: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
485: PetscFPTrapPop();
486: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
487: #if defined(PETSC_USE_COMPLEX)
488: } else if (A->hermitian) {
489: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
490: PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
491: PetscFPTrapPop();
492: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
493: #endif
494: } else { /* symmetric case */
495: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
496: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
497: PetscFPTrapPop();
498: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
499: }
500: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
502: MatDenseRestoreArrayRead(B,&b);
503: MatDenseRestoreArray(X,&x);
504: PetscLogFlops(nrhs*(2.0*m*m - m));
505: return(0);
506: }
508: static PetscErrorCode MatConjugate_SeqDense(Mat);
510: static PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
511: {
512: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
513: PetscErrorCode ierr;
514: const PetscScalar *x;
515: PetscScalar *y;
516: PetscBLASInt one = 1,info,m;
519: PetscBLASIntCast(A->rmap->n,&m);
520: VecGetArrayRead(xx,&x);
521: VecGetArray(yy,&y);
522: PetscArraycpy(y,x,A->rmap->n);
523: if (A->factortype == MAT_FACTOR_LU) {
524: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
525: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
526: PetscFPTrapPop();
527: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
528: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
529: if (A->spd) {
530: #if defined(PETSC_USE_COMPLEX)
531: MatConjugate_SeqDense(A);
532: #endif
533: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
534: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
535: PetscFPTrapPop();
536: #if defined(PETSC_USE_COMPLEX)
537: MatConjugate_SeqDense(A);
538: #endif
539: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
540: #if defined(PETSC_USE_COMPLEX)
541: } else if (A->hermitian) {
542: MatConjugate_SeqDense(A);
543: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
544: PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
545: PetscFPTrapPop();
546: MatConjugate_SeqDense(A);
547: #endif
548: } else { /* symmetric case */
549: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
550: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
551: PetscFPTrapPop();
552: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
553: }
554: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
555: VecRestoreArrayRead(xx,&x);
556: VecRestoreArray(yy,&y);
557: PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
558: return(0);
559: }
561: /* ---------------------------------------------------------------*/
562: /* COMMENT: I have chosen to hide row permutation in the pivots,
563: rather than put it in the Mat->row slot.*/
564: PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
565: {
566: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
568: PetscBLASInt n,m,info;
571: PetscBLASIntCast(A->cmap->n,&n);
572: PetscBLASIntCast(A->rmap->n,&m);
573: if (!mat->pivots) {
574: PetscMalloc1(A->rmap->n,&mat->pivots);
575: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
576: }
577: if (!A->rmap->n || !A->cmap->n) return(0);
578: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
579: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info));
580: PetscFPTrapPop();
582: if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
583: if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
585: A->ops->solve = MatSolve_SeqDense;
586: A->ops->matsolve = MatMatSolve_SeqDense;
587: A->ops->solvetranspose = MatSolveTranspose_SeqDense;
588: A->factortype = MAT_FACTOR_LU;
590: PetscFree(A->solvertype);
591: PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);
593: PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);
594: return(0);
595: }
597: /* Cholesky as L*L^T or L*D*L^T and the symmetric/hermitian complex variants */
598: PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
599: {
600: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
602: PetscBLASInt info,n;
605: PetscBLASIntCast(A->cmap->n,&n);
606: if (!A->rmap->n || !A->cmap->n) return(0);
607: if (A->spd) {
608: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
609: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info));
610: PetscFPTrapPop();
611: #if defined(PETSC_USE_COMPLEX)
612: } else if (A->hermitian) {
613: if (!mat->pivots) {
614: PetscMalloc1(A->rmap->n,&mat->pivots);
615: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
616: }
617: if (!mat->fwork) {
618: PetscScalar dummy;
620: mat->lfwork = -1;
621: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
622: PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
623: PetscFPTrapPop();
624: mat->lfwork = (PetscInt)PetscRealPart(dummy);
625: PetscMalloc1(mat->lfwork,&mat->fwork);
626: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
627: }
628: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
629: PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
630: PetscFPTrapPop();
631: #endif
632: } else { /* symmetric case */
633: if (!mat->pivots) {
634: PetscMalloc1(A->rmap->n,&mat->pivots);
635: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
636: }
637: if (!mat->fwork) {
638: PetscScalar dummy;
640: mat->lfwork = -1;
641: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
642: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
643: PetscFPTrapPop();
644: mat->lfwork = (PetscInt)PetscRealPart(dummy);
645: PetscMalloc1(mat->lfwork,&mat->fwork);
646: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
647: }
648: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
649: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
650: PetscFPTrapPop();
651: }
652: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
654: A->ops->solve = MatSolve_SeqDense;
655: A->ops->matsolve = MatMatSolve_SeqDense;
656: A->ops->solvetranspose = MatSolveTranspose_SeqDense;
657: A->factortype = MAT_FACTOR_CHOLESKY;
659: PetscFree(A->solvertype);
660: PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);
662: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
663: return(0);
664: }
667: PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
668: {
670: MatFactorInfo info;
673: info.fill = 1.0;
675: MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
676: (*fact->ops->choleskyfactor)(fact,0,&info);
677: return(0);
678: }
680: PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
681: {
683: fact->assembled = PETSC_TRUE;
684: fact->preallocated = PETSC_TRUE;
685: fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
686: fact->ops->solve = MatSolve_SeqDense;
687: fact->ops->matsolve = MatMatSolve_SeqDense;
688: fact->ops->solvetranspose = MatSolveTranspose_SeqDense;
689: return(0);
690: }
692: PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
693: {
695: fact->preallocated = PETSC_TRUE;
696: fact->assembled = PETSC_TRUE;
697: fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
698: fact->ops->solve = MatSolve_SeqDense;
699: fact->ops->matsolve = MatMatSolve_SeqDense;
700: fact->ops->solvetranspose = MatSolveTranspose_SeqDense;
701: return(0);
702: }
704: /* uses LAPACK */
705: PETSC_INTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
706: {
710: MatCreate(PetscObjectComm((PetscObject)A),fact);
711: MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
712: MatSetType(*fact,MATDENSE);
713: if (ftype == MAT_FACTOR_LU) {
714: (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
715: } else {
716: (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
717: }
718: (*fact)->factortype = ftype;
720: PetscFree((*fact)->solvertype);
721: PetscStrallocpy(MATSOLVERPETSC,&(*fact)->solvertype);
722: return(0);
723: }
725: /* ------------------------------------------------------------------*/
726: static PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
727: {
728: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
729: PetscScalar *x,*v = mat->v,zero = 0.0,xt;
730: const PetscScalar *b;
731: PetscErrorCode ierr;
732: PetscInt m = A->rmap->n,i;
733: PetscBLASInt o = 1,bm;
736: #if defined(PETSC_HAVE_CUDA)
737: if (A->offloadmask == PETSC_OFFLOAD_GPU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented");
738: #endif
739: if (shift == -1) shift = 0.0; /* negative shift indicates do not error on zero diagonal; this code never zeros on zero diagonal */
740: PetscBLASIntCast(m,&bm);
741: if (flag & SOR_ZERO_INITIAL_GUESS) {
742: /* this is a hack fix, should have another version without the second BLASdotu */
743: VecSet(xx,zero);
744: }
745: VecGetArray(xx,&x);
746: VecGetArrayRead(bb,&b);
747: its = its*lits;
748: if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
749: while (its--) {
750: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
751: for (i=0; i<m; i++) {
752: PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
753: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
754: }
755: }
756: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
757: for (i=m-1; i>=0; i--) {
758: PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
759: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
760: }
761: }
762: }
763: VecRestoreArrayRead(bb,&b);
764: VecRestoreArray(xx,&x);
765: return(0);
766: }
768: /* -----------------------------------------------------------------*/
769: PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
770: {
771: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
772: const PetscScalar *v = mat->v,*x;
773: PetscScalar *y;
774: PetscErrorCode ierr;
775: PetscBLASInt m, n,_One=1;
776: PetscScalar _DOne=1.0,_DZero=0.0;
779: PetscBLASIntCast(A->rmap->n,&m);
780: PetscBLASIntCast(A->cmap->n,&n);
781: VecGetArrayRead(xx,&x);
782: VecGetArrayWrite(yy,&y);
783: if (!A->rmap->n || !A->cmap->n) {
784: PetscBLASInt i;
785: for (i=0; i<n; i++) y[i] = 0.0;
786: } else {
787: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One));
788: PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);
789: }
790: VecRestoreArrayRead(xx,&x);
791: VecRestoreArrayWrite(yy,&y);
792: return(0);
793: }
795: PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
796: {
797: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
798: PetscScalar *y,_DOne=1.0,_DZero=0.0;
799: PetscErrorCode ierr;
800: PetscBLASInt m, n, _One=1;
801: const PetscScalar *v = mat->v,*x;
804: PetscBLASIntCast(A->rmap->n,&m);
805: PetscBLASIntCast(A->cmap->n,&n);
806: VecGetArrayRead(xx,&x);
807: VecGetArrayWrite(yy,&y);
808: if (!A->rmap->n || !A->cmap->n) {
809: PetscBLASInt i;
810: for (i=0; i<m; i++) y[i] = 0.0;
811: } else {
812: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One));
813: PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);
814: }
815: VecRestoreArrayRead(xx,&x);
816: VecRestoreArrayWrite(yy,&y);
817: return(0);
818: }
820: PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
821: {
822: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
823: const PetscScalar *v = mat->v,*x;
824: PetscScalar *y,_DOne=1.0;
825: PetscErrorCode ierr;
826: PetscBLASInt m, n, _One=1;
829: PetscBLASIntCast(A->rmap->n,&m);
830: PetscBLASIntCast(A->cmap->n,&n);
831: if (!A->rmap->n || !A->cmap->n) return(0);
832: if (zz != yy) {VecCopy(zz,yy);}
833: VecGetArrayRead(xx,&x);
834: VecGetArray(yy,&y);
835: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
836: VecRestoreArrayRead(xx,&x);
837: VecRestoreArray(yy,&y);
838: PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
839: return(0);
840: }
842: PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
843: {
844: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
845: const PetscScalar *v = mat->v,*x;
846: PetscScalar *y;
847: PetscErrorCode ierr;
848: PetscBLASInt m, n, _One=1;
849: PetscScalar _DOne=1.0;
852: PetscBLASIntCast(A->rmap->n,&m);
853: PetscBLASIntCast(A->cmap->n,&n);
854: if (!A->rmap->n || !A->cmap->n) return(0);
855: if (zz != yy) {VecCopy(zz,yy);}
856: VecGetArrayRead(xx,&x);
857: VecGetArray(yy,&y);
858: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
859: VecRestoreArrayRead(xx,&x);
860: VecRestoreArray(yy,&y);
861: PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
862: return(0);
863: }
865: /* -----------------------------------------------------------------*/
866: static PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
867: {
868: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
870: PetscInt i;
873: *ncols = A->cmap->n;
874: if (cols) {
875: PetscMalloc1(A->cmap->n+1,cols);
876: for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
877: }
878: if (vals) {
879: const PetscScalar *v;
881: MatDenseGetArrayRead(A,&v);
882: PetscMalloc1(A->cmap->n+1,vals);
883: v += row;
884: for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
885: MatDenseRestoreArrayRead(A,&v);
886: }
887: return(0);
888: }
890: static PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
891: {
895: if (cols) {PetscFree(*cols);}
896: if (vals) {PetscFree(*vals); }
897: return(0);
898: }
899: /* ----------------------------------------------------------------*/
900: static PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
901: {
902: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
903: PetscScalar *av;
904: PetscInt i,j,idx=0;
905: #if defined(PETSC_HAVE_CUDA)
906: PetscOffloadMask oldf;
907: #endif
908: PetscErrorCode ierr;
911: MatDenseGetArray(A,&av);
912: if (!mat->roworiented) {
913: if (addv == INSERT_VALUES) {
914: for (j=0; j<n; j++) {
915: if (indexn[j] < 0) {idx += m; continue;}
916: #if defined(PETSC_USE_DEBUG)
917: if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
918: #endif
919: for (i=0; i<m; i++) {
920: if (indexm[i] < 0) {idx++; continue;}
921: #if defined(PETSC_USE_DEBUG)
922: if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
923: #endif
924: av[indexn[j]*mat->lda + indexm[i]] = v[idx++];
925: }
926: }
927: } else {
928: for (j=0; j<n; j++) {
929: if (indexn[j] < 0) {idx += m; continue;}
930: #if defined(PETSC_USE_DEBUG)
931: if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
932: #endif
933: for (i=0; i<m; i++) {
934: if (indexm[i] < 0) {idx++; continue;}
935: #if defined(PETSC_USE_DEBUG)
936: if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
937: #endif
938: av[indexn[j]*mat->lda + indexm[i]] += v[idx++];
939: }
940: }
941: }
942: } else {
943: if (addv == INSERT_VALUES) {
944: for (i=0; i<m; i++) {
945: if (indexm[i] < 0) { idx += n; continue;}
946: #if defined(PETSC_USE_DEBUG)
947: if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
948: #endif
949: for (j=0; j<n; j++) {
950: if (indexn[j] < 0) { idx++; continue;}
951: #if defined(PETSC_USE_DEBUG)
952: if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
953: #endif
954: av[indexn[j]*mat->lda + indexm[i]] = v[idx++];
955: }
956: }
957: } else {
958: for (i=0; i<m; i++) {
959: if (indexm[i] < 0) { idx += n; continue;}
960: #if defined(PETSC_USE_DEBUG)
961: if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
962: #endif
963: for (j=0; j<n; j++) {
964: if (indexn[j] < 0) { idx++; continue;}
965: #if defined(PETSC_USE_DEBUG)
966: if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
967: #endif
968: av[indexn[j]*mat->lda + indexm[i]] += v[idx++];
969: }
970: }
971: }
972: }
973: /* hack to prevent unneeded copy to the GPU while returning the array */
974: #if defined(PETSC_HAVE_CUDA)
975: oldf = A->offloadmask;
976: A->offloadmask = PETSC_OFFLOAD_GPU;
977: #endif
978: MatDenseRestoreArray(A,&av);
979: #if defined(PETSC_HAVE_CUDA)
980: A->offloadmask = (oldf == PETSC_OFFLOAD_UNALLOCATED ? PETSC_OFFLOAD_UNALLOCATED : PETSC_OFFLOAD_CPU);
981: #endif
982: return(0);
983: }
985: static PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
986: {
987: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
988: const PetscScalar *vv;
989: PetscInt i,j;
990: PetscErrorCode ierr;
993: MatDenseGetArrayRead(A,&vv);
994: /* row-oriented output */
995: for (i=0; i<m; i++) {
996: if (indexm[i] < 0) {v += n;continue;}
997: if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n);
998: for (j=0; j<n; j++) {
999: if (indexn[j] < 0) {v++; continue;}
1000: if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n);
1001: *v++ = vv[indexn[j]*mat->lda + indexm[i]];
1002: }
1003: }
1004: MatDenseRestoreArrayRead(A,&vv);
1005: return(0);
1006: }
1008: /* -----------------------------------------------------------------*/
1010: PetscErrorCode MatView_Dense_Binary(Mat mat,PetscViewer viewer)
1011: {
1012: PetscErrorCode ierr;
1013: PetscBool skipHeader;
1014: PetscViewerFormat format;
1015: PetscInt header[4],M,N,m,lda,i,j,k;
1016: const PetscScalar *v;
1017: PetscScalar *vwork;
1020: PetscViewerSetUp(viewer);
1021: PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
1022: PetscViewerGetFormat(viewer,&format);
1023: if (skipHeader) format = PETSC_VIEWER_NATIVE;
1025: MatGetSize(mat,&M,&N);
1027: /* write matrix header */
1028: header[0] = MAT_FILE_CLASSID; header[1] = M; header[2] = N;
1029: header[3] = (format == PETSC_VIEWER_NATIVE) ? MATRIX_BINARY_FORMAT_DENSE : M*N;
1030: if (!skipHeader) {PetscViewerBinaryWrite(viewer,header,4,PETSC_INT);}
1032: MatGetLocalSize(mat,&m,NULL);
1033: if (format != PETSC_VIEWER_NATIVE) {
1034: PetscInt nnz = m*N, *iwork;
1035: /* store row lengths for each row */
1036: PetscMalloc1(nnz,&iwork);
1037: for (i=0; i<m; i++) iwork[i] = N;
1038: PetscViewerBinaryWriteAll(viewer,iwork,m,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1039: /* store column indices (zero start index) */
1040: for (k=0, i=0; i<m; i++)
1041: for (j=0; j<N; j++, k++)
1042: iwork[k] = j;
1043: PetscViewerBinaryWriteAll(viewer,iwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1044: PetscFree(iwork);
1045: }
1046: /* store matrix values as a dense matrix in row major order */
1047: PetscMalloc1(m*N,&vwork);
1048: MatDenseGetArrayRead(mat,&v);
1049: MatDenseGetLDA(mat,&lda);
1050: for (k=0, i=0; i<m; i++)
1051: for (j=0; j<N; j++, k++)
1052: vwork[k] = v[i+lda*j];
1053: MatDenseRestoreArrayRead(mat,&v);
1054: PetscViewerBinaryWriteAll(viewer,vwork,m*N,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR);
1055: PetscFree(vwork);
1056: return(0);
1057: }
1059: PetscErrorCode MatLoad_Dense_Binary(Mat mat,PetscViewer viewer)
1060: {
1062: PetscBool skipHeader;
1063: PetscInt header[4],M,N,m,nz,lda,i,j,k;
1064: PetscInt rows,cols;
1065: PetscScalar *v,*vwork;
1068: PetscViewerSetUp(viewer);
1069: PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
1071: if (!skipHeader) {
1072: PetscViewerBinaryRead(viewer,header,4,NULL,PETSC_INT);
1073: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Not a matrix object in file");
1074: M = header[1]; N = header[2];
1075: if (M < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix row size (%D) in file is negative",M);
1076: if (N < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix column size (%D) in file is negative",N);
1077: nz = header[3];
1078: if (nz != MATRIX_BINARY_FORMAT_DENSE && nz < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Unknown matrix format %D in file",nz);
1079: } else {
1080: MatGetSize(mat,&M,&N);
1081: if (M < 0 || N < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Matrix binary file header was skipped, thus the user must specify the global sizes of input matrix");
1082: nz = MATRIX_BINARY_FORMAT_DENSE;
1083: }
1085: /* setup global sizes if not set */
1086: if (mat->rmap->N < 0) mat->rmap->N = M;
1087: if (mat->cmap->N < 0) mat->cmap->N = N;
1088: MatSetUp(mat);
1089: /* check if global sizes are correct */
1090: MatGetSize(mat,&rows,&cols);
1091: if (M != rows || N != cols) SETERRQ4(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different sizes (%d, %d) than the input matrix (%d, %d)",M,N,rows,cols);
1093: MatGetSize(mat,NULL,&N);
1094: MatGetLocalSize(mat,&m,NULL);
1095: MatDenseGetArray(mat,&v);
1096: MatDenseGetLDA(mat,&lda);
1097: if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense format */
1098: PetscInt nnz = m*N;
1099: /* read in matrix values */
1100: PetscMalloc1(nnz,&vwork);
1101: PetscViewerBinaryReadAll(viewer,vwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR);
1102: /* store values in column major order */
1103: for (j=0; j<N; j++)
1104: for (i=0; i<m; i++)
1105: v[i+lda*j] = vwork[i*N+j];
1106: PetscFree(vwork);
1107: } else { /* matrix in file is sparse format */
1108: PetscInt nnz = 0, *rlens, *icols;
1109: /* read in row lengths */
1110: PetscMalloc1(m,&rlens);
1111: PetscViewerBinaryReadAll(viewer,rlens,m,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1112: for (i=0; i<m; i++) nnz += rlens[i];
1113: /* read in column indices and values */
1114: PetscMalloc2(nnz,&icols,nnz,&vwork);
1115: PetscViewerBinaryReadAll(viewer,icols,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1116: PetscViewerBinaryReadAll(viewer,vwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR);
1117: /* store values in column major order */
1118: for (k=0, i=0; i<m; i++)
1119: for (j=0; j<rlens[i]; j++, k++)
1120: v[i+lda*icols[k]] = vwork[k];
1121: PetscFree(rlens);
1122: PetscFree2(icols,vwork);
1123: }
1124: MatDenseRestoreArray(mat,&v);
1125: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
1126: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
1127: return(0);
1128: }
1130: PetscErrorCode MatLoad_SeqDense(Mat newMat, PetscViewer viewer)
1131: {
1132: PetscBool isbinary, ishdf5;
1138: /* force binary viewer to load .info file if it has not yet done so */
1139: PetscViewerSetUp(viewer);
1140: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1141: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5, &ishdf5);
1142: if (isbinary) {
1143: MatLoad_Dense_Binary(newMat,viewer);
1144: } else if (ishdf5) {
1145: #if defined(PETSC_HAVE_HDF5)
1146: MatLoad_Dense_HDF5(newMat,viewer);
1147: #else
1148: SETERRQ(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1149: #endif
1150: } else {
1151: SETERRQ2(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newMat)->type_name);
1152: }
1153: return(0);
1154: }
1156: static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
1157: {
1158: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1159: PetscErrorCode ierr;
1160: PetscInt i,j;
1161: const char *name;
1162: PetscScalar *v,*av;
1163: PetscViewerFormat format;
1164: #if defined(PETSC_USE_COMPLEX)
1165: PetscBool allreal = PETSC_TRUE;
1166: #endif
1169: MatDenseGetArrayRead(A,(const PetscScalar**)&av);
1170: PetscViewerGetFormat(viewer,&format);
1171: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1172: return(0); /* do nothing for now */
1173: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
1174: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1175: for (i=0; i<A->rmap->n; i++) {
1176: v = av + i;
1177: PetscViewerASCIIPrintf(viewer,"row %D:",i);
1178: for (j=0; j<A->cmap->n; j++) {
1179: #if defined(PETSC_USE_COMPLEX)
1180: if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
1181: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1182: } else if (PetscRealPart(*v)) {
1183: PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)PetscRealPart(*v));
1184: }
1185: #else
1186: if (*v) {
1187: PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)*v);
1188: }
1189: #endif
1190: v += a->lda;
1191: }
1192: PetscViewerASCIIPrintf(viewer,"\n");
1193: }
1194: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1195: } else {
1196: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1197: #if defined(PETSC_USE_COMPLEX)
1198: /* determine if matrix has all real values */
1199: v = av;
1200: for (i=0; i<A->rmap->n*A->cmap->n; i++) {
1201: if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;}
1202: }
1203: #endif
1204: if (format == PETSC_VIEWER_ASCII_MATLAB) {
1205: PetscObjectGetName((PetscObject)A,&name);
1206: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);
1207: PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);
1208: PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
1209: }
1211: for (i=0; i<A->rmap->n; i++) {
1212: v = av + i;
1213: for (j=0; j<A->cmap->n; j++) {
1214: #if defined(PETSC_USE_COMPLEX)
1215: if (allreal) {
1216: PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));
1217: } else {
1218: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1219: }
1220: #else
1221: PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);
1222: #endif
1223: v += a->lda;
1224: }
1225: PetscViewerASCIIPrintf(viewer,"\n");
1226: }
1227: if (format == PETSC_VIEWER_ASCII_MATLAB) {
1228: PetscViewerASCIIPrintf(viewer,"];\n");
1229: }
1230: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1231: }
1232: MatDenseRestoreArrayRead(A,(const PetscScalar**)&av);
1233: PetscViewerFlush(viewer);
1234: return(0);
1235: }
1237: static PetscErrorCode MatView_SeqDense_Binary(Mat mat,PetscViewer viewer)
1238: {
1239: PetscErrorCode ierr;
1240: PetscViewerFormat format;
1241: PetscInt header[4],M,N,m,lda,i,j,k;
1242: const PetscScalar *v;
1243: PetscScalar *vwork;
1246: PetscViewerSetUp(viewer);
1248: PetscViewerGetFormat(viewer,&format);
1249: MatGetSize(mat,&M,&N);
1251: /* write matrix header */
1252: header[0] = MAT_FILE_CLASSID; header[1] = M; header[2] = N;
1253: header[3] = (format == PETSC_VIEWER_NATIVE) ? MATRIX_BINARY_FORMAT_DENSE : M*N;
1254: PetscViewerBinaryWrite(viewer,header,4,PETSC_INT);
1256: MatGetLocalSize(mat,&m,NULL);
1257: if (format != PETSC_VIEWER_NATIVE) {
1258: PetscInt nnz = m*N, *iwork;
1259: /* store row lengths for each row */
1260: PetscMalloc1(nnz,&iwork);
1261: for (i=0; i<m; i++) iwork[i] = N;
1262: PetscViewerBinaryWrite(viewer,iwork,m,PETSC_INT);
1263: /* store column indices (zero start index) */
1264: for (k=0, i=0; i<m; i++)
1265: for (j=0; j<N; j++, k++)
1266: iwork[k] = j;
1267: PetscViewerBinaryWrite(viewer,iwork,nnz,PETSC_INT);
1268: PetscFree(iwork);
1269: }
1270: /* store the matrix values as a dense matrix in row major order */
1271: PetscMalloc1(m*N,&vwork);
1272: MatDenseGetArrayRead(mat,&v);
1273: MatDenseGetLDA(mat,&lda);
1274: for (k=0, i=0; i<m; i++)
1275: for (j=0; j<N; j++, k++)
1276: vwork[k] = v[i+lda*j];
1277: MatDenseRestoreArrayRead(mat,&v);
1278: PetscViewerBinaryWrite(viewer,vwork,m*N,PETSC_SCALAR);
1279: PetscFree(vwork);
1280: return(0);
1281: }
1283: #include <petscdraw.h>
1284: static PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1285: {
1286: Mat A = (Mat) Aa;
1287: PetscErrorCode ierr;
1288: PetscInt m = A->rmap->n,n = A->cmap->n,i,j;
1289: int color = PETSC_DRAW_WHITE;
1290: const PetscScalar *v;
1291: PetscViewer viewer;
1292: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r;
1293: PetscViewerFormat format;
1296: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
1297: PetscViewerGetFormat(viewer,&format);
1298: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
1300: /* Loop over matrix elements drawing boxes */
1301: MatDenseGetArrayRead(A,&v);
1302: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1303: PetscDrawCollectiveBegin(draw);
1304: /* Blue for negative and Red for positive */
1305: for (j = 0; j < n; j++) {
1306: x_l = j; x_r = x_l + 1.0;
1307: for (i = 0; i < m; i++) {
1308: y_l = m - i - 1.0;
1309: y_r = y_l + 1.0;
1310: if (PetscRealPart(v[j*m+i]) > 0.) color = PETSC_DRAW_RED;
1311: else if (PetscRealPart(v[j*m+i]) < 0.) color = PETSC_DRAW_BLUE;
1312: else continue;
1313: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1314: }
1315: }
1316: PetscDrawCollectiveEnd(draw);
1317: } else {
1318: /* use contour shading to indicate magnitude of values */
1319: /* first determine max of all nonzero values */
1320: PetscReal minv = 0.0, maxv = 0.0;
1321: PetscDraw popup;
1323: for (i=0; i < m*n; i++) {
1324: if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1325: }
1326: if (minv >= maxv) maxv = minv + PETSC_SMALL;
1327: PetscDrawGetPopup(draw,&popup);
1328: PetscDrawScalePopup(popup,minv,maxv);
1330: PetscDrawCollectiveBegin(draw);
1331: for (j=0; j<n; j++) {
1332: x_l = j;
1333: x_r = x_l + 1.0;
1334: for (i=0; i<m; i++) {
1335: y_l = m - i - 1.0;
1336: y_r = y_l + 1.0;
1337: color = PetscDrawRealToColor(PetscAbsScalar(v[j*m+i]),minv,maxv);
1338: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1339: }
1340: }
1341: PetscDrawCollectiveEnd(draw);
1342: }
1343: MatDenseRestoreArrayRead(A,&v);
1344: return(0);
1345: }
1347: static PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1348: {
1349: PetscDraw draw;
1350: PetscBool isnull;
1351: PetscReal xr,yr,xl,yl,h,w;
1355: PetscViewerDrawGetDraw(viewer,0,&draw);
1356: PetscDrawIsNull(draw,&isnull);
1357: if (isnull) return(0);
1359: xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1360: xr += w; yr += h; xl = -w; yl = -h;
1361: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
1362: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
1363: PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
1364: PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
1365: PetscDrawSave(draw);
1366: return(0);
1367: }
1369: PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1370: {
1372: PetscBool iascii,isbinary,isdraw;
1375: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1376: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1377: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1379: if (iascii) {
1380: MatView_SeqDense_ASCII(A,viewer);
1381: } else if (isbinary) {
1382: MatView_SeqDense_Binary(A,viewer);
1383: } else if (isdraw) {
1384: MatView_SeqDense_Draw(A,viewer);
1385: }
1386: return(0);
1387: }
1389: static PetscErrorCode MatDensePlaceArray_SeqDense(Mat A,const PetscScalar array[])
1390: {
1391: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1394: a->unplacedarray = a->v;
1395: a->unplaced_user_alloc = a->user_alloc;
1396: a->v = (PetscScalar*) array;
1397: #if defined(PETSC_HAVE_CUDA)
1398: A->offloadmask = PETSC_OFFLOAD_CPU;
1399: #endif
1400: return(0);
1401: }
1403: static PetscErrorCode MatDenseResetArray_SeqDense(Mat A)
1404: {
1405: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1408: a->v = a->unplacedarray;
1409: a->user_alloc = a->unplaced_user_alloc;
1410: a->unplacedarray = NULL;
1411: #if defined(PETSC_HAVE_CUDA)
1412: A->offloadmask = PETSC_OFFLOAD_CPU;
1413: #endif
1414: return(0);
1415: }
1417: PetscErrorCode MatDestroy_SeqDense(Mat mat)
1418: {
1419: Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
1423: #if defined(PETSC_USE_LOG)
1424: PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1425: #endif
1426: PetscFree(l->pivots);
1427: PetscFree(l->fwork);
1428: MatDestroy(&l->ptapwork);
1429: if (!l->user_alloc) {PetscFree(l->v);}
1430: PetscFree(mat->data);
1432: PetscObjectChangeTypeName((PetscObject)mat,0);
1433: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL);
1434: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);
1435: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);
1436: PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);
1437: PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);
1438: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",NULL);
1439: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",NULL);
1440: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqaij_C",NULL);
1441: #if defined(PETSC_HAVE_ELEMENTAL)
1442: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_elemental_C",NULL);
1443: #endif
1444: #if defined(PETSC_HAVE_CUDA)
1445: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqdensecuda_C",NULL);
1446: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdensecuda_seqdensecuda_C",NULL);
1447: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdensecuda_seqdense_C",NULL);
1448: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqaijcusparse_seqdense_C",NULL);
1449: #endif
1450: PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL);
1451: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqaij_seqdense_C",NULL);
1452: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdense_seqdense_C",NULL);
1453: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C",NULL);
1454: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C",NULL);
1455: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqbaij_seqdense_C",NULL);
1456: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqbaij_seqdense_C",NULL);
1457: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqbaij_seqdense_C",NULL);
1458: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqsbaij_seqdense_C",NULL);
1459: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqsbaij_seqdense_C",NULL);
1460: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqsbaij_seqdense_C",NULL);
1461: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaijperm_seqdense_C",NULL);
1462: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaijperm_seqdense_C",NULL);
1463: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaijsell_seqdense_C",NULL);
1464: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaijsell_seqdense_C",NULL);
1465: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaijmkl_seqdense_C",NULL);
1466: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaijmkl_seqdense_C",NULL);
1467: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_nest_seqdense_C",NULL);
1468: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_nest_seqdense_C",NULL);
1470: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",NULL);
1471: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaij_seqdense_C",NULL);
1472: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaijperm_seqdense_C",NULL);
1473: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaijperm_seqdense_C",NULL);
1474: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaijsell_seqdense_C",NULL);
1475: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaijsell_seqdense_C",NULL);
1476: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaijmkl_seqdense_C",NULL);
1477: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaijmkl_seqdense_C",NULL);
1478: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);
1479: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);
1480: return(0);
1481: }
1483: static PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1484: {
1485: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1487: PetscInt k,j,m,n,M;
1488: PetscScalar *v,tmp;
1491: m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1492: if (reuse == MAT_INPLACE_MATRIX && m == n) { /* in place transpose */
1493: MatDenseGetArray(A,&v);
1494: for (j=0; j<m; j++) {
1495: for (k=0; k<j; k++) {
1496: tmp = v[j + k*M];
1497: v[j + k*M] = v[k + j*M];
1498: v[k + j*M] = tmp;
1499: }
1500: }
1501: MatDenseRestoreArray(A,&v);
1502: } else { /* out-of-place transpose */
1503: Mat tmat;
1504: Mat_SeqDense *tmatd;
1505: PetscScalar *v2;
1506: PetscInt M2;
1508: if (reuse != MAT_REUSE_MATRIX) {
1509: MatCreate(PetscObjectComm((PetscObject)A),&tmat);
1510: MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);
1511: MatSetType(tmat,((PetscObject)A)->type_name);
1512: MatSeqDenseSetPreallocation(tmat,NULL);
1513: } else tmat = *matout;
1515: MatDenseGetArrayRead(A,(const PetscScalar**)&v);
1516: MatDenseGetArray(tmat,&v2);
1517: tmatd = (Mat_SeqDense*)tmat->data;
1518: M2 = tmatd->lda;
1519: for (j=0; j<n; j++) {
1520: for (k=0; k<m; k++) v2[j + k*M2] = v[k + j*M];
1521: }
1522: MatDenseRestoreArray(tmat,&v2);
1523: MatDenseRestoreArrayRead(A,(const PetscScalar**)&v);
1524: MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
1525: MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
1526: if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) *matout = tmat;
1527: else {
1528: MatHeaderMerge(A,&tmat);
1529: }
1530: }
1531: return(0);
1532: }
1534: static PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg)
1535: {
1536: Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1537: Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1538: PetscInt i;
1539: const PetscScalar *v1,*v2;
1540: PetscErrorCode ierr;
1543: if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; return(0);}
1544: if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; return(0);}
1545: MatDenseGetArrayRead(A1,&v1);
1546: MatDenseGetArrayRead(A2,&v2);
1547: for (i=0; i<A1->cmap->n; i++) {
1548: PetscArraycmp(v1,v2,A1->rmap->n,flg);
1549: if (*flg == PETSC_FALSE) return(0);
1550: v1 += mat1->lda;
1551: v2 += mat2->lda;
1552: }
1553: MatDenseRestoreArrayRead(A1,&v1);
1554: MatDenseRestoreArrayRead(A2,&v2);
1555: *flg = PETSC_TRUE;
1556: return(0);
1557: }
1559: static PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1560: {
1561: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1562: PetscInt i,n,len;
1563: PetscScalar *x;
1564: const PetscScalar *vv;
1565: PetscErrorCode ierr;
1568: VecGetSize(v,&n);
1569: VecGetArray(v,&x);
1570: len = PetscMin(A->rmap->n,A->cmap->n);
1571: MatDenseGetArrayRead(A,&vv);
1572: if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1573: for (i=0; i<len; i++) {
1574: x[i] = vv[i*mat->lda + i];
1575: }
1576: MatDenseRestoreArrayRead(A,&vv);
1577: VecRestoreArray(v,&x);
1578: return(0);
1579: }
1581: static PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1582: {
1583: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1584: const PetscScalar *l,*r;
1585: PetscScalar x,*v,*vv;
1586: PetscErrorCode ierr;
1587: PetscInt i,j,m = A->rmap->n,n = A->cmap->n;
1590: MatDenseGetArray(A,&vv);
1591: if (ll) {
1592: VecGetSize(ll,&m);
1593: VecGetArrayRead(ll,&l);
1594: if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1595: for (i=0; i<m; i++) {
1596: x = l[i];
1597: v = vv + i;
1598: for (j=0; j<n; j++) { (*v) *= x; v+= mat->lda;}
1599: }
1600: VecRestoreArrayRead(ll,&l);
1601: PetscLogFlops(1.0*n*m);
1602: }
1603: if (rr) {
1604: VecGetSize(rr,&n);
1605: VecGetArrayRead(rr,&r);
1606: if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1607: for (i=0; i<n; i++) {
1608: x = r[i];
1609: v = vv + i*mat->lda;
1610: for (j=0; j<m; j++) (*v++) *= x;
1611: }
1612: VecRestoreArrayRead(rr,&r);
1613: PetscLogFlops(1.0*n*m);
1614: }
1615: MatDenseRestoreArray(A,&vv);
1616: return(0);
1617: }
1619: PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1620: {
1621: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1622: PetscScalar *v,*vv;
1623: PetscReal sum = 0.0;
1624: PetscInt lda =mat->lda,m=A->rmap->n,i,j;
1625: PetscErrorCode ierr;
1628: MatDenseGetArrayRead(A,(const PetscScalar**)&vv);
1629: v = vv;
1630: if (type == NORM_FROBENIUS) {
1631: if (lda>m) {
1632: for (j=0; j<A->cmap->n; j++) {
1633: v = vv+j*lda;
1634: for (i=0; i<m; i++) {
1635: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1636: }
1637: }
1638: } else {
1639: #if defined(PETSC_USE_REAL___FP16)
1640: PetscBLASInt one = 1,cnt = A->cmap->n*A->rmap->n;
1641: *nrm = BLASnrm2_(&cnt,v,&one);
1642: }
1643: #else
1644: for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1645: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1646: }
1647: }
1648: *nrm = PetscSqrtReal(sum);
1649: #endif
1650: PetscLogFlops(2.0*A->cmap->n*A->rmap->n);
1651: } else if (type == NORM_1) {
1652: *nrm = 0.0;
1653: for (j=0; j<A->cmap->n; j++) {
1654: v = vv + j*mat->lda;
1655: sum = 0.0;
1656: for (i=0; i<A->rmap->n; i++) {
1657: sum += PetscAbsScalar(*v); v++;
1658: }
1659: if (sum > *nrm) *nrm = sum;
1660: }
1661: PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1662: } else if (type == NORM_INFINITY) {
1663: *nrm = 0.0;
1664: for (j=0; j<A->rmap->n; j++) {
1665: v = vv + j;
1666: sum = 0.0;
1667: for (i=0; i<A->cmap->n; i++) {
1668: sum += PetscAbsScalar(*v); v += mat->lda;
1669: }
1670: if (sum > *nrm) *nrm = sum;
1671: }
1672: PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1673: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
1674: MatDenseRestoreArrayRead(A,(const PetscScalar**)&vv);
1675: return(0);
1676: }
1678: static PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg)
1679: {
1680: Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
1684: switch (op) {
1685: case MAT_ROW_ORIENTED:
1686: aij->roworiented = flg;
1687: break;
1688: case MAT_NEW_NONZERO_LOCATIONS:
1689: case MAT_NEW_NONZERO_LOCATION_ERR:
1690: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1691: case MAT_NEW_DIAGONALS:
1692: case MAT_KEEP_NONZERO_PATTERN:
1693: case MAT_IGNORE_OFF_PROC_ENTRIES:
1694: case MAT_USE_HASH_TABLE:
1695: case MAT_IGNORE_ZERO_ENTRIES:
1696: case MAT_IGNORE_LOWER_TRIANGULAR:
1697: case MAT_SORTED_FULL:
1698: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1699: break;
1700: case MAT_SPD:
1701: case MAT_SYMMETRIC:
1702: case MAT_STRUCTURALLY_SYMMETRIC:
1703: case MAT_HERMITIAN:
1704: case MAT_SYMMETRY_ETERNAL:
1705: /* These options are handled directly by MatSetOption() */
1706: break;
1707: default:
1708: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
1709: }
1710: return(0);
1711: }
1713: static PetscErrorCode MatZeroEntries_SeqDense(Mat A)
1714: {
1715: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1717: PetscInt lda=l->lda,m=A->rmap->n,j;
1718: PetscScalar *v;
1721: MatDenseGetArray(A,&v);
1722: if (lda>m) {
1723: for (j=0; j<A->cmap->n; j++) {
1724: PetscArrayzero(v+j*lda,m);
1725: }
1726: } else {
1727: PetscArrayzero(v,A->rmap->n*A->cmap->n);
1728: }
1729: MatDenseRestoreArray(A,&v);
1730: return(0);
1731: }
1733: static PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1734: {
1735: PetscErrorCode ierr;
1736: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1737: PetscInt m = l->lda, n = A->cmap->n, i,j;
1738: PetscScalar *slot,*bb,*v;
1739: const PetscScalar *xx;
1742: #if defined(PETSC_USE_DEBUG)
1743: for (i=0; i<N; i++) {
1744: if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
1745: if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n);
1746: }
1747: #endif
1748: if (!N) return(0);
1750: /* fix right hand side if needed */
1751: if (x && b) {
1752: VecGetArrayRead(x,&xx);
1753: VecGetArray(b,&bb);
1754: for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
1755: VecRestoreArrayRead(x,&xx);
1756: VecRestoreArray(b,&bb);
1757: }
1759: MatDenseGetArray(A,&v);
1760: for (i=0; i<N; i++) {
1761: slot = v + rows[i];
1762: for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
1763: }
1764: if (diag != 0.0) {
1765: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
1766: for (i=0; i<N; i++) {
1767: slot = v + (m+1)*rows[i];
1768: *slot = diag;
1769: }
1770: }
1771: MatDenseRestoreArray(A,&v);
1772: return(0);
1773: }
1775: static PetscErrorCode MatDenseGetLDA_SeqDense(Mat A,PetscInt *lda)
1776: {
1777: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1780: *lda = mat->lda;
1781: return(0);
1782: }
1784: PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[])
1785: {
1786: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1789: *array = mat->v;
1790: return(0);
1791: }
1793: PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1794: {
1796: return(0);
1797: }
1799: /*@C
1800: MatDenseGetLDA - gets the leading dimension of the array returned from MatDenseGetArray()
1802: Logically Collective on Mat
1804: Input Parameter:
1805: . mat - a MATSEQDENSE or MATMPIDENSE matrix
1807: Output Parameter:
1808: . lda - the leading dimension
1810: Level: intermediate
1812: .seealso: MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatSeqDenseSetLDA()
1813: @*/
1814: PetscErrorCode MatDenseGetLDA(Mat A,PetscInt *lda)
1815: {
1819: PetscUseMethod(A,"MatDenseGetLDA_C",(Mat,PetscInt*),(A,lda));
1820: return(0);
1821: }
1823: /*@C
1824: MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored
1826: Logically Collective on Mat
1828: Input Parameter:
1829: . mat - a MATSEQDENSE or MATMPIDENSE matrix
1831: Output Parameter:
1832: . array - pointer to the data
1834: Level: intermediate
1836: .seealso: MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
1837: @*/
1838: PetscErrorCode MatDenseGetArray(Mat A,PetscScalar **array)
1839: {
1843: PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));
1844: return(0);
1845: }
1847: /*@C
1848: MatDenseRestoreArray - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()
1850: Logically Collective on Mat
1852: Input Parameters:
1853: + mat - a MATSEQDENSE or MATMPIDENSE matrix
1854: - array - pointer to the data
1856: Level: intermediate
1858: .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
1859: @*/
1860: PetscErrorCode MatDenseRestoreArray(Mat A,PetscScalar **array)
1861: {
1865: PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));
1866: if (array) *array = NULL;
1867: PetscObjectStateIncrease((PetscObject)A);
1868: return(0);
1869: }
1871: /*@C
1872: MatDenseGetArrayRead - gives access to the array where the data for a SeqDense matrix is stored
1874: Not Collective
1876: Input Parameter:
1877: . mat - a MATSEQDENSE or MATMPIDENSE matrix
1879: Output Parameter:
1880: . array - pointer to the data
1882: Level: intermediate
1884: .seealso: MatDenseRestoreArray(), MatDenseGetArray(), MatDenseRestoreArrayRead()
1885: @*/
1886: PetscErrorCode MatDenseGetArrayRead(Mat A,const PetscScalar **array)
1887: {
1891: PetscUseMethod(A,"MatDenseGetArrayRead_C",(Mat,const PetscScalar**),(A,array));
1892: return(0);
1893: }
1895: /*@C
1896: MatDenseRestoreArrayRead - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()
1898: Not Collective
1900: Input Parameters:
1901: + mat - a MATSEQDENSE or MATMPIDENSE matrix
1902: - array - pointer to the data
1904: Level: intermediate
1906: .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArray()
1907: @*/
1908: PetscErrorCode MatDenseRestoreArrayRead(Mat A,const PetscScalar **array)
1909: {
1913: PetscUseMethod(A,"MatDenseRestoreArrayRead_C",(Mat,const PetscScalar**),(A,array));
1914: if (array) *array = NULL;
1915: return(0);
1916: }
1918: static PetscErrorCode MatCreateSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
1919: {
1920: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1922: PetscInt i,j,nrows,ncols,blda;
1923: const PetscInt *irow,*icol;
1924: PetscScalar *av,*bv,*v = mat->v;
1925: Mat newmat;
1928: ISGetIndices(isrow,&irow);
1929: ISGetIndices(iscol,&icol);
1930: ISGetLocalSize(isrow,&nrows);
1931: ISGetLocalSize(iscol,&ncols);
1933: /* Check submatrixcall */
1934: if (scall == MAT_REUSE_MATRIX) {
1935: PetscInt n_cols,n_rows;
1936: MatGetSize(*B,&n_rows,&n_cols);
1937: if (n_rows != nrows || n_cols != ncols) {
1938: /* resize the result matrix to match number of requested rows/columns */
1939: MatSetSizes(*B,nrows,ncols,nrows,ncols);
1940: }
1941: newmat = *B;
1942: } else {
1943: /* Create and fill new matrix */
1944: MatCreate(PetscObjectComm((PetscObject)A),&newmat);
1945: MatSetSizes(newmat,nrows,ncols,nrows,ncols);
1946: MatSetType(newmat,((PetscObject)A)->type_name);
1947: MatSeqDenseSetPreallocation(newmat,NULL);
1948: }
1950: /* Now extract the data pointers and do the copy,column at a time */
1951: MatDenseGetArray(newmat,&bv);
1952: MatDenseGetLDA(newmat,&blda);
1953: for (i=0; i<ncols; i++) {
1954: av = v + mat->lda*icol[i];
1955: for (j=0; j<nrows; j++) bv[j] = av[irow[j]];
1956: bv += blda;
1957: }
1958: MatDenseRestoreArray(newmat,&bv);
1960: /* Assemble the matrices so that the correct flags are set */
1961: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1962: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1964: /* Free work space */
1965: ISRestoreIndices(isrow,&irow);
1966: ISRestoreIndices(iscol,&icol);
1967: *B = newmat;
1968: return(0);
1969: }
1971: static PetscErrorCode MatCreateSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1972: {
1974: PetscInt i;
1977: if (scall == MAT_INITIAL_MATRIX) {
1978: PetscCalloc1(n+1,B);
1979: }
1981: for (i=0; i<n; i++) {
1982: MatCreateSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1983: }
1984: return(0);
1985: }
1987: static PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1988: {
1990: return(0);
1991: }
1993: static PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1994: {
1996: return(0);
1997: }
1999: static PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
2000: {
2001: Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data;
2002: PetscErrorCode ierr;
2003: const PetscScalar *va;
2004: PetscScalar *vb;
2005: PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
2008: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
2009: if (A->ops->copy != B->ops->copy) {
2010: MatCopy_Basic(A,B,str);
2011: return(0);
2012: }
2013: if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
2014: MatDenseGetArrayRead(A,&va);
2015: MatDenseGetArray(B,&vb);
2016: if (lda1>m || lda2>m) {
2017: for (j=0; j<n; j++) {
2018: PetscArraycpy(vb+j*lda2,va+j*lda1,m);
2019: }
2020: } else {
2021: PetscArraycpy(vb,va,A->rmap->n*A->cmap->n);
2022: }
2023: MatDenseRestoreArray(B,&vb);
2024: MatDenseRestoreArrayRead(A,&va);
2025: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2026: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2027: return(0);
2028: }
2030: static PetscErrorCode MatSetUp_SeqDense(Mat A)
2031: {
2035: MatSeqDenseSetPreallocation(A,0);
2036: return(0);
2037: }
2039: static PetscErrorCode MatConjugate_SeqDense(Mat A)
2040: {
2041: PetscInt i,nz = A->rmap->n*A->cmap->n;
2042: PetscScalar *aa;
2046: MatDenseGetArray(A,&aa);
2047: for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
2048: MatDenseRestoreArray(A,&aa);
2049: return(0);
2050: }
2052: static PetscErrorCode MatRealPart_SeqDense(Mat A)
2053: {
2054: PetscInt i,nz = A->rmap->n*A->cmap->n;
2055: PetscScalar *aa;
2059: MatDenseGetArray(A,&aa);
2060: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
2061: MatDenseRestoreArray(A,&aa);
2062: return(0);
2063: }
2065: static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
2066: {
2067: PetscInt i,nz = A->rmap->n*A->cmap->n;
2068: PetscScalar *aa;
2072: MatDenseGetArray(A,&aa);
2073: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
2074: MatDenseRestoreArray(A,&aa);
2075: return(0);
2076: }
2078: /* ----------------------------------------------------------------*/
2079: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
2080: {
2082: PetscInt m=A->rmap->n,n=B->cmap->n;
2083: PetscBool flg;
2086: MatSetSizes(C,m,n,m,n);
2087: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2088: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
2089: MatSeqDenseSetPreallocation(C,NULL);
2090: return(0);
2091: }
2093: PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2094: {
2095: Mat_SeqDense *a,*b=(Mat_SeqDense*)B->data,*c=(Mat_SeqDense*)C->data;
2096: PetscBLASInt m,n,k;
2097: const PetscScalar *av,*bv;
2098: PetscScalar *cv;
2099: PetscScalar _DOne=1.0,_DZero=0.0;
2100: PetscBool flg;
2101: PetscErrorCode (*numeric)(Mat,Mat,Mat)=NULL;
2102: PetscErrorCode ierr;
2105: /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
2106: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
2107: if (flg) numeric = MatMatMultNumeric_SeqAIJ_SeqDense;
2108: PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&flg);
2109: if (flg) numeric = MatMatMultNumeric_SeqBAIJ_SeqDense;
2110: PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&flg);
2111: if (flg) numeric = MatMatMultNumeric_SeqSBAIJ_SeqDense;
2112: PetscObjectTypeCompare((PetscObject)A,MATNEST,&flg);
2113: if (flg) numeric = MatMatMultNumeric_Nest_Dense;
2114: if (numeric) {
2115: C->ops->matmultnumeric = numeric;
2116: (*numeric)(A,B,C);
2117: return(0);
2118: }
2119: a = (Mat_SeqDense*)A->data;
2120: PetscBLASIntCast(C->rmap->n,&m);
2121: PetscBLASIntCast(C->cmap->n,&n);
2122: PetscBLASIntCast(A->cmap->n,&k);
2123: if (!m || !n || !k) return(0);
2124: MatDenseGetArrayRead(A,&av);
2125: MatDenseGetArrayRead(B,&bv);
2126: MatDenseGetArray(C,&cv);
2127: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda));
2128: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2129: MatDenseRestoreArrayRead(A,&av);
2130: MatDenseRestoreArrayRead(B,&bv);
2131: MatDenseRestoreArray(C,&cv);
2132: return(0);
2133: }
2135: PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
2136: {
2138: PetscInt m=A->rmap->n,n=B->rmap->n;
2139: PetscBool flg;
2142: MatSetSizes(C,m,n,m,n);
2143: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2144: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
2145: MatSeqDenseSetPreallocation(C,NULL);
2146: C->ops->mattransposemultnumeric = MatMatTransposeMultNumeric_SeqDense_SeqDense;
2147: return(0);
2148: }
2150: PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2151: {
2152: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2153: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2154: Mat_SeqDense *c = (Mat_SeqDense*)C->data;
2155: PetscBLASInt m,n,k;
2156: PetscScalar _DOne=1.0,_DZero=0.0;
2160: PetscBLASIntCast(C->rmap->n,&m);
2161: PetscBLASIntCast(C->cmap->n,&n);
2162: PetscBLASIntCast(A->cmap->n,&k);
2163: if (!m || !n || !k) return(0);
2164: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2165: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2166: return(0);
2167: }
2169: PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
2170: {
2172: PetscInt m=A->cmap->n,n=B->cmap->n;
2173: PetscBool flg;
2176: MatSetSizes(C,m,n,m,n);
2177: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2178: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
2179: MatSeqDenseSetPreallocation(C,NULL);
2180: return(0);
2181: }
2183: PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2184: {
2185: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2186: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2187: Mat_SeqDense *c = (Mat_SeqDense*)C->data;
2188: PetscBLASInt m,n,k;
2189: PetscScalar _DOne=1.0,_DZero=0.0;
2193: PetscBLASIntCast(C->rmap->n,&m);
2194: PetscBLASIntCast(C->cmap->n,&n);
2195: PetscBLASIntCast(A->rmap->n,&k);
2196: if (!m || !n || !k) return(0);
2197: PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2198: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2199: return(0);
2200: }
2202: /* ----------------------------------------------- */
2203: static PetscErrorCode MatProductSetFromOptions_SeqDense_AB(Mat C)
2204: {
2206: C->ops->matmultsymbolic = MatMatMultSymbolic_SeqDense_SeqDense;
2207: C->ops->productsymbolic = MatProductSymbolic_AB;
2208: /* dense mat may not call MatProductSymbolic(), thus set C->ops->productnumeric here */
2209: C->ops->productnumeric = MatProductNumeric_AB;
2210: return(0);
2211: }
2213: static PetscErrorCode MatProductSetFromOptions_SeqDense_AtB(Mat C)
2214: {
2216: C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_SeqDense_SeqDense;
2217: C->ops->productsymbolic = MatProductSymbolic_AtB;
2218: C->ops->productnumeric = MatProductNumeric_AtB;
2219: return(0);
2220: }
2222: static PetscErrorCode MatProductSetFromOptions_SeqDense_ABt(Mat C)
2223: {
2225: C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_SeqDense_SeqDense;
2226: C->ops->productsymbolic = MatProductSymbolic_ABt;
2227: C->ops->productnumeric = MatProductNumeric_ABt;
2228: return(0);
2229: }
2231: static PetscErrorCode MatProductSetFromOptions_SeqDense_PtAP(Mat C)
2232: {
2234: C->ops->productsymbolic = MatProductSymbolic_Basic;
2235: return(0);
2236: }
2238: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat C)
2239: {
2241: Mat_Product *product = C->product;
2244: switch (product->type) {
2245: case MATPRODUCT_AB:
2246: MatProductSetFromOptions_SeqDense_AB(C);
2247: break;
2248: case MATPRODUCT_AtB:
2249: MatProductSetFromOptions_SeqDense_AtB(C);
2250: break;
2251: case MATPRODUCT_ABt:
2252: MatProductSetFromOptions_SeqDense_ABt(C);
2253: break;
2254: case MATPRODUCT_PtAP:
2255: MatProductSetFromOptions_SeqDense_PtAP(C);
2256: break;
2257: default: SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"MatProduct type is not supported");
2258: }
2259: return(0);
2260: }
2261: /* ----------------------------------------------- */
2263: static PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
2264: {
2265: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2266: PetscErrorCode ierr;
2267: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2268: PetscScalar *x;
2269: const PetscScalar *aa;
2272: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2273: VecGetArray(v,&x);
2274: VecGetLocalSize(v,&p);
2275: MatDenseGetArrayRead(A,&aa);
2276: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2277: for (i=0; i<m; i++) {
2278: x[i] = aa[i]; if (idx) idx[i] = 0;
2279: for (j=1; j<n; j++) {
2280: if (PetscRealPart(x[i]) < PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;}
2281: }
2282: }
2283: MatDenseRestoreArrayRead(A,&aa);
2284: VecRestoreArray(v,&x);
2285: return(0);
2286: }
2288: static PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
2289: {
2290: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2291: PetscErrorCode ierr;
2292: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2293: PetscScalar *x;
2294: PetscReal atmp;
2295: const PetscScalar *aa;
2298: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2299: VecGetArray(v,&x);
2300: VecGetLocalSize(v,&p);
2301: MatDenseGetArrayRead(A,&aa);
2302: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2303: for (i=0; i<m; i++) {
2304: x[i] = PetscAbsScalar(aa[i]);
2305: for (j=1; j<n; j++) {
2306: atmp = PetscAbsScalar(aa[i+a->lda*j]);
2307: if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
2308: }
2309: }
2310: MatDenseRestoreArrayRead(A,&aa);
2311: VecRestoreArray(v,&x);
2312: return(0);
2313: }
2315: static PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
2316: {
2317: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2318: PetscErrorCode ierr;
2319: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2320: PetscScalar *x;
2321: const PetscScalar *aa;
2324: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2325: MatDenseGetArrayRead(A,&aa);
2326: VecGetArray(v,&x);
2327: VecGetLocalSize(v,&p);
2328: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2329: for (i=0; i<m; i++) {
2330: x[i] = aa[i]; if (idx) idx[i] = 0;
2331: for (j=1; j<n; j++) {
2332: if (PetscRealPart(x[i]) > PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;}
2333: }
2334: }
2335: VecRestoreArray(v,&x);
2336: MatDenseRestoreArrayRead(A,&aa);
2337: return(0);
2338: }
2340: static PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
2341: {
2342: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2343: PetscErrorCode ierr;
2344: PetscScalar *x;
2345: const PetscScalar *aa;
2348: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2349: MatDenseGetArrayRead(A,&aa);
2350: VecGetArray(v,&x);
2351: PetscArraycpy(x,aa+col*a->lda,A->rmap->n);
2352: VecRestoreArray(v,&x);
2353: MatDenseRestoreArrayRead(A,&aa);
2354: return(0);
2355: }
2357: PETSC_INTERN PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
2358: {
2359: PetscErrorCode ierr;
2360: PetscInt i,j,m,n;
2361: const PetscScalar *a;
2364: MatGetSize(A,&m,&n);
2365: PetscArrayzero(norms,n);
2366: MatDenseGetArrayRead(A,&a);
2367: if (type == NORM_2) {
2368: for (i=0; i<n; i++) {
2369: for (j=0; j<m; j++) {
2370: norms[i] += PetscAbsScalar(a[j]*a[j]);
2371: }
2372: a += m;
2373: }
2374: } else if (type == NORM_1) {
2375: for (i=0; i<n; i++) {
2376: for (j=0; j<m; j++) {
2377: norms[i] += PetscAbsScalar(a[j]);
2378: }
2379: a += m;
2380: }
2381: } else if (type == NORM_INFINITY) {
2382: for (i=0; i<n; i++) {
2383: for (j=0; j<m; j++) {
2384: norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
2385: }
2386: a += m;
2387: }
2388: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
2389: MatDenseRestoreArrayRead(A,&a);
2390: if (type == NORM_2) {
2391: for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
2392: }
2393: return(0);
2394: }
2396: static PetscErrorCode MatSetRandom_SeqDense(Mat x,PetscRandom rctx)
2397: {
2399: PetscScalar *a;
2400: PetscInt m,n,i;
2403: MatGetSize(x,&m,&n);
2404: MatDenseGetArray(x,&a);
2405: for (i=0; i<m*n; i++) {
2406: PetscRandomGetValue(rctx,a+i);
2407: }
2408: MatDenseRestoreArray(x,&a);
2409: return(0);
2410: }
2412: static PetscErrorCode MatMissingDiagonal_SeqDense(Mat A,PetscBool *missing,PetscInt *d)
2413: {
2415: *missing = PETSC_FALSE;
2416: return(0);
2417: }
2419: /* vals is not const */
2420: static PetscErrorCode MatDenseGetColumn_SeqDense(Mat A,PetscInt col,PetscScalar **vals)
2421: {
2423: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2424: PetscScalar *v;
2427: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2428: MatDenseGetArray(A,&v);
2429: *vals = v+col*a->lda;
2430: MatDenseRestoreArray(A,&v);
2431: return(0);
2432: }
2434: static PetscErrorCode MatDenseRestoreColumn_SeqDense(Mat A,PetscScalar **vals)
2435: {
2437: *vals = 0; /* user cannot accidently use the array later */
2438: return(0);
2439: }
2441: /* -------------------------------------------------------------------*/
2442: static struct _MatOps MatOps_Values = { MatSetValues_SeqDense,
2443: MatGetRow_SeqDense,
2444: MatRestoreRow_SeqDense,
2445: MatMult_SeqDense,
2446: /* 4*/ MatMultAdd_SeqDense,
2447: MatMultTranspose_SeqDense,
2448: MatMultTransposeAdd_SeqDense,
2449: 0,
2450: 0,
2451: 0,
2452: /* 10*/ 0,
2453: MatLUFactor_SeqDense,
2454: MatCholeskyFactor_SeqDense,
2455: MatSOR_SeqDense,
2456: MatTranspose_SeqDense,
2457: /* 15*/ MatGetInfo_SeqDense,
2458: MatEqual_SeqDense,
2459: MatGetDiagonal_SeqDense,
2460: MatDiagonalScale_SeqDense,
2461: MatNorm_SeqDense,
2462: /* 20*/ MatAssemblyBegin_SeqDense,
2463: MatAssemblyEnd_SeqDense,
2464: MatSetOption_SeqDense,
2465: MatZeroEntries_SeqDense,
2466: /* 24*/ MatZeroRows_SeqDense,
2467: 0,
2468: 0,
2469: 0,
2470: 0,
2471: /* 29*/ MatSetUp_SeqDense,
2472: 0,
2473: 0,
2474: 0,
2475: 0,
2476: /* 34*/ MatDuplicate_SeqDense,
2477: 0,
2478: 0,
2479: 0,
2480: 0,
2481: /* 39*/ MatAXPY_SeqDense,
2482: MatCreateSubMatrices_SeqDense,
2483: 0,
2484: MatGetValues_SeqDense,
2485: MatCopy_SeqDense,
2486: /* 44*/ MatGetRowMax_SeqDense,
2487: MatScale_SeqDense,
2488: MatShift_Basic,
2489: 0,
2490: MatZeroRowsColumns_SeqDense,
2491: /* 49*/ MatSetRandom_SeqDense,
2492: 0,
2493: 0,
2494: 0,
2495: 0,
2496: /* 54*/ 0,
2497: 0,
2498: 0,
2499: 0,
2500: 0,
2501: /* 59*/ 0,
2502: MatDestroy_SeqDense,
2503: MatView_SeqDense,
2504: 0,
2505: 0,
2506: /* 64*/ 0,
2507: 0,
2508: 0,
2509: 0,
2510: 0,
2511: /* 69*/ MatGetRowMaxAbs_SeqDense,
2512: 0,
2513: 0,
2514: 0,
2515: 0,
2516: /* 74*/ 0,
2517: 0,
2518: 0,
2519: 0,
2520: 0,
2521: /* 79*/ 0,
2522: 0,
2523: 0,
2524: 0,
2525: /* 83*/ MatLoad_SeqDense,
2526: 0,
2527: MatIsHermitian_SeqDense,
2528: 0,
2529: 0,
2530: 0,
2531: /* 89*/ 0,
2532: 0,
2533: MatMatMultNumeric_SeqDense_SeqDense,
2534: 0,
2535: 0,
2536: /* 94*/ 0,
2537: 0,
2538: 0,
2539: MatMatTransposeMultNumeric_SeqDense_SeqDense,
2540: 0,
2541: /* 99*/ MatProductSetFromOptions_SeqDense,
2542: 0,
2543: 0,
2544: MatConjugate_SeqDense,
2545: 0,
2546: /*104*/ 0,
2547: MatRealPart_SeqDense,
2548: MatImaginaryPart_SeqDense,
2549: 0,
2550: 0,
2551: /*109*/ 0,
2552: 0,
2553: MatGetRowMin_SeqDense,
2554: MatGetColumnVector_SeqDense,
2555: MatMissingDiagonal_SeqDense,
2556: /*114*/ 0,
2557: 0,
2558: 0,
2559: 0,
2560: 0,
2561: /*119*/ 0,
2562: 0,
2563: 0,
2564: 0,
2565: 0,
2566: /*124*/ 0,
2567: MatGetColumnNorms_SeqDense,
2568: 0,
2569: 0,
2570: 0,
2571: /*129*/ 0,
2572: 0,
2573: 0,
2574: MatTransposeMatMultNumeric_SeqDense_SeqDense,
2575: 0,
2576: /*134*/ 0,
2577: 0,
2578: 0,
2579: 0,
2580: 0,
2581: /*139*/ 0,
2582: 0,
2583: 0,
2584: 0,
2585: 0,
2586: MatCreateMPIMatConcatenateSeqMat_SeqDense,
2587: /*145*/ 0,
2588: 0,
2589: 0
2590: };
2592: /*@C
2593: MatCreateSeqDense - Creates a sequential dense matrix that
2594: is stored in column major order (the usual Fortran 77 manner). Many
2595: of the matrix operations use the BLAS and LAPACK routines.
2597: Collective
2599: Input Parameters:
2600: + comm - MPI communicator, set to PETSC_COMM_SELF
2601: . m - number of rows
2602: . n - number of columns
2603: - data - optional location of matrix data in column major order. Set data=NULL for PETSc
2604: to control all matrix memory allocation.
2606: Output Parameter:
2607: . A - the matrix
2609: Notes:
2610: The data input variable is intended primarily for Fortran programmers
2611: who wish to allocate their own matrix memory space. Most users should
2612: set data=NULL.
2614: Level: intermediate
2616: .seealso: MatCreate(), MatCreateDense(), MatSetValues()
2617: @*/
2618: PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2619: {
2623: MatCreate(comm,A);
2624: MatSetSizes(*A,m,n,m,n);
2625: MatSetType(*A,MATSEQDENSE);
2626: MatSeqDenseSetPreallocation(*A,data);
2627: return(0);
2628: }
2630: /*@C
2631: MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
2633: Collective
2635: Input Parameters:
2636: + B - the matrix
2637: - data - the array (or NULL)
2639: Notes:
2640: The data input variable is intended primarily for Fortran programmers
2641: who wish to allocate their own matrix memory space. Most users should
2642: need not call this routine.
2644: Level: intermediate
2646: .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA()
2648: @*/
2649: PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2650: {
2654: PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));
2655: return(0);
2656: }
2658: PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2659: {
2660: Mat_SeqDense *b;
2664: B->preallocated = PETSC_TRUE;
2666: PetscLayoutSetUp(B->rmap);
2667: PetscLayoutSetUp(B->cmap);
2669: b = (Mat_SeqDense*)B->data;
2670: b->Mmax = B->rmap->n;
2671: b->Nmax = B->cmap->n;
2672: if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n;
2674: PetscIntMultError(b->lda,b->Nmax,NULL);
2675: if (!data) { /* petsc-allocated storage */
2676: if (!b->user_alloc) { PetscFree(b->v); }
2677: PetscCalloc1((size_t)b->lda*b->Nmax,&b->v);
2678: PetscLogObjectMemory((PetscObject)B,b->lda*b->Nmax*sizeof(PetscScalar));
2680: b->user_alloc = PETSC_FALSE;
2681: } else { /* user-allocated storage */
2682: if (!b->user_alloc) { PetscFree(b->v); }
2683: b->v = data;
2684: b->user_alloc = PETSC_TRUE;
2685: }
2686: B->assembled = PETSC_TRUE;
2687: return(0);
2688: }
2690: #if defined(PETSC_HAVE_ELEMENTAL)
2691: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
2692: {
2693: Mat mat_elemental;
2694: PetscErrorCode ierr;
2695: const PetscScalar *array;
2696: PetscScalar *v_colwise;
2697: PetscInt M=A->rmap->N,N=A->cmap->N,i,j,k,*rows,*cols;
2700: PetscMalloc3(M*N,&v_colwise,M,&rows,N,&cols);
2701: MatDenseGetArrayRead(A,&array);
2702: /* convert column-wise array into row-wise v_colwise, see MatSetValues_Elemental() */
2703: k = 0;
2704: for (j=0; j<N; j++) {
2705: cols[j] = j;
2706: for (i=0; i<M; i++) {
2707: v_colwise[j*M+i] = array[k++];
2708: }
2709: }
2710: for (i=0; i<M; i++) {
2711: rows[i] = i;
2712: }
2713: MatDenseRestoreArrayRead(A,&array);
2715: MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
2716: MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
2717: MatSetType(mat_elemental,MATELEMENTAL);
2718: MatSetUp(mat_elemental);
2720: /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
2721: MatSetValues(mat_elemental,M,rows,N,cols,v_colwise,ADD_VALUES);
2722: MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
2723: MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
2724: PetscFree3(v_colwise,rows,cols);
2726: if (reuse == MAT_INPLACE_MATRIX) {
2727: MatHeaderReplace(A,&mat_elemental);
2728: } else {
2729: *newmat = mat_elemental;
2730: }
2731: return(0);
2732: }
2733: #endif
2735: /*@C
2736: MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
2738: Input parameter:
2739: + A - the matrix
2740: - lda - the leading dimension
2742: Notes:
2743: This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
2744: it asserts that the preallocation has a leading dimension (the LDA parameter
2745: of Blas and Lapack fame) larger than M, the first dimension of the matrix.
2747: Level: intermediate
2749: .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2751: @*/
2752: PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda)
2753: {
2754: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2757: if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
2758: b->lda = lda;
2759: b->changelda = PETSC_FALSE;
2760: b->Mmax = PetscMax(b->Mmax,lda);
2761: return(0);
2762: }
2764: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2765: {
2767: PetscMPIInt size;
2770: MPI_Comm_size(comm,&size);
2771: if (size == 1) {
2772: if (scall == MAT_INITIAL_MATRIX) {
2773: MatDuplicate(inmat,MAT_COPY_VALUES,outmat);
2774: } else {
2775: MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);
2776: }
2777: } else {
2778: MatCreateMPIMatConcatenateSeqMat_MPIDense(comm,inmat,n,scall,outmat);
2779: }
2780: return(0);
2781: }
2783: /*MC
2784: MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
2786: Options Database Keys:
2787: . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
2789: Level: beginner
2791: .seealso: MatCreateSeqDense()
2793: M*/
2794: PetscErrorCode MatCreate_SeqDense(Mat B)
2795: {
2796: Mat_SeqDense *b;
2798: PetscMPIInt size;
2801: MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
2802: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
2804: PetscNewLog(B,&b);
2805: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2806: B->data = (void*)b;
2808: b->roworiented = PETSC_TRUE;
2810: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetLDA_C",MatDenseGetLDA_SeqDense);
2811: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense);
2812: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense);
2813: PetscObjectComposeFunction((PetscObject)B,"MatDensePlaceArray_C",MatDensePlaceArray_SeqDense);
2814: PetscObjectComposeFunction((PetscObject)B,"MatDenseResetArray_C",MatDenseResetArray_SeqDense);
2815: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayRead_C",MatDenseGetArray_SeqDense);
2816: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayRead_C",MatDenseRestoreArray_SeqDense);
2817: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ);
2818: #if defined(PETSC_HAVE_ELEMENTAL)
2819: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_elemental_C",MatConvert_SeqDense_Elemental);
2820: #endif
2821: #if defined(PETSC_HAVE_CUDA)
2822: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqdensecuda_C",MatConvert_SeqDense_SeqDenseCUDA);
2823: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdensecuda_seqdensecuda_C",MatProductSetFromOptions_SeqDense);
2824: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdensecuda_seqdense_C",MatProductSetFromOptions_SeqDense);
2825: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqaijcusparse_seqdense_C",MatProductSetFromOptions_SeqAIJ_SeqDense);
2826: #endif
2827: PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense);
2828: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqaij_seqdense_C",MatProductSetFromOptions_SeqAIJ_SeqDense);
2829: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdense_seqdense_C",MatProductSetFromOptions_SeqDense);
2830: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2831: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2832: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqbaij_seqdense_C",MatProductSetFromOptions_SeqXBAIJ_SeqDense);
2833: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqbaij_seqdense_C",MatMatMultSymbolic_SeqBAIJ_SeqDense);
2834: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqbaij_seqdense_C",MatMatMultNumeric_SeqBAIJ_SeqDense);
2835: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqsbaij_seqdense_C",MatProductSetFromOptions_SeqXBAIJ_SeqDense);
2836: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqsbaij_seqdense_C",MatMatMultSymbolic_SeqSBAIJ_SeqDense);
2837: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqsbaij_seqdense_C",MatMatMultNumeric_SeqSBAIJ_SeqDense);
2838: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijperm_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2839: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijperm_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2840: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijsell_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2841: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijsell_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2842: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijmkl_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2843: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijmkl_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2844: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_nest_seqdense_C",MatMatMultSymbolic_Nest_Dense);
2845: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_nest_seqdense_C",MatMatMultNumeric_Nest_Dense);
2847: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2848: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaij_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2849: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijperm_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2850: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijperm_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2851: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijsell_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2852: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijsell_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2854: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijmkl_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2855: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijmkl_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2856: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumn_C",MatDenseGetColumn_SeqDense);
2857: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_SeqDense);
2858: PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);
2859: return(0);
2860: }
2862: /*@C
2863: MatDenseGetColumn - gives access to a column of a dense matrix. This is only the local part of the column. You MUST call MatDenseRestoreColumn() to avoid memory bleeding.
2865: Not Collective
2867: Input Parameter:
2868: + mat - a MATSEQDENSE or MATMPIDENSE matrix
2869: - col - column index
2871: Output Parameter:
2872: . vals - pointer to the data
2874: Level: intermediate
2876: .seealso: MatDenseRestoreColumn()
2877: @*/
2878: PetscErrorCode MatDenseGetColumn(Mat A,PetscInt col,PetscScalar **vals)
2879: {
2883: PetscUseMethod(A,"MatDenseGetColumn_C",(Mat,PetscInt,PetscScalar**),(A,col,vals));
2884: return(0);
2885: }
2887: /*@C
2888: MatDenseRestoreColumn - returns access to a column of a dense matrix which is returned by MatDenseGetColumn().
2890: Not Collective
2892: Input Parameter:
2893: . mat - a MATSEQDENSE or MATMPIDENSE matrix
2895: Output Parameter:
2896: . vals - pointer to the data
2898: Level: intermediate
2900: .seealso: MatDenseGetColumn()
2901: @*/
2902: PetscErrorCode MatDenseRestoreColumn(Mat A,PetscScalar **vals)
2903: {
2907: PetscUseMethod(A,"MatDenseRestoreColumn_C",(Mat,PetscScalar**),(A,vals));
2908: return(0);
2909: }