Actual source code: ex20.c
petsc-3.8.4 2018-03-24
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: #include <petscksp.h>
9: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
10: {
11: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
12: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
13: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
14: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
15: return 0;
16: }
18: int main(int argc,char **args)
19: {
20: Mat C;
21: int i,m = 5,rank,size,N,start,end,M;
22: int ierr,idx[4];
23: PetscScalar Ke[16];
24: PetscReal h;
25: Vec u,b;
26: KSP ksp;
27: MatNullSpace nullsp;
29: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
30: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
31: N = (m+1)*(m+1); /* dimension of matrix */
32: M = m*m; /* number of elements */
33: h = 1.0/m; /* mesh width */
34: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
35: MPI_Comm_size(PETSC_COMM_WORLD,&size);
37: /* Create stiffness matrix */
38: MatCreate(PETSC_COMM_WORLD,&C);
39: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
40: MatSetFromOptions(C);
41: MatSetUp(C);
42: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
43: end = start + M/size + ((M%size) > rank);
45: /* Assemble matrix */
46: FormElementStiffness(h*h,Ke); /* element stiffness for Laplacian */
47: for (i=start; i<end; i++) {
48: /* location of lower left corner of element */
49: /* node numbers for the four corners of element */
50: idx[0] = (m+1)*(i/m) + (i % m);
51: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
52: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
53: }
54: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
57: /* Create right-hand-side and solution vectors */
58: VecCreate(PETSC_COMM_WORLD,&u);
59: VecSetSizes(u,PETSC_DECIDE,N);
60: VecSetFromOptions(u);
61: PetscObjectSetName((PetscObject)u,"Approx. Solution");
62: VecDuplicate(u,&b);
63: PetscObjectSetName((PetscObject)b,"Right hand side");
65: VecSet(b,1.0);
66: VecSetValue(b,0,1.2,ADD_VALUES);
67: VecSet(u,0.0);
69: /* Solve linear system */
70: KSPCreate(PETSC_COMM_WORLD,&ksp);
71: KSPSetOperators(ksp,C,C);
72: KSPSetFromOptions(ksp);
73: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
75: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nullsp);
76: /*
77: The KSP solver will remove this nullspace from the solution at each iteration
78: */
79: MatSetNullSpace(C,nullsp);
80: /*
81: The KSP solver will remove from the right hand side any portion in this nullspace, thus making the linear system consistent.
82: */
83: MatSetTransposeNullSpace(C,nullsp);
84: MatNullSpaceDestroy(&nullsp);
86: KSPSolve(ksp,b,u);
89: /* Free work space */
90: KSPDestroy(&ksp);
91: VecDestroy(&u);
92: VecDestroy(&b);
93: MatDestroy(&C);
94: PetscFinalize();
95: return ierr;
96: }