Actual source code: ex29.c
petsc-3.3-p7 2013-05-11
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 <petscdmda.h>
30: #include <petscksp.h>
31: #include <petscpcmg.h>
33: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,MatStructure*,void*);
34: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
36: typedef enum {DIRICHLET, NEUMANN} BCType;
38: typedef struct {
39: PetscScalar rho;
40: PetscScalar 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, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-3,-3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
61: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
62: user.rho = 1.0;
63: PetscOptionsScalar("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, PETSC_NULL);
64: user.nu = 0.1;
65: PetscOptionsScalar("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, PETSC_NULL);
66: bc = (PetscInt)DIRICHLET;
67: PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,PETSC_NULL);
68: user.bcType = (BCType)bc;
69: PetscOptionsEnd();
71: KSPSetComputeRHS(ksp,ComputeRHS,&user);
72: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
73: KSPSetDM(ksp,da);
74: KSPSetFromOptions(ksp);
75: KSPSetUp(ksp);
76: KSPSolve(ksp,PETSC_NULL,PETSC_NULL);
77: KSPGetSolution(ksp,&x);
78: KSPGetRhs(ksp,&b);
80: DMDestroy(&da);
81: KSPDestroy(&ksp);
82: PetscFinalize();
84: return 0;
85: }
89: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
90: {
91: UserContext *user = (UserContext*)ctx;
93: PetscInt i,j,mx,my,xm,ym,xs,ys;
94: PetscScalar Hx,Hy;
95: PetscScalar **array;
96: DM da;
99: KSPGetDM(ksp,&da);
100: DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
101: Hx = 1.0 / (PetscReal)(mx-1);
102: Hy = 1.0 / (PetscReal)(my-1);
103: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
104: DMDAVecGetArray(da, b, &array);
105: for (j=ys; j<ys+ym; j++){
106: for(i=xs; i<xs+xm; i++){
107: array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
108: }
109: }
110: DMDAVecRestoreArray(da, b, &array);
111: VecAssemblyBegin(b);
112: VecAssemblyEnd(b);
114: /* force right hand side to be consistent for singular matrix */
115: /* note this is really a hack, normally the model would provide you with a consistent right handside */
116: if (user->bcType == NEUMANN) {
117: MatNullSpace nullspace;
119: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
120: MatNullSpaceRemove(nullspace,b,PETSC_NULL);
121: MatNullSpaceDestroy(&nullspace);
122: }
123: return(0);
124: }
126:
129: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscScalar centerRho, PetscScalar *rho)
130: {
132: if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
133: *rho = centerRho;
134: } else {
135: *rho = 1.0;
136: }
137: return(0);
138: }
142: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,MatStructure *str,void *ctx)
143: {
144: UserContext *user = (UserContext*)ctx;
145: PetscScalar centerRho;
147: PetscInt i,j,mx,my,xm,ym,xs,ys,num;
148: PetscScalar v[5],Hx,Hy,HydHx,HxdHy,rho;
149: MatStencil row, col[5];
150: DM da;
153: KSPGetDM(ksp,&da);
154: centerRho = user->rho;
155: DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
156: Hx = 1.0 / (PetscReal)(mx-1);
157: Hy = 1.0 / (PetscReal)(my-1);
158: HxdHy = Hx/Hy;
159: HydHx = Hy/Hx;
160: DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
161: for (j=ys; j<ys+ym; j++){
162: for(i=xs; i<xs+xm; i++){
163: row.i = i; row.j = j;
164: ComputeRho(i, j, mx, my, centerRho, &rho);
165: if (i==0 || j==0 || i==mx-1 || j==my-1) {
166: if (user->bcType == DIRICHLET) {
167: v[0] = 2.0*rho*(HxdHy + HydHx);
168: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
169: } else if (user->bcType == NEUMANN) {
170: num = 0;
171: if (j!=0) {
172: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
173: num++;
174: }
175: if (i!=0) {
176: v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
177: num++;
178: }
179: if (i!=mx-1) {
180: v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
181: num++;
182: }
183: if (j!=my-1) {
184: v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
185: num++;
186: }
187: v[num] = (num/2.0)*rho*(HxdHy + HydHx); col[num].i = i; col[num].j = j;
188: num++;
189: MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
190: }
191: } else {
192: v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
193: v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
194: v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
195: v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
196: v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
197: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
198: }
199: }
200: }
201: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
202: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
203: if (user->bcType == NEUMANN) {
204: MatNullSpace nullspace;
206: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
207: MatSetNullSpace(jac,nullspace);
208: MatNullSpaceDestroy(&nullspace);
209: }
211: return(0);
212: }