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