1be1d678aSKris Buschelman 28c79f6d3SBarry Smith /* 38c79f6d3SBarry Smith Support for the parallel AIJ matrix vector multiply 48c79f6d3SBarry Smith */ 5c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 6aeda0f58SHong Zhang #include <petsc/private/vecimpl.h> 7af0996ceSBarry Smith #include <petsc/private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */ 88c79f6d3SBarry Smith 99371c9d4SSatish Balay PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat) { 1044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data; 11ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ *)(aij->B->data); 12b3c64f9dSJunchao Zhang PetscInt i, j, *aj = B->j, *garray; 13b3c64f9dSJunchao Zhang PetscInt ec = 0; /* Number of nonzero external columns */ 141eb62cbbSBarry Smith IS from, to; 151eb62cbbSBarry Smith Vec gvec; 16aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 170f5bd95cSBarry Smith PetscTable gid1_lid1; 180f5bd95cSBarry Smith PetscTablePosition tpos; 19b1d57f15SBarry Smith PetscInt gid, lid; 206f531f54SSatish Balay #else 21d0f46423SBarry Smith PetscInt N = mat->cmap->N, *indices; 222066d6f7SSatish Balay #endif 232066d6f7SSatish Balay 243a40ed3dSBarry Smith PetscFunctionBegin; 254b8d542aSHong Zhang if (!aij->garray) { 2628b400f6SJacob Faibussowitsch PetscCheck(aij->B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing B mat"); 27aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 28c5bfad50SMark F. Adams /* use a table */ 299566063dSJacob Faibussowitsch PetscCall(PetscTableCreate(aij->B->rmap->n, mat->cmap->N + 1, &gid1_lid1)); 30d0f46423SBarry Smith for (i = 0; i < aij->B->rmap->n; i++) { 312066d6f7SSatish Balay for (j = 0; j < B->ilen[i]; j++) { 32b1d57f15SBarry Smith PetscInt data, gid1 = aj[B->i[i] + j] + 1; 339566063dSJacob Faibussowitsch PetscCall(PetscTableFind(gid1_lid1, gid1, &data)); 34fa46199cSSatish Balay if (!data) { 352066d6f7SSatish Balay /* one based table */ 369566063dSJacob Faibussowitsch PetscCall(PetscTableAdd(gid1_lid1, gid1, ++ec, INSERT_VALUES)); 372066d6f7SSatish Balay } 382066d6f7SSatish Balay } 392066d6f7SSatish Balay } 402066d6f7SSatish Balay /* form array of columns we need */ 419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ec, &garray)); 429566063dSJacob Faibussowitsch PetscCall(PetscTableGetHeadPosition(gid1_lid1, &tpos)); 432066d6f7SSatish Balay while (tpos) { 449566063dSJacob Faibussowitsch PetscCall(PetscTableGetNext(gid1_lid1, &tpos, &gid, &lid)); 45b0a32e0cSBarry Smith gid--; 46b0a32e0cSBarry Smith lid--; 472066d6f7SSatish Balay garray[lid] = gid; 482066d6f7SSatish Balay } 499566063dSJacob Faibussowitsch PetscCall(PetscSortInt(ec, garray)); /* sort, and rebuild */ 509566063dSJacob Faibussowitsch PetscCall(PetscTableRemoveAll(gid1_lid1)); 5148a46eb9SPierre Jolivet for (i = 0; i < ec; i++) PetscCall(PetscTableAdd(gid1_lid1, garray[i] + 1, i + 1, INSERT_VALUES)); 522066d6f7SSatish Balay /* compact out the extra columns in B */ 53d0f46423SBarry Smith for (i = 0; i < aij->B->rmap->n; i++) { 542066d6f7SSatish Balay for (j = 0; j < B->ilen[i]; j++) { 55b1d57f15SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 569566063dSJacob Faibussowitsch PetscCall(PetscTableFind(gid1_lid1, gid1, &lid)); 57fa46199cSSatish Balay lid--; 58b3fb0a6cSMatthew Knepley aj[B->i[i] + j] = lid; 592066d6f7SSatish Balay } 602066d6f7SSatish Balay } 619566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&aij->B->cmap)); 629566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B), ec, ec, 1, &aij->B->cmap)); 639566063dSJacob Faibussowitsch PetscCall(PetscTableDestroy(&gid1_lid1)); 642066d6f7SSatish Balay #else 6511285404SBarry Smith /* Make an array as long as the number of columns */ 661eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 679566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &indices)); 68d0f46423SBarry Smith for (i = 0; i < aij->B->rmap->n; i++) { 69d6dfbf8fSBarry Smith for (j = 0; j < B->ilen[i]; j++) { 70bfec09a0SHong Zhang if (!indices[aj[B->i[i] + j]]) ec++; 71bfec09a0SHong Zhang indices[aj[B->i[i] + j]] = 1; 72416022c9SBarry Smith } 731eb62cbbSBarry Smith } 748c79f6d3SBarry Smith 751eb62cbbSBarry Smith /* form array of columns we need */ 769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ec, &garray)); 771eb62cbbSBarry Smith ec = 0; 781eb62cbbSBarry Smith for (i = 0; i < N; i++) { 791eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 801eb62cbbSBarry Smith } 811eb62cbbSBarry Smith 821eb62cbbSBarry Smith /* make indices now point into garray */ 83ad540459SPierre Jolivet for (i = 0; i < ec; i++) indices[garray[i]] = i; 841eb62cbbSBarry Smith 851eb62cbbSBarry Smith /* compact out the extra columns in B */ 86d0f46423SBarry Smith for (i = 0; i < aij->B->rmap->n; i++) { 87ad540459SPierre Jolivet for (j = 0; j < B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 88d6dfbf8fSBarry Smith } 899566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&aij->B->cmap)); 909566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B), ec, ec, 1, &aij->B->cmap)); 919566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 922066d6f7SSatish Balay #endif 934b8d542aSHong Zhang } else { 944b8d542aSHong Zhang garray = aij->garray; 954b8d542aSHong Zhang } 964b8d542aSHong Zhang 974b8d542aSHong Zhang if (!aij->lvec) { 9828b400f6SJacob Faibussowitsch PetscCheck(aij->B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing B mat"); 999566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aij->B, &aij->lvec, NULL)); 1004b8d542aSHong Zhang } 1019566063dSJacob Faibussowitsch PetscCall(VecGetSize(aij->lvec, &ec)); 1021eb62cbbSBarry Smith 103d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 1049566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ec, garray, PETSC_COPY_VALUES, &from)); 1059566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, ec, 0, 1, &to)); 1061eb62cbbSBarry Smith 1071eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 108b5eb4454SBarry Smith /* This does not allocate the array's memory so is efficient */ 1099566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat), 1, mat->cmap->n, mat->cmap->N, NULL, &gvec)); 1101eb62cbbSBarry Smith 1112d336d48SLois Curfman McInnes /* generate the scatter context */ 1129566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&aij->Mvctx)); 1139566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(gvec, from, aij->lvec, to, &aij->Mvctx)); 1149566063dSJacob Faibussowitsch PetscCall(VecScatterViewFromOptions(aij->Mvctx, (PetscObject)mat, "-matmult_vecscatter_view")); 11567bb5161SHong Zhang aij->garray = garray; 1164b8d542aSHong Zhang 1179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 1189566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 1199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gvec)); 1203a40ed3dSBarry Smith PetscFunctionReturn(0); 1218c79f6d3SBarry Smith } 1229e25ed09SBarry Smith 123fff043a9SJunchao Zhang /* Disassemble the off-diag portion of the MPIAIJXxx matrix. 124fff043a9SJunchao Zhang In other words, change the B from reduced format using local col ids 125fff043a9SJunchao Zhang to expanded format using global col ids, which will make it easier to 126fff043a9SJunchao Zhang insert new nonzeros (with global col ids) into the matrix. 127fff043a9SJunchao Zhang The off-diag B determines communication in the matrix vector multiply. 1282493cbb0SBarry Smith */ 1299371c9d4SSatish Balay PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) { 1302493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *)A->data; 1312493cbb0SBarry Smith Mat B = aij->B, Bnew; 132ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ *)B->data; 133*4dfa11a4SJacob Faibussowitsch PetscInt i, j, m = B->rmap->n, n = A->cmap->N, col, ct = 0, *garray = aij->garray, *nz; 13487828ca2SBarry Smith PetscScalar v; 135fff043a9SJunchao Zhang const PetscScalar *ba; 1362493cbb0SBarry Smith 1373a40ed3dSBarry Smith PetscFunctionBegin; 1382493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aij->lvec)); 140464493b3SBarry Smith if (aij->colmap) { 141aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1429566063dSJacob Faibussowitsch PetscCall(PetscTableDestroy(&aij->colmap)); 1432066d6f7SSatish Balay #else 1449566063dSJacob Faibussowitsch PetscCall(PetscFree(aij->colmap)); 1452066d6f7SSatish Balay #endif 146464493b3SBarry Smith } 1472493cbb0SBarry Smith 1482493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 1512493cbb0SBarry Smith 1522493cbb0SBarry Smith /* invent new B and copy stuff over */ 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m + 1, &nz)); 154ad540459SPierre Jolivet for (i = 0; i < m; i++) nz[i] = Baij->i[i + 1] - Baij->i[i]; 1559566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &Bnew)); 1569566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Bnew, m, n, m, n)); /* Bnew now uses A->cmap->N as its col size */ 1579566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(Bnew, A, A)); 1589566063dSJacob Faibussowitsch PetscCall(MatSetType(Bnew, ((PetscObject)B)->type_name)); 1599566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Bnew, 0, nz)); 1602205254eSKarl Rupp 161b38c15b3SStefano Zampini if (Baij->nonew >= 0) { /* Inherit insertion error options (if positive). */ 162b38c15b3SStefano Zampini ((Mat_SeqAIJ *)Bnew->data)->nonew = Baij->nonew; 163b38c15b3SStefano Zampini } 164b38c15b3SStefano Zampini 16577341eacSDmitry Karpeev /* 16677341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 16777341eacSDmitry Karpeev Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 16877341eacSDmitry Karpeev */ 169f69fde56SShane Stafford Bnew->nonzerostate = B->nonzerostate; 1702205254eSKarl Rupp 1719566063dSJacob Faibussowitsch PetscCall(PetscFree(nz)); 1729566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(B, &ba)); 1732493cbb0SBarry Smith for (i = 0; i < m; i++) { 174bfec09a0SHong Zhang for (j = Baij->i[i]; j < Baij->i[i + 1]; j++) { 175bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 176fff043a9SJunchao Zhang v = ba[ct++]; 1779566063dSJacob Faibussowitsch PetscCall(MatSetValues(Bnew, 1, &i, 1, &col, &v, B->insertmode)); 1782493cbb0SBarry Smith } 1792493cbb0SBarry Smith } 1809566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArrayRead(B, &ba)); 181fff043a9SJunchao Zhang 1829566063dSJacob Faibussowitsch PetscCall(PetscFree(aij->garray)); 1839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1842205254eSKarl Rupp 1852493cbb0SBarry Smith aij->B = Bnew; 186227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 1873a40ed3dSBarry Smith PetscFunctionReturn(0); 1882493cbb0SBarry Smith } 1892493cbb0SBarry Smith 1902cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 1912cd6534aSBarry Smith 192f4259b30SLisandro Dalcin static PetscInt *auglyrmapd = NULL, *auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 193f4259b30SLisandro Dalcin static Vec auglydd = NULL, auglyoo = NULL; /* work vectors used to scale the two parts of the local matrix */ 1942cd6534aSBarry Smith 1959371c9d4SSatish Balay PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA, Vec scale) { 1962cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ *)inA->data; /*access private part of matrix */ 197b1d57f15SBarry Smith PetscInt i, n, nt, cstart, cend, no, *garray = ina->garray, *lindices; 198b1d57f15SBarry Smith PetscInt *r_rmapd, *r_rmapo; 1992cd6534aSBarry Smith 2002cd6534aSBarry Smith PetscFunctionBegin; 2019566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(inA, &cstart, &cend)); 2029566063dSJacob Faibussowitsch PetscCall(MatGetSize(ina->A, NULL, &n)); 2039566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(inA->rmap->mapping->n + 1, &r_rmapd)); 2042cd6534aSBarry Smith nt = 0; 205992144d0SBarry Smith for (i = 0; i < inA->rmap->mapping->n; i++) { 206992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2072cd6534aSBarry Smith nt++; 208992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2092cd6534aSBarry Smith } 2102cd6534aSBarry Smith } 21108401ef6SPierre Jolivet PetscCheck(nt == n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Hmm nt %" PetscInt_FMT " n %" PetscInt_FMT, nt, n); 2129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &auglyrmapd)); 213992144d0SBarry Smith for (i = 0; i < inA->rmap->mapping->n; i++) { 214ad540459SPierre Jolivet if (r_rmapd[i]) auglyrmapd[(r_rmapd[i] - 1) - cstart] = i; 2152cd6534aSBarry Smith } 2169566063dSJacob Faibussowitsch PetscCall(PetscFree(r_rmapd)); 2179566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &auglydd)); 2182cd6534aSBarry Smith 2199566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(inA->cmap->N + 1, &lindices)); 220ad540459SPierre Jolivet for (i = 0; i < ina->B->cmap->n; i++) lindices[garray[i]] = i + 1; 221992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 2229566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(inA->rmap->mapping->n + 1, &r_rmapo)); 2232cd6534aSBarry Smith nt = 0; 224992144d0SBarry Smith for (i = 0; i < inA->rmap->mapping->n; i++) { 225992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2262cd6534aSBarry Smith nt++; 227992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2282cd6534aSBarry Smith } 2292cd6534aSBarry Smith } 23008401ef6SPierre Jolivet PetscCheck(nt <= no, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Hmm nt %" PetscInt_FMT " no %" PetscInt_FMT, nt, n); 2319566063dSJacob Faibussowitsch PetscCall(PetscFree(lindices)); 2329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, &auglyrmapo)); 233992144d0SBarry Smith for (i = 0; i < inA->rmap->mapping->n; i++) { 234ad540459SPierre Jolivet if (r_rmapo[i]) auglyrmapo[(r_rmapo[i] - 1)] = i; 2352cd6534aSBarry Smith } 2369566063dSJacob Faibussowitsch PetscCall(PetscFree(r_rmapo)); 2379566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, nt, &auglyoo)); 2382cd6534aSBarry Smith PetscFunctionReturn(0); 2392cd6534aSBarry Smith } 2402cd6534aSBarry Smith 2419371c9d4SSatish Balay PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A, Vec scale) { 24292b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 24392b32695SKris Buschelman 24492b32695SKris Buschelman PetscFunctionBegin; 245cac4c232SBarry Smith PetscTryMethod(A, "MatDiagonalScaleLocal_C", (Mat, Vec), (A, scale)); 24692b32695SKris Buschelman PetscFunctionReturn(0); 24792b32695SKris Buschelman } 24892b32695SKris Buschelman 2499371c9d4SSatish Balay PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A, Vec scale) { 2502cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data; /*access private part of matrix */ 251b1d57f15SBarry Smith PetscInt n, i; 252bca11509SBarry Smith PetscScalar *d, *o; 253bca11509SBarry Smith const PetscScalar *s; 2542cd6534aSBarry Smith 2552cd6534aSBarry Smith PetscFunctionBegin; 25648a46eb9SPierre Jolivet if (!auglyrmapd) PetscCall(MatMPIAIJDiagonalScaleLocalSetUp(A, scale)); 2572cd6534aSBarry Smith 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(scale, &s)); 2592cd6534aSBarry Smith 2609566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(auglydd, &n)); 2619566063dSJacob Faibussowitsch PetscCall(VecGetArray(auglydd, &d)); 2629371c9d4SSatish Balay for (i = 0; i < n; i++) { d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ } 2639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(auglydd, &d)); 2642cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 2659566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(a->A, NULL, auglydd)); 2662cd6534aSBarry Smith 2679566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(auglyoo, &n)); 2689566063dSJacob Faibussowitsch PetscCall(VecGetArray(auglyoo, &o)); 2699371c9d4SSatish Balay for (i = 0; i < n; i++) { o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ } 2709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(scale, &s)); 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(auglyoo, &o)); 2722cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 2739566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(a->B, NULL, auglyoo)); 2742cd6534aSBarry Smith PetscFunctionReturn(0); 2752cd6534aSBarry Smith } 276