libMesh
Functions | Variables
adaptivity_ex5.C File Reference

Go to the source code of this file.

Functions

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

Variables

std::unique_ptr< FunctionBase< Number > > parsed_solution
 

Function Documentation

◆ assemble_cd() [1/2]

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

Definition at line 575 of file adaptivity_ex5.C.

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

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::TypeVector< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::TransientSystem< Base >::old_solution(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseVector< T >::resize(), and libMesh::DenseMatrix< T >::resize().

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

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::TypeVector< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), exact_solution(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::libmesh_ignore(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::TransientSystem< Base >::old_solution(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and value.

Referenced by main().

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

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

◆ exact_value()

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

Definition at line 106 of file adaptivity_ex5.C.

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

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

Referenced by init_cd().

◆ exodus_filename()

std::string exodus_filename ( unsigned  number)

Definition at line 766 of file adaptivity_ex5.C.

767 {
768  std::ostringstream oss;
769 
770  oss << "out_";
771  oss << std::setw(3) << std::setfill('0') << number;
772  oss << ".e";
773 
774  return oss.str();
775 }

Referenced by main().

◆ init_cd() [1/2]

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

Definition at line 548 of file adaptivity_ex5.C.

550 {
551  // It is a good idea to make sure we are initializing
552  // the proper system.
553  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
554 
555  // Get a reference to the Convection-Diffusion system object.
557  es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
558 
559  // Project initial conditions at time 0
560  es.parameters.set<Real> ("time") = system.time = 0;
561 
562  if (parsed_solution.get())
563  system.project_solution(parsed_solution.get(), nullptr);
564  else
565  system.project_solution(exact_value, nullptr, es.parameters);
566 }

References exact_value(), libMesh::EquationSystems::get_system(), libMesh::if(), libMesh::EquationSystems::parameters, parsed_solution, 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 129 of file adaptivity_ex5.C.

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

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

Variable Documentation

◆ parsed_solution

std::unique_ptr<FunctionBase<Number> > parsed_solution

Definition at line 116 of file adaptivity_ex5.C.

Referenced by init_cd(), and main().

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::DofMap::get_periodic_boundaries
PeriodicBoundaries * get_periodic_boundaries()
Definition: dof_map.h:1295
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
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::if
if(subdm)
Definition: petsc_dm_wrapper.C:77
libMesh::GMVIO
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:54
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::RealVectorValue
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
Definition: hp_coarsentest.h:46
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::BasicOStreamProxy::flags
std::ios_base::fmtflags flags() const
Get the associated format flags.
Definition: ostream_proxy.h:158
libMesh::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
assemble_cd
void assemble_cd(EquationSystems &es, const std::string &system_name)
Definition: transient_ex1.C:293
libMesh::MeshTools::Generation::build_square
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
Definition: mesh_generation.C:1501
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
exodus_filename
std::string exodus_filename(unsigned number)
Definition: adaptivity_ex5.C:766
exact_value
Number exact_value(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
Definition: adaptivity_ex5.C:106
libMesh::DECODE
Definition: enum_xdr_mode.h:39
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::QUAD4
Definition: enum_elem_type.h:41
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
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::DofMap::add_periodic_boundary
void add_periodic_boundary(const PeriodicBoundaryBase &periodic_boundary)
Adds a copy of the specified periodic boundary to the system.
Definition: dof_map_constraints.C:4510
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::PeriodicBoundary
The definition of a periodic boundary.
Definition: periodic_boundary.h:44
libMesh::TRILINOS_SOLVERS
Definition: enum_solver_package.h:37
libMesh::ENCODE
Definition: enum_xdr_mode.h:38
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
init_cd
void init_cd(EquationSystems &es, const std::string &system_name)
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
libMesh::TransientSystem::old_local_solution
std::unique_ptr< NumericVector< Number > > old_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: transient_system.h:119
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
parsed_solution
std::unique_ptr< FunctionBase< Number > > parsed_solution
Definition: adaptivity_ex5.C:116
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 >
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557