/*@ DMDAGetAO - Gets the application ordering context for a distributed array. Collective on DMDA Input Parameter: . da - the distributed array Output Parameters: . ao - the application ordering context for DMDAs Level: intermediate Notes: In this case, the AO maps to the natural grid ordering that would be used for the DMDA if only 1 processor were employed (ordering most rapidly in the x-direction, then y, then z). Multiple degrees of freedom are numbered for each node (rather than 1 component for the whole grid, then the next component, etc.) .keywords: distributed array, get, global, indices, local-to-global .seealso: DMDACreate2d(), DMDASetAOType(), DMDAGetGhostCorners(), DMDAGetCorners(), DMDALocalToGlocal() DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetOwnershipRanges(), AO, AOPetscToApplication(), AOApplicationToPetsc() @*/ PetscErrorCode DMDAGetAO(DM da,AO *ao) { DM_DA *dd; PetscBool isdmda; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(da,DM_CLASSID,1); PetscValidPointer(ao,2); ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isdmda);CHKERRQ(ierr); if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Requires a DMDA as input"); /* now we can safely dereference */ dd = (DM_DA*)da->data; /* Build the natural ordering to PETSc ordering mappings. */ if (!dd->ao) { IS ispetsc,isnatural; PetscErrorCode ierr; PetscInt Nlocal; ierr = DMDAGetNatural_Private(da,&Nlocal,&isnatural);CHKERRQ(ierr); ierr = ISCreateStride(PetscObjectComm((PetscObject)da),Nlocal,dd->base,1,&ispetsc);CHKERRQ(ierr); ierr = AOCreate(PetscObjectComm((PetscObject)da),&dd->ao);CHKERRQ(ierr); ierr = AOSetIS(dd->ao,isnatural,ispetsc);CHKERRQ(ierr); ierr = AOSetType(dd->ao,dd->aotype);CHKERRQ(ierr); ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ao);CHKERRQ(ierr); ierr = ISDestroy(&ispetsc);CHKERRQ(ierr); ierr = ISDestroy(&isnatural);CHKERRQ(ierr); } *ao = dd->ao; PetscFunctionReturn(0); }
/*@C AOCreateBasicIS - Creates a basic application ordering using two index sets. Collective on IS Input Parameters: + isapp - index set that defines an ordering - ispetsc - index set that defines another ordering (may be NULL to use the natural ordering) Output Parameter: . aoout - the new application ordering Level: beginner Notes: the index sets isapp and ispetsc must contain the all the integers 0 to napp-1 (where napp is the length of the index sets) with no duplicates; that is there cannot be any "holes" .keywords: AO, create .seealso: AOCreateBasic(), AODestroy() @*/ PetscErrorCode AOCreateBasicIS(IS isapp,IS ispetsc,AO *aoout) { PetscErrorCode ierr; MPI_Comm comm; AO ao; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject)isapp,&comm);CHKERRQ(ierr); ierr = AOCreate(comm,&ao);CHKERRQ(ierr); ierr = AOSetIS(ao,isapp,ispetsc);CHKERRQ(ierr); ierr = AOSetType(ao,AOBASIC);CHKERRQ(ierr); *aoout = ao; PetscFunctionReturn(0); }
/*@C AOCreateMemoryScalableIS - Creates a memory scalable application ordering using two index sets. Collective on IS Input Parameters: + isapp - index set that defines an ordering - ispetsc - index set that defines another ordering (may be NULL to use the natural ordering) Output Parameter: . aoout - the new application ordering Level: beginner Notes: The index sets isapp and ispetsc must contain the all the integers 0 to napp-1 (where napp is the length of the index sets) with no duplicates; that is there cannot be any "holes". Comparing with AOCreateBasicIS(), this routine trades memory with message communication. .keywords: AO, create .seealso: AOCreateMemoryScalable(), AODestroy() @*/ PetscErrorCode AOCreateMemoryScalableIS(IS isapp,IS ispetsc,AO *aoout) { PetscErrorCode ierr; MPI_Comm comm; AO ao; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject)isapp,&comm);CHKERRQ(ierr); ierr = AOCreate(comm,&ao);CHKERRQ(ierr); ierr = AOSetIS(ao,isapp,ispetsc);CHKERRQ(ierr); ierr = AOSetType(ao,AOMEMORYSCALABLE);CHKERRQ(ierr); ierr = AOViewFromOptions(ao,NULL,"-ao_view");CHKERRQ(ierr); *aoout = ao; PetscFunctionReturn(0); }
int main(int argc,char **argv) { PetscErrorCode ierr; PetscInt i,n = 5; PetscInt getpetsc[] = {0,3,4},getapp[] = {2,1,9,7}; PetscInt getpetsc1[] = {0,3,4},getapp1[] = {2,1,9,7}; PetscInt getpetsc2[] = {0,3,4},getapp2[] = {2,1,9,7}; PetscInt getpetsc3[] = {0,3,4},getapp3[] = {2,1,9,7}; PetscInt getpetsc4[] = {0,3,4},getapp4[] = {2,1,9,7}; PetscMPIInt rank,size; IS ispetsc,isapp; AO ao; const PetscInt *app; ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); /* create the index sets */ ierr = ISCreateStride(PETSC_COMM_WORLD,n,rank,size,&isapp);CHKERRQ(ierr); ierr = ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&ispetsc);CHKERRQ(ierr); /* natural numbering */ /* create the application ordering */ ierr = AOCreateBasicIS(isapp,ispetsc,&ao);CHKERRQ(ierr); ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = AOPetscToApplication(ao,4,getapp);CHKERRQ(ierr); ierr = AOApplicationToPetsc(ao,3,getpetsc);CHKERRQ(ierr); ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 2,1,9,7 PetscToApplication %D %D %D %D\n",rank,getapp[0],getapp[1],getapp[2],getapp[3]);CHKERRQ(ierr); ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %D %D %D\n",rank,getpetsc[0],getpetsc[1],getpetsc[2]);CHKERRQ(ierr); ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); ierr = AODestroy(&ao);CHKERRQ(ierr); /* test MemoryScalable ao */ /*-------------------------*/ ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable: \n");CHKERRQ(ierr); ierr = AOCreateMemoryScalableIS(isapp,ispetsc,&ao);CHKERRQ(ierr);CHKERRQ(ierr); ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = AOPetscToApplication(ao,4,getapp1);CHKERRQ(ierr); ierr = AOApplicationToPetsc(ao,3,getpetsc1);CHKERRQ(ierr); /* Check accuracy */; for (i=0; i<4; i++) { if (getapp1[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp1 %d != getapp %d",getapp1[i],getapp[i]); } for (i=0; i<3; i++) { if (getpetsc1[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc1 %d != getpetsc %d",getpetsc1[i],getpetsc[i]); } ierr = AODestroy(&ao);CHKERRQ(ierr); /* test MemoryScalable ao: ispetsc = NULL */ /*-----------------------------------------------*/ ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest AOCreateMemoryScalable with ispetsc=NULL:\n");CHKERRQ(ierr); ierr = AOCreateMemoryScalableIS(isapp,NULL,&ao);CHKERRQ(ierr);CHKERRQ(ierr); ierr = AOView(ao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = AOPetscToApplication(ao,4,getapp2);CHKERRQ(ierr); ierr = AOApplicationToPetsc(ao,3,getpetsc2);CHKERRQ(ierr); /* Check accuracy */; for (i=0; i<4; i++) { if (getapp2[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp2 %d != getapp %d",getapp2[i],getapp[i]); } for (i=0; i<3; i++) { if (getpetsc2[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc2 %d != getpetsc %d",getpetsc2[i],getpetsc[i]); } ierr = AODestroy(&ao);CHKERRQ(ierr); /* test AOCreateMemoryScalable() ao: */ ierr = ISGetIndices(isapp,&app);CHKERRQ(ierr); ierr = AOCreateMemoryScalable(PETSC_COMM_WORLD,n,app,NULL,&ao);CHKERRQ(ierr); ierr = ISRestoreIndices(isapp,&app);CHKERRQ(ierr); ierr = AOPetscToApplication(ao,4,getapp4);CHKERRQ(ierr); ierr = AOApplicationToPetsc(ao,3,getpetsc4);CHKERRQ(ierr); /* Check accuracy */; for (i=0; i<4; i++) { if (getapp4[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp4 %d != getapp %d",getapp4[i],getapp[i]); } for (i=0; i<3; i++) { if (getpetsc4[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc4 %d != getpetsc %d",getpetsc4[i],getpetsc[i]); } ierr = AODestroy(&ao);CHKERRQ(ierr); /* test general API */ /*------------------*/ ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTest general API: \n");CHKERRQ(ierr); ierr = AOCreate(PETSC_COMM_WORLD,&ao);CHKERRQ(ierr); ierr = AOSetIS(ao,isapp,ispetsc);CHKERRQ(ierr); ierr = AOSetType(ao,AOMEMORYSCALABLE);CHKERRQ(ierr); ierr = AOSetFromOptions(ao);CHKERRQ(ierr); /* ispetsc and isapp are nolonger used. */ ierr = ISDestroy(&ispetsc);CHKERRQ(ierr); ierr = ISDestroy(&isapp);CHKERRQ(ierr); ierr = AOPetscToApplication(ao,4,getapp3);CHKERRQ(ierr); ierr = AOApplicationToPetsc(ao,3,getpetsc3);CHKERRQ(ierr); ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 2,1,9,7 PetscToApplication %D %D %D %D\n",rank,getapp3[0],getapp3[1],getapp3[2],getapp3[3]);CHKERRQ(ierr); ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] 0,3,4 ApplicationToPetsc %D %D %D\n",rank,getpetsc3[0],getpetsc3[1],getpetsc3[2]);CHKERRQ(ierr); ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); /* Check accuracy */; for (i=0; i<4; i++) { if (getapp3[i] != getapp[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getapp3 %d != getapp %d",getapp3[i],getapp[i]); } for (i=0; i<3; i++) { if (getpetsc3[i] != getpetsc[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"getpetsc3 %d != getpetsc %d",getpetsc3[i],getpetsc[i]); } ierr = AODestroy(&ao);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; }