Actual source code: ex4.c

petsc-3.3-p7 2013-05-11
  2: /* Program usage:  mpiexec ex1 [-help] [all PETSc options] */

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

  6: /*T
  7:    Concepts: mathematical functions
  8:    Processors: n
  9: T*/

 11: /* 
 12:   Include "petscpf.h" so that we can use pf functions and "petscdmda.h" so
 13:  we can use the PETSc distributed arrays
 14: */

 16: #include <petscpf.h>
 17: #include <petscdmda.h>

 21: PetscErrorCode myfunction(void *ctx,PetscInt n,const PetscScalar *xy,PetscScalar *u)
 22: {
 23:   PetscInt i;

 26:   for (i=0; i<n; i++) {
 27:     u[2*i] = xy[2*i];
 28:     u[2*i+1] = xy[2*i+1];
 29:   }
 30:   return(0);
 31: }

 35: int main(int argc,char **argv)
 36: {
 37:   Vec            u,xy;
 38:   DM             da;
 40:   PetscInt       m = 10, n = 10, dof = 2;
 41:   PF             pf;

 43:   PetscInitialize(&argc,&argv,(char*)0,help);
 44: 
 45:   DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,PETSC_DECIDE,PETSC_DECIDE,dof,1,0,0,&da);
 46:   DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);
 47:   DMCreateGlobalVector(da,&u);
 48:   DMDAGetCoordinates(da,&xy);

 50:   DMDACreatePF(da,&pf);
 51:   PFSet(pf,myfunction,0,0,0,0);
 52:   PFSetFromOptions(pf);

 54:   PFApplyVec(pf,xy,u);

 56:   VecView(u,PETSC_VIEWER_DRAW_WORLD);

 58:   /* 
 59:      Free work space.  All PETSc objects should be destroyed when they
 60:      are no longer needed.
 61:   */
 62:   PFDestroy(&pf);
 63:   DMDestroy(&da);
 64:   PetscFinalize();
 65:   return 0;
 66: }
 67: