Actual source code: ex29.c
1: /*T
2: Concepts: KSP^solving a system of linear equations
3: Concepts: KSP^Laplacian, 2d
4: Processors: n
5: T*/
7: /*
8: Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation
10: div \kappa grad u = f, 0 < x,y < 1,
12: with forcing function
14: f = e^{-(1 - x)^2/\nu} e^{-(1 - y)^2/\nu}
16: with Dirichlet boundary conditions
18: u = f(x,y) for x = 0, x = 1, y = 0, y = 1
20: or pure Neumman boundary conditions
22: This uses multigrid to solve the linear system
23: */
25: static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";
27: #include petscda.h
28: #include petscksp.h
29: #include petscmg.h
34: typedef enum {DIRICHLET, NEUMANN} BCType;
36: typedef struct {
37: PetscScalar nu;
38: BCType bcType;
39: } UserContext;
43: int main(int argc,char **argv)
44: {
45: DMMG *dmmg;
46: DA da;
47: UserContext user;
48: PetscReal norm;
49: PetscScalar mone = -1.0;
50: const char *bcTypes[2] = {"dirichlet","neumann"};
52: PetscInt l,bc;
54: PetscInitialize(&argc,&argv,(char *)0,help);
56: DMMGCreate(PETSC_COMM_WORLD,3,PETSC_NULL,&dmmg);
57: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,3,3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
58: DMMGSetDM(dmmg,(DM)da);
59: DADestroy(da);
60: for (l = 0; l < DMMGGetLevels(dmmg); l++) {
61: DMMGSetUser(dmmg,l,&user);
62: }
64: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMMG");
65: user.nu = 0.1;
66: PetscOptionsScalar("-nu", "The width of the Gaussian source", "ex29.c", 0.1, &user.nu, PETSC_NULL);
67: bc = (PetscInt)DIRICHLET;
68: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,PETSC_NULL);
69: user.bcType = (BCType)bc;
70: PetscOptionsEnd();
72: DMMGSetKSP(dmmg,ComputeRHS,ComputeJacobian);
73: if (user.bcType == NEUMANN) {
74: DMMGSetNullSpace(dmmg,PETSC_TRUE,0,PETSC_NULL);
75: }
77: DMMGSolve(dmmg);
79: MatMult(DMMGGetJ(dmmg),DMMGGetx(dmmg),DMMGGetr(dmmg));
80: VecAXPY(&mone,DMMGGetb(dmmg),DMMGGetr(dmmg));
81: VecNorm(DMMGGetr(dmmg),NORM_2,&norm);
82: /* PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",norm); */
84: DMMGDestroy(dmmg);
85: PetscFinalize();
87: return 0;
88: }
92: PetscErrorCode ComputeRHS(DMMG dmmg, Vec b)
93: {
94: DA da = (DA)dmmg->dm;
95: UserContext *user = (UserContext *) dmmg->user;
97: PetscInt i,j,mx,my,xm,ym,xs,ys;
98: PetscScalar Hx,Hy;
99: PetscScalar **array;
102: DAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0);
103: Hx = 1.0 / (PetscReal)(mx-1);
104: Hy = 1.0 / (PetscReal)(my-1);
105: DAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
106: DAVecGetArray(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: DAVecRestoreArray(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: KSPGetNullSpace(dmmg->ksp,&nullspace);
122: MatNullSpaceRemove(nullspace,b,PETSC_NULL);
123: }
124: return(0);
125: }
127:
130: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscScalar *rho)
131: {
133: if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
134: *rho = 100.0;
135: } else {
136: *rho = 1.0;
137: }
138: return(0);
139: }
143: PetscErrorCode ComputeJacobian(DMMG dmmg, Mat jac)
144: {
145: DA da = (DA) dmmg->dm;
146: UserContext *user = (UserContext *) dmmg->user;
148: PetscInt i,j,mx,my,xm,ym,xs,ys,num;
149: PetscScalar v[5],Hx,Hy,HydHx,HxdHy,rho;
150: MatStencil row, col[5];
153: DAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0);
154: Hx = 1.0 / (PetscReal)(mx-1);
155: Hy = 1.0 / (PetscReal)(my-1);
156: HxdHy = Hx/Hy;
157: HydHx = Hy/Hx;
158: DAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
159: for (j=ys; j<ys+ym; j++){
160: for(i=xs; i<xs+xm; i++){
161: row.i = i; row.j = j;
162: if (i==0 || j==0 || i==mx-1 || j==my-1) {
163: ComputeRho(i, j, mx, my, &rho);
164: if (user->bcType == DIRICHLET) {
165: v[0] = 2.0*rho*(HxdHy + HydHx);
166: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
167: } else if (user->bcType == NEUMANN) {
168: num = 0;
169: if (j!=0) {
170: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
171: num++;
172: }
173: if (i!=0) {
174: v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
175: num++;
176: }
177: if (i!=mx-1) {
178: v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
179: num++;
180: }
181: if (j!=my-1) {
182: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
183: num++;
184: }
185: v[num] = (num/2.0)*rho*(HxdHy + HydHx); col[num].i = i; col[num].j = j;
186: num++;
187: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
188: }
189: } else {
190: v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
191: v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
192: v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
193: v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
194: v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
195: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
196: }
197: }
198: }
199: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
200: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
201: return(0);
202: }