Go to the source code of this file.
◆ assemble_cd()
Definition at line 293 of file transient_ex1.C.
299 #ifdef LIBMESH_ENABLE_AMR
302 libmesh_assert_equal_to (system_name,
"Convection-Diffusion");
316 FEType fe_type = system.variable_type(0);
322 std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
323 std::unique_ptr<FEBase> fe_face (FEBase::build(
dim, fe_type));
327 QGauss qrule (
dim, fe_type.default_quadrature_order());
328 QGauss qface (
dim-1, fe_type.default_quadrature_order());
331 fe->attach_quadrature_rule (&qrule);
332 fe_face->attach_quadrature_rule (&qface);
337 const std::vector<Real> & JxW = fe->get_JxW();
338 const std::vector<Real> & JxW_face = fe_face->get_JxW();
341 const std::vector<std::vector<Real>> & phi = fe->get_phi();
342 const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
346 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
349 const std::vector<Point> & qface_points = fe_face->get_xyz();
355 const DofMap & dof_map = system.get_dof_map();
367 std::vector<dof_id_type> dof_indices;
401 Ke.
resize (dof_indices.size(),
404 Fe.
resize (dof_indices.size());
412 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
419 for (std::size_t l=0; l<phi.size(); l++)
421 u_old += phi[l][qp]*system.
old_solution (dof_indices[l]);
430 for (std::size_t i=0; i<phi.size(); i++)
440 (grad_u_old*velocity)*phi[i][qp] +
443 0.01*(grad_u_old*dphi[i][qp]))
446 for (std::size_t j=0; j<phi.size(); j++)
451 phi[i][qp]*phi[j][qp] +
455 (velocity*dphi[j][qp])*phi[i][qp] +
458 0.01*(dphi[i][qp]*dphi[j][qp]))
475 const Real penalty = 1.e10;
480 for (
auto s : elem->side_index_range())
481 if (elem->neighbor_ptr(s) ==
nullptr)
483 fe_face->reinit(elem, s);
485 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
492 for (std::size_t i=0; i<psi.size(); i++)
493 Fe(i) += penalty*JxW_face[qp]*
value*psi[i][qp];
496 for (std::size_t i=0; i<psi.size(); i++)
497 for (std::size_t j=0; j<psi.size(); j++)
498 Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
511 system.matrix->add_matrix (Ke, dof_indices);
512 system.rhs->add_vector (Fe, dof_indices);
516 #endif // #ifdef LIBMESH_ENABLE_AMR
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::TypeVector< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), exact_solution(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::libmesh_ignore(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::TransientSystem< Base >::old_solution(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and value.
Referenced by main().
◆ exact_solution()
Real exact_solution |
( |
const Real |
x, |
|
|
const Real |
y, |
|
|
const Real |
t |
|
) |
| |
◆ exact_value()
Number exact_value |
( |
const Point & |
p, |
|
|
const Parameters & |
parameters, |
|
|
const std::string & |
, |
|
|
const std::string & |
|
|
) |
| |
◆ init_cd() [1/2]
void init_cd |
( |
EquationSystems & |
es, |
|
|
const std::string & |
libmesh_dbg_varsystem_name |
|
) |
| |
◆ init_cd() [2/2]
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 104 of file transient_ex1.C.
111 "--enable-petsc, --enable-trilinos, or --enable-eigen");
116 #ifndef LIBMESH_ENABLE_AMR
117 libmesh_example_requires(
false,
"--enable-amr");
121 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
141 mesh_refinement.uniformly_refine (5);
156 system.add_variable (
"u",
FIRST);
161 system.attach_init_function (
init_cd);
164 equation_systems.init ();
167 equation_systems.print_info();
170 #ifdef LIBMESH_HAVE_EXODUS_API
190 const Real dt = 0.025;
193 for (
unsigned int t_step = 0; t_step < 50; t_step++)
199 equation_systems.parameters.set<
Real> (
"time") = system.time;
200 equation_systems.parameters.set<
Real> (
"dt") = dt;
207 std::ostringstream
out;
215 << std::setprecision(3)
235 equation_systems.get_system(
"Convection-Diffusion").solve();
238 if ((t_step+1)%10 == 0)
241 #ifdef LIBMESH_HAVE_EXODUS_API
244 exo.write_timestep (
exodus_filename, equation_systems, t_step+1, system.time);
246 std::ostringstream file_name;
261 #endif // #ifdef LIBMESH_ENABLE_AMR
References libMesh::EquationSystems::add_system(), libMesh::ExodusII_IO::append(), assemble_cd(), libMesh::default_solver_package(), exodus_filename(), libMesh::FIRST, libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), init_cd(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::TransientSystem< Base >::old_local_solution, libMesh::out, libMesh::EquationSystems::parameters, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Real, libMesh::Parameters::set(), libMesh::MeshRefinement::uniformly_refine(), libMesh::MeshOutput< MT >::write_equation_systems(), and libMesh::ExodusII_IO::write_timestep().
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
const MeshBase & get_mesh() const
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.
Manages storage and variables for transient systems.
This class implements specific orders of Gauss quadrature.
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
const T_sys & get_system(const std::string &name) const
void assemble_cd(EquationSystems &es, const std::string &system_name)
SolverPackage default_solver_package()
Number exact_value(const Point &p, const Parameters ¶meters, const std::string &, const std::string &)
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.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
This is the MeshBase class.
std::string exodus_filename(unsigned number)
void libmesh_ignore(const Args &...)
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
void init_cd(EquationSystems &es, const std::string &system_name)
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
This is the EquationSystems class.
T & set(const std::string &)
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 resize(const unsigned int n)
Resize the vector.
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.
Number old_solution(const dof_id_type global_dof_number) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::unique_ptr< NumericVector< Number > > old_local_solution
All the values I need to compute my contribution to the simulation at hand.
Real exact_solution(const Real x, const Real y, const Real t)
This is the exact solution that we are trying to obtain.
const T & get(const std::string &) const
Parameters parameters
Data structure holding arbitrary parameters.