38 #include "libmesh/libmesh.h" 39 #include "libmesh/mesh.h" 40 #include "libmesh/mesh_refinement.h" 41 #include "libmesh/gmv_io.h" 42 #include "libmesh/equation_systems.h" 43 #include "libmesh/fe.h" 44 #include "libmesh/quadrature_gauss.h" 45 #include "libmesh/dof_map.h" 46 #include "libmesh/sparse_matrix.h" 47 #include "libmesh/numeric_vector.h" 48 #include "libmesh/dense_matrix.h" 49 #include "libmesh/dense_vector.h" 50 #include "libmesh/exodusII_io.h" 51 #include "libmesh/enum_solver_package.h" 52 #include "libmesh/getpot.h" 56 #include "libmesh/linear_implicit_system.h" 57 #include "libmesh/transient_system.h" 58 #include "libmesh/vector_value.h" 61 #include "libmesh/elem.h" 74 const std::string & system_name);
82 const std::string & system_name);
105 int main (
int argc,
char ** argv)
112 "--enable-petsc, --enable-trilinos, or --enable-eigen");
117 #ifndef LIBMESH_ENABLE_AMR 118 libmesh_example_requires(
false,
"--enable-amr");
122 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
136 GetPot input(argc, argv);
137 const unsigned int n_refinements = input(
"n_refinements", 5);
159 system.add_variable (
"u",
FIRST);
164 system.attach_init_function (
init_cd);
167 equation_systems.
init ();
173 #ifdef LIBMESH_HAVE_EXODUS_API 193 const Real dt = 0.025;
196 for (
unsigned int t_step = 0; t_step < 50; t_step++)
210 std::ostringstream
out;
218 << std::setprecision(3)
238 equation_systems.
get_system(
"Convection-Diffusion").solve();
241 if ((t_step+1)%10 == 0)
244 #ifdef LIBMESH_HAVE_EXODUS_API 249 std::ostringstream file_name;
264 #endif // #ifdef LIBMESH_ENABLE_AMR 275 const std::string & libmesh_dbg_var(system_name))
279 libmesh_assert_equal_to (system_name,
"Convection-Diffusion");
297 const std::string & system_name)
302 #ifdef LIBMESH_ENABLE_AMR 305 libmesh_assert_equal_to (system_name,
"Convection-Diffusion");
319 FEType fe_type = system.variable_type(0);
330 QGauss qrule (
dim, fe_type.default_quadrature_order());
331 QGauss qface (
dim-1, fe_type.default_quadrature_order());
334 fe->attach_quadrature_rule (&qrule);
335 fe_face->attach_quadrature_rule (&qface);
340 const std::vector<Real> & JxW = fe->get_JxW();
341 const std::vector<Real> & JxW_face = fe_face->get_JxW();
344 const std::vector<std::vector<Real>> & phi = fe->get_phi();
345 const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
349 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
352 const std::vector<Point> & qface_points = fe_face->get_xyz();
358 const DofMap & dof_map = system.get_dof_map();
370 std::vector<dof_id_type> dof_indices;
386 for (
const auto & elem :
mesh.active_local_element_ptr_range())
406 Ke.
resize (dof_indices.size(),
409 Fe.
resize (dof_indices.size());
417 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
424 for (std::size_t l=0; l<phi.size(); l++)
426 u_old += phi[l][qp]*system.
old_solution (dof_indices[l]);
435 for (std::size_t i=0; i<phi.size(); i++)
445 (grad_u_old*velocity)*phi[i][qp] +
448 0.01*(grad_u_old*dphi[i][qp]))
451 for (std::size_t j=0; j<phi.size(); j++)
456 phi[i][qp]*phi[j][qp] +
460 (velocity*dphi[j][qp])*phi[i][qp] +
463 0.01*(dphi[i][qp]*dphi[j][qp]))
480 const Real penalty = 1.e10;
485 for (
auto s : elem->side_index_range())
486 if (elem->neighbor_ptr(s) ==
nullptr)
488 fe_face->reinit(elem, s);
490 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
497 for (std::size_t i=0; i<psi.size(); i++)
498 Fe(i) += penalty*JxW_face[qp]*
value*psi[i][qp];
501 for (std::size_t i=0; i<psi.size(); i++)
502 for (std::size_t j=0; j<psi.size(); j++)
503 Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
516 matrix.add_matrix (Ke, dof_indices);
517 system.rhs->add_vector (Fe, dof_indices);
521 #endif // #ifdef LIBMESH_ENABLE_AMR class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
This is the EquationSystems class.
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
This class provides the ability to map between arbitrary, user-defined strings and several data types...
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
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 dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
int main(int argc, char **argv)
void resize(const unsigned int n)
Resize the vector.
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Manages storage and variables for transient systems.
Number exact_value(const Point &p, const Parameters ¶meters, const std::string &, const std::string &)
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.
void assemble_cd(EquationSystems &es, const std::string &system_name)
const T_sys & get_system(std::string_view name) const
void append(bool val)
If true, this flag will cause the ExodusII_IO object to attempt to open an existing file for writing...
NumericVector< Number > * old_local_solution
All the values I need to compute my contribution to the simulation at hand.
void init_cd(EquationSystems &es, const std::string &system_name)
This is the MeshBase class.
SolverPackage default_solver_package()
This class handles the numbering of degrees of freedom on a mesh.
void libmesh_ignore(const Args &...)
const T & get(std::string_view) const
Implements (adaptive) mesh refinement algorithms for a MeshBase.
std::string exodus_filename(unsigned number)
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.
Real exact_solution(const Real x, const Real y, const Real t)
This is the exact solution that we are trying to obtain.
Number old_solution(const dof_id_type global_dof_number) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
void write_timestep(const std::string &fname, const EquationSystems &es, const int timestep, const Real time, const std::set< std::string > *system_names=nullptr)
Writes out the solution at a specific timestep.
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
Parameters parameters
Data structure holding arbitrary parameters.
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.
A Point defines a location in LIBMESH_DIM dimensional Real space.
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.