Actual source code: ex20.c
1: /*$Id: ex20.c,v 1.19 2001/08/07 21:30:50 bsmith Exp $*/
3: static char help[] = "This example solves a linear system in parallel with KSP. The matrix\n\
4: uses simple bilinear elements on the unit square. To test the parallel\n\
5: matrix assembly,the matrix is intentionally laid out across processors\n\
6: differently from the way it is assembled. Input arguments are:\n\
7: -m <size> : problem size\n\n";
9: #include petscksp.h
13: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
14: {
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: }
24: int main(int argc,char **args)
25: {
26: Mat C;
27: int i,m = 5,rank,size,N,start,end,M;
28: int ierr,idx[4];
29: PetscTruth flg;
30: PetscScalar zero = 0.0,Ke[16], one = 1.0;
31: PetscReal h;
32: Vec u,b;
33: KSP ksp;
34: MatNullSpace nullsp;
35: PC pc;
37: PetscInitialize(&argc,&args,(char *)0,help);
38: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
39: N = (m+1)*(m+1); /* dimension of matrix */
40: M = m*m; /* number of elements */
41: h = 1.0/m; /* mesh width */
42: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
43: MPI_Comm_size(PETSC_COMM_WORLD,&size);
45: /* Create stiffness matrix */
46: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&C);
47: MatSetFromOptions(C);
48: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
49: end = start + M/size + ((M%size) > rank);
51: /* Assemble matrix */
52: FormElementStiffness(h*h,Ke); /* element stiffness for Laplacian */
53: for (i=start; i<end; i++) {
54: /* location of lower left corner of element */
55: /* node numbers for the four corners of element */
56: idx[0] = (m+1)*(i/m) + (i % m);
57: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
58: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
59: }
60: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
61: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
63: /* Create right-hand-side and solution vectors */
64: VecCreate(PETSC_COMM_WORLD,&u);
65: VecSetSizes(u,PETSC_DECIDE,N);
66: VecSetFromOptions(u);
67: PetscObjectSetName((PetscObject)u,"Approx. Solution");
68: VecDuplicate(u,&b);
69: PetscObjectSetName((PetscObject)b,"Right hand side");
71: VecSet(&one,u);
72: MatMult(C,u,b);
73: VecSet(&zero,u);
75: /* Solve linear system */
76: KSPCreate(PETSC_COMM_WORLD,&ksp);
77: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
78: KSPSetFromOptions(ksp);
79: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
81: PetscOptionsHasName(PETSC_NULL,"-fixnullspace",&flg);
82: if (flg) {
83: KSPGetPC(ksp,&pc);
84: MatNullSpaceCreate(PETSC_COMM_WORLD,1,0,PETSC_NULL,&nullsp);
85: PCNullSpaceAttach(pc,nullsp);
86: MatNullSpaceDestroy(nullsp);
87: }
89: KSPSetRhs(ksp,b);
90: KSPSetSolution(ksp,u);
91: KSPSolve(ksp);
94: /* Free work space */
95: KSPDestroy(ksp);
96: VecDestroy(u);
97: VecDestroy(b);
98: MatDestroy(C);
99: PetscFinalize();
100: return 0;
101: }