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