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