Actual source code: ex29.c
petsc-3.5.4 2015-05-23
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;
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);
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: DMDASetUniformCoordinates(da,0,1,0,1,0,0);
61: DMDASetFieldName(da,0,"Pressure");
63: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
64: user.rho = 1.0;
65: PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL);
66: user.nu = 0.1;
67: PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL);
68: bc = (PetscInt)DIRICHLET;
69: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);
70: user.bcType = (BCType)bc;
71: PetscOptionsEnd();
73: KSPSetComputeRHS(ksp,ComputeRHS,&user);
74: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
75: KSPSetDM(ksp,da);
76: KSPSetFromOptions(ksp);
77: KSPSetUp(ksp);
78: KSPSolve(ksp,NULL,NULL);
79: KSPGetSolution(ksp,&x);
80: KSPGetRhs(ksp,&b);
82: DMDestroy(&da);
83: KSPDestroy(&ksp);
84: PetscFinalize();
86: return 0;
87: }
91: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
92: {
93: UserContext *user = (UserContext*)ctx;
95: PetscInt i,j,mx,my,xm,ym,xs,ys;
96: PetscScalar Hx,Hy;
97: PetscScalar **array;
98: DM da;
101: KSPGetDM(ksp,&da);
102: DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
103: Hx = 1.0 / (PetscReal)(mx-1);
104: Hy = 1.0 / (PetscReal)(my-1);
105: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
106: DMDAVecGetArray(da, b, &array);
107: for (j=ys; j<ys+ym; j++) {
108: for (i=xs; i<xs+xm; i++) {
109: array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
110: }
111: }
112: DMDAVecRestoreArray(da, b, &array);
113: VecAssemblyBegin(b);
114: VecAssemblyEnd(b);
116: /* force right hand side to be consistent for singular matrix */
117: /* note this is really a hack, normally the model would provide you with a consistent right handside */
118: if (user->bcType == NEUMANN) {
119: MatNullSpace nullspace;
121: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
122: MatNullSpaceRemove(nullspace,b);
123: MatNullSpaceDestroy(&nullspace);
124: }
125: return(0);
126: }
131: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
132: {
134: if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
135: *rho = centerRho;
136: } else {
137: *rho = 1.0;
138: }
139: return(0);
140: }
144: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
145: {
146: UserContext *user = (UserContext*)ctx;
147: PetscReal centerRho;
149: PetscInt i,j,mx,my,xm,ym,xs,ys;
150: PetscScalar v[5];
151: PetscReal Hx,Hy,HydHx,HxdHy,rho;
152: MatStencil row, col[5];
153: DM da;
156: KSPGetDM(ksp,&da);
157: centerRho = user->rho;
158: DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
159: Hx = 1.0 / (PetscReal)(mx-1);
160: Hy = 1.0 / (PetscReal)(my-1);
161: HxdHy = Hx/Hy;
162: HydHx = Hy/Hx;
163: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
164: for (j=ys; j<ys+ym; j++) {
165: for (i=xs; i<xs+xm; i++) {
166: row.i = i; row.j = j;
167: ComputeRho(i, j, mx, my, centerRho, &rho);
168: if (i==0 || j==0 || i==mx-1 || j==my-1) {
169: if (user->bcType == DIRICHLET) {
170: v[0] = 2.0*rho*(HxdHy + HydHx);
171: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
172: } else if (user->bcType == NEUMANN) {
173: PetscInt numx = 0, numy = 0, num = 0;
174: if (j!=0) {
175: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
176: numy++; num++;
177: }
178: if (i!=0) {
179: v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
180: numx++; num++;
181: }
182: if (i!=mx-1) {
183: v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
184: numx++; num++;
185: }
186: if (j!=my-1) {
187: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
188: numy++; num++;
189: }
190: v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i; col[num].j = j;
191: num++;
192: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
193: }
194: } else {
195: v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
196: v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
197: v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
198: v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
199: v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
200: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
201: }
202: }
203: }
204: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
205: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
206: if (user->bcType == NEUMANN) {
207: MatNullSpace nullspace;
209: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
210: MatSetNullSpace(jac,nullspace);
211: MatNullSpaceDestroy(&nullspace);
212: }
213: return(0);
214: }