xref: /petsc/src/dm/tests/ex25.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests DMLocalToGlobal() for dof > 1\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscdm.h>
5*c4762a1bSJed Brown #include <petscdmda.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown int main(int argc,char **argv)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscInt       M = 6,N = 5,P = 4, m = PETSC_DECIDE,n = PETSC_DECIDE,p = PETSC_DECIDE,i,j,k,is,js,ks,in,jen,kn;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   DM             da;
12*c4762a1bSJed Brown   Vec            local,global;
13*c4762a1bSJed Brown   PetscScalar    ****l;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16*c4762a1bSJed Brown   /* Create distributed array and get vectors */
17*c4762a1bSJed Brown   ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,P,m,n,p,2,1,NULL,NULL,NULL,&da);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
19*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
20*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = DMCreateLocalVector(da,&local);CHKERRQ(ierr);
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&is,&js,&ks,&in,&jen,&kn);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = DMDAVecGetArrayDOF(da,local,&l);CHKERRQ(ierr);
25*c4762a1bSJed Brown   for (i=is; i<is+in; i++) {
26*c4762a1bSJed Brown     for (j=js; j<js+jen; j++) {
27*c4762a1bSJed Brown       for (k=ks; k<ks+kn; k++) {
28*c4762a1bSJed Brown         l[k][j][i][0] = 2*(i + j*M + k*M*N);
29*c4762a1bSJed Brown         l[k][j][i][1] = 2*(i + j*M + k*M*N) + 1;
30*c4762a1bSJed Brown       }
31*c4762a1bSJed Brown     }
32*c4762a1bSJed Brown   }
33*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayDOF(da,local,&l);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da,local,ADD_VALUES,global);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(da,local,ADD_VALUES,global);CHKERRQ(ierr);
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown   ierr = VecView(global,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown   /* Free memory */
40*c4762a1bSJed Brown   ierr = VecDestroy(&local);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = VecDestroy(&global);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = PetscFinalize();
44*c4762a1bSJed Brown   return ierr;
45*c4762a1bSJed Brown }
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown /*TEST
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown       test:
52*c4762a1bSJed Brown          filter: grep -v -i Process
53*c4762a1bSJed Brown          output_file: output/ex25_1.out
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown       test:
56*c4762a1bSJed Brown          suffix: 2
57*c4762a1bSJed Brown          nsize: 2
58*c4762a1bSJed Brown          filter: grep -v -i Process
59*c4762a1bSJed Brown          output_file: output/ex25_2.out
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown       test:
62*c4762a1bSJed Brown          suffix: 3
63*c4762a1bSJed Brown          nsize: 3
64*c4762a1bSJed Brown          filter: grep -v -i Process
65*c4762a1bSJed Brown          output_file: output/ex25_2.out
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown       test:
68*c4762a1bSJed Brown          suffix: 4
69*c4762a1bSJed Brown          nsize: 4
70*c4762a1bSJed Brown          filter: grep -v -i Process
71*c4762a1bSJed Brown          output_file: output/ex25_2.out
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown       test:
74*c4762a1bSJed Brown          suffix: 5
75*c4762a1bSJed Brown          nsize: 5
76*c4762a1bSJed Brown          filter: grep -v -i Process
77*c4762a1bSJed Brown          output_file: output/ex25_2.out
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown       test:
80*c4762a1bSJed Brown          suffix: 6
81*c4762a1bSJed Brown          nsize: 6
82*c4762a1bSJed Brown          filter: grep -v -i Process
83*c4762a1bSJed Brown          output_file: output/ex25_2.out
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown TEST*/
86