1a12302e2SJed Brown 2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/ 3a12302e2SJed Brown 4b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Composite_Nest(DM dm,Mat *J) 5a12302e2SJed Brown { 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. */ 175f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetGlobalISs(dm,&isg)); 18b989ae6dSJed Brown 19b989ae6dSJed Brown /* Get submatrices */ 205f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 255f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(rlink->dm,&sub)); 26*28b400f6SJacob 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 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateNest(PetscObjectComm((PetscObject)dm),n,isg,n,isg,submats,J)); 32b989ae6dSJed Brown 33b989ae6dSJed Brown /* Disown references */ 345f80ce2aSJacob Faibussowitsch for (i=0; i<n; i++) CHKERRQ(ISDestroy(&isg[i])); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(isg)); 36b989ae6dSJed Brown 37b989ae6dSJed Brown for (i=0; i<n*n; i++) { 385f80ce2aSJacob Faibussowitsch if (submats[i]) CHKERRQ(MatDestroy(&submats[i])); 39b989ae6dSJed Brown } 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(submats)); 41a12302e2SJed Brown PetscFunctionReturn(0); 42a12302e2SJed Brown } 43a12302e2SJed Brown 44b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Composite_AIJ(DM dm,Mat *J) 45a12302e2SJed Brown { 46a12302e2SJed Brown PetscErrorCode ierr; 47a12302e2SJed Brown DM_Composite *com = (DM_Composite*)dm->data; 483bf036e2SBarry Smith struct DMCompositeLink *next; 49a12302e2SJed Brown PetscInt m,*dnz,*onz,i,j,mA; 50a12302e2SJed Brown Mat Atmp; 51a12302e2SJed Brown PetscMPIInt rank; 52fcfd50ebSBarry Smith PetscBool dense = PETSC_FALSE; 53a12302e2SJed Brown 54a12302e2SJed Brown PetscFunctionBegin; 55a12302e2SJed Brown /* use global vector to determine layout needed for matrix */ 56a12302e2SJed Brown m = com->n; 57a12302e2SJed Brown 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PetscObjectComm((PetscObject)dm),J)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(*J,dm->mattype)); 61a12302e2SJed Brown 62a12302e2SJed Brown /* 63a12302e2SJed Brown Extremely inefficient but will compute entire Jacobian for testing 64a12302e2SJed Brown */ 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL)); 66a12302e2SJed Brown if (dense) { 67a12302e2SJed Brown PetscInt rstart,rend,*indices; 68a12302e2SJed Brown PetscScalar *values; 69a12302e2SJed Brown 70a12302e2SJed Brown mA = com->N; 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(*J,mA,NULL,mA-m,NULL)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(*J,mA,NULL)); 73a12302e2SJed Brown 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(*J,&rstart,&rend)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(mA,&values,mA,&indices)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(values,mA)); 77a12302e2SJed Brown for (i=0; i<mA; i++) indices[i] = i; 78a12302e2SJed Brown for (i=rstart; i<rend; i++) { 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES)); 80a12302e2SJed Brown } 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(values,indices)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY)); 84a12302e2SJed Brown PetscFunctionReturn(0); 85a12302e2SJed Brown } 86a12302e2SJed Brown 875f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 88ce94432eSBarry Smith ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)dm),m,m,dnz,onz);CHKERRQ(ierr); 89a12302e2SJed Brown /* loop over packed objects, handling one at at time */ 90a12302e2SJed Brown next = com->next; 91a12302e2SJed Brown while (next) { 92a12302e2SJed Brown PetscInt nc,rstart,*ccols,maxnc; 93a12302e2SJed Brown const PetscInt *cols,*rstarts; 94a12302e2SJed Brown PetscMPIInt proc; 95a12302e2SJed Brown 965f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(next->dm,&Atmp)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(Atmp,&rstart,NULL)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRanges(Atmp,&rstarts)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(Atmp,&mA,NULL)); 100a12302e2SJed Brown 101a12302e2SJed Brown maxnc = 0; 102a12302e2SJed Brown for (i=0; i<mA; i++) { 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(Atmp,rstart+i,&nc,NULL,NULL)); 104a12302e2SJed Brown maxnc = PetscMax(nc,maxnc); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(Atmp,rstart+i,&nc,NULL,NULL)); 106a12302e2SJed Brown } 1075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(maxnc,&ccols)); 108a12302e2SJed Brown for (i=0; i<mA; i++) { 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(Atmp,rstart+i,&nc,&cols,NULL)); 110a12302e2SJed Brown /* remap the columns taking into how much they are shifted on each process */ 111a12302e2SJed Brown for (j=0; j<nc; j++) { 112a12302e2SJed Brown proc = 0; 113a12302e2SJed Brown while (cols[j] >= rstarts[proc+1]) proc++; 114a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 115a12302e2SJed Brown } 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(Atmp,rstart+i,&nc,&cols,NULL)); 118a12302e2SJed Brown } 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ccols)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Atmp)); 121a12302e2SJed Brown next = next->next; 122a12302e2SJed Brown } 123a12302e2SJed Brown if (com->FormCoupleLocations) { 1245f80ce2aSJacob Faibussowitsch CHKERRQ((*com->FormCoupleLocations)(dm,NULL,dnz,onz,__rstart,__nrows,__start,__end)); 125a12302e2SJed Brown } 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(*J,0,dnz,0,onz)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(*J,0,dnz)); 128a12302e2SJed Brown ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 129a12302e2SJed Brown 130fcfd50ebSBarry Smith if (dm->prealloc_only) PetscFunctionReturn(0); 131ff6157d0SJed Brown 132a12302e2SJed Brown next = com->next; 133a12302e2SJed Brown while (next) { 134a12302e2SJed Brown PetscInt nc,rstart,row,maxnc,*ccols; 135a12302e2SJed Brown const PetscInt *cols,*rstarts; 136a12302e2SJed Brown const PetscScalar *values; 137a12302e2SJed Brown PetscMPIInt proc; 138a12302e2SJed Brown 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(next->dm,&Atmp)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(Atmp,&rstart,NULL)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRanges(Atmp,&rstarts)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(Atmp,&mA,NULL)); 143a12302e2SJed Brown maxnc = 0; 144a12302e2SJed Brown for (i=0; i<mA; i++) { 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(Atmp,rstart+i,&nc,NULL,NULL)); 146a12302e2SJed Brown maxnc = PetscMax(nc,maxnc); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(Atmp,rstart+i,&nc,NULL,NULL)); 148a12302e2SJed Brown } 1495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(maxnc,&ccols)); 150a12302e2SJed Brown for (i=0; i<mA; i++) { 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values)); 152a12302e2SJed Brown for (j=0; j<nc; j++) { 153a12302e2SJed Brown proc = 0; 154a12302e2SJed Brown while (cols[j] >= rstarts[proc+1]) proc++; 155a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 156a12302e2SJed Brown } 157a12302e2SJed Brown row = com->rstart+next->rstart+i; 1585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values)); 160a12302e2SJed Brown } 1615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ccols)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Atmp)); 163a12302e2SJed Brown next = next->next; 164a12302e2SJed Brown } 165a12302e2SJed Brown if (com->FormCoupleLocations) { 166a12302e2SJed Brown PetscInt __rstart; 1675f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(*J,&__rstart,NULL)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ((*com->FormCoupleLocations)(dm,*J,NULL,NULL,__rstart,0,0,0)); 169a12302e2SJed Brown } 1705f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY)); 172a12302e2SJed Brown PetscFunctionReturn(0); 173a12302e2SJed Brown } 174a12302e2SJed Brown 175b412c318SBarry Smith PetscErrorCode DMCreateMatrix_Composite(DM dm,Mat *J) 176a12302e2SJed Brown { 177a12302e2SJed Brown PetscBool usenest; 17845b6f7e9SBarry Smith ISLocalToGlobalMapping ltogmap; 179a12302e2SJed Brown 180a12302e2SJed Brown PetscFunctionBegin; 1815f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(dm)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(dm->mattype,MATNEST,&usenest)); 184a12302e2SJed Brown if (usenest) { 1855f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix_Composite_Nest(dm,J)); 186a12302e2SJed Brown } else { 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix_Composite_AIJ(dm,J)); 188a12302e2SJed Brown } 189a12302e2SJed Brown 1905f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalToGlobalMapping(dm,<ogmap)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetLocalToGlobalMapping(*J,ltogmap,ltogmap)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetDM(*J,dm)); 193a12302e2SJed Brown PetscFunctionReturn(0); 194a12302e2SJed Brown } 195