Go to the documentation of this file.
37 #include "libmesh/libmesh.h"
38 #include "libmesh/mesh.h"
39 #include "libmesh/mesh_generation.h"
40 #include "libmesh/vtk_io.h"
41 #include "libmesh/linear_implicit_system.h"
42 #include "libmesh/equation_systems.h"
45 #include "libmesh/fe.h"
48 #include "libmesh/quadrature_gauss.h"
52 #include "libmesh/sparse_matrix.h"
53 #include "libmesh/numeric_vector.h"
54 #include "libmesh/dense_matrix.h"
55 #include "libmesh/dense_vector.h"
56 #include "libmesh/elem.h"
57 #include "libmesh/enum_solver_package.h"
61 #include "libmesh/dof_map.h"
73 const std::string & system_name);
80 int main (
int argc,
char ** argv)
87 "--enable-petsc, --enable-trilinos, or --enable-eigen");
93 for (
int i=1; i<argc; i++)
99 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
138 equation_systems.
init();
158 equation_systems.
get_system(
"Poisson").solve();
160 #if defined(LIBMESH_HAVE_VTK) && !defined(LIBMESH_ENABLE_PARMESH)
166 #endif // #ifdef LIBMESH_HAVE_VTK
180 const std::string & libmesh_dbg_var(system_name))
185 libmesh_assert_equal_to (system_name,
"Poisson");
204 FEType fe_type = dof_map.variable_type(0);
218 fe->attach_quadrature_rule (&qrule);
231 fe_face->attach_quadrature_rule (&qface);
237 const std::vector<Real> & JxW = fe->get_JxW();
242 const std::vector<Point> & q_point = fe->get_xyz();
245 const std::vector<std::vector<Real>> & phi = fe->get_phi();
249 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
263 std::vector<dof_id_type> dof_indices;
283 dof_map.dof_indices (elem, dof_indices);
290 const unsigned int n_dofs =
291 cast_int<unsigned int>(dof_indices.size());
301 libmesh_assert_equal_to (n_dofs, phi.size());
312 Ke.
resize (n_dofs, n_dofs);
318 for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
324 for (
unsigned int i=0; i != n_dofs; i++)
325 for (
unsigned int j=0; j != n_dofs; j++)
327 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
335 const Real x = q_point[qp](0);
336 const Real y = q_point[qp](1);
337 const Real eps = 1.e-3;
360 for (
unsigned int i=0; i != n_dofs; i++)
361 Fe(i) += JxW[qp]*fxy*phi[i][qp];
400 for (
auto side : elem->side_index_range())
401 if (elem->neighbor_ptr(side) ==
nullptr)
405 const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
409 const std::vector<Real> & JxW_face = fe_face->get_JxW();
414 const std::vector<Point> & qface_point = fe_face->get_xyz();
418 fe_face->reinit(elem, side);
423 libmesh_assert_equal_to (n_dofs, phi_face.size());
426 for (
unsigned int qp=0; qp<qface.
n_points(); qp++)
430 const Real xf = qface_point[qp](0);
431 const Real yf = qface_point[qp](1);
435 const Real penalty = 1.e10;
441 for (
unsigned int i=0; i != n_dofs; i++)
442 for (
unsigned int j=0; j != n_dofs; j++)
443 Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
447 for (
unsigned int i=0; i != n_dofs; i++)
448 Fe(i) += JxW_face[qp]*penalty*
value*phi_face[i][qp];
458 dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
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
SolverPackage default_solver_package()
unsigned int mesh_dimension() const
void assemble_poisson(EquationSystems &es, const std::string &system_name)
This class implements reading and writing meshes in the VTK format.
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
int main(int argc, char **argv)
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.
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
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.
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 ...