xref: /petsc/src/ksp/pc/tests/ex5.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests the multigrid code.  The input parameters are:\n\
3c4762a1bSJed Brown   -x N              Use a mesh in the x direction of N.  \n\
4c4762a1bSJed Brown   -c N              Use N V-cycles.  \n\
5c4762a1bSJed Brown   -l N              Use N Levels.  \n\
6c4762a1bSJed Brown   -smooths N        Use N pre smooths and N post smooths.  \n\
7c4762a1bSJed Brown   -j                Use Jacobi smoother.  \n\
8c4762a1bSJed Brown   -a use additive multigrid \n\
9c4762a1bSJed Brown   -f use full multigrid (preconditioner variant) \n\
10c4762a1bSJed Brown This example also demonstrates matrix-free methods\n\n";
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown   This is not a good example to understand the use of multigrid with PETSc.
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16c4762a1bSJed Brown #include <petscksp.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown PetscErrorCode  residual(Mat,Vec,Vec,Vec);
19c4762a1bSJed Brown PetscErrorCode  gauss_seidel(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
20c4762a1bSJed Brown PetscErrorCode  jacobi_smoother(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
21c4762a1bSJed Brown PetscErrorCode  interpolate(Mat,Vec,Vec,Vec);
22c4762a1bSJed Brown PetscErrorCode  restrct(Mat,Vec,Vec);
23c4762a1bSJed Brown PetscErrorCode  Create1dLaplacian(PetscInt,Mat*);
24c4762a1bSJed Brown PetscErrorCode  CalculateRhs(Vec);
25c4762a1bSJed Brown PetscErrorCode  CalculateError(Vec,Vec,Vec,PetscReal*);
26c4762a1bSJed Brown PetscErrorCode  CalculateSolution(PetscInt,Vec*);
27c4762a1bSJed Brown PetscErrorCode  amult(Mat,Vec,Vec);
28f2fddbb2SStefano Zampini PetscErrorCode  apply_pc(PC,Vec,Vec);
29c4762a1bSJed Brown 
30c4762a1bSJed Brown int main(int Argc,char **Args)
31c4762a1bSJed Brown {
32c4762a1bSJed Brown   PetscInt       x_mesh = 15,levels = 3,cycles = 1,use_jacobi = 0;
33c4762a1bSJed Brown   PetscInt       i,smooths = 1,*N,its;
34c4762a1bSJed Brown   PCMGType       am = PC_MG_MULTIPLICATIVE;
35c4762a1bSJed Brown   Mat            cmat,mat[20],fmat;
36c4762a1bSJed Brown   KSP            cksp,ksp[20],kspmg;
37c4762a1bSJed Brown   PetscReal      e[3];  /* l_2 error,max error, residual */
38c4762a1bSJed Brown   const char     *shellname;
39c4762a1bSJed Brown   Vec            x,solution,X[20],R[20],B[20];
40c4762a1bSJed Brown   PC             pcmg,pc;
41c4762a1bSJed Brown   PetscBool      flg;
42c4762a1bSJed Brown 
43*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&Argc,&Args,(char*)0,help));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-x",&x_mesh,NULL));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-l",&levels,NULL));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-c",&cycles,NULL));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-smooths",&smooths,NULL));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-a",&flg));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   if (flg) am = PC_MG_ADDITIVE;
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-f",&flg));
52c4762a1bSJed Brown   if (flg) am = PC_MG_FULL;
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-j",&flg));
54c4762a1bSJed Brown   if (flg) use_jacobi = 1;
55c4762a1bSJed Brown 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(levels,&N));
57c4762a1bSJed Brown   N[0] = x_mesh;
58c4762a1bSJed Brown   for (i=1; i<levels; i++) {
59c4762a1bSJed Brown     N[i] = N[i-1]/2;
602c71b3e2SJacob Faibussowitsch     PetscCheckFalse(N[i] < 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"Too many levels or N is not large enough");
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown 
635f80ce2aSJacob Faibussowitsch   CHKERRQ(Create1dLaplacian(N[levels-1],&cmat));
64c4762a1bSJed Brown 
655f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&kspmg));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(kspmg,&pcmg));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(kspmg));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pcmg,PCMG));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PCMGSetLevels(pcmg,levels,NULL));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PCMGSetType(pcmg,am));
71c4762a1bSJed Brown 
725f80ce2aSJacob Faibussowitsch   CHKERRQ(PCMGGetCoarseSolve(pcmg,&cksp));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOperators(cksp,cmat,cmat));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(cksp,&pc));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCLU));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetType(cksp,KSPPREONLY));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* zero is finest level */
79c4762a1bSJed Brown   for (i=0; i<levels-1; i++) {
80f2fddbb2SStefano Zampini     Mat dummy;
81f2fddbb2SStefano Zampini 
825f80ce2aSJacob Faibussowitsch     CHKERRQ(PCMGSetResidual(pcmg,levels - 1 - i,residual,NULL));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,N[i+1],N[i],N[i+1],N[i],NULL,&mat[i]));
845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShellSetOperation(mat[i],MATOP_MULT,(void (*)(void))restrct));
855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShellSetOperation(mat[i],MATOP_MULT_TRANSPOSE_ADD,(void (*)(void))interpolate));
865f80ce2aSJacob Faibussowitsch     CHKERRQ(PCMGSetInterpolation(pcmg,levels - 1 - i,mat[i]));
875f80ce2aSJacob Faibussowitsch     CHKERRQ(PCMGSetRestriction(pcmg,levels - 1 - i,mat[i]));
885f80ce2aSJacob Faibussowitsch     CHKERRQ(PCMGSetCycleTypeOnLevel(pcmg,levels - 1 - i,(PCMGCycleType)cycles));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown     /* set smoother */
915f80ce2aSJacob Faibussowitsch     CHKERRQ(PCMGGetSmoother(pcmg,levels - 1 - i,&ksp[i]));
925f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp[i],&pc));
935f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCSHELL));
945f80ce2aSJacob Faibussowitsch     CHKERRQ(PCShellSetName(pc,"user_precond"));
955f80ce2aSJacob Faibussowitsch     CHKERRQ(PCShellGetName(pc,&shellname));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"level=%D, PCShell name is %s\n",i,shellname));
97c4762a1bSJed Brown 
98f2fddbb2SStefano Zampini     /* this is not used unless different options are passed to the solver */
995f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,N[i],N[i],N[i],N[i],NULL,&dummy));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShellSetOperation(dummy,MATOP_MULT,(void (*)(void))amult));
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetOperators(ksp[i],dummy,dummy));
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&dummy));
103f2fddbb2SStefano Zampini 
104c4762a1bSJed Brown     /*
105c4762a1bSJed Brown         We override the matrix passed in by forcing it to use Richardson with
106c4762a1bSJed Brown         a user provided application. This is non-standard and this practice
107c4762a1bSJed Brown         should be avoided.
108c4762a1bSJed Brown     */
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(PCShellSetApply(pc,apply_pc));
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(PCShellSetApplyRichardson(pc,gauss_seidel));
111c4762a1bSJed Brown     if (use_jacobi) {
1125f80ce2aSJacob Faibussowitsch       CHKERRQ(PCShellSetApplyRichardson(pc,jacobi_smoother));
113c4762a1bSJed Brown     }
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetType(ksp[i],KSPRICHARDSON));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetInitialGuessNonzero(ksp[i],PETSC_TRUE));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(ksp[i],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,smooths));
117c4762a1bSJed Brown 
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[i],&x));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown     X[levels - 1 - i] = x;
121c4762a1bSJed Brown     if (i > 0) {
1225f80ce2aSJacob Faibussowitsch       CHKERRQ(PCMGSetX(pcmg,levels - 1 - i,x));
123c4762a1bSJed Brown     }
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[i],&x));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown     B[levels -1 - i] = x;
127c4762a1bSJed Brown     if (i > 0) {
1285f80ce2aSJacob Faibussowitsch       CHKERRQ(PCMGSetRhs(pcmg,levels - 1 - i,x));
129c4762a1bSJed Brown     }
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[i],&x));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown     R[levels - 1 - i] = x;
133c4762a1bSJed Brown 
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(PCMGSetR(pcmg,levels - 1 - i,x));
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown   /* create coarse level vectors */
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(PCMGSetX(pcmg,0,x)); X[0] = x;
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x));
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(PCMGSetRhs(pcmg,0,x)); B[0] = x;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* create matrix multiply for finest level */
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,N[0],N[0],N[0],N[0],NULL,&fmat));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(fmat,MATOP_MULT,(void (*)(void))amult));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOperators(kspmg,fmat,fmat));
146c4762a1bSJed Brown 
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(CalculateSolution(N[0],&solution));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(CalculateRhs(B[levels-1]));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(X[levels-1],0.0));
150c4762a1bSJed Brown 
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(CalculateError(solution,X[levels-1],R[levels-1],e));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"l_2 error %g max error %g resi %g\n",(double)e[0],(double)e[1],(double)e[2]));
154c4762a1bSJed Brown 
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolve(kspmg,B[levels-1],X[levels-1]));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(kspmg,&its));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(CalculateError(solution,X[levels-1],R[levels-1],e));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(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]));
160c4762a1bSJed Brown 
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(N));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&solution));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* note we have to keep a list of all vectors allocated, this is
165c4762a1bSJed Brown      not ideal, but putting it in MGDestroy is not so good either*/
166c4762a1bSJed Brown   for (i=0; i<levels; i++) {
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&X[i]));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&B[i]));
1695f80ce2aSJacob Faibussowitsch     if (i) CHKERRQ(VecDestroy(&R[i]));
170c4762a1bSJed Brown   }
171c4762a1bSJed Brown   for (i=0; i<levels-1; i++) {
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&mat[i]));
173c4762a1bSJed Brown   }
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&cmat));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&fmat));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&kspmg));
177*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
178*b122ec5aSJacob Faibussowitsch   return 0;
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown /* --------------------------------------------------------------------- */
182c4762a1bSJed Brown PetscErrorCode residual(Mat mat,Vec bb,Vec xx,Vec rr)
183c4762a1bSJed Brown {
184c4762a1bSJed Brown   PetscInt          i,n1;
185c4762a1bSJed Brown   PetscScalar       *x,*r;
186c4762a1bSJed Brown   const PetscScalar *b;
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   PetscFunctionBegin;
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(bb,&n1));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(bb,&b));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(xx,&x));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(rr,&r));
193c4762a1bSJed Brown   n1--;
194c4762a1bSJed Brown   r[0]  = b[0] + x[1] - 2.0*x[0];
195c4762a1bSJed Brown   r[n1] = b[n1] + x[n1-1] - 2.0*x[n1];
196c4762a1bSJed Brown   for (i=1; i<n1; i++) r[i] = b[i] + x[i+1] + x[i-1] - 2.0*x[i];
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(bb,&b));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(xx,&x));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(rr,&r));
200c4762a1bSJed Brown   PetscFunctionReturn(0);
201c4762a1bSJed Brown }
202f2fddbb2SStefano Zampini 
203c4762a1bSJed Brown PetscErrorCode amult(Mat mat,Vec xx,Vec yy)
204c4762a1bSJed Brown {
205c4762a1bSJed Brown   PetscInt          i,n1;
206c4762a1bSJed Brown   PetscScalar       *y;
207c4762a1bSJed Brown   const PetscScalar *x;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   PetscFunctionBegin;
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(xx,&n1));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(xx,&x));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(yy,&y));
213c4762a1bSJed Brown   n1--;
214c4762a1bSJed Brown   y[0] =  -x[1] + 2.0*x[0];
215c4762a1bSJed Brown   y[n1] = -x[n1-1] + 2.0*x[n1];
216c4762a1bSJed Brown   for (i=1; i<n1; i++) y[i] = -x[i+1] - x[i-1] + 2.0*x[i];
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(xx,&x));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(yy,&y));
219c4762a1bSJed Brown   PetscFunctionReturn(0);
220c4762a1bSJed Brown }
221f2fddbb2SStefano Zampini 
222c4762a1bSJed Brown /* --------------------------------------------------------------------- */
223f2fddbb2SStefano Zampini PetscErrorCode apply_pc(PC pc,Vec bb,Vec xx)
224f2fddbb2SStefano Zampini {
225f2fddbb2SStefano Zampini   PetscFunctionBegin;
226f2fddbb2SStefano Zampini   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not implemented");
227f2fddbb2SStefano Zampini }
228f2fddbb2SStefano Zampini 
229c4762a1bSJed 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)
230c4762a1bSJed Brown {
231c4762a1bSJed Brown   PetscInt          i,n1;
232c4762a1bSJed Brown   PetscScalar       *x;
233c4762a1bSJed Brown   const PetscScalar *b;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   PetscFunctionBegin;
236c4762a1bSJed Brown   *its    = m;
237c4762a1bSJed Brown   *reason = PCRICHARDSON_CONVERGED_ITS;
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(bb,&n1)); n1--;
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(bb,&b));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(xx,&x));
241c4762a1bSJed Brown   while (m--) {
242c4762a1bSJed Brown     x[0] =  .5*(x[1] + b[0]);
243c4762a1bSJed Brown     for (i=1; i<n1; i++) x[i] = .5*(x[i+1] + x[i-1] + b[i]);
244c4762a1bSJed Brown     x[n1] = .5*(x[n1-1] + b[n1]);
245c4762a1bSJed Brown     for (i=n1-1; i>0; i--) x[i] = .5*(x[i+1] + x[i-1] + b[i]);
246c4762a1bSJed Brown     x[0] =  .5*(x[1] + b[0]);
247c4762a1bSJed Brown   }
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(bb,&b));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(xx,&x));
250c4762a1bSJed Brown   PetscFunctionReturn(0);
251c4762a1bSJed Brown }
252c4762a1bSJed Brown /* --------------------------------------------------------------------- */
253c4762a1bSJed 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)
254c4762a1bSJed Brown {
255c4762a1bSJed Brown   PetscInt          i,n,n1;
256c4762a1bSJed Brown   PetscScalar       *r,*x;
257c4762a1bSJed Brown   const PetscScalar *b;
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   PetscFunctionBegin;
260c4762a1bSJed Brown   *its    = m;
261c4762a1bSJed Brown   *reason = PCRICHARDSON_CONVERGED_ITS;
2625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(bb,&n)); n1 = n - 1;
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(bb,&b));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(xx,&x));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(w,&r));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   while (m--) {
268c4762a1bSJed Brown     r[0] = .5*(x[1] + b[0]);
269c4762a1bSJed Brown     for (i=1; i<n1; i++) r[i] = .5*(x[i+1] + x[i-1] + b[i]);
270c4762a1bSJed Brown     r[n1] = .5*(x[n1-1] + b[n1]);
271c4762a1bSJed Brown     for (i=0; i<n; i++) x[i] = (2.0*r[i] + x[i])/3.0;
272c4762a1bSJed Brown   }
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(bb,&b));
2745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(xx,&x));
2755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(w,&r));
276c4762a1bSJed Brown   PetscFunctionReturn(0);
277c4762a1bSJed Brown }
278c4762a1bSJed Brown /*
279c4762a1bSJed Brown    We know for this application that yy  and zz are the same
280c4762a1bSJed Brown */
281c4762a1bSJed Brown /* --------------------------------------------------------------------- */
282c4762a1bSJed Brown PetscErrorCode interpolate(Mat mat,Vec xx,Vec yy,Vec zz)
283c4762a1bSJed Brown {
284c4762a1bSJed Brown   PetscInt          i,n,N,i2;
285c4762a1bSJed Brown   PetscScalar       *y;
286c4762a1bSJed Brown   const PetscScalar *x;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   PetscFunctionBegin;
2895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(yy,&N));
2905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(xx,&x));
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(yy,&y));
292c4762a1bSJed Brown   n    = N/2;
293c4762a1bSJed Brown   for (i=0; i<n; i++) {
294c4762a1bSJed Brown     i2       = 2*i;
295c4762a1bSJed Brown     y[i2]   += .5*x[i];
296c4762a1bSJed Brown     y[i2+1] +=    x[i];
297c4762a1bSJed Brown     y[i2+2] += .5*x[i];
298c4762a1bSJed Brown   }
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(xx,&x));
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(yy,&y));
301c4762a1bSJed Brown   PetscFunctionReturn(0);
302c4762a1bSJed Brown }
303c4762a1bSJed Brown /* --------------------------------------------------------------------- */
304c4762a1bSJed Brown PetscErrorCode restrct(Mat mat,Vec rr,Vec bb)
305c4762a1bSJed Brown {
306c4762a1bSJed Brown   PetscInt          i,n,N,i2;
307c4762a1bSJed Brown   PetscScalar       *b;
308c4762a1bSJed Brown   const PetscScalar *r;
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   PetscFunctionBegin;
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(rr,&N));
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(rr,&r));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(bb,&b));
314c4762a1bSJed Brown   n    = N/2;
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   for (i=0; i<n; i++) {
317c4762a1bSJed Brown     i2   = 2*i;
318c4762a1bSJed Brown     b[i] = (r[i2] + 2.0*r[i2+1] + r[i2+2]);
319c4762a1bSJed Brown   }
3205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(rr,&r));
3215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(bb,&b));
322c4762a1bSJed Brown   PetscFunctionReturn(0);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown /* --------------------------------------------------------------------- */
325c4762a1bSJed Brown PetscErrorCode Create1dLaplacian(PetscInt n,Mat *mat)
326c4762a1bSJed Brown {
327c4762a1bSJed Brown   PetscScalar    mone = -1.0,two = 2.0;
328c4762a1bSJed Brown   PetscInt       i,idx;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   PetscFunctionBegin;
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,mat));
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   idx  = n-1;
3345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(*mat,1,&idx,1,&idx,&two,INSERT_VALUES));
335c4762a1bSJed Brown   for (i=0; i<n-1; i++) {
3365f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(*mat,1,&i,1,&i,&two,INSERT_VALUES));
337c4762a1bSJed Brown     idx  = i+1;
3385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(*mat,1,&idx,1,&i,&mone,INSERT_VALUES));
3395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(*mat,1,&i,1,&idx,&mone,INSERT_VALUES));
340c4762a1bSJed Brown   }
3415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY));
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY));
343c4762a1bSJed Brown   PetscFunctionReturn(0);
344c4762a1bSJed Brown }
345c4762a1bSJed Brown /* --------------------------------------------------------------------- */
346c4762a1bSJed Brown PetscErrorCode CalculateRhs(Vec u)
347c4762a1bSJed Brown {
348c4762a1bSJed Brown   PetscInt    i,n;
3495f80ce2aSJacob Faibussowitsch   PetscReal   h;
350c4762a1bSJed Brown   PetscScalar uu;
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   PetscFunctionBegin;
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(u,&n));
354c4762a1bSJed Brown   h    = 1.0/((PetscReal)(n+1));
355c4762a1bSJed Brown   for (i=0; i<n; i++) {
3565f80ce2aSJacob Faibussowitsch     uu = 2.0*h*h;
3575f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(u,1,&i,&uu,INSERT_VALUES));
358c4762a1bSJed Brown   }
359c4762a1bSJed Brown   PetscFunctionReturn(0);
360c4762a1bSJed Brown }
361c4762a1bSJed Brown /* --------------------------------------------------------------------- */
362c4762a1bSJed Brown PetscErrorCode CalculateSolution(PetscInt n,Vec *solution)
363c4762a1bSJed Brown {
364c4762a1bSJed Brown   PetscInt       i;
365c4762a1bSJed Brown   PetscReal      h,x = 0.0;
366c4762a1bSJed Brown   PetscScalar    uu;
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   PetscFunctionBegin;
3695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,solution));
370c4762a1bSJed Brown   h    = 1.0/((PetscReal)(n+1));
371c4762a1bSJed Brown   for (i=0; i<n; i++) {
372c4762a1bSJed Brown     x   += h; uu = x*(1.-x);
3735f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(*solution,1,&i,&uu,INSERT_VALUES));
374c4762a1bSJed Brown   }
375c4762a1bSJed Brown   PetscFunctionReturn(0);
376c4762a1bSJed Brown }
377c4762a1bSJed Brown /* --------------------------------------------------------------------- */
378c4762a1bSJed Brown PetscErrorCode CalculateError(Vec solution,Vec u,Vec r,PetscReal *e)
379c4762a1bSJed Brown {
380c4762a1bSJed Brown   PetscFunctionBegin;
3815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(r,NORM_2,e+2));
3825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(r,-1.0,u,solution));
3835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(r,NORM_2,e));
3845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(r,NORM_1,e+1));
385c4762a1bSJed Brown   PetscFunctionReturn(0);
386c4762a1bSJed Brown }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown /*TEST
389c4762a1bSJed Brown 
390c4762a1bSJed Brown    test:
391c4762a1bSJed Brown 
392c4762a1bSJed Brown TEST*/
393