![]() |
|
#include <mesh_smooth.h>
Public Member Functions | |
void | operator() (SolidMesh *orig_mesh_pt, SolidMesh *copy_of_mesh_pt, const Vector< unsigned > &controlled_boundary_id, const unsigned &max_steps=100000000) |
void | operator() (SolidMesh *orig_mesh_pt, SolidMesh *copy_of_mesh_pt, const Vector< unsigned > &controlled_boundary_id, DocInfo doc_info, const unsigned &max_steps=100000000) |
~NonLinearElasticitySmoothMesh () | |
Destructor (empty) More... | |
void | actions_before_newton_solve () |
void | backup () |
void | reset () |
void | doc_solution (DocInfo &doc_info) |
Doc the solution. More... | |
![]() | |
virtual void | debug_hook_fct (const unsigned &i) |
void | set_analytic_dparameter (double *const ¶meter_pt) |
void | unset_analytic_dparameter (double *const ¶meter_pt) |
bool | is_dparameter_calculated_analytically (double *const ¶meter_pt) |
void | set_analytic_hessian_products () |
void | unset_analytic_hessian_products () |
bool | are_hessian_products_calculated_analytically () |
void | set_pinned_values_to_zero () |
bool | distributed () const |
virtual void | actions_before_adapt () |
virtual void | actions_after_adapt () |
Actions that are to be performed after a mesh adaptation. More... | |
OomphCommunicator * | communicator_pt () |
access function to the oomph-lib communicator More... | |
const OomphCommunicator * | communicator_pt () const |
access function to the oomph-lib communicator, const version More... | |
Problem () | |
Problem (const Problem &dummy)=delete | |
Broken copy constructor. More... | |
void | operator= (const Problem &)=delete |
Broken assignment operator. More... | |
virtual | ~Problem () |
Virtual destructor to clean up memory. More... | |
Mesh *& | mesh_pt () |
Return a pointer to the global mesh. More... | |
Mesh *const & | mesh_pt () const |
Return a pointer to the global mesh (const version) More... | |
Mesh *& | mesh_pt (const unsigned &imesh) |
Mesh *const & | mesh_pt (const unsigned &imesh) const |
Return a pointer to the i-th submesh (const version) More... | |
unsigned | nsub_mesh () const |
Return number of submeshes. More... | |
unsigned | add_sub_mesh (Mesh *const &mesh_pt) |
void | flush_sub_meshes () |
void | build_global_mesh () |
void | rebuild_global_mesh () |
LinearSolver *& | linear_solver_pt () |
Return a pointer to the linear solver object. More... | |
LinearSolver *const & | linear_solver_pt () const |
Return a pointer to the linear solver object (const version) More... | |
LinearSolver *& | mass_matrix_solver_for_explicit_timestepper_pt () |
LinearSolver * | mass_matrix_solver_for_explicit_timestepper_pt () const |
EigenSolver *& | eigen_solver_pt () |
Return a pointer to the eigen solver object. More... | |
EigenSolver *const & | eigen_solver_pt () const |
Return a pointer to the eigen solver object (const version) More... | |
Time *& | time_pt () |
Return a pointer to the global time object. More... | |
Time * | time_pt () const |
Return a pointer to the global time object (const version). More... | |
double & | time () |
Return the current value of continuous time. More... | |
double | time () const |
Return the current value of continuous time (const version) More... | |
TimeStepper *& | time_stepper_pt () |
const TimeStepper * | time_stepper_pt () const |
TimeStepper *& | time_stepper_pt (const unsigned &i) |
Return a pointer to the i-th timestepper. More... | |
ExplicitTimeStepper *& | explicit_time_stepper_pt () |
Return a pointer to the explicit timestepper. More... | |
unsigned long | set_timestepper_for_all_data (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data=false) |
virtual void | shift_time_values () |
Shift all values along to prepare for next timestep. More... | |
AssemblyHandler *& | assembly_handler_pt () |
Return a pointer to the assembly handler object. More... | |
AssemblyHandler *const & | assembly_handler_pt () const |
Return a pointer to the assembly handler object (const version) More... | |
double & | minimum_dt () |
Access function to min timestep in adaptive timestepping. More... | |
double & | maximum_dt () |
Access function to max timestep in adaptive timestepping. More... | |
unsigned & | max_newton_iterations () |
Access function to max Newton iterations before giving up. More... | |
void | problem_is_nonlinear (const bool &prob_lin) |
Access function to Problem_is_nonlinear. More... | |
double & | max_residuals () |
bool & | time_adaptive_newton_crash_on_solve_fail () |
Access function for Time_adaptive_newton_crash_on_solve_fail. More... | |
double & | newton_solver_tolerance () |
void | add_time_stepper_pt (TimeStepper *const &time_stepper_pt) |
void | set_explicit_time_stepper_pt (ExplicitTimeStepper *const &explicit_time_stepper_pt) |
void | initialise_dt (const double &dt) |
void | initialise_dt (const Vector< double > &dt) |
Data *& | global_data_pt (const unsigned &i) |
Return a pointer to the the i-th global data object. More... | |
void | add_global_data (Data *const &global_data_pt) |
void | flush_global_data () |
LinearAlgebraDistribution *const & | dof_distribution_pt () const |
Return the pointer to the dof distribution (read-only) More... | |
unsigned long | ndof () const |
Return the number of dofs. More... | |
unsigned | ntime_stepper () const |
Return the number of time steppers. More... | |
unsigned | nglobal_data () const |
Return the number of global data values. More... | |
unsigned | self_test () |
Self-test: Check meshes and global data. Return 0 for OK. More... | |
void | enable_store_local_dof_pt_in_elements () |
void | disable_store_local_dof_pt_in_elements () |
unsigned long | assign_eqn_numbers (const bool &assign_local_eqn_numbers=true) |
void | describe_dofs (std::ostream &out= *(oomph_info.stream_pt())) const |
void | enable_discontinuous_formulation () |
void | disable_discontinuous_formulation () |
void | get_dofs (DoubleVector &dofs) const |
void | get_dofs (const unsigned &t, DoubleVector &dofs) const |
Return vector of the t'th history value of all dofs. More... | |
void | set_dofs (const DoubleVector &dofs) |
Set the values of the dofs. More... | |
void | set_dofs (const unsigned &t, DoubleVector &dofs) |
Set the history values of the dofs. More... | |
void | set_dofs (const unsigned &t, Vector< double * > &dof_pt) |
void | add_to_dofs (const double &lambda, const DoubleVector &increment_dofs) |
Add lambda x incremenet_dofs[l] to the l-th dof. More... | |
double * | global_dof_pt (const unsigned &i) |
double & | dof (const unsigned &i) |
i-th dof in the problem More... | |
double | dof (const unsigned &i) const |
i-th dof in the problem (const version) More... | |
double *& | dof_pt (const unsigned &i) |
Pointer to i-th dof in the problem. More... | |
double * | dof_pt (const unsigned &i) const |
Pointer to i-th dof in the problem (const version) More... | |
virtual void | get_inverse_mass_matrix_times_residuals (DoubleVector &Mres) |
virtual void | get_dvaluesdt (DoubleVector &f) |
virtual void | get_residuals (DoubleVector &residuals) |
Get the total residuals Vector for the problem. More... | |
virtual void | get_jacobian (DoubleVector &residuals, DenseDoubleMatrix &jacobian) |
virtual void | get_jacobian (DoubleVector &residuals, CRDoubleMatrix &jacobian) |
virtual void | get_jacobian (DoubleVector &residuals, CCDoubleMatrix &jacobian) |
virtual void | get_jacobian (DoubleVector &residuals, SumOfMatrices &jacobian) |
void | get_fd_jacobian (DoubleVector &residuals, DenseMatrix< double > &jacobian) |
Get the full Jacobian by finite differencing. More... | |
void | get_derivative_wrt_global_parameter (double *const ¶meter_pt, DoubleVector &result) |
void | get_hessian_vector_products (DoubleVectorWithHaloEntries const &Y, Vector< DoubleVectorWithHaloEntries > const &C, Vector< DoubleVectorWithHaloEntries > &product) |
void | solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector, const bool &steady=true) |
Solve the eigenproblem. More... | |
void | solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, const bool &steady=true) |
virtual void | get_eigenproblem_matrices (CRDoubleMatrix &mass_matrix, CRDoubleMatrix &main_matrix, const double &shift=0.0) |
void | assign_eigenvector_to_dofs (DoubleVector &eigenvector) |
Assign the eigenvector passed to the function to the dofs. More... | |
void | add_eigenvector_to_dofs (const double &epsilon, const DoubleVector &eigenvector) |
void | store_current_dof_values () |
Store the current values of the degrees of freedom. More... | |
void | restore_dof_values () |
Restore the stored values of the degrees of freedom. More... | |
void | enable_jacobian_reuse () |
void | disable_jacobian_reuse () |
Disable recycling of Jacobian in Newton iteration. More... | |
bool | jacobian_reuse_is_enabled () |
Is recycling of Jacobian in Newton iteration enabled? More... | |
bool & | use_predictor_values_as_initial_guess () |
void | newton_solve () |
Use Newton method to solve the problem. More... | |
void | enable_globally_convergent_newton_method () |
enable globally convergent Newton method More... | |
void | disable_globally_convergent_newton_method () |
disable globally convergent Newton method More... | |
void | newton_solve (unsigned const &max_adapt) |
void | steady_newton_solve (unsigned const &max_adapt=0) |
void | copy (Problem *orig_problem_pt) |
virtual Problem * | make_copy () |
virtual void | read (std::ifstream &restart_file, bool &unsteady_restart) |
virtual void | read (std::ifstream &restart_file) |
virtual void | dump (std::ofstream &dump_file) const |
void | dump (const std::string &dump_file_name) const |
void | delete_all_external_storage () |
virtual void | symmetrise_eigenfunction_for_adaptive_pitchfork_tracking () |
double * | bifurcation_parameter_pt () const |
void | get_bifurcation_eigenfunction (Vector< DoubleVector > &eigenfunction) |
void | activate_fold_tracking (double *const ¶meter_pt, const bool &block_solve=true) |
void | activate_bifurcation_tracking (double *const ¶meter_pt, const DoubleVector &eigenvector, const bool &block_solve=true) |
void | activate_bifurcation_tracking (double *const ¶meter_pt, const DoubleVector &eigenvector, const DoubleVector &normalisation, const bool &block_solve=true) |
void | activate_pitchfork_tracking (double *const ¶meter_pt, const DoubleVector &symmetry_vector, const bool &block_solve=true) |
void | activate_hopf_tracking (double *const ¶meter_pt, const bool &block_solve=true) |
void | activate_hopf_tracking (double *const ¶meter_pt, const double &omega, const DoubleVector &null_real, const DoubleVector &null_imag, const bool &block_solve=true) |
void | deactivate_bifurcation_tracking () |
void | reset_assembly_handler_to_default () |
Reset the system to the standard non-augemented state. More... | |
double | arc_length_step_solve (double *const ¶meter_pt, const double &ds, const unsigned &max_adapt=0) |
double | arc_length_step_solve (Data *const &data_pt, const unsigned &data_index, const double &ds, const unsigned &max_adapt=0) |
void | reset_arc_length_parameters () |
int & | sign_of_jacobian () |
void | explicit_timestep (const double &dt, const bool &shift_values=true) |
Take an explicit timestep of size dt. More... | |
void | unsteady_newton_solve (const double &dt) |
void | unsteady_newton_solve (const double &dt, const bool &shift_values) |
void | unsteady_newton_solve (const double &dt, const unsigned &max_adapt, const bool &first, const bool &shift=true) |
double | doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const bool &first, const bool &shift=true) |
double | doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const unsigned &suppress_resolve_after_spatial_adapt_flag, const bool &first, const bool &shift=true) |
double | adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon) |
double | adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon, const bool &shift_values) |
void | assign_initial_values_impulsive () |
void | assign_initial_values_impulsive (const double &dt) |
void | calculate_predictions () |
Calculate predictions. More... | |
void | enable_mass_matrix_reuse () |
void | disable_mass_matrix_reuse () |
bool | mass_matrix_reuse_is_enabled () |
Return whether the mass matrix is being reused. More... | |
void | refine_uniformly (const Vector< unsigned > &nrefine_for_mesh) |
void | refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
void | refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh) |
void | refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
void | refine_uniformly (DocInfo &doc_info) |
void | refine_uniformly_and_prune (DocInfo &doc_info) |
void | refine_uniformly () |
void | refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info) |
Do uniform refinement for submesh i_mesh with documentation. More... | |
void | refine_uniformly (const unsigned &i_mesh) |
Do uniform refinement for submesh i_mesh without documentation. More... | |
void | p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh) |
void | p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
void | p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh) |
void | p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
void | p_refine_uniformly (DocInfo &doc_info) |
void | p_refine_uniformly_and_prune (DocInfo &doc_info) |
void | p_refine_uniformly () |
void | p_refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info) |
Do uniform p-refinement for submesh i_mesh with documentation. More... | |
void | p_refine_uniformly (const unsigned &i_mesh) |
Do uniform p-refinement for submesh i_mesh without documentation. More... | |
void | refine_selected_elements (const Vector< unsigned > &elements_to_be_refined) |
void | refine_selected_elements (const Vector< RefineableElement * > &elements_to_be_refined_pt) |
void | refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined) |
void | refine_selected_elements (const unsigned &i_mesh, const Vector< RefineableElement * > &elements_to_be_refined_pt) |
void | refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined) |
void | refine_selected_elements (const Vector< Vector< RefineableElement * >> &elements_to_be_refined_pt) |
void | p_refine_selected_elements (const Vector< unsigned > &elements_to_be_refined) |
void | p_refine_selected_elements (const Vector< PRefineableElement * > &elements_to_be_refined_pt) |
void | p_refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined) |
void | p_refine_selected_elements (const unsigned &i_mesh, const Vector< PRefineableElement * > &elements_to_be_refined_pt) |
void | p_refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined) |
void | p_refine_selected_elements (const Vector< Vector< PRefineableElement * >> &elements_to_be_refined_pt) |
unsigned | unrefine_uniformly () |
unsigned | unrefine_uniformly (const unsigned &i_mesh) |
void | p_unrefine_uniformly (DocInfo &doc_info) |
void | p_unrefine_uniformly (const unsigned &i_mesh, DocInfo &doc_info) |
Do uniform p-unrefinement for submesh i_mesh without documentation. More... | |
void | adapt (unsigned &n_refined, unsigned &n_unrefined) |
void | adapt () |
void | p_adapt (unsigned &n_refined, unsigned &n_unrefined) |
void | p_adapt () |
void | adapt_based_on_error_estimates (unsigned &n_refined, unsigned &n_unrefined, Vector< Vector< double >> &elemental_error) |
void | adapt_based_on_error_estimates (Vector< Vector< double >> &elemental_error) |
void | get_all_error_estimates (Vector< Vector< double >> &elemental_error) |
void | doc_errors (DocInfo &doc_info) |
Get max and min error for all elements in submeshes. More... | |
void | doc_errors () |
Get max and min error for all elements in submeshes. More... | |
void | enable_info_in_newton_solve () |
void | disable_info_in_newton_solve () |
Disable the output of information when in the newton solver. More... | |
![]() | |
ExplicitTimeSteppableObject () | |
Empty constructor. More... | |
ExplicitTimeSteppableObject (const ExplicitTimeSteppableObject &)=delete | |
Broken copy constructor. More... | |
void | operator= (const ExplicitTimeSteppableObject &)=delete |
Broken assignment operator. More... | |
virtual | ~ExplicitTimeSteppableObject () |
Empty destructor. More... | |
virtual void | actions_before_explicit_stage () |
virtual void | actions_after_explicit_stage () |
Private Attributes | |
Vector< Vector< double > > | Orig_node_pos |
Original nodal positions. More... | |
Vector< Vector< double > > | Backup_node_pos |
Backup nodal positions. More... | |
SolidMesh * | Orig_mesh_pt |
Bulk original mesh. More... | |
SolidMesh * | Dummy_mesh_pt |
Copy of mesh to work on. More... | |
/////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////// Auxiliary Problem to smooth a SolidMesh by adjusting the internal nodal positions via the solution of a nonlinear solid mechanics problem. The mesh will typically have been created with an unstructured mesh generator that uses a low-order (simplex) representation of the element geometry; some of the nodes, typically non-vertex nodes on the domain's curvilinear boundaries, were then moved to their new position to provide a more accurate representation of the geometry. This class should be used to deal with elements that may have become inverted during the node motion.
Important assumption:
Template argument specifies type of element. It must be a pure solid mechanics element! This shouldn't cause any problems since mesh smoothing operations tend to be performed off-line so the mesh may as well be built with pure solid elements even if it is ultimately to be used with other element types. (This restriction could easily be avoided but would require double templating and would generally be messy...)
|
inline |
|
inlinevirtual |
Update nodal positions in main mesh – also moves the nodes of the FaceElements that impose the new position
Reimplemented from oomph::Problem.
References i, j, oomph::Node::ndim(), oomph::Mesh::nnode(), oomph::SolidMesh::node_pt(), oomph::oomph_info, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Orig_mesh_pt, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Orig_node_pos, oomph::Helper_namespace_for_mesh_smoothing::Scale, oomph::Node::x(), and oomph::SolidNode::xi().
|
inline |
Backup nodal positions in dummy mesh to allow for reset after non-convergence of Newton method
References oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Backup_node_pos, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Dummy_mesh_pt, i, j, oomph::Node::ndim(), oomph::Mesh::nnode(), oomph::SolidMesh::node_pt(), and oomph::Node::x().
Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::operator()().
|
inline |
Doc the solution.
References oomph::Mesh::check_inverted_elements(), oomph::DocInfo::directory(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Dummy_mesh_pt, MergeRestartFiles::filename, oomph::DocInfo::is_doc_enabled(), oomph::DocInfo::number(), oomph::oomph_info, and oomph::Mesh::output().
Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::operator()().
|
inline |
Functor to update the nodal positions in SolidMesh pointed to by orig_mesh_pt in response to the displacement of some of its nodes relative to their original position which must still be indicated by the nodes' Lagrangian position. copy_of_mesh_pt must be a deep copy of orig_mesh_pt, with the same boundary coordinates etc. This mesh is used as workspace and can be deleted afterwards. The vector controlled_boundary_id contains the ids of the mesh boundaries in orig_mesh_pt whose position is supposed to remain fixed (while the other nodes are re-positioned to avoid the inversion of elements). The final optional argument specifies the max. number of increments in which the mesh boundary is deformed.
References oomph::DocInfo::disable_doc().
|
inline |
Functor to update the nodal positions in SolidMesh pointed to by orig_mesh_pt in response to the displacement of some of its nodes relative to their original position which must still be indicated by the nodes' Lagrangian position. copy_of_mesh_pt must be a deep copy of orig_mesh_pt, with the same boundary coordinates etc. This mesh is used as workspace and can be deleted afterwards. The vector controlled_boundary_id contains the ids of the mesh boundaries in orig_mesh_pt whose position is supposed to remain fixed (while the other nodes are re-positioned to avoid the inversion of elements). The DocInfo allows allows the output of the intermediate meshes. The final optional argument specifies the max. number of increments in which the mesh boundary is deformed.
References oomph::Problem::add_sub_mesh(), oomph::Problem::assign_eqn_numbers(), b, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::backup(), oomph::Mesh::boundary_element_pt(), oomph::Problem::build_global_mesh(), oomph::Helper_namespace_for_mesh_smoothing::Constitutive_law_pt, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::doc_solution(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Dummy_mesh_pt, e(), oomph::Mesh::element_pt(), oomph::Mesh::face_index_at_boundary(), i, j, n, oomph::Mesh::nboundary(), oomph::Mesh::nboundary_element(), oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::Problem::newton_solve(), oomph::Mesh::nnode(), oomph::SolidMesh::node_pt(), oomph::DocInfo::number(), oomph::oomph_info, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Orig_mesh_pt, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Orig_node_pos, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::reset(), oomph::Helper_namespace_for_mesh_smoothing::Scale, oomph::Helper_namespace_for_mesh_smoothing::Scale_increment, oomph::FaceElement::set_boundary_number_in_bulk_mesh(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt(), oomph::SolidMesh::set_lagrangian_nodal_coordinates(), and oomph::Node::x().
|
inline |
Reset nodal positions in dummy mesh to allow for restart of Newton method with reduced increment in Scale
References oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Backup_node_pos, oomph::NonLinearElasticitySmoothMesh< ELEMENT >::Dummy_mesh_pt, i, j, oomph::Node::ndim(), oomph::Mesh::nnode(), oomph::SolidMesh::node_pt(), and oomph::Node::x().
Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::operator()().
|
private |
Backup nodal positions.
Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::backup(), and oomph::NonLinearElasticitySmoothMesh< ELEMENT >::reset().
|
private |
|
private |
Bulk original mesh.
Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::actions_before_newton_solve(), and oomph::NonLinearElasticitySmoothMesh< ELEMENT >::operator()().
|
private |
Original nodal positions.
Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::actions_before_newton_solve(), and oomph::NonLinearElasticitySmoothMesh< ELEMENT >::operator()().