Actual source code: dense.c
petsc-3.12.2 2019-11-22
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: #if defined(PETSC_MISSING_LAPACK_POTRF)
42: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
43: #else
44: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
46: PetscBLASInt info,n;
49: if (!A->rmap->n || !A->cmap->n) return(0);
50: PetscBLASIntCast(A->cmap->n,&n);
51: if (A->factortype == MAT_FACTOR_LU) {
52: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
53: if (!mat->fwork) {
54: mat->lfwork = n;
55: PetscMalloc1(mat->lfwork,&mat->fwork);
56: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
57: }
58: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
59: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
60: PetscFPTrapPop();
61: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
62: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
63: if (A->spd) {
64: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
65: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&n,mat->v,&mat->lda,&info));
66: PetscFPTrapPop();
67: MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
68: #if defined(PETSC_USE_COMPLEX)
69: } else if (A->hermitian) {
70: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
71: if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
72: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
73: PetscStackCallBLAS("LAPACKhetri",LAPACKhetri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
74: PetscFPTrapPop();
75: MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
76: #endif
77: } else { /* symmetric case */
78: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
79: if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
80: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
81: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
82: PetscFPTrapPop();
83: MatSeqDenseSymmetrize_Private(A,PETSC_FALSE);
84: }
85: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad Inversion: zero pivot in row %D",(PetscInt)info-1);
86: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
87: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
88: #endif
90: A->ops->solve = NULL;
91: A->ops->matsolve = NULL;
92: A->ops->solvetranspose = NULL;
93: A->ops->matsolvetranspose = NULL;
94: A->ops->solveadd = NULL;
95: A->ops->solvetransposeadd = NULL;
96: A->factortype = MAT_FACTOR_NONE;
97: PetscFree(A->solvertype);
98: return(0);
99: }
101: PetscErrorCode MatZeroRowsColumns_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
102: {
103: PetscErrorCode ierr;
104: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
105: PetscInt m = l->lda, n = A->cmap->n,r = A->rmap->n, i,j;
106: PetscScalar *slot,*bb,*v;
107: const PetscScalar *xx;
110: #if defined(PETSC_USE_DEBUG)
111: for (i=0; i<N; i++) {
112: if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
113: 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);
114: 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);
115: }
116: #endif
117: if (!N) return(0);
119: /* fix right hand side if needed */
120: if (x && b) {
121: Vec xt;
123: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
124: VecDuplicate(x,&xt);
125: VecCopy(x,xt);
126: VecScale(xt,-1.0);
127: MatMultAdd(A,xt,b,b);
128: VecDestroy(&xt);
129: VecGetArrayRead(x,&xx);
130: VecGetArray(b,&bb);
131: for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
132: VecRestoreArrayRead(x,&xx);
133: VecRestoreArray(b,&bb);
134: }
136: MatDenseGetArray(A,&v);
137: for (i=0; i<N; i++) {
138: slot = v + rows[i]*m;
139: PetscArrayzero(slot,r);
140: }
141: for (i=0; i<N; i++) {
142: slot = v + rows[i];
143: for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
144: }
145: if (diag != 0.0) {
146: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
147: for (i=0; i<N; i++) {
148: slot = v + (m+1)*rows[i];
149: *slot = diag;
150: }
151: }
152: MatDenseRestoreArray(A,&v);
153: return(0);
154: }
156: PetscErrorCode MatPtAPNumeric_SeqDense_SeqDense(Mat A,Mat P,Mat C)
157: {
158: Mat_SeqDense *c = (Mat_SeqDense*)(C->data);
162: if (c->ptapwork) {
163: (*C->ops->matmultnumeric)(A,P,c->ptapwork);
164: (*C->ops->transposematmultnumeric)(P,c->ptapwork,C);
165: } else { /* first time went trough the basic. Should we add better dispatching for subclasses? */
166: MatPtAP_Basic(A,P,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);
167: }
168: return(0);
169: }
171: PetscErrorCode MatPtAPSymbolic_SeqDense_SeqDense(Mat A,Mat P,PetscReal fill,Mat *C)
172: {
173: Mat_SeqDense *c;
174: PetscBool flg;
178: PetscObjectTypeCompare((PetscObject)P,((PetscObject)A)->type_name,&flg);
179: MatCreate(PetscObjectComm((PetscObject)A),C);
180: MatSetSizes(*C,P->cmap->n,P->cmap->n,P->cmap->N,P->cmap->N);
181: MatSetType(*C,flg ? ((PetscObject)A)->type_name : MATDENSE);
182: MatSeqDenseSetPreallocation(*C,NULL);
183: c = (Mat_SeqDense*)((*C)->data);
184: MatCreate(PetscObjectComm((PetscObject)A),&c->ptapwork);
185: MatSetSizes(c->ptapwork,A->rmap->n,P->cmap->n,A->rmap->N,P->cmap->N);
186: MatSetType(c->ptapwork,flg ? ((PetscObject)A)->type_name : MATDENSE);
187: MatSeqDenseSetPreallocation(c->ptapwork,NULL);
188: return(0);
189: }
191: PETSC_INTERN PetscErrorCode MatPtAP_SeqDense_SeqDense(Mat A,Mat P,MatReuse reuse,PetscReal fill,Mat *C)
192: {
196: if (reuse == MAT_INITIAL_MATRIX) {
197: MatPtAPSymbolic_SeqDense_SeqDense(A,P,fill,C);
198: }
199: PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
200: (*(*C)->ops->ptapnumeric)(A,P,*C);
201: PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
202: return(0);
203: }
205: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
206: {
207: Mat B = NULL;
208: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
209: Mat_SeqDense *b;
211: PetscInt *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->N,i;
212: MatScalar *av=a->a;
213: PetscBool isseqdense;
216: if (reuse == MAT_REUSE_MATRIX) {
217: PetscObjectTypeCompare((PetscObject)*newmat,MATSEQDENSE,&isseqdense);
218: if (!isseqdense) SETERRQ1(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type_name);
219: }
220: if (reuse != MAT_REUSE_MATRIX) {
221: MatCreate(PetscObjectComm((PetscObject)A),&B);
222: MatSetSizes(B,m,n,m,n);
223: MatSetType(B,MATSEQDENSE);
224: MatSeqDenseSetPreallocation(B,NULL);
225: b = (Mat_SeqDense*)(B->data);
226: } else {
227: b = (Mat_SeqDense*)((*newmat)->data);
228: PetscArrayzero(b->v,m*n);
229: }
230: for (i=0; i<m; i++) {
231: PetscInt j;
232: for (j=0;j<ai[1]-ai[0];j++) {
233: b->v[*aj*m+i] = *av;
234: aj++;
235: av++;
236: }
237: ai++;
238: }
240: if (reuse == MAT_INPLACE_MATRIX) {
241: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
242: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
243: MatHeaderReplace(A,&B);
244: } else {
245: if (B) *newmat = B;
246: MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
247: MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
248: }
249: return(0);
250: }
252: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
253: {
254: Mat B;
255: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
257: PetscInt i, j;
258: PetscInt *rows, *nnz;
259: MatScalar *aa = a->v, *vals;
262: MatCreate(PetscObjectComm((PetscObject)A),&B);
263: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
264: MatSetType(B,MATSEQAIJ);
265: PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals);
266: for (j=0; j<A->cmap->n; j++) {
267: for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) ++nnz[i];
268: aa += a->lda;
269: }
270: MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz);
271: aa = a->v;
272: for (j=0; j<A->cmap->n; j++) {
273: PetscInt numRows = 0;
274: for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) {rows[numRows] = i; vals[numRows++] = aa[i];}
275: MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES);
276: aa += a->lda;
277: }
278: PetscFree3(rows,nnz,vals);
279: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
280: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
282: if (reuse == MAT_INPLACE_MATRIX) {
283: MatHeaderReplace(A,&B);
284: } else {
285: *newmat = B;
286: }
287: return(0);
288: }
290: PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
291: {
292: Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
293: const PetscScalar *xv;
294: PetscScalar *yv;
295: PetscBLASInt N,m,ldax,lday,one = 1;
296: PetscErrorCode ierr;
299: MatDenseGetArrayRead(X,&xv);
300: MatDenseGetArray(Y,&yv);
301: PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);
302: PetscBLASIntCast(X->rmap->n,&m);
303: PetscBLASIntCast(x->lda,&ldax);
304: PetscBLASIntCast(y->lda,&lday);
305: if (ldax>m || lday>m) {
306: PetscInt j;
308: for (j=0; j<X->cmap->n; j++) {
309: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&alpha,xv+j*ldax,&one,yv+j*lday,&one));
310: }
311: } else {
312: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&alpha,xv,&one,yv,&one));
313: }
314: MatDenseRestoreArrayRead(X,&xv);
315: MatDenseRestoreArray(Y,&yv);
316: PetscLogFlops(PetscMax(2*N-1,0));
317: return(0);
318: }
320: static PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
321: {
322: PetscLogDouble N = A->rmap->n*A->cmap->n;
325: info->block_size = 1.0;
326: info->nz_allocated = N;
327: info->nz_used = N;
328: info->nz_unneeded = 0;
329: info->assemblies = A->num_ass;
330: info->mallocs = 0;
331: info->memory = ((PetscObject)A)->mem;
332: info->fill_ratio_given = 0;
333: info->fill_ratio_needed = 0;
334: info->factor_mallocs = 0;
335: return(0);
336: }
338: static PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
339: {
340: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
341: PetscScalar *v;
343: PetscBLASInt one = 1,j,nz,lda;
346: MatDenseGetArray(A,&v);
347: PetscBLASIntCast(a->lda,&lda);
348: if (lda>A->rmap->n) {
349: PetscBLASIntCast(A->rmap->n,&nz);
350: for (j=0; j<A->cmap->n; j++) {
351: PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v+j*lda,&one));
352: }
353: } else {
354: PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);
355: PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v,&one));
356: }
357: PetscLogFlops(nz);
358: MatDenseRestoreArray(A,&v);
359: return(0);
360: }
362: static PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl)
363: {
364: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
365: PetscInt i,j,m = A->rmap->n,N = a->lda;
366: const PetscScalar *v;
367: PetscErrorCode ierr;
370: *fl = PETSC_FALSE;
371: if (A->rmap->n != A->cmap->n) return(0);
372: MatDenseGetArrayRead(A,&v);
373: for (i=0; i<m; i++) {
374: for (j=i; j<m; j++) {
375: if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) return(0);
376: }
377: }
378: MatDenseRestoreArrayRead(A,&v);
379: *fl = PETSC_TRUE;
380: return(0);
381: }
383: PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
384: {
385: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
387: PetscInt lda = (PetscInt)mat->lda,j,m;
390: PetscLayoutReference(A->rmap,&newi->rmap);
391: PetscLayoutReference(A->cmap,&newi->cmap);
392: MatSeqDenseSetPreallocation(newi,NULL);
393: if (cpvalues == MAT_COPY_VALUES) {
394: const PetscScalar *av;
395: PetscScalar *v;
397: MatDenseGetArrayRead(A,&av);
398: MatDenseGetArray(newi,&v);
399: if (lda>A->rmap->n) {
400: m = A->rmap->n;
401: for (j=0; j<A->cmap->n; j++) {
402: PetscArraycpy(v+j*m,av+j*lda,m);
403: }
404: } else {
405: PetscArraycpy(v,av,A->rmap->n*A->cmap->n);
406: }
407: MatDenseRestoreArray(newi,&v);
408: MatDenseRestoreArrayRead(A,&av);
409: }
410: return(0);
411: }
413: PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
414: {
418: MatCreate(PetscObjectComm((PetscObject)A),newmat);
419: MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
420: MatSetType(*newmat,((PetscObject)A)->type_name);
421: MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);
422: return(0);
423: }
425: static PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
426: {
427: MatFactorInfo info;
431: MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
432: (*fact->ops->lufactor)(fact,0,0,&info);
433: return(0);
434: }
436: static PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
437: {
438: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
439: PetscErrorCode ierr;
440: const PetscScalar *x;
441: PetscScalar *y;
442: PetscBLASInt one = 1,info,m;
445: PetscBLASIntCast(A->rmap->n,&m);
446: VecGetArrayRead(xx,&x);
447: VecGetArray(yy,&y);
448: PetscArraycpy(y,x,A->rmap->n);
449: if (A->factortype == MAT_FACTOR_LU) {
450: #if defined(PETSC_MISSING_LAPACK_GETRS)
451: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
452: #else
453: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
454: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
455: PetscFPTrapPop();
456: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
457: #endif
458: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
459: #if defined(PETSC_MISSING_LAPACK_POTRS)
460: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
461: #else
462: if (A->spd) {
463: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
464: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
465: PetscFPTrapPop();
466: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
467: #if defined(PETSC_USE_COMPLEX)
468: } else if (A->hermitian) {
469: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
470: PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
471: PetscFPTrapPop();
472: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
473: #endif
474: } else { /* symmetric case */
475: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
476: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
477: PetscFPTrapPop();
478: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
479: }
480: #endif
481: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
482: VecRestoreArrayRead(xx,&x);
483: VecRestoreArray(yy,&y);
484: PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
485: return(0);
486: }
488: static PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
489: {
490: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
491: PetscErrorCode ierr;
492: const PetscScalar *b;
493: PetscScalar *x;
494: PetscInt n;
495: PetscBLASInt nrhs,info,m;
498: PetscBLASIntCast(A->rmap->n,&m);
499: MatGetSize(B,NULL,&n);
500: PetscBLASIntCast(n,&nrhs);
501: MatDenseGetArrayRead(B,&b);
502: MatDenseGetArray(X,&x);
504: PetscArraycpy(x,b,m*nrhs);
506: if (A->factortype == MAT_FACTOR_LU) {
507: #if defined(PETSC_MISSING_LAPACK_GETRS)
508: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
509: #else
510: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
511: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
512: PetscFPTrapPop();
513: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
514: #endif
515: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
516: #if defined(PETSC_MISSING_LAPACK_POTRS)
517: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
518: #else
519: if (A->spd) {
520: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
521: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
522: PetscFPTrapPop();
523: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
524: #if defined(PETSC_USE_COMPLEX)
525: } else if (A->hermitian) {
526: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
527: PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
528: PetscFPTrapPop();
529: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
530: #endif
531: } else { /* symmetric case */
532: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
533: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
534: PetscFPTrapPop();
535: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
536: }
537: #endif
538: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
540: MatDenseRestoreArrayRead(B,&b);
541: MatDenseRestoreArray(X,&x);
542: PetscLogFlops(nrhs*(2.0*m*m - m));
543: return(0);
544: }
546: static PetscErrorCode MatConjugate_SeqDense(Mat);
548: static PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
549: {
550: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
551: PetscErrorCode ierr;
552: const PetscScalar *x;
553: PetscScalar *y;
554: PetscBLASInt one = 1,info,m;
557: PetscBLASIntCast(A->rmap->n,&m);
558: VecGetArrayRead(xx,&x);
559: VecGetArray(yy,&y);
560: PetscArraycpy(y,x,A->rmap->n);
561: if (A->factortype == MAT_FACTOR_LU) {
562: #if defined(PETSC_MISSING_LAPACK_GETRS)
563: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
564: #else
565: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
566: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
567: PetscFPTrapPop();
568: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
569: #endif
570: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
571: #if defined(PETSC_MISSING_LAPACK_POTRS)
572: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
573: #else
574: if (A->spd) {
575: #if defined(PETSC_USE_COMPLEX)
576: MatConjugate_SeqDense(A);
577: #endif
578: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
579: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
580: PetscFPTrapPop();
581: #if defined(PETSC_USE_COMPLEX)
582: MatConjugate_SeqDense(A);
583: #endif
584: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
585: #if defined(PETSC_USE_COMPLEX)
586: } else if (A->hermitian) {
587: MatConjugate_SeqDense(A);
588: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
589: PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
590: PetscFPTrapPop();
591: MatConjugate_SeqDense(A);
592: #endif
593: } else { /* symmetric case */
594: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
595: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
596: PetscFPTrapPop();
597: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
598: }
599: #endif
600: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
601: VecRestoreArrayRead(xx,&x);
602: VecRestoreArray(yy,&y);
603: PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);
604: return(0);
605: }
607: /* ---------------------------------------------------------------*/
608: /* COMMENT: I have chosen to hide row permutation in the pivots,
609: rather than put it in the Mat->row slot.*/
610: PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
611: {
612: #if defined(PETSC_MISSING_LAPACK_GETRF)
614: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
615: #else
616: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
618: PetscBLASInt n,m,info;
621: PetscBLASIntCast(A->cmap->n,&n);
622: PetscBLASIntCast(A->rmap->n,&m);
623: if (!mat->pivots) {
624: PetscMalloc1(A->rmap->n,&mat->pivots);
625: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
626: }
627: if (!A->rmap->n || !A->cmap->n) return(0);
628: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
629: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info));
630: PetscFPTrapPop();
632: if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
633: if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
635: A->ops->solve = MatSolve_SeqDense;
636: A->ops->matsolve = MatMatSolve_SeqDense;
637: A->ops->solvetranspose = MatSolveTranspose_SeqDense;
638: A->factortype = MAT_FACTOR_LU;
640: PetscFree(A->solvertype);
641: PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);
643: PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);
644: #endif
645: return(0);
646: }
648: /* Cholesky as L*L^T or L*D*L^T and the symmetric/hermitian complex variants */
649: PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
650: {
651: #if defined(PETSC_MISSING_LAPACK_POTRF)
653: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
654: #else
655: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
657: PetscBLASInt info,n;
660: PetscBLASIntCast(A->cmap->n,&n);
661: if (!A->rmap->n || !A->cmap->n) return(0);
662: if (A->spd) {
663: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
664: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info));
665: PetscFPTrapPop();
666: #if defined(PETSC_USE_COMPLEX)
667: } else if (A->hermitian) {
668: if (!mat->pivots) {
669: PetscMalloc1(A->rmap->n,&mat->pivots);
670: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
671: }
672: if (!mat->fwork) {
673: PetscScalar dummy;
675: mat->lfwork = -1;
676: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
677: PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
678: PetscFPTrapPop();
679: mat->lfwork = (PetscInt)PetscRealPart(dummy);
680: PetscMalloc1(mat->lfwork,&mat->fwork);
681: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
682: }
683: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
684: PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
685: PetscFPTrapPop();
686: #endif
687: } else { /* symmetric case */
688: if (!mat->pivots) {
689: PetscMalloc1(A->rmap->n,&mat->pivots);
690: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
691: }
692: if (!mat->fwork) {
693: PetscScalar dummy;
695: mat->lfwork = -1;
696: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
697: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
698: PetscFPTrapPop();
699: mat->lfwork = (PetscInt)PetscRealPart(dummy);
700: PetscMalloc1(mat->lfwork,&mat->fwork);
701: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
702: }
703: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
704: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
705: PetscFPTrapPop();
706: }
707: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
709: A->ops->solve = MatSolve_SeqDense;
710: A->ops->matsolve = MatMatSolve_SeqDense;
711: A->ops->solvetranspose = MatSolveTranspose_SeqDense;
712: A->factortype = MAT_FACTOR_CHOLESKY;
714: PetscFree(A->solvertype);
715: PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);
717: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
718: #endif
719: return(0);
720: }
723: PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
724: {
726: MatFactorInfo info;
729: info.fill = 1.0;
731: MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
732: (*fact->ops->choleskyfactor)(fact,0,&info);
733: return(0);
734: }
736: PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
737: {
739: fact->assembled = PETSC_TRUE;
740: fact->preallocated = PETSC_TRUE;
741: fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
742: fact->ops->solve = MatSolve_SeqDense;
743: fact->ops->matsolve = MatMatSolve_SeqDense;
744: fact->ops->solvetranspose = MatSolveTranspose_SeqDense;
745: return(0);
746: }
748: PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
749: {
751: fact->preallocated = PETSC_TRUE;
752: fact->assembled = PETSC_TRUE;
753: fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
754: fact->ops->solve = MatSolve_SeqDense;
755: fact->ops->matsolve = MatMatSolve_SeqDense;
756: fact->ops->solvetranspose = MatSolveTranspose_SeqDense;
757: return(0);
758: }
760: /* uses LAPACK */
761: PETSC_INTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
762: {
766: MatCreate(PetscObjectComm((PetscObject)A),fact);
767: MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
768: MatSetType(*fact,MATDENSE);
769: if (ftype == MAT_FACTOR_LU) {
770: (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
771: } else {
772: (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
773: }
774: (*fact)->factortype = ftype;
776: PetscFree((*fact)->solvertype);
777: PetscStrallocpy(MATSOLVERPETSC,&(*fact)->solvertype);
778: return(0);
779: }
781: /* ------------------------------------------------------------------*/
782: static PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
783: {
784: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
785: PetscScalar *x,*v = mat->v,zero = 0.0,xt;
786: const PetscScalar *b;
787: PetscErrorCode ierr;
788: PetscInt m = A->rmap->n,i;
789: PetscBLASInt o = 1,bm;
792: #if defined(PETSC_HAVE_CUDA)
793: if (A->offloadmask == PETSC_OFFLOAD_GPU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented");
794: #endif
795: if (shift == -1) shift = 0.0; /* negative shift indicates do not error on zero diagonal; this code never zeros on zero diagonal */
796: PetscBLASIntCast(m,&bm);
797: if (flag & SOR_ZERO_INITIAL_GUESS) {
798: /* this is a hack fix, should have another version without the second BLASdotu */
799: VecSet(xx,zero);
800: }
801: VecGetArray(xx,&x);
802: VecGetArrayRead(bb,&b);
803: its = its*lits;
804: if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
805: while (its--) {
806: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
807: for (i=0; i<m; i++) {
808: PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
809: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
810: }
811: }
812: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
813: for (i=m-1; i>=0; i--) {
814: PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
815: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
816: }
817: }
818: }
819: VecRestoreArrayRead(bb,&b);
820: VecRestoreArray(xx,&x);
821: return(0);
822: }
824: /* -----------------------------------------------------------------*/
825: PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
826: {
827: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
828: const PetscScalar *v = mat->v,*x;
829: PetscScalar *y;
830: PetscErrorCode ierr;
831: PetscBLASInt m, n,_One=1;
832: PetscScalar _DOne=1.0,_DZero=0.0;
835: PetscBLASIntCast(A->rmap->n,&m);
836: PetscBLASIntCast(A->cmap->n,&n);
837: VecGetArrayRead(xx,&x);
838: VecGetArrayWrite(yy,&y);
839: if (!A->rmap->n || !A->cmap->n) {
840: PetscBLASInt i;
841: for (i=0; i<n; i++) y[i] = 0.0;
842: } else {
843: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One));
844: PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);
845: }
846: VecRestoreArrayRead(xx,&x);
847: VecRestoreArrayWrite(yy,&y);
848: return(0);
849: }
851: PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
852: {
853: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
854: PetscScalar *y,_DOne=1.0,_DZero=0.0;
855: PetscErrorCode ierr;
856: PetscBLASInt m, n, _One=1;
857: const PetscScalar *v = mat->v,*x;
860: PetscBLASIntCast(A->rmap->n,&m);
861: PetscBLASIntCast(A->cmap->n,&n);
862: VecGetArrayRead(xx,&x);
863: VecGetArrayWrite(yy,&y);
864: if (!A->rmap->n || !A->cmap->n) {
865: PetscBLASInt i;
866: for (i=0; i<m; i++) y[i] = 0.0;
867: } else {
868: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One));
869: PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);
870: }
871: VecRestoreArrayRead(xx,&x);
872: VecRestoreArrayWrite(yy,&y);
873: return(0);
874: }
876: PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
877: {
878: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
879: const PetscScalar *v = mat->v,*x;
880: PetscScalar *y,_DOne=1.0;
881: PetscErrorCode ierr;
882: PetscBLASInt m, n, _One=1;
885: PetscBLASIntCast(A->rmap->n,&m);
886: PetscBLASIntCast(A->cmap->n,&n);
887: if (!A->rmap->n || !A->cmap->n) return(0);
888: if (zz != yy) {VecCopy(zz,yy);}
889: VecGetArrayRead(xx,&x);
890: VecGetArray(yy,&y);
891: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
892: VecRestoreArrayRead(xx,&x);
893: VecRestoreArray(yy,&y);
894: PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
895: return(0);
896: }
898: PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
899: {
900: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
901: const PetscScalar *v = mat->v,*x;
902: PetscScalar *y;
903: PetscErrorCode ierr;
904: PetscBLASInt m, n, _One=1;
905: PetscScalar _DOne=1.0;
908: PetscBLASIntCast(A->rmap->n,&m);
909: PetscBLASIntCast(A->cmap->n,&n);
910: if (!A->rmap->n || !A->cmap->n) return(0);
911: if (zz != yy) {VecCopy(zz,yy);}
912: VecGetArrayRead(xx,&x);
913: VecGetArray(yy,&y);
914: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
915: VecRestoreArrayRead(xx,&x);
916: VecRestoreArray(yy,&y);
917: PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
918: return(0);
919: }
921: /* -----------------------------------------------------------------*/
922: static PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
923: {
924: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
926: PetscInt i;
929: *ncols = A->cmap->n;
930: if (cols) {
931: PetscMalloc1(A->cmap->n+1,cols);
932: for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
933: }
934: if (vals) {
935: const PetscScalar *v;
937: MatDenseGetArrayRead(A,&v);
938: PetscMalloc1(A->cmap->n+1,vals);
939: v += row;
940: for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
941: MatDenseRestoreArrayRead(A,&v);
942: }
943: return(0);
944: }
946: static PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
947: {
951: if (cols) {PetscFree(*cols);}
952: if (vals) {PetscFree(*vals); }
953: return(0);
954: }
955: /* ----------------------------------------------------------------*/
956: static PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
957: {
958: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
959: PetscScalar *av;
960: PetscInt i,j,idx=0;
961: #if defined(PETSC_HAVE_CUDA)
962: PetscOffloadMask oldf;
963: #endif
964: PetscErrorCode ierr;
967: MatDenseGetArray(A,&av);
968: if (!mat->roworiented) {
969: if (addv == INSERT_VALUES) {
970: for (j=0; j<n; j++) {
971: if (indexn[j] < 0) {idx += m; continue;}
972: #if defined(PETSC_USE_DEBUG)
973: 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);
974: #endif
975: for (i=0; i<m; i++) {
976: if (indexm[i] < 0) {idx++; continue;}
977: #if defined(PETSC_USE_DEBUG)
978: 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);
979: #endif
980: av[indexn[j]*mat->lda + indexm[i]] = v[idx++];
981: }
982: }
983: } else {
984: for (j=0; j<n; j++) {
985: if (indexn[j] < 0) {idx += m; continue;}
986: #if defined(PETSC_USE_DEBUG)
987: 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);
988: #endif
989: for (i=0; i<m; i++) {
990: if (indexm[i] < 0) {idx++; continue;}
991: #if defined(PETSC_USE_DEBUG)
992: 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);
993: #endif
994: av[indexn[j]*mat->lda + indexm[i]] += v[idx++];
995: }
996: }
997: }
998: } else {
999: if (addv == INSERT_VALUES) {
1000: for (i=0; i<m; i++) {
1001: if (indexm[i] < 0) { idx += n; continue;}
1002: #if defined(PETSC_USE_DEBUG)
1003: 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);
1004: #endif
1005: for (j=0; j<n; j++) {
1006: if (indexn[j] < 0) { idx++; continue;}
1007: #if defined(PETSC_USE_DEBUG)
1008: 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);
1009: #endif
1010: av[indexn[j]*mat->lda + indexm[i]] = v[idx++];
1011: }
1012: }
1013: } else {
1014: for (i=0; i<m; i++) {
1015: if (indexm[i] < 0) { idx += n; continue;}
1016: #if defined(PETSC_USE_DEBUG)
1017: 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);
1018: #endif
1019: for (j=0; j<n; j++) {
1020: if (indexn[j] < 0) { idx++; continue;}
1021: #if defined(PETSC_USE_DEBUG)
1022: 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);
1023: #endif
1024: av[indexn[j]*mat->lda + indexm[i]] += v[idx++];
1025: }
1026: }
1027: }
1028: }
1029: /* hack to prevent unneeded copy to the GPU while returning the array */
1030: #if defined(PETSC_HAVE_CUDA)
1031: oldf = A->offloadmask;
1032: A->offloadmask = PETSC_OFFLOAD_GPU;
1033: #endif
1034: MatDenseRestoreArray(A,&av);
1035: #if defined(PETSC_HAVE_CUDA)
1036: A->offloadmask = (oldf == PETSC_OFFLOAD_UNALLOCATED ? PETSC_OFFLOAD_UNALLOCATED : PETSC_OFFLOAD_CPU);
1037: #endif
1038: return(0);
1039: }
1041: static PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
1042: {
1043: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1044: const PetscScalar *vv;
1045: PetscInt i,j;
1046: PetscErrorCode ierr;
1049: MatDenseGetArrayRead(A,&vv);
1050: /* row-oriented output */
1051: for (i=0; i<m; i++) {
1052: if (indexm[i] < 0) {v += n;continue;}
1053: 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);
1054: for (j=0; j<n; j++) {
1055: if (indexn[j] < 0) {v++; continue;}
1056: 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);
1057: *v++ = vv[indexn[j]*mat->lda + indexm[i]];
1058: }
1059: }
1060: MatDenseRestoreArrayRead(A,&vv);
1061: return(0);
1062: }
1064: /* -----------------------------------------------------------------*/
1066: static PetscErrorCode MatLoad_SeqDense_Binary(Mat newmat,PetscViewer viewer)
1067: {
1068: Mat_SeqDense *a;
1070: PetscInt *scols,i,j,nz,header[4];
1071: int fd;
1072: PetscMPIInt size;
1073: PetscInt *rowlengths = 0,M,N,*cols,grows,gcols;
1074: PetscScalar *vals,*svals,*v,*w;
1075: MPI_Comm comm;
1078: PetscObjectGetComm((PetscObject)viewer,&comm);
1079: MPI_Comm_size(comm,&size);
1080: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
1081: PetscViewerBinaryGetDescriptor(viewer,&fd);
1082: PetscBinaryRead(fd,header,4,NULL,PETSC_INT);
1083: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
1084: M = header[1]; N = header[2]; nz = header[3];
1086: /* set global size if not set already*/
1087: if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
1088: MatSetSizes(newmat,M,N,M,N);
1089: } else {
1090: /* if sizes and type are already set, check if the vector global sizes are correct */
1091: MatGetSize(newmat,&grows,&gcols);
1092: if (M != grows || N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols);
1093: }
1094: a = (Mat_SeqDense*)newmat->data;
1095: if (!a->user_alloc) {
1096: MatSeqDenseSetPreallocation(newmat,NULL);
1097: }
1099: if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
1100: a = (Mat_SeqDense*)newmat->data;
1101: v = a->v;
1102: /* Allocate some temp space to read in the values and then flip them
1103: from row major to column major */
1104: PetscMalloc1(M*N > 0 ? M*N : 1,&w);
1105: /* read in nonzero values */
1106: PetscBinaryRead(fd,w,M*N,NULL,PETSC_SCALAR);
1107: /* now flip the values and store them in the matrix*/
1108: for (j=0; j<N; j++) {
1109: for (i=0; i<M; i++) {
1110: *v++ =w[i*N+j];
1111: }
1112: }
1113: PetscFree(w);
1114: } else {
1115: /* read row lengths */
1116: PetscMalloc1(M+1,&rowlengths);
1117: PetscBinaryRead(fd,rowlengths,M,NULL,PETSC_INT);
1119: a = (Mat_SeqDense*)newmat->data;
1120: v = a->v;
1122: /* read column indices and nonzeros */
1123: PetscMalloc1(nz+1,&scols);
1124: cols = scols;
1125: PetscBinaryRead(fd,cols,nz,NULL,PETSC_INT);
1126: PetscMalloc1(nz+1,&svals);
1127: vals = svals;
1128: PetscBinaryRead(fd,vals,nz,NULL,PETSC_SCALAR);
1130: /* insert into matrix */
1131: for (i=0; i<M; i++) {
1132: for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
1133: svals += rowlengths[i]; scols += rowlengths[i];
1134: }
1135: PetscFree(vals);
1136: PetscFree(cols);
1137: PetscFree(rowlengths);
1138: }
1139: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1140: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1141: return(0);
1142: }
1144: PetscErrorCode MatLoad_SeqDense(Mat newMat, PetscViewer viewer)
1145: {
1146: PetscBool isbinary, ishdf5;
1152: /* force binary viewer to load .info file if it has not yet done so */
1153: PetscViewerSetUp(viewer);
1154: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1155: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5, &ishdf5);
1156: if (isbinary) {
1157: MatLoad_SeqDense_Binary(newMat,viewer);
1158: } else if (ishdf5) {
1159: #if defined(PETSC_HAVE_HDF5)
1160: MatLoad_Dense_HDF5(newMat,viewer);
1161: #else
1162: SETERRQ(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1163: #endif
1164: } else {
1165: 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);
1166: }
1167: return(0);
1168: }
1170: static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
1171: {
1172: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1173: PetscErrorCode ierr;
1174: PetscInt i,j;
1175: const char *name;
1176: PetscScalar *v,*av;
1177: PetscViewerFormat format;
1178: #if defined(PETSC_USE_COMPLEX)
1179: PetscBool allreal = PETSC_TRUE;
1180: #endif
1183: MatDenseGetArrayRead(A,(const PetscScalar**)&av);
1184: PetscViewerGetFormat(viewer,&format);
1185: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1186: return(0); /* do nothing for now */
1187: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
1188: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1189: for (i=0; i<A->rmap->n; i++) {
1190: v = av + i;
1191: PetscViewerASCIIPrintf(viewer,"row %D:",i);
1192: for (j=0; j<A->cmap->n; j++) {
1193: #if defined(PETSC_USE_COMPLEX)
1194: if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
1195: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1196: } else if (PetscRealPart(*v)) {
1197: PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)PetscRealPart(*v));
1198: }
1199: #else
1200: if (*v) {
1201: PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)*v);
1202: }
1203: #endif
1204: v += a->lda;
1205: }
1206: PetscViewerASCIIPrintf(viewer,"\n");
1207: }
1208: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1209: } else {
1210: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1211: #if defined(PETSC_USE_COMPLEX)
1212: /* determine if matrix has all real values */
1213: v = av;
1214: for (i=0; i<A->rmap->n*A->cmap->n; i++) {
1215: if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;}
1216: }
1217: #endif
1218: if (format == PETSC_VIEWER_ASCII_MATLAB) {
1219: PetscObjectGetName((PetscObject)A,&name);
1220: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);
1221: PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);
1222: PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
1223: }
1225: for (i=0; i<A->rmap->n; i++) {
1226: v = av + i;
1227: for (j=0; j<A->cmap->n; j++) {
1228: #if defined(PETSC_USE_COMPLEX)
1229: if (allreal) {
1230: PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));
1231: } else {
1232: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1233: }
1234: #else
1235: PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);
1236: #endif
1237: v += a->lda;
1238: }
1239: PetscViewerASCIIPrintf(viewer,"\n");
1240: }
1241: if (format == PETSC_VIEWER_ASCII_MATLAB) {
1242: PetscViewerASCIIPrintf(viewer,"];\n");
1243: }
1244: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1245: }
1246: MatDenseRestoreArrayRead(A,(const PetscScalar**)&av);
1247: PetscViewerFlush(viewer);
1248: return(0);
1249: }
1251: static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
1252: {
1253: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1254: PetscErrorCode ierr;
1255: int fd;
1256: PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
1257: PetscScalar *av,*v,*anonz,*vals;
1258: PetscViewerFormat format;
1261: PetscViewerBinaryGetDescriptor(viewer,&fd);
1262: MatDenseGetArrayRead(A,(const PetscScalar**)&av);
1263: PetscViewerGetFormat(viewer,&format);
1264: if (format == PETSC_VIEWER_NATIVE) {
1265: /* store the matrix as a dense matrix */
1266: PetscMalloc1(4,&col_lens);
1268: col_lens[0] = MAT_FILE_CLASSID;
1269: col_lens[1] = m;
1270: col_lens[2] = n;
1271: col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
1273: PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);
1274: PetscFree(col_lens);
1276: /* write out matrix, by rows */
1277: PetscMalloc1(m*n+1,&vals);
1278: v = av;
1279: for (j=0; j<n; j++) {
1280: for (i=0; i<m; i++) {
1281: vals[j + i*n] = *v++;
1282: }
1283: }
1284: PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);
1285: PetscFree(vals);
1286: } else {
1287: PetscMalloc1(4+nz,&col_lens);
1289: col_lens[0] = MAT_FILE_CLASSID;
1290: col_lens[1] = m;
1291: col_lens[2] = n;
1292: col_lens[3] = nz;
1294: /* store lengths of each row and write (including header) to file */
1295: for (i=0; i<m; i++) col_lens[4+i] = n;
1296: PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);
1298: /* Possibly should write in smaller increments, not whole matrix at once? */
1299: /* store column indices (zero start index) */
1300: ict = 0;
1301: for (i=0; i<m; i++) {
1302: for (j=0; j<n; j++) col_lens[ict++] = j;
1303: }
1304: PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);
1305: PetscFree(col_lens);
1307: /* store nonzero values */
1308: PetscMalloc1(nz+1,&anonz);
1309: ict = 0;
1310: for (i=0; i<m; i++) {
1311: v = av + i;
1312: for (j=0; j<n; j++) {
1313: anonz[ict++] = *v; v += a->lda;
1314: }
1315: }
1316: PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);
1317: PetscFree(anonz);
1318: }
1319: MatDenseRestoreArrayRead(A,(const PetscScalar**)&av);
1320: return(0);
1321: }
1323: #include <petscdraw.h>
1324: static PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1325: {
1326: Mat A = (Mat) Aa;
1327: PetscErrorCode ierr;
1328: PetscInt m = A->rmap->n,n = A->cmap->n,i,j;
1329: int color = PETSC_DRAW_WHITE;
1330: const PetscScalar *v;
1331: PetscViewer viewer;
1332: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r;
1333: PetscViewerFormat format;
1336: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
1337: PetscViewerGetFormat(viewer,&format);
1338: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
1340: /* Loop over matrix elements drawing boxes */
1341: MatDenseGetArrayRead(A,&v);
1342: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1343: PetscDrawCollectiveBegin(draw);
1344: /* Blue for negative and Red for positive */
1345: for (j = 0; j < n; j++) {
1346: x_l = j; x_r = x_l + 1.0;
1347: for (i = 0; i < m; i++) {
1348: y_l = m - i - 1.0;
1349: y_r = y_l + 1.0;
1350: if (PetscRealPart(v[j*m+i]) > 0.) color = PETSC_DRAW_RED;
1351: else if (PetscRealPart(v[j*m+i]) < 0.) color = PETSC_DRAW_BLUE;
1352: else continue;
1353: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1354: }
1355: }
1356: PetscDrawCollectiveEnd(draw);
1357: } else {
1358: /* use contour shading to indicate magnitude of values */
1359: /* first determine max of all nonzero values */
1360: PetscReal minv = 0.0, maxv = 0.0;
1361: PetscDraw popup;
1363: for (i=0; i < m*n; i++) {
1364: if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1365: }
1366: if (minv >= maxv) maxv = minv + PETSC_SMALL;
1367: PetscDrawGetPopup(draw,&popup);
1368: PetscDrawScalePopup(popup,minv,maxv);
1370: PetscDrawCollectiveBegin(draw);
1371: for (j=0; j<n; j++) {
1372: x_l = j;
1373: x_r = x_l + 1.0;
1374: for (i=0; i<m; i++) {
1375: y_l = m - i - 1.0;
1376: y_r = y_l + 1.0;
1377: color = PetscDrawRealToColor(PetscAbsScalar(v[j*m+i]),minv,maxv);
1378: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1379: }
1380: }
1381: PetscDrawCollectiveEnd(draw);
1382: }
1383: MatDenseRestoreArrayRead(A,&v);
1384: return(0);
1385: }
1387: static PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1388: {
1389: PetscDraw draw;
1390: PetscBool isnull;
1391: PetscReal xr,yr,xl,yl,h,w;
1395: PetscViewerDrawGetDraw(viewer,0,&draw);
1396: PetscDrawIsNull(draw,&isnull);
1397: if (isnull) return(0);
1399: xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1400: xr += w; yr += h; xl = -w; yl = -h;
1401: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
1402: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
1403: PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
1404: PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
1405: PetscDrawSave(draw);
1406: return(0);
1407: }
1409: PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1410: {
1412: PetscBool iascii,isbinary,isdraw;
1415: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1416: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1417: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1419: if (iascii) {
1420: MatView_SeqDense_ASCII(A,viewer);
1421: } else if (isbinary) {
1422: MatView_SeqDense_Binary(A,viewer);
1423: } else if (isdraw) {
1424: MatView_SeqDense_Draw(A,viewer);
1425: }
1426: return(0);
1427: }
1429: static PetscErrorCode MatDensePlaceArray_SeqDense(Mat A,const PetscScalar array[])
1430: {
1431: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1434: a->unplacedarray = a->v;
1435: a->unplaced_user_alloc = a->user_alloc;
1436: a->v = (PetscScalar*) array;
1437: #if defined(PETSC_HAVE_CUDA)
1438: A->offloadmask = PETSC_OFFLOAD_CPU;
1439: #endif
1440: return(0);
1441: }
1443: static PetscErrorCode MatDenseResetArray_SeqDense(Mat A)
1444: {
1445: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1448: a->v = a->unplacedarray;
1449: a->user_alloc = a->unplaced_user_alloc;
1450: a->unplacedarray = NULL;
1451: #if defined(PETSC_HAVE_CUDA)
1452: A->offloadmask = PETSC_OFFLOAD_CPU;
1453: #endif
1454: return(0);
1455: }
1457: PetscErrorCode MatDestroy_SeqDense(Mat mat)
1458: {
1459: Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
1463: #if defined(PETSC_USE_LOG)
1464: PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1465: #endif
1466: PetscFree(l->pivots);
1467: PetscFree(l->fwork);
1468: MatDestroy(&l->ptapwork);
1469: if (!l->user_alloc) {PetscFree(l->v);}
1470: PetscFree(mat->data);
1472: PetscObjectChangeTypeName((PetscObject)mat,0);
1473: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL);
1474: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);
1475: PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);
1476: PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);
1477: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);
1478: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqaij_C",NULL);
1479: #if defined(PETSC_HAVE_ELEMENTAL)
1480: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_elemental_C",NULL);
1481: #endif
1482: #if defined(PETSC_HAVE_CUDA)
1483: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqdensecuda_C",NULL);
1484: PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaijcusparse_seqdense_C",NULL);
1485: #endif
1486: #if defined(PETSC_HAVE_VIENNACL)
1487: PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaijviennacl_seqdense_C",NULL);
1488: #endif
1489: PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL);
1490: PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaij_seqdense_C",NULL);
1491: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C",NULL);
1492: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C",NULL);
1493: PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqbaij_seqdense_C",NULL);
1494: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqbaij_seqdense_C",NULL);
1495: PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqbaij_seqdense_C",NULL);
1496: PetscObjectComposeFunction((PetscObject)mat,"MatPtAP_seqaij_seqdense_C",NULL);
1497: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_seqaij_seqdense_C",NULL);
1498: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",NULL);
1499: PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaij_seqdense_C",NULL);
1500: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);
1501: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);
1502: return(0);
1503: }
1505: static PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1506: {
1507: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1509: PetscInt k,j,m,n,M;
1510: PetscScalar *v,tmp;
1513: m = A->rmap->n; M = mat->lda; n = A->cmap->n;
1514: if (reuse == MAT_INPLACE_MATRIX && m == n) { /* in place transpose */
1515: MatDenseGetArray(A,&v);
1516: for (j=0; j<m; j++) {
1517: for (k=0; k<j; k++) {
1518: tmp = v[j + k*M];
1519: v[j + k*M] = v[k + j*M];
1520: v[k + j*M] = tmp;
1521: }
1522: }
1523: MatDenseRestoreArray(A,&v);
1524: } else { /* out-of-place transpose */
1525: Mat tmat;
1526: Mat_SeqDense *tmatd;
1527: PetscScalar *v2;
1528: PetscInt M2;
1530: if (reuse != MAT_REUSE_MATRIX) {
1531: MatCreate(PetscObjectComm((PetscObject)A),&tmat);
1532: MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);
1533: MatSetType(tmat,((PetscObject)A)->type_name);
1534: MatSeqDenseSetPreallocation(tmat,NULL);
1535: } else tmat = *matout;
1537: MatDenseGetArrayRead(A,(const PetscScalar**)&v);
1538: MatDenseGetArray(tmat,&v2);
1539: tmatd = (Mat_SeqDense*)tmat->data;
1540: M2 = tmatd->lda;
1541: for (j=0; j<n; j++) {
1542: for (k=0; k<m; k++) v2[j + k*M2] = v[k + j*M];
1543: }
1544: MatDenseRestoreArray(tmat,&v2);
1545: MatDenseRestoreArrayRead(A,(const PetscScalar**)&v);
1546: MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
1547: MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
1548: if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) *matout = tmat;
1549: else {
1550: MatHeaderMerge(A,&tmat);
1551: }
1552: }
1553: return(0);
1554: }
1556: static PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg)
1557: {
1558: Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1559: Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1560: PetscInt i;
1561: const PetscScalar *v1,*v2;
1562: PetscErrorCode ierr;
1565: if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; return(0);}
1566: if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; return(0);}
1567: MatDenseGetArrayRead(A1,&v1);
1568: MatDenseGetArrayRead(A2,&v2);
1569: for (i=0; i<A1->cmap->n; i++) {
1570: PetscArraycmp(v1,v2,A1->rmap->n,flg);
1571: if (*flg == PETSC_FALSE) return(0);
1572: v1 += mat1->lda;
1573: v2 += mat2->lda;
1574: }
1575: MatDenseRestoreArrayRead(A1,&v1);
1576: MatDenseRestoreArrayRead(A2,&v2);
1577: *flg = PETSC_TRUE;
1578: return(0);
1579: }
1581: static PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1582: {
1583: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1584: PetscInt i,n,len;
1585: PetscScalar *x;
1586: const PetscScalar *vv;
1587: PetscErrorCode ierr;
1590: VecGetSize(v,&n);
1591: VecGetArray(v,&x);
1592: len = PetscMin(A->rmap->n,A->cmap->n);
1593: MatDenseGetArrayRead(A,&vv);
1594: if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1595: for (i=0; i<len; i++) {
1596: x[i] = vv[i*mat->lda + i];
1597: }
1598: MatDenseRestoreArrayRead(A,&vv);
1599: VecRestoreArray(v,&x);
1600: return(0);
1601: }
1603: static PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1604: {
1605: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1606: const PetscScalar *l,*r;
1607: PetscScalar x,*v,*vv;
1608: PetscErrorCode ierr;
1609: PetscInt i,j,m = A->rmap->n,n = A->cmap->n;
1612: MatDenseGetArray(A,&vv);
1613: if (ll) {
1614: VecGetSize(ll,&m);
1615: VecGetArrayRead(ll,&l);
1616: if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1617: for (i=0; i<m; i++) {
1618: x = l[i];
1619: v = vv + i;
1620: for (j=0; j<n; j++) { (*v) *= x; v+= mat->lda;}
1621: }
1622: VecRestoreArrayRead(ll,&l);
1623: PetscLogFlops(1.0*n*m);
1624: }
1625: if (rr) {
1626: VecGetSize(rr,&n);
1627: VecGetArrayRead(rr,&r);
1628: if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1629: for (i=0; i<n; i++) {
1630: x = r[i];
1631: v = vv + i*mat->lda;
1632: for (j=0; j<m; j++) (*v++) *= x;
1633: }
1634: VecRestoreArrayRead(rr,&r);
1635: PetscLogFlops(1.0*n*m);
1636: }
1637: MatDenseRestoreArray(A,&vv);
1638: return(0);
1639: }
1641: PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1642: {
1643: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1644: PetscScalar *v,*vv;
1645: PetscReal sum = 0.0;
1646: PetscInt lda =mat->lda,m=A->rmap->n,i,j;
1647: PetscErrorCode ierr;
1650: MatDenseGetArrayRead(A,(const PetscScalar**)&vv);
1651: v = vv;
1652: if (type == NORM_FROBENIUS) {
1653: if (lda>m) {
1654: for (j=0; j<A->cmap->n; j++) {
1655: v = vv+j*lda;
1656: for (i=0; i<m; i++) {
1657: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1658: }
1659: }
1660: } else {
1661: #if defined(PETSC_USE_REAL___FP16)
1662: PetscBLASInt one = 1,cnt = A->cmap->n*A->rmap->n;
1663: *nrm = BLASnrm2_(&cnt,v,&one);
1664: }
1665: #else
1666: for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1667: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1668: }
1669: }
1670: *nrm = PetscSqrtReal(sum);
1671: #endif
1672: PetscLogFlops(2.0*A->cmap->n*A->rmap->n);
1673: } else if (type == NORM_1) {
1674: *nrm = 0.0;
1675: for (j=0; j<A->cmap->n; j++) {
1676: v = vv + j*mat->lda;
1677: sum = 0.0;
1678: for (i=0; i<A->rmap->n; i++) {
1679: sum += PetscAbsScalar(*v); v++;
1680: }
1681: if (sum > *nrm) *nrm = sum;
1682: }
1683: PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1684: } else if (type == NORM_INFINITY) {
1685: *nrm = 0.0;
1686: for (j=0; j<A->rmap->n; j++) {
1687: v = vv + j;
1688: sum = 0.0;
1689: for (i=0; i<A->cmap->n; i++) {
1690: sum += PetscAbsScalar(*v); v += mat->lda;
1691: }
1692: if (sum > *nrm) *nrm = sum;
1693: }
1694: PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
1695: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
1696: MatDenseRestoreArrayRead(A,(const PetscScalar**)&vv);
1697: return(0);
1698: }
1700: static PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg)
1701: {
1702: Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
1706: switch (op) {
1707: case MAT_ROW_ORIENTED:
1708: aij->roworiented = flg;
1709: break;
1710: case MAT_NEW_NONZERO_LOCATIONS:
1711: case MAT_NEW_NONZERO_LOCATION_ERR:
1712: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1713: case MAT_NEW_DIAGONALS:
1714: case MAT_KEEP_NONZERO_PATTERN:
1715: case MAT_IGNORE_OFF_PROC_ENTRIES:
1716: case MAT_USE_HASH_TABLE:
1717: case MAT_IGNORE_ZERO_ENTRIES:
1718: case MAT_IGNORE_LOWER_TRIANGULAR:
1719: case MAT_SORTED_FULL:
1720: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1721: break;
1722: case MAT_SPD:
1723: case MAT_SYMMETRIC:
1724: case MAT_STRUCTURALLY_SYMMETRIC:
1725: case MAT_HERMITIAN:
1726: case MAT_SYMMETRY_ETERNAL:
1727: /* These options are handled directly by MatSetOption() */
1728: break;
1729: default:
1730: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
1731: }
1732: return(0);
1733: }
1735: static PetscErrorCode MatZeroEntries_SeqDense(Mat A)
1736: {
1737: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1739: PetscInt lda=l->lda,m=A->rmap->n,j;
1740: PetscScalar *v;
1743: MatDenseGetArray(A,&v);
1744: if (lda>m) {
1745: for (j=0; j<A->cmap->n; j++) {
1746: PetscArrayzero(v+j*lda,m);
1747: }
1748: } else {
1749: PetscArrayzero(v,A->rmap->n*A->cmap->n);
1750: }
1751: MatDenseRestoreArray(A,&v);
1752: return(0);
1753: }
1755: static PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1756: {
1757: PetscErrorCode ierr;
1758: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1759: PetscInt m = l->lda, n = A->cmap->n, i,j;
1760: PetscScalar *slot,*bb,*v;
1761: const PetscScalar *xx;
1764: #if defined(PETSC_USE_DEBUG)
1765: for (i=0; i<N; i++) {
1766: if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
1767: 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);
1768: }
1769: #endif
1770: if (!N) return(0);
1772: /* fix right hand side if needed */
1773: if (x && b) {
1774: VecGetArrayRead(x,&xx);
1775: VecGetArray(b,&bb);
1776: for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
1777: VecRestoreArrayRead(x,&xx);
1778: VecRestoreArray(b,&bb);
1779: }
1781: MatDenseGetArray(A,&v);
1782: for (i=0; i<N; i++) {
1783: slot = v + rows[i];
1784: for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
1785: }
1786: if (diag != 0.0) {
1787: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
1788: for (i=0; i<N; i++) {
1789: slot = v + (m+1)*rows[i];
1790: *slot = diag;
1791: }
1792: }
1793: MatDenseRestoreArray(A,&v);
1794: return(0);
1795: }
1797: static PetscErrorCode MatDenseGetLDA_SeqDense(Mat A,PetscInt *lda)
1798: {
1799: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1802: *lda = mat->lda;
1803: return(0);
1804: }
1806: PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[])
1807: {
1808: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1811: *array = mat->v;
1812: return(0);
1813: }
1815: PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1816: {
1818: return(0);
1819: }
1821: /*@C
1822: MatDenseGetLDA - gets the leading dimension of the array returned from MatDenseGetArray()
1824: Logically Collective on Mat
1826: Input Parameter:
1827: . mat - a MATSEQDENSE or MATMPIDENSE matrix
1829: Output Parameter:
1830: . lda - the leading dimension
1832: Level: intermediate
1834: .seealso: MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatSeqDenseSetLDA()
1835: @*/
1836: PetscErrorCode MatDenseGetLDA(Mat A,PetscInt *lda)
1837: {
1841: PetscUseMethod(A,"MatDenseGetLDA_C",(Mat,PetscInt*),(A,lda));
1842: return(0);
1843: }
1845: /*@C
1846: MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored
1848: Logically Collective on Mat
1850: Input Parameter:
1851: . mat - a MATSEQDENSE or MATMPIDENSE matrix
1853: Output Parameter:
1854: . array - pointer to the data
1856: Level: intermediate
1858: .seealso: MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
1859: @*/
1860: PetscErrorCode MatDenseGetArray(Mat A,PetscScalar **array)
1861: {
1865: PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));
1866: return(0);
1867: }
1869: /*@C
1870: MatDenseRestoreArray - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()
1872: Logically Collective on Mat
1874: Input Parameters:
1875: + mat - a MATSEQDENSE or MATMPIDENSE matrix
1876: - array - pointer to the data
1878: Level: intermediate
1880: .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
1881: @*/
1882: PetscErrorCode MatDenseRestoreArray(Mat A,PetscScalar **array)
1883: {
1887: PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));
1888: if (array) *array = NULL;
1889: PetscObjectStateIncrease((PetscObject)A);
1890: return(0);
1891: }
1893: /*@C
1894: MatDenseGetArrayRead - gives access to the array where the data for a SeqDense matrix is stored
1896: Not Collective
1898: Input Parameter:
1899: . mat - a MATSEQDENSE or MATMPIDENSE matrix
1901: Output Parameter:
1902: . array - pointer to the data
1904: Level: intermediate
1906: .seealso: MatDenseRestoreArray(), MatDenseGetArray(), MatDenseRestoreArrayRead()
1907: @*/
1908: PetscErrorCode MatDenseGetArrayRead(Mat A,const PetscScalar **array)
1909: {
1913: PetscUseMethod(A,"MatDenseGetArrayRead_C",(Mat,const PetscScalar**),(A,array));
1914: return(0);
1915: }
1917: /*@C
1918: MatDenseRestoreArrayRead - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()
1920: Not Collective
1922: Input Parameters:
1923: + mat - a MATSEQDENSE or MATMPIDENSE matrix
1924: - array - pointer to the data
1926: Level: intermediate
1928: .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArray()
1929: @*/
1930: PetscErrorCode MatDenseRestoreArrayRead(Mat A,const PetscScalar **array)
1931: {
1935: PetscUseMethod(A,"MatDenseRestoreArrayRead_C",(Mat,const PetscScalar**),(A,array));
1936: if (array) *array = NULL;
1937: return(0);
1938: }
1940: static PetscErrorCode MatCreateSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
1941: {
1942: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1944: PetscInt i,j,nrows,ncols,blda;
1945: const PetscInt *irow,*icol;
1946: PetscScalar *av,*bv,*v = mat->v;
1947: Mat newmat;
1950: ISGetIndices(isrow,&irow);
1951: ISGetIndices(iscol,&icol);
1952: ISGetLocalSize(isrow,&nrows);
1953: ISGetLocalSize(iscol,&ncols);
1955: /* Check submatrixcall */
1956: if (scall == MAT_REUSE_MATRIX) {
1957: PetscInt n_cols,n_rows;
1958: MatGetSize(*B,&n_rows,&n_cols);
1959: if (n_rows != nrows || n_cols != ncols) {
1960: /* resize the result matrix to match number of requested rows/columns */
1961: MatSetSizes(*B,nrows,ncols,nrows,ncols);
1962: }
1963: newmat = *B;
1964: } else {
1965: /* Create and fill new matrix */
1966: MatCreate(PetscObjectComm((PetscObject)A),&newmat);
1967: MatSetSizes(newmat,nrows,ncols,nrows,ncols);
1968: MatSetType(newmat,((PetscObject)A)->type_name);
1969: MatSeqDenseSetPreallocation(newmat,NULL);
1970: }
1972: /* Now extract the data pointers and do the copy,column at a time */
1973: MatDenseGetArray(newmat,&bv);
1974: MatDenseGetLDA(newmat,&blda);
1975: for (i=0; i<ncols; i++) {
1976: av = v + mat->lda*icol[i];
1977: for (j=0; j<nrows; j++) bv[j] = av[irow[j]];
1978: bv += blda;
1979: }
1980: MatDenseRestoreArray(newmat,&bv);
1982: /* Assemble the matrices so that the correct flags are set */
1983: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1984: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1986: /* Free work space */
1987: ISRestoreIndices(isrow,&irow);
1988: ISRestoreIndices(iscol,&icol);
1989: *B = newmat;
1990: return(0);
1991: }
1993: static PetscErrorCode MatCreateSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1994: {
1996: PetscInt i;
1999: if (scall == MAT_INITIAL_MATRIX) {
2000: PetscCalloc1(n+1,B);
2001: }
2003: for (i=0; i<n; i++) {
2004: MatCreateSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
2005: }
2006: return(0);
2007: }
2009: static PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
2010: {
2012: return(0);
2013: }
2015: static PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
2016: {
2018: return(0);
2019: }
2021: static PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
2022: {
2023: Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data;
2024: PetscErrorCode ierr;
2025: const PetscScalar *va;
2026: PetscScalar *vb;
2027: PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
2030: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
2031: if (A->ops->copy != B->ops->copy) {
2032: MatCopy_Basic(A,B,str);
2033: return(0);
2034: }
2035: if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
2036: MatDenseGetArrayRead(A,&va);
2037: MatDenseGetArray(B,&vb);
2038: if (lda1>m || lda2>m) {
2039: for (j=0; j<n; j++) {
2040: PetscArraycpy(vb+j*lda2,va+j*lda1,m);
2041: }
2042: } else {
2043: PetscArraycpy(vb,va,A->rmap->n*A->cmap->n);
2044: }
2045: MatDenseRestoreArray(B,&vb);
2046: MatDenseRestoreArrayRead(A,&va);
2047: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2048: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2049: return(0);
2050: }
2052: static PetscErrorCode MatSetUp_SeqDense(Mat A)
2053: {
2057: MatSeqDenseSetPreallocation(A,0);
2058: return(0);
2059: }
2061: static PetscErrorCode MatConjugate_SeqDense(Mat A)
2062: {
2063: PetscInt i,nz = A->rmap->n*A->cmap->n;
2064: PetscScalar *aa;
2068: MatDenseGetArray(A,&aa);
2069: for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
2070: MatDenseRestoreArray(A,&aa);
2071: return(0);
2072: }
2074: static PetscErrorCode MatRealPart_SeqDense(Mat A)
2075: {
2076: PetscInt i,nz = A->rmap->n*A->cmap->n;
2077: PetscScalar *aa;
2081: MatDenseGetArray(A,&aa);
2082: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
2083: MatDenseRestoreArray(A,&aa);
2084: return(0);
2085: }
2087: static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
2088: {
2089: PetscInt i,nz = A->rmap->n*A->cmap->n;
2090: PetscScalar *aa;
2094: MatDenseGetArray(A,&aa);
2095: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
2096: MatDenseRestoreArray(A,&aa);
2097: return(0);
2098: }
2100: /* ----------------------------------------------------------------*/
2101: PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2102: {
2106: if (scall == MAT_INITIAL_MATRIX) {
2107: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
2108: MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
2109: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
2110: }
2111: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
2112: if ((*C)->ops->matmultnumeric) {
2113: (*(*C)->ops->matmultnumeric)(A,B,*C);
2114: } else {
2115: MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);
2116: }
2117: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
2118: return(0);
2119: }
2121: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2122: {
2124: PetscInt m=A->rmap->n,n=B->cmap->n;
2125: Mat Cmat;
2126: PetscBool flg;
2129: MatCreate(PETSC_COMM_SELF,&Cmat);
2130: MatSetSizes(Cmat,m,n,m,n);
2131: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2132: MatSetType(Cmat,flg ? ((PetscObject)A)->type_name : MATDENSE);
2133: MatSeqDenseSetPreallocation(Cmat,NULL);
2134: *C = Cmat;
2135: return(0);
2136: }
2138: PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2139: {
2140: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2141: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2142: Mat_SeqDense *c = (Mat_SeqDense*)C->data;
2143: PetscBLASInt m,n,k;
2144: const PetscScalar *av,*bv;
2145: PetscScalar *cv;
2146: PetscScalar _DOne=1.0,_DZero=0.0;
2147: PetscErrorCode ierr;
2148: PetscBool flg;
2151: /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
2152: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
2153: if (flg) {
2154: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
2155: (*C->ops->matmultnumeric)(A,B,C);
2156: return(0);
2157: }
2158: PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&flg);
2159: if (flg) {
2160: C->ops->matmultnumeric = MatMatMultNumeric_SeqBAIJ_SeqDense;
2161: (*C->ops->matmultnumeric)(A,B,C);
2162: return(0);
2163: }
2164: PetscBLASIntCast(C->rmap->n,&m);
2165: PetscBLASIntCast(C->cmap->n,&n);
2166: PetscBLASIntCast(A->cmap->n,&k);
2167: if (!m || !n || !k) return(0);
2168: MatDenseGetArrayRead(A,&av);
2169: MatDenseGetArrayRead(B,&bv);
2170: MatDenseGetArray(C,&cv);
2171: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda));
2172: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2173: MatDenseRestoreArrayRead(A,&av);
2174: MatDenseRestoreArrayRead(B,&bv);
2175: MatDenseRestoreArray(C,&cv);
2176: return(0);
2177: }
2179: PetscErrorCode MatMatTransposeMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2180: {
2184: if (scall == MAT_INITIAL_MATRIX) {
2185: PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);
2186: MatMatTransposeMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
2187: PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);
2188: }
2189: PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);
2190: MatMatTransposeMultNumeric_SeqDense_SeqDense(A,B,*C);
2191: PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);
2192: return(0);
2193: }
2195: PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2196: {
2198: PetscInt m=A->rmap->n,n=B->rmap->n;
2199: Mat Cmat;
2200: PetscBool flg;
2203: MatCreate(PETSC_COMM_SELF,&Cmat);
2204: MatSetSizes(Cmat,m,n,m,n);
2205: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2206: MatSetType(Cmat,flg ? ((PetscObject)A)->type_name : MATDENSE);
2207: MatSeqDenseSetPreallocation(Cmat,NULL);
2208: *C = Cmat;
2209: return(0);
2210: }
2212: PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2213: {
2214: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2215: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2216: Mat_SeqDense *c = (Mat_SeqDense*)C->data;
2217: PetscBLASInt m,n,k;
2218: PetscScalar _DOne=1.0,_DZero=0.0;
2222: PetscBLASIntCast(C->rmap->n,&m);
2223: PetscBLASIntCast(C->cmap->n,&n);
2224: PetscBLASIntCast(A->cmap->n,&k);
2225: if (!m || !n || !k) return(0);
2226: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2227: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2228: return(0);
2229: }
2231: PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2232: {
2236: if (scall == MAT_INITIAL_MATRIX) {
2237: PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
2238: MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);
2239: PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
2240: }
2241: PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
2242: MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);
2243: PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
2244: return(0);
2245: }
2247: PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2248: {
2250: PetscInt m=A->cmap->n,n=B->cmap->n;
2251: Mat Cmat;
2252: PetscBool flg;
2255: MatCreate(PETSC_COMM_SELF,&Cmat);
2256: MatSetSizes(Cmat,m,n,m,n);
2257: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2258: MatSetType(Cmat,flg ? ((PetscObject)A)->type_name : MATDENSE);
2259: MatSeqDenseSetPreallocation(Cmat,NULL);
2260: *C = Cmat;
2261: return(0);
2262: }
2264: PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2265: {
2266: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2267: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2268: Mat_SeqDense *c = (Mat_SeqDense*)C->data;
2269: PetscBLASInt m,n,k;
2270: PetscScalar _DOne=1.0,_DZero=0.0;
2274: PetscBLASIntCast(C->rmap->n,&m);
2275: PetscBLASIntCast(C->cmap->n,&n);
2276: PetscBLASIntCast(A->rmap->n,&k);
2277: if (!m || !n || !k) return(0);
2278: PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2279: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2280: return(0);
2281: }
2283: static PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
2284: {
2285: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2286: PetscErrorCode ierr;
2287: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2288: PetscScalar *x;
2289: const PetscScalar *aa;
2292: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2293: VecGetArray(v,&x);
2294: VecGetLocalSize(v,&p);
2295: MatDenseGetArrayRead(A,&aa);
2296: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2297: for (i=0; i<m; i++) {
2298: x[i] = aa[i]; if (idx) idx[i] = 0;
2299: for (j=1; j<n; j++) {
2300: if (PetscRealPart(x[i]) < PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;}
2301: }
2302: }
2303: MatDenseRestoreArrayRead(A,&aa);
2304: VecRestoreArray(v,&x);
2305: return(0);
2306: }
2308: static PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
2309: {
2310: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2311: PetscErrorCode ierr;
2312: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2313: PetscScalar *x;
2314: PetscReal atmp;
2315: const PetscScalar *aa;
2318: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2319: VecGetArray(v,&x);
2320: VecGetLocalSize(v,&p);
2321: MatDenseGetArrayRead(A,&aa);
2322: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2323: for (i=0; i<m; i++) {
2324: x[i] = PetscAbsScalar(aa[i]);
2325: for (j=1; j<n; j++) {
2326: atmp = PetscAbsScalar(aa[i+a->lda*j]);
2327: if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
2328: }
2329: }
2330: MatDenseRestoreArrayRead(A,&aa);
2331: VecRestoreArray(v,&x);
2332: return(0);
2333: }
2335: static PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
2336: {
2337: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2338: PetscErrorCode ierr;
2339: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2340: PetscScalar *x;
2341: const PetscScalar *aa;
2344: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2345: MatDenseGetArrayRead(A,&aa);
2346: VecGetArray(v,&x);
2347: VecGetLocalSize(v,&p);
2348: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2349: for (i=0; i<m; i++) {
2350: x[i] = aa[i]; if (idx) idx[i] = 0;
2351: for (j=1; j<n; j++) {
2352: if (PetscRealPart(x[i]) > PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;}
2353: }
2354: }
2355: VecRestoreArray(v,&x);
2356: MatDenseRestoreArrayRead(A,&aa);
2357: return(0);
2358: }
2360: static PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
2361: {
2362: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2363: PetscErrorCode ierr;
2364: PetscScalar *x;
2365: const PetscScalar *aa;
2368: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2369: MatDenseGetArrayRead(A,&aa);
2370: VecGetArray(v,&x);
2371: PetscArraycpy(x,aa+col*a->lda,A->rmap->n);
2372: VecRestoreArray(v,&x);
2373: MatDenseRestoreArrayRead(A,&aa);
2374: return(0);
2375: }
2377: PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
2378: {
2379: PetscErrorCode ierr;
2380: PetscInt i,j,m,n;
2381: const PetscScalar *a;
2384: MatGetSize(A,&m,&n);
2385: PetscArrayzero(norms,n);
2386: MatDenseGetArrayRead(A,&a);
2387: if (type == NORM_2) {
2388: for (i=0; i<n; i++) {
2389: for (j=0; j<m; j++) {
2390: norms[i] += PetscAbsScalar(a[j]*a[j]);
2391: }
2392: a += m;
2393: }
2394: } else if (type == NORM_1) {
2395: for (i=0; i<n; i++) {
2396: for (j=0; j<m; j++) {
2397: norms[i] += PetscAbsScalar(a[j]);
2398: }
2399: a += m;
2400: }
2401: } else if (type == NORM_INFINITY) {
2402: for (i=0; i<n; i++) {
2403: for (j=0; j<m; j++) {
2404: norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
2405: }
2406: a += m;
2407: }
2408: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
2409: MatDenseRestoreArrayRead(A,&a);
2410: if (type == NORM_2) {
2411: for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
2412: }
2413: return(0);
2414: }
2416: static PetscErrorCode MatSetRandom_SeqDense(Mat x,PetscRandom rctx)
2417: {
2419: PetscScalar *a;
2420: PetscInt m,n,i;
2423: MatGetSize(x,&m,&n);
2424: MatDenseGetArray(x,&a);
2425: for (i=0; i<m*n; i++) {
2426: PetscRandomGetValue(rctx,a+i);
2427: }
2428: MatDenseRestoreArray(x,&a);
2429: return(0);
2430: }
2432: static PetscErrorCode MatMissingDiagonal_SeqDense(Mat A,PetscBool *missing,PetscInt *d)
2433: {
2435: *missing = PETSC_FALSE;
2436: return(0);
2437: }
2439: /* vals is not const */
2440: static PetscErrorCode MatDenseGetColumn_SeqDense(Mat A,PetscInt col,PetscScalar **vals)
2441: {
2443: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2444: PetscScalar *v;
2447: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2448: MatDenseGetArray(A,&v);
2449: *vals = v+col*a->lda;
2450: MatDenseRestoreArray(A,&v);
2451: return(0);
2452: }
2454: static PetscErrorCode MatDenseRestoreColumn_SeqDense(Mat A,PetscScalar **vals)
2455: {
2457: *vals = 0; /* user cannot accidently use the array later */
2458: return(0);
2459: }
2461: /* -------------------------------------------------------------------*/
2462: static struct _MatOps MatOps_Values = { MatSetValues_SeqDense,
2463: MatGetRow_SeqDense,
2464: MatRestoreRow_SeqDense,
2465: MatMult_SeqDense,
2466: /* 4*/ MatMultAdd_SeqDense,
2467: MatMultTranspose_SeqDense,
2468: MatMultTransposeAdd_SeqDense,
2469: 0,
2470: 0,
2471: 0,
2472: /* 10*/ 0,
2473: MatLUFactor_SeqDense,
2474: MatCholeskyFactor_SeqDense,
2475: MatSOR_SeqDense,
2476: MatTranspose_SeqDense,
2477: /* 15*/ MatGetInfo_SeqDense,
2478: MatEqual_SeqDense,
2479: MatGetDiagonal_SeqDense,
2480: MatDiagonalScale_SeqDense,
2481: MatNorm_SeqDense,
2482: /* 20*/ MatAssemblyBegin_SeqDense,
2483: MatAssemblyEnd_SeqDense,
2484: MatSetOption_SeqDense,
2485: MatZeroEntries_SeqDense,
2486: /* 24*/ MatZeroRows_SeqDense,
2487: 0,
2488: 0,
2489: 0,
2490: 0,
2491: /* 29*/ MatSetUp_SeqDense,
2492: 0,
2493: 0,
2494: 0,
2495: 0,
2496: /* 34*/ MatDuplicate_SeqDense,
2497: 0,
2498: 0,
2499: 0,
2500: 0,
2501: /* 39*/ MatAXPY_SeqDense,
2502: MatCreateSubMatrices_SeqDense,
2503: 0,
2504: MatGetValues_SeqDense,
2505: MatCopy_SeqDense,
2506: /* 44*/ MatGetRowMax_SeqDense,
2507: MatScale_SeqDense,
2508: MatShift_Basic,
2509: 0,
2510: MatZeroRowsColumns_SeqDense,
2511: /* 49*/ MatSetRandom_SeqDense,
2512: 0,
2513: 0,
2514: 0,
2515: 0,
2516: /* 54*/ 0,
2517: 0,
2518: 0,
2519: 0,
2520: 0,
2521: /* 59*/ 0,
2522: MatDestroy_SeqDense,
2523: MatView_SeqDense,
2524: 0,
2525: 0,
2526: /* 64*/ 0,
2527: 0,
2528: 0,
2529: 0,
2530: 0,
2531: /* 69*/ MatGetRowMaxAbs_SeqDense,
2532: 0,
2533: 0,
2534: 0,
2535: 0,
2536: /* 74*/ 0,
2537: 0,
2538: 0,
2539: 0,
2540: 0,
2541: /* 79*/ 0,
2542: 0,
2543: 0,
2544: 0,
2545: /* 83*/ MatLoad_SeqDense,
2546: 0,
2547: MatIsHermitian_SeqDense,
2548: 0,
2549: 0,
2550: 0,
2551: /* 89*/ MatMatMult_SeqDense_SeqDense,
2552: MatMatMultSymbolic_SeqDense_SeqDense,
2553: MatMatMultNumeric_SeqDense_SeqDense,
2554: MatPtAP_SeqDense_SeqDense,
2555: MatPtAPSymbolic_SeqDense_SeqDense,
2556: /* 94*/ MatPtAPNumeric_SeqDense_SeqDense,
2557: MatMatTransposeMult_SeqDense_SeqDense,
2558: MatMatTransposeMultSymbolic_SeqDense_SeqDense,
2559: MatMatTransposeMultNumeric_SeqDense_SeqDense,
2560: 0,
2561: /* 99*/ 0,
2562: 0,
2563: 0,
2564: MatConjugate_SeqDense,
2565: 0,
2566: /*104*/ 0,
2567: MatRealPart_SeqDense,
2568: MatImaginaryPart_SeqDense,
2569: 0,
2570: 0,
2571: /*109*/ 0,
2572: 0,
2573: MatGetRowMin_SeqDense,
2574: MatGetColumnVector_SeqDense,
2575: MatMissingDiagonal_SeqDense,
2576: /*114*/ 0,
2577: 0,
2578: 0,
2579: 0,
2580: 0,
2581: /*119*/ 0,
2582: 0,
2583: 0,
2584: 0,
2585: 0,
2586: /*124*/ 0,
2587: MatGetColumnNorms_SeqDense,
2588: 0,
2589: 0,
2590: 0,
2591: /*129*/ 0,
2592: MatTransposeMatMult_SeqDense_SeqDense,
2593: MatTransposeMatMultSymbolic_SeqDense_SeqDense,
2594: MatTransposeMatMultNumeric_SeqDense_SeqDense,
2595: 0,
2596: /*134*/ 0,
2597: 0,
2598: 0,
2599: 0,
2600: 0,
2601: /*139*/ 0,
2602: 0,
2603: 0,
2604: 0,
2605: 0,
2606: /*144*/ MatCreateMPIMatConcatenateSeqMat_SeqDense
2607: };
2609: /*@C
2610: MatCreateSeqDense - Creates a sequential dense matrix that
2611: is stored in column major order (the usual Fortran 77 manner). Many
2612: of the matrix operations use the BLAS and LAPACK routines.
2614: Collective
2616: Input Parameters:
2617: + comm - MPI communicator, set to PETSC_COMM_SELF
2618: . m - number of rows
2619: . n - number of columns
2620: - data - optional location of matrix data in column major order. Set data=NULL for PETSc
2621: to control all matrix memory allocation.
2623: Output Parameter:
2624: . A - the matrix
2626: Notes:
2627: The data input variable is intended primarily for Fortran programmers
2628: who wish to allocate their own matrix memory space. Most users should
2629: set data=NULL.
2631: Level: intermediate
2633: .seealso: MatCreate(), MatCreateDense(), MatSetValues()
2634: @*/
2635: PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2636: {
2640: MatCreate(comm,A);
2641: MatSetSizes(*A,m,n,m,n);
2642: MatSetType(*A,MATSEQDENSE);
2643: MatSeqDenseSetPreallocation(*A,data);
2644: return(0);
2645: }
2647: /*@C
2648: MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
2650: Collective
2652: Input Parameters:
2653: + B - the matrix
2654: - data - the array (or NULL)
2656: Notes:
2657: The data input variable is intended primarily for Fortran programmers
2658: who wish to allocate their own matrix memory space. Most users should
2659: need not call this routine.
2661: Level: intermediate
2663: .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA()
2665: @*/
2666: PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2667: {
2671: PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));
2672: return(0);
2673: }
2675: PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2676: {
2677: Mat_SeqDense *b;
2681: B->preallocated = PETSC_TRUE;
2683: PetscLayoutSetUp(B->rmap);
2684: PetscLayoutSetUp(B->cmap);
2686: b = (Mat_SeqDense*)B->data;
2687: b->Mmax = B->rmap->n;
2688: b->Nmax = B->cmap->n;
2689: if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n;
2691: PetscIntMultError(b->lda,b->Nmax,NULL);
2692: if (!data) { /* petsc-allocated storage */
2693: if (!b->user_alloc) { PetscFree(b->v); }
2694: PetscCalloc1((size_t)b->lda*b->Nmax,&b->v);
2695: PetscLogObjectMemory((PetscObject)B,b->lda*b->Nmax*sizeof(PetscScalar));
2697: b->user_alloc = PETSC_FALSE;
2698: } else { /* user-allocated storage */
2699: if (!b->user_alloc) { PetscFree(b->v); }
2700: b->v = data;
2701: b->user_alloc = PETSC_TRUE;
2702: }
2703: B->assembled = PETSC_TRUE;
2704: return(0);
2705: }
2707: #if defined(PETSC_HAVE_ELEMENTAL)
2708: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
2709: {
2710: Mat mat_elemental;
2711: PetscErrorCode ierr;
2712: const PetscScalar *array;
2713: PetscScalar *v_colwise;
2714: PetscInt M=A->rmap->N,N=A->cmap->N,i,j,k,*rows,*cols;
2717: PetscMalloc3(M*N,&v_colwise,M,&rows,N,&cols);
2718: MatDenseGetArrayRead(A,&array);
2719: /* convert column-wise array into row-wise v_colwise, see MatSetValues_Elemental() */
2720: k = 0;
2721: for (j=0; j<N; j++) {
2722: cols[j] = j;
2723: for (i=0; i<M; i++) {
2724: v_colwise[j*M+i] = array[k++];
2725: }
2726: }
2727: for (i=0; i<M; i++) {
2728: rows[i] = i;
2729: }
2730: MatDenseRestoreArrayRead(A,&array);
2732: MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
2733: MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
2734: MatSetType(mat_elemental,MATELEMENTAL);
2735: MatSetUp(mat_elemental);
2737: /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
2738: MatSetValues(mat_elemental,M,rows,N,cols,v_colwise,ADD_VALUES);
2739: MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
2740: MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
2741: PetscFree3(v_colwise,rows,cols);
2743: if (reuse == MAT_INPLACE_MATRIX) {
2744: MatHeaderReplace(A,&mat_elemental);
2745: } else {
2746: *newmat = mat_elemental;
2747: }
2748: return(0);
2749: }
2750: #endif
2752: /*@C
2753: MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
2755: Input parameter:
2756: + A - the matrix
2757: - lda - the leading dimension
2759: Notes:
2760: This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
2761: it asserts that the preallocation has a leading dimension (the LDA parameter
2762: of Blas and Lapack fame) larger than M, the first dimension of the matrix.
2764: Level: intermediate
2766: .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2768: @*/
2769: PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda)
2770: {
2771: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2774: 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);
2775: b->lda = lda;
2776: b->changelda = PETSC_FALSE;
2777: b->Mmax = PetscMax(b->Mmax,lda);
2778: return(0);
2779: }
2781: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2782: {
2784: PetscMPIInt size;
2787: MPI_Comm_size(comm,&size);
2788: if (size == 1) {
2789: if (scall == MAT_INITIAL_MATRIX) {
2790: MatDuplicate(inmat,MAT_COPY_VALUES,outmat);
2791: } else {
2792: MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);
2793: }
2794: } else {
2795: MatCreateMPIMatConcatenateSeqMat_MPIDense(comm,inmat,n,scall,outmat);
2796: }
2797: return(0);
2798: }
2800: /*MC
2801: MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
2803: Options Database Keys:
2804: . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
2806: Level: beginner
2808: .seealso: MatCreateSeqDense()
2810: M*/
2811: PetscErrorCode MatCreate_SeqDense(Mat B)
2812: {
2813: Mat_SeqDense *b;
2815: PetscMPIInt size;
2818: MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
2819: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
2821: PetscNewLog(B,&b);
2822: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2823: B->data = (void*)b;
2825: b->roworiented = PETSC_TRUE;
2827: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetLDA_C",MatDenseGetLDA_SeqDense);
2828: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense);
2829: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense);
2830: PetscObjectComposeFunction((PetscObject)B,"MatDensePlaceArray_C",MatDensePlaceArray_SeqDense);
2831: PetscObjectComposeFunction((PetscObject)B,"MatDenseResetArray_C",MatDenseResetArray_SeqDense);
2832: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayRead_C",MatDenseGetArray_SeqDense);
2833: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayRead_C",MatDenseRestoreArray_SeqDense);
2834: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ);
2835: #if defined(PETSC_HAVE_ELEMENTAL)
2836: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_elemental_C",MatConvert_SeqDense_Elemental);
2837: #endif
2838: #if defined(PETSC_HAVE_CUDA)
2839: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqdensecuda_C",MatConvert_SeqDense_SeqDenseCUDA);
2840: PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijcusparse_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2841: #endif
2842: PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense);
2843: PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaij_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2844: #if defined(PETSC_HAVE_VIENNACL)
2845: PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijviennacl_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2846: #endif
2847: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2848: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2849: PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqbaij_seqdense_C",MatMatMult_SeqBAIJ_SeqDense);
2850: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqbaij_seqdense_C",MatMatMultSymbolic_SeqBAIJ_SeqDense);
2851: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqbaij_seqdense_C",MatMatMultNumeric_SeqBAIJ_SeqDense);
2852: PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_seqdense_C",MatPtAP_SeqDense_SeqDense);
2853: PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijperm_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2854: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijperm_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2855: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijperm_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2856: PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijperm_seqdense_C",MatPtAP_SeqDense_SeqDense);
2857: PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijsell_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2858: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijsell_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2859: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijsell_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2860: PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijsell_seqdense_C",MatPtAP_SeqDense_SeqDense);
2861: PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijmkl_seqdense_C",MatMatMult_SeqAIJ_SeqDense);
2862: PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijmkl_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);
2863: PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijmkl_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);
2864: PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijmkl_seqdense_C",MatPtAP_SeqDense_SeqDense);
2866: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaij_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2867: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2868: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaij_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2869: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijperm_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2870: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijperm_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2871: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijperm_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2872: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijsell_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2873: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijsell_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2874: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijsell_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2876: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijmkl_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);
2877: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijmkl_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);
2878: PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijmkl_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);
2879: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumn_C",MatDenseGetColumn_SeqDense);
2880: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_SeqDense);
2881: PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);
2882: return(0);
2883: }
2885: /*@C
2886: 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.
2888: Not Collective
2890: Input Parameter:
2891: + mat - a MATSEQDENSE or MATMPIDENSE matrix
2892: - col - column index
2894: Output Parameter:
2895: . vals - pointer to the data
2897: Level: intermediate
2899: .seealso: MatDenseRestoreColumn()
2900: @*/
2901: PetscErrorCode MatDenseGetColumn(Mat A,PetscInt col,PetscScalar **vals)
2902: {
2906: PetscUseMethod(A,"MatDenseGetColumn_C",(Mat,PetscInt,PetscScalar**),(A,col,vals));
2907: return(0);
2908: }
2910: /*@C
2911: MatDenseRestoreColumn - returns access to a column of a dense matrix which is returned by MatDenseGetColumn().
2913: Not Collective
2915: Input Parameter:
2916: . mat - a MATSEQDENSE or MATMPIDENSE matrix
2918: Output Parameter:
2919: . vals - pointer to the data
2921: Level: intermediate
2923: .seealso: MatDenseGetColumn()
2924: @*/
2925: PetscErrorCode MatDenseRestoreColumn(Mat A,PetscScalar **vals)
2926: {
2930: PetscUseMethod(A,"MatDenseRestoreColumn_C",(Mat,PetscScalar**),(A,vals));
2931: return(0);
2932: }