1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "VecView() with a DMDA1d vector and draw viewer.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown #include <petscdm.h> 6*c4762a1bSJed Brown #include <petscdmda.h> 7*c4762a1bSJed Brown #include <petscao.h> 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown PetscErrorCode apply(void *ctx,PetscInt n,const PetscScalar *x,PetscScalar *y) 10*c4762a1bSJed Brown { 11*c4762a1bSJed Brown PetscInt i; 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown for (i=0; i<n; i++) {y[3*i] = x[i]; y[3*i+1] = x[i]*x[i]; y[3*i+2] = x[i]*x[i]*x[i];} 14*c4762a1bSJed Brown return 0; 15*c4762a1bSJed Brown } 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown int main(int argc,char **argv) 18*c4762a1bSJed Brown { 19*c4762a1bSJed Brown PetscErrorCode ierr; 20*c4762a1bSJed Brown DM da; 21*c4762a1bSJed Brown Vec global; 22*c4762a1bSJed Brown PF pf; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 25*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,10,3,1,NULL,&da);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = PFCreate(PETSC_COMM_WORLD,1,3,&pf);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PFSet(pf,apply,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = PFApplyVec(pf,NULL,global);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PFDestroy(&pf);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = VecView(global,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = VecDestroy(&global);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = PetscFinalize(); 37*c4762a1bSJed Brown return ierr; 38*c4762a1bSJed Brown } 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown /*TEST 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown test: 44*c4762a1bSJed Brown nsize: 2 45*c4762a1bSJed Brown requires: x 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown TEST*/ 48