Actual source code: ex4.c


  2: static char help[] = "Demonstrates various vector routines for DMDA.\n\n";

  4: /*
  5:   Include "petscpf.h" so that we can use pf functions and "petscdmda.h" so
  6:  we can use the PETSc distributed arrays
  7: */

  9: #include <petscpf.h>
 10: #include <petscdm.h>
 11: #include <petscdmda.h>

 13: PetscErrorCode myfunction(void *ctx,PetscInt n,const PetscScalar *xy,PetscScalar *u)
 14: {
 15:   PetscInt i;

 18:   for (i=0; i<n; i++) {
 19:     u[2*i]   = xy[2*i];
 20:     u[2*i+1] = xy[2*i+1];
 21:   }
 22:   return 0;
 23: }

 25: int main(int argc,char **argv)
 26: {
 27:   Vec            u,xy;
 28:   DM             da;
 29:   PetscInt       m = 10, n = 10, dof = 2;
 30:   PF             pf;

 32:   PetscInitialize(&argc,&argv,(char*)0,help);
 33:   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:   DMSetFromOptions(da);
 35:   DMSetUp(da);
 36:   DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);
 37:   DMCreateGlobalVector(da,&u);
 38:   DMGetCoordinates(da,&xy);

 40:   DMDACreatePF(da,&pf);
 41:   PFSet(pf,myfunction,0,0,0,0);
 42:   PFSetFromOptions(pf);

 44:   PFApplyVec(pf,xy,u);

 46:   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:   VecDestroy(&u);
 53:   PFDestroy(&pf);
 54:   DMDestroy(&da);
 55:   PetscFinalize();
 56:   return 0;
 57: }

 59: /*TEST

 61:    test:

 63: TEST*/