xref: /petsc/src/mat/impls/elemental/matelem.cxx (revision c97e14ec040eef239226df753a10ee2a903f9b90)
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 
19d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Elemental(Mat A, PetscViewer viewer)
20d71ae5a4SJacob Faibussowitsch {
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   }
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58db31f6deSJed Brown }
59db31f6deSJed Brown 
60d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetInfo_Elemental(Mat A, MatInfoType flag, MatInfo *info)
61d71ae5a4SJacob Faibussowitsch {
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) {
71462c564dSBarry Smith     //PetscCallMPI(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 */
78462c564dSBarry Smith     //PetscCallMPI(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;
854dfa11a4SJacob Faibussowitsch   info->memory            = 0; /* REVIEW ME */
8615767789SHong Zhang   info->fill_ratio_given  = 0; /* determined by Elemental */
8715767789SHong Zhang   info->fill_ratio_needed = 0;
8815767789SHong Zhang   info->factor_mallocs    = 0;
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90180a43e4SHong Zhang }
91180a43e4SHong Zhang 
9266976f2fSJacob Faibussowitsch static PetscErrorCode MatSetOption_Elemental(Mat A, MatOption op, PetscBool flg)
93d71ae5a4SJacob Faibussowitsch {
9432d7a744SHong Zhang   Mat_Elemental *a = (Mat_Elemental *)A->data;
9532d7a744SHong Zhang 
9632d7a744SHong Zhang   PetscFunctionBegin;
9732d7a744SHong Zhang   switch (op) {
98d71ae5a4SJacob Faibussowitsch   case MAT_ROW_ORIENTED:
99d71ae5a4SJacob Faibussowitsch     a->roworiented = flg;
100d71ae5a4SJacob Faibussowitsch     break;
101d71ae5a4SJacob Faibussowitsch   default:
102888c827cSStefano Zampini     break;
10332d7a744SHong Zhang   }
1043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10532d7a744SHong Zhang }
10632d7a744SHong Zhang 
107d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValues_Elemental(Mat A, PetscInt nr, const PetscInt *rows, PetscInt nc, const PetscInt *cols, const PetscScalar *vals, InsertMode imode)
108d71ae5a4SJacob Faibussowitsch {
109db31f6deSJed Brown   Mat_Elemental *a = (Mat_Elemental *)A->data;
110fbbcb730SJack Poulson   PetscInt       i, j, rrank, ridx, crank, cidx, erow, ecol, numQueues = 0;
111db31f6deSJed Brown 
112db31f6deSJed Brown   PetscFunctionBegin;
113fbbcb730SJack Poulson   // TODO: Initialize matrix to all zeros?
114fbbcb730SJack Poulson 
115fbbcb730SJack Poulson   // Count the number of queues from this process
11632d7a744SHong Zhang   if (a->roworiented) {
117db31f6deSJed Brown     for (i = 0; i < nr; i++) {
118db31f6deSJed Brown       if (rows[i] < 0) continue;
119db31f6deSJed Brown       P2RO(A, 0, rows[i], &rrank, &ridx);
120db31f6deSJed Brown       RO2E(A, 0, rrank, ridx, &erow);
121aed4548fSBarry Smith       PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect row translation");
122db31f6deSJed Brown       for (j = 0; j < nc; j++) {
123db31f6deSJed Brown         if (cols[j] < 0) continue;
124db31f6deSJed Brown         P2RO(A, 1, cols[j], &crank, &cidx);
125db31f6deSJed Brown         RO2E(A, 1, crank, cidx, &ecol);
126aed4548fSBarry Smith         PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect col translation");
127fbbcb730SJack Poulson         if (!a->emat->IsLocal(erow, ecol)) { /* off-proc entry */
128fbbcb730SJack Poulson           /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
12908401ef6SPierre Jolivet           PetscCheck(imode == ADD_VALUES, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only ADD_VALUES to off-processor entry is supported");
130fbbcb730SJack Poulson           ++numQueues;
131aae2c449SHong Zhang           continue;
132ed36708cSHong Zhang         }
133fbbcb730SJack Poulson         /* printf("Locally updating (%d,%d)\n",erow,ecol); */
134db31f6deSJed Brown         switch (imode) {
135d71ae5a4SJacob Faibussowitsch         case INSERT_VALUES:
136d71ae5a4SJacob Faibussowitsch           a->emat->Set(erow, ecol, (PetscElemScalar)vals[i * nc + j]);
137d71ae5a4SJacob Faibussowitsch           break;
138d71ae5a4SJacob Faibussowitsch         case ADD_VALUES:
139d71ae5a4SJacob Faibussowitsch           a->emat->Update(erow, ecol, (PetscElemScalar)vals[i * nc + j]);
140d71ae5a4SJacob Faibussowitsch           break;
141d71ae5a4SJacob Faibussowitsch         default:
142d71ae5a4SJacob Faibussowitsch           SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for InsertMode %d", (int)imode);
143db31f6deSJed Brown         }
144db31f6deSJed Brown       }
145db31f6deSJed Brown     }
146fbbcb730SJack Poulson 
147fbbcb730SJack Poulson     /* printf("numQueues=%d\n",numQueues); */
148fbbcb730SJack Poulson     a->emat->Reserve(numQueues);
149fbbcb730SJack Poulson     for (i = 0; i < nr; i++) {
150fbbcb730SJack Poulson       if (rows[i] < 0) continue;
151fbbcb730SJack Poulson       P2RO(A, 0, rows[i], &rrank, &ridx);
152fbbcb730SJack Poulson       RO2E(A, 0, rrank, ridx, &erow);
153fbbcb730SJack Poulson       for (j = 0; j < nc; j++) {
154fbbcb730SJack Poulson         if (cols[j] < 0) continue;
155fbbcb730SJack Poulson         P2RO(A, 1, cols[j], &crank, &cidx);
156fbbcb730SJack Poulson         RO2E(A, 1, crank, cidx, &ecol);
157fbbcb730SJack Poulson         if (!a->emat->IsLocal(erow, ecol)) { /*off-proc entry*/
158fbbcb730SJack Poulson           /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
159fbbcb730SJack Poulson           a->emat->QueueUpdate(erow, ecol, vals[i * nc + j]);
160fbbcb730SJack Poulson         }
161fbbcb730SJack Poulson       }
162fbbcb730SJack Poulson     }
1635e116b59SBarry Smith   } else { /* column-oriented */
16432d7a744SHong Zhang     for (j = 0; j < nc; j++) {
16532d7a744SHong Zhang       if (cols[j] < 0) continue;
16632d7a744SHong Zhang       P2RO(A, 1, cols[j], &crank, &cidx);
16732d7a744SHong Zhang       RO2E(A, 1, crank, cidx, &ecol);
168aed4548fSBarry Smith       PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect col translation");
16932d7a744SHong Zhang       for (i = 0; i < nr; i++) {
17032d7a744SHong Zhang         if (rows[i] < 0) continue;
17132d7a744SHong Zhang         P2RO(A, 0, rows[i], &rrank, &ridx);
17232d7a744SHong Zhang         RO2E(A, 0, rrank, ridx, &erow);
173aed4548fSBarry Smith         PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect row translation");
17432d7a744SHong Zhang         if (!a->emat->IsLocal(erow, ecol)) { /* off-proc entry */
17532d7a744SHong Zhang           /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
17608401ef6SPierre Jolivet           PetscCheck(imode == ADD_VALUES, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only ADD_VALUES to off-processor entry is supported");
17732d7a744SHong Zhang           ++numQueues;
17832d7a744SHong Zhang           continue;
17932d7a744SHong Zhang         }
18032d7a744SHong Zhang         /* printf("Locally updating (%d,%d)\n",erow,ecol); */
18132d7a744SHong Zhang         switch (imode) {
182d71ae5a4SJacob Faibussowitsch         case INSERT_VALUES:
183d71ae5a4SJacob Faibussowitsch           a->emat->Set(erow, ecol, (PetscElemScalar)vals[i + j * nr]);
184d71ae5a4SJacob Faibussowitsch           break;
185d71ae5a4SJacob Faibussowitsch         case ADD_VALUES:
186d71ae5a4SJacob Faibussowitsch           a->emat->Update(erow, ecol, (PetscElemScalar)vals[i + j * nr]);
187d71ae5a4SJacob Faibussowitsch           break;
188d71ae5a4SJacob Faibussowitsch         default:
189d71ae5a4SJacob Faibussowitsch           SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for InsertMode %d", (int)imode);
19032d7a744SHong Zhang         }
19132d7a744SHong Zhang       }
19232d7a744SHong Zhang     }
19332d7a744SHong Zhang 
19432d7a744SHong Zhang     /* printf("numQueues=%d\n",numQueues); */
19532d7a744SHong Zhang     a->emat->Reserve(numQueues);
19632d7a744SHong Zhang     for (j = 0; j < nc; j++) {
19732d7a744SHong Zhang       if (cols[j] < 0) continue;
19832d7a744SHong Zhang       P2RO(A, 1, cols[j], &crank, &cidx);
19932d7a744SHong Zhang       RO2E(A, 1, crank, cidx, &ecol);
20032d7a744SHong Zhang 
20132d7a744SHong Zhang       for (i = 0; i < nr; i++) {
20232d7a744SHong Zhang         if (rows[i] < 0) continue;
20332d7a744SHong Zhang         P2RO(A, 0, rows[i], &rrank, &ridx);
20432d7a744SHong Zhang         RO2E(A, 0, rrank, ridx, &erow);
20532d7a744SHong Zhang         if (!a->emat->IsLocal(erow, ecol)) { /*off-proc entry*/
20632d7a744SHong Zhang           /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
20732d7a744SHong Zhang           a->emat->QueueUpdate(erow, ecol, vals[i + j * nr]);
20832d7a744SHong Zhang         }
20932d7a744SHong Zhang       }
21032d7a744SHong Zhang     }
21132d7a744SHong Zhang   }
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
213db31f6deSJed Brown }
214db31f6deSJed Brown 
215d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Elemental(Mat A, Vec X, Vec Y)
216d71ae5a4SJacob Faibussowitsch {
217db31f6deSJed Brown   Mat_Elemental         *a = (Mat_Elemental *)A->data;
218e6dea9dbSXuan Zhou   const PetscElemScalar *x;
219e6dea9dbSXuan Zhou   PetscElemScalar       *y;
220df311e6cSXuan Zhou   PetscElemScalar        one = 1, zero = 0;
221db31f6deSJed Brown 
222db31f6deSJed Brown   PetscFunctionBegin;
2239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y, (PetscScalar **)&y));
225db31f6deSJed Brown   { /* Scoping so that constructor is called before pointer is returned */
226e8146892SSatish Balay     El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ye;
2270c18141cSBarry Smith     xe.LockedAttach(A->cmap->N, 1, *a->grid, 0, 0, x, A->cmap->n);
2280c18141cSBarry Smith     ye.Attach(A->rmap->N, 1, *a->grid, 0, 0, y, A->rmap->n);
229e8146892SSatish Balay     El::Gemv(El::NORMAL, one, *a->emat, xe, zero, ye);
230db31f6deSJed Brown   }
2319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x));
2329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y, (PetscScalar **)&y));
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234db31f6deSJed Brown }
235db31f6deSJed Brown 
236d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Elemental(Mat A, Vec X, Vec Y)
237d71ae5a4SJacob Faibussowitsch {
2389426833fSXuan Zhou   Mat_Elemental         *a = (Mat_Elemental *)A->data;
239df311e6cSXuan Zhou   const PetscElemScalar *x;
240df311e6cSXuan Zhou   PetscElemScalar       *y;
241e6dea9dbSXuan Zhou   PetscElemScalar        one = 1, zero = 0;
2429426833fSXuan Zhou 
2439426833fSXuan Zhou   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
2459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y, (PetscScalar **)&y));
2469426833fSXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
247e8146892SSatish Balay     El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ye;
2480c18141cSBarry Smith     xe.LockedAttach(A->rmap->N, 1, *a->grid, 0, 0, x, A->rmap->n);
2490c18141cSBarry Smith     ye.Attach(A->cmap->N, 1, *a->grid, 0, 0, y, A->cmap->n);
250e8146892SSatish Balay     El::Gemv(El::TRANSPOSE, one, *a->emat, xe, zero, ye);
2519426833fSXuan Zhou   }
2529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x));
2539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y, (PetscScalar **)&y));
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2559426833fSXuan Zhou }
2569426833fSXuan Zhou 
257d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_Elemental(Mat A, Vec X, Vec Y, Vec Z)
258d71ae5a4SJacob Faibussowitsch {
259db31f6deSJed Brown   Mat_Elemental         *a = (Mat_Elemental *)A->data;
260df311e6cSXuan Zhou   const PetscElemScalar *x;
261df311e6cSXuan Zhou   PetscElemScalar       *z;
262e6dea9dbSXuan Zhou   PetscElemScalar        one = 1;
263db31f6deSJed Brown 
264db31f6deSJed Brown   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   if (Y != Z) PetscCall(VecCopy(Y, Z));
2669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
2679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Z, (PetscScalar **)&z));
268db31f6deSJed Brown   { /* Scoping so that constructor is called before pointer is returned */
269e8146892SSatish Balay     El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ze;
2700c18141cSBarry Smith     xe.LockedAttach(A->cmap->N, 1, *a->grid, 0, 0, x, A->cmap->n);
2710c18141cSBarry Smith     ze.Attach(A->rmap->N, 1, *a->grid, 0, 0, z, A->rmap->n);
272e8146892SSatish Balay     El::Gemv(El::NORMAL, one, *a->emat, xe, one, ze);
273db31f6deSJed Brown   }
2749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x));
2759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Z, (PetscScalar **)&z));
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
277db31f6deSJed Brown }
278db31f6deSJed Brown 
279d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_Elemental(Mat A, Vec X, Vec Y, Vec Z)
280d71ae5a4SJacob Faibussowitsch {
281e883f9d5SXuan Zhou   Mat_Elemental         *a = (Mat_Elemental *)A->data;
282df311e6cSXuan Zhou   const PetscElemScalar *x;
283df311e6cSXuan Zhou   PetscElemScalar       *z;
284e6dea9dbSXuan Zhou   PetscElemScalar        one = 1;
285e883f9d5SXuan Zhou 
286e883f9d5SXuan Zhou   PetscFunctionBegin;
2879566063dSJacob Faibussowitsch   if (Y != Z) PetscCall(VecCopy(Y, Z));
2889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Z, (PetscScalar **)&z));
290e883f9d5SXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
291e8146892SSatish Balay     El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ze;
2920c18141cSBarry Smith     xe.LockedAttach(A->rmap->N, 1, *a->grid, 0, 0, x, A->rmap->n);
2930c18141cSBarry Smith     ze.Attach(A->cmap->N, 1, *a->grid, 0, 0, z, A->cmap->n);
294e8146892SSatish Balay     El::Gemv(El::TRANSPOSE, one, *a->emat, xe, one, ze);
295e883f9d5SXuan Zhou   }
2969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x));
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Z, (PetscScalar **)&z));
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
299e883f9d5SXuan Zhou }
300e883f9d5SXuan Zhou 
301d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultNumeric_Elemental(Mat A, Mat B, Mat C)
302d71ae5a4SJacob Faibussowitsch {
303c1d1b975SXuan Zhou   Mat_Elemental  *a   = (Mat_Elemental *)A->data;
304c1d1b975SXuan Zhou   Mat_Elemental  *b   = (Mat_Elemental *)B->data;
3059a9e8502SHong Zhang   Mat_Elemental  *c   = (Mat_Elemental *)C->data;
306e6dea9dbSXuan Zhou   PetscElemScalar one = 1, zero = 0;
307c1d1b975SXuan Zhou 
308c1d1b975SXuan Zhou   PetscFunctionBegin;
309aae2c449SHong Zhang   { /* Scoping so that constructor is called before pointer is returned */
310e8146892SSatish Balay     El::Gemm(El::NORMAL, El::NORMAL, one, *a->emat, *b->emat, zero, *c->emat);
311aae2c449SHong Zhang   }
3129a9e8502SHong Zhang   C->assembled = PETSC_TRUE;
3133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3149a9e8502SHong Zhang }
3159a9e8502SHong Zhang 
316cedd07caSPierre Jolivet PetscErrorCode MatMatMultSymbolic_Elemental(Mat A, Mat B, PetscReal, Mat Ce)
317d71ae5a4SJacob Faibussowitsch {
3189a9e8502SHong Zhang   PetscFunctionBegin;
3199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Ce, A->rmap->n, B->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
3209566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ce, MATELEMENTAL));
3219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Ce));
3224222ddf1SHong Zhang   Ce->ops->matmultnumeric = MatMatMultNumeric_Elemental;
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
324c1d1b975SXuan Zhou }
325c1d1b975SXuan Zhou 
326d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatTransposeMultNumeric_Elemental(Mat A, Mat B, Mat C)
327d71ae5a4SJacob Faibussowitsch {
328df311e6cSXuan Zhou   Mat_Elemental  *a   = (Mat_Elemental *)A->data;
329df311e6cSXuan Zhou   Mat_Elemental  *b   = (Mat_Elemental *)B->data;
330df311e6cSXuan Zhou   Mat_Elemental  *c   = (Mat_Elemental *)C->data;
331e6dea9dbSXuan Zhou   PetscElemScalar one = 1, zero = 0;
332df311e6cSXuan Zhou 
333df311e6cSXuan Zhou   PetscFunctionBegin;
334df311e6cSXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
335e8146892SSatish Balay     El::Gemm(El::NORMAL, El::TRANSPOSE, one, *a->emat, *b->emat, zero, *c->emat);
336df311e6cSXuan Zhou   }
337df311e6cSXuan Zhou   C->assembled = PETSC_TRUE;
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
339df311e6cSXuan Zhou }
340df311e6cSXuan Zhou 
341cedd07caSPierre Jolivet static PetscErrorCode MatMatTransposeMultSymbolic_Elemental(Mat A, Mat B, PetscReal, Mat C)
342d71ae5a4SJacob Faibussowitsch {
343df311e6cSXuan Zhou   PetscFunctionBegin;
3449566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, A->rmap->n, B->rmap->n, PETSC_DECIDE, PETSC_DECIDE));
3459566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATELEMENTAL));
3469566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
348df311e6cSXuan Zhou }
349df311e6cSXuan Zhou 
350d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Elemental_AB(Mat C)
351d71ae5a4SJacob Faibussowitsch {
3524222ddf1SHong Zhang   PetscFunctionBegin;
3534222ddf1SHong Zhang   C->ops->matmultsymbolic = MatMatMultSymbolic_Elemental;
3544222ddf1SHong Zhang   C->ops->productsymbolic = MatProductSymbolic_AB;
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3564222ddf1SHong Zhang }
3574222ddf1SHong Zhang 
358d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Elemental_ABt(Mat C)
359d71ae5a4SJacob Faibussowitsch {
3604222ddf1SHong Zhang   PetscFunctionBegin;
3614222ddf1SHong Zhang   C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_Elemental;
3624222ddf1SHong Zhang   C->ops->productsymbolic          = MatProductSymbolic_ABt;
3633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3644222ddf1SHong Zhang }
3654222ddf1SHong Zhang 
366d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Elemental(Mat C)
367d71ae5a4SJacob Faibussowitsch {
3684222ddf1SHong Zhang   Mat_Product *product = C->product;
3694222ddf1SHong Zhang 
3704222ddf1SHong Zhang   PetscFunctionBegin;
3714222ddf1SHong Zhang   switch (product->type) {
372d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AB:
373d71ae5a4SJacob Faibussowitsch     PetscCall(MatProductSetFromOptions_Elemental_AB(C));
374d71ae5a4SJacob Faibussowitsch     break;
375d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_ABt:
376d71ae5a4SJacob Faibussowitsch     PetscCall(MatProductSetFromOptions_Elemental_ABt(C));
377d71ae5a4SJacob Faibussowitsch     break;
378d71ae5a4SJacob Faibussowitsch   default:
379d71ae5a4SJacob Faibussowitsch     break;
3804222ddf1SHong Zhang   }
3813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3824222ddf1SHong Zhang }
383f6e5db1bSPierre Jolivet 
38466976f2fSJacob Faibussowitsch static PetscErrorCode MatMatMultNumeric_Elemental_MPIDense(Mat A, Mat B, Mat C)
385d71ae5a4SJacob Faibussowitsch {
386f6e5db1bSPierre Jolivet   Mat Be, Ce;
387f6e5db1bSPierre Jolivet 
388f6e5db1bSPierre Jolivet   PetscFunctionBegin;
3899566063dSJacob Faibussowitsch   PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &Be));
390fb842aefSJose E. Roman   PetscCall(MatMatMult(A, Be, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Ce));
3919566063dSJacob Faibussowitsch   PetscCall(MatConvert(Ce, MATMPIDENSE, MAT_REUSE_MATRIX, &C));
3929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Be));
3939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ce));
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
395f6e5db1bSPierre Jolivet }
396f6e5db1bSPierre Jolivet 
39766976f2fSJacob Faibussowitsch static PetscErrorCode MatMatMultSymbolic_Elemental_MPIDense(Mat A, Mat B, PetscReal, Mat C)
398d71ae5a4SJacob Faibussowitsch {
399f6e5db1bSPierre Jolivet   PetscFunctionBegin;
4009566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
4019566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATMPIDENSE));
4029566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
403f6e5db1bSPierre Jolivet   C->ops->matmultnumeric = MatMatMultNumeric_Elemental_MPIDense;
4043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
405f6e5db1bSPierre Jolivet }
406f6e5db1bSPierre Jolivet 
40766976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense_AB(Mat C)
408d71ae5a4SJacob Faibussowitsch {
409f6e5db1bSPierre Jolivet   PetscFunctionBegin;
410f6e5db1bSPierre Jolivet   C->ops->matmultsymbolic = MatMatMultSymbolic_Elemental_MPIDense;
411f6e5db1bSPierre Jolivet   C->ops->productsymbolic = MatProductSymbolic_AB;
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413f6e5db1bSPierre Jolivet }
414f6e5db1bSPierre Jolivet 
41566976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense(Mat C)
416d71ae5a4SJacob Faibussowitsch {
417f6e5db1bSPierre Jolivet   Mat_Product *product = C->product;
418f6e5db1bSPierre Jolivet 
419f6e5db1bSPierre Jolivet   PetscFunctionBegin;
42048a46eb9SPierre Jolivet   if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Elemental_MPIDense_AB(C));
4213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
422f6e5db1bSPierre Jolivet }
4234222ddf1SHong Zhang 
424d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Elemental(Mat A, Vec D)
425d71ae5a4SJacob Faibussowitsch {
426a9d89745SXuan Zhou   PetscInt        i, nrows, ncols, nD, rrank, ridx, crank, cidx;
427a9d89745SXuan Zhou   Mat_Elemental  *a = (Mat_Elemental *)A->data;
428a9d89745SXuan Zhou   PetscElemScalar v;
429ce94432eSBarry Smith   MPI_Comm        comm;
43061119200SXuan Zhou 
43161119200SXuan Zhou   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
4339566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &nrows, &ncols));
434a9d89745SXuan Zhou   nD = nrows > ncols ? ncols : nrows;
435a9d89745SXuan Zhou   for (i = 0; i < nD; i++) {
436a9d89745SXuan Zhou     PetscInt erow, ecol;
437a9d89745SXuan Zhou     P2RO(A, 0, i, &rrank, &ridx);
438a9d89745SXuan Zhou     RO2E(A, 0, rrank, ridx, &erow);
439aed4548fSBarry Smith     PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation");
440a9d89745SXuan Zhou     P2RO(A, 1, i, &crank, &cidx);
441a9d89745SXuan Zhou     RO2E(A, 1, crank, cidx, &ecol);
442aed4548fSBarry Smith     PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation");
443a9d89745SXuan Zhou     v = a->emat->Get(erow, ecol);
4449566063dSJacob Faibussowitsch     PetscCall(VecSetValues(D, 1, &i, (PetscScalar *)&v, INSERT_VALUES));
445ade3cc5eSXuan Zhou   }
4469566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(D));
4479566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(D));
4483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44961119200SXuan Zhou }
45061119200SXuan Zhou 
451d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDiagonalScale_Elemental(Mat X, Vec L, Vec R)
452d71ae5a4SJacob Faibussowitsch {
453ade3cc5eSXuan Zhou   Mat_Elemental         *x = (Mat_Elemental *)X->data;
454ade3cc5eSXuan Zhou   const PetscElemScalar *d;
455ade3cc5eSXuan Zhou 
456ade3cc5eSXuan Zhou   PetscFunctionBegin;
4579065cd98SJed Brown   if (R) {
4589566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(R, (const PetscScalar **)&d));
459e8146892SSatish Balay     El::DistMatrix<PetscElemScalar, El::VC, El::STAR> de;
4600c18141cSBarry Smith     de.LockedAttach(X->cmap->N, 1, *x->grid, 0, 0, d, X->cmap->n);
461e8146892SSatish Balay     El::DiagonalScale(El::RIGHT, El::NORMAL, de, *x->emat);
4629566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(R, (const PetscScalar **)&d));
4639065cd98SJed Brown   }
4649065cd98SJed Brown   if (L) {
4659566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(L, (const PetscScalar **)&d));
466e8146892SSatish Balay     El::DistMatrix<PetscElemScalar, El::VC, El::STAR> de;
4670c18141cSBarry Smith     de.LockedAttach(X->rmap->N, 1, *x->grid, 0, 0, d, X->rmap->n);
468e8146892SSatish Balay     El::DiagonalScale(El::LEFT, El::NORMAL, de, *x->emat);
4699566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(L, (const PetscScalar **)&d));
470ade3cc5eSXuan Zhou   }
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
472ade3cc5eSXuan Zhou }
473ade3cc5eSXuan Zhou 
474cedd07caSPierre Jolivet static PetscErrorCode MatMissingDiagonal_Elemental(Mat, PetscBool *missing, PetscInt *)
475d71ae5a4SJacob Faibussowitsch {
47634fa46e1SFrancesco Ballarin   PetscFunctionBegin;
47734fa46e1SFrancesco Ballarin   *missing = PETSC_FALSE;
4783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47934fa46e1SFrancesco Ballarin }
48034fa46e1SFrancesco Ballarin 
481d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_Elemental(Mat X, PetscScalar a)
482d71ae5a4SJacob Faibussowitsch {
48365b78793SXuan Zhou   Mat_Elemental *x = (Mat_Elemental *)X->data;
48465b78793SXuan Zhou 
48565b78793SXuan Zhou   PetscFunctionBegin;
486e8146892SSatish Balay   El::Scale((PetscElemScalar)a, *x->emat);
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48865b78793SXuan Zhou }
48965b78793SXuan Zhou 
490f82baa17SHong Zhang /*
491f82baa17SHong Zhang   MatAXPY - Computes Y = a*X + Y.
492f82baa17SHong Zhang */
493cedd07caSPierre Jolivet static PetscErrorCode MatAXPY_Elemental(Mat Y, PetscScalar a, Mat X, MatStructure)
494d71ae5a4SJacob Faibussowitsch {
495e09a3074SHong Zhang   Mat_Elemental *x = (Mat_Elemental *)X->data;
496e09a3074SHong Zhang   Mat_Elemental *y = (Mat_Elemental *)Y->data;
497e09a3074SHong Zhang 
498e09a3074SHong Zhang   PetscFunctionBegin;
499e8146892SSatish Balay   El::Axpy((PetscElemScalar)a, *x->emat, *y->emat);
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
5013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
502e09a3074SHong Zhang }
503e09a3074SHong Zhang 
504cedd07caSPierre Jolivet static PetscErrorCode MatCopy_Elemental(Mat A, Mat B, MatStructure)
505d71ae5a4SJacob Faibussowitsch {
506d6223691SXuan Zhou   Mat_Elemental *a = (Mat_Elemental *)A->data;
507d6223691SXuan Zhou   Mat_Elemental *b = (Mat_Elemental *)B->data;
508d6223691SXuan Zhou 
509d6223691SXuan Zhou   PetscFunctionBegin;
510e8146892SSatish Balay   El::Copy(*a->emat, *b->emat);
5119566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)B));
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
513d6223691SXuan Zhou }
514d6223691SXuan Zhou 
515d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_Elemental(Mat A, MatDuplicateOption op, Mat *B)
516d71ae5a4SJacob Faibussowitsch {
517df311e6cSXuan Zhou   Mat            Be;
518ce94432eSBarry Smith   MPI_Comm       comm;
519df311e6cSXuan Zhou   Mat_Elemental *a = (Mat_Elemental *)A->data;
520df311e6cSXuan Zhou 
521df311e6cSXuan Zhou   PetscFunctionBegin;
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5239566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Be));
5249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Be, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
5259566063dSJacob Faibussowitsch   PetscCall(MatSetType(Be, MATELEMENTAL));
5269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Be));
527df311e6cSXuan Zhou   *B = Be;
528df311e6cSXuan Zhou   if (op == MAT_COPY_VALUES) {
529df311e6cSXuan Zhou     Mat_Elemental *b = (Mat_Elemental *)Be->data;
530e8146892SSatish Balay     El::Copy(*a->emat, *b->emat);
531df311e6cSXuan Zhou   }
532df311e6cSXuan Zhou   Be->assembled = PETSC_TRUE;
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
534df311e6cSXuan Zhou }
535df311e6cSXuan Zhou 
536d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Elemental(Mat A, MatReuse reuse, Mat *B)
537d71ae5a4SJacob Faibussowitsch {
538ec8cb81fSBarry Smith   Mat            Be = *B;
539ce94432eSBarry Smith   MPI_Comm       comm;
5404fe7bbcaSHong Zhang   Mat_Elemental *a = (Mat_Elemental *)A->data, *b;
541d6223691SXuan Zhou 
542d6223691SXuan Zhou   PetscFunctionBegin;
5437fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5455cb544a0SHong Zhang   /* Only out-of-place supported */
5465f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX, comm, PETSC_ERR_SUP, "Only out-of-place supported");
5475262d616SXuan Zhou   if (reuse == MAT_INITIAL_MATRIX) {
5489566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &Be));
5499566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Be, A->cmap->n, A->rmap->n, PETSC_DECIDE, PETSC_DECIDE));
5509566063dSJacob Faibussowitsch     PetscCall(MatSetType(Be, MATELEMENTAL));
5519566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Be));
5525262d616SXuan Zhou     *B = Be;
5535262d616SXuan Zhou   }
5544fe7bbcaSHong Zhang   b = (Mat_Elemental *)Be->data;
555e8146892SSatish Balay   El::Transpose(*a->emat, *b->emat);
5565262d616SXuan Zhou   Be->assembled = PETSC_TRUE;
5573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
558d6223691SXuan Zhou }
559d6223691SXuan Zhou 
560d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConjugate_Elemental(Mat A)
561d71ae5a4SJacob Faibussowitsch {
562dfcb0403SXuan Zhou   Mat_Elemental *a = (Mat_Elemental *)A->data;
563dfcb0403SXuan Zhou 
564dfcb0403SXuan Zhou   PetscFunctionBegin;
565e8146892SSatish Balay   El::Conjugate(*a->emat);
5663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
567dfcb0403SXuan Zhou }
568dfcb0403SXuan Zhou 
569d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHermitianTranspose_Elemental(Mat A, MatReuse reuse, Mat *B)
570d71ae5a4SJacob Faibussowitsch {
571ec8cb81fSBarry Smith   Mat            Be = *B;
572ce94432eSBarry Smith   MPI_Comm       comm;
5734a29722dSXuan Zhou   Mat_Elemental *a = (Mat_Elemental *)A->data, *b;
5744a29722dSXuan Zhou 
5754a29722dSXuan Zhou   PetscFunctionBegin;
5769566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5775cb544a0SHong Zhang   /* Only out-of-place supported */
5784a29722dSXuan Zhou   if (reuse == MAT_INITIAL_MATRIX) {
5799566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &Be));
5809566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Be, A->cmap->n, A->rmap->n, PETSC_DECIDE, PETSC_DECIDE));
5819566063dSJacob Faibussowitsch     PetscCall(MatSetType(Be, MATELEMENTAL));
5829566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Be));
5834a29722dSXuan Zhou     *B = Be;
5844a29722dSXuan Zhou   }
5854a29722dSXuan Zhou   b = (Mat_Elemental *)Be->data;
586e8146892SSatish Balay   El::Adjoint(*a->emat, *b->emat);
5874a29722dSXuan Zhou   Be->assembled = PETSC_TRUE;
5883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5894a29722dSXuan Zhou }
5904a29722dSXuan Zhou 
591d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_Elemental(Mat A, Vec B, Vec X)
592d71ae5a4SJacob Faibussowitsch {
5931f881ff8SXuan Zhou   Mat_Elemental   *a = (Mat_Elemental *)A->data;
594df311e6cSXuan Zhou   PetscElemScalar *x;
595ea394475SSatish Balay   PetscInt         pivoting = a->pivoting;
5961f881ff8SXuan Zhou 
5971f881ff8SXuan Zhou   PetscFunctionBegin;
5989566063dSJacob Faibussowitsch   PetscCall(VecCopy(B, X));
5999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, (PetscScalar **)&x));
600ea394475SSatish Balay 
601e8146892SSatish Balay   El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe;
6020c18141cSBarry Smith   xe.Attach(A->rmap->N, 1, *a->grid, 0, 0, x, A->rmap->n);
603e8146892SSatish Balay   El::DistMatrix<PetscElemScalar, El::MC, El::MR> xer(xe);
604fc54b460SXuan Zhou   switch (A->factortype) {
605fc54b460SXuan Zhou   case MAT_FACTOR_LU:
606ea394475SSatish Balay     if (pivoting == 0) {
607e8146892SSatish Balay       El::lu::SolveAfter(El::NORMAL, *a->emat, xer);
608ea394475SSatish Balay     } else if (pivoting == 1) {
609ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, xer);
610ea394475SSatish Balay     } else { /* pivoting == 2 */
611ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *a->Q, xer);
6121f881ff8SXuan Zhou     }
613fc54b460SXuan Zhou     break;
614d71ae5a4SJacob Faibussowitsch   case MAT_FACTOR_CHOLESKY:
615d71ae5a4SJacob Faibussowitsch     El::cholesky::SolveAfter(El::UPPER, El::NORMAL, *a->emat, xer);
616d71ae5a4SJacob Faibussowitsch     break;
617d71ae5a4SJacob Faibussowitsch   default:
618d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unfactored Matrix or Unsupported MatFactorType");
619d71ae5a4SJacob Faibussowitsch     break;
6201f881ff8SXuan Zhou   }
621ea394475SSatish Balay   El::Copy(xer, xe);
622ea394475SSatish Balay 
6239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, (PetscScalar **)&x));
6243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6251f881ff8SXuan Zhou }
6261f881ff8SXuan Zhou 
627d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveAdd_Elemental(Mat A, Vec B, Vec Y, Vec X)
628d71ae5a4SJacob Faibussowitsch {
629df311e6cSXuan Zhou   PetscFunctionBegin;
6309566063dSJacob Faibussowitsch   PetscCall(MatSolve_Elemental(A, B, X));
6319566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, 1, Y));
6323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
633df311e6cSXuan Zhou }
634df311e6cSXuan Zhou 
635d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_Elemental(Mat A, Mat B, Mat X)
636d71ae5a4SJacob Faibussowitsch {
6371f0e42cfSHong Zhang   Mat_Elemental *a = (Mat_Elemental *)A->data;
63842ba530cSPierre Jolivet   Mat_Elemental *x;
63942ba530cSPierre Jolivet   Mat            C;
640ea394475SSatish Balay   PetscInt       pivoting = a->pivoting;
64142ba530cSPierre Jolivet   PetscBool      flg;
64242ba530cSPierre Jolivet   MatType        type;
6431f0e42cfSHong Zhang 
644ae844d54SHong Zhang   PetscFunctionBegin;
6459566063dSJacob Faibussowitsch   PetscCall(MatGetType(X, &type));
6469566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type, MATELEMENTAL, &flg));
64742ba530cSPierre Jolivet   if (!flg) {
6489566063dSJacob Faibussowitsch     PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &C));
64942ba530cSPierre Jolivet     x = (Mat_Elemental *)C->data;
65042ba530cSPierre Jolivet   } else {
65142ba530cSPierre Jolivet     x = (Mat_Elemental *)X->data;
65242ba530cSPierre Jolivet     El::Copy(*((Mat_Elemental *)B->data)->emat, *x->emat);
65342ba530cSPierre Jolivet   }
654fc54b460SXuan Zhou   switch (A->factortype) {
655fc54b460SXuan Zhou   case MAT_FACTOR_LU:
656ea394475SSatish Balay     if (pivoting == 0) {
657e8146892SSatish Balay       El::lu::SolveAfter(El::NORMAL, *a->emat, *x->emat);
658ea394475SSatish Balay     } else if (pivoting == 1) {
659ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *x->emat);
660ea394475SSatish Balay     } else {
661ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *a->Q, *x->emat);
662d6223691SXuan Zhou     }
663fc54b460SXuan Zhou     break;
664d71ae5a4SJacob Faibussowitsch   case MAT_FACTOR_CHOLESKY:
665d71ae5a4SJacob Faibussowitsch     El::cholesky::SolveAfter(El::UPPER, El::NORMAL, *a->emat, *x->emat);
666d71ae5a4SJacob Faibussowitsch     break;
667d71ae5a4SJacob Faibussowitsch   default:
668d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unfactored Matrix or Unsupported MatFactorType");
669d71ae5a4SJacob Faibussowitsch     break;
670fc54b460SXuan Zhou   }
67142ba530cSPierre Jolivet   if (!flg) {
6729566063dSJacob Faibussowitsch     PetscCall(MatConvert(C, type, MAT_REUSE_MATRIX, &X));
6739566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
67442ba530cSPierre Jolivet   }
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
676ae844d54SHong Zhang }
677ae844d54SHong Zhang 
678cedd07caSPierre Jolivet static PetscErrorCode MatLUFactor_Elemental(Mat A, IS, IS, const MatFactorInfo *)
679d71ae5a4SJacob Faibussowitsch {
6807c920d81SXuan Zhou   Mat_Elemental *a        = (Mat_Elemental *)A->data;
681ea394475SSatish Balay   PetscInt       pivoting = a->pivoting;
6827c920d81SXuan Zhou 
683ae844d54SHong Zhang   PetscFunctionBegin;
684ea394475SSatish Balay   if (pivoting == 0) {
685e8146892SSatish Balay     El::LU(*a->emat);
686ea394475SSatish Balay   } else if (pivoting == 1) {
687ea394475SSatish Balay     El::LU(*a->emat, *a->P);
688ea394475SSatish Balay   } else {
689ea394475SSatish Balay     El::LU(*a->emat, *a->P, *a->Q);
690d6223691SXuan Zhou   }
6911293973dSHong Zhang   A->factortype = MAT_FACTOR_LU;
692834d3fecSHong Zhang   A->assembled  = PETSC_TRUE;
693f6224b95SHong Zhang 
6949566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->solvertype));
6959566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &A->solvertype));
6963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
697ae844d54SHong Zhang }
698ae844d54SHong Zhang 
699d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_Elemental(Mat F, Mat A, const MatFactorInfo *info)
700d71ae5a4SJacob Faibussowitsch {
701d7c3f9d8SHong Zhang   PetscFunctionBegin;
7029566063dSJacob Faibussowitsch   PetscCall(MatCopy(A, F, SAME_NONZERO_PATTERN));
703f22e26b7SPierre Jolivet   PetscCall(MatLUFactor_Elemental(F, nullptr, nullptr, info));
7043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
705d7c3f9d8SHong Zhang }
706d7c3f9d8SHong Zhang 
707cedd07caSPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_Elemental(Mat, Mat, IS, IS, const MatFactorInfo *)
708d71ae5a4SJacob Faibussowitsch {
709d7c3f9d8SHong Zhang   PetscFunctionBegin;
710d8304050SJose E. Roman   /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
7113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
712d7c3f9d8SHong Zhang }
713d7c3f9d8SHong Zhang 
714cedd07caSPierre Jolivet static PetscErrorCode MatCholeskyFactor_Elemental(Mat A, IS, const MatFactorInfo *)
715d71ae5a4SJacob Faibussowitsch {
71645cf121fSXuan Zhou   Mat_Elemental                                    *a = (Mat_Elemental *)A->data;
717e8146892SSatish Balay   El::DistMatrix<PetscElemScalar, El::MC, El::STAR> d;
71845cf121fSXuan Zhou 
71945cf121fSXuan Zhou   PetscFunctionBegin;
720e8146892SSatish Balay   El::Cholesky(El::UPPER, *a->emat);
721fc54b460SXuan Zhou   A->factortype = MAT_FACTOR_CHOLESKY;
722834d3fecSHong Zhang   A->assembled  = PETSC_TRUE;
723f6224b95SHong Zhang 
7249566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->solvertype));
7259566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &A->solvertype));
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72745cf121fSXuan Zhou }
72845cf121fSXuan Zhou 
729d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F, Mat A, const MatFactorInfo *info)
730d71ae5a4SJacob Faibussowitsch {
731cb76c1d8SXuan Zhou   PetscFunctionBegin;
7329566063dSJacob Faibussowitsch   PetscCall(MatCopy(A, F, SAME_NONZERO_PATTERN));
733f22e26b7SPierre Jolivet   PetscCall(MatCholeskyFactor_Elemental(F, nullptr, info));
7343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73579673f7bSHong Zhang }
736cb76c1d8SXuan Zhou 
737cedd07caSPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat, Mat, IS, const MatFactorInfo *)
738d71ae5a4SJacob Faibussowitsch {
73979673f7bSHong Zhang   PetscFunctionBegin;
740d8304050SJose E. Roman   /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
7413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74279673f7bSHong Zhang }
74379673f7bSHong Zhang 
74466976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_elemental_elemental(Mat, MatSolverType *type)
745d71ae5a4SJacob Faibussowitsch {
74615767789SHong Zhang   PetscFunctionBegin;
74715767789SHong Zhang   *type = MATSOLVERELEMENTAL;
7483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74915767789SHong Zhang }
75015767789SHong Zhang 
751d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_elemental_elemental(Mat A, MatFactorType ftype, Mat *F)
752d71ae5a4SJacob Faibussowitsch {
7531293973dSHong Zhang   Mat B;
7541293973dSHong Zhang 
7551293973dSHong Zhang   PetscFunctionBegin;
7561293973dSHong Zhang   /* Create the factorization matrix */
7579566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
7589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
7599566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATELEMENTAL));
7609566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
7611293973dSHong Zhang   B->factortype      = ftype;
76266e17bc3SBarry Smith   B->trivialsymbolic = PETSC_TRUE;
7639566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
7649566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &B->solvertype));
76500c67f3bSHong Zhang 
7669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_elemental_elemental));
7671293973dSHong Zhang   *F = B;
7683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7691293973dSHong Zhang }
7701293973dSHong Zhang 
771d1f0640dSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Elemental(void)
772d71ae5a4SJacob Faibussowitsch {
77342c9c57cSBarry Smith   PetscFunctionBegin;
7749566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL, MATELEMENTAL, MAT_FACTOR_LU, MatGetFactor_elemental_elemental));
7759566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL, MATELEMENTAL, MAT_FACTOR_CHOLESKY, MatGetFactor_elemental_elemental));
7763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
77742c9c57cSBarry Smith }
77842c9c57cSBarry Smith 
779d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNorm_Elemental(Mat A, NormType type, PetscReal *nrm)
780d71ae5a4SJacob Faibussowitsch {
7811f881ff8SXuan Zhou   Mat_Elemental *a = (Mat_Elemental *)A->data;
7821f881ff8SXuan Zhou 
7831f881ff8SXuan Zhou   PetscFunctionBegin;
7841f881ff8SXuan Zhou   switch (type) {
785d71ae5a4SJacob Faibussowitsch   case NORM_1:
786d71ae5a4SJacob Faibussowitsch     *nrm = El::OneNorm(*a->emat);
787d71ae5a4SJacob Faibussowitsch     break;
788d71ae5a4SJacob Faibussowitsch   case NORM_FROBENIUS:
789d71ae5a4SJacob Faibussowitsch     *nrm = El::FrobeniusNorm(*a->emat);
790d71ae5a4SJacob Faibussowitsch     break;
791d71ae5a4SJacob Faibussowitsch   case NORM_INFINITY:
792d71ae5a4SJacob Faibussowitsch     *nrm = El::InfinityNorm(*a->emat);
793d71ae5a4SJacob Faibussowitsch     break;
794d71ae5a4SJacob Faibussowitsch   default:
795d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm type");
7961f881ff8SXuan Zhou   }
7973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7981f881ff8SXuan Zhou }
7991f881ff8SXuan Zhou 
800d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_Elemental(Mat A)
801d71ae5a4SJacob Faibussowitsch {
8025262d616SXuan Zhou   Mat_Elemental *a = (Mat_Elemental *)A->data;
8035262d616SXuan Zhou 
8045262d616SXuan Zhou   PetscFunctionBegin;
805e8146892SSatish Balay   El::Zero(*a->emat);
8063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8075262d616SXuan Zhou }
8085262d616SXuan Zhou 
809d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A, IS *rows, IS *cols)
810d71ae5a4SJacob Faibussowitsch {
811db31f6deSJed Brown   Mat_Elemental *a = (Mat_Elemental *)A->data;
812db31f6deSJed Brown   PetscInt       i, m, shift, stride, *idx;
813db31f6deSJed Brown 
814db31f6deSJed Brown   PetscFunctionBegin;
815db31f6deSJed Brown   if (rows) {
816db31f6deSJed Brown     m      = a->emat->LocalHeight();
817db31f6deSJed Brown     shift  = a->emat->ColShift();
818db31f6deSJed Brown     stride = a->emat->ColStride();
8199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &idx));
820db31f6deSJed Brown     for (i = 0; i < m; i++) {
821db31f6deSJed Brown       PetscInt rank, offset;
822db31f6deSJed Brown       E2RO(A, 0, shift + i * stride, &rank, &offset);
823db31f6deSJed Brown       RO2P(A, 0, rank, offset, &idx[i]);
824db31f6deSJed Brown     }
8259566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_OWN_POINTER, rows));
826db31f6deSJed Brown   }
827db31f6deSJed Brown   if (cols) {
828db31f6deSJed Brown     m      = a->emat->LocalWidth();
829db31f6deSJed Brown     shift  = a->emat->RowShift();
830db31f6deSJed Brown     stride = a->emat->RowStride();
8319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &idx));
832db31f6deSJed Brown     for (i = 0; i < m; i++) {
833db31f6deSJed Brown       PetscInt rank, offset;
834db31f6deSJed Brown       E2RO(A, 1, shift + i * stride, &rank, &offset);
835db31f6deSJed Brown       RO2P(A, 1, rank, offset, &idx[i]);
836db31f6deSJed Brown     }
8379566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_OWN_POINTER, cols));
838db31f6deSJed Brown   }
8393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
840db31f6deSJed Brown }
841db31f6deSJed Brown 
842cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Elemental_Dense(Mat A, MatType, MatReuse reuse, Mat *B)
843d71ae5a4SJacob Faibussowitsch {
8442ef0cf24SXuan Zhou   Mat             Bmpi;
845af295397SXuan Zhou   Mat_Elemental  *a = (Mat_Elemental *)A->data;
846ce94432eSBarry Smith   MPI_Comm        comm;
847fa9e26c8SHong Zhang   IS              isrows, iscols;
848fa9e26c8SHong Zhang   PetscInt        rrank, ridx, crank, cidx, nrows, ncols, i, j, erow, ecol, elrow, elcol;
849fa9e26c8SHong Zhang   const PetscInt *rows, *cols;
850df311e6cSXuan Zhou   PetscElemScalar v;
851fa9e26c8SHong Zhang   const El::Grid &grid = a->emat->Grid();
852af295397SXuan Zhou 
853af295397SXuan Zhou   PetscFunctionBegin;
8549566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
855a000e53bSHong Zhang 
856de0a22f0SHong Zhang   if (reuse == MAT_REUSE_MATRIX) {
857de0a22f0SHong Zhang     Bmpi = *B;
858de0a22f0SHong Zhang   } else {
8599566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &Bmpi));
8609566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Bmpi, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
8619566063dSJacob Faibussowitsch     PetscCall(MatSetType(Bmpi, MATDENSE));
8629566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Bmpi));
863de0a22f0SHong Zhang   }
864fa9e26c8SHong Zhang 
865fa9e26c8SHong Zhang   /* Get local entries of A */
8669566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(A, &isrows, &iscols));
8679566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows, &nrows));
8689566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows, &rows));
8699566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscols, &ncols));
8709566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscols, &cols));
871fa9e26c8SHong Zhang 
87257299f2fSHong Zhang   if (a->roworiented) {
8732ef0cf24SXuan Zhou     for (i = 0; i < nrows; i++) {
874fa9e26c8SHong Zhang       P2RO(A, 0, rows[i], &rrank, &ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
8752ef0cf24SXuan Zhou       RO2E(A, 0, rrank, ridx, &erow);
876aed4548fSBarry Smith       PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation");
8772ef0cf24SXuan Zhou       for (j = 0; j < ncols; j++) {
878fa9e26c8SHong Zhang         P2RO(A, 1, cols[j], &crank, &cidx);
8792ef0cf24SXuan Zhou         RO2E(A, 1, crank, cidx, &ecol);
880aed4548fSBarry Smith         PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation");
881fa9e26c8SHong Zhang 
882fa9e26c8SHong Zhang         elrow = erow / grid.MCSize(); /* Elemental local row index */
883fa9e26c8SHong Zhang         elcol = ecol / grid.MRSize(); /* Elemental local column index */
884fa9e26c8SHong Zhang         v     = a->emat->GetLocal(elrow, elcol);
8859566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Bmpi, 1, &rows[i], 1, &cols[j], (PetscScalar *)&v, INSERT_VALUES));
8862ef0cf24SXuan Zhou       }
8872ef0cf24SXuan Zhou     }
88857299f2fSHong Zhang   } else { /* column-oriented */
88957299f2fSHong Zhang     for (j = 0; j < ncols; j++) {
89057299f2fSHong Zhang       P2RO(A, 1, cols[j], &crank, &cidx);
89157299f2fSHong Zhang       RO2E(A, 1, crank, cidx, &ecol);
892aed4548fSBarry Smith       PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation");
89357299f2fSHong Zhang       for (i = 0; i < nrows; i++) {
89457299f2fSHong Zhang         P2RO(A, 0, rows[i], &rrank, &ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
89557299f2fSHong Zhang         RO2E(A, 0, rrank, ridx, &erow);
896aed4548fSBarry Smith         PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation");
89757299f2fSHong Zhang 
89857299f2fSHong Zhang         elrow = erow / grid.MCSize(); /* Elemental local row index */
89957299f2fSHong Zhang         elcol = ecol / grid.MRSize(); /* Elemental local column index */
90057299f2fSHong Zhang         v     = a->emat->GetLocal(elrow, elcol);
9019566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Bmpi, 1, &rows[i], 1, &cols[j], (PetscScalar *)&v, INSERT_VALUES));
90257299f2fSHong Zhang       }
90357299f2fSHong Zhang     }
90457299f2fSHong Zhang   }
9059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Bmpi, MAT_FINAL_ASSEMBLY));
9069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Bmpi, MAT_FINAL_ASSEMBLY));
907511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
9089566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &Bmpi));
909c4ad791aSXuan Zhou   } else {
910c4ad791aSXuan Zhou     *B = Bmpi;
911c4ad791aSXuan Zhou   }
9129566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
9139566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscols));
9143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
915af295397SXuan Zhou }
916af295397SXuan Zhou 
917cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
918d71ae5a4SJacob Faibussowitsch {
919af8000cdSHong Zhang   Mat                mat_elemental;
920af8000cdSHong Zhang   PetscInt           M = A->rmap->N, N = A->cmap->N, row, ncols;
921af8000cdSHong Zhang   const PetscInt    *cols;
922af8000cdSHong Zhang   const PetscScalar *vals;
923af8000cdSHong Zhang 
924af8000cdSHong Zhang   PetscFunctionBegin;
925bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
926bdeae278SStefano Zampini     mat_elemental = *newmat;
9279566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
928bdeae278SStefano Zampini   } else {
9299566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
9309566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N));
9319566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
9329566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
933bdeae278SStefano Zampini   }
934af8000cdSHong Zhang   for (row = 0; row < M; row++) {
9359566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
936bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
9379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES));
9389566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
939af8000cdSHong Zhang   }
9409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
9419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
942af8000cdSHong Zhang 
943511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
9449566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &mat_elemental));
945af8000cdSHong Zhang   } else {
946af8000cdSHong Zhang     *newmat = mat_elemental;
947af8000cdSHong Zhang   }
9483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
949af8000cdSHong Zhang }
950af8000cdSHong Zhang 
951cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
952d71ae5a4SJacob Faibussowitsch {
9535d7652ecSHong Zhang   Mat                mat_elemental;
9545d7652ecSHong Zhang   PetscInt           row, ncols, rstart = A->rmap->rstart, rend = A->rmap->rend, j;
9555d7652ecSHong Zhang   const PetscInt    *cols;
9565d7652ecSHong Zhang   const PetscScalar *vals;
9575d7652ecSHong Zhang 
9585d7652ecSHong Zhang   PetscFunctionBegin;
959bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
960bdeae278SStefano Zampini     mat_elemental = *newmat;
9619566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
962bdeae278SStefano Zampini   } else {
9639566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
9649566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, A->rmap->N, A->cmap->N));
9659566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
9669566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
967bdeae278SStefano Zampini   }
9685d7652ecSHong Zhang   for (row = rstart; row < rend; row++) {
9699566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
9705d7652ecSHong Zhang     for (j = 0; j < ncols; j++) {
971bdeae278SStefano Zampini       /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
9729566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat_elemental, 1, &row, 1, &cols[j], &vals[j], ADD_VALUES));
9735d7652ecSHong Zhang     }
9749566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
9755d7652ecSHong Zhang   }
9769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
9779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
9785d7652ecSHong Zhang 
979511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
9809566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &mat_elemental));
9815d7652ecSHong Zhang   } else {
9825d7652ecSHong Zhang     *newmat = mat_elemental;
9835d7652ecSHong Zhang   }
9843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9855d7652ecSHong Zhang }
9865d7652ecSHong Zhang 
987cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
988d71ae5a4SJacob Faibussowitsch {
9896214f412SHong Zhang   Mat                mat_elemental;
9906214f412SHong Zhang   PetscInt           M = A->rmap->N, N = A->cmap->N, row, ncols, j;
9916214f412SHong Zhang   const PetscInt    *cols;
9926214f412SHong Zhang   const PetscScalar *vals;
9936214f412SHong Zhang 
9946214f412SHong Zhang   PetscFunctionBegin;
995bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
996bdeae278SStefano Zampini     mat_elemental = *newmat;
9979566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
998bdeae278SStefano Zampini   } else {
9999566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
10009566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N));
10019566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
10029566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
1003bdeae278SStefano Zampini   }
10049566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(A));
10056214f412SHong Zhang   for (row = 0; row < M; row++) {
10069566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
1007bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10089566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES));
10096214f412SHong Zhang     for (j = 0; j < ncols; j++) { /* lower triangular part */
1010eb1ec7c1SStefano Zampini       PetscScalar v;
10116214f412SHong Zhang       if (cols[j] == row) continue;
1012b94d7dedSBarry Smith       v = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j];
10139566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat_elemental, 1, &cols[j], 1, &row, &v, ADD_VALUES));
10146214f412SHong Zhang     }
10159566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
10166214f412SHong Zhang   }
10179566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(A));
10189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
10199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
10206214f412SHong Zhang 
1021511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
10229566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &mat_elemental));
10236214f412SHong Zhang   } else {
10246214f412SHong Zhang     *newmat = mat_elemental;
10256214f412SHong Zhang   }
10263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10276214f412SHong Zhang }
10286214f412SHong Zhang 
1029cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
1030d71ae5a4SJacob Faibussowitsch {
10316214f412SHong Zhang   Mat                mat_elemental;
10326214f412SHong Zhang   PetscInt           M = A->rmap->N, N = A->cmap->N, row, ncols, j, rstart = A->rmap->rstart, rend = A->rmap->rend;
10336214f412SHong Zhang   const PetscInt    *cols;
10346214f412SHong Zhang   const PetscScalar *vals;
10356214f412SHong Zhang 
10366214f412SHong Zhang   PetscFunctionBegin;
1037bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
1038bdeae278SStefano Zampini     mat_elemental = *newmat;
10399566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(mat_elemental));
1040bdeae278SStefano Zampini   } else {
10419566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
10429566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N));
10439566063dSJacob Faibussowitsch     PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
10449566063dSJacob Faibussowitsch     PetscCall(MatSetUp(mat_elemental));
1045bdeae278SStefano Zampini   }
10469566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(A));
10476214f412SHong Zhang   for (row = rstart; row < rend; row++) {
10489566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
1049bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10509566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES));
10516214f412SHong Zhang     for (j = 0; j < ncols; j++) { /* lower triangular part */
1052eb1ec7c1SStefano Zampini       PetscScalar v;
10536214f412SHong Zhang       if (cols[j] == row) continue;
1054b94d7dedSBarry Smith       v = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j];
10559566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat_elemental, 1, &cols[j], 1, &row, &v, ADD_VALUES));
10566214f412SHong Zhang     }
10579566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
10586214f412SHong Zhang   }
10599566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(A));
10609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
10619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
10626214f412SHong Zhang 
1063511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
10649566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &mat_elemental));
10656214f412SHong Zhang   } else {
10666214f412SHong Zhang     *newmat = mat_elemental;
10676214f412SHong Zhang   }
10683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10696214f412SHong Zhang }
10706214f412SHong Zhang 
1071d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Elemental(Mat A)
1072d71ae5a4SJacob Faibussowitsch {
1073db31f6deSJed Brown   Mat_Elemental      *a = (Mat_Elemental *)A->data;
10745e9f5b67SHong Zhang   Mat_Elemental_Grid *commgrid;
10755e9f5b67SHong Zhang   PetscBool           flg;
10765e9f5b67SHong Zhang   MPI_Comm            icomm;
1077db31f6deSJed Brown 
1078db31f6deSJed Brown   PetscFunctionBegin;
1079db31f6deSJed Brown   delete a->emat;
1080ea394475SSatish Balay   delete a->P;
1081ea394475SSatish Balay   delete a->Q;
10825e9f5b67SHong Zhang 
1083e8146892SSatish Balay   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
1084f22e26b7SPierre Jolivet   PetscCall(PetscCommDuplicate(cxxcomm.comm, &icomm, nullptr));
10859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Elemental_keyval, (void **)&commgrid, (int *)&flg));
10865e9f5b67SHong Zhang   if (--commgrid->grid_refct == 0) {
10875e9f5b67SHong Zhang     delete commgrid->grid;
10889566063dSJacob Faibussowitsch     PetscCall(PetscFree(commgrid));
10899566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Elemental_keyval));
10905e9f5b67SHong Zhang   }
10919566063dSJacob Faibussowitsch   PetscCall(PetscCommDestroy(&icomm));
1092f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetOwnershipIS_C", nullptr));
1093f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", nullptr));
1094f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_elemental_mpidense_C", nullptr));
10959566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
10963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1097db31f6deSJed Brown }
1098db31f6deSJed Brown 
109966976f2fSJacob Faibussowitsch static PetscErrorCode MatSetUp_Elemental(Mat A)
1100d71ae5a4SJacob Faibussowitsch {
1101db31f6deSJed Brown   Mat_Elemental *a = (Mat_Elemental *)A->data;
1102f19600a0SHong Zhang   MPI_Comm       comm;
1103db31f6deSJed Brown   PetscMPIInt    rsize, csize;
1104f19600a0SHong Zhang   PetscInt       n;
1105db31f6deSJed Brown 
1106db31f6deSJed Brown   PetscFunctionBegin;
11079566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
11089566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
1109db31f6deSJed Brown 
1110d8304050SJose E. Roman   /* Check if local row and column sizes are equally distributed.
1111f19600a0SHong Zhang      Jed: Elemental uses "element" cyclic ordering so the sizes need to match that
1112f19600a0SHong Zhang      exactly.  The strategy in MatElemental is for PETSc to implicitly permute to block ordering (like would be returned by
1113d8304050SJose E. Roman      PetscSplitOwnership(comm,&n,&N), at which point Elemental matrices can act on PETSc vectors without redistributing the vectors. */
11149566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1115f19600a0SHong Zhang   n = PETSC_DECIDE;
11169566063dSJacob Faibussowitsch   PetscCall(PetscSplitOwnership(comm, &n, &A->rmap->N));
11175f80ce2aSJacob 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);
1118f19600a0SHong Zhang 
1119f19600a0SHong Zhang   n = PETSC_DECIDE;
11209566063dSJacob Faibussowitsch   PetscCall(PetscSplitOwnership(comm, &n, &A->cmap->N));
11215f80ce2aSJacob 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);
1122f19600a0SHong Zhang 
11232f613bf5SBarry Smith   a->emat->Resize(A->rmap->N, A->cmap->N);
1124e8146892SSatish Balay   El::Zero(*a->emat);
1125db31f6deSJed Brown 
11269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(A->rmap->comm, &rsize));
11279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(A->cmap->comm, &csize));
11285f80ce2aSJacob Faibussowitsch   PetscCheck(csize == rsize, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Cannot use row and column communicators of different sizes");
1129db31f6deSJed Brown   a->commsize = rsize;
11309371c9d4SSatish Balay   a->mr[0]    = A->rmap->N % rsize;
11319371c9d4SSatish Balay   if (!a->mr[0]) a->mr[0] = rsize;
11329371c9d4SSatish Balay   a->mr[1] = A->cmap->N % csize;
11339371c9d4SSatish Balay   if (!a->mr[1]) a->mr[1] = csize;
1134db31f6deSJed Brown   a->m[0] = A->rmap->N / rsize + (a->mr[0] != rsize);
1135db31f6deSJed Brown   a->m[1] = A->cmap->N / csize + (a->mr[1] != csize);
11363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1137db31f6deSJed Brown }
1138db31f6deSJed Brown 
113966976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType)
1140d71ae5a4SJacob Faibussowitsch {
1141aae2c449SHong Zhang   Mat_Elemental *a = (Mat_Elemental *)A->data;
1142aae2c449SHong Zhang 
1143aae2c449SHong Zhang   PetscFunctionBegin;
1144fbbcb730SJack Poulson   /* printf("Calling ProcessQueues\n"); */
1145fbbcb730SJack Poulson   a->emat->ProcessQueues();
1146fbbcb730SJack Poulson   /* printf("Finished ProcessQueues\n"); */
11473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1148aae2c449SHong Zhang }
1149aae2c449SHong Zhang 
115066976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_Elemental(Mat, MatAssemblyType)
1151d71ae5a4SJacob Faibussowitsch {
1152aae2c449SHong Zhang   PetscFunctionBegin;
1153aae2c449SHong Zhang   /* Currently does nothing */
11543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1155aae2c449SHong Zhang }
1156aae2c449SHong Zhang 
115766976f2fSJacob Faibussowitsch static PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer)
1158d71ae5a4SJacob Faibussowitsch {
11597b30e0f1SHong Zhang   Mat      Adense, Ae;
11607b30e0f1SHong Zhang   MPI_Comm comm;
11617b30e0f1SHong Zhang 
11627b30e0f1SHong Zhang   PetscFunctionBegin;
11639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)newMat, &comm));
11649566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Adense));
11659566063dSJacob Faibussowitsch   PetscCall(MatSetType(Adense, MATDENSE));
11669566063dSJacob Faibussowitsch   PetscCall(MatLoad(Adense, viewer));
11679566063dSJacob Faibussowitsch   PetscCall(MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae));
11689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Adense));
11699566063dSJacob Faibussowitsch   PetscCall(MatHeaderReplace(newMat, &Ae));
11703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11717b30e0f1SHong Zhang }
11727b30e0f1SHong Zhang 
11739371c9d4SSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_Elemental,
1174f22e26b7SPierre Jolivet                                        nullptr,
1175f22e26b7SPierre Jolivet                                        nullptr,
117640d92e34SHong Zhang                                        MatMult_Elemental,
117740d92e34SHong Zhang                                        /* 4*/ MatMultAdd_Elemental,
11789426833fSXuan Zhou                                        MatMultTranspose_Elemental,
1179e883f9d5SXuan Zhou                                        MatMultTransposeAdd_Elemental,
118040d92e34SHong Zhang                                        MatSolve_Elemental,
1181df311e6cSXuan Zhou                                        MatSolveAdd_Elemental,
1182f22e26b7SPierre Jolivet                                        nullptr,
1183f22e26b7SPierre Jolivet                                        /*10*/ nullptr,
118440d92e34SHong Zhang                                        MatLUFactor_Elemental,
118540d92e34SHong Zhang                                        MatCholeskyFactor_Elemental,
1186f22e26b7SPierre Jolivet                                        nullptr,
118740d92e34SHong Zhang                                        MatTranspose_Elemental,
118840d92e34SHong Zhang                                        /*15*/ MatGetInfo_Elemental,
1189f22e26b7SPierre Jolivet                                        nullptr,
119061119200SXuan Zhou                                        MatGetDiagonal_Elemental,
1191ade3cc5eSXuan Zhou                                        MatDiagonalScale_Elemental,
119240d92e34SHong Zhang                                        MatNorm_Elemental,
119340d92e34SHong Zhang                                        /*20*/ MatAssemblyBegin_Elemental,
119440d92e34SHong Zhang                                        MatAssemblyEnd_Elemental,
119532d7a744SHong Zhang                                        MatSetOption_Elemental,
119640d92e34SHong Zhang                                        MatZeroEntries_Elemental,
1197f22e26b7SPierre Jolivet                                        /*24*/ nullptr,
119840d92e34SHong Zhang                                        MatLUFactorSymbolic_Elemental,
119940d92e34SHong Zhang                                        MatLUFactorNumeric_Elemental,
120040d92e34SHong Zhang                                        MatCholeskyFactorSymbolic_Elemental,
120140d92e34SHong Zhang                                        MatCholeskyFactorNumeric_Elemental,
120240d92e34SHong Zhang                                        /*29*/ MatSetUp_Elemental,
1203f22e26b7SPierre Jolivet                                        nullptr,
1204f22e26b7SPierre Jolivet                                        nullptr,
1205f22e26b7SPierre Jolivet                                        nullptr,
1206f22e26b7SPierre Jolivet                                        nullptr,
1207df311e6cSXuan Zhou                                        /*34*/ MatDuplicate_Elemental,
1208f22e26b7SPierre Jolivet                                        nullptr,
1209f22e26b7SPierre Jolivet                                        nullptr,
1210f22e26b7SPierre Jolivet                                        nullptr,
1211f22e26b7SPierre Jolivet                                        nullptr,
121240d92e34SHong Zhang                                        /*39*/ MatAXPY_Elemental,
1213f22e26b7SPierre Jolivet                                        nullptr,
1214f22e26b7SPierre Jolivet                                        nullptr,
1215f22e26b7SPierre Jolivet                                        nullptr,
121640d92e34SHong Zhang                                        MatCopy_Elemental,
1217f22e26b7SPierre Jolivet                                        /*44*/ nullptr,
121840d92e34SHong Zhang                                        MatScale_Elemental,
12197d68702bSBarry Smith                                        MatShift_Basic,
1220f22e26b7SPierre Jolivet                                        nullptr,
1221f22e26b7SPierre Jolivet                                        nullptr,
1222f22e26b7SPierre Jolivet                                        /*49*/ nullptr,
1223f22e26b7SPierre Jolivet                                        nullptr,
1224f22e26b7SPierre Jolivet                                        nullptr,
1225f22e26b7SPierre Jolivet                                        nullptr,
1226f22e26b7SPierre Jolivet                                        nullptr,
1227f22e26b7SPierre Jolivet                                        /*54*/ nullptr,
1228f22e26b7SPierre Jolivet                                        nullptr,
1229f22e26b7SPierre Jolivet                                        nullptr,
1230f22e26b7SPierre Jolivet                                        nullptr,
1231f22e26b7SPierre Jolivet                                        nullptr,
1232f22e26b7SPierre Jolivet                                        /*59*/ nullptr,
123340d92e34SHong Zhang                                        MatDestroy_Elemental,
123440d92e34SHong Zhang                                        MatView_Elemental,
1235f22e26b7SPierre Jolivet                                        nullptr,
1236f22e26b7SPierre Jolivet                                        nullptr,
1237f22e26b7SPierre Jolivet                                        /*64*/ nullptr,
1238f22e26b7SPierre Jolivet                                        nullptr,
1239f22e26b7SPierre Jolivet                                        nullptr,
1240f22e26b7SPierre Jolivet                                        nullptr,
1241f22e26b7SPierre Jolivet                                        nullptr,
1242f22e26b7SPierre Jolivet                                        /*69*/ nullptr,
12432ef0cf24SXuan Zhou                                        MatConvert_Elemental_Dense,
1244f22e26b7SPierre Jolivet                                        nullptr,
1245f22e26b7SPierre Jolivet                                        nullptr,
12468bb0f5c6SPierre Jolivet                                        nullptr,
1247f22e26b7SPierre Jolivet                                        /*74*/ nullptr,
1248f22e26b7SPierre Jolivet                                        nullptr,
1249f22e26b7SPierre Jolivet                                        nullptr,
1250f22e26b7SPierre Jolivet                                        nullptr,
12518bb0f5c6SPierre Jolivet                                        MatLoad_Elemental,
1252f22e26b7SPierre Jolivet                                        /*79*/ nullptr,
1253f22e26b7SPierre Jolivet                                        nullptr,
1254f22e26b7SPierre Jolivet                                        nullptr,
1255f22e26b7SPierre Jolivet                                        nullptr,
12568bb0f5c6SPierre Jolivet                                        nullptr,
1257f22e26b7SPierre Jolivet                                        /*84*/ nullptr,
125840d92e34SHong Zhang                                        MatMatMultNumeric_Elemental,
1259f22e26b7SPierre Jolivet                                        nullptr,
1260f22e26b7SPierre Jolivet                                        nullptr,
1261df311e6cSXuan Zhou                                        MatMatTransposeMultNumeric_Elemental,
12628bb0f5c6SPierre Jolivet                                        /*89*/ nullptr,
12638bb0f5c6SPierre Jolivet                                        MatProductSetFromOptions_Elemental,
1264f22e26b7SPierre Jolivet                                        nullptr,
1265f22e26b7SPierre Jolivet                                        nullptr,
1266dfcb0403SXuan Zhou                                        MatConjugate_Elemental,
12678bb0f5c6SPierre Jolivet                                        /*94*/ nullptr,
1268f22e26b7SPierre Jolivet                                        nullptr,
1269f22e26b7SPierre Jolivet                                        nullptr,
1270f22e26b7SPierre Jolivet                                        nullptr,
1271f22e26b7SPierre Jolivet                                        nullptr,
12728bb0f5c6SPierre Jolivet                                        /*99*/ nullptr,
12738bb0f5c6SPierre Jolivet                                        MatMatSolve_Elemental,
1274f22e26b7SPierre Jolivet                                        nullptr,
1275f22e26b7SPierre Jolivet                                        nullptr,
1276f22e26b7SPierre Jolivet                                        nullptr,
12778bb0f5c6SPierre Jolivet                                        /*104*/ MatMissingDiagonal_Elemental,
12788bb0f5c6SPierre Jolivet                                        nullptr,
12798bb0f5c6SPierre Jolivet                                        nullptr,
12808bb0f5c6SPierre Jolivet                                        nullptr,
12818bb0f5c6SPierre Jolivet                                        nullptr,
12828bb0f5c6SPierre Jolivet                                        /*109*/ nullptr,
12838bb0f5c6SPierre Jolivet                                        nullptr,
12848bb0f5c6SPierre Jolivet                                        MatHermitianTranspose_Elemental,
12858bb0f5c6SPierre Jolivet                                        nullptr,
12868bb0f5c6SPierre Jolivet                                        nullptr,
1287f22e26b7SPierre Jolivet                                        /*114*/ nullptr,
1288f22e26b7SPierre Jolivet                                        nullptr,
1289f22e26b7SPierre Jolivet                                        nullptr,
1290f22e26b7SPierre Jolivet                                        nullptr,
1291f22e26b7SPierre Jolivet                                        nullptr,
1292f22e26b7SPierre Jolivet                                        /*119*/ nullptr,
12938bb0f5c6SPierre Jolivet                                        nullptr,
1294f22e26b7SPierre Jolivet                                        nullptr,
1295f22e26b7SPierre Jolivet                                        nullptr,
1296f22e26b7SPierre Jolivet                                        nullptr,
1297f22e26b7SPierre Jolivet                                        /*124*/ nullptr,
1298f22e26b7SPierre Jolivet                                        nullptr,
1299f22e26b7SPierre Jolivet                                        nullptr,
1300f22e26b7SPierre Jolivet                                        nullptr,
1301f22e26b7SPierre Jolivet                                        nullptr,
1302f22e26b7SPierre Jolivet                                        /*129*/ nullptr,
1303f22e26b7SPierre Jolivet                                        nullptr,
1304f22e26b7SPierre Jolivet                                        nullptr,
1305f22e26b7SPierre Jolivet                                        nullptr,
1306f22e26b7SPierre Jolivet                                        nullptr,
1307f22e26b7SPierre Jolivet                                        /*134*/ nullptr,
1308f22e26b7SPierre Jolivet                                        nullptr,
1309f22e26b7SPierre Jolivet                                        nullptr,
1310f22e26b7SPierre Jolivet                                        nullptr,
1311f22e26b7SPierre Jolivet                                        nullptr,
1312f22e26b7SPierre Jolivet                                        nullptr,
1313f22e26b7SPierre Jolivet                                        /*140*/ nullptr,
1314f22e26b7SPierre Jolivet                                        nullptr,
1315*c97e14ecSAlex Lindsay                                        nullptr,
1316f22e26b7SPierre Jolivet                                        nullptr};
131740d92e34SHong Zhang 
1318ed36708cSHong Zhang /*MC
1319ed36708cSHong Zhang    MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package
1320ed36708cSHong Zhang 
1321c2b89b5dSBarry Smith   Use ./configure --download-elemental to install PETSc to use Elemental
1322c2b89b5dSBarry Smith 
1323ed36708cSHong Zhang    Options Database Keys:
13245cc86fc1SJed Brown + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions()
132589bba20eSBarry Smith . -pc_factor_mat_solver_type elemental - to use this direct solver with the option -pc_type lu
13265cc86fc1SJed Brown - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix
1327ed36708cSHong Zhang 
1328ed36708cSHong Zhang   Level: beginner
1329ed36708cSHong Zhang 
133089bba20eSBarry Smith   Note:
133189bba20eSBarry Smith    Note unlike most matrix formats, this format does not store all the matrix entries for a contiguous
133289bba20eSBarry Smith    range of rows on an MPI rank. Use `MatGetOwnershipIS()` to determine what values are stored on
133389bba20eSBarry Smith    the given rank.
133489bba20eSBarry Smith 
133589bba20eSBarry Smith .seealso: `MATDENSE`, `MATSCALAPACK`, `MatGetOwnershipIS()`
1336ed36708cSHong Zhang M*/
1337f22e26b7SPierre Jolivet #if defined(__clang__)
1338f22e26b7SPierre Jolivet   #pragma clang diagnostic push
1339f22e26b7SPierre Jolivet   #pragma clang diagnostic ignored "-Wzero-as-null-pointer-constant"
1340f22e26b7SPierre Jolivet #endif
1341d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A)
1342d71ae5a4SJacob Faibussowitsch {
1343db31f6deSJed Brown   Mat_Elemental      *a;
13445682a260SJack Poulson   PetscBool           flg, flg1;
13455e9f5b67SHong Zhang   Mat_Elemental_Grid *commgrid;
13465e9f5b67SHong Zhang   MPI_Comm            icomm;
13475682a260SJack Poulson   PetscInt            optv1;
1348db31f6deSJed Brown 
1349db31f6deSJed Brown   PetscFunctionBegin;
1350aea10558SJacob Faibussowitsch   A->ops[0]     = MatOps_Values;
135140d92e34SHong Zhang   A->insertmode = NOT_SET_VALUES;
1352db31f6deSJed Brown 
13534dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
1354db31f6deSJed Brown   A->data = (void *)a;
1355db31f6deSJed Brown 
1356db31f6deSJed Brown   /* Set up the elemental matrix */
1357e8146892SSatish Balay   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
13585e9f5b67SHong Zhang 
13595e9f5b67SHong Zhang   /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */
13605e9f5b67SHong Zhang   if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) {
1361c8025a54SPierre Jolivet     PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Elemental_keyval, nullptr));
13629566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister(ElementalCitation, &ElementalCite));
13635e9f5b67SHong Zhang   }
13649566063dSJacob Faibussowitsch   PetscCall(PetscCommDuplicate(cxxcomm.comm, &icomm, NULL));
13659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Elemental_keyval, (void **)&commgrid, (int *)&flg));
13665e9f5b67SHong Zhang   if (!flg) {
13674dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&commgrid));
13685cb544a0SHong Zhang 
1369d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "Elemental Options", "Mat");
1370e8146892SSatish Balay     /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */
13719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-mat_elemental_grid_height", "Grid Height", "None", El::mpi::Size(cxxcomm), &optv1, &flg1));
13725682a260SJack Poulson     if (flg1) {
1373aed4548fSBarry 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));
1374e8146892SSatish Balay       commgrid->grid = new El::Grid(cxxcomm, optv1); /* use user-provided grid height */
13752ef0cf24SXuan Zhou     } else {
1376e8146892SSatish Balay       commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */
1377da0640a4SHong Zhang       /* printf("new commgrid->grid = %p\n",commgrid->grid);  -- memory leak revealed by valgrind? */
1378ed667823SXuan Zhou     }
13795e9f5b67SHong Zhang     commgrid->grid_refct = 1;
13809566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_set_attr(icomm, Petsc_Elemental_keyval, (void *)commgrid));
1381ea394475SSatish Balay 
1382ea394475SSatish Balay     a->pivoting = 1;
13839566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-mat_elemental_pivoting", "Pivoting", "None", a->pivoting, &a->pivoting, NULL));
1384ea394475SSatish Balay 
1385d0609cedSBarry Smith     PetscOptionsEnd();
13865e9f5b67SHong Zhang   } else {
13875e9f5b67SHong Zhang     commgrid->grid_refct++;
13885e9f5b67SHong Zhang   }
13899566063dSJacob Faibussowitsch   PetscCall(PetscCommDestroy(&icomm));
13905e9f5b67SHong Zhang   a->grid        = commgrid->grid;
1391e8146892SSatish Balay   a->emat        = new El::DistMatrix<PetscElemScalar>(*a->grid);
139232d7a744SHong Zhang   a->roworiented = PETSC_TRUE;
1393db31f6deSJed Brown 
13949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetOwnershipIS_C", MatGetOwnershipIS_Elemental));
13959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_elemental_mpidense_C", MatProductSetFromOptions_Elemental_MPIDense));
13969566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATELEMENTAL));
13973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1398db31f6deSJed Brown }
1399f22e26b7SPierre Jolivet #if defined(__clang__)
1400f22e26b7SPierre Jolivet   #pragma clang diagnostic pop
1401f22e26b7SPierre Jolivet #endif
1402