Actual source code: ex29.c
petsc-3.9.4 2018-09-11
1: /*T
2: Concepts: KSP^solving a system of linear equations
3: Concepts: KSP^Laplacian, 2d
4: Processors: n
5: T*/
9: /*
10: Added at the request of Marc Garbey.
12: Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation
14: -div \rho grad u = f, 0 < x,y < 1,
16: with forcing function
18: f = e^{-x^2/\nu} e^{-y^2/\nu}
20: with Dirichlet boundary conditions
22: u = f(x,y) for x = 0, x = 1, y = 0, y = 1
24: or pure Neumman boundary conditions
26: This uses multigrid to solve the linear system
27: */
29: static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";
31: #include <petscdm.h>
32: #include <petscdmda.h>
33: #include <petscksp.h>
35: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
36: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
38: typedef enum {DIRICHLET, NEUMANN} BCType;
40: typedef struct {
41: PetscReal rho;
42: PetscReal nu;
43: BCType bcType;
44: } UserContext;
46: int main(int argc,char **argv)
47: {
48: KSP ksp;
49: DM da;
50: UserContext user;
51: const char *bcTypes[2] = {"dirichlet","neumann"};
53: PetscInt bc;
54: Vec b,x;
56: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
57: KSPCreate(PETSC_COMM_WORLD,&ksp);
58: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,3,3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
59: DMSetFromOptions(da);
60: DMSetUp(da);
61: DMDASetUniformCoordinates(da,0,1,0,1,0,0);
62: DMDASetFieldName(da,0,"Pressure");
64: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
65: user.rho = 1.0;
66: PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL);
67: user.nu = 0.1;
68: PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL);
69: bc = (PetscInt)DIRICHLET;
70: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);
71: user.bcType = (BCType)bc;
72: PetscOptionsEnd();
74: KSPSetComputeRHS(ksp,ComputeRHS,&user);
75: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
76: KSPSetDM(ksp,da);
77: KSPSetFromOptions(ksp);
78: KSPSetUp(ksp);
79: KSPSolve(ksp,NULL,NULL);
80: KSPGetSolution(ksp,&x);
81: KSPGetRhs(ksp,&b);
83: DMDestroy(&da);
84: KSPDestroy(&ksp);
85: PetscFinalize();
86: return ierr;
87: }
89: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
90: {
91: UserContext *user = (UserContext*)ctx;
93: PetscInt i,j,mx,my,xm,ym,xs,ys;
94: PetscScalar Hx,Hy;
95: PetscScalar **array;
96: DM da;
99: KSPGetDM(ksp,&da);
100: DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
101: Hx = 1.0 / (PetscReal)(mx-1);
102: Hy = 1.0 / (PetscReal)(my-1);
103: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
104: DMDAVecGetArray(da, b, &array);
105: for (j=ys; j<ys+ym; j++) {
106: for (i=xs; i<xs+xm; i++) {
107: array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
108: }
109: }
110: DMDAVecRestoreArray(da, b, &array);
111: VecAssemblyBegin(b);
112: VecAssemblyEnd(b);
114: /* force right hand side to be consistent for singular matrix */
115: /* note this is really a hack, normally the model would provide you with a consistent right handside */
116: if (user->bcType == NEUMANN) {
117: MatNullSpace nullspace;
119: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
120: MatNullSpaceRemove(nullspace,b);
121: MatNullSpaceDestroy(&nullspace);
122: }
123: return(0);
124: }
127: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
128: {
130: if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
131: *rho = centerRho;
132: } else {
133: *rho = 1.0;
134: }
135: return(0);
136: }
138: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
139: {
140: UserContext *user = (UserContext*)ctx;
141: PetscReal centerRho;
143: PetscInt i,j,mx,my,xm,ym,xs,ys;
144: PetscScalar v[5];
145: PetscReal Hx,Hy,HydHx,HxdHy,rho;
146: MatStencil row, col[5];
147: DM da;
150: KSPGetDM(ksp,&da);
151: centerRho = user->rho;
152: DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
153: Hx = 1.0 / (PetscReal)(mx-1);
154: Hy = 1.0 / (PetscReal)(my-1);
155: HxdHy = Hx/Hy;
156: HydHx = Hy/Hx;
157: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
158: for (j=ys; j<ys+ym; j++) {
159: for (i=xs; i<xs+xm; i++) {
160: row.i = i; row.j = j;
161: ComputeRho(i, j, mx, my, centerRho, &rho);
162: if (i==0 || j==0 || i==mx-1 || j==my-1) {
163: if (user->bcType == DIRICHLET) {
164: v[0] = 2.0*rho*(HxdHy + HydHx);
165: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
166: } else if (user->bcType == NEUMANN) {
167: PetscInt numx = 0, numy = 0, num = 0;
168: if (j!=0) {
169: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
170: numy++; num++;
171: }
172: if (i!=0) {
173: v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
174: numx++; num++;
175: }
176: if (i!=mx-1) {
177: v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
178: numx++; num++;
179: }
180: if (j!=my-1) {
181: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
182: numy++; num++;
183: }
184: v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i; col[num].j = j;
185: num++;
186: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
187: }
188: } else {
189: v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
190: v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
191: v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
192: v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
193: v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
194: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
195: }
196: }
197: }
198: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
199: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
200: if (user->bcType == NEUMANN) {
201: MatNullSpace nullspace;
203: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
204: MatSetNullSpace(J,nullspace);
205: MatNullSpaceDestroy(&nullspace);
206: }
207: return(0);
208: }
211: /*TEST
213: test:
214: args: -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -ksp_rtol 1.e-3
216: test:
217: suffix: 2
218: args: -bc_type neumann -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -mg_coarse_pc_factor_shift_type nonzero
219: requires: !single
221: test:
222: suffix: telescope
223: nsize: 4
224: args: -ksp_monitor_short -da_grid_x 257 -da_grid_y 257 -pc_type mg -pc_mg_galerkin pmat -pc_mg_levels 4 -ksp_type richardson -mg_levels_ksp_type chebyshev -mg_levels_pc_type jacobi -mg_coarse_pc_type telescope -mg_coarse_pc_telescope_ignore_kspcomputeoperators -mg_coarse_telescope_pc_type mg -mg_coarse_telescope_pc_mg_galerkin pmat -mg_coarse_telescope_pc_mg_levels 3 -mg_coarse_telescope_mg_levels_ksp_type chebyshev -mg_coarse_telescope_mg_levels_pc_type jacobi -mg_coarse_pc_telescope_reduction_factor 4
226: test:
227: suffix: 3
228: args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_pc_type jacobi
230: test:
231: suffix: 4
232: args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_ksp_max_it 3 -mg_levels_ksp_max_it 4
234: TEST*/