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 SLES. The matrixn
4: uses simple bilinear elements on the unit square. To test the paralleln
5: matrix assembly,the matrix is intentionally laid out across processorsn
6: differently from the way it is assembled. Input arguments are:n
7: -m <size> : problem sizenn";
9: #include petscsles.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: }
20: int main(int argc,char **args)
21: {
22: Mat C;
23: int i,m = 5,rank,size,N,start,end,M,its;
24: int ierr,idx[4];
25: PetscTruth flg;
26: PetscScalar zero = 0.0,Ke[16];
27: PetscReal h;
28: Vec u,b;
29: SLES sles;
30: KSP ksp;
31: MatNullSpace nullsp;
32: PC pc;
33: PetscRandom rand;
35: PetscInitialize(&argc,&args,(char *)0,help);
36: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
37: N = (m+1)*(m+1); /* dimension of matrix */
38: M = m*m; /* number of elements */
39: h = 1.0/m; /* mesh width */
40: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
41: MPI_Comm_size(PETSC_COMM_WORLD,&size);
43: /* Create stiffness matrix */
44: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&C);
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: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT_REAL,&rand);
70: VecSetRandom(rand,u);
71: PetscRandomDestroy(rand);
72: MatMult(C,u,b);
73: VecSet(&zero,u);
75: /* Solve linear system */
76: SLESCreate(PETSC_COMM_WORLD,&sles);
77: SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
78: SLESSetFromOptions(sles);
79: SLESGetKSP(sles,&ksp);
80: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
82: PetscOptionsHasName(PETSC_NULL,"-fixnullspace",&flg);
83: if (flg) {
84: SLESGetPC(sles,&pc);
85: MatNullSpaceCreate(PETSC_COMM_WORLD,1,0,PETSC_NULL,&nullsp);
86: PCNullSpaceAttach(pc,nullsp);
87: MatNullSpaceDestroy(nullsp);
88: }
90: SLESSolve(sles,b,u,&its);
93: /* Free work space */
94: SLESDestroy(sles);
95: VecDestroy(u);
96: VecDestroy(b);
97: MatDestroy(C);
98: PetscFinalize();
99: return 0;
100: }