libMesh
Functions
adjoints_ex4.C File Reference

Go to the source code of this file.

Functions

void write_output (EquationSystems &es, unsigned int a_step, std::string solution_type, FEMParameters &param)
 
void set_system_parameters (LaplaceSystem &system, FEMParameters &param)
 
std::unique_ptr< MeshRefinementbuild_mesh_refinement (MeshBase &mesh, FEMParameters &param)
 
std::unique_ptr< AdjointRefinementEstimatorbuild_adjoint_refinement_error_estimator (QoISet &qois)
 
int main (int argc, char **argv)
 

Function Documentation

◆ build_adjoint_refinement_error_estimator()

std::unique_ptr<AdjointRefinementEstimator> build_adjoint_refinement_error_estimator ( QoISet qois)

Definition at line 225 of file adjoints_ex4.C.

226 {
227  libMesh::out << "Computing the error estimate using the Adjoint Refinement Error Estimator" << std::endl << std::endl;
228 
229  AdjointRefinementEstimator *adjoint_refinement_estimator = new AdjointRefinementEstimator;
230 
231  adjoint_refinement_estimator->qoi_set() = qois;
232 
233  // We enrich the FE space for the dual problem by doing 2 uniform h refinements
234  adjoint_refinement_estimator->number_h_refinements = 2;
235 
236  return std::unique_ptr<AdjointRefinementEstimator>(adjoint_refinement_estimator);
237 }

References libMesh::AdjointRefinementEstimator::number_h_refinements, libMesh::out, and libMesh::AdjointRefinementEstimator::qoi_set().

Referenced by main().

◆ build_mesh_refinement()

std::unique_ptr<MeshRefinement> build_mesh_refinement ( MeshBase mesh,
FEMParameters param 
)

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 243 of file adjoints_ex4.C.

