Actual source code: matelem.cxx
petsc-3.8.4 2018-03-24
1: #include <../src/mat/impls/elemental/matelemimpl.h>
3: /*
4: The variable Petsc_Elemental_keyval is used to indicate an MPI attribute that
5: is attached to a communicator, in this case the attribute is a Mat_Elemental_Grid
6: */
7: static PetscMPIInt Petsc_Elemental_keyval = MPI_KEYVAL_INVALID;
9: /*@C
10: PetscElementalInitializePackage - Initialize Elemental package
12: Logically Collective
14: Level: developer
16: .seealso: MATELEMENTAL, PetscElementalFinalizePackage()
17: @*/
18: PetscErrorCode PetscElementalInitializePackage(void)
19: {
23: if (El::Initialized()) return(0);
24: El::Initialize(); /* called by the 1st call of MatCreate_Elemental */
25: PetscRegisterFinalize(PetscElementalFinalizePackage);
26: return(0);
27: }
29: /*@C
30: PetscElementalFinalizePackage - Finalize Elemental package
32: Logically Collective
34: Level: developer
36: .seealso: MATELEMENTAL, PetscElementalInitializePackage()
37: @*/
38: PetscErrorCode PetscElementalFinalizePackage(void)
39: {
41: El::Finalize(); /* called by PetscFinalize() */
42: return(0);
43: }
45: static PetscErrorCode MatView_Elemental(Mat A,PetscViewer viewer)
46: {
48: Mat_Elemental *a = (Mat_Elemental*)A->data;
49: PetscBool iascii;
52: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
53: if (iascii) {
54: PetscViewerFormat format;
55: PetscViewerGetFormat(viewer,&format);
56: if (format == PETSC_VIEWER_ASCII_INFO) {
57: /* call elemental viewing function */
58: PetscViewerASCIIPrintf(viewer,"Elemental run parameters:\n");
59: PetscViewerASCIIPrintf(viewer," allocated entries=%d\n",(*a->emat).AllocatedMemory());
60: PetscViewerASCIIPrintf(viewer," grid height=%d, grid width=%d\n",(*a->emat).Grid().Height(),(*a->emat).Grid().Width());
61: if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
62: /* call elemental viewing function */
63: PetscPrintf(PetscObjectComm((PetscObject)viewer),"test matview_elemental 2\n");
64: }
66: } else if (format == PETSC_VIEWER_DEFAULT) {
67: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
68: El::Print( *a->emat, "Elemental matrix (cyclic ordering)" );
69: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
70: if (A->factortype == MAT_FACTOR_NONE){
71: Mat Adense;
72: PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");
73: MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);
74: MatView(Adense,viewer);
75: MatDestroy(&Adense);
76: }
77: } else SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Format");
78: } else {
79: /* convert to dense format and call MatView() */
80: Mat Adense;
81: PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");
82: MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);
83: MatView(Adense,viewer);
84: MatDestroy(&Adense);
85: }
86: return(0);
87: }
89: static PetscErrorCode MatGetInfo_Elemental(Mat A,MatInfoType flag,MatInfo *info)
90: {
91: Mat_Elemental *a = (Mat_Elemental*)A->data;
94: info->block_size = 1.0;
96: if (flag == MAT_LOCAL) {
97: info->nz_allocated = (double)(*a->emat).AllocatedMemory(); /* locally allocated */
98: info->nz_used = info->nz_allocated;
99: } else if (flag == MAT_GLOBAL_MAX) {
100: //MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));
101: /* see MatGetInfo_MPIAIJ() for getting global info->nz_allocated! */
102: //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_MAX not written yet");
103: } else if (flag == MAT_GLOBAL_SUM) {
104: //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_SUM not written yet");
105: info->nz_allocated = (double)(*a->emat).AllocatedMemory(); /* locally allocated */
106: info->nz_used = info->nz_allocated; /* assume Elemental does accurate allocation */
107: //MPIU_Allreduce(isend,irecv,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));
108: //PetscPrintf(PETSC_COMM_SELF," ... [%d] locally allocated %g\n",rank,info->nz_allocated);
109: }
111: info->nz_unneeded = 0.0;
112: info->assemblies = (double)A->num_ass;
113: info->mallocs = 0;
114: info->memory = ((PetscObject)A)->mem;
115: info->fill_ratio_given = 0; /* determined by Elemental */
116: info->fill_ratio_needed = 0;
117: info->factor_mallocs = 0;
118: return(0);
119: }
121: PetscErrorCode MatSetOption_Elemental(Mat A,MatOption op,PetscBool flg)
122: {
123: Mat_Elemental *a = (Mat_Elemental*)A->data;
126: switch (op) {
127: case MAT_NEW_NONZERO_LOCATIONS:
128: case MAT_NEW_NONZERO_LOCATION_ERR:
129: case MAT_NEW_NONZERO_ALLOCATION_ERR:
130: case MAT_SYMMETRIC:
131: break;
132: case MAT_ROW_ORIENTED:
133: a->roworiented = flg;
134: break;
135: default:
136: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
137: }
138: return(0);
139: }
141: static PetscErrorCode MatSetValues_Elemental(Mat A,PetscInt nr,const PetscInt *rows,PetscInt nc,const PetscInt *cols,const PetscScalar *vals,InsertMode imode)
142: {
143: Mat_Elemental *a = (Mat_Elemental*)A->data;
144: PetscInt i,j,rrank,ridx,crank,cidx,erow,ecol,numQueues=0;
147: // TODO: Initialize matrix to all zeros?
149: // Count the number of queues from this process
150: if (a->roworiented) {
151: for (i=0; i<nr; i++) {
152: if (rows[i] < 0) continue;
153: P2RO(A,0,rows[i],&rrank,&ridx);
154: RO2E(A,0,rrank,ridx,&erow);
155: if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation");
156: for (j=0; j<nc; j++) {
157: if (cols[j] < 0) continue;
158: P2RO(A,1,cols[j],&crank,&cidx);
159: RO2E(A,1,crank,cidx,&ecol);
160: if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation");
161: if (!a->emat->IsLocal(erow,ecol) ){ /* off-proc entry */
162: /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
163: if (imode != ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported");
164: ++numQueues;
165: continue;
166: }
167: /* printf("Locally updating (%d,%d)\n",erow,ecol); */
168: switch (imode) {
169: case INSERT_VALUES: a->emat->Set(erow,ecol,(PetscElemScalar)vals[i*nc+j]); break;
170: case ADD_VALUES: a->emat->Update(erow,ecol,(PetscElemScalar)vals[i*nc+j]); break;
171: default: SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode);
172: }
173: }
174: }
176: /* printf("numQueues=%d\n",numQueues); */
177: a->emat->Reserve( numQueues );
178: for (i=0; i<nr; i++) {
179: if (rows[i] < 0) continue;
180: P2RO(A,0,rows[i],&rrank,&ridx);
181: RO2E(A,0,rrank,ridx,&erow);
182: for (j=0; j<nc; j++) {
183: if (cols[j] < 0) continue;
184: P2RO(A,1,cols[j],&crank,&cidx);
185: RO2E(A,1,crank,cidx,&ecol);
186: if ( !a->emat->IsLocal(erow,ecol) ) { /*off-proc entry*/
187: /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
188: a->emat->QueueUpdate( erow, ecol, vals[i*nc+j] );
189: }
190: }
191: }
192: } else { /* columnoriented */
193: for (j=0; j<nc; j++) {
194: if (cols[j] < 0) continue;
195: P2RO(A,1,cols[j],&crank,&cidx);
196: RO2E(A,1,crank,cidx,&ecol);
197: if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation");
198: for (i=0; i<nr; i++) {
199: if (rows[i] < 0) continue;
200: P2RO(A,0,rows[i],&rrank,&ridx);
201: RO2E(A,0,rrank,ridx,&erow);
202: if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation");
203: if (!a->emat->IsLocal(erow,ecol) ){ /* off-proc entry */
204: /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
205: if (imode != ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported");
206: ++numQueues;
207: continue;
208: }
209: /* printf("Locally updating (%d,%d)\n",erow,ecol); */
210: switch (imode) {
211: case INSERT_VALUES: a->emat->Set(erow,ecol,(PetscElemScalar)vals[i+j*nr]); break;
212: case ADD_VALUES: a->emat->Update(erow,ecol,(PetscElemScalar)vals[i+j*nr]); break;
213: default: SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode);
214: }
215: }
216: }
218: /* printf("numQueues=%d\n",numQueues); */
219: a->emat->Reserve( numQueues );
220: for (j=0; j<nc; j++) {
221: if (cols[j] < 0) continue;
222: P2RO(A,1,cols[j],&crank,&cidx);
223: RO2E(A,1,crank,cidx,&ecol);
225: for (i=0; i<nr; i++) {
226: if (rows[i] < 0) continue;
227: P2RO(A,0,rows[i],&rrank,&ridx);
228: RO2E(A,0,rrank,ridx,&erow);
229: if ( !a->emat->IsLocal(erow,ecol) ) { /*off-proc entry*/
230: /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
231: a->emat->QueueUpdate( erow, ecol, vals[i+j*nr] );
232: }
233: }
234: }
235: }
236: return(0);
237: }
239: static PetscErrorCode MatMult_Elemental(Mat A,Vec X,Vec Y)
240: {
241: Mat_Elemental *a = (Mat_Elemental*)A->data;
242: PetscErrorCode ierr;
243: const PetscElemScalar *x;
244: PetscElemScalar *y;
245: PetscElemScalar one = 1,zero = 0;
248: VecGetArrayRead(X,(const PetscScalar **)&x);
249: VecGetArray(Y,(PetscScalar **)&y);
250: { /* Scoping so that constructor is called before pointer is returned */
251: El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
252: xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
253: ye.Attach(A->rmap->N,1,*a->grid,0,0,y,A->rmap->n);
254: El::Gemv(El::NORMAL,one,*a->emat,xe,zero,ye);
255: }
256: VecRestoreArrayRead(X,(const PetscScalar **)&x);
257: VecRestoreArray(Y,(PetscScalar **)&y);
258: return(0);
259: }
261: static PetscErrorCode MatMultTranspose_Elemental(Mat A,Vec X,Vec Y)
262: {
263: Mat_Elemental *a = (Mat_Elemental*)A->data;
264: PetscErrorCode ierr;
265: const PetscElemScalar *x;
266: PetscElemScalar *y;
267: PetscElemScalar one = 1,zero = 0;
270: VecGetArrayRead(X,(const PetscScalar **)&x);
271: VecGetArray(Y,(PetscScalar **)&y);
272: { /* Scoping so that constructor is called before pointer is returned */
273: El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
274: xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
275: ye.Attach(A->cmap->N,1,*a->grid,0,0,y,A->cmap->n);
276: El::Gemv(El::TRANSPOSE,one,*a->emat,xe,zero,ye);
277: }
278: VecRestoreArrayRead(X,(const PetscScalar **)&x);
279: VecRestoreArray(Y,(PetscScalar **)&y);
280: return(0);
281: }
283: static PetscErrorCode MatMultAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
284: {
285: Mat_Elemental *a = (Mat_Elemental*)A->data;
286: PetscErrorCode ierr;
287: const PetscElemScalar *x;
288: PetscElemScalar *z;
289: PetscElemScalar one = 1;
292: if (Y != Z) {VecCopy(Y,Z);}
293: VecGetArrayRead(X,(const PetscScalar **)&x);
294: VecGetArray(Z,(PetscScalar **)&z);
295: { /* Scoping so that constructor is called before pointer is returned */
296: El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
297: xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
298: ze.Attach(A->rmap->N,1,*a->grid,0,0,z,A->rmap->n);
299: El::Gemv(El::NORMAL,one,*a->emat,xe,one,ze);
300: }
301: VecRestoreArrayRead(X,(const PetscScalar **)&x);
302: VecRestoreArray(Z,(PetscScalar **)&z);
303: return(0);
304: }
306: static PetscErrorCode MatMultTransposeAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
307: {
308: Mat_Elemental *a = (Mat_Elemental*)A->data;
309: PetscErrorCode ierr;
310: const PetscElemScalar *x;
311: PetscElemScalar *z;
312: PetscElemScalar one = 1;
315: if (Y != Z) {VecCopy(Y,Z);}
316: VecGetArrayRead(X,(const PetscScalar **)&x);
317: VecGetArray(Z,(PetscScalar **)&z);
318: { /* Scoping so that constructor is called before pointer is returned */
319: El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
320: xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
321: ze.Attach(A->cmap->N,1,*a->grid,0,0,z,A->cmap->n);
322: El::Gemv(El::TRANSPOSE,one,*a->emat,xe,one,ze);
323: }
324: VecRestoreArrayRead(X,(const PetscScalar **)&x);
325: VecRestoreArray(Z,(PetscScalar **)&z);
326: return(0);
327: }
329: static PetscErrorCode MatMatMultNumeric_Elemental(Mat A,Mat B,Mat C)
330: {
331: Mat_Elemental *a = (Mat_Elemental*)A->data;
332: Mat_Elemental *b = (Mat_Elemental*)B->data;
333: Mat_Elemental *c = (Mat_Elemental*)C->data;
334: PetscElemScalar one = 1,zero = 0;
337: { /* Scoping so that constructor is called before pointer is returned */
338: El::Gemm(El::NORMAL,El::NORMAL,one,*a->emat,*b->emat,zero,*c->emat);
339: }
340: C->assembled = PETSC_TRUE;
341: return(0);
342: }
344: static PetscErrorCode MatMatMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C)
345: {
347: Mat Ce;
348: MPI_Comm comm;
351: PetscObjectGetComm((PetscObject)A,&comm);
352: MatCreate(comm,&Ce);
353: MatSetSizes(Ce,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
354: MatSetType(Ce,MATELEMENTAL);
355: MatSetUp(Ce);
356: *C = Ce;
357: return(0);
358: }
360: static PetscErrorCode MatMatMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
361: {
365: if (scall == MAT_INITIAL_MATRIX){
366: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
367: MatMatMultSymbolic_Elemental(A,B,1.0,C);
368: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
369: }
370: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
371: MatMatMultNumeric_Elemental(A,B,*C);
372: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
373: return(0);
374: }
376: static PetscErrorCode MatMatTransposeMultNumeric_Elemental(Mat A,Mat B,Mat C)
377: {
378: Mat_Elemental *a = (Mat_Elemental*)A->data;
379: Mat_Elemental *b = (Mat_Elemental*)B->data;
380: Mat_Elemental *c = (Mat_Elemental*)C->data;
381: PetscElemScalar one = 1,zero = 0;
384: { /* Scoping so that constructor is called before pointer is returned */
385: El::Gemm(El::NORMAL,El::TRANSPOSE,one,*a->emat,*b->emat,zero,*c->emat);
386: }
387: C->assembled = PETSC_TRUE;
388: return(0);
389: }
391: static PetscErrorCode MatMatTransposeMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C)
392: {
394: Mat Ce;
395: MPI_Comm comm;
398: PetscObjectGetComm((PetscObject)A,&comm);
399: MatCreate(comm,&Ce);
400: MatSetSizes(Ce,A->rmap->n,B->rmap->n,PETSC_DECIDE,PETSC_DECIDE);
401: MatSetType(Ce,MATELEMENTAL);
402: MatSetUp(Ce);
403: *C = Ce;
404: return(0);
405: }
407: static PetscErrorCode MatMatTransposeMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
408: {
412: if (scall == MAT_INITIAL_MATRIX){
413: PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);
414: MatMatMultSymbolic_Elemental(A,B,1.0,C);
415: PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);
416: }
417: PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);
418: MatMatTransposeMultNumeric_Elemental(A,B,*C);
419: PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);
420: return(0);
421: }
423: static PetscErrorCode MatGetDiagonal_Elemental(Mat A,Vec D)
424: {
425: PetscInt i,nrows,ncols,nD,rrank,ridx,crank,cidx;
426: Mat_Elemental *a = (Mat_Elemental*)A->data;
427: PetscErrorCode ierr;
428: PetscElemScalar v;
429: MPI_Comm comm;
432: PetscObjectGetComm((PetscObject)A,&comm);
433: MatGetSize(A,&nrows,&ncols);
434: nD = nrows>ncols ? ncols : nrows;
435: for (i=0; i<nD; i++) {
436: PetscInt erow,ecol;
437: P2RO(A,0,i,&rrank,&ridx);
438: RO2E(A,0,rrank,ridx,&erow);
439: if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
440: P2RO(A,1,i,&crank,&cidx);
441: RO2E(A,1,crank,cidx,&ecol);
442: if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
443: v = a->emat->Get(erow,ecol);
444: VecSetValues(D,1,&i,(PetscScalar*)&v,INSERT_VALUES);
445: }
446: VecAssemblyBegin(D);
447: VecAssemblyEnd(D);
448: return(0);
449: }
451: static PetscErrorCode MatDiagonalScale_Elemental(Mat X,Vec L,Vec R)
452: {
453: Mat_Elemental *x = (Mat_Elemental*)X->data;
454: const PetscElemScalar *d;
455: PetscErrorCode ierr;
458: if (R) {
459: VecGetArrayRead(R,(const PetscScalar **)&d);
460: El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
461: de.LockedAttach(X->cmap->N,1,*x->grid,0,0,d,X->cmap->n);
462: El::DiagonalScale(El::RIGHT,El::NORMAL,de,*x->emat);
463: VecRestoreArrayRead(R,(const PetscScalar **)&d);
464: }
465: if (L) {
466: VecGetArrayRead(L,(const PetscScalar **)&d);
467: El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
468: de.LockedAttach(X->rmap->N,1,*x->grid,0,0,d,X->rmap->n);
469: El::DiagonalScale(El::LEFT,El::NORMAL,de,*x->emat);
470: VecRestoreArrayRead(L,(const PetscScalar **)&d);
471: }
472: return(0);
473: }
475: static PetscErrorCode MatMissingDiagonal_Elemental(Mat A,PetscBool *missing,PetscInt *d)
476: {
478: *missing = PETSC_FALSE;
479: return(0);
480: }
482: static PetscErrorCode MatScale_Elemental(Mat X,PetscScalar a)
483: {
484: Mat_Elemental *x = (Mat_Elemental*)X->data;
487: El::Scale((PetscElemScalar)a,*x->emat);
488: return(0);
489: }
491: /*
492: MatAXPY - Computes Y = a*X + Y.
493: */
494: static PetscErrorCode MatAXPY_Elemental(Mat Y,PetscScalar a,Mat X,MatStructure str)
495: {
496: Mat_Elemental *x = (Mat_Elemental*)X->data;
497: Mat_Elemental *y = (Mat_Elemental*)Y->data;
501: El::Axpy((PetscElemScalar)a,*x->emat,*y->emat);
502: PetscObjectStateIncrease((PetscObject)Y);
503: return(0);
504: }
506: static PetscErrorCode MatCopy_Elemental(Mat A,Mat B,MatStructure str)
507: {
508: Mat_Elemental *a=(Mat_Elemental*)A->data;
509: Mat_Elemental *b=(Mat_Elemental*)B->data;
513: El::Copy(*a->emat,*b->emat);
514: PetscObjectStateIncrease((PetscObject)B);
515: return(0);
516: }
518: static PetscErrorCode MatDuplicate_Elemental(Mat A,MatDuplicateOption op,Mat *B)
519: {
520: Mat Be;
521: MPI_Comm comm;
522: Mat_Elemental *a=(Mat_Elemental*)A->data;
526: PetscObjectGetComm((PetscObject)A,&comm);
527: MatCreate(comm,&Be);
528: MatSetSizes(Be,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
529: MatSetType(Be,MATELEMENTAL);
530: MatSetUp(Be);
531: *B = Be;
532: if (op == MAT_COPY_VALUES) {
533: Mat_Elemental *b=(Mat_Elemental*)Be->data;
534: El::Copy(*a->emat,*b->emat);
535: }
536: Be->assembled = PETSC_TRUE;
537: return(0);
538: }
540: static PetscErrorCode MatTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
541: {
542: Mat Be = *B;
544: MPI_Comm comm;
545: Mat_Elemental *a = (Mat_Elemental*)A->data, *b;
548: PetscObjectGetComm((PetscObject)A,&comm);
549: /* Only out-of-place supported */
550: if (reuse == MAT_INITIAL_MATRIX){
551: MatCreate(comm,&Be);
552: MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);
553: MatSetType(Be,MATELEMENTAL);
554: MatSetUp(Be);
555: *B = Be;
556: }
557: b = (Mat_Elemental*)Be->data;
558: El::Transpose(*a->emat,*b->emat);
559: Be->assembled = PETSC_TRUE;
560: return(0);
561: }
563: static PetscErrorCode MatConjugate_Elemental(Mat A)
564: {
565: Mat_Elemental *a = (Mat_Elemental*)A->data;
568: El::Conjugate(*a->emat);
569: return(0);
570: }
572: static PetscErrorCode MatHermitianTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
573: {
574: Mat Be = *B;
576: MPI_Comm comm;
577: Mat_Elemental *a = (Mat_Elemental*)A->data, *b;
580: PetscObjectGetComm((PetscObject)A,&comm);
581: /* Only out-of-place supported */
582: if (reuse == MAT_INITIAL_MATRIX){
583: MatCreate(comm,&Be);
584: MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);
585: MatSetType(Be,MATELEMENTAL);
586: MatSetUp(Be);
587: *B = Be;
588: }
589: b = (Mat_Elemental*)Be->data;
590: El::Adjoint(*a->emat,*b->emat);
591: Be->assembled = PETSC_TRUE;
592: return(0);
593: }
595: static PetscErrorCode MatSolve_Elemental(Mat A,Vec B,Vec X)
596: {
597: Mat_Elemental *a = (Mat_Elemental*)A->data;
598: PetscErrorCode ierr;
599: PetscElemScalar *x;
600: PetscInt pivoting = a->pivoting;
603: VecCopy(B,X);
604: VecGetArray(X,(PetscScalar **)&x);
606: El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe;
607: xe.Attach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
608: El::DistMatrix<PetscElemScalar,El::MC,El::MR> xer(xe);
609: switch (A->factortype) {
610: case MAT_FACTOR_LU:
611: if (pivoting == 0) {
612: El::lu::SolveAfter(El::NORMAL,*a->emat,xer);
613: } else if (pivoting == 1) {
614: El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,xer);
615: } else { /* pivoting == 2 */
616: El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,xer);
617: }
618: break;
619: case MAT_FACTOR_CHOLESKY:
620: El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,xer);
621: break;
622: default:
623: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
624: break;
625: }
626: El::Copy(xer,xe);
628: VecRestoreArray(X,(PetscScalar **)&x);
629: return(0);
630: }
632: static PetscErrorCode MatSolveAdd_Elemental(Mat A,Vec B,Vec Y,Vec X)
633: {
634: PetscErrorCode ierr;
637: MatSolve_Elemental(A,B,X);
638: VecAXPY(X,1,Y);
639: return(0);
640: }
642: static PetscErrorCode MatMatSolve_Elemental(Mat A,Mat B,Mat X)
643: {
644: Mat_Elemental *a=(Mat_Elemental*)A->data;
645: Mat_Elemental *b=(Mat_Elemental*)B->data;
646: Mat_Elemental *x=(Mat_Elemental*)X->data;
647: PetscInt pivoting = a->pivoting;
650: El::Copy(*b->emat,*x->emat);
651: switch (A->factortype) {
652: case MAT_FACTOR_LU:
653: if (pivoting == 0) {
654: El::lu::SolveAfter(El::NORMAL,*a->emat,*x->emat);
655: } else if (pivoting == 1) {
656: El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*x->emat);
657: } else {
658: El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,*x->emat);
659: }
660: break;
661: case MAT_FACTOR_CHOLESKY:
662: El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,*x->emat);
663: break;
664: default:
665: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
666: break;
667: }
668: return(0);
669: }
671: static PetscErrorCode MatLUFactor_Elemental(Mat A,IS row,IS col,const MatFactorInfo *info)
672: {
673: Mat_Elemental *a = (Mat_Elemental*)A->data;
675: PetscInt pivoting = a->pivoting;
678: if (pivoting == 0) {
679: El::LU(*a->emat);
680: } else if (pivoting == 1) {
681: El::LU(*a->emat,*a->P);
682: } else {
683: El::LU(*a->emat,*a->P,*a->Q);
684: }
685: A->factortype = MAT_FACTOR_LU;
686: A->assembled = PETSC_TRUE;
688: PetscFree(A->solvertype);
689: PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype);
690: return(0);
691: }
693: static PetscErrorCode MatLUFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
694: {
698: MatCopy(A,F,SAME_NONZERO_PATTERN);
699: MatLUFactor_Elemental(F,0,0,info);
700: return(0);
701: }
703: static PetscErrorCode MatLUFactorSymbolic_Elemental(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
704: {
706: /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
707: return(0);
708: }
710: static PetscErrorCode MatCholeskyFactor_Elemental(Mat A,IS perm,const MatFactorInfo *info)
711: {
712: Mat_Elemental *a = (Mat_Elemental*)A->data;
713: El::DistMatrix<PetscElemScalar,El::MC,El::STAR> d;
717: El::Cholesky(El::UPPER,*a->emat);
718: A->factortype = MAT_FACTOR_CHOLESKY;
719: A->assembled = PETSC_TRUE;
721: PetscFree(A->solvertype);
722: PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype);
723: return(0);
724: }
726: static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
727: {
731: MatCopy(A,F,SAME_NONZERO_PATTERN);
732: MatCholeskyFactor_Elemental(F,0,info);
733: return(0);
734: }
736: static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat F,Mat A,IS perm,const MatFactorInfo *info)
737: {
739: /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
740: return(0);
741: }
743: PetscErrorCode MatFactorGetSolverPackage_elemental_elemental(Mat A,const MatSolverPackage *type)
744: {
746: *type = MATSOLVERELEMENTAL;
747: return(0);
748: }
750: static PetscErrorCode MatGetFactor_elemental_elemental(Mat A,MatFactorType ftype,Mat *F)
751: {
752: Mat B;
756: /* Create the factorization matrix */
757: MatCreate(PetscObjectComm((PetscObject)A),&B);
758: MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
759: MatSetType(B,MATELEMENTAL);
760: MatSetUp(B);
761: B->factortype = ftype;
762: PetscFree(B->solvertype);
763: PetscStrallocpy(MATSOLVERELEMENTAL,&B->solvertype);
765: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_elemental_elemental);
766: *F = B;
767: return(0);
768: }
770: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_Elemental(void)
771: {
775: MatSolverPackageRegister(MATSOLVERELEMENTAL,MATELEMENTAL, MAT_FACTOR_LU,MatGetFactor_elemental_elemental);
776: MatSolverPackageRegister(MATSOLVERELEMENTAL,MATELEMENTAL, MAT_FACTOR_CHOLESKY,MatGetFactor_elemental_elemental);
777: return(0);
778: }
780: static PetscErrorCode MatNorm_Elemental(Mat A,NormType type,PetscReal *nrm)
781: {
782: Mat_Elemental *a=(Mat_Elemental*)A->data;
785: switch (type){
786: case NORM_1:
787: *nrm = El::OneNorm(*a->emat);
788: break;
789: case NORM_FROBENIUS:
790: *nrm = El::FrobeniusNorm(*a->emat);
791: break;
792: case NORM_INFINITY:
793: *nrm = El::InfinityNorm(*a->emat);
794: break;
795: default:
796: printf("Error: unsupported norm type!\n");
797: }
798: return(0);
799: }
801: static PetscErrorCode MatZeroEntries_Elemental(Mat A)
802: {
803: Mat_Elemental *a=(Mat_Elemental*)A->data;
806: El::Zero(*a->emat);
807: return(0);
808: }
810: static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A,IS *rows,IS *cols)
811: {
812: Mat_Elemental *a = (Mat_Elemental*)A->data;
814: PetscInt i,m,shift,stride,*idx;
817: if (rows) {
818: m = a->emat->LocalHeight();
819: shift = a->emat->ColShift();
820: stride = a->emat->ColStride();
821: PetscMalloc1(m,&idx);
822: for (i=0; i<m; i++) {
823: PetscInt rank,offset;
824: E2RO(A,0,shift+i*stride,&rank,&offset);
825: RO2P(A,0,rank,offset,&idx[i]);
826: }
827: ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,rows);
828: }
829: if (cols) {
830: m = a->emat->LocalWidth();
831: shift = a->emat->RowShift();
832: stride = a->emat->RowStride();
833: PetscMalloc1(m,&idx);
834: for (i=0; i<m; i++) {
835: PetscInt rank,offset;
836: E2RO(A,1,shift+i*stride,&rank,&offset);
837: RO2P(A,1,rank,offset,&idx[i]);
838: }
839: ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,cols);
840: }
841: return(0);
842: }
844: static PetscErrorCode MatConvert_Elemental_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
845: {
846: Mat Bmpi;
847: Mat_Elemental *a = (Mat_Elemental*)A->data;
848: MPI_Comm comm;
849: PetscErrorCode ierr;
850: IS isrows,iscols;
851: PetscInt rrank,ridx,crank,cidx,nrows,ncols,i,j,erow,ecol,elrow,elcol;
852: const PetscInt *rows,*cols;
853: PetscElemScalar v;
854: const El::Grid &grid = a->emat->Grid();
857: PetscObjectGetComm((PetscObject)A,&comm);
858:
859: if (reuse == MAT_REUSE_MATRIX) {
860: Bmpi = *B;
861: } else {
862: MatCreate(comm,&Bmpi);
863: MatSetSizes(Bmpi,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
864: MatSetType(Bmpi,MATDENSE);
865: MatSetUp(Bmpi);
866: }
868: /* Get local entries of A */
869: MatGetOwnershipIS(A,&isrows,&iscols);
870: ISGetLocalSize(isrows,&nrows);
871: ISGetIndices(isrows,&rows);
872: ISGetLocalSize(iscols,&ncols);
873: ISGetIndices(iscols,&cols);
875: if (a->roworiented) {
876: for (i=0; i<nrows; i++) {
877: P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
878: RO2E(A,0,rrank,ridx,&erow);
879: if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
880: for (j=0; j<ncols; j++) {
881: P2RO(A,1,cols[j],&crank,&cidx);
882: RO2E(A,1,crank,cidx,&ecol);
883: if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
885: elrow = erow / grid.MCSize(); /* Elemental local row index */
886: elcol = ecol / grid.MRSize(); /* Elemental local column index */
887: v = a->emat->GetLocal(elrow,elcol);
888: MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES);
889: }
890: }
891: } else { /* column-oriented */
892: for (j=0; j<ncols; j++) {
893: P2RO(A,1,cols[j],&crank,&cidx);
894: RO2E(A,1,crank,cidx,&ecol);
895: if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
896: for (i=0; i<nrows; i++) {
897: P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
898: RO2E(A,0,rrank,ridx,&erow);
899: if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
901: elrow = erow / grid.MCSize(); /* Elemental local row index */
902: elcol = ecol / grid.MRSize(); /* Elemental local column index */
903: v = a->emat->GetLocal(elrow,elcol);
904: MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES);
905: }
906: }
907: }
908: MatAssemblyBegin(Bmpi,MAT_FINAL_ASSEMBLY);
909: MatAssemblyEnd(Bmpi,MAT_FINAL_ASSEMBLY);
910: if (reuse == MAT_INPLACE_MATRIX) {
911: MatHeaderReplace(A,&Bmpi);
912: } else {
913: *B = Bmpi;
914: }
915: ISDestroy(&isrows);
916: ISDestroy(&iscols);
917: return(0);
918: }
920: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
921: {
922: Mat mat_elemental;
923: PetscErrorCode ierr;
924: PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols;
925: const PetscInt *cols;
926: const PetscScalar *vals;
929: MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
930: MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
931: MatSetType(mat_elemental,MATELEMENTAL);
932: MatSetUp(mat_elemental);
933: for (row=0; row<M; row++) {
934: MatGetRow(A,row,&ncols,&cols,&vals);
935: /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
936: MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);
937: MatRestoreRow(A,row,&ncols,&cols,&vals);
938: }
939: MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
940: MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
942: if (reuse == MAT_INPLACE_MATRIX) {
943: MatHeaderReplace(A,&mat_elemental);
944: } else {
945: *newmat = mat_elemental;
946: }
947: return(0);
948: }
950: PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
951: {
952: Mat mat_elemental;
953: PetscErrorCode ierr;
954: PetscInt row,ncols,rstart=A->rmap->rstart,rend=A->rmap->rend,j;
955: const PetscInt *cols;
956: const PetscScalar *vals;
959: MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
960: MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);
961: MatSetType(mat_elemental,MATELEMENTAL);
962: MatSetUp(mat_elemental);
963: for (row=rstart; row<rend; row++) {
964: MatGetRow(A,row,&ncols,&cols,&vals);
965: for (j=0; j<ncols; j++) {
966: /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
967: MatSetValues(mat_elemental,1,&row,1,&cols[j],&vals[j],ADD_VALUES);
968: }
969: MatRestoreRow(A,row,&ncols,&cols,&vals);
970: }
971: MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
972: MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
974: if (reuse == MAT_INPLACE_MATRIX) {
975: MatHeaderReplace(A,&mat_elemental);
976: } else {
977: *newmat = mat_elemental;
978: }
979: return(0);
980: }
982: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
983: {
984: Mat mat_elemental;
985: PetscErrorCode ierr;
986: PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols,j;
987: const PetscInt *cols;
988: const PetscScalar *vals;
991: MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
992: MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
993: MatSetType(mat_elemental,MATELEMENTAL);
994: MatSetUp(mat_elemental);
995: MatGetRowUpperTriangular(A);
996: for (row=0; row<M; row++) {
997: MatGetRow(A,row,&ncols,&cols,&vals);
998: /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
999: MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);
1000: for (j=0; j<ncols; j++) { /* lower triangular part */
1001: if (cols[j] == row) continue;
1002: MatSetValues(mat_elemental,1,&cols[j],1,&row,&vals[j],ADD_VALUES);
1003: }
1004: MatRestoreRow(A,row,&ncols,&cols,&vals);
1005: }
1006: MatRestoreRowUpperTriangular(A);
1007: MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
1008: MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
1010: if (reuse == MAT_INPLACE_MATRIX) {
1011: MatHeaderReplace(A,&mat_elemental);
1012: } else {
1013: *newmat = mat_elemental;
1014: }
1015: return(0);
1016: }
1018: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
1019: {
1020: Mat mat_elemental;
1021: PetscErrorCode ierr;
1022: PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols,j,rstart=A->rmap->rstart,rend=A->rmap->rend;
1023: const PetscInt *cols;
1024: const PetscScalar *vals;
1027: MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
1028: MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);
1029: MatSetType(mat_elemental,MATELEMENTAL);
1030: MatSetUp(mat_elemental);
1031: MatGetRowUpperTriangular(A);
1032: for (row=rstart; row<rend; row++) {
1033: MatGetRow(A,row,&ncols,&cols,&vals);
1034: /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
1035: MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);
1036: for (j=0; j<ncols; j++) { /* lower triangular part */
1037: if (cols[j] == row) continue;
1038: MatSetValues(mat_elemental,1,&cols[j],1,&row,&vals[j],ADD_VALUES);
1039: }
1040: MatRestoreRow(A,row,&ncols,&cols,&vals);
1041: }
1042: MatRestoreRowUpperTriangular(A);
1043: MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
1044: MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
1046: if (reuse == MAT_INPLACE_MATRIX) {
1047: MatHeaderReplace(A,&mat_elemental);
1048: } else {
1049: *newmat = mat_elemental;
1050: }
1051: return(0);
1052: }
1054: static PetscErrorCode MatDestroy_Elemental(Mat A)
1055: {
1056: Mat_Elemental *a = (Mat_Elemental*)A->data;
1057: PetscErrorCode ierr;
1058: Mat_Elemental_Grid *commgrid;
1059: PetscBool flg;
1060: MPI_Comm icomm;
1063: delete a->emat;
1064: delete a->P;
1065: delete a->Q;
1067: El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
1068: PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);
1069: MPI_Attr_get(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);
1070: if (--commgrid->grid_refct == 0) {
1071: delete commgrid->grid;
1072: PetscFree(commgrid);
1073: MPI_Keyval_free(&Petsc_Elemental_keyval);
1074: }
1075: PetscCommDestroy(&icomm);
1076: PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",NULL);
1077: PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);
1078: PetscFree(A->data);
1079: return(0);
1080: }
1082: PetscErrorCode MatSetUp_Elemental(Mat A)
1083: {
1084: Mat_Elemental *a = (Mat_Elemental*)A->data;
1086: PetscMPIInt rsize,csize;
1089: PetscLayoutSetUp(A->rmap);
1090: PetscLayoutSetUp(A->cmap);
1092: a->emat->Resize(A->rmap->N,A->cmap->N);
1093: El::Zero(*a->emat);
1095: MPI_Comm_size(A->rmap->comm,&rsize);
1096: MPI_Comm_size(A->cmap->comm,&csize);
1097: if (csize != rsize) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot use row and column communicators of different sizes");
1098: a->commsize = rsize;
1099: a->mr[0] = A->rmap->N % rsize; if (!a->mr[0]) a->mr[0] = rsize;
1100: a->mr[1] = A->cmap->N % csize; if (!a->mr[1]) a->mr[1] = csize;
1101: a->m[0] = A->rmap->N / rsize + (a->mr[0] != rsize);
1102: a->m[1] = A->cmap->N / csize + (a->mr[1] != csize);
1103: return(0);
1104: }
1106: PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType type)
1107: {
1108: Mat_Elemental *a = (Mat_Elemental*)A->data;
1111: /* printf("Calling ProcessQueues\n"); */
1112: a->emat->ProcessQueues();
1113: /* printf("Finished ProcessQueues\n"); */
1114: return(0);
1115: }
1117: PetscErrorCode MatAssemblyEnd_Elemental(Mat A, MatAssemblyType type)
1118: {
1120: /* Currently does nothing */
1121: return(0);
1122: }
1124: PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer)
1125: {
1127: Mat Adense,Ae;
1128: MPI_Comm comm;
1131: PetscObjectGetComm((PetscObject)newMat,&comm);
1132: MatCreate(comm,&Adense);
1133: MatSetType(Adense,MATDENSE);
1134: MatLoad(Adense,viewer);
1135: MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX,&Ae);
1136: MatDestroy(&Adense);
1137: MatHeaderReplace(newMat,&Ae);
1138: return(0);
1139: }
1141: /* -------------------------------------------------------------------*/
1142: static struct _MatOps MatOps_Values = {
1143: MatSetValues_Elemental,
1144: 0,
1145: 0,
1146: MatMult_Elemental,
1147: /* 4*/ MatMultAdd_Elemental,
1148: MatMultTranspose_Elemental,
1149: MatMultTransposeAdd_Elemental,
1150: MatSolve_Elemental,
1151: MatSolveAdd_Elemental,
1152: 0,
1153: /*10*/ 0,
1154: MatLUFactor_Elemental,
1155: MatCholeskyFactor_Elemental,
1156: 0,
1157: MatTranspose_Elemental,
1158: /*15*/ MatGetInfo_Elemental,
1159: 0,
1160: MatGetDiagonal_Elemental,
1161: MatDiagonalScale_Elemental,
1162: MatNorm_Elemental,
1163: /*20*/ MatAssemblyBegin_Elemental,
1164: MatAssemblyEnd_Elemental,
1165: MatSetOption_Elemental,
1166: MatZeroEntries_Elemental,
1167: /*24*/ 0,
1168: MatLUFactorSymbolic_Elemental,
1169: MatLUFactorNumeric_Elemental,
1170: MatCholeskyFactorSymbolic_Elemental,
1171: MatCholeskyFactorNumeric_Elemental,
1172: /*29*/ MatSetUp_Elemental,
1173: 0,
1174: 0,
1175: 0,
1176: 0,
1177: /*34*/ MatDuplicate_Elemental,
1178: 0,
1179: 0,
1180: 0,
1181: 0,
1182: /*39*/ MatAXPY_Elemental,
1183: 0,
1184: 0,
1185: 0,
1186: MatCopy_Elemental,
1187: /*44*/ 0,
1188: MatScale_Elemental,
1189: MatShift_Basic,
1190: 0,
1191: 0,
1192: /*49*/ 0,
1193: 0,
1194: 0,
1195: 0,
1196: 0,
1197: /*54*/ 0,
1198: 0,
1199: 0,
1200: 0,
1201: 0,
1202: /*59*/ 0,
1203: MatDestroy_Elemental,
1204: MatView_Elemental,
1205: 0,
1206: 0,
1207: /*64*/ 0,
1208: 0,
1209: 0,
1210: 0,
1211: 0,
1212: /*69*/ 0,
1213: 0,
1214: MatConvert_Elemental_Dense,
1215: 0,
1216: 0,
1217: /*74*/ 0,
1218: 0,
1219: 0,
1220: 0,
1221: 0,
1222: /*79*/ 0,
1223: 0,
1224: 0,
1225: 0,
1226: MatLoad_Elemental,
1227: /*84*/ 0,
1228: 0,
1229: 0,
1230: 0,
1231: 0,
1232: /*89*/ MatMatMult_Elemental,
1233: MatMatMultSymbolic_Elemental,
1234: MatMatMultNumeric_Elemental,
1235: 0,
1236: 0,
1237: /*94*/ 0,
1238: MatMatTransposeMult_Elemental,
1239: MatMatTransposeMultSymbolic_Elemental,
1240: MatMatTransposeMultNumeric_Elemental,
1241: 0,
1242: /*99*/ 0,
1243: 0,
1244: 0,
1245: MatConjugate_Elemental,
1246: 0,
1247: /*104*/0,
1248: 0,
1249: 0,
1250: 0,
1251: 0,
1252: /*109*/MatMatSolve_Elemental,
1253: 0,
1254: 0,
1255: 0,
1256: MatMissingDiagonal_Elemental,
1257: /*114*/0,
1258: 0,
1259: 0,
1260: 0,
1261: 0,
1262: /*119*/0,
1263: MatHermitianTranspose_Elemental,
1264: 0,
1265: 0,
1266: 0,
1267: /*124*/0,
1268: 0,
1269: 0,
1270: 0,
1271: 0,
1272: /*129*/0,
1273: 0,
1274: 0,
1275: 0,
1276: 0,
1277: /*134*/0,
1278: 0,
1279: 0,
1280: 0,
1281: 0
1282: };
1284: /*MC
1285: MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package
1287: Use ./configure --download-elemental to install PETSc to use Elemental
1289: Use -pc_type lu -pc_factor_mat_solver_package elemental to use this direct solver
1291: Options Database Keys:
1292: + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions()
1293: - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix
1295: Level: beginner
1297: .seealso: MATDENSE
1298: M*/
1300: PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A)
1301: {
1302: Mat_Elemental *a;
1303: PetscErrorCode ierr;
1304: PetscBool flg,flg1;
1305: Mat_Elemental_Grid *commgrid;
1306: MPI_Comm icomm;
1307: PetscInt optv1;
1310: PetscElementalInitializePackage();
1311: PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));
1312: A->insertmode = NOT_SET_VALUES;
1314: PetscNewLog(A,&a);
1315: A->data = (void*)a;
1317: /* Set up the elemental matrix */
1318: El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
1320: /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */
1321: if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) {
1322: MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0);
1323: /* MPI_Comm_create_Keyval(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0); -- new version? */
1324: }
1325: PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);
1326: MPI_Attr_get(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);
1327: if (!flg) {
1328: PetscNewLog(A,&commgrid);
1330: PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Elemental Options","Mat");
1331: /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */
1332: PetscOptionsInt("-mat_elemental_grid_height","Grid Height","None",El::mpi::Size(cxxcomm),&optv1,&flg1);
1333: if (flg1) {
1334: if (El::mpi::Size(cxxcomm) % optv1 != 0) {
1335: SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Grid Height %D must evenly divide CommSize %D",optv1,(PetscInt)El::mpi::Size(cxxcomm));
1336: }
1337: commgrid->grid = new El::Grid(cxxcomm,optv1); /* use user-provided grid height */
1338: } else {
1339: commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */
1340: /* printf("new commgrid->grid = %p\n",commgrid->grid); -- memory leak revealed by valgrind? */
1341: }
1342: commgrid->grid_refct = 1;
1343: MPI_Attr_put(icomm,Petsc_Elemental_keyval,(void*)commgrid);
1345: a->pivoting = 1;
1346: PetscOptionsInt("-mat_elemental_pivoting","Pivoting","None",a->pivoting,&a->pivoting,NULL);
1348: PetscOptionsEnd();
1349: } else {
1350: commgrid->grid_refct++;
1351: }
1352: PetscCommDestroy(&icomm);
1353: a->grid = commgrid->grid;
1354: a->emat = new El::DistMatrix<PetscElemScalar>(*a->grid);
1355: a->roworiented = PETSC_TRUE;
1357: PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",MatGetOwnershipIS_Elemental);
1358: PetscObjectChangeTypeName((PetscObject)A,MATELEMENTAL);
1359: return(0);
1360: }