Go to the documentation of this file.
35 #include "libmesh/libmesh.h"
36 #include "libmesh/mesh.h"
37 #include "libmesh/mesh_generation.h"
38 #include "libmesh/exodusII_io.h"
39 #include "libmesh/gnuplot_io.h"
40 #include "libmesh/linear_implicit_system.h"
41 #include "libmesh/equation_systems.h"
44 #include "libmesh/fe.h"
47 #include "libmesh/quadrature_gauss.h"
51 #include "libmesh/dof_map.h"
55 #include "libmesh/sparse_matrix.h"
56 #include "libmesh/numeric_vector.h"
57 #include "libmesh/dense_matrix.h"
58 #include "libmesh/dense_vector.h"
63 #include "libmesh/perf_log.h"
66 #include "libmesh/elem.h"
68 #include "libmesh/string_to_enum.h"
69 #include "libmesh/getpot.h"
70 #include "libmesh/enum_solver_package.h"
84 const std::string & system_name);
92 int main (
int argc,
char ** argv)
99 "--enable-petsc, --enable-trilinos, or --enable-eigen");
105 GetPot command_line (argc, argv);
112 libmesh_error_msg(
"Usage:\n" <<
"\t " << argv[0] <<
" -d 2(3)" <<
" -n 15");
121 for (
int i=1; i<argc; i++)
132 if (command_line.search(1,
"-d"))
133 dim = command_line.next(
dim);
136 libmesh_example_requires(
dim <= LIBMESH_DIM,
"2D/3D support");
144 if (command_line.search(1,
"-n"))
145 ps = command_line.next(ps);
148 std::string order =
"SECOND";
149 if (command_line.search(2,
"-Order",
"-o"))
150 order = command_line.next(order);
153 std::string family =
"LAGRANGE";
154 if (command_line.search(2,
"-FEFamily",
"-f"))
155 family = command_line.next(family);
158 if ((family ==
"MONOMIAL") || (family ==
"XYZ"))
161 libmesh_error_msg(
"This example requires a C^0 (or higher) FE basis.");
170 Real halfwidth =
dim > 1 ? 1. : 0.;
171 Real halfheight =
dim > 2 ? 1. : 0.;
173 if ((family ==
"LAGRANGE") && (order ==
"FIRST"))
182 -halfwidth, halfwidth,
183 -halfheight, halfheight,
195 -halfwidth, halfwidth,
196 -halfheight, halfheight,
203 const Point cent = elem->centroid();
206 if ((cent(0) > 0) == (cent(1) > 0))
207 elem->subdomain_id() = 1;
209 else if (cent(0) > 0)
210 elem->subdomain_id() = 1;
225 std::set<subdomain_id_type> active_subdomains;
230 active_subdomains.
clear(); active_subdomains.insert(0);
232 Utility::string_to_enum<Order> (order),
233 Utility::string_to_enum<FEFamily>(family),
238 active_subdomains.clear(); active_subdomains.insert(1);
240 Utility::string_to_enum<Order> (order),
241 Utility::string_to_enum<FEFamily>(family),
249 equation_systems.
init();
256 equation_systems.
get_system(
"Poisson").solve();
267 #ifdef LIBMESH_HAVE_EXODUS_API
269 "out_3.e" :
"out_2.e", equation_systems);
270 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
289 const std::string & libmesh_dbg_var(system_name))
293 libmesh_assert_equal_to (system_name,
"Poisson");
299 PerfLog perf_log (
"Matrix Assembly");
318 FEType fe_type = dof_map.variable_type(0);
330 fe->attach_quadrature_rule (&qrule);
343 fe_face->attach_quadrature_rule (&qface);
349 const std::vector<Real> & JxW = fe->get_JxW();
354 const std::vector<Point> & q_point = fe->get_xyz();
357 const std::vector<std::vector<Real>> & phi = fe->get_phi();
361 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
373 std::vector<dof_id_type> dof_indices, dof_indices2;
389 perf_log.
push(
"elem init");
395 dof_map.dof_indices (elem, dof_indices, 0);
396 dof_map.dof_indices (elem, dof_indices2, 1);
416 Ke.
resize (std::max(dof_indices.size(), dof_indices2.size()),
417 std::max(dof_indices.size(), dof_indices2.size()));
419 Fe.
resize (std::max(dof_indices.size(), dof_indices2.size()));
424 perf_log.
pop(
"elem init");
435 perf_log.
push (
"Ke");
437 for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
438 for (std::size_t i=0; i<phi.size(); i++)
439 for (std::size_t j=0; j<phi.size(); j++)
440 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
451 perf_log.
push (
"Fe");
453 for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
469 const Real x = q_point[qp](0);
471 const Real y = q_point[qp](1);
476 const Real z = q_point[qp](2);
480 const Real eps = 1.e-3;
500 fxy = (0.25*
pi*
pi)*sin(.5*
pi*x);
504 fxy = - (uxx + uyy + ((
dim==2) ? 0. : uzz));
508 for (std::size_t i=0; i<phi.size(); i++)
509 Fe(i) += JxW[qp]*fxy*phi[i][qp];
524 LOG_SCOPE_WITH(
"BCs",
"", perf_log);
529 for (
auto side : elem->side_index_range())
530 if ((elem->neighbor_ptr(side) ==
nullptr) ||
531 (elem->neighbor_ptr(side)->subdomain_id() != elem->subdomain_id()))
536 const Real penalty = 1.e10;
540 const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
544 const std::vector<Real> & JxW_face = fe_face->get_JxW();
549 const std::vector<Point> & qface_point = fe_face->get_xyz();
553 fe_face->reinit(elem, side);
556 for (
unsigned int qp=0; qp<qface.
n_points(); qp++)
560 const Real xf = qface_point[qp](0);
562 const Real yf = qface_point[qp](1);
567 const Real zf = qface_point[qp](2);
577 for (std::size_t i=0; i<phi_face.size(); i++)
578 for (std::size_t j=0; j<phi_face.size(); j++)
579 Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
583 for (std::size_t i=0; i<phi_face.size(); i++)
584 Fe(i) += JxW_face[qp]*penalty*
value*phi_face[i][qp];
596 LOG_SCOPE_WITH(
"matrix insertion",
"", perf_log);
598 if (dof_indices.size())
604 if (dof_indices2.size())
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 clear() override
Clear all the data structures associated with the system.
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
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
virtual element_iterator local_elements_begin()=0
int main(int argc, char **argv)
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 SimpleRange< element_iterator > element_ptr_range()=0
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 pop(const char *label, const char *header="")
Pop the event label off the stack, resuming any lower event.
processor_id_type processor_id() const
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.
The PerfLog class allows monitoring of specific events.
virtual void init()
Initialize all the systems.
A Point defines a location in LIBMESH_DIM dimensional Real space.
SparseMatrix< Number > * matrix
The system matrix.
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
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.
void assemble_poisson(EquationSystems &es, const std::string &system_name)
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...
virtual element_iterator local_elements_end()=0
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.
void push(const char *label, const char *header="")
Push the event label onto the stack, pausing any active event.
Real exact_solution(const Real x, const Real y=0., const Real z=0.)
This is the exact solution that we are trying to obtain.
const DofMap & get_dof_map() const
This class implements writing meshes using GNUplot, designed for use only with 1D meshes.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...