1a12302e2SJed Brown 2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/ 3a12302e2SJed Brown 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Composite_Nest(DM dm, Mat *J) 5d71ae5a4SJacob Faibussowitsch { 6b989ae6dSJed Brown const DM_Composite *com = (DM_Composite *)dm->data; 7b989ae6dSJed Brown const struct DMCompositeLink *rlink, *clink; 8b989ae6dSJed Brown IS *isg; 9b989ae6dSJed Brown Mat *submats; 10b989ae6dSJed Brown PetscInt i, j, n; 11a12302e2SJed Brown 12a12302e2SJed Brown PetscFunctionBegin; 139ae5db72SJed Brown n = com->nDM; /* Total number of entries */ 14b989ae6dSJed Brown 15b989ae6dSJed Brown /* Explicit index sets are not required for MatCreateNest, but getting them here allows MatNest to do consistency 16b989ae6dSJed Brown * checking and allows ISEqual to compare by identity instead of by contents. */ 179566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, &isg)); 18b989ae6dSJed Brown 19b989ae6dSJed Brown /* Get submatrices */ 209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n * n, &submats)); 21b989ae6dSJed Brown for (i = 0, rlink = com->next; rlink; i++, rlink = rlink->next) { 22b989ae6dSJed Brown for (j = 0, clink = com->next; clink; j++, clink = clink->next) { 230298fd71SBarry Smith Mat sub = NULL; 24b989ae6dSJed Brown if (i == j) { 259566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(rlink->dm, &sub)); 2628b400f6SJacob Faibussowitsch } else PetscCheck(!com->FormCoupleLocations, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot manage off-diagonal parts yet"); 27b989ae6dSJed Brown submats[i * n + j] = sub; 28b989ae6dSJed Brown } 29b989ae6dSJed Brown } 30b989ae6dSJed Brown 319566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)dm), n, isg, n, isg, submats, J)); 32b989ae6dSJed Brown 33b989ae6dSJed Brown /* Disown references */ 349566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&isg[i])); 359566063dSJacob Faibussowitsch PetscCall(PetscFree(isg)); 36b989ae6dSJed Brown 37b989ae6dSJed Brown for (i = 0; i < n * n; i++) { 389566063dSJacob Faibussowitsch if (submats[i]) PetscCall(MatDestroy(&submats[i])); 39b989ae6dSJed Brown } 409566063dSJacob Faibussowitsch PetscCall(PetscFree(submats)); 41*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42a12302e2SJed Brown } 43a12302e2SJed Brown 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Composite_AIJ(DM dm, Mat *J) 45d71ae5a4SJacob Faibussowitsch { 46a12302e2SJed Brown DM_Composite *com = (DM_Composite *)dm->data; 473bf036e2SBarry Smith struct DMCompositeLink *next; 48a12302e2SJed Brown PetscInt m, *dnz, *onz, i, j, mA; 49a12302e2SJed Brown Mat Atmp; 50a12302e2SJed Brown PetscMPIInt rank; 51fcfd50ebSBarry Smith PetscBool dense = PETSC_FALSE; 52a12302e2SJed Brown 53a12302e2SJed Brown PetscFunctionBegin; 54a12302e2SJed Brown /* use global vector to determine layout needed for matrix */ 55a12302e2SJed Brown m = com->n; 56a12302e2SJed Brown 579566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, m, PETSC_DETERMINE, PETSC_DETERMINE)); 599566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, dm->mattype)); 60a12302e2SJed Brown 61a12302e2SJed Brown /* 62a12302e2SJed Brown Extremely inefficient but will compute entire Jacobian for testing 63a12302e2SJed Brown */ 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL)); 65a12302e2SJed Brown if (dense) { 66a12302e2SJed Brown PetscInt rstart, rend, *indices; 67a12302e2SJed Brown PetscScalar *values; 68a12302e2SJed Brown 69a12302e2SJed Brown mA = com->N; 709566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, mA, NULL, mA - m, NULL)); 719566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, mA, NULL)); 72a12302e2SJed Brown 739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*J, &rstart, &rend)); 749566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(mA, &values, mA, &indices)); 759566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, mA)); 76a12302e2SJed Brown for (i = 0; i < mA; i++) indices[i] = i; 7748a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(*J, 1, &i, mA, indices, values, INSERT_VALUES)); 789566063dSJacob Faibussowitsch PetscCall(PetscFree2(values, indices)); 799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY)); 809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY)); 81*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82a12302e2SJed Brown } 83a12302e2SJed Brown 849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 85d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)dm), m, m, dnz, onz); 86a12302e2SJed Brown /* loop over packed objects, handling one at at time */ 87a12302e2SJed Brown next = com->next; 88a12302e2SJed Brown while (next) { 89a12302e2SJed Brown PetscInt nc, rstart, *ccols, maxnc; 90a12302e2SJed Brown const PetscInt *cols, *rstarts; 91a12302e2SJed Brown PetscMPIInt proc; 92a12302e2SJed Brown 939566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(next->dm, &Atmp)); 949566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Atmp, &rstart, NULL)); 959566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Atmp, &rstarts)); 969566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Atmp, &mA, NULL)); 97a12302e2SJed Brown 98a12302e2SJed Brown maxnc = 0; 99a12302e2SJed Brown for (i = 0; i < mA; i++) { 1009566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, NULL, NULL)); 101a12302e2SJed Brown maxnc = PetscMax(nc, maxnc); 1029566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, NULL, NULL)); 103a12302e2SJed Brown } 1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxnc, &ccols)); 105a12302e2SJed Brown for (i = 0; i < mA; i++) { 1069566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, &cols, NULL)); 107a12302e2SJed Brown /* remap the columns taking into how much they are shifted on each process */ 108a12302e2SJed Brown for (j = 0; j < nc; j++) { 109a12302e2SJed Brown proc = 0; 110a12302e2SJed Brown while (cols[j] >= rstarts[proc + 1]) proc++; 111a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 112a12302e2SJed Brown } 1139566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(com->rstart + next->rstart + i, nc, ccols, dnz, onz)); 1149566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, &cols, NULL)); 115a12302e2SJed Brown } 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(ccols)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atmp)); 118a12302e2SJed Brown next = next->next; 119a12302e2SJed Brown } 1201baa6e33SBarry Smith if (com->FormCoupleLocations) PetscCall((*com->FormCoupleLocations)(dm, NULL, dnz, onz, __rstart, __nrows, __start, __end)); 1219566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, 0, dnz, 0, onz)); 1229566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, 0, dnz)); 123d0609cedSBarry Smith MatPreallocateEnd(dnz, onz); 124a12302e2SJed Brown 125*3ba16761SJacob Faibussowitsch if (dm->prealloc_only) PetscFunctionReturn(PETSC_SUCCESS); 126ff6157d0SJed Brown 127a12302e2SJed Brown next = com->next; 128a12302e2SJed Brown while (next) { 129a12302e2SJed Brown PetscInt nc, rstart, row, maxnc, *ccols; 130a12302e2SJed Brown const PetscInt *cols, *rstarts; 131a12302e2SJed Brown const PetscScalar *values; 132a12302e2SJed Brown PetscMPIInt proc; 133a12302e2SJed Brown 1349566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(next->dm, &Atmp)); 1359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Atmp, &rstart, NULL)); 1369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Atmp, &rstarts)); 1379566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Atmp, &mA, NULL)); 138a12302e2SJed Brown maxnc = 0; 139a12302e2SJed Brown for (i = 0; i < mA; i++) { 1409566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, NULL, NULL)); 141a12302e2SJed Brown maxnc = PetscMax(nc, maxnc); 1429566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, NULL, NULL)); 143a12302e2SJed Brown } 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxnc, &ccols)); 145a12302e2SJed Brown for (i = 0; i < mA; i++) { 1469566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, (const PetscInt **)&cols, &values)); 147a12302e2SJed Brown for (j = 0; j < nc; j++) { 148a12302e2SJed Brown proc = 0; 149a12302e2SJed Brown while (cols[j] >= rstarts[proc + 1]) proc++; 150a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 151a12302e2SJed Brown } 152a12302e2SJed Brown row = com->rstart + next->rstart + i; 1539566063dSJacob Faibussowitsch PetscCall(MatSetValues(*J, 1, &row, nc, ccols, values, INSERT_VALUES)); 1549566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, (const PetscInt **)&cols, &values)); 155a12302e2SJed Brown } 1569566063dSJacob Faibussowitsch PetscCall(PetscFree(ccols)); 1579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atmp)); 158a12302e2SJed Brown next = next->next; 159a12302e2SJed Brown } 160a12302e2SJed Brown if (com->FormCoupleLocations) { 161a12302e2SJed Brown PetscInt __rstart; 1629566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*J, &__rstart, NULL)); 1639566063dSJacob Faibussowitsch PetscCall((*com->FormCoupleLocations)(dm, *J, NULL, NULL, __rstart, 0, 0, 0)); 164a12302e2SJed Brown } 1659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY)); 1669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY)); 167*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 168a12302e2SJed Brown } 169a12302e2SJed Brown 170d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateMatrix_Composite(DM dm, Mat *J) 171d71ae5a4SJacob Faibussowitsch { 172a12302e2SJed Brown PetscBool usenest; 17345b6f7e9SBarry Smith ISLocalToGlobalMapping ltogmap; 174a12302e2SJed Brown 175a12302e2SJed Brown PetscFunctionBegin; 1769566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 1779566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 1789566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(dm->mattype, MATNEST, &usenest)); 1791baa6e33SBarry Smith if (usenest) PetscCall(DMCreateMatrix_Composite_Nest(dm, J)); 1801baa6e33SBarry Smith else PetscCall(DMCreateMatrix_Composite_AIJ(dm, J)); 181a12302e2SJed Brown 1829566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dm, <ogmap)); 1839566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(*J, ltogmap, ltogmap)); 1849566063dSJacob Faibussowitsch PetscCall(MatSetDM(*J, dm)); 185*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 186a12302e2SJed Brown } 187