Go to the documentation of this file.
34 #include "libmesh/libmesh.h"
35 #include "libmesh/mesh.h"
36 #include "libmesh/mesh_generation.h"
37 #include "libmesh/linear_implicit_system.h"
38 #include "libmesh/equation_systems.h"
39 #include "libmesh/exodusII_io.h"
40 #include "libmesh/gmv_io.h"
41 #include "libmesh/enum_solver_package.h"
44 #include "libmesh/fe.h"
47 #include "libmesh/quadrature_gauss.h"
51 #include "libmesh/sparse_matrix.h"
52 #include "libmesh/numeric_vector.h"
53 #include "libmesh/dense_matrix.h"
54 #include "libmesh/dense_vector.h"
55 #include "libmesh/elem.h"
59 #include "libmesh/dof_map.h"
71 const std::string & system_name);
79 int main (
int argc,
char ** argv)
86 "--enable-petsc, --enable-trilinos, or --enable-eigen");
92 for (
int i=1; i<argc; i++)
98 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
135 equation_systems.
init();
155 equation_systems.
get_system(
"Poisson").solve();
157 #ifdef LIBMESH_HAVE_EXODUS_API
161 #ifdef LIBMESH_HAVE_GMV
177 const std::string & libmesh_dbg_var(system_name))
182 libmesh_assert_equal_to (system_name,
"Poisson");
201 FEType fe_type = dof_map.variable_type(0);
212 fe->attach_quadrature_rule (&qrule);
225 fe_face->attach_quadrature_rule (&qface);
231 const std::vector<Real> & JxW = fe->get_JxW();
236 const std::vector<Point> & q_point = fe->get_xyz();
240 const std::vector<std::vector<RealGradient>> & phi = fe->get_phi();
244 const std::vector<std::vector<RealTensor>> & dphi = fe->get_dphi();
258 std::vector<dof_id_type> dof_indices;
278 dof_map.dof_indices (elem, dof_indices);
295 Ke.
resize (dof_indices.size(),
298 Fe.
resize (dof_indices.size());
302 for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
307 for (std::size_t i=0; i<phi.size(); i++)
308 for (std::size_t j=0; j<phi.size(); j++)
309 Ke(i,j) += JxW[qp] * dphi[i][qp].contract(dphi[j][qp]);
316 const Real x = q_point[qp](0);
317 const Real y = q_point[qp](1);
318 const Real eps = 1.e-3;
348 for (std::size_t i=0; i<phi.size(); i++)
349 Fe(i) += JxW[qp]*f*phi[i][qp];
387 for (
auto side : elem->side_index_range())
388 if (elem->neighbor_ptr(side) ==
nullptr)
392 const std::vector<std::vector<RealGradient>> & phi_face = fe_face->get_phi();
396 const std::vector<Real> & JxW_face = fe_face->get_JxW();
401 const std::vector<Point> & qface_point = fe_face->get_xyz();
405 fe_face->reinit(elem, side);
408 for (
unsigned int qp=0; qp<qface.
n_points(); qp++)
412 const Real xf = qface_point[qp](0);
413 const Real yf = qface_point[qp](1);
417 const Real penalty = 1.e10;
424 for (std::size_t i=0; i<phi_face.size(); i++)
425 for (std::size_t j=0; j<phi_face.size(); j++)
426 Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
430 for (std::size_t i=0; i<phi_face.size(); i++)
431 Fe(i) += JxW_face[qp]*penalty*f*phi_face[i][qp];
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 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
void assemble_poisson(EquationSystems &es, const std::string &system_name)
SolverPackage default_solver_package()
unsigned int mesh_dimension() const
This class implements writing meshes in the GMV format.
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
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.
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.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
virtual void init()
Initialize all the systems.
SparseMatrix< Number > * matrix
The system matrix.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
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 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.
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
int main(int argc, char **argv)
const DofMap & get_dof_map() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
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.