1ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/ 2a12302e2SJed Brown 3d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Composite_Nest(DM dm, Mat *J) 4d71ae5a4SJacob Faibussowitsch { 5b989ae6dSJed Brown const DM_Composite *com = (DM_Composite *)dm->data; 6b989ae6dSJed Brown const struct DMCompositeLink *rlink, *clink; 7b989ae6dSJed Brown IS *isg; 8b989ae6dSJed Brown Mat *submats; 9b989ae6dSJed Brown PetscInt i, j, n; 10a12302e2SJed Brown 11a12302e2SJed Brown PetscFunctionBegin; 129ae5db72SJed Brown n = com->nDM; /* Total number of entries */ 13b989ae6dSJed Brown 14b989ae6dSJed Brown /* Explicit index sets are not required for MatCreateNest, but getting them here allows MatNest to do consistency 15b989ae6dSJed Brown * checking and allows ISEqual to compare by identity instead of by contents. */ 169566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, &isg)); 17b989ae6dSJed Brown 18b989ae6dSJed Brown /* Get submatrices */ 199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n * n, &submats)); 20b989ae6dSJed Brown for (i = 0, rlink = com->next; rlink; i++, rlink = rlink->next) { 21b989ae6dSJed Brown for (j = 0, clink = com->next; clink; j++, clink = clink->next) { 220298fd71SBarry Smith Mat sub = NULL; 23b989ae6dSJed Brown if (i == j) { 249566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(rlink->dm, &sub)); 2528b400f6SJacob Faibussowitsch } else PetscCheck(!com->FormCoupleLocations, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot manage off-diagonal parts yet"); 26b989ae6dSJed Brown submats[i * n + j] = sub; 27b989ae6dSJed Brown } 28b989ae6dSJed Brown } 29b989ae6dSJed Brown 309566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)dm), n, isg, n, isg, submats, J)); 31b989ae6dSJed Brown 32b989ae6dSJed Brown /* Disown references */ 339566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&isg[i])); 349566063dSJacob Faibussowitsch PetscCall(PetscFree(isg)); 35b989ae6dSJed Brown 36b989ae6dSJed Brown for (i = 0; i < n * n; i++) { 379566063dSJacob Faibussowitsch if (submats[i]) PetscCall(MatDestroy(&submats[i])); 38b989ae6dSJed Brown } 399566063dSJacob Faibussowitsch PetscCall(PetscFree(submats)); 403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41a12302e2SJed Brown } 42a12302e2SJed Brown 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Composite_AIJ(DM dm, Mat *J) 44d71ae5a4SJacob Faibussowitsch { 45a12302e2SJed Brown DM_Composite *com = (DM_Composite *)dm->data; 463bf036e2SBarry Smith struct DMCompositeLink *next; 47a12302e2SJed Brown PetscInt m, *dnz, *onz, i, j, mA; 48a12302e2SJed Brown Mat Atmp; 49a12302e2SJed Brown PetscMPIInt rank; 50fcfd50ebSBarry Smith PetscBool dense = PETSC_FALSE; 51a12302e2SJed Brown 52a12302e2SJed Brown PetscFunctionBegin; 53a12302e2SJed Brown /* use global vector to determine layout needed for matrix */ 54a12302e2SJed Brown m = com->n; 55a12302e2SJed Brown 569566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 579566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, m, PETSC_DETERMINE, PETSC_DETERMINE)); 589566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, dm->mattype)); 59a12302e2SJed Brown 60a12302e2SJed Brown /* 61a12302e2SJed Brown Extremely inefficient but will compute entire Jacobian for testing 62a12302e2SJed Brown */ 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL)); 64a12302e2SJed Brown if (dense) { 65a12302e2SJed Brown PetscInt rstart, rend, *indices; 66a12302e2SJed Brown PetscScalar *values; 67a12302e2SJed Brown 68a12302e2SJed Brown mA = com->N; 699566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, mA, NULL, mA - m, NULL)); 709566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, mA, NULL)); 71a12302e2SJed Brown 729566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*J, &rstart, &rend)); 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(mA, &values, mA, &indices)); 749566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, mA)); 75a12302e2SJed Brown for (i = 0; i < mA; i++) indices[i] = i; 7648a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(*J, 1, &i, mA, indices, values, INSERT_VALUES)); 779566063dSJacob Faibussowitsch PetscCall(PetscFree2(values, indices)); 789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY)); 799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY)); 803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81a12302e2SJed Brown } 82a12302e2SJed Brown 839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 84d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)dm), m, m, dnz, onz); 85*15229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 86a12302e2SJed Brown next = com->next; 87a12302e2SJed Brown while (next) { 88a12302e2SJed Brown PetscInt nc, rstart, *ccols, maxnc; 89a12302e2SJed Brown const PetscInt *cols, *rstarts; 90a12302e2SJed Brown PetscMPIInt proc; 91a12302e2SJed Brown 929566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(next->dm, &Atmp)); 939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Atmp, &rstart, NULL)); 949566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Atmp, &rstarts)); 959566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Atmp, &mA, NULL)); 96a12302e2SJed Brown 97a12302e2SJed Brown maxnc = 0; 98a12302e2SJed Brown for (i = 0; i < mA; i++) { 999566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, NULL, NULL)); 100a12302e2SJed Brown maxnc = PetscMax(nc, maxnc); 1019566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, NULL, NULL)); 102a12302e2SJed Brown } 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxnc, &ccols)); 104a12302e2SJed Brown for (i = 0; i < mA; i++) { 1059566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, &cols, NULL)); 106a12302e2SJed Brown /* remap the columns taking into how much they are shifted on each process */ 107a12302e2SJed Brown for (j = 0; j < nc; j++) { 108a12302e2SJed Brown proc = 0; 109a12302e2SJed Brown while (cols[j] >= rstarts[proc + 1]) proc++; 110a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 111a12302e2SJed Brown } 1129566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(com->rstart + next->rstart + i, nc, ccols, dnz, onz)); 1139566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, &cols, NULL)); 114a12302e2SJed Brown } 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(ccols)); 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atmp)); 117a12302e2SJed Brown next = next->next; 118a12302e2SJed Brown } 1191baa6e33SBarry Smith if (com->FormCoupleLocations) PetscCall((*com->FormCoupleLocations)(dm, NULL, dnz, onz, __rstart, __nrows, __start, __end)); 1209566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, 0, dnz, 0, onz)); 1219566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, 0, dnz)); 122d0609cedSBarry Smith MatPreallocateEnd(dnz, onz); 123a12302e2SJed Brown 1243ba16761SJacob Faibussowitsch if (dm->prealloc_only) PetscFunctionReturn(PETSC_SUCCESS); 125ff6157d0SJed Brown 126a12302e2SJed Brown next = com->next; 127a12302e2SJed Brown while (next) { 128a12302e2SJed Brown PetscInt nc, rstart, row, maxnc, *ccols; 129a12302e2SJed Brown const PetscInt *cols, *rstarts; 130a12302e2SJed Brown const PetscScalar *values; 131a12302e2SJed Brown PetscMPIInt proc; 132a12302e2SJed Brown 1339566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(next->dm, &Atmp)); 1349566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Atmp, &rstart, NULL)); 1359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Atmp, &rstarts)); 1369566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Atmp, &mA, NULL)); 137a12302e2SJed Brown maxnc = 0; 138a12302e2SJed Brown for (i = 0; i < mA; i++) { 1399566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, NULL, NULL)); 140a12302e2SJed Brown maxnc = PetscMax(nc, maxnc); 1419566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, NULL, NULL)); 142a12302e2SJed Brown } 1439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxnc, &ccols)); 144a12302e2SJed Brown for (i = 0; i < mA; i++) { 1459566063dSJacob Faibussowitsch PetscCall(MatGetRow(Atmp, rstart + i, &nc, (const PetscInt **)&cols, &values)); 146a12302e2SJed Brown for (j = 0; j < nc; j++) { 147a12302e2SJed Brown proc = 0; 148a12302e2SJed Brown while (cols[j] >= rstarts[proc + 1]) proc++; 149a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 150a12302e2SJed Brown } 151a12302e2SJed Brown row = com->rstart + next->rstart + i; 1529566063dSJacob Faibussowitsch PetscCall(MatSetValues(*J, 1, &row, nc, ccols, values, INSERT_VALUES)); 1539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Atmp, rstart + i, &nc, (const PetscInt **)&cols, &values)); 154a12302e2SJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(PetscFree(ccols)); 1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atmp)); 157a12302e2SJed Brown next = next->next; 158a12302e2SJed Brown } 159a12302e2SJed Brown if (com->FormCoupleLocations) { 160a12302e2SJed Brown PetscInt __rstart; 1619566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*J, &__rstart, NULL)); 1629566063dSJacob Faibussowitsch PetscCall((*com->FormCoupleLocations)(dm, *J, NULL, NULL, __rstart, 0, 0, 0)); 163a12302e2SJed Brown } 1649566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY)); 1659566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY)); 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167a12302e2SJed Brown } 168a12302e2SJed Brown 169d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateMatrix_Composite(DM dm, Mat *J) 170d71ae5a4SJacob Faibussowitsch { 171a12302e2SJed Brown PetscBool usenest; 17245b6f7e9SBarry Smith ISLocalToGlobalMapping ltogmap; 173a12302e2SJed Brown 174a12302e2SJed Brown PetscFunctionBegin; 1759566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 1769566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 1779566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(dm->mattype, MATNEST, &usenest)); 1781baa6e33SBarry Smith if (usenest) PetscCall(DMCreateMatrix_Composite_Nest(dm, J)); 1791baa6e33SBarry Smith else PetscCall(DMCreateMatrix_Composite_AIJ(dm, J)); 180a12302e2SJed Brown 1819566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dm, <ogmap)); 1829566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(*J, ltogmap, ltogmap)); 1839566063dSJacob Faibussowitsch PetscCall(MatSetDM(*J, dm)); 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 185a12302e2SJed Brown } 186