| 
    libMesh
    
   | 
 
This class uses the error estimate given by different types of derefinement (h coarsening or p reduction) to choose between h refining and p elevation. More...
#include <hp_coarsentest.h>
Public Member Functions | |
| HPCoarsenTest () | |
| Constructor.  More... | |
| HPCoarsenTest (const HPCoarsenTest &)=delete | |
| This class cannot be (default) copy constructed/assigned because it has unique_ptr members.  More... | |
| HPCoarsenTest & | operator= (const HPCoarsenTest &)=delete | 
| HPCoarsenTest (HPCoarsenTest &&)=default | |
| Defaulted move ctor, move assignment operator, and destructor.  More... | |
| HPCoarsenTest & | operator= (HPCoarsenTest &&)=default | 
| virtual | ~HPCoarsenTest ()=default | 
| virtual void | select_refinement (System &system) override | 
| This pure virtual function must be redefined in derived classes to take a mesh flagged for h refinement and potentially change the desired refinement type.  More... | |
Public Attributes | |
| Real | p_weight | 
| Because the coarsening test seems to always choose p refinement, we're providing an option to make h refinement more likely.  More... | |
| std::vector< float > | component_scale | 
| This vector can be used to "scale" certain variables in a system.  More... | |
Protected Member Functions | |
| void | add_projection (const System &, const Elem *, unsigned int var) | 
| The helper function which adds individual fine element data to the coarse element projection.  More... | |
Protected Attributes | |
| Elem * | coarse | 
| The coarse element on which a solution projection is cached.  More... | |
| std::vector< dof_id_type > | dof_indices | 
| Global DOF indices for fine elements.  More... | |
| std::unique_ptr< FEBase > | fe | 
| The finite element objects for fine and coarse elements.  More... | |
| std::unique_ptr< FEBase > | fe_coarse | 
| const std::vector< std::vector< Real > > * | phi | 
| The shape functions and their derivatives.  More... | |
| const std::vector< std::vector< Real > > * | phi_coarse | 
| const std::vector< std::vector< RealGradient > > * | dphi | 
| const std::vector< std::vector< RealGradient > > * | dphi_coarse | 
| const std::vector< std::vector< RealTensor > > * | d2phi | 
| const std::vector< std::vector< RealTensor > > * | d2phi_coarse | 
| const std::vector< Real > * | JxW | 
| Mapping jacobians.  More... | |
| const std::vector< Point > * | xyz_values | 
| Quadrature locations.  More... | |
| std::vector< Point > | coarse_qpoints | 
| std::unique_ptr< QBase > | qrule | 
| The quadrature rule for the fine element.  More... | |
| DenseMatrix< Number > | Ke | 
| Linear system for projections.  More... | |
| DenseVector< Number > | Fe | 
| DenseVector< Number > | Uc | 
| Coefficients for projected coarse and projected p-derefined solutions.  More... | |
| DenseVector< Number > | Up | 
This class uses the error estimate given by different types of derefinement (h coarsening or p reduction) to choose between h refining and p elevation.
Currently we assume that a set of elements has already been flagged for h refinement, and we may want to change some of those elements to be flagged for p refinement.
This code is currently experimental and will not produce optimal hp meshes without significant improvement.
Definition at line 67 of file hp_coarsentest.h.
      
  | 
  inline | 
Constructor.
Definition at line 74 of file hp_coarsentest.h.
      
  | 
  delete | 
This class cannot be (default) copy constructed/assigned because it has unique_ptr members.
Explicitly deleting these functions is the best way to document this fact.
      
  | 
  default | 
Defaulted move ctor, move assignment operator, and destructor.
      
  | 
  virtualdefault | 
      
  | 
  protected | 
The helper function which adds individual fine element data to the coarse element projection.
Definition at line 47 of file hp_coarsentest.C.
References libMesh::Elem::active(), libMesh::TypeVector< T >::add_scaled(), libMesh::TypeTensor< T >::add_scaled(), libMesh::C_ONE, libMesh::C_ZERO, libMesh::Elem::child_ref_range(), coarse, coarse_qpoints, libMesh::TypeTensor< T >::contract(), libMesh::System::current_solution(), d2phi, d2phi_coarse, dof_indices, libMesh::DofMap::dof_indices(), dphi, fe, Fe, fe_coarse, libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMesh::index_range(), libMesh::FEMap::inverse_map(), Ke, libMesh::libmesh_assert(), libMesh::MeshBase::mesh_dimension(), phi_coarse, qrule, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::DenseVector< T >::size(), libMesh::Elem::subactive(), Uc, xyz_values, libMesh::DenseMatrix< T >::zero(), libMesh::DenseVector< T >::zero(), and libMesh::zero.
Referenced by select_refinement().
      
  | 
  delete | 
      
  | 
  default | 
      
  | 
  overridevirtual | 
