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