Go to the documentation of this file.
   35 #include "libmesh/mesh.h" 
   36 #include "libmesh/mesh_generation.h" 
   37 #include "libmesh/edge_edge3.h" 
   38 #include "libmesh/gnuplot_io.h" 
   39 #include "libmesh/equation_systems.h" 
   40 #include "libmesh/linear_implicit_system.h" 
   41 #include "libmesh/fe.h" 
   42 #include "libmesh/getpot.h" 
   43 #include "libmesh/quadrature_gauss.h" 
   44 #include "libmesh/sparse_matrix.h" 
   45 #include "libmesh/dof_map.h" 
   46 #include "libmesh/numeric_vector.h" 
   47 #include "libmesh/dense_matrix.h" 
   48 #include "libmesh/dense_vector.h" 
   49 #include "libmesh/error_vector.h" 
   50 #include "libmesh/kelly_error_estimator.h" 
   51 #include "libmesh/mesh_refinement.h" 
   52 #include "libmesh/enum_solver_package.h" 
   58                  const std::string & system_name);
 
   60 int main(
int argc, 
char ** argv)
 
   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);
 
  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++)
 
  134       if (r_step != max_r_steps)
 
  160           equation_systems.
reinit();
 
  172 #endif // #ifndef LIBMESH_ENABLE_AMR 
  186                  const std::string & system_name)
 
  191 #ifdef LIBMESH_ENABLE_AMR 
  194   libmesh_assert_equal_to (system_name, 
"1D");
 
  212   FEType fe_type = dof_map.variable_type(0);
 
  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 
  
bool refine_and_coarsen_elements()
Refines and coarsens user-requested elements.
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
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
 
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
 
unsigned int n_points() const
 
The libMesh namespace provides an interface to certain functionality in the library.
 
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 &...)
 
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
 
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.
 
virtual void init()
Initialize all the systems.
 
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.
 
int main(int argc, char **argv)
 
unsigned int & max_h_level()
The max_h_level is the greatest refinement level an element should reach.
 
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 ...
 
void resize(const unsigned int n)
Resize the vector.
 
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...
 
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.
 
Real & refine_fraction()
The refine_fraction sets either a desired target or a desired maximum number of elements to flag for ...
 
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 ...
 
void flag_elements_by_error_fraction(const ErrorVector &error_per_cell, const Real refine_fraction=0.3, const Real coarsen_fraction=0.0, const unsigned int max_level=libMesh::invalid_uint)
Flags elements for coarsening and refinement based on the computed error passed in error_per_cell.
 
virtual T maximum() const
 
void assemble_1D(EquationSystems &es, const std::string &system_name)