Actual source code: ex20.c
petsc-3.4.2 2013-07-02
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>
11: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
12: {
13: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
14: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
15: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
16: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
17: return 0;
18: }
22: int main(int argc,char **args)
23: {
24: Mat C;
25: int i,m = 5,rank,size,N,start,end,M;
26: int ierr,idx[4];
27: PetscBool flg;
28: PetscScalar Ke[16];
29: PetscReal h;
30: Vec u,b;
31: KSP ksp;
32: MatNullSpace nullsp;
34: PetscInitialize(&argc,&args,(char*)0,help);
35: PetscOptionsGetInt(NULL,"-m",&m,NULL);
36: N = (m+1)*(m+1); /* dimension of matrix */
37: M = m*m; /* number of elements */
38: h = 1.0/m; /* mesh width */
39: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
40: MPI_Comm_size(PETSC_COMM_WORLD,&size);
42: /* Create stiffness matrix */
43: MatCreate(PETSC_COMM_WORLD,&C);
44: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
45: MatSetFromOptions(C);
46: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
47: end = start + M/size + ((M%size) > rank);
49: /* Assemble matrix */
50: FormElementStiffness(h*h,Ke); /* element stiffness for Laplacian */
51: for (i=start; i<end; i++) {
52: /* location of lower left corner of element */
53: /* node numbers for the four corners of element */
54: idx[0] = (m+1)*(i/m) + (i % m);
55: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
56: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
57: }
58: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
61: /* Create right-hand-side and solution vectors */
62: VecCreate(PETSC_COMM_WORLD,&u);
63: VecSetSizes(u,PETSC_DECIDE,N);
64: VecSetFromOptions(u);
65: PetscObjectSetName((PetscObject)u,"Approx. Solution");
66: VecDuplicate(u,&b);
67: PetscObjectSetName((PetscObject)b,"Right hand side");
69: VecSet(u,1.0);
70: MatMult(C,u,b);
71: VecSet(u,0.0);
73: /* Solve linear system */
74: KSPCreate(PETSC_COMM_WORLD,&ksp);
75: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
76: KSPSetFromOptions(ksp);
77: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
79: flg = PETSC_FALSE;
80: PetscOptionsGetBool(NULL,"-fixnullspace",&flg,NULL);
81: if (flg) {
82: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nullsp);
83: KSPSetNullSpace(ksp,nullsp);
84: MatNullSpaceDestroy(&&nullsp);
85: }
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 0;
96: }