244 {
245  // Initialize libMesh.
246  LibMeshInit init (argc, argv);
247 
248  // This example requires a linear solver package.
249  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
250  "--enable-petsc, --enable-trilinos, or --enable-eigen");
251 
252  // Skip adaptive examples on a non-adaptive libMesh build
253 #ifndef LIBMESH_ENABLE_AMR
254  libmesh_example_requires(false, "--enable-amr");
255 #else
256 
257  // This doesn't converge with Eigen BICGSTAB for some reason...
258  libmesh_example_requires(libMesh::default_solver_package() != EIGEN_SOLVERS, "--enable-petsc");
259 
260  libMesh::out << "Started " << argv[0] << std::endl;
261 
262  // Make sure the general input file exists, and parse it
263  {
264  std::ifstream i("general.in");
265  if (!i)
266  libmesh_error_msg('[' << init.comm().rank() << "] Can't find general.in; exiting early.");
267  }
268  GetPot infile("general.in");
269 
270  // Read in parameters from the input file
271  FEMParameters param(init.comm());
272  param.read(infile);
273 
274  // Skip this default-2D example if libMesh was compiled as 1D-only.
275  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
276 
277  // Create a mesh, with dimension to be overridden later, distributed
278  // across the default MPI communicator.
279  Mesh mesh(init.comm());
280 
281  // And an object to refine it
282  std::unique_ptr<MeshRefinement> mesh_refinement =
283  build_mesh_refinement(mesh, param);
284 
285  // And an EquationSystems to run on it
286  EquationSystems equation_systems (mesh);
287 
288  libMesh::out << "Reading in and building the mesh" << std::endl;
289 
290  // Read in the mesh
291  mesh.read(param.domainfile.c_str());
292  // Make all the elements of the mesh second order so we can compute
293  // with a higher order basis
295 
296  // Create a mesh refinement object to do the initial uniform refinements
297  // on the coarse grid read in from lshaped.xda
298  MeshRefinement initial_uniform_refinements(mesh);
299  initial_uniform_refinements.uniformly_refine(param.coarserefinements);
300 
301  libMesh::out << "Building system" << std::endl;
302 
303  // Build the FEMSystem
304  LaplaceSystem & system = equation_systems.add_system<LaplaceSystem> ("LaplaceSystem");
305 
306  // Set its parameters
307  set_system_parameters(system, param);
308 
309  libMesh::out << "Initializing systems" << std::endl;
310 
311  equation_systems.init ();
312 
313  // Print information about the mesh and system to the screen.
314  mesh.print_info();
315  equation_systems.print_info();
316  LinearSolver<Number> *linear_solver = system.get_linear_solver();
317 
318  {
319  // Adaptively solve the timestep
320  unsigned int a_step = 0;
321  for (; a_step != param.max_adaptivesteps; ++a_step)
322  {
323  // We can't adapt to both a tolerance and a
324  // target mesh size
325  if (param.global_tolerance != 0.)
326  libmesh_assert_equal_to (param.nelem_target, 0);
327  // If we aren't adapting to a tolerance we need a
328  // target mesh size
329  else
330  libmesh_assert_greater (param.nelem_target, 0);
331 
332  linear_solver->reuse_preconditioner(false);
333 
334  // Solve the forward problem
335  system.solve();
336 
337  // Write out the computed primal solution
338  write_output(equation_systems, a_step, "primal", param);
339 
340  // Get a pointer to the primal solution vector
341  NumericVector<Number> & primal_solution = *system.solution;
342 
343  // Declare a QoISet object, we need this object to set weights for our QoI error contributions
344  QoISet qois;
345 
346  // Declare a qoi_indices vector, each index will correspond to a QoI
347  std::vector<unsigned int> qoi_indices;
348  qoi_indices.push_back(0);
349  qoi_indices.push_back(1);
350  qois.add_indices(qoi_indices);
351 
352  // Set weights for each index, these will weight the contribution of each QoI in the final error
353  // estimate to be used for flagging elements for refinement
354  qois.set_weight(0, 0.5);
355  qois.set_weight(1, 0.5);
356 
357  // Make sure we get the contributions to the adjoint RHS from the sides
358  system.assemble_qoi_sides = true;
359 
360  // We are about to solve the adjoint system, but before we do this we see the same preconditioner
361  // flag to reuse the preconditioner from the forward solver
362  linear_solver->reuse_preconditioner(param.reuse_preconditioner);
363 
364  // Solve the adjoint system. This takes the transpose of the stiffness matrix and then
365  // solves the resulting system
366  system.adjoint_solve();
367 
368  // Now that we have solved the adjoint, set the adjoint_already_solved boolean to true, so we dont solve unnecessarily in the error estimator
369  system.set_adjoint_already_solved(true);
370 
371  // Get a pointer to the solution vector of the adjoint problem for QoI 0
372  NumericVector<Number> & dual_solution_0 = system.get_adjoint_solution(0);
373 
374  // Swap the primal and dual solutions so we can write out the adjoint solution
375  primal_solution.swap(dual_solution_0);
376  write_output(equation_systems, a_step, "adjoint_0", param);
377 
378  // Swap back
379  primal_solution.swap(dual_solution_0);
380 
381  // Get a pointer to the solution vector of the adjoint problem for QoI 0
382  NumericVector<Number> & dual_solution_1 = system.get_adjoint_solution(1);
383 
384  // Swap again
385  primal_solution.swap(dual_solution_1);
386  write_output(equation_systems, a_step, "adjoint_1", param);
387 
388  // Swap back again
389  primal_solution.swap(dual_solution_1);
390 
391  libMesh::out << "Adaptive step "
392  << a_step
393  << ", we have "
394  << mesh.n_active_elem()
395  << " active elements and "
396  << equation_systems.n_active_dofs()
397  << " active dofs."
398  << std::endl;
399 
400  // Postprocess, compute the approximate QoIs and write them out to the console
401  libMesh::out << "Postprocessing: " << std::endl;
402  system.postprocess_sides = true;
403  system.postprocess();
404  Number QoI_0_computed = system.get_QoI_value("computed", 0);
405  Number QoI_0_exact = system.get_QoI_value("exact", 0);
406  Number QoI_1_computed = system.get_QoI_value("computed", 1);
407  Number QoI_1_exact = system.get_QoI_value("exact", 1);
408 
409  libMesh::out << "The relative error in QoI 0 is "
410  << std::setprecision(17)
411  << std::abs(QoI_0_computed - QoI_0_exact) / std::abs(QoI_0_exact)
412  << std::endl;
413 
414  libMesh::out << "The relative error in QoI 1 is "
415  << std::setprecision(17)
416  << std::abs(QoI_1_computed - QoI_1_exact) / std::abs(QoI_1_exact)
417  << std::endl
418  << std::endl;
419 
420  // We will declare an error vector for passing to the adjoint refinement error estimator
421  ErrorVector QoI_elementwise_error;
422 
423  // Build an adjoint refinement error estimator object
424  std::unique_ptr<AdjointRefinementEstimator> adjoint_refinement_error_estimator =
426 
427  // Estimate the error in each element using the Adjoint Refinement estimator
428  adjoint_refinement_error_estimator->estimate_error(system, QoI_elementwise_error);
429 
430  // Print out the computed error estimate, note that we access the global error estimates
431  // using an accessor function, right now sum(QoI_elementwise_error) != global_QoI_error_estimate
432  libMesh::out << "The computed relative error in QoI 0 is "
433  << std::setprecision(17)
434  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_exact)
435  << std::endl;
436 
437  libMesh::out << "The computed relative error in QoI 1 is "
438  << std::setprecision(17)
439  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_exact)
440  << std::endl
441  << std::endl;
442 
443  // Also print out effectivity indices (estimated error/true error)
444  libMesh::out << "The effectivity index for the computed error in QoI 0 is "
445  << std::setprecision(17)
446  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact)
447  << std::endl;
448 
449  libMesh::out << "The effectivity index for the computed error in QoI 1 is "
450  << std::setprecision(17)
451  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact)
452  << std::endl
453  << std::endl;
454 
455  // For refinement purposes we need to sort by error
456  // *magnitudes*, but AdjointRefinement gives us signed errors.
457  if (!param.refine_uniformly)
458  for (std::size_t i=0; i<QoI_elementwise_error.size(); i++)
459  if (QoI_elementwise_error[i] != 0.)
460  QoI_elementwise_error[i] = std::abs(QoI_elementwise_error[i]);
461 
462  // We have to refine either based on reaching an error tolerance or
463  // a number of elements target, which should be verified above
464  // Otherwise we flag elements by error tolerance or nelem target
465 
466  // Uniform refinement
467  if (param.refine_uniformly)
468  {
469  mesh_refinement->uniformly_refine(1);
470  }
471  // Adaptively refine based on reaching an error tolerance
472  else if (param.global_tolerance >= 0. && param.nelem_target == 0.)
473  {
474  mesh_refinement->flag_elements_by_error_tolerance (QoI_elementwise_error);
475 
476  mesh_refinement->refine_and_coarsen_elements();
477  }
478  // Adaptively refine based on reaching a target number of elements
479  else
480  {
481  if (mesh.n_active_elem() >= param.nelem_target)
482  {
483  libMesh::out << "We reached the target number of elements." << std::endl << std::endl;
484  break;
485  }
486 
487  mesh_refinement->flag_elements_by_nelem_target (QoI_elementwise_error);
488 
489  mesh_refinement->refine_and_coarsen_elements();
490  }
491 
492  // Dont forget to reinit the system after each adaptive refinement !
493  equation_systems.reinit();
494 
495  libMesh::out << "Refined mesh to "
496  << mesh.n_active_elem()
497  << " active elements and "
498  << equation_systems.n_active_dofs()
499  << " active dofs."
500  << std::endl;
501  }
502 
503  // Do one last solve if necessary
504  if (a_step == param.max_adaptivesteps)
505  {
506  linear_solver->reuse_preconditioner(false);
507  system.solve();
508 
509  write_output(equation_systems, a_step, "primal", param);
510 
511  NumericVector<Number> & primal_solution = *system.solution;
512 
513  QoISet qois;
514  std::vector<unsigned int> qoi_indices;
515 
516  qoi_indices.push_back(0);
517  qoi_indices.push_back(1);
518  qois.add_indices(qoi_indices);
519 
520  qois.set_weight(0, 0.5);
521  qois.set_weight(1, 0.5);
522 
523  system.assemble_qoi_sides = true;
524  linear_solver->reuse_preconditioner(param.reuse_preconditioner);
525  system.adjoint_solve();
526 
527  // Now that we have solved the adjoint, set the adjoint_already_solved boolean to true, so we dont solve unnecessarily in the error estimator
528  system.set_adjoint_already_solved(true);
529 
530  NumericVector<Number> & dual_solution_0 = system.get_adjoint_solution(0);
531 
532  primal_solution.swap(dual_solution_0);
533  write_output(equation_systems, a_step, "adjoint_0", param);
534 
535  primal_solution.swap(dual_solution_0);
536 
537  NumericVector<Number> & dual_solution_1 = system.get_adjoint_solution(1);
538 
539  primal_solution.swap(dual_solution_1);
540  write_output(equation_systems, a_step, "adjoint_1", param);
541 
542  primal_solution.swap(dual_solution_1);
543 
544  libMesh::out << "Adaptive step "
545  << a_step
546  << ", we have "
547  << mesh.n_active_elem()
548  << " active elements and "
549  << equation_systems.n_active_dofs()
550  << " active dofs."
551  << std::endl;
552 
553  libMesh::out << "Postprocessing: " << std::endl;
554  system.postprocess_sides = true;
555  system.postprocess();
556 
557  Number QoI_0_computed = system.get_QoI_value("computed", 0);
558  Number QoI_0_exact = system.get_QoI_value("exact", 0);
559  Number QoI_1_computed = system.get_QoI_value("computed", 1);
560  Number QoI_1_exact = system.get_QoI_value("exact", 1);
561 
562  libMesh::out << "The relative error in QoI 0 is "
563  << std::setprecision(17)
564  << std::abs(QoI_0_computed - QoI_0_exact) / std::abs(QoI_0_exact)
565  << std::endl;
566 
567  libMesh::out << "The relative error in QoI 1 is "
568  << std::setprecision(17)
569  << std::abs(QoI_1_computed - QoI_1_exact) / std::abs(QoI_1_exact)
570  << std::endl
571  << std::endl;
572 
573  // We will declare an error vector for passing to the adjoint refinement error estimator
574  // Right now, only the first entry of this vector will be filled (with the global QoI error estimate)
575  // Later, each entry of the vector will contain elementwise error that the user can sum to get the total error
576  ErrorVector QoI_elementwise_error;
577 
578  // Build an adjoint refinement error estimator object
579  std::unique_ptr<AdjointRefinementEstimator> adjoint_refinement_error_estimator =
581 
582  // Estimate the error in each element using the Adjoint Refinement estimator
583  adjoint_refinement_error_estimator->estimate_error(system, QoI_elementwise_error);
584 
585  // Print out the computed error estimate, note that we access the global error estimates
586  // using an accessor function, right now sum(QoI_elementwise_error) != global_QoI_error_estimate
587  libMesh::out << "The computed relative error in QoI 0 is "
588  << std::setprecision(17)
589  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_exact)
590  << std::endl;
591 
592  libMesh::out << "The computed relative error in QoI 1 is "
593  << std::setprecision(17)
594  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_exact)
595  << std::endl
596  << std::endl;
597 
598  // Also print out effectivity indices (estimated error/true error)
599  libMesh::out << "The effectivity index for the computed error in QoI 0 is "
600  << std::setprecision(17)
601  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact)
602  << std::endl;
603 
604  libMesh::out << "The effectivity index for the computed error in QoI 1 is "
605  << std::setprecision(17)
606  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact)
607  << std::endl
608  << std::endl;
609 
610  // Hard coded assert to ensure that the actual numbers we are getting are what they should be
611 
612  // The effectivity index isn't exactly reproducible at single precision
613  // libmesh_assert_less(std::abs(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact) - 0.84010976704434637), 1.e-5);
614  // libmesh_assert_less(std::abs(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact) - 0.48294428289950514), 1.e-5);
615 
616  // But the effectivity indices should always be sane
617  libmesh_assert_less(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact), 2.5);
618  libmesh_assert_greater(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact), .4);
619  libmesh_assert_less(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact), 2.5);
620  libmesh_assert_greater(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact), .4);
621 
622  // And the computed errors should still be low
623  libmesh_assert_less(std::abs(QoI_0_computed - QoI_0_exact), 2e-4);
624  libmesh_assert_less(std::abs(QoI_1_computed - QoI_1_exact), 2e-4);
625  }
626  }
627 
628  libMesh::err << '[' << mesh.processor_id()
629  << "] Completing output."
630  << std::endl;
631 
632 #endif // #ifndef LIBMESH_ENABLE_AMR
633 
634  // All done.
635  return 0;
636 }

