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

Function Documentation

◆ assemble_cd()

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

Definition at line 296 of file transient_ex1.C.

References libMesh::TypeVector< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), exact_solution(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), 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().

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

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

References libMesh::pi, libMesh::Utility::pow(), and libMesh::Real.

Referenced by assemble_cd(), and exact_value().

46 {
47  static const Real pi = acos(-1.);
48 
49  return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
50 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real pi
.
Definition: libmesh.h:274

◆ exact_value()

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

Definition at line 92 of file transient_ex1.C.

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

Referenced by init_cd().

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

◆ init_cd() [1/2]

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

Referenced by main().

◆ init_cd() [2/2]

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

Definition at line 274 of file transient_ex1.C.

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

276 {
277  // It is a good idea to make sure we are initializing
278  // the proper system.
279  libmesh_assert_equal_to (system_name, "Convection-Diffusion");
280 
281  // Get a reference to the Convection-Diffusion system object.
283  es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
284 
285  // Project initial conditions at time 0
286  es.parameters.set<Real> ("time") = system.time = 0;
287 
288  system.project_solution(exact_value, nullptr, es.parameters);
289 }
Manages storage and variables for transient systems.
Number exact_value(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
Definition: transient_ex1.C:92
const T_sys & get_system(std::string_view name) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Definition: parameters.h:469
Parameters parameters
Data structure holding arbitrary parameters.

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 105 of file transient_ex1.C.

References libMesh::EquationSystems::add_system(), libMesh::ExodusII_IO::append(), assemble_cd(), libMesh::default_solver_package(), exodus_filename(), libMesh::FIRST, libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), init_cd(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::TransientSystem< Base >::old_local_solution, libMesh::out, libMesh::EquationSystems::parameters, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Real, libMesh::Parameters::set(), libMesh::MeshRefinement::uniformly_refine(), libMesh::MeshOutput< MT >::write_equation_systems(), and libMesh::ExodusII_IO::write_timestep().

106 {
107  // Initialize libMesh.
108  LibMeshInit init (argc, argv);
109 
110  // This example requires a linear solver package.
111  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
112  "--enable-petsc, --enable-trilinos, or --enable-eigen");
113 
114  // This example requires Adaptive Mesh Refinement support - although
115  // it only refines uniformly, the refinement code used is the same
116  // underneath
117 #ifndef LIBMESH_ENABLE_AMR
118  libmesh_example_requires(false, "--enable-amr");
119 #else
120 
121  // Skip this 2D example if libMesh was compiled as 1D-only.
122  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
123 
124  // Read the mesh from file. This is the coarse mesh that will be used
125  // in example 10 to demonstrate adaptive mesh refinement. Here we will
126  // simply read it in and uniformly refine it before we compute with
127  // it.
128  //
129  // Create a mesh object, with dimension to be overridden later,
130  // distributed across the default MPI communicator.
131  Mesh mesh(init.comm());
132 
133  mesh.read ("mesh.xda");
134 
135  // Query the command line for the number of mesh refinements to use
136  GetPot input(argc, argv);
137  const unsigned int n_refinements = input("n_refinements", 5);
138 
139  // Create a MeshRefinement object to handle refinement of our mesh.
140  // This class handles all the details of mesh refinement and coarsening.
141  MeshRefinement mesh_refinement (mesh);
142 
143  // Uniformly refine the mesh as requested.
144  mesh_refinement.uniformly_refine (n_refinements);
145 
146  // Print information about the mesh to the screen.
147  mesh.print_info();
148 
149  // Create an equation systems object.
150  EquationSystems equation_systems (mesh);
151 
152  // Add a transient system to the EquationSystems
153  // object named "Convection-Diffusion".
155  equation_systems.add_system<TransientLinearImplicitSystem> ("Convection-Diffusion");
156 
157  // Adds the variable "u" to "Convection-Diffusion". "u"
158  // will be approximated using first-order approximation.
159  system.add_variable ("u", FIRST);
160 
161  // Give the system a pointer to the matrix assembly
162  // and initialization functions.
163  system.attach_assemble_function (assemble_cd);
164  system.attach_init_function (init_cd);
165 
166  // Initialize the data structures for the equation system.
167  equation_systems.init ();
168 
169  // Prints information about the system to the screen.
170  equation_systems.print_info();
171 
172  // Write out the initial conditions.
173 #ifdef LIBMESH_HAVE_EXODUS_API
174  // If Exodus is available, we'll write all timesteps to the same file
175  // rather than one file per timestep.
176  std::string exodus_filename = "transient_ex1.e";
178 #else
179  GMVIO(mesh).write_equation_systems ("out_000.gmv", equation_systems);
180 #endif
181 
182  // The Convection-Diffusion system requires that we specify
183  // the flow velocity. We will specify it as a RealVectorValue
184  // data type and then use the Parameters object to pass it to
185  // the assemble function.
186  equation_systems.parameters.set<RealVectorValue>("velocity") =
187  RealVectorValue (0.8, 0.8);
188 
189  // Solve the system "Convection-Diffusion". This will be done by
190  // looping over the specified time interval and calling the
191  // solve() member at each time step. This will assemble the
192  // system and call the linear solver.
193  const Real dt = 0.025;
194  system.time = 0.;
195 
196  for (unsigned int t_step = 0; t_step < 50; t_step++)
197  {
198  // Increment the time counter, set the time and the
199  // time step size as parameters in the EquationSystem.
200  system.time += dt;
201 
202  equation_systems.parameters.set<Real> ("time") = system.time;
203  equation_systems.parameters.set<Real> ("dt") = dt;
204 
205  // A pretty update message
206  libMesh::out << " Solving time step ";
207 
208  // Do fancy zero-padded formatting of the current time.
209  {
210  std::ostringstream out;
211 
212  out << std::setw(2)
213  << std::right
214  << t_step
215  << ", time="
216  << std::fixed
217  << std::setw(6)
218  << std::setprecision(3)
219  << std::setfill('0')
220  << std::left
221  << system.time
222  << "...";
223 
224  libMesh::out << out.str() << std::endl;
225  }
226 
227  // At this point we need to update the old
228  // solution vector. The old solution vector
229  // will be the current solution vector from the
230  // previous time step. We will do this by extracting the
231  // system from the EquationSystems object and using
232  // vector assignment. Since only TransientSystems
233  // (and systems derived from them) contain old solutions
234  // we need to specify the system type when we ask for it.
235  *system.old_local_solution = *system.current_local_solution;
236 
237  // Assemble & solve the linear system
238  equation_systems.get_system("Convection-Diffusion").solve();
239 
240  // Output every 10 timesteps to file.
241  if ((t_step+1)%10 == 0)
242  {
243 
244 #ifdef LIBMESH_HAVE_EXODUS_API
245  ExodusII_IO exo(mesh);
246  exo.append(true);
247  exo.write_timestep (exodus_filename, equation_systems, t_step+1, system.time);
248 #else
249  std::ostringstream file_name;
250 
251  file_name << "out_"
252  << std::setw(3)
253  << std::setfill('0')
254  << std::right
255  << t_step+1
256  << ".gmv";
257 
258 
259  GMVIO(mesh).write_equation_systems (file_name.str(),
260  equation_systems);
261 #endif
262  }
263  }
264 #endif // #ifdef LIBMESH_ENABLE_AMR
265 
266  // All done.
267  return 0;
268 }
This is the EquationSystems class.
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
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
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
MeshBase & mesh
Manages storage and variables for transient systems.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:46
void assemble_cd(EquationSystems &es, const std::string &system_name)
NumericVector< Number > * old_local_solution
All the values I need to compute my contribution to the simulation at hand.
void init_cd(EquationSystems &es, const std::string &system_name)
SolverPackage default_solver_package()
Definition: libmesh.C:1050
Implements (adaptive) mesh refinement algorithms for a MeshBase.
std::string exodus_filename(unsigned number)
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1489
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50