34 #include "libmesh/enum_solver_type.h" 35 #include "libmesh/libmesh.h" 36 #include "libmesh/mesh.h" 37 #include "libmesh/mesh_generation.h" 38 #include "libmesh/linear_implicit_system.h" 39 #include "libmesh/linear_solver.h" 40 #include "libmesh/equation_systems.h" 41 #include "libmesh/exodusII_io.h" 42 #include "libmesh/gmv_io.h" 45 #include "libmesh/enum_norm_type.h" 46 #include "libmesh/enum_solver_package.h" 47 #include "libmesh/getpot.h" 48 #include "libmesh/string_to_enum.h" 51 #include "libmesh/fe.h" 52 #include "libmesh/fe_interface.h" 55 #include "libmesh/quadrature_gauss.h" 59 #include "libmesh/sparse_matrix.h" 60 #include "libmesh/numeric_vector.h" 61 #include "libmesh/dense_matrix.h" 62 #include "libmesh/dense_vector.h" 63 #include "libmesh/elem.h" 67 #include "libmesh/dof_map.h" 79 const std::string & system_name);
87 int main (
int argc,
char ** argv)
94 "--enable-petsc, --enable-trilinos, or --enable-eigen");
100 for (
int i=1; i<argc; i++)
106 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
136 std::string order_str =
"SECOND";
139 const Order order = Utility::string_to_enum<Order>(order_str);
142 std::string family_str =
"LAGRANGE_VEC";
145 const FEFamily family = Utility::string_to_enum<FEFamily>(family_str);
148 "FE family " + family_str +
" isn't vector-valued");
161 equation_systems.
init();
192 libMesh::out <<
"L2 norm of solution = " << std::setprecision(17) <<
193 l2_norm << std::endl;
196 "Failed to calculate solution");
198 const Real error_in_norm = std::abs(l2_norm - sqrt(
Real(2)));
200 libMesh::out <<
"error in L2 norm = " << std::setprecision(17) <<
201 error_in_norm << std::endl;
206 const int n = std::min(nx, ny);
207 const int p =
static_cast<int>(order);
209 libMesh::out <<
"error bound = " << std::setprecision(17) <<
210 expected_error_bound << std::endl;
211 libMesh::out <<
"error ratio = " << std::setprecision(17) <<
212 error_in_norm / expected_error_bound << std::endl;
213 libmesh_error_msg_if (error_in_norm > expected_error_bound,
214 "Error exceeds expected bound of " <<
215 expected_error_bound);
217 #ifdef LIBMESH_HAVE_EXODUS_API 221 #ifdef LIBMESH_HAVE_GMV 237 const std::string & libmesh_dbg_var(system_name))
242 libmesh_assert_equal_to (system_name,
"Poisson");
261 FEType fe_type = dof_map.variable_type(0);
272 fe->attach_quadrature_rule (&qrule);
285 fe_face->attach_quadrature_rule (&qface);
291 const std::vector<Real> & JxW = fe->get_JxW();
296 const std::vector<Point> & q_point = fe->get_xyz();
300 const std::vector<std::vector<RealGradient>> & phi = fe->get_phi();
304 const std::vector<std::vector<RealTensor>> & dphi = fe->get_dphi();
318 std::vector<dof_id_type> dof_indices;
335 for (
const auto & elem :
mesh.active_local_element_ptr_range())
341 dof_map.dof_indices (elem, dof_indices);
358 Ke.
resize (dof_indices.size(),
361 Fe.
resize (dof_indices.size());
365 const Real eps = 1.e-3 * elem->hmin();
369 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
374 for (std::size_t i=0; i<phi.size(); i++)
375 for (std::size_t j=0; j<phi.size(); j++)
376 Ke(i,j) += JxW[qp] * dphi[i][qp].contract(dphi[j][qp]);
383 const Real x = q_point[qp](0);
384 const Real y = q_point[qp](1);
414 for (std::size_t i=0; i<phi.size(); i++)
415 Fe(i) += JxW[qp]*f*phi[i][qp];
453 for (
auto side : elem->side_index_range())
454 if (elem->neighbor_ptr(side) ==
nullptr)
458 const std::vector<std::vector<RealGradient>> & phi_face = fe_face->get_phi();
462 const std::vector<Real> & JxW_face = fe_face->get_JxW();
467 const std::vector<Point> & qface_point = fe_face->get_xyz();
471 fe_face->reinit(elem, side);
474 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
478 const Real xf = qface_point[qp](0);
479 const Real yf = qface_point[qp](1);
483 const Real penalty = 1.e10;
490 for (std::size_t i=0; i<phi_face.size(); i++)
491 for (std::size_t j=0; j<phi_face.size(); j++)
492 Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
496 for (std::size_t i=0; i<phi_face.size(); i++)
497 Fe(i) += JxW_face[qp]*penalty*f*phi_face[i][qp];
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
T command_line_next(std::string name, T default_value)
Use GetPot's search()/next() functions to get following arguments from the command line...
This is the EquationSystems class.
Order
defines an enum for polynomial orders.
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 ...
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void resize(const unsigned int n)
Resize the vector.
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]...
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
static FEFieldType field_type(const FEType &fe_type)
NumericVector< Number > * rhs
The system matrix.
Order default_quadrature_order() const
virtual LinearSolver< Number > * get_linear_solver() const override
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
The libMesh namespace provides an interface to certain functionality in the library.
This class implements writing meshes in the GMV format.
virtual void solve() override
Assembles & solves the linear system A*x=b.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
SolverPackage default_solver_package()
This class handles the numbering of degrees of freedom on a mesh.
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.
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
unsigned int add_variable(std::string_view 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.
void set_solver_type(const SolverType st)
Sets the type of solver to use.
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.
int main(int argc, char **argv)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MeshBase & get_mesh() const
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Real exact_solution(const int component, const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
virtual void init()
Initialize all the systems.
FEFamily
defines an enum for finite element families.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
const DofMap & get_dof_map() const
const SparseMatrix< Number > & get_system_matrix() const
void assemble_poisson(EquationSystems &es, const std::string &system_name)