libMesh
Functions
systems_of_equations_ex1.C File Reference

Go to the source code of this file.

Functions

void assemble_stokes (EquationSystems &es, const std::string &system_name)
 
int main (int argc, char **argv)
 
void assemble_stokes (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble_stokes() [1/2]

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

Definition at line 148 of file systems_of_equations_ex1.C.

150 {
151  // It is a good idea to make sure we are assembling
152  // the proper system.
153  libmesh_assert_equal_to (system_name, "Stokes");
154 
155  // Get a constant reference to the mesh object.
156  const MeshBase & mesh = es.get_mesh();
157 
158  // The dimension that we are running
159  const unsigned int dim = mesh.mesh_dimension();
160 
161  // Get a reference to the Convection-Diffusion system object.
162  LinearImplicitSystem & system =
163  es.get_system<LinearImplicitSystem> ("Stokes");
164 
165  // Numeric ids corresponding to each variable in the system
166  const unsigned int u_var = system.variable_number ("u");
167  const unsigned int v_var = system.variable_number ("v");
168  const unsigned int p_var = system.variable_number ("p");
169 
170  // Get the Finite Element type for "u". Note this will be
171  // the same as the type for "v".
172  FEType fe_vel_type = system.variable_type(u_var);
173 
174  // Get the Finite Element type for "p".
175  FEType fe_pres_type = system.variable_type(p_var);
176 
177  // Build a Finite Element object of the specified type for
178  // the velocity variables.
179  std::unique_ptr<FEBase> fe_vel (FEBase::build(dim, fe_vel_type));
180 
181  // Build a Finite Element object of the specified type for
182  // the pressure variables.
183  std::unique_ptr<FEBase> fe_pres (FEBase::build(dim, fe_pres_type));
184 
185  // A Gauss quadrature rule for numerical integration.
186  // Let the FEType object decide what order rule is appropriate.
187  QGauss qrule (dim, fe_vel_type.default_quadrature_order());
188 
189  // Tell the finite element objects to use our quadrature rule.
190  fe_vel->attach_quadrature_rule (&qrule);
191  fe_pres->attach_quadrature_rule (&qrule);
192 
193  // Here we define some references to cell-specific data that
194  // will be used to assemble the linear system.
195  //
196  // The element Jacobian * quadrature weight at each integration point.
197  const std::vector<Real> & JxW = fe_vel->get_JxW();
198 
199  // The element shape function gradients for the velocity
200  // variables evaluated at the quadrature points.
201  const std::vector<std::vector<RealGradient>> & dphi = fe_vel->get_dphi();
202 
203  // The element shape functions for the pressure variable
204  // evaluated at the quadrature points.
205  const std::vector<std::vector<Real>> & psi = fe_pres->get_phi();
206 
207  // A reference to the DofMap object for this system. The DofMap
208  // object handles the index translation from node and element numbers
209  // to degree of freedom numbers. We will talk more about the DofMap
210  // in future examples.
211  const DofMap & dof_map = system.get_dof_map();
212 
213  // Define data structures to contain the element matrix
214  // and right-hand-side vector contribution. Following
215  // basic finite element terminology we will denote these
216  // "Ke" and "Fe".
219 
221  Kuu(Ke), Kuv(Ke), Kup(Ke),
222  Kvu(Ke), Kvv(Ke), Kvp(Ke),
223  Kpu(Ke), Kpv(Ke), Kpp(Ke);
224 
226  Fu(Fe),
227  Fv(Fe),
228  Fp(Fe);
229 
230  // This vector will hold the degree of freedom indices for
231  // the element. These define where in the global system
232  // the element degrees of freedom get mapped.
233  std::vector<dof_id_type> dof_indices;
234  std::vector<dof_id_type> dof_indices_u;
235  std::vector<dof_id_type> dof_indices_v;
236  std::vector<dof_id_type> dof_indices_p;
237 
238  // Now we will loop over all the elements in the mesh that
239  // live on the local processor. We will compute the element
240  // matrix and right-hand-side contribution. In case users later
241  // modify this program to include refinement, we will be safe and
242  // will only consider the active elements; hence we use a variant of
243  // the active_elem_iterator.
244  for (const auto & elem : mesh.active_local_element_ptr_range())
245  {
246  // Get the degree of freedom indices for the
247  // current element. These define where in the global
248  // matrix and right-hand-side this element will
249  // contribute to.
250  dof_map.dof_indices (elem, dof_indices);
251  dof_map.dof_indices (elem, dof_indices_u, u_var);
252  dof_map.dof_indices (elem, dof_indices_v, v_var);
253  dof_map.dof_indices (elem, dof_indices_p, p_var);
254 
255  const unsigned int n_dofs = dof_indices.size();
256  const unsigned int n_u_dofs = dof_indices_u.size();
257  const unsigned int n_v_dofs = dof_indices_v.size();
258  const unsigned int n_p_dofs = dof_indices_p.size();
259 
260  // Compute the element-specific data for the current
261  // element. This involves computing the location of the
262  // quadrature points (q_point) and the shape functions
263  // (phi, dphi) for the current element.
264  fe_vel->reinit (elem);
265  fe_pres->reinit (elem);
266 
267  // Zero the element matrix and right-hand side before
268  // summing them. We use the resize member here because
269  // the number of degrees of freedom might have changed from
270  // the last element. Note that this will be the case if the
271  // element type is different (i.e. the last element was a
272  // triangle, now we are on a quadrilateral).
273  Ke.resize (n_dofs, n_dofs);
274  Fe.resize (n_dofs);
275 
276  // Reposition the submatrices... The idea is this:
277  //
278  // - - - -
279  // | Kuu Kuv Kup | | Fu |
280  // Ke = | Kvu Kvv Kvp |; Fe = | Fv |
281  // | Kpu Kpv Kpp | | Fp |
282  // - - - -
283  //
284  // The DenseSubMatrix.reposition () member takes the
285  // (row_offset, column_offset, row_size, column_size).
286  //
287  // Similarly, the DenseSubVector.reposition () member
288  // takes the (row_offset, row_size)
289  Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
290  Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
291  Kup.reposition (u_var*n_u_dofs, p_var*n_u_dofs, n_u_dofs, n_p_dofs);
292 
293  Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
294  Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
295  Kvp.reposition (v_var*n_v_dofs, p_var*n_v_dofs, n_v_dofs, n_p_dofs);
296 
297  Kpu.reposition (p_var*n_u_dofs, u_var*n_u_dofs, n_p_dofs, n_u_dofs);
298  Kpv.reposition (p_var*n_u_dofs, v_var*n_u_dofs, n_p_dofs, n_v_dofs);
299  Kpp.reposition (p_var*n_u_dofs, p_var*n_u_dofs, n_p_dofs, n_p_dofs);
300 
301  Fu.reposition (u_var*n_u_dofs, n_u_dofs);
302  Fv.reposition (v_var*n_u_dofs, n_v_dofs);
303  Fp.reposition (p_var*n_u_dofs, n_p_dofs);
304 
305  // Now we will build the element matrix.
306  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
307  {
308  // Assemble the u-velocity row
309  // uu coupling
310  for (unsigned int i=0; i<n_u_dofs; i++)
311  for (unsigned int j=0; j<n_u_dofs; j++)
312  Kuu(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
313 
314  // up coupling
315  for (unsigned int i=0; i<n_u_dofs; i++)
316  for (unsigned int j=0; j<n_p_dofs; j++)
317  Kup(i,j) += -JxW[qp]*psi[j][qp]*dphi[i][qp](0);
318 
319 
320  // Assemble the v-velocity row
321  // vv coupling
322  for (unsigned int i=0; i<n_v_dofs; i++)
323  for (unsigned int j=0; j<n_v_dofs; j++)
324  Kvv(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
325 
326  // vp coupling
327  for (unsigned int i=0; i<n_v_dofs; i++)
328  for (unsigned int j=0; j<n_p_dofs; j++)
329  Kvp(i,j) += -JxW[qp]*psi[j][qp]*dphi[i][qp](1);
330 
331 
332  // Assemble the pressure row
333  // pu coupling
334  for (unsigned int i=0; i<n_p_dofs; i++)
335  for (unsigned int j=0; j<n_u_dofs; j++)
336  Kpu(i,j) += -JxW[qp]*psi[i][qp]*dphi[j][qp](0);
337 
338  // pv coupling
339  for (unsigned int i=0; i<n_p_dofs; i++)
340  for (unsigned int j=0; j<n_v_dofs; j++)
341  Kpv(i,j) += -JxW[qp]*psi[i][qp]*dphi[j][qp](1);
342 
343  } // end of the quadrature point qp-loop
344 
345  // At this point the interior element integration has
346  // been completed. However, we have not yet addressed
347  // boundary conditions. For this example we will only
348  // consider simple Dirichlet boundary conditions imposed
349  // via the penalty method. The penalty method used here
350  // is equivalent (for Lagrange basis functions) to lumping
351  // the matrix resulting from the L2 projection penalty
352  // approach introduced in example 3.
353  {
354  // The following loops over the sides of the element.
355  // If the element has no neighbor on a side then that
356  // side MUST live on a boundary of the domain.
357  for (auto s : elem->side_index_range())
358  if (elem->neighbor_ptr(s) == nullptr)
359  {
360  std::unique_ptr<const Elem> side (elem->build_side_ptr(s));
361 
362  // Loop over the nodes on the side.
363  for (auto ns : side->node_index_range())
364  {
365  // The location on the boundary of the current
366  // node.
367 
368  // const Real xf = side->point(ns)(0);
369  const Real yf = side->point(ns)(1);
370 
371  // The penalty value. \f$ \frac{1}{\epsilon \f$
372  const Real penalty = 1.e10;
373 
374  // The boundary values.
375 
376  // Set u = 1 on the top boundary, 0 everywhere else
377  const Real u_value = (yf > .99) ? 1. : 0.;
378 
379  // Set v = 0 everywhere
380  const Real v_value = 0.;
381 
382  // Find the node on the element matching this node on
383  // the side. That defined where in the element matrix
384  // the boundary condition will be applied.
385  for (auto n : elem->node_index_range())
386  if (elem->node_id(n) == side->node_id(ns))
387  {
388  // Matrix contribution.
389  Kuu(n,n) += penalty;
390  Kvv(n,n) += penalty;
391 
392  // Right-hand-side contribution.
393  Fu(n) += penalty*u_value;
394  Fv(n) += penalty*v_value;
395  }
396  } // end face node loop
397  } // end if (elem->neighbor(side) == nullptr)
398  } // end boundary condition section
399 
400  // If this assembly program were to be used on an adaptive mesh,
401  // we would have to apply any hanging node constraint equations.
402  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
403 
404  // The element matrix and right-hand-side are now built
405  // for this element. Add them to the global matrix and
406  // right-hand-side vector. The NumericMatrix::add_matrix()
407  // and NumericVector::add_vector() members do this for us.
408  system.matrix->add_matrix (Ke, dof_indices);
409  system.rhs->add_vector (Fe, dof_indices);
410  } // end of element loop
411 }

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::FEType::default_quadrature_order(), dim, libMesh::DofMap::dof_indices(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::Real, libMesh::DenseSubVector< T >::reposition(), libMesh::DenseSubMatrix< T >::reposition(), libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, libMesh::System::variable_number(), and libMesh::System::variable_type().

◆ assemble_stokes() [2/2]

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

Referenced by main().

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 68 of file systems_of_equations_ex1.C.

69 {
70  // Initialize libMesh.
71  LibMeshInit init (argc, argv);
72 
73  // This example requires a linear solver package.
74  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
75  "--enable-petsc, --enable-trilinos, or --enable-eigen");
76 
77  // Skip this 2D example if libMesh was compiled as 1D-only.
78  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
79 
80  // This example NaNs with the Eigen sparse linear solvers and
81  // Trilinos solvers, but should work OK with either PETSc or
82  // Laspack.
83  libmesh_example_requires(libMesh::default_solver_package() != EIGEN_SOLVERS, "--enable-petsc or --enable-laspack");
84  libmesh_example_requires(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc or --enable-laspack");
85 
86  // Create a mesh, with dimension to be overridden later, distributed
87  // across the default MPI communicator.
88  Mesh mesh(init.comm());
89 
90  // Use the MeshTools::Generation mesh generator to create a uniform
91  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
92  // to build a mesh of 8x8 Quad9 elements. Building these
93  // higher-order elements allows us to use higher-order
94  // approximation, as in example 3.
96  15, 15,
97  0., 1.,
98  0., 1.,
99  QUAD9);
100 
101  // Print information about the mesh to the screen.
102  mesh.print_info();
103 
104  // Create an equation systems object.
105  EquationSystems equation_systems (mesh);
106 
107  // Declare the system and its variables.
108  // Create a transient system named "Stokes"
109  LinearImplicitSystem & system =
110  equation_systems.add_system<LinearImplicitSystem> ("Stokes");
111 
112  // Add the variables "u" & "v" to "Stokes". They
113  // will be approximated using second-order approximation.
114  system.add_variable ("u", SECOND);
115  system.add_variable ("v", SECOND);
116 
117  // Add the variable "p" to "Stokes". This will
118  // be approximated with a first-order basis,
119  // providing an LBB-stable pressure-velocity pair.
120  system.add_variable ("p", FIRST);
121 
122  // Give the system a pointer to the matrix assembly
123  // function.
125 
126  // Initialize the data structures for the equation system.
127  equation_systems.init ();
128 
129  equation_systems.parameters.set<unsigned int>("linear solver maximum iterations") = 250;
130  equation_systems.parameters.set<Real> ("linear solver tolerance") = TOLERANCE;
131 
132  // Prints information about the system to the screen.
133  equation_systems.print_info();
134 
135  // Assemble & solve the linear system,
136  // then write the solution.
137  equation_systems.get_system("Stokes").solve();
138 
139 #ifdef LIBMESH_HAVE_EXODUS_API
141  equation_systems);
142 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
143 
144  // All done.
145  return 0;
146 }

References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_stokes(), libMesh::System::attach_assemble_function(), libMesh::MeshTools::Generation::build_square(), libMesh::default_solver_package(), libMesh::EIGEN_SOLVERS, libMesh::FIRST, libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::EquationSystems::parameters, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::QUAD9, libMesh::Real, libMesh::SECOND, libMesh::Parameters::set(), libMesh::TOLERANCE, libMesh::TRILINOS_SOLVERS, and libMesh::MeshOutput< MT >::write_equation_systems().

libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
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::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::DenseSubMatrix
Defines a dense submatrix for use in Finite Element-type computations.
Definition: dense_submatrix.h:45
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::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::SECOND
Definition: enum_order.h:43
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::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
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::System::attach_assemble_function
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
Definition: system.C:1755
libMesh::NumericVector::add_vector
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
Definition: numeric_vector.C:363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::System::add_variable
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1069
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
assemble_stokes
void assemble_stokes(EquationSystems &es, const std::string &system_name)
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
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::DenseSubVector< Number >
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::System::variable_type
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2233
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::SparseMatrix::add_matrix
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
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::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEType::default_quadrature_order
Order default_quadrature_order() const
Definition: fe_type.h:353
libMesh::LinearImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
Definition: linear_implicit_system.h:55
libMesh::FIRST
Definition: enum_order.h:42
libMesh::EIGEN_SOLVERS
Definition: enum_solver_package.h:40
libMesh::DenseVector< Number >