Actual source code: ex3.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: /*
  8:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
  9:   automatically includes:
 10:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 11:      petscmat.h - matrices
 12:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 13:      petscviewer.h - viewers               petscpc.h  - preconditioners
 14: */
 15: #include <petscksp.h>

 17: /* Declare user-defined routines */
 18: extern PetscErrorCode FormElementStiffness(PetscReal,PetscScalar*);
 19: extern PetscErrorCode FormElementRhs(PetscScalar,PetscScalar,PetscReal,PetscScalar*);

 21: int main(int argc,char **args)
 22: {
 23:   Vec            u,b,ustar; /* approx solution, RHS, exact solution */
 24:   Mat            A;           /* linear system matrix */
 25:   KSP            ksp;         /* Krylov subspace method context */
 26:   PetscInt       N;           /* dimension of system (global) */
 27:   PetscInt       M;           /* number of elements (global) */
 28:   PetscMPIInt    rank;        /* processor rank */
 29:   PetscMPIInt    size;        /* size of communicator */
 30:   PetscScalar    Ke[16];      /* element matrix */
 31:   PetscScalar    r[4];        /* element vector */
 32:   PetscReal      h;           /* mesh width */
 33:   PetscReal      norm;        /* norm of solution error */
 34:   PetscScalar    x,y;
 35:   PetscInt       idx[4],count,*rows,i,m = 5,start,end,its;

 37:   PetscInitialize(&argc,&args,(char*)0,help);
 38:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 39:   N    = (m+1)*(m+1);
 40:   M    = m*m;
 41:   h    = 1.0/m;
 42:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 43:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 45:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 46:          Compute the matrix and right-hand-side vector that define
 47:          the linear system, Au = b.
 48:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 50:   /*
 51:      Create stiffness matrix
 52:   */
 53:   MatCreate(PETSC_COMM_WORLD,&A);
 54:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 55:   MatSetFromOptions(A);
 56:   MatSeqAIJSetPreallocation(A,9,NULL);
 57:   MatMPIAIJSetPreallocation(A,9,NULL,8,NULL);
 58:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 59:   end   = start + M/size + ((M%size) > rank);

 61:   /*
 62:      Assemble matrix
 63:   */
 64:   FormElementStiffness(h*h,Ke);
 65:   for (i=start; i<end; i++) {
 66:     /* node numbers for the four corners of element */
 67:     idx[0] = (m+1)*(i/m) + (i % m);
 68:     idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 69:     MatSetValues(A,4,idx,4,idx,Ke,ADD_VALUES);
 70:   }
 71:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 72:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 74:   /*
 75:      Create right-hand-side and solution vectors
 76:   */
 77:   VecCreate(PETSC_COMM_WORLD,&u);
 78:   VecSetSizes(u,PETSC_DECIDE,N);
 79:   VecSetFromOptions(u);
 80:   PetscObjectSetName((PetscObject)u,"Approx. Solution");
 81:   VecDuplicate(u,&b);
 82:   PetscObjectSetName((PetscObject)b,"Right hand side");
 83:   VecDuplicate(b,&ustar);
 84:   VecSet(u,0.0);
 85:   VecSet(b,0.0);

 87:   /*
 88:      Assemble right-hand-side vector
 89:   */
 90:   for (i=start; i<end; i++) {
 91:     /* location of lower left corner of element */
 92:     x = h*(i % m); y = h*(i/m);
 93:     /* node numbers for the four corners of element */
 94:     idx[0] = (m+1)*(i/m) + (i % m);
 95:     idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 96:     FormElementRhs(x,y,h*h,r);
 97:     VecSetValues(b,4,idx,r,ADD_VALUES);
 98:   }
 99:   VecAssemblyBegin(b);
100:   VecAssemblyEnd(b);

102:   /*
103:      Modify matrix and right-hand-side for Dirichlet boundary conditions
104:   */
105:   PetscMalloc1(4*m,&rows);
106:   for (i=0; i<m+1; i++) {
107:     rows[i] = i; /* bottom */
108:     rows[3*m - 1 +i] = m*(m+1) + i; /* top */
109:   }
110:   count = m+1; /* left side */
111:   for (i=m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
112:   count = 2*m; /* left side */
113:   for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
114:   for (i=0; i<4*m; i++) {
115:     y = h*(rows[i]/(m+1));
116:     VecSetValues(u,1,&rows[i],&y,INSERT_VALUES);
117:     VecSetValues(b,1,&rows[i],&y,INSERT_VALUES);
118:   }
119:   MatZeroRows(A,4*m,rows,1.0,0,0);
120:   PetscFree(rows);

122:   VecAssemblyBegin(u);
123:   VecAssemblyEnd(u);
124:   VecAssemblyBegin(b);
125:   VecAssemblyEnd(b);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:                 Create the linear solver and set various options
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

131:   KSPCreate(PETSC_COMM_WORLD,&ksp);
132:   KSPSetOperators(ksp,A,A);
133:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
134:   KSPSetFromOptions(ksp);

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:                       Solve the linear system
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

140:   KSPSolve(ksp,b,u);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:                       Check solution and clean up
144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

146:   /* Check error */
147:   VecGetOwnershipRange(ustar,&start,&end);
148:   for (i=start; i<end; i++) {
149:     y = h*(i/(m+1));
150:     VecSetValues(ustar,1,&i,&y,INSERT_VALUES);
151:   }
152:   VecAssemblyBegin(ustar);
153:   VecAssemblyEnd(ustar);
154:   VecAXPY(u,-1.0,ustar);
155:   VecNorm(u,NORM_2,&norm);
156:   KSPGetIterationNumber(ksp,&its);
157:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)(norm*h),its);

159:   /*
160:      Free work space.  All PETSc objects should be destroyed when they
161:      are no longer needed.
162:   */
163:   KSPDestroy(&ksp)); PetscCall(VecDestroy(&u);
164:   VecDestroy(&ustar)); PetscCall(VecDestroy(&b);
165:   MatDestroy(&A);

167:   /*
168:      Always call PetscFinalize() before exiting a program.  This routine
169:        - finalizes the PETSc libraries as well as MPI
170:        - provides summary and diagnostic information if certain runtime
171:          options are chosen (e.g., -log_view).
172:   */
173:   PetscFinalize();
174:   return 0;
175: }

177: /* --------------------------------------------------------------------- */
178: /* element stiffness for Laplacian */
179: PetscErrorCode FormElementStiffness(PetscReal H,PetscScalar *Ke)
180: {
182:   Ke[0]  = H/6.0;    Ke[1]  = -.125*H; Ke[2]  = H/12.0;   Ke[3]  = -.125*H;
183:   Ke[4]  = -.125*H;  Ke[5]  = H/6.0;   Ke[6]  = -.125*H;  Ke[7]  = H/12.0;
184:   Ke[8]  = H/12.0;   Ke[9]  = -.125*H; Ke[10] = H/6.0;    Ke[11] = -.125*H;
185:   Ke[12] = -.125*H;  Ke[13] = H/12.0;  Ke[14] = -.125*H;  Ke[15] = H/6.0;
186:   return 0;
187: }
188: /* --------------------------------------------------------------------- */
189: PetscErrorCode FormElementRhs(PetscScalar x,PetscScalar y,PetscReal H,PetscScalar *r)
190: {
192:   r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
193:   return 0;
194: }

196: /*TEST

198:    test:
199:       suffix: 1
200:       nsize: 2
201:       args: -ksp_monitor_short

203: TEST*/