Actual source code: ex29.c

petsc-master 2020-12-02
Report Typos and Errors
  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;
 53:   PetscBool      testsolver = PETSC_FALSE;

 55:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 56:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 57:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,3,3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
 58:   DMSetFromOptions(da);
 59:   DMSetUp(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:   PetscOptionsBool("-testsolver", "Run solver multiple times, useful for performance studies of solver", "ex29.c", testsolver, &testsolver, NULL);
 72:   PetscOptionsEnd();

 74:   KSPSetComputeRHS(ksp,ComputeRHS,&user);
 75:   KSPSetComputeOperators(ksp,ComputeMatrix,&user);
 76:   KSPSetDM(ksp,da);
 77:   KSPSetFromOptions(ksp);
 78:   KSPSetUp(ksp);
 79:   KSPSolve(ksp,NULL,NULL);

 81:   if (testsolver) {
 82:     KSPGetSolution(ksp,&x);
 83:     KSPGetRhs(ksp,&b);
 84:     KSPSetDMActive(ksp,PETSC_FALSE);
 85:     KSPSolve(ksp,b,x);
 86:     {
 87: #if defined(PETSC_USE_LOG)
 88:       PetscLogStage stage;
 89: #endif
 90:       PetscInt      i,n = 20;

 92:       PetscLogStageRegister("Solve only",&stage);
 93:       PetscLogStagePush(stage);
 94:       for (i=0; i<n; i++) {
 95:         KSPSolve(ksp,b,x);
 96:       }
 97:       PetscLogStagePop();
 98:     }
 99:   }

101:   DMDestroy(&da);
102:   KSPDestroy(&ksp);
103:   PetscFinalize();
104:   return ierr;
105: }

107: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
108: {
109:   UserContext    *user = (UserContext*)ctx;
111:   PetscInt       i,j,mx,my,xm,ym,xs,ys;
112:   PetscScalar    Hx,Hy;
113:   PetscScalar    **array;
114:   DM             da;

117:   KSPGetDM(ksp,&da);
118:   DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
119:   Hx   = 1.0 / (PetscReal)(mx-1);
120:   Hy   = 1.0 / (PetscReal)(my-1);
121:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
122:   DMDAVecGetArray(da, b, &array);
123:   for (j=ys; j<ys+ym; j++) {
124:     for (i=xs; i<xs+xm; i++) {
125:       array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
126:     }
127:   }
128:   DMDAVecRestoreArray(da, b, &array);
129:   VecAssemblyBegin(b);
130:   VecAssemblyEnd(b);

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

137:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
138:     MatNullSpaceRemove(nullspace,b);
139:     MatNullSpaceDestroy(&nullspace);
140:   }
141:   return(0);
142: }


145: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
146: {
148:   if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
149:     *rho = centerRho;
150:   } else {
151:     *rho = 1.0;
152:   }
153:   return(0);
154: }

