1a12302e2SJed Brown 2*3c48a1e8SJed Brown #include "packimpl.h" /*I "petscdmcomposite.h" I*/ 3a12302e2SJed Brown 4a12302e2SJed Brown #undef __FUNCT__ 5a12302e2SJed Brown #define __FUNCT__ "DMGetMatrix_Composite_Nest" 6a12302e2SJed Brown static PetscErrorCode DMGetMatrix_Composite_Nest(DM dm,const MatType mtype,Mat *J) 7a12302e2SJed Brown { 8b989ae6dSJed Brown const DM_Composite *com = (DM_Composite*)dm->data; 9b989ae6dSJed Brown const struct DMCompositeLink *rlink,*clink; 10b989ae6dSJed Brown PetscErrorCode ierr; 11b989ae6dSJed Brown IS *isg; 12b989ae6dSJed Brown Mat *submats; 13b989ae6dSJed Brown PetscInt i,j,n; 14a12302e2SJed Brown 15a12302e2SJed Brown PetscFunctionBegin; 16b989ae6dSJed Brown n = com->nDM + com->nredundant; /* Total number of entries */ 17b989ae6dSJed Brown 18b989ae6dSJed Brown /* Explicit index sets are not required for MatCreateNest, but getting them here allows MatNest to do consistency 19b989ae6dSJed Brown * checking and allows ISEqual to compare by identity instead of by contents. */ 20b989ae6dSJed Brown ierr = DMCompositeGetGlobalISs(dm,&isg);CHKERRQ(ierr); 21b989ae6dSJed Brown 22b989ae6dSJed Brown /* Get submatrices */ 23b989ae6dSJed Brown ierr = PetscMalloc(n*n*sizeof(Mat),&submats);CHKERRQ(ierr); 24b989ae6dSJed Brown for (i=0,rlink=com->next; rlink; i++,rlink=rlink->next) { 25b989ae6dSJed Brown for (j=0,clink=com->next; clink; j++,clink=clink->next) { 26b989ae6dSJed Brown Mat sub = PETSC_NULL; 27b989ae6dSJed Brown if (i == j) { 28b989ae6dSJed Brown switch (rlink->type) { 29b989ae6dSJed Brown case DMCOMPOSITE_ARRAY: 30b989ae6dSJed Brown ierr = MatCreateMPIDense(((PetscObject)dm)->comm,rlink->n,clink->n,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_NULL,&sub);CHKERRQ(ierr); 31b989ae6dSJed Brown break; 32b989ae6dSJed Brown case DMCOMPOSITE_DM: 33b989ae6dSJed Brown ierr = DMGetMatrix(rlink->dm,PETSC_NULL,&sub);CHKERRQ(ierr); 34b989ae6dSJed Brown break; 35b989ae6dSJed Brown default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"unknown type"); 36b989ae6dSJed Brown } 37b989ae6dSJed Brown } else if (com->FormCoupleLocations) { 38b989ae6dSJed Brown SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot manage off-diagonal parts yet"); 39b989ae6dSJed Brown } 40b989ae6dSJed Brown submats[i*n+j] = sub; 41b989ae6dSJed Brown } 42b989ae6dSJed Brown } 43b989ae6dSJed Brown 44b989ae6dSJed Brown ierr = MatCreateNest(((PetscObject)dm)->comm,n,isg,n,isg,submats,J);CHKERRQ(ierr); 45b989ae6dSJed Brown 46b989ae6dSJed Brown /* Disown references */ 47b989ae6dSJed Brown for (i=0; i<n; i++) {ierr = ISDestroy(isg[i]);CHKERRQ(ierr);} 48b989ae6dSJed Brown ierr = PetscFree(isg);CHKERRQ(ierr); 49b989ae6dSJed Brown 50b989ae6dSJed Brown for (i=0; i<n*n; i++) { 51b989ae6dSJed Brown if (submats[i]) {ierr = MatDestroy(submats[i]);CHKERRQ(ierr);} 52b989ae6dSJed Brown } 53b989ae6dSJed Brown ierr = PetscFree(submats);CHKERRQ(ierr); 54a12302e2SJed Brown PetscFunctionReturn(0); 55a12302e2SJed Brown } 56a12302e2SJed Brown 57a12302e2SJed Brown #undef __FUNCT__ 58a12302e2SJed Brown #define __FUNCT__ "DMGetMatrix_Composite_AIJ" 59a12302e2SJed Brown static PetscErrorCode DMGetMatrix_Composite_AIJ(DM dm,const MatType mtype,Mat *J) 60a12302e2SJed Brown { 61a12302e2SJed Brown PetscErrorCode ierr; 62a12302e2SJed Brown DM_Composite *com = (DM_Composite*)dm->data; 63a12302e2SJed Brown struct DMCompositeLink *next = com->next; 64a12302e2SJed Brown PetscInt m,*dnz,*onz,i,j,mA; 65a12302e2SJed Brown Mat Atmp; 66a12302e2SJed Brown PetscMPIInt rank; 67330e01e6SJed Brown PetscBool dense = PETSC_FALSE,prealloc_only = PETSC_FALSE; 68a12302e2SJed Brown 69a12302e2SJed Brown PetscFunctionBegin; 70a12302e2SJed Brown /* use global vector to determine layout needed for matrix */ 71a12302e2SJed Brown m = com->n; 72a12302e2SJed Brown 73a12302e2SJed Brown ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr); 74a12302e2SJed Brown ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 75a12302e2SJed Brown ierr = MatSetType(*J,mtype);CHKERRQ(ierr); 76a12302e2SJed Brown 77a12302e2SJed Brown /* 78a12302e2SJed Brown Extremely inefficient but will compute entire Jacobian for testing 79a12302e2SJed Brown */ 80ff6157d0SJed Brown ierr = PetscOptionsGetBool(((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr); 81a12302e2SJed Brown if (dense) { 82a12302e2SJed Brown PetscInt rstart,rend,*indices; 83a12302e2SJed Brown PetscScalar *values; 84a12302e2SJed Brown 85a12302e2SJed Brown mA = com->N; 86a12302e2SJed Brown ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr); 87a12302e2SJed Brown ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr); 88a12302e2SJed Brown 89a12302e2SJed Brown ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr); 90a12302e2SJed Brown ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr); 91a12302e2SJed Brown ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr); 92a12302e2SJed Brown for (i=0; i<mA; i++) indices[i] = i; 93a12302e2SJed Brown for (i=rstart; i<rend; i++) { 94a12302e2SJed Brown ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr); 95a12302e2SJed Brown } 96a12302e2SJed Brown ierr = PetscFree2(values,indices);CHKERRQ(ierr); 97a12302e2SJed Brown ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98a12302e2SJed Brown ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99a12302e2SJed Brown PetscFunctionReturn(0); 100a12302e2SJed Brown } 101a12302e2SJed Brown 102a12302e2SJed Brown ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 103a12302e2SJed Brown ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr); 104a12302e2SJed Brown /* loop over packed objects, handling one at at time */ 105a12302e2SJed Brown next = com->next; 106a12302e2SJed Brown while (next) { 107a12302e2SJed Brown if (next->type == DMCOMPOSITE_ARRAY) { 108a12302e2SJed Brown if (rank == next->rank) { /* zero the "little" block */ 109a12302e2SJed Brown for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) { 110a12302e2SJed Brown for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 111a12302e2SJed Brown ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr); 112a12302e2SJed Brown } 113a12302e2SJed Brown } 114a12302e2SJed Brown } 115a12302e2SJed Brown } else if (next->type == DMCOMPOSITE_DM) { 116a12302e2SJed Brown PetscInt nc,rstart,*ccols,maxnc; 117a12302e2SJed Brown const PetscInt *cols,*rstarts; 118a12302e2SJed Brown PetscMPIInt proc; 119a12302e2SJed Brown 120a12302e2SJed Brown ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr); 121a12302e2SJed Brown ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr); 122a12302e2SJed Brown ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 123a12302e2SJed Brown ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr); 124a12302e2SJed Brown 125a12302e2SJed Brown maxnc = 0; 126a12302e2SJed Brown for (i=0; i<mA; i++) { 127a12302e2SJed Brown ierr = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 128a12302e2SJed Brown ierr = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 129a12302e2SJed Brown maxnc = PetscMax(nc,maxnc); 130a12302e2SJed Brown } 131a12302e2SJed Brown ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr); 132a12302e2SJed Brown for (i=0; i<mA; i++) { 133a12302e2SJed Brown ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr); 134a12302e2SJed Brown /* remap the columns taking into how much they are shifted on each process */ 135a12302e2SJed Brown for (j=0; j<nc; j++) { 136a12302e2SJed Brown proc = 0; 137a12302e2SJed Brown while (cols[j] >= rstarts[proc+1]) proc++; 138a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 139a12302e2SJed Brown } 140a12302e2SJed Brown ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr); 141a12302e2SJed Brown ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr); 142a12302e2SJed Brown } 143a12302e2SJed Brown ierr = PetscFree(ccols);CHKERRQ(ierr); 144a12302e2SJed Brown ierr = MatDestroy(Atmp);CHKERRQ(ierr); 145a12302e2SJed Brown } else { 146a12302e2SJed Brown SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 147a12302e2SJed Brown } 148a12302e2SJed Brown next = next->next; 149a12302e2SJed Brown } 150a12302e2SJed Brown if (com->FormCoupleLocations) { 151a12302e2SJed Brown ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr); 152a12302e2SJed Brown } 153a12302e2SJed Brown ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr); 154a12302e2SJed Brown ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr); 155a12302e2SJed Brown ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 156a12302e2SJed Brown 157ff6157d0SJed Brown ierr = PetscOptionsGetBool(((PetscObject)dm)->prefix,"-dmcomposite_preallocate_only",&prealloc_only,PETSC_NULL);CHKERRQ(ierr); 158ff6157d0SJed Brown if (prealloc_only) PetscFunctionReturn(0); 159ff6157d0SJed Brown 160a12302e2SJed Brown next = com->next; 161a12302e2SJed Brown while (next) { 162a12302e2SJed Brown if (next->type == DMCOMPOSITE_ARRAY) { 163a12302e2SJed Brown if (rank == next->rank) { 164a12302e2SJed Brown for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) { 165a12302e2SJed Brown for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 166a12302e2SJed Brown ierr = MatSetValue(*J,j,i,0.0,INSERT_VALUES);CHKERRQ(ierr); 167a12302e2SJed Brown } 168a12302e2SJed Brown } 169a12302e2SJed Brown } 170a12302e2SJed Brown } else if (next->type == DMCOMPOSITE_DM) { 171a12302e2SJed Brown PetscInt nc,rstart,row,maxnc,*ccols; 172a12302e2SJed Brown const PetscInt *cols,*rstarts; 173a12302e2SJed Brown const PetscScalar *values; 174a12302e2SJed Brown PetscMPIInt proc; 175a12302e2SJed Brown 176a12302e2SJed Brown ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr); 177a12302e2SJed Brown ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr); 178a12302e2SJed Brown ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 179a12302e2SJed Brown ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr); 180a12302e2SJed Brown maxnc = 0; 181a12302e2SJed Brown for (i=0; i<mA; i++) { 182a12302e2SJed Brown ierr = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 183a12302e2SJed Brown ierr = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 184a12302e2SJed Brown maxnc = PetscMax(nc,maxnc); 185a12302e2SJed Brown } 186a12302e2SJed Brown ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr); 187a12302e2SJed Brown for (i=0; i<mA; i++) { 188a12302e2SJed Brown ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr); 189a12302e2SJed Brown for (j=0; j<nc; j++) { 190a12302e2SJed Brown proc = 0; 191a12302e2SJed Brown while (cols[j] >= rstarts[proc+1]) proc++; 192a12302e2SJed Brown ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 193a12302e2SJed Brown } 194a12302e2SJed Brown row = com->rstart+next->rstart+i; 195a12302e2SJed Brown ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr); 196a12302e2SJed Brown ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr); 197a12302e2SJed Brown } 198a12302e2SJed Brown ierr = PetscFree(ccols);CHKERRQ(ierr); 199a12302e2SJed Brown ierr = MatDestroy(Atmp);CHKERRQ(ierr); 200a12302e2SJed Brown } else { 201a12302e2SJed Brown SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 202a12302e2SJed Brown } 203a12302e2SJed Brown next = next->next; 204a12302e2SJed Brown } 205a12302e2SJed Brown if (com->FormCoupleLocations) { 206a12302e2SJed Brown PetscInt __rstart; 207a12302e2SJed Brown ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr); 208a12302e2SJed Brown ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr); 209a12302e2SJed Brown } 210a12302e2SJed Brown ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 211a12302e2SJed Brown ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 212a12302e2SJed Brown PetscFunctionReturn(0); 213a12302e2SJed Brown } 214a12302e2SJed Brown 215a12302e2SJed Brown #undef __FUNCT__ 216a12302e2SJed Brown #define __FUNCT__ "DMGetMatrix_Composite" 2177087cfbeSBarry Smith PetscErrorCode DMGetMatrix_Composite(DM dm,const MatType mtype,Mat *J) 218a12302e2SJed Brown { 219a12302e2SJed Brown PetscErrorCode ierr; 220a12302e2SJed Brown PetscBool usenest; 221a12302e2SJed Brown ISLocalToGlobalMapping ltogmap,ltogmapb; 222a12302e2SJed Brown 223a12302e2SJed Brown PetscFunctionBegin; 224a12302e2SJed Brown ierr = PetscStrcmp(mtype,MATNEST,&usenest);CHKERRQ(ierr); 225a12302e2SJed Brown if (usenest) { 226a12302e2SJed Brown ierr = DMGetMatrix_Composite_Nest(dm,mtype,J);CHKERRQ(ierr); 227a12302e2SJed Brown } else { 228a12302e2SJed Brown ierr = DMGetMatrix_Composite_AIJ(dm,mtype?mtype:MATAIJ,J);CHKERRQ(ierr); 229a12302e2SJed Brown } 230a12302e2SJed Brown 231a12302e2SJed Brown ierr = DMGetLocalToGlobalMapping(dm,<ogmap);CHKERRQ(ierr); 232a12302e2SJed Brown ierr = DMGetLocalToGlobalMappingBlock(dm,<ogmapb);CHKERRQ(ierr); 233a12302e2SJed Brown ierr = MatSetLocalToGlobalMapping(*J,ltogmap,ltogmap);CHKERRQ(ierr); 234a12302e2SJed Brown ierr = MatSetLocalToGlobalMappingBlock(*J,ltogmapb,ltogmapb);CHKERRQ(ierr); 235a12302e2SJed Brown PetscFunctionReturn(0); 236a12302e2SJed Brown } 237