Actual source code: ex4.c
1: static char help[] = "Demonstrates various vector routines for DMDA.\n\n";
3: /*
4: Include "petscpf.h" so that we can use pf functions and "petscdmda.h" so
5: we can use the PETSc distributed arrays
6: */
8: #include <petscpf.h>
9: #include <petscdm.h>
10: #include <petscdmda.h>
12: PetscErrorCode myfunction(void *ctx, PetscInt n, const PetscScalar *xy, PetscScalar *u)
13: {
14: PetscInt i;
16: PetscFunctionBeginUser;
17: for (i = 0; i < n; i++) {
18: u[2 * i] = xy[2 * i];
19: u[2 * i + 1] = xy[2 * i + 1];
20: }
21: PetscFunctionReturn(PETSC_SUCCESS);
22: }
24: int main(int argc, char **argv)
25: {
26: Vec u, xy;
27: DM da;
28: PetscInt m = 10, n = 10, dof = 2;
29: PF pf;
31: PetscFunctionBeginUser;
32: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
33: PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, m, n, PETSC_DECIDE, PETSC_DECIDE, dof, 1, 0, 0, &da));
34: PetscCall(DMSetFromOptions(da));
35: PetscCall(DMSetUp(da));
36: PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
37: PetscCall(DMCreateGlobalVector(da, &u));
38: PetscCall(DMGetCoordinates(da, &xy));
40: PetscCall(DMDACreatePF(da, &pf));
41: PetscCall(PFSet(pf, myfunction, 0, 0, 0, 0));
42: PetscCall(PFSetFromOptions(pf));
44: PetscCall(PFApplyVec(pf, xy, u));
46: PetscCall(VecView(u, PETSC_VIEWER_DRAW_WORLD));
48: /*
49: Free work space. All PETSc objects should be destroyed when they
50: are no longer needed.
51: */
52: PetscCall(VecDestroy(&u));
53: PetscCall(PFDestroy(&pf));
54: PetscCall(DMDestroy(&da));
55: PetscCall(PetscFinalize());
56: return 0;
57: }
59: /*TEST
61: test:
63: TEST*/