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

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

Definition at line 292 of file transient_ex1.C.

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 479 of file adaptivity_ex2.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_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.

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

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

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

106 {
107  return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
108 }
const T & get(const std::string &) const
Definition: parameters.h:430
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
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 455 of file adaptivity_ex2.C.

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

457 {
458  // It is a good idea to make sure we are initializing
459  // the proper system.
460  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
461 
462  // Get a reference to the Convection-Diffusion system object.
464  es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
465 
466  // Project initial conditions at time 0
467  es.parameters.set<Real> ("time") = system.time = 0;
468 
469  system.project_solution(exact_value, libmesh_nullptr, es.parameters);
470 }
const class libmesh_nullptr_t libmesh_nullptr
This class provides a specific system class.
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.
const T_sys & get_system(const std::string &name) const
int main ( int  argc,
char **  argv 
)

Definition at line 117 of file adaptivity_ex2.C.

References assemble_cd(), libMesh::MeshRefinement::coarsen_fraction(), libMesh::LibMeshInit::comm(), 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::MeshRefinement::max_h_level(), mesh, 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::WRITE, libMesh::MeshBase::write(), and libMesh::MeshOutput< MT >::write_equation_systems().

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