libMesh
Functions | Variables
adaptivity_ex5.C File Reference

Go to the source code of this file.

Functions

void assemble_cd (EquationSystems &es, const std::string &system_name)
 
void init_cd (EquationSystems &es, const std::string &system_name)
 
Real exact_solution (const Real x, const Real y, const Real t)
 This is the exact solution that we are trying to obtain. More...
 
Number exact_value (const Point &p, const Parameters &parameters, const std::string &, const std::string &)
 
std::string exodus_filename (unsigned number)
 
int main (int argc, char **argv)
 
void init_cd (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 
void assemble_cd (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Variables

UniquePtr< FunctionBase< Number > > parsed_solution
 

Function Documentation

void assemble_cd ( EquationSystems es,
const std::string &  system_name 
)

Definition at line 292 of file transient_ex1.C.

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(), libmesh_nullptr, 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().

294 {
295  // Ignore unused parameter warnings when !LIBMESH_ENABLE_AMR.
296  libmesh_ignore(es);
297  libmesh_ignore(system_name);
298 
299 #ifdef LIBMESH_ENABLE_AMR
300  // It is a good idea to make sure we are assembling
301  // the proper system.
302  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
303 
304  // Get a constant reference to the mesh object.
305  const MeshBase & mesh = es.get_mesh();
306 
307  // The dimension that we are running
308  const unsigned int dim = mesh.mesh_dimension();
309 
310  // Get a reference to the Convection-Diffusion system object.
312  es.get_system<TransientLinearImplicitSystem> ("Convection-Diffusion");
313 
314  // Get a constant reference to the Finite Element type
315  // for the first (and only) variable in the system.
316  FEType fe_type = system.variable_type(0);
317 
318  // Build a Finite Element object of the specified type. Since the
319  // FEBase::build() member dynamically creates memory we will
320  // store the object as a UniquePtr<FEBase>. This can be thought
321  // of as a pointer that will clean up after itself.
322  UniquePtr<FEBase> fe (FEBase::build(dim, fe_type));
323  UniquePtr<FEBase> fe_face (FEBase::build(dim, fe_type));
324 
325  // A Gauss quadrature rule for numerical integration.
326  // Let the FEType object decide what order rule is appropriate.
327  QGauss qrule (dim, fe_type.default_quadrature_order());
328  QGauss qface (dim-1, fe_type.default_quadrature_order());
329 
330  // Tell the finite element object to use our quadrature rule.
331  fe->attach_quadrature_rule (&qrule);
332  fe_face->attach_quadrature_rule (&qface);
333 
334  // Here we define some references to cell-specific data that
335  // will be used to assemble the linear system. We will start
336  // with the element Jacobian * quadrature weight at each integration point.
337  const std::vector<Real> & JxW = fe->get_JxW();
338  const std::vector<Real> & JxW_face = fe_face->get_JxW();
339 
340  // The element shape functions evaluated at the quadrature points.
341  const std::vector<std::vector<Real>> & phi = fe->get_phi();
342  const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
343 
344  // The element shape function gradients evaluated at the quadrature
345  // points.
346  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
347 
348  // The XY locations of the quadrature points used for face integration
349  const std::vector<Point> & qface_points = fe_face->get_xyz();
350 
351  // A reference to the DofMap object for this system. The DofMap
352  // object handles the index translation from node and element numbers
353  // to degree of freedom numbers. We will talk more about the DofMap
354  // in future examples.
355  const DofMap & dof_map = system.get_dof_map();
356 
357  // Define data structures to contain the element matrix
358  // and right-hand-side vector contribution. Following
359  // basic finite element terminology we will denote these
360  // "Ke" and "Fe".
363 
364  // This vector will hold the degree of freedom indices for
365  // the element. These define where in the global system
366  // the element degrees of freedom get mapped.
367  std::vector<dof_id_type> dof_indices;
368 
369  // Here we extract the velocity & parameters that we put in the
370  // EquationSystems object.
371  const RealVectorValue velocity =
372  es.parameters.get<RealVectorValue> ("velocity");
373 
374  const Real dt = es.parameters.get<Real> ("dt");
375 
376  // Now we will loop over all the elements in the mesh that
377  // live on the local processor. We will compute the element
378  // matrix and right-hand-side contribution. Since the mesh
379  // will be refined we want to only consider the ACTIVE elements,
380  // hence we use a variant of the active_elem_iterator.
381  for (const auto & elem : mesh.active_local_element_ptr_range())
382  {
383  // Get the degree of freedom indices for the
384  // current element. These define where in the global
385  // matrix and right-hand-side this element will
386  // contribute to.
387  dof_map.dof_indices (elem, dof_indices);
388 
389  // Compute the element-specific data for the current
390  // element. This involves computing the location of the
391  // quadrature points (q_point) and the shape functions
392  // (phi, dphi) for the current element.
393  fe->reinit (elem);
394 
395  // Zero the element matrix and right-hand side before
396  // summing them. We use the resize member here because
397  // the number of degrees of freedom might have changed from
398  // the last element. Note that this will be the case if the
399  // element type is different (i.e. the last element was a
400  // triangle, now we are on a quadrilateral).
401  Ke.resize (dof_indices.size(),
402  dof_indices.size());
403 
404  Fe.resize (dof_indices.size());
405 
406  // Now we will build the element matrix and right-hand-side.
407  // Constructing the RHS requires the solution and its
408  // gradient from the previous timestep. This myst be
409  // calculated at each quadrature point by summing the
410  // solution degree-of-freedom values by the appropriate
411  // weight functions.
412  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
413  {
414  // Values to hold the old solution & its gradient.
415  Number u_old = 0.;
416  Gradient grad_u_old;
417 
418  // Compute the old solution & its gradient.
419  for (std::size_t l=0; l<phi.size(); l++)
420  {
421  u_old += phi[l][qp]*system.old_solution (dof_indices[l]);
422 
423  // This will work,
424  // grad_u_old += dphi[l][qp]*system.old_solution (dof_indices[l]);
425  // but we can do it without creating a temporary like this:
426  grad_u_old.add_scaled (dphi[l][qp], system.old_solution (dof_indices[l]));
427  }
428 
429  // Now compute the element matrix and RHS contributions.
430  for (std::size_t i=0; i<phi.size(); i++)
431  {
432  // The RHS contribution
433  Fe(i) += JxW[qp]*(
434  // Mass matrix term
435  u_old*phi[i][qp] +
436  -.5*dt*(
437  // Convection term
438  // (grad_u_old may be complex, so the
439  // order here is important!)
440  (grad_u_old*velocity)*phi[i][qp] +
441 
442  // Diffusion term
443  0.01*(grad_u_old*dphi[i][qp]))
444  );
445 
446  for (std::size_t j=0; j<phi.size(); j++)
447  {
448  // The matrix contribution
449  Ke(i,j) += JxW[qp]*(
450  // Mass-matrix
451  phi[i][qp]*phi[j][qp] +
452 
453  .5*dt*(
454  // Convection term
455  (velocity*dphi[j][qp])*phi[i][qp] +
456 
457  // Diffusion term
458  0.01*(dphi[i][qp]*dphi[j][qp]))
459  );
460  }
461  }
462  }
463 
464  // At this point the interior element integration has
465  // been completed. However, we have not yet addressed
466  // boundary conditions. For this example we will only
467  // consider simple Dirichlet boundary conditions imposed
468  // via the penalty method.
469  //
470  // The following loops over the sides of the element.
471  // If the element has no neighbor on a side then that
472  // side MUST live on a boundary of the domain.
473  {
474  // The penalty value.
475  const Real penalty = 1.e10;
476 
477  // The following loops over the sides of the element.
478  // If the element has no neighbor on a side then that
479  // side MUST live on a boundary of the domain.
480  for (auto s : elem->side_index_range())
481  if (elem->neighbor_ptr(s) == libmesh_nullptr)
482  {
483  fe_face->reinit(elem, s);
484 
485  for (unsigned int qp=0; qp<qface.n_points(); qp++)
486  {
487  const Number value = exact_solution (qface_points[qp](0),
488  qface_points[qp](1),
489  system.time);
490 
491  // RHS contribution
492  for (std::size_t i=0; i<psi.size(); i++)
493  Fe(i) += penalty*JxW_face[qp]*value*psi[i][qp];
494 
495  // Matrix contribution
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];
499  }
500  }
501  }
502 
503  // If this assembly program were to be used on an adaptive mesh,
504  // we would have to apply any hanging node constraint equations
505  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
506 
507  // The element matrix and right-hand-side are now built
508  // for this element. Add them to the global matrix and
509  // right-hand-side vector. The SparseMatrix::add_matrix()
510  // and NumericVector::add_vector() members do this for us.
511  system.matrix->add_matrix (Ke, dof_indices);
512  system.rhs->add_vector (Fe, dof_indices);
513  }
514 
515  // That concludes the system matrix assembly routine.
516 #endif // #ifdef LIBMESH_ENABLE_AMR
517 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
void add_scaled(const TypeVector< T2 > &, const T)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:624
unsigned int dim
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:350
MeshBase & mesh
const class libmesh_nullptr_t libmesh_nullptr
This class provides a specific system class.
const T & get(const std::string &) const
Definition: parameters.h:430
This is the MeshBase class.
Definition: mesh_base.h:68
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:167
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.
Definition: dof_map.h:1806
Real exact_solution(const Real x, const Real y, const Real t)
This is the exact solution that we are trying to obtain.
void libmesh_ignore(const T &)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
unsigned int mesh_dimension() const
Definition: mesh_base.C:148
static const bool value
Definition: xdr_io.C:108
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:776
This class implements specific orders of Gauss quadrature.
Parameters parameters
Data structure holding arbitrary parameters.
const MeshBase & get_mesh() const
Number old_solution(const dof_id_type global_dof_number) const
const T_sys & get_system(const std::string &name) const
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.
Definition: dof_map.C:1917
void assemble_cd ( EquationSystems es,
const std::string &  libmesh_dbg_varsystem_name 
)

Definition at line 568 of file adaptivity_ex5.C.

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(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::TransientSystem< Base >::old_solution(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseVector< T >::resize(), and libMesh::DenseMatrix< T >::resize().

570 {
571  // It is a good idea to make sure we are assembling
572  // the proper system.
573  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
574 
575  // Get a constant reference to the mesh object.
576  const MeshBase & mesh = es.get_mesh();
577 
578  // The dimension that we are running
579  const unsigned int dim = mesh.mesh_dimension();
580 
581  // Get a reference to the Convection-Diffusion system object.
583  es.get_system<TransientLinearImplicitSystem> ("Convection-Diffusion");
584 
585  // Get the Finite Element type for the first (and only)
586  // variable in the system.
587  FEType fe_type = system.variable_type(0);
588 
589  // Build a Finite Element object of the specified type. Since the
590  // FEBase::build() member dynamically creates memory we will
591  // store the object as a UniquePtr<FEBase>. This can be thought
592  // of as a pointer that will clean up after itself.
593  UniquePtr<FEBase> fe (FEBase::build(dim, fe_type));
594  UniquePtr<FEBase> fe_face (FEBase::build(dim, fe_type));
595 
596  // A Gauss quadrature rule for numerical integration.
597  // Let the FEType object decide what order rule is appropriate.
598  QGauss qrule (dim, fe_type.default_quadrature_order());
599  QGauss qface (dim-1, fe_type.default_quadrature_order());
600 
601  // Tell the finite element object to use our quadrature rule.
602  fe->attach_quadrature_rule (&qrule);
603  fe_face->attach_quadrature_rule (&qface);
604 
605  // Here we define some references to cell-specific data that
606  // will be used to assemble the linear system. We will start
607  // with the element Jacobian * quadrature weight at each integration point.
608  const std::vector<Real> & JxW = fe->get_JxW();
609 
610  // The element shape functions evaluated at the quadrature points.
611  const std::vector<std::vector<Real>> & phi = fe->get_phi();
612 
613  // The element shape function gradients evaluated at the quadrature
614  // points.
615  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
616 
617  // A reference to the DofMap object for this system. The DofMap
618  // object handles the index translation from node and element numbers
619  // to degree of freedom numbers. We will talk more about the DofMap
620  // in future examples.
621  const DofMap & dof_map = system.get_dof_map();
622 
623  // Define data structures to contain the element matrix
624  // and right-hand-side vector contribution. Following
625  // basic finite element terminology we will denote these
626  // "Ke" and "Fe".
629 
630  // This vector will hold the degree of freedom indices for
631  // the element. These define where in the global system
632  // the element degrees of freedom get mapped.
633  std::vector<dof_id_type> dof_indices;
634 
635  // Here we extract the velocity & parameters that we put in the
636  // EquationSystems object.
637  const RealVectorValue velocity =
638  es.parameters.get<RealVectorValue> ("velocity");
639 
640  const Real diffusivity =
641  es.parameters.get<Real> ("diffusivity");
642 
643  const Real dt = es.parameters.get<Real> ("dt");
644 
645  // Now we will loop over all the elements in the mesh that
646  // live on the local processor. We will compute the element
647  // matrix and right-hand-side contribution. Since the mesh
648  // will be refined we want to only consider the ACTIVE elements,
649  // hence we use a variant of the active_elem_iterator.
650  for (const auto & elem : mesh.active_local_element_ptr_range())
651  {
652  // Get the degree of freedom indices for the
653  // current element. These define where in the global
654  // matrix and right-hand-side this element will
655  // contribute to.
656  dof_map.dof_indices (elem, dof_indices);
657 
658  // Compute the element-specific data for the current
659  // element. This involves computing the location of the
660  // quadrature points (q_point) and the shape functions
661  // (phi, dphi) for the current element.
662  fe->reinit (elem);
663 
664  // Zero the element matrix and right-hand side before
665  // summing them. We use the resize member here because
666  // the number of degrees of freedom might have changed from
667  // the last element. Note that this will be the case if the
668  // element type is different (i.e. the last element was a
669  // triangle, now we are on a quadrilateral).
670  Ke.resize (dof_indices.size(),
671  dof_indices.size());
672 
673  Fe.resize (dof_indices.size());
674 
675  // Now we will build the element matrix and right-hand-side.
676  // Constructing the RHS requires the solution and its
677  // gradient from the previous timestep. This myst be
678  // calculated at each quadrature point by summing the
679  // solution degree-of-freedom values by the appropriate
680  // weight functions.
681  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
682  {
683  // Values to hold the old solution & its gradient.
684  Number u_old = 0.;
685  Gradient grad_u_old;
686 
687  // Compute the old solution & its gradient.
688  for (std::size_t l=0; l<phi.size(); l++)
689  {
690  u_old += phi[l][qp]*system.old_solution (dof_indices[l]);
691 
692  // This will work,
693  // grad_u_old += dphi[l][qp]*system.old_solution (dof_indices[l]);
694  // but we can do it without creating a temporary like this:
695  grad_u_old.add_scaled (dphi[l][qp], system.old_solution (dof_indices[l]));
696  }
697 
698  // Now compute the element matrix and RHS contributions.
699  for (std::size_t i=0; i<phi.size(); i++)
700  {
701  // The RHS contribution
702  Fe(i) += JxW[qp]*(
703  // Mass matrix term
704  u_old*phi[i][qp] +
705  -.5*dt*(
706  // Convection term
707  // (grad_u_old may be complex, so the
708  // order here is important!)
709  (grad_u_old*velocity)*phi[i][qp] +
710 
711  // Diffusion term
712  diffusivity*(grad_u_old*dphi[i][qp]))
713  );
714 
715  for (std::size_t j=0; j<phi.size(); j++)
716  {
717  // The matrix contribution
718  Ke(i,j) += JxW[qp]*(
719  // Mass-matrix
720  phi[i][qp]*phi[j][qp] +
721  .5*dt*(
722  // Convection term
723  (velocity*dphi[j][qp])*phi[i][qp] +
724  // Diffusion term
725  diffusivity*(dphi[i][qp]*dphi[j][qp]))
726  );
727  }
728  }
729  }
730 
731  // We have now built the element matrix and RHS vector in terms
732  // of the element degrees of freedom. However, it is possible
733  // that some of the element DOFs are constrained to enforce
734  // solution continuity, i.e. they are not really "free". We need
735  // to constrain those DOFs in terms of non-constrained DOFs to
736  // ensure a continuous solution. The
737  // DofMap::constrain_element_matrix_and_vector() method does
738  // just that.
739  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
740 
741  // The element matrix and right-hand-side are now built
742  // for this element. Add them to the global matrix and
743  // right-hand-side vector. The SparseMatrix::add_matrix()
744  // and NumericVector::add_vector() members do this for us.
745  system.matrix->add_matrix (Ke, dof_indices);
746  system.rhs->add_vector (Fe, dof_indices);
747 
748  }
749  // Finished computing the system matrix and right-hand side.
750 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
void add_scaled(const TypeVector< T2 > &, const T)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:624
unsigned int dim
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:350
MeshBase & mesh
This class provides a specific system class.
const T & get(const std::string &) const
Definition: parameters.h:430
This is the MeshBase class.
Definition: mesh_base.h:68
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:167
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.
Definition: dof_map.h:1806
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
unsigned int mesh_dimension() const
Definition: mesh_base.C:148
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:776
This class implements specific orders of Gauss quadrature.
Parameters parameters
Data structure holding arbitrary parameters.
const MeshBase & get_mesh() const
Number old_solution(const dof_id_type global_dof_number) const
const T_sys & get_system(const std::string &name) const
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.
Definition: dof_map.C:1917
Real exact_solution ( const Real  x,
const Real  y,
const Real  t 
)

This is the exact solution that we are trying to obtain.

We will solve

  • (u_xx + u_yy) = f

and take a finite difference approximation using this function to get f. This is the well-known "method of manufactured solutions".

Definition at line 43 of file exact_solution.C.

Referenced by exact_value().

46 {
47  static const Real pi = acos(-1.);
48 
49  return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
50 }
PetscErrorCode Vec x
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real pi
.
Definition: libmesh.h:172
Number exact_value ( const Point p,
const Parameters parameters,
const std::string &  ,
const std::string &   
)

Definition at line 104 of file adaptivity_ex5.C.

References exact_solution(), libMesh::Parameters::get(), and libMesh::Real.

Referenced by init_cd().

108 {
109  return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
110 }
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
Definition: parameters.h:430
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::string exodus_filename ( unsigned  number)

Definition at line 756 of file adaptivity_ex5.C.

Referenced by main().

757 {
758  std::ostringstream oss;
759 
760  oss << "out_";
761  oss << std::setw(3) << std::setfill('0') << number;
762  oss << ".e";
763 
764  return oss.str();
765 }
void init_cd ( EquationSystems es,
const std::string &  system_name 
)

Referenced by main().

void init_cd ( EquationSystems es,
const std::string &  libmesh_dbg_varsystem_name 
)

Definition at line 541 of file adaptivity_ex5.C.

References exact_value(), libMesh::EquationSystems::get_system(), if(), libmesh_nullptr, libMesh::EquationSystems::parameters, parsed_solution, libMesh::Real, and libMesh::Parameters::set().

543 {
544  // It is a good idea to make sure we are initializing
545  // the proper system.
546  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
547 
548  // Get a reference to the Convection-Diffusion system object.
550  es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
551 
552  // Project initial conditions at time 0
553  es.parameters.set<Real> ("time") = system.time = 0;
554 
555  if (parsed_solution.get())
556  system.project_solution(parsed_solution.get(), libmesh_nullptr);
557  else
558  system.project_solution(exact_value, libmesh_nullptr, es.parameters);
559 }
UniquePtr< FunctionBase< Number > > parsed_solution
const class libmesh_nullptr_t libmesh_nullptr
This class provides a specific system class.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Definition: parameters.h:469
Number exact_value(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
Parameters parameters
Data structure holding arbitrary parameters.
const T_sys & get_system(const std::string &name) const
if(!eq) SETERRQ2(((PetscObject) dm) -> comm, PETSC_ERR_ARG_WRONG,"DM of type %s, not of type %s",((PetscObject) dm) ->type, DMLIBMESH)
int main ( int  argc,
char **  argv 
)

Definition at line 127 of file adaptivity_ex5.C.

References libMesh::DofMap::add_periodic_boundary(), assemble_cd(), libMesh::MeshTools::Generation::build_square(), libMesh::MeshRefinement::coarsen_fraction(), libMesh::LibMeshInit::comm(), libMesh::DECODE, libMesh::default_solver_package(), libMesh::ENCODE, libMesh::JumpErrorEstimator::estimate_error(), exodus_filename(), libMesh::FIRST, libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::BasicOStreamProxy< charT, traits >::flags(), libMesh::DofMap::get_periodic_boundaries(), libMesh::H1, libMesh::TriangleWrapper::init(), init_cd(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::MeshRefinement::max_h_level(), mesh, libMesh::PeriodicBoundaryBase::myboundary, libMesh::TransientSystem< Base >::old_local_solution, libMesh::out, libMesh::PeriodicBoundaryBase::pairedboundary, parsed_solution, libMesh::MeshBase::print_info(), libMesh::QUAD4, libMesh::MeshBase::read(), libMesh::Real, libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::MeshRefinement::refine_fraction(), libMesh::MeshRefinement::set_periodic_boundaries_ptr(), libMesh::TOLERANCE, libMesh::TRILINOS_SOLVERS, libMesh::MeshRefinement::uniformly_refine(), libMesh::MeshBase::write(), and libMesh::MeshOutput< MT >::write_equation_systems().

128 {
129  // Initialize libMesh.
130  LibMeshInit init (argc, argv);
131 
132  // This example requires a linear solver package.
133  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
134  "--enable-petsc, --enable-trilinos, or --enable-eigen");
135 
136  // Skip this 2D example if libMesh was compiled as 1D-only.
137  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
138 
139 #if !defined(LIBMESH_ENABLE_AMR)
140  libmesh_example_requires(false, "--enable-amr");
141 #elif !defined(LIBMESH_HAVE_XDR)
142  // We use XDR support in our output here
143  libmesh_example_requires(false, "--enable-xdr");
144 #elif !defined(LIBMESH_ENABLE_PERIODIC)
145  libmesh_example_requires(false, "--enable-periodic");
146 #elif (LIBMESH_DOF_ID_BYTES == 8)
147  libmesh_example_requires(false, "--with-dof-id-bytes=4");
148 #else
149 
150  // Our Trilinos interface does not yet support adaptive transient
151  // problems
152  libmesh_example_requires(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc");
153 
154  // Brief message to the user regarding the program name
155  // and command line arguments.
156 
157  // Use commandline parameter to specify if we are to
158  // read in an initial solution or generate it ourself
159  libMesh::out << "Usage:\n"
160  <<"\t " << argv[0] << " -init_timestep 0\n"
161  << "OR\n"
162  <<"\t " << argv[0] << " -read_solution -init_timestep 26\n"
163  << std::endl;
164 
165  libMesh::out << "Running: " << argv[0];
166 
167  for (int i=1; i<argc; i++)
168  libMesh::out << " " << argv[i];
169 
170  libMesh::out << std::endl << std::endl;
171 
172  // Create a GetPot object to parse the command line
173  GetPot command_line (argc, argv);
174 
175 
176  // This boolean value is obtained from the command line, it is true
177  // if the flag "-read_solution" is present, false otherwise.
178  // It indicates whether we are going to read in
179  // the mesh and solution files "saved_mesh.xda" and "saved_solution.xda"
180  // or whether we are going to start from scratch by just reading
181  // "mesh.xda"
182  const bool read_solution = command_line.search("-read_solution");
183 
184  // This value is also obtained from the commandline and it specifies the
185  // initial value for the t_step looping variable. We must
186  // distinguish between the two cases here, whether we read in the
187  // solution or we started from scratch, so that we do not overwrite the
188  // gmv output files.
189  unsigned int init_timestep = 0;
190 
191  // Search the command line for the "init_timestep" flag and if it is
192  // present, set init_timestep accordingly.
193  if (command_line.search("-init_timestep"))
194  init_timestep = command_line.next(0);
195  else
196  {
197  // This handy function will print the file name, line number,
198  // specified message, and then throw an exception.
199  libmesh_error_msg("ERROR: Initial timestep not specified!");
200  }
201 
202  // This value is also obtained from the command line, and specifies
203  // the number of time steps to take.
204  unsigned int n_timesteps = 0;
205 
206  // Again do a search on the command line for the argument
207  if (command_line.search("-n_timesteps"))
208  n_timesteps = command_line.next(0);
209  else
210  libmesh_error_msg("ERROR: Number of timesteps not specified");
211 
212  // The user can specify a different exact solution on the command
213  // line, if we have an expression parser compiled in
214 #ifdef LIBMESH_HAVE_FPARSER
215  const bool have_expression = command_line.search("-exact_solution");
216 #else
217  const bool have_expression = false;
218 #endif
219  if (have_expression)
220  parsed_solution.reset
221  (new ParsedFunction<Number>(command_line.next(std::string())));
222 
223  // Skip this 2D example if libMesh was compiled as 1D-only.
224  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
225 
226  // Create a new mesh on the default MPI communicator.
227  // DistributedMesh currently has a bug which is triggered by this
228  // example.
229  ReplicatedMesh mesh(init.comm());
230 
231  // Create an equation systems object.
232  EquationSystems equation_systems (mesh);
233  MeshRefinement mesh_refinement (mesh);
234 
235  // Declare the system and its variables.
236  // Begin by creating a transient system
237  // named "Convection-Diffusion".
239  equation_systems.add_system<TransientLinearImplicitSystem>
240  ("Convection-Diffusion");
241 
242  // Give the system a pointer to the assembly function.
243  system.attach_assemble_function (assemble_cd);
244 
245  // Creating and attaching Periodic Boundaries
246  DofMap & dof_map = system.get_dof_map();
247 
248  // Create a boundary periodic with one displaced 2.0 in the x
249  // direction
250  PeriodicBoundary horz(RealVectorValue(2.0, 0., 0.));
251 
252  // Connect boundary ids 3 and 1 with it
253  horz.myboundary = 3;
254  horz.pairedboundary = 1;
255 
256  // Add it to the PeriodicBoundaries
257  dof_map.add_periodic_boundary(horz);
258 
259  // Create a boundary periodic with one displaced 2.0 in the y
260  // direction
261  PeriodicBoundary vert(RealVectorValue(0., 2.0, 0.));
262 
263  // Connect boundary ids 0 and 2 with it
264  vert.myboundary = 0;
265  vert.pairedboundary = 2;
266 
267  // Add it to the PeriodicBoundaries
268  dof_map.add_periodic_boundary(vert);
269 
270  // Next build or read the mesh. We do this only *after* generating
271  // periodic boundaries; otherwise a DistributedMesh won't know to
272  // retain periodic neighbor elements.
273 
274  // First we process the case where we do not read in the solution
275  if (!read_solution)
276  {
277  MeshTools::Generation::build_square(mesh, 2, 2, 0., 2., 0., 2., QUAD4);
278 
279  // Again do a search on the command line for an argument
280  unsigned int n_refinements = 5;
281  if (command_line.search("-n_refinements"))
282  n_refinements = command_line.next(0);
283 
284  // Uniformly refine the mesh 5 times
285  if (!read_solution)
286  mesh_refinement.uniformly_refine (n_refinements);
287 
288  // Print information about the mesh to the screen.
289  mesh.print_info();
290 
291  // Adds the variable "u" to "Convection-Diffusion". "u"
292  // will be approximated using first-order approximation.
293  system.add_variable ("u", FIRST);
294 
295  // Give the system a pointer to the initialization function.
296  system.attach_init_function (init_cd);
297  }
298  // Otherwise we read in the solution and mesh
299  else
300  {
301  // Read in the mesh stored in "saved_mesh.xda"
302  mesh.read("saved_mesh.xdr");
303 
304  // Print information about the mesh to the screen.
305  mesh.print_info();
306 
307  // Read in the solution stored in "saved_solution.xda"
308  equation_systems.read("saved_solution.xdr", DECODE);
309  }
310 
311  // Initialize the data structures for the equation system.
312  if (!read_solution)
313  equation_systems.init ();
314  else
315  equation_systems.reinit ();
316 
317  // Print out the H1 norm of the initialized or saved solution, for
318  // verification purposes:
319  Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
320 
321  libMesh::out << "Initial H1 norm = " << H1norm << std::endl << std::endl;
322 
323  // Prints information about the system to the screen.
324  equation_systems.print_info();
325 
326  equation_systems.parameters.set<unsigned int>
327  ("linear solver maximum iterations") = 250;
328  equation_systems.parameters.set<Real>
329  ("linear solver tolerance") = TOLERANCE;
330 
331  if (!read_solution)
332  {
333  // Write out the initial condition
334 #ifdef LIBMESH_HAVE_GMV
335  GMVIO(mesh).write_equation_systems ("out.gmv.000",
336  equation_systems);
337 #endif
338 #ifdef LIBMESH_HAVE_EXODUS_API
340  equation_systems);
341 #endif
342  }
343  else
344  {
345  // Write out the solution that was read in
346 #ifdef LIBMESH_HAVE_GMV
347  GMVIO(mesh).write_equation_systems ("solution_read_in.gmv",
348  equation_systems);
349 #endif
350 #ifdef LIBMESH_HAVE_EXODUS_API
351  ExodusII_IO(mesh).write_equation_systems ("solution_read_in.e",
352  equation_systems);
353 #endif
354  }
355 
356 
357  // The Convection-Diffusion system requires that we specify
358  // the flow velocity. We will specify it as a RealVectorValue
359  // data type and then use the Parameters object to pass it to
360  // the assemble function.
361  equation_systems.parameters.set<RealVectorValue>("velocity") =
362  RealVectorValue (0.8, 0.8);
363 
364  // The Convection-Diffusion system also requires a specified
365  // diffusivity. We use an isotropic (hence Real) value.
366  equation_systems.parameters.set<Real>("diffusivity") = 0.01;
367 
368  // Solve the system "Convection-Diffusion". This will be done by
369  // looping over the specified time interval and calling the
370  // solve() member at each time step. This will assemble the
371  // system and call the linear solver.
372  const Real dt = 0.025;
373  system.time = init_timestep*dt;
374 
375  // Tell the MeshRefinement object about the periodic boundaries so
376  // that it can get heuristics like level-one conformity and
377  // unrefined island elimination right.
378  mesh_refinement.set_periodic_boundaries_ptr(dof_map.get_periodic_boundaries());
379 
380  // We do 25 timesteps both before and after writing out the
381  // intermediate solution
382  for (unsigned int t_step=init_timestep; t_step<(init_timestep+n_timesteps); t_step++)
383  {
384  // Increment the time counter, set the time and the
385  // time step size as parameters in the EquationSystem.
386  system.time += dt;
387 
388  equation_systems.parameters.set<Real> ("time") = system.time;
389  equation_systems.parameters.set<Real> ("dt") = dt;
390 
391  // A pretty update message
392  libMesh::out << " Solving time step ";
393 
394  {
395  // Save flags to avoid polluting cout with custom precision values, etc.
396  std::ios_base::fmtflags os_flags = libMesh::out.flags();
397 
398  libMesh::out << t_step
399  << ", time="
400  << std::setw(6)
401  << std::setprecision(3)
402  << std::setfill('0')
403  << std::left
404  << system.time
405  << "..."
406  << std::endl;
407 
408  // Restore flags
409  libMesh::out.flags(os_flags);
410  }
411 
412  // At this point we need to update the old
413  // solution vector. The old solution vector
414  // will be the current solution vector from the
415  // previous time step. We will do this by extracting the
416  // system from the EquationSystems object and using
417  // vector assignment. Since only TransientLinearImplicitSystems
418  // (and systems derived from them) contain old solutions
419  // we need to specify the system type when we ask for it.
420  *system.old_local_solution = *system.current_local_solution;
421 
422  // The number of refinement steps per time step.
423  unsigned int max_r_steps = 1;
424  if (command_line.search("-max_r_steps"))
425  max_r_steps = command_line.next(0);
426 
427  // A refinement loop.
428  for (unsigned int r_step=0; r_step<max_r_steps+1; r_step++)
429  {
430  // Assemble & solve the linear system
431  system.solve();
432 
433  // Print out the H1 norm, for verification purposes:
434  H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
435  libMesh::out << "H1 norm = " << H1norm << std::endl;
436 
437  // Possibly refine the mesh
438  if (r_step+1 <= max_r_steps)
439  {
440  libMesh::out << " Refining the mesh..." << std::endl;
441 
442  // The ErrorVector is a particular StatisticsVector
443  // for computing error information on a finite element mesh.
444  ErrorVector error;
445 
446  // The ErrorEstimator class interrogates a finite element
447  // solution and assigns to each element a positive error value.
448  // This value is used for deciding which elements to refine
449  // and which to coarsen.
450  //ErrorEstimator* error_estimator = new KellyErrorEstimator;
451  KellyErrorEstimator error_estimator;
452 
453  // Compute the error for each active element using the provided
454  // flux_jump indicator. Note in general you will need to
455  // provide an error estimator specifically designed for your
456  // application.
457  error_estimator.estimate_error (system,
458  error);
459 
460  // This takes the error in error and decides which elements
461  // will be coarsened or refined. Any element within 20% of the
462  // maximum error on any element will be refined, and any
463  // element within 7% of the minimum error on any element might
464  // be coarsened. Note that the elements flagged for refinement
465  // will be refined, but those flagged for coarsening _might_ be
466  // coarsened.
467  mesh_refinement.refine_fraction() = 0.80;
468  mesh_refinement.coarsen_fraction() = 0.07;
469  mesh_refinement.max_h_level() = 5;
470  mesh_refinement.flag_elements_by_error_fraction (error);
471 
472  // This call actually refines and coarsens the flagged
473  // elements.
474  mesh_refinement.refine_and_coarsen_elements();
475 
476  // This call reinitializes the EquationSystems object for
477  // the newly refined mesh. One of the steps in the
478  // reinitialization is projecting the solution,
479  // old_solution, etc... vectors from the old mesh to
480  // the current one.
481  equation_systems.reinit ();
482  }
483  }
484 
485  // Again do a search on the command line for an argument
486  unsigned int output_freq = 10;
487  if (command_line.search("-output_freq"))
488  output_freq = command_line.next(0);
489 
490  // Output every 10 timesteps to file.
491  if ((t_step+1)%output_freq == 0)
492  {
493 #ifdef LIBMESH_HAVE_GMV
494  // std::ostringstream file_name;
495  // out << "out.gmv."
496  // << std::setw(3)
497  // << std::setfill('0')
498  // << std::right
499  // << t_step+1;
500  // GMVIO(mesh).write_equation_systems (file_name.str(),
501  // equation_systems);
502 #endif
503 #ifdef LIBMESH_HAVE_EXODUS_API
504  // So... if paraview is told to open a file called out.e.{N}, it automatically tries to
505  // open out.e.{N-1}, out.e.{N-2}, etc. If we name the file something else, we can work
506  // around that issue, but the right thing to do (for adaptive meshes) is to write a filename
507  // with the adaptation step into a separate file.
509  equation_systems);
510 #endif
511  }
512  }
513 
514  if (!read_solution)
515  {
516  // Print out the H1 norm of the saved solution, for verification purposes:
517  H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
518 
519  libMesh::out << "Final H1 norm = " << H1norm << std::endl << std::endl;
520 
521  mesh.write("saved_mesh.xdr");
522  equation_systems.write("saved_solution.xdr", ENCODE);
523 #ifdef LIBMESH_HAVE_GMV
524  GMVIO(mesh).write_equation_systems ("saved_solution.gmv",
525  equation_systems);
526 #endif
527 #ifdef LIBMESH_HAVE_EXODUS_API
528  ExodusII_IO(mesh).write_equation_systems ("saved_solution.e",
529  equation_systems);
530 #endif
531  }
532 #endif // #ifndef LIBMESH_ENABLE_AMR
533 
534  return 0;
535 }
void assemble_cd(EquationSystems &es, const std::string &system_name)
This is the EquationSystems class.
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
UniquePtr< NumericVector< Number > > old_local_solution
All the values I need to compute my contribution to the simulation at hand.
void add_periodic_boundary(const PeriodicBoundaryBase &periodic_boundary)
Adds a copy of the specified periodic boundary to the system.
A Function generated (via FParser) by parsing a mathematical expression.
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
UniquePtr< FunctionBase< Number > > parsed_solution
MeshBase & mesh
This class provides a specific system class.
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
Definition: system_norm.h:43
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
static const Real TOLERANCE
void init_cd(EquationSystems &es, const std::string &system_name)
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:62
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:46
The definition of a periodic boundary.
SolverPackage default_solver_package()
Definition: libmesh.C:995
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:167
This is the MeshRefinement class.
std::string exodus_filename(unsigned number)
PeriodicBoundaries * get_periodic_boundaries()
Definition: dof_map.h:1107
void init(triangulateio &t)
Initializes the fields of t to NULL/0 as necessary.
virtual void write(const std::string &name)=0
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=libmesh_nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
This class implements the Kelly error indicator which is based on the flux jumps between elements...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=libmesh_nullptr, bool estimate_parent_error=false) libmesh_override
This function uses the derived class&#39;s jump error estimate formula to estimate the error on each cell...
std::ios_base::fmtflags flags() const
Get the associated format flags.
virtual void read(const std::string &name, void *mesh_data=libmesh_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 print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:448

Variable Documentation

UniquePtr<FunctionBase<Number> > parsed_solution

Definition at line 114 of file adaptivity_ex5.C.

Referenced by init_cd(), and main().