Go to the source code of this file.
◆ assemble_1D()
Definition at line 185 of file adaptivity_ex1.C.
191 #ifdef LIBMESH_ENABLE_AMR
194 libmesh_assert_equal_to (system_name,
"1D");
212 FEType fe_type = dof_map.variable_type(0);
218 std::unique_ptr<FEBase> fe(FEBase::build(
dim, fe_type));
222 fe->attach_quadrature_rule(&qrule);
228 const std::vector<Real> & JxW = fe->get_JxW();
231 const std::vector<std::vector<Real>> & phi = fe->get_phi();
234 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
244 std::vector<dof_id_type> dof_indices;
255 dof_map.dof_indices(elem, dof_indices);
263 const unsigned int n_dofs =
264 cast_int<unsigned int>(dof_indices.size());
265 libmesh_assert_equal_to (n_dofs, phi.size());
271 Ke.
resize(n_dofs, n_dofs);
276 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
280 for (
unsigned int i=0; i != n_dofs; i++)
282 Fe(i) += JxW[qp]*phi[i][qp];
284 for (
unsigned int j=0; j != n_dofs; j++)
286 Ke(i,j) += JxW[qp]*(1.e-3*dphi[i][qp]*dphi[j][qp] +
287 phi[i][qp]*phi[j][qp]);
298 double penalty = 1.e10;
303 for (
auto s : elem->side_index_range())
308 if (elem->neighbor_ptr(s) ==
nullptr)
317 dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
323 #endif // #ifdef LIBMESH_ENABLE_AMR
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), dim, libMesh::FIFTH, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::libmesh_ignore(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and libMesh::ExplicitSystem::rhs.
Referenced by main().
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 60 of file adaptivity_ex1.C.
71 "--enable-petsc, --enable-trilinos, or --enable-eigen");
74 #ifndef LIBMESH_ENABLE_AMR
75 libmesh_example_requires(
false,
"--enable-amr");
82 GetPot command_line (argc, argv);
85 if (command_line.search(1,
"-n"))
86 n = command_line.next(n);
114 mesh_refinement.refine_fraction() = 0.7;
115 mesh_refinement.coarsen_fraction() = 0.3;
117 mesh_refinement.max_h_level() = 5;
120 equation_systems.init();
123 const unsigned int max_r_steps = 5;
126 for (
unsigned int r_step=0; r_step<=max_r_steps; r_step++)
129 equation_systems.get_system(
"1D").solve();
134 if (r_step != max_r_steps)
152 mesh_refinement.flag_elements_by_error_fraction (error);
155 mesh_refinement.refine_and_coarsen_elements();
160 equation_systems.reinit();
167 GnuPlotIO plot(
mesh,
"Adaptivity Example 1", GnuPlotIO::GRID_ON);
171 plot.write_equation_systems(
"gnuplot_script", equation_systems);
172 #endif // #ifndef LIBMESH_ENABLE_AMR
References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_1D(), libMesh::System::attach_assemble_function(), libMesh::MeshTools::Generation::build_line(), libMesh::MeshRefinement::coarsen_fraction(), libMesh::default_solver_package(), libMesh::EDGE3, libMesh::JumpErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::EquationSystems::get_system(), libMesh::GnuPlotIO::GRID_ON, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::StatisticsVector< T >::l2_norm(), libMesh::MeshRefinement::max_h_level(), libMesh::StatisticsVector< T >::maximum(), mesh, libMesh::out, libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::MeshRefinement::refine_fraction(), libMesh::EquationSystems::reinit(), libMesh::SECOND, libMesh::JumpErrorEstimator::use_unweighted_quadrature_rules, and libMesh::MeshOutput< MT >::write_equation_systems().
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
const MeshBase & get_mesh() const
NumericVector< Number > * rhs
The system matrix.
This class implements specific orders of Gauss quadrature.
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
This function uses the derived class's jump error estimate formula to estimate the error on each cell...
virtual Real l2_norm() const
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
const T_sys & get_system(const std::string &name) const
SolverPackage default_solver_package()
unsigned int mesh_dimension() const
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
This class implements the Kelly error indicator which is based on the flux jumps between elements.
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
This is the MeshBase class.
void libmesh_ignore(const Args &...)
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
SparseMatrix< Number > * matrix
The system matrix.
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
This is the EquationSystems class.
void resize(const unsigned int n)
Resize the vector.
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
This class handles the numbering of degrees of freedom on a mesh.
bool use_unweighted_quadrature_rules
This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweigh...
const DofMap & get_dof_map() const
This class implements writing meshes using GNUplot, designed for use only with 1D meshes.
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
virtual T maximum() const
void assemble_1D(EquationSystems &es, const std::string &system_name)