xref: /petsc/src/ksp/pc/tests/ex5.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests the multigrid code.  The input parameters are:\n\
3*c4762a1bSJed Brown   -x N              Use a mesh in the x direction of N.  \n\
4*c4762a1bSJed Brown   -c N              Use N V-cycles.  \n\
5*c4762a1bSJed Brown   -l N              Use N Levels.  \n\
6*c4762a1bSJed Brown   -smooths N        Use N pre smooths and N post smooths.  \n\
7*c4762a1bSJed Brown   -j                Use Jacobi smoother.  \n\
8*c4762a1bSJed Brown   -a use additive multigrid \n\
9*c4762a1bSJed Brown   -f use full multigrid (preconditioner variant) \n\
10*c4762a1bSJed Brown This example also demonstrates matrix-free methods\n\n";
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown /*
13*c4762a1bSJed Brown   This is not a good example to understand the use of multigrid with PETSc.
14*c4762a1bSJed Brown */
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown #include <petscksp.h>
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown PetscErrorCode  residual(Mat,Vec,Vec,Vec);
19*c4762a1bSJed Brown PetscErrorCode  gauss_seidel(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
20*c4762a1bSJed Brown PetscErrorCode  jacobi_smoother(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
21*c4762a1bSJed Brown PetscErrorCode  interpolate(Mat,Vec,Vec,Vec);
22*c4762a1bSJed Brown PetscErrorCode  restrct(Mat,Vec,Vec);
23*c4762a1bSJed Brown PetscErrorCode  Create1dLaplacian(PetscInt,Mat*);
24*c4762a1bSJed Brown PetscErrorCode  CalculateRhs(Vec);
25*c4762a1bSJed Brown PetscErrorCode  CalculateError(Vec,Vec,Vec,PetscReal*);
26*c4762a1bSJed Brown PetscErrorCode  CalculateSolution(PetscInt,Vec*);
27*c4762a1bSJed Brown PetscErrorCode  amult(Mat,Vec,Vec);
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown int main(int Argc,char **Args)
30*c4762a1bSJed Brown {
31*c4762a1bSJed Brown   PetscInt       x_mesh = 15,levels = 3,cycles = 1,use_jacobi = 0;
32*c4762a1bSJed Brown   PetscInt       i,smooths = 1,*N,its;
33*c4762a1bSJed Brown   PetscErrorCode ierr;
34*c4762a1bSJed Brown   PCMGType       am = PC_MG_MULTIPLICATIVE;
35*c4762a1bSJed Brown   Mat            cmat,mat[20],fmat;
36*c4762a1bSJed Brown   KSP            cksp,ksp[20],kspmg;
37*c4762a1bSJed Brown   PetscReal      e[3];  /* l_2 error,max error, residual */
38*c4762a1bSJed Brown   const char     *shellname;
39*c4762a1bSJed Brown   Vec            x,solution,X[20],R[20],B[20];
40*c4762a1bSJed Brown   PC             pcmg,pc;
41*c4762a1bSJed Brown   PetscBool      flg;
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   ierr = PetscInitialize(&Argc,&Args,(char*)0,help);if (ierr) return ierr;
44*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-x",&x_mesh,NULL);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-l",&levels,NULL);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-c",&cycles,NULL);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-smooths",&smooths,NULL);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-a",&flg);CHKERRQ(ierr);
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   if (flg) am = PC_MG_ADDITIVE;
51*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-f",&flg);CHKERRQ(ierr);
52*c4762a1bSJed Brown   if (flg) am = PC_MG_FULL;
53*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-j",&flg);CHKERRQ(ierr);
54*c4762a1bSJed Brown   if (flg) use_jacobi = 1;
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   ierr = PetscMalloc1(levels,&N);CHKERRQ(ierr);
57*c4762a1bSJed Brown   N[0] = x_mesh;
58*c4762a1bSJed Brown   for (i=1; i<levels; i++) {
59*c4762a1bSJed Brown     N[i] = N[i-1]/2;
60*c4762a1bSJed Brown     if (N[i] < 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Too many levels or N is not large enough");
61*c4762a1bSJed Brown   }
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown   ierr = Create1dLaplacian(N[levels-1],&cmat);CHKERRQ(ierr);
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   ierr = KSPCreate(PETSC_COMM_WORLD,&kspmg);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = KSPGetPC(kspmg,&pcmg);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = KSPSetFromOptions(kspmg);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = PCSetType(pcmg,PCMG);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = PCMGSetLevels(pcmg,levels,NULL);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = PCMGSetType(pcmg,am);CHKERRQ(ierr);
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown   ierr = PCMGGetCoarseSolve(pcmg,&cksp);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = KSPSetOperators(cksp,cmat,cmat);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = KSPGetPC(cksp,&pc);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = KSPSetType(cksp,KSPPREONLY);CHKERRQ(ierr);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   /* zero is finest level */
79*c4762a1bSJed Brown   for (i=0; i<levels-1; i++) {
80*c4762a1bSJed Brown     ierr = PCMGSetResidual(pcmg,levels - 1 - i,residual,(Mat)0);CHKERRQ(ierr);
81*c4762a1bSJed Brown     ierr = MatCreateShell(PETSC_COMM_WORLD,N[i+1],N[i],N[i+1],N[i],(void*)0,&mat[i]);CHKERRQ(ierr);
82*c4762a1bSJed Brown     ierr = MatShellSetOperation(mat[i],MATOP_MULT,(void (*)(void))restrct);CHKERRQ(ierr);
83*c4762a1bSJed Brown     ierr = MatShellSetOperation(mat[i],MATOP_MULT_TRANSPOSE_ADD,(void (*)(void))interpolate);CHKERRQ(ierr);
84*c4762a1bSJed Brown     ierr = PCMGSetInterpolation(pcmg,levels - 1 - i,mat[i]);CHKERRQ(ierr);
85*c4762a1bSJed Brown     ierr = PCMGSetRestriction(pcmg,levels - 1 - i,mat[i]);CHKERRQ(ierr);
86*c4762a1bSJed Brown     ierr = PCMGSetCycleTypeOnLevel(pcmg,levels - 1 - i,(PCMGCycleType)cycles);CHKERRQ(ierr);
87*c4762a1bSJed Brown 
88*c4762a1bSJed Brown     /* set smoother */
89*c4762a1bSJed Brown     ierr = PCMGGetSmoother(pcmg,levels - 1 - i,&ksp[i]);CHKERRQ(ierr);
90*c4762a1bSJed Brown     ierr = KSPGetPC(ksp[i],&pc);CHKERRQ(ierr);
91*c4762a1bSJed Brown     ierr = PCSetType(pc,PCSHELL);CHKERRQ(ierr);
92*c4762a1bSJed Brown     ierr = PCShellSetName(pc,"user_precond");CHKERRQ(ierr);
93*c4762a1bSJed Brown     ierr = PCShellGetName(pc,&shellname);CHKERRQ(ierr);
94*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"level=%D, PCShell name is %s\n",i,shellname);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown     /* this is a dummy! since KSP requires a matrix passed in  */
97*c4762a1bSJed Brown     ierr = KSPSetOperators(ksp[i],mat[i],mat[i]);CHKERRQ(ierr);
98*c4762a1bSJed Brown     /*
99*c4762a1bSJed Brown         We override the matrix passed in by forcing it to use Richardson with
100*c4762a1bSJed Brown         a user provided application. This is non-standard and this practice
101*c4762a1bSJed Brown         should be avoided.
102*c4762a1bSJed Brown     */
103*c4762a1bSJed Brown     ierr = PCShellSetApplyRichardson(pc,gauss_seidel);CHKERRQ(ierr);
104*c4762a1bSJed Brown     if (use_jacobi) {
105*c4762a1bSJed Brown       ierr = PCShellSetApplyRichardson(pc,jacobi_smoother);CHKERRQ(ierr);
106*c4762a1bSJed Brown     }
107*c4762a1bSJed Brown     ierr = KSPSetType(ksp[i],KSPRICHARDSON);CHKERRQ(ierr);
108*c4762a1bSJed Brown     ierr = KSPSetInitialGuessNonzero(ksp[i],PETSC_TRUE);CHKERRQ(ierr);
109*c4762a1bSJed Brown     ierr = KSPSetTolerances(ksp[i],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,smooths);CHKERRQ(ierr);
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,N[i],&x);CHKERRQ(ierr);
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown     X[levels - 1 - i] = x;
114*c4762a1bSJed Brown     if (i > 0) {
115*c4762a1bSJed Brown       ierr = PCMGSetX(pcmg,levels - 1 - i,x);CHKERRQ(ierr);
116*c4762a1bSJed Brown     }
117*c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,N[i],&x);CHKERRQ(ierr);
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown     B[levels -1 - i] = x;
120*c4762a1bSJed Brown     if (i > 0) {
121*c4762a1bSJed Brown       ierr = PCMGSetRhs(pcmg,levels - 1 - i,x);CHKERRQ(ierr);
122*c4762a1bSJed Brown     }
123*c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,N[i],&x);CHKERRQ(ierr);
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown     R[levels - 1 - i] = x;
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown     ierr = PCMGSetR(pcmg,levels - 1 - i,x);CHKERRQ(ierr);
128*c4762a1bSJed Brown   }
129*c4762a1bSJed Brown   /* create coarse level vectors */
130*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = PCMGSetX(pcmg,0,x);CHKERRQ(ierr); X[0] = x;
132*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = PCMGSetRhs(pcmg,0,x);CHKERRQ(ierr); B[0] = x;
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown   /* create matrix multiply for finest level */
136*c4762a1bSJed Brown   ierr = MatCreateShell(PETSC_COMM_WORLD,N[0],N[0],N[0],N[0],(void*)0,&fmat);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = MatShellSetOperation(fmat,MATOP_MULT,(void (*)(void))amult);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = KSPSetOperators(kspmg,fmat,fmat);CHKERRQ(ierr);
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown   ierr = CalculateSolution(N[0],&solution);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = CalculateRhs(B[levels-1]);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = VecSet(X[levels-1],0.0);CHKERRQ(ierr);
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown   ierr = residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = CalculateError(solution,X[levels-1],R[levels-1],e);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"l_2 error %g max error %g resi %g\n",(double)e[0],(double)e[1],(double)e[2]);CHKERRQ(ierr);
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown   ierr = KSPSolve(kspmg,B[levels-1],X[levels-1]);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = KSPGetIterationNumber(kspmg,&its);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = CalculateError(solution,X[levels-1],R[levels-1],e);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"its %D l_2 error %g max error %g resi %g\n",its,(double)e[0],(double)e[1],(double)e[2]);CHKERRQ(ierr);
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   ierr = PetscFree(N);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr = VecDestroy(&solution);CHKERRQ(ierr);
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown   /* note we have to keep a list of all vectors allocated, this is
158*c4762a1bSJed Brown      not ideal, but putting it in MGDestroy is not so good either*/
159*c4762a1bSJed Brown   for (i=0; i<levels; i++) {
160*c4762a1bSJed Brown     ierr = VecDestroy(&X[i]);CHKERRQ(ierr);
161*c4762a1bSJed Brown     ierr = VecDestroy(&B[i]);CHKERRQ(ierr);
162*c4762a1bSJed Brown     if (i) {ierr = VecDestroy(&R[i]);CHKERRQ(ierr);}
163*c4762a1bSJed Brown   }
164*c4762a1bSJed Brown   for (i=0; i<levels-1; i++) {
165*c4762a1bSJed Brown     ierr = MatDestroy(&mat[i]);CHKERRQ(ierr);
166*c4762a1bSJed Brown   }
167*c4762a1bSJed Brown   ierr = MatDestroy(&cmat);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = MatDestroy(&fmat);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = KSPDestroy(&kspmg);CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = PetscFinalize();
171*c4762a1bSJed Brown   return ierr;
172*c4762a1bSJed Brown }
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
175*c4762a1bSJed Brown PetscErrorCode residual(Mat mat,Vec bb,Vec xx,Vec rr)
176*c4762a1bSJed Brown {
177*c4762a1bSJed Brown   PetscInt          i,n1;
178*c4762a1bSJed Brown   PetscErrorCode    ierr;
179*c4762a1bSJed Brown   PetscScalar       *x,*r;
180*c4762a1bSJed Brown   const PetscScalar *b;
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   PetscFunctionBegin;
183*c4762a1bSJed Brown   ierr = VecGetSize(bb,&n1);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
187*c4762a1bSJed Brown   n1--;
188*c4762a1bSJed Brown   r[0]  = b[0] + x[1] - 2.0*x[0];
189*c4762a1bSJed Brown   r[n1] = b[n1] + x[n1-1] - 2.0*x[n1];
190*c4762a1bSJed Brown   for (i=1; i<n1; i++) r[i] = b[i] + x[i+1] + x[i-1] - 2.0*x[i];
191*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
192*c4762a1bSJed Brown   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
193*c4762a1bSJed Brown   ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
194*c4762a1bSJed Brown   PetscFunctionReturn(0);
195*c4762a1bSJed Brown }
196*c4762a1bSJed Brown PetscErrorCode amult(Mat mat,Vec xx,Vec yy)
197*c4762a1bSJed Brown {
198*c4762a1bSJed Brown   PetscInt          i,n1;
199*c4762a1bSJed Brown   PetscErrorCode    ierr;
200*c4762a1bSJed Brown   PetscScalar       *y;
201*c4762a1bSJed Brown   const PetscScalar *x;
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown   PetscFunctionBegin;
204*c4762a1bSJed Brown   ierr = VecGetSize(xx,&n1);CHKERRQ(ierr);
205*c4762a1bSJed Brown   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
206*c4762a1bSJed Brown   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
207*c4762a1bSJed Brown   n1--;
208*c4762a1bSJed Brown   y[0] =  -x[1] + 2.0*x[0];
209*c4762a1bSJed Brown   y[n1] = -x[n1-1] + 2.0*x[n1];
210*c4762a1bSJed Brown   for (i=1; i<n1; i++) y[i] = -x[i+1] - x[i-1] + 2.0*x[i];
211*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
212*c4762a1bSJed Brown   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
213*c4762a1bSJed Brown   PetscFunctionReturn(0);
214*c4762a1bSJed Brown }
215*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
216*c4762a1bSJed Brown PetscErrorCode gauss_seidel(PC pc,Vec bb,Vec xx,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt m,PetscBool guesszero,PetscInt *its,PCRichardsonConvergedReason *reason)
217*c4762a1bSJed Brown {
218*c4762a1bSJed Brown   PetscInt          i,n1;
219*c4762a1bSJed Brown   PetscErrorCode    ierr;
220*c4762a1bSJed Brown   PetscScalar       *x;
221*c4762a1bSJed Brown   const PetscScalar *b;
222*c4762a1bSJed Brown 
223*c4762a1bSJed Brown   PetscFunctionBegin;
224*c4762a1bSJed Brown   *its    = m;
225*c4762a1bSJed Brown   *reason = PCRICHARDSON_CONVERGED_ITS;
226*c4762a1bSJed Brown   ierr    = VecGetSize(bb,&n1);CHKERRQ(ierr); n1--;
227*c4762a1bSJed Brown   ierr    = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr    = VecGetArray(xx,&x);CHKERRQ(ierr);
229*c4762a1bSJed Brown   while (m--) {
230*c4762a1bSJed Brown     x[0] =  .5*(x[1] + b[0]);
231*c4762a1bSJed Brown     for (i=1; i<n1; i++) x[i] = .5*(x[i+1] + x[i-1] + b[i]);
232*c4762a1bSJed Brown     x[n1] = .5*(x[n1-1] + b[n1]);
233*c4762a1bSJed Brown     for (i=n1-1; i>0; i--) x[i] = .5*(x[i+1] + x[i-1] + b[i]);
234*c4762a1bSJed Brown     x[0] =  .5*(x[1] + b[0]);
235*c4762a1bSJed Brown   }
236*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
237*c4762a1bSJed Brown   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
238*c4762a1bSJed Brown   PetscFunctionReturn(0);
239*c4762a1bSJed Brown }
240*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
241*c4762a1bSJed Brown PetscErrorCode jacobi_smoother(PC pc,Vec bb,Vec xx,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt m,PetscBool guesszero,PetscInt *its,PCRichardsonConvergedReason *reason)
242*c4762a1bSJed Brown {
243*c4762a1bSJed Brown   PetscInt          i,n,n1;
244*c4762a1bSJed Brown   PetscErrorCode    ierr;
245*c4762a1bSJed Brown   PetscScalar       *r,*x;
246*c4762a1bSJed Brown   const PetscScalar *b;
247*c4762a1bSJed Brown 
248*c4762a1bSJed Brown   PetscFunctionBegin;
249*c4762a1bSJed Brown   *its    = m;
250*c4762a1bSJed Brown   *reason = PCRICHARDSON_CONVERGED_ITS;
251*c4762a1bSJed Brown   ierr    = VecGetSize(bb,&n);CHKERRQ(ierr); n1 = n - 1;
252*c4762a1bSJed Brown   ierr    = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
253*c4762a1bSJed Brown   ierr    = VecGetArray(xx,&x);CHKERRQ(ierr);
254*c4762a1bSJed Brown   ierr    = VecGetArray(w,&r);CHKERRQ(ierr);
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown   while (m--) {
257*c4762a1bSJed Brown     r[0] = .5*(x[1] + b[0]);
258*c4762a1bSJed Brown     for (i=1; i<n1; i++) r[i] = .5*(x[i+1] + x[i-1] + b[i]);
259*c4762a1bSJed Brown     r[n1] = .5*(x[n1-1] + b[n1]);
260*c4762a1bSJed Brown     for (i=0; i<n; i++) x[i] = (2.0*r[i] + x[i])/3.0;
261*c4762a1bSJed Brown   }
262*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
263*c4762a1bSJed Brown   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
264*c4762a1bSJed Brown   ierr = VecRestoreArray(w,&r);CHKERRQ(ierr);
265*c4762a1bSJed Brown   PetscFunctionReturn(0);
266*c4762a1bSJed Brown }
267*c4762a1bSJed Brown /*
268*c4762a1bSJed Brown    We know for this application that yy  and zz are the same
269*c4762a1bSJed Brown */
270*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
271*c4762a1bSJed Brown PetscErrorCode interpolate(Mat mat,Vec xx,Vec yy,Vec zz)
272*c4762a1bSJed Brown {
273*c4762a1bSJed Brown   PetscInt          i,n,N,i2;
274*c4762a1bSJed Brown   PetscErrorCode    ierr;
275*c4762a1bSJed Brown   PetscScalar       *y;
276*c4762a1bSJed Brown   const PetscScalar *x;
277*c4762a1bSJed Brown 
278*c4762a1bSJed Brown   PetscFunctionBegin;
279*c4762a1bSJed Brown   ierr = VecGetSize(yy,&N);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
282*c4762a1bSJed Brown   n    = N/2;
283*c4762a1bSJed Brown   for (i=0; i<n; i++) {
284*c4762a1bSJed Brown     i2       = 2*i;
285*c4762a1bSJed Brown     y[i2]   += .5*x[i];
286*c4762a1bSJed Brown     y[i2+1] +=    x[i];
287*c4762a1bSJed Brown     y[i2+2] += .5*x[i];
288*c4762a1bSJed Brown   }
289*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
290*c4762a1bSJed Brown   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
291*c4762a1bSJed Brown   PetscFunctionReturn(0);
292*c4762a1bSJed Brown }
293*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
294*c4762a1bSJed Brown PetscErrorCode restrct(Mat mat,Vec rr,Vec bb)
295*c4762a1bSJed Brown {
296*c4762a1bSJed Brown   PetscInt          i,n,N,i2;
297*c4762a1bSJed Brown   PetscErrorCode    ierr;
298*c4762a1bSJed Brown   PetscScalar       *b;
299*c4762a1bSJed Brown   const PetscScalar *r;
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown   PetscFunctionBegin;
302*c4762a1bSJed Brown   ierr = VecGetSize(rr,&N);CHKERRQ(ierr);
303*c4762a1bSJed Brown   ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr);
304*c4762a1bSJed Brown   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
305*c4762a1bSJed Brown   n    = N/2;
306*c4762a1bSJed Brown 
307*c4762a1bSJed Brown   for (i=0; i<n; i++) {
308*c4762a1bSJed Brown     i2   = 2*i;
309*c4762a1bSJed Brown     b[i] = (r[i2] + 2.0*r[i2+1] + r[i2+2]);
310*c4762a1bSJed Brown   }
311*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr);
312*c4762a1bSJed Brown   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
313*c4762a1bSJed Brown   PetscFunctionReturn(0);
314*c4762a1bSJed Brown }
315*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
316*c4762a1bSJed Brown PetscErrorCode Create1dLaplacian(PetscInt n,Mat *mat)
317*c4762a1bSJed Brown {
318*c4762a1bSJed Brown   PetscScalar    mone = -1.0,two = 2.0;
319*c4762a1bSJed Brown   PetscInt       i,idx;
320*c4762a1bSJed Brown   PetscErrorCode ierr;
321*c4762a1bSJed Brown 
322*c4762a1bSJed Brown   PetscFunctionBegin;
323*c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,mat);CHKERRQ(ierr);
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown   idx  = n-1;
326*c4762a1bSJed Brown   ierr = MatSetValues(*mat,1,&idx,1,&idx,&two,INSERT_VALUES);CHKERRQ(ierr);
327*c4762a1bSJed Brown   for (i=0; i<n-1; i++) {
328*c4762a1bSJed Brown     ierr = MatSetValues(*mat,1,&i,1,&i,&two,INSERT_VALUES);CHKERRQ(ierr);
329*c4762a1bSJed Brown     idx  = i+1;
330*c4762a1bSJed Brown     ierr = MatSetValues(*mat,1,&idx,1,&i,&mone,INSERT_VALUES);CHKERRQ(ierr);
331*c4762a1bSJed Brown     ierr = MatSetValues(*mat,1,&i,1,&idx,&mone,INSERT_VALUES);CHKERRQ(ierr);
332*c4762a1bSJed Brown   }
333*c4762a1bSJed Brown   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
334*c4762a1bSJed Brown   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
335*c4762a1bSJed Brown   PetscFunctionReturn(0);
336*c4762a1bSJed Brown }
337*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
338*c4762a1bSJed Brown PetscErrorCode CalculateRhs(Vec u)
339*c4762a1bSJed Brown {
340*c4762a1bSJed Brown   PetscErrorCode ierr;
341*c4762a1bSJed Brown   PetscInt       i,n;
342*c4762a1bSJed Brown   PetscReal      h,x = 0.0;
343*c4762a1bSJed Brown   PetscScalar    uu;
344*c4762a1bSJed Brown 
345*c4762a1bSJed Brown   PetscFunctionBegin;
346*c4762a1bSJed Brown   ierr = VecGetSize(u,&n);CHKERRQ(ierr);
347*c4762a1bSJed Brown   h    = 1.0/((PetscReal)(n+1));
348*c4762a1bSJed Brown   for (i=0; i<n; i++) {
349*c4762a1bSJed Brown     x   += h; uu = 2.0*h*h;
350*c4762a1bSJed Brown     ierr = VecSetValues(u,1,&i,&uu,INSERT_VALUES);CHKERRQ(ierr);
351*c4762a1bSJed Brown   }
352*c4762a1bSJed Brown   PetscFunctionReturn(0);
353*c4762a1bSJed Brown }
354*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
355*c4762a1bSJed Brown PetscErrorCode CalculateSolution(PetscInt n,Vec *solution)
356*c4762a1bSJed Brown {
357*c4762a1bSJed Brown   PetscErrorCode ierr;
358*c4762a1bSJed Brown   PetscInt       i;
359*c4762a1bSJed Brown   PetscReal      h,x = 0.0;
360*c4762a1bSJed Brown   PetscScalar    uu;
361*c4762a1bSJed Brown 
362*c4762a1bSJed Brown   PetscFunctionBegin;
363*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,n,solution);CHKERRQ(ierr);
364*c4762a1bSJed Brown   h    = 1.0/((PetscReal)(n+1));
365*c4762a1bSJed Brown   for (i=0; i<n; i++) {
366*c4762a1bSJed Brown     x   += h; uu = x*(1.-x);
367*c4762a1bSJed Brown     ierr = VecSetValues(*solution,1,&i,&uu,INSERT_VALUES);CHKERRQ(ierr);
368*c4762a1bSJed Brown   }
369*c4762a1bSJed Brown   PetscFunctionReturn(0);
370*c4762a1bSJed Brown }
371*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
372*c4762a1bSJed Brown PetscErrorCode CalculateError(Vec solution,Vec u,Vec r,PetscReal *e)
373*c4762a1bSJed Brown {
374*c4762a1bSJed Brown   PetscErrorCode ierr;
375*c4762a1bSJed Brown 
376*c4762a1bSJed Brown   PetscFunctionBegin;
377*c4762a1bSJed Brown   ierr = VecNorm(r,NORM_2,e+2);CHKERRQ(ierr);
378*c4762a1bSJed Brown   ierr = VecWAXPY(r,-1.0,u,solution);CHKERRQ(ierr);
379*c4762a1bSJed Brown   ierr = VecNorm(r,NORM_2,e);CHKERRQ(ierr);
380*c4762a1bSJed Brown   ierr = VecNorm(r,NORM_1,e+1);CHKERRQ(ierr);
381*c4762a1bSJed Brown   PetscFunctionReturn(0);
382*c4762a1bSJed Brown }
383*c4762a1bSJed Brown 
384*c4762a1bSJed Brown 
385*c4762a1bSJed Brown 
386*c4762a1bSJed Brown 
387*c4762a1bSJed Brown /*TEST
388*c4762a1bSJed Brown 
389*c4762a1bSJed Brown    test:
390*c4762a1bSJed Brown 
391*c4762a1bSJed Brown TEST*/
392