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 &  libmesh_dbg_varsystem_name 
)

Definition at line 487 of file adaptivity_ex2.C.

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

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

◆ assemble_cd() [2/2]

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

Definition at line 293 of file transient_ex1.C.

295 {
296  // Ignore unused parameter warnings when !LIBMESH_ENABLE_AMR.
297  libmesh_ignore(es, 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 std::unique_ptr<FEBase>. This can be thought
321  // of as a pointer that will clean up after itself.
322  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
323  std::unique_ptr<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) == 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 }

Referenced by main().

◆ 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.

46 {
47  static const Real pi = acos(-1.);
48 
49  return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
50 }

Referenced by assemble_cd(), and exact_value().

◆ 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.

108 {
109  return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
110 }

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

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

◆ init_cd() [1/2]

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

Definition at line 463 of file adaptivity_ex2.C.

465 {
466  // It is a good idea to make sure we are initializing
467  // the proper system.
468  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
469 
470  // Get a reference to the Convection-Diffusion system object.
472  es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
473 
474  // Project initial conditions at time 0
475  es.parameters.set<Real> ("time") = system.time = 0;
476 
477  system.project_solution(exact_value, nullptr, es.parameters);
478 }

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

◆ init_cd() [2/2]

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

Referenced by main().

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 119 of file adaptivity_ex2.C.

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

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

exact_solution
Real exact_solution(const Real x, const Real y, const Real t)
This is the exact solution that we are trying to obtain.
Definition: exact_solution.C:43
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::pi
const Real pi
.
Definition: libmesh.h:237
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::MeshBase::read
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.
libMesh::TransientSystem
Manages storage and variables for transient systems.
Definition: transient_system.h:57
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::DofMap::dof_indices
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:1967
libMesh::JumpErrorEstimator::estimate_error
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
This function uses the derived class's jump error estimate formula to estimate the error on each cell...
Definition: jump_error_estimator.C:53
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::MeshBase::write
virtual void write(const std::string &name)=0
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::DenseMatrix< Number >
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
libMesh::WRITE
Definition: enum_xdr_mode.h:40
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::GMVIO
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:54
init_cd
void init_cd(EquationSystems &es, const std::string &system_name)
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::VectorValue< Real >
libMesh::ReplicatedMesh
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
Definition: replicated_mesh.h:47
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
libMesh::KellyErrorEstimator
This class implements the Kelly error indicator which is based on the flux jumps between elements.
Definition: kelly_error_estimator.h:59
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::TypeVector::add_scaled
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:665
libMesh::ErrorVector
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::READ
Definition: enum_xdr_mode.h:41
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::Parameters::set
T & set(const std::string &)
Definition: parameters.h:460
libMesh::DofMap::constrain_element_matrix_and_vector
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:2034
libMesh::MeshOutput::write_equation_systems
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
libMesh::TRILINOS_SOLVERS
Definition: enum_solver_package.h:37
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
value
static const bool value
Definition: xdr_io.C:56
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::TransientSystem::old_solution
Number old_solution(const dof_id_type global_dof_number) const
Definition: transient_system.C:122
libMesh::JumpErrorEstimator::use_unweighted_quadrature_rules
bool use_unweighted_quadrature_rules
This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweigh...
Definition: jump_error_estimator.h:113
libMesh::SystemNorm
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
Definition: system_norm.h:51
assemble_cd
void assemble_cd(EquationSystems &es, const std::string &system_name)
Definition: transient_ex1.C:293
libMesh::H1
Definition: enum_norm_type.h:37
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
exact_solution
Real exact_solution(const Real x, const Real y, const Real t)
This is the exact solution that we are trying to obtain.
Definition: exact_solution.C:43
libMesh::Parameters::get
const T & get(const std::string &) const
Definition: parameters.h:421
libMesh::out
OStreamProxy out
libMesh::FIRST
Definition: enum_order.h:42
libMesh::DenseVector< Number >
exact_value
Number exact_value(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
Definition: adaptivity_ex2.C:104
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557