libMesh
Functions
adjoints_ex2.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< ErrorEstimatorbuild_error_estimator (FEMParameters &param)
 
int main (int argc, char **argv)
 

Function Documentation

◆ build_error_estimator()

std::unique_ptr<ErrorEstimator> build_error_estimator ( FEMParameters param)

Definition at line 241 of file adjoints_ex2.C.

242 {
243  if (param.indicator_type == "kelly")
244  {
245  libMesh::out << "Using Kelly Error Estimator" << std::endl;
246 
247  return libmesh_make_unique<KellyErrorEstimator>();
248  }
249  else if (param.indicator_type == "adjoint_residual")
250  {
251  libMesh::out << "Using Adjoint Residual Error Estimator with Patch Recovery Weights" << std::endl << std::endl;
252 
253  AdjointResidualErrorEstimator * adjoint_residual_estimator = new AdjointResidualErrorEstimator;
254 
255  adjoint_residual_estimator->error_plot_suffix = "error.gmv";
256 
258  adjoint_residual_estimator->primal_error_estimator().reset(p1);
259 
261  adjoint_residual_estimator->dual_error_estimator().reset(p2);
262 
263  adjoint_residual_estimator->primal_error_estimator()->error_norm.set_type(0, H1_SEMINORM);
264 
265  adjoint_residual_estimator->dual_error_estimator()->error_norm.set_type(0, H1_SEMINORM);
266 
267  return std::unique_ptr<ErrorEstimator>(adjoint_residual_estimator);
268  }
269  else
270  libmesh_error_msg("Unknown indicator_type = " << param.indicator_type);
271 }

References libMesh::AdjointResidualErrorEstimator::dual_error_estimator(), libMesh::AdjointResidualErrorEstimator::error_plot_suffix, libMesh::H1_SEMINORM, FEMParameters::indicator_type, libMesh::out, and libMesh::AdjointResidualErrorEstimator::primal_error_estimator().

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 274 of file adjoints_ex2.C.

