libMesh
Functions
adaptivity_ex2.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 &)
 
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))
 

Function Documentation

◆ assemble_cd() [1/2]

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

Definition at line 296 of file transient_ex1.C.

Referenced by main().

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

◆ assemble_cd() [2/2]

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

Definition at line 480 of file adaptivity_ex2.C.

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

482 {
483  // It is a good idea to make sure we are assembling
484  // the proper system.
485  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
486 
487  // Get a constant reference to the mesh object.
488  const MeshBase & mesh = es.get_mesh();
489 
490  // The dimension that we are running
491  const unsigned int dim = mesh.mesh_dimension();
492 
493  // Get a reference to the Convection-Diffusion system object.
495  es.get_system<TransientLinearImplicitSystem> ("Convection-Diffusion");
496 
497  // Get the Finite Element type for the first (and only)
498  // variable in the system.
499  FEType fe_type = system.variable_type(0);
500 
501  // Build a Finite Element object of the specified type. Since the
502  // FEBase::build() member dynamically creates memory we will
503  // store the object as a std::unique_ptr<FEBase>. This can be thought
504  // of as a pointer that will clean up after itself.
505  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
506  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
507 
508  // A Gauss quadrature rule for numerical integration.
509  // Let the FEType object decide what order rule is appropriate.
510  QGauss qrule (dim, fe_type.default_quadrature_order());
511  QGauss qface (dim-1, fe_type.default_quadrature_order());
512 
513  // Tell the finite element object to use our quadrature rule.
514  fe->attach_quadrature_rule (&qrule);
515  fe_face->attach_quadrature_rule (&qface);
516 
517  // Here we define some references to cell-specific data that
518  // will be used to assemble the linear system. We will start
519  // with the element Jacobian * quadrature weight at each integration point.
520  const std::vector<Real> & JxW = fe->get_JxW();
521  const std::vector<Real> & JxW_face = fe_face->get_JxW();
522 
523  // The element shape functions evaluated at the quadrature points.
524  const std::vector<std::vector<Real>> & phi = fe->get_phi();
525  const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
526 
527  // The element shape function gradients evaluated at the quadrature
528  // points.
529  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
530 
531  // The XY locations of the quadrature points used for face integration
532  const std::vector<Point> & qface_points = fe_face->get_xyz();
533 
534  // A reference to the DofMap object for this system. The DofMap
535  // object handles the index translation from node and element numbers
536  // to degree of freedom numbers. We will talk more about the DofMap
537  // in future examples.
538  const DofMap & dof_map = system.get_dof_map();
539 
540  // Define data structures to contain the element matrix
541  // and right-hand-side vector contribution. Following
542  // basic finite element terminology we will denote these
543  // "Ke" and "Fe".
546 
547  // This vector will hold the degree of freedom indices for
548  // the element. These define where in the global system
549  // the element degrees of freedom get mapped.
550  std::vector<dof_id_type> dof_indices;
551 
552  // Here we extract the velocity & parameters that we put in the
553  // EquationSystems object.
554  const RealVectorValue velocity =
555  es.parameters.get<RealVectorValue> ("velocity");
556 
557  const Real diffusivity =
558  es.parameters.get<Real> ("diffusivity");
559 
560  const Real dt = es.parameters.get<Real> ("dt");
561 
562  // The global system matrix
563  SparseMatrix<Number> & matrix = system.get_system_matrix();
564 
565  // Now we will loop over all the elements in the mesh that
566  // live on the local processor. We will compute the element
567  // matrix and right-hand-side contribution. Since the mesh
568  // will be refined we want to only consider the ACTIVE elements,
569  // hence we use a variant of the active_elem_iterator.
570  for (const auto & elem : mesh.active_local_element_ptr_range())
571  {
572  // Get the degree of freedom indices for the
573  // current element. These define where in the global
574  // matrix and right-hand-side this element will
575  // contribute to.
576  dof_map.dof_indices (elem, dof_indices);
577 
578  // Compute the element-specific data for the current
579  // element. This involves computing the location of the
580  // quadrature points (q_point) and the shape functions
581  // (phi, dphi) for the current element.
582  fe->reinit (elem);
583 
584  const unsigned int n_dofs =
585  cast_int<unsigned int>(dof_indices.size());
586  libmesh_assert_equal_to (n_dofs, phi.size());
587 
588  // Zero the element matrix and right-hand side before
589  // summing them. We use the resize member here because
590  // the number of degrees of freedom might have changed from
591  // the last element. Note that this will be the case if the
592  // element type is different (i.e. the last element was a
593  // triangle, now we are on a quadrilateral).
594  Ke.resize (n_dofs, n_dofs);
595 
596  Fe.resize (n_dofs);
597 
598  // Now we will build the element matrix and right-hand-side.
599  // Constructing the RHS requires the solution and its
600  // gradient from the previous timestep. This myst be
601  // calculated at each quadrature point by summing the
602  // solution degree-of-freedom values by the appropriate
603  // weight functions.
604  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
605  {
606  // Values to hold the old solution & its gradient.
607  Number u_old = 0.;
608  Gradient grad_u_old;
609 
610  // Compute the old solution & its gradient.
611  for (unsigned int l=0; l != n_dofs; l++)
612  {
613  u_old += phi[l][qp]*system.old_solution (dof_indices[l]);
614 
615  // This will work,
616  // grad_u_old += dphi[l][qp]*system.old_solution (dof_indices[l]);
617  // but we can do it without creating a temporary like this:
618  grad_u_old.add_scaled (dphi[l][qp], system.old_solution (dof_indices[l]));
619  }
620 
621  // Now compute the element matrix and RHS contributions.
622  for (unsigned int i=0; i != n_dofs; i++)
623  {
624  // The RHS contribution
625  Fe(i) += JxW[qp]*(
626  // Mass matrix term
627  u_old*phi[i][qp] +
628  -.5*dt*(
629  // Convection term
630  // (grad_u_old may be complex, so the
631  // order here is important!)
632  (grad_u_old*velocity)*phi[i][qp] +
633 
634  // Diffusion term
635  diffusivity*(grad_u_old*dphi[i][qp]))
636  );
637 
638  for (unsigned int j=0; j != n_dofs; j++)
639  {
640  // The matrix contribution
641  Ke(i,j) += JxW[qp]*(
642  // Mass-matrix
643  phi[i][qp]*phi[j][qp] +
644  .5*dt*(
645  // Convection term
646  (velocity*dphi[j][qp])*phi[i][qp] +
647  // Diffusion term
648  diffusivity*(dphi[i][qp]*dphi[j][qp]))
649  );
650  }
651  }
652  }
653 
654  // At this point the interior element integration has
655  // been completed. However, we have not yet addressed
656  // boundary conditions. For this example we will only
657  // consider simple Dirichlet boundary conditions imposed
658  // via the penalty method.
659  //
660  // The following loops over the sides of the element.
661  // If the element has no neighbor on a side then that
662  // side MUST live on a boundary of the domain.
663  {
664  // The penalty value.
665  const Real penalty = 1.e10;
666 
667  // The following loops over the sides of the element.
668  // If the element has no neighbor on a side then that
669  // side MUST live on a boundary of the domain.
670  for (auto s : elem->side_index_range())
671  if (elem->neighbor_ptr(s) == nullptr)
672  {
673  fe_face->reinit(elem, s);
674 
675  libmesh_assert_equal_to (n_dofs, psi.size());
676 
677  for (unsigned int qp=0; qp<qface.n_points(); qp++)
678  {
679  const Number value = exact_solution (qface_points[qp](0),
680  qface_points[qp](1),
681  system.time);
682 
683  // RHS contribution
684  for (unsigned int i=0; i != n_dofs; i++)
685  Fe(i) += penalty*JxW_face[qp]*value*psi[i][qp];
686 
687  // Matrix contribution
688  for (unsigned int i=0; i != n_dofs; i++)
689  for (unsigned int j=0; j != n_dofs; j++)
690  Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
691  }
692  }
693  }
694 
695 
696  // We have now built the element matrix and RHS vector in terms
697  // of the element degrees of freedom. However, it is possible
698  // that some of the element DOFs are constrained to enforce
699  // solution continuity, i.e. they are not really "free". We need
700  // to constrain those DOFs in terms of non-constrained DOFs to
701  // ensure a continuous solution. The
702  // DofMap::constrain_element_matrix_and_vector() method does
703  // just that.
704  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
705 
706  // The element matrix and right-hand-side are now built
707  // for this element. Add them to the global matrix and
708  // right-hand-side vector. The SparseMatrix::add_matrix()
709  // and NumericVector::add_vector() members do this for us.
710  matrix.add_matrix (Ke, dof_indices);
711  system.rhs->add_vector (Fe, dof_indices);
712 
713  }
714  // Finished computing the system matrix and right-hand side.
715 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:651
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:2254
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:1992
unsigned int dim
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:374
MeshBase & mesh
Manages storage and variables for transient systems.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:74
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
const T & get(std::string_view) const
Definition: parameters.h:426
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
const MeshBase & get_mesh() const
static const bool value
Definition: xdr_io.C:54
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:324
Parameters parameters
Data structure holding arbitrary parameters.

◆ exact_solution()

Real exact_solution ( const Real  x,
const Real  y,
const Real  z = 0. 
)

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 assemble_cd(), and 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 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real pi
.
Definition: libmesh.h:274

◆ exact_value()

Number exact_value ( const Point p,
const Parameters parameters,
const std::string &  ,
const std::string &   
)

Definition at line 104 of file adaptivity_ex2.C.

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

Referenced by init_cd(), and libMesh::DofMap::max_constraint_error().

108 {
109  return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
110 }
const T & get(std::string_view) const
Definition: parameters.h:426
Real exact_solution(const Real x, const Real y, const Real t)
This is the exact solution that we are trying to obtain.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ init_cd() [1/2]

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

Referenced by main().

◆ init_cd() [2/2]

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

Definition at line 456 of file adaptivity_ex2.C.

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

458 {
459  // It is a good idea to make sure we are initializing
460  // the proper system.
461  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
462 
463  // Get a reference to the Convection-Diffusion system object.
465  es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
466 
467  // Project initial conditions at time 0
468  es.parameters.set<Real> ("time") = system.time = 0;
469 
470  system.project_solution(exact_value, nullptr, es.parameters);
471 }
Manages storage and variables for transient systems.
const T_sys & get_system(std::string_view name) const
Number exact_value(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Definition: parameters.h:469
Parameters parameters
Data structure holding arbitrary parameters.

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 119 of file adaptivity_ex2.C.

References assemble_cd(), libMesh::MeshRefinement::coarsen_fraction(), libMesh::command_line_next(), libMesh::default_solver_package(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::FIRST, libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::H1, libMesh::TriangleWrapper::init(), init_cd(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::invalid_uint, libMesh::MeshRefinement::max_h_level(), mesh, libMesh::on_command_line(), libMesh::out, libMesh::MeshBase::print_info(), libMesh::READ, libMesh::MeshBase::read(), libMesh::Real, libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::MeshRefinement::refine_fraction(), libMesh::TOLERANCE, libMesh::TRILINOS_SOLVERS, libMesh::MeshRefinement::uniformly_refine(), libMesh::JumpErrorEstimator::use_unweighted_quadrature_rules, libMesh::WRITE, libMesh::MeshBase::write(), and libMesh::MeshOutput< MT >::write_equation_systems().

120 {
121  // Initialize libMesh.
122  LibMeshInit init (argc, argv);
123 
124  // This example requires a linear solver package.
125  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
126  "--enable-petsc, --enable-trilinos, or --enable-eigen");
127 
128 #ifndef LIBMESH_ENABLE_AMR
129  libmesh_example_requires(false, "--enable-amr");
130 #else
131 
132  // Our Trilinos interface does not yet support adaptive transient
133  // problems
134  libmesh_example_requires(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc");
135 
136  // Brief message to the user regarding the program name
137  // and command line arguments.
138 
139  // Use commandline parameter to specify if we are to
140  // read in an initial solution or generate it ourself
141  libMesh::out << "Usage:\n"
142  <<"\t " << argv[0] << " -init_timestep 0 -n_timesteps 25 [-n_refinements 5]\n"
143  << "OR\n"
144  <<"\t " << argv[0] << " -read_solution -init_timestep 26 -n_timesteps 25\n"
145  << std::endl;
146 
147  libMesh::out << "Running: " << argv[0];
148 
149  for (int i=1; i<argc; i++)
150  libMesh::out << " " << argv[i];
151 
152  libMesh::out << std::endl << std::endl;
153 
154  // This boolean value is obtained from the command line, it is true
155  // if the flag "-read_solution" is present, false otherwise.
156  // It indicates whether we are going to read in
157  // the mesh and solution files "saved_mesh.xda" and "saved_solution.xda"
158  // or whether we are going to start from scratch by just reading
159  // "mesh.xda"
160  const bool read_solution = libMesh::on_command_line("-read_solution");
161 
162  // This value is also obtained from the commandline and it specifies the
163  // initial value for the t_step looping variable. We must
164  // distinguish between the two cases here, whether we read in the
165  // solution or we started from scratch, so that we do not overwrite the
166  // gmv output files.
167  const unsigned int init_timestep =
168  libMesh::command_line_next("-init_timestep",
170 
171  if (init_timestep == libMesh::invalid_uint)
172  {
173  // This handy function will print the file name, line number,
174  // specified message, and then throw an exception.
175  libmesh_error_msg("ERROR: Initial timestep not specified!");
176  }
177 
178  // This command line value specifies how many time steps to take.
179  const unsigned int n_timesteps =
180  libMesh::command_line_next("-n_timesteps",
182 
183  if (n_timesteps == libMesh::invalid_uint)
184  libmesh_error_msg("ERROR: Number of timesteps not specified");
185 
186  // This command line value specifies how far to allow refinement
187  const unsigned int max_h_level =
188  libMesh::command_line_next("-max_h_level", 5);
189 
190  // Skip this 2D example if libMesh was compiled as 1D-only.
191  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
192 
193  // Create a new mesh on the default MPI communicator.
194  // We still need some work on automatic parallel restarts with
195  // DistributedMesh
196  ReplicatedMesh mesh(init.comm());
197 
198  // Create an equation systems object.
199  EquationSystems equation_systems (mesh);
200  MeshRefinement mesh_refinement (mesh);
201 
202  // First we process the case where we do not read in the solution
203  if (!read_solution)
204  {
205  // Read the mesh from file.
206  mesh.read ("mesh.xda");
207 
208  // Again do a search on the command line for an argument
209  const unsigned int n_refinements = 5;
210  libMesh::command_line_next("-n_refinements", 5);
211 
212  // Uniformly refine the mesh 5 times
213  if (!read_solution)
214  mesh_refinement.uniformly_refine (n_refinements);
215 
216  // Print information about the mesh to the screen.
217  mesh.print_info();
218 
219  // Declare the system and its variables.
220  // Begin by creating a transient system
221  // named "Convection-Diffusion".
223  equation_systems.add_system<TransientLinearImplicitSystem>("Convection-Diffusion");
224 
225  // Adds the variable "u" to "Convection-Diffusion". "u"
226  // will be approximated using first-order approximation.
227  system.add_variable ("u", FIRST);
228 
229  // Give the system a pointer to the matrix assembly
230  // and initialization functions.
231  system.attach_assemble_function (assemble_cd);
232  system.attach_init_function (init_cd);
233 
234  // Initialize the data structures for the equation system.
235  equation_systems.init ();
236  }
237  // Otherwise we read in the solution and mesh
238  else
239  {
240  // Read in the mesh stored in "saved_mesh.xda"
241  mesh.read("saved_mesh.xda");
242 
243  // Print information about the mesh to the screen.
244  mesh.print_info();
245 
246  // Read in the solution stored in "saved_solution.xda"
247  equation_systems.read("saved_solution.xda", READ);
248 
249  // Get a reference to the system so that we can call update() on it
251  equation_systems.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
252 
253  // We need to call update to put system in a consistent state
254  // with the solution that was read in
255  system.update();
256 
257  // Attach the same matrix assembly function as above. Note, we do not
258  // have to attach an init() function since we are initializing the
259  // system by reading in "saved_solution.xda"
260  system.attach_assemble_function (assemble_cd);
261 
262  // Print out the H1 norm of the saved solution, for verification purposes:
263  Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
264 
265  libMesh::out << "Initial H1 norm = " << H1norm << std::endl << std::endl;
266  }
267 
268  // Prints information about the system to the screen.
269  equation_systems.print_info();
270 
271  equation_systems.parameters.set<unsigned int>("linear solver maximum iterations") = 250;
272  equation_systems.parameters.set<Real>("linear solver tolerance") = TOLERANCE;
273 
274  if (!read_solution)
275  // Write out the initial condition
276  GMVIO(mesh).write_equation_systems ("out.gmv.000", equation_systems);
277  else
278  // Write out the solution that was read in
279  GMVIO(mesh).write_equation_systems ("solution_read_in.gmv", equation_systems);
280 
281  // The Convection-Diffusion system requires that we specify
282  // the flow velocity. We will specify it as a RealVectorValue
283  // data type and then use the Parameters object to pass it to
284  // the assemble function.
285  equation_systems.parameters.set<RealVectorValue>("velocity") =
286  RealVectorValue (0.8, 0.8);
287 
288  // The Convection-Diffusion system also requires a specified
289  // diffusivity. We use an isotropic (hence Real) value.
290  equation_systems.parameters.set<Real>("diffusivity") = 0.01;
291 
292  // Solve the system "Convection-Diffusion". This will be done by
293  // looping over the specified time interval and calling the
294  // solve() member at each time step. This will assemble the
295  // system and call the linear solver.
296 
297  // Since only TransientLinearImplicitSystems (and systems
298  // derived from them) contain old solutions, to use the
299  // old_local_solution later we now need to specify the system
300  // type when we ask for it.
302  equation_systems.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
303 
304  const Real dt = 0.025;
305  system.time = init_timestep*dt;
306 
307  // We do 25 timesteps both before and after writing out the
308  // intermediate solution
309  for (unsigned int t_step=init_timestep; t_step<(init_timestep+n_timesteps); t_step++)
310  {
311  // Increment the time counter, set the time and the
312  // time step size as parameters in the EquationSystem.
313  system.time += dt;
314 
315  equation_systems.parameters.set<Real> ("time") = system.time;
316  equation_systems.parameters.set<Real> ("dt") = dt;
317 
318  // A pretty update message
319  libMesh::out << " Solving time step ";
320 
321  // Add a set of scope braces to enforce data locality.
322  {
323  std::ostringstream out;
324 
325  out << std::setw(2)
326  << std::right
327  << t_step
328  << ", time="
329  << std::fixed
330  << std::setw(6)
331  << std::setprecision(3)
332  << std::setfill('0')
333  << std::left
334  << system.time
335  << "...";
336 
337  libMesh::out << out.str() << std::endl;
338  }
339 
340  // At this point we need to update the old
341  // solution vector. The old solution vector
342  // will be the current solution vector from the
343  // previous time step.
344 
345  *system.old_local_solution = *system.current_local_solution;
346 
347  // The number of refinement steps per time step.
348  const unsigned int max_r_steps = 2;
349 
350  // A refinement loop.
351  for (unsigned int r_step=0; r_step<max_r_steps; r_step++)
352  {
353  // Assemble & solve the linear system
354  system.solve();
355 
356  // Print out the H1 norm, for verification purposes:
357  Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
358 
359  libMesh::out << "H1 norm = " << H1norm << std::endl;
360 
361  // Possibly refine the mesh
362  if (r_step+1 != max_r_steps)
363  {
364  libMesh::out << " Refining the mesh..." << std::endl;
365 
366  // The ErrorVector is a particular StatisticsVector
367  // for computing error information on a finite element mesh.
368  ErrorVector error;
369 
370  // The ErrorEstimator class interrogates a finite element
371  // solution and assigns to each element a positive error value.
372  // This value is used for deciding which elements to refine
373  // and which to coarsen.
374  KellyErrorEstimator error_estimator;
375 
376  // This is a subclass of JumpErrorEstimator, based on
377  // measuring discontinuities across sides between
378  // elements, and we can tell it to use a cheaper
379  // "unweighted" quadrature rule when numerically
380  // integrating those discontinuities.
381  error_estimator.use_unweighted_quadrature_rules = true;
382 
383  // Compute the error for each active element using the provided
384  // flux_jump indicator. Note in general you will need to
385  // provide an error estimator specifically designed for your
386  // application.
387  error_estimator.estimate_error (system,
388  error);
389 
390  // This takes the error in error and decides which elements
391  // will be coarsened or refined. Any element within 20% of the
392  // maximum error on any element will be refined, and any
393  // element within 7% of the minimum error on any element might
394  // be coarsened. Note that the elements flagged for refinement
395  // will be refined, but those flagged for coarsening _might_ be
396  // coarsened.
397  mesh_refinement.refine_fraction() = 0.80;
398  mesh_refinement.coarsen_fraction() = 0.07;
399  mesh_refinement.max_h_level() = max_h_level;
400  mesh_refinement.flag_elements_by_error_fraction (error);
401 
402  // This call actually refines and coarsens the flagged
403  // elements.
404  mesh_refinement.refine_and_coarsen_elements();
405 
406  // This call reinitializes the EquationSystems object for
407  // the newly refined mesh. One of the steps in the
408  // reinitialization is projecting the solution,
409  // old_solution, etc... vectors from the old mesh to
410  // the current one.
411  equation_systems.reinit ();
412  }
413  }
414 
415  // Again do a search on the command line for an argument
416  const unsigned int output_freq =
417  libMesh::command_line_next("-output_freq", 10);
418 
419  // Output every 10 timesteps to file.
420  if ((t_step+1)%output_freq == 0)
421  {
422  std::ostringstream file_name;
423 
424  file_name << "out.gmv."
425  << std::setw(3)
426  << std::setfill('0')
427  << std::right
428  << t_step+1;
429 
430  GMVIO(mesh).write_equation_systems (file_name.str(),
431  equation_systems);
432  }
433  }
434 
435  if (!read_solution)
436  {
437  // Print out the H1 norm of the saved solution, for verification purposes:
438  Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
439 
440  libMesh::out << "Final H1 norm = " << H1norm << std::endl << std::endl;
441 
442  mesh.write("saved_mesh.xda");
443  equation_systems.write("saved_solution.xda", WRITE);
444  GMVIO(mesh).write_equation_systems ("saved_solution.gmv",
445  equation_systems);
446  }
447 #endif // #ifndef LIBMESH_ENABLE_AMR
448 
449  return 0;
450 }
T command_line_next(std::string name, T default_value)
Use GetPot&#39;s search()/next() functions to get following arguments from the command line...
Definition: libmesh.C:1011
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.
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:286
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 ...
Definition: mesh_output.C:31
static constexpr Real TOLERANCE
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
MeshBase & mesh
Manages storage and variables for transient systems.
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
Definition: system_norm.h:45
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:46
SolverPackage default_solver_package()
Definition: libmesh.C:1050
Implements (adaptive) mesh refinement algorithms for a MeshBase.
void assemble_cd(EquationSystems &es, const std::string &system_name)
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1489
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
virtual void write(const std::string &name)=0
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
void init_cd(EquationSystems &es, const std::string &system_name)
OStreamProxy out
bool use_unweighted_quadrature_rules
This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweigh...
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&#39;s jump error estimate formula to estimate the error on each cell...
bool on_command_line(std::string arg)
Definition: libmesh.C:924