Actual source code: ex20.c


  2: static char help[] = "Bilinear elements on the unit square for Laplacian.  To test the parallel\n\
  3: matrix assembly,the matrix is intentionally laid out across processors\n\
  4: differently from the way it is assembled.  Input arguments are:\n\
  5:   -m <size> : problem size\n\n";

  7: #include <petscksp.h>

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

 18: int main(int argc,char **args)
 19: {
 20:   Mat            C;
 21:   PetscMPIInt    rank,size;
 22:   PetscInt       i,m = 5,N,start,end,M;
 23:   PetscInt       idx[4];
 24:   PetscScalar    Ke[16];
 25:   PetscReal      h;
 26:   Vec            u,b;
 27:   KSP            ksp;
 28:   MatNullSpace   nullsp;

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

 38:   /* Create stiffness matrix */
 39:   MatCreate(PETSC_COMM_WORLD,&C);
 40:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 41:   MatSetFromOptions(C);
 42:   MatSetUp(C);
 43:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 44:   end   = start + M/size + ((M%size) > rank);

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

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

 66:   VecSet(b,1.0);
 67:   VecSetValue(b,0,1.2,ADD_VALUES);
 68:   VecSet(u,0.0);

 70:   /* Solve linear system */
 71:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 72:   KSPSetOperators(ksp,C,C);
 73:   KSPSetFromOptions(ksp);
 74:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);

 76:   MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nullsp);
 77:   /*
 78:      The KSP solver will remove this nullspace from the solution at each iteration
 79:   */
 80:   MatSetNullSpace(C,nullsp);
 81:   /*
 82:      The KSP solver will remove from the right hand side any portion in this nullspace, thus making the linear system consistent.
 83:   */
 84:   MatSetTransposeNullSpace(C,nullsp);
 85:   MatNullSpaceDestroy(&nullsp);

 87:   KSPSolve(ksp,b,u);

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

 98: /*TEST

100:     test:
101:       args: -ksp_monitor_short

103: TEST*/