Go to the source code of this file.
◆ build_adjoint_refinement_error_estimator()
◆ build_mesh_refinement()
Definition at line 208 of file adjoints_ex4.C.
219 return std::unique_ptr<MeshRefinement>(mesh_refinement);
References libMesh::MeshRefinement::absolute_global_tolerance(), libMesh::MeshRefinement::coarsen_by_parents(), FEMParameters::coarsen_fraction, libMesh::MeshRefinement::coarsen_fraction(), FEMParameters::coarsen_threshold, libMesh::MeshRefinement::coarsen_threshold(), FEMParameters::global_tolerance, mesh, FEMParameters::nelem_target, libMesh::MeshRefinement::nelem_target(), FEMParameters::refine_fraction, and libMesh::MeshRefinement::refine_fraction().
Referenced by main().
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 243 of file adjoints_ex4.C.
250 "--enable-petsc, --enable-trilinos, or --enable-eigen");
253 #ifndef LIBMESH_ENABLE_AMR
254 libmesh_example_requires(
false,
"--enable-amr");
264 std::ifstream i(
"general.in");
266 libmesh_error_msg(
'[' <<
init.comm().rank() <<
"] Can't find general.in; exiting early.");
268 GetPot infile(
"general.in");
275 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
282 std::unique_ptr<MeshRefinement> mesh_refinement =
288 libMesh::out <<
"Reading in and building the mesh" << std::endl;
291 mesh.
read(param.domainfile.c_str());
299 initial_uniform_refinements.uniformly_refine(param.coarserefinements);
311 equation_systems.init ();
315 equation_systems.print_info();
320 unsigned int a_step = 0;
321 for (; a_step != param.max_adaptivesteps; ++a_step)
325 if (param.global_tolerance != 0.)
326 libmesh_assert_equal_to (param.nelem_target, 0);
330 libmesh_assert_greater (param.nelem_target, 0);
338 write_output(equation_systems, a_step,
"primal", param);
347 std::vector<unsigned int> qoi_indices;
348 qoi_indices.push_back(0);
349 qoi_indices.push_back(1);
375 primal_solution.
swap(dual_solution_0);
376 write_output(equation_systems, a_step,
"adjoint_0", param);
379 primal_solution.
swap(dual_solution_0);
385 primal_solution.
swap(dual_solution_1);
386 write_output(equation_systems, a_step,
"adjoint_1", param);
389 primal_solution.
swap(dual_solution_1);
395 <<
" active elements and "
396 << equation_systems.n_active_dofs()
410 << std::setprecision(17)
415 << std::setprecision(17)
424 std::unique_ptr<AdjointRefinementEstimator> adjoint_refinement_error_estimator =
428 adjoint_refinement_error_estimator->estimate_error(system, QoI_elementwise_error);
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)
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)
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)
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)
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]);
467 if (param.refine_uniformly)
469 mesh_refinement->uniformly_refine(1);
472 else if (param.global_tolerance >= 0. && param.nelem_target == 0.)
474 mesh_refinement->flag_elements_by_error_tolerance (QoI_elementwise_error);
476 mesh_refinement->refine_and_coarsen_elements();
483 libMesh::out <<
"We reached the target number of elements." << std::endl << std::endl;
487 mesh_refinement->flag_elements_by_nelem_target (QoI_elementwise_error);
489 mesh_refinement->refine_and_coarsen_elements();
493 equation_systems.reinit();
497 <<
" active elements and "
498 << equation_systems.n_active_dofs()
504 if (a_step == param.max_adaptivesteps)
509 write_output(equation_systems, a_step,
"primal", param);
514 std::vector<unsigned int> qoi_indices;
516 qoi_indices.push_back(0);
517 qoi_indices.push_back(1);
532 primal_solution.
swap(dual_solution_0);
533 write_output(equation_systems, a_step,
"adjoint_0", param);
535 primal_solution.
swap(dual_solution_0);
539 primal_solution.
swap(dual_solution_1);
540 write_output(equation_systems, a_step,
"adjoint_1", param);
542 primal_solution.
swap(dual_solution_1);
548 <<
" active elements and "
549 << equation_systems.n_active_dofs()
563 << std::setprecision(17)
568 << std::setprecision(17)
579 std::unique_ptr<AdjointRefinementEstimator> adjoint_refinement_error_estimator =
583 adjoint_refinement_error_estimator->estimate_error(system, QoI_elementwise_error);
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)
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)
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)
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)
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);
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);
629 <<
"] Completing output."
632 #endif // #ifndef LIBMESH_ENABLE_AMR
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()
Definition at line 155 of file adjoints_ex4.C.
177 system.
time_solver = libmesh_make_unique<SteadySolver>(system);
182 system.
time_solver->diff_solver() = std::unique_ptr<DiffSolver>(solver);
191 if (system.
time_solver->reduce_deltat_on_diffsolver_failure)
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()
virtual void reuse_preconditioner(bool)
Set the same_preconditioner flag, which indicates if we reuse the same preconditioner for subsequent ...
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.
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.
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.
void set_system_parameters(LaplaceSystem &system, FEMParameters ¶m)
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
std::unique_ptr< MeshRefinement > build_mesh_refinement(MeshBase &mesh, FEMParameters ¶m)
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.
std::unique_ptr< AdjointRefinementEstimator > build_adjoint_refinement_error_estimator(QoISet &qois)
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.
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)
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 &...)
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
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...
unsigned char number_h_refinements
How many h refinements to perform to get the fine grid.
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 ...
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 ...
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.
QoISet & qoi_set()
Access to the QoISet (default: weight all QoIs equally) to use when computing errors.
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)
void write_output(EquationSystems &es, unsigned int a_step, std::string solution_type, FEMParameters ¶m)
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.
This class implements a "brute force" goal-oriented error estimator which computes an estimate of err...
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)
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.