xref: /petsc/src/mat/impls/elemental/matelem.cxx (revision f6e5db1b5bb81357eff7eb8aee412106c0d7d08b)
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 PetscErrorCode MatMatMultNumeric_Elemental_MPIDense(Mat A,Mat B,Mat C)
383 {
384   Mat            Be,Ce;
385   PetscErrorCode ierr;
386 
387   PetscFunctionBegin;
388   ierr = MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX,&Be);CHKERRQ(ierr);
389   ierr = MatMatMult(A,Be,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Ce);CHKERRQ(ierr);
390   ierr = MatConvert(Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
391   ierr = MatDestroy(&Be);CHKERRQ(ierr);
392   ierr = MatDestroy(&Ce);CHKERRQ(ierr);
393   PetscFunctionReturn(0);
394 }
395 
396 PetscErrorCode MatMatMultSymbolic_Elemental_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
397 {
398   PetscErrorCode ierr;
399 
400   PetscFunctionBegin;
401   ierr = MatSetSizes(C,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
402   ierr = MatSetType(C,MATMPIDENSE);CHKERRQ(ierr);
403   ierr = MatSetUp(C);CHKERRQ(ierr);
404   C->ops->matmultnumeric = MatMatMultNumeric_Elemental_MPIDense;
405   PetscFunctionReturn(0);
406 }
407 
408 PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense_AB(Mat C)
409 {
410   PetscFunctionBegin;
411   C->ops->matmultsymbolic = MatMatMultSymbolic_Elemental_MPIDense;
412   C->ops->productsymbolic = MatProductSymbolic_AB;
413   PetscFunctionReturn(0);
414 }
415 
416 PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense(Mat C)
417 {
418   PetscErrorCode ierr;
419   Mat_Product    *product = C->product;
420 
421   PetscFunctionBegin;
422   if (product->type == MATPRODUCT_AB) {
423     ierr = MatProductSetFromOptions_Elemental_MPIDense_AB(C);CHKERRQ(ierr);
424   }
425   PetscFunctionReturn(0);
426 }
427 /* --------------------------------------- */
428 
429 static PetscErrorCode MatGetDiagonal_Elemental(Mat A,Vec D)
430 {
431   PetscInt        i,nrows,ncols,nD,rrank,ridx,crank,cidx;
432   Mat_Elemental   *a = (Mat_Elemental*)A->data;
433   PetscErrorCode  ierr;
434   PetscElemScalar v;
435   MPI_Comm        comm;
436 
437   PetscFunctionBegin;
438   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
439   ierr = MatGetSize(A,&nrows,&ncols);CHKERRQ(ierr);
440   nD = nrows>ncols ? ncols : nrows;
441   for (i=0; i<nD; i++) {
442     PetscInt erow,ecol;
443     P2RO(A,0,i,&rrank,&ridx);
444     RO2E(A,0,rrank,ridx,&erow);
445     if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
446     P2RO(A,1,i,&crank,&cidx);
447     RO2E(A,1,crank,cidx,&ecol);
448     if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
449     v = a->emat->Get(erow,ecol);
450     ierr = VecSetValues(D,1,&i,(PetscScalar*)&v,INSERT_VALUES);CHKERRQ(ierr);
451   }
452   ierr = VecAssemblyBegin(D);CHKERRQ(ierr);
453   ierr = VecAssemblyEnd(D);CHKERRQ(ierr);
454   PetscFunctionReturn(0);
455 }
456 
457 static PetscErrorCode MatDiagonalScale_Elemental(Mat X,Vec L,Vec R)
458 {
459   Mat_Elemental         *x = (Mat_Elemental*)X->data;
460   const PetscElemScalar *d;
461   PetscErrorCode        ierr;
462 
463   PetscFunctionBegin;
464   if (R) {
465     ierr = VecGetArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr);
466     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
467     de.LockedAttach(X->cmap->N,1,*x->grid,0,0,d,X->cmap->n);
468     El::DiagonalScale(El::RIGHT,El::NORMAL,de,*x->emat);
469     ierr = VecRestoreArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr);
470   }
471   if (L) {
472     ierr = VecGetArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr);
473     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
474     de.LockedAttach(X->rmap->N,1,*x->grid,0,0,d,X->rmap->n);
475     El::DiagonalScale(El::LEFT,El::NORMAL,de,*x->emat);
476     ierr = VecRestoreArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr);
477   }
478   PetscFunctionReturn(0);
479 }
480 
481 static PetscErrorCode MatMissingDiagonal_Elemental(Mat A,PetscBool *missing,PetscInt *d)
482 {
483   PetscFunctionBegin;
484   *missing = PETSC_FALSE;
485   PetscFunctionReturn(0);
486 }
487 
488 static PetscErrorCode MatScale_Elemental(Mat X,PetscScalar a)
489 {
490   Mat_Elemental  *x = (Mat_Elemental*)X->data;
491 
492   PetscFunctionBegin;
493   El::Scale((PetscElemScalar)a,*x->emat);
494   PetscFunctionReturn(0);
495 }
496 
497 /*
498   MatAXPY - Computes Y = a*X + Y.
499 */
500 static PetscErrorCode MatAXPY_Elemental(Mat Y,PetscScalar a,Mat X,MatStructure str)
501 {
502   Mat_Elemental  *x = (Mat_Elemental*)X->data;
503   Mat_Elemental  *y = (Mat_Elemental*)Y->data;
504   PetscErrorCode ierr;
505 
506   PetscFunctionBegin;
507   El::Axpy((PetscElemScalar)a,*x->emat,*y->emat);
508   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
509   PetscFunctionReturn(0);
510 }
511 
512 static PetscErrorCode MatCopy_Elemental(Mat A,Mat B,MatStructure str)
513 {
514   Mat_Elemental *a=(Mat_Elemental*)A->data;
515   Mat_Elemental *b=(Mat_Elemental*)B->data;
516   PetscErrorCode ierr;
517 
518   PetscFunctionBegin;
519   El::Copy(*a->emat,*b->emat);
520   ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 static PetscErrorCode MatDuplicate_Elemental(Mat A,MatDuplicateOption op,Mat *B)
525 {
526   Mat            Be;
527   MPI_Comm       comm;
528   Mat_Elemental  *a=(Mat_Elemental*)A->data;
529   PetscErrorCode ierr;
530 
531   PetscFunctionBegin;
532   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
533   ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
534   ierr = MatSetSizes(Be,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
535   ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
536   ierr = MatSetUp(Be);CHKERRQ(ierr);
537   *B = Be;
538   if (op == MAT_COPY_VALUES) {
539     Mat_Elemental *b=(Mat_Elemental*)Be->data;
540     El::Copy(*a->emat,*b->emat);
541   }
542   Be->assembled = PETSC_TRUE;
543   PetscFunctionReturn(0);
544 }
545 
546 static PetscErrorCode MatTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
547 {
548   Mat            Be = *B;
549   PetscErrorCode ierr;
550   MPI_Comm       comm;
551   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
552 
553   PetscFunctionBegin;
554   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
555   /* Only out-of-place supported */
556   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Only out-of-place supported");
557   if (reuse == MAT_INITIAL_MATRIX) {
558     ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
559     ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
560     ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
561     ierr = MatSetUp(Be);CHKERRQ(ierr);
562     *B = Be;
563   }
564   b = (Mat_Elemental*)Be->data;
565   El::Transpose(*a->emat,*b->emat);
566   Be->assembled = PETSC_TRUE;
567   PetscFunctionReturn(0);
568 }
569 
570 static PetscErrorCode MatConjugate_Elemental(Mat A)
571 {
572   Mat_Elemental  *a = (Mat_Elemental*)A->data;
573 
574   PetscFunctionBegin;
575   El::Conjugate(*a->emat);
576   PetscFunctionReturn(0);
577 }
578 
579 static PetscErrorCode MatHermitianTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
580 {
581   Mat            Be = *B;
582   PetscErrorCode ierr;
583   MPI_Comm       comm;
584   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
585 
586   PetscFunctionBegin;
587   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
588   /* Only out-of-place supported */
589   if (reuse == MAT_INITIAL_MATRIX){
590     ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
591     ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
592     ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
593     ierr = MatSetUp(Be);CHKERRQ(ierr);
594     *B = Be;
595   }
596   b = (Mat_Elemental*)Be->data;
597   El::Adjoint(*a->emat,*b->emat);
598   Be->assembled = PETSC_TRUE;
599   PetscFunctionReturn(0);
600 }
601 
602 static PetscErrorCode MatSolve_Elemental(Mat A,Vec B,Vec X)
603 {
604   Mat_Elemental     *a = (Mat_Elemental*)A->data;
605   PetscErrorCode    ierr;
606   PetscElemScalar   *x;
607   PetscInt          pivoting = a->pivoting;
608 
609   PetscFunctionBegin;
610   ierr = VecCopy(B,X);CHKERRQ(ierr);
611   ierr = VecGetArray(X,(PetscScalar **)&x);CHKERRQ(ierr);
612 
613   El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe;
614   xe.Attach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
615   El::DistMatrix<PetscElemScalar,El::MC,El::MR> xer(xe);
616   switch (A->factortype) {
617   case MAT_FACTOR_LU:
618     if (pivoting == 0) {
619       El::lu::SolveAfter(El::NORMAL,*a->emat,xer);
620     } else if (pivoting == 1) {
621       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,xer);
622     } else { /* pivoting == 2 */
623       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,xer);
624     }
625     break;
626   case MAT_FACTOR_CHOLESKY:
627     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,xer);
628     break;
629   default:
630     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
631     break;
632   }
633   El::Copy(xer,xe);
634 
635   ierr = VecRestoreArray(X,(PetscScalar **)&x);CHKERRQ(ierr);
636   PetscFunctionReturn(0);
637 }
638 
639 static PetscErrorCode MatSolveAdd_Elemental(Mat A,Vec B,Vec Y,Vec X)
640 {
641   PetscErrorCode    ierr;
642 
643   PetscFunctionBegin;
644   ierr = MatSolve_Elemental(A,B,X);CHKERRQ(ierr);
645   ierr = VecAXPY(X,1,Y);CHKERRQ(ierr);
646   PetscFunctionReturn(0);
647 }
648 
649 static PetscErrorCode MatMatSolve_Elemental(Mat A,Mat B,Mat X)
650 {
651   Mat_Elemental *a=(Mat_Elemental*)A->data;
652   Mat_Elemental *b=(Mat_Elemental*)B->data;
653   Mat_Elemental *x=(Mat_Elemental*)X->data;
654   PetscInt      pivoting = a->pivoting;
655 
656   PetscFunctionBegin;
657   El::Copy(*b->emat,*x->emat);
658   switch (A->factortype) {
659   case MAT_FACTOR_LU:
660     if (pivoting == 0) {
661       El::lu::SolveAfter(El::NORMAL,*a->emat,*x->emat);
662     } else if (pivoting == 1) {
663       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*x->emat);
664     } else {
665       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,*x->emat);
666     }
667     break;
668   case MAT_FACTOR_CHOLESKY:
669     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,*x->emat);
670     break;
671   default:
672     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
673     break;
674   }
675   PetscFunctionReturn(0);
676 }
677 
678 static PetscErrorCode MatLUFactor_Elemental(Mat A,IS row,IS col,const MatFactorInfo *info)
679 {
680   Mat_Elemental  *a = (Mat_Elemental*)A->data;
681   PetscErrorCode ierr;
682   PetscInt       pivoting = a->pivoting;
683 
684   PetscFunctionBegin;
685   if (pivoting == 0) {
686     El::LU(*a->emat);
687   } else if (pivoting == 1) {
688     El::LU(*a->emat,*a->P);
689   } else {
690     El::LU(*a->emat,*a->P,*a->Q);
691   }
692   A->factortype = MAT_FACTOR_LU;
693   A->assembled  = PETSC_TRUE;
694 
695   ierr = PetscFree(A->solvertype);CHKERRQ(ierr);
696   ierr = PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype);CHKERRQ(ierr);
697   PetscFunctionReturn(0);
698 }
699 
700 static PetscErrorCode  MatLUFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
701 {
702   PetscErrorCode ierr;
703 
704   PetscFunctionBegin;
705   ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
706   ierr = MatLUFactor_Elemental(F,0,0,info);CHKERRQ(ierr);
707   PetscFunctionReturn(0);
708 }
709 
710 static PetscErrorCode  MatLUFactorSymbolic_Elemental(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
711 {
712   PetscFunctionBegin;
713   /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
714   PetscFunctionReturn(0);
715 }
716 
717 static PetscErrorCode MatCholeskyFactor_Elemental(Mat A,IS perm,const MatFactorInfo *info)
718 {
719   Mat_Elemental  *a = (Mat_Elemental*)A->data;
720   El::DistMatrix<PetscElemScalar,El::MC,El::STAR> d;
721   PetscErrorCode ierr;
722 
723   PetscFunctionBegin;
724   El::Cholesky(El::UPPER,*a->emat);
725   A->factortype = MAT_FACTOR_CHOLESKY;
726   A->assembled  = PETSC_TRUE;
727 
728   ierr = PetscFree(A->solvertype);CHKERRQ(ierr);
729   ierr = PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype);CHKERRQ(ierr);
730   PetscFunctionReturn(0);
731 }
732 
733 static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
734 {
735   PetscErrorCode ierr;
736 
737   PetscFunctionBegin;
738   ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
739   ierr = MatCholeskyFactor_Elemental(F,0,info);CHKERRQ(ierr);
740   PetscFunctionReturn(0);
741 }
742 
743 static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat F,Mat A,IS perm,const MatFactorInfo *info)
744 {
745   PetscFunctionBegin;
746   /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
747   PetscFunctionReturn(0);
748 }
749 
750 PetscErrorCode MatFactorGetSolverType_elemental_elemental(Mat A,MatSolverType *type)
751 {
752   PetscFunctionBegin;
753   *type = MATSOLVERELEMENTAL;
754   PetscFunctionReturn(0);
755 }
756 
757 static PetscErrorCode MatGetFactor_elemental_elemental(Mat A,MatFactorType ftype,Mat *F)
758 {
759   Mat            B;
760   PetscErrorCode ierr;
761 
762   PetscFunctionBegin;
763   /* Create the factorization matrix */
764   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
765   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
766   ierr = MatSetType(B,MATELEMENTAL);CHKERRQ(ierr);
767   ierr = MatSetUp(B);CHKERRQ(ierr);
768   B->factortype = ftype;
769   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
770   ierr = PetscStrallocpy(MATSOLVERELEMENTAL,&B->solvertype);CHKERRQ(ierr);
771 
772   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_elemental_elemental);CHKERRQ(ierr);
773   *F            = B;
774   PetscFunctionReturn(0);
775 }
776 
777 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_Elemental(void)
778 {
779   PetscErrorCode ierr;
780 
781   PetscFunctionBegin;
782   ierr = MatSolverTypeRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_LU,MatGetFactor_elemental_elemental);CHKERRQ(ierr);
783   ierr = MatSolverTypeRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_CHOLESKY,MatGetFactor_elemental_elemental);CHKERRQ(ierr);
784   PetscFunctionReturn(0);
785 }
786 
787 static PetscErrorCode MatNorm_Elemental(Mat A,NormType type,PetscReal *nrm)
788 {
789   Mat_Elemental *a=(Mat_Elemental*)A->data;
790 
791   PetscFunctionBegin;
792   switch (type){
793   case NORM_1:
794     *nrm = El::OneNorm(*a->emat);
795     break;
796   case NORM_FROBENIUS:
797     *nrm = El::FrobeniusNorm(*a->emat);
798     break;
799   case NORM_INFINITY:
800     *nrm = El::InfinityNorm(*a->emat);
801     break;
802   default:
803     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported norm type");
804   }
805   PetscFunctionReturn(0);
806 }
807 
808 static PetscErrorCode MatZeroEntries_Elemental(Mat A)
809 {
810   Mat_Elemental *a=(Mat_Elemental*)A->data;
811 
812   PetscFunctionBegin;
813   El::Zero(*a->emat);
814   PetscFunctionReturn(0);
815 }
816 
817 static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A,IS *rows,IS *cols)
818 {
819   Mat_Elemental  *a = (Mat_Elemental*)A->data;
820   PetscErrorCode ierr;
821   PetscInt       i,m,shift,stride,*idx;
822 
823   PetscFunctionBegin;
824   if (rows) {
825     m = a->emat->LocalHeight();
826     shift = a->emat->ColShift();
827     stride = a->emat->ColStride();
828     ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr);
829     for (i=0; i<m; i++) {
830       PetscInt rank,offset;
831       E2RO(A,0,shift+i*stride,&rank,&offset);
832       RO2P(A,0,rank,offset,&idx[i]);
833     }
834     ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,rows);CHKERRQ(ierr);
835   }
836   if (cols) {
837     m = a->emat->LocalWidth();
838     shift = a->emat->RowShift();
839     stride = a->emat->RowStride();
840     ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr);
841     for (i=0; i<m; i++) {
842       PetscInt rank,offset;
843       E2RO(A,1,shift+i*stride,&rank,&offset);
844       RO2P(A,1,rank,offset,&idx[i]);
845     }
846     ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,cols);CHKERRQ(ierr);
847   }
848   PetscFunctionReturn(0);
849 }
850 
851 static PetscErrorCode MatConvert_Elemental_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
852 {
853   Mat                Bmpi;
854   Mat_Elemental      *a = (Mat_Elemental*)A->data;
855   MPI_Comm           comm;
856   PetscErrorCode     ierr;
857   IS                 isrows,iscols;
858   PetscInt           rrank,ridx,crank,cidx,nrows,ncols,i,j,erow,ecol,elrow,elcol;
859   const PetscInt     *rows,*cols;
860   PetscElemScalar    v;
861   const El::Grid     &grid = a->emat->Grid();
862 
863   PetscFunctionBegin;
864   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
865 
866   if (reuse == MAT_REUSE_MATRIX) {
867     Bmpi = *B;
868   } else {
869     ierr = MatCreate(comm,&Bmpi);CHKERRQ(ierr);
870     ierr = MatSetSizes(Bmpi,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
871     ierr = MatSetType(Bmpi,MATDENSE);CHKERRQ(ierr);
872     ierr = MatSetUp(Bmpi);CHKERRQ(ierr);
873   }
874 
875   /* Get local entries of A */
876   ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr);
877   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
878   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
879   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
880   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
881 
882   if (a->roworiented) {
883     for (i=0; i<nrows; i++) {
884       P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
885       RO2E(A,0,rrank,ridx,&erow);
886       if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
887       for (j=0; j<ncols; j++) {
888         P2RO(A,1,cols[j],&crank,&cidx);
889         RO2E(A,1,crank,cidx,&ecol);
890         if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
891 
892         elrow = erow / grid.MCSize(); /* Elemental local row index */
893         elcol = ecol / grid.MRSize(); /* Elemental local column index */
894         v = a->emat->GetLocal(elrow,elcol);
895         ierr = MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES);CHKERRQ(ierr);
896       }
897     }
898   } else { /* column-oriented */
899     for (j=0; j<ncols; j++) {
900       P2RO(A,1,cols[j],&crank,&cidx);
901       RO2E(A,1,crank,cidx,&ecol);
902       if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
903       for (i=0; i<nrows; i++) {
904         P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
905         RO2E(A,0,rrank,ridx,&erow);
906         if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
907 
908         elrow = erow / grid.MCSize(); /* Elemental local row index */
909         elcol = ecol / grid.MRSize(); /* Elemental local column index */
910         v = a->emat->GetLocal(elrow,elcol);
911         ierr = MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES);CHKERRQ(ierr);
912       }
913     }
914   }
915   ierr = MatAssemblyBegin(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
916   ierr = MatAssemblyEnd(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
917   if (reuse == MAT_INPLACE_MATRIX) {
918     ierr = MatHeaderReplace(A,&Bmpi);CHKERRQ(ierr);
919   } else {
920     *B = Bmpi;
921   }
922   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
923   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
928 {
929   Mat               mat_elemental;
930   PetscErrorCode    ierr;
931   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols;
932   const PetscInt    *cols;
933   const PetscScalar *vals;
934 
935   PetscFunctionBegin;
936   if (reuse == MAT_REUSE_MATRIX) {
937     mat_elemental = *newmat;
938     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
939   } else {
940     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
941     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
942     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
943     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
944   }
945   for (row=0; row<M; row++) {
946     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
947     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
948     ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
949     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
950   }
951   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
952   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
953 
954   if (reuse == MAT_INPLACE_MATRIX) {
955     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
956   } else {
957     *newmat = mat_elemental;
958   }
959   PetscFunctionReturn(0);
960 }
961 
962 PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
963 {
964   Mat               mat_elemental;
965   PetscErrorCode    ierr;
966   PetscInt          row,ncols,rstart=A->rmap->rstart,rend=A->rmap->rend,j;
967   const PetscInt    *cols;
968   const PetscScalar *vals;
969 
970   PetscFunctionBegin;
971   if (reuse == MAT_REUSE_MATRIX) {
972     mat_elemental = *newmat;
973     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
974   } else {
975     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
976     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
977     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
978     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
979   }
980   for (row=rstart; row<rend; row++) {
981     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
982     for (j=0; j<ncols; j++) {
983       /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
984       ierr = MatSetValues(mat_elemental,1,&row,1,&cols[j],&vals[j],ADD_VALUES);CHKERRQ(ierr);
985     }
986     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
987   }
988   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
989   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
990 
991   if (reuse == MAT_INPLACE_MATRIX) {
992     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
993   } else {
994     *newmat = mat_elemental;
995   }
996   PetscFunctionReturn(0);
997 }
998 
999 PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
1000 {
1001   Mat               mat_elemental;
1002   PetscErrorCode    ierr;
1003   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols,j;
1004   const PetscInt    *cols;
1005   const PetscScalar *vals;
1006 
1007   PetscFunctionBegin;
1008   if (reuse == MAT_REUSE_MATRIX) {
1009     mat_elemental = *newmat;
1010     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
1011   } else {
1012     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
1013     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1014     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
1015     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
1016   }
1017   ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
1018   for (row=0; row<M; row++) {
1019     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1020     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
1021     ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1022     for (j=0; j<ncols; j++) { /* lower triangular part */
1023       PetscScalar v;
1024       if (cols[j] == row) continue;
1025       v    = A->hermitian ? PetscConj(vals[j]) : vals[j];
1026       ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&v,ADD_VALUES);CHKERRQ(ierr);
1027     }
1028     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1029   }
1030   ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
1031   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1032   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1033 
1034   if (reuse == MAT_INPLACE_MATRIX) {
1035     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
1036   } else {
1037     *newmat = mat_elemental;
1038   }
1039   PetscFunctionReturn(0);
1040 }
1041 
1042 PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
1043 {
1044   Mat               mat_elemental;
1045   PetscErrorCode    ierr;
1046   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols,j,rstart=A->rmap->rstart,rend=A->rmap->rend;
1047   const PetscInt    *cols;
1048   const PetscScalar *vals;
1049 
1050   PetscFunctionBegin;
1051   if (reuse == MAT_REUSE_MATRIX) {
1052     mat_elemental = *newmat;
1053     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
1054   } else {
1055     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
1056     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1057     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
1058     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
1059   }
1060   ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
1061   for (row=rstart; row<rend; row++) {
1062     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1063     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
1064     ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1065     for (j=0; j<ncols; j++) { /* lower triangular part */
1066       PetscScalar v;
1067       if (cols[j] == row) continue;
1068       v    = A->hermitian ? PetscConj(vals[j]) : vals[j];
1069       ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&v,ADD_VALUES);CHKERRQ(ierr);
1070     }
1071     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1072   }
1073   ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
1074   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1075   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1076 
1077   if (reuse == MAT_INPLACE_MATRIX) {
1078     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
1079   } else {
1080     *newmat = mat_elemental;
1081   }
1082   PetscFunctionReturn(0);
1083 }
1084 
1085 static PetscErrorCode MatDestroy_Elemental(Mat A)
1086 {
1087   Mat_Elemental      *a = (Mat_Elemental*)A->data;
1088   PetscErrorCode     ierr;
1089   Mat_Elemental_Grid *commgrid;
1090   PetscBool          flg;
1091   MPI_Comm           icomm;
1092 
1093   PetscFunctionBegin;
1094   delete a->emat;
1095   delete a->P;
1096   delete a->Q;
1097 
1098   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
1099   ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr);
1100   ierr = MPI_Comm_get_attr(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr);
1101   if (--commgrid->grid_refct == 0) {
1102     delete commgrid->grid;
1103     ierr = PetscFree(commgrid);CHKERRQ(ierr);
1104     ierr = MPI_Comm_free_keyval(&Petsc_Elemental_keyval);CHKERRQ(ierr);
1105   }
1106   ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr);
1107   ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",NULL);CHKERRQ(ierr);
1108   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr);
1109   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_elemental_mpidense_C",NULL);CHKERRQ(ierr);
1110   ierr = PetscFree(A->data);CHKERRQ(ierr);
1111   PetscFunctionReturn(0);
1112 }
1113 
1114 PetscErrorCode MatSetUp_Elemental(Mat A)
1115 {
1116   Mat_Elemental  *a = (Mat_Elemental*)A->data;
1117   PetscErrorCode ierr;
1118   MPI_Comm       comm;
1119   PetscMPIInt    rsize,csize;
1120   PetscInt       n;
1121 
1122   PetscFunctionBegin;
1123   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1124   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1125 
1126   /* Check if local row and column sizes are equally distributed.
1127      Jed: Elemental uses "element" cyclic ordering so the sizes need to match that
1128      exactly.  The strategy in MatElemental is for PETSc to implicitly permute to block ordering (like would be returned by
1129      PetscSplitOwnership(comm,&n,&N), at which point Elemental matrices can act on PETSc vectors without redistributing the vectors. */
1130   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1131   n = PETSC_DECIDE;
1132   ierr = PetscSplitOwnership(comm,&n,&A->rmap->N);CHKERRQ(ierr);
1133   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);
1134 
1135   n = PETSC_DECIDE;
1136   ierr = PetscSplitOwnership(comm,&n,&A->cmap->N);CHKERRQ(ierr);
1137   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);
1138 
1139   a->emat->Resize(A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1140   El::Zero(*a->emat);
1141 
1142   ierr = MPI_Comm_size(A->rmap->comm,&rsize);CHKERRQ(ierr);
1143   ierr = MPI_Comm_size(A->cmap->comm,&csize);CHKERRQ(ierr);
1144   if (csize != rsize) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot use row and column communicators of different sizes");
1145   a->commsize = rsize;
1146   a->mr[0] = A->rmap->N % rsize; if (!a->mr[0]) a->mr[0] = rsize;
1147   a->mr[1] = A->cmap->N % csize; if (!a->mr[1]) a->mr[1] = csize;
1148   a->m[0]  = A->rmap->N / rsize + (a->mr[0] != rsize);
1149   a->m[1]  = A->cmap->N / csize + (a->mr[1] != csize);
1150   PetscFunctionReturn(0);
1151 }
1152 
1153 PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType type)
1154 {
1155   Mat_Elemental  *a = (Mat_Elemental*)A->data;
1156 
1157   PetscFunctionBegin;
1158   /* printf("Calling ProcessQueues\n"); */
1159   a->emat->ProcessQueues();
1160   /* printf("Finished ProcessQueues\n"); */
1161   PetscFunctionReturn(0);
1162 }
1163 
1164 PetscErrorCode MatAssemblyEnd_Elemental(Mat A, MatAssemblyType type)
1165 {
1166   PetscFunctionBegin;
1167   /* Currently does nothing */
1168   PetscFunctionReturn(0);
1169 }
1170 
1171 PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer)
1172 {
1173   PetscErrorCode ierr;
1174   Mat            Adense,Ae;
1175   MPI_Comm       comm;
1176 
1177   PetscFunctionBegin;
1178   ierr = PetscObjectGetComm((PetscObject)newMat,&comm);CHKERRQ(ierr);
1179   ierr = MatCreate(comm,&Adense);CHKERRQ(ierr);
1180   ierr = MatSetType(Adense,MATDENSE);CHKERRQ(ierr);
1181   ierr = MatLoad(Adense,viewer);CHKERRQ(ierr);
1182   ierr = MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX,&Ae);CHKERRQ(ierr);
1183   ierr = MatDestroy(&Adense);CHKERRQ(ierr);
1184   ierr = MatHeaderReplace(newMat,&Ae);CHKERRQ(ierr);
1185   PetscFunctionReturn(0);
1186 }
1187 
1188 /* -------------------------------------------------------------------*/
1189 static struct _MatOps MatOps_Values = {
1190        MatSetValues_Elemental,
1191        0,
1192        0,
1193        MatMult_Elemental,
1194 /* 4*/ MatMultAdd_Elemental,
1195        MatMultTranspose_Elemental,
1196        MatMultTransposeAdd_Elemental,
1197        MatSolve_Elemental,
1198        MatSolveAdd_Elemental,
1199        0,
1200 /*10*/ 0,
1201        MatLUFactor_Elemental,
1202        MatCholeskyFactor_Elemental,
1203        0,
1204        MatTranspose_Elemental,
1205 /*15*/ MatGetInfo_Elemental,
1206        0,
1207        MatGetDiagonal_Elemental,
1208        MatDiagonalScale_Elemental,
1209        MatNorm_Elemental,
1210 /*20*/ MatAssemblyBegin_Elemental,
1211        MatAssemblyEnd_Elemental,
1212        MatSetOption_Elemental,
1213        MatZeroEntries_Elemental,
1214 /*24*/ 0,
1215        MatLUFactorSymbolic_Elemental,
1216        MatLUFactorNumeric_Elemental,
1217        MatCholeskyFactorSymbolic_Elemental,
1218        MatCholeskyFactorNumeric_Elemental,
1219 /*29*/ MatSetUp_Elemental,
1220        0,
1221        0,
1222        0,
1223        0,
1224 /*34*/ MatDuplicate_Elemental,
1225        0,
1226        0,
1227        0,
1228        0,
1229 /*39*/ MatAXPY_Elemental,
1230        0,
1231        0,
1232        0,
1233        MatCopy_Elemental,
1234 /*44*/ 0,
1235        MatScale_Elemental,
1236        MatShift_Basic,
1237        0,
1238        0,
1239 /*49*/ 0,
1240        0,
1241        0,
1242        0,
1243        0,
1244 /*54*/ 0,
1245        0,
1246        0,
1247        0,
1248        0,
1249 /*59*/ 0,
1250        MatDestroy_Elemental,
1251        MatView_Elemental,
1252        0,
1253        0,
1254 /*64*/ 0,
1255        0,
1256        0,
1257        0,
1258        0,
1259 /*69*/ 0,
1260        0,
1261        MatConvert_Elemental_Dense,
1262        0,
1263        0,
1264 /*74*/ 0,
1265        0,
1266        0,
1267        0,
1268        0,
1269 /*79*/ 0,
1270        0,
1271        0,
1272        0,
1273        MatLoad_Elemental,
1274 /*84*/ 0,
1275        0,
1276        0,
1277        0,
1278        0,
1279 /*89*/ 0,
1280        0,
1281        MatMatMultNumeric_Elemental,
1282        0,
1283        0,
1284 /*94*/ 0,
1285        0,
1286        0,
1287        MatMatTransposeMultNumeric_Elemental,
1288        0,
1289 /*99*/ MatProductSetFromOptions_Elemental,
1290        0,
1291        0,
1292        MatConjugate_Elemental,
1293        0,
1294 /*104*/0,
1295        0,
1296        0,
1297        0,
1298        0,
1299 /*109*/MatMatSolve_Elemental,
1300        0,
1301        0,
1302        0,
1303        MatMissingDiagonal_Elemental,
1304 /*114*/0,
1305        0,
1306        0,
1307        0,
1308        0,
1309 /*119*/0,
1310        MatHermitianTranspose_Elemental,
1311        0,
1312        0,
1313        0,
1314 /*124*/0,
1315        0,
1316        0,
1317        0,
1318        0,
1319 /*129*/0,
1320        0,
1321        0,
1322        0,
1323        0,
1324 /*134*/0,
1325        0,
1326        0,
1327        0,
1328        0,
1329        0,
1330 /*140*/0,
1331        0,
1332        0,
1333        0,
1334        0,
1335 /*145*/0,
1336        0,
1337        0
1338 };
1339 
1340 /*MC
1341    MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package
1342 
1343   Use ./configure --download-elemental to install PETSc to use Elemental
1344 
1345   Use -pc_type lu -pc_factor_mat_solver_type elemental to use this direct solver
1346 
1347    Options Database Keys:
1348 + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions()
1349 - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix
1350 
1351   Level: beginner
1352 
1353 .seealso: MATDENSE
1354 M*/
1355 
1356 PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A)
1357 {
1358   Mat_Elemental      *a;
1359   PetscErrorCode     ierr;
1360   PetscBool          flg,flg1;
1361   Mat_Elemental_Grid *commgrid;
1362   MPI_Comm           icomm;
1363   PetscInt           optv1;
1364 
1365   PetscFunctionBegin;
1366   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1367   A->insertmode = NOT_SET_VALUES;
1368 
1369   ierr = PetscNewLog(A,&a);CHKERRQ(ierr);
1370   A->data = (void*)a;
1371 
1372   /* Set up the elemental matrix */
1373   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
1374 
1375   /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */
1376   if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) {
1377     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0);CHKERRQ(ierr);
1378   }
1379   ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr);
1380   ierr = MPI_Comm_get_attr(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr);
1381   if (!flg) {
1382     ierr = PetscNewLog(A,&commgrid);CHKERRQ(ierr);
1383 
1384     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Elemental Options","Mat");CHKERRQ(ierr);
1385     /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */
1386     ierr = PetscOptionsInt("-mat_elemental_grid_height","Grid Height","None",El::mpi::Size(cxxcomm),&optv1,&flg1);CHKERRQ(ierr);
1387     if (flg1) {
1388       if (El::mpi::Size(cxxcomm) % optv1) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Grid Height %D must evenly divide CommSize %D",optv1,(PetscInt)El::mpi::Size(cxxcomm));
1389       commgrid->grid = new El::Grid(cxxcomm,optv1); /* use user-provided grid height */
1390     } else {
1391       commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */
1392       /* printf("new commgrid->grid = %p\n",commgrid->grid);  -- memory leak revealed by valgrind? */
1393     }
1394     commgrid->grid_refct = 1;
1395     ierr = MPI_Comm_set_attr(icomm,Petsc_Elemental_keyval,(void*)commgrid);CHKERRQ(ierr);
1396 
1397     a->pivoting    = 1;
1398     ierr = PetscOptionsInt("-mat_elemental_pivoting","Pivoting","None",a->pivoting,&a->pivoting,NULL);CHKERRQ(ierr);
1399 
1400     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1401   } else {
1402     commgrid->grid_refct++;
1403   }
1404   ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr);
1405   a->grid        = commgrid->grid;
1406   a->emat        = new El::DistMatrix<PetscElemScalar>(*a->grid);
1407   a->roworiented = PETSC_TRUE;
1408 
1409   ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",MatGetOwnershipIS_Elemental);CHKERRQ(ierr);
1410   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_elemental_mpidense_C",MatProductSetFromOptions_Elemental_MPIDense);CHKERRQ(ierr);
1411   ierr = PetscObjectChangeTypeName((PetscObject)A,MATELEMENTAL);CHKERRQ(ierr);
1412   PetscFunctionReturn(0);
1413 }
1414