Actual source code: ex14.c

petsc-3.3-p7 2013-05-11
  2: static char help[] = "Tests basic assembly of MATIJ using the edges of a 2d parallel DA.\n";

  4: #include <petscmat.h> 
  5: #include <petscdmda.h>
  8: int main(int argc,char **args)
  9: {
 11:   DM       da;
 12:   Mat      A;
 13:   PetscInt i,j,k,i0,j0,m,n,gi0,gj0,gm,gn,M = 4,N = 4, e0[2],e1[2];
 14: #if 0
 15:   PetscBool preallocate,flag;
 16: #endif

 18:   PetscInitialize(&argc,&args,(char *)0,help);

 20:   /* Create the DA and extract its size info. */
 21:   DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_STENCIL_STAR,-M,-N,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,&da);
 22:   DMDAGetInfo(da, PETSC_NULL, &M,&N,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 23:   DMDAGetCorners(da,&i0,&j0,PETSC_NULL,&m,&n,PETSC_NULL);
 24:   DMDAGetGhostCorners(da,&gi0,&gj0,PETSC_NULL,&gm,&gn,PETSC_NULL);

 26:   /* Create MATIJ with m*n local rows (and columns). */
 27:   MatCreate(PETSC_COMM_WORLD,&A);
 28:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 29:   MatSetType(A,MATIJ);
 30: #if 0
 31:   PetscOptionsGetBool(PETSC_NULL, "--preallocate", &preallocate, &flag);
 32:   if(preallocate) {
 33:   }
 34: #endif
 35: 
 36:   /* 
 37:    Add local and ghosted edges to A: grid points are indexed by i first, 
 38:    so that points with the same i-index differ by a multiple of M.
 39:    */
 40:   for(j = j0; j < j0+n; ++j) {
 41:     for(i = i0; i < i0+m; ++i) {
 42:       k = 0;
 43:       if(i+1 < gi0+gm) {/* there is a point to the right, so draw an edge to it.*/
 44:         e0[k] = i*M+j; e1[k] = (i+i)*M+j;
 45:         ++k;
 46:       }
 47:       if(j+1 < gj0+gn) {/* there is a point above, so draw an edge to it.*/
 48:         e0[k] = i*M+j; e1[k] = (i)*M+j+1;
 49:         ++k;
 50:       }
 51:       MatIJSetEdges(A,k,e0,e1);
 52: 
 53:     }
 54:   }
 55:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);

 58:   /* View A. */
 59:   MatView(A, PETSC_VIEWER_STDOUT_WORLD);
 60: 
 61:   MatDestroy(&A);
 62:   DMDestroy(&da);
 63:   PetscFinalize();
 64:   return 0;
 65: }