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) { 711c2dc1cbSBarry Smith //PetscCall(MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin))); 7215767789SHong Zhang /* see MatGetInfo_MPIAIJ() for getting global info->nz_allocated! */ 7315767789SHong Zhang //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_MAX not written yet"); 7415767789SHong Zhang } else if (flag == MAT_GLOBAL_SUM) { 7515767789SHong Zhang //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_SUM not written yet"); 763966268fSBarry Smith info->nz_allocated = (*a->emat).AllocatedMemory(); /* locally allocated */ 7715767789SHong Zhang info->nz_used = info->nz_allocated; /* assume Elemental does accurate allocation */ 781c2dc1cbSBarry Smith //PetscCall(MPIU_Allreduce(isend,irecv,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A))); 7915767789SHong Zhang //PetscPrintf(PETSC_COMM_SELF," ... [%d] locally allocated %g\n",rank,info->nz_allocated); 8015767789SHong Zhang } 8115767789SHong Zhang 8215767789SHong Zhang info->nz_unneeded = 0.0; 833966268fSBarry Smith info->assemblies = A->num_ass; 8415767789SHong Zhang info->mallocs = 0; 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 92d71ae5a4SJacob Faibussowitsch 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) { 9832d7a744SHong Zhang case MAT_NEW_NONZERO_LOCATIONS: 9932d7a744SHong Zhang case MAT_NEW_NONZERO_LOCATION_ERR: 10032d7a744SHong Zhang case MAT_NEW_NONZERO_ALLOCATION_ERR: 101891a0571SHong Zhang case MAT_SYMMETRIC: 102071fcb05SBarry Smith case MAT_SORTED_FULL: 103d71ae5a4SJacob Faibussowitsch case MAT_HERMITIAN: 104d71ae5a4SJacob Faibussowitsch break; 105d71ae5a4SJacob Faibussowitsch case MAT_ROW_ORIENTED: 106d71ae5a4SJacob Faibussowitsch a->roworiented = flg; 107d71ae5a4SJacob Faibussowitsch break; 108d71ae5a4SJacob Faibussowitsch default: 109d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %s", MatOptions[op]); 11032d7a744SHong Zhang } 1113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11232d7a744SHong Zhang } 11332d7a744SHong Zhang 114d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValues_Elemental(Mat A, PetscInt nr, const PetscInt *rows, PetscInt nc, const PetscInt *cols, const PetscScalar *vals, InsertMode imode) 115d71ae5a4SJacob Faibussowitsch { 116db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data; 117fbbcb730SJack Poulson PetscInt i, j, rrank, ridx, crank, cidx, erow, ecol, numQueues = 0; 118db31f6deSJed Brown 119db31f6deSJed Brown PetscFunctionBegin; 120fbbcb730SJack Poulson // TODO: Initialize matrix to all zeros? 121fbbcb730SJack Poulson 122fbbcb730SJack Poulson // Count the number of queues from this process 12332d7a744SHong Zhang if (a->roworiented) { 124db31f6deSJed Brown for (i = 0; i < nr; i++) { 125db31f6deSJed Brown if (rows[i] < 0) continue; 126db31f6deSJed Brown P2RO(A, 0, rows[i], &rrank, &ridx); 127db31f6deSJed Brown RO2E(A, 0, rrank, ridx, &erow); 128aed4548fSBarry Smith PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect row translation"); 129db31f6deSJed Brown for (j = 0; j < nc; j++) { 130db31f6deSJed Brown if (cols[j] < 0) continue; 131db31f6deSJed Brown P2RO(A, 1, cols[j], &crank, &cidx); 132db31f6deSJed Brown RO2E(A, 1, crank, cidx, &ecol); 133aed4548fSBarry Smith PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect col translation"); 134fbbcb730SJack Poulson if (!a->emat->IsLocal(erow, ecol)) { /* off-proc entry */ 135fbbcb730SJack Poulson /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */ 13608401ef6SPierre Jolivet PetscCheck(imode == ADD_VALUES, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only ADD_VALUES to off-processor entry is supported"); 137fbbcb730SJack Poulson ++numQueues; 138aae2c449SHong Zhang continue; 139ed36708cSHong Zhang } 140fbbcb730SJack Poulson /* printf("Locally updating (%d,%d)\n",erow,ecol); */ 141db31f6deSJed Brown switch (imode) { 142d71ae5a4SJacob Faibussowitsch case INSERT_VALUES: 143d71ae5a4SJacob Faibussowitsch a->emat->Set(erow, ecol, (PetscElemScalar)vals[i * nc + j]); 144d71ae5a4SJacob Faibussowitsch break; 145d71ae5a4SJacob Faibussowitsch case ADD_VALUES: 146d71ae5a4SJacob Faibussowitsch a->emat->Update(erow, ecol, (PetscElemScalar)vals[i * nc + j]); 147d71ae5a4SJacob Faibussowitsch break; 148d71ae5a4SJacob Faibussowitsch default: 149d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for InsertMode %d", (int)imode); 150db31f6deSJed Brown } 151db31f6deSJed Brown } 152db31f6deSJed Brown } 153fbbcb730SJack Poulson 154fbbcb730SJack Poulson /* printf("numQueues=%d\n",numQueues); */ 155fbbcb730SJack Poulson a->emat->Reserve(numQueues); 156fbbcb730SJack Poulson for (i = 0; i < nr; i++) { 157fbbcb730SJack Poulson if (rows[i] < 0) continue; 158fbbcb730SJack Poulson P2RO(A, 0, rows[i], &rrank, &ridx); 159fbbcb730SJack Poulson RO2E(A, 0, rrank, ridx, &erow); 160fbbcb730SJack Poulson for (j = 0; j < nc; j++) { 161fbbcb730SJack Poulson if (cols[j] < 0) continue; 162fbbcb730SJack Poulson P2RO(A, 1, cols[j], &crank, &cidx); 163fbbcb730SJack Poulson RO2E(A, 1, crank, cidx, &ecol); 164fbbcb730SJack Poulson if (!a->emat->IsLocal(erow, ecol)) { /*off-proc entry*/ 165fbbcb730SJack Poulson /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */ 166fbbcb730SJack Poulson a->emat->QueueUpdate(erow, ecol, vals[i * nc + j]); 167fbbcb730SJack Poulson } 168fbbcb730SJack Poulson } 169fbbcb730SJack Poulson } 17032d7a744SHong Zhang } else { /* columnoriented */ 17132d7a744SHong Zhang for (j = 0; j < nc; j++) { 17232d7a744SHong Zhang if (cols[j] < 0) continue; 17332d7a744SHong Zhang P2RO(A, 1, cols[j], &crank, &cidx); 17432d7a744SHong Zhang RO2E(A, 1, crank, cidx, &ecol); 175aed4548fSBarry Smith PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect col translation"); 17632d7a744SHong Zhang for (i = 0; i < nr; i++) { 17732d7a744SHong Zhang if (rows[i] < 0) continue; 17832d7a744SHong Zhang P2RO(A, 0, rows[i], &rrank, &ridx); 17932d7a744SHong Zhang RO2E(A, 0, rrank, ridx, &erow); 180aed4548fSBarry Smith PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Incorrect row translation"); 18132d7a744SHong Zhang if (!a->emat->IsLocal(erow, ecol)) { /* off-proc entry */ 18232d7a744SHong Zhang /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */ 18308401ef6SPierre Jolivet PetscCheck(imode == ADD_VALUES, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only ADD_VALUES to off-processor entry is supported"); 18432d7a744SHong Zhang ++numQueues; 18532d7a744SHong Zhang continue; 18632d7a744SHong Zhang } 18732d7a744SHong Zhang /* printf("Locally updating (%d,%d)\n",erow,ecol); */ 18832d7a744SHong Zhang switch (imode) { 189d71ae5a4SJacob Faibussowitsch case INSERT_VALUES: 190d71ae5a4SJacob Faibussowitsch a->emat->Set(erow, ecol, (PetscElemScalar)vals[i + j * nr]); 191d71ae5a4SJacob Faibussowitsch break; 192d71ae5a4SJacob Faibussowitsch case ADD_VALUES: 193d71ae5a4SJacob Faibussowitsch a->emat->Update(erow, ecol, (PetscElemScalar)vals[i + j * nr]); 194d71ae5a4SJacob Faibussowitsch break; 195d71ae5a4SJacob Faibussowitsch default: 196d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for InsertMode %d", (int)imode); 19732d7a744SHong Zhang } 19832d7a744SHong Zhang } 19932d7a744SHong Zhang } 20032d7a744SHong Zhang 20132d7a744SHong Zhang /* printf("numQueues=%d\n",numQueues); */ 20232d7a744SHong Zhang a->emat->Reserve(numQueues); 20332d7a744SHong Zhang for (j = 0; j < nc; j++) { 20432d7a744SHong Zhang if (cols[j] < 0) continue; 20532d7a744SHong Zhang P2RO(A, 1, cols[j], &crank, &cidx); 20632d7a744SHong Zhang RO2E(A, 1, crank, cidx, &ecol); 20732d7a744SHong Zhang 20832d7a744SHong Zhang for (i = 0; i < nr; i++) { 20932d7a744SHong Zhang if (rows[i] < 0) continue; 21032d7a744SHong Zhang P2RO(A, 0, rows[i], &rrank, &ridx); 21132d7a744SHong Zhang RO2E(A, 0, rrank, ridx, &erow); 21232d7a744SHong Zhang if (!a->emat->IsLocal(erow, ecol)) { /*off-proc entry*/ 21332d7a744SHong Zhang /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */ 21432d7a744SHong Zhang a->emat->QueueUpdate(erow, ecol, vals[i + j * nr]); 21532d7a744SHong Zhang } 21632d7a744SHong Zhang } 21732d7a744SHong Zhang } 21832d7a744SHong Zhang } 2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 220db31f6deSJed Brown } 221db31f6deSJed Brown 222d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Elemental(Mat A, Vec X, Vec Y) 223d71ae5a4SJacob Faibussowitsch { 224db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data; 225e6dea9dbSXuan Zhou const PetscElemScalar *x; 226e6dea9dbSXuan Zhou PetscElemScalar *y; 227df311e6cSXuan Zhou PetscElemScalar one = 1, zero = 0; 228db31f6deSJed Brown 229db31f6deSJed Brown PetscFunctionBegin; 2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); 2319566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, (PetscScalar **)&y)); 232db31f6deSJed Brown { /* Scoping so that constructor is called before pointer is returned */ 233e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ye; 2340c18141cSBarry Smith xe.LockedAttach(A->cmap->N, 1, *a->grid, 0, 0, x, A->cmap->n); 2350c18141cSBarry Smith ye.Attach(A->rmap->N, 1, *a->grid, 0, 0, y, A->rmap->n); 236e8146892SSatish Balay El::Gemv(El::NORMAL, one, *a->emat, xe, zero, ye); 237db31f6deSJed Brown } 2389566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x)); 2399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, (PetscScalar **)&y)); 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 241db31f6deSJed Brown } 242db31f6deSJed Brown 243d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Elemental(Mat A, Vec X, Vec Y) 244d71ae5a4SJacob Faibussowitsch { 2459426833fSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 246df311e6cSXuan Zhou const PetscElemScalar *x; 247df311e6cSXuan Zhou PetscElemScalar *y; 248e6dea9dbSXuan Zhou PetscElemScalar one = 1, zero = 0; 2499426833fSXuan Zhou 2509426833fSXuan Zhou PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); 2529566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, (PetscScalar **)&y)); 2539426833fSXuan Zhou { /* Scoping so that constructor is called before pointer is returned */ 254e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ye; 2550c18141cSBarry Smith xe.LockedAttach(A->rmap->N, 1, *a->grid, 0, 0, x, A->rmap->n); 2560c18141cSBarry Smith ye.Attach(A->cmap->N, 1, *a->grid, 0, 0, y, A->cmap->n); 257e8146892SSatish Balay El::Gemv(El::TRANSPOSE, one, *a->emat, xe, zero, ye); 2589426833fSXuan Zhou } 2599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x)); 2609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, (PetscScalar **)&y)); 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2629426833fSXuan Zhou } 2639426833fSXuan Zhou 264d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_Elemental(Mat A, Vec X, Vec Y, Vec Z) 265d71ae5a4SJacob Faibussowitsch { 266db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data; 267df311e6cSXuan Zhou const PetscElemScalar *x; 268df311e6cSXuan Zhou PetscElemScalar *z; 269e6dea9dbSXuan Zhou PetscElemScalar one = 1; 270db31f6deSJed Brown 271db31f6deSJed Brown PetscFunctionBegin; 2729566063dSJacob Faibussowitsch if (Y != Z) PetscCall(VecCopy(Y, Z)); 2739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); 2749566063dSJacob Faibussowitsch PetscCall(VecGetArray(Z, (PetscScalar **)&z)); 275db31f6deSJed Brown { /* Scoping so that constructor is called before pointer is returned */ 276e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ze; 2770c18141cSBarry Smith xe.LockedAttach(A->cmap->N, 1, *a->grid, 0, 0, x, A->cmap->n); 2780c18141cSBarry Smith ze.Attach(A->rmap->N, 1, *a->grid, 0, 0, z, A->rmap->n); 279e8146892SSatish Balay El::Gemv(El::NORMAL, one, *a->emat, xe, one, ze); 280db31f6deSJed Brown } 2819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x)); 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Z, (PetscScalar **)&z)); 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284db31f6deSJed Brown } 285db31f6deSJed Brown 286d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_Elemental(Mat A, Vec X, Vec Y, Vec Z) 287d71ae5a4SJacob Faibussowitsch { 288e883f9d5SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 289df311e6cSXuan Zhou const PetscElemScalar *x; 290df311e6cSXuan Zhou PetscElemScalar *z; 291e6dea9dbSXuan Zhou PetscElemScalar one = 1; 292e883f9d5SXuan Zhou 293e883f9d5SXuan Zhou PetscFunctionBegin; 2949566063dSJacob Faibussowitsch if (Y != Z) PetscCall(VecCopy(Y, Z)); 2959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); 2969566063dSJacob Faibussowitsch PetscCall(VecGetArray(Z, (PetscScalar **)&z)); 297e883f9d5SXuan Zhou { /* Scoping so that constructor is called before pointer is returned */ 298e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe, ze; 2990c18141cSBarry Smith xe.LockedAttach(A->rmap->N, 1, *a->grid, 0, 0, x, A->rmap->n); 3000c18141cSBarry Smith ze.Attach(A->cmap->N, 1, *a->grid, 0, 0, z, A->cmap->n); 301e8146892SSatish Balay El::Gemv(El::TRANSPOSE, one, *a->emat, xe, one, ze); 302e883f9d5SXuan Zhou } 3039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x)); 3049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Z, (PetscScalar **)&z)); 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 306e883f9d5SXuan Zhou } 307e883f9d5SXuan Zhou 308d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultNumeric_Elemental(Mat A, Mat B, Mat C) 309d71ae5a4SJacob Faibussowitsch { 310c1d1b975SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 311c1d1b975SXuan Zhou Mat_Elemental *b = (Mat_Elemental *)B->data; 3129a9e8502SHong Zhang Mat_Elemental *c = (Mat_Elemental *)C->data; 313e6dea9dbSXuan Zhou PetscElemScalar one = 1, zero = 0; 314c1d1b975SXuan Zhou 315c1d1b975SXuan Zhou PetscFunctionBegin; 316aae2c449SHong Zhang { /* Scoping so that constructor is called before pointer is returned */ 317e8146892SSatish Balay El::Gemm(El::NORMAL, El::NORMAL, one, *a->emat, *b->emat, zero, *c->emat); 318aae2c449SHong Zhang } 3199a9e8502SHong Zhang C->assembled = PETSC_TRUE; 3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3219a9e8502SHong Zhang } 3229a9e8502SHong Zhang 323cedd07caSPierre Jolivet PetscErrorCode MatMatMultSymbolic_Elemental(Mat A, Mat B, PetscReal, Mat Ce) 324d71ae5a4SJacob Faibussowitsch { 3259a9e8502SHong Zhang PetscFunctionBegin; 3269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Ce, A->rmap->n, B->cmap->n, PETSC_DECIDE, PETSC_DECIDE)); 3279566063dSJacob Faibussowitsch PetscCall(MatSetType(Ce, MATELEMENTAL)); 3289566063dSJacob Faibussowitsch PetscCall(MatSetUp(Ce)); 3294222ddf1SHong Zhang Ce->ops->matmultnumeric = MatMatMultNumeric_Elemental; 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 331c1d1b975SXuan Zhou } 332c1d1b975SXuan Zhou 333d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatTransposeMultNumeric_Elemental(Mat A, Mat B, Mat C) 334d71ae5a4SJacob Faibussowitsch { 335df311e6cSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 336df311e6cSXuan Zhou Mat_Elemental *b = (Mat_Elemental *)B->data; 337df311e6cSXuan Zhou Mat_Elemental *c = (Mat_Elemental *)C->data; 338e6dea9dbSXuan Zhou PetscElemScalar one = 1, zero = 0; 339df311e6cSXuan Zhou 340df311e6cSXuan Zhou PetscFunctionBegin; 341df311e6cSXuan Zhou { /* Scoping so that constructor is called before pointer is returned */ 342e8146892SSatish Balay El::Gemm(El::NORMAL, El::TRANSPOSE, one, *a->emat, *b->emat, zero, *c->emat); 343df311e6cSXuan Zhou } 344df311e6cSXuan Zhou C->assembled = PETSC_TRUE; 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346df311e6cSXuan Zhou } 347df311e6cSXuan Zhou 348cedd07caSPierre Jolivet static PetscErrorCode MatMatTransposeMultSymbolic_Elemental(Mat A, Mat B, PetscReal, Mat C) 349d71ae5a4SJacob Faibussowitsch { 350df311e6cSXuan Zhou PetscFunctionBegin; 3519566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->rmap->n, B->rmap->n, PETSC_DECIDE, PETSC_DECIDE)); 3529566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATELEMENTAL)); 3539566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355df311e6cSXuan Zhou } 356df311e6cSXuan Zhou 357d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Elemental_AB(Mat C) 358d71ae5a4SJacob Faibussowitsch { 3594222ddf1SHong Zhang PetscFunctionBegin; 3604222ddf1SHong Zhang C->ops->matmultsymbolic = MatMatMultSymbolic_Elemental; 3614222ddf1SHong Zhang C->ops->productsymbolic = MatProductSymbolic_AB; 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3634222ddf1SHong Zhang } 3644222ddf1SHong Zhang 365d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Elemental_ABt(Mat C) 366d71ae5a4SJacob Faibussowitsch { 3674222ddf1SHong Zhang PetscFunctionBegin; 3684222ddf1SHong Zhang C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_Elemental; 3694222ddf1SHong Zhang C->ops->productsymbolic = MatProductSymbolic_ABt; 3703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3714222ddf1SHong Zhang } 3724222ddf1SHong Zhang 373d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Elemental(Mat C) 374d71ae5a4SJacob Faibussowitsch { 3754222ddf1SHong Zhang Mat_Product *product = C->product; 3764222ddf1SHong Zhang 3774222ddf1SHong Zhang PetscFunctionBegin; 3784222ddf1SHong Zhang switch (product->type) { 379d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AB: 380d71ae5a4SJacob Faibussowitsch PetscCall(MatProductSetFromOptions_Elemental_AB(C)); 381d71ae5a4SJacob Faibussowitsch break; 382d71ae5a4SJacob Faibussowitsch case MATPRODUCT_ABt: 383d71ae5a4SJacob Faibussowitsch PetscCall(MatProductSetFromOptions_Elemental_ABt(C)); 384d71ae5a4SJacob Faibussowitsch break; 385d71ae5a4SJacob Faibussowitsch default: 386d71ae5a4SJacob Faibussowitsch break; 3874222ddf1SHong Zhang } 3883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3894222ddf1SHong Zhang } 390f6e5db1bSPierre Jolivet 391d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultNumeric_Elemental_MPIDense(Mat A, Mat B, Mat C) 392d71ae5a4SJacob Faibussowitsch { 393f6e5db1bSPierre Jolivet Mat Be, Ce; 394f6e5db1bSPierre Jolivet 395f6e5db1bSPierre Jolivet PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &Be)); 3979566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, Be, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Ce)); 3989566063dSJacob Faibussowitsch PetscCall(MatConvert(Ce, MATMPIDENSE, MAT_REUSE_MATRIX, &C)); 3999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Be)); 4009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ce)); 4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 402f6e5db1bSPierre Jolivet } 403f6e5db1bSPierre Jolivet 404cedd07caSPierre Jolivet PetscErrorCode MatMatMultSymbolic_Elemental_MPIDense(Mat A, Mat B, PetscReal, Mat C) 405d71ae5a4SJacob Faibussowitsch { 406f6e5db1bSPierre Jolivet PetscFunctionBegin; 4079566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, PETSC_DECIDE, PETSC_DECIDE)); 4089566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATMPIDENSE)); 4099566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 410f6e5db1bSPierre Jolivet C->ops->matmultnumeric = MatMatMultNumeric_Elemental_MPIDense; 4113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 412f6e5db1bSPierre Jolivet } 413f6e5db1bSPierre Jolivet 414d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense_AB(Mat C) 415d71ae5a4SJacob Faibussowitsch { 416f6e5db1bSPierre Jolivet PetscFunctionBegin; 417f6e5db1bSPierre Jolivet C->ops->matmultsymbolic = MatMatMultSymbolic_Elemental_MPIDense; 418f6e5db1bSPierre Jolivet C->ops->productsymbolic = MatProductSymbolic_AB; 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 420f6e5db1bSPierre Jolivet } 421f6e5db1bSPierre Jolivet 422d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductSetFromOptions_Elemental_MPIDense(Mat C) 423d71ae5a4SJacob Faibussowitsch { 424f6e5db1bSPierre Jolivet Mat_Product *product = C->product; 425f6e5db1bSPierre Jolivet 426f6e5db1bSPierre Jolivet PetscFunctionBegin; 42748a46eb9SPierre Jolivet if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Elemental_MPIDense_AB(C)); 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 429f6e5db1bSPierre Jolivet } 4304222ddf1SHong Zhang 431d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Elemental(Mat A, Vec D) 432d71ae5a4SJacob Faibussowitsch { 433a9d89745SXuan Zhou PetscInt i, nrows, ncols, nD, rrank, ridx, crank, cidx; 434a9d89745SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 435a9d89745SXuan Zhou PetscElemScalar v; 436ce94432eSBarry Smith MPI_Comm comm; 43761119200SXuan Zhou 43861119200SXuan Zhou PetscFunctionBegin; 4399566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 4409566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &nrows, &ncols)); 441a9d89745SXuan Zhou nD = nrows > ncols ? ncols : nrows; 442a9d89745SXuan Zhou for (i = 0; i < nD; i++) { 443a9d89745SXuan Zhou PetscInt erow, ecol; 444a9d89745SXuan Zhou P2RO(A, 0, i, &rrank, &ridx); 445a9d89745SXuan Zhou RO2E(A, 0, rrank, ridx, &erow); 446aed4548fSBarry Smith PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation"); 447a9d89745SXuan Zhou P2RO(A, 1, i, &crank, &cidx); 448a9d89745SXuan Zhou RO2E(A, 1, crank, cidx, &ecol); 449aed4548fSBarry Smith PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation"); 450a9d89745SXuan Zhou v = a->emat->Get(erow, ecol); 4519566063dSJacob Faibussowitsch PetscCall(VecSetValues(D, 1, &i, (PetscScalar *)&v, INSERT_VALUES)); 452ade3cc5eSXuan Zhou } 4539566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(D)); 4549566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(D)); 4553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45661119200SXuan Zhou } 45761119200SXuan Zhou 458d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDiagonalScale_Elemental(Mat X, Vec L, Vec R) 459d71ae5a4SJacob Faibussowitsch { 460ade3cc5eSXuan Zhou Mat_Elemental *x = (Mat_Elemental *)X->data; 461ade3cc5eSXuan Zhou const PetscElemScalar *d; 462ade3cc5eSXuan Zhou 463ade3cc5eSXuan Zhou PetscFunctionBegin; 4649065cd98SJed Brown if (R) { 4659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(R, (const PetscScalar **)&d)); 466e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> de; 4670c18141cSBarry Smith de.LockedAttach(X->cmap->N, 1, *x->grid, 0, 0, d, X->cmap->n); 468e8146892SSatish Balay El::DiagonalScale(El::RIGHT, El::NORMAL, de, *x->emat); 4699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(R, (const PetscScalar **)&d)); 4709065cd98SJed Brown } 4719065cd98SJed Brown if (L) { 4729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(L, (const PetscScalar **)&d)); 473e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> de; 4740c18141cSBarry Smith de.LockedAttach(X->rmap->N, 1, *x->grid, 0, 0, d, X->rmap->n); 475e8146892SSatish Balay El::DiagonalScale(El::LEFT, El::NORMAL, de, *x->emat); 4769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(L, (const PetscScalar **)&d)); 477ade3cc5eSXuan Zhou } 4783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 479ade3cc5eSXuan Zhou } 480ade3cc5eSXuan Zhou 481cedd07caSPierre Jolivet static PetscErrorCode MatMissingDiagonal_Elemental(Mat, PetscBool *missing, PetscInt *) 482d71ae5a4SJacob Faibussowitsch { 48334fa46e1SFrancesco Ballarin PetscFunctionBegin; 48434fa46e1SFrancesco Ballarin *missing = PETSC_FALSE; 4853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48634fa46e1SFrancesco Ballarin } 48734fa46e1SFrancesco Ballarin 488d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_Elemental(Mat X, PetscScalar a) 489d71ae5a4SJacob Faibussowitsch { 49065b78793SXuan Zhou Mat_Elemental *x = (Mat_Elemental *)X->data; 49165b78793SXuan Zhou 49265b78793SXuan Zhou PetscFunctionBegin; 493e8146892SSatish Balay El::Scale((PetscElemScalar)a, *x->emat); 4943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49565b78793SXuan Zhou } 49665b78793SXuan Zhou 497f82baa17SHong Zhang /* 498f82baa17SHong Zhang MatAXPY - Computes Y = a*X + Y. 499f82baa17SHong Zhang */ 500cedd07caSPierre Jolivet static PetscErrorCode MatAXPY_Elemental(Mat Y, PetscScalar a, Mat X, MatStructure) 501d71ae5a4SJacob Faibussowitsch { 502e09a3074SHong Zhang Mat_Elemental *x = (Mat_Elemental *)X->data; 503e09a3074SHong Zhang Mat_Elemental *y = (Mat_Elemental *)Y->data; 504e09a3074SHong Zhang 505e09a3074SHong Zhang PetscFunctionBegin; 506e8146892SSatish Balay El::Axpy((PetscElemScalar)a, *x->emat, *y->emat); 5079566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)Y)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 509e09a3074SHong Zhang } 510e09a3074SHong Zhang 511cedd07caSPierre Jolivet static PetscErrorCode MatCopy_Elemental(Mat A, Mat B, MatStructure) 512d71ae5a4SJacob Faibussowitsch { 513d6223691SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 514d6223691SXuan Zhou Mat_Elemental *b = (Mat_Elemental *)B->data; 515d6223691SXuan Zhou 516d6223691SXuan Zhou PetscFunctionBegin; 517e8146892SSatish Balay El::Copy(*a->emat, *b->emat); 5189566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)B)); 5193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 520d6223691SXuan Zhou } 521d6223691SXuan Zhou 522d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_Elemental(Mat A, MatDuplicateOption op, Mat *B) 523d71ae5a4SJacob Faibussowitsch { 524df311e6cSXuan Zhou Mat Be; 525ce94432eSBarry Smith MPI_Comm comm; 526df311e6cSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 527df311e6cSXuan Zhou 528df311e6cSXuan Zhou PetscFunctionBegin; 5299566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 5309566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Be)); 5319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Be, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE)); 5329566063dSJacob Faibussowitsch PetscCall(MatSetType(Be, MATELEMENTAL)); 5339566063dSJacob Faibussowitsch PetscCall(MatSetUp(Be)); 534df311e6cSXuan Zhou *B = Be; 535df311e6cSXuan Zhou if (op == MAT_COPY_VALUES) { 536df311e6cSXuan Zhou Mat_Elemental *b = (Mat_Elemental *)Be->data; 537e8146892SSatish Balay El::Copy(*a->emat, *b->emat); 538df311e6cSXuan Zhou } 539df311e6cSXuan Zhou Be->assembled = PETSC_TRUE; 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 541df311e6cSXuan Zhou } 542df311e6cSXuan Zhou 543d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Elemental(Mat A, MatReuse reuse, Mat *B) 544d71ae5a4SJacob Faibussowitsch { 545ec8cb81fSBarry Smith Mat Be = *B; 546ce94432eSBarry Smith MPI_Comm comm; 5474fe7bbcaSHong Zhang Mat_Elemental *a = (Mat_Elemental *)A->data, *b; 548d6223691SXuan Zhou 549d6223691SXuan Zhou PetscFunctionBegin; 5507fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B)); 5519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 5525cb544a0SHong Zhang /* Only out-of-place supported */ 5535f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX, comm, PETSC_ERR_SUP, "Only out-of-place supported"); 5545262d616SXuan Zhou if (reuse == MAT_INITIAL_MATRIX) { 5559566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Be)); 5569566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Be, A->cmap->n, A->rmap->n, PETSC_DECIDE, PETSC_DECIDE)); 5579566063dSJacob Faibussowitsch PetscCall(MatSetType(Be, MATELEMENTAL)); 5589566063dSJacob Faibussowitsch PetscCall(MatSetUp(Be)); 5595262d616SXuan Zhou *B = Be; 5605262d616SXuan Zhou } 5614fe7bbcaSHong Zhang b = (Mat_Elemental *)Be->data; 562e8146892SSatish Balay El::Transpose(*a->emat, *b->emat); 5635262d616SXuan Zhou Be->assembled = PETSC_TRUE; 5643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 565d6223691SXuan Zhou } 566d6223691SXuan Zhou 567d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConjugate_Elemental(Mat A) 568d71ae5a4SJacob Faibussowitsch { 569dfcb0403SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 570dfcb0403SXuan Zhou 571dfcb0403SXuan Zhou PetscFunctionBegin; 572e8146892SSatish Balay El::Conjugate(*a->emat); 5733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 574dfcb0403SXuan Zhou } 575dfcb0403SXuan Zhou 576d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHermitianTranspose_Elemental(Mat A, MatReuse reuse, Mat *B) 577d71ae5a4SJacob Faibussowitsch { 578ec8cb81fSBarry Smith Mat Be = *B; 579ce94432eSBarry Smith MPI_Comm comm; 5804a29722dSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data, *b; 5814a29722dSXuan Zhou 5824a29722dSXuan Zhou PetscFunctionBegin; 5839566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 5845cb544a0SHong Zhang /* Only out-of-place supported */ 5854a29722dSXuan Zhou if (reuse == MAT_INITIAL_MATRIX) { 5869566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Be)); 5879566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Be, A->cmap->n, A->rmap->n, PETSC_DECIDE, PETSC_DECIDE)); 5889566063dSJacob Faibussowitsch PetscCall(MatSetType(Be, MATELEMENTAL)); 5899566063dSJacob Faibussowitsch PetscCall(MatSetUp(Be)); 5904a29722dSXuan Zhou *B = Be; 5914a29722dSXuan Zhou } 5924a29722dSXuan Zhou b = (Mat_Elemental *)Be->data; 593e8146892SSatish Balay El::Adjoint(*a->emat, *b->emat); 5944a29722dSXuan Zhou Be->assembled = PETSC_TRUE; 5953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5964a29722dSXuan Zhou } 5974a29722dSXuan Zhou 598d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_Elemental(Mat A, Vec B, Vec X) 599d71ae5a4SJacob Faibussowitsch { 6001f881ff8SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 601df311e6cSXuan Zhou PetscElemScalar *x; 602ea394475SSatish Balay PetscInt pivoting = a->pivoting; 6031f881ff8SXuan Zhou 6041f881ff8SXuan Zhou PetscFunctionBegin; 6059566063dSJacob Faibussowitsch PetscCall(VecCopy(B, X)); 6069566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, (PetscScalar **)&x)); 607ea394475SSatish Balay 608e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::VC, El::STAR> xe; 6090c18141cSBarry Smith xe.Attach(A->rmap->N, 1, *a->grid, 0, 0, x, A->rmap->n); 610e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::MC, El::MR> xer(xe); 611fc54b460SXuan Zhou switch (A->factortype) { 612fc54b460SXuan Zhou case MAT_FACTOR_LU: 613ea394475SSatish Balay if (pivoting == 0) { 614e8146892SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, xer); 615ea394475SSatish Balay } else if (pivoting == 1) { 616ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, xer); 617ea394475SSatish Balay } else { /* pivoting == 2 */ 618ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *a->Q, xer); 6191f881ff8SXuan Zhou } 620fc54b460SXuan Zhou break; 621d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY: 622d71ae5a4SJacob Faibussowitsch El::cholesky::SolveAfter(El::UPPER, El::NORMAL, *a->emat, xer); 623d71ae5a4SJacob Faibussowitsch break; 624d71ae5a4SJacob Faibussowitsch default: 625d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unfactored Matrix or Unsupported MatFactorType"); 626d71ae5a4SJacob Faibussowitsch break; 6271f881ff8SXuan Zhou } 628ea394475SSatish Balay El::Copy(xer, xe); 629ea394475SSatish Balay 6309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, (PetscScalar **)&x)); 6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6321f881ff8SXuan Zhou } 6331f881ff8SXuan Zhou 634d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveAdd_Elemental(Mat A, Vec B, Vec Y, Vec X) 635d71ae5a4SJacob Faibussowitsch { 636df311e6cSXuan Zhou PetscFunctionBegin; 6379566063dSJacob Faibussowitsch PetscCall(MatSolve_Elemental(A, B, X)); 6389566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, 1, Y)); 6393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 640df311e6cSXuan Zhou } 641df311e6cSXuan Zhou 642d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_Elemental(Mat A, Mat B, Mat X) 643d71ae5a4SJacob Faibussowitsch { 6441f0e42cfSHong Zhang Mat_Elemental *a = (Mat_Elemental *)A->data; 64542ba530cSPierre Jolivet Mat_Elemental *x; 64642ba530cSPierre Jolivet Mat C; 647ea394475SSatish Balay PetscInt pivoting = a->pivoting; 64842ba530cSPierre Jolivet PetscBool flg; 64942ba530cSPierre Jolivet MatType type; 6501f0e42cfSHong Zhang 651ae844d54SHong Zhang PetscFunctionBegin; 6529566063dSJacob Faibussowitsch PetscCall(MatGetType(X, &type)); 6539566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type, MATELEMENTAL, &flg)); 65442ba530cSPierre Jolivet if (!flg) { 6559566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &C)); 65642ba530cSPierre Jolivet x = (Mat_Elemental *)C->data; 65742ba530cSPierre Jolivet } else { 65842ba530cSPierre Jolivet x = (Mat_Elemental *)X->data; 65942ba530cSPierre Jolivet El::Copy(*((Mat_Elemental *)B->data)->emat, *x->emat); 66042ba530cSPierre Jolivet } 661fc54b460SXuan Zhou switch (A->factortype) { 662fc54b460SXuan Zhou case MAT_FACTOR_LU: 663ea394475SSatish Balay if (pivoting == 0) { 664e8146892SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *x->emat); 665ea394475SSatish Balay } else if (pivoting == 1) { 666ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *x->emat); 667ea394475SSatish Balay } else { 668ea394475SSatish Balay El::lu::SolveAfter(El::NORMAL, *a->emat, *a->P, *a->Q, *x->emat); 669d6223691SXuan Zhou } 670fc54b460SXuan Zhou break; 671d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY: 672d71ae5a4SJacob Faibussowitsch El::cholesky::SolveAfter(El::UPPER, El::NORMAL, *a->emat, *x->emat); 673d71ae5a4SJacob Faibussowitsch break; 674d71ae5a4SJacob Faibussowitsch default: 675d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unfactored Matrix or Unsupported MatFactorType"); 676d71ae5a4SJacob Faibussowitsch break; 677fc54b460SXuan Zhou } 67842ba530cSPierre Jolivet if (!flg) { 6799566063dSJacob Faibussowitsch PetscCall(MatConvert(C, type, MAT_REUSE_MATRIX, &X)); 6809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 68142ba530cSPierre Jolivet } 6823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 683ae844d54SHong Zhang } 684ae844d54SHong Zhang 685cedd07caSPierre Jolivet static PetscErrorCode MatLUFactor_Elemental(Mat A, IS, IS, const MatFactorInfo *) 686d71ae5a4SJacob Faibussowitsch { 6877c920d81SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 688ea394475SSatish Balay PetscInt pivoting = a->pivoting; 6897c920d81SXuan Zhou 690ae844d54SHong Zhang PetscFunctionBegin; 691ea394475SSatish Balay if (pivoting == 0) { 692e8146892SSatish Balay El::LU(*a->emat); 693ea394475SSatish Balay } else if (pivoting == 1) { 694ea394475SSatish Balay El::LU(*a->emat, *a->P); 695ea394475SSatish Balay } else { 696ea394475SSatish Balay El::LU(*a->emat, *a->P, *a->Q); 697d6223691SXuan Zhou } 6981293973dSHong Zhang A->factortype = MAT_FACTOR_LU; 699834d3fecSHong Zhang A->assembled = PETSC_TRUE; 700f6224b95SHong Zhang 7019566063dSJacob Faibussowitsch PetscCall(PetscFree(A->solvertype)); 7029566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &A->solvertype)); 7033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 704ae844d54SHong Zhang } 705ae844d54SHong Zhang 706d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_Elemental(Mat F, Mat A, const MatFactorInfo *info) 707d71ae5a4SJacob Faibussowitsch { 708d7c3f9d8SHong Zhang PetscFunctionBegin; 7099566063dSJacob Faibussowitsch PetscCall(MatCopy(A, F, SAME_NONZERO_PATTERN)); 710*f22e26b7SPierre Jolivet PetscCall(MatLUFactor_Elemental(F, nullptr, nullptr, info)); 7113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 712d7c3f9d8SHong Zhang } 713d7c3f9d8SHong Zhang 714cedd07caSPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_Elemental(Mat, Mat, IS, IS, const MatFactorInfo *) 715d71ae5a4SJacob Faibussowitsch { 716d7c3f9d8SHong Zhang PetscFunctionBegin; 717d8304050SJose E. Roman /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */ 7183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 719d7c3f9d8SHong Zhang } 720d7c3f9d8SHong Zhang 721cedd07caSPierre Jolivet static PetscErrorCode MatCholeskyFactor_Elemental(Mat A, IS, const MatFactorInfo *) 722d71ae5a4SJacob Faibussowitsch { 72345cf121fSXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 724e8146892SSatish Balay El::DistMatrix<PetscElemScalar, El::MC, El::STAR> d; 72545cf121fSXuan Zhou 72645cf121fSXuan Zhou PetscFunctionBegin; 727e8146892SSatish Balay El::Cholesky(El::UPPER, *a->emat); 728fc54b460SXuan Zhou A->factortype = MAT_FACTOR_CHOLESKY; 729834d3fecSHong Zhang A->assembled = PETSC_TRUE; 730f6224b95SHong Zhang 7319566063dSJacob Faibussowitsch PetscCall(PetscFree(A->solvertype)); 7329566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &A->solvertype)); 7333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 73445cf121fSXuan Zhou } 73545cf121fSXuan Zhou 736d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F, Mat A, const MatFactorInfo *info) 737d71ae5a4SJacob Faibussowitsch { 738cb76c1d8SXuan Zhou PetscFunctionBegin; 7399566063dSJacob Faibussowitsch PetscCall(MatCopy(A, F, SAME_NONZERO_PATTERN)); 740*f22e26b7SPierre Jolivet PetscCall(MatCholeskyFactor_Elemental(F, nullptr, info)); 7413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74279673f7bSHong Zhang } 743cb76c1d8SXuan Zhou 744cedd07caSPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat, Mat, IS, const MatFactorInfo *) 745d71ae5a4SJacob Faibussowitsch { 74679673f7bSHong Zhang PetscFunctionBegin; 747d8304050SJose E. Roman /* F is created and allocated by MatGetFactor_elemental_petsc(), skip this routine. */ 7483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74979673f7bSHong Zhang } 75079673f7bSHong Zhang 751cedd07caSPierre Jolivet PetscErrorCode MatFactorGetSolverType_elemental_elemental(Mat, MatSolverType *type) 752d71ae5a4SJacob Faibussowitsch { 75315767789SHong Zhang PetscFunctionBegin; 75415767789SHong Zhang *type = MATSOLVERELEMENTAL; 7553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75615767789SHong Zhang } 75715767789SHong Zhang 758d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_elemental_elemental(Mat A, MatFactorType ftype, Mat *F) 759d71ae5a4SJacob Faibussowitsch { 7601293973dSHong Zhang Mat B; 7611293973dSHong Zhang 7621293973dSHong Zhang PetscFunctionBegin; 7631293973dSHong Zhang /* Create the factorization matrix */ 7649566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 7659566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE)); 7669566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATELEMENTAL)); 7679566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 7681293973dSHong Zhang B->factortype = ftype; 76966e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE; 7709566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 7719566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERELEMENTAL, &B->solvertype)); 77200c67f3bSHong Zhang 7739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_elemental_elemental)); 7741293973dSHong Zhang *F = B; 7753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7761293973dSHong Zhang } 7771293973dSHong Zhang 778d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_Elemental(void) 779d71ae5a4SJacob Faibussowitsch { 78042c9c57cSBarry Smith PetscFunctionBegin; 7819566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL, MATELEMENTAL, MAT_FACTOR_LU, MatGetFactor_elemental_elemental)); 7829566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERELEMENTAL, MATELEMENTAL, MAT_FACTOR_CHOLESKY, MatGetFactor_elemental_elemental)); 7833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78442c9c57cSBarry Smith } 78542c9c57cSBarry Smith 786d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNorm_Elemental(Mat A, NormType type, PetscReal *nrm) 787d71ae5a4SJacob Faibussowitsch { 7881f881ff8SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 7891f881ff8SXuan Zhou 7901f881ff8SXuan Zhou PetscFunctionBegin; 7911f881ff8SXuan Zhou switch (type) { 792d71ae5a4SJacob Faibussowitsch case NORM_1: 793d71ae5a4SJacob Faibussowitsch *nrm = El::OneNorm(*a->emat); 794d71ae5a4SJacob Faibussowitsch break; 795d71ae5a4SJacob Faibussowitsch case NORM_FROBENIUS: 796d71ae5a4SJacob Faibussowitsch *nrm = El::FrobeniusNorm(*a->emat); 797d71ae5a4SJacob Faibussowitsch break; 798d71ae5a4SJacob Faibussowitsch case NORM_INFINITY: 799d71ae5a4SJacob Faibussowitsch *nrm = El::InfinityNorm(*a->emat); 800d71ae5a4SJacob Faibussowitsch break; 801d71ae5a4SJacob Faibussowitsch default: 802d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm type"); 8031f881ff8SXuan Zhou } 8043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8051f881ff8SXuan Zhou } 8061f881ff8SXuan Zhou 807d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_Elemental(Mat A) 808d71ae5a4SJacob Faibussowitsch { 8095262d616SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 8105262d616SXuan Zhou 8115262d616SXuan Zhou PetscFunctionBegin; 812e8146892SSatish Balay El::Zero(*a->emat); 8133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8145262d616SXuan Zhou } 8155262d616SXuan Zhou 816d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A, IS *rows, IS *cols) 817d71ae5a4SJacob Faibussowitsch { 818db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data; 819db31f6deSJed Brown PetscInt i, m, shift, stride, *idx; 820db31f6deSJed Brown 821db31f6deSJed Brown PetscFunctionBegin; 822db31f6deSJed Brown if (rows) { 823db31f6deSJed Brown m = a->emat->LocalHeight(); 824db31f6deSJed Brown shift = a->emat->ColShift(); 825db31f6deSJed Brown stride = a->emat->ColStride(); 8269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx)); 827db31f6deSJed Brown for (i = 0; i < m; i++) { 828db31f6deSJed Brown PetscInt rank, offset; 829db31f6deSJed Brown E2RO(A, 0, shift + i * stride, &rank, &offset); 830db31f6deSJed Brown RO2P(A, 0, rank, offset, &idx[i]); 831db31f6deSJed Brown } 8329566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_OWN_POINTER, rows)); 833db31f6deSJed Brown } 834db31f6deSJed Brown if (cols) { 835db31f6deSJed Brown m = a->emat->LocalWidth(); 836db31f6deSJed Brown shift = a->emat->RowShift(); 837db31f6deSJed Brown stride = a->emat->RowStride(); 8389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx)); 839db31f6deSJed Brown for (i = 0; i < m; i++) { 840db31f6deSJed Brown PetscInt rank, offset; 841db31f6deSJed Brown E2RO(A, 1, shift + i * stride, &rank, &offset); 842db31f6deSJed Brown RO2P(A, 1, rank, offset, &idx[i]); 843db31f6deSJed Brown } 8449566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_OWN_POINTER, cols)); 845db31f6deSJed Brown } 8463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 847db31f6deSJed Brown } 848db31f6deSJed Brown 849cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Elemental_Dense(Mat A, MatType, MatReuse reuse, Mat *B) 850d71ae5a4SJacob Faibussowitsch { 8512ef0cf24SXuan Zhou Mat Bmpi; 852af295397SXuan Zhou Mat_Elemental *a = (Mat_Elemental *)A->data; 853ce94432eSBarry Smith MPI_Comm comm; 854fa9e26c8SHong Zhang IS isrows, iscols; 855fa9e26c8SHong Zhang PetscInt rrank, ridx, crank, cidx, nrows, ncols, i, j, erow, ecol, elrow, elcol; 856fa9e26c8SHong Zhang const PetscInt *rows, *cols; 857df311e6cSXuan Zhou PetscElemScalar v; 858fa9e26c8SHong Zhang const El::Grid &grid = a->emat->Grid(); 859af295397SXuan Zhou 860af295397SXuan Zhou PetscFunctionBegin; 8619566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 862a000e53bSHong Zhang 863de0a22f0SHong Zhang if (reuse == MAT_REUSE_MATRIX) { 864de0a22f0SHong Zhang Bmpi = *B; 865de0a22f0SHong Zhang } else { 8669566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Bmpi)); 8679566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Bmpi, A->rmap->n, A->cmap->n, PETSC_DECIDE, PETSC_DECIDE)); 8689566063dSJacob Faibussowitsch PetscCall(MatSetType(Bmpi, MATDENSE)); 8699566063dSJacob Faibussowitsch PetscCall(MatSetUp(Bmpi)); 870de0a22f0SHong Zhang } 871fa9e26c8SHong Zhang 872fa9e26c8SHong Zhang /* Get local entries of A */ 8739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A, &isrows, &iscols)); 8749566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows, &nrows)); 8759566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows, &rows)); 8769566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols, &ncols)); 8779566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols, &cols)); 878fa9e26c8SHong Zhang 87957299f2fSHong Zhang if (a->roworiented) { 8802ef0cf24SXuan Zhou for (i = 0; i < nrows; i++) { 881fa9e26c8SHong Zhang P2RO(A, 0, rows[i], &rrank, &ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */ 8822ef0cf24SXuan Zhou RO2E(A, 0, rrank, ridx, &erow); 883aed4548fSBarry Smith PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation"); 8842ef0cf24SXuan Zhou for (j = 0; j < ncols; j++) { 885fa9e26c8SHong Zhang P2RO(A, 1, cols[j], &crank, &cidx); 8862ef0cf24SXuan Zhou RO2E(A, 1, crank, cidx, &ecol); 887aed4548fSBarry Smith PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation"); 888fa9e26c8SHong Zhang 889fa9e26c8SHong Zhang elrow = erow / grid.MCSize(); /* Elemental local row index */ 890fa9e26c8SHong Zhang elcol = ecol / grid.MRSize(); /* Elemental local column index */ 891fa9e26c8SHong Zhang v = a->emat->GetLocal(elrow, elcol); 8929566063dSJacob Faibussowitsch PetscCall(MatSetValues(Bmpi, 1, &rows[i], 1, &cols[j], (PetscScalar *)&v, INSERT_VALUES)); 8932ef0cf24SXuan Zhou } 8942ef0cf24SXuan Zhou } 89557299f2fSHong Zhang } else { /* column-oriented */ 89657299f2fSHong Zhang for (j = 0; j < ncols; j++) { 89757299f2fSHong Zhang P2RO(A, 1, cols[j], &crank, &cidx); 89857299f2fSHong Zhang RO2E(A, 1, crank, cidx, &ecol); 899aed4548fSBarry Smith PetscCheck(crank >= 0 && cidx >= 0 && ecol >= 0, comm, PETSC_ERR_PLIB, "Incorrect col translation"); 90057299f2fSHong Zhang for (i = 0; i < nrows; i++) { 90157299f2fSHong Zhang P2RO(A, 0, rows[i], &rrank, &ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */ 90257299f2fSHong Zhang RO2E(A, 0, rrank, ridx, &erow); 903aed4548fSBarry Smith PetscCheck(rrank >= 0 && ridx >= 0 && erow >= 0, comm, PETSC_ERR_PLIB, "Incorrect row translation"); 90457299f2fSHong Zhang 90557299f2fSHong Zhang elrow = erow / grid.MCSize(); /* Elemental local row index */ 90657299f2fSHong Zhang elcol = ecol / grid.MRSize(); /* Elemental local column index */ 90757299f2fSHong Zhang v = a->emat->GetLocal(elrow, elcol); 9089566063dSJacob Faibussowitsch PetscCall(MatSetValues(Bmpi, 1, &rows[i], 1, &cols[j], (PetscScalar *)&v, INSERT_VALUES)); 90957299f2fSHong Zhang } 91057299f2fSHong Zhang } 91157299f2fSHong Zhang } 9129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Bmpi, MAT_FINAL_ASSEMBLY)); 9139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Bmpi, MAT_FINAL_ASSEMBLY)); 914511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 9159566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &Bmpi)); 916c4ad791aSXuan Zhou } else { 917c4ad791aSXuan Zhou *B = Bmpi; 918c4ad791aSXuan Zhou } 9199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 9209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 9213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 922af295397SXuan Zhou } 923af295397SXuan Zhou 924cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat) 925d71ae5a4SJacob Faibussowitsch { 926af8000cdSHong Zhang Mat mat_elemental; 927af8000cdSHong Zhang PetscInt M = A->rmap->N, N = A->cmap->N, row, ncols; 928af8000cdSHong Zhang const PetscInt *cols; 929af8000cdSHong Zhang const PetscScalar *vals; 930af8000cdSHong Zhang 931af8000cdSHong Zhang PetscFunctionBegin; 932bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 933bdeae278SStefano Zampini mat_elemental = *newmat; 9349566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental)); 935bdeae278SStefano Zampini } else { 9369566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental)); 9379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N)); 9389566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL)); 9399566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental)); 940bdeae278SStefano Zampini } 941af8000cdSHong Zhang for (row = 0; row < M; row++) { 9429566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals)); 943bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 9449566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES)); 9459566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals)); 946af8000cdSHong Zhang } 9479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY)); 9489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY)); 949af8000cdSHong Zhang 950511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 9519566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental)); 952af8000cdSHong Zhang } else { 953af8000cdSHong Zhang *newmat = mat_elemental; 954af8000cdSHong Zhang } 9553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 956af8000cdSHong Zhang } 957af8000cdSHong Zhang 958cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat) 959d71ae5a4SJacob Faibussowitsch { 9605d7652ecSHong Zhang Mat mat_elemental; 9615d7652ecSHong Zhang PetscInt row, ncols, rstart = A->rmap->rstart, rend = A->rmap->rend, j; 9625d7652ecSHong Zhang const PetscInt *cols; 9635d7652ecSHong Zhang const PetscScalar *vals; 9645d7652ecSHong Zhang 9655d7652ecSHong Zhang PetscFunctionBegin; 966bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 967bdeae278SStefano Zampini mat_elemental = *newmat; 9689566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental)); 969bdeae278SStefano Zampini } else { 9709566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental)); 9719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, A->rmap->N, A->cmap->N)); 9729566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL)); 9739566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental)); 974bdeae278SStefano Zampini } 9755d7652ecSHong Zhang for (row = rstart; row < rend; row++) { 9769566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals)); 9775d7652ecSHong Zhang for (j = 0; j < ncols; j++) { 978bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 9799566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, 1, &cols[j], &vals[j], ADD_VALUES)); 9805d7652ecSHong Zhang } 9819566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals)); 9825d7652ecSHong Zhang } 9839566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY)); 9849566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY)); 9855d7652ecSHong Zhang 986511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 9879566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental)); 9885d7652ecSHong Zhang } else { 9895d7652ecSHong Zhang *newmat = mat_elemental; 9905d7652ecSHong Zhang } 9913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9925d7652ecSHong Zhang } 9935d7652ecSHong Zhang 994cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat) 995d71ae5a4SJacob Faibussowitsch { 9966214f412SHong Zhang Mat mat_elemental; 9976214f412SHong Zhang PetscInt M = A->rmap->N, N = A->cmap->N, row, ncols, j; 9986214f412SHong Zhang const PetscInt *cols; 9996214f412SHong Zhang const PetscScalar *vals; 10006214f412SHong Zhang 10016214f412SHong Zhang PetscFunctionBegin; 1002bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 1003bdeae278SStefano Zampini mat_elemental = *newmat; 10049566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental)); 1005bdeae278SStefano Zampini } else { 10069566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental)); 10079566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N)); 10089566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL)); 10099566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental)); 1010bdeae278SStefano Zampini } 10119566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(A)); 10126214f412SHong Zhang for (row = 0; row < M; row++) { 10139566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals)); 1014bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 10159566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES)); 10166214f412SHong Zhang for (j = 0; j < ncols; j++) { /* lower triangular part */ 1017eb1ec7c1SStefano Zampini PetscScalar v; 10186214f412SHong Zhang if (cols[j] == row) continue; 1019b94d7dedSBarry Smith v = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j]; 10209566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &cols[j], 1, &row, &v, ADD_VALUES)); 10216214f412SHong Zhang } 10229566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals)); 10236214f412SHong Zhang } 10249566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(A)); 10259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY)); 10269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY)); 10276214f412SHong Zhang 1028511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 10299566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental)); 10306214f412SHong Zhang } else { 10316214f412SHong Zhang *newmat = mat_elemental; 10326214f412SHong Zhang } 10333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10346214f412SHong Zhang } 10356214f412SHong Zhang 1036cedd07caSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType, MatReuse reuse, Mat *newmat) 1037d71ae5a4SJacob Faibussowitsch { 10386214f412SHong Zhang Mat mat_elemental; 10396214f412SHong Zhang PetscInt M = A->rmap->N, N = A->cmap->N, row, ncols, j, rstart = A->rmap->rstart, rend = A->rmap->rend; 10406214f412SHong Zhang const PetscInt *cols; 10416214f412SHong Zhang const PetscScalar *vals; 10426214f412SHong Zhang 10436214f412SHong Zhang PetscFunctionBegin; 1044bdeae278SStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 1045bdeae278SStefano Zampini mat_elemental = *newmat; 10469566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat_elemental)); 1047bdeae278SStefano Zampini } else { 10489566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental)); 10499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, M, N)); 10509566063dSJacob Faibussowitsch PetscCall(MatSetType(mat_elemental, MATELEMENTAL)); 10519566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat_elemental)); 1052bdeae278SStefano Zampini } 10539566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(A)); 10546214f412SHong Zhang for (row = rstart; row < rend; row++) { 10559566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals)); 1056bdeae278SStefano Zampini /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 10579566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &row, ncols, cols, vals, ADD_VALUES)); 10586214f412SHong Zhang for (j = 0; j < ncols; j++) { /* lower triangular part */ 1059eb1ec7c1SStefano Zampini PetscScalar v; 10606214f412SHong Zhang if (cols[j] == row) continue; 1061b94d7dedSBarry Smith v = A->hermitian == PETSC_BOOL3_TRUE ? PetscConj(vals[j]) : vals[j]; 10629566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat_elemental, 1, &cols[j], 1, &row, &v, ADD_VALUES)); 10636214f412SHong Zhang } 10649566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals)); 10656214f412SHong Zhang } 10669566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(A)); 10679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY)); 10689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY)); 10696214f412SHong Zhang 1070511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 10719566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &mat_elemental)); 10726214f412SHong Zhang } else { 10736214f412SHong Zhang *newmat = mat_elemental; 10746214f412SHong Zhang } 10753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10766214f412SHong Zhang } 10776214f412SHong Zhang 1078d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Elemental(Mat A) 1079d71ae5a4SJacob Faibussowitsch { 1080db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data; 10815e9f5b67SHong Zhang Mat_Elemental_Grid *commgrid; 10825e9f5b67SHong Zhang PetscBool flg; 10835e9f5b67SHong Zhang MPI_Comm icomm; 1084db31f6deSJed Brown 1085db31f6deSJed Brown PetscFunctionBegin; 1086db31f6deSJed Brown delete a->emat; 1087ea394475SSatish Balay delete a->P; 1088ea394475SSatish Balay delete a->Q; 10895e9f5b67SHong Zhang 1090e8146892SSatish Balay El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A)); 1091*f22e26b7SPierre Jolivet PetscCall(PetscCommDuplicate(cxxcomm.comm, &icomm, nullptr)); 10929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Elemental_keyval, (void **)&commgrid, (int *)&flg)); 10935e9f5b67SHong Zhang if (--commgrid->grid_refct == 0) { 10945e9f5b67SHong Zhang delete commgrid->grid; 10959566063dSJacob Faibussowitsch PetscCall(PetscFree(commgrid)); 10969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Elemental_keyval)); 10975e9f5b67SHong Zhang } 10989566063dSJacob Faibussowitsch PetscCall(PetscCommDestroy(&icomm)); 1099*f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetOwnershipIS_C", nullptr)); 1100*f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", nullptr)); 1101*f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_elemental_mpidense_C", nullptr)); 11029566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 11033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1104db31f6deSJed Brown } 1105db31f6deSJed Brown 1106d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUp_Elemental(Mat A) 1107d71ae5a4SJacob Faibussowitsch { 1108db31f6deSJed Brown Mat_Elemental *a = (Mat_Elemental *)A->data; 1109f19600a0SHong Zhang MPI_Comm comm; 1110db31f6deSJed Brown PetscMPIInt rsize, csize; 1111f19600a0SHong Zhang PetscInt n; 1112db31f6deSJed Brown 1113db31f6deSJed Brown PetscFunctionBegin; 11149566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 11159566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 1116db31f6deSJed Brown 1117d8304050SJose E. Roman /* Check if local row and column sizes are equally distributed. 1118f19600a0SHong Zhang Jed: Elemental uses "element" cyclic ordering so the sizes need to match that 1119f19600a0SHong Zhang exactly. The strategy in MatElemental is for PETSc to implicitly permute to block ordering (like would be returned by 1120d8304050SJose E. Roman PetscSplitOwnership(comm,&n,&N), at which point Elemental matrices can act on PETSc vectors without redistributing the vectors. */ 11219566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 1122f19600a0SHong Zhang n = PETSC_DECIDE; 11239566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(comm, &n, &A->rmap->N)); 11245f80ce2aSJacob 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); 1125f19600a0SHong Zhang 1126f19600a0SHong Zhang n = PETSC_DECIDE; 11279566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(comm, &n, &A->cmap->N)); 11285f80ce2aSJacob 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); 1129f19600a0SHong Zhang 11302f613bf5SBarry Smith a->emat->Resize(A->rmap->N, A->cmap->N); 1131e8146892SSatish Balay El::Zero(*a->emat); 1132db31f6deSJed Brown 11339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(A->rmap->comm, &rsize)); 11349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(A->cmap->comm, &csize)); 11355f80ce2aSJacob Faibussowitsch PetscCheck(csize == rsize, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Cannot use row and column communicators of different sizes"); 1136db31f6deSJed Brown a->commsize = rsize; 11379371c9d4SSatish Balay a->mr[0] = A->rmap->N % rsize; 11389371c9d4SSatish Balay if (!a->mr[0]) a->mr[0] = rsize; 11399371c9d4SSatish Balay a->mr[1] = A->cmap->N % csize; 11409371c9d4SSatish Balay if (!a->mr[1]) a->mr[1] = csize; 1141db31f6deSJed Brown a->m[0] = A->rmap->N / rsize + (a->mr[0] != rsize); 1142db31f6deSJed Brown a->m[1] = A->cmap->N / csize + (a->mr[1] != csize); 11433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1144db31f6deSJed Brown } 1145db31f6deSJed Brown 1146cedd07caSPierre Jolivet PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType) 1147d71ae5a4SJacob Faibussowitsch { 1148aae2c449SHong Zhang Mat_Elemental *a = (Mat_Elemental *)A->data; 1149aae2c449SHong Zhang 1150aae2c449SHong Zhang PetscFunctionBegin; 1151fbbcb730SJack Poulson /* printf("Calling ProcessQueues\n"); */ 1152fbbcb730SJack Poulson a->emat->ProcessQueues(); 1153fbbcb730SJack Poulson /* printf("Finished ProcessQueues\n"); */ 11543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1155aae2c449SHong Zhang } 1156aae2c449SHong Zhang 1157cedd07caSPierre Jolivet PetscErrorCode MatAssemblyEnd_Elemental(Mat, MatAssemblyType) 1158d71ae5a4SJacob Faibussowitsch { 1159aae2c449SHong Zhang PetscFunctionBegin; 1160aae2c449SHong Zhang /* Currently does nothing */ 11613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1162aae2c449SHong Zhang } 1163aae2c449SHong Zhang 1164d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer) 1165d71ae5a4SJacob Faibussowitsch { 11667b30e0f1SHong Zhang Mat Adense, Ae; 11677b30e0f1SHong Zhang MPI_Comm comm; 11687b30e0f1SHong Zhang 11697b30e0f1SHong Zhang PetscFunctionBegin; 11709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)newMat, &comm)); 11719566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &Adense)); 11729566063dSJacob Faibussowitsch PetscCall(MatSetType(Adense, MATDENSE)); 11739566063dSJacob Faibussowitsch PetscCall(MatLoad(Adense, viewer)); 11749566063dSJacob Faibussowitsch PetscCall(MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae)); 11759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adense)); 11769566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(newMat, &Ae)); 11773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11787b30e0f1SHong Zhang } 11797b30e0f1SHong Zhang 11809371c9d4SSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_Elemental, 1181*f22e26b7SPierre Jolivet nullptr, 1182*f22e26b7SPierre Jolivet nullptr, 118340d92e34SHong Zhang MatMult_Elemental, 118440d92e34SHong Zhang /* 4*/ MatMultAdd_Elemental, 11859426833fSXuan Zhou MatMultTranspose_Elemental, 1186e883f9d5SXuan Zhou MatMultTransposeAdd_Elemental, 118740d92e34SHong Zhang MatSolve_Elemental, 1188df311e6cSXuan Zhou MatSolveAdd_Elemental, 1189*f22e26b7SPierre Jolivet nullptr, 1190*f22e26b7SPierre Jolivet /*10*/ nullptr, 119140d92e34SHong Zhang MatLUFactor_Elemental, 119240d92e34SHong Zhang MatCholeskyFactor_Elemental, 1193*f22e26b7SPierre Jolivet nullptr, 119440d92e34SHong Zhang MatTranspose_Elemental, 119540d92e34SHong Zhang /*15*/ MatGetInfo_Elemental, 1196*f22e26b7SPierre Jolivet nullptr, 119761119200SXuan Zhou MatGetDiagonal_Elemental, 1198ade3cc5eSXuan Zhou MatDiagonalScale_Elemental, 119940d92e34SHong Zhang MatNorm_Elemental, 120040d92e34SHong Zhang /*20*/ MatAssemblyBegin_Elemental, 120140d92e34SHong Zhang MatAssemblyEnd_Elemental, 120232d7a744SHong Zhang MatSetOption_Elemental, 120340d92e34SHong Zhang MatZeroEntries_Elemental, 1204*f22e26b7SPierre Jolivet /*24*/ nullptr, 120540d92e34SHong Zhang MatLUFactorSymbolic_Elemental, 120640d92e34SHong Zhang MatLUFactorNumeric_Elemental, 120740d92e34SHong Zhang MatCholeskyFactorSymbolic_Elemental, 120840d92e34SHong Zhang MatCholeskyFactorNumeric_Elemental, 120940d92e34SHong Zhang /*29*/ MatSetUp_Elemental, 1210*f22e26b7SPierre Jolivet nullptr, 1211*f22e26b7SPierre Jolivet nullptr, 1212*f22e26b7SPierre Jolivet nullptr, 1213*f22e26b7SPierre Jolivet nullptr, 1214df311e6cSXuan Zhou /*34*/ MatDuplicate_Elemental, 1215*f22e26b7SPierre Jolivet nullptr, 1216*f22e26b7SPierre Jolivet nullptr, 1217*f22e26b7SPierre Jolivet nullptr, 1218*f22e26b7SPierre Jolivet nullptr, 121940d92e34SHong Zhang /*39*/ MatAXPY_Elemental, 1220*f22e26b7SPierre Jolivet nullptr, 1221*f22e26b7SPierre Jolivet nullptr, 1222*f22e26b7SPierre Jolivet nullptr, 122340d92e34SHong Zhang MatCopy_Elemental, 1224*f22e26b7SPierre Jolivet /*44*/ nullptr, 122540d92e34SHong Zhang MatScale_Elemental, 12267d68702bSBarry Smith MatShift_Basic, 1227*f22e26b7SPierre Jolivet nullptr, 1228*f22e26b7SPierre Jolivet nullptr, 1229*f22e26b7SPierre Jolivet /*49*/ nullptr, 1230*f22e26b7SPierre Jolivet nullptr, 1231*f22e26b7SPierre Jolivet nullptr, 1232*f22e26b7SPierre Jolivet nullptr, 1233*f22e26b7SPierre Jolivet nullptr, 1234*f22e26b7SPierre Jolivet /*54*/ nullptr, 1235*f22e26b7SPierre Jolivet nullptr, 1236*f22e26b7SPierre Jolivet nullptr, 1237*f22e26b7SPierre Jolivet nullptr, 1238*f22e26b7SPierre Jolivet nullptr, 1239*f22e26b7SPierre Jolivet /*59*/ nullptr, 124040d92e34SHong Zhang MatDestroy_Elemental, 124140d92e34SHong Zhang MatView_Elemental, 1242*f22e26b7SPierre Jolivet nullptr, 1243*f22e26b7SPierre Jolivet nullptr, 1244*f22e26b7SPierre Jolivet /*64*/ nullptr, 1245*f22e26b7SPierre Jolivet nullptr, 1246*f22e26b7SPierre Jolivet nullptr, 1247*f22e26b7SPierre Jolivet nullptr, 1248*f22e26b7SPierre Jolivet nullptr, 1249*f22e26b7SPierre Jolivet /*69*/ nullptr, 1250*f22e26b7SPierre Jolivet nullptr, 12512ef0cf24SXuan Zhou MatConvert_Elemental_Dense, 1252*f22e26b7SPierre Jolivet nullptr, 1253*f22e26b7SPierre Jolivet nullptr, 1254*f22e26b7SPierre Jolivet /*74*/ nullptr, 1255*f22e26b7SPierre Jolivet nullptr, 1256*f22e26b7SPierre Jolivet nullptr, 1257*f22e26b7SPierre Jolivet nullptr, 1258*f22e26b7SPierre Jolivet nullptr, 1259*f22e26b7SPierre Jolivet /*79*/ nullptr, 1260*f22e26b7SPierre Jolivet nullptr, 1261*f22e26b7SPierre Jolivet nullptr, 1262*f22e26b7SPierre Jolivet nullptr, 12637b30e0f1SHong Zhang MatLoad_Elemental, 1264*f22e26b7SPierre Jolivet /*84*/ nullptr, 1265*f22e26b7SPierre Jolivet nullptr, 1266*f22e26b7SPierre Jolivet nullptr, 1267*f22e26b7SPierre Jolivet nullptr, 1268*f22e26b7SPierre Jolivet nullptr, 1269*f22e26b7SPierre Jolivet /*89*/ nullptr, 1270*f22e26b7SPierre Jolivet nullptr, 127140d92e34SHong Zhang MatMatMultNumeric_Elemental, 1272*f22e26b7SPierre Jolivet nullptr, 1273*f22e26b7SPierre Jolivet nullptr, 1274*f22e26b7SPierre Jolivet /*94*/ nullptr, 1275*f22e26b7SPierre Jolivet nullptr, 1276*f22e26b7SPierre Jolivet nullptr, 1277df311e6cSXuan Zhou MatMatTransposeMultNumeric_Elemental, 1278*f22e26b7SPierre Jolivet nullptr, 12794222ddf1SHong Zhang /*99*/ MatProductSetFromOptions_Elemental, 1280*f22e26b7SPierre Jolivet nullptr, 1281*f22e26b7SPierre Jolivet nullptr, 1282dfcb0403SXuan Zhou MatConjugate_Elemental, 1283*f22e26b7SPierre Jolivet nullptr, 1284*f22e26b7SPierre Jolivet /*104*/ nullptr, 1285*f22e26b7SPierre Jolivet nullptr, 1286*f22e26b7SPierre Jolivet nullptr, 1287*f22e26b7SPierre Jolivet nullptr, 1288*f22e26b7SPierre Jolivet nullptr, 128940d92e34SHong Zhang /*109*/ MatMatSolve_Elemental, 1290*f22e26b7SPierre Jolivet nullptr, 1291*f22e26b7SPierre Jolivet nullptr, 1292*f22e26b7SPierre Jolivet nullptr, 129334fa46e1SFrancesco Ballarin MatMissingDiagonal_Elemental, 1294*f22e26b7SPierre Jolivet /*114*/ nullptr, 1295*f22e26b7SPierre Jolivet nullptr, 1296*f22e26b7SPierre Jolivet nullptr, 1297*f22e26b7SPierre Jolivet nullptr, 1298*f22e26b7SPierre Jolivet nullptr, 1299*f22e26b7SPierre Jolivet /*119*/ nullptr, 13004a29722dSXuan Zhou MatHermitianTranspose_Elemental, 1301*f22e26b7SPierre Jolivet nullptr, 1302*f22e26b7SPierre Jolivet nullptr, 1303*f22e26b7SPierre Jolivet nullptr, 1304*f22e26b7SPierre Jolivet /*124*/ nullptr, 1305*f22e26b7SPierre Jolivet nullptr, 1306*f22e26b7SPierre Jolivet nullptr, 1307*f22e26b7SPierre Jolivet nullptr, 1308*f22e26b7SPierre Jolivet nullptr, 1309*f22e26b7SPierre Jolivet /*129*/ nullptr, 1310*f22e26b7SPierre Jolivet nullptr, 1311*f22e26b7SPierre Jolivet nullptr, 1312*f22e26b7SPierre Jolivet nullptr, 1313*f22e26b7SPierre Jolivet nullptr, 1314*f22e26b7SPierre Jolivet /*134*/ nullptr, 1315*f22e26b7SPierre Jolivet nullptr, 1316*f22e26b7SPierre Jolivet nullptr, 1317*f22e26b7SPierre Jolivet nullptr, 1318*f22e26b7SPierre Jolivet nullptr, 1319*f22e26b7SPierre Jolivet nullptr, 1320*f22e26b7SPierre Jolivet /*140*/ nullptr, 1321*f22e26b7SPierre Jolivet nullptr, 1322*f22e26b7SPierre Jolivet nullptr, 1323*f22e26b7SPierre Jolivet nullptr, 1324*f22e26b7SPierre Jolivet nullptr, 1325*f22e26b7SPierre Jolivet /*145*/ nullptr, 1326*f22e26b7SPierre Jolivet nullptr, 1327*f22e26b7SPierre Jolivet nullptr, 1328*f22e26b7SPierre Jolivet nullptr, 1329*f22e26b7SPierre Jolivet nullptr, 1330*f22e26b7SPierre Jolivet /*150*/ nullptr, 1331*f22e26b7SPierre Jolivet nullptr}; 133240d92e34SHong Zhang 1333ed36708cSHong Zhang /*MC 1334ed36708cSHong Zhang MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package 1335ed36708cSHong Zhang 1336c2b89b5dSBarry Smith Use ./configure --download-elemental to install PETSc to use Elemental 1337c2b89b5dSBarry Smith 1338ed36708cSHong Zhang Options Database Keys: 13395cc86fc1SJed Brown + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions() 134089bba20eSBarry Smith . -pc_factor_mat_solver_type elemental - to use this direct solver with the option -pc_type lu 13415cc86fc1SJed Brown - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix 1342ed36708cSHong Zhang 1343ed36708cSHong Zhang Level: beginner 1344ed36708cSHong Zhang 134589bba20eSBarry Smith Note: 134689bba20eSBarry Smith Note unlike most matrix formats, this format does not store all the matrix entries for a contiguous 134789bba20eSBarry Smith range of rows on an MPI rank. Use `MatGetOwnershipIS()` to determine what values are stored on 134889bba20eSBarry Smith the given rank. 134989bba20eSBarry Smith 135089bba20eSBarry Smith .seealso: `MATDENSE`, `MATSCALAPACK`, `MatGetOwnershipIS()` 1351ed36708cSHong Zhang M*/ 1352*f22e26b7SPierre Jolivet #if defined(__clang__) 1353*f22e26b7SPierre Jolivet #pragma clang diagnostic push 1354*f22e26b7SPierre Jolivet #pragma clang diagnostic ignored "-Wzero-as-null-pointer-constant" 1355*f22e26b7SPierre Jolivet #endif 1356d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A) 1357d71ae5a4SJacob Faibussowitsch { 1358db31f6deSJed Brown Mat_Elemental *a; 13595682a260SJack Poulson PetscBool flg, flg1; 13605e9f5b67SHong Zhang Mat_Elemental_Grid *commgrid; 13615e9f5b67SHong Zhang MPI_Comm icomm; 13625682a260SJack Poulson PetscInt optv1; 1363db31f6deSJed Brown 1364db31f6deSJed Brown PetscFunctionBegin; 13659566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops, &MatOps_Values, sizeof(struct _MatOps))); 136640d92e34SHong Zhang A->insertmode = NOT_SET_VALUES; 1367db31f6deSJed Brown 13684dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 1369db31f6deSJed Brown A->data = (void *)a; 1370db31f6deSJed Brown 1371db31f6deSJed Brown /* Set up the elemental matrix */ 1372e8146892SSatish Balay El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A)); 13735e9f5b67SHong Zhang 13745e9f5b67SHong Zhang /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */ 13755e9f5b67SHong Zhang if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) { 13769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Elemental_keyval, (void *)0)); 13779566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(ElementalCitation, &ElementalCite)); 13785e9f5b67SHong Zhang } 13799566063dSJacob Faibussowitsch PetscCall(PetscCommDuplicate(cxxcomm.comm, &icomm, NULL)); 13809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Elemental_keyval, (void **)&commgrid, (int *)&flg)); 13815e9f5b67SHong Zhang if (!flg) { 13824dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&commgrid)); 13835cb544a0SHong Zhang 1384d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "Elemental Options", "Mat"); 1385e8146892SSatish Balay /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */ 13869566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_elemental_grid_height", "Grid Height", "None", El::mpi::Size(cxxcomm), &optv1, &flg1)); 13875682a260SJack Poulson if (flg1) { 1388aed4548fSBarry 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)); 1389e8146892SSatish Balay commgrid->grid = new El::Grid(cxxcomm, optv1); /* use user-provided grid height */ 13902ef0cf24SXuan Zhou } else { 1391e8146892SSatish Balay commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */ 1392da0640a4SHong Zhang /* printf("new commgrid->grid = %p\n",commgrid->grid); -- memory leak revealed by valgrind? */ 1393ed667823SXuan Zhou } 13945e9f5b67SHong Zhang commgrid->grid_refct = 1; 13959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_attr(icomm, Petsc_Elemental_keyval, (void *)commgrid)); 1396ea394475SSatish Balay 1397ea394475SSatish Balay a->pivoting = 1; 13989566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_elemental_pivoting", "Pivoting", "None", a->pivoting, &a->pivoting, NULL)); 1399ea394475SSatish Balay 1400d0609cedSBarry Smith PetscOptionsEnd(); 14015e9f5b67SHong Zhang } else { 14025e9f5b67SHong Zhang commgrid->grid_refct++; 14035e9f5b67SHong Zhang } 14049566063dSJacob Faibussowitsch PetscCall(PetscCommDestroy(&icomm)); 14055e9f5b67SHong Zhang a->grid = commgrid->grid; 1406e8146892SSatish Balay a->emat = new El::DistMatrix<PetscElemScalar>(*a->grid); 140732d7a744SHong Zhang a->roworiented = PETSC_TRUE; 1408db31f6deSJed Brown 14099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetOwnershipIS_C", MatGetOwnershipIS_Elemental)); 14109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_elemental_mpidense_C", MatProductSetFromOptions_Elemental_MPIDense)); 14119566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATELEMENTAL)); 14123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1413db31f6deSJed Brown } 1414*f22e26b7SPierre Jolivet #if defined(__clang__) 1415*f22e26b7SPierre Jolivet #pragma clang diagnostic pop 1416*f22e26b7SPierre Jolivet #endif 1417