275 {
276  // Initialize libMesh.
277  LibMeshInit init (argc, argv);
278 
279  // This example requires a linear solver package.
280  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
281  "--enable-petsc, --enable-trilinos, or --enable-eigen");
282 
283  // Skip adaptive examples on a non-adaptive libMesh build
284 #ifndef LIBMESH_ENABLE_AMR
285  libmesh_example_requires(false, "--enable-amr");
286 #else
287 
288  libMesh::out << "Started " << argv[0] << std::endl;
289 
290  // Make sure the general input file exists, and parse it
291  {
292  std::ifstream i("general.in");
293  if (!i)
294  libmesh_error_msg('[' << init.comm().rank() << "] Can't find general.in; exiting early.");
295  }
296  GetPot infile("general.in");
297 
298  // Read in parameters from the input file
299  FEMParameters param(init.comm());
300  param.read(infile);
301 
302  // Skip this default-2D example if libMesh was compiled as 1D-only.
303  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
304 
305  // Create a mesh, with dimension to be overridden later, distributed
306  // across the default MPI communicator.
307  Mesh mesh(init.comm());
308 
309  // And an object to refine it
310  std::unique_ptr<MeshRefinement> mesh_refinement =
311  build_mesh_refinement(mesh, param);
312 
313  // And an EquationSystems to run on it
314  EquationSystems equation_systems (mesh);
315 
316  libMesh::out << "Reading in and building the mesh" << std::endl;
317 
318  // Read in the mesh
319  mesh.read(param.domainfile.c_str());
320  // Make all the elements of the mesh second order so we can compute
321  // with a higher order basis
323 
324  // Create a mesh refinement object to do the initial uniform refinements
325  // on the coarse grid read in from lshaped.xda
326  MeshRefinement initial_uniform_refinements(mesh);
327  initial_uniform_refinements.uniformly_refine(param.coarserefinements);
328 
329  libMesh::out << "Building system" << std::endl;
330 
331  // Build the FEMSystem
332  LaplaceSystem & system = equation_systems.add_system<LaplaceSystem> ("LaplaceSystem");
333 
334  QoISet qois;
335 
336  std::vector<unsigned int> qoi_indices;
337  qoi_indices.push_back(0);
338  qois.add_indices(qoi_indices);
339 
340  qois.set_weight(0, 0.5);
341 
342  // Put some scope here to test that the cloning is working right
343  {
344  LaplaceQoI qoi;
345  system.attach_qoi(&qoi);
346  }
347 
348  // Set its parameters
349  set_system_parameters(system, param);
350 
351  libMesh::out << "Initializing systems" << std::endl;
352 
353  equation_systems.init ();
354 
355  // Print information about the mesh and system to the screen.
356  mesh.print_info();
357  equation_systems.print_info();
358 
359  {
360  // Adaptively solve the timestep
361  unsigned int a_step = 0;
362  for (; a_step != param.max_adaptivesteps; ++a_step)
363  {
364  // We can't adapt to both a tolerance and a
365  // target mesh size
366  if (param.global_tolerance != 0.)
367  libmesh_assert_equal_to (param.nelem_target, 0);
368  // If we aren't adapting to a tolerance we need a
369  // target mesh size
370  else
371  libmesh_assert_greater (param.nelem_target, 0);
372 
373  // Solve the forward problem
374  system.solve();
375 
376  // Write out the computed primal solution
377  write_output(equation_systems, a_step, "primal", param);
378 
379  // Get a pointer to the primal solution vector
380  NumericVector<Number> & primal_solution = *system.solution;
381 
382  // A SensitivityData object to hold the qois and parameters
383  SensitivityData sensitivities(qois, system, system.get_parameter_vector());
384 
385  // Make sure we get the contributions to the adjoint RHS from the sides
386  system.assemble_qoi_sides = true;
387 
388  // Here we solve the adjoint problem inside the adjoint_qoi_parameter_sensitivity
389  // function, so we have to set the adjoint_already_solved boolean to false
390  system.set_adjoint_already_solved(false);
391 
392  // Compute the sensitivities
393  system.adjoint_qoi_parameter_sensitivity(qois, system.get_parameter_vector(), sensitivities);
394 
395  // Now that we have solved the adjoint, set the adjoint_already_solved boolean to true, so we dont solve unnecessarily in the error estimator
396  system.set_adjoint_already_solved(true);
397 
398  GetPot infile_l_shaped("l-shaped.in");
399 
400  Number sensitivity_QoI_0_0_computed = sensitivities[0][0];
401  Number sensitivity_QoI_0_0_exact = infile_l_shaped("sensitivity_0_0", 0.0);
402  Number sensitivity_QoI_0_1_computed = sensitivities[0][1];
403  Number sensitivity_QoI_0_1_exact = infile_l_shaped("sensitivity_0_1", 0.0);
404 
405  libMesh::out << "Adaptive step "
406  << a_step
407  << ", we have "
408  << mesh.n_active_elem()
409  << " active elements and "
410  << equation_systems.n_active_dofs()
411  << " active dofs."
412  << std::endl;
413 
414  libMesh::out << "Sensitivity of QoI one to Parameter one is "
415  << sensitivity_QoI_0_0_computed
416  << std::endl;
417  libMesh::out << "Sensitivity of QoI one to Parameter two is "
418  << sensitivity_QoI_0_1_computed
419  << std::endl;
420 
421  libMesh::out << "The relative error in sensitivity QoI_0_0 is "
422  << std::setprecision(17)
423  << std::abs(sensitivity_QoI_0_0_computed - sensitivity_QoI_0_0_exact) / std::abs(sensitivity_QoI_0_0_exact)
424  << std::endl;
425 
426  libMesh::out << "The relative error in sensitivity QoI_0_1 is "
427  << std::setprecision(17)
428  << std::abs(sensitivity_QoI_0_1_computed - sensitivity_QoI_0_1_exact) / std::abs(sensitivity_QoI_0_1_exact)
429  << std::endl
430  << std::endl;
431 
432  // Get a pointer to the solution vector of the adjoint problem for QoI 0
433  NumericVector<Number> & dual_solution_0 = system.get_adjoint_solution(0);
434 
435  // Swap the primal and dual solutions so we can write out the adjoint solution
436  primal_solution.swap(dual_solution_0);
437  write_output(equation_systems, a_step, "adjoint_0", param);
438 
439  // Swap back
440  primal_solution.swap(dual_solution_0);
441 
442  // We have to refine either based on reaching an error tolerance or
443  // a number of elements target, which should be verified above
444  // Otherwise we flag elements by error tolerance or nelem target
445 
446  // Uniform refinement
447  if (param.refine_uniformly)
448  {
449  libMesh::out << "Refining Uniformly" << std::endl << std::endl;
450 
451  mesh_refinement->uniformly_refine(1);
452  }
453  // Adaptively refine based on reaching an error tolerance
454  else if (param.global_tolerance >= 0. && param.nelem_target == 0.)
455  {
456  // Now we construct the data structures for the mesh refinement process
457  ErrorVector error;
458 
459  // Build an error estimator object
460  std::unique_ptr<ErrorEstimator> error_estimator =
461  build_error_estimator(param);
462 
463  // Estimate the error in each element using the Adjoint Residual or Kelly error estimator
464  error_estimator->estimate_error(system, error);
465 
466  mesh_refinement->flag_elements_by_error_tolerance (error);
467 
468  mesh_refinement->refine_and_coarsen_elements();
469  }
470  // Adaptively refine based on reaching a target number of elements
471  else
472  {
473  // Now we construct the data structures for the mesh refinement process
474  ErrorVector error;
475 
476  // Build an error estimator object
477  std::unique_ptr<ErrorEstimator> error_estimator =
478  build_error_estimator(param);
479 
480  // Estimate the error in each element using the Adjoint Residual or Kelly error estimator
481  error_estimator->estimate_error(system, error);
482 
483  if (mesh.n_active_elem() >= param.nelem_target)
484  {
485  libMesh::out<<"We reached the target number of elements."<<std::endl <<std::endl;
486  break;
487  }
488 
489  mesh_refinement->flag_elements_by_nelem_target (error);
490 
491  mesh_refinement->refine_and_coarsen_elements();
492  }
493 
494  // Dont forget to reinit the system after each adaptive refinement !
495  equation_systems.reinit();
496 
497  libMesh::out << "Refined mesh to "
498  << mesh.n_active_elem()
499  << " active elements and "
500  << equation_systems.n_active_dofs()
501  << " active dofs."
502  << std::endl;
503  }
504 
505  // Do one last solve if necessary
506  if (a_step == param.max_adaptivesteps)
507  {
508  system.solve();
509 
510  write_output(equation_systems, a_step, "primal", param);
511 
512  NumericVector<Number> & primal_solution = *system.solution;
513 
514  SensitivityData sensitivities(qois, system, system.get_parameter_vector());
515 
516  system.assemble_qoi_sides = true;
517 
518  // Here we solve the adjoint problem inside the adjoint_qoi_parameter_sensitivity
519  // function, so we have to set the adjoint_already_solved boolean to false
520  system.set_adjoint_already_solved(false);
521 
522  system.adjoint_qoi_parameter_sensitivity(qois, system.get_parameter_vector(), sensitivities);
523 
524  // Now that we have solved the adjoint, set the adjoint_already_solved boolean to true, so we dont solve unnecessarily in the error estimator
525  system.set_adjoint_already_solved(true);
526 
527  GetPot infile_l_shaped("l-shaped.in");
528 
529  Number sensitivity_QoI_0_0_computed = sensitivities[0][0];
530  Number sensitivity_QoI_0_0_exact = infile_l_shaped("sensitivity_0_0", 0.0);
531  Number sensitivity_QoI_0_1_computed = sensitivities[0][1];
532  Number sensitivity_QoI_0_1_exact = infile_l_shaped("sensitivity_0_1", 0.0);
533 
534  libMesh::out << "Adaptive step "
535  << a_step
536  << ", we have "
537  << mesh.n_active_elem()
538  << " active elements and "
539  << equation_systems.n_active_dofs()
540  << " active dofs."
541  << std::endl;
542 
543  libMesh::out << "Sensitivity of QoI one to Parameter one is "
544  << sensitivity_QoI_0_0_computed
545  << std::endl;
546 
547  libMesh::out << "Sensitivity of QoI one to Parameter two is "
548  << sensitivity_QoI_0_1_computed
549  << std::endl;
550 
551  libMesh::out << "The error in sensitivity QoI_0_0 is "
552  << std::setprecision(17)
553  << std::abs(sensitivity_QoI_0_0_computed - sensitivity_QoI_0_0_exact)/sensitivity_QoI_0_0_exact
554  << std::endl;
555 
556  libMesh::out << "The error in sensitivity QoI_0_1 is "
557  << std::setprecision(17)
558  << std::abs(sensitivity_QoI_0_1_computed - sensitivity_QoI_0_1_exact)/sensitivity_QoI_0_1_exact
559  << std::endl
560  << std::endl;
561 
562  // Hard coded asserts to ensure that the actual numbers we are getting are what they should be
563  libmesh_assert_less(std::abs((sensitivity_QoI_0_0_computed - sensitivity_QoI_0_0_exact)/sensitivity_QoI_0_0_exact), 2.e-4);
564  libmesh_assert_less(std::abs((sensitivity_QoI_0_1_computed - sensitivity_QoI_0_1_exact)/sensitivity_QoI_0_1_exact), 2.e-4);
565 
566  NumericVector<Number> & dual_solution_0 = system.get_adjoint_solution(0);
567 
568  primal_solution.swap(dual_solution_0);
569  write_output(equation_systems, a_step, "adjoint_0", param);
570 
571  primal_solution.swap(dual_solution_0);
572  }
573  }
574 
575  libMesh::err << '[' << system.processor_id()
576  << "] Completing output."
577  << std::endl;
578 
579 #endif // #ifndef LIBMESH_ENABLE_AMR
580 
581  // All done.
582  return 0;
583 }

