Actual source code: dense.c
petsc-3.15.0 2021-03-30
2: /*
3: Defines the basic matrix operations for sequential dense.
4: */
6: #include <../src/mat/impls/dense/seq/dense.h>
7: #include <petscblaslapack.h>
9: #include <../src/mat/impls/aij/seq/aij.h>
11: PetscErrorCode MatSeqDenseSymmetrize_Private(Mat A, PetscBool hermitian)
12: {
13: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
14: PetscInt j, k, n = A->rmap->n;
15: PetscScalar *v;
19: if (A->rmap->n != A->cmap->n) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot symmetrize a rectangular matrix");
20: MatDenseGetArray(A,&v);
21: if (!hermitian) {
22: for (k=0;k<n;k++) {
23: for (j=k;j<n;j++) {
24: v[j*mat->lda + k] = v[k*mat->lda + j];
25: }
26: }
27: } else {
28: for (k=0;k<n;k++) {
29: for (j=k;j<n;j++) {
30: v[j*mat->lda + k] = PetscConj(v[k*mat->lda + j]);
31: }
32: }
33: }
34: MatDenseRestoreArray(A,&v);
35: return(0);
36: }
38: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat A)
39: {
40: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
42: PetscBLASInt info,n;
45: if (!A->rmap->n || !A->cmap->n) return(0);
46: PetscBLASIntCast(A->cmap->n,&n);
47: if (A->factortype == MAT_FACTOR_LU) {
48: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
49: if (!mat->fwork) {
50: mat->lfwork = n;
51: PetscMalloc1(mat->lfwork,&mat->fwork);
52: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
53: }
54: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
55: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
56: PetscFPTrapPop();
57: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
58: } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
59: if (A->spd) {
60: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
61: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&n,mat->v,&mat->lda,&info));
62: PetscFPTrapPop();
63: MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
64: #if defined(PETSC_USE_COMPLEX)
65: } else if (A->hermitian) {
66: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
67: if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
68: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
69: PetscStackCallBLAS("LAPACKhetri",LAPACKhetri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
70: PetscFPTrapPop();
71: MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);
72: #endif
73: } else { /* symmetric case */
74: if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
75: if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
76: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
77: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
78: PetscFPTrapPop();
79: MatSeqDenseSymmetrize_Private(A,PETSC_FALSE);
80: }
81: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad Inversion: zero pivot in row %D",(PetscInt)info-1);
82: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
83: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
85: A->ops->solve = NULL;
86: A->ops->matsolve = NULL;
87: A->ops->solvetranspose = NULL;
88: A->ops->matsolvetranspose = NULL;
89: A->ops->solveadd = NULL;
90: A->ops->solvetransposeadd = NULL;
91: A->factortype = MAT_FACTOR_NONE;
92: PetscFree(A->solvertype);
93: return(0);
94: }
96: PetscErrorCode MatZeroRowsColumns_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
97: {
98: PetscErrorCode ierr;
99: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
100: PetscInt m = l->lda, n = A->cmap->n,r = A->rmap->n, i,j;
101: PetscScalar *slot,*bb,*v;
102: const PetscScalar *xx;
105: if (PetscDefined(USE_DEBUG)) {
106: for (i=0; i<N; i++) {
107: if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
108: if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n);
109: if (rows[i] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Col %D requested to be zeroed greater than or equal number of cols %D",rows[i],A->cmap->n);
110: }
111: }
112: if (!N) return(0);
114: /* fix right hand side if needed */
115: if (x && b) {
116: Vec xt;
118: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
119: VecDuplicate(x,&xt);
120: VecCopy(x,xt);
121: VecScale(xt,-1.0);
122: MatMultAdd(A,xt,b,b);
123: VecDestroy(&xt);
124: VecGetArrayRead(x,&xx);
125: VecGetArray(b,&bb);
126: for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
127: VecRestoreArrayRead(x,&xx);
128: VecRestoreArray(b,&bb);
129: }
131: MatDenseGetArray(A,&v);
132: for (i=0; i<N; i++) {
133: slot = v + rows[i]*m;
134: PetscArrayzero(slot,r);
135: }
136: for (i=0; i<N; i++) {
137: slot = v + rows[i];
138: for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
139: }
140: if (diag != 0.0) {
141: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
142: for (i=0; i<N; i++) {
143: slot = v + (m+1)*rows[i];
144: *slot = diag;
145: }
146: }
147: MatDenseRestoreArray(A,&v);
148: return(0);
149: }
151: PetscErrorCode MatPtAPNumeric_SeqDense_SeqDense(Mat A,Mat P,Mat C)
152: {
153: Mat_SeqDense *c = (Mat_SeqDense*)(C->data);
157: if (c->ptapwork) {
158: (*C->ops->matmultnumeric)(A,P,c->ptapwork);
159: (*C->ops->transposematmultnumeric)(P,c->ptapwork,C);
160: } else SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Must call MatPtAPSymbolic_SeqDense_SeqDense() first");
161: return(0);
162: }
164: PetscErrorCode MatPtAPSymbolic_SeqDense_SeqDense(Mat A,Mat P,PetscReal fill,Mat C)
165: {
166: Mat_SeqDense *c;
167: PetscBool cisdense;
171: MatSetSizes(C,P->cmap->n,P->cmap->n,P->cmap->N,P->cmap->N);
172: PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"");
173: if (!cisdense) {
174: PetscBool flg;
176: PetscObjectTypeCompare((PetscObject)P,((PetscObject)A)->type_name,&flg);
177: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
178: }
179: MatSetUp(C);
180: c = (Mat_SeqDense*)C->data;
181: MatCreate(PetscObjectComm((PetscObject)A),&c->ptapwork);
182: MatSetSizes(c->ptapwork,A->rmap->n,P->cmap->n,A->rmap->N,P->cmap->N);
183: MatSetType(c->ptapwork,((PetscObject)C)->type_name);
184: MatSetUp(c->ptapwork);
185: return(0);
186: }
188: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
189: {
190: Mat B = NULL;
191: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
192: Mat_SeqDense *b;
193: PetscErrorCode ierr;
194: PetscInt *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->N,i;
195: const MatScalar *av;
196: PetscBool isseqdense;
199: if (reuse == MAT_REUSE_MATRIX) {
200: PetscObjectTypeCompare((PetscObject)*newmat,MATSEQDENSE,&isseqdense);
201: if (!isseqdense) SETERRQ1(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type_name);
202: }
203: if (reuse != MAT_REUSE_MATRIX) {
204: MatCreate(PetscObjectComm((PetscObject)A),&B);
205: MatSetSizes(B,m,n,m,n);
206: MatSetType(B,MATSEQDENSE);
207: MatSeqDenseSetPreallocation(B,NULL);
208: b = (Mat_SeqDense*)(B->data);
209: } else {
210: b = (Mat_SeqDense*)((*newmat)->data);
211: PetscArrayzero(b->v,m*n);
212: }
213: MatSeqAIJGetArrayRead(A,&av);
214: for (i=0; i<m; i++) {
215: PetscInt j;
216: for (j=0;j<ai[1]-ai[0];j++) {
217: b->v[*aj*m+i] = *av;
218: aj++;
219: av++;
220: }
221: ai++;
222: }
223: MatSeqAIJRestoreArrayRead(A,&av);
225: if (reuse == MAT_INPLACE_MATRIX) {
226: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
227: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
228: MatHeaderReplace(A,&B);
229: } else {
230: if (B) *newmat = B;
231: MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
232: MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
233: }
234: return(0);
235: }
237: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
238: {
239: Mat B = NULL;
240: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
242: PetscInt i, j;
243: PetscInt *rows, *nnz;
244: MatScalar *aa = a->v, *vals;
247: PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals);
248: if (reuse != MAT_REUSE_MATRIX) {
249: MatCreate(PetscObjectComm((PetscObject)A),&B);
250: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
251: MatSetType(B,MATSEQAIJ);
252: for (j=0; j<A->cmap->n; j++) {
253: for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || (i == j && A->cmap->n == A->rmap->n)) ++nnz[i];
254: aa += a->lda;
255: }
256: MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz);
257: } else B = *newmat;
258: aa = a->v;
259: for (j=0; j<A->cmap->n; j++) {
260: PetscInt numRows = 0;
261: for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || (i == j && A->cmap->n == A->rmap->n)) {rows[numRows] = i; vals[numRows++] = aa[i];}
262: MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES);
263: aa += a->lda;
264: }
265: PetscFree3(rows,nnz,vals);
266: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
267: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
269: if (reuse == MAT_INPLACE_MATRIX) {
270: MatHeaderReplace(A,&B);
271: } else if (reuse != MAT_REUSE_MATRIX) *newmat = B;
272: return(0);
273: }
275: PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
276: {
277: Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
278: const PetscScalar *xv;
279: PetscScalar *yv;
280: PetscBLASInt N,m,ldax = 0,lday = 0,one = 1;
281: PetscErrorCode ierr;
284: MatDenseGetArrayRead(X,&xv);
285: MatDenseGetArray(Y,&yv);
286: PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);
287: PetscBLASIntCast(X->rmap->n,&m);
288: PetscBLASIntCast(x->lda,&ldax);
289: PetscBLASIntCast(y->lda,&lday);
290: if (ldax>m || lday>m) {
291: PetscInt j;
293: for (j=0; j<X->cmap->n; j++) {
294: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&alpha,xv+j*ldax,&one,yv+j*lday,&one));
295: }
296: } else {
297: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&alpha,xv,&one,yv,&one));
298: }
299: MatDenseRestoreArrayRead(X,&xv);
300: MatDenseRestoreArray(Y,&yv);
301: PetscLogFlops(PetscMax(2.0*N-1,0));
302: return(0);
303: }
305: static PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
306: {
307: PetscLogDouble N = A->rmap->n*A->cmap->n;
310: info->block_size = 1.0;
311: info->nz_allocated = N;
312: info->nz_used = N;
313: info->nz_unneeded = 0;
314: info->assemblies = A->num_ass;
315: info->mallocs = 0;
316: info->memory = ((PetscObject)A)->mem;
317: info->fill_ratio_given = 0;
318: info->fill_ratio_needed = 0;
319: info->factor_mallocs = 0;
320: return(0);
321: }
323: PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
324: {
325: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
326: PetscScalar *v;
328: PetscBLASInt one = 1,j,nz,lda = 0;
331: MatDenseGetArray(A,&v);
332: PetscBLASIntCast(a->lda,&lda);
333: if (lda>A->rmap->n) {
334: PetscBLASIntCast(A->rmap->n,&nz);
335: for (j=0; j<A->cmap->n; j++) {
336: PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v+j*lda,&one));
337: }
338: } else {
339: PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);
340: PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v,&one));
341: }
342: PetscLogFlops(nz);
343: MatDenseRestoreArray(A,&v);
344: return(0);
345: }
347: static PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl)
348: {
349: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
350: PetscInt i,j,m = A->rmap->n,N = a->lda;
351: const PetscScalar *v;
352: PetscErrorCode ierr;
355: *fl = PETSC_FALSE;
356: if (A->rmap->n != A->cmap->n) return(0);
357: MatDenseGetArrayRead(A,&v);
358: for (i=0; i<m; i++) {
359: for (j=i; j<m; j++) {
360: if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) {
361: goto restore;
362: }
363: }
364: }
365: *fl = PETSC_TRUE;
366: restore:
367: MatDenseRestoreArrayRead(A,&v);
368: return(0);
369: }
371: static PetscErrorCode MatIsSymmetric_SeqDense(Mat A,PetscReal rtol,PetscBool *fl)
372: {
373: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
374: PetscInt i,j,m = A->rmap->n,N = a->lda;
375: const PetscScalar *v;
376: PetscErrorCode ierr;
379: *fl = PETSC_FALSE;
380: if (A->rmap->n != A->cmap->n) return(0);
381: MatDenseGetArrayRead(A,&v);
382: for (i=0; i<m; i++) {
383: for (j=i; j<m; j++) {
384: if (PetscAbsScalar(v[i+j*N] - v[j+i*N]) > rtol) {
385: goto restore;
386: }
387: }
388: }
389: *fl = PETSC_TRUE;
390: restore:
391: MatDenseRestoreArrayRead(A,&v);
392: return(0);
393: }
395: PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
396: {
397: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
399: PetscInt lda = (PetscInt)mat->lda,j,m,nlda = lda;
402: PetscLayoutReference(A->rmap,&newi->rmap);
403: PetscLayoutReference(A->cmap,&newi->cmap);
404: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { /* propagate LDA */
405: MatDenseSetLDA(newi,lda);
406: }
407: MatSeqDenseSetPreallocation(newi,NULL);
408: if (cpvalues == MAT_COPY_VALUES) {
409: const PetscScalar *av;
410: PetscScalar *v;
412: MatDenseGetArrayRead(A,&av);
413: MatDenseGetArray(newi,&v);
414: MatDenseGetLDA(newi,&nlda);
415: m = A->rmap->n;
416: if (lda>m || nlda>m) {
417: for (j=0; j<A->cmap->n; j++) {
418: PetscArraycpy(v+j*nlda,av+j*lda,m);
419: }
420: } else {
421: PetscArraycpy(v,av,A->rmap->n*A->cmap->n);
422: }
423: MatDenseRestoreArray(newi,&v);
424: MatDenseRestoreArrayRead(A,&av);
425: }
426: return(0);
427: }
429: PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
430: {
434: MatCreate(PetscObjectComm((PetscObject)A),newmat);
435: MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
436: MatSetType(*newmat,((PetscObject)A)->type_name);
437: MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);
438: return(0);
439: }
441: static PetscErrorCode MatSolve_SeqDense_Internal_LU(Mat A, PetscScalar *x, PetscBLASInt xlda, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k, PetscBool T)
442: {
443: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
444: PetscBLASInt info;
445: PetscErrorCode ierr;
448: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
449: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_(T ? "T" : "N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
450: PetscFPTrapPop();
451: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
452: PetscLogFlops(nrhs*(2.0*m*m - m));
453: return(0);
454: }
456: static PetscErrorCode MatConjugate_SeqDense(Mat);
458: static PetscErrorCode MatSolve_SeqDense_Internal_Cholesky(Mat A, PetscScalar *x, PetscBLASInt xlda, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k, PetscBool T)
459: {
460: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
461: PetscBLASInt info;
462: PetscErrorCode ierr;
465: if (A->spd) {
466: if (PetscDefined(USE_COMPLEX) && T) {MatConjugate_SeqDense(A);}
467: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
468: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
469: PetscFPTrapPop();
470: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
471: if (PetscDefined(USE_COMPLEX) && T) {MatConjugate_SeqDense(A);}
472: #if defined(PETSC_USE_COMPLEX)
473: } else if (A->hermitian) {
474: if (T) {MatConjugate_SeqDense(A);}
475: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
476: PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
477: PetscFPTrapPop();
478: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
479: if (T) {MatConjugate_SeqDense(A);}
480: #endif
481: } else { /* symmetric case */
482: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
483: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
484: PetscFPTrapPop();
485: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
486: }
487: PetscLogFlops(nrhs*(2.0*m*m - m));
488: return(0);
489: }
491: static PetscErrorCode MatSolve_SeqDense_Internal_QR(Mat A, PetscScalar *x, PetscBLASInt xlda, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k)
492: {
493: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
494: PetscBLASInt info;
495: char trans;
496: PetscErrorCode ierr;
499: if (PetscDefined(USE_COMPLEX)) {
500: trans = 'C';
501: } else {
502: trans = 'T';
503: }
504: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
505: PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("L", &trans, &m,&nrhs,&mat->rank,mat->v,&mat->lda,mat->tau,x,&xlda,mat->fwork,&mat->lfwork,&info));
506: PetscFPTrapPop();
507: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ORMQR - Bad orthogonal transform");
508: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
509: PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U", "N", "N", &mat->rank,&nrhs,mat->v,&mat->lda,x,&xlda,&info));
510: PetscFPTrapPop();
511: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"TRTRS - Bad triangular solve");
512: for (PetscInt j = 0; j < nrhs; j++) {
513: for (PetscInt i = mat->rank; i < k; i++) {
514: x[j*xlda + i] = 0.;
515: }
516: }
517: PetscLogFlops(nrhs*(4.0*m*mat->rank - PetscSqr(mat->rank)));
518: return(0);
519: }
521: static PetscErrorCode MatSolveTranspose_SeqDense_Internal_QR(Mat A, PetscScalar *x, PetscBLASInt xlda, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k)
522: {
523: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
524: PetscBLASInt info;
525: PetscErrorCode ierr;
528: if (A->rmap->n == A->cmap->n && mat->rank == A->rmap->n) {
529: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
530: PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U", "T", "N", &m,&nrhs,mat->v,&mat->lda,x,&xlda,&info));
531: PetscFPTrapPop();
532: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"TRTRS - Bad triangular solve");
533: if (PetscDefined(USE_COMPLEX)) {MatConjugate_SeqDense(A);}
534: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
535: PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("L", "N", &m,&nrhs,&mat->rank,mat->v,&mat->lda,mat->tau,x,&xlda,mat->fwork,&mat->lfwork,&info));
536: PetscFPTrapPop();
537: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ORMQR - Bad orthogonal transform");
538: if (PetscDefined(USE_COMPLEX)) {MatConjugate_SeqDense(A);}
539: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"QR factored matrix cannot be used for transpose solve");
540: PetscLogFlops(nrhs*(4.0*m*mat->rank - PetscSqr(mat->rank)));
541: return(0);
542: }
544: static PetscErrorCode MatSolve_SeqDense_SetUp(Mat A, Vec xx, Vec yy, PetscScalar **_y, PetscBLASInt *_m, PetscBLASInt *_k)
545: {
546: Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
547: PetscScalar *y;
548: PetscBLASInt m=0, k=0;
549: PetscErrorCode ierr;
552: PetscBLASIntCast(A->rmap->n,&m);
553: PetscBLASIntCast(A->cmap->n,&k);
554: if (k < m) {
555: VecCopy(xx, mat->qrrhs);
556: VecGetArray(mat->qrrhs,&y);
557: } else {
558: VecCopy(xx, yy);
559: VecGetArray(yy,&y);
560: }
561: *_y = y;
562: *_k = k;
563: *_m = m;
564: return(0);
565: }
567: static PetscErrorCode MatSolve_SeqDense_TearDown(Mat A, Vec xx, Vec yy, PetscScalar **_y, PetscBLASInt *_m, PetscBLASInt *_k)
568: {
569: Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
570: PetscScalar *y = NULL;
571: PetscBLASInt m, k;
575: y = *_y;
576: *_y = NULL;
577: k = *_k;
578: m = *_m;
579: if (k < m) {
580: PetscScalar *yv;
581: VecGetArray(yy,&yv);
582: PetscArraycpy(yv, y, k);
583: VecRestoreArray(yy,&yv);
584: VecRestoreArray(mat->qrrhs, &y);
585: } else {
586: VecRestoreArray(yy,&y);
587: }
588: return(0);
589: }
591: static PetscErrorCode MatSolve_SeqDense_LU(Mat A, Vec xx, Vec yy)
592: {
593: PetscScalar *y = NULL;
594: PetscBLASInt m = 0, k = 0;
598: MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k);
599: MatSolve_SeqDense_Internal_LU(A, y, m, m, 1, k, PETSC_FALSE);
600: MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k);
601: return(0);
602: }
604: static PetscErrorCode MatSolveTranspose_SeqDense_LU(Mat A, Vec xx, Vec yy)
605: {
606: PetscScalar *y = NULL;
607: PetscBLASInt m = 0, k = 0;
611: MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k);
612: MatSolve_SeqDense_Internal_LU(A, y, m, m, 1, k, PETSC_TRUE);
613: MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k);
614: return(0);
615: }
617: static PetscErrorCode MatSolve_SeqDense_Cholesky(Mat A, Vec xx, Vec yy)
618: {
619: PetscScalar *y = NULL;
620: PetscBLASInt m = 0, k = 0;
624: MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k);
625: MatSolve_SeqDense_Internal_Cholesky(A, y, m, m, 1, k, PETSC_FALSE);
626: MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k);
627: return(0);
628: }
630: static PetscErrorCode MatSolveTranspose_SeqDense_Cholesky(Mat A, Vec xx, Vec yy)
631: {
632: PetscScalar *y = NULL;
633: PetscBLASInt m = 0, k = 0;
637: MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k);
638: MatSolve_SeqDense_Internal_Cholesky(A, y, m, m, 1, k, PETSC_TRUE);
639: MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k);
640: return(0);
641: }
643: static PetscErrorCode MatSolve_SeqDense_QR(Mat A, Vec xx, Vec yy)
644: {
645: PetscScalar *y = NULL;
646: PetscBLASInt m = 0, k = 0;
650: MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k);
651: MatSolve_SeqDense_Internal_QR(A, y, PetscMax(m,k), m, 1, k);
652: MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k);
653: return(0);
654: }
656: static PetscErrorCode MatSolveTranspose_SeqDense_QR(Mat A, Vec xx, Vec yy)
657: {
658: PetscScalar *y = NULL;
659: PetscBLASInt m = 0, k = 0;
663: MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k);
664: MatSolveTranspose_SeqDense_Internal_QR(A, y, PetscMax(m,k), m, 1, k);
665: MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k);
666: return(0);
667: }
669: static PetscErrorCode MatMatSolve_SeqDense_SetUp(Mat A, Mat B, Mat X, PetscScalar **_y, PetscBLASInt *_ylda, PetscBLASInt *_m, PetscBLASInt *_nrhs, PetscBLASInt *_k)
670: {
671: PetscErrorCode ierr;
672: const PetscScalar *b;
673: PetscScalar *y;
674: PetscInt n, _blda, _xlda;
675: PetscBLASInt nrhs=0,m=0,k=0,blda=0,xlda=0,ylda=0;
678: *_ylda=0; *_m=0; *_nrhs=0; *_k=0;
679: PetscBLASIntCast(A->rmap->n,&m);
680: PetscBLASIntCast(A->cmap->n,&k);
681: MatGetSize(B,NULL,&n);
682: PetscBLASIntCast(n,&nrhs);
683: MatDenseGetLDA(B,&_blda);
684: PetscBLASIntCast(_blda, &blda);
685: MatDenseGetLDA(X,&_xlda);
686: PetscBLASIntCast(_xlda, &xlda);
687: if (xlda < m) {
688: MatDenseGetArrayRead(B,&b);
689: PetscMalloc1(nrhs * m, &y);
690: if (blda == m) {
691: PetscArraycpy(y,b,blda*nrhs);
692: } else {
693: for (PetscInt j = 0; j < nrhs; j++) {
694: PetscArraycpy(&y[j*m],&b[j*blda],m);
695: }
696: }
697: ylda = m;
698: MatDenseRestoreArrayRead(B,&b);
699: } else {
700: if (blda == xlda) {
701: MatCopy(B, X, SAME_NONZERO_PATTERN);
702: MatDenseGetArray(X,&y);
703: } else {
704: MatDenseGetArray(X,&y);
705: MatDenseGetArrayRead(B,&b);
706: for (PetscInt j = 0; j < nrhs; j++) {
707: PetscArraycpy(&y[j*xlda],&b[j*blda],m);
708: }
709: MatDenseRestoreArrayRead(B,&b);
710: }
711: ylda = xlda;
712: }
713: *_y = y;
714: *_ylda = ylda;
715: *_k = k;
716: *_m = m;
717: *_nrhs = nrhs;
718: return(0);
719: }
721: static PetscErrorCode MatMatSolve_SeqDense_TearDown(Mat A, Mat B, Mat X, PetscScalar **_y, PetscBLASInt *_ylda, PetscBLASInt *_m, PetscBLASInt *_nrhs, PetscBLASInt *_k)
722: {
723: PetscScalar *y;
724: PetscInt _xlda;
725: PetscBLASInt k,ylda,nrhs,xlda=0;
726: PetscErrorCode ierr;
729: y = *_y;
730: *_y = NULL;
731: k = *_k;
732: ylda = *_ylda;
733: nrhs = *_nrhs;
734: MatDenseGetLDA(X,&_xlda);
735: PetscBLASIntCast(_xlda, &xlda);
736: if (xlda != ylda) {
737: PetscScalar *xv;
738: MatDenseGetArray(X,&xv);
739: for (PetscInt j = 0; j < nrhs; j++) {
740: PetscArraycpy(&xv[j*xlda],&y[j*ylda],k);
741: }
742: MatDenseRestoreArray(X,&xv);
743: PetscFree(y);
744: } else {
745: MatDenseRestoreArray(X,&y);
746: }
747: return(0);
748: }
750: static PetscErrorCode MatMatSolve_SeqDense_LU(Mat A, Mat B, Mat X)
751: {
752: PetscScalar *y;
753: PetscBLASInt m, k, ylda, nrhs;
757: MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ylda, &m, &nrhs, &k);
758: MatSolve_SeqDense_Internal_LU(A, y, ylda, m, nrhs, k, PETSC_FALSE);
759: MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ylda, &m, &nrhs, &k);
760: return(0);
761: }
763: static PetscErrorCode MatMatSolveTranspose_SeqDense_LU(Mat A, Mat B, Mat X)
764: {
765: PetscScalar *y;
766: PetscBLASInt m, k, ylda, nrhs;
770: MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ylda, &m, &nrhs, &k);
771: MatSolve_SeqDense_Internal_LU(A, y, ylda, m, nrhs, k, PETSC_TRUE);
772: MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ylda, &m, &nrhs, &k);
773: return(0);
774: }
776: static PetscErrorCode MatMatSolve_SeqDense_Cholesky(Mat A, Mat B, Mat X)
777: {
778: PetscScalar *y;
779: PetscBLASInt m, k, ylda, nrhs;
783: MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ylda, &m, &nrhs, &k);
784: MatSolve_SeqDense_Internal_Cholesky(A, y, ylda, m, nrhs, k, PETSC_FALSE);
785: MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ylda, &m, &nrhs, &k);
786: return(0);
787: }
789: static PetscErrorCode MatMatSolveTranspose_SeqDense_Cholesky(Mat A, Mat B, Mat X)
790: {
791: PetscScalar *y;
792: PetscBLASInt m, k, ylda, nrhs;
796: MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ylda, &m, &nrhs, &k);
797: MatSolve_SeqDense_Internal_Cholesky(A, y, ylda, m, nrhs, k, PETSC_TRUE);
798: MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ylda, &m, &nrhs, &k);
799: return(0);
800: }
802: static PetscErrorCode MatMatSolve_SeqDense_QR(Mat A, Mat B, Mat X)
803: {
804: PetscScalar *y;
805: PetscBLASInt m, k, ylda, nrhs;
809: MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ylda, &m, &nrhs, &k);
810: MatSolve_SeqDense_Internal_QR(A, y, ylda, m, nrhs, k);
811: MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ylda, &m, &nrhs, &k);
812: return(0);
813: }
815: static PetscErrorCode MatMatSolveTranspose_SeqDense_QR(Mat A, Mat B, Mat X)
816: {
817: PetscScalar *y;
818: PetscBLASInt m, k, ylda, nrhs;
822: MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ylda, &m, &nrhs, &k);
823: MatSolveTranspose_SeqDense_Internal_QR(A, y, ylda, m, nrhs, k);
824: MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ylda, &m, &nrhs, &k);
825: return(0);
826: }
828: static PetscErrorCode MatConjugate_SeqDense(Mat);
830: /* ---------------------------------------------------------------*/
831: /* COMMENT: I have chosen to hide row permutation in the pivots,
832: rather than put it in the Mat->row slot.*/
833: PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
834: {
835: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
837: PetscBLASInt n,m,info;
840: PetscBLASIntCast(A->cmap->n,&n);
841: PetscBLASIntCast(A->rmap->n,&m);
842: if (!mat->pivots) {
843: PetscMalloc1(A->rmap->n,&mat->pivots);
844: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
845: }
846: if (!A->rmap->n || !A->cmap->n) return(0);
847: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
848: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info));
849: PetscFPTrapPop();
851: if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
852: if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
854: A->ops->solve = MatSolve_SeqDense_LU;
855: A->ops->matsolve = MatMatSolve_SeqDense_LU;
856: A->ops->solvetranspose = MatSolveTranspose_SeqDense_LU;
857: A->ops->matsolvetranspose = MatMatSolveTranspose_SeqDense_LU;
858: A->factortype = MAT_FACTOR_LU;
860: PetscFree(A->solvertype);
861: PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);
863: PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);
864: return(0);
865: }
867: static PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
868: {
869: MatFactorInfo info;
873: MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
874: (*fact->ops->lufactor)(fact,NULL,NULL,&info);
875: return(0);
876: }
878: PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
879: {
881: fact->preallocated = PETSC_TRUE;
882: fact->assembled = PETSC_TRUE;
883: fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense;
884: fact->ops->solve = MatSolve_SeqDense_LU;
885: fact->ops->matsolve = MatMatSolve_SeqDense_LU;
886: fact->ops->solvetranspose = MatSolveTranspose_SeqDense_LU;
887: return(0);
888: }
890: /* Cholesky as L*L^T or L*D*L^T and the symmetric/hermitian complex variants */
891: PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
892: {
893: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
895: PetscBLASInt info,n;
898: PetscBLASIntCast(A->cmap->n,&n);
899: if (!A->rmap->n || !A->cmap->n) return(0);
900: if (A->spd) {
901: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
902: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info));
903: PetscFPTrapPop();
904: #if defined(PETSC_USE_COMPLEX)
905: } else if (A->hermitian) {
906: if (!mat->pivots) {
907: PetscMalloc1(A->rmap->n,&mat->pivots);
908: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
909: }
910: if (!mat->fwork) {
911: PetscScalar dummy;
913: mat->lfwork = -1;
914: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
915: PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
916: PetscFPTrapPop();
917: mat->lfwork = (PetscInt)PetscRealPart(dummy);
918: PetscMalloc1(mat->lfwork,&mat->fwork);
919: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
920: }
921: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
922: PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
923: PetscFPTrapPop();
924: #endif
925: } else { /* symmetric case */
926: if (!mat->pivots) {
927: PetscMalloc1(A->rmap->n,&mat->pivots);
928: PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));
929: }
930: if (!mat->fwork) {
931: PetscScalar dummy;
933: mat->lfwork = -1;
934: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
935: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
936: PetscFPTrapPop();
937: mat->lfwork = (PetscInt)PetscRealPart(dummy);
938: PetscMalloc1(mat->lfwork,&mat->fwork);
939: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
940: }
941: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
942: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
943: PetscFPTrapPop();
944: }
945: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
947: A->ops->solve = MatSolve_SeqDense_Cholesky;
948: A->ops->matsolve = MatMatSolve_SeqDense_Cholesky;
949: A->ops->solvetranspose = MatSolveTranspose_SeqDense_Cholesky;
950: A->ops->matsolvetranspose = MatMatSolveTranspose_SeqDense_Cholesky;
951: A->factortype = MAT_FACTOR_CHOLESKY;
953: PetscFree(A->solvertype);
954: PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);
956: PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);
957: return(0);
958: }
960: static PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
961: {
963: MatFactorInfo info;
966: info.fill = 1.0;
968: MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
969: (*fact->ops->choleskyfactor)(fact,NULL,&info);
970: return(0);
971: }
973: PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
974: {
976: fact->assembled = PETSC_TRUE;
977: fact->preallocated = PETSC_TRUE;
978: fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
979: fact->ops->solve = MatSolve_SeqDense_Cholesky;
980: fact->ops->matsolve = MatMatSolve_SeqDense_Cholesky;
981: fact->ops->solvetranspose = MatSolve_SeqDense_Cholesky;
982: return(0);
983: }
985: static PetscErrorCode MatQRFactor_SeqDense(Mat A,IS col,const MatFactorInfo *minfo)
986: {
987: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
989: PetscBLASInt n,m,info, min, max;
992: PetscBLASIntCast(A->cmap->n,&n);
993: PetscBLASIntCast(A->rmap->n,&m);
994: max = PetscMax(m, n);
995: min = PetscMin(m, n);
996: if (!mat->tau) {
997: PetscMalloc1(min,&mat->tau);
998: PetscLogObjectMemory((PetscObject)A,min*sizeof(PetscScalar));
999: }
1000: if (!mat->pivots) {
1001: PetscMalloc1(m,&mat->pivots);
1002: PetscLogObjectMemory((PetscObject)A,m*sizeof(PetscScalar));
1003: }
1004: if (!mat->qrrhs) {
1005: MatCreateVecs(A, NULL, &(mat->qrrhs));
1006: }
1007: if (!A->rmap->n || !A->cmap->n) return(0);
1008: if (!mat->fwork) {
1009: PetscScalar dummy;
1011: mat->lfwork = -1;
1012: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1013: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&m,&n,mat->v,&mat->lda,mat->tau,&dummy,&mat->lfwork,&info));
1014: PetscFPTrapPop();
1015: mat->lfwork = (PetscInt)PetscRealPart(dummy);
1016: PetscMalloc1(mat->lfwork,&mat->fwork);
1017: PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));
1018: }
1019: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1020: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&m,&n,mat->v,&mat->lda,mat->tau,mat->fwork,&mat->lfwork,&info));
1021: PetscFPTrapPop();
1022: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to QR factorization");
1023: // TODO: try to estimate rank or test for and use geqp3 for rank revealing QR. For now just say rank is min of m and n
1024: mat->rank = min;
1026: A->ops->solve = MatSolve_SeqDense_QR;
1027: A->ops->matsolve = MatMatSolve_SeqDense_QR;
1028: A->factortype = MAT_FACTOR_QR;
1029: if (m == n) {
1030: A->ops->solvetranspose = MatSolveTranspose_SeqDense_QR;
1031: A->ops->matsolvetranspose = MatMatSolveTranspose_SeqDense_QR;
1032: }
1034: PetscFree(A->solvertype);
1035: PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);
1037: PetscLogFlops(2.0*min*min*(max-min/3.0));
1038: return(0);
1039: }
1041: static PetscErrorCode MatQRFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
1042: {
1044: MatFactorInfo info;
1047: info.fill = 1.0;
1049: MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);
1050: MatQRFactor_SeqDense(fact,NULL,&info);
1051: return(0);
1052: }
1054: static PetscErrorCode MatQRFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
1055: {
1059: fact->assembled = PETSC_TRUE;
1060: fact->preallocated = PETSC_TRUE;
1061: fact->ops->solve = MatSolve_SeqDense_QR;
1062: fact->ops->matsolve = MatMatSolve_SeqDense_QR;
1063: if (A->cmap->n == A->rmap->n) {
1064: fact->ops->solvetranspose = MatSolveTranspose_SeqDense_QR;
1065: fact->ops->matsolvetranspose = MatMatSolveTranspose_SeqDense_QR;
1066: }
1067: PetscObjectComposeFunction((PetscObject)fact,"MatQRFactorNumeric_C",MatQRFactorNumeric_SeqDense);
1068: return(0);
1069: }
1072: /* uses LAPACK */
1073: PETSC_INTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
1074: {
1078: MatCreate(PetscObjectComm((PetscObject)A),fact);
1079: MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
1080: MatSetType(*fact,MATDENSE);
1081: if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
1082: (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
1083: (*fact)->ops->ilufactorsymbolic = MatLUFactorSymbolic_SeqDense;
1084: } else {
1085: (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
1086: }
1087: (*fact)->factortype = ftype;
1089: PetscFree((*fact)->solvertype);
1090: PetscStrallocpy(MATSOLVERPETSC,&(*fact)->solvertype);
1091: return(0);
1092: }
1094: /* ------------------------------------------------------------------*/
1095: static PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
1096: {
1097: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1098: PetscScalar *x,*v = mat->v,zero = 0.0,xt;
1099: const PetscScalar *b;
1100: PetscErrorCode ierr;
1101: PetscInt m = A->rmap->n,i;
1102: PetscBLASInt o = 1,bm = 0;
1105: #if defined(PETSC_HAVE_CUDA)
1106: if (A->offloadmask == PETSC_OFFLOAD_GPU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented");
1107: #endif
1108: if (shift == -1) shift = 0.0; /* negative shift indicates do not error on zero diagonal; this code never zeros on zero diagonal */
1109: PetscBLASIntCast(m,&bm);
1110: if (flag & SOR_ZERO_INITIAL_GUESS) {
1111: /* this is a hack fix, should have another version without the second BLASdotu */
1112: VecSet(xx,zero);
1113: }
1114: VecGetArray(xx,&x);
1115: VecGetArrayRead(bb,&b);
1116: its = its*lits;
1117: if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
1118: while (its--) {
1119: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
1120: for (i=0; i<m; i++) {
1121: PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
1122: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
1123: }
1124: }
1125: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
1126: for (i=m-1; i>=0; i--) {
1127: PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
1128: x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
1129: }
1130: }
1131: }
1132: VecRestoreArrayRead(bb,&b);
1133: VecRestoreArray(xx,&x);
1134: return(0);
1135: }
1137: /* -----------------------------------------------------------------*/
1138: PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
1139: {
1140: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1141: const PetscScalar *v = mat->v,*x;
1142: PetscScalar *y;
1143: PetscErrorCode ierr;
1144: PetscBLASInt m, n,_One=1;
1145: PetscScalar _DOne=1.0,_DZero=0.0;
1148: PetscBLASIntCast(A->rmap->n,&m);
1149: PetscBLASIntCast(A->cmap->n,&n);
1150: VecGetArrayRead(xx,&x);
1151: VecGetArrayWrite(yy,&y);
1152: if (!A->rmap->n || !A->cmap->n) {
1153: PetscBLASInt i;
1154: for (i=0; i<n; i++) y[i] = 0.0;
1155: } else {
1156: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One));
1157: PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);
1158: }
1159: VecRestoreArrayRead(xx,&x);
1160: VecRestoreArrayWrite(yy,&y);
1161: return(0);
1162: }
1164: PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
1165: {
1166: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1167: PetscScalar *y,_DOne=1.0,_DZero=0.0;
1168: PetscErrorCode ierr;
1169: PetscBLASInt m, n, _One=1;
1170: const PetscScalar *v = mat->v,*x;
1173: PetscBLASIntCast(A->rmap->n,&m);
1174: PetscBLASIntCast(A->cmap->n,&n);
1175: VecGetArrayRead(xx,&x);
1176: VecGetArrayWrite(yy,&y);
1177: if (!A->rmap->n || !A->cmap->n) {
1178: PetscBLASInt i;
1179: for (i=0; i<m; i++) y[i] = 0.0;
1180: } else {
1181: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One));
1182: PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);
1183: }
1184: VecRestoreArrayRead(xx,&x);
1185: VecRestoreArrayWrite(yy,&y);
1186: return(0);
1187: }
1189: PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
1190: {
1191: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1192: const PetscScalar *v = mat->v,*x;
1193: PetscScalar *y,_DOne=1.0;
1194: PetscErrorCode ierr;
1195: PetscBLASInt m, n, _One=1;
1198: PetscBLASIntCast(A->rmap->n,&m);
1199: PetscBLASIntCast(A->cmap->n,&n);
1200: VecCopy(zz,yy);
1201: if (!A->rmap->n || !A->cmap->n) return(0);
1202: VecGetArrayRead(xx,&x);
1203: VecGetArray(yy,&y);
1204: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
1205: VecRestoreArrayRead(xx,&x);
1206: VecRestoreArray(yy,&y);
1207: PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
1208: return(0);
1209: }
1211: PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
1212: {
1213: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1214: const PetscScalar *v = mat->v,*x;
1215: PetscScalar *y;
1216: PetscErrorCode ierr;
1217: PetscBLASInt m, n, _One=1;
1218: PetscScalar _DOne=1.0;
1221: PetscBLASIntCast(A->rmap->n,&m);
1222: PetscBLASIntCast(A->cmap->n,&n);
1223: VecCopy(zz,yy);
1224: if (!A->rmap->n || !A->cmap->n) return(0);
1225: VecGetArrayRead(xx,&x);
1226: VecGetArray(yy,&y);
1227: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
1228: VecRestoreArrayRead(xx,&x);
1229: VecRestoreArray(yy,&y);
1230: PetscLogFlops(2.0*A->rmap->n*A->cmap->n);
1231: return(0);
1232: }
1234: /* -----------------------------------------------------------------*/
1235: static PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
1236: {
1237: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1239: PetscInt i;
1242: *ncols = A->cmap->n;
1243: if (cols) {
1244: PetscMalloc1(A->cmap->n+1,cols);
1245: for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
1246: }
1247: if (vals) {
1248: const PetscScalar *v;
1250: MatDenseGetArrayRead(A,&v);
1251: PetscMalloc1(A->cmap->n+1,vals);
1252: v += row;
1253: for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
1254: MatDenseRestoreArrayRead(A,&v);
1255: }
1256: return(0);
1257: }
1259: static PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
1260: {
1264: if (cols) {PetscFree(*cols);}
1265: if (vals) {PetscFree(*vals); }
1266: return(0);
1267: }
1268: /* ----------------------------------------------------------------*/
1269: static PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
1270: {
1271: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1272: PetscScalar *av;
1273: PetscInt i,j,idx=0;
1274: #if defined(PETSC_HAVE_CUDA)
1275: PetscOffloadMask oldf;
1276: #endif
1277: PetscErrorCode ierr;
1280: MatDenseGetArray(A,&av);
1281: if (!mat->roworiented) {
1282: if (addv == INSERT_VALUES) {
1283: for (j=0; j<n; j++) {
1284: if (indexn[j] < 0) {idx += m; continue;}
1285: if (PetscUnlikelyDebug(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);
1286: for (i=0; i<m; i++) {
1287: if (indexm[i] < 0) {idx++; continue;}
1288: if (PetscUnlikelyDebug(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);
1289: av[indexn[j]*mat->lda + indexm[i]] = v[idx++];
1290: }
1291: }
1292: } else {
1293: for (j=0; j<n; j++) {
1294: if (indexn[j] < 0) {idx += m; continue;}
1295: if (PetscUnlikelyDebug(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);
1296: for (i=0; i<m; i++) {
1297: if (indexm[i] < 0) {idx++; continue;}
1298: if (PetscUnlikelyDebug(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);
1299: av[indexn[j]*mat->lda + indexm[i]] += v[idx++];
1300: }
1301: }
1302: }
1303: } else {
1304: if (addv == INSERT_VALUES) {
1305: for (i=0; i<m; i++) {
1306: if (indexm[i] < 0) { idx += n; continue;}
1307: if (PetscUnlikelyDebug(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);
1308: for (j=0; j<n; j++) {
1309: if (indexn[j] < 0) { idx++; continue;}
1310: if (PetscUnlikelyDebug(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);
1311: av[indexn[j]*mat->lda + indexm[i]] = v[idx++];
1312: }
1313: }
1314: } else {
1315: for (i=0; i<m; i++) {
1316: if (indexm[i] < 0) { idx += n; continue;}
1317: if (PetscUnlikelyDebug(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);
1318: for (j=0; j<n; j++) {
1319: if (indexn[j] < 0) { idx++; continue;}
1320: if (PetscUnlikelyDebug(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);
1321: av[indexn[j]*mat->lda + indexm[i]] += v[idx++];
1322: }
1323: }
1324: }
1325: }
1326: /* hack to prevent unneeded copy to the GPU while returning the array */
1327: #if defined(PETSC_HAVE_CUDA)
1328: oldf = A->offloadmask;
1329: A->offloadmask = PETSC_OFFLOAD_GPU;
1330: #endif
1331: MatDenseRestoreArray(A,&av);
1332: #if defined(PETSC_HAVE_CUDA)
1333: A->offloadmask = (oldf == PETSC_OFFLOAD_UNALLOCATED ? PETSC_OFFLOAD_UNALLOCATED : PETSC_OFFLOAD_CPU);
1334: #endif
1335: return(0);
1336: }
1338: static PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
1339: {
1340: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1341: const PetscScalar *vv;
1342: PetscInt i,j;
1343: PetscErrorCode ierr;
1346: MatDenseGetArrayRead(A,&vv);
1347: /* row-oriented output */
1348: for (i=0; i<m; i++) {
1349: if (indexm[i] < 0) {v += n;continue;}
1350: 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);
1351: for (j=0; j<n; j++) {
1352: if (indexn[j] < 0) {v++; continue;}
1353: 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);
1354: *v++ = vv[indexn[j]*mat->lda + indexm[i]];
1355: }
1356: }
1357: MatDenseRestoreArrayRead(A,&vv);
1358: return(0);
1359: }
1361: /* -----------------------------------------------------------------*/
1363: PetscErrorCode MatView_Dense_Binary(Mat mat,PetscViewer viewer)
1364: {
1365: PetscErrorCode ierr;
1366: PetscBool skipHeader;
1367: PetscViewerFormat format;
1368: PetscInt header[4],M,N,m,lda,i,j,k;
1369: const PetscScalar *v;
1370: PetscScalar *vwork;
1373: PetscViewerSetUp(viewer);
1374: PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
1375: PetscViewerGetFormat(viewer,&format);
1376: if (skipHeader) format = PETSC_VIEWER_NATIVE;
1378: MatGetSize(mat,&M,&N);
1380: /* write matrix header */
1381: header[0] = MAT_FILE_CLASSID; header[1] = M; header[2] = N;
1382: header[3] = (format == PETSC_VIEWER_NATIVE) ? MATRIX_BINARY_FORMAT_DENSE : M*N;
1383: if (!skipHeader) {PetscViewerBinaryWrite(viewer,header,4,PETSC_INT);}
1385: MatGetLocalSize(mat,&m,NULL);
1386: if (format != PETSC_VIEWER_NATIVE) {
1387: PetscInt nnz = m*N, *iwork;
1388: /* store row lengths for each row */
1389: PetscMalloc1(nnz,&iwork);
1390: for (i=0; i<m; i++) iwork[i] = N;
1391: PetscViewerBinaryWriteAll(viewer,iwork,m,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1392: /* store column indices (zero start index) */
1393: for (k=0, i=0; i<m; i++)
1394: for (j=0; j<N; j++, k++)
1395: iwork[k] = j;
1396: PetscViewerBinaryWriteAll(viewer,iwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1397: PetscFree(iwork);
1398: }
1399: /* store matrix values as a dense matrix in row major order */
1400: PetscMalloc1(m*N,&vwork);
1401: MatDenseGetArrayRead(mat,&v);
1402: MatDenseGetLDA(mat,&lda);
1403: for (k=0, i=0; i<m; i++)
1404: for (j=0; j<N; j++, k++)
1405: vwork[k] = v[i+lda*j];
1406: MatDenseRestoreArrayRead(mat,&v);
1407: PetscViewerBinaryWriteAll(viewer,vwork,m*N,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR);
1408: PetscFree(vwork);
1409: return(0);
1410: }
1412: PetscErrorCode MatLoad_Dense_Binary(Mat mat,PetscViewer viewer)
1413: {
1415: PetscBool skipHeader;
1416: PetscInt header[4],M,N,m,nz,lda,i,j,k;
1417: PetscInt rows,cols;
1418: PetscScalar *v,*vwork;
1421: PetscViewerSetUp(viewer);
1422: PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
1424: if (!skipHeader) {
1425: PetscViewerBinaryRead(viewer,header,4,NULL,PETSC_INT);
1426: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Not a matrix object in file");
1427: M = header[1]; N = header[2];
1428: if (M < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix row size (%D) in file is negative",M);
1429: if (N < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix column size (%D) in file is negative",N);
1430: nz = header[3];
1431: if (nz != MATRIX_BINARY_FORMAT_DENSE && nz < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Unknown matrix format %D in file",nz);
1432: } else {
1433: MatGetSize(mat,&M,&N);
1434: if (M < 0 || N < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Matrix binary file header was skipped, thus the user must specify the global sizes of input matrix");
1435: nz = MATRIX_BINARY_FORMAT_DENSE;
1436: }
1438: /* setup global sizes if not set */
1439: if (mat->rmap->N < 0) mat->rmap->N = M;
1440: if (mat->cmap->N < 0) mat->cmap->N = N;
1441: MatSetUp(mat);
1442: /* check if global sizes are correct */
1443: MatGetSize(mat,&rows,&cols);
1444: if (M != rows || N != cols) SETERRQ4(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different sizes (%d, %d) than the input matrix (%d, %d)",M,N,rows,cols);
1446: MatGetSize(mat,NULL,&N);
1447: MatGetLocalSize(mat,&m,NULL);
1448: MatDenseGetArray(mat,&v);
1449: MatDenseGetLDA(mat,&lda);
1450: if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense format */
1451: PetscInt nnz = m*N;
1452: /* read in matrix values */
1453: PetscMalloc1(nnz,&vwork);
1454: PetscViewerBinaryReadAll(viewer,vwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR);
1455: /* store values in column major order */
1456: for (j=0; j<N; j++)
1457: for (i=0; i<m; i++)
1458: v[i+lda*j] = vwork[i*N+j];
1459: PetscFree(vwork);
1460: } else { /* matrix in file is sparse format */
1461: PetscInt nnz = 0, *rlens, *icols;
1462: /* read in row lengths */
1463: PetscMalloc1(m,&rlens);
1464: PetscViewerBinaryReadAll(viewer,rlens,m,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1465: for (i=0; i<m; i++) nnz += rlens[i];
1466: /* read in column indices and values */
1467: PetscMalloc2(nnz,&icols,nnz,&vwork);
1468: PetscViewerBinaryReadAll(viewer,icols,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
1469: PetscViewerBinaryReadAll(viewer,vwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR);
1470: /* store values in column major order */
1471: for (k=0, i=0; i<m; i++)
1472: for (j=0; j<rlens[i]; j++, k++)
1473: v[i+lda*icols[k]] = vwork[k];
1474: PetscFree(rlens);
1475: PetscFree2(icols,vwork);
1476: }
1477: MatDenseRestoreArray(mat,&v);
1478: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
1479: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
1480: return(0);
1481: }
1483: PetscErrorCode MatLoad_SeqDense(Mat newMat, PetscViewer viewer)
1484: {
1485: PetscBool isbinary, ishdf5;
1491: /* force binary viewer to load .info file if it has not yet done so */
1492: PetscViewerSetUp(viewer);
1493: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1494: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5, &ishdf5);
1495: if (isbinary) {
1496: MatLoad_Dense_Binary(newMat,viewer);
1497: } else if (ishdf5) {
1498: #if defined(PETSC_HAVE_HDF5)
1499: MatLoad_Dense_HDF5(newMat,viewer);
1500: #else
1501: SETERRQ(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1502: #endif
1503: } else {
1504: 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);
1505: }
1506: return(0);
1507: }
1509: static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
1510: {
1511: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1512: PetscErrorCode ierr;
1513: PetscInt i,j;
1514: const char *name;
1515: PetscScalar *v,*av;
1516: PetscViewerFormat format;
1517: #if defined(PETSC_USE_COMPLEX)
1518: PetscBool allreal = PETSC_TRUE;
1519: #endif
1522: MatDenseGetArrayRead(A,(const PetscScalar**)&av);
1523: PetscViewerGetFormat(viewer,&format);
1524: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1525: return(0); /* do nothing for now */
1526: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
1527: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1528: for (i=0; i<A->rmap->n; i++) {
1529: v = av + i;
1530: PetscViewerASCIIPrintf(viewer,"row %D:",i);
1531: for (j=0; j<A->cmap->n; j++) {
1532: #if defined(PETSC_USE_COMPLEX)
1533: if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
1534: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1535: } else if (PetscRealPart(*v)) {
1536: PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)PetscRealPart(*v));
1537: }
1538: #else
1539: if (*v) {
1540: PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)*v);
1541: }
1542: #endif
1543: v += a->lda;
1544: }
1545: PetscViewerASCIIPrintf(viewer,"\n");
1546: }
1547: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1548: } else {
1549: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1550: #if defined(PETSC_USE_COMPLEX)
1551: /* determine if matrix has all real values */
1552: v = av;
1553: for (i=0; i<A->rmap->n*A->cmap->n; i++) {
1554: if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;}
1555: }
1556: #endif
1557: if (format == PETSC_VIEWER_ASCII_MATLAB) {
1558: PetscObjectGetName((PetscObject)A,&name);
1559: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);
1560: PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);
1561: PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
1562: }
1564: for (i=0; i<A->rmap->n; i++) {
1565: v = av + i;
1566: for (j=0; j<A->cmap->n; j++) {
1567: #if defined(PETSC_USE_COMPLEX)
1568: if (allreal) {
1569: PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));
1570: } else {
1571: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
1572: }
1573: #else
1574: PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);
1575: #endif
1576: v += a->lda;
1577: }
1578: PetscViewerASCIIPrintf(viewer,"\n");
1579: }
1580: if (format == PETSC_VIEWER_ASCII_MATLAB) {
1581: PetscViewerASCIIPrintf(viewer,"];\n");
1582: }
1583: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1584: }
1585: MatDenseRestoreArrayRead(A,(const PetscScalar**)&av);
1586: PetscViewerFlush(viewer);
1587: return(0);
1588: }
1590: #include <petscdraw.h>
1591: static PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1592: {
1593: Mat A = (Mat) Aa;
1594: PetscErrorCode ierr;
1595: PetscInt m = A->rmap->n,n = A->cmap->n,i,j;
1596: int color = PETSC_DRAW_WHITE;
1597: const PetscScalar *v;
1598: PetscViewer viewer;
1599: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r;
1600: PetscViewerFormat format;
1603: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
1604: PetscViewerGetFormat(viewer,&format);
1605: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
1607: /* Loop over matrix elements drawing boxes */
1608: MatDenseGetArrayRead(A,&v);
1609: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1610: PetscDrawCollectiveBegin(draw);
1611: /* Blue for negative and Red for positive */
1612: for (j = 0; j < n; j++) {
1613: x_l = j; x_r = x_l + 1.0;
1614: for (i = 0; i < m; i++) {
1615: y_l = m - i - 1.0;
1616: y_r = y_l + 1.0;
1617: if (PetscRealPart(v[j*m+i]) > 0.) color = PETSC_DRAW_RED;
1618: else if (PetscRealPart(v[j*m+i]) < 0.) color = PETSC_DRAW_BLUE;
1619: else continue;
1620: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1621: }
1622: }
1623: PetscDrawCollectiveEnd(draw);
1624: } else {
1625: /* use contour shading to indicate magnitude of values */
1626: /* first determine max of all nonzero values */
1627: PetscReal minv = 0.0, maxv = 0.0;
1628: PetscDraw popup;
1630: for (i=0; i < m*n; i++) {
1631: if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1632: }
1633: if (minv >= maxv) maxv = minv + PETSC_SMALL;
1634: PetscDrawGetPopup(draw,&popup);
1635: PetscDrawScalePopup(popup,minv,maxv);
1637: PetscDrawCollectiveBegin(draw);
1638: for (j=0; j<n; j++) {
1639: x_l = j;
1640: x_r = x_l + 1.0;
1641: for (i=0; i<m; i++) {
1642: y_l = m - i - 1.0;
1643: y_r = y_l + 1.0;
1644: color = PetscDrawRealToColor(PetscAbsScalar(v[j*m+i]),minv,maxv);
1645: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
1646: }
1647: }
1648: PetscDrawCollectiveEnd(draw);
1649: }
1650: MatDenseRestoreArrayRead(A,&v);
1651: return(0);
1652: }
1654: static PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1655: {
1656: PetscDraw draw;
1657: PetscBool isnull;
1658: PetscReal xr,yr,xl,yl,h,w;
1662: PetscViewerDrawGetDraw(viewer,0,&draw);
1663: PetscDrawIsNull(draw,&isnull);
1664: if (isnull) return(0);
1666: xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1667: xr += w; yr += h; xl = -w; yl = -h;
1668: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
1669: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
1670: PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);
1671: PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
1672: PetscDrawSave(draw);
1673: return(0);
1674: }
1676: PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1677: {
1679: PetscBool iascii,isbinary,isdraw;
1682: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1683: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1684: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1686: if (iascii) {
1687: MatView_SeqDense_ASCII(A,viewer);
1688: } else if (isbinary) {
1689: MatView_Dense_Binary(A,viewer);
1690: } else if (isdraw) {
1691: MatView_SeqDense_Draw(A,viewer);
1692: }
1693: return(0);
1694: }
1696: static PetscErrorCode MatDensePlaceArray_SeqDense(Mat A,const PetscScalar *array)
1697: {
1698: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1701: if (a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1702: if (a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1703: if (a->unplacedarray) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreArray() first");
1704: a->unplacedarray = a->v;
1705: a->unplaced_user_alloc = a->user_alloc;
1706: a->v = (PetscScalar*) array;
1707: a->user_alloc = PETSC_TRUE;
1708: #if defined(PETSC_HAVE_CUDA)
1709: A->offloadmask = PETSC_OFFLOAD_CPU;
1710: #endif
1711: return(0);
1712: }
1714: static PetscErrorCode MatDenseResetArray_SeqDense(Mat A)
1715: {
1716: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1719: if (a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1720: if (a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1721: a->v = a->unplacedarray;
1722: a->user_alloc = a->unplaced_user_alloc;
1723: a->unplacedarray = NULL;
1724: #if defined(PETSC_HAVE_CUDA)
1725: A->offloadmask = PETSC_OFFLOAD_CPU;
1726: #endif
1727: return(0);
1728: }
1730: static PetscErrorCode MatDenseReplaceArray_SeqDense(Mat A,const PetscScalar *array)
1731: {
1732: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
1736: if (a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1737: if (a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1738: if (!a->user_alloc) { PetscFree(a->v); }
1739: a->v = (PetscScalar*) array;
1740: a->user_alloc = PETSC_FALSE;
1741: #if defined(PETSC_HAVE_CUDA)
1742: A->offloadmask = PETSC_OFFLOAD_CPU;
1743: #endif
1744: return(0);
1745: }
1747: PetscErrorCode MatDestroy_SeqDense(Mat mat)
1748: {
1749: Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
1753: #if defined(PETSC_USE_LOG)
1754: PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1755: #endif
1756: VecDestroy(&(l->qrrhs));
1757: PetscFree(l->tau);
1758: PetscFree(l->pivots);
1759: PetscFree(l->fwork);
1760: MatDestroy(&l->ptapwork);
1761: if (!l->user_alloc) {PetscFree(l->v);}
1762: if (!l->unplaced_user_alloc) {PetscFree(l->unplacedarray);}
1763: if (l->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1764: if (l->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1765: VecDestroy(&l->cvec);
1766: MatDestroy(&l->cmat);
1767: PetscFree(mat->data);
1769: PetscObjectChangeTypeName((PetscObject)mat,NULL);
1770: PetscObjectComposeFunction((PetscObject)mat,"MatQRFactor_C",NULL);
1771: PetscObjectComposeFunction((PetscObject)mat,"MatQRFactorNumeric_C",NULL);
1772: PetscObjectComposeFunction((PetscObject)mat,"MatQRFactorSymbolic_C",NULL);
1773: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL);
1774: PetscObjectComposeFunction((PetscObject)mat,"MatDenseSetLDA_C",NULL);
1775: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);
1776: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);
1777: PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);
1778: PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);
1779: PetscObjectComposeFunction((PetscObject)mat,"MatDenseReplaceArray_C",NULL);
1780: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",NULL);
1781: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",NULL);
1782: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayWrite_C",NULL);
1783: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayWrite_C",NULL);
1784: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqaij_C",NULL);
1785: #if defined(PETSC_HAVE_ELEMENTAL)
1786: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_elemental_C",NULL);
1787: #endif
1788: #if defined(PETSC_HAVE_SCALAPACK)
1789: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_scalapack_C",NULL);
1790: #endif
1791: #if defined(PETSC_HAVE_CUDA)
1792: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqdensecuda_C",NULL);
1793: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdensecuda_seqdensecuda_C",NULL);
1794: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdensecuda_seqdense_C",NULL);
1795: #endif
1796: PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL);
1797: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqaij_seqdense_C",NULL);
1798: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdense_seqdense_C",NULL);
1799: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqbaij_seqdense_C",NULL);
1800: PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqsbaij_seqdense_C",NULL);
1802: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);
1803: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);
1804: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVec_C",NULL);
1805: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVec_C",NULL);
1806: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecRead_C",NULL);
1807: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecRead_C",NULL);
1808: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecWrite_C",NULL);
1809: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecWrite_C",NULL);
1810: PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetSubMatrix_C",NULL);
1811: PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreSubMatrix_C",NULL);
1812: return(0);
1813: }
1815: static PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1816: {
1817: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1819: PetscInt k,j,m = A->rmap->n, M = mat->lda, n = A->cmap->n;
1820: PetscScalar *v,tmp;
1823: if (reuse == MAT_INPLACE_MATRIX) {
1824: if (m == n) { /* in place transpose */
1825: MatDenseGetArray(A,&v);
1826: for (j=0; j<m; j++) {
1827: for (k=0; k<j; k++) {
1828: tmp = v[j + k*M];
1829: v[j + k*M] = v[k + j*M];
1830: v[k + j*M] = tmp;
1831: }
1832: }
1833: MatDenseRestoreArray(A,&v);
1834: } else { /* reuse memory, temporary allocates new memory */
1835: PetscScalar *v2;
1836: PetscLayout tmplayout;
1838: PetscMalloc1((size_t)m*n,&v2);
1839: MatDenseGetArray(A,&v);
1840: for (j=0; j<n; j++) {
1841: for (k=0; k<m; k++) v2[j + (size_t)k*n] = v[k + (size_t)j*M];
1842: }
1843: PetscArraycpy(v,v2,(size_t)m*n);
1844: PetscFree(v2);
1845: MatDenseRestoreArray(A,&v);
1846: /* cleanup size dependent quantities */
1847: VecDestroy(&mat->cvec);
1848: MatDestroy(&mat->cmat);
1849: PetscFree(mat->pivots);
1850: PetscFree(mat->fwork);
1851: MatDestroy(&mat->ptapwork);
1852: /* swap row/col layouts */
1853: mat->lda = n;
1854: tmplayout = A->rmap;
1855: A->rmap = A->cmap;
1856: A->cmap = tmplayout;
1857: }
1858: } else { /* out-of-place transpose */
1859: Mat tmat;
1860: Mat_SeqDense *tmatd;
1861: PetscScalar *v2;
1862: PetscInt M2;
1864: if (reuse == MAT_INITIAL_MATRIX) {
1865: MatCreate(PetscObjectComm((PetscObject)A),&tmat);
1866: MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);
1867: MatSetType(tmat,((PetscObject)A)->type_name);
1868: MatSeqDenseSetPreallocation(tmat,NULL);
1869: } else tmat = *matout;
1871: MatDenseGetArrayRead(A,(const PetscScalar**)&v);
1872: MatDenseGetArray(tmat,&v2);
1873: tmatd = (Mat_SeqDense*)tmat->data;
1874: M2 = tmatd->lda;
1875: for (j=0; j<n; j++) {
1876: for (k=0; k<m; k++) v2[j + k*M2] = v[k + j*M];
1877: }
1878: MatDenseRestoreArray(tmat,&v2);
1879: MatDenseRestoreArrayRead(A,(const PetscScalar**)&v);
1880: MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
1881: MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
1882: *matout = tmat;
1883: }
1884: return(0);
1885: }
1887: static PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg)
1888: {
1889: Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1890: Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1891: PetscInt i;
1892: const PetscScalar *v1,*v2;
1893: PetscErrorCode ierr;
1896: if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; return(0);}
1897: if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; return(0);}
1898: MatDenseGetArrayRead(A1,&v1);
1899: MatDenseGetArrayRead(A2,&v2);
1900: for (i=0; i<A1->cmap->n; i++) {
1901: PetscArraycmp(v1,v2,A1->rmap->n,flg);
1902: if (*flg == PETSC_FALSE) return(0);
1903: v1 += mat1->lda;
1904: v2 += mat2->lda;
1905: }
1906: MatDenseRestoreArrayRead(A1,&v1);
1907: MatDenseRestoreArrayRead(A2,&v2);
1908: *flg = PETSC_TRUE;
1909: return(0);
1910: }
1912: static PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1913: {
1914: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1915: PetscInt i,n,len;
1916: PetscScalar *x;
1917: const PetscScalar *vv;
1918: PetscErrorCode ierr;
1921: VecGetSize(v,&n);
1922: VecGetArray(v,&x);
1923: len = PetscMin(A->rmap->n,A->cmap->n);
1924: MatDenseGetArrayRead(A,&vv);
1925: if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1926: for (i=0; i<len; i++) {
1927: x[i] = vv[i*mat->lda + i];
1928: }
1929: MatDenseRestoreArrayRead(A,&vv);
1930: VecRestoreArray(v,&x);
1931: return(0);
1932: }
1934: static PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1935: {
1936: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1937: const PetscScalar *l,*r;
1938: PetscScalar x,*v,*vv;
1939: PetscErrorCode ierr;
1940: PetscInt i,j,m = A->rmap->n,n = A->cmap->n;
1943: MatDenseGetArray(A,&vv);
1944: if (ll) {
1945: VecGetSize(ll,&m);
1946: VecGetArrayRead(ll,&l);
1947: if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1948: for (i=0; i<m; i++) {
1949: x = l[i];
1950: v = vv + i;
1951: for (j=0; j<n; j++) { (*v) *= x; v+= mat->lda;}
1952: }
1953: VecRestoreArrayRead(ll,&l);
1954: PetscLogFlops(1.0*n*m);
1955: }
1956: if (rr) {
1957: VecGetSize(rr,&n);
1958: VecGetArrayRead(rr,&r);
1959: if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1960: for (i=0; i<n; i++) {
1961: x = r[i];
1962: v = vv + i*mat->lda;
1963: for (j=0; j<m; j++) (*v++) *= x;
1964: }
1965: VecRestoreArrayRead(rr,&r);
1966: PetscLogFlops(1.0*n*m);
1967: }
1968: MatDenseRestoreArray(A,&vv);
1969: return(0);
1970: }
1972: PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1973: {
1974: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1975: PetscScalar *v,*vv;
1976: PetscReal sum = 0.0;
1977: PetscInt lda =mat->lda,m=A->rmap->n,i,j;
1978: PetscErrorCode ierr;
1981: MatDenseGetArrayRead(A,(const PetscScalar**)&vv);
1982: v = vv;
1983: if (type == NORM_FROBENIUS) {
1984: if (lda>m) {
1985: for (j=0; j<A->cmap->n; j++) {
1986: v = vv+j*lda;
1987: for (i=0; i<m; i++) {
1988: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1989: }
1990: }
1991: } else {
1992: #if defined(PETSC_USE_REAL___FP16)
1993: PetscBLASInt one = 1,cnt = A->cmap->n*A->rmap->n;
1994: PetscStackCallBLAS("BLASnrm2",*nrm = BLASnrm2_(&cnt,v,&one));
1995: }
1996: #else
1997: for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1998: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1999: }
2000: }
2001: *nrm = PetscSqrtReal(sum);
2002: #endif
2003: PetscLogFlops(2.0*A->cmap->n*A->rmap->n);
2004: } else if (type == NORM_1) {
2005: *nrm = 0.0;
2006: for (j=0; j<A->cmap->n; j++) {
2007: v = vv + j*mat->lda;
2008: sum = 0.0;
2009: for (i=0; i<A->rmap->n; i++) {
2010: sum += PetscAbsScalar(*v); v++;
2011: }
2012: if (sum > *nrm) *nrm = sum;
2013: }
2014: PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
2015: } else if (type == NORM_INFINITY) {
2016: *nrm = 0.0;
2017: for (j=0; j<A->rmap->n; j++) {
2018: v = vv + j;
2019: sum = 0.0;
2020: for (i=0; i<A->cmap->n; i++) {
2021: sum += PetscAbsScalar(*v); v += mat->lda;
2022: }
2023: if (sum > *nrm) *nrm = sum;
2024: }
2025: PetscLogFlops(1.0*A->cmap->n*A->rmap->n);
2026: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
2027: MatDenseRestoreArrayRead(A,(const PetscScalar**)&vv);
2028: return(0);
2029: }
2031: static PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg)
2032: {
2033: Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
2037: switch (op) {
2038: case MAT_ROW_ORIENTED:
2039: aij->roworiented = flg;
2040: break;
2041: case MAT_NEW_NONZERO_LOCATIONS:
2042: case MAT_NEW_NONZERO_LOCATION_ERR:
2043: case MAT_NEW_NONZERO_ALLOCATION_ERR:
2044: case MAT_FORCE_DIAGONAL_ENTRIES:
2045: case MAT_KEEP_NONZERO_PATTERN:
2046: case MAT_IGNORE_OFF_PROC_ENTRIES:
2047: case MAT_USE_HASH_TABLE:
2048: case MAT_IGNORE_ZERO_ENTRIES:
2049: case MAT_IGNORE_LOWER_TRIANGULAR:
2050: case MAT_SORTED_FULL:
2051: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
2052: break;
2053: case MAT_SPD:
2054: case MAT_SYMMETRIC:
2055: case MAT_STRUCTURALLY_SYMMETRIC:
2056: case MAT_HERMITIAN:
2057: case MAT_SYMMETRY_ETERNAL:
2058: /* These options are handled directly by MatSetOption() */
2059: break;
2060: default:
2061: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
2062: }
2063: return(0);
2064: }
2066: static PetscErrorCode MatZeroEntries_SeqDense(Mat A)
2067: {
2068: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
2070: PetscInt lda=l->lda,m=A->rmap->n,j;
2071: PetscScalar *v;
2074: MatDenseGetArray(A,&v);
2075: if (lda>m) {
2076: for (j=0; j<A->cmap->n; j++) {
2077: PetscArrayzero(v+j*lda,m);
2078: }
2079: } else {
2080: PetscArrayzero(v,A->rmap->n*A->cmap->n);
2081: }
2082: MatDenseRestoreArray(A,&v);
2083: return(0);
2084: }
2086: static PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2087: {
2088: PetscErrorCode ierr;
2089: Mat_SeqDense *l = (Mat_SeqDense*)A->data;
2090: PetscInt m = l->lda, n = A->cmap->n, i,j;
2091: PetscScalar *slot,*bb,*v;
2092: const PetscScalar *xx;
2095: if (PetscDefined(USE_DEBUG)) {
2096: for (i=0; i<N; i++) {
2097: if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
2098: 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);
2099: }
2100: }
2101: if (!N) return(0);
2103: /* fix right hand side if needed */
2104: if (x && b) {
2105: VecGetArrayRead(x,&xx);
2106: VecGetArray(b,&bb);
2107: for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
2108: VecRestoreArrayRead(x,&xx);
2109: VecRestoreArray(b,&bb);
2110: }
2112: MatDenseGetArray(A,&v);
2113: for (i=0; i<N; i++) {
2114: slot = v + rows[i];
2115: for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
2116: }
2117: if (diag != 0.0) {
2118: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
2119: for (i=0; i<N; i++) {
2120: slot = v + (m+1)*rows[i];
2121: *slot = diag;
2122: }
2123: }
2124: MatDenseRestoreArray(A,&v);
2125: return(0);
2126: }
2128: static PetscErrorCode MatDenseGetLDA_SeqDense(Mat A,PetscInt *lda)
2129: {
2130: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2133: *lda = mat->lda;
2134: return(0);
2135: }
2137: PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar **array)
2138: {
2139: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2142: if (mat->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
2143: *array = mat->v;
2144: return(0);
2145: }
2147: PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar **array)
2148: {
2150: *array = NULL;
2151: return(0);
2152: }
2154: /*@C
2155: MatDenseGetLDA - gets the leading dimension of the array returned from MatDenseGetArray()
2157: Not collective
2159: Input Parameter:
2160: . mat - a MATSEQDENSE or MATMPIDENSE matrix
2162: Output Parameter:
2163: . lda - the leading dimension
2165: Level: intermediate
2167: .seealso: MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatDenseSetLDA()
2168: @*/
2169: PetscErrorCode MatDenseGetLDA(Mat A,PetscInt *lda)
2170: {
2176: PetscUseMethod(A,"MatDenseGetLDA_C",(Mat,PetscInt*),(A,lda));
2177: return(0);
2178: }
2180: /*@C
2181: MatDenseSetLDA - Sets the leading dimension of the array used by the dense matrix
2183: Not collective
2185: Input Parameter:
2186: + mat - a MATSEQDENSE or MATMPIDENSE matrix
2187: - lda - the leading dimension
2189: Level: intermediate
2191: .seealso: MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatDenseGetLDA()
2192: @*/
2193: PetscErrorCode MatDenseSetLDA(Mat A,PetscInt lda)
2194: {
2199: PetscTryMethod(A,"MatDenseSetLDA_C",(Mat,PetscInt),(A,lda));
2200: return(0);
2201: }
2203: /*@C
2204: MatDenseGetArray - gives read-write access to the array where the data for a dense matrix is stored
2206: Logically Collective on Mat
2208: Input Parameter:
2209: . mat - a dense matrix
2211: Output Parameter:
2212: . array - pointer to the data
2214: Level: intermediate
2216: .seealso: MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatDenseGetArrayWrite(), MatDenseRestoreArrayWrite()
2217: @*/
2218: PetscErrorCode MatDenseGetArray(Mat A,PetscScalar **array)
2219: {
2225: PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));
2226: return(0);
2227: }
2229: /*@C
2230: MatDenseRestoreArray - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()
2232: Logically Collective on Mat
2234: Input Parameters:
2235: + mat - a dense matrix
2236: - array - pointer to the data
2238: Level: intermediate
2240: .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatDenseGetArrayWrite(), MatDenseRestoreArrayWrite()
2241: @*/
2242: PetscErrorCode MatDenseRestoreArray(Mat A,PetscScalar **array)
2243: {
2249: PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));
2250: PetscObjectStateIncrease((PetscObject)A);
2251: #if defined(PETSC_HAVE_CUDA)
2252: A->offloadmask = PETSC_OFFLOAD_CPU;
2253: #endif
2254: return(0);
2255: }
2257: /*@C
2258: MatDenseGetArrayRead - gives read-only access to the array where the data for a dense matrix is stored
2260: Not Collective
2262: Input Parameter:
2263: . mat - a dense matrix
2265: Output Parameter:
2266: . array - pointer to the data
2268: Level: intermediate
2270: .seealso: MatDenseRestoreArrayRead(), MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayWrite(), MatDenseRestoreArrayWrite()
2271: @*/
2272: PetscErrorCode MatDenseGetArrayRead(Mat A,const PetscScalar **array)
2273: {
2279: PetscUseMethod(A,"MatDenseGetArrayRead_C",(Mat,const PetscScalar**),(A,array));
2280: return(0);
2281: }
2283: /*@C
2284: MatDenseRestoreArrayRead - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArrayRead()
2286: Not Collective
2288: Input Parameters:
2289: + mat - a dense matrix
2290: - array - pointer to the data
2292: Level: intermediate
2294: .seealso: MatDenseGetArrayRead(), MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayWrite(), MatDenseRestoreArrayWrite()
2295: @*/
2296: PetscErrorCode MatDenseRestoreArrayRead(Mat A,const PetscScalar **array)
2297: {
2303: PetscUseMethod(A,"MatDenseRestoreArrayRead_C",(Mat,const PetscScalar**),(A,array));
2304: return(0);
2305: }
2307: /*@C
2308: MatDenseGetArrayWrite - gives write-only access to the array where the data for a dense matrix is stored
2310: Not Collective
2312: Input Parameter:
2313: . mat - a dense matrix
2315: Output Parameter:
2316: . array - pointer to the data
2318: Level: intermediate
2320: .seealso: MatDenseRestoreArrayWrite(), MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
2321: @*/
2322: PetscErrorCode MatDenseGetArrayWrite(Mat A,PetscScalar **array)
2323: {
2329: PetscUseMethod(A,"MatDenseGetArrayWrite_C",(Mat,PetscScalar**),(A,array));
2330: return(0);
2331: }
2333: /*@C
2334: MatDenseRestoreArrayWrite - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArrayWrite()
2336: Not Collective
2338: Input Parameters:
2339: + mat - a dense matrix
2340: - array - pointer to the data
2342: Level: intermediate
2344: .seealso: MatDenseGetArrayWrite(), MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
2345: @*/
2346: PetscErrorCode MatDenseRestoreArrayWrite(Mat A,PetscScalar **array)
2347: {
2353: PetscUseMethod(A,"MatDenseRestoreArrayWrite_C",(Mat,PetscScalar**),(A,array));
2354: PetscObjectStateIncrease((PetscObject)A);
2355: #if defined(PETSC_HAVE_CUDA)
2356: A->offloadmask = PETSC_OFFLOAD_CPU;
2357: #endif
2358: return(0);
2359: }
2361: static PetscErrorCode MatCreateSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
2362: {
2363: Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2365: PetscInt i,j,nrows,ncols,blda;
2366: const PetscInt *irow,*icol;
2367: PetscScalar *av,*bv,*v = mat->v;
2368: Mat newmat;
2371: ISGetIndices(isrow,&irow);
2372: ISGetIndices(iscol,&icol);
2373: ISGetLocalSize(isrow,&nrows);
2374: ISGetLocalSize(iscol,&ncols);
2376: /* Check submatrixcall */
2377: if (scall == MAT_REUSE_MATRIX) {
2378: PetscInt n_cols,n_rows;
2379: MatGetSize(*B,&n_rows,&n_cols);
2380: if (n_rows != nrows || n_cols != ncols) {
2381: /* resize the result matrix to match number of requested rows/columns */
2382: MatSetSizes(*B,nrows,ncols,nrows,ncols);
2383: }
2384: newmat = *B;
2385: } else {
2386: /* Create and fill new matrix */
2387: MatCreate(PetscObjectComm((PetscObject)A),&newmat);
2388: MatSetSizes(newmat,nrows,ncols,nrows,ncols);
2389: MatSetType(newmat,((PetscObject)A)->type_name);
2390: MatSeqDenseSetPreallocation(newmat,NULL);
2391: }
2393: /* Now extract the data pointers and do the copy,column at a time */
2394: MatDenseGetArray(newmat,&bv);
2395: MatDenseGetLDA(newmat,&blda);
2396: for (i=0; i<ncols; i++) {
2397: av = v + mat->lda*icol[i];
2398: for (j=0; j<nrows; j++) bv[j] = av[irow[j]];
2399: bv += blda;
2400: }
2401: MatDenseRestoreArray(newmat,&bv);
2403: /* Assemble the matrices so that the correct flags are set */
2404: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2405: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2407: /* Free work space */
2408: ISRestoreIndices(isrow,&irow);
2409: ISRestoreIndices(iscol,&icol);
2410: *B = newmat;
2411: return(0);
2412: }
2414: static PetscErrorCode MatCreateSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
2415: {
2417: PetscInt i;
2420: if (scall == MAT_INITIAL_MATRIX) {
2421: PetscCalloc1(n+1,B);
2422: }
2424: for (i=0; i<n; i++) {
2425: MatCreateSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);
2426: }
2427: return(0);
2428: }
2430: static PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
2431: {
2433: return(0);
2434: }
2436: static PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
2437: {
2439: return(0);
2440: }
2442: PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
2443: {
2444: Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data;
2445: PetscErrorCode ierr;
2446: const PetscScalar *va;
2447: PetscScalar *vb;
2448: PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
2451: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
2452: if (A->ops->copy != B->ops->copy) {
2453: MatCopy_Basic(A,B,str);
2454: return(0);
2455: }
2456: if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
2457: MatDenseGetArrayRead(A,&va);
2458: MatDenseGetArray(B,&vb);
2459: if (lda1>m || lda2>m) {
2460: for (j=0; j<n; j++) {
2461: PetscArraycpy(vb+j*lda2,va+j*lda1,m);
2462: }
2463: } else {
2464: PetscArraycpy(vb,va,A->rmap->n*A->cmap->n);
2465: }
2466: MatDenseRestoreArray(B,&vb);
2467: MatDenseRestoreArrayRead(A,&va);
2468: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2469: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2470: return(0);
2471: }
2473: static PetscErrorCode MatSetUp_SeqDense(Mat A)
2474: {
2478: PetscLayoutSetUp(A->rmap);
2479: PetscLayoutSetUp(A->cmap);
2480: if (!A->preallocated) {
2481: MatSeqDenseSetPreallocation(A,NULL);
2482: }
2483: return(0);
2484: }
2486: static PetscErrorCode MatConjugate_SeqDense(Mat A)
2487: {
2488: Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
2489: PetscInt i,nz = A->rmap->n*A->cmap->n;
2490: PetscInt min = PetscMin(A->rmap->n,A->cmap->n);
2491: PetscScalar *aa;
2495: MatDenseGetArray(A,&aa);
2496: for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
2497: MatDenseRestoreArray(A,&aa);
2498: if (mat->tau) for (i = 0; i < min; i++) mat->tau[i] = PetscConj(mat->tau[i]);
2499: return(0);
2500: }
2502: static PetscErrorCode MatRealPart_SeqDense(Mat A)
2503: {
2504: PetscInt i,nz = A->rmap->n*A->cmap->n;
2505: PetscScalar *aa;
2509: MatDenseGetArray(A,&aa);
2510: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
2511: MatDenseRestoreArray(A,&aa);
2512: return(0);
2513: }
2515: static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
2516: {
2517: PetscInt i,nz = A->rmap->n*A->cmap->n;
2518: PetscScalar *aa;
2522: MatDenseGetArray(A,&aa);
2523: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
2524: MatDenseRestoreArray(A,&aa);
2525: return(0);
2526: }
2528: /* ----------------------------------------------------------------*/
2529: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
2530: {
2532: PetscInt m=A->rmap->n,n=B->cmap->n;
2533: PetscBool cisdense;
2536: MatSetSizes(C,m,n,m,n);
2537: PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"");
2538: if (!cisdense) {
2539: PetscBool flg;
2541: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2542: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
2543: }
2544: MatSetUp(C);
2545: return(0);
2546: }
2548: PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2549: {
2550: Mat_SeqDense *a=(Mat_SeqDense*)A->data,*b=(Mat_SeqDense*)B->data,*c=(Mat_SeqDense*)C->data;
2551: PetscBLASInt m,n,k;
2552: const PetscScalar *av,*bv;
2553: PetscScalar *cv;
2554: PetscScalar _DOne=1.0,_DZero=0.0;
2555: PetscErrorCode ierr;
2558: PetscBLASIntCast(C->rmap->n,&m);
2559: PetscBLASIntCast(C->cmap->n,&n);
2560: PetscBLASIntCast(A->cmap->n,&k);
2561: if (!m || !n || !k) return(0);
2562: MatDenseGetArrayRead(A,&av);
2563: MatDenseGetArrayRead(B,&bv);
2564: MatDenseGetArrayWrite(C,&cv);
2565: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda));
2566: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2567: MatDenseRestoreArrayRead(A,&av);
2568: MatDenseRestoreArrayRead(B,&bv);
2569: MatDenseRestoreArrayWrite(C,&cv);
2570: return(0);
2571: }
2573: PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
2574: {
2576: PetscInt m=A->rmap->n,n=B->rmap->n;
2577: PetscBool cisdense;
2580: MatSetSizes(C,m,n,m,n);
2581: PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"");
2582: if (!cisdense) {
2583: PetscBool flg;
2585: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2586: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
2587: }
2588: MatSetUp(C);
2589: return(0);
2590: }
2592: PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2593: {
2594: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2595: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2596: Mat_SeqDense *c = (Mat_SeqDense*)C->data;
2597: const PetscScalar *av,*bv;
2598: PetscScalar *cv;
2599: PetscBLASInt m,n,k;
2600: PetscScalar _DOne=1.0,_DZero=0.0;
2601: PetscErrorCode ierr;
2604: PetscBLASIntCast(C->rmap->n,&m);
2605: PetscBLASIntCast(C->cmap->n,&n);
2606: PetscBLASIntCast(A->cmap->n,&k);
2607: if (!m || !n || !k) return(0);
2608: MatDenseGetArrayRead(A,&av);
2609: MatDenseGetArrayRead(B,&bv);
2610: MatDenseGetArrayWrite(C,&cv);
2611: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda));
2612: MatDenseRestoreArrayRead(A,&av);
2613: MatDenseRestoreArrayRead(B,&bv);
2614: MatDenseRestoreArrayWrite(C,&cv);
2615: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2616: return(0);
2617: }
2619: PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
2620: {
2622: PetscInt m=A->cmap->n,n=B->cmap->n;
2623: PetscBool cisdense;
2626: MatSetSizes(C,m,n,m,n);
2627: PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"");
2628: if (!cisdense) {
2629: PetscBool flg;
2631: PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg);
2632: MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE);
2633: }
2634: MatSetUp(C);
2635: return(0);
2636: }
2638: PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2639: {
2640: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2641: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
2642: Mat_SeqDense *c = (Mat_SeqDense*)C->data;
2643: const PetscScalar *av,*bv;
2644: PetscScalar *cv;
2645: PetscBLASInt m,n,k;
2646: PetscScalar _DOne=1.0,_DZero=0.0;
2647: PetscErrorCode ierr;
2650: PetscBLASIntCast(C->rmap->n,&m);
2651: PetscBLASIntCast(C->cmap->n,&n);
2652: PetscBLASIntCast(A->rmap->n,&k);
2653: if (!m || !n || !k) return(0);
2654: MatDenseGetArrayRead(A,&av);
2655: MatDenseGetArrayRead(B,&bv);
2656: MatDenseGetArrayWrite(C,&cv);
2657: PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda));
2658: MatDenseRestoreArrayRead(A,&av);
2659: MatDenseRestoreArrayRead(B,&bv);
2660: MatDenseRestoreArrayWrite(C,&cv);
2661: PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1));
2662: return(0);
2663: }
2665: /* ----------------------------------------------- */
2666: static PetscErrorCode MatProductSetFromOptions_SeqDense_AB(Mat C)
2667: {
2669: C->ops->matmultsymbolic = MatMatMultSymbolic_SeqDense_SeqDense;
2670: C->ops->productsymbolic = MatProductSymbolic_AB;
2671: return(0);
2672: }
2674: static PetscErrorCode MatProductSetFromOptions_SeqDense_AtB(Mat C)
2675: {
2677: C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_SeqDense_SeqDense;
2678: C->ops->productsymbolic = MatProductSymbolic_AtB;
2679: return(0);
2680: }
2682: static PetscErrorCode MatProductSetFromOptions_SeqDense_ABt(Mat C)
2683: {
2685: C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_SeqDense_SeqDense;
2686: C->ops->productsymbolic = MatProductSymbolic_ABt;
2687: return(0);
2688: }
2690: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat C)
2691: {
2693: Mat_Product *product = C->product;
2696: switch (product->type) {
2697: case MATPRODUCT_AB:
2698: MatProductSetFromOptions_SeqDense_AB(C);
2699: break;
2700: case MATPRODUCT_AtB:
2701: MatProductSetFromOptions_SeqDense_AtB(C);
2702: break;
2703: case MATPRODUCT_ABt:
2704: MatProductSetFromOptions_SeqDense_ABt(C);
2705: break;
2706: default:
2707: break;
2708: }
2709: return(0);
2710: }
2711: /* ----------------------------------------------- */
2713: static PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
2714: {
2715: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2716: PetscErrorCode ierr;
2717: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2718: PetscScalar *x;
2719: const PetscScalar *aa;
2722: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2723: VecGetArray(v,&x);
2724: VecGetLocalSize(v,&p);
2725: MatDenseGetArrayRead(A,&aa);
2726: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2727: for (i=0; i<m; i++) {
2728: x[i] = aa[i]; if (idx) idx[i] = 0;
2729: for (j=1; j<n; j++) {
2730: if (PetscRealPart(x[i]) < PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;}
2731: }
2732: }
2733: MatDenseRestoreArrayRead(A,&aa);
2734: VecRestoreArray(v,&x);
2735: return(0);
2736: }
2738: static PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
2739: {
2740: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2741: PetscErrorCode ierr;
2742: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2743: PetscScalar *x;
2744: PetscReal atmp;
2745: const PetscScalar *aa;
2748: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2749: VecGetArray(v,&x);
2750: VecGetLocalSize(v,&p);
2751: MatDenseGetArrayRead(A,&aa);
2752: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2753: for (i=0; i<m; i++) {
2754: x[i] = PetscAbsScalar(aa[i]);
2755: for (j=1; j<n; j++) {
2756: atmp = PetscAbsScalar(aa[i+a->lda*j]);
2757: if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
2758: }
2759: }
2760: MatDenseRestoreArrayRead(A,&aa);
2761: VecRestoreArray(v,&x);
2762: return(0);
2763: }
2765: static PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
2766: {
2767: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2768: PetscErrorCode ierr;
2769: PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p;
2770: PetscScalar *x;
2771: const PetscScalar *aa;
2774: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2775: MatDenseGetArrayRead(A,&aa);
2776: VecGetArray(v,&x);
2777: VecGetLocalSize(v,&p);
2778: if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2779: for (i=0; i<m; i++) {
2780: x[i] = aa[i]; if (idx) idx[i] = 0;
2781: for (j=1; j<n; j++) {
2782: if (PetscRealPart(x[i]) > PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;}
2783: }
2784: }
2785: VecRestoreArray(v,&x);
2786: MatDenseRestoreArrayRead(A,&aa);
2787: return(0);
2788: }
2790: PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
2791: {
2792: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2793: PetscErrorCode ierr;
2794: PetscScalar *x;
2795: const PetscScalar *aa;
2798: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2799: MatDenseGetArrayRead(A,&aa);
2800: VecGetArray(v,&x);
2801: PetscArraycpy(x,aa+col*a->lda,A->rmap->n);
2802: VecRestoreArray(v,&x);
2803: MatDenseRestoreArrayRead(A,&aa);
2804: return(0);
2805: }
2807: PETSC_INTERN PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
2808: {
2809: PetscErrorCode ierr;
2810: PetscInt i,j,m,n;
2811: const PetscScalar *a;
2814: MatGetSize(A,&m,&n);
2815: PetscArrayzero(norms,n);
2816: MatDenseGetArrayRead(A,&a);
2817: if (type == NORM_2) {
2818: for (i=0; i<n; i++) {
2819: for (j=0; j<m; j++) {
2820: norms[i] += PetscAbsScalar(a[j]*a[j]);
2821: }
2822: a += m;
2823: }
2824: } else if (type == NORM_1) {
2825: for (i=0; i<n; i++) {
2826: for (j=0; j<m; j++) {
2827: norms[i] += PetscAbsScalar(a[j]);
2828: }
2829: a += m;
2830: }
2831: } else if (type == NORM_INFINITY) {
2832: for (i=0; i<n; i++) {
2833: for (j=0; j<m; j++) {
2834: norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
2835: }
2836: a += m;
2837: }
2838: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
2839: MatDenseRestoreArrayRead(A,&a);
2840: if (type == NORM_2) {
2841: for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
2842: }
2843: return(0);
2844: }
2846: static PetscErrorCode MatSetRandom_SeqDense(Mat x,PetscRandom rctx)
2847: {
2849: PetscScalar *a;
2850: PetscInt lda,m,n,i,j;
2853: MatGetSize(x,&m,&n);
2854: MatDenseGetLDA(x,&lda);
2855: MatDenseGetArray(x,&a);
2856: for (j=0; j<n; j++) {
2857: for (i=0; i<m; i++) {
2858: PetscRandomGetValue(rctx,a+j*lda+i);
2859: }
2860: }
2861: MatDenseRestoreArray(x,&a);
2862: return(0);
2863: }
2865: static PetscErrorCode MatMissingDiagonal_SeqDense(Mat A,PetscBool *missing,PetscInt *d)
2866: {
2868: *missing = PETSC_FALSE;
2869: return(0);
2870: }
2872: /* vals is not const */
2873: static PetscErrorCode MatDenseGetColumn_SeqDense(Mat A,PetscInt col,PetscScalar **vals)
2874: {
2876: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2877: PetscScalar *v;
2880: if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2881: MatDenseGetArray(A,&v);
2882: *vals = v+col*a->lda;
2883: MatDenseRestoreArray(A,&v);
2884: return(0);
2885: }
2887: static PetscErrorCode MatDenseRestoreColumn_SeqDense(Mat A,PetscScalar **vals)
2888: {
2890: *vals = NULL; /* user cannot accidently use the array later */
2891: return(0);
2892: }
2894: /* -------------------------------------------------------------------*/
2895: static struct _MatOps MatOps_Values = { MatSetValues_SeqDense,
2896: MatGetRow_SeqDense,
2897: MatRestoreRow_SeqDense,
2898: MatMult_SeqDense,
2899: /* 4*/ MatMultAdd_SeqDense,
2900: MatMultTranspose_SeqDense,
2901: MatMultTransposeAdd_SeqDense,
2902: NULL,
2903: NULL,
2904: NULL,
2905: /* 10*/ NULL,
2906: MatLUFactor_SeqDense,
2907: MatCholeskyFactor_SeqDense,
2908: MatSOR_SeqDense,
2909: MatTranspose_SeqDense,
2910: /* 15*/ MatGetInfo_SeqDense,
2911: MatEqual_SeqDense,
2912: MatGetDiagonal_SeqDense,
2913: MatDiagonalScale_SeqDense,
2914: MatNorm_SeqDense,
2915: /* 20*/ MatAssemblyBegin_SeqDense,
2916: MatAssemblyEnd_SeqDense,
2917: MatSetOption_SeqDense,
2918: MatZeroEntries_SeqDense,
2919: /* 24*/ MatZeroRows_SeqDense,
2920: NULL,
2921: NULL,
2922: NULL,
2923: NULL,
2924: /* 29*/ MatSetUp_SeqDense,
2925: NULL,
2926: NULL,
2927: NULL,
2928: NULL,
2929: /* 34*/ MatDuplicate_SeqDense,
2930: NULL,
2931: NULL,
2932: NULL,
2933: NULL,
2934: /* 39*/ MatAXPY_SeqDense,
2935: MatCreateSubMatrices_SeqDense,
2936: NULL,
2937: MatGetValues_SeqDense,
2938: MatCopy_SeqDense,
2939: /* 44*/ MatGetRowMax_SeqDense,
2940: MatScale_SeqDense,
2941: MatShift_Basic,
2942: NULL,
2943: MatZeroRowsColumns_SeqDense,
2944: /* 49*/ MatSetRandom_SeqDense,
2945: NULL,
2946: NULL,
2947: NULL,
2948: NULL,
2949: /* 54*/ NULL,
2950: NULL,
2951: NULL,
2952: NULL,
2953: NULL,
2954: /* 59*/ MatCreateSubMatrix_SeqDense,
2955: MatDestroy_SeqDense,
2956: MatView_SeqDense,
2957: NULL,
2958: NULL,
2959: /* 64*/ NULL,
2960: NULL,
2961: NULL,
2962: NULL,
2963: NULL,
2964: /* 69*/ MatGetRowMaxAbs_SeqDense,
2965: NULL,
2966: NULL,
2967: NULL,
2968: NULL,
2969: /* 74*/ NULL,
2970: NULL,
2971: NULL,
2972: NULL,
2973: NULL,
2974: /* 79*/ NULL,
2975: NULL,
2976: NULL,
2977: NULL,
2978: /* 83*/ MatLoad_SeqDense,
2979: MatIsSymmetric_SeqDense,
2980: MatIsHermitian_SeqDense,
2981: NULL,
2982: NULL,
2983: NULL,
2984: /* 89*/ NULL,
2985: NULL,
2986: MatMatMultNumeric_SeqDense_SeqDense,
2987: NULL,
2988: NULL,
2989: /* 94*/ NULL,
2990: NULL,
2991: NULL,
2992: MatMatTransposeMultNumeric_SeqDense_SeqDense,
2993: NULL,
2994: /* 99*/ MatProductSetFromOptions_SeqDense,
2995: NULL,
2996: NULL,
2997: MatConjugate_SeqDense,
2998: NULL,
2999: /*104*/ NULL,
3000: MatRealPart_SeqDense,
3001: MatImaginaryPart_SeqDense,
3002: NULL,
3003: NULL,
3004: /*109*/ NULL,
3005: NULL,
3006: MatGetRowMin_SeqDense,
3007: MatGetColumnVector_SeqDense,
3008: MatMissingDiagonal_SeqDense,
3009: /*114*/ NULL,
3010: NULL,
3011: NULL,
3012: NULL,
3013: NULL,
3014: /*119*/ NULL,
3015: NULL,
3016: NULL,
3017: NULL,
3018: NULL,
3019: /*124*/ NULL,
3020: MatGetColumnNorms_SeqDense,
3021: NULL,
3022: NULL,
3023: NULL,
3024: /*129*/ NULL,
3025: NULL,
3026: NULL,
3027: MatTransposeMatMultNumeric_SeqDense_SeqDense,
3028: NULL,
3029: /*134*/ NULL,
3030: NULL,
3031: NULL,
3032: NULL,
3033: NULL,
3034: /*139*/ NULL,
3035: NULL,
3036: NULL,
3037: NULL,
3038: NULL,
3039: MatCreateMPIMatConcatenateSeqMat_SeqDense,
3040: /*145*/ NULL,
3041: NULL,
3042: NULL
3043: };
3045: /*@C
3046: MatCreateSeqDense - Creates a sequential dense matrix that
3047: is stored in column major order (the usual Fortran 77 manner). Many
3048: of the matrix operations use the BLAS and LAPACK routines.
3050: Collective
3052: Input Parameters:
3053: + comm - MPI communicator, set to PETSC_COMM_SELF
3054: . m - number of rows
3055: . n - number of columns
3056: - data - optional location of matrix data in column major order. Set data=NULL for PETSc
3057: to control all matrix memory allocation.
3059: Output Parameter:
3060: . A - the matrix
3062: Notes:
3063: The data input variable is intended primarily for Fortran programmers
3064: who wish to allocate their own matrix memory space. Most users should
3065: set data=NULL.
3067: Level: intermediate
3069: .seealso: MatCreate(), MatCreateDense(), MatSetValues()
3070: @*/
3071: PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
3072: {
3076: MatCreate(comm,A);
3077: MatSetSizes(*A,m,n,m,n);
3078: MatSetType(*A,MATSEQDENSE);
3079: MatSeqDenseSetPreallocation(*A,data);
3080: return(0);
3081: }
3083: /*@C
3084: MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
3086: Collective
3088: Input Parameters:
3089: + B - the matrix
3090: - data - the array (or NULL)
3092: Notes:
3093: The data input variable is intended primarily for Fortran programmers
3094: who wish to allocate their own matrix memory space. Most users should
3095: need not call this routine.
3097: Level: intermediate
3099: .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatDenseSetLDA()
3101: @*/
3102: PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
3103: {
3108: PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));
3109: return(0);
3110: }
3112: PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
3113: {
3114: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
3118: if (b->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
3119: B->preallocated = PETSC_TRUE;
3121: PetscLayoutSetUp(B->rmap);
3122: PetscLayoutSetUp(B->cmap);
3124: if (b->lda <= 0) b->lda = B->rmap->n;
3126: PetscIntMultError(b->lda,B->cmap->n,NULL);
3127: if (!data) { /* petsc-allocated storage */
3128: if (!b->user_alloc) { PetscFree(b->v); }
3129: PetscCalloc1((size_t)b->lda*B->cmap->n,&b->v);
3130: PetscLogObjectMemory((PetscObject)B,b->lda*B->cmap->n*sizeof(PetscScalar));
3132: b->user_alloc = PETSC_FALSE;
3133: } else { /* user-allocated storage */
3134: if (!b->user_alloc) { PetscFree(b->v); }
3135: b->v = data;
3136: b->user_alloc = PETSC_TRUE;
3137: }
3138: B->assembled = PETSC_TRUE;
3139: return(0);
3140: }
3142: #if defined(PETSC_HAVE_ELEMENTAL)
3143: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
3144: {
3145: Mat mat_elemental;
3146: PetscErrorCode ierr;
3147: const PetscScalar *array;
3148: PetscScalar *v_colwise;
3149: PetscInt M=A->rmap->N,N=A->cmap->N,i,j,k,*rows,*cols;
3152: PetscMalloc3(M*N,&v_colwise,M,&rows,N,&cols);
3153: MatDenseGetArrayRead(A,&array);
3154: /* convert column-wise array into row-wise v_colwise, see MatSetValues_Elemental() */
3155: k = 0;
3156: for (j=0; j<N; j++) {
3157: cols[j] = j;
3158: for (i=0; i<M; i++) {
3159: v_colwise[j*M+i] = array[k++];
3160: }
3161: }
3162: for (i=0; i<M; i++) {
3163: rows[i] = i;
3164: }
3165: MatDenseRestoreArrayRead(A,&array);
3167: MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
3168: MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
3169: MatSetType(mat_elemental,MATELEMENTAL);
3170: MatSetUp(mat_elemental);
3172: /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
3173: MatSetValues(mat_elemental,M,rows,N,cols,v_colwise,ADD_VALUES);
3174: MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
3175: MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
3176: PetscFree3(v_colwise,rows,cols);
3178: if (reuse == MAT_INPLACE_MATRIX) {
3179: MatHeaderReplace(A,&mat_elemental);
3180: } else {
3181: *newmat = mat_elemental;
3182: }
3183: return(0);
3184: }
3185: #endif
3187: PetscErrorCode MatDenseSetLDA_SeqDense(Mat B,PetscInt lda)
3188: {
3189: Mat_SeqDense *b = (Mat_SeqDense*)B->data;
3190: PetscBool data;
3193: data = (PetscBool)((B->rmap->n > 0 && B->cmap->n > 0) ? (b->v ? PETSC_TRUE : PETSC_FALSE) : PETSC_FALSE);
3194: if (!b->user_alloc && data && b->lda!=lda) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"LDA cannot be changed after allocation of internal storage");
3195: 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);
3196: b->lda = lda;
3197: return(0);
3198: }
3200: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
3201: {
3203: PetscMPIInt size;
3206: MPI_Comm_size(comm,&size);
3207: if (size == 1) {
3208: if (scall == MAT_INITIAL_MATRIX) {
3209: MatDuplicate(inmat,MAT_COPY_VALUES,outmat);
3210: } else {
3211: MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);
3212: }
3213: } else {
3214: MatCreateMPIMatConcatenateSeqMat_MPIDense(comm,inmat,n,scall,outmat);
3215: }
3216: return(0);
3217: }
3219: PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat A,PetscInt col,Vec *v)
3220: {
3221: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3225: if (a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
3226: if (a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
3227: if (!a->cvec) {
3228: VecCreateSeqWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,NULL,&a->cvec);
3229: PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec);
3230: }
3231: a->vecinuse = col + 1;
3232: MatDenseGetArray(A,(PetscScalar**)&a->ptrinuse);
3233: VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)a->lda);
3234: *v = a->cvec;
3235: return(0);
3236: }
3238: PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat A,PetscInt col,Vec *v)
3239: {
3240: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3244: if (!a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
3245: if (!a->cvec) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector");
3246: a->vecinuse = 0;
3247: MatDenseRestoreArray(A,(PetscScalar**)&a->ptrinuse);
3248: VecResetArray(a->cvec);
3249: *v = NULL;
3250: return(0);
3251: }
3253: PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat A,PetscInt col,Vec *v)
3254: {
3255: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3259: if (a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
3260: if (a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
3261: if (!a->cvec) {
3262: VecCreateSeqWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,NULL,&a->cvec);
3263: PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec);
3264: }
3265: a->vecinuse = col + 1;
3266: MatDenseGetArrayRead(A,&a->ptrinuse);
3267: VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)a->lda);
3268: VecLockReadPush(a->cvec);
3269: *v = a->cvec;
3270: return(0);
3271: }
3273: PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat A,PetscInt col,Vec *v)
3274: {
3275: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3279: if (!a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
3280: if (!a->cvec) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector");
3281: a->vecinuse = 0;
3282: MatDenseRestoreArrayRead(A,&a->ptrinuse);
3283: VecLockReadPop(a->cvec);
3284: VecResetArray(a->cvec);
3285: *v = NULL;
3286: return(0);
3287: }
3289: PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat A,PetscInt col,Vec *v)
3290: {
3291: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3295: if (a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
3296: if (a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
3297: if (!a->cvec) {
3298: VecCreateSeqWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,NULL,&a->cvec);
3299: PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec);
3300: }
3301: a->vecinuse = col + 1;
3302: MatDenseGetArrayWrite(A,(PetscScalar**)&a->ptrinuse);
3303: VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)a->lda);
3304: *v = a->cvec;
3305: return(0);
3306: }
3308: PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat A,PetscInt col,Vec *v)
3309: {
3310: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3314: if (!a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
3315: if (!a->cvec) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector");
3316: a->vecinuse = 0;
3317: MatDenseRestoreArrayWrite(A,(PetscScalar**)&a->ptrinuse);
3318: VecResetArray(a->cvec);
3319: *v = NULL;
3320: return(0);
3321: }
3323: PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat A,PetscInt cbegin,PetscInt cend,Mat *v)
3324: {
3325: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3329: if (a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
3330: if (a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
3331: if (a->cmat && cend-cbegin != a->cmat->cmap->N) {
3332: MatDestroy(&a->cmat);
3333: }
3334: if (!a->cmat) {
3335: MatCreateDense(PetscObjectComm((PetscObject)A),A->rmap->n,PETSC_DECIDE,A->rmap->N,cend-cbegin,(PetscScalar*)a->v + (size_t)cbegin * (size_t)a->lda,&a->cmat);
3336: PetscLogObjectParent((PetscObject)A,(PetscObject)a->cmat);
3337: } else {
3338: MatDensePlaceArray(a->cmat,a->v + (size_t)cbegin * (size_t)a->lda);
3339: }
3340: MatDenseSetLDA(a->cmat,a->lda);
3341: a->matinuse = cbegin + 1;
3342: *v = a->cmat;
3343: return(0);
3344: }
3346: PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat A,Mat *v)
3347: {
3348: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3352: if (!a->matinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetSubMatrix() first");
3353: if (!a->cmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column matrix");
3354: if (*v != a->cmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not the matrix obtained from MatDenseGetSubMatrix()");
3355: a->matinuse = 0;
3356: MatDenseResetArray(a->cmat);
3357: *v = NULL;
3358: return(0);
3359: }
3361: /*MC
3362: MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
3364: Options Database Keys:
3365: . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
3367: Level: beginner
3369: .seealso: MatCreateSeqDense()
3371: M*/
3372: PetscErrorCode MatCreate_SeqDense(Mat B)
3373: {
3374: Mat_SeqDense *b;
3376: PetscMPIInt size;
3379: MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
3380: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
3382: PetscNewLog(B,&b);
3383: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
3384: B->data = (void*)b;
3386: b->roworiented = PETSC_TRUE;
3388: PetscObjectComposeFunction((PetscObject)B,"MatQRFactor_C",MatQRFactor_SeqDense);
3389: PetscObjectComposeFunction((PetscObject)B,"MatQRFactorNumeric_C",MatQRFactorNumeric_SeqDense);
3390: PetscObjectComposeFunction((PetscObject)B,"MatQRFactorSymbolic_C",MatQRFactorSymbolic_SeqDense);
3391: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetLDA_C",MatDenseGetLDA_SeqDense);
3392: PetscObjectComposeFunction((PetscObject)B,"MatDenseSetLDA_C",MatDenseSetLDA_SeqDense);
3393: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense);
3394: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense);
3395: PetscObjectComposeFunction((PetscObject)B,"MatDensePlaceArray_C",MatDensePlaceArray_SeqDense);
3396: PetscObjectComposeFunction((PetscObject)B,"MatDenseResetArray_C",MatDenseResetArray_SeqDense);
3397: PetscObjectComposeFunction((PetscObject)B,"MatDenseReplaceArray_C",MatDenseReplaceArray_SeqDense);
3398: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayRead_C",MatDenseGetArray_SeqDense);
3399: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayRead_C",MatDenseRestoreArray_SeqDense);
3400: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayWrite_C",MatDenseGetArray_SeqDense);
3401: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayWrite_C",MatDenseRestoreArray_SeqDense);
3402: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ);
3403: #if defined(PETSC_HAVE_ELEMENTAL)
3404: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_elemental_C",MatConvert_SeqDense_Elemental);
3405: #endif
3406: #if defined(PETSC_HAVE_SCALAPACK)
3407: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_scalapack_C",MatConvert_Dense_ScaLAPACK);
3408: #endif
3409: #if defined(PETSC_HAVE_CUDA)
3410: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqdensecuda_C",MatConvert_SeqDense_SeqDenseCUDA);
3411: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdensecuda_seqdensecuda_C",MatProductSetFromOptions_SeqDense);
3412: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdensecuda_seqdense_C",MatProductSetFromOptions_SeqDense);
3413: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdense_seqdensecuda_C",MatProductSetFromOptions_SeqDense);
3414: #endif
3415: PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense);
3416: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqaij_seqdense_C",MatProductSetFromOptions_SeqAIJ_SeqDense);
3417: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdense_seqdense_C",MatProductSetFromOptions_SeqDense);
3418: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqbaij_seqdense_C",MatProductSetFromOptions_SeqXBAIJ_SeqDense);
3419: PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqsbaij_seqdense_C",MatProductSetFromOptions_SeqXBAIJ_SeqDense);
3421: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumn_C",MatDenseGetColumn_SeqDense);
3422: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_SeqDense);
3423: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumnVec_C",MatDenseGetColumnVec_SeqDense);
3424: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumnVec_C",MatDenseRestoreColumnVec_SeqDense);
3425: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumnVecRead_C",MatDenseGetColumnVecRead_SeqDense);
3426: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumnVecRead_C",MatDenseRestoreColumnVecRead_SeqDense);
3427: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumnVecWrite_C",MatDenseGetColumnVecWrite_SeqDense);
3428: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumnVecWrite_C",MatDenseRestoreColumnVecWrite_SeqDense);
3429: PetscObjectComposeFunction((PetscObject)B,"MatDenseGetSubMatrix_C",MatDenseGetSubMatrix_SeqDense);
3430: PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreSubMatrix_C",MatDenseRestoreSubMatrix_SeqDense);
3431: PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);
3432: return(0);
3433: }
3435: /*@C
3436: 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.
3438: Not Collective
3440: Input Parameters:
3441: + mat - a MATSEQDENSE or MATMPIDENSE matrix
3442: - col - column index
3444: Output Parameter:
3445: . vals - pointer to the data
3447: Level: intermediate
3449: .seealso: MatDenseRestoreColumn()
3450: @*/
3451: PetscErrorCode MatDenseGetColumn(Mat A,PetscInt col,PetscScalar **vals)
3452: {
3459: PetscUseMethod(A,"MatDenseGetColumn_C",(Mat,PetscInt,PetscScalar**),(A,col,vals));
3460: return(0);
3461: }
3463: /*@C
3464: MatDenseRestoreColumn - returns access to a column of a dense matrix which is returned by MatDenseGetColumn().
3466: Not Collective
3468: Input Parameter:
3469: . mat - a MATSEQDENSE or MATMPIDENSE matrix
3471: Output Parameter:
3472: . vals - pointer to the data
3474: Level: intermediate
3476: .seealso: MatDenseGetColumn()
3477: @*/
3478: PetscErrorCode MatDenseRestoreColumn(Mat A,PetscScalar **vals)
3479: {
3485: PetscUseMethod(A,"MatDenseRestoreColumn_C",(Mat,PetscScalar**),(A,vals));
3486: return(0);
3487: }
3489: /*@C
3490: MatDenseGetColumnVec - Gives read-write access to a column of a dense matrix, represented as a Vec.
3492: Collective
3494: Input Parameters:
3495: + mat - the Mat object
3496: - col - the column index
3498: Output Parameter:
3499: . v - the vector
3501: Notes:
3502: The vector is owned by PETSc. Users need to call MatDenseRestoreColumnVec() when the vector is no longer needed.
3503: Use MatDenseGetColumnVecRead() to obtain read-only access or MatDenseGetColumnVecWrite() for write-only access.
3505: Level: intermediate
3507: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVecRead(), MatDenseGetColumnVecWrite(), MatDenseRestoreColumnVec(), MatDenseRestoreColumnVecRead(), MatDenseRestoreColumnVecWrite()
3508: @*/
3509: PetscErrorCode MatDenseGetColumnVec(Mat A,PetscInt col,Vec *v)
3510: {
3518: if (!A->preallocated) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated");
3519: if (col < 0 || col > A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %D, should be in [0,%D)",col,A->cmap->N);
3520: PetscUseMethod(A,"MatDenseGetColumnVec_C",(Mat,PetscInt,Vec*),(A,col,v));
3521: return(0);
3522: }
3524: /*@C
3525: MatDenseRestoreColumnVec - Returns access to a column of a dense matrix obtained from MatDenseGetColumnVec().
3527: Collective
3529: Input Parameters:
3530: + mat - the Mat object
3531: . col - the column index
3532: - v - the Vec object
3534: Level: intermediate
3536: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVec(), MatDenseGetColumnVecRead(), MatDenseGetColumnVecWrite(), MatDenseRestoreColumnVecRead(), MatDenseRestoreColumnVecWrite()
3537: @*/
3538: PetscErrorCode MatDenseRestoreColumnVec(Mat A,PetscInt col,Vec *v)
3539: {
3547: if (!A->preallocated) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated");
3548: if (col < 0 || col > A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %D, should be in [0,%D)",col,A->cmap->N);
3549: PetscUseMethod(A,"MatDenseRestoreColumnVec_C",(Mat,PetscInt,Vec*),(A,col,v));
3550: return(0);
3551: }
3553: /*@C
3554: MatDenseGetColumnVecRead - Gives read-only access to a column of a dense matrix, represented as a Vec.
3556: Collective
3558: Input Parameters:
3559: + mat - the Mat object
3560: - col - the column index
3562: Output Parameter:
3563: . v - the vector
3565: Notes:
3566: The vector is owned by PETSc and users cannot modify it.
3567: Users need to call MatDenseRestoreColumnVecRead() when the vector is no longer needed.
3568: Use MatDenseGetColumnVec() to obtain read-write access or MatDenseGetColumnVecWrite() for write-only access.
3570: Level: intermediate
3572: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVec(), MatDenseGetColumnVecWrite(), MatDenseRestoreColumnVec(), MatDenseRestoreColumnVecRead(), MatDenseRestoreColumnVecWrite()
3573: @*/
3574: PetscErrorCode MatDenseGetColumnVecRead(Mat A,PetscInt col,Vec *v)
3575: {
3583: if (!A->preallocated) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated");
3584: if (col < 0 || col > A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %D, should be in [0,%D)",col,A->cmap->N);
3585: PetscUseMethod(A,"MatDenseGetColumnVecRead_C",(Mat,PetscInt,Vec*),(A,col,v));
3586: return(0);
3587: }
3589: /*@C
3590: MatDenseRestoreColumnVecRead - Returns access to a column of a dense matrix obtained from MatDenseGetColumnVecRead().
3592: Collective
3594: Input Parameters:
3595: + mat - the Mat object
3596: . col - the column index
3597: - v - the Vec object
3599: Level: intermediate
3601: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVec(), MatDenseGetColumnVecRead(), MatDenseGetColumnVecWrite(), MatDenseRestoreColumnVec(), MatDenseRestoreColumnVecWrite()
3602: @*/
3603: PetscErrorCode MatDenseRestoreColumnVecRead(Mat A,PetscInt col,Vec *v)
3604: {
3612: if (!A->preallocated) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated");
3613: if (col < 0 || col > A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %D, should be in [0,%D)",col,A->cmap->N);
3614: PetscUseMethod(A,"MatDenseRestoreColumnVecRead_C",(Mat,PetscInt,Vec*),(A,col,v));
3615: return(0);
3616: }
3618: /*@C
3619: MatDenseGetColumnVecWrite - Gives write-only access to a column of a dense matrix, represented as a Vec.
3621: Collective
3623: Input Parameters:
3624: + mat - the Mat object
3625: - col - the column index
3627: Output Parameter:
3628: . v - the vector
3630: Notes:
3631: The vector is owned by PETSc. Users need to call MatDenseRestoreColumnVecWrite() when the vector is no longer needed.
3632: Use MatDenseGetColumnVec() to obtain read-write access or MatDenseGetColumnVecRead() for read-only access.
3634: Level: intermediate
3636: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVec(), MatDenseGetColumnVecRead(), MatDenseRestoreColumnVec(), MatDenseRestoreColumnVecRead(), MatDenseRestoreColumnVecWrite()
3637: @*/
3638: PetscErrorCode MatDenseGetColumnVecWrite(Mat A,PetscInt col,Vec *v)
3639: {
3647: if (!A->preallocated) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated");
3648: if (col < 0 || col > A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %D, should be in [0,%D)",col,A->cmap->N);
3649: PetscUseMethod(A,"MatDenseGetColumnVecWrite_C",(Mat,PetscInt,Vec*),(A,col,v));
3650: return(0);
3651: }
3653: /*@C
3654: MatDenseRestoreColumnVecWrite - Returns access to a column of a dense matrix obtained from MatDenseGetColumnVecWrite().
3656: Collective
3658: Input Parameters:
3659: + mat - the Mat object
3660: . col - the column index
3661: - v - the Vec object
3663: Level: intermediate
3665: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVec(), MatDenseGetColumnVecRead(), MatDenseGetColumnVecWrite(), MatDenseRestoreColumnVec(), MatDenseRestoreColumnVecRead()
3666: @*/
3667: PetscErrorCode MatDenseRestoreColumnVecWrite(Mat A,PetscInt col,Vec *v)
3668: {
3676: if (!A->preallocated) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated");
3677: if (col < 0 || col > A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %D, should be in [0,%D)",col,A->cmap->N);
3678: PetscUseMethod(A,"MatDenseRestoreColumnVecWrite_C",(Mat,PetscInt,Vec*),(A,col,v));
3679: return(0);
3680: }
3682: /*@C
3683: MatDenseGetSubMatrix - Gives access to a block of columns of a dense matrix, represented as a Mat.
3685: Collective
3687: Input Parameters:
3688: + mat - the Mat object
3689: . cbegin - the first index in the block
3690: - cend - the last index in the block
3692: Output Parameter:
3693: . v - the matrix
3695: Notes:
3696: The matrix is owned by PETSc. Users need to call MatDenseRestoreSubMatrix() when the matrix is no longer needed.
3698: Level: intermediate
3700: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVec(), MatDenseRestoreColumnVec(), MatDenseRestoreSubMatrix()
3701: @*/
3702: PetscErrorCode MatDenseGetSubMatrix(Mat A,PetscInt cbegin,PetscInt cend,Mat *v)
3703: {
3712: if (!A->preallocated) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated");
3713: if (cbegin < 0 || cbegin > A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid cbegin %D, should be in [0,%D)",cbegin,A->cmap->N);
3714: if (cend < cbegin || cend > A->cmap->N) SETERRQ3(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid cend %D, should be in [%D,%D)",cend,cbegin,A->cmap->N);
3715: PetscUseMethod(A,"MatDenseGetSubMatrix_C",(Mat,PetscInt,PetscInt,Mat*),(A,cbegin,cend,v));
3716: return(0);
3717: }
3719: /*@C
3720: MatDenseRestoreSubMatrix - Returns access to a block of columns of a dense matrix obtained from MatDenseGetSubMatrix().
3722: Collective
3724: Input Parameters:
3725: + mat - the Mat object
3726: - v - the Mat object
3728: Level: intermediate
3730: .seealso: MATDENSE, MATDENSECUDA, MatDenseGetColumnVec(), MatDenseRestoreColumnVec(), MatDenseGetSubMatrix()
3731: @*/
3732: PetscErrorCode MatDenseRestoreSubMatrix(Mat A,Mat *v)
3733: {
3740: PetscUseMethod(A,"MatDenseRestoreSubMatrix_C",(Mat,Mat*),(A,v));
3741: return(0);
3742: }