Actual source code: ex3.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  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: /*T
  8:    Concepts: KSP^basic parallel example
  9:    Concepts: Matrices^inserting elements by blocks
 10:    Processors: n
 11: T*/

 13: /*
 14:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 15:   automatically includes:
 16:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 17:      petscmat.h - matrices
 18:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 19:      petscviewer.h - viewers               petscpc.h  - preconditioners
 20: */
 21: #include <petscksp.h>

 23: /* Declare user-defined routines */
 24: extern PetscErrorCode FormElementStiffness(PetscReal,PetscScalar*);
 25: extern PetscErrorCode FormElementRhs(PetscReal,PetscReal,PetscReal,PetscScalar*);

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

 47:   PetscInitialize(&argc,&args,(char*)0,help);
 48:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 49:   N    = (m+1)*(m+1);
 50:   M    = m*m;
 51:   h    = 1.0/m;
 52:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 53:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 55:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 56:          Compute the matrix and right-hand-side vector that define
 57:          the linear system, Au = b.
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 60:   /*
 61:      Create stiffness matrix
 62:   */
 63:   MatCreate(PETSC_COMM_WORLD,&A);
 64:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 65:   MatSetFromOptions(A);
 66:   MatSeqAIJSetPreallocation(A,9,NULL);
 67:   MatMPIAIJSetPreallocation(A,9,NULL,8,NULL);
 68:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 69:   end   = start + M/size + ((M%size) > rank);

 71:   /*
 72:      Assemble matrix
 73:   */
 74:   FormElementStiffness(h*h,Ke);
 75:   for (i=start; i<end; i++) {
 76:     /* location of lower left corner of element */
 77:     x = h*(i % m); y = h*(i/m);
 78:     /* node numbers for the four corners of element */
 79:     idx[0] = (m+1)*(i/m) + (i % m);
 80:     idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 81:     MatSetValues(A,4,idx,4,idx,Ke,ADD_VALUES);
 82:   }
 83:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 86:   /*
 87:      Create right-hand-side and solution vectors
 88:   */
 89:   VecCreate(PETSC_COMM_WORLD,&u);
 90:   VecSetSizes(u,PETSC_DECIDE,N);
 91:   VecSetFromOptions(u);
 92:   PetscObjectSetName((PetscObject)u,"Approx. Solution");
 93:   VecDuplicate(u,&b);
 94:   PetscObjectSetName((PetscObject)b,"Right hand side");
 95:   VecDuplicate(b,&ustar);
 96:   VecSet(u,0.0);
 97:   VecSet(b,0.0);

 99:   /*
100:      Assemble right-hand-side vector
101:   */
102:   for (i=start; i<end; i++) {
103:     /* location of lower left corner of element */
104:     x = h*(i % m); y = h*(i/m);
105:     /* node numbers for the four corners of element */
106:     idx[0] = (m+1)*(i/m) + (i % m);
107:     idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
108:     FormElementRhs(x,y,h*h,r);
109:     VecSetValues(b,4,idx,r,ADD_VALUES);
110:   }
111:   VecAssemblyBegin(b);
112:   VecAssemblyEnd(b);

114:   /*
115:      Modify matrix and right-hand-side for Dirichlet boundary conditions
116:   */
117:   PetscMalloc1(4*m,&rows);
118:   for (i=0; i<m+1; i++) {
119:     rows[i] = i; /* bottom */
120:     rows[3*m - 1 +i] = m*(m+1) + i; /* top */
121:   }
122:   count = m+1; /* left side */
123:   for (i=m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
124:   count = 2*m; /* left side */
125:   for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
126:   for (i=0; i<4*m; i++) {
127:     x    = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
128:     val  = y;
129:     VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
130:     VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
131:   }
132:   MatZeroRows(A,4*m,rows,1.0,0,0);
133:   PetscFree(rows);

135:   VecAssemblyBegin(u);
136:   VecAssemblyEnd(u);
137:   VecAssemblyBegin(b);
138:   VecAssemblyEnd(b);

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:                 Create the linear solver and set various options
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

144:   KSPCreate(PETSC_COMM_WORLD,&ksp);
145:   KSPSetOperators(ksp,A,A);
146:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
147:   KSPSetFromOptions(ksp);

149:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150:                       Solve the linear system
151:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

153:   KSPSolve(ksp,b,u);

155:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156:                       Check solution and clean up
157:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

159:   /* Check error */
160:   VecGetOwnershipRange(ustar,&start,&end);
161:   for (i=start; i<end; i++) {
162:     x    = h*(i % (m+1)); y = h*(i/(m+1));
163:     val  = y;
164:     VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
165:   }
166:   VecAssemblyBegin(ustar);
167:   VecAssemblyEnd(ustar);
168:   VecAXPY(u,-1.0,ustar);
169:   VecNorm(u,NORM_2,&norm);
170:   KSPGetIterationNumber(ksp,&its);
171:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)(norm*h),its);

173:   /*
174:      Free work space.  All PETSc objects should be destroyed when they
175:      are no longer needed.
176:   */
177:   KSPDestroy(&ksp); VecDestroy(&u);
178:   VecDestroy(&ustar); VecDestroy(&b);
179:   MatDestroy(&A);

181:   /*
182:      Always call PetscFinalize() before exiting a program.  This routine
183:        - finalizes the PETSc libraries as well as MPI
184:        - provides summary and diagnostic information if certain runtime
185:          options are chosen (e.g., -log_summary).
186:   */
187:   PetscFinalize();
188:   return 0;
189: }

191: /* --------------------------------------------------------------------- */
194: /* element stiffness for Laplacian */
195: PetscErrorCode FormElementStiffness(PetscReal H,PetscScalar *Ke)
196: {
198:   Ke[0]  = H/6.0;    Ke[1]  = -.125*H; Ke[2]  = H/12.0;   Ke[3]  = -.125*H;
199:   Ke[4]  = -.125*H;  Ke[5]  = H/6.0;   Ke[6]  = -.125*H;  Ke[7]  = H/12.0;
200:   Ke[8]  = H/12.0;   Ke[9]  = -.125*H; Ke[10] = H/6.0;    Ke[11] = -.125*H;
201:   Ke[12] = -.125*H;  Ke[13] = H/12.0;  Ke[14] = -.125*H;  Ke[15] = H/6.0;
202:   return(0);
203: }
204: /* --------------------------------------------------------------------- */
207: PetscErrorCode FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
208: {
210:   r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
211:   return(0);
212: }