Actual source code: ex32.c
petsc-3.8.4 2018-03-24
1: /*T
2: Concepts: KSP^solving a system of linear equations
3: Concepts: KSP^Laplacian, 2d
4: Processors: n
5: T*/
7: /*
8: Laplacian in 2D. Modeled by the partial differential equation
10: div grad u = f, 0 < x,y < 1,
12: with forcing function
14: f = e^{-(1 - x)^2/\nu} e^{-(1 - y)^2/\nu}
16: with pure Neumann boundary conditions
18: The functions are cell-centered
20: This uses multigrid to solve the linear system
22: Contributed by Andrei Draganescu <aidraga@sandia.gov>
24: Note the nice multigrid convergence despite the fact it is only using
25: peicewise constant interpolation/restriction. This is because cell-centered multigrid
26: does not need the same rule:
28: polynomial degree(interpolation) + polynomial degree(restriction) + 2 > degree of PDE
30: that vertex based multigrid needs.
31: */
33: static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";
35: #include <petscdm.h>
36: #include <petscdmda.h>
37: #include <petscksp.h>
39: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
40: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
42: typedef enum {DIRICHLET, NEUMANN} BCType;
44: typedef struct {
45: PetscScalar nu;
46: BCType bcType;
47: } UserContext;
49: int main(int argc,char **argv)
50: {
51: KSP ksp;
52: DM da;
53: UserContext user;
54: const char *bcTypes[2] = {"dirichlet","neumann"};
56: PetscInt bc;
58: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
59: KSPCreate(PETSC_COMM_WORLD,&ksp);
60: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,12,12,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
61: DMSetFromOptions(da);
62: DMSetUp(da);
63: DMDASetInterpolationType(da, DMDA_Q0);
65: KSPSetDM(ksp,da);
68: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DM");
69: user.nu = 0.1;
70: PetscOptionsScalar("-nu", "The width of the Gaussian source", "ex29.c", 0.1, &user.nu, NULL);
71: bc = (PetscInt)NEUMANN;
72: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);
73: user.bcType = (BCType)bc;
74: PetscOptionsEnd();
76: KSPSetComputeRHS(ksp,ComputeRHS,&user);
77: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
78: KSPSetFromOptions(ksp);
79: KSPSolve(ksp,NULL,NULL);
80: KSPDestroy(&ksp);
81: DMDestroy(&da);
82: PetscFinalize();
83: return ierr;
84: }
86: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
87: {
88: UserContext *user = (UserContext*)ctx;
90: PetscInt i,j,mx,my,xm,ym,xs,ys;
91: PetscScalar Hx,Hy;
92: PetscScalar **array;
93: DM da;
96: KSPGetDM(ksp,&da);
97: DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
98: Hx = 1.0 / (PetscReal)(mx);
99: Hy = 1.0 / (PetscReal)(my);
100: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
101: DMDAVecGetArray(da, b, &array);
102: for (j=ys; j<ys+ym; j++) {
103: for (i=xs; i<xs+xm; i++) {
104: array[j][i] = PetscExpScalar(-(((PetscReal)i+0.5)*Hx)*(((PetscReal)i+0.5)*Hx)/user->nu)*PetscExpScalar(-(((PetscReal)j+0.5)*Hy)*(((PetscReal)j+0.5)*Hy)/user->nu)*Hx*Hy;
105: }
106: }
107: DMDAVecRestoreArray(da, b, &array);
108: VecAssemblyBegin(b);
109: VecAssemblyEnd(b);
111: /* force right hand side to be consistent for singular matrix */
112: /* note this is really a hack, normally the model would provide you with a consistent right handside */
113: if (user->bcType == NEUMANN) {
114: MatNullSpace nullspace;
116: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
117: MatNullSpaceRemove(nullspace,b);
118: MatNullSpaceDestroy(&nullspace);
119: }
120: return(0);
121: }
124: PetscErrorCode ComputeMatrix(KSP ksp, Mat J,Mat jac, void *ctx)
125: {
126: UserContext *user = (UserContext*)ctx;
128: PetscInt i,j,mx,my,xm,ym,xs,ys,num, numi, numj;
129: PetscScalar v[5],Hx,Hy,HydHx,HxdHy;
130: MatStencil row, col[5];
131: DM da;
134: KSPGetDM(ksp,&da);
135: DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
136: Hx = 1.0 / (PetscReal)(mx);
137: Hy = 1.0 / (PetscReal)(my);
138: HxdHy = Hx/Hy;
139: HydHx = Hy/Hx;
140: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
141: for (j=ys; j<ys+ym; j++) {
142: for (i=xs; i<xs+xm; i++) {
143: row.i = i; row.j = j;
144: if (i==0 || j==0 || i==mx-1 || j==my-1) {
145: if (user->bcType == DIRICHLET) {
146: v[0] = 2.0*(HxdHy + HydHx);
147: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
148: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Dirichlet boundary conditions not supported !\n");
149: } else if (user->bcType == NEUMANN) {
150: num = 0; numi=0; numj=0;
151: if (j!=0) {
152: v[num] = -HxdHy;
153: col[num].i = i;
154: col[num].j = j-1;
155: num++; numj++;
156: }
157: if (i!=0) {
158: v[num] = -HydHx;
159: col[num].i = i-1;
160: col[num].j = j;
161: num++; numi++;
162: }
163: if (i!=mx-1) {
164: v[num] = -HydHx;
165: col[num].i = i+1;
166: col[num].j = j;
167: num++; numi++;
168: }
169: if (j!=my-1) {
170: v[num] = -HxdHy;
171: col[num].i = i;
172: col[num].j = j+1;
173: num++; numj++;
174: }
175: v[num] = (PetscReal)(numj)*HxdHy + (PetscReal)(numi)*HydHx; col[num].i = i; col[num].j = j;
176: num++;
177: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
178: }
179: } else {
180: v[0] = -HxdHy; col[0].i = i; col[0].j = j-1;
181: v[1] = -HydHx; col[1].i = i-1; col[1].j = j;
182: v[2] = 2.0*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
183: v[3] = -HydHx; col[3].i = i+1; col[3].j = j;
184: v[4] = -HxdHy; col[4].i = i; col[4].j = j+1;
185: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
186: }
187: }
188: }
189: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
190: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
191: if (user->bcType == NEUMANN) {
192: MatNullSpace nullspace;
194: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
195: MatSetNullSpace(J,nullspace);
196: MatNullSpaceDestroy(&nullspace);
197: }
198: return(0);
199: }