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