xref: /petsc/src/mat/impls/elemental/matelem.cxx (revision 320f2790d574cc2255cb4d463df0df828ded99b1)
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   { /* We have already initialized MPI, so this song and dance is just to pass these variables (which won't be used by Elemental) through the interface that needs references */
27     int zero = 0;
28     char **nothing = 0;
29     El::Initialize(zero,nothing);   /* called by the 1st call of MatCreate_Elemental */
30   }
31   ierr = PetscRegisterFinalize(PetscElementalFinalizePackage);CHKERRQ(ierr);
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNCT__
36 #define __FUNCT__ "PetscElementalFinalizePackage"
37 /*@C
38    PetscElementalFinalizePackage - Finalize Elemental package
39 
40    Logically Collective
41 
42    Level: developer
43 
44 .seealso: MATELEMENTAL, PetscElementalInitializePackage()
45 @*/
46 PetscErrorCode PetscElementalFinalizePackage(void)
47 {
48   PetscFunctionBegin;
49   El::Finalize();  /* called by PetscFinalize() */
50   PetscFunctionReturn(0);
51 }
52 
53 #undef __FUNCT__
54 #define __FUNCT__ "MatView_Elemental"
55 static PetscErrorCode MatView_Elemental(Mat A,PetscViewer viewer)
56 {
57   PetscErrorCode ierr;
58   Mat_Elemental  *a = (Mat_Elemental*)A->data;
59   PetscBool      iascii;
60 
61   PetscFunctionBegin;
62   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
63   if (iascii) {
64     PetscViewerFormat format;
65     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
66     if (format == PETSC_VIEWER_ASCII_INFO) {
67       /* call elemental viewing function */
68       ierr = PetscViewerASCIIPrintf(viewer,"Elemental run parameters:\n");CHKERRQ(ierr);
69       ierr = PetscViewerASCIIPrintf(viewer,"  allocated entries=%d\n",(*a->emat).AllocatedMemory());CHKERRQ(ierr);
70       ierr = PetscViewerASCIIPrintf(viewer,"  grid height=%d, grid width=%d\n",(*a->emat).Grid().Height(),(*a->emat).Grid().Width());CHKERRQ(ierr);
71       if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
72         /* call elemental viewing function */
73         ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"test matview_elemental 2\n");CHKERRQ(ierr);
74       }
75 
76     } else if (format == PETSC_VIEWER_DEFAULT) {
77       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
78       El::Print( *a->emat, "Elemental matrix (cyclic ordering)" );
79       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
80       if (A->factortype == MAT_FACTOR_NONE){
81         Mat Adense;
82         ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");CHKERRQ(ierr);
83         ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr);
84         ierr = MatView(Adense,viewer);CHKERRQ(ierr);
85         ierr = MatDestroy(&Adense);CHKERRQ(ierr);
86       }
87     } else SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Format");
88   } else {
89     /* convert to dense format and call MatView() */
90     Mat Adense;
91     ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");CHKERRQ(ierr);
92     ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr);
93     ierr = MatView(Adense,viewer);CHKERRQ(ierr);
94     ierr = MatDestroy(&Adense);CHKERRQ(ierr);
95   }
96   PetscFunctionReturn(0);
97 }
98 
99 #undef __FUNCT__
100 #define __FUNCT__ "MatGetInfo_Elemental"
101 static PetscErrorCode MatGetInfo_Elemental(Mat A,MatInfoType flag,MatInfo *info)
102 {
103   Mat_Elemental  *a = (Mat_Elemental*)A->data;
104 
105   PetscFunctionBegin;
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 = MPI_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 = MPI_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   Mat_Elemental  *a = (Mat_Elemental*)A->data;
138   PetscInt       i,j,rrank,ridx,crank,cidx,erow,ecol,elrow,elcol;
139   const El::Grid &grid = a->emat->Grid();
140 
141   PetscFunctionBegin;
142   for (i=0; i<nr; i++) {
143     if (rows[i] < 0) continue;
144     P2RO(A,0,rows[i],&rrank,&ridx);
145     RO2E(A,0,rrank,ridx,&erow);
146     if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation");
147     for (j=0; j<nc; j++) {
148       if (cols[j] < 0) continue;
149       P2RO(A,1,cols[j],&crank,&cidx);
150       RO2E(A,1,crank,cidx,&ecol);
151       if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation");
152       if (erow % grid.MCSize() != grid.MCRank() || ecol % grid.MRSize() != grid.MRRank()){ /* off-proc entry */
153         if (imode != ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported");
154         /* 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); */
155         a->esubmat->Set(0,0, (PetscElemScalar)vals[i*nc+j]);
156         a->interface->Axpy(1.0,*(a->esubmat),erow,ecol);
157         continue;
158       }
159       elrow = erow / grid.MCSize();
160       elcol = ecol / grid.MRSize();
161       switch (imode) {
162       case INSERT_VALUES: a->emat->SetLocal(elrow,elcol,(PetscElemScalar)vals[i*nc+j]); break;
163       case ADD_VALUES: a->emat->UpdateLocal(elrow,elcol,(PetscElemScalar)vals[i*nc+j]); break;
164       default: SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode);
165       }
166     }
167   }
168   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "MatMult_Elemental"
173 static PetscErrorCode MatMult_Elemental(Mat A,Vec X,Vec Y)
174 {
175   Mat_Elemental         *a = (Mat_Elemental*)A->data;
176   PetscErrorCode        ierr;
177   const PetscElemScalar *x;
178   PetscElemScalar       *y;
179   PetscElemScalar       one = 1,zero = 0;
180 
181   PetscFunctionBegin;
182   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
183   ierr = VecGetArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
184   { /* Scoping so that constructor is called before pointer is returned */
185     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
186     xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
187     ye.Attach(A->rmap->N,1,*a->grid,0,0,y,A->rmap->n);
188     El::Gemv(El::NORMAL,one,*a->emat,xe,zero,ye);
189   }
190   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
191   ierr = VecRestoreArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
192   PetscFunctionReturn(0);
193 }
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "MatMultTranspose_Elemental"
197 static PetscErrorCode MatMultTranspose_Elemental(Mat A,Vec X,Vec Y)
198 {
199   Mat_Elemental         *a = (Mat_Elemental*)A->data;
200   PetscErrorCode        ierr;
201   const PetscElemScalar *x;
202   PetscElemScalar       *y;
203   PetscElemScalar       one = 1,zero = 0;
204 
205   PetscFunctionBegin;
206   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
207   ierr = VecGetArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
208   { /* Scoping so that constructor is called before pointer is returned */
209     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
210     xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
211     ye.Attach(A->cmap->N,1,*a->grid,0,0,y,A->cmap->n);
212     El::Gemv(El::TRANSPOSE,one,*a->emat,xe,zero,ye);
213   }
214   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
215   ierr = VecRestoreArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "MatMultAdd_Elemental"
221 static PetscErrorCode MatMultAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
222 {
223   Mat_Elemental         *a = (Mat_Elemental*)A->data;
224   PetscErrorCode        ierr;
225   const PetscElemScalar *x;
226   PetscElemScalar       *z;
227   PetscElemScalar       one = 1;
228 
229   PetscFunctionBegin;
230   if (Y != Z) {ierr = VecCopy(Y,Z);CHKERRQ(ierr);}
231   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
232   ierr = VecGetArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
233   { /* Scoping so that constructor is called before pointer is returned */
234     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
235     xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
236     ze.Attach(A->rmap->N,1,*a->grid,0,0,z,A->rmap->n);
237     El::Gemv(El::NORMAL,one,*a->emat,xe,one,ze);
238   }
239   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
240   ierr = VecRestoreArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "MatMultTransposeAdd_Elemental"
246 static PetscErrorCode MatMultTransposeAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
247 {
248   Mat_Elemental         *a = (Mat_Elemental*)A->data;
249   PetscErrorCode        ierr;
250   const PetscElemScalar *x;
251   PetscElemScalar       *z;
252   PetscElemScalar       one = 1;
253 
254   PetscFunctionBegin;
255   if (Y != Z) {ierr = VecCopy(Y,Z);CHKERRQ(ierr);}
256   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
257   ierr = VecGetArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
258   { /* Scoping so that constructor is called before pointer is returned */
259     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
260     xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
261     ze.Attach(A->cmap->N,1,*a->grid,0,0,z,A->cmap->n);
262     El::Gemv(El::TRANSPOSE,one,*a->emat,xe,one,ze);
263   }
264   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
265   ierr = VecRestoreArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "MatMatMultNumeric_Elemental"
271 static PetscErrorCode MatMatMultNumeric_Elemental(Mat A,Mat B,Mat C)
272 {
273   Mat_Elemental    *a = (Mat_Elemental*)A->data;
274   Mat_Elemental    *b = (Mat_Elemental*)B->data;
275   Mat_Elemental    *c = (Mat_Elemental*)C->data;
276   PetscElemScalar  one = 1,zero = 0;
277 
278   PetscFunctionBegin;
279   { /* Scoping so that constructor is called before pointer is returned */
280     El::Gemm(El::NORMAL,El::NORMAL,one,*a->emat,*b->emat,zero,*c->emat);
281   }
282   C->assembled = PETSC_TRUE;
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "MatMatMultSymbolic_Elemental"
288 static PetscErrorCode MatMatMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C)
289 {
290   PetscErrorCode ierr;
291   Mat            Ce;
292   MPI_Comm       comm;
293 
294   PetscFunctionBegin;
295   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
296   ierr = MatCreate(comm,&Ce);CHKERRQ(ierr);
297   ierr = MatSetSizes(Ce,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
298   ierr = MatSetType(Ce,MATELEMENTAL);CHKERRQ(ierr);
299   ierr = MatSetUp(Ce);CHKERRQ(ierr);
300   *C = Ce;
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "MatMatMult_Elemental"
306 static PetscErrorCode MatMatMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
307 {
308   PetscErrorCode ierr;
309 
310   PetscFunctionBegin;
311   if (scall == MAT_INITIAL_MATRIX){
312     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
313     ierr = MatMatMultSymbolic_Elemental(A,B,1.0,C);CHKERRQ(ierr);
314     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
315   }
316   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
317   ierr = MatMatMultNumeric_Elemental(A,B,*C);CHKERRQ(ierr);
318   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
319   PetscFunctionReturn(0);
320 }
321 
322 #undef __FUNCT__
323 #define __FUNCT__ "MatMatTransposeMultNumeric_Elemental"
324 static PetscErrorCode MatMatTransposeMultNumeric_Elemental(Mat A,Mat B,Mat C)
325 {
326   Mat_Elemental      *a = (Mat_Elemental*)A->data;
327   Mat_Elemental      *b = (Mat_Elemental*)B->data;
328   Mat_Elemental      *c = (Mat_Elemental*)C->data;
329   PetscElemScalar    one = 1,zero = 0;
330 
331   PetscFunctionBegin;
332   { /* Scoping so that constructor is called before pointer is returned */
333     El::Gemm(El::NORMAL,El::TRANSPOSE,one,*a->emat,*b->emat,zero,*c->emat);
334   }
335   C->assembled = PETSC_TRUE;
336   PetscFunctionReturn(0);
337 }
338 
339 #undef __FUNCT__
340 #define __FUNCT__ "MatMatTransposeMultSymbolic_Elemental"
341 static PetscErrorCode MatMatTransposeMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C)
342 {
343   PetscErrorCode ierr;
344   Mat            Ce;
345   MPI_Comm       comm;
346 
347   PetscFunctionBegin;
348   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
349   ierr = MatCreate(comm,&Ce);CHKERRQ(ierr);
350   ierr = MatSetSizes(Ce,A->rmap->n,B->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
351   ierr = MatSetType(Ce,MATELEMENTAL);CHKERRQ(ierr);
352   ierr = MatSetUp(Ce);CHKERRQ(ierr);
353   *C = Ce;
354   PetscFunctionReturn(0);
355 }
356 
357 #undef __FUNCT__
358 #define __FUNCT__ "MatMatTransposeMult_Elemental"
359 static PetscErrorCode MatMatTransposeMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
360 {
361   PetscErrorCode ierr;
362 
363   PetscFunctionBegin;
364   if (scall == MAT_INITIAL_MATRIX){
365     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
366     ierr = MatMatMultSymbolic_Elemental(A,B,1.0,C);CHKERRQ(ierr);
367     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
368   }
369   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
370   ierr = MatMatTransposeMultNumeric_Elemental(A,B,*C);CHKERRQ(ierr);
371   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "MatGetDiagonal_Elemental"
377 static PetscErrorCode MatGetDiagonal_Elemental(Mat A,Vec D)
378 {
379   PetscInt        i,nrows,ncols,nD,rrank,ridx,crank,cidx;
380   Mat_Elemental   *a = (Mat_Elemental*)A->data;
381   PetscErrorCode  ierr;
382   PetscElemScalar v;
383   MPI_Comm        comm;
384 
385   PetscFunctionBegin;
386   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
387   ierr = MatGetSize(A,&nrows,&ncols);CHKERRQ(ierr);
388   nD = nrows>ncols ? ncols : nrows;
389   for (i=0; i<nD; i++) {
390     PetscInt erow,ecol;
391     P2RO(A,0,i,&rrank,&ridx);
392     RO2E(A,0,rrank,ridx,&erow);
393     if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
394     P2RO(A,1,i,&crank,&cidx);
395     RO2E(A,1,crank,cidx,&ecol);
396     if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
397     v = a->emat->Get(erow,ecol);
398     ierr = VecSetValues(D,1,&i,(PetscScalar*)&v,INSERT_VALUES);CHKERRQ(ierr);
399   }
400   ierr = VecAssemblyBegin(D);CHKERRQ(ierr);
401   ierr = VecAssemblyEnd(D);CHKERRQ(ierr);
402   PetscFunctionReturn(0);
403 }
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "MatDiagonalScale_Elemental"
407 static PetscErrorCode MatDiagonalScale_Elemental(Mat X,Vec L,Vec R)
408 {
409   Mat_Elemental         *x = (Mat_Elemental*)X->data;
410   const PetscElemScalar *d;
411   PetscErrorCode        ierr;
412 
413   PetscFunctionBegin;
414   if (R) {
415     ierr = VecGetArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr);
416     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
417     de.LockedAttach(X->cmap->N,1,*x->grid,0,0,d,X->cmap->n);
418     El::DiagonalScale(El::RIGHT,El::NORMAL,de,*x->emat);
419     ierr = VecRestoreArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr);
420   }
421   if (L) {
422     ierr = VecGetArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr);
423     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
424     de.LockedAttach(X->rmap->N,1,*x->grid,0,0,d,X->rmap->n);
425     El::DiagonalScale(El::LEFT,El::NORMAL,de,*x->emat);
426     ierr = VecRestoreArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr);
427   }
428   PetscFunctionReturn(0);
429 }
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "MatScale_Elemental"
433 static PetscErrorCode MatScale_Elemental(Mat X,PetscScalar a)
434 {
435   Mat_Elemental  *x = (Mat_Elemental*)X->data;
436 
437   PetscFunctionBegin;
438   El::Scale((PetscElemScalar)a,*x->emat);
439   PetscFunctionReturn(0);
440 }
441 
442 /*
443   MatAXPY - Computes Y = a*X + Y.
444 */
445 #undef __FUNCT__
446 #define __FUNCT__ "MatAXPY_Elemental"
447 static PetscErrorCode MatAXPY_Elemental(Mat Y,PetscScalar a,Mat X,MatStructure str)
448 {
449   Mat_Elemental  *x = (Mat_Elemental*)X->data;
450   Mat_Elemental  *y = (Mat_Elemental*)Y->data;
451   PetscErrorCode ierr;
452 
453   PetscFunctionBegin;
454   El::Axpy((PetscElemScalar)a,*x->emat,*y->emat);
455   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
456   PetscFunctionReturn(0);
457 }
458 
459 #undef __FUNCT__
460 #define __FUNCT__ "MatCopy_Elemental"
461 static PetscErrorCode MatCopy_Elemental(Mat A,Mat B,MatStructure str)
462 {
463   Mat_Elemental *a=(Mat_Elemental*)A->data;
464   Mat_Elemental *b=(Mat_Elemental*)B->data;
465 
466   PetscFunctionBegin;
467   El::Copy(*a->emat,*b->emat);
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "MatDuplicate_Elemental"
473 static PetscErrorCode MatDuplicate_Elemental(Mat A,MatDuplicateOption op,Mat *B)
474 {
475   Mat            Be;
476   MPI_Comm       comm;
477   Mat_Elemental  *a=(Mat_Elemental*)A->data;
478   PetscErrorCode ierr;
479 
480   PetscFunctionBegin;
481   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
482   ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
483   ierr = MatSetSizes(Be,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
484   ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
485   ierr = MatSetUp(Be);CHKERRQ(ierr);
486   *B = Be;
487   if (op == MAT_COPY_VALUES) {
488     Mat_Elemental *b=(Mat_Elemental*)Be->data;
489     El::Copy(*a->emat,*b->emat);
490   }
491   Be->assembled = PETSC_TRUE;
492   PetscFunctionReturn(0);
493 }
494 
495 #undef __FUNCT__
496 #define __FUNCT__ "MatTranspose_Elemental"
497 static PetscErrorCode MatTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
498 {
499   Mat            Be = *B;
500   PetscErrorCode ierr;
501   MPI_Comm       comm;
502   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
503 
504   PetscFunctionBegin;
505   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
506   /* Only out-of-place supported */
507   if (reuse == MAT_INITIAL_MATRIX){
508     ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
509     ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
510     ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
511     ierr = MatSetUp(Be);CHKERRQ(ierr);
512     *B = Be;
513   }
514   b = (Mat_Elemental*)Be->data;
515   El::Transpose(*a->emat,*b->emat);
516   Be->assembled = PETSC_TRUE;
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "MatConjugate_Elemental"
522 static PetscErrorCode MatConjugate_Elemental(Mat A)
523 {
524   Mat_Elemental  *a = (Mat_Elemental*)A->data;
525 
526   PetscFunctionBegin;
527   El::Conjugate(*a->emat);
528   PetscFunctionReturn(0);
529 }
530 
531 #undef __FUNCT__
532 #define __FUNCT__ "MatHermitianTranspose_Elemental"
533 static PetscErrorCode MatHermitianTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
534 {
535   Mat            Be = *B;
536   PetscErrorCode ierr;
537   MPI_Comm       comm;
538   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
539 
540   PetscFunctionBegin;
541   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
542   /* Only out-of-place supported */
543   if (reuse == MAT_INITIAL_MATRIX){
544     ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
545     ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
546     ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
547     ierr = MatSetUp(Be);CHKERRQ(ierr);
548     *B = Be;
549   }
550   b = (Mat_Elemental*)Be->data;
551   El::Adjoint(*a->emat,*b->emat);
552   Be->assembled = PETSC_TRUE;
553   PetscFunctionReturn(0);
554 }
555 
556 #undef __FUNCT__
557 #define __FUNCT__ "MatSolve_Elemental"
558 static PetscErrorCode MatSolve_Elemental(Mat A,Vec B,Vec X)
559 {
560   Mat_Elemental     *a = (Mat_Elemental*)A->data;
561   PetscErrorCode    ierr;
562   PetscElemScalar   *x;
563 
564   PetscFunctionBegin;
565   ierr = VecCopy(B,X);CHKERRQ(ierr);
566   ierr = VecGetArray(X,(PetscScalar **)&x);CHKERRQ(ierr);
567   El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe;
568   xe.Attach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
569   El::DistMatrix<PetscElemScalar,El::MC,El::MR> xer(xe);
570   switch (A->factortype) {
571   case MAT_FACTOR_LU:
572     if ((*a->pivot).AllocatedMemory()) {
573       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->pivot,xer);
574       El::Copy(xer,xe);
575     } else {
576       El::lu::SolveAfter(El::NORMAL,*a->emat,xer);
577       El::Copy(xer,xe);
578     }
579     break;
580   case MAT_FACTOR_CHOLESKY:
581     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,xer);
582     El::Copy(xer,xe);
583     break;
584   default:
585     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
586     break;
587   }
588   ierr = VecRestoreArray(X,(PetscScalar **)&x);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591 
592 #undef __FUNCT__
593 #define __FUNCT__ "MatSolveAdd_Elemental"
594 static PetscErrorCode MatSolveAdd_Elemental(Mat A,Vec B,Vec Y,Vec X)
595 {
596   PetscErrorCode    ierr;
597 
598   PetscFunctionBegin;
599   ierr = MatSolve_Elemental(A,B,X);CHKERRQ(ierr);
600   ierr = VecAXPY(X,1,Y);CHKERRQ(ierr);
601   PetscFunctionReturn(0);
602 }
603 
604 #undef __FUNCT__
605 #define __FUNCT__ "MatMatSolve_Elemental"
606 static PetscErrorCode MatMatSolve_Elemental(Mat A,Mat B,Mat X)
607 {
608   Mat_Elemental *a=(Mat_Elemental*)A->data;
609   Mat_Elemental *b=(Mat_Elemental*)B->data;
610   Mat_Elemental *x=(Mat_Elemental*)X->data;
611 
612   PetscFunctionBegin;
613   El::Copy(*b->emat,*x->emat);
614   switch (A->factortype) {
615   case MAT_FACTOR_LU:
616     if ((*a->pivot).AllocatedMemory()) {
617       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->pivot,*x->emat);
618     } else {
619       El::lu::SolveAfter(El::NORMAL,*a->emat,*x->emat);
620     }
621     break;
622   case MAT_FACTOR_CHOLESKY:
623     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,*x->emat);
624     break;
625   default:
626     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
627     break;
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "MatLUFactor_Elemental"
634 static PetscErrorCode MatLUFactor_Elemental(Mat A,IS row,IS col,const MatFactorInfo *info)
635 {
636   Mat_Elemental  *a = (Mat_Elemental*)A->data;
637 
638   PetscFunctionBegin;
639   if (info->dtcol){
640     El::LU(*a->emat,*a->pivot);
641   } else {
642     El::LU(*a->emat);
643   }
644   A->factortype = MAT_FACTOR_LU;
645   A->assembled  = PETSC_TRUE;
646   PetscFunctionReturn(0);
647 }
648 
649 #undef __FUNCT__
650 #define __FUNCT__ "MatLUFactorNumeric_Elemental"
651 static PetscErrorCode  MatLUFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
652 {
653   PetscErrorCode ierr;
654 
655   PetscFunctionBegin;
656   ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
657   ierr = MatLUFactor_Elemental(F,0,0,info);CHKERRQ(ierr);
658   PetscFunctionReturn(0);
659 }
660 
661 #undef __FUNCT__
662 #define __FUNCT__ "MatLUFactorSymbolic_Elemental"
663 static PetscErrorCode  MatLUFactorSymbolic_Elemental(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
664 {
665   PetscFunctionBegin;
666   /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
667   PetscFunctionReturn(0);
668 }
669 
670 #undef __FUNCT__
671 #define __FUNCT__ "MatCholeskyFactor_Elemental"
672 static PetscErrorCode MatCholeskyFactor_Elemental(Mat A,IS perm,const MatFactorInfo *info)
673 {
674   Mat_Elemental  *a = (Mat_Elemental*)A->data;
675   El::DistMatrix<PetscElemScalar,El::MC,El::STAR> d;
676 
677   PetscFunctionBegin;
678   El::Cholesky(El::UPPER,*a->emat);
679   A->factortype = MAT_FACTOR_CHOLESKY;
680   A->assembled  = PETSC_TRUE;
681   PetscFunctionReturn(0);
682 }
683 
684 #undef __FUNCT__
685 #define __FUNCT__ "MatCholeskyFactorNumeric_Elemental"
686 static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
687 {
688   PetscErrorCode ierr;
689 
690   PetscFunctionBegin;
691   ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
692   ierr = MatCholeskyFactor_Elemental(F,0,info);CHKERRQ(ierr);
693   PetscFunctionReturn(0);
694 }
695 
696 #undef __FUNCT__
697 #define __FUNCT__ "MatCholeskyFactorSymbolic_Elemental"
698 static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat F,Mat A,IS perm,const MatFactorInfo *info)
699 {
700   PetscFunctionBegin;
701   /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
702   PetscFunctionReturn(0);
703 }
704 
705 #undef __FUNCT__
706 #define __FUNCT__ "MatFactorGetSolverPackage_elemental_elemental"
707 PetscErrorCode MatFactorGetSolverPackage_elemental_elemental(Mat A,const MatSolverPackage *type)
708 {
709   PetscFunctionBegin;
710   *type = MATSOLVERELEMENTAL;
711   PetscFunctionReturn(0);
712 }
713 
714 #undef __FUNCT__
715 #define __FUNCT__ "MatGetFactor_elemental_elemental"
716 static PetscErrorCode MatGetFactor_elemental_elemental(Mat A,MatFactorType ftype,Mat *F)
717 {
718   Mat            B;
719   PetscErrorCode ierr;
720 
721   PetscFunctionBegin;
722   /* Create the factorization matrix */
723   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
724   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
725   ierr = MatSetType(B,MATELEMENTAL);CHKERRQ(ierr);
726   ierr = MatSetUp(B);CHKERRQ(ierr);
727   B->factortype = ftype;
728   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_elemental_elemental);CHKERRQ(ierr);
729   *F            = B;
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNCT__
734 #define __FUNCT__ "MatSolverPackageRegister_Elemental"
735 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_Elemental(void)
736 {
737   PetscErrorCode ierr;
738 
739   PetscFunctionBegin;
740   ierr = MatSolverPackageRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_LU,MatGetFactor_elemental_elemental);CHKERRQ(ierr);
741   ierr = MatSolverPackageRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_CHOLESKY,MatGetFactor_elemental_elemental);CHKERRQ(ierr);
742   PetscFunctionReturn(0);
743 }
744 
745 #undef __FUNCT__
746 #define __FUNCT__ "MatNorm_Elemental"
747 static PetscErrorCode MatNorm_Elemental(Mat A,NormType type,PetscReal *nrm)
748 {
749   Mat_Elemental *a=(Mat_Elemental*)A->data;
750 
751   PetscFunctionBegin;
752   switch (type){
753   case NORM_1:
754     *nrm = El::OneNorm(*a->emat);
755     break;
756   case NORM_FROBENIUS:
757     *nrm = El::FrobeniusNorm(*a->emat);
758     break;
759   case NORM_INFINITY:
760     *nrm = El::InfinityNorm(*a->emat);
761     break;
762   default:
763     printf("Error: unsupported norm type!\n");
764   }
765   PetscFunctionReturn(0);
766 }
767 
768 #undef __FUNCT__
769 #define __FUNCT__ "MatZeroEntries_Elemental"
770 static PetscErrorCode MatZeroEntries_Elemental(Mat A)
771 {
772   Mat_Elemental *a=(Mat_Elemental*)A->data;
773 
774   PetscFunctionBegin;
775   El::Zero(*a->emat);
776   PetscFunctionReturn(0);
777 }
778 
779 #undef __FUNCT__
780 #define __FUNCT__ "MatGetOwnershipIS_Elemental"
781 static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A,IS *rows,IS *cols)
782 {
783   Mat_Elemental  *a = (Mat_Elemental*)A->data;
784   PetscErrorCode ierr;
785   PetscInt       i,m,shift,stride,*idx;
786 
787   PetscFunctionBegin;
788   if (rows) {
789     m = a->emat->LocalHeight();
790     shift = a->emat->ColShift();
791     stride = a->emat->ColStride();
792     ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr);
793     for (i=0; i<m; i++) {
794       PetscInt rank,offset;
795       E2RO(A,0,shift+i*stride,&rank,&offset);
796       RO2P(A,0,rank,offset,&idx[i]);
797     }
798     ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,rows);CHKERRQ(ierr);
799   }
800   if (cols) {
801     m = a->emat->LocalWidth();
802     shift = a->emat->RowShift();
803     stride = a->emat->RowStride();
804     ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr);
805     for (i=0; i<m; i++) {
806       PetscInt rank,offset;
807       E2RO(A,1,shift+i*stride,&rank,&offset);
808       RO2P(A,1,rank,offset,&idx[i]);
809     }
810     ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,cols);CHKERRQ(ierr);
811   }
812   PetscFunctionReturn(0);
813 }
814 
815 #undef __FUNCT__
816 #define __FUNCT__ "MatConvert_Elemental_Dense"
817 static PetscErrorCode MatConvert_Elemental_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
818 {
819   Mat                Bmpi;
820   Mat_Elemental      *a = (Mat_Elemental*)A->data;
821   MPI_Comm           comm;
822   PetscErrorCode     ierr;
823   PetscInt           rrank,ridx,crank,cidx,nrows,ncols,i,j;
824   PetscElemScalar    v;
825   PetscBool          s1,s2,s3;
826 
827   PetscFunctionBegin;
828   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
829 
830   ierr = PetscStrcmp(newtype,MATDENSE,&s1);CHKERRQ(ierr);
831   ierr = PetscStrcmp(newtype,MATSEQDENSE,&s2);CHKERRQ(ierr);
832   ierr = PetscStrcmp(newtype,MATMPIDENSE,&s3);CHKERRQ(ierr);
833   if (!s1 && !s2 && !s3) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported New MatType: must be MATDENSE, MATSEQDENSE or MATMPIDENSE");
834   ierr = MatCreate(comm,&Bmpi);CHKERRQ(ierr);
835   ierr = MatSetSizes(Bmpi,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
836   ierr = MatSetType(Bmpi,MATDENSE);CHKERRQ(ierr);
837   ierr = MatSetUp(Bmpi);CHKERRQ(ierr);
838   ierr = MatGetSize(A,&nrows,&ncols);CHKERRQ(ierr); /* global nrows and ncols */
839   for (i=0; i<nrows; i++) {
840     PetscInt erow,ecol;
841     P2RO(A,0,i,&rrank,&ridx);
842     RO2E(A,0,rrank,ridx,&erow);
843     if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
844     for (j=0; j<ncols; j++) {
845       P2RO(A,1,j,&crank,&cidx);
846       RO2E(A,1,crank,cidx,&ecol);
847       if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
848       v = a->emat->Get(erow,ecol);
849 #if 0
850       PetscMPIInt rank;
851       printf("[%d] Ae(%d, %d, %g) -> B(%d, %d)\n",rank,erow,ecol,v,i,j);
852 #endif
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   if (--commgrid->grid_refct == 0) {
1029     delete commgrid->grid;
1030     ierr = PetscFree(commgrid);CHKERRQ(ierr);
1031     ierr = MPI_Keyval_free(&Petsc_Elemental_keyval);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     /* ierr = MPI_Comm_create_Keyval(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0); -- new version? */
1378   }
1379   ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr);
1380   ierr = MPI_Attr_get(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr);
1381   if (!flg) {
1382     ierr = PetscNewLog(A,&commgrid);CHKERRQ(ierr);
1383 
1384     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Elemental Options","Mat");CHKERRQ(ierr);
1385     /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */
1386     ierr = PetscOptionsInt("-mat_elemental_grid_height","Grid Height","None",El::mpi::Size(cxxcomm),&optv1,&flg1);CHKERRQ(ierr);
1387     if (flg1) {
1388       if (El::mpi::Size(cxxcomm) % optv1 != 0) {
1389         SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Grid Height %D must evenly divide CommSize %D",optv1,(PetscInt)El::mpi::Size(cxxcomm));
1390       }
1391       commgrid->grid = new El::Grid(cxxcomm,optv1); /* use user-provided grid height */
1392     } else {
1393       commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */
1394       /* printf("new commgrid->grid = %p\n",commgrid->grid);  -- memory leak revealed by valgrind? */
1395     }
1396     commgrid->grid_refct = 1;
1397     ierr = MPI_Attr_put(icomm,Petsc_Elemental_keyval,(void*)commgrid);CHKERRQ(ierr);
1398     PetscOptionsEnd();
1399   } else {
1400     commgrid->grid_refct++;
1401   }
1402   ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr);
1403   a->grid      = commgrid->grid;
1404   a->emat      = new El::DistMatrix<PetscElemScalar>(*a->grid);
1405   a->esubmat   = new El::Matrix<PetscElemScalar>(1,1);
1406   a->interface = new El::AxpyInterface<PetscElemScalar>;
1407   a->pivot     = new El::DistMatrix<PetscInt,El::VC,El::STAR>;
1408 
1409   /* build cache for off array entries formed */
1410   a->interface->Attach(El::LOCAL_TO_GLOBAL,*(a->emat));
1411 
1412   ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",MatGetOwnershipIS_Elemental);CHKERRQ(ierr);
1413   ierr = PetscObjectComposeFunction((PetscObject)A,"MatElementalHermitianGenDefEig_C",MatElementalHermitianGenDefEig_Elemental);CHKERRQ(ierr);
1414 
1415   ierr = PetscObjectChangeTypeName((PetscObject)A,MATELEMENTAL);CHKERRQ(ierr);
1416   PetscFunctionReturn(0);
1417 }
1418