Actual source code: ex20.c

  1: /*$Id: ex20.c,v 1.19 2001/08/07 21:30:50 bsmith Exp $*/

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

 9:  #include petscsles.h

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

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

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

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

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

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

 68:   VecSet(&one,u);
 69:   MatMult(C,u,b);
 70:   VecSet(&zero,u);

 72:   /* Solve linear system */
 73:   SLESCreate(PETSC_COMM_WORLD,&sles);
 74:   SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
 75:   SLESSetFromOptions(sles);
 76:   SLESGetKSP(sles,&ksp);
 77:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);

 79:   PetscOptionsHasName(PETSC_NULL,"-fixnullspace",&flg);
 80:   if (flg) {
 81:     SLESGetPC(sles,&pc);
 82:     MatNullSpaceCreate(PETSC_COMM_WORLD,1,0,PETSC_NULL,&nullsp);
 83:     PCNullSpaceAttach(pc,nullsp);
 84:     MatNullSpaceDestroy(nullsp);
 85:   }

 87:   SLESSolve(sles,b,u,&its);


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