xref: /petsc/src/mat/impls/elemental/matelem.cxx (revision 7fb60732ccf1806988e9909cd147cdb09c395756)
11673f470SJose E. Roman #include <petsc/private/petscelemental.h>
2db31f6deSJed Brown 
327e75052SPierre Jolivet const char ElementalCitation[] = "@Article{Elemental2012,\n"
427e75052SPierre Jolivet "  author  = {Jack Poulson and Bryan Marker and Jeff R. Hammond and Nichols A. Romero and Robert {v}an~{d}e~{G}eijn},\n"
527e75052SPierre Jolivet "  title   = {Elemental: A New Framework for Distributed Memory Dense Matrix Computations},\n"
627e75052SPierre Jolivet "  journal = {{ACM} Transactions on Mathematical Software},\n"
727e75052SPierre Jolivet "  volume  = {39},\n"
827e75052SPierre Jolivet "  number  = {2},\n"
927e75052SPierre Jolivet "  year    = {2013}\n"
1027e75052SPierre Jolivet "}\n";
1127e75052SPierre Jolivet static PetscBool ElementalCite = PETSC_FALSE;
1227e75052SPierre Jolivet 
135e9f5b67SHong Zhang /*
145e9f5b67SHong Zhang     The variable Petsc_Elemental_keyval is used to indicate an MPI attribute that
155e9f5b67SHong Zhang   is attached to a communicator, in this case the attribute is a Mat_Elemental_Grid
165e9f5b67SHong Zhang */
175e9f5b67SHong Zhang static PetscMPIInt Petsc_Elemental_keyval = MPI_KEYVAL_INVALID;
185e9f5b67SHong Zhang 
19db31f6deSJed Brown static PetscErrorCode MatView_Elemental(Mat A,PetscViewer viewer)
20db31f6deSJed Brown {
21db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
22db31f6deSJed Brown   PetscBool      iascii;
23db31f6deSJed Brown 
24db31f6deSJed Brown   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
26db31f6deSJed Brown   if (iascii) {
27db31f6deSJed Brown     PetscViewerFormat format;
289566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
29db31f6deSJed Brown     if (format == PETSC_VIEWER_ASCII_INFO) {
3079673f7bSHong Zhang       /* call elemental viewing function */
319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Elemental run parameters:\n"));
3263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  allocated entries=%zu\n",(*a->emat).AllocatedMemory()));
339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  grid height=%d, grid width=%d\n",(*a->emat).Grid().Height(),(*a->emat).Grid().Width()));
344fe7bbcaSHong Zhang       if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
3579673f7bSHong Zhang         /* call elemental viewing function */
369566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)viewer),"test matview_elemental 2\n"));
374fe7bbcaSHong Zhang       }
3879673f7bSHong Zhang 
39db31f6deSJed Brown     } else if (format == PETSC_VIEWER_DEFAULT) {
409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
41e8146892SSatish Balay       El::Print( *a->emat, "Elemental matrix (cyclic ordering)");
429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
43834d3fecSHong Zhang       if (A->factortype == MAT_FACTOR_NONE) {
4461119200SXuan Zhou         Mat Adense;
459566063dSJacob Faibussowitsch         PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense));
469566063dSJacob Faibussowitsch         PetscCall(MatView(Adense,viewer));
479566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&Adense));
48834d3fecSHong Zhang       }
49ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Format");
50d2daa67eSHong Zhang   } else {
515cb544a0SHong Zhang     /* convert to dense format and call MatView() */
5261119200SXuan Zhou     Mat Adense;
539566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense));
549566063dSJacob Faibussowitsch     PetscCall(MatView(Adense,viewer));
559566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Adense));
56d2daa67eSHong Zhang   }
57db31f6deSJed Brown   PetscFunctionReturn(0);
58db31f6deSJed Brown }
59db31f6deSJed Brown 
6015767789SHong Zhang static PetscErrorCode MatGetInfo_Elemental(Mat A,MatInfoType flag,MatInfo *info)
61180a43e4SHong Zhang {
6215767789SHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data;
6315767789SHong Zhang 
64180a43e4SHong Zhang   PetscFunctionBegin;
655cb544a0SHong Zhang   info->block_size = 1.0;
6615767789SHong Zhang 
6715767789SHong Zhang   if (flag == MAT_LOCAL) {
683966268fSBarry Smith     info->nz_allocated   = (*a->emat).AllocatedMemory(); /* locally allocated */
6915767789SHong Zhang     info->nz_used        = info->nz_allocated;
7015767789SHong Zhang   } else if (flag == MAT_GLOBAL_MAX) {
711c2dc1cbSBarry Smith     //PetscCall(MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin)));
7215767789SHong Zhang     /* see MatGetInfo_MPIAIJ() for getting global info->nz_allocated! */
7315767789SHong Zhang     //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_MAX not written yet");
7415767789SHong Zhang   } else if (flag == MAT_GLOBAL_SUM) {
7515767789SHong Zhang     //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_SUM not written yet");
763966268fSBarry Smith     info->nz_allocated   = (*a->emat).AllocatedMemory(); /* locally allocated */
7715767789SHong Zhang     info->nz_used        = info->nz_allocated; /* assume Elemental does accurate allocation */
781c2dc1cbSBarry Smith     //PetscCall(MPIU_Allreduce(isend,irecv,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A)));
7915767789SHong Zhang     //PetscPrintf(PETSC_COMM_SELF,"    ... [%d] locally allocated %g\n",rank,info->nz_allocated);
8015767789SHong Zhang   }
8115767789SHong Zhang 
8215767789SHong Zhang   info->nz_unneeded       = 0.0;
833966268fSBarry Smith   info->assemblies        = A->num_ass;
8415767789SHong Zhang   info->mallocs           = 0;
8515767789SHong Zhang   info->memory            = ((PetscObject)A)->mem;
8615767789SHong Zhang   info->fill_ratio_given  = 0; /* determined by Elemental */
8715767789SHong Zhang   info->fill_ratio_needed = 0;
8815767789SHong Zhang   info->factor_mallocs    = 0;
89180a43e4SHong Zhang   PetscFunctionReturn(0);
90180a43e4SHong Zhang }
91180a43e4SHong Zhang 
9232d7a744SHong Zhang PetscErrorCode MatSetOption_Elemental(Mat A,MatOption op,PetscBool flg)
9332d7a744SHong Zhang {
9432d7a744SHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data;
9532d7a744SHong Zhang 
9632d7a744SHong Zhang   PetscFunctionBegin;
9732d7a744SHong Zhang   switch (op) {
9832d7a744SHong Zhang   case MAT_NEW_NONZERO_LOCATIONS:
9932d7a744SHong Zhang   case MAT_NEW_NONZERO_LOCATION_ERR:
10032d7a744SHong Zhang   case MAT_NEW_NONZERO_ALLOCATION_ERR:
101891a0571SHong Zhang   case MAT_SYMMETRIC:
102071fcb05SBarry Smith   case MAT_SORTED_FULL:
103eb1ec7c1SStefano Zampini   case MAT_HERMITIAN:
104891a0571SHong Zhang     break;
10532d7a744SHong Zhang   case MAT_ROW_ORIENTED:
10632d7a744SHong Zhang     a->roworiented = flg;
10732d7a744SHong Zhang     break;
10832d7a744SHong Zhang   default:
10998921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
11032d7a744SHong Zhang   }
11132d7a744SHong Zhang   PetscFunctionReturn(0);
11232d7a744SHong Zhang }
11332d7a744SHong Zhang 
114e6dea9dbSXuan Zhou static PetscErrorCode MatSetValues_Elemental(Mat A,PetscInt nr,const PetscInt *rows,PetscInt nc,const PetscInt *cols,const PetscScalar *vals,InsertMode imode)
115db31f6deSJed Brown {
116db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
117fbbcb730SJack Poulson   PetscInt       i,j,rrank,ridx,crank,cidx,erow,ecol,numQueues=0;
118db31f6deSJed Brown 
119db31f6deSJed Brown   PetscFunctionBegin;
120fbbcb730SJack Poulson   // TODO: Initialize matrix to all zeros?
121fbbcb730SJack Poulson 
122fbbcb730SJack Poulson   // Count the number of queues from this process
12332d7a744SHong Zhang   if (a->roworiented) {
124db31f6deSJed Brown     for (i=0; i<nr; i++) {
125db31f6deSJed Brown       if (rows[i] < 0) continue;
126db31f6deSJed Brown       P2RO(A,0,rows[i],&rrank,&ridx);
127db31f6deSJed Brown       RO2E(A,0,rrank,ridx,&erow);
128aed4548fSBarry Smith       PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation");
129db31f6deSJed Brown       for (j=0; j<nc; j++) {
130db31f6deSJed Brown         if (cols[j] < 0) continue;
131db31f6deSJed Brown         P2RO(A,1,cols[j],&crank,&cidx);
132db31f6deSJed Brown         RO2E(A,1,crank,cidx,&ecol);
133aed4548fSBarry Smith         PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation");
134fbbcb730SJack Poulson         if (!a->emat->IsLocal(erow,ecol)) { /* off-proc entry */
135fbbcb730SJack Poulson           /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
13608401ef6SPierre Jolivet           PetscCheck(imode == ADD_VALUES,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported");
137fbbcb730SJack Poulson           ++numQueues;
138aae2c449SHong Zhang           continue;
139ed36708cSHong Zhang         }
140fbbcb730SJack Poulson         /* printf("Locally updating (%d,%d)\n",erow,ecol); */
141db31f6deSJed Brown         switch (imode) {
142fbbcb730SJack Poulson         case INSERT_VALUES: a->emat->Set(erow,ecol,(PetscElemScalar)vals[i*nc+j]); break;
143fbbcb730SJack Poulson         case ADD_VALUES: a->emat->Update(erow,ecol,(PetscElemScalar)vals[i*nc+j]); break;
14498921bdaSJacob Faibussowitsch         default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode);
145db31f6deSJed Brown         }
146db31f6deSJed Brown       }
147db31f6deSJed Brown     }
148fbbcb730SJack Poulson 
149fbbcb730SJack Poulson     /* printf("numQueues=%d\n",numQueues); */
150fbbcb730SJack Poulson     a->emat->Reserve( numQueues);
151fbbcb730SJack Poulson     for (i=0; i<nr; i++) {
152fbbcb730SJack Poulson       if (rows[i] < 0) continue;
153fbbcb730SJack Poulson       P2RO(A,0,rows[i],&rrank,&ridx);
154fbbcb730SJack Poulson       RO2E(A,0,rrank,ridx,&erow);
155fbbcb730SJack Poulson       for (j=0; j<nc; j++) {
156fbbcb730SJack Poulson         if (cols[j] < 0) continue;
157fbbcb730SJack Poulson         P2RO(A,1,cols[j],&crank,&cidx);
158fbbcb730SJack Poulson         RO2E(A,1,crank,cidx,&ecol);
159fbbcb730SJack Poulson         if (!a->emat->IsLocal(erow,ecol)) { /*off-proc entry*/
160fbbcb730SJack Poulson           /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
161fbbcb730SJack Poulson           a->emat->QueueUpdate( erow, ecol, vals[i*nc+j]);
162fbbcb730SJack Poulson         }
163fbbcb730SJack Poulson       }
164fbbcb730SJack Poulson     }
16532d7a744SHong Zhang   } else { /* columnoriented */
16632d7a744SHong Zhang     for (j=0; j<nc; j++) {
16732d7a744SHong Zhang       if (cols[j] < 0) continue;
16832d7a744SHong Zhang       P2RO(A,1,cols[j],&crank,&cidx);
16932d7a744SHong Zhang       RO2E(A,1,crank,cidx,&ecol);
170aed4548fSBarry Smith       PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation");
17132d7a744SHong Zhang       for (i=0; i<nr; i++) {
17232d7a744SHong Zhang         if (rows[i] < 0) continue;
17332d7a744SHong Zhang         P2RO(A,0,rows[i],&rrank,&ridx);
17432d7a744SHong Zhang         RO2E(A,0,rrank,ridx,&erow);
175aed4548fSBarry Smith         PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation");
17632d7a744SHong Zhang         if (!a->emat->IsLocal(erow,ecol)) { /* off-proc entry */
17732d7a744SHong Zhang           /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
17808401ef6SPierre Jolivet           PetscCheck(imode == ADD_VALUES,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported");
17932d7a744SHong Zhang           ++numQueues;
18032d7a744SHong Zhang           continue;
18132d7a744SHong Zhang         }
18232d7a744SHong Zhang         /* printf("Locally updating (%d,%d)\n",erow,ecol); */
18332d7a744SHong Zhang         switch (imode) {
18432d7a744SHong Zhang         case INSERT_VALUES: a->emat->Set(erow,ecol,(PetscElemScalar)vals[i+j*nr]); break;
18532d7a744SHong Zhang         case ADD_VALUES: a->emat->Update(erow,ecol,(PetscElemScalar)vals[i+j*nr]); break;
18698921bdaSJacob Faibussowitsch         default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode);
18732d7a744SHong Zhang         }
18832d7a744SHong Zhang       }
18932d7a744SHong Zhang     }
19032d7a744SHong Zhang 
19132d7a744SHong Zhang     /* printf("numQueues=%d\n",numQueues); */
19232d7a744SHong Zhang     a->emat->Reserve( numQueues);
19332d7a744SHong Zhang     for (j=0; j<nc; j++) {
19432d7a744SHong Zhang       if (cols[j] < 0) continue;
19532d7a744SHong Zhang       P2RO(A,1,cols[j],&crank,&cidx);
19632d7a744SHong Zhang       RO2E(A,1,crank,cidx,&ecol);
19732d7a744SHong Zhang 
19832d7a744SHong Zhang       for (i=0; i<nr; i++) {
19932d7a744SHong Zhang         if (rows[i] < 0) continue;
20032d7a744SHong Zhang         P2RO(A,0,rows[i],&rrank,&ridx);
20132d7a744SHong Zhang         RO2E(A,0,rrank,ridx,&erow);
20232d7a744SHong Zhang         if (!a->emat->IsLocal(erow,ecol)) { /*off-proc entry*/
20332d7a744SHong Zhang           /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
20432d7a744SHong Zhang           a->emat->QueueUpdate( erow, ecol, vals[i+j*nr]);
20532d7a744SHong Zhang         }
20632d7a744SHong Zhang       }
20732d7a744SHong Zhang     }
20832d7a744SHong Zhang   }
209db31f6deSJed Brown   PetscFunctionReturn(0);
210db31f6deSJed Brown }
211db31f6deSJed Brown 
212db31f6deSJed Brown static PetscErrorCode MatMult_Elemental(Mat A,Vec X,Vec Y)
213db31f6deSJed Brown {
214db31f6deSJed Brown   Mat_Elemental         *a = (Mat_Elemental*)A->data;
215e6dea9dbSXuan Zhou   const PetscElemScalar *x;
216e6dea9dbSXuan Zhou   PetscElemScalar       *y;
217df311e6cSXuan Zhou   PetscElemScalar       one = 1,zero = 0;
218db31f6deSJed Brown 
219db31f6deSJed Brown   PetscFunctionBegin;
2209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,(const PetscScalar **)&x));
2219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y,(PetscScalar **)&y));
222db31f6deSJed Brown   { /* Scoping so that constructor is called before pointer is returned */
223e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
2240c18141cSBarry Smith     xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
2250c18141cSBarry Smith     ye.Attach(A->rmap->N,1,*a->grid,0,0,y,A->rmap->n);
226e8146892SSatish Balay     El::Gemv(El::NORMAL,one,*a->emat,xe,zero,ye);
227db31f6deSJed Brown   }
2289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,(const PetscScalar **)&x));
2299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y,(PetscScalar **)&y));
230db31f6deSJed Brown   PetscFunctionReturn(0);
231db31f6deSJed Brown }
232db31f6deSJed Brown 
2339426833fSXuan Zhou static PetscErrorCode MatMultTranspose_Elemental(Mat A,Vec X,Vec Y)
2349426833fSXuan Zhou {
2359426833fSXuan Zhou   Mat_Elemental         *a = (Mat_Elemental*)A->data;
236df311e6cSXuan Zhou   const PetscElemScalar *x;
237df311e6cSXuan Zhou   PetscElemScalar       *y;
238e6dea9dbSXuan Zhou   PetscElemScalar       one = 1,zero = 0;
2399426833fSXuan Zhou 
2409426833fSXuan Zhou   PetscFunctionBegin;
2419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,(const PetscScalar **)&x));
2429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y,(PetscScalar **)&y));
2439426833fSXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
244e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
2450c18141cSBarry Smith     xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
2460c18141cSBarry Smith     ye.Attach(A->cmap->N,1,*a->grid,0,0,y,A->cmap->n);
247e8146892SSatish Balay     El::Gemv(El::TRANSPOSE,one,*a->emat,xe,zero,ye);
2489426833fSXuan Zhou   }
2499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,(const PetscScalar **)&x));
2509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y,(PetscScalar **)&y));
2519426833fSXuan Zhou   PetscFunctionReturn(0);
2529426833fSXuan Zhou }
2539426833fSXuan Zhou 
254db31f6deSJed Brown static PetscErrorCode MatMultAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
255db31f6deSJed Brown {
256db31f6deSJed Brown   Mat_Elemental         *a = (Mat_Elemental*)A->data;
257df311e6cSXuan Zhou   const PetscElemScalar *x;
258df311e6cSXuan Zhou   PetscElemScalar       *z;
259e6dea9dbSXuan Zhou   PetscElemScalar       one = 1;
260db31f6deSJed Brown 
261db31f6deSJed Brown   PetscFunctionBegin;
2629566063dSJacob Faibussowitsch   if (Y != Z) PetscCall(VecCopy(Y,Z));
2639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,(const PetscScalar **)&x));
2649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Z,(PetscScalar **)&z));
265db31f6deSJed Brown   { /* Scoping so that constructor is called before pointer is returned */
266e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
2670c18141cSBarry Smith     xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
2680c18141cSBarry Smith     ze.Attach(A->rmap->N,1,*a->grid,0,0,z,A->rmap->n);
269e8146892SSatish Balay     El::Gemv(El::NORMAL,one,*a->emat,xe,one,ze);
270db31f6deSJed Brown   }
2719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,(const PetscScalar **)&x));
2729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Z,(PetscScalar **)&z));
273db31f6deSJed Brown   PetscFunctionReturn(0);
274db31f6deSJed Brown }
275db31f6deSJed Brown 
276e883f9d5SXuan Zhou static PetscErrorCode MatMultTransposeAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
277e883f9d5SXuan Zhou {
278e883f9d5SXuan Zhou   Mat_Elemental         *a = (Mat_Elemental*)A->data;
279df311e6cSXuan Zhou   const PetscElemScalar *x;
280df311e6cSXuan Zhou   PetscElemScalar       *z;
281e6dea9dbSXuan Zhou   PetscElemScalar       one = 1;
282e883f9d5SXuan Zhou 
283e883f9d5SXuan Zhou   PetscFunctionBegin;
2849566063dSJacob Faibussowitsch   if (Y != Z) PetscCall(VecCopy(Y,Z));
2859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,(const PetscScalar **)&x));
2869566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Z,(PetscScalar **)&z));
287e883f9d5SXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
288e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
2890c18141cSBarry Smith     xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
2900c18141cSBarry Smith     ze.Attach(A->cmap->N,1,*a->grid,0,0,z,A->cmap->n);
291e8146892SSatish Balay     El::Gemv(El::TRANSPOSE,one,*a->emat,xe,one,ze);
292e883f9d5SXuan Zhou   }
2939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,(const PetscScalar **)&x));
2949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Z,(PetscScalar **)&z));
295e883f9d5SXuan Zhou   PetscFunctionReturn(0);
296e883f9d5SXuan Zhou }
297e883f9d5SXuan Zhou 
2984222ddf1SHong Zhang PetscErrorCode MatMatMultNumeric_Elemental(Mat A,Mat B,Mat C)
299c1d1b975SXuan Zhou {
300c1d1b975SXuan Zhou   Mat_Elemental    *a = (Mat_Elemental*)A->data;
301c1d1b975SXuan Zhou   Mat_Elemental    *b = (Mat_Elemental*)B->data;
3029a9e8502SHong Zhang   Mat_Elemental    *c = (Mat_Elemental*)C->data;
303e6dea9dbSXuan Zhou   PetscElemScalar  one = 1,zero = 0;
304c1d1b975SXuan Zhou 
305c1d1b975SXuan Zhou   PetscFunctionBegin;
306aae2c449SHong Zhang   { /* Scoping so that constructor is called before pointer is returned */
307e8146892SSatish Balay     El::Gemm(El::NORMAL,El::NORMAL,one,*a->emat,*b->emat,zero,*c->emat);
308aae2c449SHong Zhang   }
3099a9e8502SHong Zhang   C->assembled = PETSC_TRUE;
3109a9e8502SHong Zhang   PetscFunctionReturn(0);
3119a9e8502SHong Zhang }
3129a9e8502SHong Zhang 
3134222ddf1SHong Zhang PetscErrorCode MatMatMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat Ce)
3149a9e8502SHong Zhang {
3159a9e8502SHong Zhang   PetscFunctionBegin;
3169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Ce,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE));
3179566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ce,MATELEMENTAL));
3189566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Ce));
3194222ddf1SHong Zhang   Ce->ops->matmultnumeric = MatMatMultNumeric_Elemental;
320c1d1b975SXuan Zhou   PetscFunctionReturn(0);
321c1d1b975SXuan Zhou }
322c1d1b975SXuan Zhou 
323df311e6cSXuan Zhou static PetscErrorCode MatMatTransposeMultNumeric_Elemental(Mat A,Mat B,Mat C)
324df311e6cSXuan Zhou {
325df311e6cSXuan Zhou   Mat_Elemental      *a = (Mat_Elemental*)A->data;
326df311e6cSXuan Zhou   Mat_Elemental      *b = (Mat_Elemental*)B->data;
327df311e6cSXuan Zhou   Mat_Elemental      *c = (Mat_Elemental*)C->data;
328e6dea9dbSXuan Zhou   PetscElemScalar    one = 1,zero = 0;
329df311e6cSXuan Zhou 
330df311e6cSXuan Zhou   PetscFunctionBegin;
331df311e6cSXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
332e8146892SSatish Balay     El::Gemm(El::NORMAL,El::TRANSPOSE,one,*a->emat,*b->emat,zero,*c->emat);
333df311e6cSXuan Zhou   }
334df311e6cSXuan Zhou   C->assembled = PETSC_TRUE;
335df311e6cSXuan Zhou   PetscFunctionReturn(0);
336df311e6cSXuan Zhou }
337df311e6cSXuan Zhou 
3384222ddf1SHong Zhang static PetscErrorCode MatMatTransposeMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat C)
339df311e6cSXuan Zhou {
340df311e6cSXuan Zhou   PetscFunctionBegin;
3419566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,A->rmap->n,B->rmap->n,PETSC_DECIDE,PETSC_DECIDE));
3429566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATELEMENTAL));
3439566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
344df311e6cSXuan Zhou   PetscFunctionReturn(0);
345df311e6cSXuan Zhou }
346df311e6cSXuan Zhou 
3474222ddf1SHong Zhang /* --------------------------------------- */
3484222ddf1SHong Zhang static PetscErrorCode MatProductSetFromOptions_Elemental_AB(Mat C)
3494222ddf1SHong Zhang {
3504222ddf1SHong Zhang   PetscFunctionBegin;
3514222ddf1SHong Zhang   C->ops->matmultsymbolic = MatMatMultSymbolic_Elemental;
3524222ddf1SHong Zhang   C->ops->productsymbolic = MatProductSymbolic_AB;
3534222ddf1SHong Zhang   PetscFunctionReturn(0);
3544222ddf1SHong Zhang }
3554222ddf1SHong Zhang 
3564222ddf1SHong Zhang static PetscErrorCode MatProductSetFromOptions_Elemental_ABt(Mat C)
3574222ddf1SHong Zhang {
3584222ddf1SHong Zhang   PetscFunctionBegin;
3594222ddf1SHong Zhang   C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_Elemental;
3604222ddf1SHong Zhang   C->ops->productsymbolic          = MatProductSymbolic_ABt;
3614222ddf1SHong Zhang   PetscFunctionReturn(0);
3624222ddf1SHong Zhang }
3634222ddf1SHong Zhang 
3644222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Elemental(Mat C)
3654222ddf1SHong Zhang {
3664222ddf1SHong Zhang   Mat_Product    *product = C->product;
3674222ddf1SHong Zhang 
3684222ddf1SHong Zhang   PetscFunctionBegin;
3694222ddf1SHong Zhang   switch (product->type) {
3704222ddf1SHong Zhang   case MATPRODUCT_AB:
3719566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions_Elemental_AB(C));
3724222ddf1SHong Zhang     break;
3734222ddf1SHong Zhang   case MATPRODUCT_ABt:
3749566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions_Elemental_ABt(C));
3754222ddf1SHong Zhang     break;
3766718818eSStefano Zampini   default:
3776718818eSStefano Zampini     break;
3784222ddf1SHong Zhang   }
3794222ddf1SHong Zhang   PetscFunctionReturn(0);
3804222ddf1SHong Zhang }
381f6e5db1bSPierre Jolivet 
382f6e5db1bSPierre Jolivet PetscErrorCode MatMatMultNumeric_Elemental_MPIDense(Mat A,Mat B,Mat C)
383f6e5db1bSPierre Jolivet {
384f6e5db1bSPierre Jolivet   Mat            Be,Ce;
385f6e5db1bSPierre Jolivet 
386f6e5db1bSPierre Jolivet   PetscFunctionBegin;
3879566063dSJacob Faibussowitsch   PetscCall(MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX,&Be));
3889566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,Be,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Ce));
3899566063dSJacob Faibussowitsch   PetscCall(MatConvert(Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C));
3909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Be));
3919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ce));
392f6e5db1bSPierre Jolivet   PetscFunctionReturn(0);
393f6e5db1bSPierre Jolivet }
394f6e5db1bSPierre Jolivet 
395f6e5db1bSPierre Jolivet PetscErrorCode MatMatMultSymbolic_Elemental_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
396f6e5db1bSPierre Jolivet {
397f6e5db1bSPierre Jolivet   PetscFunctionBegin;
3989566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE));
3999566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATMPIDENSE));
4009566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
401f6e5db1bSPierre Jolivet   C->ops->matmultnumeric = MatMatMultNumeric_Elemental_MPIDense;
402f6e5db1bSPierre Jolivet   PetscFunctionReturn(0);
403f6e5db1bSPierre Jolivet }
404f6e5db1bSPierre Jolivet 
405f6e5db1bSPierre Jolivet PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense_AB(Mat C)
406f6e5db1bSPierre Jolivet {
407f6e5db1bSPierre Jolivet   PetscFunctionBegin;
408f6e5db1bSPierre Jolivet   C->ops->matmultsymbolic = MatMatMultSymbolic_Elemental_MPIDense;
409f6e5db1bSPierre Jolivet   C->ops->productsymbolic = MatProductSymbolic_AB;
410f6e5db1bSPierre Jolivet   PetscFunctionReturn(0);
411f6e5db1bSPierre Jolivet }
412f6e5db1bSPierre Jolivet 
413f6e5db1bSPierre Jolivet PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense(Mat C)
414f6e5db1bSPierre Jolivet {
415f6e5db1bSPierre Jolivet   Mat_Product    *product = C->product;
416f6e5db1bSPierre Jolivet 
417f6e5db1bSPierre Jolivet   PetscFunctionBegin;
418f6e5db1bSPierre Jolivet   if (product->type == MATPRODUCT_AB) {
4199566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions_Elemental_MPIDense_AB(C));
420f6e5db1bSPierre Jolivet   }
421f6e5db1bSPierre Jolivet   PetscFunctionReturn(0);
422f6e5db1bSPierre Jolivet }
4234222ddf1SHong Zhang /* --------------------------------------- */
4244222ddf1SHong Zhang 
425a9d89745SXuan Zhou static PetscErrorCode MatGetDiagonal_Elemental(Mat A,Vec D)
42661119200SXuan Zhou {
427a9d89745SXuan Zhou   PetscInt        i,nrows,ncols,nD,rrank,ridx,crank,cidx;
428a9d89745SXuan Zhou   Mat_Elemental   *a = (Mat_Elemental*)A->data;
429a9d89745SXuan Zhou   PetscElemScalar v;
430ce94432eSBarry Smith   MPI_Comm        comm;
43161119200SXuan Zhou 
43261119200SXuan Zhou   PetscFunctionBegin;
4339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
4349566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&nrows,&ncols));
435a9d89745SXuan Zhou   nD = nrows>ncols ? ncols : nrows;
436a9d89745SXuan Zhou   for (i=0; i<nD; i++) {
437a9d89745SXuan Zhou     PetscInt erow,ecol;
438a9d89745SXuan Zhou     P2RO(A,0,i,&rrank,&ridx);
439a9d89745SXuan Zhou     RO2E(A,0,rrank,ridx,&erow);
440aed4548fSBarry Smith     PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0,comm,PETSC_ERR_PLIB,"Incorrect row translation");
441a9d89745SXuan Zhou     P2RO(A,1,i,&crank,&cidx);
442a9d89745SXuan Zhou     RO2E(A,1,crank,cidx,&ecol);
443aed4548fSBarry Smith     PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0,comm,PETSC_ERR_PLIB,"Incorrect col translation");
444a9d89745SXuan Zhou     v = a->emat->Get(erow,ecol);
4459566063dSJacob Faibussowitsch     PetscCall(VecSetValues(D,1,&i,(PetscScalar*)&v,INSERT_VALUES));
446ade3cc5eSXuan Zhou   }
4479566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(D));
4489566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(D));
44961119200SXuan Zhou   PetscFunctionReturn(0);
45061119200SXuan Zhou }
45161119200SXuan Zhou 
452ade3cc5eSXuan Zhou static PetscErrorCode MatDiagonalScale_Elemental(Mat X,Vec L,Vec R)
453ade3cc5eSXuan Zhou {
454ade3cc5eSXuan Zhou   Mat_Elemental         *x = (Mat_Elemental*)X->data;
455ade3cc5eSXuan Zhou   const PetscElemScalar *d;
456ade3cc5eSXuan Zhou 
457ade3cc5eSXuan Zhou   PetscFunctionBegin;
4589065cd98SJed Brown   if (R) {
4599566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(R,(const PetscScalar **)&d));
460e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
4610c18141cSBarry Smith     de.LockedAttach(X->cmap->N,1,*x->grid,0,0,d,X->cmap->n);
462e8146892SSatish Balay     El::DiagonalScale(El::RIGHT,El::NORMAL,de,*x->emat);
4639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(R,(const PetscScalar **)&d));
4649065cd98SJed Brown   }
4659065cd98SJed Brown   if (L) {
4669566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(L,(const PetscScalar **)&d));
467e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
4680c18141cSBarry Smith     de.LockedAttach(X->rmap->N,1,*x->grid,0,0,d,X->rmap->n);
469e8146892SSatish Balay     El::DiagonalScale(El::LEFT,El::NORMAL,de,*x->emat);
4709566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(L,(const PetscScalar **)&d));
471ade3cc5eSXuan Zhou   }
472ade3cc5eSXuan Zhou   PetscFunctionReturn(0);
473ade3cc5eSXuan Zhou }
474ade3cc5eSXuan Zhou 
47534fa46e1SFrancesco Ballarin static PetscErrorCode MatMissingDiagonal_Elemental(Mat A,PetscBool *missing,PetscInt *d)
47634fa46e1SFrancesco Ballarin {
47734fa46e1SFrancesco Ballarin   PetscFunctionBegin;
47834fa46e1SFrancesco Ballarin   *missing = PETSC_FALSE;
47934fa46e1SFrancesco Ballarin   PetscFunctionReturn(0);
48034fa46e1SFrancesco Ballarin }
48134fa46e1SFrancesco Ballarin 
482e6dea9dbSXuan Zhou static PetscErrorCode MatScale_Elemental(Mat X,PetscScalar a)
48365b78793SXuan Zhou {
48465b78793SXuan Zhou   Mat_Elemental  *x = (Mat_Elemental*)X->data;
48565b78793SXuan Zhou 
48665b78793SXuan Zhou   PetscFunctionBegin;
487e8146892SSatish Balay   El::Scale((PetscElemScalar)a,*x->emat);
48865b78793SXuan Zhou   PetscFunctionReturn(0);
48965b78793SXuan Zhou }
49065b78793SXuan Zhou 
491f82baa17SHong Zhang /*
492f82baa17SHong Zhang   MatAXPY - Computes Y = a*X + Y.
493f82baa17SHong Zhang */
494e6dea9dbSXuan Zhou static PetscErrorCode MatAXPY_Elemental(Mat Y,PetscScalar a,Mat X,MatStructure str)
495e09a3074SHong Zhang {
496e09a3074SHong Zhang   Mat_Elemental  *x = (Mat_Elemental*)X->data;
497e09a3074SHong Zhang   Mat_Elemental  *y = (Mat_Elemental*)Y->data;
498e09a3074SHong Zhang 
499e09a3074SHong Zhang   PetscFunctionBegin;
500e8146892SSatish Balay   El::Axpy((PetscElemScalar)a,*x->emat,*y->emat);
5019566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
502e09a3074SHong Zhang   PetscFunctionReturn(0);
503e09a3074SHong Zhang }
504e09a3074SHong Zhang 
505d6223691SXuan Zhou static PetscErrorCode MatCopy_Elemental(Mat A,Mat B,MatStructure str)
506d6223691SXuan Zhou {
507d6223691SXuan Zhou   Mat_Elemental *a=(Mat_Elemental*)A->data;
508d6223691SXuan Zhou   Mat_Elemental *b=(Mat_Elemental*)B->data;
509d6223691SXuan Zhou 
510d6223691SXuan Zhou   PetscFunctionBegin;
511e8146892SSatish Balay   El::Copy(*a->emat,*b->emat);
5129566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)B));
513d6223691SXuan Zhou   PetscFunctionReturn(0);
514d6223691SXuan Zhou }
515d6223691SXuan Zhou 
516df311e6cSXuan Zhou static PetscErrorCode MatDuplicate_Elemental(Mat A,MatDuplicateOption op,Mat *B)
517df311e6cSXuan Zhou {
518df311e6cSXuan Zhou   Mat            Be;
519ce94432eSBarry Smith   MPI_Comm       comm;
520df311e6cSXuan Zhou   Mat_Elemental  *a=(Mat_Elemental*)A->data;
521df311e6cSXuan Zhou 
522df311e6cSXuan Zhou   PetscFunctionBegin;
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
5249566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&Be));
5259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Be,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE));
5269566063dSJacob Faibussowitsch   PetscCall(MatSetType(Be,MATELEMENTAL));
5279566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Be));
528df311e6cSXuan Zhou   *B = Be;
529df311e6cSXuan Zhou   if (op == MAT_COPY_VALUES) {
530df311e6cSXuan Zhou     Mat_Elemental *b=(Mat_Elemental*)Be->data;
531e8146892SSatish Balay     El::Copy(*a->emat,*b->emat);
532df311e6cSXuan Zhou   }
533df311e6cSXuan Zhou   Be->assembled = PETSC_TRUE;
534df311e6cSXuan Zhou   PetscFunctionReturn(0);
535df311e6cSXuan Zhou }
536df311e6cSXuan Zhou 
537d6223691SXuan Zhou static PetscErrorCode MatTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
538d6223691SXuan Zhou {
539ec8cb81fSBarry Smith   Mat            Be = *B;
540ce94432eSBarry Smith   MPI_Comm       comm;
5414fe7bbcaSHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
542d6223691SXuan Zhou 
543d6223691SXuan Zhou   PetscFunctionBegin;
544*7fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A,*B));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
5465cb544a0SHong Zhang   /* Only out-of-place supported */
5475f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX,comm,PETSC_ERR_SUP,"Only out-of-place supported");
5485262d616SXuan Zhou   if (reuse == MAT_INITIAL_MATRIX) {
5499566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm,&Be));
5509566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE));
5519566063dSJacob Faibussowitsch     PetscCall(MatSetType(Be,MATELEMENTAL));
5529566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Be));
5535262d616SXuan Zhou     *B = Be;
5545262d616SXuan Zhou   }
5554fe7bbcaSHong Zhang   b = (Mat_Elemental*)Be->data;
556e8146892SSatish Balay   El::Transpose(*a->emat,*b->emat);
5575262d616SXuan Zhou   Be->assembled = PETSC_TRUE;
558d6223691SXuan Zhou   PetscFunctionReturn(0);
559d6223691SXuan Zhou }
560d6223691SXuan Zhou 
561dfcb0403SXuan Zhou static PetscErrorCode MatConjugate_Elemental(Mat A)
562dfcb0403SXuan Zhou {
563dfcb0403SXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data;
564dfcb0403SXuan Zhou 
565dfcb0403SXuan Zhou   PetscFunctionBegin;
566e8146892SSatish Balay   El::Conjugate(*a->emat);
567dfcb0403SXuan Zhou   PetscFunctionReturn(0);
568dfcb0403SXuan Zhou }
569dfcb0403SXuan Zhou 
5704a29722dSXuan Zhou static PetscErrorCode MatHermitianTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
5714a29722dSXuan Zhou {
572ec8cb81fSBarry Smith   Mat            Be = *B;
573ce94432eSBarry Smith   MPI_Comm       comm;
5744a29722dSXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
5754a29722dSXuan Zhou 
5764a29722dSXuan Zhou   PetscFunctionBegin;
5779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
5785cb544a0SHong Zhang   /* Only out-of-place supported */
5794a29722dSXuan Zhou   if (reuse == MAT_INITIAL_MATRIX) {
5809566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm,&Be));
5819566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE));
5829566063dSJacob Faibussowitsch     PetscCall(MatSetType(Be,MATELEMENTAL));
5839566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Be));
5844a29722dSXuan Zhou     *B = Be;
5854a29722dSXuan Zhou   }
5864a29722dSXuan Zhou   b = (Mat_Elemental*)Be->data;
587e8146892SSatish Balay   El::Adjoint(*a->emat,*b->emat);
5884a29722dSXuan Zhou   Be->assembled = PETSC_TRUE;
5894a29722dSXuan Zhou   PetscFunctionReturn(0);
5904a29722dSXuan Zhou }
5914a29722dSXuan Zhou 
5921f881ff8SXuan Zhou static PetscErrorCode MatSolve_Elemental(Mat A,Vec B,Vec X)
5931f881ff8SXuan Zhou {
5941f881ff8SXuan Zhou   Mat_Elemental     *a = (Mat_Elemental*)A->data;
595df311e6cSXuan Zhou   PetscElemScalar   *x;
596ea394475SSatish Balay   PetscInt          pivoting = a->pivoting;
5971f881ff8SXuan Zhou 
5981f881ff8SXuan Zhou   PetscFunctionBegin;
5999566063dSJacob Faibussowitsch   PetscCall(VecCopy(B,X));
6009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,(PetscScalar **)&x));
601ea394475SSatish Balay 
602e8146892SSatish Balay   El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe;
6030c18141cSBarry Smith   xe.Attach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
604e8146892SSatish Balay   El::DistMatrix<PetscElemScalar,El::MC,El::MR> xer(xe);
605fc54b460SXuan Zhou   switch (A->factortype) {
606fc54b460SXuan Zhou   case MAT_FACTOR_LU:
607ea394475SSatish Balay     if (pivoting == 0) {
608e8146892SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,xer);
609ea394475SSatish Balay     } else if (pivoting == 1) {
610ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,xer);
611ea394475SSatish Balay     } else { /* pivoting == 2 */
612ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,xer);
6131f881ff8SXuan Zhou     }
614fc54b460SXuan Zhou     break;
615fc54b460SXuan Zhou   case MAT_FACTOR_CHOLESKY:
616e8146892SSatish Balay     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,xer);
617fc54b460SXuan Zhou     break;
618fc54b460SXuan Zhou   default:
6194fe7bbcaSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
620fc54b460SXuan Zhou     break;
6211f881ff8SXuan Zhou   }
622ea394475SSatish Balay   El::Copy(xer,xe);
623ea394475SSatish Balay 
6249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,(PetscScalar **)&x));
6251f881ff8SXuan Zhou   PetscFunctionReturn(0);
6261f881ff8SXuan Zhou }
6271f881ff8SXuan Zhou 
628df311e6cSXuan Zhou static PetscErrorCode MatSolveAdd_Elemental(Mat A,Vec B,Vec Y,Vec X)
629df311e6cSXuan Zhou {
630df311e6cSXuan Zhou   PetscFunctionBegin;
6319566063dSJacob Faibussowitsch   PetscCall(MatSolve_Elemental(A,B,X));
6329566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,1,Y));
633df311e6cSXuan Zhou   PetscFunctionReturn(0);
634df311e6cSXuan Zhou }
635df311e6cSXuan Zhou 
636ae844d54SHong Zhang static PetscErrorCode MatMatSolve_Elemental(Mat A,Mat B,Mat X)
637ae844d54SHong Zhang {
6381f0e42cfSHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data;
63942ba530cSPierre Jolivet   Mat_Elemental  *x;
64042ba530cSPierre Jolivet   Mat            C;
641ea394475SSatish Balay   PetscInt       pivoting = a->pivoting;
64242ba530cSPierre Jolivet   PetscBool      flg;
64342ba530cSPierre Jolivet   MatType        type;
6441f0e42cfSHong Zhang 
645ae844d54SHong Zhang   PetscFunctionBegin;
6469566063dSJacob Faibussowitsch   PetscCall(MatGetType(X,&type));
6479566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type,MATELEMENTAL,&flg));
64842ba530cSPierre Jolivet   if (!flg) {
6499566063dSJacob Faibussowitsch     PetscCall(MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX,&C));
65042ba530cSPierre Jolivet     x = (Mat_Elemental*)C->data;
65142ba530cSPierre Jolivet   } else {
65242ba530cSPierre Jolivet     x = (Mat_Elemental*)X->data;
65342ba530cSPierre Jolivet     El::Copy(*((Mat_Elemental*)B->data)->emat,*x->emat);
65442ba530cSPierre Jolivet   }
655fc54b460SXuan Zhou   switch (A->factortype) {
656fc54b460SXuan Zhou   case MAT_FACTOR_LU:
657ea394475SSatish Balay     if (pivoting == 0) {
658e8146892SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*x->emat);
659ea394475SSatish Balay     } else if (pivoting == 1) {
660ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*x->emat);
661ea394475SSatish Balay     } else {
662ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,*x->emat);
663d6223691SXuan Zhou     }
664fc54b460SXuan Zhou     break;
665fc54b460SXuan Zhou   case MAT_FACTOR_CHOLESKY:
666e8146892SSatish Balay     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,*x->emat);
667fc54b460SXuan Zhou     break;
668fc54b460SXuan Zhou   default:
6694fe7bbcaSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
670fc54b460SXuan Zhou     break;
671fc54b460SXuan Zhou   }
67242ba530cSPierre Jolivet   if (!flg) {
6739566063dSJacob Faibussowitsch     PetscCall(MatConvert(C,type,MAT_REUSE_MATRIX,&X));
6749566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
67542ba530cSPierre Jolivet   }
676ae844d54SHong Zhang   PetscFunctionReturn(0);
677ae844d54SHong Zhang }
678ae844d54SHong Zhang 
679ae844d54SHong Zhang static PetscErrorCode MatLUFactor_Elemental(Mat A,IS row,IS col,const MatFactorInfo *info)
680ae844d54SHong Zhang {
6817c920d81SXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data;
682ea394475SSatish Balay   PetscInt       pivoting = a->pivoting;
6837c920d81SXuan Zhou 
684ae844d54SHong Zhang   PetscFunctionBegin;
685ea394475SSatish Balay   if (pivoting == 0) {
686e8146892SSatish Balay     El::LU(*a->emat);
687ea394475SSatish Balay   } else if (pivoting == 1) {
688ea394475SSatish Balay     El::LU(*a->emat,*a->P);
689ea394475SSatish Balay   } else {
690ea394475SSatish Balay     El::LU(*a->emat,*a->P,*a->Q);
691d6223691SXuan Zhou   }
6921293973dSHong Zhang   A->factortype = MAT_FACTOR_LU;
693834d3fecSHong Zhang   A->assembled  = PETSC_TRUE;
694f6224b95SHong Zhang 
6959566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->solvertype));
6969566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype));
697ae844d54SHong Zhang   PetscFunctionReturn(0);
698ae844d54SHong Zhang }
699ae844d54SHong Zhang 
700d7c3f9d8SHong Zhang static PetscErrorCode  MatLUFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
701d7c3f9d8SHong Zhang {
702d7c3f9d8SHong Zhang   PetscFunctionBegin;
7039566063dSJacob Faibussowitsch   PetscCall(MatCopy(A,F,SAME_NONZERO_PATTERN));
7049566063dSJacob Faibussowitsch   PetscCall(MatLUFactor_Elemental(F,0,0,info));
705d7c3f9d8SHong Zhang   PetscFunctionReturn(0);
706d7c3f9d8SHong Zhang }
707d7c3f9d8SHong Zhang 
708d7c3f9d8SHong Zhang static PetscErrorCode  MatLUFactorSymbolic_Elemental(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
709d7c3f9d8SHong Zhang {
710d7c3f9d8SHong Zhang   PetscFunctionBegin;
711d8304050SJose E. Roman   /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
712d7c3f9d8SHong Zhang   PetscFunctionReturn(0);
713d7c3f9d8SHong Zhang }
714d7c3f9d8SHong Zhang 
71545cf121fSXuan Zhou static PetscErrorCode MatCholeskyFactor_Elemental(Mat A,IS perm,const MatFactorInfo *info)
71645cf121fSXuan Zhou {
71745cf121fSXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data;
718e8146892SSatish Balay   El::DistMatrix<PetscElemScalar,El::MC,El::STAR> d;
71945cf121fSXuan Zhou 
72045cf121fSXuan Zhou   PetscFunctionBegin;
721e8146892SSatish Balay   El::Cholesky(El::UPPER,*a->emat);
722fc54b460SXuan Zhou   A->factortype = MAT_FACTOR_CHOLESKY;
723834d3fecSHong Zhang   A->assembled  = PETSC_TRUE;
724f6224b95SHong Zhang 
7259566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->solvertype));
7269566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype));
72745cf121fSXuan Zhou   PetscFunctionReturn(0);
72845cf121fSXuan Zhou }
72945cf121fSXuan Zhou 
73079673f7bSHong Zhang static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
73179673f7bSHong Zhang {
732cb76c1d8SXuan Zhou   PetscFunctionBegin;
7339566063dSJacob Faibussowitsch   PetscCall(MatCopy(A,F,SAME_NONZERO_PATTERN));
7349566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactor_Elemental(F,0,info));
735cb76c1d8SXuan Zhou   PetscFunctionReturn(0);
73679673f7bSHong Zhang }
737cb76c1d8SXuan Zhou 
73879673f7bSHong Zhang static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat F,Mat A,IS perm,const MatFactorInfo *info)
73979673f7bSHong Zhang {
74079673f7bSHong Zhang   PetscFunctionBegin;
741d8304050SJose E. Roman   /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
74279673f7bSHong Zhang   PetscFunctionReturn(0);
74379673f7bSHong Zhang }
74479673f7bSHong Zhang 
745ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_elemental_elemental(Mat A,MatSolverType *type)
74615767789SHong Zhang {
74715767789SHong Zhang   PetscFunctionBegin;
74815767789SHong Zhang   *type = MATSOLVERELEMENTAL;
74915767789SHong Zhang   PetscFunctionReturn(0);
75015767789SHong Zhang }
75115767789SHong Zhang 
75215767789SHong Zhang static PetscErrorCode MatGetFactor_elemental_elemental(Mat A,MatFactorType ftype,Mat *F)
7531293973dSHong Zhang {
7541293973dSHong Zhang   Mat            B;
7551293973dSHong Zhang 
7561293973dSHong Zhang   PetscFunctionBegin;
7571293973dSHong Zhang   /* Create the factorization matrix */
7589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
7599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE));
7609566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,MATELEMENTAL));
7619566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
7621293973dSHong Zhang   B->factortype = ftype;
76366e17bc3SBarry Smith   B->trivialsymbolic = PETSC_TRUE;
7649566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
7659566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL,&B->solvertype));
76600c67f3bSHong Zhang 
7679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_elemental_elemental));
7681293973dSHong Zhang   *F            = B;
7691293973dSHong Zhang   PetscFunctionReturn(0);
7701293973dSHong Zhang }
7711293973dSHong Zhang 
7723ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_Elemental(void)
77342c9c57cSBarry Smith {
77442c9c57cSBarry Smith   PetscFunctionBegin;
7759566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_LU,MatGetFactor_elemental_elemental));
7769566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_CHOLESKY,MatGetFactor_elemental_elemental));
77742c9c57cSBarry Smith   PetscFunctionReturn(0);
77842c9c57cSBarry Smith }
77942c9c57cSBarry Smith 
7801f881ff8SXuan Zhou static PetscErrorCode MatNorm_Elemental(Mat A,NormType type,PetscReal *nrm)
7811f881ff8SXuan Zhou {
7821f881ff8SXuan Zhou   Mat_Elemental *a=(Mat_Elemental*)A->data;
7831f881ff8SXuan Zhou 
7841f881ff8SXuan Zhou   PetscFunctionBegin;
7851f881ff8SXuan Zhou   switch (type) {
7861f881ff8SXuan Zhou   case NORM_1:
787e8146892SSatish Balay     *nrm = El::OneNorm(*a->emat);
7881f881ff8SXuan Zhou     break;
7891f881ff8SXuan Zhou   case NORM_FROBENIUS:
790e8146892SSatish Balay     *nrm = El::FrobeniusNorm(*a->emat);
7911f881ff8SXuan Zhou     break;
7921f881ff8SXuan Zhou   case NORM_INFINITY:
793e8146892SSatish Balay     *nrm = El::InfinityNorm(*a->emat);
7941f881ff8SXuan Zhou     break;
7951f881ff8SXuan Zhou   default:
796d8304050SJose E. Roman     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported norm type");
7971f881ff8SXuan Zhou   }
7981f881ff8SXuan Zhou   PetscFunctionReturn(0);
7991f881ff8SXuan Zhou }
8001f881ff8SXuan Zhou 
8015262d616SXuan Zhou static PetscErrorCode MatZeroEntries_Elemental(Mat A)
8025262d616SXuan Zhou {
8035262d616SXuan Zhou   Mat_Elemental *a=(Mat_Elemental*)A->data;
8045262d616SXuan Zhou 
8055262d616SXuan Zhou   PetscFunctionBegin;
806e8146892SSatish Balay   El::Zero(*a->emat);
8075262d616SXuan Zhou   PetscFunctionReturn(0);
8085262d616SXuan Zhou }
8095262d616SXuan Zhou 
810db31f6deSJed Brown static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A,IS *rows,IS *cols)
811db31f6deSJed Brown {
812db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
813db31f6deSJed Brown   PetscInt       i,m,shift,stride,*idx;
814db31f6deSJed Brown 
815db31f6deSJed Brown   PetscFunctionBegin;
816db31f6deSJed Brown   if (rows) {
817db31f6deSJed Brown     m = a->emat->LocalHeight();
818db31f6deSJed Brown     shift = a->emat->ColShift();
819db31f6deSJed Brown     stride = a->emat->ColStride();
8209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m,&idx));
821db31f6deSJed Brown     for (i=0; i<m; i++) {
822db31f6deSJed Brown       PetscInt rank,offset;
823db31f6deSJed Brown       E2RO(A,0,shift+i*stride,&rank,&offset);
824db31f6deSJed Brown       RO2P(A,0,rank,offset,&idx[i]);
825db31f6deSJed Brown     }
8269566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,rows));
827db31f6deSJed Brown   }
828db31f6deSJed Brown   if (cols) {
829db31f6deSJed Brown     m = a->emat->LocalWidth();
830db31f6deSJed Brown     shift = a->emat->RowShift();
831db31f6deSJed Brown     stride = a->emat->RowStride();
8329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m,&idx));
833db31f6deSJed Brown     for (i=0; i<m; i++) {
834db31f6deSJed Brown       PetscInt rank,offset;
835db31f6deSJed Brown       E2RO(A,1,shift+i*stride,&rank,&offset);
836db31f6deSJed Brown       RO2P(A,1,rank,offset,&idx[i]);
837db31f6deSJed Brown     }
8389566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,cols));
839db31f6deSJed Brown   }
840db31f6deSJed Brown   PetscFunctionReturn(0);
841db31f6deSJed Brown }
842db31f6deSJed Brown 
84319fd82e9SBarry Smith static PetscErrorCode MatConvert_Elemental_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
844af295397SXuan Zhou {
8452ef0cf24SXuan Zhou   Mat                Bmpi;
846af295397SXuan Zhou   Mat_Elemental      *a = (Mat_Elemental*)A->data;
847ce94432eSBarry Smith   MPI_Comm           comm;
848fa9e26c8SHong Zhang   IS                 isrows,iscols;
849fa9e26c8SHong Zhang   PetscInt           rrank,ridx,crank,cidx,nrows,ncols,i,j,erow,ecol,elrow,elcol;
850fa9e26c8SHong Zhang   const PetscInt     *rows,*cols;
851df311e6cSXuan Zhou   PetscElemScalar    v;
852fa9e26c8SHong Zhang   const El::Grid     &grid = a->emat->Grid();
853af295397SXuan Zhou 
854af295397SXuan Zhou   PetscFunctionBegin;
8559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
856a000e53bSHong Zhang 
857de0a22f0SHong Zhang   if (reuse == MAT_REUSE_MATRIX) {
858de0a22f0SHong Zhang     Bmpi = *B;
859de0a22f0SHong Zhang   } else {
8609566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm,&Bmpi));
8619566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Bmpi,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE));
8629566063dSJacob Faibussowitsch     PetscCall(MatSetType(Bmpi,MATDENSE));
8639566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Bmpi));
864de0a22f0SHong Zhang   }
865fa9e26c8SHong Zhang 
866fa9e26c8SHong Zhang   /* Get local entries of A */
8679566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(A,&isrows,&iscols));
8689566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows,&nrows));
8699566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows,&rows));
8709566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscols,&ncols));
8719566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscols,&cols));
872fa9e26c8SHong Zhang 
87357299f2fSHong Zhang   if (a->roworiented) {
8742ef0cf24SXuan Zhou     for (i=0; i<nrows; i++) {
875fa9e26c8SHong Zhang       P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
8762ef0cf24SXuan Zhou       RO2E(A,0,rrank,ridx,&erow);
877aed4548fSBarry Smith       PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0,comm,PETSC_ERR_PLIB,"Incorrect row translation");
8782ef0cf24SXuan Zhou       for (j=0; j<ncols; j++) {
879fa9e26c8SHong Zhang         P2RO(A,1,cols[j],&crank,&cidx);
8802ef0cf24SXuan Zhou         RO2E(A,1,crank,cidx,&ecol);
881aed4548fSBarry Smith         PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0,comm,PETSC_ERR_PLIB,"Incorrect col translation");
882fa9e26c8SHong Zhang 
883fa9e26c8SHong Zhang         elrow = erow / grid.MCSize(); /* Elemental local row index */
884fa9e26c8SHong Zhang         elcol = ecol / grid.MRSize(); /* Elemental local column index */
885fa9e26c8SHong Zhang         v = a->emat->GetLocal(elrow,elcol);
8869566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES));
8872ef0cf24SXuan Zhou       }
8882ef0cf24SXuan Zhou     }
88957299f2fSHong Zhang   } else { /* column-oriented */
89057299f2fSHong Zhang     for (j=0; j<ncols; j++) {
89157299f2fSHong Zhang       P2RO(A,1,cols[j],&crank,&cidx);
89257299f2fSHong Zhang       RO2E(A,1,crank,cidx,&ecol);
893aed4548fSBarry Smith       PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0,comm,PETSC_ERR_PLIB,"Incorrect col translation");
89457299f2fSHong Zhang       for (i=0; i<nrows; i++) {
89557299f2fSHong Zhang         P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
89657299f2fSHong Zhang         RO2E(A,0,rrank,ridx,&erow);
897aed4548fSBarry Smith         PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0,comm,PETSC_ERR_PLIB,"Incorrect row translation");
89857299f2fSHong Zhang 
89957299f2fSHong Zhang         elrow = erow / grid.MCSize(); /* Elemental local row index */
90057299f2fSHong Zhang         elcol = ecol / grid.MRSize(); /* Elemental local column index */
90157299f2fSHong Zhang         v = a->emat->GetLocal(elrow,elcol);
9029566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES));
90357299f2fSHong Zhang       }
90457299f2fSHong Zhang     }
90557299f2fSHong Zhang   }
9069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Bmpi,MAT_FINAL_ASSEMBLY));
9079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Bmpi,MAT_FINAL_ASSEMBLY));
908511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
9099566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&Bmpi));
910c4ad791aSXuan Zhou   } else {
911c4ad791aSXuan Zhou     *B = Bmpi;
912c4ad791aSXuan Zhou   }
9139566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
9149566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscols));
915af295397SXuan Zhou   PetscFunctionReturn(0);
916af295397SXuan Zhou }
917af295397SXuan Zhou 
918cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
919af8000cdSHong Zhang {
920af8000cdSHong Zhang   Mat               mat_elemental;
921af8000cdSHong Zhang   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols;
922af8000cdSHong Zhang   const PetscInt    *cols;
923af8000cdSHong Zhang   const PetscScalar *vals;
924af8000cdSHong Zhang 
925af8000cdSHong Zhang   PetscFunctionBegin;
926bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
927bdeae278SStefano Zampini     mat_elemental = *newmat;
9289566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
929bdeae278SStefano Zampini   } else {
9309566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
9319566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N));
9329566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental,MATELEMENTAL));
9339566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
934bdeae278SStefano Zampini   }
935af8000cdSHong Zhang   for (row=0; row<M; row++) {
9369566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A,row,&ncols,&cols,&vals));
937bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
9389566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES));
9399566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A,row,&ncols,&cols,&vals));
940af8000cdSHong Zhang   }
9419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
9429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
943af8000cdSHong Zhang 
944511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
9459566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&mat_elemental));
946af8000cdSHong Zhang   } else {
947af8000cdSHong Zhang     *newmat = mat_elemental;
948af8000cdSHong Zhang   }
949af8000cdSHong Zhang   PetscFunctionReturn(0);
950af8000cdSHong Zhang }
951af8000cdSHong Zhang 
952cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
9535d7652ecSHong Zhang {
9545d7652ecSHong Zhang   Mat               mat_elemental;
9555d7652ecSHong Zhang   PetscInt          row,ncols,rstart=A->rmap->rstart,rend=A->rmap->rend,j;
9565d7652ecSHong Zhang   const PetscInt    *cols;
9575d7652ecSHong Zhang   const PetscScalar *vals;
9585d7652ecSHong Zhang 
9595d7652ecSHong Zhang   PetscFunctionBegin;
960bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
961bdeae278SStefano Zampini     mat_elemental = *newmat;
9629566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
963bdeae278SStefano Zampini   } else {
9649566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
9659566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N));
9669566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental,MATELEMENTAL));
9679566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
968bdeae278SStefano Zampini   }
9695d7652ecSHong Zhang   for (row=rstart; row<rend; row++) {
9709566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A,row,&ncols,&cols,&vals));
9715d7652ecSHong Zhang     for (j=0; j<ncols; j++) {
972bdeae278SStefano Zampini       /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
9739566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat_elemental,1,&row,1,&cols[j],&vals[j],ADD_VALUES));
9745d7652ecSHong Zhang     }
9759566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A,row,&ncols,&cols,&vals));
9765d7652ecSHong Zhang   }
9779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
9789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
9795d7652ecSHong Zhang 
980511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
9819566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&mat_elemental));
9825d7652ecSHong Zhang   } else {
9835d7652ecSHong Zhang     *newmat = mat_elemental;
9845d7652ecSHong Zhang   }
9855d7652ecSHong Zhang   PetscFunctionReturn(0);
9865d7652ecSHong Zhang }
9875d7652ecSHong Zhang 
988cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
9896214f412SHong Zhang {
9906214f412SHong Zhang   Mat               mat_elemental;
9916214f412SHong Zhang   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols,j;
9926214f412SHong Zhang   const PetscInt    *cols;
9936214f412SHong Zhang   const PetscScalar *vals;
9946214f412SHong Zhang 
9956214f412SHong Zhang   PetscFunctionBegin;
996bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
997bdeae278SStefano Zampini     mat_elemental = *newmat;
9989566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
999bdeae278SStefano Zampini   } else {
10009566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
10019566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N));
10029566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental,MATELEMENTAL));
10039566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
1004bdeae278SStefano Zampini   }
10059566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(A));
10066214f412SHong Zhang   for (row=0; row<M; row++) {
10079566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A,row,&ncols,&cols,&vals));
1008bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10099566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES));
10106214f412SHong Zhang     for (j=0; j<ncols; j++) { /* lower triangular part */
1011eb1ec7c1SStefano Zampini       PetscScalar v;
10126214f412SHong Zhang       if (cols[j] == row) continue;
1013b94d7dedSBarry Smith       v    = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j];
10149566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat_elemental,1,&cols[j],1,&row,&v,ADD_VALUES));
10156214f412SHong Zhang     }
10169566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A,row,&ncols,&cols,&vals));
10176214f412SHong Zhang   }
10189566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(A));
10199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
10209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
10216214f412SHong Zhang 
1022511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
10239566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&mat_elemental));
10246214f412SHong Zhang   } else {
10256214f412SHong Zhang     *newmat = mat_elemental;
10266214f412SHong Zhang   }
10276214f412SHong Zhang   PetscFunctionReturn(0);
10286214f412SHong Zhang }
10296214f412SHong Zhang 
1030cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
10316214f412SHong Zhang {
10326214f412SHong Zhang   Mat               mat_elemental;
10336214f412SHong Zhang   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols,j,rstart=A->rmap->rstart,rend=A->rmap->rend;
10346214f412SHong Zhang   const PetscInt    *cols;
10356214f412SHong Zhang   const PetscScalar *vals;
10366214f412SHong Zhang 
10376214f412SHong Zhang   PetscFunctionBegin;
1038bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
1039bdeae278SStefano Zampini     mat_elemental = *newmat;
10409566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
1041bdeae278SStefano Zampini   } else {
10429566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
10439566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N));
10449566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental,MATELEMENTAL));
10459566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
1046bdeae278SStefano Zampini   }
10479566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(A));
10486214f412SHong Zhang   for (row=rstart; row<rend; row++) {
10499566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A,row,&ncols,&cols,&vals));
1050bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10519566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES));
10526214f412SHong Zhang     for (j=0; j<ncols; j++) { /* lower triangular part */
1053eb1ec7c1SStefano Zampini       PetscScalar v;
10546214f412SHong Zhang       if (cols[j] == row) continue;
1055b94d7dedSBarry Smith       v    = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j];
10569566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat_elemental,1,&cols[j],1,&row,&v,ADD_VALUES));
10576214f412SHong Zhang     }
10589566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A,row,&ncols,&cols,&vals));
10596214f412SHong Zhang   }
10609566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(A));
10619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
10629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
10636214f412SHong Zhang 
1064511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
10659566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&mat_elemental));
10666214f412SHong Zhang   } else {
10676214f412SHong Zhang     *newmat = mat_elemental;
10686214f412SHong Zhang   }
10696214f412SHong Zhang   PetscFunctionReturn(0);
10706214f412SHong Zhang }
10716214f412SHong Zhang 
1072db31f6deSJed Brown static PetscErrorCode MatDestroy_Elemental(Mat A)
1073db31f6deSJed Brown {
1074db31f6deSJed Brown   Mat_Elemental      *a = (Mat_Elemental*)A->data;
10755e9f5b67SHong Zhang   Mat_Elemental_Grid *commgrid;
10765e9f5b67SHong Zhang   PetscBool          flg;
10775e9f5b67SHong Zhang   MPI_Comm           icomm;
1078db31f6deSJed Brown 
1079db31f6deSJed Brown   PetscFunctionBegin;
1080db31f6deSJed Brown   delete a->emat;
1081ea394475SSatish Balay   delete a->P;
1082ea394475SSatish Balay   delete a->Q;
10835e9f5b67SHong Zhang 
1084e8146892SSatish Balay   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
10859566063dSJacob Faibussowitsch   PetscCall(PetscCommDuplicate(cxxcomm.comm,&icomm,NULL));
10869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_get_attr(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg));
10875e9f5b67SHong Zhang   if (--commgrid->grid_refct == 0) {
10885e9f5b67SHong Zhang     delete commgrid->grid;
10899566063dSJacob Faibussowitsch     PetscCall(PetscFree(commgrid));
10909566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Elemental_keyval));
10915e9f5b67SHong Zhang   }
10929566063dSJacob Faibussowitsch   PetscCall(PetscCommDestroy(&icomm));
10939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",NULL));
10949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL));
10959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_elemental_mpidense_C",NULL));
10969566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1097db31f6deSJed Brown   PetscFunctionReturn(0);
1098db31f6deSJed Brown }
1099db31f6deSJed Brown 
1100db31f6deSJed Brown PetscErrorCode MatSetUp_Elemental(Mat A)
1101db31f6deSJed Brown {
1102db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
1103f19600a0SHong Zhang   MPI_Comm       comm;
1104db31f6deSJed Brown   PetscMPIInt    rsize,csize;
1105f19600a0SHong Zhang   PetscInt       n;
1106db31f6deSJed Brown 
1107db31f6deSJed Brown   PetscFunctionBegin;
11089566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
11099566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
1110db31f6deSJed Brown 
1111d8304050SJose E. Roman   /* Check if local row and column sizes are equally distributed.
1112f19600a0SHong Zhang      Jed: Elemental uses "element" cyclic ordering so the sizes need to match that
1113f19600a0SHong Zhang      exactly.  The strategy in MatElemental is for PETSc to implicitly permute to block ordering (like would be returned by
1114d8304050SJose E. Roman      PetscSplitOwnership(comm,&n,&N), at which point Elemental matrices can act on PETSc vectors without redistributing the vectors. */
11159566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
1116f19600a0SHong Zhang   n = PETSC_DECIDE;
11179566063dSJacob Faibussowitsch   PetscCall(PetscSplitOwnership(comm,&n,&A->rmap->N));
11185f80ce2aSJacob Faibussowitsch   PetscCheck(n == A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %" PetscInt_FMT " of ELEMENTAL matrix must be equally distributed",A->rmap->n);
1119f19600a0SHong Zhang 
1120f19600a0SHong Zhang   n = PETSC_DECIDE;
11219566063dSJacob Faibussowitsch   PetscCall(PetscSplitOwnership(comm,&n,&A->cmap->N));
11225f80ce2aSJacob Faibussowitsch   PetscCheck(n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %" PetscInt_FMT " of ELEMENTAL matrix must be equally distributed",A->cmap->n);
1123f19600a0SHong Zhang 
11242f613bf5SBarry Smith   a->emat->Resize(A->rmap->N,A->cmap->N);
1125e8146892SSatish Balay   El::Zero(*a->emat);
1126db31f6deSJed Brown 
11279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(A->rmap->comm,&rsize));
11289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(A->cmap->comm,&csize));
11295f80ce2aSJacob Faibussowitsch   PetscCheck(csize == rsize,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot use row and column communicators of different sizes");
1130db31f6deSJed Brown   a->commsize = rsize;
1131db31f6deSJed Brown   a->mr[0] = A->rmap->N % rsize; if (!a->mr[0]) a->mr[0] = rsize;
1132db31f6deSJed Brown   a->mr[1] = A->cmap->N % csize; if (!a->mr[1]) a->mr[1] = csize;
1133db31f6deSJed Brown   a->m[0]  = A->rmap->N / rsize + (a->mr[0] != rsize);
1134db31f6deSJed Brown   a->m[1]  = A->cmap->N / csize + (a->mr[1] != csize);
1135db31f6deSJed Brown   PetscFunctionReturn(0);
1136db31f6deSJed Brown }
1137db31f6deSJed Brown 
1138aae2c449SHong Zhang PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType type)
1139aae2c449SHong Zhang {
1140aae2c449SHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data;
1141aae2c449SHong Zhang 
1142aae2c449SHong Zhang   PetscFunctionBegin;
1143fbbcb730SJack Poulson   /* printf("Calling ProcessQueues\n"); */
1144fbbcb730SJack Poulson   a->emat->ProcessQueues();
1145fbbcb730SJack Poulson   /* printf("Finished ProcessQueues\n"); */
1146aae2c449SHong Zhang   PetscFunctionReturn(0);
1147aae2c449SHong Zhang }
1148aae2c449SHong Zhang 
1149aae2c449SHong Zhang PetscErrorCode MatAssemblyEnd_Elemental(Mat A, MatAssemblyType type)
1150aae2c449SHong Zhang {
1151aae2c449SHong Zhang   PetscFunctionBegin;
1152aae2c449SHong Zhang   /* Currently does nothing */
1153aae2c449SHong Zhang   PetscFunctionReturn(0);
1154aae2c449SHong Zhang }
1155aae2c449SHong Zhang 
11567b30e0f1SHong Zhang PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer)
11577b30e0f1SHong Zhang {
11587b30e0f1SHong Zhang   Mat            Adense,Ae;
11597b30e0f1SHong Zhang   MPI_Comm       comm;
11607b30e0f1SHong Zhang 
11617b30e0f1SHong Zhang   PetscFunctionBegin;
11629566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)newMat,&comm));
11639566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&Adense));
11649566063dSJacob Faibussowitsch   PetscCall(MatSetType(Adense,MATDENSE));
11659566063dSJacob Faibussowitsch   PetscCall(MatLoad(Adense,viewer));
11669566063dSJacob Faibussowitsch   PetscCall(MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX,&Ae));
11679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Adense));
11689566063dSJacob Faibussowitsch   PetscCall(MatHeaderReplace(newMat,&Ae));
11697b30e0f1SHong Zhang   PetscFunctionReturn(0);
11707b30e0f1SHong Zhang }
11717b30e0f1SHong Zhang 
117240d92e34SHong Zhang /* -------------------------------------------------------------------*/
117340d92e34SHong Zhang static struct _MatOps MatOps_Values = {
117440d92e34SHong Zhang        MatSetValues_Elemental,
117540d92e34SHong Zhang        0,
117640d92e34SHong Zhang        0,
117740d92e34SHong Zhang        MatMult_Elemental,
117840d92e34SHong Zhang /* 4*/ MatMultAdd_Elemental,
11799426833fSXuan Zhou        MatMultTranspose_Elemental,
1180e883f9d5SXuan Zhou        MatMultTransposeAdd_Elemental,
118140d92e34SHong Zhang        MatSolve_Elemental,
1182df311e6cSXuan Zhou        MatSolveAdd_Elemental,
11832ef15aa2SHong Zhang        0,
11842ef15aa2SHong Zhang /*10*/ 0,
118540d92e34SHong Zhang        MatLUFactor_Elemental,
118640d92e34SHong Zhang        MatCholeskyFactor_Elemental,
118740d92e34SHong Zhang        0,
118840d92e34SHong Zhang        MatTranspose_Elemental,
118940d92e34SHong Zhang /*15*/ MatGetInfo_Elemental,
119040d92e34SHong Zhang        0,
119161119200SXuan Zhou        MatGetDiagonal_Elemental,
1192ade3cc5eSXuan Zhou        MatDiagonalScale_Elemental,
119340d92e34SHong Zhang        MatNorm_Elemental,
119440d92e34SHong Zhang /*20*/ MatAssemblyBegin_Elemental,
119540d92e34SHong Zhang        MatAssemblyEnd_Elemental,
119632d7a744SHong Zhang        MatSetOption_Elemental,
119740d92e34SHong Zhang        MatZeroEntries_Elemental,
119840d92e34SHong Zhang /*24*/ 0,
119940d92e34SHong Zhang        MatLUFactorSymbolic_Elemental,
120040d92e34SHong Zhang        MatLUFactorNumeric_Elemental,
120140d92e34SHong Zhang        MatCholeskyFactorSymbolic_Elemental,
120240d92e34SHong Zhang        MatCholeskyFactorNumeric_Elemental,
120340d92e34SHong Zhang /*29*/ MatSetUp_Elemental,
120440d92e34SHong Zhang        0,
120540d92e34SHong Zhang        0,
120640d92e34SHong Zhang        0,
120740d92e34SHong Zhang        0,
1208df311e6cSXuan Zhou /*34*/ MatDuplicate_Elemental,
120940d92e34SHong Zhang        0,
121040d92e34SHong Zhang        0,
121140d92e34SHong Zhang        0,
121240d92e34SHong Zhang        0,
121340d92e34SHong Zhang /*39*/ MatAXPY_Elemental,
121440d92e34SHong Zhang        0,
121540d92e34SHong Zhang        0,
121640d92e34SHong Zhang        0,
121740d92e34SHong Zhang        MatCopy_Elemental,
121840d92e34SHong Zhang /*44*/ 0,
121940d92e34SHong Zhang        MatScale_Elemental,
12207d68702bSBarry Smith        MatShift_Basic,
122140d92e34SHong Zhang        0,
122240d92e34SHong Zhang        0,
122340d92e34SHong Zhang /*49*/ 0,
122440d92e34SHong Zhang        0,
122540d92e34SHong Zhang        0,
122640d92e34SHong Zhang        0,
122740d92e34SHong Zhang        0,
122840d92e34SHong Zhang /*54*/ 0,
122940d92e34SHong Zhang        0,
123040d92e34SHong Zhang        0,
123140d92e34SHong Zhang        0,
123240d92e34SHong Zhang        0,
123340d92e34SHong Zhang /*59*/ 0,
123440d92e34SHong Zhang        MatDestroy_Elemental,
123540d92e34SHong Zhang        MatView_Elemental,
123640d92e34SHong Zhang        0,
123740d92e34SHong Zhang        0,
123840d92e34SHong Zhang /*64*/ 0,
123940d92e34SHong Zhang        0,
124040d92e34SHong Zhang        0,
124140d92e34SHong Zhang        0,
124240d92e34SHong Zhang        0,
124340d92e34SHong Zhang /*69*/ 0,
124440d92e34SHong Zhang        0,
12452ef0cf24SXuan Zhou        MatConvert_Elemental_Dense,
124640d92e34SHong Zhang        0,
124740d92e34SHong Zhang        0,
124840d92e34SHong Zhang /*74*/ 0,
124940d92e34SHong Zhang        0,
125040d92e34SHong Zhang        0,
125140d92e34SHong Zhang        0,
125240d92e34SHong Zhang        0,
125340d92e34SHong Zhang /*79*/ 0,
125440d92e34SHong Zhang        0,
125540d92e34SHong Zhang        0,
125640d92e34SHong Zhang        0,
12577b30e0f1SHong Zhang        MatLoad_Elemental,
125840d92e34SHong Zhang /*84*/ 0,
125940d92e34SHong Zhang        0,
126040d92e34SHong Zhang        0,
126140d92e34SHong Zhang        0,
126240d92e34SHong Zhang        0,
12634222ddf1SHong Zhang /*89*/ 0,
12644222ddf1SHong Zhang        0,
126540d92e34SHong Zhang        MatMatMultNumeric_Elemental,
126640d92e34SHong Zhang        0,
126740d92e34SHong Zhang        0,
126840d92e34SHong Zhang /*94*/ 0,
12694222ddf1SHong Zhang        0,
12704222ddf1SHong Zhang        0,
1271df311e6cSXuan Zhou        MatMatTransposeMultNumeric_Elemental,
127240d92e34SHong Zhang        0,
12734222ddf1SHong Zhang /*99*/ MatProductSetFromOptions_Elemental,
127440d92e34SHong Zhang        0,
127540d92e34SHong Zhang        0,
1276dfcb0403SXuan Zhou        MatConjugate_Elemental,
127740d92e34SHong Zhang        0,
127840d92e34SHong Zhang /*104*/0,
127940d92e34SHong Zhang        0,
128040d92e34SHong Zhang        0,
128140d92e34SHong Zhang        0,
128240d92e34SHong Zhang        0,
128340d92e34SHong Zhang /*109*/MatMatSolve_Elemental,
128440d92e34SHong Zhang        0,
128540d92e34SHong Zhang        0,
128640d92e34SHong Zhang        0,
128734fa46e1SFrancesco Ballarin        MatMissingDiagonal_Elemental,
128840d92e34SHong Zhang /*114*/0,
128940d92e34SHong Zhang        0,
129040d92e34SHong Zhang        0,
129140d92e34SHong Zhang        0,
129240d92e34SHong Zhang        0,
129340d92e34SHong Zhang /*119*/0,
12944a29722dSXuan Zhou        MatHermitianTranspose_Elemental,
129540d92e34SHong Zhang        0,
129640d92e34SHong Zhang        0,
129740d92e34SHong Zhang        0,
129840d92e34SHong Zhang /*124*/0,
129940d92e34SHong Zhang        0,
130040d92e34SHong Zhang        0,
130140d92e34SHong Zhang        0,
130240d92e34SHong Zhang        0,
130340d92e34SHong Zhang /*129*/0,
130440d92e34SHong Zhang        0,
130540d92e34SHong Zhang        0,
130640d92e34SHong Zhang        0,
130740d92e34SHong Zhang        0,
130840d92e34SHong Zhang /*134*/0,
130940d92e34SHong Zhang        0,
131040d92e34SHong Zhang        0,
131140d92e34SHong Zhang        0,
13124222ddf1SHong Zhang        0,
13134222ddf1SHong Zhang        0,
13144222ddf1SHong Zhang /*140*/0,
13154222ddf1SHong Zhang        0,
13164222ddf1SHong Zhang        0,
13174222ddf1SHong Zhang        0,
13184222ddf1SHong Zhang        0,
13194222ddf1SHong Zhang /*145*/0,
13204222ddf1SHong Zhang        0,
132199a7f59eSMark Adams        0,
132299a7f59eSMark Adams        0,
1323*7fb60732SBarry Smith        0,
1324*7fb60732SBarry Smith /*150*/0
132540d92e34SHong Zhang };
132640d92e34SHong Zhang 
1327ed36708cSHong Zhang /*MC
1328ed36708cSHong Zhang    MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package
1329ed36708cSHong Zhang 
1330c2b89b5dSBarry Smith   Use ./configure --download-elemental to install PETSc to use Elemental
1331c2b89b5dSBarry Smith 
133289bba20eSBarry Smith   Option Database Keys:
1333c2b89b5dSBarry Smith 
1334ed36708cSHong Zhang    Options Database Keys:
13355cc86fc1SJed Brown + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions()
133689bba20eSBarry Smith . -pc_factor_mat_solver_type elemental - to use this direct solver with the option -pc_type lu
13375cc86fc1SJed Brown - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix
1338ed36708cSHong Zhang 
1339ed36708cSHong Zhang   Level: beginner
1340ed36708cSHong Zhang 
134189bba20eSBarry Smith   Note:
134289bba20eSBarry Smith    Note unlike most matrix formats, this format does not store all the matrix entries for a contiguous
134389bba20eSBarry Smith    range of rows on an MPI rank. Use `MatGetOwnershipIS()` to determine what values are stored on
134489bba20eSBarry Smith    the given rank.
134589bba20eSBarry Smith 
134689bba20eSBarry Smith .seealso: `MATDENSE`, `MATSCALAPACK`, `MatGetOwnershipIS()`
1347ed36708cSHong Zhang M*/
13484a29722dSXuan Zhou 
13498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A)
1350db31f6deSJed Brown {
1351db31f6deSJed Brown   Mat_Elemental      *a;
13525682a260SJack Poulson   PetscBool          flg,flg1;
13535e9f5b67SHong Zhang   Mat_Elemental_Grid *commgrid;
13545e9f5b67SHong Zhang   MPI_Comm           icomm;
13555682a260SJack Poulson   PetscInt           optv1;
1356db31f6deSJed Brown 
1357db31f6deSJed Brown   PetscFunctionBegin;
13589566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps)));
135940d92e34SHong Zhang   A->insertmode = NOT_SET_VALUES;
1360db31f6deSJed Brown 
13619566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(A,&a));
1362db31f6deSJed Brown   A->data = (void*)a;
1363db31f6deSJed Brown 
1364db31f6deSJed Brown   /* Set up the elemental matrix */
1365e8146892SSatish Balay   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
13665e9f5b67SHong Zhang 
13675e9f5b67SHong Zhang   /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */
13685e9f5b67SHong Zhang   if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) {
13699566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0));
13709566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister(ElementalCitation,&ElementalCite));
13715e9f5b67SHong Zhang   }
13729566063dSJacob Faibussowitsch   PetscCall(PetscCommDuplicate(cxxcomm.comm,&icomm,NULL));
13739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_get_attr(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg));
13745e9f5b67SHong Zhang   if (!flg) {
13759566063dSJacob Faibussowitsch     PetscCall(PetscNewLog(A,&commgrid));
13765cb544a0SHong Zhang 
1377d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Elemental Options","Mat");
1378e8146892SSatish Balay     /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */
13799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-mat_elemental_grid_height","Grid Height","None",El::mpi::Size(cxxcomm),&optv1,&flg1));
13805682a260SJack Poulson     if (flg1) {
1381aed4548fSBarry Smith       PetscCheck((El::mpi::Size(cxxcomm) % optv1) == 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Grid Height %" PetscInt_FMT " must evenly divide CommSize %" PetscInt_FMT,optv1,(PetscInt)El::mpi::Size(cxxcomm));
1382e8146892SSatish Balay       commgrid->grid = new El::Grid(cxxcomm,optv1); /* use user-provided grid height */
13832ef0cf24SXuan Zhou     } else {
1384e8146892SSatish Balay       commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */
1385da0640a4SHong Zhang       /* printf("new commgrid->grid = %p\n",commgrid->grid);  -- memory leak revealed by valgrind? */
1386ed667823SXuan Zhou     }
13875e9f5b67SHong Zhang     commgrid->grid_refct = 1;
13889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_set_attr(icomm,Petsc_Elemental_keyval,(void*)commgrid));
1389ea394475SSatish Balay 
1390ea394475SSatish Balay     a->pivoting    = 1;
13919566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-mat_elemental_pivoting","Pivoting","None",a->pivoting,&a->pivoting,NULL));
1392ea394475SSatish Balay 
1393d0609cedSBarry Smith     PetscOptionsEnd();
13945e9f5b67SHong Zhang   } else {
13955e9f5b67SHong Zhang     commgrid->grid_refct++;
13965e9f5b67SHong Zhang   }
13979566063dSJacob Faibussowitsch   PetscCall(PetscCommDestroy(&icomm));
13985e9f5b67SHong Zhang   a->grid        = commgrid->grid;
1399e8146892SSatish Balay   a->emat        = new El::DistMatrix<PetscElemScalar>(*a->grid);
140032d7a744SHong Zhang   a->roworiented = PETSC_TRUE;
1401db31f6deSJed Brown 
14029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",MatGetOwnershipIS_Elemental));
14039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_elemental_mpidense_C",MatProductSetFromOptions_Elemental_MPIDense));
14049566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATELEMENTAL));
1405db31f6deSJed Brown   PetscFunctionReturn(0);
1406db31f6deSJed Brown }
1407