This pure virtual function must be redefined in derived classes to take a mesh flagged for h refinement and potentially change the desired refinement type.
Implements libMesh::HPSelector.
Definition at line 137 of file hp_coarsentest.C.
References add_projection(), libMesh::TypeVector< T >::add_scaled(), libMesh::TypeTensor< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), libMesh::C_ONE, libMesh::C_ZERO, libMesh::DenseMatrix< T >::cholesky_solve(), coarse, coarse_qpoints, libMesh::HPSelector::component_scale, libMesh::TypeTensor< T >::contract(), libMesh::System::current_solution(), d2phi, d2phi_coarse, libMesh::FEType::default_quadrature_rule(), dim, libMesh::DISCONTINUOUS, libMesh::Elem::DO_NOTHING, dof_indices, libMesh::DofMap::dof_indices(), libMesh::DofObject::dof_number(), dphi, dphi_coarse, fe, Fe, fe_coarse, libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMesh::Elem::hack_p_level(), libMesh::index_range(), libMesh::FEMap::inverse_map(), JxW, Ke, libMesh::libmesh_assert(), libMesh::MeshRefinement::make_flags_parallel_consistent(), mesh, libMesh::FEInterface::n_dofs(), n_nodes, n_vars, libMesh::System::n_vars(), libMesh::TensorTools::norm_sq(), libMesh::TypeVector< T >::norm_sq(), libMesh::TypeTensor< T >::norm_sq(), libMesh::System::number(), libMesh::FEType::order, libMesh::Elem::p_level(), p_weight, libMesh::Elem::parent(), phi, phi_coarse, qrule, libMesh::Real, libMesh::Elem::REFINE, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), std::sqrt(), libMesh::TypeVector< T >::subtract_scaled(), libMesh::TypeTensor< T >::subtract_scaled(), Uc, Up, libMesh::DofMap::variable_type(), xyz_values, libMesh::DenseMatrix< T >::zero(), libMesh::DenseVector< T >::zero(), and libMesh::zero.
Referenced by main().
      
  | 
  protected | 
The coarse element on which a solution projection is cached.
Definition at line 118 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Definition at line 146 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  inherited | 
This vector can be used to "scale" certain variables in a system.
If the mask is not empty, the consideration given to each component's h and p error estimates will be scaled by component_scale[c].
Definition at line 82 of file hp_selector.h.
Referenced by select_refinement().
      
  | 
  protected | 
Definition at line 135 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Definition at line 135 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Global DOF indices for fine elements.
Definition at line 123 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Definition at line 134 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Definition at line 134 of file hp_coarsentest.h.
Referenced by select_refinement().
      
  | 
  protected | 
The finite element objects for fine and coarse elements.
Definition at line 128 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Definition at line 157 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Definition at line 128 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Mapping jacobians.
Definition at line 140 of file hp_coarsentest.h.
Referenced by select_refinement().
      
  | 
  protected | 
Linear system for projections.
Definition at line 156 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
| Real libMesh::HPCoarsenTest::p_weight | 
Because the coarsening test seems to always choose p refinement, we're providing an option to make h refinement more likely.
Definition at line 106 of file hp_coarsentest.h.
Referenced by select_refinement().
      
  | 
  protected | 
The shape functions and their derivatives.
Definition at line 133 of file hp_coarsentest.h.
Referenced by select_refinement().
      
  | 
  protected | 
Definition at line 133 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
The quadrature rule for the fine element.
Definition at line 151 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Coefficients for projected coarse and projected p-derefined solutions.
Definition at line 162 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
      
  | 
  protected | 
Definition at line 163 of file hp_coarsentest.h.
Referenced by select_refinement().
      
  | 
  protected | 
Quadrature locations.
Definition at line 145 of file hp_coarsentest.h.
Referenced by add_projection(), and select_refinement().
 1.8.16