Go to the documentation of this file.
   60 #include "libmesh/equation_systems.h" 
   61 #include "libmesh/error_vector.h" 
   62 #include "libmesh/factory.h" 
   63 #include "libmesh/linear_solver.h" 
   64 #include "libmesh/mesh.h" 
   65 #include "libmesh/mesh_refinement.h" 
   66 #include "libmesh/newton_solver.h" 
   67 #include "libmesh/numeric_vector.h" 
   68 #include "libmesh/partitioner.h" 
   69 #include "libmesh/steady_solver.h" 
   70 #include "libmesh/system_norm.h" 
   71 #include "libmesh/auto_ptr.h"  
   72 #include "libmesh/enum_solver_package.h" 
   75 #include "libmesh/kelly_error_estimator.h" 
   76 #include "libmesh/patch_recovery_error_estimator.h" 
   79 #include "libmesh/adjoint_residual_error_estimator.h" 
   80 #include "libmesh/qoi_set.h" 
   83 #include "libmesh/getpot.h" 
   84 #include "libmesh/gmv_io.h" 
   85 #include "libmesh/exodusII_io.h" 
   88 #include "femparameters.h" 
  103                   std::string solution_type,  
 
  109 #ifdef LIBMESH_HAVE_GMV 
  114       std::ostringstream file_name_gmv;
 
  115       file_name_gmv << solution_type
 
  123         (file_name_gmv.str(), es);
 
  127 #ifdef LIBMESH_HAVE_EXODUS_API 
  140       std::ostringstream file_name_exodus;
 
  142       file_name_exodus << solution_type << 
".e";
 
  144         file_name_exodus << 
"-s" 
  184   system.
time_solver = libmesh_make_unique<SteadySolver>(system);
 
  189     system.
time_solver->diff_solver() = std::unique_ptr<DiffSolver>(solver);
 
  199     if (system.
time_solver->reduce_deltat_on_diffsolver_failure)
 
  214 #ifdef LIBMESH_ENABLE_AMR 
  227   return std::unique_ptr<MeshRefinement>(mesh_refinement);
 
  230 #endif // LIBMESH_ENABLE_AMR 
  244       libMesh::out << 
"Using Kelly Error Estimator" << std::endl;
 
  246       return libmesh_make_unique<KellyErrorEstimator>();
 
  250       libMesh::out << 
"Using Adjoint Residual Error Estimator with Patch Recovery Weights" << std::endl;
 
  254       adjoint_residual_estimator->
qoi_set() = qois;
 
  270       return std::unique_ptr<ErrorEstimator>(adjoint_residual_estimator);
 
  273     libmesh_error_msg(
"Unknown indicator_type = " << param.
indicator_type);
 
  277 int main (
int argc, 
char ** argv)
 
  284                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
  287 #ifndef LIBMESH_ENABLE_EXCEPTIONS 
  288   libmesh_example_requires(
false, 
"--enable-exceptions");
 
  292 #ifndef LIBMESH_ENABLE_AMR 
  293   libmesh_example_requires(
false, 
"--enable-amr");
 
  300     std::ifstream i(
"general.in");
 
  302       libmesh_error_msg(
'[' << 
init.comm().rank() << 
"] Can't find general.in; exiting early.");
 
  306   GetPot infile(
"general.in");
 
  309   infile.parse_command_line(argc, argv);
 
  315   libmesh_example_requires(2 <= LIBMESH_DIM, 
"2D support");
 
  323   if (param.mesh_partitioner_type != 
"Default")
 
  325 #ifndef LIBMESH_HAVE_RTTI 
  326       libmesh_example_requires(
false, 
"RTTI support");
 
  338               (param.mesh_partitioner_type == 
"Centroid" ||
 
  339                param.mesh_partitioner_type == 
"Hilbert" ||
 
  340                param.mesh_partitioner_type == 
"Morton" ||
 
  341                param.mesh_partitioner_type == 
"SFCurves" ||
 
  342                param.mesh_partitioner_type == 
"Metis"))
 
  343             libmesh_example_requires(
false, 
"--disable-parmesh");
 
  350           libmesh_example_requires(
false, param.mesh_partitioner_type + 
" partitioner support");
 
  353 #endif // LIBMESH_HAVE_RTI 
  357   std::unique_ptr<MeshRefinement> mesh_refinement =
 
  363   libMesh::out << 
"Reading in and building the mesh" << std::endl;
 
  366   mesh.
read(param.domainfile.c_str());
 
  386   equation_systems.
init ();
 
  395     unsigned int a_step = 0;
 
  396     for (; a_step != param.max_adaptivesteps; ++a_step)
 
  400         if (param.global_tolerance != 0.)
 
  401           libmesh_assert_equal_to (param.nelem_target, 0);
 
  405           libmesh_assert_greater (param.nelem_target, 0);
 
  413         write_output(equation_systems, a_step, 
"primal", param);
 
  422         std::vector<unsigned int> qoi_indices;
 
  423         qoi_indices.push_back(0);
 
  424         qoi_indices.push_back(1);
 
  450         primal_solution.
swap(dual_solution_0);
 
  451         write_output(equation_systems, a_step, 
"adjoint_0", param);
 
  454         primal_solution.
swap(dual_solution_0);
 
  460         primal_solution.
