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