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