Actual source code: ex3.c
2: static char help[] = "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: /*T
9: Concepts: KSP^basic parallel example
10: Concepts: Matrices^inserting elements by blocks
11: Processors: n
12: T*/
14: /*
15: Include "petscksp.h" so that we can use KSP solvers. Note that this file
16: automatically includes:
17: petsc.h - base PETSc routines petscvec.h - vectors
18: petscsys.h - system routines petscmat.h - matrices
19: petscis.h - index sets petscksp.h - Krylov subspace methods
20: petscviewer.h - viewers petscpc.h - preconditioners
21: */
22: #include petscksp.h
24: /* Declare user-defined routines */
30: int main(int argc,char **args)
31: {
32: Vec u,b,ustar; /* approx solution, RHS, exact solution */
33: Mat A; /* linear system matrix */
34: KSP ksp; /* Krylov subspace method context */
35: IS is; /* index set - used for boundary conditions */
36: PetscInt N; /* dimension of system (global) */
37: PetscInt M; /* number of elements (global) */
38: PetscMPIInt rank; /* processor rank */
39: PetscMPIInt size; /* size of communicator */
40: PetscScalar Ke[16]; /* element matrix */
41: PetscScalar r[4]; /* element vector */
42: PetscReal h; /* mesh width */
43: PetscReal norm; /* norm of solution error */
44: PetscReal x,y;
45: PetscScalar val,zero = 0.0,one = 1.0,none = -1.0;
47: PetscInt idx[4],count,*rows,i,m = 5,start,end,its;
49: PetscInitialize(&argc,&args,(char *)0,help);
50: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
51: N = (m+1)*(m+1);
52: M = m*m;
53: h = 1.0/m;
54: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
55: MPI_Comm_size(PETSC_COMM_WORLD,&size);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Compute the matrix and right-hand-side vector that define
59: the linear system, Au = b.
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: /*
63: Create stiffness matrix
64: */
65: MatCreate(PETSC_COMM_WORLD,&A);
66: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
67: MatSetFromOptions(A);
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,zero);
97: VecSet(b,zero);
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: PetscMalloc(4*m*sizeof(PetscInt),&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) {
124: rows[count++] = i;
125: }
126: count = 2*m; /* left side */
127: for (i=2*m+1; i<m*(m+1); i+= m+1) {
128: rows[count++] = i;
129: }
130: for (i=0; i<4*m; i++) {
131: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
132: val = y;
133: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
134: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
135: }
136: PetscFree(rows);
137: VecAssemblyBegin(u);
138: VecAssemblyEnd(u);
139: VecAssemblyBegin(b);
140: VecAssemblyEnd(b);
142: MatZeroRows(A,4*m,rows,one);
144: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145: Create the linear solver and set various options
146: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148: KSPCreate(PETSC_COMM_WORLD,&ksp);
149: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
150: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
151: KSPSetFromOptions(ksp);
153: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: Solve the linear system
155: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157: KSPSolve(ksp,b,u);
159: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: Check solution and clean up
161: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: /* Check error */
164: VecGetOwnershipRange(ustar,&start,&end);
165: for (i=start; i<end; i++) {
166: x = h*(i % (m+1)); y = h*(i/(m+1));
167: val = y;
168: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
169: }
170: VecAssemblyBegin(ustar);
171: VecAssemblyEnd(ustar);
172: VecAXPY(u,none,ustar);
173: VecNorm(u,NORM_2,&norm);
174: KSPGetIterationNumber(ksp,&its);
175: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %D\n",norm*h,its);
177: /*
178: Free work space. All PETSc objects should be destroyed when they
179: are no longer needed.
180: */
181: KSPDestroy(ksp); VecDestroy(u);
182: VecDestroy(ustar); VecDestroy(b);
183: MatDestroy(A);
185: /*
186: Always call PetscFinalize() before exiting a program. This routine
187: - finalizes the PETSc libraries as well as MPI
188: - provides summary and diagnostic information if certain runtime
189: options are chosen (e.g., -log_summary).
190: */
191: PetscFinalize();
192: return 0;
193: }
195: /* --------------------------------------------------------------------- */
198: /* element stiffness for Laplacian */
199: PetscErrorCode FormElementStiffness(PetscReal H,PetscScalar *Ke)
200: {
202: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
203: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
204: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
205: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
206: return(0);
207: }
208: /* --------------------------------------------------------------------- */
211: PetscErrorCode FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
212: {
214: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
215: return(0);
216: }