swap(dual_solution_1);
 
  461         write_output(equation_systems, a_step, 
"adjoint_1", param);
 
  464         primal_solution.
swap(dual_solution_1);
 
  470                      << 
" active elements and " 
  485                      << std::setprecision(17)
 
  490                      << std::setprecision(17)
 
  499         std::unique_ptr<ErrorEstimator> error_estimator =
 
  503         error_estimator->estimate_error(system, error);
 
  510         if (param.refine_uniformly)
 
  512             mesh_refinement->uniformly_refine(1);
 
  515         else if (param.global_tolerance >= 0. && param.nelem_target == 0.)
 
  517             mesh_refinement->flag_elements_by_error_tolerance (error);
 
  519             mesh_refinement->refine_and_coarsen_elements();
 
  526                 libMesh::out << 
"We reached the target number of elements." << std::endl << std::endl;
 
  530             mesh_refinement->flag_elements_by_nelem_target (error);
 
  532             mesh_refinement->refine_and_coarsen_elements();
 
  536         equation_systems.
reinit();
 
  540                      << 
" active elements and " 
  547     if (a_step == param.max_adaptivesteps)
 
  552         write_output(equation_systems, a_step, 
"primal", param);
 
  557         std::vector<unsigned int> qoi_indices;
 
  559         qoi_indices.push_back(0);
 
  560         qoi_indices.push_back(1);
 
  575         primal_solution.
swap(dual_solution_0);
 
  576         write_output(equation_systems, a_step, 
"adjoint_0", param);
 
  578         primal_solution.
swap(dual_solution_0);
 
  582         primal_solution.
swap(dual_solution_1);
 
  583         write_output(equation_systems, a_step, 
"adjoint_1", param);
 
  585         primal_solution.
swap(dual_solution_1);
 
  591                      << 
" active elements and " 
  606                      << std::setprecision(17)
 
  611                      << std::setprecision(17)
 
  617         if (param.max_adaptivesteps > 5 && param.coarserefinements > 2)
 
  619             libmesh_assert_less(
std::abs(QoI_0_computed - QoI_0_exact)/
std::abs(QoI_0_exact), 4.e-5);
 
  620             libmesh_assert_less(
std::abs(QoI_1_computed - QoI_1_exact)/
std::abs(QoI_1_exact), 1.e-4);
 
  626             libmesh_assert_less(
std::abs(QoI_0_computed - QoI_0_exact)/
std::abs(QoI_0_exact), 4.e-4);
 
  627             libmesh_assert_less(
std::abs(QoI_1_computed - QoI_1_exact)/
std::abs(QoI_1_exact), 2.e-3);
 
  633                << 
"] Completing output." 
  636 #endif // #ifndef LIBMESH_ENABLE_AMR 
  
virtual void reuse_preconditioner(bool)
Set the same_preconditioner flag, which indicates if we reuse the same preconditioner for subsequent ...
 
void set_system_parameters(LaplaceSystem &system, FEMParameters ¶m)
 
virtual void solve() override
Invokes the solver associated with the system.
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
virtual void postprocess()
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
 
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.
 
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
 
streambufT * rdbuf() const
Get the associated stream buffer.
 
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...
 
libMesh::Real refine_fraction
 
bool print_solution_norms
Set print_residual_norms to true to print |U| whenever it is used in an assembly() call.
 
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...
 
static std::unique_ptr< Base > build(const std::string &name)
Builds an object of type Base identified by name.
 
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)
 
int main(int argc, char **argv)
 
std::size_t n_active_dofs() const
 
double minimum_linear_tolerance
The tolerance for linear solves is kept above this minimum.
 
QoISet & qoi_set()
Access to the QoISet (default: weight all QoIs equally) to use when computing errors.
 
Real relative_step_tolerance
 
bool verbose
The DiffSolver may print a lot more to libMesh::out if verbose is set to true; default is false.
 
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 &...)
 
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
 
virtual void init()
Initialize all the systems.
 
std::unique_ptr< ErrorEstimator > build_error_estimator(FEMParameters ¶m, QoISet &qois)
 
bool continue_after_max_iterations
Defaults to true, telling the DiffSolver to continue rather than exit when a solve has reached its ma...
 
virtual LinearSolver< Number > * get_linear_solver() const override
 
libMesh::Real relative_residual_tolerance
 
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
 
void set_weight(std::size_t, Real)
Set the weight for this index.
 
libMesh::Real relative_step_tolerance
 
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...
 
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 write_output(EquationSystems &es, unsigned int a_step, std::string solution_type, FEMParameters ¶m)
 
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...
 
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.
 
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.
 
void set_patch_reuse(bool)
 
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...
 
bool & coarsen_by_parents()
If coarsen_by_parents is true, complete groups of sibling elements (elements with the same parent) wi...
 
bool postprocess_sides
If postprocess_sides is true (it is false by default), the postprocessing loop will loop over all sid...
 
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,...
 
Number & get_QoI_value(std::string type, unsigned int QoI_index)
 
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
 
std::unique_ptr< MeshRefinement > build_mesh_refinement(MeshBase &mesh, FEMParameters ¶m)
 
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
 
virtual bool is_replicated() const
 
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.
 
virtual std::unique_ptr< Partitioner > & partitioner()
A partitioner to use at each prepare_for_use()