Actual source code: ex29.c
petsc-3.13.6 2020-09-29
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;
55: PetscBool testsolver = PETSC_FALSE;
57: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
58: KSPCreate(PETSC_COMM_WORLD,&ksp);
59: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,3,3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
60: DMSetFromOptions(da);
61: DMSetUp(da);
62: DMDASetUniformCoordinates(da,0,1,0,1,0,0);
63: DMDASetFieldName(da,0,"Pressure");
65: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
66: user.rho = 1.0;
67: PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL);
68: user.nu = 0.1;
69: PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL);
70: bc = (PetscInt)DIRICHLET;
71: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);
72: user.bcType = (BCType)bc;
73: PetscOptionsBool("-testsolver", "Run solver multiple times, useful for performance studies of solver", "ex29.c", testsolver, &testsolver, NULL);
74: PetscOptionsEnd();
76: KSPSetComputeRHS(ksp,ComputeRHS,&user);
77: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
78: KSPSetDM(ksp,da);
79: KSPSetFromOptions(ksp);
80: KSPSetUp(ksp);
81: KSPSolve(ksp,NULL,NULL);
83: if (testsolver) {
84: KSPGetSolution(ksp,&x);
85: KSPGetRhs(ksp,&b);
86: KSPSetDMActive(ksp,PETSC_FALSE);
87: KSPSolve(ksp,b,x);
88: {
89: PetscLogStage stage;
90: PetscInt i,n = 20;
91: PetscLogStageRegister("Solve only",&stage);
92: PetscLogStagePush(stage);
93: for (i=0; i<n; i++) {
94: KSPSolve(ksp,b,x);
95: }
96: PetscLogStagePop();
97: }
98: }
100: DMDestroy(&da);
101: KSPDestroy(&ksp);
102: PetscFinalize();
103: return ierr;
104: }
106: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
107: {
108: UserContext *user = (UserContext*)ctx;
110: PetscInt i,j,mx,my,xm,ym,xs,ys;
111: PetscScalar Hx,Hy;
112: PetscScalar **array;
113: DM da;
116: KSPGetDM(ksp,&da);
117: DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
118: Hx = 1.0 / (PetscReal)(mx-1);
119: Hy = 1.0 / (PetscReal)(my-1);
120: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
121: DMDAVecGetArray(da, b, &array);
122: for (j=ys; j<ys+ym; j++) {
123: for (i=xs; i<xs+xm; i++) {
124: array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
125: }
126: }
127: DMDAVecRestoreArray(da, b, &array);
128: VecAssemblyBegin(b);
129: VecAssemblyEnd(b);
131: /* force right hand side to be consistent for singular matrix */
132: /* note this is really a hack, normally the model would provide you with a consistent right handside */
133: if (user->bcType == NEUMANN) {
134: MatNullSpace nullspace;
136: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
137: MatNullSpaceRemove(nullspace,b);
138: MatNullSpaceDestroy(&nullspace);
139: }
140: return(0);
141: }
144: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
145: {
147: if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
148: *rho = centerRho;
149: } else {
150: *rho = 1.0;
151: }
152: return(0);
153: }
155: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
156: {
157: UserContext *user = (UserContext*)ctx;
158: PetscReal centerRho;
160: PetscInt i,j,mx,my,xm,ym,xs,ys;
161: PetscScalar v[5];
162: PetscReal Hx,Hy,HydHx,HxdHy,rho;
163: MatStencil row, col[5];
164: DM da;
165: PetscBool check_matis = PETSC_FALSE;
168: KSPGetDM(ksp,&da);
169: centerRho = user->rho;
170: DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
171: Hx = 1.0 / (PetscReal)(mx-1);
172: Hy = 1.0 / (PetscReal)(my-1);
173: HxdHy = Hx/Hy;
174: HydHx = Hy/Hx;
175: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
176: for (j=ys; j<ys+ym; j++) {
177: for (i=xs; i<xs+xm; i++) {
178: row.i = i; row.j = j;
179: ComputeRho(i, j, mx, my, centerRho, &rho);
180: if (i==0 || j==0 || i==mx-1 || j==my-1) {
181: if (user->bcType == DIRICHLET) {
182: v[0] = 2.0*rho*(HxdHy + HydHx);
183: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
184: } else if (user->bcType == NEUMANN) {
185: PetscInt numx = 0, numy = 0, num = 0;
186: if (j!=0) {
187: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
188: numy++; num++;
189: }
190: if (i!=0) {
191: v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
192: numx++; num++;
193: }
194: if (i!=mx-1) {
195: v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
196: numx++; num++;
197: }
198: if (j!=my-1) {
199: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
200: numy++; num++;
201: }
202: v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i; col[num].j = j;
203: num++;
204: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
205: }
206: } else {
207: v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
208: v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
209: v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
210: v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
211: v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
212: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
213: }
214: }
215: }
216: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
217: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
218: MatViewFromOptions(jac,NULL,"-view_mat");
219: PetscOptionsGetBool(NULL,NULL,"-check_matis",&check_matis,NULL);
220: if (check_matis) {
221: void (*f)(void);
222: Mat J2;
223: MatType jtype;
224: PetscReal nrm;
226: MatGetType(jac,&jtype);
227: MatConvert(jac,MATIS,MAT_INITIAL_MATRIX,&J2);
228: MatViewFromOptions(J2,NULL,"-view_conv");
229: MatConvert(J2,jtype,MAT_INPLACE_MATRIX,&J2);
230: MatGetOperation(jac,MATOP_VIEW,&f);
231: MatSetOperation(J2,MATOP_VIEW,f);
232: MatSetDM(J2,da);
233: MatViewFromOptions(J2,NULL,"-view_conv_assembled");
234: MatAXPY(J2,-1.,jac,DIFFERENT_NONZERO_PATTERN);
235: MatNorm(J2,NORM_FROBENIUS,&nrm);
236: PetscPrintf(PETSC_COMM_WORLD,"Error MATIS %g\n",(double)nrm);
237: MatViewFromOptions(J2,NULL,"-view_conv_err");
238: MatDestroy(&J2);
239: }
240: if (user->bcType == NEUMANN) {
241: MatNullSpace nullspace;
243: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
244: MatSetNullSpace(J,nullspace);
245: MatNullSpaceDestroy(&nullspace);
246: }
247: return(0);
248: }
251: /*TEST
253: test:
254: args: -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -ksp_rtol 1.e-3
256: test:
257: suffix: 2
258: 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
259: requires: !single
261: test:
262: suffix: telescope
263: nsize: 4
264: 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
266: test:
267: suffix: 3
268: args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_pc_type jacobi
270: test:
271: suffix: 4
272: 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
274: test:
275: suffix: 5
276: nsize: 2
277: requires: hypre !complex
278: args: -pc_type mg -da_refine 2 -ksp_monitor -matptap_via hypre -pc_mg_galerkin both
280: TEST*/