Go to the documentation of this file.
65 #include "libmesh/equation_systems.h"
66 #include "libmesh/error_vector.h"
67 #include "libmesh/mesh.h"
68 #include "libmesh/mesh_refinement.h"
69 #include "libmesh/newton_solver.h"
70 #include "libmesh/numeric_vector.h"
71 #include "libmesh/steady_solver.h"
72 #include "libmesh/system_norm.h"
73 #include "libmesh/auto_ptr.h"
74 #include "libmesh/enum_solver_package.h"
77 #include "libmesh/parameter_vector.h"
78 #include "libmesh/sensitivity_data.h"
81 #include "libmesh/kelly_error_estimator.h"
82 #include "libmesh/patch_recovery_error_estimator.h"
85 #include "libmesh/adjoint_residual_error_estimator.h"
86 #include "libmesh/qoi_set.h"
89 #include "libmesh/getpot.h"
90 #include "libmesh/gmv_io.h"
91 #include "libmesh/exodusII_io.h"
94 #include "femparameters.h"
109 std::string solution_type,
115 #ifdef LIBMESH_HAVE_GMV
120 std::ostringstream file_name_gmv;
121 file_name_gmv << solution_type
132 #ifdef LIBMESH_HAVE_EXODUS_API
145 std::ostringstream file_name_exodus;
147 file_name_exodus << solution_type <<
".e";
149 file_name_exodus <<
"-s"
187 system.
time_solver = libmesh_make_unique<SteadySolver>(system);
192 system.
time_solver->diff_solver() = std::unique_ptr<DiffSolver>(solver);
201 if (system.
time_solver->reduce_deltat_on_diffsolver_failure)
216 #ifdef LIBMESH_ENABLE_AMR
229 return std::unique_ptr<MeshRefinement>(mesh_refinement);
232 #endif // LIBMESH_ENABLE_AMR
245 libMesh::out <<
"Using Kelly Error Estimator" << std::endl;
247 return libmesh_make_unique<KellyErrorEstimator>();
251 libMesh::out <<
"Using Adjoint Residual Error Estimator with Patch Recovery Weights" << std::endl << std::endl;
267 return std::unique_ptr<ErrorEstimator>(adjoint_residual_estimator);
270 libmesh_error_msg(
"Unknown indicator_type = " << param.
indicator_type);
274 int main (
int argc,
char ** argv)
281 "--enable-petsc, --enable-trilinos, or --enable-eigen");
284 #ifndef LIBMESH_ENABLE_AMR
285 libmesh_example_requires(
false,
"--enable-amr");
292 std::ifstream i(
"general.in");
294 libmesh_error_msg(
'[' <<
init.comm().rank() <<
"] Can't find general.in; exiting early.");
296 GetPot infile(
"general.in");
303 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
310 std::unique_ptr<MeshRefinement> mesh_refinement =
316 libMesh::out <<
"Reading in and building the mesh" << std::endl;
319 mesh.
read(param.domainfile.c_str());
336 std::vector<unsigned int> qoi_indices;
337 qoi_indices.push_back(0);
340 qois.set_weight(0, 0.5);
353 equation_systems.
init ();
361 unsigned int a_step = 0;
362 for (; a_step != param.max_adaptivesteps; ++a_step)
366 if (param.global_tolerance != 0.)
367 libmesh_assert_equal_to (param.nelem_target, 0);
371 libmesh_assert_greater (param.nelem_target, 0);
377 write_output(equation_systems, a_step,
"primal", param);
398 GetPot infile_l_shaped(
"l-shaped.in");
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);
409 <<
" active elements and "
414 libMesh::out <<
"Sensitivity of QoI one to Parameter one is "
415 << sensitivity_QoI_0_0_computed
417 libMesh::out <<
"Sensitivity of QoI one to Parameter two is "
418 << sensitivity_QoI_0_1_computed
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)
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)
436 primal_solution.
swap(dual_solution_0);
437 write_output(equation_systems, a_step,
"adjoint_0", param);
440 primal_solution.
swap(dual_solution_0);
447 if (param.refine_uniformly)
449 libMesh::out <<
"Refining Uniformly" << std::endl << std::endl;
451 mesh_refinement->uniformly_refine(1);
454 else if (param.global_tolerance >= 0. && param.nelem_target == 0.)
460 std::unique_ptr<ErrorEstimator> error_estimator =
464 error_estimator->estimate_error(system, error);
466 mesh_refinement->flag_elements_by_error_tolerance (error);
468 mesh_refinement->refine_and_coarsen_elements();
477 std::unique_ptr<ErrorEstimator> error_estimator =
481 error_estimator->estimate_error(system, error);
485 libMesh::out<<
"We reached the target number of elements."<<std::endl <<std::endl;
489 mesh_refinement->flag_elements_by_nelem_target (error);
491 mesh_refinement->refine_and_coarsen_elements();
495 equation_systems.
reinit();
499 <<
" active elements and "
506 if (a_step == param.max_adaptivesteps)
510 write_output(equation_systems, a_step,
"primal", param);
527 GetPot infile_l_shaped(
"l-shaped.in");
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);
538 <<
" active elements and "
543 libMesh::out <<
"Sensitivity of QoI one to Parameter one is "
544 << sensitivity_QoI_0_0_computed
547 libMesh::out <<
"Sensitivity of QoI one to Parameter two is "
548 << sensitivity_QoI_0_1_computed
552 << std::setprecision(17)
553 <<
std::abs(sensitivity_QoI_0_0_computed - sensitivity_QoI_0_0_exact)/sensitivity_QoI_0_0_exact
557 << std::setprecision(17)
558 <<
std::abs(sensitivity_QoI_0_1_computed - sensitivity_QoI_0_1_exact)/sensitivity_QoI_0_1_exact
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);
568 primal_solution.
swap(dual_solution_0);
569 write_output(equation_systems, a_step,
"adjoint_0", param);
571 primal_solution.
swap(dual_solution_0);
576 <<
"] Completing output."
579 #endif // #ifndef LIBMESH_ENABLE_AMR
virtual void solve() override
Invokes the solver associated with the system.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
bool print_jacobian_norms
Set print_jacobian_norms to true to print |J| whenever it is assembled.
Real & coarsen_threshold()
The coarsen_threshold provides hysteresis in AMR/C strategies.
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
int main(int argc, char **argv)
const MeshBase & get_mesh() const
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.
libMesh::Real coarsen_fraction
bool print_residual_norms
Set print_residual_norms to true to print |F| whenever it is assembled.
std::string error_plot_suffix
To aid in investigating error estimator behavior, set this string to a suffix with which to plot (pre...
This class implements the Patch Recovery error indicator.
unsigned int max_nonlinear_iterations
bool require_residual_reduction
If this is set to true, the solver is forced to test the residual after each Newton step,...
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...
virtual void swap(NumericVector< T > &v)
Swaps the contents of this with v.
bool print_solution_norms
bool print_residual_norms
bool print_residuals
Set print_residuals to true to print F whenever it is assembled.
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
The libMesh namespace provides an interface to certain functionality in the library.
bool continue_after_backtrack_failure
Defaults to false, telling the DiffSolver to throw an error when the backtracking scheme fails to fin...
std::unique_ptr< ErrorEstimator > build_error_estimator(FEMParameters ¶m)
void attach_qoi(DifferentiableQoI *qoi_in)
Attach external QoI object.
libMesh::Real refine_fraction
void write_output(EquationSystems &es, unsigned int a_step, std::string solution_type, FEMParameters ¶m)
bool print_solution_norms
Set print_residual_norms to true to print |U| whenever it is used in an assembly() call.
Data structure for holding completed parameter sensitivity calculations.
SolverPackage default_solver_package()
double initial_linear_tolerance
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.
This class implements writing meshes in the GMV format.
double linear_tolerance_multiplier
The tolerance for linear solves is kept below this multiplier (which defaults to 1e-3) times the norm...
std::vector< unsigned int > fe_order
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
bool require_residual_reduction
libMesh::Real coarsen_threshold
unsigned int max_linear_iterations
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
std::string indicator_type
bool & analytic_jacobians()
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
std::vector< std::string > fe_family
This is the MeshBase class.
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
std::size_t n_active_dofs() const
double minimum_linear_tolerance
The tolerance for linear solves is kept above this minimum.
Real relative_step_tolerance
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
processor_id_type processor_id() const
void libmesh_ignore(const Args &...)
virtual void adjoint_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector ¶meters, SensitivityData &sensitivities) override
Solves for the derivative of each of the system's quantities of interest q in qoi[qoi_indices] with r...
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
virtual void init()
Initialize all the systems.
bool continue_after_max_iterations
Defaults to true, telling the DiffSolver to continue rather than exit when a solve has reached its ma...
libMesh::Real relative_residual_tolerance
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
libMesh::Real relative_step_tolerance
double minimum_linear_tolerance
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
unsigned int & fe_order()
dof_id_type & nelem_target()
If nelem_target is set to a nonzero value, methods like flag_elements_by_nelem_target() will attempt ...
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
This is the EquationSystems class.
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 ...
bool quiet
The DiffSolver should not print anything to libMesh::out unless quiet is set to false; default is tru...
void set_system_parameters(LaplaceSystem &system, FEMParameters ¶m)
bool assemble_qoi_sides
If assemble_qoi_sides is true (it is false by default), the assembly loop for a quantity of interest ...
std::unique_ptr< ErrorEstimator > & dual_error_estimator()
Access to the "subestimator" (default: PatchRecovery) to use on the dual/adjoint solution.
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
std::unique_ptr< MeshRefinement > build_mesh_refinement(MeshBase &mesh, FEMParameters ¶m)
Real & coarsen_fraction()
The coarsen_fraction sets either a desired target or a desired maximum number of elements to flag for...
Real relative_residual_tolerance
libMesh::Real global_tolerance
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Real & refine_fraction()
The refine_fraction sets either a desired target or a desired maximum number of elements to flag for ...
libMesh::Real min_step_length
void read(GetPot &input, const std::vector< std::string > *other_variable_names=nullptr)
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
void add_indices(const std::vector< unsigned int > &indices)
Add this indices to the set to be calculated.
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
Real & absolute_global_tolerance()
If absolute_global_tolerance is set to a nonzero value, methods like flag_elements_by_global_toleranc...
ParameterVector & get_parameter_vector()
bool & coarsen_by_parents()
If coarsen_by_parents is true, complete groups of sibling elements (elements with the same parent) wi...
bool print_solutions
Set print_solutions to true to print U whenever it is used in an assembly() call.
libMesh::Real verify_analytic_jacobians
bool print_jacobian_norms
unsigned int max_nonlinear_iterations
The DiffSolver should exit in failure if max_nonlinear_iterations is exceeded and continue_after_max_...
unsigned int nelem_target
double linear_tolerance_multiplier
std::string & fe_family()
Real minsteplength
If the quasi-Newton step length must be reduced to below this factor to give a residual reduction,...
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
This class implements a goal oriented error indicator, by weighting residual-based estimates from the...
std::unique_ptr< ErrorEstimator > & primal_error_estimator()
Access to the "subestimator" (default: PatchRecovery) to use on the primal/forward solution.
virtual dof_id_type n_active_elem() const =0
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the DiffSolver may tighten the ...
void set_adjoint_already_solved(bool setting)
Setter for the adjoint_already_solved boolean.