156: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
157: {
158:   UserContext    *user = (UserContext*)ctx;
159:   PetscReal      centerRho;
161:   PetscInt       i,j,mx,my,xm,ym,xs,ys;
162:   PetscScalar    v[5];
163:   PetscReal      Hx,Hy,HydHx,HxdHy,rho;
164:   MatStencil     row, col[5];
165:   DM             da;
166:   PetscBool      check_matis = PETSC_FALSE;

169:   KSPGetDM(ksp,&da);
170:   centerRho = user->rho;
171:   DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
172:   Hx        = 1.0 / (PetscReal)(mx-1);
173:   Hy        = 1.0 / (PetscReal)(my-1);
174:   HxdHy     = Hx/Hy;
175:   HydHx     = Hy/Hx;
176:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
177:   for (j=ys; j<ys+ym; j++) {
178:     for (i=xs; i<xs+xm; i++) {
179:       row.i = i; row.j = j;
180:       ComputeRho(i, j, mx, my, centerRho, &rho);
181:       if (i==0 || j==0 || i==mx-1 || j==my-1) {
182:         if (user->bcType == DIRICHLET) {
183:           v[0] = 2.0*rho*(HxdHy + HydHx);
184:           MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
185:         } else if (user->bcType == NEUMANN) {
186:           PetscInt numx = 0, numy = 0, num = 0;
187:           if (j!=0) {
188:             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j-1;
189:             numy++; num++;
190:           }
191:           if (i!=0) {
192:             v[num] = -rho*HydHx;              col[num].i = i-1; col[num].j = j;
193:             numx++; num++;
194:           }
195:           if (i!=mx-1) {
196:             v[num] = -rho*HydHx;              col[num].i = i+1; col[num].j = j;
197:             numx++; num++;
198:           }
199:           if (j!=my-1) {
200:             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j+1;
201:             numy++; num++;
202:           }
203:           v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i;   col[num].j = j;
204:           num++;
205:           MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
206:         }
207:       } else {
208:         v[0] = -rho*HxdHy;              col[0].i = i;   col[0].j = j-1;
209:         v[1] = -rho*HydHx;              col[1].i = i-1; col[1].j = j;
210:         v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
211:         v[3] = -rho*HydHx;              col[3].i = i+1; col[3].j = j;
212:         v[4] = -rho*HxdHy;              col[4].i = i;   col[4].j = j+1;
213:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
214:       }
215:     }
216:   }
217:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
218:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
219:   MatViewFromOptions(jac,NULL,"-view_mat");
220:   PetscOptionsGetBool(NULL,NULL,"-check_matis",&check_matis,NULL);
221:   if (check_matis) {
222:     void      (*f)(void);
223:     Mat       J2;
224:     MatType   jtype;
225:     PetscReal nrm;

227:     MatGetType(jac,&jtype);
228:     MatConvert(jac,MATIS,MAT_INITIAL_MATRIX,&J2);
229:     MatViewFromOptions(J2,NULL,"-view_conv");
230:     MatConvert(J2,jtype,MAT_INPLACE_MATRIX,&J2);
231:     MatGetOperation(jac,MATOP_VIEW,&f);
232:     MatSetOperation(J2,MATOP_VIEW,f);
233:     MatSetDM(J2,da);
234:     MatViewFromOptions(J2,NULL,"-view_conv_assembled");
235:     MatAXPY(J2,-1.,jac,DIFFERENT_NONZERO_PATTERN);
236:     MatNorm(J2,NORM_FROBENIUS,&nrm);
237:     PetscPrintf(PETSC_COMM_WORLD,"Error MATIS %g\n",(double)nrm);
238:     MatViewFromOptions(J2,NULL,"-view_conv_err");
239:     MatDestroy(&J2);
240:   }
241:   if (user->bcType == NEUMANN) {
242:     MatNullSpace nullspace;

244:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
245:     MatSetNullSpace(J,nullspace);
246:     MatNullSpaceDestroy(&nullspace);
247:   }
248:   return(0);
249: }


252: /*TEST

254:    test:
255:       args: -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -ksp_rtol 1.e-3

257:    test:
258:       suffix: 2
259:       args: -bc_type neumann -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -mg_coarse_pc_factor_shift_type nonzero
260:       requires: !single

262:    test:
263:       suffix: telescope
264:       nsize: 4
265:       args: -ksp_monitor_short -da_grid_x 257 -da_grid_y 257 -pc_type mg -pc_mg_galerkin pmat -pc_mg_levels 4 -ksp_type richardson -mg_levels_ksp_type chebyshev -mg_levels_pc_type jacobi -mg_coarse_pc_type telescope -mg_coarse_pc_telescope_ignore_kspcomputeoperators -mg_coarse_telescope_pc_type mg -mg_coarse_telescope_pc_mg_galerkin pmat -mg_coarse_telescope_pc_mg_levels 3 -mg_coarse_telescope_mg_levels_ksp_type chebyshev -mg_coarse_telescope_mg_levels_pc_type jacobi -mg_coarse_pc_telescope_reduction_factor 4

267:    test:
268:       suffix: 3
269:       args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_pc_type jacobi

271:    test:
272:       suffix: 4
273:       args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_ksp_max_it 3 -mg_levels_ksp_max_it 4

275:    test:
276:       suffix: 5
277:       nsize: 2
278:       requires: hypre !complex
279:       args: -pc_type mg  -da_refine 2 -ksp_monitor  -matptap_via hypre -pc_mg_galerkin both

281: TEST*/