Actual source code: ex29.c
2: /*
3: Added at the request of Marc Garbey.
5: Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation
7: -div \rho grad u = f, 0 < x,y < 1,
9: with forcing function
11: f = e^{-x^2/\nu} e^{-y^2/\nu}
13: with Dirichlet boundary conditions
15: u = f(x,y) for x = 0, x = 1, y = 0, y = 1
17: or pure Neumman boundary conditions
19: This uses multigrid to solve the linear system
20: */
22: static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";
24: #include <petscdm.h>
25: #include <petscdmda.h>
26: #include <petscksp.h>
28: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
29: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
31: typedef enum {DIRICHLET, NEUMANN} BCType;
33: typedef struct {
34: PetscReal rho;
35: PetscReal nu;
36: BCType bcType;
37: } UserContext;
39: int main(int argc,char **argv)
40: {
41: KSP ksp;
42: DM da;
43: UserContext user;
44: const char *bcTypes[2] = {"dirichlet","neumann"};
46: PetscInt bc;
47: Vec b,x;
48: PetscBool testsolver = PETSC_FALSE;
50: PetscInitialize(&argc,&argv,(char*)0,help);
51: KSPCreate(PETSC_COMM_WORLD,&ksp);
52: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,3,3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
53: DMSetFromOptions(da);
54: DMSetUp(da);
55: DMDASetUniformCoordinates(da,0,1,0,1,0,0);
56: DMDASetFieldName(da,0,"Pressure");
58: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
59: user.rho = 1.0;
60: PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL);
61: user.nu = 0.1;
62: PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL);
63: bc = (PetscInt)DIRICHLET;
64: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);
65: user.bcType = (BCType)bc;
66: PetscOptionsBool("-testsolver", "Run solver multiple times, useful for performance studies of solver", "ex29.c", testsolver, &testsolver, NULL);
67: PetscOptionsEnd();
69: KSPSetComputeRHS(ksp,ComputeRHS,&user);
70: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
71: KSPSetDM(ksp,da);
72: KSPSetFromOptions(ksp);
73: KSPSetUp(ksp);
74: KSPSolve(ksp,NULL,NULL);
76: if (testsolver) {
77: KSPGetSolution(ksp,&x);
78: KSPGetRhs(ksp,&b);
79: KSPSetDMActive(ksp,PETSC_FALSE);
80: KSPSolve(ksp,b,x);
81: {
82: #if defined(PETSC_USE_LOG)
83: PetscLogStage stage;
84: #endif
85: PetscInt i,n = 20;
87: PetscLogStageRegister("Solve only",&stage);
88: PetscLogStagePush(stage);
89: for (i=0; i<n; i++) {
90: KSPSolve(ksp,b,x);
91: }
92: PetscLogStagePop();
93: }
94: }
96: DMDestroy(&da);
97: KSPDestroy(&ksp);
98: PetscFinalize();
99: return 0;
100: }
102: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
103: {
104: UserContext *user = (UserContext*)ctx;
105: PetscInt i,j,mx,my,xm,ym,xs,ys;
106: PetscScalar Hx,Hy;
107: PetscScalar **array;
108: DM da;
111: KSPGetDM(ksp,&da);
112: DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
113: Hx = 1.0 / (PetscReal)(mx-1);
114: Hy = 1.0 / (PetscReal)(my-1);
115: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
116: DMDAVecGetArray(da, b, &array);
117: for (j=ys; j<ys+ym; j++) {
118: for (i=xs; i<xs+xm; i++) {
119: array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
120: }
121: }
122: DMDAVecRestoreArray(da, b, &array);
123: VecAssemblyBegin(b);
124: VecAssemblyEnd(b);
126: /* force right hand side to be consistent for singular matrix */
127: /* note this is really a hack, normally the model would provide you with a consistent right handside */
128: if (user->bcType == NEUMANN) {
129: MatNullSpace nullspace;
131: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
132: MatNullSpaceRemove(nullspace,b);
133: MatNullSpaceDestroy(&nullspace);
134: }
135: return 0;
136: }
138: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
139: {
141: if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
142: *rho = centerRho;
143: } else {
144: *rho = 1.0;
145: }
146: return 0;
147: }
149: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
150: {
151: UserContext *user = (UserContext*)ctx;
152: PetscReal centerRho;
153: PetscInt i,j,mx,my,xm,ym,xs,ys;
154: PetscScalar v[5];
155: PetscReal Hx,Hy,HydHx,HxdHy,rho;
156: MatStencil row, col[5];
157: DM da;
158: PetscBool check_matis = PETSC_FALSE;
161: KSPGetDM(ksp,&da);
162: centerRho = user->rho;
163: DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
164: Hx = 1.0 / (PetscReal)(mx-1);
165: Hy = 1.0 / (PetscReal)(my-1);
166: HxdHy = Hx/Hy;
167: HydHx = Hy/Hx;
168: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
169: for (j=ys; j<ys+ym; j++) {
170: for (i=xs; i<xs+xm; i++) {
171: row.i = i; row.j = j;
172: ComputeRho(i, j, mx, my, centerRho, &rho);
173: if (i==0 || j==0 || i==mx-1 || j==my-1) {
174: if (user->bcType == DIRICHLET) {
175: v[0] = 2.0*rho*(HxdHy + HydHx);
176: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
177: } else if (user->bcType == NEUMANN) {
178: PetscInt numx = 0, numy = 0, num = 0;
179: if (j!=0) {
180: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
181: numy++; num++;
182: }
183: if (i!=0) {
184: v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
185: numx++; num++;
186: }
187: if (i!=mx-1) {
188: v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
189: numx++; num++;
190: }
191: if (j!=my-1) {
192: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
193: numy++; num++;
194: }
195: v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i; col[num].j = j;
196: num++;
197: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
198: }
199: } else {
200: v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
201: v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
202: v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
203: v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
204: v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
205: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
206: }
207: }
208: }
209: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
210: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
211: MatViewFromOptions(jac,NULL,"-view_mat");
212: PetscOptionsGetBool(NULL,NULL,"-check_matis",&check_matis,NULL);
213: if (check_matis) {
214: void (*f)(void);
215: Mat J2;
216: MatType jtype;
217: PetscReal nrm;
219: MatGetType(jac,&jtype);
220: MatConvert(jac,MATIS,MAT_INITIAL_MATRIX,&J2);
221: MatViewFromOptions(J2,NULL,"-view_conv");
222: MatConvert(J2,jtype,MAT_INPLACE_MATRIX,&J2);
223: MatGetOperation(jac,MATOP_VIEW,&f);
224: MatSetOperation(J2,MATOP_VIEW,f);
225: MatSetDM(J2,da);
226: MatViewFromOptions(J2,NULL,"-view_conv_assembled");
227: MatAXPY(J2,-1.,jac,DIFFERENT_NONZERO_PATTERN);
228: MatNorm(J2,NORM_FROBENIUS,&nrm);
229: PetscPrintf(PETSC_COMM_WORLD,"Error MATIS %g\n",(double)nrm);
230: MatViewFromOptions(J2,NULL,"-view_conv_err");
231: MatDestroy(&J2);
232: }
233: if (user->bcType == NEUMANN) {
234: MatNullSpace nullspace;
236: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
237: MatSetNullSpace(J,nullspace);
238: MatNullSpaceDestroy(&nullspace);
239: }
240: return 0;
241: }
243: /*TEST
245: test:
246: args: -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -ksp_rtol 1.e-3
248: test:
249: suffix: 2
250: 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
251: requires: !single
253: test:
254: suffix: telescope
255: nsize: 4
256: 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
258: test:
259: suffix: 3
260: args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_pc_type jacobi
262: test:
263: suffix: 4
264: 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
266: test:
267: suffix: 5
268: nsize: 2
269: requires: hypre !complex
270: args: -pc_type mg -da_refine 2 -ksp_monitor -matptap_via hypre -pc_mg_galerkin both
272: test:
273: suffix: 6
274: args: -pc_type svd -pc_svd_monitor ::all
276: TEST*/