Actual source code: ex29.c

petsc-3.4.5 2014-06-29
  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>

 32: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,MatStructure*,void*);
 33: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);

 35: typedef enum {DIRICHLET, NEUMANN} BCType;

 37: typedef struct {
 38:   PetscReal rho;
 39:   PetscReal nu;
 40:   BCType    bcType;
 41: } UserContext;

 45: int main(int argc,char **argv)
 46: {
 47:   KSP            ksp;
 48:   DM             da;
 49:   UserContext    user;
 50:   const char     *bcTypes[2] = {"dirichlet","neumann"};
 52:   PetscInt       bc;
 53:   Vec            b,x;

 55:   PetscInitialize(&argc,&argv,(char*)0,help);

 57:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 58:   DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-3,-3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&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();

 85:   return 0;
 86: }

 90: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
 91: {
 92:   UserContext    *user = (UserContext*)ctx;
 94:   PetscInt       i,j,mx,my,xm,ym,xs,ys;
 95:   PetscScalar    Hx,Hy;
 96:   PetscScalar    **array;
 97:   DM             da;

100:   KSPGetDM(ksp,&da);
101:   DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
102:   Hx   = 1.0 / (PetscReal)(mx-1);
103:   Hy   = 1.0 / (PetscReal)(my-1);
104:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
105:   DMDAVecGetArray(da, b, &array);
106:   for (j=ys; j<ys+ym; j++) {
107:     for (i=xs; i<xs+xm; i++) {
108:       array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
109:     }
110:   }
111:   DMDAVecRestoreArray(da, b, &array);
112:   VecAssemblyBegin(b);
113:   VecAssemblyEnd(b);

115:   /* force right hand side to be consistent for singular matrix */
116:   /* note this is really a hack, normally the model would provide you with a consistent right handside */
117:   if (user->bcType == NEUMANN) {
118:     MatNullSpace nullspace;

120:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
121:     MatNullSpaceRemove(nullspace,b,NULL);
122:     MatNullSpaceDestroy(&nullspace);
123:   }
124:   return(0);
125: }


130: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *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 = centerRho;
135:   } else {
136:     *rho = 1.0;
137:   }
138:   return(0);
139: }

143: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,MatStructure *str,void *ctx)
144: {
145:   UserContext    *user = (UserContext*)ctx;
146:   PetscReal      centerRho;
148:   PetscInt       i,j,mx,my,xm,ym,xs,ys;
149:   PetscScalar    v[5];
150:   PetscReal      Hx,Hy,HydHx,HxdHy,rho;
151:   MatStencil     row, col[5];
152:   DM             da;

155:   KSPGetDM(ksp,&da);
156:   centerRho = user->rho;
157:   DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
158:   Hx        = 1.0 / (PetscReal)(mx-1);
159:   Hy        = 1.0 / (PetscReal)(my-1);
160:   HxdHy     = Hx/Hy;
161:   HydHx     = Hy/Hx;
162:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
163:   for (j=ys; j<ys+ym; j++) {
164:     for (i=xs; i<xs+xm; i++) {
165:       row.i = i; row.j = j;
166:       ComputeRho(i, j, mx, my, centerRho, &rho);
167:       if (i==0 || j==0 || i==mx-1 || j==my-1) {
168:         if (user->bcType == DIRICHLET) {
169:           v[0] = 2.0*rho*(HxdHy + HydHx);
170:           MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
171:         } else if (user->bcType == NEUMANN) {
172:           PetscInt numx = 0, numy = 0, num = 0;
173:           if (j!=0) {
174:             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j-1;
175:             numy++; num++;
176:           }
177:           if (i!=0) {
178:             v[num] = -rho*HydHx;              col[num].i = i-1; col[num].j = j;
179:             numx++; num++;
180:           }
181:           if (i!=mx-1) {
182:             v[num] = -rho*HydHx;              col[num].i = i+1; col[num].j = j;
183:             numx++; num++;
184:           }
185:           if (j!=my-1) {
186:             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j+1;
187:             numy++; num++;
188:           }
189:           v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i;   col[num].j = j;
190:           num++;
191:           MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
192:         }
193:       } else {
194:         v[0] = -rho*HxdHy;              col[0].i = i;   col[0].j = j-1;
195:         v[1] = -rho*HydHx;              col[1].i = i-1; col[1].j = j;
196:         v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
197:         v[3] = -rho*HydHx;              col[3].i = i+1; col[3].j = j;
198:         v[4] = -rho*HxdHy;              col[4].i = i;   col[4].j = j+1;
199:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
200:       }
201:     }
202:   }
203:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
204:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
205:   if (user->bcType == NEUMANN) {
206:     MatNullSpace nullspace;

208:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
209:     MatSetNullSpace(jac,nullspace);
210:     MatNullSpaceDestroy(&nullspace);
211:   }
212:   return(0);
213: }