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