Actual source code: ex38.c
2: static char help[] = "Tests the aSA multigrid code.\n"
3: "Parameters:\n"
4: "-n n to use a matrix size of n\n";
6: #include petscda.h
7: #include petscksp.h
8: #include petscasa.h
10: PetscErrorCode Create1dLaplacian(PetscInt,Mat*);
11: PetscErrorCode CalculateRhs(Vec);
15: int main(int Argc,char **Args)
16: {
17: PetscInt n = 60;
18: PetscErrorCode ierr;
19: Mat cmat;
20: Vec b,x;
21: KSP kspmg;
22: PC pcmg;
23: DA da;
25: PetscInitialize(&Argc,&Args,(char *)0,help);
26:
27: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
28: Create1dLaplacian(n,&cmat);
29: MatGetVecs(cmat,&x,0);
30: MatGetVecs(cmat,&b,0);
31: CalculateRhs(b);
32: VecSet(x,0.0);
34: KSPCreate(PETSC_COMM_WORLD,&kspmg);
35: KSPSetType(kspmg, KSPCG);
37: KSPGetPC(kspmg,&pcmg);
38: PCSetType(pcmg,PCASA);
40: /* maybe user wants to override some of the choices */
41: KSPSetFromOptions(kspmg);
43: KSPSetOperators(kspmg,cmat,cmat,DIFFERENT_NONZERO_PATTERN);
45: DACreate1d(PETSC_COMM_WORLD, DA_NONPERIODIC, n, 1, 1, 0, &da);
46: DASetRefinementFactor(da, 3, 3, 3);
47: PCASASetDM(pcmg, (DM) da);
49: PCASASetTolerances(pcmg, 1.e-10, 1.e-10, PETSC_DEFAULT,
50: PETSC_DEFAULT);
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Solve the linear system
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: KSPSolve(kspmg,b,x);
57: KSPDestroy(kspmg);
58: VecDestroy(x);
59: VecDestroy(b);
60: MatDestroy(cmat);
61: DADestroy(da);
62: PetscFinalize();
63: return 0;
64: }
66: /* --------------------------------------------------------------------- */
69: PetscErrorCode Create1dLaplacian(PetscInt n,Mat *mat)
70: {
71: PetscScalar mone = -1.0,two = 2.0;
72: PetscInt i,j,loc_start,loc_end;
76: MatCreateMPIAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, n, n,
77: 1, PETSC_NULL, 2, PETSC_NULL, mat);
78:
79: MatGetOwnershipRange(*mat,&loc_start,&loc_end);
80: for (i=loc_start; i<loc_end; i++) {
81: if(i>0) { j=i-1; MatSetValues(*mat,1,&i,1,&j,&mone,INSERT_VALUES); }
82: MatSetValues(*mat,1,&i,1,&i,&two,INSERT_VALUES);
83: if(i<n-1) { j=i+1; MatSetValues(*mat,1,&i,1,&j,&mone,INSERT_VALUES); }
84: }
85: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
86: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
87: return(0);
88: }
89: /* --------------------------------------------------------------------- */
92: PetscErrorCode CalculateRhs(Vec u)
93: {
95: PetscInt i,n,loc_start,loc_end;
96: PetscReal h;
97: PetscScalar uu;
100: VecGetSize(u,&n);
101: VecGetOwnershipRange(u,&loc_start,&loc_end);
102: h = 1.0/((PetscReal)(n+1));
103: uu = 2.0*h*h;
104: for (i=loc_start; i<loc_end; i++) {
105: VecSetValues(u,1,&i,&uu,INSERT_VALUES);
106: }
107: VecAssemblyBegin(u);
108: VecAssemblyEnd(u);
110: return(0);
111: }