1 static char help[] = "Tests various DMComposite routines.\n\n"; 2 3 #include <petscdm.h> 4 #include <petscdmda.h> 5 #include <petscdmcomposite.h> 6 7 int main(int argc, char **argv) 8 { 9 PetscMPIInt rank; 10 DM da1, da2, packer; 11 Vec local, global, globals[2], buffer; 12 PetscScalar value; 13 PetscViewer viewer; 14 15 PetscFunctionBeginUser; 16 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 17 18 PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &packer)); 19 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 8, 1, 1, NULL, &da1)); 20 PetscCall(DMSetFromOptions(da1)); 21 PetscCall(DMSetUp(da1)); 22 PetscCall(DMCompositeAddDM(packer, da1)); 23 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 6, 1, 1, NULL, &da2)); 24 PetscCall(DMSetFromOptions(da2)); 25 PetscCall(DMSetUp(da2)); 26 PetscCall(DMCompositeAddDM(packer, da2)); 27 28 PetscCall(DMCreateGlobalVector(packer, &global)); 29 PetscCall(DMCreateLocalVector(packer, &local)); 30 PetscCall(DMCreateLocalVector(packer, &buffer)); 31 32 PetscCall(DMCompositeGetAccessArray(packer, global, 2, NULL, globals)); 33 value = 1; 34 PetscCall(VecSet(globals[0], value)); 35 value = -1; 36 PetscCall(VecSet(globals[1], value)); 37 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 38 value = rank + 1; 39 PetscCall(VecScale(globals[0], value)); 40 PetscCall(VecScale(globals[1], value)); 41 PetscCall(DMCompositeRestoreAccessArray(packer, global, 2, NULL, globals)); 42 43 /* Test GlobalToLocal in insert mode */ 44 PetscCall(DMGlobalToLocalBegin(packer, global, INSERT_VALUES, local)); 45 PetscCall(DMGlobalToLocalEnd(packer, global, INSERT_VALUES, local)); 46 47 PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 48 PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "\nLocal Vector: processor %d\n", rank)); 49 PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer)); 50 PetscCall(VecView(local, viewer)); 51 PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer)); 52 PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 53 54 /* Test LocalToGlobal in insert mode */ 55 PetscCall(DMLocalToGlobalBegin(packer, local, INSERT_VALUES, global)); 56 PetscCall(DMLocalToGlobalEnd(packer, local, INSERT_VALUES, global)); 57 58 PetscCall(VecView(global, PETSC_VIEWER_STDOUT_WORLD)); 59 60 /* Test LocalToLocal in insert mode */ 61 PetscCall(DMLocalToLocalBegin(packer, local, INSERT_VALUES, buffer)); 62 PetscCall(DMLocalToLocalEnd(packer, local, INSERT_VALUES, buffer)); 63 64 PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 65 PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "\nLocal Vector: processor %d\n", rank)); 66 PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer)); 67 PetscCall(VecView(buffer, viewer)); 68 PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer)); 69 PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 70 71 PetscCall(VecDestroy(&buffer)); 72 PetscCall(VecDestroy(&local)); 73 PetscCall(VecDestroy(&global)); 74 PetscCall(DMDestroy(&packer)); 75 PetscCall(DMDestroy(&da2)); 76 PetscCall(DMDestroy(&da1)); 77 78 PetscCall(PetscFinalize()); 79 return 0; 80 } 81 82 /*TEST 83 84 test: 85 nsize: 3 86 87 TEST*/ 88