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