1c4762a1bSJed Brown static char help[] = "Tests the multigrid code. The input parameters are:\n\
2c4762a1bSJed Brown -x N Use a mesh in the x direction of N. \n\
3c4762a1bSJed Brown -c N Use N V-cycles. \n\
4c4762a1bSJed Brown -l N Use N Levels. \n\
5c4762a1bSJed Brown -smooths N Use N pre smooths and N post smooths. \n\
6c4762a1bSJed Brown -j Use Jacobi smoother. \n\
7c4762a1bSJed Brown -a use additive multigrid \n\
8c4762a1bSJed Brown -f use full multigrid (preconditioner variant) \n\
9c4762a1bSJed Brown This example also demonstrates matrix-free methods\n\n";
10c4762a1bSJed Brown
11c4762a1bSJed Brown /*
12c4762a1bSJed Brown This is not a good example to understand the use of multigrid with PETSc.
13c4762a1bSJed Brown */
14c4762a1bSJed Brown
15c4762a1bSJed Brown #include <petscksp.h>
16c4762a1bSJed Brown
17c4762a1bSJed Brown PetscErrorCode residual(Mat, Vec, Vec, Vec);
18c4762a1bSJed Brown PetscErrorCode gauss_seidel(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *);
19c4762a1bSJed Brown PetscErrorCode jacobi_smoother(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *);
20c4762a1bSJed Brown PetscErrorCode interpolate(Mat, Vec, Vec, Vec);
21c4762a1bSJed Brown PetscErrorCode restrct(Mat, Vec, Vec);
22c4762a1bSJed Brown PetscErrorCode Create1dLaplacian(PetscInt, Mat *);
23c4762a1bSJed Brown PetscErrorCode CalculateRhs(Vec);
24c4762a1bSJed Brown PetscErrorCode CalculateError(Vec, Vec, Vec, PetscReal *);
25c4762a1bSJed Brown PetscErrorCode CalculateSolution(PetscInt, Vec *);
26c4762a1bSJed Brown PetscErrorCode amult(Mat, Vec, Vec);
27f2fddbb2SStefano Zampini PetscErrorCode apply_pc(PC, Vec, Vec);
28c4762a1bSJed Brown
main(int Argc,char ** Args)29d71ae5a4SJacob Faibussowitsch int main(int Argc, char **Args)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown PetscInt x_mesh = 15, levels = 3, cycles = 1, use_jacobi = 0;
32c4762a1bSJed Brown PetscInt i, smooths = 1, *N, its;
33c4762a1bSJed Brown PCMGType am = PC_MG_MULTIPLICATIVE;
34c4762a1bSJed Brown Mat cmat, mat[20], fmat;
35c4762a1bSJed Brown KSP cksp, ksp[20], kspmg;
36c4762a1bSJed Brown PetscReal e[3]; /* l_2 error,max error, residual */
37c4762a1bSJed Brown const char *shellname;
38c4762a1bSJed Brown Vec x, solution, X[20], R[20], B[20];
39c4762a1bSJed Brown PC pcmg, pc;
40c4762a1bSJed Brown PetscBool flg;
41c4762a1bSJed Brown
42c8025a54SPierre Jolivet PetscCall(PetscInitialize(&Argc, &Args, NULL, help));
439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-x", &x_mesh, NULL));
449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-l", &levels, NULL));
459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-c", &cycles, NULL));
469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-smooths", &smooths, NULL));
479566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-a", &flg));
48c4762a1bSJed Brown
49c4762a1bSJed Brown if (flg) am = PC_MG_ADDITIVE;
509566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-f", &flg));
51c4762a1bSJed Brown if (flg) am = PC_MG_FULL;
529566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-j", &flg));
53c4762a1bSJed Brown if (flg) use_jacobi = 1;
54c4762a1bSJed Brown
559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(levels, &N));
56c4762a1bSJed Brown N[0] = x_mesh;
57c4762a1bSJed Brown for (i = 1; i < levels; i++) {
58c4762a1bSJed Brown N[i] = N[i - 1] / 2;
597827d75bSBarry Smith PetscCheck(N[i] >= 1, PETSC_COMM_WORLD, PETSC_ERR_USER, "Too many levels or N is not large enough");
60c4762a1bSJed Brown }
61c4762a1bSJed Brown
629566063dSJacob Faibussowitsch PetscCall(Create1dLaplacian(N[levels - 1], &cmat));
63c4762a1bSJed Brown
649566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &kspmg));
659566063dSJacob Faibussowitsch PetscCall(KSPGetPC(kspmg, &pcmg));
669566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(kspmg));
679566063dSJacob Faibussowitsch PetscCall(PCSetType(pcmg, PCMG));
689566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pcmg, levels, NULL));
699566063dSJacob Faibussowitsch PetscCall(PCMGSetType(pcmg, am));
70c4762a1bSJed Brown
719566063dSJacob Faibussowitsch PetscCall(PCMGGetCoarseSolve(pcmg, &cksp));
729566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(cksp, cmat, cmat));
739566063dSJacob Faibussowitsch PetscCall(KSPGetPC(cksp, &pc));
749566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCLU));
759566063dSJacob Faibussowitsch PetscCall(KSPSetType(cksp, KSPPREONLY));
76c4762a1bSJed Brown
77c4762a1bSJed Brown /* zero is finest level */
78c4762a1bSJed Brown for (i = 0; i < levels - 1; i++) {
79f2fddbb2SStefano Zampini Mat dummy;
80f2fddbb2SStefano Zampini
819566063dSJacob Faibussowitsch PetscCall(PCMGSetResidual(pcmg, levels - 1 - i, residual, NULL));
829566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, N[i + 1], N[i], N[i + 1], N[i], NULL, &mat[i]));
83*57d50842SBarry Smith PetscCall(MatShellSetOperation(mat[i], MATOP_MULT, (PetscErrorCodeFn *)restrct));
84*57d50842SBarry Smith PetscCall(MatShellSetOperation(mat[i], MATOP_MULT_TRANSPOSE_ADD, (PetscErrorCodeFn *)interpolate));
859566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pcmg, levels - 1 - i, mat[i]));
869566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pcmg, levels - 1 - i, mat[i]));
879566063dSJacob Faibussowitsch PetscCall(PCMGSetCycleTypeOnLevel(pcmg, levels - 1 - i, (PCMGCycleType)cycles));
88c4762a1bSJed Brown
89c4762a1bSJed Brown /* set smoother */
909566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pcmg, levels - 1 - i, &ksp[i]));
919566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp[i], &pc));
929566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCSHELL));
939566063dSJacob Faibussowitsch PetscCall(PCShellSetName(pc, "user_precond"));
949566063dSJacob Faibussowitsch PetscCall(PCShellGetName(pc, &shellname));
9563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "level=%" PetscInt_FMT ", PCShell name is %s\n", i, shellname));
96c4762a1bSJed Brown
97f2fddbb2SStefano Zampini /* this is not used unless different options are passed to the solver */
989566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, N[i], N[i], N[i], N[i], NULL, &dummy));
99*57d50842SBarry Smith PetscCall(MatShellSetOperation(dummy, MATOP_MULT, (PetscErrorCodeFn *)amult));
1009566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp[i], dummy, dummy));
1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&dummy));
102f2fddbb2SStefano Zampini
103c4762a1bSJed Brown /*
104c4762a1bSJed Brown We override the matrix passed in by forcing it to use Richardson with
105c4762a1bSJed Brown a user provided application. This is non-standard and this practice
106c4762a1bSJed Brown should be avoided.
107c4762a1bSJed Brown */
1089566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(pc, apply_pc));
1099566063dSJacob Faibussowitsch PetscCall(PCShellSetApplyRichardson(pc, gauss_seidel));
1101baa6e33SBarry Smith if (use_jacobi) PetscCall(PCShellSetApplyRichardson(pc, jacobi_smoother));
1119566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp[i], KSPRICHARDSON));
1129566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(ksp[i], PETSC_TRUE));
113fb842aefSJose E. Roman PetscCall(KSPSetTolerances(ksp[i], PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, smooths));
114c4762a1bSJed Brown
1159566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N[i], &x));
116c4762a1bSJed Brown
117c4762a1bSJed Brown X[levels - 1 - i] = x;
11848a46eb9SPierre Jolivet if (i > 0) PetscCall(PCMGSetX(pcmg, levels - 1 - i, x));
1199566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N[i], &x));
120c4762a1bSJed Brown
121c4762a1bSJed Brown B[levels - 1 - i] = x;
12248a46eb9SPierre Jolivet if (i > 0) PetscCall(PCMGSetRhs(pcmg, levels - 1 - i, x));
1239566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N[i], &x));
124c4762a1bSJed Brown
125c4762a1bSJed Brown R[levels - 1 - i] = x;
126c4762a1bSJed Brown
1279566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pcmg, levels - 1 - i, x));
128c4762a1bSJed Brown }
129c4762a1bSJed Brown /* create coarse level vectors */
1309566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N[levels - 1], &x));
1319371c9d4SSatish Balay PetscCall(PCMGSetX(pcmg, 0, x));
1329371c9d4SSatish Balay X[0] = x;
1339566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N[levels - 1], &x));
1349371c9d4SSatish Balay PetscCall(PCMGSetRhs(pcmg, 0, x));
1359371c9d4SSatish Balay B[0] = x;
136c4762a1bSJed Brown
137c4762a1bSJed Brown /* create matrix multiply for finest level */
1389566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, N[0], N[0], N[0], N[0], NULL, &fmat));
139*57d50842SBarry Smith PetscCall(MatShellSetOperation(fmat, MATOP_MULT, (PetscErrorCodeFn *)amult));
1409566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(kspmg, fmat, fmat));
141c4762a1bSJed Brown
1429566063dSJacob Faibussowitsch PetscCall(CalculateSolution(N[0], &solution));
1439566063dSJacob Faibussowitsch PetscCall(CalculateRhs(B[levels - 1]));
1449566063dSJacob Faibussowitsch PetscCall(VecSet(X[levels - 1], 0.0));
145c4762a1bSJed Brown
1469566063dSJacob Faibussowitsch PetscCall(residual((Mat)0, B[levels - 1], X[levels - 1], R[levels - 1]));
1479566063dSJacob Faibussowitsch PetscCall(CalculateError(solution, X[levels - 1], R[levels - 1], e));
1489566063dSJacob 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]));
149c4762a1bSJed Brown
1509566063dSJacob Faibussowitsch PetscCall(KSPSolve(kspmg, B[levels - 1], X[levels - 1]));
1519566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(kspmg, &its));
1529566063dSJacob Faibussowitsch PetscCall(residual((Mat)0, B[levels - 1], X[levels - 1], R[levels - 1]));
1539566063dSJacob Faibussowitsch PetscCall(CalculateError(solution, X[levels - 1], R[levels - 1], e));
15463a3b9bcSJacob 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]));
155c4762a1bSJed Brown
1569566063dSJacob Faibussowitsch PetscCall(PetscFree(N));
1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&solution));
158c4762a1bSJed Brown
159c4762a1bSJed Brown /* note we have to keep a list of all vectors allocated, this is
160c4762a1bSJed Brown not ideal, but putting it in MGDestroy is not so good either*/
161c4762a1bSJed Brown for (i = 0; i < levels; i++) {
1629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X[i]));
1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&B[i]));
1649566063dSJacob Faibussowitsch if (i) PetscCall(VecDestroy(&R[i]));
165c4762a1bSJed Brown }
16648a46eb9SPierre Jolivet for (i = 0; i < levels - 1; i++) PetscCall(MatDestroy(&mat[i]));
1679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&cmat));
1689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fmat));
1699566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&kspmg));
1709566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
171b122ec5aSJacob Faibussowitsch return 0;
172c4762a1bSJed Brown }
173c4762a1bSJed Brown
residual(Mat mat,Vec bb,Vec xx,Vec rr)174d71ae5a4SJacob Faibussowitsch PetscErrorCode residual(Mat mat, Vec bb, Vec xx, Vec rr)
175d71ae5a4SJacob Faibussowitsch {
176c4762a1bSJed Brown PetscInt i, n1;
177c4762a1bSJed Brown PetscScalar *x, *r;
178c4762a1bSJed Brown const PetscScalar *b;
179c4762a1bSJed Brown
180c4762a1bSJed Brown PetscFunctionBegin;
1819566063dSJacob Faibussowitsch PetscCall(VecGetSize(bb, &n1));
1829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
1839566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
1849566063dSJacob Faibussowitsch PetscCall(VecGetArray(rr, &r));
185c4762a1bSJed Brown n1--;
186c4762a1bSJed Brown r[0] = b[0] + x[1] - 2.0 * x[0];
187c4762a1bSJed Brown r[n1] = b[n1] + x[n1 - 1] - 2.0 * x[n1];
188c4762a1bSJed Brown for (i = 1; i < n1; i++) r[i] = b[i] + x[i + 1] + x[i - 1] - 2.0 * x[i];
1899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
1909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
1919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(rr, &r));
1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
193c4762a1bSJed Brown }
194f2fddbb2SStefano Zampini
amult(Mat mat,Vec xx,Vec yy)195d71ae5a4SJacob Faibussowitsch PetscErrorCode amult(Mat mat, Vec xx, Vec yy)
196d71ae5a4SJacob Faibussowitsch {
197c4762a1bSJed Brown PetscInt i, n1;
198c4762a1bSJed Brown PetscScalar *y;
199c4762a1bSJed Brown const PetscScalar *x;
200c4762a1bSJed Brown
201c4762a1bSJed Brown PetscFunctionBegin;
2029566063dSJacob Faibussowitsch PetscCall(VecGetSize(xx, &n1));
2039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xx, &x));
2049566063dSJacob Faibussowitsch PetscCall(VecGetArray(yy, &y));
205c4762a1bSJed Brown n1--;
206c4762a1bSJed Brown y[0] = -x[1] + 2.0 * x[0];
207c4762a1bSJed Brown y[n1] = -x[n1 - 1] + 2.0 * x[n1];
208c4762a1bSJed Brown for (i = 1; i < n1; i++) y[i] = -x[i + 1] - x[i - 1] + 2.0 * x[i];
2099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xx, &x));
2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(yy, &y));
2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
212c4762a1bSJed Brown }
213f2fddbb2SStefano Zampini
apply_pc(PC pc,Vec bb,Vec xx)214d71ae5a4SJacob Faibussowitsch PetscErrorCode apply_pc(PC pc, Vec bb, Vec xx)
215d71ae5a4SJacob Faibussowitsch {
216f2fddbb2SStefano Zampini PetscFunctionBegin;
217f2fddbb2SStefano Zampini SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not implemented");
218f2fddbb2SStefano Zampini }
219f2fddbb2SStefano Zampini
gauss_seidel(PC pc,Vec bb,Vec xx,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt m,PetscBool guesszero,PetscInt * its,PCRichardsonConvergedReason * reason)220d71ae5a4SJacob Faibussowitsch 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)
221d71ae5a4SJacob Faibussowitsch {
222c4762a1bSJed Brown PetscInt i, n1;
223c4762a1bSJed Brown PetscScalar *x;
224c4762a1bSJed Brown const PetscScalar *b;
225c4762a1bSJed Brown
226c4762a1bSJed Brown PetscFunctionBegin;
227c4762a1bSJed Brown *its = m;
228c4762a1bSJed Brown *reason = PCRICHARDSON_CONVERGED_ITS;
2299371c9d4SSatish Balay PetscCall(VecGetSize(bb, &n1));
2309371c9d4SSatish Balay n1--;
2319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
2329566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
233c4762a1bSJed Brown while (m--) {
234c4762a1bSJed Brown x[0] = .5 * (x[1] + b[0]);
235c4762a1bSJed Brown for (i = 1; i < n1; i++) x[i] = .5 * (x[i + 1] + x[i - 1] + b[i]);
236c4762a1bSJed Brown x[n1] = .5 * (x[n1 - 1] + b[n1]);
237c4762a1bSJed Brown for (i = n1 - 1; i > 0; i--) x[i] = .5 * (x[i + 1] + x[i - 1] + b[i]);
238c4762a1bSJed Brown x[0] = .5 * (x[1] + b[0]);
239c4762a1bSJed Brown }
2409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
2419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
243c4762a1bSJed Brown }
244f1580f4eSBarry Smith
jacobi_smoother(PC pc,Vec bb,Vec xx,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt m,PetscBool guesszero,PetscInt * its,PCRichardsonConvergedReason * reason)245d71ae5a4SJacob Faibussowitsch 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)
246d71ae5a4SJacob Faibussowitsch {
247c4762a1bSJed Brown PetscInt i, n, n1;
248c4762a1bSJed Brown PetscScalar *r, *x;
249c4762a1bSJed Brown const PetscScalar *b;
250c4762a1bSJed Brown
251c4762a1bSJed Brown PetscFunctionBegin;
252c4762a1bSJed Brown *its = m;
253c4762a1bSJed Brown *reason = PCRICHARDSON_CONVERGED_ITS;
2549371c9d4SSatish Balay PetscCall(VecGetSize(bb, &n));
2559371c9d4SSatish Balay n1 = n - 1;
2569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
2579566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
2589566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &r));
259c4762a1bSJed Brown
260c4762a1bSJed Brown while (m--) {
261c4762a1bSJed Brown r[0] = .5 * (x[1] + b[0]);
262c4762a1bSJed Brown for (i = 1; i < n1; i++) r[i] = .5 * (x[i + 1] + x[i - 1] + b[i]);
263c4762a1bSJed Brown r[n1] = .5 * (x[n1 - 1] + b[n1]);
264c4762a1bSJed Brown for (i = 0; i < n; i++) x[i] = (2.0 * r[i] + x[i]) / 3.0;
265c4762a1bSJed Brown }
2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &r));
2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
270c4762a1bSJed Brown }
271c4762a1bSJed Brown /*
272c4762a1bSJed Brown We know for this application that yy and zz are the same
273c4762a1bSJed Brown */
274f1580f4eSBarry Smith
interpolate(Mat mat,Vec xx,Vec yy,Vec zz)275d71ae5a4SJacob Faibussowitsch PetscErrorCode interpolate(Mat mat, Vec xx, Vec yy, Vec zz)
276d71ae5a4SJacob Faibussowitsch {
277c4762a1bSJed Brown PetscInt i, n, N, i2;
278c4762a1bSJed Brown PetscScalar *y;
279c4762a1bSJed Brown const PetscScalar *x;
280c4762a1bSJed Brown
281c4762a1bSJed Brown PetscFunctionBegin;
2829566063dSJacob Faibussowitsch PetscCall(VecGetSize(yy, &N));
2839566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xx, &x));
2849566063dSJacob Faibussowitsch PetscCall(VecGetArray(yy, &y));
285c4762a1bSJed Brown n = N / 2;
286c4762a1bSJed Brown for (i = 0; i < n; i++) {
287c4762a1bSJed Brown i2 = 2 * i;
288c4762a1bSJed Brown y[i2] += .5 * x[i];
289c4762a1bSJed Brown y[i2 + 1] += x[i];
290c4762a1bSJed Brown y[i2 + 2] += .5 * x[i];
291c4762a1bSJed Brown }
2929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xx, &x));
2939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(yy, &y));
2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
295c4762a1bSJed Brown }
296f1580f4eSBarry Smith
restrct(Mat mat,Vec rr,Vec bb)297d71ae5a4SJacob Faibussowitsch PetscErrorCode restrct(Mat mat, Vec rr, Vec bb)
298d71ae5a4SJacob Faibussowitsch {
299c4762a1bSJed Brown PetscInt i, n, N, i2;
300c4762a1bSJed Brown PetscScalar *b;
301c4762a1bSJed Brown const PetscScalar *r;
302c4762a1bSJed Brown
303c4762a1bSJed Brown PetscFunctionBegin;
3049566063dSJacob Faibussowitsch PetscCall(VecGetSize(rr, &N));
3059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(rr, &r));
3069566063dSJacob Faibussowitsch PetscCall(VecGetArray(bb, &b));
307c4762a1bSJed Brown n = N / 2;
308c4762a1bSJed Brown
309c4762a1bSJed Brown for (i = 0; i < n; i++) {
310c4762a1bSJed Brown i2 = 2 * i;
311c4762a1bSJed Brown b[i] = (r[i2] + 2.0 * r[i2 + 1] + r[i2 + 2]);
312c4762a1bSJed Brown }
3139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(rr, &r));
3149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(bb, &b));
3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
316c4762a1bSJed Brown }
317f1580f4eSBarry Smith
Create1dLaplacian(PetscInt n,Mat * mat)318d71ae5a4SJacob Faibussowitsch PetscErrorCode Create1dLaplacian(PetscInt n, Mat *mat)
319d71ae5a4SJacob Faibussowitsch {
32042e5ec60SJeff-Hadley PetscFunctionBeginUser;
3219566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n, n, 3, NULL, mat));
32242e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, n - 1, n - 1, 2.0, INSERT_VALUES));
32342e5ec60SJeff-Hadley for (PetscInt i = 0; i < n - 1; i++) {
32442e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, i, i, 2.0, INSERT_VALUES));
32542e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, i + 1, i, -1.0, INSERT_VALUES));
32642e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, i, i + 1, -1.0, INSERT_VALUES));
327c4762a1bSJed Brown }
3289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
3299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
331c4762a1bSJed Brown }
332f1580f4eSBarry Smith
CalculateRhs(Vec u)333d71ae5a4SJacob Faibussowitsch PetscErrorCode CalculateRhs(Vec u)
334d71ae5a4SJacob Faibussowitsch {
335c4762a1bSJed Brown PetscInt i, n;
3365f80ce2aSJacob Faibussowitsch PetscReal h;
337c4762a1bSJed Brown PetscScalar uu;
338c4762a1bSJed Brown
339c4762a1bSJed Brown PetscFunctionBegin;
3409566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n));
341c4762a1bSJed Brown h = 1.0 / ((PetscReal)(n + 1));
342c4762a1bSJed Brown for (i = 0; i < n; i++) {
3435f80ce2aSJacob Faibussowitsch uu = 2.0 * h * h;
3449566063dSJacob Faibussowitsch PetscCall(VecSetValues(u, 1, &i, &uu, INSERT_VALUES));
345c4762a1bSJed Brown }
3463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
347c4762a1bSJed Brown }
348f1580f4eSBarry Smith
CalculateSolution(PetscInt n,Vec * solution)349d71ae5a4SJacob Faibussowitsch PetscErrorCode CalculateSolution(PetscInt n, Vec *solution)
350d71ae5a4SJacob Faibussowitsch {
351c4762a1bSJed Brown PetscInt i;
352c4762a1bSJed Brown PetscReal h, x = 0.0;
353c4762a1bSJed Brown PetscScalar uu;
354c4762a1bSJed Brown
355c4762a1bSJed Brown PetscFunctionBegin;
3569566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, solution));
357c4762a1bSJed Brown h = 1.0 / ((PetscReal)(n + 1));
358c4762a1bSJed Brown for (i = 0; i < n; i++) {
3599371c9d4SSatish Balay x += h;
3609371c9d4SSatish Balay uu = x * (1. - x);
3619566063dSJacob Faibussowitsch PetscCall(VecSetValues(*solution, 1, &i, &uu, INSERT_VALUES));
362c4762a1bSJed Brown }
3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
364c4762a1bSJed Brown }
365f1580f4eSBarry Smith
CalculateError(Vec solution,Vec u,Vec r,PetscReal * e)366d71ae5a4SJacob Faibussowitsch PetscErrorCode CalculateError(Vec solution, Vec u, Vec r, PetscReal *e)
367d71ae5a4SJacob Faibussowitsch {
368c4762a1bSJed Brown PetscFunctionBegin;
3699566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, e + 2));
3709566063dSJacob Faibussowitsch PetscCall(VecWAXPY(r, -1.0, u, solution));
3719566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, e));
3729566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_1, e + 1));
3733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown
376c4762a1bSJed Brown /*TEST
377c4762a1bSJed Brown
378c4762a1bSJed Brown test:
379c4762a1bSJed Brown
380c4762a1bSJed Brown TEST*/
381