Go to the documentation of this file.
38 #include "libmesh/libmesh.h"
39 #include "libmesh/replicated_mesh.h"
40 #include "libmesh/mesh_refinement.h"
41 #include "libmesh/gmv_io.h"
42 #include "libmesh/exodusII_io.h"
43 #include "libmesh/equation_systems.h"
44 #include "libmesh/fe.h"
45 #include "libmesh/quadrature_gauss.h"
46 #include "libmesh/dof_map.h"
47 #include "libmesh/sparse_matrix.h"
48 #include "libmesh/numeric_vector.h"
49 #include "libmesh/dense_matrix.h"
50 #include "libmesh/dense_vector.h"
51 #include "libmesh/periodic_boundaries.h"
52 #include "libmesh/periodic_boundary.h"
53 #include "libmesh/mesh_generation.h"
54 #include "libmesh/parsed_function.h"
55 #include "libmesh/getpot.h"
56 #include "libmesh/auto_ptr.h"
57 #include "libmesh/enum_solver_package.h"
58 #include "libmesh/enum_xdr_mode.h"
59 #include "libmesh/enum_norm_type.h"
63 #include "libmesh/transient_system.h"
64 #include "libmesh/linear_implicit_system.h"
65 #include "libmesh/vector_value.h"
69 #include "libmesh/error_vector.h"
70 #include "libmesh/kelly_error_estimator.h"
73 #include "libmesh/elem.h"
85 #ifdef LIBMESH_ENABLE_AMR
87 const std::string & system_name);
96 const std::string & system_name);
129 int main (
int argc,
char ** argv)
136 "--enable-petsc, --enable-trilinos, or --enable-eigen");
139 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
141 #if !defined(LIBMESH_ENABLE_AMR)
142 libmesh_example_requires(
false,
"--enable-amr");
143 #elif !defined(LIBMESH_HAVE_XDR)
145 libmesh_example_requires(
false,
"--enable-xdr");
146 #elif !defined(LIBMESH_ENABLE_PERIODIC)
147 libmesh_example_requires(
false,
"--enable-periodic");
148 #elif (LIBMESH_DOF_ID_BYTES == 8)
149 libmesh_example_requires(
false,
"--with-dof-id-bytes=4");
162 <<
"\t " << argv[0] <<
" -init_timestep 0\n"
164 <<
"\t " << argv[0] <<
" -read_solution -init_timestep 26\n"
169 for (
int i=1; i<argc; i++)
175 GetPot command_line (argc, argv);
184 const bool read_solution = command_line.search(
"-read_solution");
191 unsigned int init_timestep = 0;
195 if (command_line.search(
"-init_timestep"))
196 init_timestep = command_line.next(0);
201 libmesh_error_msg(
"ERROR: Initial timestep not specified!");
206 unsigned int n_timesteps = 0;
209 if (command_line.search(
"-n_timesteps"))
210 n_timesteps = command_line.next(0);
212 libmesh_error_msg(
"ERROR: Number of timesteps not specified");
216 #ifdef LIBMESH_HAVE_FPARSER
217 const bool have_expression = command_line.search(
"-exact_solution");
219 const bool have_expression =
false;
222 parsed_solution = libmesh_make_unique<ParsedFunction<Number>>(command_line.next(std::string()));
225 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
241 (
"Convection-Diffusion");
247 DofMap & dof_map = system.get_dof_map();
281 unsigned int n_refinements = 5;
282 if (command_line.search(
"-n_refinements"))
283 n_refinements = command_line.next(0);
294 system.add_variable (
"u",
FIRST);
297 system.attach_init_function (
init_cd);
309 equation_systems.read(
"saved_solution.xdr",
DECODE);
314 equation_systems.init ();
316 equation_systems.reinit ();
322 libMesh::out <<
"Initial H1 norm = " << H1norm << std::endl << std::endl;
325 equation_systems.print_info();
327 equation_systems.parameters.set<
unsigned int>
328 (
"linear solver maximum iterations") = 250;
329 equation_systems.parameters.set<
Real>
335 #ifdef LIBMESH_HAVE_GMV
339 #ifdef LIBMESH_HAVE_EXODUS_API
347 #ifdef LIBMESH_HAVE_GMV
351 #ifdef LIBMESH_HAVE_EXODUS_API
367 equation_systems.parameters.set<
Real>(
"diffusivity") = 0.01;
373 const Real dt = 0.025;
374 system.time = init_timestep*dt;
383 for (
unsigned int t_step=init_timestep; t_step<(init_timestep+n_timesteps); t_step++)
389 equation_systems.parameters.set<
Real> (
"time") = system.time;
390 equation_systems.parameters.set<
Real> (
"dt") = dt;
402 << std::setprecision(3)
424 unsigned int max_r_steps = 1;
425 if (command_line.search(
"-max_r_steps"))
426 max_r_steps = command_line.next(0);
429 for (
unsigned int r_step=0; r_step<max_r_steps+1; r_step++)
435 H1norm = system.calculate_norm(*system.solution,
SystemNorm(
H1));
439 if (r_step+1 <= max_r_steps)
488 equation_systems.reinit ();
493 unsigned int output_freq = 10;
494 if (command_line.search(
"-output_freq"))
495 output_freq = command_line.next(0);
498 if ((t_step+1)%output_freq == 0)
500 #ifdef LIBMESH_HAVE_GMV
510 #ifdef LIBMESH_HAVE_EXODUS_API
524 H1norm = system.calculate_norm(*system.solution,
SystemNorm(
H1));
526 libMesh::out <<
"Final H1 norm = " << H1norm << std::endl << std::endl;
529 equation_systems.write(
"saved_solution.xdr",
ENCODE);
530 #ifdef LIBMESH_HAVE_GMV
534 #ifdef LIBMESH_HAVE_EXODUS_API
539 #endif // #ifndef LIBMESH_ENABLE_AMR
549 const std::string & libmesh_dbg_var(system_name))
553 libmesh_assert_equal_to (system_name,
"Convection-Diffusion");
574 #ifdef LIBMESH_ENABLE_AMR
576 const std::string & libmesh_dbg_var(system_name))
580 libmesh_assert_equal_to (system_name,
"Convection-Diffusion");
594 FEType fe_type = system.variable_type(0);
605 QGauss qrule (
dim, fe_type.default_quadrature_order());
606 QGauss qface (
dim-1, fe_type.default_quadrature_order());
609 fe->attach_quadrature_rule (&qrule);
610 fe_face->attach_quadrature_rule (&qface);
615 const std::vector<Real> & JxW = fe->get_JxW();
618 const std::vector<std::vector<Real>> & phi = fe->get_phi();
622 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
628 const DofMap & dof_map = system.get_dof_map();
640 std::vector<dof_id_type> dof_indices;
647 const Real diffusivity =
671 const unsigned int n_dofs =
672 cast_int<unsigned int>(dof_indices.size());
673 libmesh_assert_equal_to (n_dofs, phi.size());
681 Ke.
resize (n_dofs, n_dofs);
691 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
698 for (
unsigned int l=0; l != n_dofs; l++)
700 u_old += phi[l][qp]*system.
old_solution (dof_indices[l]);
709 for (
unsigned int i=0; i != n_dofs; i++)
719 (grad_u_old*velocity)*phi[i][qp] +
722 diffusivity*(grad_u_old*dphi[i][qp]))
725 for (
unsigned int j=0; j != n_dofs; j++)
730 phi[i][qp]*phi[j][qp] +
733 (velocity*dphi[j][qp])*phi[i][qp] +
735 diffusivity*(dphi[i][qp]*dphi[j][qp]))
755 system.matrix->add_matrix (Ke, dof_indices);
756 system.rhs->add_vector (Fe, dof_indices);
761 #endif // #ifdef LIBMESH_ENABLE_AMR
768 std::ostringstream oss;
771 oss << std::setw(3) << std::setfill(
'0') << number;
bool refine_and_coarsen_elements()
Refines and coarsens user-requested elements.
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 void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
This function uses the derived class's jump error estimate formula to estimate the error on each cell...
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
virtual void write(const std::string &name)=0
PeriodicBoundaries * get_periodic_boundaries()
The libMesh namespace provides an interface to certain functionality in the library.
const T_sys & get_system(const std::string &name) const
static const Real TOLERANCE
SolverPackage default_solver_package()
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...
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
std::ios_base::fmtflags flags() const
Get the associated format flags.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
void assemble_cd(EquationSystems &es, const std::string &system_name)
This class implements the Kelly error indicator which is based on the flux jumps between elements.
This is the MeshBase class.
std::string exodus_filename(unsigned number)
Number exact_value(const Point &p, const Parameters ¶meters, const std::string &, const std::string &)
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
A Point defines a location in LIBMESH_DIM dimensional Real space.
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
boundary_id_type myboundary
The boundary ID of this boundary and its counterpart.
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Real exact_solution(const Real x, const Real y, const Real t)
This is the exact solution that we are trying to obtain.
void add_periodic_boundary(const PeriodicBoundaryBase &periodic_boundary)
Adds a copy of the specified periodic boundary to the system.
unsigned int & max_h_level()
The max_h_level is the greatest refinement level an element should reach.
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 ...
The definition of a periodic boundary.
void resize(const unsigned int n)
Resize the vector.
Real & coarsen_fraction()
The coarsen_fraction sets either a desired target or a desired maximum number of elements to flag for...
Real & refine_fraction()
The refine_fraction sets either a desired target or a desired maximum number of elements to flag for ...
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.
int main(int argc, char **argv)
Number old_solution(const dof_id_type global_dof_number) const
void set_periodic_boundaries_ptr(PeriodicBoundaries *pb_ptr)
Sets the PeriodicBoundaries pointer.
bool use_unweighted_quadrature_rules
This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweigh...
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
void init_cd(EquationSystems &es, const std::string &system_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
boundary_id_type pairedboundary
std::unique_ptr< NumericVector< Number > > old_local_solution
All the values I need to compute my contribution to the simulation at hand.
void flag_elements_by_error_fraction(const ErrorVector &error_per_cell, const Real refine_fraction=0.3, const Real coarsen_fraction=0.0, const unsigned int max_level=libMesh::invalid_uint)
Flags elements for coarsening and refinement based on the computed error passed in error_per_cell.
std::unique_ptr< FunctionBase< Number > > parsed_solution
const T & get(const std::string &) const
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
This class provides the ability to map between arbitrary, user-defined strings and several data types...
Parameters parameters
Data structure holding arbitrary parameters.