Actual source code: ex29.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: Added at the request of Marc Garbey.
10: Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation
12: -div \rho grad u = f, 0 < x,y < 1,
14: with forcing function
16: f = e^{-x^2/\nu} e^{-y^2/\nu}
18: with Dirichlet boundary conditions
20: u = f(x,y) for x = 0, x = 1, y = 0, y = 1
22: or pure Neumman boundary conditions
24: This uses multigrid to solve the linear system
25: */
27: static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";
29: #include <petscdm.h>
30: #include <petscdmda.h>
31: #include <petscksp.h>
33: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
34: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
36: typedef enum {DIRICHLET, NEUMANN} BCType;
38: typedef struct {
39: PetscReal rho;
40: PetscReal nu;
41: BCType bcType;
42: } UserContext;
44: int main(int argc,char **argv)
45: {
46: KSP ksp;
47: DM da;
48: UserContext user;
49: const char *bcTypes[2] = {"dirichlet","neumann"};
51: PetscInt bc;
52: Vec b,x;
54: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
55: KSPCreate(PETSC_COMM_WORLD,&ksp);
56: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,3,3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
57: DMSetFromOptions(da);
58: DMSetUp(da);
59: DMDASetUniformCoordinates(da,0,1,0,1,0,0);
60: DMDASetFieldName(da,0,"Pressure");
62: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
63: user.rho = 1.0;
64: PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL);
65: user.nu = 0.1;
66: PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL);
67: bc = (PetscInt)DIRICHLET;
68: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);
69: user.bcType = (BCType)bc;
70: PetscOptionsEnd();
72: KSPSetComputeRHS(ksp,ComputeRHS,&user);
73: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
74: KSPSetDM(ksp,da);
75: KSPSetFromOptions(ksp);
76: KSPSetUp(ksp);
77: KSPSolve(ksp,NULL,NULL);
78: KSPGetSolution(ksp,&x);
79: KSPGetRhs(ksp,&b);
81: DMDestroy(&da);
82: KSPDestroy(&ksp);
83: PetscFinalize();
84: return ierr;
85: }
87: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
88: {
89: UserContext *user = (UserContext*)ctx;
91: PetscInt i,j,mx,my,xm,ym,xs,ys;
92: PetscScalar Hx,Hy;
93: PetscScalar **array;
94: DM da;
97: KSPGetDM(ksp,&da);
98: DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
99: Hx = 1.0 / (PetscReal)(mx-1);
100: Hy = 1.0 / (PetscReal)(my-1);
101: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
102: DMDAVecGetArray(da, b, &array);
103: for (j=ys; j<ys+ym; j++) {
104: for (i=xs; i<xs+xm; i++) {
105: array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
106: }
107: }
108: DMDAVecRestoreArray(da, b, &array);
109: VecAssemblyBegin(b);
110: VecAssemblyEnd(b);
112: /* force right hand side to be consistent for singular matrix */
113: /* note this is really a hack, normally the model would provide you with a consistent right handside */
114: if (user->bcType == NEUMANN) {
115: MatNullSpace nullspace;
117: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
118: MatNullSpaceRemove(nullspace,b);
119: MatNullSpaceDestroy(&nullspace);
120: }
121: return(0);
122: }
125: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
126: {
128: if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
129: *rho = centerRho;
130: } else {
131: *rho = 1.0;
132: }
133: return(0);
134: }
136: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
137: {
138: UserContext *user = (UserContext*)ctx;
139: PetscReal centerRho;
141: PetscInt i,j,mx,my,xm,ym,xs,ys;
142: PetscScalar v[5];
143: PetscReal Hx,Hy,HydHx,HxdHy,rho;
144: MatStencil row, col[5];
145: DM da;
148: KSPGetDM(ksp,&da);
149: centerRho = user->rho;
150: DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
151: Hx = 1.0 / (PetscReal)(mx-1);
152: Hy = 1.0 / (PetscReal)(my-1);
153: HxdHy = Hx/Hy;
154: HydHx = Hy/Hx;
155: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
156: for (j=ys; j<ys+ym; j++) {
157: for (i=xs; i<xs+xm; i++) {
158: row.i = i; row.j = j;
159: ComputeRho(i, j, mx, my, centerRho, &rho);
160: if (i==0 || j==0 || i==mx-1 || j==my-1) {
161: if (user->bcType == DIRICHLET) {
162: v[0] = 2.0*rho*(HxdHy + HydHx);
163: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
164: } else if (user->bcType == NEUMANN) {
165: PetscInt numx = 0, numy = 0, num = 0;
166: if (j!=0) {
167: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
168: numy++; num++;
169: }
170: if (i!=0) {
171: v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
172: numx++; num++;
173: }
174: if (i!=mx-1) {
175: v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
176: numx++; num++;
177: }
178: if (j!=my-1) {
179: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
180: numy++; num++;
181: }
182: v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i; col[num].j = j;
183: num++;
184: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
185: }
186: } else {
187: v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
188: v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
189: v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
190: v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
191: v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
192: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
193: }
194: }
195: }
196: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
197: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
198: if (user->bcType == NEUMANN) {
199: MatNullSpace nullspace;
201: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
202: MatSetNullSpace(J,nullspace);
203: MatNullSpaceDestroy(&nullspace);
204: }
205: return(0);
206: }