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
MatView_Elemental(Mat A,PetscViewer viewer)19d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Elemental(Mat A, PetscViewer viewer)
20d71ae5a4SJacob Faibussowitsch {
21db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data;
229f196a02SMartin Diehl PetscBool isascii;
23db31f6deSJed Brown
24db31f6deSJed Brown PetscFunctionBegin;
259f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
269f196a02SMartin Diehl if (isascii) {
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
MatGetInfo_Elemental(Mat A,MatInfoType flag,MatInfo * info)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
MatSetOption_Elemental(Mat A,MatOption op,PetscBool flg)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
MatSetValues_Elemental(Mat A,PetscInt nr,const PetscInt * rows,PetscInt nc,const PetscInt * cols,const PetscScalar * vals,InsertMode imode)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
MatMult_Elemental(Mat A,Vec X,Vec Y)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
MatMultTranspose_Elemental(Mat A,Vec X,Vec Y)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
MatMultAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)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
MatMultTransposeAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)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
MatMatMultNumeric_Elemental(Mat A,Mat B,Mat C)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
MatMatMultSymbolic_Elemental(Mat A,Mat B,PetscReal,Mat Ce)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
MatMatTransposeMultNumeric_Elemental(Mat A,Mat B,Mat C)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
MatMatTransposeMultSymbolic_Elemental(Mat A,Mat B,PetscReal,Mat C)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
MatProductSetFromOptions_Elemental_AB(Mat C)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
MatProductSetFromOptions_Elemental_ABt(Mat C)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
MatProductSetFromOptions_Elemental(Mat C)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
MatMatMultNumeric_Elemental_MPIDense(Mat A,Mat B,Mat C)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
MatMatMultSymbolic_Elemental_MPIDense(Mat A,Mat B,PetscReal,Mat C)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
MatProductSetFromOptions_Elemental_MPIDense_AB(Mat C)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
MatProductSetFromOptions_Elemental_MPIDense(Mat C)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
MatGetDiagonal_Elemental(Mat A,Vec D)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
MatDiagonalScale_Elemental(Mat X,Vec L,Vec R)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
MatScale_Elemental(Mat X,PetscScalar a)474d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_Elemental(Mat X, PetscScalar a)
475d71ae5a4SJacob Faibussowitsch {
47665b78793SXuan Zhou Mat_Elemental *x = (Mat_Elemental *)X->data;
47765b78793SXuan Zhou
47865b78793SXuan Zhou PetscFunctionBegin;
479e8146892SSatish Balay El::Scale((PetscElemScalar)a, *x->emat);
4803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
48165b78793SXuan Zhou }
48265b78793SXuan Zhou
483f82baa17SHong Zhang /*
484f82baa17SHong Zhang MatAXPY - Computes Y = a*X + Y.
485f82baa17SHong Zhang */
MatAXPY_Elemental(Mat Y,PetscScalar a,Mat X,MatStructure)486cedd07caSPierre Jolivet static PetscErrorCode MatAXPY_Elemental(Mat Y, PetscScalar a, Mat X, MatStructure)
487d71ae5a4SJacob Faibussowitsch {
488e09a3074SHong Zhang Mat_Elemental *x = (Mat_Elemental *)X->data;
489e09a3074SHong Zhang Mat_Elemental *y = (Mat_Elemental *)Y->data;
490e09a3074SHong Zhang
491e09a3074SHong Zhang PetscFunctionBegin;
492e8146892SSatish Balay El::Axpy((PetscElemScalar)a, *x->emat, *y->emat);
4939566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)Y));
4943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
495e09a3074SHong Zhang }
496e09a3074SHong Zhang
MatCopy_Elemental(Mat A,Mat B,MatStructure)497cedd07caSPierre Jolivet static PetscErrorCode MatCopy_Elemental(Mat A, Mat B, MatStructure)
498d71ae5a4SJacob Faibussowitsch {
499d6223691SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
500d6223691SXuan Zhou Mat_Elemental *b = (Mat_Elemental *)B->data;
501d6223691SXuan Zhou
502d6223691SXuan Zhou PetscFunctionBegin;
503e8146892SSatish Balay El::Copy(*a->emat, *b->emat);
5049566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)B));
5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
506d6223691SXuan Zhou }
507d6223691SXuan Zhou
MatDuplicate_Elemental(Mat A,MatDuplicateOption op,Mat * B)508d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_Elemental(Mat A, MatDuplicateOption op, Mat *B)
509d71ae5a4SJacob Faibussowitsch {
510df311e6cSXuan Zhou Mat Be;
511ce94432eSBarry Smith MPI_Comm comm;
512df311e6cSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
513df311e6cSXuan Zhou
514df311e6cSXuan Zhou PetscFunctionBegin;
5159566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5169566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Be));
5179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Be, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
5189566063dSJacob Faibussowitsch PetscCall(MatSetType(Be, MATELEMENTAL));
5199566063dSJacob Faibussowitsch PetscCall(MatSetUp(Be));
520df311e6cSXuan Zhou *B = Be;
521df311e6cSXuan Zhou if (op == MAT_COPY_VALUES) {
522df311e6cSXuan Zhou Mat_Elemental *b = (Mat_Elemental *)Be->data;
523e8146892SSatish Balay El::Copy(*a->emat, *b->emat);
524df311e6cSXuan Zhou }
525df311e6cSXuan Zhou Be->assembled = PETSC_TRUE;
5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
527df311e6cSXuan Zhou }
528df311e6cSXuan Zhou
MatTranspose_Elemental(Mat A,MatReuse reuse,Mat * B)529d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Elemental(Mat A, MatReuse reuse, Mat *B)
530d71ae5a4SJacob Faibussowitsch {
531ec8cb81fSBarry Smith Mat Be = *B;
532ce94432eSBarry Smith MPI_Comm comm;
5334fe7bbcaSHong Zhang Mat_Elemental *a = (Mat_Elemental *)A->data, *b;
534d6223691SXuan Zhou
535d6223691SXuan Zhou PetscFunctionBegin;
5367fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
5379566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5385cb544a0SHong Zhang /* Only out-of-place supported */
5395f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX, comm, PETSC_ERR_SUP, "Only out-of-place supported");
5405262d616SXuan Zhou if (reuse == MAT_INITIAL_MATRIX) {
5419566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Be));
5429566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Be, A->cmap->n, A->rmap->n, PETSC_DECIDE, PETSC_DECIDE));
5439566063dSJacob Faibussowitsch PetscCall(MatSetType(Be, MATELEMENTAL));
5449566063dSJacob Faibussowitsch PetscCall(MatSetUp(Be));
5455262d616SXuan Zhou *B = Be;
5465262d616SXuan Zhou }
5474fe7bbcaSHong Zhang b = (Mat_Elemental *)Be->data;
548e8146892SSatish Balay El::Transpose(*a->emat, *b->emat);
5495262d616SXuan Zhou Be->assembled = PETSC_TRUE;
5503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
551d6223691SXuan Zhou }
552d6223691SXuan Zhou
MatConjugate_Elemental(Mat A)553d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConjugate_Elemental(Mat A)
554d71ae5a4SJacob Faibussowitsch {
555dfcb0403SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
556dfcb0403SXuan Zhou
557dfcb0403SXuan Zhou PetscFunctionBegin;
558e8146892SSatish Balay El::Conjugate(*a->emat);
5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
560dfcb0403SXuan Zhou }
561dfcb0403SXuan Zhou
MatHermitianTranspose_Elemental(Mat A,MatReuse reuse,Mat * B)562d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHermitianTranspose_Elemental(Mat A, MatReuse reuse, Mat *B)
563d71ae5a4SJacob Faibussowitsch {
564ec8cb81fSBarry Smith Mat Be = *B;
565ce94432eSBarry Smith MPI_Comm comm;
5664a29722dSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data, *b;
5674a29722dSXuan Zhou
5684a29722dSXuan Zhou PetscFunctionBegin;
5699566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5705cb544a0SHong Zhang /* Only out-of-place supported */
5714a29722dSXuan Zhou if (reuse == MAT_INITIAL_MATRIX) {
5729566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Be));
5739566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Be, A->cmap->n, A->rmap->n, PETSC_DECIDE, PETSC_DECIDE));
5749566063dSJacob Faibussowitsch PetscCall(MatSetType(Be, MATELEMENTAL));
5759566063dSJacob Faibussowitsch PetscCall(MatSetUp(Be));
5764a29722dSXuan Zhou *B = Be;
5774a29722dSXuan Zhou }
5784a29722dSXuan Zhou b = (Mat_Elemental *)Be->data;
579e8146892SSatish Balay El::Adjoint(*a->emat, *b->emat);
5804a29722dSXuan Zhou Be->assembled = PETSC_TRUE;
5813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5824a29722dSXuan Zhou }
5834a29722dSXuan Zhou
MatSolve_Elemental(Mat A,Vec B,Vec X)584d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_Elemental(Mat A, Vec B, Vec X)
585d71ae5a4SJacob Faibussowitsch {
5861f881ff8SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
587df311e6cSXuan Zhou PetscElemScalar *x;
588ea394475SSatish Balay PetscInt pivoting = a->pivoting;
5891f881ff8SXuan Zhou
5901f881ff8SXuan Zhou PetscFunctionBegin;
5919566063dSJacob Faibussowitsch PetscCall(VecCopy(B, X));
5929566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, (PetscScalar **)&x));
593ea394475SSatish Balay
594e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe;
5950c18141cSBarry Smith xe.Attach(A->rmap->N, 1, *a->grid, 0, 0, x, A->rmap->n);
596e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::MC, El::MR> xer(xe);
597fc54b460SXuan Zhou switch (A->factortype) {
598fc54b460SXuan Zhou case MAT_FACTOR_LU:
599ea394475SSatish Balay if (pivoting == 0) {
600e8146892SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, xer);
601ea394475SSatish Balay } else if (pivoting == 1) {
602ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, xer);
603ea394475SSatish Balay } else { /* pivoting == 2 */
604ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *a->Q, xer);
6051f881ff8SXuan Zhou }
606fc54b460SXuan Zhou break;
607d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY:
608d71ae5a4SJacob Faibussowitsch El::cholesky::SolveAfter(El::UPPER, El::NORMAL, *a->emat, xer);
609d71ae5a4SJacob Faibussowitsch break;
610d71ae5a4SJacob Faibussowitsch default:
611d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unfactored Matrix or Unsupported MatFactorType");
612d71ae5a4SJacob Faibussowitsch break;
6131f881ff8SXuan Zhou }
614ea394475SSatish Balay El::Copy(xer, xe);
615ea394475SSatish Balay
6169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, (PetscScalar **)&x));
6173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6181f881ff8SXuan Zhou }
6191f881ff8SXuan Zhou
MatSolveAdd_Elemental(Mat A,Vec B,Vec Y,Vec X)620d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveAdd_Elemental(Mat A, Vec B, Vec Y, Vec X)
621d71ae5a4SJacob Faibussowitsch {
622df311e6cSXuan Zhou PetscFunctionBegin;
6239566063dSJacob Faibussowitsch PetscCall(MatSolve_Elemental(A, B, X));
6249566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, 1, Y));
6253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
626df311e6cSXuan Zhou }
627df311e6cSXuan Zhou
MatMatSolve_Elemental(Mat A,Mat B,Mat X)628d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_Elemental(Mat A, Mat B, Mat X)
629d71ae5a4SJacob Faibussowitsch {
6301f0e42cfSHong Zhang Mat_Elemental *a = (Mat_Elemental *)A->data;
63142ba530cSPierre Jolivet Mat_Elemental *x;
63242ba530cSPierre Jolivet Mat C;
633ea394475SSatish Balay PetscInt pivoting = a->pivoting;
63442ba530cSPierre Jolivet PetscBool flg;
63542ba530cSPierre Jolivet MatType type;
6361f0e42cfSHong Zhang
637ae844d54SHong Zhang PetscFunctionBegin;
6389566063dSJacob Faibussowitsch PetscCall(MatGetType(X, &type));
6399566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type, MATELEMENTAL, &flg));
64042ba530cSPierre Jolivet if (!flg) {
6419566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &C));
64242ba530cSPierre Jolivet x = (Mat_Elemental *)C->data;
64342ba530cSPierre Jolivet } else {
64442ba530cSPierre Jolivet x = (Mat_Elemental *)X->data;
64542ba530cSPierre Jolivet El::Copy(*((Mat_Elemental *)B->data)->emat, *x->emat);
64642ba530cSPierre Jolivet }
647fc54b460SXuan Zhou switch (A->factortype) {
648fc54b460SXuan Zhou case MAT_FACTOR_LU:
649ea394475SSatish Balay if (pivoting == 0) {
650e8146892SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *x->emat);
651ea394475SSatish Balay } else if (pivoting == 1) {
652ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *x->emat);
653ea394475SSatish Balay } else {
654ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *a->Q, *x->emat);
655d6223691SXuan Zhou }
656fc54b460SXuan Zhou break;
657d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY:
658d71ae5a4SJacob Faibussowitsch El::cholesky::SolveAfter(El::UPPER, El::NORMAL, *a->emat, *x->emat);
659d71ae5a4SJacob Faibussowitsch break;
660d71ae5a4SJacob Faibussowitsch default:
661d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unfactored Matrix or Unsupported MatFactorType");
662d71ae5a4SJacob Faibussowitsch break;
663fc54b460SXuan Zhou }
66442ba530cSPierre Jolivet if (!flg) {
6659566063dSJacob Faibussowitsch PetscCall(MatConvert(C, type, MAT_REUSE_MATRIX, &X));
6669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
66742ba530cSPierre Jolivet }
6683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
669ae844d54SHong Zhang }
670ae844d54SHong Zhang
MatLUFactor_Elemental(Mat A,IS,IS,const MatFactorInfo *)671cedd07caSPierre Jolivet static PetscErrorCode MatLUFactor_Elemental(Mat A, IS, IS, const MatFactorInfo *)
672d71ae5a4SJacob Faibussowitsch {
6737c920d81SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
674ea394475SSatish Balay PetscInt pivoting = a->pivoting;
6757c920d81SXuan Zhou
676ae844d54SHong Zhang PetscFunctionBegin;
677ea394475SSatish Balay if (pivoting == 0) {
678e8146892SSatish Balay El::LU(*a->emat);
679ea394475SSatish Balay } else if (pivoting == 1) {
680ea394475SSatish Balay El::LU(*a->emat, *a->P);
681ea394475SSatish Balay } else {
682ea394475SSatish Balay El::LU(*a->emat, *a->P, *a->Q);
683d6223691SXuan Zhou }
6841293973dSHong Zhang A->factortype = MAT_FACTOR_LU;
685834d3fecSHong Zhang A->assembled = PETSC_TRUE;
686f6224b95SHong Zhang
6879566063dSJacob Faibussowitsch PetscCall(PetscFree(A->solvertype));
6889566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &A->solvertype));
6893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
690ae844d54SHong Zhang }
691ae844d54SHong Zhang
MatLUFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo * info)692d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_Elemental(Mat F, Mat A, const MatFactorInfo *info)
693d71ae5a4SJacob Faibussowitsch {
694d7c3f9d8SHong Zhang PetscFunctionBegin;
6959566063dSJacob Faibussowitsch PetscCall(MatCopy(A, F, SAME_NONZERO_PATTERN));
696f22e26b7SPierre Jolivet PetscCall(MatLUFactor_Elemental(F, nullptr, nullptr, info));
6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
698d7c3f9d8SHong Zhang }
699d7c3f9d8SHong Zhang
MatLUFactorSymbolic_Elemental(Mat,Mat,IS,IS,const MatFactorInfo *)700cedd07caSPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_Elemental(Mat, Mat, IS, IS, const MatFactorInfo *)
701d71ae5a4SJacob Faibussowitsch {
702d7c3f9d8SHong Zhang PetscFunctionBegin;
703d8304050SJose E. Roman /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
7043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
705d7c3f9d8SHong Zhang }
706d7c3f9d8SHong Zhang
MatCholeskyFactor_Elemental(Mat A,IS,const MatFactorInfo *)707cedd07caSPierre Jolivet static PetscErrorCode MatCholeskyFactor_Elemental(Mat A, IS, const MatFactorInfo *)
708d71ae5a4SJacob Faibussowitsch {
70945cf121fSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
710e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::MC, El::STAR> d;
71145cf121fSXuan Zhou
71245cf121fSXuan Zhou PetscFunctionBegin;
713e8146892SSatish Balay El::Cholesky(El::UPPER, *a->emat);
714fc54b460SXuan Zhou A->factortype = MAT_FACTOR_CHOLESKY;
715834d3fecSHong Zhang A->assembled = PETSC_TRUE;
716f6224b95SHong Zhang
7179566063dSJacob Faibussowitsch PetscCall(PetscFree(A->solvertype));
7189566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &A->solvertype));
7193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72045cf121fSXuan Zhou }
72145cf121fSXuan Zhou
MatCholeskyFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo * info)722d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F, Mat A, const MatFactorInfo *info)
723d71ae5a4SJacob Faibussowitsch {
724cb76c1d8SXuan Zhou PetscFunctionBegin;
7259566063dSJacob Faibussowitsch PetscCall(MatCopy(A, F, SAME_NONZERO_PATTERN));
726f22e26b7SPierre Jolivet PetscCall(MatCholeskyFactor_Elemental(F, nullptr, info));
7273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72879673f7bSHong Zhang }
729cb76c1d8SXuan Zhou
MatCholeskyFactorSymbolic_Elemental(Mat,Mat,IS,const MatFactorInfo *)730cedd07caSPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat, Mat, IS, const MatFactorInfo *)
731d71ae5a4SJacob Faibussowitsch {
73279673f7bSHong Zhang PetscFunctionBegin;
733d8304050SJose E. Roman /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
7343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
73579673f7bSHong Zhang }
73679673f7bSHong Zhang
MatFactorGetSolverType_elemental_elemental(Mat,MatSolverType * type)73766976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_elemental_elemental(Mat, MatSolverType *type)
738d71ae5a4SJacob Faibussowitsch {
73915767789SHong Zhang PetscFunctionBegin;
74015767789SHong Zhang *type = MATSOLVERELEMENTAL;
7413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
74215767789SHong Zhang }
74315767789SHong Zhang
MatGetFactor_elemental_elemental(Mat A,MatFactorType ftype,Mat * F)744d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_elemental_elemental(Mat A, MatFactorType ftype, Mat *F)
745d71ae5a4SJacob Faibussowitsch {
7461293973dSHong Zhang Mat B;
7471293973dSHong Zhang
7481293973dSHong Zhang PetscFunctionBegin;
7491293973dSHong Zhang /* Create the factorization matrix */
7509566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
7519566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
7529566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATELEMENTAL));
7539566063dSJacob Faibussowitsch PetscCall(MatSetUp(B));
7541293973dSHong Zhang B->factortype = ftype;
75566e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE;
7569566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype));
7579566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &B->solvertype));
75800c67f3bSHong Zhang
7599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_elemental_elemental));
7601293973dSHong Zhang *F = B;
7613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7621293973dSHong Zhang }
7631293973dSHong Zhang
MatSolverTypeRegister_Elemental(void)764d1f0640dSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Elemental(void)
765d71ae5a4SJacob Faibussowitsch {
76642c9c57cSBarry Smith PetscFunctionBegin;
7679566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL, MATELEMENTAL, MAT_FACTOR_LU, MatGetFactor_elemental_elemental));
7689566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL, MATELEMENTAL, MAT_FACTOR_CHOLESKY, MatGetFactor_elemental_elemental));
7693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
77042c9c57cSBarry Smith }
77142c9c57cSBarry Smith
MatNorm_Elemental(Mat A,NormType type,PetscReal * nrm)772d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNorm_Elemental(Mat A, NormType type, PetscReal *nrm)
773d71ae5a4SJacob Faibussowitsch {
7741f881ff8SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
7751f881ff8SXuan Zhou
7761f881ff8SXuan Zhou PetscFunctionBegin;
7771f881ff8SXuan Zhou switch (type) {
778d71ae5a4SJacob Faibussowitsch case NORM_1:
779d71ae5a4SJacob Faibussowitsch *nrm = El::OneNorm(*a->emat);
780d71ae5a4SJacob Faibussowitsch break;
781d71ae5a4SJacob Faibussowitsch case NORM_FROBENIUS:
782d71ae5a4SJacob Faibussowitsch *nrm = El::FrobeniusNorm(*a->emat);
783d71ae5a4SJacob Faibussowitsch break;
784d71ae5a4SJacob Faibussowitsch case NORM_INFINITY:
785d71ae5a4SJacob Faibussowitsch *nrm = El::InfinityNorm(*a->emat);
786d71ae5a4SJacob Faibussowitsch break;
787d71ae5a4SJacob Faibussowitsch default:
788d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm type");
7891f881ff8SXuan Zhou }
7903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7911f881ff8SXuan Zhou }
7921f881ff8SXuan Zhou
MatZeroEntries_Elemental(Mat A)793d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_Elemental(Mat A)
794d71ae5a4SJacob Faibussowitsch {
7955262d616SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
7965262d616SXuan Zhou
7975262d616SXuan Zhou PetscFunctionBegin;
798e8146892SSatish Balay El::Zero(*a->emat);
7993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8005262d616SXuan Zhou }
8015262d616SXuan Zhou
MatGetOwnershipIS_Elemental(Mat A,IS * rows,IS * cols)802d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A, IS *rows, IS *cols)
803d71ae5a4SJacob Faibussowitsch {
804db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data;
805db31f6deSJed Brown PetscInt i, m, shift, stride, *idx;
806db31f6deSJed Brown
807db31f6deSJed Brown PetscFunctionBegin;
808db31f6deSJed Brown if (rows) {
809db31f6deSJed Brown m = a->emat->LocalHeight();
810db31f6deSJed Brown shift = a->emat->ColShift();
811db31f6deSJed Brown stride = a->emat->ColStride();
8129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx));
813db31f6deSJed Brown for (i = 0; i < m; i++) {
814db31f6deSJed Brown PetscInt rank, offset;
815db31f6deSJed Brown E2RO(A, 0, shift + i * stride, &rank, &offset);
816db31f6deSJed Brown RO2P(A, 0, rank, offset, &idx[i]);
817db31f6deSJed Brown }
8189566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_OWN_POINTER, rows));
819db31f6deSJed Brown }
820db31f6deSJed Brown if (cols) {
821db31f6deSJed Brown m = a->emat->LocalWidth();
822db31f6deSJed Brown shift = a->emat->RowShift();
823db31f6deSJed Brown stride = a->emat->RowStride();
8249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx));
825db31f6deSJed Brown for (i = 0; i < m; i++) {
826db31f6deSJed Brown PetscInt rank, offset;
827db31f6deSJed Brown E2RO(A, 1, shift + i * stride, &rank, &offset);
828db31f6deSJed Brown RO2P(A, 1, rank, offset, &idx[i]);
829db31f6deSJed Brown }
8309566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_OWN_POINTER, cols));
831db31f6deSJed Brown }
8323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
833db31f6deSJed Brown }
834db31f6deSJed Brown
MatConvert_Elemental_Dense(Mat A,MatType,MatReuse reuse,Mat * B)835cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Elemental_Dense(Mat A, MatType, MatReuse reuse, Mat *B)
836d71ae5a4SJacob Faibussowitsch {
8372ef0cf24SXuan Zhou Mat Bmpi;
838af295397SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data;
839ce94432eSBarry Smith MPI_Comm comm;
840fa9e26c8SHong Zhang IS isrows, iscols;
841fa9e26c8SHong Zhang PetscInt rrank, ridx, crank, cidx, nrows, ncols, i, j, erow, ecol, elrow, elcol;
842fa9e26c8SHong Zhang const PetscInt *rows, *cols;
843df311e6cSXuan Zhou PetscElemScalar v;
844fa9e26c8SHong Zhang const El::Grid &grid = a->emat->Grid();
845af295397SXuan Zhou
846af295397SXuan Zhou PetscFunctionBegin;
8479566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
848a000e53bSHong Zhang
849de0a22f0SHong Zhang if (reuse == MAT_REUSE_MATRIX) {
850de0a22f0SHong Zhang Bmpi = *B;
851de0a22f0SHong Zhang } else {
8529566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Bmpi));
8539566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Bmpi, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
8549566063dSJacob Faibussowitsch PetscCall(MatSetType(Bmpi, MATDENSE));
8559566063dSJacob Faibussowitsch PetscCall(MatSetUp(Bmpi));
856de0a22f0SHong Zhang }
857fa9e26c8SHong Zhang
858fa9e26c8SHong Zhang /* Get local entries of A */
8599566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A, &isrows, &iscols));
8609566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows, &nrows));
8619566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows, &rows));
8629566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols, &ncols));
8639566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols, &cols));
864fa9e26c8SHong Zhang
86557299f2fSHong Zhang if (a->roworiented) {
8662ef0cf24SXuan Zhou for (i = 0; i < nrows; i++) {
867fa9e26c8SHong Zhang P2RO(A, 0, rows[i], &rrank, &ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
8682ef0cf24SXuan Zhou RO2E(A, 0, rrank, ridx, &erow);
869aed4548fSBarry Smith PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation");
8702ef0cf24SXuan Zhou for (j = 0; j < ncols; j++) {
871fa9e26c8SHong Zhang P2RO(A, 1, cols[j], &crank, &cidx);
8722ef0cf24SXuan Zhou RO2E(A, 1, crank, cidx, &ecol);
873aed4548fSBarry Smith PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation");
874fa9e26c8SHong Zhang
875fa9e26c8SHong Zhang elrow = erow / grid.MCSize(); /* Elemental local row index */
876fa9e26c8SHong Zhang elcol = ecol / grid.MRSize(); /* Elemental local column index */
877fa9e26c8SHong Zhang v = a->emat->GetLocal(elrow, elcol);
8789566063dSJacob Faibussowitsch PetscCall(MatSetValues(Bmpi, 1, &rows[i], 1, &cols[j], (PetscScalar *)&v, INSERT_VALUES));
8792ef0cf24SXuan Zhou }
8802ef0cf24SXuan Zhou }
88157299f2fSHong Zhang } else { /* column-oriented */
88257299f2fSHong Zhang for (j = 0; j < ncols; j++) {
88357299f2fSHong Zhang P2RO(A, 1, cols[j], &crank, &cidx);
88457299f2fSHong Zhang RO2E(A, 1, crank, cidx, &ecol);
885aed4548fSBarry Smith PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation");
88657299f2fSHong Zhang for (i = 0; i < nrows; i++) {
88757299f2fSHong Zhang P2RO(A, 0, rows[i], &rrank, &ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
88857299f2fSHong Zhang RO2E(A, 0, rrank, ridx, &erow);
889aed4548fSBarry Smith PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation");
89057299f2fSHong Zhang
89157299f2fSHong Zhang elrow = erow / grid.MCSize(); /* Elemental local row index */
89257299f2fSHong Zhang elcol = ecol / grid.MRSize(); /* Elemental local column index */
89357299f2fSHong Zhang v = a->emat->GetLocal(elrow, elcol);
8949566063dSJacob Faibussowitsch PetscCall(MatSetValues(Bmpi, 1, &rows[i], 1, &cols[j], (PetscScalar *)&v, INSERT_VALUES));
89557299f2fSHong Zhang }
89657299f2fSHong Zhang }
89757299f2fSHong Zhang }
8989566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Bmpi, MAT_FINAL_ASSEMBLY));
8999566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Bmpi, MAT_FINAL_ASSEMBLY));
900511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) {
9019566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &Bmpi));
902c4ad791aSXuan Zhou } else {
903c4ad791aSXuan Zhou *B = Bmpi;
904c4ad791aSXuan Zhou }
9059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows));
9069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols));
9073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
908af295397SXuan Zhou }
909af295397SXuan Zhou
MatConvert_SeqAIJ_Elemental(Mat A,MatType,MatReuse reuse,Mat * newmat)910cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
911d71ae5a4SJacob Faibussowitsch {
912af8000cdSHong Zhang Mat mat_elemental;
913af8000cdSHong Zhang PetscInt M = A->rmap->N, N = A->cmap->N, row, ncols;
914af8000cdSHong Zhang const PetscInt *cols;
915af8000cdSHong Zhang const PetscScalar *vals;
916af8000cdSHong Zhang
917af8000cdSHong Zhang PetscFunctionBegin;
918bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) {
919bdeae278SStefano Zampini mat_elemental = *newmat;
9209566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental));
921bdeae278SStefano Zampini } else {
9229566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
9239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N));
9249566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
9259566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental));
926bdeae278SStefano Zampini }
927af8000cdSHong Zhang for (row = 0; row < M; row++) {
9289566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
929bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
9309566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES));
9319566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
932af8000cdSHong Zhang }
9339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
9349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
935af8000cdSHong Zhang
936511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) {
9379566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental));
938af8000cdSHong Zhang } else {
939af8000cdSHong Zhang *newmat = mat_elemental;
940af8000cdSHong Zhang }
9413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
942af8000cdSHong Zhang }
943af8000cdSHong Zhang
MatConvert_MPIAIJ_Elemental(Mat A,MatType,MatReuse reuse,Mat * newmat)944cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
945d71ae5a4SJacob Faibussowitsch {
9465d7652ecSHong Zhang Mat mat_elemental;
9475d7652ecSHong Zhang PetscInt row, ncols, rstart = A->rmap->rstart, rend = A->rmap->rend, j;
9485d7652ecSHong Zhang const PetscInt *cols;
9495d7652ecSHong Zhang const PetscScalar *vals;
9505d7652ecSHong Zhang
9515d7652ecSHong Zhang PetscFunctionBegin;
952bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) {
953bdeae278SStefano Zampini mat_elemental = *newmat;
9549566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental));
955bdeae278SStefano Zampini } else {
9569566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
9579566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, A->rmap->N, A->cmap->N));
9589566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
9599566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental));
960bdeae278SStefano Zampini }
9615d7652ecSHong Zhang for (row = rstart; row < rend; row++) {
9629566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
9635d7652ecSHong Zhang for (j = 0; j < ncols; j++) {
964bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
9659566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, 1, &cols[j], &vals[j], ADD_VALUES));
9665d7652ecSHong Zhang }
9679566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
9685d7652ecSHong Zhang }
9699566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
9709566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
9715d7652ecSHong Zhang
972511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) {
9739566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental));
9745d7652ecSHong Zhang } else {
9755d7652ecSHong Zhang *newmat = mat_elemental;
9765d7652ecSHong Zhang }
9773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9785d7652ecSHong Zhang }
9795d7652ecSHong Zhang
MatConvert_SeqSBAIJ_Elemental(Mat A,MatType,MatReuse reuse,Mat * newmat)980cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
981d71ae5a4SJacob Faibussowitsch {
9826214f412SHong Zhang Mat mat_elemental;
9836214f412SHong Zhang PetscInt M = A->rmap->N, N = A->cmap->N, row, ncols, j;
9846214f412SHong Zhang const PetscInt *cols;
9856214f412SHong Zhang const PetscScalar *vals;
9866214f412SHong Zhang
9876214f412SHong Zhang PetscFunctionBegin;
988bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) {
989bdeae278SStefano Zampini mat_elemental = *newmat;
9909566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental));
991bdeae278SStefano Zampini } else {
9929566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
9939566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N));
9949566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
9959566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental));
996bdeae278SStefano Zampini }
9979566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(A));
9986214f412SHong Zhang for (row = 0; row < M; row++) {
9999566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
1000bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10019566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES));
10026214f412SHong Zhang for (j = 0; j < ncols; j++) { /* lower triangular part */
1003eb1ec7c1SStefano Zampini PetscScalar v;
10046214f412SHong Zhang if (cols[j] == row) continue;
1005b94d7dedSBarry Smith v = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j];
10069566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &cols[j], 1, &row, &v, ADD_VALUES));
10076214f412SHong Zhang }
10089566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
10096214f412SHong Zhang }
10109566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(A));
10119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
10129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
10136214f412SHong Zhang
1014511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) {
10159566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental));
10166214f412SHong Zhang } else {
10176214f412SHong Zhang *newmat = mat_elemental;
10186214f412SHong Zhang }
10193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10206214f412SHong Zhang }
10216214f412SHong Zhang
MatConvert_MPISBAIJ_Elemental(Mat A,MatType,MatReuse reuse,Mat * newmat)1022cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat)
1023d71ae5a4SJacob Faibussowitsch {
10246214f412SHong Zhang Mat mat_elemental;
10256214f412SHong Zhang PetscInt M = A->rmap->N, N = A->cmap->N, row, ncols, j, rstart = A->rmap->rstart, rend = A->rmap->rend;
10266214f412SHong Zhang const PetscInt *cols;
10276214f412SHong Zhang const PetscScalar *vals;
10286214f412SHong Zhang
10296214f412SHong Zhang PetscFunctionBegin;
1030bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) {
1031bdeae278SStefano Zampini mat_elemental = *newmat;
10329566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental));
1033bdeae278SStefano Zampini } else {
10349566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
10359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N));
10369566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
10379566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental));
1038bdeae278SStefano Zampini }
10399566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(A));
10406214f412SHong Zhang for (row = rstart; row < rend; row++) {
10419566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals));
1042bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10439566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES));
10446214f412SHong Zhang for (j = 0; j < ncols; j++) { /* lower triangular part */
1045eb1ec7c1SStefano Zampini PetscScalar v;
10466214f412SHong Zhang if (cols[j] == row) continue;
1047b94d7dedSBarry Smith v = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j];
10489566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &cols[j], 1, &row, &v, ADD_VALUES));
10496214f412SHong Zhang }
10509566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals));
10516214f412SHong Zhang }
10529566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(A));
10539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
10549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
10556214f412SHong Zhang
1056511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) {
10579566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental));
10586214f412SHong Zhang } else {
10596214f412SHong Zhang *newmat = mat_elemental;
10606214f412SHong Zhang }
10613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10626214f412SHong Zhang }
10636214f412SHong Zhang
MatDestroy_Elemental(Mat A)1064d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Elemental(Mat A)
1065d71ae5a4SJacob Faibussowitsch {
1066db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data;
10675e9f5b67SHong Zhang Mat_Elemental_Grid *commgrid;
1068b8b5be36SMartin Diehl PetscMPIInt iflg;
10695e9f5b67SHong Zhang MPI_Comm icomm;
1070db31f6deSJed Brown
1071db31f6deSJed Brown PetscFunctionBegin;
1072db31f6deSJed Brown delete a->emat;
1073ea394475SSatish Balay delete a->P;
1074ea394475SSatish Balay delete a->Q;
10755e9f5b67SHong Zhang
1076e8146892SSatish Balay El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
1077f22e26b7SPierre Jolivet PetscCall(PetscCommDuplicate(cxxcomm.comm, &icomm, nullptr));
1078b8b5be36SMartin Diehl PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Elemental_keyval, (void **)&commgrid, &iflg));
10795e9f5b67SHong Zhang if (--commgrid->grid_refct == 0) {
10805e9f5b67SHong Zhang delete commgrid->grid;
10819566063dSJacob Faibussowitsch PetscCall(PetscFree(commgrid));
10829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Elemental_keyval));
10835e9f5b67SHong Zhang }
10849566063dSJacob Faibussowitsch PetscCall(PetscCommDestroy(&icomm));
1085f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetOwnershipIS_C", nullptr));
1086f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", nullptr));
1087f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_elemental_mpidense_C", nullptr));
10889566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data));
10893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1090db31f6deSJed Brown }
1091db31f6deSJed Brown
MatSetUp_Elemental(Mat A)109266976f2fSJacob Faibussowitsch static PetscErrorCode MatSetUp_Elemental(Mat A)
1093d71ae5a4SJacob Faibussowitsch {
1094db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data;
1095f19600a0SHong Zhang MPI_Comm comm;
1096db31f6deSJed Brown PetscMPIInt rsize, csize;
1097f19600a0SHong Zhang PetscInt n;
1098db31f6deSJed Brown
1099db31f6deSJed Brown PetscFunctionBegin;
11009566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap));
11019566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap));
1102db31f6deSJed Brown
1103d8304050SJose E. Roman /* Check if local row and column sizes are equally distributed.
1104f19600a0SHong Zhang Jed: Elemental uses "element" cyclic ordering so the sizes need to match that
1105f19600a0SHong Zhang exactly. The strategy in MatElemental is for PETSc to implicitly permute to block ordering (like would be returned by
1106d8304050SJose E. Roman PetscSplitOwnership(comm,&n,&N), at which point Elemental matrices can act on PETSc vectors without redistributing the vectors. */
11079566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1108f19600a0SHong Zhang n = PETSC_DECIDE;
11099566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(comm, &n, &A->rmap->N));
11105f80ce2aSJacob 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);
1111f19600a0SHong Zhang
1112f19600a0SHong Zhang n = PETSC_DECIDE;
11139566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(comm, &n, &A->cmap->N));
11145f80ce2aSJacob 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);
1115f19600a0SHong Zhang
11162f613bf5SBarry Smith a->emat->Resize(A->rmap->N, A->cmap->N);
1117e8146892SSatish Balay El::Zero(*a->emat);
1118db31f6deSJed Brown
11199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(A->rmap->comm, &rsize));
11209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(A->cmap->comm, &csize));
11215f80ce2aSJacob Faibussowitsch PetscCheck(csize == rsize, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Cannot use row and column communicators of different sizes");
1122db31f6deSJed Brown a->commsize = rsize;
11239371c9d4SSatish Balay a->mr[0] = A->rmap->N % rsize;
11249371c9d4SSatish Balay if (!a->mr[0]) a->mr[0] = rsize;
11259371c9d4SSatish Balay a->mr[1] = A->cmap->N % csize;
11269371c9d4SSatish Balay if (!a->mr[1]) a->mr[1] = csize;
1127db31f6deSJed Brown a->m[0] = A->rmap->N / rsize + (a->mr[0] != rsize);
1128db31f6deSJed Brown a->m[1] = A->cmap->N / csize + (a->mr[1] != csize);
11293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1130db31f6deSJed Brown }
1131db31f6deSJed Brown
MatAssemblyBegin_Elemental(Mat A,MatAssemblyType)113266976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType)
1133d71ae5a4SJacob Faibussowitsch {
1134aae2c449SHong Zhang Mat_Elemental *a = (Mat_Elemental *)A->data;
1135aae2c449SHong Zhang
1136aae2c449SHong Zhang PetscFunctionBegin;
1137fbbcb730SJack Poulson /* printf("Calling ProcessQueues\n"); */
1138fbbcb730SJack Poulson a->emat->ProcessQueues();
1139fbbcb730SJack Poulson /* printf("Finished ProcessQueues\n"); */
11403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1141aae2c449SHong Zhang }
1142aae2c449SHong Zhang
MatAssemblyEnd_Elemental(Mat,MatAssemblyType)114366976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_Elemental(Mat, MatAssemblyType)
1144d71ae5a4SJacob Faibussowitsch {
1145aae2c449SHong Zhang PetscFunctionBegin;
1146aae2c449SHong Zhang /* Currently does nothing */
11473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1148aae2c449SHong Zhang }
1149aae2c449SHong Zhang
MatLoad_Elemental(Mat newMat,PetscViewer viewer)115066976f2fSJacob Faibussowitsch static PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer)
1151d71ae5a4SJacob Faibussowitsch {
11527b30e0f1SHong Zhang Mat Adense, Ae;
11537b30e0f1SHong Zhang MPI_Comm comm;
11547b30e0f1SHong Zhang
11557b30e0f1SHong Zhang PetscFunctionBegin;
11569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)newMat, &comm));
11579566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Adense));
11589566063dSJacob Faibussowitsch PetscCall(MatSetType(Adense, MATDENSE));
11599566063dSJacob Faibussowitsch PetscCall(MatLoad(Adense, viewer));
11609566063dSJacob Faibussowitsch PetscCall(MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae));
11619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adense));
11629566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(newMat, &Ae));
11633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11647b30e0f1SHong Zhang }
11657b30e0f1SHong Zhang
11669371c9d4SSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_Elemental,
1167f22e26b7SPierre Jolivet nullptr,
1168f22e26b7SPierre Jolivet nullptr,
116940d92e34SHong Zhang MatMult_Elemental,
117040d92e34SHong Zhang /* 4*/ MatMultAdd_Elemental,
11719426833fSXuan Zhou MatMultTranspose_Elemental,
1172e883f9d5SXuan Zhou MatMultTransposeAdd_Elemental,
117340d92e34SHong Zhang MatSolve_Elemental,
1174df311e6cSXuan Zhou MatSolveAdd_Elemental,
1175f22e26b7SPierre Jolivet nullptr,
1176f22e26b7SPierre Jolivet /*10*/ nullptr,
117740d92e34SHong Zhang MatLUFactor_Elemental,
117840d92e34SHong Zhang MatCholeskyFactor_Elemental,
1179f22e26b7SPierre Jolivet nullptr,
118040d92e34SHong Zhang MatTranspose_Elemental,
118140d92e34SHong Zhang /*15*/ MatGetInfo_Elemental,
1182f22e26b7SPierre Jolivet nullptr,
118361119200SXuan Zhou MatGetDiagonal_Elemental,
1184ade3cc5eSXuan Zhou MatDiagonalScale_Elemental,
118540d92e34SHong Zhang MatNorm_Elemental,
118640d92e34SHong Zhang /*20*/ MatAssemblyBegin_Elemental,
118740d92e34SHong Zhang MatAssemblyEnd_Elemental,
118832d7a744SHong Zhang MatSetOption_Elemental,
118940d92e34SHong Zhang MatZeroEntries_Elemental,
1190f22e26b7SPierre Jolivet /*24*/ nullptr,
119140d92e34SHong Zhang MatLUFactorSymbolic_Elemental,
119240d92e34SHong Zhang MatLUFactorNumeric_Elemental,
119340d92e34SHong Zhang MatCholeskyFactorSymbolic_Elemental,
119440d92e34SHong Zhang MatCholeskyFactorNumeric_Elemental,
119540d92e34SHong Zhang /*29*/ MatSetUp_Elemental,
1196f22e26b7SPierre Jolivet nullptr,
1197f22e26b7SPierre Jolivet nullptr,
1198f22e26b7SPierre Jolivet nullptr,
1199f22e26b7SPierre Jolivet nullptr,
1200df311e6cSXuan Zhou /*34*/ MatDuplicate_Elemental,
1201f22e26b7SPierre Jolivet nullptr,
1202f22e26b7SPierre Jolivet nullptr,
1203f22e26b7SPierre Jolivet nullptr,
1204f22e26b7SPierre Jolivet nullptr,
120540d92e34SHong Zhang /*39*/ MatAXPY_Elemental,
1206f22e26b7SPierre Jolivet nullptr,
1207f22e26b7SPierre Jolivet nullptr,
1208f22e26b7SPierre Jolivet nullptr,
120940d92e34SHong Zhang MatCopy_Elemental,
1210f22e26b7SPierre Jolivet /*44*/ nullptr,
121140d92e34SHong Zhang MatScale_Elemental,
12127d68702bSBarry Smith MatShift_Basic,
1213f22e26b7SPierre Jolivet nullptr,
1214f22e26b7SPierre Jolivet nullptr,
1215f22e26b7SPierre Jolivet /*49*/ nullptr,
1216f22e26b7SPierre Jolivet nullptr,
1217f22e26b7SPierre Jolivet nullptr,
1218f22e26b7SPierre Jolivet nullptr,
1219f22e26b7SPierre Jolivet nullptr,
1220f22e26b7SPierre Jolivet /*54*/ nullptr,
1221f22e26b7SPierre Jolivet nullptr,
1222f22e26b7SPierre Jolivet nullptr,
1223f22e26b7SPierre Jolivet nullptr,
1224f22e26b7SPierre Jolivet nullptr,
1225f22e26b7SPierre Jolivet /*59*/ nullptr,
122640d92e34SHong Zhang MatDestroy_Elemental,
122740d92e34SHong Zhang MatView_Elemental,
1228f22e26b7SPierre Jolivet nullptr,
1229f22e26b7SPierre Jolivet nullptr,
1230f22e26b7SPierre Jolivet /*64*/ nullptr,
1231f22e26b7SPierre Jolivet nullptr,
1232f22e26b7SPierre Jolivet nullptr,
1233f22e26b7SPierre Jolivet nullptr,
1234f22e26b7SPierre Jolivet nullptr,
1235f22e26b7SPierre Jolivet /*69*/ nullptr,
12362ef0cf24SXuan Zhou MatConvert_Elemental_Dense,
1237f22e26b7SPierre Jolivet nullptr,
1238f22e26b7SPierre Jolivet nullptr,
12398bb0f5c6SPierre Jolivet nullptr,
1240f22e26b7SPierre Jolivet /*74*/ nullptr,
1241f22e26b7SPierre Jolivet nullptr,
1242f22e26b7SPierre Jolivet nullptr,
1243f22e26b7SPierre Jolivet nullptr,
12448bb0f5c6SPierre Jolivet MatLoad_Elemental,
1245f22e26b7SPierre Jolivet /*79*/ nullptr,
1246f22e26b7SPierre Jolivet nullptr,
1247f22e26b7SPierre Jolivet nullptr,
1248f22e26b7SPierre Jolivet nullptr,
12498bb0f5c6SPierre Jolivet nullptr,
1250f22e26b7SPierre Jolivet /*84*/ nullptr,
125140d92e34SHong Zhang MatMatMultNumeric_Elemental,
1252f22e26b7SPierre Jolivet nullptr,
1253f22e26b7SPierre Jolivet nullptr,
1254df311e6cSXuan Zhou MatMatTransposeMultNumeric_Elemental,
12558bb0f5c6SPierre Jolivet /*89*/ nullptr,
12568bb0f5c6SPierre Jolivet MatProductSetFromOptions_Elemental,
1257f22e26b7SPierre Jolivet nullptr,
1258f22e26b7SPierre Jolivet nullptr,
1259dfcb0403SXuan Zhou MatConjugate_Elemental,
12608bb0f5c6SPierre Jolivet /*94*/ nullptr,
1261f22e26b7SPierre Jolivet nullptr,
1262f22e26b7SPierre Jolivet nullptr,
1263f22e26b7SPierre Jolivet nullptr,
1264f22e26b7SPierre Jolivet nullptr,
12658bb0f5c6SPierre Jolivet /*99*/ nullptr,
12668bb0f5c6SPierre Jolivet MatMatSolve_Elemental,
1267f22e26b7SPierre Jolivet nullptr,
1268f22e26b7SPierre Jolivet nullptr,
1269f22e26b7SPierre Jolivet nullptr,
1270*421480d9SBarry Smith /*104*/ nullptr,
12718bb0f5c6SPierre Jolivet nullptr,
12728bb0f5c6SPierre Jolivet nullptr,
12738bb0f5c6SPierre Jolivet nullptr,
12748bb0f5c6SPierre Jolivet nullptr,
12758bb0f5c6SPierre Jolivet /*109*/ nullptr,
12768bb0f5c6SPierre Jolivet MatHermitianTranspose_Elemental,
12778bb0f5c6SPierre Jolivet nullptr,
12788bb0f5c6SPierre Jolivet nullptr,
1279f22e26b7SPierre Jolivet /*114*/ nullptr,
1280f22e26b7SPierre Jolivet nullptr,
1281f22e26b7SPierre Jolivet nullptr,
1282f22e26b7SPierre Jolivet nullptr,
1283f22e26b7SPierre Jolivet nullptr,
1284f22e26b7SPierre Jolivet /*119*/ nullptr,
12858bb0f5c6SPierre Jolivet nullptr,
1286f22e26b7SPierre Jolivet nullptr,
1287f22e26b7SPierre Jolivet nullptr,
1288f22e26b7SPierre Jolivet nullptr,
1289f22e26b7SPierre Jolivet /*124*/ nullptr,
1290f22e26b7SPierre Jolivet nullptr,
1291f22e26b7SPierre Jolivet nullptr,
1292f22e26b7SPierre Jolivet nullptr,
1293f22e26b7SPierre Jolivet nullptr,
1294f22e26b7SPierre Jolivet /*129*/ nullptr,
1295f22e26b7SPierre Jolivet nullptr,
1296f22e26b7SPierre Jolivet nullptr,
1297f22e26b7SPierre Jolivet nullptr,
1298f22e26b7SPierre Jolivet nullptr,
1299f22e26b7SPierre Jolivet /*134*/ nullptr,
1300f22e26b7SPierre Jolivet nullptr,
1301f22e26b7SPierre Jolivet nullptr,
1302f22e26b7SPierre Jolivet nullptr,
1303f22e26b7SPierre Jolivet nullptr,
1304f22e26b7SPierre Jolivet nullptr,
1305f22e26b7SPierre Jolivet /*140*/ nullptr,
1306f22e26b7SPierre Jolivet nullptr,
1307c97e14ecSAlex Lindsay nullptr,
1308c2be7ffeSStefano Zampini nullptr,
1309f22e26b7SPierre Jolivet nullptr};
131040d92e34SHong Zhang
1311ed36708cSHong Zhang /*MC
1312ed36708cSHong Zhang MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package
1313ed36708cSHong Zhang
1314c2b89b5dSBarry Smith Use ./configure --download-elemental to install PETSc to use Elemental
1315c2b89b5dSBarry Smith
1316ed36708cSHong Zhang Options Database Keys:
13175cc86fc1SJed Brown + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions()
131889bba20eSBarry Smith . -pc_factor_mat_solver_type elemental - to use this direct solver with the option -pc_type lu
13195cc86fc1SJed Brown - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix
1320ed36708cSHong Zhang
1321ed36708cSHong Zhang Level: beginner
1322ed36708cSHong Zhang
132389bba20eSBarry Smith Note:
132489bba20eSBarry Smith Note unlike most matrix formats, this format does not store all the matrix entries for a contiguous
132589bba20eSBarry Smith range of rows on an MPI rank. Use `MatGetOwnershipIS()` to determine what values are stored on
132689bba20eSBarry Smith the given rank.
132789bba20eSBarry Smith
132889bba20eSBarry Smith .seealso: `MATDENSE`, `MATSCALAPACK`, `MatGetOwnershipIS()`
1329ed36708cSHong Zhang M*/
1330f22e26b7SPierre Jolivet #if defined(__clang__)
1331f22e26b7SPierre Jolivet #pragma clang diagnostic push
1332f22e26b7SPierre Jolivet #pragma clang diagnostic ignored "-Wzero-as-null-pointer-constant"
1333f22e26b7SPierre Jolivet #endif
MatCreate_Elemental(Mat A)1334d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A)
1335d71ae5a4SJacob Faibussowitsch {
1336db31f6deSJed Brown Mat_Elemental *a;
1337b8b5be36SMartin Diehl PetscBool flg;
1338b8b5be36SMartin Diehl PetscMPIInt iflg;
13395e9f5b67SHong Zhang Mat_Elemental_Grid *commgrid;
13405e9f5b67SHong Zhang MPI_Comm icomm;
13415682a260SJack Poulson PetscInt optv1;
1342db31f6deSJed Brown
1343db31f6deSJed Brown PetscFunctionBegin;
1344aea10558SJacob Faibussowitsch A->ops[0] = MatOps_Values;
134540d92e34SHong Zhang A->insertmode = NOT_SET_VALUES;
1346db31f6deSJed Brown
13474dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a));
1348db31f6deSJed Brown A->data = (void *)a;
1349db31f6deSJed Brown
1350db31f6deSJed Brown /* Set up the elemental matrix */
1351e8146892SSatish Balay El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
13525e9f5b67SHong Zhang
13535e9f5b67SHong Zhang /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */
13545e9f5b67SHong Zhang if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) {
1355c8025a54SPierre Jolivet PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Elemental_keyval, nullptr));
13569566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(ElementalCitation, &ElementalCite));
13575e9f5b67SHong Zhang }
13589566063dSJacob Faibussowitsch PetscCall(PetscCommDuplicate(cxxcomm.comm, &icomm, NULL));
1359b8b5be36SMartin Diehl PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Elemental_keyval, (void **)&commgrid, &iflg));
1360b8b5be36SMartin Diehl if (!iflg) {
13614dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&commgrid));
13625cb544a0SHong Zhang
1363d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "Elemental Options", "Mat");
1364e8146892SSatish Balay /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */
1365b8b5be36SMartin Diehl PetscCall(PetscOptionsInt("-mat_elemental_grid_height", "Grid Height", "None", El::mpi::Size(cxxcomm), &optv1, &flg));
1366b8b5be36SMartin Diehl if (flg) {
1367aed4548fSBarry 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));
1368e8146892SSatish Balay commgrid->grid = new El::Grid(cxxcomm, optv1); /* use user-provided grid height */
13692ef0cf24SXuan Zhou } else {
1370e8146892SSatish Balay commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */
1371da0640a4SHong Zhang /* printf("new commgrid->grid = %p\n",commgrid->grid); -- memory leak revealed by valgrind? */
1372ed667823SXuan Zhou }
13735e9f5b67SHong Zhang commgrid->grid_refct = 1;
13749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_attr(icomm, Petsc_Elemental_keyval, (void *)commgrid));
1375ea394475SSatish Balay
1376ea394475SSatish Balay a->pivoting = 1;
13779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_elemental_pivoting", "Pivoting", "None", a->pivoting, &a->pivoting, NULL));
1378ea394475SSatish Balay
1379d0609cedSBarry Smith PetscOptionsEnd();
13805e9f5b67SHong Zhang } else {
13815e9f5b67SHong Zhang commgrid->grid_refct++;
13825e9f5b67SHong Zhang }
13839566063dSJacob Faibussowitsch PetscCall(PetscCommDestroy(&icomm));
13845e9f5b67SHong Zhang a->grid = commgrid->grid;
1385e8146892SSatish Balay a->emat = new El::DistMatrix<PetscElemScalar>(*a->grid);
138632d7a744SHong Zhang a->roworiented = PETSC_TRUE;
1387db31f6deSJed Brown
13889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetOwnershipIS_C", MatGetOwnershipIS_Elemental));
13899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_elemental_mpidense_C", MatProductSetFromOptions_Elemental_MPIDense));
13909566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATELEMENTAL));
13913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1392db31f6deSJed Brown }
1393f22e26b7SPierre Jolivet #if defined(__clang__)
1394f22e26b7SPierre Jolivet #pragma clang diagnostic pop
1395f22e26b7SPierre Jolivet #endif
1396