1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
17 
18 
19 
20 // Adjoints Example 5 - SolutionHistory (both Memory and File), General Localized Vectors and Unsteady Adjoints.
21 // \author Vikram Garg
22 // \date 2020
23 //
24 // This example showcases the solution storage and retrival capabilities of libMesh, in the context of
25 // unsteady adjoints. The primary motivation is adjoint sensitivity analysis for unsteady
26 // problems. The PDE we are interested in is the simple 2-d heat
27 // equation:
28 // partial(T)/partial(t) - K Laplacian(T) = 0
29 // with initial condition:
30 // T(x,y;0) = sin(pi*x) sin(pi*y) and boundary conditions:
31 // T(boundary;t) = 0
32 
33 // For these initial and boundary conditions, the exact solution
34 // u = exp(-K 2*pi^2 t) * sin(pi*x) * sin(pi*y)
35 
36 // We specify our Quantity of Interest (QoI) as
37 // Q(u) = int_{domain} u(x,y;1) sin(pi*x) sin(pi*y) dx dy, and
38 // are interested in computing the sensitivity dQ/dK
39 
40 // The exact value of this sensitivity is:
41 // dQ/dK = int_{domain} du/dK sin(pi*x) sin(pi*y) dx dy
42 // = int_{domain} (-2*pi^2 * exp(-K pi^2) ) sin^2(pi*x) sin^2(pi*y) dx dy
43 // = (-2*pi^2 * exp(-K 2*pi^2) )/4 = -4.9022 (K = 1.0e-3)
44 
45 // For this QoI, the continuous adjoint problem reads,
46 // -partial(z)/partial(t) - K Laplacian(z) = 0
47 // with initial condition:
48 // T(x,y;1) = sin(pi*x) sin(pi*y)
49 // and boundary condition:
50 // T(boundary;t) = 0
51 
52 // which has the exact solution,
53 // z = exp(-K 2*pi^2 (1 - t)) * sin(pi*x) * sin(pi*y)
54 // which is the mirror image in time of the forward solution
55 
56 // For an adjoint consistent space-time formulation, the discrete
57 // adjoint can be obtained by marching backwards from the adjoint
58 // initial condition and solving the transpose of the discrete primal
59 // problem at the last nonlinear solve of the corresponding primal
60 // timestep. This necessitates the storage of the primal solution at
61 // all timesteps, which is accomplished here using either a
62 // MemorySolutionHistory or a FileSolutionHistory object. As the name suggests, this object
63 // simply stores the primal solution (and other vectors we may choose
64 // to save) in memory or disk, so that we can retrieve them later, whenever necessary.
65 
66 // The discrete adjoint system for implicit time steppers requires the
67 // localization of vectors other than system.solution, which is
68 // accomplished using the localize_vectors method. In this particular
69 // example, we use the localized adjoint solution to assemble the
70 // residual contribution for the current adjoint timestep from the last
71 // computed adjoint timestep.
72 
73 // Finally, The adjoint_advance_timestep method, the backwards time
74 // analog of advance_timestep prepares the time solver for solving the
75 // adjoint system, while the retrieve_timestep method retrieves the
76 // saved solutions at the current system.time, so that the adjoint
77 // sensitivity contribution for the current time can be computed.
78 
79 // Local includes
80 #include "initial.h"
81 #include "adjoint_initial.h"
82 #include "femparameters.h"
83 #include "heatsystem.h"
84 
85 // Libmesh includes
86 #include "libmesh/equation_systems.h"
87 #include "libmesh/dirichlet_boundaries.h"
88 #include "libmesh/system_norm.h"
89 #include "libmesh/numeric_vector.h"
90 #include "libmesh/auto_ptr.h" // libmesh_make_unique
91 #include "libmesh/enum_solver_package.h"
92 
93 #include "libmesh/mesh.h"
94 #include "libmesh/mesh_generation.h"
95 #include "libmesh/mesh_refinement.h"
96 
97 #include "libmesh/petsc_diff_solver.h"
98 #include "libmesh/steady_solver.h"
99 #include "libmesh/euler_solver.h"
100 #include "libmesh/euler2_solver.h"
101 #include "libmesh/newton_solver.h"
102 #include "libmesh/twostep_time_solver.h"
103 
104 #include "libmesh/getpot.h"
105 #include "libmesh/tecplot_io.h"
106 #include "libmesh/gmv_io.h"
107 #include "libmesh/exodusII_io.h"
108 
109 #include "libmesh/sensitivity_data.h"
110 
111 // SolutionHistory Includes
112 #include "libmesh/solution_history.h"
113 #include "libmesh/memory_solution_history.h"
114 #include "libmesh/file_solution_history.h"
115 
116 // C++ includes
117 #include <iostream>
118 #include <sys/time.h>
119 #include <iomanip>
120 
write_output(EquationSystems & es,unsigned int t_step,std::string solution_type,FEMParameters & param)121 void write_output(EquationSystems & es,
122                   unsigned int t_step,       // The current time step count
123                   std::string solution_type, // primal or adjoint solve
124                   FEMParameters & param)
125 {
126   // Ignore parameters when there are no output formats available.
127   libmesh_ignore(es, t_step, solution_type, param);
128 
129 #ifdef LIBMESH_HAVE_GMV
130   if (param.output_gmv)
131     {
132       MeshBase & mesh = es.get_mesh();
133 
134       std::ostringstream file_name_gmv;
135       file_name_gmv << solution_type
136                     << ".out.gmv."
137                     << std::setw(2)
138                     << std::setfill('0')
139                     << std::right
140                     << t_step;
141 
142       GMVIO(mesh).write_equation_systems(file_name_gmv.str(), es);
143     }
144 #endif
145 
146 #ifdef LIBMESH_HAVE_EXODUS_API
147   if (param.output_exodus)
148     {
149       MeshBase & mesh = es.get_mesh();
150 
151       // We write out one file per timestep. The files are named in
152       // the following way:
153       // foo.e
154       // foo.e-s002
155       // foo.e-s003
156       // ...
157       // so that, if you open the first one with Paraview, it actually
158       // opens the entire sequence of adapted files.
159       std::ostringstream file_name_exodus;
160 
161       file_name_exodus << solution_type << ".e";
162       if (t_step > 0)
163         file_name_exodus << "-s"
164                          << std::setw(3)
165                          << std::setfill('0')
166                          << std::right
167                          << t_step + 1;
168 
169       // TODO: Get the current time from the System...
170       ExodusII_IO(mesh).write_timestep(file_name_exodus.str(),
171                                        es,
172                                        1,
173                                        /*time=*/t_step + 1);
174     }
175 #endif
176 }
177 
set_system_parameters(HeatSystem & system,FEMParameters & param)178 void set_system_parameters(HeatSystem &system, FEMParameters &param)
179 {
180   // Use the prescribed FE type
181   system.fe_family() = param.fe_family[0];
182   system.fe_order() = param.fe_order[0];
183 
184   // Use analytical jacobians?
185   system.analytic_jacobians() = param.analytic_jacobians;
186 
187   // Verify analytic jacobians against numerical ones?
188   system.verify_analytic_jacobians = param.verify_analytic_jacobians;
189   system.numerical_jacobian_h = param.numerical_jacobian_h;
190 
191   // More desperate debugging options
192   system.print_solution_norms    = param.print_solution_norms;
193   system.print_solutions         = param.print_solutions;
194   system.print_residual_norms    = param.print_residual_norms;
195   system.print_residuals         = param.print_residuals;
196   system.print_jacobian_norms    = param.print_jacobian_norms;
197   system.print_jacobians         = param.print_jacobians;
198   system.print_element_jacobians = param.print_element_jacobians;
199   system.print_element_residuals = param.print_element_residuals;
200 
201   // Solve this as a time-dependent or steady system
202   if (param.transient)
203     {
204       UnsteadySolver *innersolver;
205       if (param.timesolver_core == "euler2")
206         {
207           Euler2Solver *euler2solver =
208             new Euler2Solver(system);
209 
210           euler2solver->theta = param.timesolver_theta;
211           innersolver = euler2solver;
212         }
213       else if (param.timesolver_core == "euler")
214         {
215           EulerSolver *eulersolver =
216             new EulerSolver(system);
217 
218           eulersolver->theta = param.timesolver_theta;
219           innersolver = eulersolver;
220         }
221       else
222         libmesh_error_msg("This example (and unsteady adjoints in libMesh) only support Backward Euler and explicit methods.");
223 
224       if (param.timesolver_tolerance)
225         {
226           TwostepTimeSolver *timesolver =
227             new TwostepTimeSolver(system);
228 
229           timesolver->max_growth       = param.timesolver_maxgrowth;
230           timesolver->target_tolerance = param.timesolver_tolerance;
231           timesolver->upper_tolerance  = param.timesolver_upper_tolerance;
232           timesolver->component_norm   = SystemNorm(param.timesolver_norm);
233 
234           timesolver->core_time_solver =
235             std::unique_ptr<UnsteadySolver>(innersolver);
236 
237           system.time_solver =
238             std::unique_ptr<UnsteadySolver>(timesolver);
239 
240           // The Memory/File Solution History object we will set the system SolutionHistory object to
241           if(param.solution_history_type == "memory")
242           {
243             MemorySolutionHistory heatsystem_solution_history(system);
244             system.time_solver->set_solution_history(heatsystem_solution_history);
245           }
246           else if (param.solution_history_type == "file")
247           {
248             FileSolutionHistory heatsystem_solution_history(system);
249             system.time_solver->set_solution_history(heatsystem_solution_history);
250           }
251           else
252           libmesh_error_msg("Unrecognized solution history type: " << param.solution_history_type);
253 
254         }
255       else
256       {
257         system.time_solver =
258           std::unique_ptr<TimeSolver>(innersolver);
259 
260         // The Memory/File Solution History object we will set the system SolutionHistory object to
261         if(param.solution_history_type == "memory")
262         {
263           MemorySolutionHistory heatsystem_solution_history(system);
264           system.time_solver->set_solution_history(heatsystem_solution_history);
265         }
266         else if (param.solution_history_type == "file")
267         {
268           FileSolutionHistory heatsystem_solution_history(system);
269           system.time_solver->set_solution_history(heatsystem_solution_history);
270         }
271         else
272         libmesh_error_msg("Unrecognized solution history type: " << param.solution_history_type);
273 
274         // The Memory/File Solution History object we will set the system SolutionHistory object to
275         FileSolutionHistory heatsystem_solution_history(system);
276         system.time_solver->set_solution_history(heatsystem_solution_history);
277 
278       }
279     }
280   else
281     system.time_solver =
282       std::unique_ptr<TimeSolver>(new SteadySolver(system));
283 
284   system.time_solver->reduce_deltat_on_diffsolver_failure =
285                                         param.deltat_reductions;
286   system.time_solver->quiet           = param.time_solver_quiet;
287 
288  #ifdef LIBMESH_ENABLE_DIRICHLET
289   // Create any Dirichlet boundary conditions
290   typedef
291     std::map<boundary_id_type, FunctionBase<Number> *>::
292     const_iterator Iter;
293 
294   for (Iter i = param.dirichlet_conditions.begin();
295        i != param.dirichlet_conditions.end(); ++i)
296     {
297       boundary_id_type b = i->first;
298       FunctionBase<Number> *f = i->second;
299       std::set<boundary_id_type> bdys; bdys.insert(b);
300 
301       system.get_dof_map().add_dirichlet_boundary(DirichletBoundary(bdys,
302                                                                     param.dirichlet_condition_variables[b],
303                                                                     f));
304 
305       libMesh::out << "Added Dirichlet boundary " << b << " for variables ";
306       for (std::size_t vi=0; vi != param.dirichlet_condition_variables[b].size(); ++vi)
307         libMesh::out << param.dirichlet_condition_variables[b][vi];
308       libMesh::out << std::endl;
309     }
310  #endif // LIBMESH_ENABLE_DIRICHLET
311 
312   // Set the time stepping options
313   system.deltat = param.deltat;
314 
315   // And the integration options
316   system.extra_quadrature_order = param.extra_quadrature_order;
317 
318   // And the nonlinear solver options
319   if (param.use_petsc_snes)
320     {
321 #ifdef LIBMESH_HAVE_PETSC
322       PetscDiffSolver *solver = new PetscDiffSolver(system);
323       system.time_solver->diff_solver() = std::unique_ptr<DiffSolver>(solver);
324 #else
325       libmesh_error_msg("This example requires libMesh to be compiled with PETSc support.");
326 #endif
327     }
328   else
329     {
330       NewtonSolver *solver = new NewtonSolver(system);
331       system.time_solver->diff_solver() = std::unique_ptr<DiffSolver>(solver);
332 
333       solver->quiet                       = param.solver_quiet;
334       solver->verbose                     = param.solver_verbose;
335       solver->max_nonlinear_iterations    = param.max_nonlinear_iterations;
336       solver->minsteplength               = param.min_step_length;
337       solver->relative_step_tolerance     = param.relative_step_tolerance;
338       solver->relative_residual_tolerance = param.relative_residual_tolerance;
339       solver->require_residual_reduction  = param.require_residual_reduction;
340       solver->linear_tolerance_multiplier = param.linear_tolerance_multiplier;
341       if (system.time_solver->reduce_deltat_on_diffsolver_failure)
342         {
343           solver->continue_after_max_iterations = true;
344           solver->continue_after_backtrack_failure = true;
345         }
346 
347       // And the linear solver options
348       solver->max_linear_iterations       = param.max_linear_iterations;
349       solver->initial_linear_tolerance    = param.initial_linear_tolerance;
350       solver->minimum_linear_tolerance    = param.minimum_linear_tolerance;
351     }
352 }
353 
354 // The main program.
main(int argc,char ** argv)355 int main (int argc, char ** argv)
356 {
357   // Skip adaptive examples on a non-adaptive libMesh build
358 #ifndef LIBMESH_ENABLE_AMR
359   libmesh_ignore(argc, argv);
360   libmesh_example_requires(false, "--enable-amr");
361 #else
362   // Skip this 2D example if libMesh was compiled as 1D-only.
363   libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
364 
365   // We use Dirichlet boundary conditions here
366 #ifndef LIBMESH_ENABLE_DIRICHLET
367   libmesh_example_requires(false, "--enable-dirichlet");
368 #endif
369 
370   // Initialize libMesh.
371   LibMeshInit init (argc, argv);
372 
373   // This doesn't converge with Trilinos for some reason...
374   libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
375 
376   libMesh::out << "Started " << argv[0] << std::endl;
377 
378   // Make sure the general input file exists, and parse it
379   {
380     std::ifstream i("general.in");
381     libmesh_error_msg_if(!i, '[' << init.comm().rank() << "] Can't find general.in; exiting early.");
382   }
383   GetPot infile("general.in");
384 
385   // But allow the command line to override it.
386   infile.parse_command_line(argc, argv);
387 
388   // Read in parameters from the input file
389   FEMParameters param(init.comm());
390   param.read(infile);
391 
392   // Create a mesh with the given dimension, distributed
393   // across the default MPI communicator.
394   Mesh mesh(init.comm(), cast_int<unsigned char>(param.dimension));
395 
396   // And an object to refine it
397   auto mesh_refinement = libmesh_make_unique<MeshRefinement>(mesh);
398 
399   // And an EquationSystems to run on it
400   EquationSystems equation_systems (mesh);
401 
402   libMesh::out << "Building mesh" << std::endl;
403 
404   // Build a unit square
405   ElemType elemtype;
406 
407   if (param.elementtype == "tri" ||
408       param.elementtype == "unstructured")
409     elemtype = TRI3;
410   else
411     elemtype = QUAD4;
412 
413   MeshTools::Generation::build_square (mesh, param.coarsegridx, param.coarsegridy,
414                                        param.domain_xmin, param.domain_xmin + param.domain_edge_width,
415                                        param.domain_ymin, param.domain_ymin + param.domain_edge_length,
416                                        elemtype);
417 
418   libMesh::out << "Building system" << std::endl;
419 
420   HeatSystem & system = equation_systems.add_system<HeatSystem> ("HeatSystem");
421 
422   set_system_parameters(system, param);
423 
424   libMesh::out << "Initializing systems" << std::endl;
425 
426   // Initialize the system
427   equation_systems.init ();
428 
429   // Refine the grid again if requested
430   for (unsigned int i=0; i != param.extrarefinements; ++i)
431     {
432       mesh_refinement->uniformly_refine(1);
433       equation_systems.reinit();
434     }
435 
436   libMesh::out << "Setting primal initial conditions" << std::endl;
437 
438   read_initial_parameters();
439 
440   system.project_solution(initial_value, initial_grad,
441                           equation_systems.parameters);
442 
443   // Output the H1 norm of the initial conditions
444   libMesh::out << "|U("
445                << system.time
446                << ")|= "
447                << system.calculate_norm(*system.solution, 0, H1)
448                << std::endl
449                << std::endl;
450 
451   // Add an adjoint vector, this will be computed after the forward
452   // time stepping is complete
453   //
454   // Tell the library not to save adjoint solutions during the forward
455   // solve
456   //
457   // Tell the library not to project this vector, and hence, memory/file
458   // solution history to not save it.
459   //
460   // Make this vector ghosted so we can localize it to each element
461   // later.
462   const std::string & adjoint_solution_name = "adjoint_solution0";
463   system.add_vector("adjoint_solution0", false, GHOSTED);
464 
465   // To keep the number of vectors consistent between the primal and adjoint
466   // loops, we will also pre-add the adjoint rhs vector
467   system.add_vector("adjoint_rhs0", false, GHOSTED);
468 
469   // Close up any resources initial.C needed
470   finish_initialization();
471 
472   // Plot the initial conditions
473   write_output(equation_systems, 0, "primal", param);
474 
475   // Print information about the mesh and system to the screen.
476   mesh.print_info();
477   equation_systems.print_info();
478 
479   // In optimized mode we catch any solver errors, so that we can
480   // write the proper footers before closing.  In debug mode we just
481   // let the exception throw so that gdb can grab it.
482 #ifdef NDEBUG
483   try
484     {
485 #endif
486       // Now we begin the timestep loop to compute the time-accurate
487       // solution of the equations.
488       for (unsigned int t_step=param.initial_timestep;
489            t_step != param.initial_timestep + param.n_timesteps; ++t_step)
490         {
491           // A pretty update message
492           libMesh::out << " Solving time step "
493                        << t_step
494                        << ", time = "
495                        << system.time
496                        << std::endl;
497 
498           // Solve the forward problem at time t, to obtain the solution at time t + dt
499           system.solve();
500 
501 
502           // Output the H1 norm of the computed solution
503           libMesh::out << "|U("
504                       << system.time + system.time_solver->last_complete_deltat()
505                       << ")|= "
506                       << system.calculate_norm(*system.solution, 0, H1)
507                       << std::endl;
508 
509           // Advance to the next timestep in a transient problem
510           libMesh::out << "Advancing timestep" << std::endl << std::endl;
511           system.time_solver->advance_timestep();
512 
513           // Write out this timestep
514           write_output(equation_systems, t_step+1, "primal", param);
515         }
516       // End timestep loop
517 
518       ///////////////// Now for the Adjoint Solution //////////////////////////////////////
519 
520       // Now we will solve the backwards in time adjoint problem
521       libMesh::out << std::endl << "Solving the adjoint problem" << std::endl;
522 
523       // We need to tell the library that it needs to project the adjoint, so
524       // MemorySolutionHistory knows it has to save it
525 
526       // Tell the library to project the adjoint vector, and hence, memory solution history to
527       // save it
528       system.set_vector_preservation(adjoint_solution_name, true);
529 
530       libMesh::out << "Setting adjoint initial conditions Z("
531                    << system.time
532                    << ")"
533                    <<std::endl;
534 
535       // The first thing we have to do is to apply the adjoint initial
536       // condition. The user should supply these. Here they are specified
537       // in the functions adjoint_initial_value and adjoint_initial_gradient
538       system.project_vector(adjoint_initial_value,
539                             adjoint_initial_grad,
540                             equation_systems.parameters,
541                             system.get_adjoint_solution(0));
542 
543       // Since we have specified an adjoint solution for the current
544       // time (T), set the adjoint_already_solved boolean to true, so
545       // we dont solve unnecessarily in the adjoint sensitivity method
546       system.set_adjoint_already_solved(true);
547 
548       libMesh::out << "|Z("
549                    << system.time
550                    << ")|= "
551                    << system.calculate_norm(system.get_adjoint_solution(), 0, H1)
552                    << std::endl
553                    << std::endl;
554 
555       // Call the adjoint advance timestep here to save the adjoint
556       // initial conditions to disk
557       system.time_solver->adjoint_advance_timestep();
558 
559       write_output(equation_systems, param.n_timesteps, "dual", param);
560 
561       // Now that the adjoint initial condition is set, we will start the
562       // backwards in time adjoint integration
563 
564       // For loop stepping backwards in time
565       for (unsigned int t_step=param.initial_timestep;
566            t_step != param.initial_timestep + param.n_timesteps; ++t_step)
567         {
568           //A pretty update message
569           libMesh::out << " Solving adjoint time step "
570                        << t_step
571                        << ", time = "
572                        << system.time
573                        << std::endl;
574 
575           // Output the H1 norm of the retrieved primal solution from the last call
576           // to adjoint_advance_timestep
577           libMesh::out << "|U("
578                        << system.time
579                        << ")|= "
580                        << system.calculate_norm(*system.solution, 0, H1)
581                        << std::endl;
582 
583           libMesh::out << "|U("
584                          << system.time - (system.time_solver->last_complete_deltat())/((param.timesolver_tolerance) ? 2.0 : 1.0)
585                          << ")|= "
586                          << system.calculate_norm(system.get_vector("_old_nonlinear_solution"), 0, H1)
587                          << std::endl;
588 
589           system.set_adjoint_already_solved(false);
590 
591           system.adjoint_solve();
592 
593           // Now that we have solved the adjoint, set the
594           // adjoint_already_solved boolean to true, so we dont solve
595           // unnecessarily in the error estimator
596           system.set_adjoint_already_solved(true);
597 
598           libMesh::out << "Saving adjoint and retrieving primal solutions at time t=" << system.time - system.deltat << std::endl;
599 
600           // The adjoint_advance_timestep function calls the retrieve and store
601           // function of the memory_solution_history class via the
602           // memory_solution_history object we declared earlier.  The
603           // retrieve function sets the system primal vectors to their
604           // values at the current time instance.
605           system.time_solver->adjoint_advance_timestep();
606 
607           libMesh::out << "|Z("
608                        << system.time
609                        << ")|= "
610                        << system.calculate_norm(system.get_adjoint_solution(), 0, H1)
611                        << std::endl
612                        << std::endl;
613 
614           // Get a pointer to the primal solution vector
615           NumericVector<Number> & primal_solution = *system.solution;
616 
617           // Get a pointer to the solution vector of the adjoint problem for QoI 0
618           NumericVector<Number> & dual_solution_0 = system.get_adjoint_solution(0);
619 
620           // Swap the primal and dual solutions so we can write out the adjoint solution
621           primal_solution.swap(dual_solution_0);
622 
623           write_output(equation_systems, param.n_timesteps - (t_step + 1), "dual", param);
624 
625           // Swap back
626           primal_solution.swap(dual_solution_0);
627         }
628       // End adjoint timestep loop
629 
630       // Now that we have computed both the primal and adjoint solutions, we compute the sensitivities to the parameter p
631       // dQ/dp = int_{0}^{T} partialQ/partialp - partialR/partialp(u,z;p) dt
632       // The quantity partialQ/partialp - partialR/partialp(u,z;p) is evaluated internally by the ImplicitSystem::adjoint_qoi_parameter_sensitivity function.
633       // This sensitivity evaluation is called internally by an overloaded TimeSolver::integrate_adjoint_sensitivity method which we call below.
634 
635       // Prepare the quantities we need to pass to TimeSolver::integrate_adjoint_sensitivity
636       QoISet qois;
637 
638       std::vector<unsigned int> qoi_indices;
639       qoi_indices.push_back(0);
640       qois.add_indices(qoi_indices);
641       qois.set_weight(0, 1.0);
642 
643       // A SensitivityData object
644       SensitivityData sensitivities(qois, system, system.get_parameter_vector());
645 
646       // Accumulator for the time integrated total sensitivity
647       Number total_sensitivity = 0.0;
648 
649       // Retrieve the primal and adjoint solutions at the current timestep
650       system.time_solver->retrieve_timestep();
651 
652       // A pretty update message
653       libMesh::out << "Retrieved, "
654                << "time = "
655                << system.time
656                << std::endl;
657 
658       libMesh::out << "|U("
659                    << system.time
660                    << ")|= "
661                    << system.calculate_norm(*system.solution, 0, H1)
662                    << std::endl;
663 
664       libMesh::out << "|U_old("
665                    << system.time
666                    << ")|= "
667                    << system.calculate_norm(system.get_vector("_old_nonlinear_solution"), 0, H1)
668                    << std::endl;
669 
670       libMesh::out << "|Z("
671                    << system.time
672                    << ")|= "
673                    << system.calculate_norm(system.get_adjoint_solution(0), 0, H1)
674                    << std::endl
675                    << std::endl;
676 
677       // Now we begin the timestep loop to compute the time-accurate
678       // adjoint sensitivities
679       for (unsigned int t_step=param.initial_timestep;
680            t_step != param.initial_timestep + param.n_timesteps; ++t_step)
681         {
682           // Call the postprocess function which we have overloaded to compute
683           // accumulate the perturbed residuals
684           system.time_solver->integrate_adjoint_sensitivity(qois, system.get_parameter_vector(), sensitivities);
685 
686           // A pretty update message
687           libMesh::out << "Retrieved, "
688                << "time = "
689                << system.time
690                << std::endl;
691 
692           libMesh::out << "|U("
693                    << system.time
694                    << ")|= "
695                    << system.calculate_norm(*system.solution, 0, H1)
696                    << std::endl;
697 
698           libMesh::out << "|U_old("
699                    << system.time
700                    << ")|= "
701                    << system.calculate_norm(system.get_vector("_old_nonlinear_solution"), 0, H1)
702                    << std::endl;
703 
704           libMesh::out << "|Z("
705                    << system.time
706                    << ")|= "
707                    << system.calculate_norm(system.get_adjoint_solution(0), 0, H1)
708                    << std::endl
709                    << std::endl;
710 
711           // Add the contribution of the current timestep to the total sensitivity
712           total_sensitivity += sensitivities[0][0];
713         }
714 
715       // Print it out
716       libMesh::out << "Sensitivity of QoI 0 w.r.t parameter 0 is: "
717                    << total_sensitivity
718                    << std::endl;
719 
720       // Hard coded test to ensure that the actual numbers we are
721       // getting are what they should be
722       // The 2e-4 tolerance is chosen to ensure success even with
723       // 32-bit floats
724       if(param.timesolver_tolerance)
725       {
726         libmesh_error_msg_if(std::abs(system.time - (1.0089)) >= 2.e-4,
727                              "Mismatch in end time reached by adaptive timestepper!");
728 
729         libmesh_error_msg_if(std::abs(total_sensitivity - 4.87767) >= 2.e-4,
730                              "Mismatch in sensitivity gold value!");
731       }
732       else
733       {
734         libmesh_error_msg_if(std::abs(total_sensitivity - 4.83551) >= 2.e-4,
735                              "Mismatch in sensitivity gold value!");
736       }
737 #ifdef NDEBUG
738     }
739   catch (...)
740     {
741       libMesh::err << '[' << mesh.processor_id()
742                    << "] Caught exception; exiting early." << std::endl;
743     }
744 #endif
745 
746   libMesh::err << '[' << mesh.processor_id()
747                << "] Completing output."
748                << std::endl;
749 
750   // All done.
751   return 0;
752 
753 #endif // LIBMESH_ENABLE_AMR
754 }
755