Actual source code: ex21.c
petsc-3.5.4 2015-05-23
1: static const char help[] = "Test DMCreateInjection() for mapping coordinates in 3D";
3: #include <petscvec.h>
4: #include <petscmat.h>
5: #include <petscdm.h>
6: #include <petscdmda.h>
10: PetscErrorCode test1_DAInjection3d(PetscInt mx, PetscInt my, PetscInt mz)
11: {
12: PetscErrorCode ierr;
13: DM dac,daf;
14: PetscViewer vv;
15: Vec ac,af;
16: PetscInt periodicity;
17: DMBoundaryType bx,by,bz;
20: bx = DM_BOUNDARY_NONE;
21: by = DM_BOUNDARY_NONE;
22: bz = DM_BOUNDARY_NONE;
24: periodicity = 0;
26: PetscOptionsGetInt(NULL,"-periodic", &periodicity, NULL);
27: if (periodicity==1) {
28: bx = DM_BOUNDARY_PERIODIC;
29: } else if (periodicity==2) {
30: by = DM_BOUNDARY_PERIODIC;
31: } else if (periodicity==3) {
32: bz = DM_BOUNDARY_PERIODIC;
33: }
35: DMDACreate3d(PETSC_COMM_WORLD, bx,by,bz, DMDA_STENCIL_BOX,
36: mx+1, my+1,mz+1,
37: PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE,
38: 1, /* 1 dof */
39: 1, /* stencil = 1 */
40: NULL,NULL,NULL,
41: &daf);
43: DMSetFromOptions(daf);
45: DMCoarsen(daf,MPI_COMM_NULL,&dac);
47: DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, -1.0,1.0);
48: DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, -1.0,1.0);
50: {
51: DM cdaf,cdac;
52: Vec coordsc,coordsf,coordsf2;
53: VecScatter inject;
54: Mat interp;
55: PetscReal norm;
57: DMGetCoordinateDM(dac,&cdac);
58: DMGetCoordinateDM(daf,&cdaf);
60: DMGetCoordinates(dac,&coordsc);
61: DMGetCoordinates(daf,&coordsf);
63: DMCreateInjection(cdac,cdaf,&inject);
65: VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);
66: VecScatterEnd(inject ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);
67: VecScatterDestroy(&inject);
69: DMCreateInterpolation(cdac,cdaf,&interp,NULL);
70: VecDuplicate(coordsf,&coordsf2);
71: MatInterpolate(interp,coordsc,coordsf2);
72: VecAXPY(coordsf2,-1.0,coordsf);
73: VecNorm(coordsf2,NORM_MAX,&norm);
74: /* The fine coordinates are only reproduced in certain cases */
75: if (!bx && !by && !bz && norm > 1.e-10) {PetscPrintf(PETSC_COMM_WORLD,"Norm %g\n",(double)norm);}
76: VecDestroy(&coordsf2);
77: MatDestroy(&interp);
78: }
80: if (0) {
81: DMCreateGlobalVector(dac,&ac);
82: VecZeroEntries(ac);
84: DMCreateGlobalVector(daf,&af);
85: VecZeroEntries(af);
87: PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_7.vtk", &vv);
88: PetscViewerSetFormat(vv, PETSC_VIEWER_ASCII_VTK);
89: DMView(dac, vv);
90: VecView(ac, vv);
91: PetscViewerDestroy(&vv);
93: PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_7.vtk", &vv);
94: PetscViewerSetFormat(vv, PETSC_VIEWER_ASCII_VTK);
95: DMView(daf, vv);
96: VecView(af, vv);
97: PetscViewerDestroy(&vv);
98: VecDestroy(&ac);
99: VecDestroy(&af);
100: }
101: DMDestroy(&dac);
102: DMDestroy(&daf);
103: return(0);
104: }
108: int main(int argc,char **argv)
109: {
111: PetscInt mx,my,mz;
113: PetscInitialize(&argc,&argv,0,help);
114: mx = 2;
115: my = 2;
116: mz = 2;
117: PetscOptionsGetInt(NULL,"-mx", &mx, 0);
118: PetscOptionsGetInt(NULL,"-my", &my, 0);
119: PetscOptionsGetInt(NULL,"-mz", &mz, 0);
121: test1_DAInjection3d(mx,my,mz);
123: PetscFinalize();
124: return 0;
125: }