References std::abs(), libMesh::QoISet::add_indices(), libMesh::EquationSystems::add_system(), libMesh::ImplicitSystem::adjoint_qoi_parameter_sensitivity(), libMesh::MeshBase::all_second_order(), libMesh::DifferentiableQoI::assemble_qoi_sides, libMesh::DifferentiableSystem::attach_qoi(), build_error_estimator(), build_mesh_refinement(), libMesh::default_solver_package(), libMesh::err, libMesh::System::get_adjoint_solution(), LaplaceSystem::get_parameter_vector(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::EquationSystems::n_active_dofs(), libMesh::MeshBase::n_active_elem(), libMesh::out, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::ParallelObject::processor_id(), FEMParameters::read(), libMesh::MeshBase::read(), libMesh::EquationSystems::reinit(), libMesh::System::set_adjoint_already_solved(), set_system_parameters(), 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 166 of file adjoints_ex2.C.

167 {
168  // Use analytical jacobians?
169  system.analytic_jacobians() = param.analytic_jacobians;
170 
171  // Verify analytic jacobians against numerical ones?
173 
174  // Use the prescribed FE type
175  system.fe_family() = param.fe_family[0];
176  system.fe_order() = param.fe_order[0];
177 
178  // More desperate debugging options
180  system.print_solutions = param.print_solutions;
182  system.print_residuals = param.print_residuals;
184  system.print_jacobians = param.print_jacobians;
185 
186  // No transient time solver
187  system.time_solver = libmesh_make_unique<SteadySolver>(system);
188 
189  // Nonlinear solver options
190  {
191  NewtonSolver * solver = new NewtonSolver(system);
192  system.time_solver->diff_solver() = std::unique_ptr<DiffSolver>(solver);
193 
194  solver->quiet = param.solver_quiet;
196  solver->minsteplength = param.min_step_length;
201  if (system.time_solver->reduce_deltat_on_diffsolver_failure)
202  {
203  solver->continue_after_max_iterations = true;
204  solver->continue_after_backtrack_failure = true;
205  }
206 
207  // And the linear solver options
211  }
212 }

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 107 of file adjoints_ex2.C.

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

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().

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
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
libMesh::AdjointResidualErrorEstimator::error_plot_suffix
std::string error_plot_suffix
To aid in investigating error estimator behavior, set this string to a suffix with which to plot (pre...
Definition: adjoint_residual_error_estimator.h:104
libMesh::PatchRecoveryErrorEstimator
This class implements the Patch Recovery error indicator.
Definition: patch_recovery_error_estimator.h:55
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
libMesh::System::get_adjoint_solution
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
Definition: system.C:957
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
libMesh::H1_SEMINORM
Definition: enum_norm_type.h:43
build_error_estimator
std::unique_ptr< ErrorEstimator > build_error_estimator(FEMParameters &param)
Definition: adjoints_ex2.C:241
FEMParameters::output_gmv
bool output_gmv
Definition: femparameters.h:68
libMesh::DifferentiableSystem::attach_qoi
void attach_qoi(DifferentiableQoI *qoi_in)
Attach external QoI object.
Definition: diff_system.h:224
FEMParameters::refine_fraction
libMesh::Real refine_fraction
Definition: femparameters.h:61
FEMParameters::output_exodus
bool output_exodus
Definition: femparameters.h:68
write_output
void write_output(EquationSystems &es, unsigned int a_step, std::string solution_type, FEMParameters &param)
Definition: adjoints_ex2.C:107
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::SensitivityData
Data structure for holding completed parameter sensitivity calculations.
Definition: sensitivity_data.h:46
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
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 >
FEMParameters::indicator_type
std::string indicator_type
Definition: femparameters.h:145
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
LaplaceQoI
Definition: L-qoi.h:17
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::ImplicitSystem::adjoint_qoi_parameter_sensitivity
virtual void adjoint_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector &parameters, SensitivityData &sensitivities) override
Solves for the derivative of each of the system's quantities of interest q in qoi[qoi_indices] with r...
Definition: implicit_system.C:686
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
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
FEMParameters::relative_step_tolerance
libMesh::Real relative_step_tolerance
Definition: femparameters.h:132
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
set_system_parameters
void set_system_parameters(LaplaceSystem &system, FEMParameters &param)
Definition: adjoints_ex2.C:166
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::AdjointResidualErrorEstimator::dual_error_estimator
std::unique_ptr< ErrorEstimator > & dual_error_estimator()
Access to the "subestimator" (default: PatchRecovery) to use on the dual/adjoint solution.
Definition: adjoint_residual_error_estimator.h:83
build_mesh_refinement
std::unique_ptr< MeshRefinement > build_mesh_refinement(MeshBase &mesh, FEMParameters &param)
Definition: adjoints_ex2.C:218
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::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
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::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
LaplaceSystem::get_parameter_vector
ParameterVector & get_parameter_vector()
Definition: L-shaped.h:36
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
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
libMesh::out
OStreamProxy out
libMesh::AdjointResidualErrorEstimator
This class implements a goal oriented error indicator, by weighting residual-based estimates from the...
Definition: adjoint_residual_error_estimator.h:49
libMesh::AdjointResidualErrorEstimator::primal_error_estimator
std::unique_ptr< ErrorEstimator > & primal_error_estimator()
Access to the "subestimator" (default: PatchRecovery) to use on the primal/forward solution.
Definition: adjoint_residual_error_estimator.h:77
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