xref: /petsc/src/tao/tutorials/ex3.c (revision 2fa40bb9206b96114faa7cb222621ec184d31cd2)
1 static char help[] = "Reduced formulation of the mother problem of PDE-constrained optimisation.\n";
2 
3 /*F
4   We solve the mother problem
5 
6   min 1/2 || y - y_d ||^2_2 + \alpha/2 || u ||^2_2
7 
8   subject to
9 
10           - \laplace y = u          on \Omega
11                      y = 0          on \Gamma
12 
13   where u is in L^2 and y is in H^1_0.
14 
15   We formulate the reduced problem solely in terms of the control
16   by using the state equation to express y in terms of u, and then
17   apply LMVM/BLMVM to the resulting reduced problem.
18 
19   Mesh independence is achieved by configuring the Riesz map for the control
20   space.
21 
22   Example meshes where the Riesz map is crucial can be downloaded from the
23   http://people.maths.ox.ac.uk/~farrellp/meshes.tar.gz
24 
25   Contributed by: Patrick Farrell <patrick.farrell@maths.ox.ac.uk>
26 
27   Run with e.g.:
28   ./ex3 -laplace_ksp_type cg -laplace_pc_type hypre -mat_lmvm_ksp_type cg -mat_lmvm_pc_type gamg -laplace_ksp_monitor_true_residual -tao_monitor -petscspace_degree 1 -tao_converged_reason -tao_gatol 1.0e-9 -dm_view hdf5:solution.h5 -sol_view hdf5:solution.h5::append -use_riesz 1 -f $DATAFILESPATH/meshes/mesh-1.h5
29 
30   and visualise in paraview with ../../../../petsc_gen_xdmf.py solution.h5.
31 
32   Toggle the Riesz map (-use_riesz 0) to see the difference setting the Riesz maps makes.
33 
34   TODO: broken for parallel runs
35 F*/
36 
37 #include <petsc.h>
38 #include <petscfe.h>
39 #include <petscviewerhdf5.h>
40 
41 typedef struct {
42   DM  dm;
43   Mat mass;
44   Vec data;
45   Vec state;
46   Vec tmp1;
47   Vec tmp2;
48   Vec adjoint;
49   Mat laplace;
50   KSP ksp_laplace;
51   PetscInt  num_bc_dofs;
52   PetscInt* bc_indices;
53   PetscScalar* bc_values;
54   PetscBool use_riesz;
55 } AppCtx;
56 
57 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
58 {
59   PetscBool      flg;
60   char           filename[2048];
61   PetscErrorCode ierr;
62 
63   PetscFunctionBeginUser;
64   filename[0] = '\0';
65   user->use_riesz = PETSC_TRUE;
66 
67   ierr = PetscOptionsBegin(comm, "", "Poisson mother problem options", "DMPLEX");CHKERRQ(ierr);
68   ierr = PetscOptionsBool("-use_riesz", "Use the Riesz map to achieve mesh independence", "ex3.c", user->use_riesz, &user->use_riesz, NULL);CHKERRQ(ierr);
69   ierr = PetscOptionsString("-f", "filename to read", "ex3.c", filename, filename, sizeof(filename), &flg);CHKERRQ(ierr);
70   ierr = PetscOptionsEnd();CHKERRQ(ierr);
71 
72   if (!flg) {
73     ierr = DMCreate(comm, dm);CHKERRQ(ierr);
74     ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
75   } else {
76     /* TODO Eliminate this in favor of DMLoad() in new code */
77 #if defined(PETSC_HAVE_HDF5)
78     const PetscInt vertices_per_cell = 3;
79     PetscViewer    viewer;
80     Vec            coordinates;
81     Vec            topology;
82     PetscInt       dim = 2, numCells;
83     PetscInt       numVertices;
84     PetscScalar*   coords;
85     PetscScalar*   topo_f;
86     PetscInt*      cells;
87     PetscInt       j;
88     DMLabel        label;
89 
90     /* Read in FEniCS HDF5 output */
91     ierr = PetscViewerHDF5Open(comm, filename, FILE_MODE_READ, &viewer);CHKERRQ(ierr);
92 
93     /* create Vecs to read in the data from H5 */
94     ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
95     ierr = PetscObjectSetName((PetscObject)coordinates, "coordinates");CHKERRQ(ierr);
96     ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
97     ierr = VecCreate(comm, &topology);CHKERRQ(ierr);
98     ierr = PetscObjectSetName((PetscObject)topology, "topology");CHKERRQ(ierr);
99     ierr = VecSetBlockSize(topology, vertices_per_cell);CHKERRQ(ierr);
100 
101     /* navigate to the right group */
102     ierr = PetscViewerHDF5PushGroup(viewer, "/Mesh/mesh");CHKERRQ(ierr);
103 
104     /* Read the Vecs */
105     ierr = VecLoad(coordinates, viewer);CHKERRQ(ierr);
106     ierr = VecLoad(topology, viewer);CHKERRQ(ierr);
107 
108     /* do some ugly calculations */
109     ierr = VecGetSize(topology, &numCells);CHKERRQ(ierr);
110     numCells = numCells / vertices_per_cell;
111     ierr = VecGetSize(coordinates, &numVertices);CHKERRQ(ierr);
112     numVertices = numVertices / dim;
113 
114     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
115     ierr = VecGetArray(topology, &topo_f);CHKERRQ(ierr);
116     /* and now we have to convert the double representation to integers to pass over, argh */
117     ierr = PetscMalloc1(numCells*vertices_per_cell, &cells);CHKERRQ(ierr);
118     for (j = 0; j < numCells*vertices_per_cell; j++) {
119       cells[j] = (PetscInt) topo_f[j];
120     }
121 
122     /* Now create the DM */
123     ierr = DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, vertices_per_cell, PETSC_TRUE, cells, dim, coords, dm);CHKERRQ(ierr);
124     /* Check for flipped first cell */
125     {
126       PetscReal v0[3], J[9], invJ[9], detJ;
127 
128       ierr = DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
129       if (detJ < 0) {
130         ierr = DMPlexOrientPoint(*dm, 0, -1);CHKERRQ(ierr);
131         ierr = DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
132         if (detJ < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Something is wrong");
133       }
134     }
135     ierr = DMPlexOrient(*dm);CHKERRQ(ierr);
136     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
137     ierr = DMGetLabel(*dm, "marker", &label);CHKERRQ(ierr);
138     ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr);
139     ierr = DMPlexLabelComplete(*dm, label);CHKERRQ(ierr);
140 
141     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
142     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
143     ierr = VecRestoreArray(topology, &topo_f);CHKERRQ(ierr);
144     ierr = PetscFree(cells);CHKERRQ(ierr);
145     ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
146     ierr = VecDestroy(&topology);CHKERRQ(ierr);
147 #else
148     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Reconfigure PETSc with --download-hdf5");
149 #endif
150   }
151   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
152   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
153   PetscFunctionReturn(0);
154 }
155 
156 void mass_kernel(PetscInt dim, PetscInt Nf, PetscInt NfAux,
157            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
158            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
159            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
160 {
161   g0[0] = 1.0;
162 }
163 
164 void laplace_kernel(PetscInt dim, PetscInt Nf, PetscInt NfAux,
165            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
166            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
167            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
168 {
169   PetscInt d;
170   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
171 }
172 
173 /* data we seek to match */
174 PetscErrorCode data_kernel(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *y, void *ctx)
175 {
176   *y = 1.0/(2*PETSC_PI*PETSC_PI) * PetscSinReal(PETSC_PI*x[0]) * PetscSinReal(PETSC_PI*x[1]);
177   /* the associated control is sin(pi*x[0])*sin(pi*x[1]) */
178   return 0;
179 }
180 PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
181 {
182   *u = 0.0;
183   return 0;
184 }
185 
186 PetscErrorCode CreateCtx(DM dm, AppCtx* user)
187 {
188   PetscErrorCode ierr;
189 
190   DM             dm_mass;
191   DM             dm_laplace;
192   PetscDS        prob_mass;
193   PetscDS        prob_laplace;
194   PetscFE        fe;
195   DMLabel        label;
196   PetscSection   section;
197   PetscInt       n, k, p, d;
198   PetscInt       dof, off;
199   IS             is;
200   const PetscInt* points;
201   const PetscInt dim = 2;
202   const PetscInt id  = 1;
203   PetscErrorCode (**wtf)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
204 
205   PetscFunctionBeginUser;
206 
207   /* make the data we seek to match */
208   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, PETSC_TRUE, NULL, 4, &fe);CHKERRQ(ierr);
209 
210   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
211   ierr = DMCreateDS(dm);CHKERRQ(ierr);
212   ierr = DMCreateGlobalVector(dm, &user->data);CHKERRQ(ierr);
213 
214   /* ugh, this is hideous */
215   /* y_d = interpolate(Expression("sin(x[0]) + .."), V) */
216   ierr = PetscMalloc(1 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &wtf);CHKERRQ(ierr);
217   wtf[0] = data_kernel;
218   ierr = DMProjectFunction(dm, 0.0, wtf, NULL, INSERT_VALUES, user->data);CHKERRQ(ierr);
219   ierr = PetscFree(wtf);CHKERRQ(ierr);
220 
221   /* assemble(inner(u, v)*dx), almost */
222   ierr = DMClone(dm, &dm_mass);CHKERRQ(ierr);
223   ierr = DMCopyDisc(dm, dm_mass);CHKERRQ(ierr);
224   ierr = DMSetNumFields(dm_mass, 1);CHKERRQ(ierr);
225   ierr = DMPlexCopyCoordinates(dm, dm_mass);CHKERRQ(ierr); /* why do I have to do this separately? */
226   ierr = DMGetDS(dm_mass, &prob_mass);CHKERRQ(ierr);
227   ierr = PetscDSSetJacobian(prob_mass, 0, 0, mass_kernel, NULL, NULL, NULL);CHKERRQ(ierr);
228   ierr = PetscDSSetDiscretization(prob_mass, 0, (PetscObject) fe);CHKERRQ(ierr);
229   ierr = DMCreateMatrix(dm_mass, &user->mass);CHKERRQ(ierr);
230   ierr = DMPlexSNESComputeJacobianFEM(dm_mass, user->data, user->mass, user->mass, NULL);CHKERRQ(ierr);
231   ierr = MatSetOption(user->mass, MAT_SYMMETRIC, PETSC_TRUE);CHKERRQ(ierr);
232   ierr = DMDestroy(&dm_mass);CHKERRQ(ierr);
233 
234   /* inner(grad(u), grad(v))*dx with homogeneous Dirichlet boundary conditions */
235   ierr = DMClone(dm, &dm_laplace);CHKERRQ(ierr);
236   ierr = DMCopyDisc(dm, dm_laplace);CHKERRQ(ierr);
237   ierr = DMSetNumFields(dm_laplace, 1);CHKERRQ(ierr);
238   ierr = DMPlexCopyCoordinates(dm, dm_laplace);CHKERRQ(ierr);
239   ierr = DMGetDS(dm_laplace, &prob_laplace);CHKERRQ(ierr);
240   ierr = PetscDSSetJacobian(prob_laplace, 0, 0, NULL, NULL, NULL, laplace_kernel);CHKERRQ(ierr);
241   ierr = PetscDSSetDiscretization(prob_laplace, 0, (PetscObject) fe);CHKERRQ(ierr);
242   ierr = DMCreateMatrix(dm_laplace, &user->laplace);CHKERRQ(ierr);
243   ierr = DMPlexSNESComputeJacobianFEM(dm_laplace, user->data, user->laplace, user->laplace, NULL);CHKERRQ(ierr);
244 
245   /* Code from Matt to get the indices associated with the boundary dofs */
246   ierr = DMGetLabel(dm_laplace, "marker", &label);CHKERRQ(ierr);
247   ierr = DMAddBoundary(dm_laplace, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, NULL, NULL);CHKERRQ(ierr);
248   ierr = DMGetLocalSection(dm_laplace, &section);CHKERRQ(ierr);
249   ierr = DMLabelGetStratumSize(label, 1, &n);CHKERRQ(ierr);
250   ierr = DMLabelGetStratumIS(label, 1, &is);CHKERRQ(ierr);
251   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
252   user->num_bc_dofs = 0;
253   for (p = 0; p < n; ++p) {
254     ierr = PetscSectionGetDof(section, points[p], &dof);CHKERRQ(ierr);
255     user->num_bc_dofs += dof;
256   }
257   ierr = PetscMalloc1(user->num_bc_dofs, &user->bc_indices);CHKERRQ(ierr);
258   for (p = 0, k = 0; p < n; ++p) {
259     ierr = PetscSectionGetDof(section, points[p], &dof);CHKERRQ(ierr);
260     ierr = PetscSectionGetOffset(section, points[p], &off);CHKERRQ(ierr);
261     for (d = 0; d < dof; ++d) user->bc_indices[k++] = off+d;
262   }
263   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
264   ierr = ISDestroy(&is);CHKERRQ(ierr);
265   ierr = DMDestroy(&dm_laplace);CHKERRQ(ierr);
266 
267   /* This is how I handle boundary conditions. I can't figure out how to get
268      plex to play with the way I want to impose the BCs. This loses symmetry,
269      but not in a disastrous way. If someone can improve it, please do! */
270   ierr = MatZeroRows(user->laplace, user->num_bc_dofs, user->bc_indices, 1.0, NULL, NULL);CHKERRQ(ierr);
271   ierr = PetscCalloc1(user->num_bc_dofs, &user->bc_values);CHKERRQ(ierr);
272 
273   /* also create the KSP for solving the Laplace system */
274   ierr = KSPCreate(PETSC_COMM_WORLD, &user->ksp_laplace);CHKERRQ(ierr);
275   ierr = KSPSetOperators(user->ksp_laplace, user->laplace, user->laplace);CHKERRQ(ierr);
276   ierr = KSPSetOptionsPrefix(user->ksp_laplace, "laplace_");CHKERRQ(ierr);
277   ierr = KSPSetFromOptions(user->ksp_laplace);CHKERRQ(ierr);
278 
279   /* A bit of setting up the user context */
280   user->dm = dm;
281   ierr = VecDuplicate(user->data, &user->state);CHKERRQ(ierr);
282   ierr = VecDuplicate(user->data, &user->adjoint);CHKERRQ(ierr);
283   ierr = VecDuplicate(user->data, &user->tmp1);CHKERRQ(ierr);
284   ierr = VecDuplicate(user->data, &user->tmp2);CHKERRQ(ierr);
285 
286   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
287 
288   PetscFunctionReturn(0);
289 }
290 
291 PetscErrorCode DestroyCtx(AppCtx* user)
292 {
293   PetscErrorCode ierr;
294 
295   PetscFunctionBeginUser;
296 
297   ierr = MatDestroy(&user->mass);CHKERRQ(ierr);
298   ierr = MatDestroy(&user->laplace);CHKERRQ(ierr);
299   ierr = KSPDestroy(&user->ksp_laplace);CHKERRQ(ierr);
300   ierr = VecDestroy(&user->data);CHKERRQ(ierr);
301   ierr = VecDestroy(&user->state);CHKERRQ(ierr);
302   ierr = VecDestroy(&user->adjoint);CHKERRQ(ierr);
303   ierr = VecDestroy(&user->tmp1);CHKERRQ(ierr);
304   ierr = VecDestroy(&user->tmp2);CHKERRQ(ierr);
305   ierr = PetscFree(user->bc_indices);CHKERRQ(ierr);
306   ierr = PetscFree(user->bc_values);CHKERRQ(ierr);
307 
308   PetscFunctionReturn(0);
309 }
310 
311 PetscErrorCode ReducedFunctionGradient(Tao tao, Vec u, PetscReal* func, Vec g, void* userv)
312 {
313   PetscErrorCode ierr;
314   AppCtx* user = (AppCtx*) userv;
315   const PetscReal alpha = 1.0e-6; /* regularisation parameter */
316   PetscReal inner;
317 
318   PetscFunctionBeginUser;
319 
320   ierr = MatMult(user->mass, u, user->tmp1);CHKERRQ(ierr);
321   ierr = VecDot(u, user->tmp1, &inner);CHKERRQ(ierr);               /* regularisation contribution to */
322   *func = alpha * 0.5 * inner;                                      /* the functional                 */
323 
324   ierr = VecSet(g, 0.0);CHKERRQ(ierr);
325   ierr = VecAXPY(g, alpha, user->tmp1);CHKERRQ(ierr);               /* regularisation contribution to the gradient */
326 
327   /* Now compute the forward state. */
328   ierr = VecSetValues(user->tmp1, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES);CHKERRQ(ierr);
329   ierr = VecAssemblyBegin(user->tmp1);CHKERRQ(ierr);
330   ierr = VecAssemblyEnd(user->tmp1);CHKERRQ(ierr);
331   ierr = KSPSolve(user->ksp_laplace, user->tmp1, user->state);CHKERRQ(ierr); /* forward solve */
332 
333   /* Now compute the adjoint state also. */
334   ierr = VecCopy(user->state, user->tmp1);CHKERRQ(ierr);
335   ierr = VecAXPY(user->tmp1, -1.0, user->data);CHKERRQ(ierr);
336   ierr = MatMult(user->mass, user->tmp1, user->tmp2);CHKERRQ(ierr);
337   ierr = VecDot(user->tmp1, user->tmp2, &inner);CHKERRQ(ierr);      /* misfit contribution to */
338   *func += 0.5 * inner;                                             /* the functional         */
339 
340   ierr = VecSetValues(user->tmp2, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES);CHKERRQ(ierr);
341   ierr = VecAssemblyBegin(user->tmp2);CHKERRQ(ierr);
342   ierr = VecAssemblyEnd(user->tmp2);CHKERRQ(ierr);
343   ierr = KSPSolve(user->ksp_laplace, user->tmp2, user->adjoint);CHKERRQ(ierr); /* adjoint solve */
344 
345   /* And bring it home with the gradient. */
346   ierr = MatMult(user->mass, user->adjoint, user->tmp1);CHKERRQ(ierr);
347   ierr = VecAXPY(g, 1.0, user->tmp1);CHKERRQ(ierr);                 /* adjoint contribution to the gradient */
348 
349   PetscFunctionReturn(0);
350 }
351 
352 int main(int argc, char **argv)
353 {
354   DM             dm;
355   Tao            tao;
356   Vec            u, lb, ub;
357   AppCtx         user;
358   PetscErrorCode ierr;
359 
360   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
361   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
362   ierr = CreateCtx(dm, &user);CHKERRQ(ierr);
363 
364   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
365   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
366   ierr = VecDuplicate(u, &lb);CHKERRQ(ierr);
367   ierr = VecDuplicate(u, &ub);CHKERRQ(ierr);
368   ierr = VecSet(lb, 0.0);CHKERRQ(ierr); /* satisfied at the minimum anyway */
369   ierr = VecSet(ub, 0.8);CHKERRQ(ierr); /* a nontrivial upper bound */
370 
371   ierr = TaoCreate(PETSC_COMM_WORLD, &tao);CHKERRQ(ierr);
372   ierr = TaoSetInitialVector(tao, u);CHKERRQ(ierr);
373   ierr = TaoSetObjectiveAndGradientRoutine(tao, ReducedFunctionGradient, &user);CHKERRQ(ierr);
374   ierr = TaoSetVariableBounds(tao, lb, ub);CHKERRQ(ierr);
375   ierr = TaoSetType(tao, TAOBLMVM);CHKERRQ(ierr);
376   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
377 
378   if (user.use_riesz) {
379     ierr = TaoLMVMSetH0(tao, user.mass);CHKERRQ(ierr);       /* crucial for mesh independence */
380     ierr = TaoSetGradientNorm(tao, user.mass);CHKERRQ(ierr);
381   }
382 
383   ierr = TaoSolve(tao);CHKERRQ(ierr);
384 
385   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
386   ierr = VecViewFromOptions(u, NULL, "-sol_view");CHKERRQ(ierr);
387 
388   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
389   ierr = DMDestroy(&dm);CHKERRQ(ierr);
390   ierr = VecDestroy(&u);CHKERRQ(ierr);
391   ierr = VecDestroy(&lb);CHKERRQ(ierr);
392   ierr = VecDestroy(&ub);CHKERRQ(ierr);
393   ierr = DestroyCtx(&user);CHKERRQ(ierr);
394   ierr = PetscFinalize();
395   return ierr;
396 }
397 
398 /*TEST
399 
400     build:
401       requires: !complex !single
402 
403     test:
404       requires: hdf5 double datafilespath !defined(PETSC_USE_64BIT_INDICES) hypre
405       args: -laplace_ksp_type cg -laplace_pc_type hypre -laplace_ksp_monitor_true_residual -laplace_ksp_max_it 4 -mat_lmvm_ksp_type cg -mat_lmvm_pc_type gamg -tao_monitor -petscspace_degree 1 -tao_converged_reason -tao_gatol 1.0e-9 -dm_view hdf5:solution.h5 -sol_view hdf5:solution.h5::append -use_riesz 1 -f $DATAFILESPATH/meshes/mesh-1.h5
406       filter: sed -e "s/-nan/nan/g"
407 
408     test:
409       suffix: guess_pod
410       requires: double triangle
411       args: -laplace_ksp_type cg -laplace_pc_type gamg  -laplace_pc_gamg_esteig_ksp_type cg -laplace_ksp_converged_reason -mat_lmvm_ksp_type cg -mat_lmvm_pc_type gamg -tao_monitor -petscspace_degree 1 -tao_converged_reason -dm_refine 3 -laplace_ksp_guess_type pod -tao_gatol 1e-6 -mat_lmvm_pc_gamg_esteig_ksp_type cg
412       filter: sed -e "s/-nan/nan/g" -e "s/CONVERGED_RTOL iterations 9/CONVERGED_RTOL iterations 8/g" -e "s/CONVERGED_RTOL iterations 4/CONVERGED_RTOL iterations 3/g"
413 
414 TEST*/
415