Actual source code: ex25.c
petsc-3.8.4 2018-03-24
2: /*
3: Partial differential equation
5: d (1 + e*sine(2*pi*k*x)) d u = 1, 0 < x < 1,
6: -- ---
7: dx dx
8: with boundary conditions
10: u = 0 for x = 0, x = 1
12: This uses multigrid to solve the linear system
14: */
16: static char help[] = "Solves 1D variable coefficient Laplacian using multigrid.\n\n";
18: #include <petscdm.h>
19: #include <petscdmda.h>
20: #include <petscksp.h>
22: static PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
23: static PetscErrorCode ComputeRHS(KSP,Vec,void*);
25: typedef struct {
26: PetscInt k;
27: PetscScalar e;
28: } AppCtx;
30: int main(int argc,char **argv)
31: {
33: KSP ksp;
34: DM da;
35: AppCtx user;
36: Mat A;
37: Vec b,b2;
38: Vec x;
39: PetscReal nrm;
41: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
42: user.k = 1;
43: user.e = .99;
44: PetscOptionsGetInt(NULL,0,"-k",&user.k,0);
45: PetscOptionsGetScalar(NULL,0,"-e",&user.e,0);
47: KSPCreate(PETSC_COMM_WORLD,&ksp);
48: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,128,1,1,0,&da);
49: DMSetFromOptions(da);
50: DMSetUp(da);
51: KSPSetDM(ksp,da);
52: KSPSetComputeRHS(ksp,ComputeRHS,&user);
53: KSPSetComputeOperators(ksp,ComputeMatrix,&user);
54: KSPSetFromOptions(ksp);
55: KSPSolve(ksp,NULL,NULL);
57: KSPGetOperators(ksp,&A,NULL);
58: KSPGetSolution(ksp,&x);
59: KSPGetRhs(ksp,&b);
60: VecDuplicate(b,&b2);
61: MatMult(A,x,b2);
62: VecAXPY(b2,-1.0,b);
63: VecNorm(b2,NORM_MAX,&nrm);
64: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)nrm);
66: VecDestroy(&b2);
67: KSPDestroy(&ksp);
68: DMDestroy(&da);
69: PetscFinalize();
70: return ierr;
71: }
73: static PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
74: {
76: PetscInt mx,idx[2];
77: PetscScalar h,v[2];
78: DM da;
81: KSPGetDM(ksp,&da);
82: DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
83: h = 1.0/((mx-1));
84: VecSet(b,h);
85: idx[0] = 0; idx[1] = mx -1;
86: v[0] = v[1] = 0.0;
87: VecSetValues(b,2,idx,v,INSERT_VALUES);
88: VecAssemblyBegin(b);
89: VecAssemblyEnd(b);
90: return(0);
91: }
93: static PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
94: {
95: AppCtx *user = (AppCtx*)ctx;
97: PetscInt i,mx,xm,xs;
98: PetscScalar v[3],h,xlow,xhigh;
99: MatStencil row,col[3];
100: DM da;
103: KSPGetDM(ksp,&da);
104: DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
105: DMDAGetCorners(da,&xs,0,0,&xm,0,0);
106: h = 1.0/(mx-1);
108: for (i=xs; i<xs+xm; i++) {
109: row.i = i;
110: if (i==0 || i==mx-1) {
111: v[0] = 2.0/h;
112: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
113: } else {
114: xlow = h*(PetscReal)i - .5*h;
115: xhigh = xlow + h;
116: v[0] = (-1.0 - user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xlow))/h;col[0].i = i-1;
117: v[1] = (2.0 + user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xlow) + user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xhigh))/h;col[1].i = row.i;
118: v[2] = (-1.0 - user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xhigh))/h;col[2].i = i+1;
119: MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);
120: }
121: }
122: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
123: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
124: return(0);
125: }