Actual source code: ex20.c

  2: static char help[] = "This example solves a linear system in parallel with KSP.  The matrix\n\
  3: uses simple bilinear elements on the unit square.  To test the parallel\n\
  4: matrix assembly,the matrix is intentionally laid out across processors\n\
  5: differently from the way it is assembled.  Input arguments are:\n\
  6:   -m <size> : problem size\n\n";

 8:  #include petscksp.h

 12: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
 13: {
 14:   Ke[0]  = H/6.0;    Ke[1]  = -.125*H; Ke[2]  = H/12.0;   Ke[3]  = -.125*H;
 15:   Ke[4]  = -.125*H;  Ke[5]  = H/6.0;   Ke[6]  = -.125*H;  Ke[7]  = H/12.0;
 16:   Ke[8]  = H/12.0;   Ke[9]  = -.125*H; Ke[10] = H/6.0;    Ke[11] = -.125*H;
 17:   Ke[12] = -.125*H;  Ke[13] = H/12.0;  Ke[14] = -.125*H;  Ke[15] = H/6.0;
 18:   return 0;
 19: }

 23: int main(int argc,char **args)
 24: {
 25:   Mat          C;
 26:   int          i,m = 5,rank,size,N,start,end,M;
 27:   int          ierr,idx[4];
 28:   PetscTruth   flg;
 29:   PetscScalar  Ke[16];
 30:   PetscReal    h;
 31:   Vec          u,b;
 32:   KSP          ksp;
 33:   MatNullSpace nullsp;

 35:   PetscInitialize(&argc,&args,(char *)0,help);
 36:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 37:   N = (m+1)*(m+1); /* dimension of matrix */
 38:   M = m*m; /* number of elements */
 39:   h = 1.0/m;       /* mesh width */
 40:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 41:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 43:   /* Create stiffness matrix */
 44:   MatCreate(PETSC_COMM_WORLD,&C);
 45:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 46:   MatSetFromOptions(C);
 47:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 48:   end   = start + M/size + ((M%size) > rank);

 50:   /* Assemble matrix */
 51:   FormElementStiffness(h*h,Ke);   /* element stiffness for Laplacian */
 52:   for (i=start; i<end; i++) {
 53:      /* location of lower left corner of element */
 54:      /* node numbers for the four corners of element */
 55:      idx[0] = (m+1)*(i/m) + (i % m);
 56:      idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 57:      MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
 58:   }
 59:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 62:   /* Create right-hand-side and solution vectors */
 63:   VecCreate(PETSC_COMM_WORLD,&u);
 64:   VecSetSizes(u,PETSC_DECIDE,N);
 65:   VecSetFromOptions(u);
 66:   PetscObjectSetName((PetscObject)u,"Approx. Solution");
 67:   VecDuplicate(u,&b);
 68:   PetscObjectSetName((PetscObject)b,"Right hand side");

 70:   VecSet(u,1.0);
 71:   MatMult(C,u,b);
 72:   VecSet(u,0.0);

 74:   /* Solve linear system */
 75:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 76:   KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
 77:   KSPSetFromOptions(ksp);
 78:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);

 80:   flg  = PETSC_FALSE;
 81:   PetscOptionsGetTruth(PETSC_NULL,"-fixnullspace",&flg,PETSC_NULL);
 82:   if (flg) {
 83:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,PETSC_NULL,&nullsp);
 84:     KSPSetNullSpace(ksp,nullsp);
 85:     MatNullSpaceDestroy(nullsp);
 86:   }
 87:   KSPSolve(ksp,b,u);


 90:   /* Free work space */
 91:   KSPDestroy(ksp);
 92:   VecDestroy(u);
 93:   VecDestroy(b);
 94:   MatDestroy(C);
 95:   PetscFinalize();
 96:   return 0;
 97: }