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 PetscErrorCode ierr; 9b989ae6dSJed Brown IS *isg; 10b989ae6dSJed Brown Mat *submats; 11b989ae6dSJed Brown PetscInt i,j,n; 12a12302e2SJed Brown 13a12302e2SJed Brown PetscFunctionBegin; 149ae5db72SJed Brown n = com->nDM; /* Total number of entries */ 15b989ae6dSJed Brown 16b989ae6dSJed Brown /* Explicit index sets are not required for MatCreateNest, but getting them here allows MatNest to do consistency 17b989ae6dSJed Brown * checking and allows ISEqual to compare by identity instead of by contents. */ 18b989ae6dSJed Brown ierr = DMCompositeGetGlobalISs(dm,&isg);CHKERRQ(ierr); 19b989ae6dSJed Brown 20b989ae6dSJed Brown /* Get submatrices */ 21785e854fSJed Brown ierr = PetscMalloc1(n*n,&submats);CHKERRQ(ierr); 22b989ae6dSJed Brown for (i=0,rlink=com->next; rlink; i++,rlink=rlink->next) { 23b989ae6dSJed Brown for (j=0,clink=com->next; clink; j++,clink=clink->next) { 240298fd71SBarry Smith Mat sub = NULL; 25b989ae6dSJed Brown if (i == j) { 26b412c318SBarry Smith ierr = DMCreateMatrix(rlink->dm,&sub);CHKERRQ(ierr); 27*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(com->FormCoupleLocations,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot manage off-diagonal parts yet"); 28b989ae6dSJed Brown submats[i*n+j] = sub; 29b989ae6dSJed Brown } 30b989ae6dSJed Brown } 31b989ae6dSJed Brown 32ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)dm),n,isg,n,isg,submats,J);CHKERRQ(ierr); 33b989ae6dSJed Brown 34b989ae6dSJed Brown /* Disown references */ 35fcfd50ebSBarry Smith for (i=0; i<n; i++) {ierr = ISDestroy(&isg[i]);CHKERRQ(ierr);} 36b989ae6dSJed Brown ierr = PetscFree(isg);CHKERRQ(ierr); 37b989ae6dSJed Brown 38b989ae6dSJed Brown for (i=0; i<n*n; i++) { 39fcfd50ebSBarry Smith if (submats[i]) {ierr = MatDestroy(&submats[i]);CHKERRQ(ierr);} 40b989ae6dSJed Brown } 41b989ae6dSJed Brown ierr = PetscFree(submats);CHKERRQ(ierr); 42a12302e2SJed Brown PetscFunctionReturn(0); 43a12302e2SJed Brown } 44a12302e2SJed Brown 45b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Composite_AIJ(DM dm,Mat *J) 46a12302e2SJed Brown { 47a12302e2SJed Brown PetscErrorCode ierr; 48a12302e2SJed Brown DM_Composite *com = (DM_Composite*)dm->data; 493bf036e2SBarry Smith struct DMCompositeLink *next; 50a12302e2SJed Brown PetscInt m,*dnz,*onz,i,j,mA; 51a12302e2SJed Brown Mat Atmp; 52a12302e2SJed Brown PetscMPIInt rank; 53fcfd50ebSBarry Smith PetscBool dense = PETSC_FALSE; 54a12302e2SJed Brown 55a12302e2SJed Brown PetscFunctionBegin; 56a12302e2SJed Brown /* use global vector to determine layout needed for matrix */ 57a12302e2SJed Brown m = com->n; 58a12302e2SJed Brown 59ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 60a12302e2SJed Brown ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 61b412c318SBarry Smith ierr = MatSetType(*J,dm->mattype);CHKERRQ(ierr); 62a12302e2SJed Brown 63a12302e2SJed Brown /* 64a12302e2SJed Brown Extremely inefficient but will compute entire Jacobian for testing 65a12302e2SJed Brown */ 66c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr); 67a12302e2SJed Brown if (dense) { 68a12302e2SJed Brown PetscInt rstart,rend,*indices; 69a12302e2SJed Brown PetscScalar *values; 70a12302e2SJed Brown 71a12302e2SJed Brown mA = com->N; 720298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,mA,NULL,mA-m,NULL);CHKERRQ(ierr); 730298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,mA,NULL);CHKERRQ(ierr); 74a12302e2SJed Brown 75a12302e2SJed Brown ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr); 76dcca6d9dSJed Brown ierr = PetscMalloc2(mA,&values,mA,&indices);CHKERRQ(ierr); 77580bdb30SBarry Smith ierr = PetscArrayzero(values,mA);CHKERRQ(ierr); 78a12302e2SJed Brown for (i=0; i<mA; i++) indices[i] = i; 79a12302e2SJed Brown for (i=rstart; i<rend; i++) { 80a12302e2SJed Brown ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr); 81a12302e2SJed Brown } 82a12302e2SJed Brown ierr = PetscFree2(values,indices);CHKERRQ(ierr); 83a12302e2SJed Brown ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 84a12302e2SJed Brown ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 85a12302e2SJed Brown PetscFunctionReturn(0); 86a12302e2SJed Brown } 87a12302e2SJed Brown 88ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr); 89ce94432eSBarry Smith ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)dm),m,m,dnz,onz);CHKERRQ(ierr); 90a12302e2SJed Brown /* loop over packed objects, handling one at at time */ 91a12302e2SJed Brown next = com->next; 92a12302e2SJed Brown while (next) { 93a12302e2SJed Brown PetscInt nc,rstart,*ccols,maxnc; 94a12302e2SJed Brown const PetscInt *cols,*rstarts; 95a12302e2SJed Brown PetscMPIInt proc; 96a12302e2SJed Brown 97b412c318SBarry Smith ierr = DMCreateMatrix(next->dm,&Atmp);CHKERRQ(ierr); 980298fd71SBarry Smith ierr = MatGetOwnershipRange(Atmp,&rstart,NULL);CHKERRQ(ierr); 99a12302e2SJed Brown ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 1000298fd71SBarry Smith ierr = MatGetLocalSize(Atmp,&mA,NULL);CHKERRQ(ierr); 101a12302e2SJed Brown 102a12302e2SJed Brown maxnc = 0; 103a12302e2SJed Brown for (i=0; i<mA; i++) { 1040298fd71SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,NULL,NULL);CHKERRQ(ierr); 105a12302e2SJed Brown maxnc = PetscMax(nc,maxnc); 1062e5835c6SStefano Zampini ierr = MatRestoreRow(Atmp,rstart+i,&nc,NULL,NULL);CHKERRQ(ierr); 107a12302e2SJed Brown } 108785e854fSJed Brown ierr = PetscMalloc1(maxnc,&ccols);CHKERRQ(ierr); 109a12302e2SJed Brown for (i=0; i<mA; i++) { 1100298fd71SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,NULL);CHKERRQ(ierr); 111a12302e2SJed Brown /* remap the columns taking into how much they are shifted on each process */ 112a12302e2SJed Brown for (j=0; j<nc; j++) { 113a12302e2SJed Brown proc = 0; 114a12302e2SJed Brown while (cols[j] >= rstarts[proc+1]) proc++; 115a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 116a12302e2SJed Brown } 117a12302e2SJed Brown ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr); 1180298fd71SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,NULL);CHKERRQ(ierr); 119a12302e2SJed Brown } 120a12302e2SJed Brown ierr = PetscFree(ccols);CHKERRQ(ierr); 121fcfd50ebSBarry Smith ierr = MatDestroy(&Atmp);CHKERRQ(ierr); 122a12302e2SJed Brown next = next->next; 123a12302e2SJed Brown } 124a12302e2SJed Brown if (com->FormCoupleLocations) { 1250298fd71SBarry Smith ierr = (*com->FormCoupleLocations)(dm,NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr); 126a12302e2SJed Brown } 127a12302e2SJed Brown ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr); 128a12302e2SJed Brown ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr); 129a12302e2SJed Brown ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 130a12302e2SJed Brown 131fcfd50ebSBarry Smith if (dm->prealloc_only) PetscFunctionReturn(0); 132ff6157d0SJed Brown 133a12302e2SJed Brown next = com->next; 134a12302e2SJed Brown while (next) { 135a12302e2SJed Brown PetscInt nc,rstart,row,maxnc,*ccols; 136a12302e2SJed Brown const PetscInt *cols,*rstarts; 137a12302e2SJed Brown const PetscScalar *values; 138a12302e2SJed Brown PetscMPIInt proc; 139a12302e2SJed Brown 140b412c318SBarry Smith ierr = DMCreateMatrix(next->dm,&Atmp);CHKERRQ(ierr); 1410298fd71SBarry Smith ierr = MatGetOwnershipRange(Atmp,&rstart,NULL);CHKERRQ(ierr); 142a12302e2SJed Brown ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 1430298fd71SBarry Smith ierr = MatGetLocalSize(Atmp,&mA,NULL);CHKERRQ(ierr); 144a12302e2SJed Brown maxnc = 0; 145a12302e2SJed Brown for (i=0; i<mA; i++) { 1460298fd71SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,NULL,NULL);CHKERRQ(ierr); 147a12302e2SJed Brown maxnc = PetscMax(nc,maxnc); 14818ff422cSJed Brown ierr = MatRestoreRow(Atmp,rstart+i,&nc,NULL,NULL);CHKERRQ(ierr); 149a12302e2SJed Brown } 150785e854fSJed Brown ierr = PetscMalloc1(maxnc,&ccols);CHKERRQ(ierr); 151a12302e2SJed Brown for (i=0; i<mA; i++) { 152a12302e2SJed Brown ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values);CHKERRQ(ierr); 153a12302e2SJed Brown for (j=0; j<nc; j++) { 154a12302e2SJed Brown proc = 0; 155a12302e2SJed Brown while (cols[j] >= rstarts[proc+1]) proc++; 156a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 157a12302e2SJed Brown } 158a12302e2SJed Brown row = com->rstart+next->rstart+i; 159a12302e2SJed Brown ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr); 160a12302e2SJed Brown ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values);CHKERRQ(ierr); 161a12302e2SJed Brown } 162a12302e2SJed Brown ierr = PetscFree(ccols);CHKERRQ(ierr); 163fcfd50ebSBarry Smith ierr = MatDestroy(&Atmp);CHKERRQ(ierr); 164a12302e2SJed Brown next = next->next; 165a12302e2SJed Brown } 166a12302e2SJed Brown if (com->FormCoupleLocations) { 167a12302e2SJed Brown PetscInt __rstart; 1680298fd71SBarry Smith ierr = MatGetOwnershipRange(*J,&__rstart,NULL);CHKERRQ(ierr); 1690298fd71SBarry Smith ierr = (*com->FormCoupleLocations)(dm,*J,NULL,NULL,__rstart,0,0,0);CHKERRQ(ierr); 170a12302e2SJed Brown } 171a12302e2SJed Brown ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172a12302e2SJed Brown ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173a12302e2SJed Brown PetscFunctionReturn(0); 174a12302e2SJed Brown } 175a12302e2SJed Brown 176b412c318SBarry Smith PetscErrorCode DMCreateMatrix_Composite(DM dm,Mat *J) 177a12302e2SJed Brown { 178a12302e2SJed Brown PetscErrorCode ierr; 179a12302e2SJed Brown PetscBool usenest; 18045b6f7e9SBarry Smith ISLocalToGlobalMapping ltogmap; 181a12302e2SJed Brown 182a12302e2SJed Brown PetscFunctionBegin; 18305ec3129SRichard Tran Mills ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 184f692024eSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 185b412c318SBarry Smith ierr = PetscStrcmp(dm->mattype,MATNEST,&usenest);CHKERRQ(ierr); 186a12302e2SJed Brown if (usenest) { 187b412c318SBarry Smith ierr = DMCreateMatrix_Composite_Nest(dm,J);CHKERRQ(ierr); 188a12302e2SJed Brown } else { 189b412c318SBarry Smith ierr = DMCreateMatrix_Composite_AIJ(dm,J);CHKERRQ(ierr); 190a12302e2SJed Brown } 191a12302e2SJed Brown 192a12302e2SJed Brown ierr = DMGetLocalToGlobalMapping(dm,<ogmap);CHKERRQ(ierr); 193a12302e2SJed Brown ierr = MatSetLocalToGlobalMapping(*J,ltogmap,ltogmap);CHKERRQ(ierr); 194c6b011d8SStefano Zampini ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 195a12302e2SJed Brown PetscFunctionReturn(0); 196a12302e2SJed Brown } 197