Actual source code: ex3.c
petsc-3.3-p6 2013-02-11
2: static char help[] = "This example 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: #include <petscksp.h>
12: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
13: {
15: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
16: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
17: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
18: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
19: return(0);
20: }
23: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
24: {
26: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
27: return(0);
28: }
32: int main(int argc,char **args)
33: {
34: Mat C;
35: PetscMPIInt rank,size;
36: PetscInt i,m = 5,N,start,end,M,its;
37: PetscScalar val,Ke[16],r[4];
38: PetscReal x,y,h,norm,tol=1.e-14;
40: PetscInt idx[4],count,*rows;
41: Vec u,ustar,b;
42: KSP ksp;
44: PetscInitialize(&argc,&args,(char *)0,help);
45: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
46: N = (m+1)*(m+1); /* dimension of matrix */
47: M = m*m; /* number of elements */
48: h = 1.0/m; /* mesh width */
49: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
50: MPI_Comm_size(PETSC_COMM_WORLD,&size);
52: /* Create stiffness matrix */
53: MatCreate(PETSC_COMM_WORLD,&C);
54: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
55: MatSetFromOptions(C);
56: MatSetUp(C);
57: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
58: end = start + M/size + ((M%size) > rank);
60: /* Assemble matrix */
61: FormElementStiffness(h*h,Ke); /* element stiffness for Laplacian */
62: for (i=start; i<end; i++) {
63: /* location of lower left corner of element */
64: x = h*(i % m); y = h*(i/m);
65: /* node numbers for the four corners of element */
66: idx[0] = (m+1)*(i/m) + (i % m);
67: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
68: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
69: }
70: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
71: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
73: /* Create right-hand-side and solution vectors */
74: VecCreate(PETSC_COMM_WORLD,&u);
75: VecSetSizes(u,PETSC_DECIDE,N);
76: VecSetFromOptions(u);
77: PetscObjectSetName((PetscObject)u,"Approx. Solution");
78: VecDuplicate(u,&b);
79: PetscObjectSetName((PetscObject)b,"Right hand side");
80: VecDuplicate(b,&ustar);
81: VecSet(u,0.0);
82: VecSet(b,0.0);
84: /* Assemble right-hand-side vector */
85: for (i=start; i<end; i++) {
86: /* location of lower left corner of element */
87: x = h*(i % m); y = h*(i/m);
88: /* node numbers for the four corners of element */
89: idx[0] = (m+1)*(i/m) + (i % m);
90: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
91: FormElementRhs(x,y,h*h,r);
92: VecSetValues(b,4,idx,r,ADD_VALUES);
93: }
94: VecAssemblyBegin(b);
95: VecAssemblyEnd(b);
97: /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
98: PetscMalloc(4*m*sizeof(PetscInt),&rows);
99: for (i=0; i<m+1; i++) {
100: rows[i] = i; /* bottom */
101: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
102: }
103: count = m+1; /* left side */
104: for (i=m+1; i<m*(m+1); i+= m+1) {
105: rows[count++] = i;
106: }
107: count = 2*m; /* left side */
108: for (i=2*m+1; i<m*(m+1); i+= m+1) {
109: rows[count++] = i;
110: }
111: for (i=0; i<4*m; i++) {
112: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
113: val = y;
114: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
115: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
116: }
117: MatZeroRows(C,4*m,rows,1.0,0,0);
119: PetscFree(rows);
120: VecAssemblyBegin(u);
121: VecAssemblyEnd(u);
122: VecAssemblyBegin(b);
123: VecAssemblyEnd(b);
125: { Mat A;
126: MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
127: MatDestroy(&C);
128: MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);
129: MatDestroy(&A);
130: }
132: /* Solve linear system */
133: KSPCreate(PETSC_COMM_WORLD,&ksp);
134: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
135: KSPSetFromOptions(ksp);
136: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
137: KSPSolve(ksp,b,u);
139: /* Check error */
140: VecGetOwnershipRange(ustar,&start,&end);
141: for (i=start; i<end; i++) {
142: x = h*(i % (m+1)); y = h*(i/(m+1));
143: val = y;
144: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
145: }
146: VecAssemblyBegin(ustar);
147: VecAssemblyEnd(ustar);
148: VecAXPY(u,-1.0,ustar);
149: VecNorm(u,NORM_2,&norm);
150: KSPGetIterationNumber(ksp,&its);
151: if (norm > tol){
152: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %G Iterations %D\n",norm*h,its);
153: }
155: /* Free work space */
156: KSPDestroy(&ksp);
157: VecDestroy(&ustar);
158: VecDestroy(&u);
159: VecDestroy(&b);
160: MatDestroy(&C);
161: PetscFinalize();
162: return 0;
163: }