xref: /petsc/src/dm/tutorials/ex4.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Demonstrates various vector routines for DMDA.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*T
5c4762a1bSJed Brown    Concepts: mathematical functions
6c4762a1bSJed Brown    Processors: n
7c4762a1bSJed Brown T*/
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /*
10c4762a1bSJed Brown   Include "petscpf.h" so that we can use pf functions and "petscdmda.h" so
11c4762a1bSJed Brown  we can use the PETSc distributed arrays
12c4762a1bSJed Brown */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscpf.h>
15c4762a1bSJed Brown #include <petscdm.h>
16c4762a1bSJed Brown #include <petscdmda.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown PetscErrorCode myfunction(void *ctx,PetscInt n,const PetscScalar *xy,PetscScalar *u)
19c4762a1bSJed Brown {
20c4762a1bSJed Brown   PetscInt i;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   PetscFunctionBeginUser;
23c4762a1bSJed Brown   for (i=0; i<n; i++) {
24c4762a1bSJed Brown     u[2*i]   = xy[2*i];
25c4762a1bSJed Brown     u[2*i+1] = xy[2*i+1];
26c4762a1bSJed Brown   }
27c4762a1bSJed Brown   PetscFunctionReturn(0);
28c4762a1bSJed Brown }
29c4762a1bSJed Brown 
30c4762a1bSJed Brown int main(int argc,char **argv)
31c4762a1bSJed Brown {
32c4762a1bSJed Brown   Vec            u,xy;
33c4762a1bSJed Brown   DM             da;
34c4762a1bSJed Brown   PetscInt       m = 10, n = 10, dof = 2;
35c4762a1bSJed Brown   PF             pf;
36c4762a1bSJed Brown 
37*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,PETSC_DECIDE,PETSC_DECIDE,dof,1,0,0,&da));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&u));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(da,&xy));
44c4762a1bSJed Brown 
455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreatePF(da,&pf));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PFSet(pf,myfunction,0,0,0,0));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PFSetFromOptions(pf));
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PFApplyVec(pf,xy,u));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(u,PETSC_VIEWER_DRAW_WORLD));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /*
54c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
55c4762a1bSJed Brown      are no longer needed.
56c4762a1bSJed Brown   */
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PFDestroy(&pf));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
60*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
61*b122ec5aSJacob Faibussowitsch   return 0;
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown /*TEST
65c4762a1bSJed Brown 
66c4762a1bSJed Brown    test:
67c4762a1bSJed Brown 
68c4762a1bSJed Brown TEST*/
69