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: PetscInt N; /* dimension of system (global) */
36: PetscInt M; /* number of elements (global) */
37: PetscMPIInt rank; /* processor rank */
38: PetscMPIInt size; /* size of communicator */
39: PetscScalar Ke[16]; /* element matrix */
40: PetscScalar r[4]; /* element vector */
41: PetscReal h; /* mesh width */
42: PetscReal norm; /* norm of solution error */
43: PetscReal x,y;
44: PetscScalar val;
46: PetscInt idx[4],count,*rows,i,m = 5,start,end,its;
48: PetscInitialize(&argc,&args,(char *)0,help);
49: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
50: N = (m+1)*(m+1);
51: M = m*m;
52: h = 1.0/m;
53: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
54: MPI_Comm_size(PETSC_COMM_WORLD,&size);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Compute the matrix and right-hand-side vector that define
58: the linear system, Au = b.
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: /*
62: Create stiffness matrix
63: */
64: MatCreate(PETSC_COMM_WORLD,&A);
65: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
66: MatSetFromOptions(A);
67: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
68: end = start + M/size + ((M%size) > rank);
70: /*
71: Assemble matrix
72: */
73: FormElementStiffness(h*h,Ke);
74: for (i=start; i<end; i++) {
75: /* location of lower left corner of element */
76: x = h*(i % m); y = h*(i/m);
77: /* node numbers for the four corners of element */
78: idx[0] = (m+1)*(i/m) + (i % m);
79: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
80: MatSetValues(A,4,idx,4,idx,Ke,ADD_VALUES);
81: }
82: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
85: /*
86: Create right-hand-side and solution vectors
87: */
88: VecCreate(PETSC_COMM_WORLD,&u);
89: VecSetSizes(u,PETSC_DECIDE,N);
90: VecSetFromOptions(u);
91: PetscObjectSetName((PetscObject)u,"Approx. Solution");
92: VecDuplicate(u,&b);
93: PetscObjectSetName((PetscObject)b,"Right hand side");
94: VecDuplicate(b,&ustar);
95: VecSet(u,0.0);
96: VecSet(b,0.0);
98: /*
99: Assemble right-hand-side vector
100: */
101: for (i=start; i<end; i++) {
102: /* location of lower left corner of element */
103: x = h*(i % m); y = h*(i/m);
104: /* node numbers for the four corners of element */
105: idx[0] = (m+1)*(i/m) + (i % m);
106: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
107: FormElementRhs(x,y,h*h,r);
108: VecSetValues(b,4,idx,r,ADD_VALUES);
109: }
110: VecAssemblyBegin(b);
111: VecAssemblyEnd(b);
113: /*
114: Modify matrix and right-hand-side for Dirichlet boundary conditions
115: */
116: PetscMalloc(4*m*sizeof(PetscInt),&rows);
117: for (i=0; i<m+1; i++) {
118: rows[i] = i; /* bottom */
119: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
120: }
121: count = m+1; /* left side */
122: for (i=m+1; i<m*(m+1); i+= m+1) {
123: rows[count++] = i;
124: }
125: count = 2*m; /* left side */
126: for (i=2*m+1; i<m*(m+1); i+= m+1) {
127: rows[count++] = i;
128: }
129: for (i=0; i<4*m; i++) {
130: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
131: val = y;
132: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
133: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
134: }
135: MatZeroRows(A,4*m,rows,1.0);
136: PetscFree(rows);
138: VecAssemblyBegin(u);
139: VecAssemblyEnd(u);
140: VecAssemblyBegin(b);
141: VecAssemblyEnd(b);
143: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: Create the linear solver and set various options
145: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147: KSPCreate(PETSC_COMM_WORLD,&ksp);
148: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
149: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
150: KSPSetFromOptions(ksp);
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Solve the linear system
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: KSPSolve(ksp,b,u);
158: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: Check solution and clean up
160: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162: /* Check error */
163: VecGetOwnershipRange(ustar,&start,&end);
164: for (i=start; i<end; i++) {
165: x = h*(i % (m+1)); y = h*(i/(m+1));
166: val = y;
167: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
168: }
169: VecAssemblyBegin(ustar);
170: VecAssemblyEnd(ustar);
171: VecAXPY(u,-1.0,ustar);
172: VecNorm(u,NORM_2,&norm);
173: KSPGetIterationNumber(ksp,&its);
174: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %D\n",norm*h,its);
176: /*
177: Free work space. All PETSc objects should be destroyed when they
178: are no longer needed.
179: */
180: KSPDestroy(ksp); VecDestroy(u);
181: VecDestroy(ustar); VecDestroy(b);
182: MatDestroy(A);
184: /*
185: Always call PetscFinalize() before exiting a program. This routine
186: - finalizes the PETSc libraries as well as MPI
187: - provides summary and diagnostic information if certain runtime
188: options are chosen (e.g., -log_summary).
189: */
190: PetscFinalize();
191: return 0;
192: }
194: /* --------------------------------------------------------------------- */
197: /* element stiffness for Laplacian */
198: PetscErrorCode FormElementStiffness(PetscReal H,PetscScalar *Ke)
199: {
201: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
202: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
203: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
204: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
205: return(0);
206: }
207: /* --------------------------------------------------------------------- */
210: PetscErrorCode FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
211: {
213: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
214: return(0);
215: }