References std::abs(), libMesh::QoISet::add_indices(), libMesh::EquationSystems::add_system(), libMesh::DifferentiableSystem::adjoint_solve(), libMesh::MeshBase::all_second_order(), libMesh::DifferentiableQoI::assemble_qoi_sides, build_adjoint_refinement_error_estimator(), build_mesh_refinement(), libMesh::default_solver_package(), libMesh::EIGEN_SOLVERS, libMesh::err, libMesh::System::get_adjoint_solution(), libMesh::DifferentiableSystem::get_linear_solver(), LaplaceSystem::get_QoI_value(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::EquationSystems::n_active_dofs(), libMesh::MeshBase::n_active_elem(), libMesh::out, LaplaceSystem::postprocess(), libMesh::DifferentiableSystem::postprocess_sides, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::ParallelObject::processor_id(), FEMParameters::read(), libMesh::MeshBase::read(), libMesh::EquationSystems::reinit(), libMesh::LinearSolver< T >::reuse_preconditioner(), libMesh::System::set_adjoint_already_solved(), set_system_parameters(), libMesh::QoISet::set_weight(), libMesh::System::solution, libMesh::FEMSystem::solve(), libMesh::NumericVector< T >::swap(), libMesh::MeshRefinement::uniformly_refine(), and write_output().

◆ set_system_parameters()

void set_system_parameters ( LaplaceSystem system,
FEMParameters param 
)

Definition at line 155 of file adjoints_ex4.C.

157 {
158  // Use analytical jacobians?
159  system.analytic_jacobians() = param.analytic_jacobians;
160 
161  // Verify analytic jacobians against numerical ones?
163 
164  // Use the prescribed FE type
165  system.fe_family() = param.fe_family[0];
166  system.fe_order() = param.fe_order[0];
167 
168  // More desperate debugging options
170  system.print_solutions = param.print_solutions;
172  system.print_residuals = param.print_residuals;
174  system.print_jacobians = param.print_jacobians;
175 
176  // No transient time solver
177  system.time_solver = libmesh_make_unique<SteadySolver>(system);
178 
179  // Nonlinear solver options
180  {
181  NewtonSolver * solver = new NewtonSolver(system);
182  system.time_solver->diff_solver() = std::unique_ptr<DiffSolver>(solver);
183 
184  solver->quiet = param.solver_quiet;
186  solver->minsteplength = param.min_step_length;
191  if (system.time_solver->reduce_deltat_on_diffsolver_failure)
192  {
193  solver->continue_after_max_iterations = true;
194  solver->continue_after_backtrack_failure = true;
195  }
196 
197  // And the linear solver options
201  }
202 }

References LaplaceSystem::analytic_jacobians(), FEMParameters::analytic_jacobians, libMesh::DiffSolver::continue_after_backtrack_failure, libMesh::DiffSolver::continue_after_max_iterations, LaplaceSystem::fe_family(), FEMParameters::fe_family, LaplaceSystem::fe_order(), FEMParameters::fe_order, FEMParameters::initial_linear_tolerance, libMesh::DiffSolver::initial_linear_tolerance, FEMParameters::linear_tolerance_multiplier, libMesh::NewtonSolver::linear_tolerance_multiplier, FEMParameters::max_linear_iterations, libMesh::DiffSolver::max_linear_iterations, FEMParameters::max_nonlinear_iterations, libMesh::DiffSolver::max_nonlinear_iterations, FEMParameters::min_step_length, FEMParameters::minimum_linear_tolerance, libMesh::DiffSolver::minimum_linear_tolerance, libMesh::NewtonSolver::minsteplength, FEMParameters::print_jacobian_norms, libMesh::DifferentiableSystem::print_jacobian_norms, FEMParameters::print_jacobians, libMesh::DifferentiableSystem::print_jacobians, FEMParameters::print_residual_norms, libMesh::DifferentiableSystem::print_residual_norms, FEMParameters::print_residuals, libMesh::DifferentiableSystem::print_residuals, FEMParameters::print_solution_norms, libMesh::DifferentiableSystem::print_solution_norms, FEMParameters::print_solutions, libMesh::DifferentiableSystem::print_solutions, libMesh::DiffSolver::quiet, FEMParameters::relative_residual_tolerance, libMesh::DiffSolver::relative_residual_tolerance, FEMParameters::relative_step_tolerance, libMesh::DiffSolver::relative_step_tolerance, libMesh::NewtonSolver::require_residual_reduction, FEMParameters::require_residual_reduction, FEMParameters::solver_quiet, libMesh::DifferentiableSystem::time_solver, FEMParameters::verify_analytic_jacobians, and libMesh::FEMSystem::verify_analytic_jacobians.

Referenced by main().

◆ write_output()

void write_output ( EquationSystems es,
unsigned int  a_step,
std::string  solution_type,
FEMParameters param 
)

Definition at line 94 of file adjoints_ex4.C.

98 {
99  // Ignore parameters when there are no output formats available.
100  libmesh_ignore(es, a_step, solution_type, param);
101 
102 #ifdef LIBMESH_HAVE_GMV
103  if (param.output_gmv)
104  {
105  MeshBase & mesh = es.get_mesh();
106 
107  std::ostringstream file_name_gmv;
108  file_name_gmv << solution_type
109  << ".out.gmv."
110  << std::setw(2)
111  << std::setfill('0')
112  << std::right
113  << a_step;
114 
116  (file_name_gmv.str(), es);
117  }
118 #endif
119 
120 #ifdef LIBMESH_HAVE_EXODUS_API
121  if (param.output_exodus)
122  {
123  MeshBase & mesh = es.get_mesh();
124 
125  // We write out one file per adaptive step. The files are named in
126  // the following way:
127  // foo.e
128  // foo.e-s002
129  // foo.e-s003
130  // ...
131  // so that, if you open the first one with Paraview, it actually
132  // opens the entire sequence of adapted files.
133  std::ostringstream file_name_exodus;
134 
135  file_name_exodus << solution_type << ".e";
136  if (a_step > 0)
137  file_name_exodus << "-s"
138  << std::setw(3)
139  << std::setfill('0')
140  << std::right
141  << a_step + 1;
142 
143  // We write each adaptive step as a pseudo "time" step, where the
144  // time simply matches the (1-based) adaptive step we are on.
145  ExodusII_IO(mesh).write_timestep(file_name_exodus.str(),
146  es,
147  1,
148  /*time=*/a_step + 1);
149  }
150 #endif
151 }

References libMesh::EquationSystems::get_mesh(), libMesh::libmesh_ignore(), mesh, FEMParameters::output_exodus, FEMParameters::output_gmv, libMesh::MeshOutput< MT >::write_equation_systems(), and libMesh::ExodusII_IO::write_timestep().

Referenced by main().

libMesh::LinearSolver::reuse_preconditioner
virtual void reuse_preconditioner(bool)
Set the same_preconditioner flag, which indicates if we reuse the same preconditioner for subsequent ...
Definition: linear_solver.C:129
FEMParameters::print_jacobians
bool print_jacobians
Definition: femparameters.h:118
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::FEMSystem::solve
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1046
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
LaplaceSystem::postprocess
virtual void postprocess()
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
Definition: L-shaped.C:168
libMesh::DifferentiableSystem::print_jacobian_norms
bool print_jacobian_norms
Set print_jacobian_norms to true to print |J| whenever it is assembled.
Definition: diff_system.h:341
libMesh::MeshRefinement::coarsen_threshold
Real & coarsen_threshold()
The coarsen_threshold provides hysteresis in AMR/C strategies.
Definition: mesh_refinement.h:897
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::MeshBase::read
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
FEMParameters::coarsen_fraction
libMesh::Real coarsen_fraction
Definition: femparameters.h:61
libMesh::DifferentiableSystem::print_residual_norms
bool print_residual_norms
Set print_residual_norms to true to print |F| whenever it is assembled.
Definition: diff_system.h:331
FEMParameters::max_nonlinear_iterations
unsigned int max_nonlinear_iterations
Definition: femparameters.h:131
libMesh::NewtonSolver::require_residual_reduction
bool require_residual_reduction
If this is set to true, the solver is forced to test the residual after each Newton step,...
Definition: newton_solver.h:98
libMesh::MeshBase::all_second_order
virtual void all_second_order(const bool full_ordered=true)=0
Converts a (conforming, non-refined) mesh with linear elements into a mesh with second-order elements...
libMesh::NumericVector::swap
virtual void swap(NumericVector< T > &v)
Swaps the contents of this with v.
Definition: numeric_vector.h:989
FEMParameters::print_residuals
bool print_residuals
Definition: femparameters.h:118
FEMParameters::print_solution_norms
bool print_solution_norms
Definition: femparameters.h:118
FEMParameters::print_residual_norms
bool print_residual_norms
Definition: femparameters.h:118
libMesh::DifferentiableSystem::print_residuals
bool print_residuals
Set print_residuals to true to print F whenever it is assembled.
Definition: diff_system.h:336
set_system_parameters
void set_system_parameters(LaplaceSystem &system, FEMParameters &param)
Definition: adjoints_ex4.C:155
libMesh::System::get_adjoint_solution
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
Definition: system.C:957
build_mesh_refinement
std::unique_ptr< MeshRefinement > build_mesh_refinement(MeshBase &mesh, FEMParameters &param)
Definition: adjoints_ex4.C:208
libMesh::DiffSolver::continue_after_backtrack_failure
bool continue_after_backtrack_failure
Defaults to false, telling the DiffSolver to throw an error when the backtracking scheme fails to fin...
Definition: diff_solver.h:180
FEMParameters::output_gmv
bool output_gmv
Definition: femparameters.h:68
FEMParameters::refine_fraction
libMesh::Real refine_fraction
Definition: femparameters.h:61
FEMParameters::output_exodus
bool output_exodus
Definition: femparameters.h:68
libMesh::DifferentiableSystem::print_solution_norms
bool print_solution_norms
Set print_residual_norms to true to print |U| whenever it is used in an assembly() call.
Definition: diff_system.h:320
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
FEMParameters::initial_linear_tolerance
double initial_linear_tolerance
Definition: femparameters.h:134
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::ExodusII_IO::write_timestep
void write_timestep(const std::string &fname, const EquationSystems &es, const int timestep, const Real time, const std::set< std::string > *system_names=nullptr)
Writes out the solution at a specific timestep.
Definition: exodusII_io.C:1286
build_adjoint_refinement_error_estimator
std::unique_ptr< AdjointRefinementEstimator > build_adjoint_refinement_error_estimator(QoISet &qois)
Definition: adjoints_ex4.C:225
libMesh::GMVIO
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:54
libMesh::NewtonSolver::linear_tolerance_multiplier
double linear_tolerance_multiplier
The tolerance for linear solves is kept below this multiplier (which defaults to 1e-3) times the norm...
Definition: newton_solver.h:143
FEMParameters::fe_order
std::vector< unsigned int > fe_order
Definition: femparameters.h:109
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
FEMParameters::solver_quiet
bool solver_quiet
Definition: femparameters.h:128
FEMParameters::require_residual_reduction
bool require_residual_reduction
Definition: femparameters.h:128
FEMParameters::coarsen_threshold
libMesh::Real coarsen_threshold
Definition: femparameters.h:61
FEMParameters::max_linear_iterations
unsigned int max_linear_iterations
Definition: femparameters.h:131
libMesh::DifferentiableSystem::time_solver
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
Definition: diff_system.h:233
libMesh::NumericVector< Number >
LaplaceSystem::analytic_jacobians
bool & analytic_jacobians()
Definition: L-shaped.h:26
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
FEMParameters::print_solutions
bool print_solutions
Definition: femparameters.h:118
FEMParameters::fe_family
std::vector< std::string > fe_family
Definition: femparameters.h:108
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::DiffSolver::minimum_linear_tolerance
double minimum_linear_tolerance
The tolerance for linear solves is kept above this minimum.
Definition: diff_solver.h:215
libMesh::DiffSolver::relative_step_tolerance
Real relative_step_tolerance
Definition: diff_solver.h:204
libMesh::QoISet
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
LaplaceSystem
Definition: L-shaped.h:11
libMesh::DifferentiableSystem::print_jacobians
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
Definition: diff_system.h:346
libMesh::DiffSolver::continue_after_max_iterations
bool continue_after_max_iterations
Defaults to true, telling the DiffSolver to continue rather than exit when a solve has reached its ma...
Definition: diff_solver.h:174
libMesh::DifferentiableSystem::get_linear_solver
virtual LinearSolver< Number > * get_linear_solver() const override
Definition: diff_system.C:175
FEMParameters::relative_residual_tolerance
libMesh::Real relative_residual_tolerance
Definition: femparameters.h:132
libMesh::ErrorVector
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
libMesh::AdjointRefinementEstimator::number_h_refinements
unsigned char number_h_refinements
How many h refinements to perform to get the fine grid.
Definition: adjoint_refinement_estimator.h:116
libMesh::QoISet::set_weight
void set_weight(std::size_t, Real)
Set the weight for this index.
Definition: qoi_set.h:229
FEMParameters::relative_step_tolerance
libMesh::Real relative_step_tolerance
Definition: femparameters.h:132
libMesh::DifferentiableSystem::adjoint_solve
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices=QoISet()) override
This function sets the _is_adjoint boolean member of TimeSolver to true and then calls the adjoint_so...
Definition: diff_system.C:164
FEMParameters::minimum_linear_tolerance
double minimum_linear_tolerance
Definition: femparameters.h:134
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::NewtonSolver
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
Definition: newton_solver.h:46
LaplaceSystem::fe_order
unsigned int & fe_order()
Definition: L-shaped.h:25
libMesh::MeshRefinement::nelem_target
dof_id_type & nelem_target()
If nelem_target is set to a nonzero value, methods like flag_elements_by_nelem_target() will attempt ...
Definition: mesh_refinement.h:903
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::MeshOutput::write_equation_systems
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
libMesh::DiffSolver::quiet
bool quiet
The DiffSolver should not print anything to libMesh::out unless quiet is set to false; default is tru...
Definition: diff_solver.h:162
libMesh::DifferentiableQoI::assemble_qoi_sides
bool assemble_qoi_sides
If assemble_qoi_sides is true (it is false by default), the assembly loop for a quantity of interest ...
Definition: diff_qoi.h:85
libMesh::MeshRefinement::coarsen_fraction
Real & coarsen_fraction()
The coarsen_fraction sets either a desired target or a desired maximum number of elements to flag for...
Definition: mesh_refinement.h:885
libMesh::DiffSolver::relative_residual_tolerance
Real relative_residual_tolerance
Definition: diff_solver.h:192
FEMParameters::global_tolerance
libMesh::Real global_tolerance
Definition: femparameters.h:60
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::AdjointRefinementEstimator::qoi_set
QoISet & qoi_set()
Access to the QoISet (default: weight all QoIs equally) to use when computing errors.
Definition: adjoint_refinement_estimator.h:73
libMesh::MeshRefinement::refine_fraction
Real & refine_fraction()
The refine_fraction sets either a desired target or a desired maximum number of elements to flag for ...
Definition: mesh_refinement.h:879
FEMParameters::min_step_length
libMesh::Real min_step_length
Definition: femparameters.h:130
FEMParameters::read
void read(GetPot &input, const std::vector< std::string > *other_variable_names=nullptr)
Definition: femparameters.C:154
write_output
void write_output(EquationSystems &es, unsigned int a_step, std::string solution_type, FEMParameters &param)
Definition: adjoints_ex4.C:94
libMesh::FEMSystem::verify_analytic_jacobians
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
Definition: fem_system.h:209
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::AdjointRefinementEstimator
This class implements a "brute force" goal-oriented error estimator which computes an estimate of err...
Definition: adjoint_refinement_estimator.h:50
libMesh::LinearSolver< Number >
libMesh::QoISet::add_indices
void add_indices(const std::vector< unsigned int > &indices)
Add this indices to the set to be calculated.
Definition: qoi_set.C:46
libMesh::DiffSolver::max_linear_iterations
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
Definition: diff_solver.h:148
libMesh::MeshRefinement::absolute_global_tolerance
Real & absolute_global_tolerance()
If absolute_global_tolerance is set to a nonzero value, methods like flag_elements_by_global_toleranc...
Definition: mesh_refinement.h:909
libMesh::MeshRefinement::coarsen_by_parents
bool & coarsen_by_parents()
If coarsen_by_parents is true, complete groups of sibling elements (elements with the same parent) wi...
Definition: mesh_refinement.h:873
libMesh::DifferentiableSystem::postprocess_sides
bool postprocess_sides
If postprocess_sides is true (it is false by default), the postprocessing loop will loop over all sid...
Definition: diff_system.h:314
FEMParameters::analytic_jacobians
bool analytic_jacobians
Definition: femparameters.h:114
libMesh::DifferentiableSystem::print_solutions
bool print_solutions
Set print_solutions to true to print U whenever it is used in an assembly() call.
Definition: diff_system.h:326
FEMParameters::verify_analytic_jacobians
libMesh::Real verify_analytic_jacobians
Definition: femparameters.h:115
FEMParameters::print_jacobian_norms
bool print_jacobian_norms
Definition: femparameters.h:118
libMesh::err
OStreamProxy err
FEMParameters
Definition: femparameters.h:22
libMesh::DiffSolver::max_nonlinear_iterations
unsigned int max_nonlinear_iterations
The DiffSolver should exit in failure if max_nonlinear_iterations is exceeded and continue_after_max_...
Definition: diff_solver.h:156
FEMParameters::nelem_target
unsigned int nelem_target
Definition: femparameters.h:59
FEMParameters::linear_tolerance_multiplier
double linear_tolerance_multiplier
Definition: femparameters.h:134
LaplaceSystem::fe_family
std::string & fe_family()
Definition: L-shaped.h:24
libMesh::NewtonSolver::minsteplength
Real minsteplength
If the quasi-Newton step length must be reduced to below this factor to give a residual reduction,...
Definition: newton_solver.h:137
LaplaceSystem::get_QoI_value
Number & get_QoI_value(std::string type, unsigned int QoI_index)
Definition: L-shaped.h:32
libMesh::out
OStreamProxy out
libMesh::EIGEN_SOLVERS
Definition: enum_solver_package.h:40
libMesh::MeshBase::n_active_elem
virtual dof_id_type n_active_elem() const =0
libMesh::DiffSolver::initial_linear_tolerance
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the DiffSolver may tighten the ...
Definition: diff_solver.h:210
libMesh::System::set_adjoint_already_solved
void set_adjoint_already_solved(bool setting)
Setter for the adjoint_already_solved boolean.
Definition: system.h:402