xref: /petsc/src/ksp/pc/tests/ex5.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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 
439566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&Argc,&Args,(char*)0,help));
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-x",&x_mesh,NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&levels,NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-c",&cycles,NULL));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-smooths",&smooths,NULL));
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-a",&flg));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   if (flg) am = PC_MG_ADDITIVE;
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-f",&flg));
52c4762a1bSJed Brown   if (flg) am = PC_MG_FULL;
539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-j",&flg));
54c4762a1bSJed Brown   if (flg) use_jacobi = 1;
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(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;
607827d75bSBarry Smith     PetscCheck(N[i] >= 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"Too many levels or N is not large enough");
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(Create1dLaplacian(N[levels-1],&cmat));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD,&kspmg));
669566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(kspmg,&pcmg));
679566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(kspmg));
689566063dSJacob Faibussowitsch   PetscCall(PCSetType(pcmg,PCMG));
699566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pcmg,levels,NULL));
709566063dSJacob Faibussowitsch   PetscCall(PCMGSetType(pcmg,am));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(PCMGGetCoarseSolve(pcmg,&cksp));
739566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(cksp,cmat,cmat));
749566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(cksp,&pc));
759566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCLU));
769566063dSJacob Faibussowitsch   PetscCall(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 
829566063dSJacob Faibussowitsch     PetscCall(PCMGSetResidual(pcmg,levels - 1 - i,residual,NULL));
839566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PETSC_COMM_WORLD,N[i+1],N[i],N[i+1],N[i],NULL,&mat[i]));
849566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(mat[i],MATOP_MULT,(void (*)(void))restrct));
859566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(mat[i],MATOP_MULT_TRANSPOSE_ADD,(void (*)(void))interpolate));
869566063dSJacob Faibussowitsch     PetscCall(PCMGSetInterpolation(pcmg,levels - 1 - i,mat[i]));
879566063dSJacob Faibussowitsch     PetscCall(PCMGSetRestriction(pcmg,levels - 1 - i,mat[i]));
889566063dSJacob Faibussowitsch     PetscCall(PCMGSetCycleTypeOnLevel(pcmg,levels - 1 - i,(PCMGCycleType)cycles));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown     /* set smoother */
919566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pcmg,levels - 1 - i,&ksp[i]));
929566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp[i],&pc));
939566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc,PCSHELL));
949566063dSJacob Faibussowitsch     PetscCall(PCShellSetName(pc,"user_precond"));
959566063dSJacob Faibussowitsch     PetscCall(PCShellGetName(pc,&shellname));
9663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"level=%" PetscInt_FMT ", PCShell name is %s\n",i,shellname));
97c4762a1bSJed Brown 
98f2fddbb2SStefano Zampini     /* this is not used unless different options are passed to the solver */
999566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PETSC_COMM_WORLD,N[i],N[i],N[i],N[i],NULL,&dummy));
1009566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(dummy,MATOP_MULT,(void (*)(void))amult));
1019566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp[i],dummy,dummy));
1029566063dSJacob Faibussowitsch     PetscCall(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     */
1099566063dSJacob Faibussowitsch     PetscCall(PCShellSetApply(pc,apply_pc));
1109566063dSJacob Faibussowitsch     PetscCall(PCShellSetApplyRichardson(pc,gauss_seidel));
111*1baa6e33SBarry Smith     if (use_jacobi) PetscCall(PCShellSetApplyRichardson(pc,jacobi_smoother));
1129566063dSJacob Faibussowitsch     PetscCall(KSPSetType(ksp[i],KSPRICHARDSON));
1139566063dSJacob Faibussowitsch     PetscCall(KSPSetInitialGuessNonzero(ksp[i],PETSC_TRUE));
1149566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(ksp[i],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,smooths));
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,N[i],&x));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown     X[levels - 1 - i] = x;
119c4762a1bSJed Brown     if (i > 0) {
1209566063dSJacob Faibussowitsch       PetscCall(PCMGSetX(pcmg,levels - 1 - i,x));
121c4762a1bSJed Brown     }
1229566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,N[i],&x));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown     B[levels -1 - i] = x;
125c4762a1bSJed Brown     if (i > 0) {
1269566063dSJacob Faibussowitsch       PetscCall(PCMGSetRhs(pcmg,levels - 1 - i,x));
127c4762a1bSJed Brown     }
1289566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,N[i],&x));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown     R[levels - 1 - i] = x;
131c4762a1bSJed Brown 
1329566063dSJacob Faibussowitsch     PetscCall(PCMGSetR(pcmg,levels - 1 - i,x));
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown   /* create coarse level vectors */
1359566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x));
1369566063dSJacob Faibussowitsch   PetscCall(PCMGSetX(pcmg,0,x)); X[0] = x;
1379566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x));
1389566063dSJacob Faibussowitsch   PetscCall(PCMGSetRhs(pcmg,0,x)); B[0] = x;
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* create matrix multiply for finest level */
1419566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,N[0],N[0],N[0],N[0],NULL,&fmat));
1429566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(fmat,MATOP_MULT,(void (*)(void))amult));
1439566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(kspmg,fmat,fmat));
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(CalculateSolution(N[0],&solution));
1469566063dSJacob Faibussowitsch   PetscCall(CalculateRhs(B[levels-1]));
1479566063dSJacob Faibussowitsch   PetscCall(VecSet(X[levels-1],0.0));
148c4762a1bSJed Brown 
1499566063dSJacob Faibussowitsch   PetscCall(residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]));
1509566063dSJacob Faibussowitsch   PetscCall(CalculateError(solution,X[levels-1],R[levels-1],e));
1519566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"l_2 error %g max error %g resi %g\n",(double)e[0],(double)e[1],(double)e[2]));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(KSPSolve(kspmg,B[levels-1],X[levels-1]));
1549566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(kspmg,&its));
1559566063dSJacob Faibussowitsch   PetscCall(residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]));
1569566063dSJacob Faibussowitsch   PetscCall(CalculateError(solution,X[levels-1],R[levels-1],e));
15763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"its %" PetscInt_FMT " l_2 error %g max error %g resi %g\n",its,(double)e[0],(double)e[1],(double)e[2]));
158c4762a1bSJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(PetscFree(N));
1609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&solution));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* note we have to keep a list of all vectors allocated, this is
163c4762a1bSJed Brown      not ideal, but putting it in MGDestroy is not so good either*/
164c4762a1bSJed Brown   for (i=0; i<levels; i++) {
1659566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&X[i]));
1669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&B[i]));
1679566063dSJacob Faibussowitsch     if (i) PetscCall(VecDestroy(&R[i]));
168c4762a1bSJed Brown   }
169c4762a1bSJed Brown   for (i=0; i<levels-1; i++) {
1709566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&mat[i]));
171c4762a1bSJed Brown   }
1729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&cmat));
1739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fmat));
1749566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&kspmg));
1759566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
176b122ec5aSJacob Faibussowitsch   return 0;
177c4762a1bSJed Brown }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown /* --------------------------------------------------------------------- */
180c4762a1bSJed Brown PetscErrorCode residual(Mat mat,Vec bb,Vec xx,Vec rr)
181c4762a1bSJed Brown {
182c4762a1bSJed Brown   PetscInt          i,n1;
183c4762a1bSJed Brown   PetscScalar       *x,*r;
184c4762a1bSJed Brown   const PetscScalar *b;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(VecGetSize(bb,&n1));
1889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb,&b));
1899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx,&x));
1909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(rr,&r));
191c4762a1bSJed Brown   n1--;
192c4762a1bSJed Brown   r[0]  = b[0] + x[1] - 2.0*x[0];
193c4762a1bSJed Brown   r[n1] = b[n1] + x[n1-1] - 2.0*x[n1];
194c4762a1bSJed Brown   for (i=1; i<n1; i++) r[i] = b[i] + x[i+1] + x[i-1] - 2.0*x[i];
1959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb,&b));
1969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx,&x));
1979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(rr,&r));
198c4762a1bSJed Brown   PetscFunctionReturn(0);
199c4762a1bSJed Brown }
200f2fddbb2SStefano Zampini 
201c4762a1bSJed Brown PetscErrorCode amult(Mat mat,Vec xx,Vec yy)
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   PetscInt          i,n1;
204c4762a1bSJed Brown   PetscScalar       *y;
205c4762a1bSJed Brown   const PetscScalar *x;
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   PetscFunctionBegin;
2089566063dSJacob Faibussowitsch   PetscCall(VecGetSize(xx,&n1));
2099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx,&x));
2109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(yy,&y));
211c4762a1bSJed Brown   n1--;
212c4762a1bSJed Brown   y[0] =  -x[1] + 2.0*x[0];
213c4762a1bSJed Brown   y[n1] = -x[n1-1] + 2.0*x[n1];
214c4762a1bSJed Brown   for (i=1; i<n1; i++) y[i] = -x[i+1] - x[i-1] + 2.0*x[i];
2159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx,&x));
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(yy,&y));
217c4762a1bSJed Brown   PetscFunctionReturn(0);
218c4762a1bSJed Brown }
219f2fddbb2SStefano Zampini 
220c4762a1bSJed Brown /* --------------------------------------------------------------------- */
221f2fddbb2SStefano Zampini PetscErrorCode apply_pc(PC pc,Vec bb,Vec xx)
222f2fddbb2SStefano Zampini {
223f2fddbb2SStefano Zampini   PetscFunctionBegin;
224f2fddbb2SStefano Zampini   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not implemented");
225f2fddbb2SStefano Zampini }
226f2fddbb2SStefano Zampini 
227c4762a1bSJed 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)
228c4762a1bSJed Brown {
229c4762a1bSJed Brown   PetscInt          i,n1;
230c4762a1bSJed Brown   PetscScalar       *x;
231c4762a1bSJed Brown   const PetscScalar *b;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   PetscFunctionBegin;
234c4762a1bSJed Brown   *its    = m;
235c4762a1bSJed Brown   *reason = PCRICHARDSON_CONVERGED_ITS;
2369566063dSJacob Faibussowitsch   PetscCall(VecGetSize(bb,&n1)); n1--;
2379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb,&b));
2389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx,&x));
239c4762a1bSJed Brown   while (m--) {
240c4762a1bSJed Brown     x[0] =  .5*(x[1] + b[0]);
241c4762a1bSJed Brown     for (i=1; i<n1; i++) x[i] = .5*(x[i+1] + x[i-1] + b[i]);
242c4762a1bSJed Brown     x[n1] = .5*(x[n1-1] + b[n1]);
243c4762a1bSJed Brown     for (i=n1-1; i>0; i--) x[i] = .5*(x[i+1] + x[i-1] + b[i]);
244c4762a1bSJed Brown     x[0] =  .5*(x[1] + b[0]);
245c4762a1bSJed Brown   }
2469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb,&b));
2479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx,&x));
248c4762a1bSJed Brown   PetscFunctionReturn(0);
249c4762a1bSJed Brown }
250c4762a1bSJed Brown /* --------------------------------------------------------------------- */
251c4762a1bSJed 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)
252c4762a1bSJed Brown {
253c4762a1bSJed Brown   PetscInt          i,n,n1;
254c4762a1bSJed Brown   PetscScalar       *r,*x;
255c4762a1bSJed Brown   const PetscScalar *b;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   PetscFunctionBegin;
258c4762a1bSJed Brown   *its    = m;
259c4762a1bSJed Brown   *reason = PCRICHARDSON_CONVERGED_ITS;
2609566063dSJacob Faibussowitsch   PetscCall(VecGetSize(bb,&n)); n1 = n - 1;
2619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb,&b));
2629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx,&x));
2639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(w,&r));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   while (m--) {
266c4762a1bSJed Brown     r[0] = .5*(x[1] + b[0]);
267c4762a1bSJed Brown     for (i=1; i<n1; i++) r[i] = .5*(x[i+1] + x[i-1] + b[i]);
268c4762a1bSJed Brown     r[n1] = .5*(x[n1-1] + b[n1]);
269c4762a1bSJed Brown     for (i=0; i<n; i++) x[i] = (2.0*r[i] + x[i])/3.0;
270c4762a1bSJed Brown   }
2719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb,&b));
2729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx,&x));
2739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(w,&r));
274c4762a1bSJed Brown   PetscFunctionReturn(0);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown /*
277c4762a1bSJed Brown    We know for this application that yy  and zz are the same
278c4762a1bSJed Brown */
279c4762a1bSJed Brown /* --------------------------------------------------------------------- */
280c4762a1bSJed Brown PetscErrorCode interpolate(Mat mat,Vec xx,Vec yy,Vec zz)
281c4762a1bSJed Brown {
282c4762a1bSJed Brown   PetscInt          i,n,N,i2;
283c4762a1bSJed Brown   PetscScalar       *y;
284c4762a1bSJed Brown   const PetscScalar *x;
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   PetscFunctionBegin;
2879566063dSJacob Faibussowitsch   PetscCall(VecGetSize(yy,&N));
2889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx,&x));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(yy,&y));
290c4762a1bSJed Brown   n    = N/2;
291c4762a1bSJed Brown   for (i=0; i<n; i++) {
292c4762a1bSJed Brown     i2       = 2*i;
293c4762a1bSJed Brown     y[i2]   += .5*x[i];
294c4762a1bSJed Brown     y[i2+1] +=    x[i];
295c4762a1bSJed Brown     y[i2+2] += .5*x[i];
296c4762a1bSJed Brown   }
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx,&x));
2989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(yy,&y));
299c4762a1bSJed Brown   PetscFunctionReturn(0);
300c4762a1bSJed Brown }
301c4762a1bSJed Brown /* --------------------------------------------------------------------- */
302c4762a1bSJed Brown PetscErrorCode restrct(Mat mat,Vec rr,Vec bb)
303c4762a1bSJed Brown {
304c4762a1bSJed Brown   PetscInt          i,n,N,i2;
305c4762a1bSJed Brown   PetscScalar       *b;
306c4762a1bSJed Brown   const PetscScalar *r;
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(VecGetSize(rr,&N));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(rr,&r));
3119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(bb,&b));
312c4762a1bSJed Brown   n    = N/2;
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   for (i=0; i<n; i++) {
315c4762a1bSJed Brown     i2   = 2*i;
316c4762a1bSJed Brown     b[i] = (r[i2] + 2.0*r[i2+1] + r[i2+2]);
317c4762a1bSJed Brown   }
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(rr,&r));
3199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(bb,&b));
320c4762a1bSJed Brown   PetscFunctionReturn(0);
321c4762a1bSJed Brown }
322c4762a1bSJed Brown /* --------------------------------------------------------------------- */
323c4762a1bSJed Brown PetscErrorCode Create1dLaplacian(PetscInt n,Mat *mat)
324c4762a1bSJed Brown {
325c4762a1bSJed Brown   PetscScalar    mone = -1.0,two = 2.0;
326c4762a1bSJed Brown   PetscInt       i,idx;
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   PetscFunctionBegin;
3299566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,mat));
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   idx  = n-1;
3329566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*mat,1,&idx,1,&idx,&two,INSERT_VALUES));
333c4762a1bSJed Brown   for (i=0; i<n-1; i++) {
3349566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*mat,1,&i,1,&i,&two,INSERT_VALUES));
335c4762a1bSJed Brown     idx  = i+1;
3369566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*mat,1,&idx,1,&i,&mone,INSERT_VALUES));
3379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*mat,1,&i,1,&idx,&mone,INSERT_VALUES));
338c4762a1bSJed Brown   }
3399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY));
3409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY));
341c4762a1bSJed Brown   PetscFunctionReturn(0);
342c4762a1bSJed Brown }
343c4762a1bSJed Brown /* --------------------------------------------------------------------- */
344c4762a1bSJed Brown PetscErrorCode CalculateRhs(Vec u)
345c4762a1bSJed Brown {
346c4762a1bSJed Brown   PetscInt    i,n;
3475f80ce2aSJacob Faibussowitsch   PetscReal   h;
348c4762a1bSJed Brown   PetscScalar uu;
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(VecGetSize(u,&n));
352c4762a1bSJed Brown   h    = 1.0/((PetscReal)(n+1));
353c4762a1bSJed Brown   for (i=0; i<n; i++) {
3545f80ce2aSJacob Faibussowitsch     uu = 2.0*h*h;
3559566063dSJacob Faibussowitsch     PetscCall(VecSetValues(u,1,&i,&uu,INSERT_VALUES));
356c4762a1bSJed Brown   }
357c4762a1bSJed Brown   PetscFunctionReturn(0);
358c4762a1bSJed Brown }
359c4762a1bSJed Brown /* --------------------------------------------------------------------- */
360c4762a1bSJed Brown PetscErrorCode CalculateSolution(PetscInt n,Vec *solution)
361c4762a1bSJed Brown {
362c4762a1bSJed Brown   PetscInt       i;
363c4762a1bSJed Brown   PetscReal      h,x = 0.0;
364c4762a1bSJed Brown   PetscScalar    uu;
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   PetscFunctionBegin;
3679566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,solution));
368c4762a1bSJed Brown   h    = 1.0/((PetscReal)(n+1));
369c4762a1bSJed Brown   for (i=0; i<n; i++) {
370c4762a1bSJed Brown     x   += h; uu = x*(1.-x);
3719566063dSJacob Faibussowitsch     PetscCall(VecSetValues(*solution,1,&i,&uu,INSERT_VALUES));
372c4762a1bSJed Brown   }
373c4762a1bSJed Brown   PetscFunctionReturn(0);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown /* --------------------------------------------------------------------- */
376c4762a1bSJed Brown PetscErrorCode CalculateError(Vec solution,Vec u,Vec r,PetscReal *e)
377c4762a1bSJed Brown {
378c4762a1bSJed Brown   PetscFunctionBegin;
3799566063dSJacob Faibussowitsch   PetscCall(VecNorm(r,NORM_2,e+2));
3809566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(r,-1.0,u,solution));
3819566063dSJacob Faibussowitsch   PetscCall(VecNorm(r,NORM_2,e));
3829566063dSJacob Faibussowitsch   PetscCall(VecNorm(r,NORM_1,e+1));
383c4762a1bSJed Brown   PetscFunctionReturn(0);
384c4762a1bSJed Brown }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown /*TEST
387c4762a1bSJed Brown 
388c4762a1bSJed Brown    test:
389c4762a1bSJed Brown 
390c4762a1bSJed Brown TEST*/
391