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