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