34 #include "libmesh/libmesh.h" 35 #include "libmesh/mesh.h" 36 #include "libmesh/mesh_generation.h" 37 #include "libmesh/exodusII_io.h" 38 #include "libmesh/gmv_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/mesh_refinement.h" 71 #include "libmesh/system_subset_by_subdomain.h" 73 #include "libmesh/string_to_enum.h" 74 #include "libmesh/getpot.h" 75 #include "libmesh/enum_solver_package.h" 89 const std::string & system_name);
97 int main (
int argc,
char ** argv)
107 #ifndef LIBMESH_ENABLE_AMR 108 libmesh_example_requires(
false,
"--enable-amr");
115 GetPot command_line (argc, argv);
118 libmesh_error_msg_if(argc < 3,
"Usage:\n" <<
"\t " << argv[0] <<
" -d 2(3)" <<
" -n 15");
124 for (
int i=1; i<argc; i++)
135 libmesh_example_requires(
dim <= LIBMESH_DIM,
"2D/3D support");
142 std::string order =
"FIRST";
147 std::string family =
"LAGRANGE";
152 libmesh_error_msg_if((family ==
"MONOMIAL") || (family ==
"XYZ"),
153 "This example requires a C^0 (or higher) FE basis.");
165 Real halfwidth =
dim > 1 ? 1. : 0.;
166 Real halfheight =
dim > 2 ? 1. : 0.;
168 if ((family ==
"LAGRANGE") && (order ==
"FIRST"))
177 -halfwidth, halfwidth,
178 -halfheight, halfheight,
190 -halfwidth, halfwidth,
191 -halfheight, halfheight,
207 for (
auto & elem :
mesh.element_ptr_range())
212 bool node_in =
false;
213 bool node_out =
false;
214 for (
auto & n : elem->node_ref_range())
222 if (node_in && node_out)
239 for (
auto elem :
mesh.element_ptr_range())
241 Real d = elem->vertex_average().norm();
243 elem->subdomain_id() = 1;
258 Utility::string_to_enum<Order> (order),
259 Utility::string_to_enum<FEFamily>(family));
266 equation_systems.
init();
273 std::set<subdomain_id_type> id_list;
284 equation_systems.
get_system(
"Poisson").solve();
296 "out_3.gmv" :
"out_2.gmv", equation_systems);
297 #ifdef LIBMESH_HAVE_EXODUS_API 299 "out_3.e" :
"out_2.e", equation_systems);
300 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 303 #endif // #ifndef LIBMESH_ENABLE_AMR 318 const std::string & libmesh_dbg_var(system_name))
322 libmesh_assert_equal_to (system_name,
"Poisson");
328 PerfLog perf_log (
"Matrix Assembly");
347 FEType fe_type = dof_map.variable_type(0);
359 fe->attach_quadrature_rule (&qrule);
372 fe_face->attach_quadrature_rule (&qface);
378 const std::vector<Real> & JxW = fe->get_JxW();
383 const std::vector<Point> & q_point = fe->get_xyz();
386 const std::vector<std::vector<Real>> & phi = fe->get_phi();
390 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
402 std::vector<dof_id_type> dof_indices;
411 for (
const auto & elem :
mesh.active_local_element_ptr_range())
415 if (elem->subdomain_id()==1)
420 perf_log.
push(
"elem init");
426 dof_map.dof_indices (elem, dof_indices);
440 Ke.
resize (dof_indices.size(),
443 Fe.
resize (dof_indices.size());
448 perf_log.
pop(
"elem init");
459 perf_log.
push (
"Ke");
461 for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
462 for (std::size_t i=0; i<phi.size(); i++)
463 for (std::size_t j=0; j<phi.size(); j++)
464 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
475 perf_log.
push (
"Fe");
477 for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
493 const Real x = q_point[qp](0);
495 const Real y = q_point[qp](1);
500 const Real z = q_point[qp](2);
504 const Real eps = 1.e-3;
524 fxy = (0.25*
pi*
pi)*sin(.5*
pi*x);
528 fxy = - (uxx + uyy + ((
dim==2) ? 0. : uzz));
532 for (std::size_t i=0; i<phi.size(); i++)
533 Fe(i) += JxW[qp]*fxy*phi[i][qp];
548 LOG_SCOPE_WITH(
"BCs",
"", perf_log);
556 for (
auto side : elem->side_index_range())
557 if ((elem->neighbor_ptr(side) ==
nullptr) ||
558 (elem->neighbor_ptr(side)->subdomain_id()!=1))
563 const Real penalty = 1.e10;
567 const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
571 const std::vector<Real> & JxW_face = fe_face->get_JxW();
576 const std::vector<Point> & qface_point = fe_face->get_xyz();
580 fe_face->reinit(elem, side);
583 for (
unsigned int qp=0; qp<qface.
n_points(); qp++)
587 const Real xf = qface_point[qp](0);
589 const Real yf = qface_point[qp](1);
594 const Real zf = qface_point[qp](2);
604 for (std::size_t i=0; i<phi_face.size(); i++)
605 for (std::size_t j=0; j<phi_face.size(); j++)
606 Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
610 for (std::size_t i=0; i<phi_face.size(); i++)
611 Fe(i) += JxW_face[qp]*penalty*
value*phi_face[i][qp];
618 dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
626 LOG_SCOPE_WITH(
"matrix insertion",
"", perf_log);
virtual void restrict_solve_to(const SystemSubset *subset, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override
After calling this method, any solve will be limited to the given subset.
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.
void pop(const char *label, const char *header="")
Pop the event label off the stack, resuming any lower event.
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...
bool refine_elements()
Only refines the user-requested elements.
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.
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.
NumericVector< Number > * rhs
The system matrix.
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.
const T_sys & get_system(std::string_view name) const
This class represents a subset of the dofs of a System, selected by the subdomain_id and possible the...
This is the MeshBase class.
int main(int argc, char **argv)
SolverPackage default_solver_package()
The PerfLog class allows monitoring of specific events.
This class handles the numbering of degrees of freedom on a mesh.
This class implements writing meshes using GNUplot, designed for use only with 1D meshes...
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.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
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.
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 push(const char *label, const char *header="")
Push the event label onto the stack, pausing any active event.
unsigned int n_points() const
Selection of subdomain ids by a list.
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.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void assemble_poisson(EquationSystems &es, const std::string &system_name)
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
virtual void init()
Initialize all the systems.
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